ATLAS Offline Software
Loading...
Searching...
No Matches
TCTDecorCheckInTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include <cassert>
9
10//Initialize---------------------------------------------------------------
12
13
14 ATH_CHECK( m_particlesKey.initialize() );
15 ATH_CHECK( m_verticesKey.initialize() );
16 ATH_CHECK( m_jetsKey.initialize() );
17
18 //from https://acode-browser1.usatlas.bnl.gov/lxr/source/athena/Event/xAOD/xAODTrackingCnv/src/TrackParticleCnvAlg.cxx
19 m_trackReadDecorKeyTCTScore = "InDetTrackParticles.TCTScore_"+m_jetsKey.key();
20 m_trackReadDecorKeyJetLink = "InDetTrackParticles.TCTJetLink_"+m_jetsKey.key();
21
22 m_jetReadDecorKeyTCTScore = m_jetsKey.key()+ ".TCTScore";
23 m_jetReadDecorKeyTrackLink = m_jetsKey.key() +".TCTTrackLink";
24
27
30
31 //-------
32 //check that the TrackClassificationTool can be accessed
33 if (m_trackClassificationTool.retrieve().isFailure()) {
34 ATH_MSG_DEBUG("Could not find InDet::InDetTrkInJetType");
35 return StatusCode::SUCCESS;
36 } else {
37 ATH_MSG_DEBUG("InDet::InDetTrkInJetType found");
38 }
39
40 return StatusCode::SUCCESS;
41 }
42
43 StatusCode TCTDecorCheckInTool::execute(const EventContext& ctx) const
44 {
45 ATH_MSG_DEBUG( "Executing..." );
48
49 //JetRead handles
52
53 // Retrieve the track particles:
55 if ( !trackTES.isValid() ) {
56 ATH_MSG_WARNING( "No TrackParticle container found in TDS" );
57 return StatusCode::SUCCESS; }
58 ATH_MSG_DEBUG( "TrackParticleContainer successfully retrieved" );
59
61 if ( !pvTES.isValid() ) {
62 ATH_MSG_WARNING( "No Primary Vertices container found in TDS" );
63 return StatusCode::SUCCESS; }
64 ATH_MSG_DEBUG( "Primary Vertices container successfully retrieved" );
65 const xAOD::Vertex *primVertex=*(pvTES->begin());
66
67
68
69 //==========================================================================
71 if ( !jetTES.isValid() ) {
72 ATH_MSG_WARNING( "No AntiKt4EMPflow jet container found in TDS" );
73 return StatusCode::SUCCESS; }
74 ATH_MSG_DEBUG( "AntiKt4EMPflow jet container successfully retrieved" );
75
76 xAOD::TrackParticleContainer::const_iterator trackItr = trackTES->begin();
77 xAOD::TrackParticleContainer::const_iterator trackItrE = trackTES->end();
78
79 xAOD::JetContainer::const_iterator jetItr = jetTES->begin();
80 xAOD::JetContainer::const_iterator jetItrE = jetTES->end();
81
82 //decorator methods for decorating the tracks
83 if(m_decoratorMethod=="decorateTrack")
84 {
85 ATH_MSG_DEBUG("Using decorateTrack method ");
86 for(trackItr = trackTES->begin(); trackItr != trackItrE; ++trackItr){
87 const xAOD::TrackParticle* itrk = (*trackItr);
88 //first decorate the track with the InDetTrkInJetType method
89 //find closest jet to the track in Delta R
90 float minDeltaR=1.0;
91 const xAOD::Jet* closestJet = *(jetTES->begin());
92 for(jetItr = jetTES->begin(); jetItr != jetItrE; ++jetItr){
93 const xAOD::Jet* curJet = (*jetItr);
94 float curDeltaR = (itrk)->p4().DeltaR(curJet->p4());
95 if(curDeltaR < minDeltaR) {minDeltaR = curDeltaR; closestJet = curJet;}
96 }
97 m_trackClassificationTool->decorateTrack(ctx, itrk, *primVertex, *jetTES, closestJet);
98 }
99
100 //loop over tracks and check if decoration was correctly added (using either decorateTrack)
101 for(trackItr = trackTES->begin(); trackItr != trackItrE; ++trackItr){
102 const xAOD::TrackParticle* itrk = (*trackItr);
103 const std::vector<float>& v_tctScoresDeco = trackReadDecorHandleTCTScore(*itrk);
104 const ElementLink<xAOD::JetContainer>& v_jetLinks = trackReadDecorHandleJetLink(*itrk);
105
106 ATH_MSG_DEBUG("TCT score from decoration: " << v_tctScoresDeco.at(0) << ", " << v_tctScoresDeco.at(1) << ", "<< v_tctScoresDeco.at(2));
107 std::vector<float> v_tctScore = m_trackClassificationTool->trkTypeWgts(ctx,itrk,*primVertex,(*v_jetLinks)->p4());
108 ATH_MSG_DEBUG("Calculated TCT score: " << v_tctScore.at(0) << ", " << v_tctScore.at(1) << ", " << v_tctScore.at(2));
109
110 for(int j=0; j<=2 ; j++) {assert(Athena_test::isEqual(v_tctScore.at(j),v_tctScoresDeco.at(j)));}
111
112 }//end track loop
113 }//end if decorator methods for decorating tracks
114 else if(m_decoratorMethod=="decorateJet")
115 {
116 ATH_MSG_DEBUG("Using decorateJets method ");
117 for(jetItr = jetTES->begin(); jetItr != jetItrE; ++jetItr){
118 const xAOD::Jet* ijet = (*jetItr);
119 std::vector<const xAOD::TrackParticle*> trkparticles(0);
120 for(trackItr = trackTES->begin(); trackItr != trackItrE; ++trackItr){
121 const xAOD::TrackParticle* itrk = (*trackItr);
122 if((itrk)->p4().DeltaR(ijet->p4()) < 0.4) {trkparticles.push_back(itrk); }
123 }
124 m_trackClassificationTool->decorateJet(ctx,trkparticles,*trackTES,*primVertex, ijet);
125 }
126
127 //loop over jets and check if decoration was correctly added
128 for(jetItr = jetTES->begin(); jetItr != jetItrE; ++jetItr){
129 const xAOD::Jet* ijet = (*jetItr);
130 const std::vector<std::vector<float>>& v_tctScoresDeco = jetReadDecorHandleTCTScore(*ijet);
131 const std::vector<ElementLink<xAOD::TrackParticleContainer>>& v_trackLinks = jetReadDecorHandleTrackLink(*ijet);
132
133 for(unsigned int i=0; i<v_tctScoresDeco.size(); i++)
134 {
135 ATH_MSG_DEBUG("TCT score from decoration: " << v_tctScoresDeco.at(i).at(0) << ", " << v_tctScoresDeco.at(i).at(1) << ", "<< v_tctScoresDeco.at(i).at(2));
136 std::vector<float> v_tctScore = m_trackClassificationTool->trkTypeWgts(ctx,*v_trackLinks.at(i),*primVertex,ijet->p4());
137 ATH_MSG_DEBUG("Calculated TCT score: " << v_tctScore.at(0) << ", " << v_tctScore.at(1) << ", " << v_tctScore.at(2));
138
139 for(int j=0; j<=2 ; j++) {assert(Athena_test::isEqual(v_tctScore.at(j),v_tctScoresDeco.at(i).at(j)));}
140 }//end of tct vector loop
141 }//end of jet loop
142 }
143 else
144 {
145 ATH_MSG_ERROR("Specified decorator method not implemented in InDetTrkInJetType: " << m_decoratorMethod );
146 }
147
148
149 return StatusCode::SUCCESS;
150 }
151
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
Handle class for reading a decoration on an object.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
ToolHandle< InDet::IInDetTrkInJetType > m_trackClassificationTool
StringProperty m_decoratorMethod
SG::ReadDecorHandleKey< xAOD::JetContainer > m_jetReadDecorKeyTCTScore
The write key for adding TCT score as decoration to Jet objects.
SG::ReadHandleKey< xAOD::JetContainer > m_jetsKey
SG::ReadHandleKey< xAOD::VertexContainer > m_verticesKey
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_particlesKey
virtual StatusCode execute(const EventContext &ctx) const override
SG::ReadDecorHandleKey< xAOD::JetContainer > m_jetReadDecorKeyTrackLink
SG::ReadDecorHandleKey< xAOD::TrackParticleContainer > m_trackReadDecorKeyTCTScore
The read key for adding TCT score as decoration to TrackParticle objects.
virtual StatusCode initialize() override
SG::ReadDecorHandleKey< xAOD::TrackParticleContainer > m_trackReadDecorKeyJetLink
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition Jet_v1.cxx:71
bool isEqual(double x1, double x2, double thresh=1e-6)
Definition FLOATassert.h:26
Jet_v1 Jet
Definition of the current "jet version".
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Vertex_v1 Vertex
Define the latest version of the vertex class.