ATLAS Offline Software
GenParticleGenericFilter.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 /*
6  * @author Elmar Ritsch
7  * @date October 2016
8  * @brief A generic particle filter tool for HepMC::GenParticle types
9  */
10 
11 
12 // class header include
14 
15 // HepMC includes
16 #include "AtlasHepMC/GenParticle.h"
17 #include "AtlasHepMC/GenVertex.h"
19 
20 // STL includes
21 #include <limits>
22 #include <algorithm>
23 
26  const std::string& n,
27  const IInterface* p )
28  : base_class(t,n,p),
29  m_minEta(std::numeric_limits<decltype(m_minEta)>::lowest()),
30  m_maxEta(std::numeric_limits<decltype(m_maxEta)>::max()),
31  m_minPhi(-M_PI),
32  m_maxPhi(M_PI),
33  m_minMom(std::numeric_limits<decltype(m_minMom)>::lowest()),
34  m_maxMom(std::numeric_limits<decltype(m_maxMom)>::max()),
35  m_pdgs(),
36  m_maxApplicableRadius(std::numeric_limits<decltype(m_maxApplicableRadius)>::max())
37 {
38  // different cut parameters
39  declareProperty("MinEta",
40  m_minEta,
41  "Minimum Particle Pseudorapidity");
42  declareProperty("MaxEta",
43  m_maxEta,
44  "Maximum Particle Pseudorapidity");
45  declareProperty("MinPhi",
46  m_minPhi,
47  "Minimum Particle Phi");
48  declareProperty("MaxPhi",
49  m_maxPhi,
50  "Maximum Particle Phi");
51  declareProperty("MinMom",
52  m_minMom,
53  "Minimum Particle Momentum");
54  declareProperty("MaxMom",
55  m_maxMom,
56  "Maximum Particle Momentum");
57  declareProperty("ParticlePDG",
58  m_pdgs,
59  "List of accepted particle PDG IDs (any accepted if empty)");
60  declareProperty("MaxApplicableRadius",
62  "Only particles with ProductionVertexRadius<MaxApplicableRadius may get filtered out");
63 }
64 
65 
68 {
69  ATH_MSG_VERBOSE("initialize() ...");
70  ATH_MSG_VERBOSE("initialize() successful");
71  return StatusCode::SUCCESS;
72 }
73 
74 
77 {
78  ATH_MSG_VERBOSE("finalize() ...");
79  ATH_MSG_VERBOSE("finalize() successful");
80  return StatusCode::SUCCESS;
81 }
82 
83 
85 #ifdef HEPMC3
87 {
88  bool pass = true;
89  const HepMC::ConstGenVertexPtr productionVertex = particle?particle->production_vertex():nullptr;
90  if (!productionVertex || productionVertex->position().perp()<=m_maxApplicableRadius) {
91  pass = check_cuts_passed(particle);
92  }
93  const auto momentum = particle->momentum();
94  ATH_MSG_VERBOSE( "GenParticle '" << particle << "' with "
95  << (productionVertex ? "pos: r=" + std::to_string(productionVertex->position().perp()) : "")
96  << ", mom: eta=" << momentum.eta() << " phi=" << momentum.phi()
97  << " did " << (pass ? "" : "NOT ")
98  << "pass the cuts.");
99  return pass;
100 }
101 #else
103 {
104  bool pass = true;
105  HepMC::ConstGenVertexPtr productionVertex = particle.production_vertex();
106  if (!productionVertex || productionVertex->position().perp()<=m_maxApplicableRadius) {
107  pass = check_cuts_passed(particle);
108  }
109  const auto momentum = particle.momentum();
110  ATH_MSG_VERBOSE( "GenParticle '" << particle << "' with "
111  << (productionVertex ? "pos: r=" + std::to_string(productionVertex->position().perp()) : "")
112  << ", mom: eta=" << momentum.eta() << " phi=" << momentum.phi()
113  << " did " << (pass ? "" : "NOT ")
114  << "pass the cuts.");
115  return pass;
116 }
117 #endif
118 
119 
121 #ifdef HEPMC3
123  const auto momentum = particle?particle->momentum():HepMC::FourVector(0,0,0,0);
124  int pdg = particle?particle->pdg_id():0;
125 #else
127  const auto& momentum = particle.momentum();
128  int pdg = particle.pdg_id();
129 #endif
130  double mom = std::sqrt(momentum.x()*momentum.x()+momentum.y()*momentum.y()+momentum.z()*momentum.z());
131  double eta = momentum.eta();
132  double phi = momentum.phi();
133 
134  // check the particle pdg code
135  if( m_pdgs.size() && std::find(std::begin(m_pdgs), std::end(m_pdgs), pdg) == std::end(m_pdgs) ) {
136  return false;
137  }
138 
139  // check the momentum cuts
140  if (mom<m_minMom) {
141  return false;
142  }
143  if (mom>m_maxMom) {
144  return false;
145  }
146 
147  // check the eta cuts
148  if (eta<m_minEta || eta>m_maxEta) {
149  return false;
150  }
151 
152  // check the phi cuts
153  if (phi<m_minPhi || phi>m_maxPhi) {
154  return false;
155  }
156 
157  // all cuts passed
158  return true;
159 }
Trk::ParticleSwitcher::particle
constexpr ParticleHypothesis particle[PARTICLEHYPOTHESES]
the array of masses
Definition: ParticleHypothesis.h:76
GenParticleGenericFilter.h
python.PerfMonSerializer.p
def p
Definition: PerfMonSerializer.py:743
ISF::GenParticleGenericFilter::m_minEta
double m_minEta
the cuts defined by the use
Definition: GenParticleGenericFilter.h:70
max
#define max(a, b)
Definition: cfImp.cxx:41
phi
Scalar phi() const
phi method
Definition: AmgMatrixBasePlugin.h:64
find
std::string find(const std::string &s)
return a remapped string
Definition: hcg.cxx:135
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:79
GenVertex.h
PlotCalibFromCool.begin
begin
Definition: PlotCalibFromCool.py:94
M_PI
#define M_PI
Definition: ActiveFraction.h:11
ISF::GenParticleGenericFilter::m_maxPhi
double m_maxPhi
max phi cut
Definition: GenParticleGenericFilter.h:73
ISF::GenParticleGenericFilter::m_maxMom
double m_maxMom
max momentum cut
Definition: GenParticleGenericFilter.h:75
ISF::GenParticleGenericFilter::m_maxApplicableRadius
double m_maxApplicableRadius
Geometrical region (=cylindrical volume around z-axis) within which this filter is applicable.
Definition: GenParticleGenericFilter.h:79
read_hist_ntuple.t
t
Definition: read_hist_ntuple.py:5
ATH_MSG_VERBOSE
#define ATH_MSG_VERBOSE(x)
Definition: AthMsgStreamMacros.h:28
ISF::GenParticleGenericFilter::m_minMom
double m_minMom
min momentum cut
Definition: GenParticleGenericFilter.h:74
GenParticle.h
mergePhysValFiles.end
end
Definition: DataQuality/DataQualityUtils/scripts/mergePhysValFiles.py:93
ISF::GenParticleGenericFilter::m_minPhi
double m_minPhi
min phi cut
Definition: GenParticleGenericFilter.h:72
ParticleGun_EoverP_Config.mom
mom
Definition: ParticleGun_EoverP_Config.py:63
SimpleVector.h
ISF::GenParticleGenericFilter::pass
bool pass(const HepMC::GenParticle &particle) const
Interface method that returns whether the given particle passes all cuts or not.
Definition: GenParticleGenericFilter.cxx:102
ParticleGun_EoverP_Config.momentum
momentum
Definition: ParticleGun_EoverP_Config.py:63
beamspotman.n
n
Definition: beamspotman.py:731
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
ISF::GenParticleGenericFilter::m_maxEta
double m_maxEta
max pseudorapidity cut
Definition: GenParticleGenericFilter.h:71
HepMC::ConstGenParticlePtr
const GenParticle * ConstGenParticlePtr
Definition: GenParticle.h:38
ActsTrk::to_string
std::string to_string(const DetectorType &type)
Definition: GeometryDefs.h:34
ISF::GenParticleGenericFilter::GenParticleGenericFilter
GenParticleGenericFilter(const std::string &t, const std::string &n, const IInterface *p)
Constructor with framework parameters.
Definition: GenParticleGenericFilter.cxx:25
ISF::GenParticleGenericFilter::m_pdgs
PDGCodes m_pdgs
list of accepted particle PDG IDs (any accepted if empty)
Definition: GenParticleGenericFilter.h:76
ISF::GenParticleGenericFilter::initialize
StatusCode initialize()
Athena algtool's Hooks.
Definition: GenParticleGenericFilter.cxx:67
declareProperty
#define declareProperty(n, p, h)
Definition: BaseFakeBkgTool.cxx:15
HepMC::ConstGenVertexPtr
const HepMC::GenVertex * ConstGenVertexPtr
Definition: GenVertex.h:60
ISF::GenParticleGenericFilter::check_cuts_passed
bool check_cuts_passed(const HepMC::GenParticle &particle) const
Check whether the given particle passes all configure cuts or not.
Definition: GenParticleGenericFilter.cxx:126
ISF::GenParticleGenericFilter::finalize
StatusCode finalize()
Athena algtool's Hooks.
Definition: GenParticleGenericFilter.cxx:76
GenParticle
@ GenParticle
Definition: TruthClasses.h:30