ATLAS Offline Software
Loading...
Searching...
No Matches
FsrPhotonTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include <utility>
6
7// Local include(s):
12#include "xAODEgamma/Photon.h"
17
18
19namespace FSR {
20
21 FsrPhotonTool::FsrPhotonTool( const std::string& name )
22 : asg::AsgTool( name ),
23 m_fsr_type(FsrCandidate::FsrUnknown),
24 m_isoSelTool("", this),
25 m_energyRescaler("", this),
28 {}
29
32
33
35
36 // Greet the user:
37 ATH_MSG_INFO( "Initialising..." );
38
39 ATH_MSG_INFO( "initialize: high_et_min " << m_high_et_min);
40 ATH_MSG_INFO( "initialize: overlap_el_ph " << m_overlap_el_ph);
41 ATH_MSG_INFO( "initialize: overlap_el_mu " << m_overlap_el_mu);
42 ATH_MSG_INFO( "initialize: far_fsr_drcut " << m_far_fsr_drcut);
43 ATH_MSG_INFO( "initialize: far_fsr_etcut " << m_far_fsr_etcut);
44 ATH_MSG_INFO( "initialize: far_fsr_isoWorkingPoint " << m_far_fsr_isoWorkingPoint);
45 ATH_MSG_INFO( "initialize: drcut " << m_drcut);
46 ATH_MSG_INFO( "initialize: etcut " << m_etcut);
47 ATH_MSG_INFO( "initialize: f1cut " << m_f1cut);
48 ATH_MSG_INFO( "initialize: topo_drcut " << m_topo_drcut);
49 ATH_MSG_INFO( "initialize: topo_f1cut " << m_topo_f1cut);
50 ATH_MSG_INFO( "initialize: egCalibToolName " << m_energyRescalerName);
51
52 // Create IsolationSelectionTool for far fsr selection.
53 // Will use the closeBy corrected isolation variables if DoCloseByCorrection is set to true, which is default.
54 asg::AsgToolConfig config ("CP::IsolationSelectionTool/isolationSelectionTool");
55 RETURN_CHECK("initialize", config.setProperty( "PhotonWP", m_far_fsr_isoWorkingPoint));
56
57 if (m_doCloseByIso) {
58 RETURN_CHECK("initialize", config.setProperty("IsoDecSuffix", "CloseByCorr"));
59 if (msg().level() <= MSG::DEBUG) {
60 RETURN_CHECK( "initialize", config.setProperty( "OutputLevel", MSG::DEBUG) );
61 }
62 }
63 RETURN_CHECK("initialize", config.makePrivateTool (m_isoSelTool));
64 ATH_MSG_INFO("initialize - photon IsolationSelectionTool initialized " << m_isoSelTool->name());
65
66 RETURN_CHECK("initialize", m_isoSelTool.retrieve());
67
68 // get egamma calibration tool - needed to recalibrate electron as a photon
69 if (m_energyRescalerName.size()) {
70 m_energyRescaler = ToolHandle<CP::IEgammaCalibrationAndSmearingTool>(m_energyRescalerName);
71 if (m_energyRescaler.empty()) {
72 ATH_MSG_ERROR("initialize - unable to get IEgammaCalibrationAndSmearingTool with name: " << m_energyRescalerName);
73 return StatusCode::FAILURE;
74 }
75
76 ATH_MSG_INFO("initialize - found IEgammaCalibrationAndSmearingTool with name: " << m_energyRescalerName);
77 }
78 else {
79 ATH_MSG_ERROR("initialize - CANNOT retrieve IEgammaCalibrationAndSmearingTool. Please set the property 'egCalibToolName' to be able to recalibrate fsr photon found as an electron ");
80 return StatusCode::FAILURE;
81 }
82
83 // setup electron and photon object quality tools - checks object quality and Dead HV cell cuts
84
85 // electron OQ
86 asg::AsgToolConfig config1 ("CP::EgammaIsGoodOQSelectionTool/electronIsGoodOQSelectionTool");
89 RETURN_CHECK("initialize", m_elIsGoodOQSelectionTool.retrieve());
90 ATH_MSG_INFO("initialize - electronIsGoodOQSelectionTool initialized " << m_elIsGoodOQSelectionTool->name());
91
92 // photon OQ
93 asg::AsgToolConfig config2 ("CP::EgammaIsGoodOQSelectionTool/photonIsGoodOQSelectionTool");
94 RETURN_CHECK("initialize", config2.setProperty( "Mask", xAOD::EgammaParameters::BADCLUSPHOTON));
96 RETURN_CHECK("initialize", m_phIsGoodOQSelectionTool.retrieve());
97 ATH_MSG_INFO("initialize - electronIsGoodOQSelectionTool initialized " << m_phIsGoodOQSelectionTool->name());
98
99 // Return gracefully:
100 return StatusCode::SUCCESS;
101 }
102
103
105 xAOD::PhotonContainer* photons,
106 const xAOD::ElectronContainer* electrons) {
107
108 if (photons == 0 ||electrons == 0) {
109 ATH_MSG_ERROR( "getFsrPhoton: You MUST provide photon AND electron containers" );
110 candidate = FsrCandidate();
112 }
113
114
115 std::vector<FsrCandidate>* cands = getFsrCandidateList(part, photons, electrons);
116
117 ATH_MSG_DEBUG( "Fsr candidate size = " << cands->size() );
118
119 // Return gracefully:
120 if (cands->size() > 0) {
121 candidate = cands->at(0);
123 }
124 candidate = FsrCandidate();
125 ATH_MSG_DEBUG( "Fsr candidate f1 = " << candidate.f1
126 << ", deltaR = " << candidate.deltaR
127 << ", Et = " << candidate.Et
128 << ", Eta = " << candidate.eta
129 << ", Phi = " << candidate.phi
130 << ", type = " <<candidate.type);
132
133 }
134
135
136 std::vector<FsrCandidate>*
138 xAOD::PhotonContainer* photons,
139 const xAOD::ElectronContainer* electrons) {
140 xAOD::PhotonContainer* photons_cont = NULL;
141 const xAOD::ElectronContainer* electrons_cont = NULL;
142
143 // check for photon container
144 if (0 == photons) {
145 ATH_MSG_ERROR( "getFsrCandidateList: No Photon container provided" );
146 return NULL;
147 }
148 else {
149 // set the containers
150 photons_cont = photons;
151 }
152
153 // check for electron container
154 if (0 == electrons) {
155 ATH_MSG_ERROR( "getFsrCandidateList: No Electron container provided" );
156 return NULL;
157 }
158 else {
159 // set the containers
160 electrons_cont = electrons;
161 }
162
166
169
170 const xAOD::Electron* electron = dynamic_cast<const xAOD::Electron*>(part);
171 if (electron) {
172 return getFarFsrCandidateList(part, photons_cont);
173 }
174
175 const xAOD::Muon* muon = dynamic_cast<const xAOD::Muon*>(part);
176 if (muon) {
177 std::vector<FsrCandidate>* nearlist = getNearFsrCandidateList(muon, photons_cont, electrons_cont);
178 return (nearlist->size() > 0) ? nearlist : getFarFsrCandidateList(part, photons_cont);
179 }
180
181 ATH_MSG_WARNING( "The FSR search not done !" );
182 return NULL;
183 }
184
185 std::vector<FsrCandidate>* FsrPhotonTool::getFarFsrCandidateList(const xAOD::IParticle* part,
186 xAOD::PhotonContainer* photons_cont) {
187
188
189 static const SG::AuxElement::Accessor<char> DFCommonPhotonsIsEMTight ("DFCommonPhotonsIsEMTight");
190 static const SG::AuxElement::Accessor<float> topoetcone20 ("topoetcone20");
191
194 std::vector< std::pair <const xAOD::IParticle*, double> > farFsrCandList;
195 farFsrCandList.clear();
196 farFsrCandList.reserve(photons_cont->size());
197
198 // Save muon or electron of overlap removal for isolation
200 parts.push_back(part);
201
202 ATH_MSG_DEBUG( "In getFarFsrCandidateList function : photon size = " << photons_cont->size());
203
204 for (auto ph : *photons_cont) {
205
206 bool oqIsOK = (bool) m_phIsGoodOQSelectionTool->accept(ph);
207
208 bool is_tight_photon = DFCommonPhotonsIsEMTight(*ph);
209 if ( oqIsOK && (ph->p4().Et() > m_far_fsr_etcut) && is_tight_photon) {
210 // correct isolation leakage
211
212 ATH_MSG_VERBOSE( "Far Fsr ph bef : pt " << ph->pt() << " topoetcone20 = " << topoetcone20(*ph));
213
214 // apply isolation selection
215 bool farPhIsoOK = (bool) m_isoSelTool->accept(*ph);
216
217 ATH_MSG_VERBOSE( "Far Fsr ph aft: pt " << ph->pt() << " topoetcone20 = " << topoetcone20(*ph));
218
219 bool far_fsr_drcut_isOK = false;
220 if (farPhIsoOK) {
221 double dr = deltaR(part->eta(), part->phi(), ph->eta(), ph->phi());
222 far_fsr_drcut_isOK = (dr > m_far_fsr_drcut);
223
224 if (far_fsr_drcut_isOK && dr < 0.2) {
225 ATH_MSG_VERBOSE( "Far Fsr candidate kinematics : author " << ph->author()
226 << " Et = " << ph->p4().Et()
227 << " dr = " << dr
228 << " isoIsOK = " << farPhIsoOK );
229 }
230
231 if(far_fsr_drcut_isOK) farFsrCandList.push_back(std::make_pair(ph, dr));
232 }
233
234 ATH_MSG_DEBUG( "Far Fsr candidate kinematics : author " << ph->author()
235 << " Et = " << ph->p4().Et()
236 << " tight = " << is_tight_photon
237 << " farPhIsoOK = " << farPhIsoOK
238 << " far_fsr_drcut_isOK = " << far_fsr_drcut_isOK);
239 }
240 }
241
243 return sortFsrCandidates(farFsrCandList);
244 }
245
246
247 std::vector<FsrCandidate>* FsrPhotonTool::getNearFsrCandidateList(const xAOD::Muon* muon,
248 const xAOD::PhotonContainer* photons_cont,
249 const xAOD::ElectronContainer* electrons_cont) {
254 std::vector< std::pair <const xAOD::IParticle*, double> > nearFsrCandList;
255
256 nearFsrCandList.reserve(photons_cont->size()+electrons_cont->size());
257 ATH_MSG_DEBUG( "In getNearFsrCandidateList function : photon size = " << photons_cont->size()
258 << ", electron size = " << electrons_cont->size());
259
260 for (auto photon : *photons_cont) {
261 float photon_f1;
262 photon->showerShapeValue(photon_f1, xAOD::EgammaParameters::f1);
263
264 bool oqIsOK = (bool) m_phIsGoodOQSelectionTool->accept(photon);
265
266 // Selection is tighter for photons below high_et_min
267 bool high_et_photon = (photon->p4().Et() > m_high_et_min);
268 if( ( oqIsOK && high_et_photon && (photon->p4().Et() > m_etcut) && (photon_f1 > m_f1cut) )
269 ||
270 ( oqIsOK && !high_et_photon && (photon_f1 > m_topo_f1cut)
271 && (m_etcut < photon->p4().Et())) ) {
272
273 double dr = deltaR(muon->eta(), muon->phi(), photon->eta(), photon->phi());
274
275 // ph_cl_eta/phi should be used in duplicate
276 if ( (!high_et_photon && dr < m_topo_drcut)
277 || (high_et_photon && dr < m_drcut) ) {
278 nearFsrCandList.push_back(std::make_pair(photon, dr));
279 ATH_MSG_DEBUG( "Near Fsr candidates ( photon ) kinematics ; author "
280 << photon->author()
281 << " Et = " << photon->p4().Et()
282 << " f1 = " << photon_f1
283 << " dr = " << dr );
284 }
285 }
286 }
287
288 unsigned int nofPhFsr = nearFsrCandList.size();
289
290 // look for fsr photons in electron container
291 for (auto electron : *electrons_cont) {
292 // remove overlap with fsr already found
293 if ( (nearFsrCandList.size() > 0) && isOverlap(electron, nearFsrCandList, nofPhFsr) ) continue;
294
295 // Check object quality
296 // bool oqIsOK = electron->isGoodOQ(xAOD::EgammaParameters::BADCLUSELECTRON);
297bool oqIsOK = (bool) m_elIsGoodOQSelectionTool->accept(electron);
298
299 const xAOD::TrackParticle* electron_track = electron->trackParticle();
300 const xAOD::TrackParticle* muon_track = muon->primaryTrackParticle();
301
302 bool elmutrackmatch = ( (fabs(electron_track->theta()- muon_track->theta()) < m_overlap_el_mu) &&
303 (deltaPhi(electron_track->phi(), muon_track->phi()) < m_overlap_el_mu) );
304 float electron_f1;
305 electron->showerShapeValue(electron_f1, xAOD::EgammaParameters::f1);
306
307 // Get the corrected cluster energy:
308 // if the electron energy has been calibrated, we apply the change applied to the
309 // electron four vector assuming that the uncorrected energy is
310 // Eclus/cosh(track_eta). This will not work if an e-p combination has also been
311 // applied, and so the following will have to change. RDS 04/2015.
312 float eCorr = electron->p4().Et()/(electron->caloCluster()->e()/cosh(electron->trackParticle()->eta()));
313 float clEt = eCorr*electron->caloCluster()->et();
314
315 ATH_MSG_VERBOSE( "Near Fsr candidate ( electron ) Et = " << clEt << " eCorr " << eCorr << " OQ is ok: " << (int) oqIsOK << ", pt,eta,phi " << electron->pt()/1000 << ", " << electron->eta() << ", " << electron->phi());
316
317 // Allow topo-clusters to come in as electrons - apply f1 sliding window cut for Et
318 // > 3.5 GeV and f1 topo cut for Et 1-3.5 GeV
319 if( elmutrackmatch && oqIsOK &&
322
323 double dr = deltaR(muon->eta(), muon->phi(), electron->caloCluster()->eta(), electron->caloCluster()->phi());
324
325 if ( dr < m_drcut ) {
326 nearFsrCandList.push_back(std::make_pair(electron, dr));
327 ATH_MSG_DEBUG( "Near Fsr candidates ( electron ) kinematics : author "
328 << electron->author()
329 << " Et = " << clEt
330 << " f1 = " << electron_f1
331 << " dr = " << dr );
332 }
333 else ATH_MSG_VERBOSE( "FAILED Near Fsr candidates ( electron ) kinematics : author "
334 << electron->author()
335 << " Et = " << clEt
336 << " f1 = " << electron_f1
337 << " dr = " << dr );
338 }
339 else {
340 double dr = deltaR(muon->eta(), muon->phi(), electron->caloCluster()->eta(), electron->caloCluster()->phi());
341
342 if (elmutrackmatch && dr < m_drcut ) {
343 ATH_MSG_DEBUG( "FAILED Near Fsr candidates ( electron ) kinematics : author "
344 << electron->author()
345 << " Et = " << clEt
346 << " f1 = " << electron_f1
347 << " dr = " << deltaR(muon->eta(), muon->phi(), electron->caloCluster()->eta(), electron->caloCluster()->phi())
348 << " theta/phi el/mu " << electron_track->theta() << "/" << muon->p4().Theta()
349 << "/" << electron_track->phi() << "/" << muon->phi()
350 << " theta/phi mu trk " << muon_track->theta() << "/" << muon_track->phi()
351 );
352 }
353 }
354 }
355
356 return sortFsrCandidates(nearFsrCandList);
357
358 } // end of getFsrPhotons
359
360
361 std::vector<FsrCandidate>* FsrPhotonTool::sortFsrCandidates(
362 const std::vector< std::pair <const xAOD::IParticle*, double> >& FsrCandList,
363 const std::string& option) {
364
365 m_fsrPhotons.clear();
366 m_fsrPhotons.reserve(FsrCandList.size());
367
368 for (unsigned int i=0; i < FsrCandList.size(); i++ ) {
369 FsrCandidate c;
370 const xAOD::IParticle* particle = FsrCandList.at(i).first;
371 const xAOD::Electron* electron = dynamic_cast<const xAOD::Electron*>(particle);
372 const xAOD::Photon* photon = dynamic_cast<const xAOD::Photon*>(particle);
373 float part_f1(-999);
374 if ( electron ){
375 c.container = "electron" ;
376
377 // If we have an energy rescaler, calibrate correctly the electron as a photon
378 if (!m_energyRescaler.empty()) {
379 xAOD::Photon photon;
380 photon.Egamma_v1::operator=(*electron);
381 if (m_energyRescaler->applyCorrection(photon) != CP::CorrectionCode::Ok) {
382 ATH_MSG_ERROR("FsrPhotonTool::sortFsrCandidates: Unable to applyCorrection to photon ");
383 }
384 if ( photon.showerShapeValue(part_f1, xAOD::EgammaParameters::f1) ) c.f1 = part_f1;
385 c.Et = photon.pt();
386 c.eta = photon.caloCluster()->eta();
387 c.phi = photon.caloCluster()->phi();
388
389 ATH_MSG_DEBUG("FsrPhotonTool::sortFsrCandidates: el/ph et " << electron->pt() << "/"
390 << electron->caloCluster()->eta() << "/" << electron->caloCluster()->phi() << " "
391 << photon.pt() << "/"
392 << photon.caloCluster()->eta() << "/" << photon.caloCluster()->phi());
393 }
394 else {
395 // Use the cluster pt/eta/phi when considering the electron to be an FSR photon
396 // setup accessor for electron cluster corrected energy
397
398 // Get the corrected cluster energy:
399 // if the electron energy has been calibrated, we apply the change applied to the
400 // electron four vector assuming that the uncorrected energy is
401 // Eclus/cosh(track_eta). This will not work if an e-p combination has also been
402 // applied, and so the following will have to change. RDS 04/2015.
403 float eCorr = electron->p4().Et()/(electron->caloCluster()->e()/cosh(electron->trackParticle()->eta()));
404 float clEt = eCorr*electron->caloCluster()->et();
405 if ( electron->showerShapeValue(part_f1, xAOD::EgammaParameters::f1) ) c.f1 = part_f1;
406 c.Et = clEt;
407 c.eta = electron->caloCluster()->eta();
408 c.phi = electron->caloCluster()->phi();
409 }
410
411 }
412 else if ( photon ) {
413 c.container = "photon";
414 if ( photon->showerShapeValue(part_f1, xAOD::EgammaParameters::f1) ) c.f1 = part_f1;
415 c.Et = photon->p4().Et();
416 c.eta = photon->eta();
417 c.phi = photon->phi();
418 }
419 else {
420 ATH_MSG_WARNING( "sortFsrCandidates: undefined particle - NOT electron nor photon. Should never get here!" );
421 c.container = "";
422 c.Et = 0;
423 c.eta = 0;
424 c.phi = 0;
425 }
426 c.particle = particle;
427 c.deltaR = FsrCandList.at(i).second;
428 if(c.deltaR < 0.05) c.Et -= 400./cosh(particle->eta());
430 m_fsrPhotons.push_back(c);
431
432 ATH_MSG_DEBUG( "sortFsrCandidates: save fsr candidate f1 = " << c.f1
433 << ", deltaR = " << c.deltaR
434 << ", Et = " << c.Et
435 << ", Eta = " << c.eta
436 << ", Phi = " << c.phi
437 << ", type = " <<c.type);
438
439 }
440
441 ATH_MSG_DEBUG( "sortFsrCandidates: found " << m_fsrPhotons.size() << " FSR photons" );
442
443 // sort FSR candidates w.r.t min dR or ET
444 if (option=="ET") std::sort(m_fsrPhotons.begin(), m_fsrPhotons.end(), compareEt);
445 else std::sort(m_fsrPhotons.begin(), m_fsrPhotons.end());
446 return &m_fsrPhotons;
447 }
448
450 const std::vector< std::pair <const xAOD::IParticle*, double> >& phfsr,
451 unsigned int nofPhFsr) {
452 for (unsigned int indx=0; indx < nofPhFsr; indx++ ) {
453 const xAOD::Photon* ph = dynamic_cast<const xAOD::Photon*>(phfsr.at(indx).first);
454 const xAOD::CaloCluster* ph_cl = ph->caloCluster();
455 const xAOD::CaloCluster* el_cl = electron->caloCluster();
456 double dr = deltaR(el_cl->eta(), el_cl->phi(),
457 ph_cl->eta(), ph_cl->phi());
458 ATH_MSG_VERBOSE( "isOverlap dr, ets " << (dr < m_overlap_el_ph) << "/" << dr << "/"
459 << ph_cl->et() << "/" << el_cl->et());
460 if ( dr < m_overlap_el_ph ) return true;
461 }
462 return false;
463 }
464
465 double FsrPhotonTool::deltaPhi(float phi1, float phi2) const {
466 double dphi= fabs(phi1 - phi2);
467 if (dphi > M_PI) dphi=2*M_PI-dphi;
468 return dphi;
469 }
470
471
472 double FsrPhotonTool::deltaR(float meta, float mphi, float peta, float pphi) const {
473
474 double deta = fabs(meta - peta);
475 double dphi = deltaPhi(mphi, pphi);
476 double dR = sqrt((dphi*dphi)+(deta*deta));
477
478 return dR;
479 }
480
481} // namespace FSR
#define M_PI
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
DataVector adapter that acts like it holds const pointers.
#define RETURN_CHECK(CONTEXT, EXP)
Helper macro for checking return codes in a compact form in the code.
Definition ReturnCheck.h:26
Return value from object correction CP tools.
@ Error
Some error happened during the object correction.
@ Ok
The correction was done successfully.
DataVector adapter that acts like it holds const pointers.
size_type size() const noexcept
Returns the number of elements in the collection.
Simple interface for searching the FSR candidate.
virtual CP::CorrectionCode getFsrPhoton(const xAOD::IParticle *part, FsrCandidate &candidate, xAOD::PhotonContainer *photons, const xAOD::ElectronContainer *electrons)
Get the "FSR candidate" as a return value for a muon (collinar and far FSR)
Gaudi::Property< bool > m_doCloseByIso
Gaudi::Property< double > m_drcut
Gaudi::Property< std::string > m_energyRescalerName
Gaudi::Property< double > m_far_fsr_etcut
ToolHandle< IAsgSelectionTool > m_phIsGoodOQSelectionTool
Gaudi::Property< double > m_overlap_el_ph
static bool compareEt(const FsrCandidate &c1, const FsrCandidate &c2)
Gaudi::Property< double > m_topo_drcut
Gaudi::Property< double > m_far_fsr_drcut
Gaudi::Property< std::string > m_far_fsr_isoWorkingPoint
FsrCandidate::FsrType m_fsr_type
Gaudi::Property< double > m_etcut
FsrPhotonTool(const std::string &name)
Create a proper constructor for Athena.
std::vector< FsrCandidate > * sortFsrCandidates(const std::vector< std::pair< const xAOD::IParticle *, double > > &FsrCandList, const std::string &option="ET")
Need for the FSR search.
Gaudi::Property< double > m_f1cut
double deltaR(float muonEta, float muonPhi, float phEta, float phPhi) const
virtual StatusCode initialize()
Function initialising the tool.
virtual std::vector< FsrCandidate > * getNearFsrCandidateList(const xAOD::Muon *part, const xAOD::PhotonContainer *photons_cont, const xAOD::ElectronContainer *electrons_cont)
Find and Return ALL NEAR FSR candidates.
bool isOverlap(const xAOD::Electron_v1 *electron, const std::vector< std::pair< const xAOD::IParticle *, double > > &phfsr, unsigned int nofPhFsr)
double deltaPhi(float phi1, float phi2) const
virtual std::vector< FsrCandidate > * getFsrCandidateList(const xAOD::IParticle *part, xAOD::PhotonContainer *photons, const xAOD::ElectronContainer *electrons)
Find ALL FSR candidates for a given particle (electron or muon) providing newly calibrated photon and...
std::vector< FsrCandidate > m_fsrPhotons
ToolHandle< CP::IEgammaCalibrationAndSmearingTool > m_energyRescaler
ToolHandle< IAsgSelectionTool > m_elIsGoodOQSelectionTool
virtual std::vector< FsrCandidate > * getFarFsrCandidateList(const xAOD::IParticle *part, xAOD::PhotonContainer *photons_cont)
Find and Return ALL FAR FSR candidates.
Gaudi::Property< double > m_topo_f1cut
Gaudi::Property< double > m_overlap_el_mu
Gaudi::Property< double > m_high_et_min
ToolHandle< CP::IIsolationSelectionTool > m_isoSelTool
~FsrPhotonTool()
Create a destructor.
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
StatusCode setProperty(const std::string &name, const T &value)
set the given property
an object that can create a AsgTool
::StatusCode makePrivateTool(ToolHandle< T > &toolHandle) const
make a private tool with the given configuration
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
const xAOD::CaloCluster * caloCluster(size_t index=0) const
Pointer to the xAOD::CaloCluster/s that define the electron candidate.
Class providing the definition of the 4-vector interface.
float theta() const
Returns the parameter, which has range 0 to .
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .)
@ VIEW_ELEMENTS
this data object is a view, it does not own its elmts
-diff
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
const uint32_t BADCLUSELECTRON
Definition EgammaDefs.h:116
const uint32_t BADCLUSPHOTON
Definition EgammaDefs.h:124
@ f1
E1/E = fraction of energy reconstructed in the first sampling, where E1 is energy in all strips belon...
Definition EgammaEnums.h:53
PhotonContainer_v1 PhotonContainer
Definition of the current "photon container version".
ElectronContainer_v1 ElectronContainer
Definition of the current "electron container version".
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Muon_v1 Muon
Reference the current persistent version:
Photon_v1 Photon
Definition of the current "egamma version".
Electron_v1 Electron
Definition of the current "egamma version".
MsgStream & msg
Definition testRead.cxx:32