ATLAS Offline Software
Loading...
Searching...
No Matches
LVL1::JGTowerBuilder Class Reference

#include <JGTowerBuilder.h>

Inheritance diagram for LVL1::JGTowerBuilder:
Collaboration diagram for LVL1::JGTowerBuilder:

Public Member Functions

 JGTowerBuilder (const std::string &name, ISvcLocator *pSvcLocator)
virtual ~JGTowerBuilder () override
virtual StatusCode initialize () override
virtual StatusCode execute (const EventContext &ctx) const override
virtual StatusCode sysInitialize () override
 Override sysInitialize.
virtual bool isClonable () const override
 Specify if the algorithm is clonable.
virtual unsigned int cardinality () const override
 Cardinality (Maximum number of clones that can exist) special value 0 means that algorithm is reentrant.
virtual StatusCode sysExecute (const EventContext &ctx) override
 Execute an algorithm.
virtual const DataObjIDColl & extraOutputDeps () const override
 Return the list of extra output dependencies.
virtual bool filterPassed (const EventContext &ctx) const
virtual void setFilterPassed (bool state, const EventContext &ctx) 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 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

SG::ReadHandleKey< CaloCellContainerm_cellsKey
SG::ReadHandleKey< CaloCellContainerm_superCellsKey
SG::ReadHandleKey< xAOD::TriggerTowerContainerm_triggerTowersKey
SG::WriteHandleKey< xAOD::JGTowerContainerm_towersKey
SG::ReadCondHandleKey< JGTowerMappingDatam_mappingKey
bool m_useAllCalo
bool m_useAllRun2TT
bool m_emulateSCTiming
bool m_useSCQuality
uint16_t m_scQuality
float m_minSCETp
float m_maxSCETm
float m_minTowerEt
const CaloCell_SuperCell_IDm_scid = nullptr
DataObjIDColl m_extendedExtraObjects
 Extra output dependency collection, extended by AthAlgorithmDHUpdate to add symlinks.
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 20 of file JGTowerBuilder.h.

Member Typedef Documentation

◆ StoreGateSvc_t

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

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ JGTowerBuilder()

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

Definition at line 23 of file JGTowerBuilder.cxx.

24 : AthReentrantAlgorithm(name, pSvcLocator)
25 {
26 declareProperty("InputCells", m_cellsKey = "AllCalo");
27 declareProperty("InputSuperCells", m_superCellsKey = "SCell");
28 declareProperty("InputTriggerTowers", m_triggerTowersKey = "xAODTriggerTowers");
29 declareProperty("UseAllCalo", m_useAllCalo = false);
30 declareProperty("UseAllRun2TT", m_useAllRun2TT = false);
31 declareProperty("OutputTowers", m_towersKey);
32 declareProperty("EmulateSuperCellTiming", m_emulateSCTiming = false);
33 declareProperty("UseSCQuality", m_useSCQuality = true);
34 declareProperty("SuperCellQuality", m_scQuality = 0x200);
35 declareProperty("MinSCETp", m_minSCETp = -1, "Minimum Et for positive energy scells");
36 declareProperty("MaxSCETm", m_maxSCETm = 1, "Maximum Et for negative energy scells");
37 declareProperty("MinTowerET", m_minTowerEt = -9e9);
38 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
SG::ReadHandleKey< CaloCellContainer > m_cellsKey
SG::ReadHandleKey< CaloCellContainer > m_superCellsKey
SG::ReadHandleKey< xAOD::TriggerTowerContainer > m_triggerTowersKey
SG::WriteHandleKey< xAOD::JGTowerContainer > m_towersKey

◆ ~JGTowerBuilder()

LVL1::JGTowerBuilder::~JGTowerBuilder ( )
overridevirtual

Definition at line 40 of file JGTowerBuilder.cxx.

40{}

Member Function Documentation

◆ cardinality()

unsigned int AthCommonReentrantAlgorithm< Gaudi::Algorithm >::cardinality ( ) const
overridevirtualinherited

Cardinality (Maximum number of clones that can exist) special value 0 means that algorithm is reentrant.

Override this to return 0 for reentrant algorithms.

Definition at line 75 of file AthCommonReentrantAlgorithm.cxx.

64{
65 return 0;
66}

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Gaudi::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< Gaudi::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< Gaudi::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< Gaudi::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 LVL1::JGTowerBuilder::execute ( const EventContext & ctx) const
overridevirtual

Definition at line 54 of file JGTowerBuilder.cxx.

55 {
56 ATH_MSG_DEBUG("Executing " << name());
57 auto superCells = SG::makeHandle(m_superCellsKey, ctx);
58 ATH_CHECK(superCells.isValid());
59 const CaloCellContainer *cells = nullptr;
60 const xAOD::TriggerTowerContainer *triggerTowers = nullptr;
61 if (m_useAllCalo)
62 {
63 auto cellsHandle = SG::makeHandle(m_cellsKey, ctx);
64 ATH_CHECK(cellsHandle.isValid());
65 cells = cellsHandle.ptr();
66 }
67 else
68 {
69 auto triggerTowersHandle = SG::makeHandle(m_triggerTowersKey, ctx);
70 ATH_CHECK(triggerTowersHandle.isValid());
71 triggerTowers = triggerTowersHandle.ptr();
72 }
73
74 auto towerHandle = SG::makeHandle(m_towersKey, ctx);
75 auto towers = std::make_unique<xAOD::JGTowerContainer>();
76 auto towersAux = std::make_unique<xAOD::JGTowerAuxContainer>();
77 towers->setStore(towersAux.get());
78
79 SG::ReadCondHandle<JGTowerMappingData> mapping(m_mappingKey, ctx);
80 if (!mapping.isValid())
81 {
82 ATH_MSG_ERROR("Failed to retrieve mapping " << m_mappingKey);
83 return StatusCode::FAILURE;
84 }
85 for (const JGTowerHelper &towerHelper : **mapping)
86 {
87 float LArEt = 0;
88 float tileEt = 0;
90 {
91 for (const xAOD::TriggerTower *triggerTower : *triggerTowers)
92 {
93 if (triggerTower->sampling() != towerHelper.sampling() && std::abs(towerHelper.Eta()) < 3.1)
94 continue;
95 if (!towerHelper.inBox(triggerTower->eta(), triggerTower->phi()))
96 continue;
97 tileEt += 500 * triggerTower->cpET();
98 }
99 }
100 else
101 {
102 // Fill towers from SC
103 for (unsigned int scIdx : towerHelper.GetSCIndices())
104 {
105 Identifier scID = m_scid->cell_id(scIdx);
106 IdentifierHash scHash = m_scid->calo_cell_hash(scID);
107 const CaloCell *superCell = superCells->findCell(scHash);
108 if (!superCell)
109 // TODO: Should this at least be a warning?
110 continue;
111 float et = superCell->et();
112 float time = superCell->time();
113 if (std::isnan(superCell->et()))
114 // TODO: Should this only be a warning, or else a job failure?
115 ATH_MSG_ERROR("Supercell ET is nan. Likely due to BCID correction or something else upstream");
117 {
118 if (
119 (et < 10e3 && std::abs(time) > 8) ||
120 (et >= 10e3 && (time > 16 || time < -8)))
121 continue;
122 }
123 else
124 {
125 if (
126 (m_useSCQuality && !(superCell->provenance() & m_scQuality)) ||
127 (et > 0 && et < m_minSCETp) ||
129 continue;
130 }
131 LArEt += et;
132 }
133
134 if (towerHelper.sampling() == 1)
135 {
136 if (m_useAllCalo)
137 for (unsigned int tileIdx : towerHelper.GetTileIndices())
138 {
139 const CaloCell *cell = cells->findCell(tileIdx);
140 if (!cell)
141 // TODO: Should this be at least a warning?
142 continue;
143 tileEt += cell->e() * cell->sinTh();
144 }
145 else
146 for (const xAOD::TriggerTower *triggerTower : *triggerTowers) {
147 if (std::abs(triggerTower->eta()) < 1.5) {
148 if (triggerTower->sampling() == 1 &&
149 towerHelper.inBox(triggerTower->eta(), triggerTower->phi()))
150 tileEt += 500 * triggerTower->cpET();
151 }
152 }
153 }
154 }
155 float towerEt = LArEt + tileEt;
156 if (towerEt < m_minTowerEt)
157 towerEt = 0;
158 xAOD::JGTower *tower = new xAOD::JGTower();
159 towers->push_back(tower);
160 tower->initialize(towers->size() - 1, towerHelper.Eta(), towerHelper.Phi());
161 tower->setdEta(towerHelper.dEta());
162 tower->setdPhi(towerHelper.dPhi());
163 tower->setEt(towerEt);
164 tower->setSCIndex(towerHelper.GetSCIndices());
165 tower->setTileIndex(towerHelper.GetTileIndices());
166 tower->setSampling(towerHelper.sampling());
167 decArea(*tower) = towerHelper.area();
168 // Decorate the relations
169 decNextEtaIndex(*tower) = towerHelper.nextEtaIndex();
170 decPreviousEtaIndex(*tower) = towerHelper.previousEtaIndex();
171 decNextPhiIndex(*tower) = towerHelper.nextPhiIndex();
172 decPreviousPhiIndex(*tower) = towerHelper.previousPhiIndex();
173 decIndexInFront(*tower) = towerHelper.indexInFront();
174 decIndexBehind(*tower) = towerHelper.indexBehind();
175 } //> end loop over tower helpers
176
177 ATH_CHECK(towerHandle.record(std::move(towers), std::move(towersAux)));
178 return StatusCode::SUCCESS;
179 }
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
float et(const xAOD::jFexSRJetRoI *j)
float time() const
get time (data member)
Definition CaloCell.h:368
uint16_t provenance() const
get provenance (data member)
Definition CaloCell.h:354
virtual double et() const override final
get et
Definition CaloCell.h:423
const CaloCell_SuperCell_ID * m_scid
SG::ReadCondHandleKey< JGTowerMappingData > m_mappingKey
void setdEta(float)
void setSCIndex(const std::vector< int > &)
set SCIndex
void setEt(float)
void setTileIndex(const std::vector< int > &)
set TileIndex
void setdPhi(float)
void setSampling(int)
virtual void initialize(const int Id, const float Eta, const float Phi, const float Et)
time(flags, cells_name, *args, **kw)
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
JGTower_v1 JGTower
Define the latest version of the JGTower class.
Definition JGTower.h:15
TriggerTowerContainer_v2 TriggerTowerContainer
Define the latest version of the TriggerTower container.
TriggerTower_v2 TriggerTower
Define the latest version of the TriggerTower class.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< Gaudi::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 & AthCommonReentrantAlgorithm< Gaudi::Algorithm >::extraOutputDeps ( ) const
overridevirtualinherited

Return the list of extra output dependencies.

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

Definition at line 94 of file AthCommonReentrantAlgorithm.cxx.

90{
91 // If we didn't find any symlinks to add, just return the collection
92 // from the base class. Otherwise, return the extended collection.
93 if (!m_extendedExtraObjects.empty()) {
95 }
97}
An algorithm that can be simultaneously executed in multiple threads.

◆ filterPassed()

virtual bool AthCommonReentrantAlgorithm< Gaudi::Algorithm >::filterPassed ( const EventContext & ctx) const
inlinevirtualinherited

Definition at line 96 of file AthCommonReentrantAlgorithm.h.

96 {
97 return execState( ctx ).filterPassed();
98 }
virtual bool filterPassed(const EventContext &ctx) const

◆ initialize()

StatusCode LVL1::JGTowerBuilder::initialize ( )
overridevirtual

Definition at line 42 of file JGTowerBuilder.cxx.

43 {
44 ATH_CHECK(m_superCellsKey.initialize());
47 ATH_CHECK(m_towersKey.initialize());
48 ATH_CHECK(m_mappingKey.initialize());
49
51 return StatusCode::SUCCESS;
52 }
const ServiceHandle< StoreGateSvc > & detStore() const
retrieve(aClass, aKey=None)
Definition PyKernel.py:110

◆ inputHandles()

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

◆ isClonable()

◆ msg()

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

Definition at line 24 of file AthCommonMsg.h.

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

◆ msgLvl()

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

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ setFilterPassed()

virtual void AthCommonReentrantAlgorithm< Gaudi::Algorithm >::setFilterPassed ( bool state,
const EventContext & ctx ) const
inlinevirtualinherited

Definition at line 100 of file AthCommonReentrantAlgorithm.h.

100 {
102 }
virtual void setFilterPassed(bool state, const EventContext &ctx) const

◆ sysExecute()

StatusCode AthCommonReentrantAlgorithm< Gaudi::Algorithm >::sysExecute ( const EventContext & ctx)
overridevirtualinherited

Execute an algorithm.

We override this in order to work around an issue with the Algorithm base class storing the event context in a member variable that can cause crashes in MT jobs.

Definition at line 85 of file AthCommonReentrantAlgorithm.cxx.

77{
78 return BaseAlg::sysExecute (ctx);
79}

◆ sysInitialize()

StatusCode AthCommonReentrantAlgorithm< Gaudi::Algorithm >::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< Gaudi::Algorithm > >.

Reimplemented in HypoBase, and InputMakerBase.

Definition at line 61 of file AthCommonReentrantAlgorithm.cxx.

107 {
109
110 if (sc.isFailure()) {
111 return sc;
112 }
113
114 ServiceHandle<ICondSvc> cs("CondSvc",name());
115 for (auto h : outputHandles()) {
116 if (h->isCondition() && h->mode() == Gaudi::DataHandle::Writer) {
117 // do this inside the loop so we don't create the CondSvc until needed
118 if ( cs.retrieve().isFailure() ) {
119 ATH_MSG_WARNING("no CondSvc found: won't autoreg WriteCondHandles");
120 return StatusCode::SUCCESS;
121 }
122 if (cs->regHandle(this,*h).isFailure()) {
124 ATH_MSG_ERROR("unable to register WriteCondHandle " << h->fullKey()
125 << " with CondSvc");
126 }
127 }
128 }
129 return sc;
130}
#define ATH_MSG_WARNING(x)
virtual std::vector< Gaudi::DataHandle * > outputHandles() const override

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< Gaudi::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< Gaudi::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 }

Member Data Documentation

◆ m_cellsKey

SG::ReadHandleKey<CaloCellContainer> LVL1::JGTowerBuilder::m_cellsKey
private

Definition at line 30 of file JGTowerBuilder.h.

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_emulateSCTiming

bool LVL1::JGTowerBuilder::m_emulateSCTiming
private

Definition at line 38 of file JGTowerBuilder.h.

◆ m_evtStore

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

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_extendedExtraObjects

DataObjIDColl AthCommonReentrantAlgorithm< Gaudi::Algorithm >::m_extendedExtraObjects
privateinherited

Extra output dependency collection, extended by AthAlgorithmDHUpdate to add symlinks.

Empty if no symlinks were found.

Definition at line 114 of file AthCommonReentrantAlgorithm.h.

◆ m_mappingKey

SG::ReadCondHandleKey<JGTowerMappingData> LVL1::JGTowerBuilder::m_mappingKey
private
Initial value:
{
this, "MappingData", ""}

Definition at line 34 of file JGTowerBuilder.h.

34 {
35 this, "MappingData", ""};

◆ m_maxSCETm

float LVL1::JGTowerBuilder::m_maxSCETm
private

Definition at line 42 of file JGTowerBuilder.h.

◆ m_minSCETp

float LVL1::JGTowerBuilder::m_minSCETp
private

Definition at line 41 of file JGTowerBuilder.h.

◆ m_minTowerEt

float LVL1::JGTowerBuilder::m_minTowerEt
private

Definition at line 43 of file JGTowerBuilder.h.

◆ m_scid

const CaloCell_SuperCell_ID* LVL1::JGTowerBuilder::m_scid = nullptr
private

Definition at line 46 of file JGTowerBuilder.h.

◆ m_scQuality

uint16_t LVL1::JGTowerBuilder::m_scQuality
private

Definition at line 40 of file JGTowerBuilder.h.

◆ m_superCellsKey

SG::ReadHandleKey<CaloCellContainer> LVL1::JGTowerBuilder::m_superCellsKey
private

Definition at line 31 of file JGTowerBuilder.h.

◆ m_towersKey

SG::WriteHandleKey<xAOD::JGTowerContainer> LVL1::JGTowerBuilder::m_towersKey
private

Definition at line 33 of file JGTowerBuilder.h.

◆ m_triggerTowersKey

SG::ReadHandleKey<xAOD::TriggerTowerContainer> LVL1::JGTowerBuilder::m_triggerTowersKey
private

Definition at line 32 of file JGTowerBuilder.h.

◆ m_useAllCalo

bool LVL1::JGTowerBuilder::m_useAllCalo
private

Definition at line 36 of file JGTowerBuilder.h.

◆ m_useAllRun2TT

bool LVL1::JGTowerBuilder::m_useAllRun2TT
private

Definition at line 37 of file JGTowerBuilder.h.

◆ m_useSCQuality

bool LVL1::JGTowerBuilder::m_useSCQuality
private

Definition at line 39 of file JGTowerBuilder.h.

◆ m_varHandleArraysDeclared

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

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

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

Definition at line 398 of file AthCommonDataStore.h.


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