ATLAS Offline Software
Loading...
Searching...
No Matches
MCTruthClassifier.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef MCTRUTHCLASSIFIER_MCTRUTHCLASSIFIER_H
6#define MCTRUTHCLASSIFIER_MCTRUTHCLASSIFIER_H
7/********************************************************************
8NAME: MCTruthClassifier.h
9PACKAGE: atlasoff/PhysicsAnalysis/MCTruthClassifier
10AUTHORS: O. Fedin
11CREATED: Sep 2007
12 ********************************************************************/
13
15#include "AsgTools/AsgTool.h"
18// EDM includes
22// For making PID selections easier
26
27#ifndef XAOD_ANALYSIS
28#include "GaudiKernel/ToolHandle.h"
31#endif
32
33#ifndef GENERATIONBASE
34//EDM includes
38#include "xAODEgamma/Electron.h"
39#include "xAODEgamma/Photon.h"
40#include "xAODMuon/Muon.h"
41#include "xAODJet/Jet.h"
42#endif
43
44#if !defined(XAOD_ANALYSIS) && !defined(GENERATIONBASE)
50#include "AthenaKernel/Units.h"
51#endif
52
53#include <cmath>
54#include <utility>
55class MCTruthClassifier : virtual public IMCTruthClassifier , public asg::AsgTool
56{
58public:
59 // constructor
60 MCTruthClassifier(const std::string& type) : asg::AsgTool(type) {
61#if !defined(XAOD_ANALYSIS) && !defined(GENERATIONBASE)
62 declareProperty("FwdElectronUseG4Sel", m_FwdElectronUseG4Sel = true,
63 "Use Geant4 selection for forward electrons calo clusters");
64 declareProperty("FwdElectronTruthExtrEtaCut",
66 "Cut on the eta of the truth Particles to be extrapolated "
67 "for Fwd electrons");
69 "FwdElectronTruthExtrEtaWindowCut",
71 "Cut on the delta eta of the truth Particles to be extrapolated for "
72 "Fwd electrons and the current FwdElectron");
73 declareProperty("partExtrConePhi", m_partExtrConePhi = 0.4);
74 declareProperty("partExtrConeEta", m_partExtrConeEta = 0.2);
75 declareProperty("phtClasConePhi", m_phtClasConePhi = 0.05);
76 declareProperty("phtClasConeEta", m_phtClasConeEta = 0.025);
77 declareProperty("useCaching", m_useCaching = true);
78 declareProperty("phtdRtoTrCut", m_phtdRtoTrCut = 0.1);
79 declareProperty("fwrdEledRtoTrCut", m_fwrdEledRtoTrCut = 0.15);
80 declareProperty("ROICone", m_ROICone = false);
81 //AV: those below are needed in egammaClusMatch
82 declareProperty("pTChargePartCut", m_pTChargePartCut = 1.0);
83 declareProperty("pTNeutralPartCut", m_pTNeutralPartCut = 0.);
84 declareProperty("inclG4part", m_inclG4part = false);
85#endif
86#ifndef GENERATIONBASE
87 declareProperty("deltaRMatchCut", m_deltaRMatchCut = 0.2);
88 declareProperty("deltaPhiMatchCut", m_deltaPhiMatchCut = 0.2);
89 declareProperty("NumOfSiHitsCut", m_NumOfSiHitsCut = 3);
90 declareProperty("jetPartDRMatch", m_jetPartDRMatch = 0.4);
91#endif
92 }
93 virtual ~MCTruthClassifier() = default ;
94
95 // Gaudi algorithm hooks
96 virtual StatusCode initialize() override {
97 ATH_MSG_INFO(" Initializing MCTruthClassifier");
98#ifndef XAOD_ANALYSIS
99 // Only needed for GenParticle interface
101#endif
103
104#if !defined(XAOD_ANALYSIS) && !defined(GENERATIONBASE)
105 if (!m_caloExtensionTool.empty()) {
106 ATH_CHECK(m_caloExtensionTool.retrieve());
107 } else {
108 m_caloExtensionTool.disable();
109 }
110
112
113 if (!m_truthInConeTool.empty()) {
114 ATH_CHECK(m_truthInConeTool.retrieve());
115 } else {
116 m_truthInConeTool.disable();
117 }
118#endif
119 return StatusCode::SUCCESS;
120 }
121 virtual std::pair<MCTruthPartClassifier::ParticleType, MCTruthPartClassifier::ParticleOrigin>
123
124#ifndef XAOD_ANALYSIS /*These can not run in Analysis Base*/
125 virtual std::pair<MCTruthPartClassifier::ParticleType, MCTruthPartClassifier::ParticleOrigin>
126 particleHepMCTruthClassifier(const HepMcParticleLink& theLink,MCTruthPartClassifier::Info* info = nullptr) const override final;
127
128#endif
129
130#ifndef GENERATIONBASE /*These can not run in Generation only release*/
131 //Main method for Track to Truth association
132 virtual const xAOD::TruthParticle* getGenPart(const xAOD::TrackParticle*, MCTruthPartClassifier::Info* info = nullptr) const override final;
133
134#ifndef XAOD_ANALYSIS
135 //Main method for egamma clusters to Truth Particle association
136 virtual const xAOD::TruthParticle* egammaClusMatch(const xAOD::CaloCluster*,bool,MCTruthPartClassifier::Info* info) const override final;
137#endif
138
139 virtual std::pair<MCTruthPartClassifier::ParticleType, MCTruthPartClassifier::ParticleOrigin>
141
143 particleTruthClassifier(const xAOD::Electron*,MCTruthPartClassifier::Info* info = nullptr) const override final;
144
146 particleTruthClassifier(const xAOD::Photon*,MCTruthPartClassifier::Info* info = nullptr) const override final;
147
149 particleTruthClassifier(const xAOD::Muon*,MCTruthPartClassifier::Info* info = nullptr) const override final;
150
153
155 particleTruthClassifier(const xAOD::Jet*, bool DR, MCTruthPartClassifier::Info* info = nullptr) const override final;
156#endif
157
158private:
159 inline double detEta(double x, double y) const { return std::abs(x - y); }
160 inline double detPhi(double x, double y) const {
161 double det = x - y;
162 if (det > M_PI) det = det - 2. * M_PI;
163 if (det < -M_PI) det = det + 2. * M_PI;
164 return std::abs(det);
165 }
166
167 // Temporary helper methods for detecting loops in the truth record
168 // Method1: Returns true if the parent particle is in the list of
169 // children of its decay vertex. Otherwise, returns the result of
170 // Method3.
171 bool TruthLoopDetectionMethod1(const xAOD::TruthVertex* childOrigVtx, const xAOD::TruthParticle* parent) const;
172 // Method2: Returns true if the parent production vertex is the
173 // child decay vertex and the child production vertex is the parent
174 // decay vertex.
175 bool TruthLoopDetectionMethod2(const xAOD::TruthParticle* child, const xAOD::TruthParticle* parent) const;
176 // Method3: Returns true if the parent and child production vertices
177 // are the same.
178 bool TruthLoopDetectionMethod3(const xAOD::TruthVertex* childOrigVtx, const xAOD::TruthParticle* parent) const;
179
181 const xAOD::TruthParticle*,
182 bool& isPrompt,
183 MCTruthPartClassifier::Info& info) const;
184
186 const xAOD::TruthParticle*,
187 bool& isPrompt,
188 MCTruthPartClassifier::Info& info) const;
189
191 const xAOD::TruthParticle*,
192 int motherPDG,
193 MCTruthPartClassifier::Info& info) const;
194
196 const xAOD::TruthParticle*,
197 bool& isPrompt,
198 MCTruthPartClassifier::Info& info) const;
199
201 const xAOD::TruthParticle*,
202 bool& isPrompt,
203 MCTruthPartClassifier::Info& info) const;
204
205#if !defined(XAOD_ANALYSIS) && !defined(GENERATIONBASE)
206 bool genPartToCalo(const EventContext& ctx,
207 const xAOD::CaloCluster* clus,
208 const xAOD::TruthParticle* thePart,
209 bool isFwrdEle,
210 double& dRmatch,
211 bool& isNarrowCone,
212 const CaloDetDescrManager& caloDDMgr) const;
213
214#endif
215
216#ifndef GENERATIONBASE
217 double fracParticleInJet(const xAOD::TruthParticle*, const xAOD::Jet*, bool DR, bool nparts) const;
218 void findJetConstituents(const xAOD::Jet*, std::set<const xAOD::TruthParticle*>& constituents, bool DR) const;
219#endif
220
221 /* Data members*/
223 m_truthParticleContainerKey{this,"xAODTruthParticleContainerName","TruthParticles","ReadHandleKey for xAOD::TruthParticleContainer"};
224
225#if !defined(XAOD_ANALYSIS) && !defined(GENERATIONBASE)
226 ToolHandle<Trk::IParticleCaloExtensionTool> m_caloExtensionTool{this,"ParticleCaloExtensionTool",""};
228 ToolHandle<xAOD::ITruthParticlesInConeTool>
229 m_truthInConeTool{this,"TruthInConeTool","xAOD::TruthParticlesInConeTool/TruthParticlesInConeTool"};
230
242
246#endif
247
248#ifndef XAOD_ANALYSIS
250 m_truthLinkVecReadHandleKey{this,"xAODTruthLinkVector","xAODTruthLinks", "ReadHandleKey for xAODTruthParticleLinkVector"};
251#endif
252#ifndef GENERATIONBASE
257#endif
258};
259#endif // MCTRUTHCLASSIFIER_MCTRUTHCLASSIFIER_H
#define M_PI
#define ASG_TOOL_CLASS(CLASSNAME, INT1)
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)
Property holding a SG store/key/clid from which a ReadHandle is made.
Definition of CaloDetDescrManager.
ATLAS-specific HepMC functions.
ParticleType
Definition TruthClasses.h:8
ParticleOrigin
Wrapper to avoid constant divisions when using units.
#define y
#define x
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
Principal data class for CaloCell clusters.
This class provides the client interface for accessing the detector description information common to...
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
bool TruthLoopDetectionMethod3(const xAOD::TruthVertex *childOrigVtx, const xAOD::TruthParticle *parent) const
virtual std::pair< MCTruthPartClassifier::ParticleType, MCTruthPartClassifier::ParticleOrigin > particleHepMCTruthClassifier(const HepMcParticleLink &theLink, MCTruthPartClassifier::Info *info=nullptr) const override final
MCTruthPartClassifier::ParticleOrigin defOrigOfElectron(const xAOD::TruthParticleContainer &xTruthParticleContainer, const xAOD::TruthParticle *, bool &isPrompt, MCTruthPartClassifier::Info &info) const
SG::ReadHandleKey< xAODTruthParticleLinkVector > m_truthLinkVecReadHandleKey
bool TruthLoopDetectionMethod1(const xAOD::TruthVertex *childOrigVtx, const xAOD::TruthParticle *parent) const
ToolHandle< xAOD::ITruthParticlesInConeTool > m_truthInConeTool
MCTruthPartClassifier::ParticleOrigin defOrigOfNeutrino(const xAOD::TruthParticleContainer &xTruthParticleContainer, const xAOD::TruthParticle *, bool &isPrompt, MCTruthPartClassifier::Info &info) const
double fracParticleInJet(const xAOD::TruthParticle *, const xAOD::Jet *, bool DR, bool nparts) const
SG::ReadCondHandleKey< CaloDetDescrManager > m_caloMgrKey
bool genPartToCalo(const EventContext &ctx, const xAOD::CaloCluster *clus, const xAOD::TruthParticle *thePart, bool isFwrdEle, double &dRmatch, bool &isNarrowCone, const CaloDetDescrManager &caloDDMgr) const
virtual const xAOD::TruthParticle * egammaClusMatch(const xAOD::CaloCluster *, bool, MCTruthPartClassifier::Info *info) const override final
MCTruthPartClassifier::ParticleOrigin defOrigOfMuon(const xAOD::TruthParticleContainer &xTruthParticleContainer, const xAOD::TruthParticle *, bool &isPrompt, MCTruthPartClassifier::Info &info) const
MCTruthPartClassifier::ParticleOrigin defOrigOfPhoton(const xAOD::TruthParticleContainer &xTruthParticleContainer, const xAOD::TruthParticle *, bool &isPrompt, MCTruthPartClassifier::Info &info) const
virtual ~MCTruthClassifier()=default
float m_FwdElectronTruthExtrEtaWindowCut
double detEta(double x, double y) const
virtual std::pair< MCTruthPartClassifier::ParticleType, MCTruthPartClassifier::ParticleOrigin > particleTruthClassifier(const xAOD::TruthParticle *, MCTruthPartClassifier::Info *info=nullptr) const override final
void findJetConstituents(const xAOD::Jet *, std::set< const xAOD::TruthParticle * > &constituents, bool DR) const
MCTruthClassifier(const std::string &type)
ToolHandle< Trk::IParticleCaloExtensionTool > m_caloExtensionTool
bool TruthLoopDetectionMethod2(const xAOD::TruthParticle *child, const xAOD::TruthParticle *parent) const
double detPhi(double x, double y) const
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_truthParticleContainerKey
virtual const xAOD::TruthParticle * getGenPart(const xAOD::TrackParticle *, MCTruthPartClassifier::Info *info=nullptr) const override final
MCTruthPartClassifier::ParticleOrigin defOrigOfTau(const xAOD::TruthParticleContainer &xTruthParticleContainer, const xAOD::TruthParticle *, int motherPDG, MCTruthPartClassifier::Info &info) const
Property holding a SG store/key/clid from which a ReadHandle is made.
Base class for the dual-use tool implementation classes.
Definition AsgTool.h:47
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
STL class.
Class describing an electron.
Class describing an photon.
NRpcCablingAlg reads raw condition data and writes derived condition data to the condition store.
STL namespace.
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.
Jet_v1 Jet
Definition of the current "jet version".
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition TruthVertex.h:15
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TruthParticle_v1 TruthParticle
Typedef to implementation.
TruthParticleContainer_v1 TruthParticleContainer
Declare the latest version of the truth particle container.