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

#include <AODReader.h>

Inheritance diagram for AODReader:
Collaboration diagram for AODReader:

Public Member Functions

 AODReader (const std::string &name, ISvcLocator *pSvcLocator)
virtual ~AODReader ()
virtual StatusCode initialize ()
virtual StatusCode execute ()
virtual StatusCode finalize ()
virtual StatusCode sysInitialize () override
 Override sysInitialize.
virtual const DataObjIDColl & extraOutputDeps () const override
 Return the list of extra output dependencies.
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 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

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

NTuple::Tuple * m_nt
NTuple::Item< float > m_energy
NTuple::Item< float > m_eta
NTuple::Item< float > m_pt
NTuple::Item< float > m_e237
NTuple::Item< float > m_e277
NTuple::Item< float > m_weta1
NTuple::Item< float > m_weta2
NTuple::Item< float > m_e2tsts1
NTuple::Item< float > m_fracs1
NTuple::Item< float > m_wtots1
NTuple::Item< float > m_f1
NTuple::Item< float > m_f1core
NTuple::Item< float > m_et
NTuple::Item< float > m_ethad1
NTuple::Item< float > m_emins1
NTuple::Item< float > m_truth_energy
NTuple::Item< float > m_truth_px
NTuple::Item< float > m_truth_py
NTuple::Item< float > m_truth_pz
DataObjIDColl m_extendedExtraObjects
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 14 of file AODReader.h.

Member Typedef Documentation

◆ StoreGateSvc_t

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

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ AODReader()

AODReader::AODReader ( const std::string & name,
ISvcLocator * pSvcLocator )

Definition at line 19 of file AODReader.cxx.

20 : AthAlgorithm(name, pSvcLocator),
21 m_nt(nullptr)
22{
23
24}
NTuple::Tuple * m_nt
Definition AODReader.h:27
AthAlgorithm()
Default constructor:

◆ ~AODReader()

AODReader::~AODReader ( )
virtual

Definition at line 26 of file AODReader.cxx.

26{ }

Member Function Documentation

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::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< Algorithm > >::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< Algorithm > >::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< Algorithm > >::evtStore ( )
inlineinherited

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

Definition at line 85 of file AthCommonDataStore.h.

◆ execute()

StatusCode AODReader::execute ( )
virtual

read the AOD electron container from persistecy storage

loop over the AOD electron container, take the first one (only one?) with the correct author

Definition at line 105 of file AODReader.cxx.

106{
107 StatusCode sc = StatusCode::SUCCESS;
108
110 const ElectronContainer* elecTES = nullptr;
111 sc=evtStore()->retrieve( elecTES, "ElectronAODCollection");
112 if( sc.isFailure() || !elecTES ) {
113 msg(MSG::WARNING) << "No AOD electron container found in TDS" << endmsg;
114 return sc;
115 }
116 msg(MSG::INFO) << "ElectronContainer successfully retrieved" << endmsg;
117
118 const McEventCollection * mcEvtColl = nullptr;
119
120 if (( mcEvtColl = evtStore()->retrieve<const McEventCollection>() )) {
121 msg(MSG::INFO) << "TruthParticleContainer successfully retrieved" << endmsg;
122 } else {
123 msg(MSG::WARNING) << "Could not retrieve TruthParticleContainer" << endmsg;
124 }
125
126 ElectronContainer::const_iterator elecItr = elecTES->begin();
127 ElectronContainer::const_iterator elecItrE = elecTES->end();
128
129 HepMC::ConstGenParticlePtr trPart{nullptr};
130 if (mcEvtColl) {
131 McEventCollection::const_iterator mcTrPart = mcEvtColl->begin();
132 if (mcTrPart != mcEvtColl->end()) {
133 trPart = (*mcTrPart)->particles_size()?*HepMC::begin(**mcTrPart):nullptr;
134 if (!trPart) {
135 msg(MSG::WARNING) << "Not a single particle event. Truth information won't be available" << endmsg;
136 }
137 } else {
138 msg(MSG::WARNING) << "Empty TruthParticleContainer. Truth information won't be available" << endmsg;
139 }
140 }
141
142 const Analysis::Electron * primElec = nullptr;
143
144
147 for (; elecItr != elecItrE; ++elecItr) {
148 if (((*elecItr)->author() & (8 + 1)) > 0) {
149 primElec = (*elecItr);
150
151 const EMShower * emshower = primElec->detail<EMShower> ("egDetailAOD");
152 if (emshower) { //exist for non-forward electrons
153 m_e237 = emshower->e237();
154 m_e277 = emshower->e277();
155 m_weta1 = emshower->weta1();
156 m_weta2 = emshower->weta2();
157 m_e2tsts1 = emshower->e2tsts1();
158 m_emins1 = emshower->emins1();
159 m_fracs1 = emshower->fracs1();
160 m_wtots1 = emshower->wtots1();
161 m_f1 = emshower->f1();
162 m_f1core = emshower->f1core();
163 m_ethad1 = emshower->ethad1();
164 } else {
165 m_e237 = 0;
166 m_e277 = 0;
167 m_weta1 = 0;
168 m_weta2 = 0;
169 m_e2tsts1 = 0;
170 m_emins1 = 0;
171 m_fracs1 = 0;
172 m_wtots1 = 0;
173 m_f1 = 0;
174 m_f1core = 0;
175 m_ethad1 = 0;
176 m_et = 0;
177 }
178 if (primElec->cluster()) { // don't know if it's possible for this not to exist, but just in case
179 m_et = primElec->cluster()->et();
180 } else {
181 m_et = 0;
182 }
183 m_energy = primElec->e();
184 m_eta = primElec->eta();
185 m_pt = primElec->pt();
186 if (trPart) {
187 m_truth_energy = trPart->momentum().e();
188 m_truth_px = trPart->momentum().px();
189 m_truth_py = trPart->momentum().py();
190 m_truth_pz = trPart->momentum().pz();
191 } else {
192 m_truth_energy = 0;
193 m_truth_px = 0;
194 m_truth_py = 0;
195 m_truth_pz = 0;
196 }
197 // writing a record
198 ATH_CHECK(ntupleSvc()->writeRecord(m_nt));
199 msg(MSG::VERBOSE) << "n tuple entry writed" << endmsg;
200 }
201 }
202
203 return StatusCode::SUCCESS;
204}
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
static Double_t sc
INTupleSvc * ntupleSvc()
NTuple::Item< float > m_truth_pz
Definition AODReader.h:47
NTuple::Item< float > m_e237
Definition AODReader.h:31
NTuple::Item< float > m_et
Definition AODReader.h:40
NTuple::Item< float > m_eta
Definition AODReader.h:29
NTuple::Item< float > m_weta1
Definition AODReader.h:33
NTuple::Item< float > m_f1
Definition AODReader.h:38
NTuple::Item< float > m_truth_px
Definition AODReader.h:45
NTuple::Item< float > m_truth_energy
Definition AODReader.h:44
NTuple::Item< float > m_e277
Definition AODReader.h:32
NTuple::Item< float > m_wtots1
Definition AODReader.h:37
NTuple::Item< float > m_truth_py
Definition AODReader.h:46
NTuple::Item< float > m_weta2
Definition AODReader.h:34
NTuple::Item< float > m_emins1
Definition AODReader.h:42
NTuple::Item< float > m_f1core
Definition AODReader.h:39
NTuple::Item< float > m_energy
Definition AODReader.h:28
NTuple::Item< float > m_fracs1
Definition AODReader.h:36
NTuple::Item< float > m_ethad1
Definition AODReader.h:41
NTuple::Item< float > m_pt
Definition AODReader.h:30
NTuple::Item< float > m_e2tsts1
Definition AODReader.h:35
MsgStream & msg() const
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
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.
double f1core() const
E1(3x1)/E = fraction of the energy reconstructed in the first longitudinal compartment of the electro...
double weta1() const
shower width using +/-3 strips around the one with the maximal energy deposit: w3 strips = sqrt{sum(E...
double ethad1() const
transverse energy in the first sampling of the hadronic calorimeters behind the cluster calculated fr...
double e237() const
uncalibrated energy (sum of cells) of the middle sampling in a rectangle of size 3x7
double e277() const
uncalibrated energy (sum of cells) of the middle sampling in a rectangle of size 7x7
double fracs1() const
shower shape in the shower core : [E(+/-3)-E(+/-1)]/E(+/-1), where E(+/-n) is the energy in +- n stri...
double weta2() const
the lateral width is calculated with a window of 3x5 cells using the energy weighted sum over all cel...
double e2tsts1() const
energy of the cell corresponding to second energy maximum in the first sampling
double emins1() const
energy reconstructed in the strip with the minimal value between the first and second maximum
double wtots1() const
shower width is determined in a window detaxdphi = 0,0625 x~0,2, corresponding typically to 20 strips...
double f1() const
E1/E = fraction of energy reconstructed in the first sampling, where E1 is energy in all strips belon...
virtual double et() const
transverse energy defined to be e*sin(theta)
virtual double e() const
energy
virtual double pt() const
transverse momentum
virtual double eta() const
pseudo rapidity
const T * detail(const std::string &name="", unsigned int index=0) const
retrieve eg-detail objects:
Definition egamma.h:318
const CaloCluster * cluster() const
pointer to CaloCluster
Definition egamma.cxx:360
ElectronContainer
::StatusCode StatusCode
StatusCode definition for legacy code.
GenEvent::particle_iterator begin(HepMC::GenEvent &e)
Definition GenEvent.h:622
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::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

◆ extraOutputDeps()

const DataObjIDColl & AthAlgorithm::extraOutputDeps ( ) const
overridevirtualinherited

Return the list of extra output dependencies.

This list is extended to include symlinks implied by inheritance relations.

Definition at line 50 of file AthAlgorithm.cxx.

51{
52 // If we didn't find any symlinks to add, just return the collection
53 // from the base class. Otherwise, return the extended collection.
54 if (!m_extendedExtraObjects.empty()) {
56 }
57 return Algorithm::extraOutputDeps();
58}
DataObjIDColl m_extendedExtraObjects

◆ finalize()

StatusCode AODReader::finalize ( )
virtual

Definition at line 99 of file AODReader.cxx.

100{
101 msg(MSG::INFO) << "Finalised " << name() << endmsg;
102 return StatusCode::SUCCESS;
103}

◆ initialize()

StatusCode AODReader::initialize ( )
virtual

Definition at line 28 of file AODReader.cxx.

29{
30 msg(MSG::INFO) << "Initializing " << name() << endmsg;
31
32
33 msg(MSG::INFO) << "creating n tuple" << endmsg;
34 NTupleFilePtr file1(ntupleSvc(), "/NTUPLES/FILE");
35
36 //Checking if the n tuple is already booked
37 NTuplePtr nt(ntupleSvc(),"/NTUPLES/FILE/AOD");
38
39 if ( !nt ) { // Normally, this should always be true, but who knows...
40 nt = ntupleSvc()->book ("/NTUPLES/FILE/AOD", CLID_ColumnWiseTuple, "Readed info from aod container");
41
42 if ( nt ) {
43 msg(MSG::INFO) << "booked ntuple " << endmsg;
44
45 StatusCode sc = StatusCode::SUCCESS;
46
47 sc = nt->addItem("eta" , this->m_eta);
48 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for eta" << endmsg;
49 sc = nt->addItem("energy" , this->m_energy);
50 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for erec" << endmsg;
51 sc = nt->addItem("pt" , this->m_pt);
52 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for pt" << endmsg;
53 sc = nt->addItem("e237" , this->m_e237);
54 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for e237" << endmsg;
55 sc = nt->addItem("e277" , this->m_e277);
56 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for e277" << endmsg;
57 sc = nt->addItem("weta1" , this->m_weta1);
58 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for weta1" << endmsg;
59 sc = nt->addItem("weta2" , this->m_weta2);
60 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for weta2" << endmsg;
61 sc = nt->addItem("e2tsts1" , this->m_e2tsts1);
62 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for e2tsts1" << endmsg;
63 sc = nt->addItem("emins1" , this->m_emins1);
64 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for emins1" << endmsg;
65 sc = nt->addItem("fracs1" , this->m_fracs1);
66 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for fracs1" << endmsg;
67 sc = nt->addItem("wtots1" , this->m_wtots1);
68 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for wtots1" << endmsg;
69 sc = nt->addItem("f1" , this->m_f1);
70 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for f1" << endmsg;
71 sc = nt->addItem("f1core" , this->m_f1core);
72 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for f1core" << endmsg;
73 sc = nt->addItem("ethad1" , this->m_ethad1);
74 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for ethad1" << endmsg;
75 sc = nt->addItem("et" , this->m_et);
76 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for et" << endmsg;
77
78 sc = nt->addItem("truth_energy" , this->m_truth_energy);
79 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for egen" << endmsg;
80 sc = nt->addItem("truth_px" , this->m_truth_px);
81 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for tr_px" << endmsg;
82 sc = nt->addItem("truth_py" , this->m_truth_py);
83 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for tr_py" << endmsg;
84 sc = nt->addItem("truth_pz" , this->m_truth_pz);
85 if( sc.isFailure() ) msg(MSG::ERROR) << "Can not create ntuple entry for tr_pz" << endmsg;
86 }
87
88 } else {
89 msg(MSG::ERROR) << "didn't manage to book ntuple" << endmsg;
90 return StatusCode::FAILURE;
91 }
92
93
94 m_nt = nt;
95
96 return StatusCode::SUCCESS;
97}

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::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.

◆ msg()

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

Definition at line 24 of file AthCommonMsg.h.

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

◆ msgLvl()

bool AthCommonMsg< Algorithm >::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< Algorithm > >::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.

◆ 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< Algorithm > >::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< Algorithm > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sysInitialize()

StatusCode AthAlgorithm::sysInitialize ( )
overridevirtualinherited

Override sysInitialize.

Override sysInitialize from the base class.

Loop through all output handles, and if they're WriteCondHandles, automatically register them and this Algorithm with the CondSvc

Scan through all outputHandles, and if they're WriteCondHandles, register them with the CondSvc

Reimplemented from AthCommonDataStore< AthCommonMsg< Algorithm > >.

Reimplemented in AthAnalysisAlgorithm, AthFilterAlgorithm, AthHistogramAlgorithm, and PyAthena::Alg.

Definition at line 66 of file AthAlgorithm.cxx.

66 {
68
69 if (sc.isFailure()) {
70 return sc;
71 }
72 ServiceHandle<ICondSvc> cs("CondSvc",name());
73 for (auto h : outputHandles()) {
74 if (h->isCondition() && h->mode() == Gaudi::DataHandle::Writer) {
75 // do this inside the loop so we don't create the CondSvc until needed
76 if ( cs.retrieve().isFailure() ) {
77 ATH_MSG_WARNING("no CondSvc found: won't autoreg WriteCondHandles");
78 return StatusCode::SUCCESS;
79 }
80 if (cs->regHandle(this,*h).isFailure()) {
81 sc = StatusCode::FAILURE;
82 ATH_MSG_ERROR("unable to register WriteCondHandle " << h->fullKey()
83 << " with CondSvc");
84 }
85 }
86 }
87 return sc;
88}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
virtual StatusCode sysInitialize() override
Override sysInitialize.
AthCommonDataStore(const std::string &name, T... args)
virtual std::vector< Gaudi::DataHandle * > outputHandles() const override

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< Algorithm > >::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< Algorithm > >::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< Algorithm > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_e237

NTuple::Item<float> AODReader::m_e237
private

Definition at line 31 of file AODReader.h.

◆ m_e277

NTuple::Item<float> AODReader::m_e277
private

Definition at line 32 of file AODReader.h.

◆ m_e2tsts1

NTuple::Item<float> AODReader::m_e2tsts1
private

Definition at line 35 of file AODReader.h.

◆ m_emins1

NTuple::Item<float> AODReader::m_emins1
private

Definition at line 42 of file AODReader.h.

◆ m_energy

NTuple::Item<float> AODReader::m_energy
private

Definition at line 28 of file AODReader.h.

◆ m_et

NTuple::Item<float> AODReader::m_et
private

Definition at line 40 of file AODReader.h.

◆ m_eta

NTuple::Item<float> AODReader::m_eta
private

Definition at line 29 of file AODReader.h.

◆ m_ethad1

NTuple::Item<float> AODReader::m_ethad1
private

Definition at line 41 of file AODReader.h.

◆ m_evtStore

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

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_extendedExtraObjects

DataObjIDColl AthAlgorithm::m_extendedExtraObjects
privateinherited

Definition at line 79 of file AthAlgorithm.h.

◆ m_f1

NTuple::Item<float> AODReader::m_f1
private

Definition at line 38 of file AODReader.h.

◆ m_f1core

NTuple::Item<float> AODReader::m_f1core
private

Definition at line 39 of file AODReader.h.

◆ m_fracs1

NTuple::Item<float> AODReader::m_fracs1
private

Definition at line 36 of file AODReader.h.

◆ m_nt

NTuple::Tuple* AODReader::m_nt
private

Definition at line 27 of file AODReader.h.

◆ m_pt

NTuple::Item<float> AODReader::m_pt
private

Definition at line 30 of file AODReader.h.

◆ m_truth_energy

NTuple::Item<float> AODReader::m_truth_energy
private

Definition at line 44 of file AODReader.h.

◆ m_truth_px

NTuple::Item<float> AODReader::m_truth_px
private

Definition at line 45 of file AODReader.h.

◆ m_truth_py

NTuple::Item<float> AODReader::m_truth_py
private

Definition at line 46 of file AODReader.h.

◆ m_truth_pz

NTuple::Item<float> AODReader::m_truth_pz
private

Definition at line 47 of file AODReader.h.

◆ m_varHandleArraysDeclared

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

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

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

Definition at line 398 of file AthCommonDataStore.h.

◆ m_weta1

NTuple::Item<float> AODReader::m_weta1
private

Definition at line 33 of file AODReader.h.

◆ m_weta2

NTuple::Item<float> AODReader::m_weta2
private

Definition at line 34 of file AODReader.h.

◆ m_wtots1

NTuple::Item<float> AODReader::m_wtots1
private

Definition at line 37 of file AODReader.h.


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