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

This algorithm allows to read wave forms from ntuples and builds a LArPhysWaveContainer containing the corresponding PhysWave. More...

#include <LArShapeFromStdNtuple.h>

Inheritance diagram for LArShapeFromStdNtuple:
Collaboration diagram for LArShapeFromStdNtuple:

Public Member Functions

 LArShapeFromStdNtuple (const std::string &name, ISvcLocator *pSvcLocator)
virtual ~LArShapeFromStdNtuple ()
virtual StatusCode initialize () override
 implements IAlgorithm::initialize()
virtual StatusCode execute () override
 implements IAlgorithm::execute() : Does nothing
virtual StatusCode finalize () override
virtual StatusCode stop () override
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

unsigned int m_skipPoints
 the first m_skipPoints points of the waveform in the ntuple are skipped
unsigned int m_prefixPoints
 make a Shape with the first m_prefixPoints as zeros
std::vector< std::string > m_root_file_names
 list of input ntuple file names
std::string m_ntuple_name
 ntuple name
std::string m_store_key
 key of the LArShape collection in Storegate
std::string m_groupingType
 Grouping type.
bool m_isComplete
 Shape type.
bool m_done
SG::ReadCondHandleKey< LArMCSymm_mcSymKey {this, "MCSymKey", "LArMCSym", "SG Key of LArMCSym object"}
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

This algorithm allows to read wave forms from ntuples and builds a LArPhysWaveContainer containing the corresponding PhysWave.

Version for standard Ntuple, produced by LArCalibTools algos....

Definition at line 22 of file LArShapeFromStdNtuple.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

◆ LArShapeFromStdNtuple()

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

Definition at line 26 of file LArShapeFromStdNtuple.cxx.

26 : AthAlgorithm(name, pSvcLocator)
27{
28 declareProperty("SkipPoints", m_skipPoints = 0);
29 declareProperty("PrefixPoints", m_prefixPoints = 0);
31 declareProperty("NtupleName", m_ntuple_name="SHAPE");
32 declareProperty("StoreKey", m_store_key="FromStdNtuple");
33 declareProperty("GroupingType", m_groupingType="ExtendedSubDetector");
34 declareProperty("isComplete", m_isComplete=false);
35
36 m_done=false;
37}
AthAlgorithm()
Default constructor:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
unsigned int m_skipPoints
the first m_skipPoints points of the waveform in the ntuple are skipped
unsigned int m_prefixPoints
make a Shape with the first m_prefixPoints as zeros
std::string m_groupingType
Grouping type.
std::vector< std::string > m_root_file_names
list of input ntuple file names
std::string m_store_key
key of the LArShape collection in Storegate
std::string m_ntuple_name
ntuple name

◆ ~LArShapeFromStdNtuple()

LArShapeFromStdNtuple::~LArShapeFromStdNtuple ( )
virtualdefault

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 }

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

virtual StatusCode LArShapeFromStdNtuple::execute ( )
inlineoverridevirtual

implements IAlgorithm::execute() : Does nothing

Definition at line 34 of file LArShapeFromStdNtuple.h.

34{return StatusCode::SUCCESS;}

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

virtual StatusCode LArShapeFromStdNtuple::finalize ( )
inlineoverridevirtual

Definition at line 36 of file LArShapeFromStdNtuple.h.

36{return StatusCode::SUCCESS;}

◆ initialize()

StatusCode LArShapeFromStdNtuple::initialize ( )
overridevirtual

implements IAlgorithm::initialize()

Definition at line 42 of file LArShapeFromStdNtuple.cxx.

43{
44 ATH_CHECK ( m_mcSymKey.initialize() );
45 return StatusCode::SUCCESS ;
46}
#define ATH_CHECK
Evaluate an expression and check for errors.
SG::ReadCondHandleKey< LArMCSym > m_mcSymKey

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

◆ stop()

StatusCode LArShapeFromStdNtuple::stop ( )
overridevirtual

Definition at line 49 of file LArShapeFromStdNtuple.cxx.

50{
51 if(m_done) return StatusCode::SUCCESS;
52
53 ATH_MSG_INFO ( "... in stop()" );
54
55 const EventContext& ctx = Gaudi::Hive::currentContext();
56 SG::ReadCondHandle<LArMCSym> mcsym (m_mcSymKey, ctx);
57
58 // get LArOnlineID helper
59 const LArOnlineID* onlineHelper = nullptr;
60 ATH_CHECK( detStore()->retrieve(onlineHelper, "LArOnlineID") );
61
62 TChain* outfit = new TChain(m_ntuple_name.c_str());
63 for (const std::string& s : m_root_file_names) {
64 outfit->Add(s.c_str());
65 }
66
67
68 Int_t phase=0;
69 Int_t det=0;
70 Float_t phaseTime=0;
71 Float_t timeOffset=0;
72 Int_t channelId=0;
73 Int_t FT=0, slot=0, channel=0;
74
75 Float_t Shape[2000]; // The function
76 Float_t ShapeDer[2000]; // The function
77 Int_t gain = 0; // LARHIGHGAIN = 0, LARMEDIUMGAIN = 1, LARLOWGAIN = 2,
78 Int_t nsamples=0;
79
80 outfit->SetBranchAddress("channelId", &channelId);
81 outfit->SetBranchAddress("FT", &FT);
82 outfit->SetBranchAddress("slot", &slot);
83 outfit->SetBranchAddress("channel", &channel);
84 outfit->SetBranchAddress("detector", &det);
85 if(m_isComplete) {
86 outfit->SetBranchAddress("PhaseTime", &phaseTime);
87 outfit->SetBranchAddress("Phase", &phase);
88 outfit->SetBranchAddress("timeOffset", &timeOffset);
89 }
90 outfit->SetBranchAddress("nSamples", &nsamples);
91 outfit->SetBranchAddress("Gain", &gain);
92 outfit->SetBranchAddress("Shape", Shape);
93 outfit->SetBranchAddress("ShapeDer", ShapeDer);
94
95 Float_t timeBinWidth=-1, prevTime=-1, timeOff=-1;
96 Int_t prevPhase=-1;
97
98 // Create new LArShapeContainer
99 std::unique_ptr<LArShapeComplete> larShapeComplete;
100 std::unique_ptr<LArShape32MC> larShapeMC;
101
102 if(m_isComplete) {
103 larShapeComplete = std::make_unique<LArShapeComplete>();
104 ATH_CHECK( larShapeComplete->setGroupingType(m_groupingType, msg()) );
105 ATH_CHECK( larShapeComplete->initialize() );
106 } else {
107 larShapeMC = std::make_unique<LArShape32MC>();
108 ATH_CHECK( larShapeMC->setGroupingType(m_groupingType, msg()) );
109 ATH_CHECK( larShapeMC->initialize() );
110 }
111
112 std::vector<float> shapemc;
113 std::vector<float> shape_dermc;
114 typedef std::vector<std::vector<float> > wave2d;
115 std::map<std::pair<unsigned int,int>, wave2d> shape;
116 std::map<std::pair<unsigned int,int>, wave2d> shape_der;
117 unsigned int hwid;
118 // loop over entries in the Tuple, one entry = one channel
119 Long64_t nentries = outfit->GetEntries();
120 for ( Long64_t iev = 0; iev < nentries; ++iev )
121 {
122 outfit->GetEvent(iev);
123 ATH_MSG_INFO ( " Chan " << std::hex << channelId << " det. "<< det << std::dec );
124
125 hwid = channelId;
126 //if (det != 4) hwid = (channelId<<4);
127 HWIdentifier id(hwid);
128
129 /*
130 if(FT != onlineHelper->feedthrough(id) || slot != onlineHelper->slot(id) || channel != onlineHelper->channel(id)) {
131 ATH_MSG_ERROR ( "Inconsistency in decoding HWID !!!!" );
132 ATH_MSG_INFO ( "Trying to fix..." );
133 hwid = (channelId<<4);
134 id=HWIdentifier(hwid);
135 */
136 if(FT != onlineHelper->feedthrough(id) || slot != onlineHelper->slot(id) || channel != onlineHelper->channel(id)) {
137 ATH_MSG_ERROR ( "Inconsistency in decoding HWID !!!!" );
138 ATH_MSG_ERROR ( FT << " - " << onlineHelper->feedthrough(id) );
139 ATH_MSG_ERROR ( slot << " - " << onlineHelper->slot(id) );
140 ATH_MSG_ERROR ( channel << " - " << onlineHelper->channel(id) );
141 //if(det == 4) {
142 // ATH_MSG_ERROR ( "Ignoring for sFcal" );
143 //} else {
144 // ATH_MSG_ERROR ( "Not creating Shape !!!!" );
145 // continue;
146 // }
147 continue;
148 /* }
149 ATH_MSG_INFO ( "Fixed....." );
150 */
151 }
152
153 // Catch potential array index out of range error.
154 if ( nsamples >= 2000 ) {
155 ATH_MSG_ERROR ( " Too many points specified vs the expected content of the ntuple ! " );
156 ATH_MSG_ERROR ( "Not creating Shape !!!!" );
157 continue;
158 } else {
159 ATH_MSG_DEBUG ( "Working with " << nsamples << " samples" );
160 }
161 if(timeBinWidth < 0) {
162 if(prevTime < 0) {
163 prevTime=phaseTime;
164 prevPhase=phase;
165 } else {
166 if(abs(phase-prevPhase) == 1) {
167 timeBinWidth=std::fabs(phaseTime - prevTime);
168 }
169 }
170 }
171
172 if(timeOff < 0) timeOff=timeOffset;
173
174 if(m_isComplete) {
175 if(shape[std::make_pair(hwid,gain)].empty()) shape[std::make_pair(hwid,gain)].reserve(50);
176 if(shape_der[std::make_pair(hwid,gain)].empty()) shape_der[std::make_pair(hwid,gain)].reserve(50);
177 shape[std::make_pair(hwid,gain)][phase].reserve(nsamples);
178 shape_der[std::make_pair(hwid,gain)][phase].reserve(nsamples);
179 for(int i=0;i<nsamples; ++i) {shape[std::make_pair(hwid,gain)][phase][i]=0.;
180 shape_der[std::make_pair(hwid,gain)][phase][i]=0.;}
181 } else {
182 shapemc.resize(nsamples);
183 shape_dermc.resize(nsamples);
184 for(int i=0;i<nsamples; ++i) {shapemc[i]=0.; shape_dermc[i]=0.;}
185 }
186 unsigned int skipped = 0;
187 unsigned int limit = nsamples;
189 for ( unsigned int i = 0; i < limit; i++ ) {
190 if ( skipped >= m_skipPoints ) {
191 if(m_isComplete) { // accumulate into map
192 shape[std::make_pair(hwid,gain)][phase][i-m_skipPoints+m_prefixPoints]=Shape[i];
193 shape_der[std::make_pair(hwid,gain)][phase][i-m_skipPoints+m_prefixPoints]=ShapeDer[i];
194
195 } else {
196 shapemc[i-m_skipPoints+m_prefixPoints]=Shape[i];
197 shape_dermc[i-m_skipPoints+m_prefixPoints]=ShapeDer[i];
198 }
199 } else skipped++;
200 }
201
202 if(!m_isComplete) {
203 if (id != mcsym->ZPhiSymOnl(id) ) {
204 ATH_MSG_INFO( "Symmetrized, not stored" );
205 } else {
206
207 ATH_MSG_INFO( "Storing shape with length: " << shapemc.size() );
208
209 //LArShapeP1 t(shapemc,shape_dermc);
210
211 //larShapeMC->setPdata(id, LArShapeP1(shapemc,shape_dermc),gain);
212 larShapeMC->set(id, gain, shapemc, shape_dermc);
213 ATH_MSG_INFO( larShapeMC->Shape(id,gain).size() << " " << larShapeMC->ShapeDer(id,gain).size() );
214 ATH_MSG_INFO( "Shape[2] =" << larShapeMC->Shape(id,gain)[2] << "shapemc[2] =" << shapemc[2] );
215 }
216 }
217
218 }
219 // for complete, could fill only now
220 if(m_isComplete){
221 std::map<std::pair<unsigned int,int>, wave2d>::iterator ibeg = shape.begin();
222 std::map<std::pair<unsigned int,int>, wave2d>::iterator iend = shape.end();
223 for(;ibeg != iend; ++ibeg) {
224 larShapeComplete->set(HWIdentifier((ibeg->first).first),(ibeg->first).second,
225 ibeg->second, shape_der[std::make_pair((ibeg->first).first, (ibeg->first).second)],
226 timeOff,timeBinWidth);
227 }
228 }
229
230 if(m_isComplete) {
231 ATH_CHECK( detStore()->record(std::move(larShapeComplete),m_store_key) );
232 } else {
233 ATH_CHECK( detStore()->record(std::move(larShapeMC),m_store_key) );
234 }
235 return StatusCode::SUCCESS;
236}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
static const Attributes_t empty
const ServiceHandle< StoreGateSvc > & detStore() const
MsgStream & msg() const
int feedthrough(const HWIdentifier id) const
Return the feedthrough of a hardware cell identifier : feedthrough = [0,31] Barrel - A/C side or H/...
int slot(const HWIdentifier id) const
Return the slot number of a hardware cell identifier: slot = [1,15] Slot-ID in top part of the crat...
int channel(const HWIdentifier id) const
Return the channel number of a hardware cell identifier channel = [0,127] in all FEB.
retrieve(aClass, aKey=None)
Definition PyKernel.py:110

◆ 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_WARNING(x)
static Double_t sc
virtual StatusCode sysInitialize() override
Override sysInitialize.
AthCommonDataStore(const std::string &name, T... args)
virtual std::vector< Gaudi::DataHandle * > outputHandles() const override
::StatusCode StatusCode
StatusCode definition for legacy code.

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

bool LArShapeFromStdNtuple::m_done
private

Definition at line 55 of file LArShapeFromStdNtuple.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_groupingType

std::string LArShapeFromStdNtuple::m_groupingType
private

Grouping type.

Definition at line 51 of file LArShapeFromStdNtuple.h.

◆ m_isComplete

bool LArShapeFromStdNtuple::m_isComplete
private

Shape type.

Definition at line 53 of file LArShapeFromStdNtuple.h.

◆ m_mcSymKey

SG::ReadCondHandleKey<LArMCSym> LArShapeFromStdNtuple::m_mcSymKey {this, "MCSymKey", "LArMCSym", "SG Key of LArMCSym object"}
private

Definition at line 57 of file LArShapeFromStdNtuple.h.

58{this, "MCSymKey", "LArMCSym", "SG Key of LArMCSym object"};

◆ m_ntuple_name

std::string LArShapeFromStdNtuple::m_ntuple_name
private

ntuple name

Definition at line 47 of file LArShapeFromStdNtuple.h.

◆ m_prefixPoints

unsigned int LArShapeFromStdNtuple::m_prefixPoints
private

make a Shape with the first m_prefixPoints as zeros

Definition at line 43 of file LArShapeFromStdNtuple.h.

◆ m_root_file_names

std::vector<std::string> LArShapeFromStdNtuple::m_root_file_names
private

list of input ntuple file names

Definition at line 45 of file LArShapeFromStdNtuple.h.

◆ m_skipPoints

unsigned int LArShapeFromStdNtuple::m_skipPoints
private

the first m_skipPoints points of the waveform in the ntuple are skipped

Definition at line 41 of file LArShapeFromStdNtuple.h.

◆ m_store_key

std::string LArShapeFromStdNtuple::m_store_key
private

key of the LArShape collection in Storegate

Definition at line 49 of file LArShapeFromStdNtuple.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.


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