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

TrigEgammaFastElectronReAlgo is a Trigger Fex Algorithm that retrieves the TrigEMCluster container and the TrackCollection containers and then creates a TrigElectron Container with a subset of calorimeter-ID selection variables that are calculated (eg. More...

#include <TrigEgammaFastElectronReAlgo.h>

Inheritance diagram for TrigEgammaFastElectronReAlgo:
Collaboration diagram for TrigEgammaFastElectronReAlgo:

Public Member Functions

 TrigEgammaFastElectronReAlgo (const std::string &name, ISvcLocator *pSvcLocator)
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

bool extrapolate (const EventContext &ctx, const CaloDetDescrManager &caloDD, const xAOD::TrigEMCluster *, const xAOD::TrackParticle *, double &, double &) const
Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Static Private Member Functions

static double getCaloPt (const xAOD::TrigElectron *el)
static double getTkPt (const xAOD::TrigElectron *el)

Private Attributes

Gaudi::Property< bool > m_acceptAll {this, "AcceptAll", false, "Build electrons for all tracks"}
Gaudi::Property< float > m_clusEtthr {this, "ClusEt", 20.0*Gaudi::Units::GeV , " lower limit on cluster Et"}
Gaudi::Property< float > m_trackPtthr {this, "TrackPt", 5.0*Gaudi::Units::GeV , "lower limit on TrackPt cut" }
Gaudi::Property< float > m_calotrkdeta_noextrap {this, "CaloTrackdEtaNoExtrap", 0.5, "Upper limit on DEta between Calo cluster and Track for track preselection before extrapolation"}
Gaudi::Property< float > m_trackPtthr_highet {this, "TrackPtHighEt", 2.0*Gaudi::Units::GeV , "lower limit on TrackPt cut High Et Cluster (20GeV)"}
Gaudi::Property< float > m_calotrkdeta_noextrap_highet {this, "CaloTrackdEtaNoExtrapHighEt", 0, "upper limit on DEta between Calo cluster and Track for track preselection before extrapolation for High Et cluster (20GeV)"}
Gaudi::Property< float > m_calotrackdeta {this, "CaloTrackdETA", 0, "Upper limit on DEta between Calo cluster and Track"}
Gaudi::Property< float > m_calotrackdphi {this, "CaloTrackdPHI", 0, "Upper limit on DPhi between Calo cluster and Track"}
Gaudi::Property< float > m_calotrackdeoverp_low {this, "CaloTrackdEoverPLow", 0, "lower limit on E(calo)/p(track)"}
Gaudi::Property< float > m_calotrackdeoverp_high {this, "CaloTrackdEoverPHigh", 0, "upper limit on track E(calo)/p(track)"}
Gaudi::Property< float > m_RCAL {this, "RCalBarrelFace", 1470.0*Gaudi::Units::mm , "Radius of inner face of the barrel calorimeter"}
Gaudi::Property< float > m_ZCAL {this, "ZCalEndcapFace", 3800.0*Gaudi::Units::mm, "z of the inner face of endcap calorimeter"}
ToolHandle< Trk::IParticleCaloExtensionToolm_caloExtensionTool
SG::ReadHandleKey< TrigRoiDescriptorCollectionm_roiCollectionKey
SG::ReadHandleKey< xAOD::TrigEMClusterContainerm_TrigEMClusterContainerKey
SG::ReadHandleKey< xAOD::TrackParticleContainerm_TrackParticleContainerKey
SG::ReadCondHandleKey< CaloDetDescrManagerm_caloDetDescrMgrKey
SG::WriteHandleKey< xAOD::TrigElectronContainerm_outputElectronsKey
SG::WriteHandleKey< xAOD::TrigElectronContainerm_outputDummyElectronsKey
ToolHandle< GenericMonitoringToolm_monTool { this, "MonTool", "", "Monitoring 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

TrigEgammaFastElectronReAlgo is a Trigger Fex Algorithm that retrieves the TrigEMCluster container and the TrackCollection containers and then creates a TrigElectron Container with a subset of calorimeter-ID selection variables that are calculated (eg.

E/p) The algo will apply some very loose selection cuts to the TrigElectronContainer created which is of TrigParticle type TrigElectron The TrigElectron conatiner will then be retrieved by the hypothesis algorithm TrigEgammaFastElectronHypo that will perform the corresponding fast electron selection

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

◆ TrigEgammaFastElectronReAlgo()

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

Definition at line 28 of file TrigEgammaFastElectronReAlgo.cxx.

29 : AthReentrantAlgorithm(name, pSvcLocator)
30{}

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 TrigEgammaFastElectronReAlgo::execute ( const EventContext & ctx) const
overridevirtual

< monitor preselection between track eta and cluster before extrapolation

Create a TrigElectron corresponding to this candidate assume cluster quantities give better estimate of transverse energy (probably a safe assumption for large pT) and that track parameters at perigee give better estimates of angular quantities

Definition at line 74 of file TrigEgammaFastElectronReAlgo.cxx.

75{
76 using namespace xAOD;
77
78 ATH_MSG_DEBUG( "Executing " <<name());
79
80 auto trigElecColl = SG::makeHandle (m_outputElectronsKey, ctx);
81 ATH_CHECK( trigElecColl.record (std::make_unique<xAOD::TrigElectronContainer>(),
82 std::make_unique<xAOD::TrigEMClusterAuxContainer>()) );
83
84 ATH_MSG_DEBUG( "Made WriteHandle " << m_outputElectronsKey );
85
86 auto trigDummyElecColl = SG::makeHandle (m_outputDummyElectronsKey, ctx);
87
88 ATH_MSG_DEBUG( "Made Dummy WriteHandle " << m_outputDummyElectronsKey );
89 ATH_CHECK( trigDummyElecColl.record (std::make_unique<xAOD::TrigElectronContainer>(),
90 std::make_unique<xAOD::TrigEMClusterAuxContainer>()) );
91
92 auto roiCollection = SG::makeHandle(m_roiCollectionKey, ctx);
93 ATH_MSG_DEBUG( "Made handle " << m_roiCollectionKey );
94
95
96 if (roiCollection->size()==0) {
97 ATH_MSG_DEBUG(" RoI collection size = 0");
98 return StatusCode::SUCCESS;
99 }
100
101 const TrigRoiDescriptor* roiDescriptor = *(roiCollection->begin());
102
103 ATH_MSG_DEBUG(" RoI ID = " << (roiDescriptor)->roiId() << ": Eta = " << (roiDescriptor)->eta()
104 << ", Phi = " << (roiDescriptor)->phi());
105
106 float calo_eta(999), calo_phi(999), calo_et(-1);
107
108 auto clusterCont = SG::makeHandle (m_TrigEMClusterContainerKey, ctx);
109 ATH_MSG_DEBUG( "Made handle " << m_TrigEMClusterContainerKey );
110
111 //JTB Should only be 1 cluster in each RoI
112
113 const xAOD::TrigEMCluster *emCluster=(*clusterCont->begin());
114 ElementLink<xAOD::TrigEMClusterContainer> clusEL=ElementLink<xAOD::TrigEMClusterContainer> (*clusterCont,0);
115
116 // copy relevant quantities to local variables
117 calo_eta = emCluster->eta();
118 calo_phi = emCluster->phi();
119 calo_et = emCluster->et();
120
121
122 // Transverse em energy
123 ATH_MSG_DEBUG("Cluster: ET=" << calo_et);
124 ATH_MSG_DEBUG("searching a matching track: loop over tracks");
125
126 SG::ReadHandle<xAOD::TrackParticleContainer> tracks(m_TrackParticleContainerKey, ctx);
127 ATH_MSG_DEBUG( "Made handle " << m_TrackParticleContainerKey );
128
129
130 if (tracks->size() == 0){
131 return StatusCode::SUCCESS; // Exit early if there are no tracks
132 }
133
134 size_t coll_size = tracks->size();
135 trigElecColl->reserve(coll_size);
136
137 // monitoring
138 std::vector<float> calotrkdeta_noextrap_mon;
139
140 auto caloPtMon = Monitored::Collection("PtCalo", *trigElecColl, getCaloPt );
141 auto trackPtMon = Monitored::Collection("PtTrack", *trigElecColl, getTkPt );
142 auto caloTrackDEtaMon = Monitored::Collection("CaloTrackdEta", *trigElecColl, &xAOD::TrigElectron::trkClusDeta );
143 auto caloTrackDPhiMon = Monitored::Collection("CaloTrackdPhi", *trigElecColl, &xAOD::TrigElectron::trkClusDphi );
144 auto etOverPtMon = Monitored::Collection("CaloTrackEoverP", *trigElecColl, &xAOD::TrigElectron::etOverPt );
145 auto caloTrackDEtaNoExtrapMon = Monitored::Collection("CaloTrackdEtaNoExtrapMon", calotrkdeta_noextrap_mon );
146
147 auto mon = Monitored::Group(m_monTool, caloPtMon, trackPtMon, caloTrackDEtaMon, caloTrackDPhiMon, etOverPtMon, caloTrackDEtaNoExtrapMon );
148
149 // Make Dummy Electron
150 xAOD::TrigElectron* trigDummyElec = new xAOD::TrigElectron();
151
152 ElementLink<xAOD::TrackParticleContainer> trackDummyEL = ElementLink<xAOD::TrackParticleContainer> (*tracks, 0);
153
154 trigDummyElecColl->push_back(trigDummyElec);
155 trigDummyElec->init( 0, 0, 0, 0, clusEL, trackDummyEL);
156
157 //Pick the Calo Det Descr
158 SG::ReadCondHandle<CaloDetDescrManager> caloDetDescrMgrHandle{
160 ATH_CHECK(caloDetDescrMgrHandle.isValid());
161 const CaloDetDescrManager* caloDD = *caloDetDescrMgrHandle;
162
163 // loop over tracks
164 for(unsigned int track_index = 0; track_index < tracks->size(); track_index++)
165 {
166 const xAOD::TrackParticle_v1* trkIter = (*tracks)[track_index];
167 ATH_MSG_DEBUG("Track loop starts");
168 ATH_MSG_VERBOSE("AlgoId = " << (trkIter)->patternRecoInfo());
169 ATH_MSG_VERBOSE("At perigee:");
170 ATH_MSG_DEBUG(" Pt = " << std::abs((trkIter)->pt()));
171 ATH_MSG_DEBUG(" phi = " << std::abs((trkIter)->phi0()));
172 ATH_MSG_DEBUG(" eta = " << std::abs((trkIter)->eta()));
173 ATH_MSG_DEBUG(" z0 = " << std::abs((trkIter)->z0()));
174 ATH_MSG_DEBUG(" d0 = " << std::abs((trkIter)->d0()));
175 // ============================================= //
176 // Pt cut
177 float trkPt = std::abs((trkIter)->pt());
178 float etoverpt = std::abs(calo_et/trkPt);
179 float calotrkdeta_noextrap = (trkIter)->eta() - calo_eta;
180
181 double etaAtCalo=999.;
182 double phiAtCalo=999.;
183 if(m_acceptAll){
184 if(!extrapolate(ctx,*caloDD,emCluster,trkIter,etaAtCalo,phiAtCalo)){
185 ATH_MSG_VERBOSE("extrapolator failed");
186 continue;
187 }
188 else{
189 ATH_MSG_VERBOSE("REGTEST: TrigElectron: cluster index = " << clusEL.index() <<
190 " eta = " << etaAtCalo << " phi = " << phiAtCalo);
191
192 xAOD::TrigElectron* trigElec = new xAOD::TrigElectron();
193 trigElecColl->push_back(trigElec);
194 ElementLink<xAOD::TrackParticleContainer> trackEL = ElementLink<xAOD::TrackParticleContainer> (*tracks, track_index);
195 trigElec->init( roiDescriptor->roiWord(),
196 etaAtCalo, phiAtCalo, etoverpt,
197 clusEL,
198 trackEL);
199 calotrkdeta_noextrap_mon.push_back(calotrkdeta_noextrap);
200 }
201 }else {
202
203 ATH_MSG_DEBUG("Apply cuts on track with index: "<<track_index);
204 if(trkPt < m_trackPtthr){
205 ATH_MSG_DEBUG("Failed track pt cut " << trkPt);
206 continue;
207 }
208 if(std::abs(calotrkdeta_noextrap) > m_calotrkdeta_noextrap){
209 ATH_MSG_DEBUG("Failed pre extrapolation calo track deta " << calotrkdeta_noextrap);
210 continue;
211 }
212 if(calo_et > m_clusEtthr){
213 if(trkPt < m_trackPtthr_highet){
214 ATH_MSG_DEBUG("Failed track pt cut for high et cluster");
215 continue;
216 }
217 if(calotrkdeta_noextrap > m_calotrkdeta_noextrap_highet){
218 ATH_MSG_DEBUG("Failed pre extrapolation calo track deta for high et");
219 continue;
220 }
221 }
222 if (etoverpt < m_calotrackdeoverp_low){
223 ATH_MSG_DEBUG("failed low cut on ET/PT");
224 continue;
225 }
226 if (etoverpt > m_calotrackdeoverp_high){
227 ATH_MSG_DEBUG("failed high cut on ET/PT");
228 continue;
229 }
230 if(!extrapolate(ctx,*caloDD,emCluster,trkIter,etaAtCalo,phiAtCalo)){
231 ATH_MSG_DEBUG("extrapolator failed 1");
232 continue;
233 }
234 // all ok: do track-matching cuts
235 ATH_MSG_DEBUG("extrapolated eta/phi=" << etaAtCalo << "/" << phiAtCalo);
236 // match in eta
237 float dEtaCalo = std::abs(etaAtCalo - calo_eta);
238 ATH_MSG_DEBUG("deta = " << dEtaCalo);
239 if ( dEtaCalo > m_calotrackdeta){
240 ATH_MSG_DEBUG("failed eta match cut " << dEtaCalo);
241 continue;
242 }
243 // match in phi: deal with differences larger than Pi
244 float dPhiCalo = std::abs(phiAtCalo - calo_phi);
245 dPhiCalo = ( dPhiCalo < M_PI ? dPhiCalo : 2*M_PI - dPhiCalo );
246 ATH_MSG_DEBUG("dphi = " << dPhiCalo);
247 if ( dPhiCalo > m_calotrackdphi) {
248 ATH_MSG_DEBUG("failed phi match cut " << dPhiCalo);
249 continue;
250 }
251 // all cuts passed
256 ElementLink<xAOD::TrackParticleContainer> trackEL = ElementLink<xAOD::TrackParticleContainer> (*tracks, track_index);
257 ATH_MSG_DEBUG("pT of trackEL: "<<(*trackEL)->pt());
258
259 ATH_MSG_DEBUG("REGTEST: TrigElectron: cluster index = " << clusEL.index() <<
260 " track index= " << track_index << " eta = " << etaAtCalo << " phi = " << phiAtCalo <<
261 " deta = " << dEtaCalo << "dphi = " << dPhiCalo);
262 xAOD::TrigElectron* trigElec = new xAOD::TrigElectron();
263 trigElecColl->push_back(trigElec);
264 trigElec->init( roiDescriptor->roiWord(),
265 etaAtCalo, phiAtCalo, etoverpt,
266 clusEL,
267 trackEL);
268
269 ATH_MSG_DEBUG(" deta = " << dEtaCalo << " deta = " << trigElec->trkClusDeta()
270 << " dphi = " << dPhiCalo << " dphi = " << trigElec->trkClusDphi()
271 << " caloEta = " << calo_eta << " caloEta = " << trigElec->caloEta()
272 << " caloPhi = " << calo_phi << " calophi = " << trigElec->caloPhi()
273 << " etaAtCalo = " << etaAtCalo << " etaAtCalo = " << trigElec->trkEtaAtCalo()
274 << " phiAtCalo = " << phiAtCalo << " phiAtCalo = " << trigElec->trkPhiAtCalo()
275 );
276 calotrkdeta_noextrap_mon.push_back(calotrkdeta_noextrap);
277
278 }
279
280 }
281
282
283 ATH_MSG_DEBUG("REGTEST: returning an xAOD::TrigElectronContainer with size "<< trigElecColl->size() << ".");
284
285 return StatusCode::SUCCESS;
286}
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
Athena::TPCnvVers::Current TrigRoiDescriptor
Gaudi::Property< float > m_calotrackdeoverp_high
SG::WriteHandleKey< xAOD::TrigElectronContainer > m_outputElectronsKey
static double getCaloPt(const xAOD::TrigElectron *el)
Gaudi::Property< float > m_calotrkdeta_noextrap_highet
SG::ReadCondHandleKey< CaloDetDescrManager > m_caloDetDescrMgrKey
Gaudi::Property< float > m_trackPtthr_highet
SG::ReadHandleKey< xAOD::TrigEMClusterContainer > m_TrigEMClusterContainerKey
bool extrapolate(const EventContext &ctx, const CaloDetDescrManager &caloDD, const xAOD::TrigEMCluster *, const xAOD::TrackParticle *, double &, double &) const
Gaudi::Property< float > m_calotrkdeta_noextrap
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_TrackParticleContainerKey
static double getTkPt(const xAOD::TrigElectron *el)
Gaudi::Property< float > m_calotrackdeoverp_low
SG::ReadHandleKey< TrigRoiDescriptorCollection > m_roiCollectionKey
SG::WriteHandleKey< xAOD::TrigElectronContainer > m_outputDummyElectronsKey
ToolHandle< GenericMonitoringTool > m_monTool
virtual unsigned int roiWord() const override final
float et() const
get Et (calibrated)
float eta() const
get Eta (calibrated)
float phi() const
get Phi (calibrated)
float trkClusDeta() const
The absolute value of the track-calo eta measurement difference.
float trkEtaAtCalo() const
Get the track's pseudorapidity extrapolated to the calorimeter.
float trkPhiAtCalo() const
Get the track's azimuthal angle extrapolated to the calorimeter.
float etOverPt() const
Get for the electron.
float caloEta() const
Pseudorapidity ( ) of the electron in the calorimeter.
void init(uint32_t roi, float trkEtaAtCalo, float trkPhiAtCalo, float etOverPt, const EMClusterLink_t &clLink, const TrackParticleLink_t &tpLink)
Initialisation function, setting most properties of the object.
float caloPhi() const
Azimuthal angle ( ) of the electron in the calorimeter.
float trkClusDphi() const
The absolute value of the track-calo phi measurement difference.
ValuesCollection< T > Collection(std::string name, const T &collection)
Declare a monitored (double-convertible) collection.
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
TrigElectron_v1 TrigElectron
Declare the latest version of the class.
TrigEMCluster_v1 TrigEMCluster
Define the latest version of the trigger EM cluster class.
setTeId setLumiBlock roiId

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

◆ extrapolate()

bool TrigEgammaFastElectronReAlgo::extrapolate ( const EventContext & ctx,
const CaloDetDescrManager & caloDD,
const xAOD::TrigEMCluster * clus,
const xAOD::TrackParticle * trk,
double & etaAtCalo,
double & phiAtCalo ) const
private

Definition at line 288 of file TrigEgammaFastElectronReAlgo.cxx.

291 {
292
293 // use the provided EM cluster to "guide" the extrapolation
294 // 1st figure which layer we want to shoot at
295 // in this case we chose between EMB2 or EME2
296 // given the cluster.
297 std::vector<CaloSampling::CaloSample> samples;
298 if (clus->energy(CaloSampling::CaloSample::EME2) >
299 clus->energy(CaloSampling::CaloSample::EMB2)) {
300 samples.push_back(CaloSampling::CaloSample::EME2);
301 } else {
302 samples.push_back(CaloSampling::CaloSample::EMB2);
303 }
304
305 // create the surfaces we want to reach.
306 // Aka either a cylinder for EME2 or a disc for EMB2
307 std::vector<std::unique_ptr<Trk::Surface>> caloSurfaces =
308 m_caloExtensionTool->caloSurfacesFromLayers(samples, clus->eta(), caloDD);
309
310 // And then try to reach them
311 const auto extension = m_caloExtensionTool->surfaceCaloExtension(
312 ctx, trk->perigeeParameters(), samples, caloSurfaces,
314
315 if (extension.empty()) {
316 ATH_MSG_VERBOSE("extrapolator failed 1");
317 return false;
318 }
319 // We target exactly one EMB2 or EME2 (the vector has 1 entry if not empty)
320 etaAtCalo = extension[0].second->position().eta();
321 phiAtCalo = extension[0].second->position().phi();
322
323 ATH_MSG_VERBOSE("Hit sampling :" << extension.at(0).first << " at eta : "
324 << etaAtCalo << " at phi : " << phiAtCalo);
325
326 return true;
327}
ToolHandle< Trk::IParticleCaloExtensionTool > m_caloExtensionTool
const Trk::Perigee & perigeeParameters() const
Returns the Trk::MeasuredPerigee track parameters.
float energy() const
get Energy (calibrated)

◆ 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

◆ getCaloPt()

double TrigEgammaFastElectronReAlgo::getCaloPt ( const xAOD::TrigElectron * el)
inlinestaticprivate

Definition at line 63 of file TrigEgammaFastElectronReAlgo.h.

63 {
64 if(!el) return -999.0;
65 return el->pt();
66 }

◆ getTkPt()

double TrigEgammaFastElectronReAlgo::getTkPt ( const xAOD::TrigElectron * el)
inlinestaticprivate

Definition at line 68 of file TrigEgammaFastElectronReAlgo.h.

68 {
69 if(!el) return -999.0;
70 return (el->trackParticle() ? el->trackParticle()->pt() : -999);
71 }

◆ initialize()

StatusCode TrigEgammaFastElectronReAlgo::initialize ( )
overridevirtual

Definition at line 33 of file TrigEgammaFastElectronReAlgo.cxx.

34{
35 ATH_MSG_DEBUG("Initialization:");
36
37 if (!m_monTool.empty()) {
38 ATH_MSG_DEBUG("Retrieving monTool");
39 CHECK(m_monTool.retrieve());
40 } else {
41 ATH_MSG_DEBUG("No monTool configured => NO MONITOING");
42 }
43
44 if ( m_caloExtensionTool.retrieve().isFailure() ) {
45 ATH_MSG_ERROR("Unable to locate TrackExtrapolator tool ");
46 return StatusCode::FAILURE;
47 }
48
49 // print out settings
50 ATH_MSG_DEBUG("Initialization completed successfully:");
51 ATH_MSG_DEBUG( "AcceptAll = "<< (m_acceptAll==true ? "True" : "False"));
52 ATH_MSG_DEBUG("TrackPt = " << m_trackPtthr);
53 ATH_MSG_DEBUG("TrackPtHighEt = " << m_trackPtthr_highet);
54 ATH_MSG_DEBUG("ClusEt = " << m_clusEtthr);
55 ATH_MSG_DEBUG("CaloTrackdEtaNoExtrap= " << m_calotrkdeta_noextrap);
56 ATH_MSG_DEBUG("CaloTrackdEtaNoExtrapHighEt= " << m_calotrkdeta_noextrap_highet);
57 ATH_MSG_DEBUG("CaloTrackdETA = " << m_calotrackdeta);
58 ATH_MSG_DEBUG("CaloTrackdPHI = " << m_calotrackdphi);
59 ATH_MSG_DEBUG("CaloTrackdEoverPLow = " << m_calotrackdeoverp_low);
60 ATH_MSG_DEBUG("CaloTrackdEoverPHigh = " << m_calotrackdeoverp_high);
61
62 ATH_CHECK( m_roiCollectionKey.initialize() );
65 ATH_CHECK( m_caloDetDescrMgrKey.initialize() );
66 ATH_CHECK( m_outputElectronsKey.initialize() );
68
69 return StatusCode::SUCCESS;
70}
#define ATH_MSG_ERROR(x)
#define CHECK(...)
Evaluate an expression and check for errors.

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

Gaudi::Property<bool> TrigEgammaFastElectronReAlgo::m_acceptAll {this, "AcceptAll", false, "Build electrons for all tracks"}
private

Definition at line 80 of file TrigEgammaFastElectronReAlgo.h.

80{this, "AcceptAll", false, "Build electrons for all tracks"};

◆ m_caloDetDescrMgrKey

SG::ReadCondHandleKey<CaloDetDescrManager> TrigEgammaFastElectronReAlgo::m_caloDetDescrMgrKey
private
Initial value:
{
this,
"CaloDetDescrManager",
"CaloDetDescrManager",
"SG Key for CaloDetDescrManager in the Condition Store"}

Definition at line 114 of file TrigEgammaFastElectronReAlgo.h.

114 {
115 this,
116 "CaloDetDescrManager",
117 "CaloDetDescrManager",
118 "SG Key for CaloDetDescrManager in the Condition Store"};

◆ m_caloExtensionTool

ToolHandle<Trk::IParticleCaloExtensionTool> TrigEgammaFastElectronReAlgo::m_caloExtensionTool
private
Initial value:
{
this, "ParticleCaloExtensionTool",
"Trk::ParticleCaloExtensionTool/ParticleCaloExtensionTool",
"Tool to extrapolate Track to Calo inner surface"}

Definition at line 94 of file TrigEgammaFastElectronReAlgo.h.

94 {
95 this, "ParticleCaloExtensionTool",
96 "Trk::ParticleCaloExtensionTool/ParticleCaloExtensionTool",
97 "Tool to extrapolate Track to Calo inner surface"};

◆ m_calotrackdeoverp_high

Gaudi::Property<float> TrigEgammaFastElectronReAlgo::m_calotrackdeoverp_high {this, "CaloTrackdEoverPHigh", 0, "upper limit on track E(calo)/p(track)"}
private

Definition at line 89 of file TrigEgammaFastElectronReAlgo.h.

89{this, "CaloTrackdEoverPHigh", 0, "upper limit on track E(calo)/p(track)"};

◆ m_calotrackdeoverp_low

Gaudi::Property<float> TrigEgammaFastElectronReAlgo::m_calotrackdeoverp_low {this, "CaloTrackdEoverPLow", 0, "lower limit on E(calo)/p(track)"}
private

Definition at line 88 of file TrigEgammaFastElectronReAlgo.h.

88{this, "CaloTrackdEoverPLow", 0, "lower limit on E(calo)/p(track)"};

◆ m_calotrackdeta

Gaudi::Property<float> TrigEgammaFastElectronReAlgo::m_calotrackdeta {this, "CaloTrackdETA", 0, "Upper limit on DEta between Calo cluster and Track"}
private

Definition at line 86 of file TrigEgammaFastElectronReAlgo.h.

86{this, "CaloTrackdETA", 0, "Upper limit on DEta between Calo cluster and Track"};

◆ m_calotrackdphi

Gaudi::Property<float> TrigEgammaFastElectronReAlgo::m_calotrackdphi {this, "CaloTrackdPHI", 0, "Upper limit on DPhi between Calo cluster and Track"}
private

Definition at line 87 of file TrigEgammaFastElectronReAlgo.h.

87{this, "CaloTrackdPHI", 0, "Upper limit on DPhi between Calo cluster and Track"};

◆ m_calotrkdeta_noextrap

Gaudi::Property<float> TrigEgammaFastElectronReAlgo::m_calotrkdeta_noextrap {this, "CaloTrackdEtaNoExtrap", 0.5, "Upper limit on DEta between Calo cluster and Track for track preselection before extrapolation"}
private

Definition at line 83 of file TrigEgammaFastElectronReAlgo.h.

83{this, "CaloTrackdEtaNoExtrap", 0.5, "Upper limit on DEta between Calo cluster and Track for track preselection before extrapolation"};

◆ m_calotrkdeta_noextrap_highet

Gaudi::Property<float> TrigEgammaFastElectronReAlgo::m_calotrkdeta_noextrap_highet {this, "CaloTrackdEtaNoExtrapHighEt", 0, "upper limit on DEta between Calo cluster and Track for track preselection before extrapolation for High Et cluster (20GeV)"}
private

Definition at line 85 of file TrigEgammaFastElectronReAlgo.h.

85{this, "CaloTrackdEtaNoExtrapHighEt", 0, "upper limit on DEta between Calo cluster and Track for track preselection before extrapolation for High Et cluster (20GeV)"};

◆ m_clusEtthr

Gaudi::Property<float> TrigEgammaFastElectronReAlgo::m_clusEtthr {this, "ClusEt", 20.0*Gaudi::Units::GeV , " lower limit on cluster Et"}
private

Definition at line 81 of file TrigEgammaFastElectronReAlgo.h.

81{this, "ClusEt", 20.0*Gaudi::Units::GeV , " lower limit on cluster Et"};

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

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

Definition at line 131 of file TrigEgammaFastElectronReAlgo.h.

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

◆ m_outputDummyElectronsKey

SG::WriteHandleKey<xAOD::TrigElectronContainer> TrigEgammaFastElectronReAlgo::m_outputDummyElectronsKey
private
Initial value:
{ this,
"DummyElectronsName",
"Electrons",
"output Dummy Electron container name "}

Definition at line 126 of file TrigEgammaFastElectronReAlgo.h.

126 { this,
127 "DummyElectronsName", // property name
128 "Electrons", // default value of StoreGate key
129 "output Dummy Electron container name "};

◆ m_outputElectronsKey

SG::WriteHandleKey<xAOD::TrigElectronContainer> TrigEgammaFastElectronReAlgo::m_outputElectronsKey
private
Initial value:
{
this,
"ElectronsName",
"Electrons",
"output Electron container name "}

Definition at line 120 of file TrigEgammaFastElectronReAlgo.h.

120 {
121 this,
122 "ElectronsName", // property name
123 "Electrons", // default value of StoreGate key
124 "output Electron container name "};

◆ m_RCAL

Gaudi::Property<float> TrigEgammaFastElectronReAlgo::m_RCAL {this, "RCalBarrelFace", 1470.0*Gaudi::Units::mm , "Radius of inner face of the barrel calorimeter"}
private

Definition at line 90 of file TrigEgammaFastElectronReAlgo.h.

90{this, "RCalBarrelFace", 1470.0*Gaudi::Units::mm , "Radius of inner face of the barrel calorimeter"};

◆ m_roiCollectionKey

SG::ReadHandleKey<TrigRoiDescriptorCollection> TrigEgammaFastElectronReAlgo::m_roiCollectionKey
private
Initial value:
{ this,
"RoIs",
"rois",
"input RoI Collection name"}

Definition at line 99 of file TrigEgammaFastElectronReAlgo.h.

99 { this,
100 "RoIs", // property name
101 "rois", // default value of StoreGate key
102 "input RoI Collection name"};

◆ m_TrackParticleContainerKey

SG::ReadHandleKey<xAOD::TrackParticleContainer> TrigEgammaFastElectronReAlgo::m_TrackParticleContainerKey
private
Initial value:
{ this,
"TrackParticlesName",
"Tracks",
"input TrackParticle container name"}

Definition at line 109 of file TrigEgammaFastElectronReAlgo.h.

109 { this,
110 "TrackParticlesName", // property name
111 "Tracks", // default value of StoreGate key
112 "input TrackParticle container name"};

◆ m_trackPtthr

Gaudi::Property<float> TrigEgammaFastElectronReAlgo::m_trackPtthr {this, "TrackPt", 5.0*Gaudi::Units::GeV , "lower limit on TrackPt cut" }
private

Definition at line 82 of file TrigEgammaFastElectronReAlgo.h.

82{this, "TrackPt", 5.0*Gaudi::Units::GeV , "lower limit on TrackPt cut" };

◆ m_trackPtthr_highet

Gaudi::Property<float> TrigEgammaFastElectronReAlgo::m_trackPtthr_highet {this, "TrackPtHighEt", 2.0*Gaudi::Units::GeV , "lower limit on TrackPt cut High Et Cluster (20GeV)"}
private

Definition at line 84 of file TrigEgammaFastElectronReAlgo.h.

84{this, "TrackPtHighEt", 2.0*Gaudi::Units::GeV , "lower limit on TrackPt cut High Et Cluster (20GeV)"};

◆ m_TrigEMClusterContainerKey

SG::ReadHandleKey<xAOD::TrigEMClusterContainer> TrigEgammaFastElectronReAlgo::m_TrigEMClusterContainerKey
private
Initial value:
{ this,
"TrigEMClusterName",
"clusters",
"input TrigEMCluster Container name"}

Definition at line 104 of file TrigEgammaFastElectronReAlgo.h.

104 { this,
105 "TrigEMClusterName", // property name
106 "clusters", // default value of StoreGate key
107 "input TrigEMCluster Container name"};

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

Gaudi::Property<float> TrigEgammaFastElectronReAlgo::m_ZCAL {this, "ZCalEndcapFace", 3800.0*Gaudi::Units::mm, "z of the inner face of endcap calorimeter"}
private

Definition at line 91 of file TrigEgammaFastElectronReAlgo.h.

91{this, "ZCalEndcapFace", 3800.0*Gaudi::Units::mm, "z of the inner face of endcap calorimeter"};

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