ATLAS Offline Software
JetTruthParticleSelectorTool.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
4 
7 
8 #include <string>
9 
10 #include <map>
11 
12 #include <functional>
13 
14 
15 #include "xAODTruth/TruthVertex.h"
16 
18 
21 
22 namespace {
23 
24  bool isWZDecay(const xAOD::TruthParticle* p) {
25  int pdg_id = std::abs(p->pdgId() );
26  int mom_pdg_id = pdg_id;
27  const xAOD::TruthVertex* vprod = p->prodVtx();
28  const xAOD::TruthVertex*oldVprod = vprod;
29  // Ascend decay chain looking for when actual decay occurs (not jsut evolution of particle)
30  while (pdg_id == mom_pdg_id) {
31  const xAOD::TruthParticle* mother = vprod->incomingParticle(0) ;
32  if (mother) {
33  mom_pdg_id = abs(mother->pdgId());
34  } else break;
35  if (pdg_id != mom_pdg_id) break;
36 
37  // Protect against strange infinite reference in sherpa samples
38  if(!mother->hasProdVtx()) break;
39  if(oldVprod==mother->prodVtx()) break;
40  oldVprod = vprod;
41  vprod = mother->prodVtx();
42  if (!vprod || vprod->nIncomingParticles() == 0) break;
43  } // End while loop
44 
45  // The following vertex-based identification of W/Z's is needed for SHERPA
46  // samples where the W/Z particle is not explicitly in the particle record.
47  // At this point if we have a valid vertex, it should be a true decay vertex.
48  // If it is a W or Z then two of those decay products should be lepton/neutrino
49  int nDecay=0;
50  // Prompt W/Z's should come from a vertex with more than one incoming particle
51  // This suppresses fave Z's from conversions
52  if (vprod && vprod->nIncomingParticles() >1)
53  {
54  for( const auto& elTruth : vprod->outgoingParticleLinks() ){
55  if (MC::isSMLepton(*elTruth)) nDecay++;
56  }
57  }
58 
59  bool isWZ = (nDecay==2);
60  return (mom_pdg_id == 23 || mom_pdg_id == 24 || isWZ);
61  }
62 
63 
64  // This part is copied and adapted from PhysicsAnalysis/D3PDMaker/TruthD3PDMaker/trunk/src/TruthJetFilterTool.cxx
65  bool isWZDecay(const xAOD::TruthParticle* p, std::vector<const xAOD::TruthParticle*>& wzLeptons, double photonCone2) {
66  int pdg_id = abs(p->pdgId() );
67 
68  if(! p->hasProdVtx() ) return false;
69  const xAOD::TruthVertex* vprod = p->prodVtx();
70 
71  if (vprod && vprod->nIncomingParticles() > 0) {
72  if (isWZDecay(p)) { // paricle decends from W or Z
73  // Save lepton reference for comparison to FSR photons (only for muons and electrons)
74  if(MC::isStable(p) && ( (pdg_id==11) || (pdg_id==13) ) ) {
75  wzLeptons.push_back(p);
76  }
77 
78  // Only exclude photons within deltaR of leptons (if m_photonCone<0, exclude all photons)
79  if( pdg_id == 22 && photonCone2>0)
80  {
81  for( const auto *lep: wzLeptons) {
82  double deltaR2 = xAOD::P4Helpers::deltaR2(*p, *lep, false);
83  if( deltaR2 < photonCone2 ) {
84  return true;
85  }
86  }
87  }
88  else // not a photon so exclude
89  return true;
90  }
91  }
92  return false;
93  }
94 
95 
96  // This part is copied and adapted from PhysicsAnalysis/D3PDMaker/TruthD3PDMaker/trunk/src/TruthJetFilterTool.cxx
97  bool isLeptonFromTau(const xAOD::TruthParticle* part){
98 
99  int pdg = part->pdgId();
100 
101  if(!MC::isSMLepton(pdg)) return false; // all leptons including tau.
102  if(!part->hasProdVtx()) return false;
103  const xAOD::TruthVertex* prod = part->prodVtx();
104  if(!prod) return false; // no parent.
105 
106  for( const auto& elParent : prod->incomingParticleLinks() ){
107  int parentId = (*elParent)->pdgId();
108  if(abs(parentId) == 15 && isWZDecay(*elParent)) {
109  return true; // Has tau parent
110  }
111 
112  if(parentId == pdg) { // Same particle just a different MC status
113  // Go up the generator record until a tau is found or not.
114  // Note that this requires a connected *lepton* chain, while calling
115  // isFromTau would allow leptons from hadrons from taus
116  if(isLeptonFromTau(*elParent)) {
117  return true;
118  }
119  }
120  }
121  return false;
122  }
123 
124 
125 }
126 
127 
129  asg::AsgTool(n),
130  m_min_pt(-1),
131  m_max_pt(-1),
132  m_maxAbsEta(-1),
133  m_minEta(-5),
134  m_maxEta( 5),
135  m_wzPhotonCone(0.1),
136  m_useOnlyCharged(false),
137  m_listPDGofStables(false),
138  m_selectionModeName("StableNoMuonNoNu"),
139  m_selectionMode()
140 {
141 
142  declareInterface<JetTruthParticleSelectorTool>(this);
143  declareProperty("min_pT", m_min_pt);
144  declareProperty("max_pT", m_max_pt);
145  declareProperty("max_absEta", m_maxAbsEta);
146  declareProperty("minEta", m_minEta);
147  declareProperty("maxEta", m_maxEta);
148  declareProperty("WZModePhotonCone", m_wzPhotonCone);
149 
150  declareProperty("useChargedParticlesOnly", m_useOnlyCharged);
151  declareProperty("MakeListOfStableParticles", m_listPDGofStables);
152  declareProperty("SelectionMode", m_selectionModeName);
153 }
154 
156 = default;
157 
159  m_wzLeptons.clear();
160 }
161 
164 {
165 
166 
167  if( m_selectionModeName == "StableNoMuonNoNu") m_selectionMode = StableNoMuonNoNu ;
168  else if(m_selectionModeName == "MuonOnly") m_selectionMode = MuonOnly ;
169  else if(m_selectionModeName == "NuOnly") m_selectionMode = NuOnly ;
170  else if(m_selectionModeName == "NoWZDecay") m_selectionMode = NoWZDecay ;
171  else {
172  ATH_MSG_ERROR(" Unknown truth selection mode "<< m_selectionMode);
173  return StatusCode::FAILURE;
174  }
175 
176  m_pdgList.clear();
177  m_avP.clear();
178  m_av2P.clear();
179  m_avPt.clear();
180  m_av2Pt.clear();
181  m_avEta.clear();
182  m_av2Eta.clear();
183  m_avPhi.clear();
184  m_av2Phi.clear();
185  m_avM.clear();
186  m_av2M.clear();
187 
188  return StatusCode::SUCCESS;
189 }
190 
191 
193  auto truth4Mom = truthPart->p4();
194  // some safety checks ... both will crash the eta calculation
195  if ( truth4Mom.E() <= 0 ) return false;
196  if ( truth4Mom.P() - truth4Mom.Pz() <= 0 ) return false;
197 
198  double pt = truth4Mom.Pt();
199  if ( m_min_pt > 0 && pt < m_min_pt ) return false;
200  if ( m_max_pt > 0 && pt > m_max_pt ) return false;
201  double eta = truth4Mom.Eta();
202  if ( m_maxAbsEta > 0 && fabs(eta) > m_maxAbsEta ) return false;
203  if ( eta < m_minEta ) return false;
204  if ( eta > m_maxEta ) return false;
205  return true;
206 }
207 
209 
210  bool result = false;
211  if (MC::isGeantino(truthPart)) return false; // protect against unexpected geantinos
212 
213  switch( m_selectionMode) { // For now only 4 modes used in practice for jets...
214  // a switch statement is probably not optimal here...
215  case StableNoMuonNoNu:
216  result =( MC::isStableOrSimDecayed(truthPart) && !HepMC::is_simulation_particle(truthPart) && !MC::isZeroEnergyPhoton(truthPart) &&
217  passKinematics(truthPart) && !MC::isMuon(truthPart) && (!( MC::isStable(truthPart) && MC::isSpecialNonInteracting(truthPart))));
218  break;
219  case NuOnly:
220  result= ( MC::isStableOrSimDecayed(truthPart) && !HepMC::is_simulation_particle(truthPart) && !MC::isZeroEnergyPhoton(truthPart) &&
221  passKinematics(truthPart) && MC::isStable(truthPart) && MC::isSpecialNonInteracting(truthPart));
222  break;
223  case MuonOnly:
224  result = MC::isMuon(truthPart);
225  break;
226  case NoWZDecay:
227  result = (MC::isStableOrSimDecayed(truthPart) && !HepMC::is_simulation_particle(truthPart) && !MC::isZeroEnergyPhoton(truthPart) &&
228  passKinematics(truthPart) &&
229  !isWZDecay(truthPart, m_wzLeptons, m_wzPhotonCone*m_wzPhotonCone) &&
230  !isLeptonFromTau(truthPart) );
231  break;
232  default:
233  break;
234  }
235 
236  if (m_listPDGofStables && result ) {
237  int apId = abs(truthPart->pdgId());
238  auto truth4Mom = truthPart->p4();
239  double p = truth4Mom.P();
240  double phi = truth4Mom.Phi();
241  double pt = truth4Mom.Pt();
242  double eta = truth4Mom.Eta();
243  PDGList::iterator it = m_pdgList.find(apId);
244  if ( it == m_pdgList.end() ) {
245  m_pdgList[apId]=1;
246  m_avP[apId]=p;
247  m_av2P[apId]=p * p;
248  m_avPt[apId]=pt;
249  m_av2Pt[apId]=pt * pt;
250  m_avEta[apId]=eta;
251  m_av2Eta[apId]=eta * eta;
252  m_avPhi[apId]=phi;
253  m_av2Phi[apId]=phi * phi;
254  m_avM[apId]=truthPart->m();
255  m_av2M[apId]=truthPart->m() * truthPart->m();
256  } else {
257  it->second++;
258  m_avP[apId]+=p;
259  m_av2P[apId]+=p * p;
260  m_avPt[apId]+=pt;
261  m_av2Pt[apId]+=pt * pt;
262  m_avEta[apId]+=eta;
263  m_av2Eta[apId]+=eta * eta;
264  m_avPhi[apId]+=phi;
265  m_av2Phi[apId]+=phi * phi;
266  m_avM[apId]+=truthPart->m();
267  m_av2M[apId]+=truthPart->m() * truthPart->m();
268  }
269  }
270 
271  ATH_MSG_VERBOSE("Truth selection result: " << truthPart<<", result=" << result);
272 
273  return result;
274 }
275 
278 {
279  if (m_listPDGofStables)
280  {
281  msg(MSG::INFO) << "Counts of PDGs of all stable particles :" << endmsg;
282  msg(MSG::INFO) << "========================================" << endmsg;
283  msg(MSG::INFO) << "| PDG | # particles |" << endmsg;
284  for ( PDGList::iterator it = m_pdgList.begin(); it != m_pdgList.end(); ++it )
285  {
286  msg(MSG::INFO) << "|"
287  << std::setw(10) << it->first << " |"
288  << std::setw(10) << it->second << " |"
289  << endmsg;
290  }
291  msg(MSG::INFO) << "| PDG | <p> | rms | <pt> | rms | <eta> | rms |" << endmsg;
292  for ( PDGList::iterator it = m_pdgList.begin(); it != m_pdgList.end(); ++it )
293  {
294  int n=it->second;
295  double p=m_avP[it->first]/n;
296  double p2=std::sqrt(m_av2P[it->first]/n-
297  m_avP[it->first]*m_avP[it->first]/n/n);
298  double pt=m_avPt[it->first]/n;
299  double pt2=std::sqrt(m_av2Pt[it->first]/n-
300  m_avPt[it->first]*m_avPt[it->first]/n/n);
301  double eta=m_avEta[it->first]/n;
302  double eta2=std::sqrt(m_av2Eta[it->first]/n-
303  m_avEta[it->first]*m_avEta[it->first]/n/n);
304  msg(MSG::INFO) << "|"
305  << std::setw(15) << it->first << " |"
306  << std::setw(15) << p << " |"
307  << std::setw(15) << p2 << " |"
308  << std::setw(15) << pt << " |"
309  << std::setw(15) << pt2 << " |"
310  << std::setw(15) << eta << " |"
311  << std::setw(15) << eta2 << " |"
312  << endmsg;
313  }
314  msg(MSG::INFO) << "| PDG | <phi> | rms | <m> | rms |" << endmsg;
315  for ( PDGList::iterator it = m_pdgList.begin(); it != m_pdgList.end(); ++it )
316  {
317  int n=it->second;
318  double phi=m_avPhi[it->first]/n;
319  double phi2=std::sqrt(m_av2Phi[it->first]/n-
320  m_avPhi[it->first]*m_avPhi[it->first]/n/n);
321  double m=m_avM[it->first]/n;
322  double m2=0;
323  if ( m_av2M[it->first]/n > m_avM[it->first]*m_avM[it->first]/n/n)
324  m2=std::sqrt(m_av2M[it->first]/n-
325  m_avM[it->first]*m_avM[it->first]/n/n);
326  msg(MSG::INFO) << "|"
327  << std::setw(15) << it->first << " |"
328  << std::setw(15) << phi << " |"
329  << std::setw(15) << phi2 << " |"
330  << std::setw(15) << m << " |"
331  << std::setw(15) << m2 << " |"
332  << endmsg;
333  }
334  }
335  return StatusCode::SUCCESS;
336 }
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
xAOD::iterator
JetConstituentVector::iterator iterator
Definition: JetConstituentVector.cxx:68
JetTruthParticleSelectorTool::JetTruthParticleSelectorTool
JetTruthParticleSelectorTool(const std::string &s)
Definition: JetTruthParticleSelectorTool.cxx:128
get_generator_info.result
result
Definition: get_generator_info.py:21
python.SystemOfUnits.m2
int m2
Definition: SystemOfUnits.py:92
python.SystemOfUnits.m
int m
Definition: SystemOfUnits.py:91
ParticleGun_SamplingFraction.eta2
eta2
Definition: ParticleGun_SamplingFraction.py:96
phi
Scalar phi() const
phi method
Definition: AmgMatrixBasePlugin.h:67
JetTruthParticleSelectorTool::StableNoMuonNoNu
@ StableNoMuonNoNu
Definition: JetTruthParticleSelectorTool.h:52
JetTruthParticleSelectorTool::setupEvent
void setupEvent()
Definition: JetTruthParticleSelectorTool.cxx:158
MC::isSpecialNonInteracting
bool isSpecialNonInteracting(const T &p)
Identify a special non-interacting particles.
Definition: HepMCHelpers.h:83
JetTruthParticleSelectorTool::passKinematics
bool passKinematics(const xAOD::TruthParticle *truth) const
Definition: JetTruthParticleSelectorTool.cxx:192
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:83
xAODP4Helpers.h
AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T > &t)
Definition: AthCommonDataStore.h:145
JetTruthParticleSelectorTool::NoWZDecay
@ NoWZDecay
Definition: JetTruthParticleSelectorTool.h:55
TruthParticleContainer.h
skel.it
it
Definition: skel.GENtoEVGEN.py:396
JetTruthParticleSelectorTool::m_avPt
Average m_avPt
Definition: JetTruthParticleSelectorTool.h:83
asg
Definition: DataHandleTestTool.h:28
test_pyathena.pt
pt
Definition: test_pyathena.py:11
JetTruthParticleSelectorTool::MuonOnly
@ MuonOnly
Definition: JetTruthParticleSelectorTool.h:53
JetTruthParticleSelectorTool::m_min_pt
double m_min_pt
Definition: JetTruthParticleSelectorTool.h:61
JetTruthParticleSelectorTool::m_max_pt
double m_max_pt
Definition: JetTruthParticleSelectorTool.h:62
JetTruthParticleSelectorTool.h
ATH_MSG_VERBOSE
#define ATH_MSG_VERBOSE(x)
Definition: AthMsgStreamMacros.h:28
xAOD::P4Helpers::deltaR2
double deltaR2(double rapidity1, double phi1, double rapidity2, double phi2)
from bare rapidity,phi
Definition: xAODP4Helpers.h:111
JetTruthParticleSelectorTool::NuOnly
@ NuOnly
Definition: JetTruthParticleSelectorTool.h:54
JetTruthParticleSelectorTool::m_avPhi
Average m_avPhi
Definition: JetTruthParticleSelectorTool.h:87
JetTruthParticleSelectorTool::m_av2Pt
Average m_av2Pt
Definition: JetTruthParticleSelectorTool.h:84
xAOD::TruthVertex_v1::outgoingParticleLinks
const TPLinks_t & outgoingParticleLinks() const
Get all the outgoing particles.
JetTruthParticleSelectorTool::m_wzLeptons
std::vector< const xAOD::TruthParticle * > m_wzLeptons
Definition: JetTruthParticleSelectorTool.h:92
isSMLepton
bool isSMLepton(const T &p)
Definition: AtlasPID.h:143
JetTruthParticleSelectorTool::m_wzPhotonCone
double m_wzPhotonCone
Definition: JetTruthParticleSelectorTool.h:66
TRTCalib_cfilter.p2
p2
Definition: TRTCalib_cfilter.py:131
JetTruthParticleSelectorTool::finalize
virtual StatusCode finalize() override
Definition: JetTruthParticleSelectorTool.cxx:277
python.utils.AtlRunQueryDQUtils.p
p
Definition: AtlRunQueryDQUtils.py:210
JetTruthParticleSelectorTool::m_av2Phi
Average m_av2Phi
Definition: JetTruthParticleSelectorTool.h:88
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
JetTruthParticleSelectorTool::m_selectionModeName
std::string m_selectionModeName
Definition: JetTruthParticleSelectorTool.h:74
HepMC::is_simulation_particle
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
Definition: MagicNumbers.h:342
JetTruthParticleSelectorTool::m_av2M
Average m_av2M
Definition: JetTruthParticleSelectorTool.h:90
beamspotman.n
n
Definition: beamspotman.py:731
endmsg
#define endmsg
Definition: AnalysisConfig_Ntuple.cxx:63
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:37
JetTruthParticleSelectorTool::m_av2P
Average m_av2P
Definition: JetTruthParticleSelectorTool.h:82
xAOD::TruthParticle_v1::hasProdVtx
bool hasProdVtx() const
Check for a production vertex on this particle.
Definition: TruthParticle_v1.cxx:74
xAOD::TruthVertex_v1::incomingParticle
const TruthParticle_v1 * incomingParticle(size_t index) const
Get one of the incoming particles.
Definition: TruthVertex_v1.cxx:69
JetTruthParticleSelectorTool::m_listPDGofStables
bool m_listPDGofStables
Definition: JetTruthParticleSelectorTool.h:72
TruthVertex.h
xAOD::TruthParticle_v1::prodVtx
const TruthVertex_v1 * prodVtx() const
The production vertex of this particle.
Definition: TruthParticle_v1.cxx:80
JetTruthParticleSelectorTool::m_useOnlyCharged
bool m_useOnlyCharged
Definition: JetTruthParticleSelectorTool.h:71
xAOD::TruthVertex_v1
Class describing a truth vertex in the MC record.
Definition: TruthVertex_v1.h:37
JetTruthParticleSelectorTool::selector
bool selector(const xAOD::TruthParticle *truth)
Definition: JetTruthParticleSelectorTool.cxx:208
JetTruthParticleSelectorTool::m_maxEta
double m_maxEta
Definition: JetTruthParticleSelectorTool.h:65
JetTruthParticleSelectorTool::m_avP
Average m_avP
Definition: JetTruthParticleSelectorTool.h:81
MagicNumbers.h
JetTruthParticleSelectorTool::~JetTruthParticleSelectorTool
virtual ~JetTruthParticleSelectorTool()
isGeantino
bool isGeantino(const T &p)
Definition: AtlasPID.h:229
JetTruthParticleSelectorTool::m_maxAbsEta
double m_maxAbsEta
Definition: JetTruthParticleSelectorTool.h:63
MC::isStableOrSimDecayed
bool isStableOrSimDecayed(const T &p)
Identify if particle is satble or decayed in simulation.
Definition: HepMCHelpers.h:65
JetTruthParticleSelectorTool::m_av2Eta
Average m_av2Eta
Definition: JetTruthParticleSelectorTool.h:86
MC::isStable
bool isStable(const T &p)
Identify if the particle is stable, i.e. has not decayed.
Definition: HepMCHelpers.h:45
JetTruthParticleSelectorTool::m_avM
Average m_avM
Definition: JetTruthParticleSelectorTool.h:89
xAOD::parentId
@ parentId
Definition: TrackingPrimitives.h:516
xAOD::TruthVertex_v1::nIncomingParticles
size_t nIncomingParticles() const
Get the number of incoming particles.
Definition: TruthVertex_v1.cxx:47
AthCommonMsg< AlgTool >::msg
MsgStream & msg() const
Definition: AthCommonMsg.h:24
JetTruthParticleSelectorTool::m_pdgList
PDGList m_pdgList
Definition: JetTruthParticleSelectorTool.h:78
JetTruthParticleSelectorTool::initialize
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
Definition: JetTruthParticleSelectorTool.cxx:163
xAOD::TruthParticle_v1::p4
virtual FourMom_t p4() const override final
The full 4-momentum of the particle.
Definition: TruthParticle_v1.cxx:196
MC::isZeroEnergyPhoton
bool isZeroEnergyPhoton(const T &p)
Identify a photon with zero energy. Probably a workaround for a generator bug.
Definition: HepMCHelpers.h:71
JetTruthParticleSelectorTool::m_avEta
Average m_avEta
Definition: JetTruthParticleSelectorTool.h:85
xAOD::TruthParticle_v1::pdgId
int pdgId() const
PDG ID code.
JetTruthParticleSelectorTool::m_minEta
double m_minEta
Definition: JetTruthParticleSelectorTool.h:64
JetTruthParticleSelectorTool::m_selectionMode
SelectionMode m_selectionMode
Definition: JetTruthParticleSelectorTool.h:75
xAOD::TruthParticle_v1::m
virtual double m() const override final
The mass of the particle.
HepMCHelpers.h
isMuon
bool isMuon(const T &p)
Definition: AtlasPID.h:154
xAOD::TruthVertex_v1::incomingParticleLinks
const TPLinks_t & incomingParticleLinks() const
Get all the incoming particles.