ATLAS Offline Software
Loading...
Searching...
No Matches
Analysis::JpsiFinder_ee Class Reference

#include <JpsiFinder_ee.h>

Inheritance diagram for Analysis::JpsiFinder_ee:
Collaboration diagram for Analysis::JpsiFinder_ee:

Public Member Functions

 JpsiFinder_ee (const std::string &t, const std::string &n, const IInterface *p)
 ~JpsiFinder_ee ()
virtual StatusCode initialize () override
virtual StatusCode performSearch (const EventContext &ctx, xAOD::VertexContainer &vxContainer) const override
std::vector< JpsiEECandidategetPairs (const std::vector< const xAOD::TrackParticle * > &) const
std::vector< JpsiEECandidategetPairs (const std::vector< const xAOD::Electron * > &) const
std::vector< JpsiEECandidategetPairs2Colls (const std::vector< const xAOD::TrackParticle * > &, const std::vector< const xAOD::Electron * > &, bool) const
double getInvariantMass (const JpsiEECandidate &, const std::vector< double > &) const
std::vector< JpsiEECandidateselectCharges (const std::vector< JpsiEECandidate > &, const std::string &) const
std::unique_ptr< xAOD::Vertexfit (const EventContext &ctx, const std::vector< const xAOD::TrackParticle * > &, const xAOD::TrackParticleContainer *importedTrackCollection) const
bool passesEgammaCuts (const xAOD::Electron *) const
bool isContainedIn (const xAOD::TrackParticle *, const xAOD::TrackParticleContainer *) const
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

Static Public Member Functions

static const InterfaceID & interfaceID ()

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

bool m_elel
bool m_eltrk
bool m_trktrk
bool m_allElectrons
bool m_useTrackMeasurement
bool m_diElectrons
double m_trk1M
double m_trk2M
double m_thresholdPt
double m_higherPt
double m_trkThresholdPt
double m_invMassUpper
double m_invMassLower
double m_collAngleTheta
double m_collAnglePhi
double m_Chi2Cut
bool m_oppChOnly
bool m_sameChOnly
bool m_allChCombs
SG::ReadHandleKey< xAOD::ElectronContainerm_electronCollectionKey
SG::ReadHandleKey< xAOD::TrackParticleContainerm_TrkParticleCollection
SG::ReadDecorHandleKey< xAOD::ElectronContainerm_gsfCaloLinkKey
ToolHandle< Trk::IVertexFitterm_iVertexFitter
ToolHandle< Trk::ITrackSelectorToolm_trkSelector
ToolHandle< InDet::VertexPointEstimatorm_vertexEstimator
ServiceHandle< IPartPropSvc > m_partPropSvc {this, "PartPropSvc", "PartPropSvc"}
bool m_egammaCuts
std::string m_elSelection
bool m_doTagAndProbe
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

Definition at line 58 of file JpsiFinder_ee.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

◆ JpsiFinder_ee()

Analysis::JpsiFinder_ee::JpsiFinder_ee ( const std::string & t,
const std::string & n,
const IInterface * p )

Definition at line 92 of file JpsiFinder_ee.cxx.

92 : AthAlgTool(t,n,p),
93 m_elel(true),
94 m_eltrk(false),
95 m_trktrk(false),
96 m_allElectrons(false),
98 m_diElectrons(true),
101 m_thresholdPt(0.0),
102 m_higherPt(0.0),
103 m_trkThresholdPt(0.0),
104 m_invMassUpper(100000.0),
105 m_invMassLower(0.0),
106 m_collAngleTheta(0.0),
107 m_collAnglePhi(0.0),
108 m_Chi2Cut(50.),
109 m_oppChOnly(true),
110 m_sameChOnly(false),
111 m_allChCombs(false),
112 m_electronCollectionKey("Electrons"),
113 m_TrkParticleCollection("InDetTrackParticles"),
114 m_iVertexFitter("Trk::TrkVKalVrtFitter"),
115 m_trkSelector("InDet::TrackSelectorTool"),
116 m_vertexEstimator("InDet::VertexPointEstimator"),
117 m_egammaCuts(true),
118 m_elSelection("d0_or_nod0"),
119 m_doTagAndProbe(false)
120
121 {
122 declareInterface<JpsiFinder_ee>(this);
123 declareProperty("elAndEl",m_elel);
124 declareProperty("elAndTrack",m_eltrk);
125 declareProperty("TrackAndTrack",m_trktrk);
126 declareProperty("allElectrons",m_allElectrons);
127 declareProperty("useElectronTrackMeasurement",m_useTrackMeasurement);
128 declareProperty("assumeDiElectrons",m_diElectrons);
129// declareProperty("electronLHValue",m_electronLHValue);
130 declareProperty("track1Mass",m_trk1M);
131 declareProperty("track2Mass",m_trk2M);
132 declareProperty("elThresholdPt",m_thresholdPt);
133 declareProperty("higherPt",m_higherPt);
134 declareProperty("trackThresholdPt",m_trkThresholdPt);
135 declareProperty("invMassUpper",m_invMassUpper);
136 declareProperty("invMassLower",m_invMassLower);
137 declareProperty("collAngleTheta",m_collAngleTheta);
138 declareProperty("collAnglePhi",m_collAnglePhi);
139 declareProperty("Chi2Cut",m_Chi2Cut);
140 declareProperty("oppChargesOnly",m_oppChOnly);
141 declareProperty("sameChargesOnly",m_sameChOnly);
142 declareProperty("allChargeCombinations",m_allChCombs);
143 declareProperty("electronCollectionKey",m_electronCollectionKey);
144 declareProperty("TrackParticleCollection",m_TrkParticleCollection);
145 declareProperty("TrkVertexFitterTool",m_iVertexFitter);
146 declareProperty("TrackSelectorTool",m_trkSelector);
147 declareProperty("VertexPointEstimator",m_vertexEstimator);
148 declareProperty("useEgammaCuts",m_egammaCuts);
149 declareProperty("doTagAndProbe",m_doTagAndProbe);
150 declareProperty("ElectronSelection",m_elSelection);
151 }
SG::ReadHandleKey< xAOD::ElectronContainer > m_electronCollectionKey
ToolHandle< Trk::ITrackSelectorTool > m_trkSelector
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_TrkParticleCollection
ToolHandle< InDet::VertexPointEstimator > m_vertexEstimator
ToolHandle< Trk::IVertexFitter > m_iVertexFitter
AthAlgTool()
Default constructor:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
constexpr double electronMassInMeV
the mass of the electron (in MeV)

◆ ~JpsiFinder_ee()

Analysis::JpsiFinder_ee::~JpsiFinder_ee ( )

Definition at line 153 of file JpsiFinder_ee.cxx.

153{ }

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 }

◆ 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>

◆ 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

◆ fit()

std::unique_ptr< xAOD::Vertex > Analysis::JpsiFinder_ee::fit ( const EventContext & ctx,
const std::vector< const xAOD::TrackParticle * > & inputTracks,
const xAOD::TrackParticleContainer * importedTrackCollection ) const

Definition at line 389 of file JpsiFinder_ee.cxx.

389 {
390 ATH_MSG_DEBUG("inside JpsiFinder_ee::fit");
391
392 const Trk::Perigee& aPerigee1 = inputTracks[0]->perigeeParameters();
393 const Trk::Perigee& aPerigee2 = inputTracks[1]->perigeeParameters();
394 int sflag = 0;
395 int errorcode = 0;
396 Amg::Vector3D startingPoint = m_vertexEstimator->getCirclesIntersectionPoint(&aPerigee1,&aPerigee2,sflag,errorcode);
397 if (errorcode != 0) {startingPoint(0) = 0.0; startingPoint(1) = 0.0; startingPoint(2) = 0.0;}
398
399 auto myVxCandidate = m_iVertexFitter->fit(ctx, inputTracks, startingPoint);
400 ATH_MSG_DEBUG("Initial fit was a success! " << myVxCandidate);
401 // Added by ASC
402 if(myVxCandidate != 0){
403 std::vector<ElementLink<DataVector<xAOD::TrackParticle> > > newLinkVector;
404 for(unsigned int i=0; i< myVxCandidate->trackParticleLinks().size(); i++){
405 ElementLink<DataVector<xAOD::TrackParticle> > mylink=myVxCandidate->trackParticleLinks()[i]; //makes a copy (non-const)
406 mylink.setStorableObject(*importedTrackCollection, true);
407 newLinkVector.push_back( mylink );
408 ATH_MSG_DEBUG("Set a link!");
409 }
410 myVxCandidate->clearTracks();
411 myVxCandidate->setTrackParticleLinks( newLinkVector );
412 ATH_MSG_DEBUG("Set all links");
413 }
414
415 return myVxCandidate;
416
417 } // End of fit method
#define ATH_MSG_DEBUG(x)
size_t size() const
Number of registered mappings.
Eigen::Matrix< double, 3, 1 > Vector3D
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee

◆ getInvariantMass()

double Analysis::JpsiFinder_ee::getInvariantMass ( const JpsiEECandidate & jpsiIn,
const std::vector< double > & massHypotheses ) const

Definition at line 530 of file JpsiFinder_ee.cxx.

530 {
531
532 // construct 4-vectors from track perigee parameters using given mass hypotheses.
533 // NOTE: in new data model (xAOD) the defining parameters are expressed as perigee parameters w.r.t. the beamspot
534 // NOTE2: TrackParticle::p4() method already returns TLorentzVector, however, we want to enforce our own mass hypothesis
535 auto mu1 = jpsiIn.trackParticle1->genvecP4();
536 auto mu2 = jpsiIn.trackParticle2->genvecP4();
537 mu1.SetM(massHypotheses[0]);
538 mu2.SetM(massHypotheses[1]);
539
540 return (mu1+mu2).M();
541
542 }

◆ getPairs() [1/2]

std::vector< JpsiEECandidate > Analysis::JpsiFinder_ee::getPairs ( const std::vector< const xAOD::Electron * > & electronsIn) const

Definition at line 453 of file JpsiFinder_ee.cxx.

453 {
454
455 std::vector<JpsiEECandidate> myPairs;
456 JpsiEECandidate pair;
457 std::vector<const xAOD::Electron*>::const_iterator outerItr;
458 std::vector<const xAOD::Electron*>::const_iterator innerItr;
459
460 if(electronsIn.size()>=2){
461 for(outerItr=electronsIn.begin();outerItr<electronsIn.end();++outerItr){
462 for(innerItr=(outerItr+1);innerItr!=electronsIn.end();++innerItr){
463 pair.el1 = *innerItr;
464 pair.el2 = *outerItr;
465 pair.pairType = ELEL;
466 myPairs.push_back(pair);
467 }
468 }
469 }
470
471 return(myPairs);
472 }

◆ getPairs() [2/2]

std::vector< JpsiEECandidate > Analysis::JpsiFinder_ee::getPairs ( const std::vector< const xAOD::TrackParticle * > & TracksIn) const

Definition at line 426 of file JpsiFinder_ee.cxx.

426 {
427
428 std::vector<JpsiEECandidate> myPairs;
429 JpsiEECandidate pair;
430 std::vector<const xAOD::TrackParticle*>::const_iterator outerItr;
431 std::vector<const xAOD::TrackParticle*>::const_iterator innerItr;
432
433 if(TracksIn.size()>=2){
434 for(outerItr=TracksIn.begin();outerItr<TracksIn.end();++outerItr){
435 for(innerItr=(outerItr+1);innerItr!=TracksIn.end();++innerItr){
436 pair.trackParticle1 = *innerItr;
437 pair.trackParticle2 = *outerItr;
438 pair.pairType = TRK2;
439 myPairs.push_back(pair);
440 }
441 }
442 }
443
444 return(myPairs);
445 }

◆ getPairs2Colls()

std::vector< JpsiEECandidate > Analysis::JpsiFinder_ee::getPairs2Colls ( const std::vector< const xAOD::TrackParticle * > & tracks,
const std::vector< const xAOD::Electron * > & electrons,
bool tagAndProbe ) const

Definition at line 480 of file JpsiFinder_ee.cxx.

480 {
481
482 std::vector<JpsiEECandidate> myPairs;
483 if (m_TrkParticleCollection.key() == "GSFCaloContainer") {
484 ATH_MSG_FATAL("GSFCaloContainer mode not implemented in getPairs2Colls.");
485 return myPairs;
486 }
487 JpsiEECandidate pair;
488
489 // Unless user is running in tag and probe mode, remove tracks which are also identified as muons
490 std::vector<const xAOD::TrackParticle*> tracksToKeep;
491 if (!tagAndProbe) {
492 if(tracks.size()>=1 && electrons.size()>=1){
493 for (const xAOD::TrackParticle* trk : tracks) {
494 bool trackIsElectron(false);
495 for (const xAOD::Electron* ele : electrons) {
496 if ( ele->trackParticleLink().cachedElement() == trk ) {
497 trackIsElectron=true;
498 break;
499 }
500 }
501 if (!trackIsElectron) tracksToKeep.push_back(trk);
502 }
503 }
504 } else {tracksToKeep = tracks;}
505
506 if(tracksToKeep.size()>=1 && electrons.size()>=1){
507 for (const xAOD::TrackParticle* trk : tracks) {
508 for (const xAOD::Electron* ele : electrons) {
509 pair.el1 = ele;
510 // Muon track 1st
511 pair.trackParticle1 = ele->trackParticleLink().cachedElement();
512 pair.trackParticle2 = trk;
513 pair.pairType = ELTRK;
514 myPairs.push_back(pair);
515 }
516 }
517 }
518
519 return(myPairs);
520 }
#define ATH_MSG_FATAL(x)
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Electron_v1 Electron
Definition of the current "egamma version".

◆ initialize()

StatusCode Analysis::JpsiFinder_ee::initialize ( )
overridevirtual

Definition at line 31 of file JpsiFinder_ee.cxx.

31 {
32
33 // Initialize the ReadHandles
36
37 // Initialize ReadDecorHandles
38 ATH_CHECK(m_gsfCaloLinkKey.initialize());
39
40 // retrieving vertex Fitter
41 ATH_CHECK(m_iVertexFitter.retrieve());
42
43 // Get the track selector tool from ToolSvc
44 ATH_CHECK(m_trkSelector.retrieve());
45
46
47 // Get the vertex point estimator tool from ToolSvc
48 ATH_CHECK(m_vertexEstimator.retrieve());
49
50
51 if (m_diElectrons) {
52 // Get the Particle Properties Service
53 ATH_CHECK(m_partPropSvc.retrieve());
54 auto particleDataTable = m_partPropSvc->PDT();
55 const HepPDT::ParticleData* pd_el = particleDataTable->particle(MC::ELECTRON);
56 m_trk1M = pd_el->mass();
57 m_trk2M = pd_el->mass();
58 }
59
60 if (m_doTagAndProbe) ATH_MSG_WARNING("You have requested tag and probe mode. Duplicate mu+trk pairs WILL be allowed, charge ordering WILL NOT be done. Tag track will be first in each candidate");
61
62
63// // Check that the user's settings are sensible
64 bool illogicalOptions(false);
65 if ( (m_elel && m_eltrk) || (m_elel && m_trktrk) || (m_eltrk && m_trktrk) ) {
66 ATH_MSG_WARNING("You are requesting incompatible combinations of muons and tracks in the pairs. JpsiEECandidates will be EMPTY!");
67 illogicalOptions=true;
68 };
70 ATH_MSG_WARNING("You are requesting Tag and Probe analysis but have not requested mu+trk mode. This is impossible. JpsiEECandidates will be EMPTY!");
71 illogicalOptions=true;
72 };
74 ATH_MSG_WARNING("You are requesting incompatible combinations of charges in the pairs. JpsiEECandidates will be EMPTY!");
75 illogicalOptions=true;
76 };
78 ATH_MSG_WARNING("You are requesting same-sign or all-sign combinations in a tag and probe analysis. This doesn't make sense. JpsiEECandidates will be EMPTY!");
79 illogicalOptions=true;
80 }
81 if (illogicalOptions) return StatusCode::FAILURE;
82
83
84
85 ATH_MSG_DEBUG("Initialize successful");
86
87 return StatusCode::SUCCESS;
88
89 }
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_WARNING(x)
ServiceHandle< IPartPropSvc > m_partPropSvc
SG::ReadDecorHandleKey< xAOD::ElectronContainer > m_gsfCaloLinkKey
static const int ELECTRON

◆ 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.

◆ interfaceID()

const InterfaceID & Analysis::JpsiFinder_ee::interfaceID ( )
inlinestatic

Definition at line 65 of file JpsiFinder_ee.h.

65{ return IID_JpsiFinder_ee;}
static const InterfaceID IID_JpsiFinder_ee("JpsiFinder_ee", 1, 0)

◆ isContainedIn()

bool Analysis::JpsiFinder_ee::isContainedIn ( const xAOD::TrackParticle * theTrack,
const xAOD::TrackParticleContainer * theCollection ) const

Definition at line 614 of file JpsiFinder_ee.cxx.

614 {
615 return std::find(theCollection->begin(), theCollection->end(), theTrack) != theCollection->end();
616 }
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.

◆ msg()

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

Definition at line 24 of file AthCommonMsg.h.

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

◆ 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.

◆ passesEgammaCuts()

bool Analysis::JpsiFinder_ee::passesEgammaCuts ( const xAOD::Electron * electron) const

Definition at line 590 of file JpsiFinder_ee.cxx.

590 {
591
592 static const SG::AuxElement::ConstAccessor<char> isLHVeryLoosenod0("DFCommonElectronsLHVeryLoosenod0");
593 static const SG::AuxElement::ConstAccessor<char> isLHVeryLoose("DFCommonElectronsLHVeryLoose");
594
595 bool passesSelection = false;
596 bool passesLHVLoose = isLHVeryLoose(*electron);
597 bool passesLHVLoosenod0 = isLHVeryLoosenod0(*electron);
598
599 if(m_elSelection == "d0") passesSelection = passesLHVLoose;
600 else if(m_elSelection == "nod0") passesSelection = passesLHVLoosenod0;
601 else if(m_elSelection == "d0_or_nod0") passesSelection = passesLHVLoose || passesLHVLoosenod0;
602 else ATH_MSG_ERROR("Invalid electron selection " << m_elSelection);
603
604 ATH_MSG_DEBUG("Electron with pT, eta: " << electron->pt() << " " << electron->eta() << " passes " << m_elSelection << " " << passesSelection);
605 return passesSelection;
606
607 }
#define ATH_MSG_ERROR(x)
passesSelection(gain, offset)
Definition PlotRamps.py:8

◆ performSearch()

StatusCode Analysis::JpsiFinder_ee::performSearch ( const EventContext & ctx,
xAOD::VertexContainer & vxContainer ) const
overridevirtual

Implements Analysis::ICandidateSearch.

Definition at line 160 of file JpsiFinder_ee.cxx.

161 {
162 ATH_MSG_DEBUG( "JpsiFinder_ee::performSearch" );
163
164
165 // Get the electrons from StoreGate
166 const xAOD::ElectronContainer* importedElectronCollection=nullptr;
167 SG::ReadHandle<xAOD::ElectronContainer> ehandle(m_electronCollectionKey,ctx);
168 if(!ehandle.isValid()){
169 ATH_MSG_WARNING("No electron collection with key " << m_electronCollectionKey.key() << " found in StoreGate. JpsiEECandidates will be EMPTY!");
170 return StatusCode::SUCCESS;;
171 }else{
172 importedElectronCollection = ehandle.cptr();
173 ATH_MSG_DEBUG("Found electron collections with key "<<m_electronCollectionKey.key());
174 }
175 ATH_MSG_DEBUG("Electron container size "<<importedElectronCollection->size());
176
177 // Get ID tracks
178 SG::ReadHandle<xAOD::TrackParticleContainer> thandle(m_TrkParticleCollection,ctx);
179 const xAOD::TrackParticleContainer* importedTrackCollection(0);
180 if(!thandle.isValid()){
181 ATH_MSG_WARNING("No TrackParticle collection with name " << m_TrkParticleCollection << " found in StoreGate!");
182 return StatusCode::SUCCESS;;
183 } else {
184 importedTrackCollection = thandle.cptr();
185 }
186
187 // Typedef for vectors of tracks and muons
188 typedef std::vector<const xAOD::TrackParticle*> TrackBag;
189 typedef std::vector<const xAOD::Electron*> ElectronBag;
190
191 // Select the inner detector tracks
192 const xAOD::Vertex* vx = 0;
193 TrackBag theIDTracksAfterSelection;
194 if (m_trktrk || m_eltrk) {
196 for (trkCItr=importedTrackCollection->begin(); trkCItr!=importedTrackCollection->end(); ++trkCItr) {
197 const xAOD::TrackParticle* TP = (*trkCItr);
198 if ( fabs(TP->pt())<m_trkThresholdPt ) continue;
199 if ( !m_trkSelector->decision(*TP, vx) ) continue;
200 theIDTracksAfterSelection.push_back(TP);
201 }
202 if (theIDTracksAfterSelection.size() == 0) return StatusCode::SUCCESS;;
203 ATH_MSG_DEBUG("Number of tracks after ID track selection: " << theIDTracksAfterSelection.size());
204 }
205
206 // Select the muons
207 ElectronBag theElectronsAfterSelection;
209 if (m_elel || m_eltrk) {
210 for (elItr=importedElectronCollection->begin(); elItr!=importedElectronCollection->end(); ++elItr) {
211 if ( *elItr == NULL ) continue;
212 if (!(*elItr)->trackParticleLink().isValid()) continue; // No electrons without ID tracks
213 const xAOD::TrackParticle* elTrk(0);
214 if (m_TrkParticleCollection.key() == "GSFCaloContainer") {
215 SG::ReadDecorHandle<xAOD::ElectronContainer, ElementLink<xAOD::TrackParticleContainer>>
216 refittedTrackParticleLink(m_gsfCaloLinkKey, ctx);
217 const ElementLink<xAOD::TrackParticleContainer>& refittedTrackLink = refittedTrackParticleLink(*(*elItr));
218 if (!refittedTrackLink.isValid()) continue;
219 elTrk = *refittedTrackLink;
220 } else {
221 if (!(*elItr)->trackParticleLink().isValid()) continue;
222 elTrk = (*elItr)->trackParticleLink().cachedElement();
223 }
224
225 if ( elTrk==NULL) continue;
226 if ( !m_trkSelector->decision(*elTrk, vx) ) continue; // all ID tracks must pass basic tracking cuts
227 if ( fabs(elTrk->pt())<m_thresholdPt ) continue; // higher pt cut if needed
228
229 if ( m_egammaCuts && !passesEgammaCuts(*elItr)) continue; // egamma cuts
230 theElectronsAfterSelection.push_back(*elItr);
231 }
232 if (theElectronsAfterSelection.size() == 0) return StatusCode::SUCCESS;;
233 ATH_MSG_DEBUG("Number of electrons after selection: " << theElectronsAfterSelection.size());
234 }
235
236 // Sort into pairs - end result will be a vector of JpsiEECandidate structs
237 std::vector<JpsiEECandidate> jpsiCandidates;
238 if (m_elel) jpsiCandidates = getPairs(theElectronsAfterSelection);
239 if (m_trktrk) jpsiCandidates = getPairs(theIDTracksAfterSelection);
240 if (m_eltrk) jpsiCandidates = getPairs2Colls(theIDTracksAfterSelection,theElectronsAfterSelection,m_doTagAndProbe);
241
242 ATH_MSG_DEBUG("Number of pairs with ee from a B decay: " << jpsiCandidates.size() );
243
244 // Pair-wise selections
245 std::vector<JpsiEECandidate>::iterator jpsiItr;
246
247 // (1) Enforce one combined muon [deleted, no electron equivalent]
248
249 // (2) Establish track content for candidates
250 // and set the appropriate track collections for the combined muon tracks where appropriate (for saving to persistency later)
251
252 // el+trk or trk+trk - always ID track collection
253 if (m_eltrk || m_trktrk) {
254 for (jpsiItr=jpsiCandidates.begin(); jpsiItr!=jpsiCandidates.end(); ++jpsiItr) {
255 (*jpsiItr).collection1 = importedTrackCollection;
256 (*jpsiItr).collection2 = importedTrackCollection;
257 }
258 }
259
260 if (m_elel) {
261 for (jpsiItr=jpsiCandidates.begin(); jpsiItr!=jpsiCandidates.end(); ++jpsiItr) {
262 if ( m_useTrackMeasurement ) {
263 if (m_TrkParticleCollection.key() == "GSFCaloContainer") {
264 SG::ReadDecorHandle<xAOD::ElectronContainer, ElementLink<xAOD::TrackParticleContainer>>
265 refittedTrackParticleLink(m_gsfCaloLinkKey, ctx);
266 const ElementLink<xAOD::TrackParticleContainer>& refittedTrackLink1 = refittedTrackParticleLink(*((*jpsiItr).el1));
267 const ElementLink<xAOD::TrackParticleContainer>& refittedTrackLink2 = refittedTrackParticleLink(*((*jpsiItr).el2));
268 (*jpsiItr).trackParticle1 = *refittedTrackLink1;
269 (*jpsiItr).trackParticle2 = *refittedTrackLink2;
270 } else {
271 (*jpsiItr).trackParticle1 = (*jpsiItr).el1->trackParticleLink().cachedElement();
272 (*jpsiItr).trackParticle2 = (*jpsiItr).el2->trackParticleLink().cachedElement();
273 }
274 (*jpsiItr).collection1 = importedTrackCollection;
275 (*jpsiItr).collection2 = importedTrackCollection;
276 } else {
277 ATH_MSG_WARNING("Not setup for non-track electron measurements yet....");
278 }
279 } // iteration over candidates
280 }
281
282
283 // (3) Enforce higher track pt if requested
284 if (m_higherPt>0.0) {
285 int index(0);
286 std::vector<int> listToDelete;
287 std::vector<int>::reverse_iterator ii;
288 for(jpsiItr=jpsiCandidates.begin(); jpsiItr!=jpsiCandidates.end();++jpsiItr,++index) {
289 if( (fabs((*jpsiItr).trackParticle1->pt()) < m_higherPt) && (fabs((*jpsiItr).trackParticle2->pt()) < m_higherPt) ) listToDelete.push_back(index);
290 }
291 for (ii=listToDelete.rbegin(); ii!=listToDelete.rend(); ++ii) {
292 jpsiCandidates.erase(jpsiCandidates.begin() + (*ii) );
293 }
294 ATH_MSG_DEBUG("Number of candidates after higherPt cut: " << jpsiCandidates.size() );
295 }
296
297 // (4) Select all opp/same charged track pairs
298 std::vector<JpsiEECandidate> sortedJpsiEECandidates;
299 if (m_oppChOnly) sortedJpsiEECandidates = selectCharges(jpsiCandidates,"OPPOSITE");
300 if (m_sameChOnly) sortedJpsiEECandidates = selectCharges(jpsiCandidates,"SAME");
301 if (m_allChCombs) sortedJpsiEECandidates = selectCharges(jpsiCandidates,"ALL");
302 ATH_MSG_DEBUG("Number of candidates after charge selection: " << sortedJpsiEECandidates.size() );
303
304 // (5) Select for decay angle, if requested
305 if (m_collAnglePhi>0.0 && m_collAngleTheta>0.0) {
306 int index(0);
307 std::vector<int> listToDelete;
308 std::vector<int>::reverse_iterator ii;
309 for(jpsiItr=sortedJpsiEECandidates.begin(); jpsiItr!=sortedJpsiEECandidates.end();++jpsiItr,++index) {
310 double deltatheta = fabs( (*jpsiItr).trackParticle1->theta() - (*jpsiItr).trackParticle2->theta() );
311 // -3.14 < phi < +3.14 ==> correction
312 double deltaphi = std::abs(xAOD::P4Helpers::deltaPhi((*jpsiItr).trackParticle1->phi0() , (*jpsiItr).trackParticle2->phi0()));
313 // perform the angle cuts
314 if ((deltatheta > m_collAngleTheta) || (deltaphi > m_collAnglePhi)) listToDelete.push_back(index);
315 }
316 for (ii=listToDelete.rbegin(); ii!=listToDelete.rend(); ++ii) {
317 sortedJpsiEECandidates.erase(sortedJpsiEECandidates.begin() + (*ii) );
318 }
319 ATH_MSG_DEBUG("Number of collimated candidates: " << sortedJpsiEECandidates.size() );
320 }
321
322 // (6) Select for invariant mass, if requested
323 std::vector<double> trkMasses;
324 trkMasses.push_back(m_trk1M);
325 trkMasses.push_back(m_trk2M);
326 if ( (m_invMassLower > 0.0) || (m_invMassUpper > 0.0) ) {
327 int index(0);
328 std::vector<int> listToDelete;
329 std::vector<int>::reverse_iterator ii;
330 for(jpsiItr=sortedJpsiEECandidates.begin(); jpsiItr!=sortedJpsiEECandidates.end(); ++jpsiItr,++index) {
331 double invMass = getInvariantMass(*jpsiItr,trkMasses);
332 // std::cout << "inv. mass: " << invMass << std::endl;
333 if ( invMass < m_invMassLower || invMass > m_invMassUpper ) {
334 listToDelete.push_back(index);
335 }
336 }
337 for (ii=listToDelete.rbegin(); ii!=listToDelete.rend(); ++ii) {
338 sortedJpsiEECandidates.erase(sortedJpsiEECandidates.begin() + (*ii) );
339 }
340 ATH_MSG_DEBUG("Number of candidates passing invariant mass selection: " << sortedJpsiEECandidates.size() );
341 }
342
343 if (sortedJpsiEECandidates.size() == 0) return StatusCode::SUCCESS;;
344
345 // Fit each pair of tracks to a vertex
346 int itritn = 0;
347 for(jpsiItr=sortedJpsiEECandidates.begin(); jpsiItr!=sortedJpsiEECandidates.end(); ++jpsiItr) {
348 ATH_MSG_DEBUG("jpsiItr: " << itritn); itritn++;
349 std::vector<const xAOD::TrackParticle*> theTracks; theTracks.clear();
350 theTracks.push_back((*jpsiItr).trackParticle1);
351 theTracks.push_back((*jpsiItr).trackParticle2);
352 ATH_MSG_DEBUG("theTracks size (should be two!) " << theTracks.size() << " being vertexed with tracks " << importedTrackCollection);
353 std::unique_ptr<xAOD::Vertex> myVxCandidate{fit(ctx, theTracks,importedTrackCollection)}; // This line actually does the fitting and object making
354 if (myVxCandidate != 0) {
355 // Chi2 cut if requested
356 double chi2 = myVxCandidate->chiSquared();
357 ATH_MSG_DEBUG("chi2 is: " << chi2);
358 if (m_Chi2Cut == 0.0 || chi2 <= m_Chi2Cut) {
359 // decorate the candidate with refitted tracks and muons via the BPhysHelper
360 xAOD::BPhysHelper jpsiHelper(myVxCandidate.get());
361 jpsiHelper.setRefTrks();
362 if (m_elel || m_eltrk) {
363 std::vector<const xAOD::Electron*> theStoredElectrons;
364 theStoredElectrons.push_back((*jpsiItr).el1);
365 if (m_elel) theStoredElectrons.push_back((*jpsiItr).el2);
366 jpsiHelper.setElectrons(theStoredElectrons,importedElectronCollection);
367 }
368 // Retain the vertex
369 vxContainer.push_back(std::move(myVxCandidate));
370 }
371 } else { // fit failed
372 ATH_MSG_DEBUG("Fitter failed!");
373 // Don't try to delete the object, since we arrived here,
374 // because this pointer is null...
375 //delete myVxCandidate;
376 }
377 }
378 ATH_MSG_DEBUG("vxContainer size " << vxContainer.size());
379
380 return StatusCode::SUCCESS;;
381 }
bool passesEgammaCuts(const xAOD::Electron *) const
std::vector< JpsiEECandidate > getPairs(const std::vector< const xAOD::TrackParticle * > &) const
std::unique_ptr< xAOD::Vertex > fit(const EventContext &ctx, const std::vector< const xAOD::TrackParticle * > &, const xAOD::TrackParticleContainer *importedTrackCollection) const
std::vector< JpsiEECandidate > selectCharges(const std::vector< JpsiEECandidate > &, const std::string &) const
double getInvariantMass(const JpsiEECandidate &, const std::vector< double > &) const
std::vector< JpsiEECandidate > getPairs2Colls(const std::vector< const xAOD::TrackParticle * > &, const std::vector< const xAOD::Electron * > &, bool) const
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
value_type push_back(value_type pElem)
Add an element to the end of the collection.
size_type size() const noexcept
Returns the number of elements in the collection.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
double chi2(TH1 *h0, TH1 *h1)
str index
Definition DeMoScan.py:362
std::vector< const xAOD::TrackParticle * > TrackBag
double invMass(const I4Momentum &pA, const I4Momentum &pB)
invariant mass from two I4momentum references
Definition P4Helpers.h:252
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[
ElectronContainer_v1 ElectronContainer
Definition of the current "electron container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".

◆ 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 }

◆ selectCharges()

std::vector< JpsiEECandidate > Analysis::JpsiFinder_ee::selectCharges ( const std::vector< JpsiEECandidate > & jpsisIn,
const std::string & selection ) const

Definition at line 550 of file JpsiFinder_ee.cxx.

550 {
551
552 bool opposite(false),same(false),all(false);
553 if (selection=="OPPOSITE") opposite=true;
554 if (selection=="SAME") same=true;
555 if (selection=="ALL") all=true;
556
557 JpsiEECandidate tmpJpsi;
558 std::vector<JpsiEECandidate> jpsis;
559 double qOverP1=0.;
560 double qOverP2=0.;
561 for(auto jpsiItr=jpsisIn.cbegin();jpsiItr!=jpsisIn.cend();jpsiItr++){
562 bool oppCh(false),sameCh(false);
563 tmpJpsi = *jpsiItr;
564 qOverP1=(*jpsiItr).trackParticle1->qOverP();
565 qOverP2=(*jpsiItr).trackParticle2->qOverP();
566 if(qOverP1*qOverP2<0.0) oppCh=true; // product charge < 0
567 if(qOverP1*qOverP2>0.0) sameCh=true; // product charge > 0
568 // +ve should be first so swap
569 // Don't do it for tag and probe analyses (because tag muon must not change position)
570 if (oppCh && qOverP1<0.0 && !m_doTagAndProbe) {
571 tmpJpsi.trackParticle1 = (*jpsiItr).trackParticle2;
572 tmpJpsi.trackParticle2 = (*jpsiItr).trackParticle1;
573 tmpJpsi.el1 = (*jpsiItr).el2;
574 tmpJpsi.el2 = (*jpsiItr).el1;
575 tmpJpsi.collection1 = (*jpsiItr).collection2;
576 tmpJpsi.collection2 = (*jpsiItr).collection1;
577 }
578 if (oppCh && (opposite || all) ) jpsis.push_back(tmpJpsi);
579 if (sameCh && (same || all) ) jpsis.push_back(tmpJpsi);
580
581 } // end of for loop
582
583 return(jpsis);
584 }
const std::string selection

◆ 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_allChCombs

bool Analysis::JpsiFinder_ee::m_allChCombs
private

Definition at line 99 of file JpsiFinder_ee.h.

◆ m_allElectrons

bool Analysis::JpsiFinder_ee::m_allElectrons
private

Definition at line 84 of file JpsiFinder_ee.h.

◆ m_Chi2Cut

double Analysis::JpsiFinder_ee::m_Chi2Cut
private

Definition at line 96 of file JpsiFinder_ee.h.

◆ m_collAnglePhi

double Analysis::JpsiFinder_ee::m_collAnglePhi
private

Definition at line 95 of file JpsiFinder_ee.h.

◆ m_collAngleTheta

double Analysis::JpsiFinder_ee::m_collAngleTheta
private

Definition at line 94 of file JpsiFinder_ee.h.

◆ 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_diElectrons

bool Analysis::JpsiFinder_ee::m_diElectrons
private

Definition at line 86 of file JpsiFinder_ee.h.

◆ m_doTagAndProbe

bool Analysis::JpsiFinder_ee::m_doTagAndProbe
private

Definition at line 113 of file JpsiFinder_ee.h.

◆ m_egammaCuts

bool Analysis::JpsiFinder_ee::m_egammaCuts
private

Definition at line 111 of file JpsiFinder_ee.h.

◆ m_electronCollectionKey

SG::ReadHandleKey<xAOD::ElectronContainer> Analysis::JpsiFinder_ee::m_electronCollectionKey
private

Definition at line 100 of file JpsiFinder_ee.h.

◆ m_elel

bool Analysis::JpsiFinder_ee::m_elel
private

Definition at line 81 of file JpsiFinder_ee.h.

◆ m_elSelection

std::string Analysis::JpsiFinder_ee::m_elSelection
private

Definition at line 112 of file JpsiFinder_ee.h.

◆ m_eltrk

bool Analysis::JpsiFinder_ee::m_eltrk
private

Definition at line 82 of file JpsiFinder_ee.h.

◆ 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_gsfCaloLinkKey

SG::ReadDecorHandleKey<xAOD::ElectronContainer> Analysis::JpsiFinder_ee::m_gsfCaloLinkKey
private
Initial value:
{
this, "GSFCaloLink", "Electrons.gsfCaloTrackParticleLink",
"ReadHandleKey for electron link to GSFCalo refitted TrackParticle"
}

Definition at line 102 of file JpsiFinder_ee.h.

102 {
103 this, "GSFCaloLink", "Electrons.gsfCaloTrackParticleLink",
104 "ReadHandleKey for electron link to GSFCalo refitted TrackParticle"
105 };

◆ m_higherPt

double Analysis::JpsiFinder_ee::m_higherPt
private

Definition at line 90 of file JpsiFinder_ee.h.

◆ m_invMassLower

double Analysis::JpsiFinder_ee::m_invMassLower
private

Definition at line 93 of file JpsiFinder_ee.h.

◆ m_invMassUpper

double Analysis::JpsiFinder_ee::m_invMassUpper
private

Definition at line 92 of file JpsiFinder_ee.h.

◆ m_iVertexFitter

ToolHandle< Trk::IVertexFitter > Analysis::JpsiFinder_ee::m_iVertexFitter
private

Definition at line 107 of file JpsiFinder_ee.h.

◆ m_oppChOnly

bool Analysis::JpsiFinder_ee::m_oppChOnly
private

Definition at line 97 of file JpsiFinder_ee.h.

◆ m_partPropSvc

ServiceHandle<IPartPropSvc> Analysis::JpsiFinder_ee::m_partPropSvc {this, "PartPropSvc", "PartPropSvc"}
private

Definition at line 110 of file JpsiFinder_ee.h.

110{this, "PartPropSvc", "PartPropSvc"};

◆ m_sameChOnly

bool Analysis::JpsiFinder_ee::m_sameChOnly
private

Definition at line 98 of file JpsiFinder_ee.h.

◆ m_thresholdPt

double Analysis::JpsiFinder_ee::m_thresholdPt
private

Definition at line 89 of file JpsiFinder_ee.h.

◆ m_trk1M

double Analysis::JpsiFinder_ee::m_trk1M
private

Definition at line 87 of file JpsiFinder_ee.h.

◆ m_trk2M

double Analysis::JpsiFinder_ee::m_trk2M
private

Definition at line 88 of file JpsiFinder_ee.h.

◆ m_TrkParticleCollection

SG::ReadHandleKey<xAOD::TrackParticleContainer> Analysis::JpsiFinder_ee::m_TrkParticleCollection
private

Definition at line 101 of file JpsiFinder_ee.h.

◆ m_trkSelector

ToolHandle< Trk::ITrackSelectorTool > Analysis::JpsiFinder_ee::m_trkSelector
private

Definition at line 108 of file JpsiFinder_ee.h.

◆ m_trkThresholdPt

double Analysis::JpsiFinder_ee::m_trkThresholdPt
private

Definition at line 91 of file JpsiFinder_ee.h.

◆ m_trktrk

bool Analysis::JpsiFinder_ee::m_trktrk
private

Definition at line 83 of file JpsiFinder_ee.h.

◆ m_useTrackMeasurement

bool Analysis::JpsiFinder_ee::m_useTrackMeasurement
private

Definition at line 85 of file JpsiFinder_ee.h.

◆ m_varHandleArraysDeclared

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

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vertexEstimator

ToolHandle< InDet::VertexPointEstimator > Analysis::JpsiFinder_ee::m_vertexEstimator
private

Definition at line 109 of file JpsiFinder_ee.h.

◆ 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: