ATLAS Offline Software
Loading...
Searching...
No Matches
JetTrackSumMomentsTool.cxx
Go to the documentation of this file.
1
2
3/*
4 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
5*/
6
7// JetTrackSumMomentsTool.cxx
8// Implementation file for class JetTrackSumMomentsTool
9// Author: James Frost <james.frost@cern.ch>
11
13
16
18//using AODTracking::FourMom_t;
19//using xAOD::TrackParticle::FourMom_t;
20
22 : asg::AsgTool(name) {
23}
24
25//**********************************************************************
26
28 ATH_MSG_INFO("Initializing JetTrackSumMomentsTool " << name());
29 if ( m_htsel.empty() ) {
30 ATH_MSG_INFO(" No track selector.");
31 } else {
32 ATH_MSG_INFO(" Track selector: " << m_htsel->name());
33 }
34
35 ATH_CHECK(m_vertexContainer_key.initialize());
36 ATH_CHECK(m_tva_key.initialize());
37
38 if(m_jetContainerName.empty()) {
39 ATH_MSG_ERROR("JetTrackSumMomentsTool needs to have its input jet container configured!");
40 return StatusCode::FAILURE;
41 }
44
45 ATH_CHECK(m_trackSumPtKey.initialize());
46 ATH_CHECK(m_trackSumMassKey.initialize());
47
48 return StatusCode::SUCCESS;
49}
50
51//**********************************************************************
52
54 // Get input vertex collection
55
57 if (!handle_v.isValid()){
58 ATH_MSG_ERROR("Could not retrieve the VertexContainer: "
59 << m_vertexContainer_key.key());
60 return StatusCode::FAILURE;
61 }
62
63 const auto *vertexContainer = handle_v.cptr();
64
65 // Get the track-vertex association
66 auto handle_tva = SG::makeHandle(m_tva_key);
67 if (!handle_tva.isValid()){
68 ATH_MSG_ERROR("Could not retrieve the TrackVertexAssociation: "
69 << m_tva_key.key());
70 return StatusCode::FAILURE;
71 }
72
73 const auto *tva = handle_tva.cptr();
74
77
78 // Get the tracks associated to the jet
79 // Note that there may be no tracks - this is both normal and an error case
80 std::vector<const xAOD::TrackParticle*> tracks;
81 for(const xAOD::Jet* jet : jets) {
82 if ( ! jet->getAssociatedObjects(m_assocTracksName, tracks) ) {
83 ATH_MSG_DEBUG("Associated tracks not found.");
84 }
85
86 if (vertexContainer->empty() ) {
87 ATH_MSG_WARNING("There are no vertices in the container. Exiting");
88 return StatusCode::FAILURE;
89 }
90
91 const xAOD::Vertex* HSvertex = nullptr;
93 HSvertex = findHSVertex(vertexContainer);
94 else
95 {
96 HSvertex = jet->getAssociatedObject<xAOD::Vertex>("OriginVertex");
97 if (!HSvertex) // nullptr if the attribute doesn't exist
98 {
99 ATH_MSG_ERROR("OriginVertex was requested, but the jet does not contain an OriginVertex");
100 return StatusCode::FAILURE;
101 }
102 else
103 ATH_MSG_VERBOSE("JetTrackSumMomentsTool " << name() << " is using OriginVertex at index: " << HSvertex->index());
104 }
105 const std::pair<float,float> tracksums = getJetTrackSums(HSvertex, tracks, tva);
106
107 trackSumPtHandle(*jet) = tracksums.first;
108 trackSumMassHandle(*jet) = tracksums.second;
109
110 ATH_MSG_VERBOSE("JetTrackSumMomentsTool " << name()
111 << ": TrackSumPt=" << tracksums.first
112 << ", TrackSumMass=" << tracksums.second );
113 }
114 return StatusCode::SUCCESS;
115}
116
117//**********************************************************************
118
120{
121 for ( size_t iVertex = 0; iVertex < vertices->size(); ++iVertex ) {
122 if(vertices->at(iVertex)->vertexType() == xAOD::VxType::PriVtx) {
123
124 ATH_MSG_VERBOSE("JetVertexTaggerTool " << name() << " Found HS vertex at index: "<< iVertex);
125 return vertices->at(iVertex);
126 }
127 }
128 ATH_MSG_VERBOSE("There is no vertex of type PriVx. Taking default vertex.");
129 return vertices->at(0);
130}
131
132//**********************************************************************
133
134std::pair<float,float> JetTrackSumMomentsTool::getJetTrackSums(const xAOD::Vertex* vertex,
135 const std::vector<const xAOD::TrackParticle*>& tracks,
136 const jet::TrackVertexAssociation* tva) const {
137
138 bool notsel = m_htsel.empty();
139 unsigned int nkeep = 0;
140 unsigned int nskip = 0;
141 xAOD::TrackParticle::FourMom_t tracksum(0,0,0,0);
142 for (size_t iTrack = 0; iTrack < tracks.size(); ++iTrack) {
143 const xAOD::TrackParticle* track = tracks.at(iTrack);
144 if ( notsel || m_htsel->accept(*track) ) {
145 const xAOD::Vertex* ptvtx = tva->associatedVertex(track);
146 if (!m_requireTrackPV || ptvtx != nullptr ) { // Track has vertex assigned
147 // Add to sums
148 if (!m_requireTrackPV || ptvtx->index() == vertex->index() ) { tracksum += track->p4(); }
149 }
150 ++nkeep;
151 }
152 else { ++nskip; }
153 }
154
155 ATH_MSG_VERBOSE("JetTrackSumMomentsTool " << name()
156 << ": nsel=" << nkeep
157 << ", nrej=" << nskip );
158
159 return std::make_pair(tracksum.Pt(),tracksum.M());
160
161}
162
163//**********************************************************************
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Handle class for reading from StoreGate.
Handle class for adding a decoration to an object.
const T * at(size_type n) const
Access an element, as an rvalue.
size_type size() const noexcept
Returns the number of elements in the collection.
SG::ReadHandleKey< jet::TrackVertexAssociation > m_tva_key
ToolHandle< InDet::IInDetTrackSelectionTool > m_htsel
const xAOD::Vertex * findHSVertex(const xAOD::VertexContainer *&) const
virtual StatusCode decorate(const xAOD::JetContainer &jets) const override
Decorate a jet collection without otherwise modifying it.
SG::WriteDecorHandleKey< xAOD::JetContainer > m_trackSumMassKey
Gaudi::Property< bool > m_requireTrackPV
Gaudi::Property< bool > m_useOriginVertex
Gaudi::Property< std::string > m_jetContainerName
SG::WriteDecorHandleKey< xAOD::JetContainer > m_trackSumPtKey
std::pair< float, float > getJetTrackSums(const xAOD::Vertex *, const std::vector< const xAOD::TrackParticle * > &, const jet::TrackVertexAssociation *) const
Gaudi::Property< std::string > m_assocTracksName
JetTrackSumMomentsTool(const std::string &name)
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexContainer_key
size_t index() const
Return the index of this element within its container.
Handle class for adding a decoration to an object.
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
Class to hold N-to-one aassociations between tracks and vertices.
const xAOD::Vertex * associatedVertex(const xAOD::TrackParticle *trk) const
IParticle::FourMom_t FourMom_t
Definition of the 4-momentum type.
VxType::VertexType vertexType() const
The type of the vertex.
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
@ PriVtx
Primary vertex.
Jet_v1 Jet
Definition of the current "jet version".
TrackParticle_v1 TrackParticle
Reference the current persistent version:
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
JetContainer_v1 JetContainer
Definition of the current "jet container version".
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< double > > JetFourMom_t
Base 4 Momentum type for Jet.
Definition JetTypes.h:17