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 30 of file TrigEgammaFastElectronReAlgo.cxx.

31 : AthReentrantAlgorithm(name, pSvcLocator)
32{}

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.

62{
63 return 0;
64}

◆ 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 76 of file TrigEgammaFastElectronReAlgo.cxx.

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

88{
89 // If we didn't find any symlinks to add, just return the collection
90 // from the base class. Otherwise, return the extended collection.
91 if (!m_extendedExtraObjects.empty()) {
93 }
95}
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 290 of file TrigEgammaFastElectronReAlgo.cxx.

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

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

75{
76 return BaseAlg::sysExecute (ctx);
77}

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

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