ATLAS Offline Software
Loading...
Searching...
No Matches
JetQGTaggerVariableTool Class Reference

Properties: VertexContainer - name of the vertex container EventInfo - name of EventInfo container AssociatedTracks - name for attribute holding the list of associated tracks TVATool - tool to do track-vertex association TrkSelTool - tool to select tracks (none ==> no selection) More...

#include <JetQGTaggerVariableTool.h>

Inheritance diagram for JetQGTaggerVariableTool:
Collaboration diagram for JetQGTaggerVariableTool:

Public Member Functions

 JetQGTaggerVariableTool (const std::string &name)
StatusCode initialize () override
 Dummy implementation of the initialisation function.
virtual StatusCode decorate (const xAOD::JetContainer &jetCont) const override
 Decorate a jet collection without otherwise modifying it.
const xAOD::VertexfindHSVertex (const xAOD::VertexContainer *&) const
virtual void print () const
 Print the state of the tool.
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysInitialize () override
 Perform system initialization for an algorithm.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const
virtual StatusCode modify (xAOD::JetContainer &jets) const override final
 Concrete implementation of the function inherited from IJetModifier.
Additional helper functions, not directly mimicking Athena
template<class T>
const T * getProperty (const std::string &name) const
 Get one of the tool's properties.
const std::string & msg_level_name () const __attribute__((deprecated))
 A deprecated function for getting the message level's name.
const std::string & getName (const void *ptr) const
 Get the name of an object that is / should be in the event store.
SG::sgkey_t getKey (const void *ptr) const
 Get the (hashed) key of an object that is in the event store.

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

Gaudi::Property< std::string > m_jetContainerName {this,"JetContainer", "", "SG key for the input jet container"}
SG::ReadHandleKey< xAOD::VertexContainerm_vertexContainer_key {this, "VertexContainer", "PrimaryVertices", "SG key for input vertex container"}
SG::ReadHandleKey< xAOD::EventInfom_eventInfo_key {this, "EventInfo", "EventInfo", "SG key for input EventInfo"}
ToolHandle< InDet::IInDetTrackSelectionToolm_trkSelectionTool {this, "TrackSelector", "", "Track selector"}
SG::ReadHandleKey< jet::TrackVertexAssociationm_tva_key {this, "TrackVertexAssociation", "", "Track vertex association key"}
SG::WriteDecorHandleKey< xAOD::JetContainerm_nTrkKey {this, "NTrksDecorName", "DFCommonJets_QGTagger_NTracks", "SG key for output NTracks decoration"}
SG::WriteDecorHandleKey< xAOD::JetContainerm_trkWidthKey {this, "TracksWidthDecorName", "DFCommonJets_QGTagger_TracksWidth", "SG key for output TracksWidth decoration"}
SG::WriteDecorHandleKey< xAOD::JetContainerm_trkC1Key {this, "TracksC1DecorName", "DFCommonJets_QGTagger_TracksC1", "SG key for output TracksC1 decoration"}
SG::WriteDecorHandleKey< xAOD::JetContainerm_nChargedKey {this, "nChargedDecorName", "DFCommonJets_QGTagger_truthjet_nCharged", "SG key for output truthjet_nCharged decoration"}
SG::WriteDecorHandleKey< xAOD::JetContainerm_truthPtKey {this, "truthPtDecorName", "DFCommonJets_QGTagger_truthjet_pt", "SG key for output truthjet_pt decoration"}
SG::WriteDecorHandleKey< xAOD::JetContainerm_truthEtaKey {this, "truthEtaDecorName", "DFCommonJets_QGTagger_truthjet_eta", "SG key for output truthjet_eta decoration"}
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

Properties: VertexContainer - name of the vertex container EventInfo - name of EventInfo container AssociatedTracks - name for attribute holding the list of associated tracks TVATool - tool to do track-vertex association TrkSelTool - tool to select tracks (none ==> no selection)

Definition at line 41 of file JetQGTaggerVariableTool.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< AlgTool > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ JetQGTaggerVariableTool()

JetQGTaggerVariableTool::JetQGTaggerVariableTool ( const std::string & name)

Definition at line 23 of file JetQGTaggerVariableTool.cxx.

24: asg::AsgTool(name){
25}

Member Function Documentation

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ decorate()

StatusCode JetQGTaggerVariableTool::decorate ( const xAOD::JetContainer & jets) const
overridevirtual

Decorate a jet collection without otherwise modifying it.

Implements IJetDecorator.

Definition at line 61 of file JetQGTaggerVariableTool.cxx.

61 {
62
63
64 SG::WriteDecorHandle<xAOD::JetContainer, int> nTrkHandle(m_nTrkKey);
65 SG::WriteDecorHandle<xAOD::JetContainer, float> trkWidthHandle(m_trkWidthKey);
66 SG::WriteDecorHandle<xAOD::JetContainer, float> trkC1Handle(m_trkC1Key);
67 SG::WriteDecorHandle<xAOD::JetContainer, int> nChargedHandle(m_nChargedKey);
68 SG::WriteDecorHandle<xAOD::JetContainer, float> truthPtHandle(m_truthPtKey);
69 SG::WriteDecorHandle<xAOD::JetContainer, float> truthEtaHandle(m_truthEtaKey);
70
71 // Get input vertex collection
72 auto vertexContainer = SG::makeHandle (m_vertexContainer_key);
73 if (!vertexContainer.isValid()){
74 ATH_MSG_ERROR("Invalid VertexContainer datahandle: " << m_vertexContainer_key.key());
75 return StatusCode::FAILURE;
76 }
77
78 const auto *vertices = vertexContainer.cptr();
79
80 ATH_MSG_DEBUG("Successfully retrieved VertexContainer: " << m_vertexContainer_key.key());
81
82 if (vertices->empty() ) {
83 ATH_MSG_WARNING("There are no vertices in the container. Exiting");
84
85 for(const xAOD::Jet* jet : jetCont){
86 nTrkHandle(*jet) = -1;
87 trkWidthHandle(*jet) = -1.;
88 trkC1Handle(*jet) = -1.;
89 nChargedHandle(*jet) = -1;
90 truthPtHandle(*jet) = -1.;
91 truthEtaHandle(*jet) = -1.;
92 }
93
94 return StatusCode::SUCCESS;
95 }
96
97 //Find the HS vertex
98 const xAOD::Vertex* HSvertex = findHSVertex(vertices);
99
100 // Get the track-vertex association
101 auto handle_tva = SG::makeHandle (m_tva_key);
102 if (!handle_tva.isValid()){
103 ATH_MSG_ERROR("Could not retrieve the TrackVertexAssociation: " << m_tva_key.key());
104 return StatusCode::FAILURE;
105 }
106
107 const auto *tva = handle_tva.cptr();
108
109 //Loop over the jets
110 for(const xAOD::Jet* jet : jetCont){
111
112 int nTracksCount = 0;
113 float TracksWidth = 0., SumTracks_pTs = 0., TracksC1 = 0., beta = 0.2;
114
115 //Get ghost-associated tracks
116 std::vector<const xAOD::IParticle*> jettracks;
117 jet->getAssociatedObjects<xAOD::IParticle>(xAOD::JetAttribute::GhostTrack,jettracks);
118
119 //Loop over the tracks
120 std::vector<bool> IsGoodTrack;
121 TLorentzVector tracki_TLV, trackj_TLV;
122 TLorentzVector jet_TLV = jet -> p4();
123 for (size_t i = 0; i < jettracks.size(); i++) {
124
125 if(!jettracks[i]){
126 IsGoodTrack.push_back(false);
127 continue;
128 }
129
130 const xAOD::TrackParticle* trk = static_cast<const xAOD::TrackParticle*>(jettracks[i]);
131
132 // only count tracks with selections
133 // 1) pt>500 MeV
134 // 2) accepted track from InDetTrackSelectionTool with CutLevel==Loose
135 // 3) associated to primary vertex OR within 3mm of the primary vertex
136
137 bool accept = (trk->pt()>500 &&
138 m_trkSelectionTool->accept(*trk) &&
139 (HSvertex == tva->associatedVertex(trk) || (HSvertex != tva->associatedVertex(trk) && std::abs((trk->z0()+trk->vz()-HSvertex->z())*sin(trk->theta()))<3.))
140 );
141
142 IsGoodTrack.push_back(accept);
143 if (!accept) continue;
144
145 nTracksCount++;
146
147 tracki_TLV = trk -> p4();
148 double DR_tracki_jet = tracki_TLV.DeltaR(jet_TLV);
149 TracksWidth += trk -> pt() * DR_tracki_jet;
150 SumTracks_pTs += trk -> pt();
151
152 }// end loop over jettracks
153
154 if(SumTracks_pTs>0.) TracksWidth = TracksWidth / SumTracks_pTs;
155 else TracksWidth = -1.;
156
157 //Calculate C1 from tracks
158 for(size_t i = 0; i < jettracks.size(); i++) {
159 const xAOD::TrackParticle* trki = static_cast<const xAOD::TrackParticle*>(jettracks[i]);
160 if( !( IsGoodTrack.at(i) ) ) continue;
161
162 for(size_t j = i+1; j < jettracks.size(); j++) {
163 const xAOD::TrackParticle* trkj = static_cast<const xAOD::TrackParticle*>(jettracks[j]);
164 if( !( IsGoodTrack.at(j) ) ) continue;
165
166 tracki_TLV = trki -> p4();
167 trackj_TLV = trkj -> p4();
168 double DR_tracki_trackj = tracki_TLV.DeltaR(trackj_TLV);
169 TracksC1 += trki -> pt() * trkj -> pt() * pow( DR_tracki_trackj, beta) ;
170
171 }//end loop over j
172 }//end double loop over ij
173
174 if(SumTracks_pTs>0.) TracksC1 = TracksC1 / ( pow(SumTracks_pTs, 2.) );
175 else TracksC1 = -1.;
176
177 // Add truth variables for QG tagging
178
179 auto eventInfoContainer = SG::makeHandle (m_eventInfo_key);
180 if (!eventInfoContainer.isValid()){
181 ATH_MSG_ERROR("Invalid EventInfo datahandle: " << m_eventInfo_key.key());
182 return StatusCode::FAILURE;
183 }
184
185 const auto *eventInfo = eventInfoContainer.cptr();
186
187 bool isMC = eventInfo->eventType(xAOD::EventInfo::IS_SIMULATION);
188 int tntrk = 0;
189 float truthjet_pt = -999.0;
190 float truthjet_eta = -999.0;
191
192 if(isMC){
193 const xAOD::Jet* tjet=nullptr;
194 static const SG::ConstAccessor< ElementLink<xAOD::JetContainer> > GhostTruthAssociationLinkAcc ("GhostTruthAssociationLink");
195 if(GhostTruthAssociationLinkAcc.isAvailable(*jet) ){
196 ATH_MSG_DEBUG("Accessing GhostTruthAssociationLink: is available");
197 ElementLink<xAOD::JetContainer> truthlink = GhostTruthAssociationLinkAcc(*jet);
198 if(truthlink.isValid() ){
199 ATH_MSG_DEBUG("Accessing GhostTruthAssociationLink: is valid");
200 if(truthlink)
201 tjet = * truthlink;
202 else{
203 ATH_MSG_DEBUG("Skipping...truth link is broken");
204 }//endelse NULL pointer
205 }
206 else {
207 ATH_MSG_DEBUG("Invalid truth link: setting weight to 1");
208 } //endelse isValid
209 } //endif isAvailable
210 else {
211 ATH_MSG_DEBUG("Cannot access truth Link: setting weight to 1");
212 }//endelse isAvailable
213
214 if(tjet){
215 ATH_MSG_DEBUG("Truth Jet: " << tjet->numConstituents());
216
217 truthjet_pt = tjet->pt();
218 truthjet_eta = tjet->eta();
219
220 //Calculate the number of charged final state particles with pT > 500 MeV
221 for (size_t ind = 0; ind < tjet->numConstituents(); ind++) {
222 const xAOD::TruthParticle *part = static_cast<const xAOD::TruthParticle*>(tjet->rawConstituent(ind));
223 ATH_MSG_DEBUG("part: " << part );
224 // dont count invalid truth particles
225 if (!part || ! MC::isStable(part) || HepMC::is_simulation_particle(part)) continue;
226 // pt>500 MeV
227 ATH_MSG_DEBUG("pt: " << (part->pt()) );
228 if( ! (part->pt()>500.) ) continue;
229 // charged
230 ATH_MSG_DEBUG("isCharged: " << (part->isCharged()) );
231 if( !(part->isCharged()) ) continue;
232 tntrk++;
233 }
234 }
235 }
236
237 nTrkHandle(*jet) = nTracksCount;
238 trkWidthHandle(*jet) = TracksWidth;
239 trkC1Handle(*jet) = TracksC1;
240 nChargedHandle(*jet) = tntrk;
241 truthPtHandle(*jet) = truthjet_pt;
242 truthEtaHandle(*jet) = truthjet_eta;
243
244 }
245
246 return StatusCode::SUCCESS;
247}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
constexpr int pow(int base, int exp) noexcept
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexContainer_key
SG::ReadHandleKey< jet::TrackVertexAssociation > m_tva_key
const xAOD::Vertex * findHSVertex(const xAOD::VertexContainer *&) const
SG::WriteDecorHandleKey< xAOD::JetContainer > m_truthEtaKey
SG::WriteDecorHandleKey< xAOD::JetContainer > m_trkC1Key
SG::WriteDecorHandleKey< xAOD::JetContainer > m_trkWidthKey
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfo_key
SG::WriteDecorHandleKey< xAOD::JetContainer > m_nTrkKey
ToolHandle< InDet::IInDetTrackSelectionTool > m_trkSelectionTool
SG::WriteDecorHandleKey< xAOD::JetContainer > m_nChargedKey
SG::WriteDecorHandleKey< xAOD::JetContainer > m_truthPtKey
@ IS_SIMULATION
true: simulation, false: data
size_t numConstituents() const
Number of constituents in this jets (this is valid even when reading a file where the constituents ha...
Definition Jet_v1.cxx:153
virtual double pt() const
The transverse momentum ( ) of the particle.
Definition Jet_v1.cxx:44
const IParticle * rawConstituent(size_t i) const
Direct access to constituents. WARNING expert use only.
Definition Jet_v1.cxx:158
virtual double eta() const
The pseudorapidity ( ) of the particle.
Definition Jet_v1.cxx:49
float z0() const
Returns the parameter.
float theta() const
Returns the parameter, which has range 0 to .
float vz() const
The z origin for the parameters.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
float z() const
Returns the z position.
StatusCode accept(const xAOD::Muon *mu)
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
bool isStable(const T &p)
Identify if the particle is stable, i.e. has not decayed.
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
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.
TruthParticle_v1 TruthParticle
Typedef to implementation.

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ findHSVertex()

const xAOD::Vertex * JetQGTaggerVariableTool::findHSVertex ( const xAOD::VertexContainer *& vertices) const

Definition at line 251 of file JetQGTaggerVariableTool.cxx.

252{
253 for ( size_t iVertex = 0; iVertex < vertices->size(); ++iVertex ) {
254 if(vertices->at(iVertex)->vertexType() == xAOD::VxType::PriVtx) {
255
256 ATH_MSG_VERBOSE("JetQGTaggerVariableTool " << name() << " Found HS vertex at index: "<< iVertex);
257 return vertices->at(iVertex);
258 }
259 }
260 ATH_MSG_VERBOSE("There is no vertex of type PriVx. Taking default vertex.");
261 return vertices->at(0);
262}
#define ATH_MSG_VERBOSE(x)
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.
VxType::VertexType vertexType() const
The type of the vertex.
@ PriVtx
Primary vertex.

◆ getKey()

SG::sgkey_t asg::AsgTool::getKey ( const void * ptr) const
inherited

Get the (hashed) key of an object that is in the event store.

This is a bit of a special one. StoreGateSvc and xAOD::TEvent both provide ways for getting the SG::sgkey_t key for an object that is in the store, based on a bare pointer. But they provide different interfaces for doing so.

In order to allow tools to efficiently perform this operation, they can use this helper function.

See also
asg::AsgTool::getName
Parameters
ptrThe bare pointer to the object that the event store should know about
Returns
The hashed key of the object in the store. If not found, an invalid (zero) key.

Definition at line 119 of file AsgTool.cxx.

119 {
120
121#ifdef XAOD_STANDALONE
122 // In case we use @c xAOD::TEvent, we have a direct function call
123 // for this.
124 return evtStore()->event()->getKey( ptr );
125#else
126 const SG::DataProxy* proxy = evtStore()->proxy( ptr );
127 return ( proxy == nullptr ? 0 : proxy->sgkey() );
128#endif // XAOD_STANDALONE
129 }
ServiceHandle< StoreGateSvc > & evtStore()

◆ getName()

const std::string & asg::AsgTool::getName ( const void * ptr) const
inherited

Get the name of an object that is / should be in the event store.

This is a bit of a special one. StoreGateSvc and xAOD::TEvent both provide ways for getting the std::string name for an object that is in the store, based on a bare pointer. But they provide different interfaces for doing so.

In order to allow tools to efficiently perform this operation, they can use this helper function.

See also
asg::AsgTool::getKey
Parameters
ptrThe bare pointer to the object that the event store should know about
Returns
The string name of the object in the store. If not found, an empty string.

Definition at line 106 of file AsgTool.cxx.

106 {
107
108#ifdef XAOD_STANDALONE
109 // In case we use @c xAOD::TEvent, we have a direct function call
110 // for this.
111 return evtStore()->event()->getName( ptr );
112#else
113 const SG::DataProxy* proxy = evtStore()->proxy( ptr );
114 static const std::string dummy = "";
115 return ( proxy == nullptr ? dummy : proxy->name() );
116#endif // XAOD_STANDALONE
117 }

◆ getProperty()

template<class T>
const T * asg::AsgTool::getProperty ( const std::string & name) const
inherited

Get one of the tool's properties.

◆ initialize()

StatusCode JetQGTaggerVariableTool::initialize ( void )
overridevirtual

Dummy implementation of the initialisation function.

It's here to allow the dual-use tools to skip defining an initialisation function. Since many are doing so...

Reimplemented from asg::AsgTool.

Definition at line 29 of file JetQGTaggerVariableTool.cxx.

29 {
30 ATH_MSG_INFO("Initializing JetQGTaggerVariableTool " << name());
31
32 if(m_jetContainerName.empty()){
33 ATH_MSG_ERROR("JetQGTaggerVariableTool needs to have its input jet container configured!");
34 return StatusCode::FAILURE;
35 }
36
43
44 ATH_CHECK(m_vertexContainer_key.initialize());
45 ATH_CHECK(m_eventInfo_key.initialize());
46
47 ATH_CHECK(m_tva_key.initialize());
48
49 ATH_CHECK(m_nTrkKey.initialize());
50 ATH_CHECK(m_trkWidthKey.initialize());
51 ATH_CHECK(m_trkC1Key.initialize());
52 ATH_CHECK(m_nChargedKey.initialize());
53 ATH_CHECK(m_truthPtKey.initialize());
54 ATH_CHECK(m_truthEtaKey.initialize());
55
56 return StatusCode::SUCCESS;
57}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)
Gaudi::Property< std::string > m_jetContainerName

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ modify()

virtual StatusCode IJetDecorator::modify ( xAOD::JetContainer & jets) const
inlinefinaloverridevirtualinherited

Concrete implementation of the function inherited from IJetModifier.

Implements IJetModifier.

Definition at line 32 of file IJetDecorator.h.

32{return decorate(jets);};
virtual StatusCode decorate(const xAOD::JetContainer &jets) const =0
Decorate a jet collection without otherwise modifying it.

◆ msg()

MsgStream & AthCommonMsg< AlgTool >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msg_level_name()

const std::string & asg::AsgTool::msg_level_name ( ) const
inherited

A deprecated function for getting the message level's name.

Instead of using this, weirdly named function, user code should get the string name of the current minimum message level (in case they really need it...), with:

MSG::name( msg().level() )

This function's name doesn't follow the ATLAS coding rules, and as such will be removed in the not too distant future.

Returns
The string name of the current minimum message level that's printed

Definition at line 101 of file AsgTool.cxx.

101 {
102
103 return MSG::name( msg().level() );
104 }
MsgStream & msg() const
const std::string & name(Level lvl)
Convenience function for translating message levels to strings.
Definition MsgLevel.cxx:19

◆ msgLvl()

bool AthCommonMsg< AlgTool >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ print()

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< AlgTool > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sysInitialize()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysInitialize ( )
overridevirtualinherited

Perform system initialization for an algorithm.

We override this to declare all the elements of handle key arrays at the end of initialization. See comments on updateVHKA.

Reimplemented in asg::AsgMetadataTool, AthCheckedComponent< AthAlgTool >, AthCheckedComponent<::AthAlgTool >, and DerivationFramework::CfAthAlgTool.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_eventInfo_key

SG::ReadHandleKey<xAOD::EventInfo> JetQGTaggerVariableTool::m_eventInfo_key {this, "EventInfo", "EventInfo", "SG key for input EventInfo"}
private

Definition at line 63 of file JetQGTaggerVariableTool.h.

63{this, "EventInfo", "EventInfo", "SG key for input EventInfo"};

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_jetContainerName

Gaudi::Property<std::string> JetQGTaggerVariableTool::m_jetContainerName {this,"JetContainer", "", "SG key for the input jet container"}
private

Definition at line 61 of file JetQGTaggerVariableTool.h.

61{this,"JetContainer", "", "SG key for the input jet container"};

◆ m_nChargedKey

SG::WriteDecorHandleKey<xAOD::JetContainer> JetQGTaggerVariableTool::m_nChargedKey {this, "nChargedDecorName", "DFCommonJets_QGTagger_truthjet_nCharged", "SG key for output truthjet_nCharged decoration"}
private

Definition at line 71 of file JetQGTaggerVariableTool.h.

71{this, "nChargedDecorName", "DFCommonJets_QGTagger_truthjet_nCharged", "SG key for output truthjet_nCharged decoration"};

◆ m_nTrkKey

SG::WriteDecorHandleKey<xAOD::JetContainer> JetQGTaggerVariableTool::m_nTrkKey {this, "NTrksDecorName", "DFCommonJets_QGTagger_NTracks", "SG key for output NTracks decoration"}
private

Definition at line 68 of file JetQGTaggerVariableTool.h.

68{this, "NTrksDecorName", "DFCommonJets_QGTagger_NTracks", "SG key for output NTracks decoration"};

◆ m_trkC1Key

SG::WriteDecorHandleKey<xAOD::JetContainer> JetQGTaggerVariableTool::m_trkC1Key {this, "TracksC1DecorName", "DFCommonJets_QGTagger_TracksC1", "SG key for output TracksC1 decoration"}
private

Definition at line 70 of file JetQGTaggerVariableTool.h.

70{this, "TracksC1DecorName", "DFCommonJets_QGTagger_TracksC1", "SG key for output TracksC1 decoration"};

◆ m_trkSelectionTool

ToolHandle<InDet::IInDetTrackSelectionTool> JetQGTaggerVariableTool::m_trkSelectionTool {this, "TrackSelector", "", "Track selector"}
private

Definition at line 65 of file JetQGTaggerVariableTool.h.

65{this, "TrackSelector", "", "Track selector"};

◆ m_trkWidthKey

SG::WriteDecorHandleKey<xAOD::JetContainer> JetQGTaggerVariableTool::m_trkWidthKey {this, "TracksWidthDecorName", "DFCommonJets_QGTagger_TracksWidth", "SG key for output TracksWidth decoration"}
private

Definition at line 69 of file JetQGTaggerVariableTool.h.

69{this, "TracksWidthDecorName", "DFCommonJets_QGTagger_TracksWidth", "SG key for output TracksWidth decoration"};

◆ m_truthEtaKey

SG::WriteDecorHandleKey<xAOD::JetContainer> JetQGTaggerVariableTool::m_truthEtaKey {this, "truthEtaDecorName", "DFCommonJets_QGTagger_truthjet_eta", "SG key for output truthjet_eta decoration"}
private

Definition at line 73 of file JetQGTaggerVariableTool.h.

73{this, "truthEtaDecorName", "DFCommonJets_QGTagger_truthjet_eta", "SG key for output truthjet_eta decoration"};

◆ m_truthPtKey

SG::WriteDecorHandleKey<xAOD::JetContainer> JetQGTaggerVariableTool::m_truthPtKey {this, "truthPtDecorName", "DFCommonJets_QGTagger_truthjet_pt", "SG key for output truthjet_pt decoration"}
private

Definition at line 72 of file JetQGTaggerVariableTool.h.

72{this, "truthPtDecorName", "DFCommonJets_QGTagger_truthjet_pt", "SG key for output truthjet_pt decoration"};

◆ m_tva_key

SG::ReadHandleKey<jet::TrackVertexAssociation> JetQGTaggerVariableTool::m_tva_key {this, "TrackVertexAssociation", "", "Track vertex association key"}
private

Definition at line 66 of file JetQGTaggerVariableTool.h.

66{this, "TrackVertexAssociation", "", "Track vertex association key"};

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< AlgTool > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vertexContainer_key

SG::ReadHandleKey<xAOD::VertexContainer> JetQGTaggerVariableTool::m_vertexContainer_key {this, "VertexContainer", "PrimaryVertices", "SG key for input vertex container"}
private

Definition at line 62 of file JetQGTaggerVariableTool.h.

62{this, "VertexContainer", "PrimaryVertices", "SG key for input vertex container"};

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< AlgTool > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.


The documentation for this class was generated from the following files: