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

Main LVL2 Algorithm. More...

#include <T2CaloEgammaReFastAlgo.h>

Inheritance diagram for T2CaloEgammaReFastAlgo:
Collaboration diagram for T2CaloEgammaReFastAlgo:

Public Member Functions

 T2CaloEgammaReFastAlgo (const std::string &name, ISvcLocator *pSvcLocator)
virtual StatusCode initialize () override
virtual StatusCode execute (const EventContext &context) 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

Gaudi::Property< bool > m_doForward
 calculate zo mass
Gaudi::Property< bool > m_doCalibWithRings
Gaudi::Property< bool > m_useRings
Gaudi::Property< float > m_l1eta {this, "L1ForceEta", -10.0, "Forced LVL1 eta"}
Gaudi::Property< float > m_l1phi {this, "L1ForcePhi", -10.0, "Forced LVL1 phi"}
Gaudi::Property< double > m_etaWidth {this, "EtaWidth", 0.1, "Eta Width of the Region of Interest"}
Gaudi::Property< double > m_phiWidth {this, "PhiWidth", 0.1, "Phi Width of the Region of Interest"}
Gaudi::Property< bool > m_storeCells
ToolHandleArray< IEgammaCalibrationm_calibsBarrel
ToolHandleArray< IEgammaCalibrationm_calibsEndcap
Gaudi::Property< std::vector< float > > m_rhoEta
Gaudi::Property< std::vector< float > > m_zEta
ToolHandleArray< IReAlgToolCalom_emAlgTools
SG::ReadHandleKey< CaloBCIDAveragem_bcidAvgKey
SG::ReadHandleKey< TrigRoiDescriptorCollectionm_roiCollectionKey
SG::WriteHandleKey< xAOD::TrigEMClusterContainerm_clusterContainerKey
ToolHandle< GenericMonitoringToolm_monTool { this, "MonTool", "", "Monitoring tool" }
ToolHandle< TrigFastCalibWithRingsm_calibWRingsTool { this, "CalibWRingsTool", "", "FastCalo Calib with Rings tool" }
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

Main LVL2 Algorithm.

Processes LVL1 information, call FEX IReAlgToolCalos and produces the TrigEMCluster output.

Definition at line 39 of file T2CaloEgammaReFastAlgo.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

◆ T2CaloEgammaReFastAlgo()

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

Definition at line 26 of file T2CaloEgammaReFastAlgo.cxx.

26 :
27 AthReentrantAlgorithm(name, pSvcLocator)
28{}

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 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ 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 T2CaloEgammaReFastAlgo::execute ( const EventContext & context) const
overridevirtual

Definition at line 42 of file T2CaloEgammaReFastAlgo.cxx.

43{
44
45 auto timer = Monitored::Timer("TIME_exec");
46 auto clET = Monitored::Scalar("TrigEMCluster_eT",-999.0);
47 auto clHET = Monitored::Scalar("TrigEMCluster_had1",-999.0);
48 auto clEta = Monitored::Scalar("TrigEMCluster_eta",-999.0);
49 auto clPhi = Monitored::Scalar("TrigEMCluster_phi",-999.0);
50 auto clReta = Monitored::Scalar("TrigEMCluster_rEta",-999.0);
51 auto clETrings = Monitored::Scalar("TrigEMCluster_et_rings",-999.0);
52 auto res_et = Monitored::Scalar("Resolution_et",-999.0);
53 auto monitoring = Monitored::Group( m_monTool, timer, clET, clHET, clEta, clPhi, clReta,clETrings,res_et);
54
55
56 SG::WriteHandle<xAOD::TrigEMClusterContainer> trigEmClusterCollection(m_clusterContainerKey, context);
57 ATH_CHECK( trigEmClusterCollection.record(std::make_unique<xAOD::TrigEMClusterContainer>(),
58 std::make_unique<xAOD::TrigEMClusterAuxContainer>()) );
59
60
61 auto roisHandle = SG::makeHandle(m_roiCollectionKey, context);
62 if (!roisHandle.isValid()) {
63 ATH_MSG_DEBUG("no RoI");
64 return StatusCode::SUCCESS;
65 }
66
67
68 trigEmClusterCollection->reserve(roisHandle->size());
69
70 ATH_MSG_DEBUG("RoI descriptor size is " << roisHandle->size() );
71
72 for (const TrigRoiDescriptor* roiDescriptor : *roisHandle) {
73 float etaL1, phiL1;
74 double etamin, etamax, phimin, phimax;
75
76 ATH_MSG_DEBUG( "RoI eta = " << roiDescriptor->eta() << " RoI phi = " << roiDescriptor->phi());
77
78 if ((m_l1eta < -9.9) && (m_l1phi < -9.9)) {
79
80 etamin = std::max(-2.5, roiDescriptor->eta() - m_etaWidth);
81 etamax = std::min( 2.5, roiDescriptor->eta() + m_etaWidth);
82
83 phimin = CxxUtils::wrapToPi(roiDescriptor->phi() - m_phiWidth);
84 phimax = CxxUtils::wrapToPi(roiDescriptor->phi() + m_phiWidth);
85
86 etaL1 = roiDescriptor->eta();
87 phiL1 = roiDescriptor->phi();
88 }
89 else {
90 etamin = std::max(-2.5, m_l1eta - m_etaWidth);
91 etamax = std::min( 2.5, static_cast<double>(m_l1eta) + m_etaWidth);
92
94 phimax = CxxUtils::wrapToPi(static_cast<double>(m_l1phi) + m_phiWidth);
95
96 etaL1 = m_l1eta;
97 phiL1 = m_l1phi;
98 }
99
100 TrigRoiDescriptor newroi(roiDescriptor->eta(), etamin, etamax,
101 roiDescriptor->phi(), phimin, phimax);
102
103
104 ATH_MSG_DEBUG(" etamin = " << etamin << " etamax = " << etamax <<
105 " phimin = " << phimin << " phimax = " << phimax);
106
107 xAOD::TrigEMCluster* ptrigEmCluster = new xAOD::TrigEMCluster();
108 trigEmClusterCollection->push_back(ptrigEmCluster);
109 ptrigEmCluster->setEnergy(0.0);
110 ptrigEmCluster->setEt(0.0);
111 ptrigEmCluster->setRawEnergy(0.0);
112 ptrigEmCluster->setRawEt(0.0);
113 ptrigEmCluster->setE277(0);
114 ptrigEmCluster->setEmaxs1(0);
115 ptrigEmCluster->setE2tsts1(0);
116 ptrigEmCluster->setEhad1(-999);
117 ptrigEmCluster->setWeta2(-999);
118 ptrigEmCluster->setFracs1(-999);
119 ptrigEmCluster->setE233(-999);
120 ptrigEmCluster->setE237(-999);
121 ptrigEmCluster->setWstot(-999);
122 ptrigEmCluster->setEta1(-999);
123 ptrigEmCluster->setNCells(0);
124 ptrigEmCluster->setRawEta(-999);
125 ptrigEmCluster->setRawPhi(-999);
126
127 // It is a good idea to clear the energies
128 for (int i = 0; i < CaloSampling::CaloSample::MINIFCAL0; i++) {
129 ptrigEmCluster->setEnergy((CaloSampling::CaloSample)i, 0.);
130 ptrigEmCluster->setRawEnergy((CaloSampling::CaloSample)i, 0.);
131 }
132 // Initial cluster position is the LVL1 position
133 ptrigEmCluster->setEta(etaL1);
134 ptrigEmCluster->setPhi(phiL1);
135
136 // Add RoI word to TrigEMCluster
137 // Dangerous !!!! we need to define a *new* roiDescriptor if we want to
138 // change the size, so we should be careful about *which* roi descriptor
139 // we save, and *which* "roiWord" (if any) we store if we need to use it
140 // again
141 ptrigEmCluster->setRoIword(roiDescriptor->roiWord());
142 const CaloDetDescrElement* caloDDE = nullptr;
143
144 uint32_t error = 0;
145 for (const auto& tool : m_emAlgTools) {
146 ATH_CHECK( tool->execute(*ptrigEmCluster, newroi, caloDDE, context) );
147 }
148 // // support to new monitoring
149 // Cluster quality is a collection of possible errors
150 // No error quality=0
151 ptrigEmCluster->setClusterQuality(error);
152 if ((error & 0xC0000000) || ptrigEmCluster->phi() < -M_PI || ptrigEmCluster->phi() > +M_PI ||
153 std::abs(ptrigEmCluster->eta()) > 10.0) {
154 // Clustering failed. Transmit ahead L1
155 ptrigEmCluster->setEta(etaL1);
156 ptrigEmCluster->setPhi(phiL1);
157 ptrigEmCluster->setEnergy(0.0);
158 ptrigEmCluster->setEt(0.0);
159 }
160
161
162 if ( caloDDE != nullptr ){
163 if ( caloDDE->is_lar_em_barrel() ){
164 for( ToolHandleArray<IEgammaCalibration>::const_iterator
165 ical=m_calibsBarrel.begin();
166 ical != m_calibsBarrel.end(); ++ical )
167 (*ical)->makeCorrection(ptrigEmCluster,caloDDE);
168 }else{
169 for( ToolHandleArray<IEgammaCalibration>::const_iterator
170 ical=m_calibsEndcap.begin();
171 ical != m_calibsEndcap.end(); ++ical )
172 (*ical)->makeCorrection(ptrigEmCluster,caloDDE);
173 }
174 }
175
176 float et_calib = -999.0;
177 float et_uncalib = ptrigEmCluster->et();
178 if ( m_doCalibWithRings ){
179 if (m_calibWRingsTool->checkRings(context)){
180 et_calib = m_calibWRingsTool->makeCalibWRings(context);
181 ptrigEmCluster->setEt(et_calib);
182 }
183 }
184
185
186 float calZ0 = 0;
187
188 // Print out Cluster produced
189 if (msgLvl(MSG::DEBUG)) {
190 ATH_MSG_DEBUG(" Values of Cluster produced: ");
191 ATH_MSG_DEBUG(" REGTEST: emEnergy = " << ptrigEmCluster->energy());
192 ATH_MSG_DEBUG(" REGTEST: hadEnergy = " << ptrigEmCluster->ehad1());
193 ATH_MSG_DEBUG(" REGTEST: e237= " << ptrigEmCluster->e237());
194 ATH_MSG_DEBUG(" REGTEST: e277= " << ptrigEmCluster->e277());
195 ATH_MSG_DEBUG(" REGTEST: clusterWidth = " << ptrigEmCluster->weta2());
196 ATH_MSG_DEBUG(" REGTEST: frac73 = " << ptrigEmCluster->fracs1());
197 ATH_MSG_DEBUG(" REGTEST: e233 = " << ptrigEmCluster->e233());
198 ATH_MSG_DEBUG(" REGTEST: wstot = " << ptrigEmCluster->wstot());
199 ATH_MSG_DEBUG(" REGTEST: eta = " << ptrigEmCluster->eta());
200 ATH_MSG_DEBUG(" REGTEST: phi = " << ptrigEmCluster->phi());
201 ATH_MSG_DEBUG(" REGTEST: Eta1 = " << ptrigEmCluster->eta1());
202 ATH_MSG_DEBUG(" REGTEST: calZ0 = " << calZ0);
203 ATH_MSG_DEBUG(" REGTEST: quality = " << ptrigEmCluster->clusterQuality());
204 ATH_MSG_DEBUG(std::hex << " REGTEST: roiWord = 0x" << ptrigEmCluster->RoIword()
205 << std::dec);
206 }
207 // my monitoring
208 clET = et_uncalib*1e-3;
209 clHET = ptrigEmCluster->ehad1()*1e-3;
210 clEta = ptrigEmCluster->eta();
211 clPhi = ptrigEmCluster->phi();
212 clETrings = et_calib*1e-3;
213
214 if (et_uncalib != 0.) {
215 res_et = (et_calib - et_uncalib)/et_uncalib;
216 }
217
218 if ( ptrigEmCluster->e277() > 0.01 ) clReta = ptrigEmCluster->e237()/ptrigEmCluster->e277();
219
220 } // end of roiCollection iterator
221
222 return StatusCode::SUCCESS;
223}
#define M_PI
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_DEBUG(x)
Athena::TPCnvVers::Current TrigRoiDescriptor
bool msgLvl(const MSG::Level lvl) const
bool is_lar_em_barrel() const
cell belongs to EM barrel
Gaudi::Property< double > m_phiWidth
ToolHandleArray< IReAlgToolCalo > m_emAlgTools
Gaudi::Property< double > m_etaWidth
ToolHandle< TrigFastCalibWithRings > m_calibWRingsTool
SG::ReadHandleKey< TrigRoiDescriptorCollection > m_roiCollectionKey
SG::WriteHandleKey< xAOD::TrigEMClusterContainer > m_clusterContainerKey
ToolHandleArray< IEgammaCalibration > m_calibsBarrel
ToolHandle< GenericMonitoringTool > m_monTool
ToolHandleArray< IEgammaCalibration > m_calibsEndcap
Gaudi::Property< float > m_l1eta
Gaudi::Property< bool > m_doCalibWithRings
Gaudi::Property< float > m_l1phi
unsigned int clusterQuality() const
get quality of cluster built (to be defined)
void setClusterQuality(unsigned int)
set quality of cluster built (to be defined)
void setNCells(int)
set number of cells used from RoI
void setRoIword(long)
set RoI Word
long RoIword() const
get RoI Word
void setRawPhi(float)
set Raw Phi (no calibration)
void setRawEnergy(float)
set Raw Energy (no calibration)
void setRawEta(float)
set Raw Eta (no calibration)
void setRawEt(float)
set Raw Et (no calibration)
void setPhi(float)
set Phi (calibrated)
float wstot() const
get width in first layer
void setEta(float)
set Eta (calibrated)
void setEt(float)
set Et (calibrated)
void setEhad1(float)
set hadronic Energy (first hadronic layer)
float et() const
get Et (calibrated)
float eta() const
get Eta (calibrated)
void setWeta2(float)
set cluster width (based on a 3x5 cluster - 2nd layer)
void setEnergy(float energy)
set Energy (calibrated)
void setE2tsts1(float)
set second maximum energy in sampling 1 (strip layer)
float phi() const
get Phi (calibrated)
float e237() const
get Energy in a 3x7 cluster (no calibration) around hottest cell
float e233() const
get Energy in a 3x3 cluster (no calibration) around hottest cell
void setE237(float)
set Energy in a 3x7 cluster (no calibration) around hottest cell
void setE277(float)
set Energy in a 7x7 cluster (no calibration) around hottest cell
void setEmaxs1(float)
set maximum energy in sampling 1 (strip layer)
float eta1() const
get Eta sampling 1 (strip layer)
void setFracs1(float)
set Energy in a 7 strips (around hottest strip) minus energy in 3 strips divided by energy in 3 strip...
float e277() const
get Energy in a 7x7 cluster (no calibration) around hottest cell
float ehad1() const
get hadronic Energy (first hadronic layer)
float weta2() const
get cluster width (based on a 3x5 cluster - 2nd layer)
void setWstot(float)
set width in first layer
float energy() const
get Energy (calibrated)
float fracs1() const
get Energy in a 7 strips (around hottest strip) minus energy in 3 strips divided by energy in 3 strip...
void setE233(float)
set Energy in a 3x3 cluster (no calibration) around hottest cell
void setEta1(float)
set Eta sampling 1 (strip layer)
T wrapToPi(T phi)
Wrap angle in radians to [-pi, pi].
Definition phihelper.h:24
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
timer(name, disabled=False)
TrigEMCluster_v1 TrigEMCluster
Define the latest version of the trigger EM cluster class.
setEventNumber uint32_t

◆ 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 T2CaloEgammaReFastAlgo::initialize ( )
overridevirtual

Definition at line 30 of file T2CaloEgammaReFastAlgo.cxx.

31{
32 m_emAlgTools.retrieve().ignore();
33 ATH_CHECK(m_clusterContainerKey.initialize());
34 ATH_CHECK(m_roiCollectionKey.initialize());
35 ATH_CHECK( m_bcidAvgKey.initialize() );
36 if (! m_monTool.empty() ) ATH_CHECK( m_monTool.retrieve() );
38
39 return StatusCode::SUCCESS;
40}
SG::ReadHandleKey< CaloBCIDAverage > m_bcidAvgKey

◆ 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_ERROR(x)
#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_bcidAvgKey

SG::ReadHandleKey<CaloBCIDAverage> T2CaloEgammaReFastAlgo::m_bcidAvgKey
private
Initial value:
{
this, "BCIDAvgKey", "CaloBCIDAverage", "" }

Definition at line 82 of file T2CaloEgammaReFastAlgo.h.

82 {
83 this, "BCIDAvgKey", "CaloBCIDAverage", "" };

◆ m_calibsBarrel

ToolHandleArray<IEgammaCalibration> T2CaloEgammaReFastAlgo::m_calibsBarrel
private
Initial value:
{
this, "CalibListBarrel", {}, "list of calib tools for the Barrel clusters"}

Definition at line 70 of file T2CaloEgammaReFastAlgo.h.

70 {
71 this, "CalibListBarrel", {}, "list of calib tools for the Barrel clusters"};

◆ m_calibsEndcap

ToolHandleArray<IEgammaCalibration> T2CaloEgammaReFastAlgo::m_calibsEndcap
private
Initial value:
{
this, "CalibListEndcap", {}, "list of calib tools for the EndCap clusters"}

Definition at line 72 of file T2CaloEgammaReFastAlgo.h.

72 {
73 this, "CalibListEndcap", {}, "list of calib tools for the EndCap clusters"};

◆ m_calibWRingsTool

ToolHandle< TrigFastCalibWithRings > T2CaloEgammaReFastAlgo::m_calibWRingsTool { this, "CalibWRingsTool", "", "FastCalo Calib with Rings tool" }
private

Definition at line 91 of file T2CaloEgammaReFastAlgo.h.

91{ this, "CalibWRingsTool", "", "FastCalo Calib with Rings tool" };

◆ m_clusterContainerKey

SG::WriteHandleKey<xAOD::TrigEMClusterContainer> T2CaloEgammaReFastAlgo::m_clusterContainerKey
private
Initial value:
{
this, "ClustersName", "CaloClusters", "Calo cluster container"}

Definition at line 88 of file T2CaloEgammaReFastAlgo.h.

88 {
89 this, "ClustersName", "CaloClusters", "Calo cluster container"};

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

Gaudi::Property<bool> T2CaloEgammaReFastAlgo::m_doCalibWithRings
private
Initial value:
{this, "DoCalibWithRings", false,
"FastCalo Et Calibration using Rings"}

Definition at line 55 of file T2CaloEgammaReFastAlgo.h.

55 {this, "DoCalibWithRings", false,
56 "FastCalo Et Calibration using Rings"};

◆ m_doForward

Gaudi::Property<bool> T2CaloEgammaReFastAlgo::m_doForward
private
Initial value:
{this, "DoForward", false,
"Do Forward clusters"}

calculate zo mass

Definition at line 52 of file T2CaloEgammaReFastAlgo.h.

52 {this, "DoForward", false,
53 "Do Forward clusters"};

◆ m_emAlgTools

ToolHandleArray<IReAlgToolCalo> T2CaloEgammaReFastAlgo::m_emAlgTools
private
Initial value:
{
this, "IReAlgToolList", {}, "list of ReAlgToolCalos for feature extraction"}

Definition at line 80 of file T2CaloEgammaReFastAlgo.h.

80 {
81 this, "IReAlgToolList", {}, "list of ReAlgToolCalos for feature extraction"};

◆ m_etaWidth

Gaudi::Property<double> T2CaloEgammaReFastAlgo::m_etaWidth {this, "EtaWidth", 0.1, "Eta Width of the Region of Interest"}
private

Definition at line 64 of file T2CaloEgammaReFastAlgo.h.

64{this, "EtaWidth", 0.1, "Eta Width of the Region of Interest"};

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

Gaudi::Property<float> T2CaloEgammaReFastAlgo::m_l1eta {this, "L1ForceEta", -10.0, "Forced LVL1 eta"}
private

Definition at line 62 of file T2CaloEgammaReFastAlgo.h.

62{this, "L1ForceEta", -10.0, "Forced LVL1 eta"};

◆ m_l1phi

Gaudi::Property<float> T2CaloEgammaReFastAlgo::m_l1phi {this, "L1ForcePhi", -10.0, "Forced LVL1 phi"}
private

Definition at line 63 of file T2CaloEgammaReFastAlgo.h.

63{this, "L1ForcePhi", -10.0, "Forced LVL1 phi"};

◆ m_monTool

ToolHandle< GenericMonitoringTool > T2CaloEgammaReFastAlgo::m_monTool { this, "MonTool", "", "Monitoring tool" }
private

Definition at line 90 of file T2CaloEgammaReFastAlgo.h.

90{ this, "MonTool", "", "Monitoring tool" };

◆ m_phiWidth

Gaudi::Property<double> T2CaloEgammaReFastAlgo::m_phiWidth {this, "PhiWidth", 0.1, "Phi Width of the Region of Interest"}
private

Definition at line 65 of file T2CaloEgammaReFastAlgo.h.

65{this, "PhiWidth", 0.1, "Phi Width of the Region of Interest"};

◆ m_rhoEta

Gaudi::Property<std::vector<float> > T2CaloEgammaReFastAlgo::m_rhoEta
private
Initial value:
{
this, "RhoEta", {}, "Variables to calculate Z0 position"}

Definition at line 75 of file T2CaloEgammaReFastAlgo.h.

75 {
76 this, "RhoEta", {}, "Variables to calculate Z0 position"};

◆ m_roiCollectionKey

SG::ReadHandleKey<TrigRoiDescriptorCollection> T2CaloEgammaReFastAlgo::m_roiCollectionKey
private
Initial value:
{
this, "RoIs", "OutputRoIs", "input RoIs"}

Definition at line 85 of file T2CaloEgammaReFastAlgo.h.

85 {
86 this, "RoIs", "OutputRoIs", "input RoIs"};

◆ m_storeCells

Gaudi::Property<bool> T2CaloEgammaReFastAlgo::m_storeCells
private
Initial value:
{this, "StoreCells", false,
"store cells in container attached to RoI"}

Definition at line 67 of file T2CaloEgammaReFastAlgo.h.

67 {this, "StoreCells", false,
68 "store cells in container attached to RoI"};

◆ m_useRings

Gaudi::Property<bool> T2CaloEgammaReFastAlgo::m_useRings
private
Initial value:
{this, "UseRings", false,
"Rings are Used in FastCaloFex"}

Definition at line 58 of file T2CaloEgammaReFastAlgo.h.

58 {this, "UseRings", false,
59 "Rings are Used in FastCaloFex"};

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

◆ m_zEta

Gaudi::Property<std::vector<float> > T2CaloEgammaReFastAlgo::m_zEta
private
Initial value:
{
this, "ZEta", {}, "Variables to calculate Z0 position"}

Definition at line 77 of file T2CaloEgammaReFastAlgo.h.

77 {
78 this, "ZEta", {}, "Variables to calculate Z0 position"};

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