ATLAS Offline Software
Loading...
Searching...
No Matches
DerivationFramework::MenuTruthThinning Class Reference

#include <MenuTruthThinning.h>

Inheritance diagram for DerivationFramework::MenuTruthThinning:
Collaboration diagram for DerivationFramework::MenuTruthThinning:

Public Member Functions

 MenuTruthThinning (const std::string &t, const std::string &n, const IInterface *p)
virtual ~MenuTruthThinning ()
virtual StatusCode initialize () override
virtual StatusCode finalize () override
virtual StatusCode doThinning (const EventContext &ctx) const override
bool isAccepted (const xAOD::TruthParticle *) const
bool matchHadronIncTau (const xAOD::TruthParticle *part) const
bool matchQuarkIncTau (const xAOD::TruthParticle *part) const
bool isOrphanIncTau (const xAOD::TruthParticle *part) const
bool matchGenParticle (const xAOD::TruthParticle *part, std::vector< int > &targetIDs, std::vector< int > &intermediateIDs, bool targetsAreRange) const
bool isLeptonFromTau (const xAOD::TruthParticle *) const
bool isFromTau (const xAOD::TruthParticle *, std::unordered_set< int > &barcode_trace) const
bool isFsrFromLepton (const xAOD::TruthParticle *) const
bool parentIsLongLived (const xAOD::TruthParticle *) const

Static Public Member Functions

static bool isttHFHadron (const xAOD::TruthParticle *)

Private Attributes

StringProperty m_streamName { this, "StreamName", "", "Name of the stream being thinned" }
SG::ThinningHandleKey< xAOD::TruthParticleContainerm_particlesKey { this, "ParticlesKey", "TruthParticles", "TruthParticle container name" }
SG::ThinningHandleKey< xAOD::TruthVertexContainerm_verticesKey { this, "VerticesKey", "TruthVertices", "TruthVertex container name" }
bool m_writePartons
 Parameter: Keep partons?
bool m_writeHadrons
 Parameter: Keep hadrons?
bool m_writeBHadrons
 Parameter: Keep b-hadrons?
bool m_writeCHadrons
 Parameter: Keep c-hadrons?
bool m_writeGeant
 Parameter: Keep geant particles?
float m_geantPhotonPtThresh
 Parameter: Write Geant photons with Pt above this threshold.
bool m_writeTauHad
 Parameter: Keep hadronic tau decays?
bool m_writeBSM
 Parameter: Keep BSM particles?
bool m_writeBosons
 Parameter: Keep bosons?
float m_photonPtCut
float m_partonPtThresh
 Parameter: Write partons with Pt above this threshold.
bool m_writeEverything
 Parameter: Write absolutely everything.
bool m_writeBosonProducts
 Parameter: Write boson decay products.
bool m_writeBSMProducts
 Parameter: Write BSM decay products.
bool m_writeTopAndDecays
 Parameter: Write top and decay products.
bool m_writeAllLeptons
 Parameter: Write all leptons.
bool m_writeLeptonsNotFromHadrons
 Parameter: Write all leptons.
bool m_writeAllStable
bool m_writeNotPhysical
 Parameter: Write particles with status code 3.
bool m_writettHFHadrons
 Parameter: Write particles for tt+HF classification.
int m_writeFirstN
 Parameter: First N particles to write.
bool m_preserveDescendants
 Parameter: preserve descendant/ancestor graph completeness.
bool m_preserveGeneratorDescendants
bool m_preserveAncestors
bool m_preserveImmediate
bool m_preserveHadVtx
std::vector< int > m_pdgIdsToKeep
 Parameter: List of PDG IDs to always keep.
std::vector< int > m_longLivedPdgIds
 Parameter: List of PDG IDs of long lived particles so that one can keep their children.
std::atomic< unsigned int > m_totpart
std::atomic< unsigned int > m_removedpart
std::atomic< int > m_eventCount {}

Detailed Description

Definition at line 27 of file MenuTruthThinning.h.

Constructor & Destructor Documentation

◆ MenuTruthThinning()

DerivationFramework::MenuTruthThinning::MenuTruthThinning ( const std::string & t,
const std::string & n,
const IInterface * p )

Definition at line 29 of file MenuTruthThinning.cxx.

31 :
32base_class(t,n,p),
34m_totpart(0),
36{
37 declareProperty ("WritePartons",
38 m_writePartons = true,
39 "Keep partons?");
40
41 declareProperty ("WriteHadrons",
42 m_writeHadrons = false,
43 "Keep hadrons?");
44
45 declareProperty ("WriteBHadrons",
46 m_writeBHadrons = true,
47 "Keep b-hadrons?");
48
49 declareProperty ("WriteCHadrons",
50 m_writeCHadrons = true,
51 "Keep c-hadrons?");
52
53 declareProperty ("WriteGeant",
54 m_writeGeant = false,
55 "Keep Geant particles?");
56
57 declareProperty ("GeantPhotonPtThresh",
59 "Write Geant photons with Pt above this threshold. "
60 "Set to < 0 to not write any.");
61
62 declareProperty ("WriteTauHad",
63 m_writeTauHad = false,
64 "Keep hadronic taus?");
65
66 declareProperty ("PartonPtThresh",
68 "Write partons with Pt above this threshold.");
69
70 declareProperty ("WriteBSM",
71 m_writeBSM = false,
72 "Keep BSM particles?");
73
74 declareProperty ("WriteBosons",
75 m_writeBosons = false,
76 "Keep bosons?");
77
78 declareProperty ("PhotonPtThresh",
79 m_photonPtCut = 3000.,
80 "Photon pt cut with WriteBosons");
81
82 declareProperty ("WriteBSMProducts",
83 m_writeBSMProducts = false,
84 "Keep BSM particle decay products?");
85
86 declareProperty ("WriteBosonProducts",
88 "Keep boson decay products?");
89
90 declareProperty ("WriteTopAndDecays",
91 m_writeTopAndDecays = false,
92 "Keep top partons and immediate decay products?");
93
94 declareProperty("WriteEverything",
95 m_writeEverything = false,
96 "Keep absolutely everything (overrides all other flags)");
97
98 declareProperty("WriteAllLeptons",
99 m_writeAllLeptons = false,
100 "Keep absolutely all leptons");
101
102 declareProperty("WriteLeptonsNotFromHadrons",
104 "Keep leptons not from hadron decays");
105
106 declareProperty("WriteAllStable",
107 m_writeAllStable = false,
108 "Keep all stable particles");
109
110 declareProperty("WriteNotPhysical",
111 m_writeNotPhysical = false,
112 "Save also non-physical particles");
113
114 declareProperty("WriteFirstN",
115 m_writeFirstN = -1,
116 "Keep first N particles in record");
117
118 declareProperty("PreserveDescendants",
119 m_preserveDescendants = false,
120 "Preserve descendants of retained particles");
121
122 declareProperty("PreserveGeneratorDescendants",
124 "Preserve descendants of retained particles excluding Geant particles");
125
126 declareProperty("PreserveAncestors",
127 m_preserveAncestors = false,
128 "Preserve ancestors of retained particles");
129
130 declareProperty ("PreserveParentsSiblingsChildren",
131 m_preserveImmediate = false,
132 "Preserve the parents, siblings and children of retained particles");
133
134 declareProperty ("PreserveHadronizationVertices",
135 m_preserveHadVtx = false,
136 "Preserve hadronization vertices with parents/children.");
137
138 declareProperty ("WritettHFHadrons",
139 m_writettHFHadrons = false,
140 "Keep tt+HF hadrons?");
141
142 declareProperty ("PDGIDsToKeep",
144 "List of PDG IDs to always keep");
145
146 declareProperty ("LongLivedPDGIDs",
148 "List of PDG IDs of long lived particles so that one can store their children");
149}
std::atomic< unsigned int > m_removedpart
std::atomic< unsigned int > m_totpart
float m_geantPhotonPtThresh
Parameter: Write Geant photons with Pt above this threshold.
bool m_writeGeant
Parameter: Keep geant particles?
bool m_writeTauHad
Parameter: Keep hadronic tau decays?
std::vector< int > m_longLivedPdgIds
Parameter: List of PDG IDs of long lived particles so that one can keep their children.
bool m_writeNotPhysical
Parameter: Write particles with status code 3.
bool m_writeTopAndDecays
Parameter: Write top and decay products.
bool m_writeBHadrons
Parameter: Keep b-hadrons?
bool m_writeLeptonsNotFromHadrons
Parameter: Write all leptons.
bool m_writePartons
Parameter: Keep partons?
bool m_preserveDescendants
Parameter: preserve descendant/ancestor graph completeness.
float m_partonPtThresh
Parameter: Write partons with Pt above this threshold.
bool m_writeEverything
Parameter: Write absolutely everything.
int m_writeFirstN
Parameter: First N particles to write.
bool m_writettHFHadrons
Parameter: Write particles for tt+HF classification.
bool m_writeHadrons
Parameter: Keep hadrons?
std::vector< int > m_pdgIdsToKeep
Parameter: List of PDG IDs to always keep.
bool m_writeBSMProducts
Parameter: Write BSM decay products.
bool m_writeBSM
Parameter: Keep BSM particles?
bool m_writeBosons
Parameter: Keep bosons?
bool m_writeAllLeptons
Parameter: Write all leptons.
bool m_writeBosonProducts
Parameter: Write boson decay products.
bool m_writeCHadrons
Parameter: Keep c-hadrons?

◆ ~MenuTruthThinning()

DerivationFramework::MenuTruthThinning::~MenuTruthThinning ( )
virtual

Definition at line 152 of file MenuTruthThinning.cxx.

152 {
153}

Member Function Documentation

◆ doThinning()

StatusCode DerivationFramework::MenuTruthThinning::doThinning ( const EventContext & ctx) const
overridevirtual

Definition at line 184 of file MenuTruthThinning.cxx.

185{
186
187 // Retrieve the truth collections
188 SG::ThinningHandle<xAOD::TruthParticleContainer> importedTruthParticles
189 (m_particlesKey, ctx);
190 m_totpart += importedTruthParticles->size();
191
192 SG::ThinningHandle<xAOD::TruthVertexContainer> importedTruthVertices
193 (m_verticesKey, ctx);
194
195 // Print events
196 //if( m_eventCount < 20 ){
197 // printxAODTruth(m_eventCount, importedTruthParticles);
198 //}
199 ++m_eventCount;
200
201 // Set up a mask with the same number of entries as the full TruthParticle collection
202 std::vector<bool> particleMask, vertexMask;
203 int nTruthParticles = importedTruthParticles->size();
204 int nTruthVertices = importedTruthVertices->size();
205 particleMask.assign(nTruthParticles,true); // default: keep all particles
206 vertexMask.assign(nTruthVertices,false); // throw all vertices: to be discussed
207
208 // Standard particle loop
209 for (int particleCounter = 0; particleCounter < nTruthParticles; ++particleCounter) {
210 const xAOD::TruthParticle* particle = (*importedTruthParticles)[particleCounter];
211 // Keep first N particles
212 if ( (particleCounter > m_writeFirstN) && (!isAccepted(particle)) ){
213 // Particle removal - note mask is true by default
214 particleMask[particleCounter] = false;
216 }
217 }
218
219 // If user requested preservation of descendants/ancestors:
220 // - loop over the masks and work out which particles need to be descended/ascended from
221 // - do the recursive loop
222 // - update the masks including the descendants/ancestors
223 // To ensure graph completeness, this over-rides anything set by the special treatment
224 // of taus in the section above
225 DerivationFramework::DecayGraphHelper decayHelper;
226 std::unordered_set<int> encounteredUniqueIDs; // For loop handling
228 for (int i=0; i<nTruthParticles; ++i) {
229 bool toKeep = particleMask[i];
230 if (!toKeep) continue;
231 const xAOD::TruthParticle* particle = (*importedTruthParticles)[i];
232 encounteredUniqueIDs.clear();
233 if (m_preserveDescendants) decayHelper.descendants(particle,particleMask,vertexMask,encounteredUniqueIDs,true);
234 encounteredUniqueIDs.clear();
235 if (m_preserveGeneratorDescendants) decayHelper.descendants(particle,particleMask,vertexMask,encounteredUniqueIDs,false);
236 encounteredUniqueIDs.clear();
237 if (m_preserveAncestors) decayHelper.ancestors(particle,particleMask,vertexMask,encounteredUniqueIDs);
238 encounteredUniqueIDs.clear();
239 }
240 }
241 // User only wants to keep parents, siblings and children of retained states
242 // Much simpler case - no recursion so no need for uniqueIDs etc
244 // Make a copy of the particle mask to avoid changes being propagated
245 // down the chain
246 std::vector<bool> particleMaskCopy = particleMask;
247 for (int i=0; i<nTruthParticles; ++i) {
248 bool toKeep = particleMask[i]; // still loop over the orginal list, not the copy
249 if (!toKeep) continue;
250 const xAOD::TruthParticle* particle = (*importedTruthParticles)[i];
251 decayHelper.immediateRelatives(particle,particleMaskCopy,vertexMask,
252 m_preserveHadVtx); // but only update the copy
253 }
254 particleMask = particleMaskCopy; // Now update the original list in one go
255 }
256
257 // Execute the thinning service based on the mask. Finish.
258 importedTruthParticles.keep (particleMask);
259 importedTruthVertices.keep (vertexMask);
260
261 return StatusCode::SUCCESS;
262}
bool isAccepted(const xAOD::TruthParticle *) const
SG::ThinningHandleKey< xAOD::TruthParticleContainer > m_particlesKey
SG::ThinningHandleKey< xAOD::TruthVertexContainer > m_verticesKey
constexpr ParticleHypothesis particle[PARTICLEHYPOTHESES]
the array of masses
TruthParticle_v1 TruthParticle
Typedef to implementation.
void descendants(const xAOD::TruthParticle *pHead, std::vector< int > &particleList, std::unordered_set< int > &encounteredUniqueIDs)
void ancestors(const xAOD::TruthParticle *pHead, std::vector< bool > &particleMask, std::vector< bool > &vertexMask, std::unordered_set< int > &encounteredUniqueIDs)
void immediateRelatives(const xAOD::TruthParticle *pHead, std::vector< bool > &particleMask, std::vector< bool > &vertexMask, bool keepHadVtx)

◆ finalize()

StatusCode DerivationFramework::MenuTruthThinning::finalize ( )
overridevirtual

Definition at line 175 of file MenuTruthThinning.cxx.

176{
177 ATH_MSG_VERBOSE("finalize() ...");
178 ATH_MSG_INFO("Total of " << m_totpart << " truth particles");
179 ATH_MSG_INFO("Removed "<< m_removedpart << " truth particles");
180 return StatusCode::SUCCESS;
181}
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)

◆ initialize()

StatusCode DerivationFramework::MenuTruthThinning::initialize ( )
overridevirtual

Definition at line 156 of file MenuTruthThinning.cxx.

157{
159 ATH_MSG_FATAL("You are asking to keep both all descendants, and only those from the event generator. Please check your job options.");
160 return StatusCode::FAILURE;
161 }
163 ATH_MSG_FATAL("You are asking to preserve only parents/children/siblings of retained states, and also more distant relatives. Please check your job options.");
164 return StatusCode::FAILURE;
165 }
166 m_totpart = 0;
167 m_removedpart = 0;
168 m_eventCount = 0;
169
170 ATH_CHECK( m_particlesKey.initialize (m_streamName) );
171 ATH_CHECK( m_verticesKey.initialize (m_streamName) );
172 return StatusCode::SUCCESS;
173}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_FATAL(x)

◆ isAccepted()

bool DerivationFramework::MenuTruthThinning::isAccepted ( const xAOD::TruthParticle * p) const

Definition at line 265 of file MenuTruthThinning.cxx.

266{
267 bool ok = false;
268
269 int pdg_id = p->absPdgId();
270
271 // All explicitly requested PDG IDs of long lived particles, this is needed
272 // because their childrens uniqueIDs can be above the cut off m_geantOffset
273 if(!m_longLivedPdgIds.empty() && parentIsLongLived(p)) ok = true;
274
275
277 if (!(MC::isPhoton(pdg_id) && m_geantPhotonPtThresh >= 0 && p->pt() > m_geantPhotonPtThresh) )
278 return false;
279 }
280
281 // Do we want to save everything?
282 if (m_writeEverything) ok = true;
283
284 // Save NotPhysical particles
285 if (m_writeNotPhysical && !MC::isPhysical(p)) ok = true;
286
287 // OK if we select partons and are at beginning of event record
288 if( m_writePartons && !MC::isHadron(pdg_id) &&
290 ok = true;
291
292 // OK if we should select hadrons and are in hadron range
293 // JRC: cut changed from PHOTOSMIN to m_geantOffset
294 if( m_writeHadrons && MC::isHadron (pdg_id) && !HepMC::is_simulation_particle(p) ) ok = true;
295
296 // OK if we should select b hadrons and are in hadron range
297 // JRC: cut changed from PHOTOSMIN to m_geantOffset
298 if( m_writeBHadrons && !HepMC::is_simulation_particle(p) && MC::isBottomHadron (pdg_id)) ok = true;
299
300 // OK if we should select c hadrons and are in hadron range
301 // JRC: cut changed from PHOTOSMIN to m_geantOffset
302 if( m_writeCHadrons && !HepMC::is_simulation_particle(p) && MC::isCharmHadron (pdg_id)) ok = true;
303
304 // PHOTOS range: check whether photons come from parton range or
305 // hadron range
306 int motherPDGID = 999999999;
307 // AV: dropped the HepMC::PHOTOSMIN condition
308 if ( !HepMC::is_simulation_particle(p) && p->hasProdVtx() )
309 {
310 const xAOD::TruthVertex* vprod = p->prodVtx();
311 if (vprod->nIncomingParticles() > 0) {
312 const xAOD::TruthParticle* mother = vprod->incomingParticle(0);
313 if (mother) motherPDGID = mother->pdgId();
314 }
315 if( m_writePartons && !MC::isHadron( motherPDGID ) ) ok = true;
316 if( m_writeHadrons && MC::isHadron( motherPDGID ) ) ok = true;
317 }
318
319 // OK if we should select G4 particles and are in G4 range
320 if( m_writeGeant && HepMC::is_simulation_particle(p) ) ok = true;
321
322 if(isLeptonFromTau(p)) ok = true;
323
324 if(isFsrFromLepton(p)) ok = true;
325
326 // Hadronic tau decays
327 if(m_writeTauHad){
328 std::unordered_set<int> uniqueID_trace;
329 if (isFromTau(p, uniqueID_trace))
330 ok = true;
331 }
332
333 // Bosons
334 if(m_writeBosons && (MC::isZ(p)||MC::isW(p)||MC::isHiggs(p) || (MC::isPhoton(p) && p->pt()>m_photonPtCut) )) ok = true;
335
336 // BSM particles
337 if(m_writeBSM && MC::isBSM(p)) ok = true;
338
339 // tt+HF hadrons
340 if (m_writettHFHadrons && isttHFHadron(p)) ok = true;
341
342 // Top quarks
343 if(m_writeTopAndDecays && MC::isTop(p)) ok = true;
344
345 // All leptons
346 if(m_writeAllLeptons && MC::isLepton(p))ok = true; // Include 4th generation...
347
348 // All stable
350
351 // All leptons not from hadron decays
352 if(!ok && m_writeLeptonsNotFromHadrons && MC::isLepton(p) && MC::isStable(p)) {// Include 4th generation...
354 }
355
356 if ((m_writeBSMProducts || m_writeBosonProducts || m_writeTopAndDecays) && p->hasProdVtx()){ // Either way we have to go through parents
357 const xAOD::TruthVertex* prodVtx = p->prodVtx();
358 unsigned int nIncoming = prodVtx->nIncomingParticles();
359 for(unsigned int itr=0; itr<nIncoming; ++itr) {
360 const xAOD::TruthParticle* incomingParticle = prodVtx->incomingParticle(itr);
361 if (!incomingParticle) continue;
362 if ((m_writeBSMProducts && MC::isBSM(incomingParticle)) ||
363 (m_writeBosonProducts && (MC::isZ(incomingParticle)||MC::isW(incomingParticle)||MC::isHiggs(incomingParticle) || (MC::isPhoton(incomingParticle) && incomingParticle->pt()>m_photonPtCut) )) ||
364 (m_writeTopAndDecays && MC::isTop(incomingParticle)) ){
365 ok = true;
366 break;
367 }
368 }
369
370 }
371
372 // All explicitly requested PDG IDs
373 for (const auto id : m_pdgIdsToKeep){
374 if (pdg_id==id) ok=true;
375 } // Loop over PDG IDs
376
377 return ok;
378}
static bool isttHFHadron(const xAOD::TruthParticle *)
bool matchQuarkIncTau(const xAOD::TruthParticle *part) const
bool parentIsLongLived(const xAOD::TruthParticle *) const
bool matchHadronIncTau(const xAOD::TruthParticle *part) const
bool isOrphanIncTau(const xAOD::TruthParticle *part) const
bool isFsrFromLepton(const xAOD::TruthParticle *) const
bool isLeptonFromTau(const xAOD::TruthParticle *) const
bool isFromTau(const xAOD::TruthParticle *, std::unordered_set< int > &barcode_trace) const
int pdgId() const
PDG ID code.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
const TruthParticle_v1 * incomingParticle(size_t index) const
Get one of the incoming particles.
size_t nIncomingParticles() const
Get the number of incoming particles.
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...
bool isZ(const T &p)
bool isCharmHadron(const T &p)
bool isPhoton(const T &p)
bool isW(const T &p)
bool isBottomHadron(const T &p)
bool isStable(const T &p)
Identify if the particle is stable, i.e. has not decayed.
bool isTop(const T &p)
bool isHiggs(const T &p)
APID: HIGGS boson is only one particle.
bool isHadron(const T &p)
bool isLepton(const T &p)
APID: the fourth generation leptons are leptons.
bool isPhysical(const T &p)
Identify if the particle is physical, i.e. is stable or decayed.
bool isBSM(const T &p)
APID: graviton and all Higgs extensions are BSM.
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition TruthVertex.h:15

◆ isFromTau()

bool DerivationFramework::MenuTruthThinning::isFromTau ( const xAOD::TruthParticle * part,
std::unordered_set< int > & barcode_trace ) const

Definition at line 538 of file MenuTruthThinning.cxx.

539 {
540
541 int pdg = part->pdgId();
542
543 const xAOD::TruthVertex* prod = part->prodVtx();
544 if(!prod) return false; // no parent.
545
546 // Simple loop catch
547 if (prod==part->decayVtx()) return false;
548
549 // More complex loop catch
550 std::unordered_set<int>::const_iterator foundVtx = uniqueID_trace.find( HepMC::uniqueID(prod) );
551 if( foundVtx != uniqueID_trace.end() ) {
552 ATH_MSG_DEBUG( "Found a loop (a la Sherpa sample). Backing out." );
553 return false;
554 }
555 uniqueID_trace.insert(HepMC::uniqueID(prod));
556
557 // Loop over the parents of this particle.
558 unsigned int nIncoming = prod->nIncomingParticles();
559 for (unsigned int pitr = 0; pitr<nIncoming; ++pitr){
560 const xAOD::TruthParticle* itrParent = prod->incomingParticle(pitr);
561 if (!itrParent) continue;
562 if (pitr>2) break; // No point in trying - this vertex does not have a quantum meaning...
563
564 if( MC::isTau(itrParent) ) {
565 // Check if one of the children of this parent was a tau - if it is, then it is
566 // photon radiation, and we already cover that under FSR
567 bool has_fsr = false;
568 if ( itrParent->hasDecayVtx() ){
569 unsigned int nChildren = itrParent->decayVtx()->nOutgoingParticles();
570 for (unsigned int citr = 0; citr<nChildren; ++citr) {
571 const xAOD::TruthParticle* itrChild = itrParent->decayVtx()->outgoingParticle(citr);
572 if (!itrChild) continue;
573 if (MC::isTau(itrChild)){
574 has_fsr = true;
575 break;
576 } // Caught FSR check
577 } // loop over immediate children
578 } // Parent has an end vertex
579 if (has_fsr) return false;
580 ATH_MSG_DEBUG("Particle with pdgId = " << pdg << ", matched to tau");
581 return true; // Has tau parent
582 }
583
584 // Go up the generator record until a tau is found or not.
585 if(isFromTau(itrParent, uniqueID_trace)) {
586 return true;
587 }
588 }
589
590 return false;
591}
#define ATH_MSG_DEBUG(x)
const TruthVertex_v1 * decayVtx() const
The decay vertex of this particle.
bool hasDecayVtx() const
Check for a decay vertex on this particle.
const TruthParticle_v1 * outgoingParticle(size_t index) const
Get one of the outgoing particles.
size_t nOutgoingParticles() const
Get the number of outgoing particles.
int uniqueID(const T &p)
bool isTau(const T &p)

◆ isFsrFromLepton()

bool DerivationFramework::MenuTruthThinning::isFsrFromLepton ( const xAOD::TruthParticle * part) const

Definition at line 614 of file MenuTruthThinning.cxx.

614 {
615 int pdg = part->pdgId();
616 if(!MC::isPhoton(part)) return false;
617 if(HepMC::is_simulation_particle(part) ) return false; // Geant photon
618 const xAOD::TruthVertex* prod = part->prodVtx();
619 if(!prod) return false; // no parent.
620 // Simple loop check
621 if (prod==part->decayVtx()) return false;
622 // Loop over the parents of this particle.
623 unsigned int nIncoming = prod->nIncomingParticles();
624 for(unsigned int pitr=0; pitr<nIncoming; ++pitr){
625 const xAOD::TruthParticle* itrParent = prod->incomingParticle(pitr);
626 int parentId = itrParent->pdgId();
627 if(MC::isElectron(parentId) ||
628 MC::isMuon(parentId) ||
629 MC::isTau(parentId)) {
630 ATH_MSG_DEBUG("Photon with uniqueID " << HepMC::uniqueID(part) << " matched to particle with pdgId = " << parentId );
631 return true; // Has lepton parent
632 }
633 if(parentId == pdg) { // Same particle just a different MC status
634 // Go up the generator record until a tau is found or not.
635 if(isFsrFromLepton(itrParent)) return true;
636 }
637 }
638 return false;
639}
bool isElectron(const T &p)
bool isMuon(const T &p)

◆ isLeptonFromTau()

bool DerivationFramework::MenuTruthThinning::isLeptonFromTau ( const xAOD::TruthParticle * part) const

Definition at line 504 of file MenuTruthThinning.cxx.

504 {
505
506 int pdg = part->pdgId();
507
508 if(!MC::isSMLepton(part)) return false; // all leptons including tau.
509
510 const xAOD::TruthVertex* prod = part->prodVtx();
511 if(!prod) return false; // no parent.
512
513 // Simple loop catch
514 if (prod==part->decayVtx()) return false;
515
516 // Loop over the parents of this particle.
517 unsigned int nIncoming = prod->nIncomingParticles();
518 for(unsigned int itr = 0; itr<nIncoming; ++itr){
519 int parentId = prod->incomingParticle(itr)->pdgId();
520 if( MC::isTau(parentId) ) {
521 ATH_MSG_DEBUG("Particle with pdgId = " << pdg << ", matched to tau");
522 return true; // Has tau parent
523 }
524
525 if(parentId == pdg) { // Same particle just a different MC status
526 // Go up the generator record until a tau is found or not.
527 // Note that this requires a connected *lepton* chain, while calling
528 // isFromTau would allow leptons from hadrons from taus
529 if(isLeptonFromTau(prod->incomingParticle(itr))) {
530 return true;
531 }
532 }
533 }
534
535 return false;
536}
bool isSMLepton(const T &p)
APID: the fourth generation leptons are not standard model leptons.

◆ isOrphanIncTau()

bool DerivationFramework::MenuTruthThinning::isOrphanIncTau ( const xAOD::TruthParticle * part) const

Definition at line 396 of file MenuTruthThinning.cxx.

396 {
397 //check if part descends from nothing, possibly passing through a tau decay
398 //(indicating a broken truth record)
399 //ported from TopFiducial
400
401 int pdgId = part->pdgId();
402
403 if (!part->hasProdVtx()) return true;
404
405 const xAOD::TruthVertex* prodVtx = part->prodVtx();
406 unsigned int nIncoming = prodVtx->nIncomingParticles();
407 unsigned int itrParent = 0;
408
409 // Simple loop catch
410 if (prodVtx==part->decayVtx()) return false;
411
412 //0 incoming-> orphan
413 if(nIncoming == 0) return true;
414
415 unsigned int nParents = 0;
416 while(nParents==0 && itrParent<nIncoming) {
417 const xAOD::TruthParticle* incomingParticle = prodVtx->incomingParticle(itrParent);
418
419 // If the parent is itself migrate up the generator record
420 if(incomingParticle->pdgId() == pdgId) {
421 if(!isOrphanIncTau(incomingParticle)) nParents++;
422 }
423 else if(MC::isTau(incomingParticle)) {
424 if(!isOrphanIncTau(incomingParticle)) nParents++;
425 }
426 else {
427 nParents++;
428 }
429 ++itrParent;
430 }
431 return (nParents==0);
432}

◆ isttHFHadron()

bool DerivationFramework::MenuTruthThinning::isttHFHadron ( const xAOD::TruthParticle * part)
static

Definition at line 595 of file MenuTruthThinning.cxx.

595 {
596
597 int ttHFClassification=6;
598
599 static const SG::ConstAccessor<int> TopHadronOriginFlagAcc("TopHadronOriginFlag");
600 if (TopHadronOriginFlagAcc.isAvailable(*part)){
601 ttHFClassification = TopHadronOriginFlagAcc(*part);
602 }
603 else{
604 return false;
605 }
606
607 if (ttHFClassification < 6 ) // useful Hadrons
608 return true;
609
610 return false;
611}

◆ matchGenParticle()

bool DerivationFramework::MenuTruthThinning::matchGenParticle ( const xAOD::TruthParticle * part,
std::vector< int > & targetIDs,
std::vector< int > & intermediateIDs,
bool targetsAreRange ) const

Definition at line 434 of file MenuTruthThinning.cxx.

435 {
436 //check if part descends from a particle, possibly passing through an intermediate decay
437 //ported from TopFiducial
438
439 int pdgId = part->pdgId();
440
441 bool found = false;
442
443 // Iterators for the target pdg or the intermediate pdg ids
444 std::vector<int>::const_iterator itrPdgId, itrPdgIdEnd;
445
446 if (!part->hasProdVtx()) return false;
447
448 // Loop over the parents
449 const xAOD::TruthVertex* prodVtx = part->prodVtx();
450 unsigned int nIncoming = prodVtx->nIncomingParticles();
451
452 // Simple loop catch
453 if (prodVtx==part->decayVtx()) return false;
454
455 unsigned int itrParent = 0;
456 while(!found && itrParent<nIncoming) {
457 const xAOD::TruthParticle* incomingParticle = prodVtx->incomingParticle(itrParent);
458
459 if(!targetsAreRange) {
460
461 // If the parent is one of the particles in the target pdg list
462 itrPdgId = targetPdgIds.begin();
463 itrPdgIdEnd = targetPdgIds.end();
464 for(;itrPdgId != itrPdgIdEnd; ++itrPdgId) {
465 if(incomingParticle->pdgId() == (*itrPdgId)) return true;
466 }
467 }
468 else {
469 int absPdgId = incomingParticle->absPdgId();
470
471 // If the parent is within the range given in the target pdg list
472 if(targetPdgIds.size() == 1) {
473 if(absPdgId >= targetPdgIds.at(0)) return true;
474 }
475 else if(targetPdgIds.size() >= 2) {
476 if(absPdgId >= targetPdgIds.at(0) &&
477 absPdgId <= targetPdgIds.at(1)) return true;
478 }
479 }
480
481 // If the parent is itself migrate up the generator record
482 if(incomingParticle->pdgId() == pdgId) {
483 found = matchGenParticle(incomingParticle, targetPdgIds, intermediatePdgIds, targetsAreRange);
484 }
485 else {
486 // If the parent is in the intermediatePdgIds list migrate up the generator record
487 itrPdgId = intermediatePdgIds.begin();
488 itrPdgIdEnd = intermediatePdgIds.end();
489 bool foundIntermediate = false;
490 while(!foundIntermediate && itrPdgId != itrPdgIdEnd) {
491 if(incomingParticle->pdgId() == (*itrPdgId)) foundIntermediate = true;
492 else ++itrPdgId;
493 }
494 if(foundIntermediate) {
495 found = matchGenParticle(incomingParticle, targetPdgIds, intermediatePdgIds, targetsAreRange);
496 }
497 }
498
499 ++itrParent;
500 }
501 return found;
502}
bool matchGenParticle(const xAOD::TruthParticle *part, std::vector< int > &targetIDs, std::vector< int > &intermediateIDs, bool targetsAreRange) const
int absPdgId() const
Absolute PDG ID code (often useful).

◆ matchHadronIncTau()

bool DerivationFramework::MenuTruthThinning::matchHadronIncTau ( const xAOD::TruthParticle * part) const

Definition at line 380 of file MenuTruthThinning.cxx.

380 {
381 //check if part descends from a hadron, possibly passing through a tau decay
382 //ported from TopFiducial
383 std::vector<int> hadrons = std::vector<int>{111};
384 std::vector<int> taus = std::vector<int>{15,-15};
385 return matchGenParticle(part, hadrons, taus, true);
386}
static Double_t taus

◆ matchQuarkIncTau()

bool DerivationFramework::MenuTruthThinning::matchQuarkIncTau ( const xAOD::TruthParticle * part) const

Definition at line 388 of file MenuTruthThinning.cxx.

388 {
389 //check if part descends from a hadron, possibly passing through a tau decay
390 //ported from TopFiducial
391 std::vector<int> quarks = std::vector<int>{1,6};
392 std::vector<int> taus = std::vector<int>{15,-15};
393 return matchGenParticle(part, quarks, taus, true);
394}

◆ parentIsLongLived()

bool DerivationFramework::MenuTruthThinning::parentIsLongLived ( const xAOD::TruthParticle * part) const

Definition at line 641 of file MenuTruthThinning.cxx.

641 {
642 //loop over the parents of the truth particle, if the parent is in the list of long lived particles
643 //store this truth particle.
644 for(size_t parent_itr = 0; parent_itr < part->nParents(); parent_itr++){
645 if(!part->parent(parent_itr)) continue;
646 const xAOD::TruthParticle* parent = part->parent(parent_itr);
647 const int parent_abs_pdgid = abs(parent->pdgId());
648 if(std::find(m_longLivedPdgIds.begin(), m_longLivedPdgIds.end(), parent_abs_pdgid) != m_longLivedPdgIds.end() ){
649 return true;
650 }
651 }
652 return false;
653}

Member Data Documentation

◆ m_eventCount

std::atomic<int> DerivationFramework::MenuTruthThinning::m_eventCount {}
mutableprivate

Definition at line 138 of file MenuTruthThinning.h.

138{};

◆ m_geantPhotonPtThresh

float DerivationFramework::MenuTruthThinning::m_geantPhotonPtThresh
private

Parameter: Write Geant photons with Pt above this threshold.

Set to < 0 to not write any.

Definition at line 77 of file MenuTruthThinning.h.

◆ m_longLivedPdgIds

std::vector<int> DerivationFramework::MenuTruthThinning::m_longLivedPdgIds
private

Parameter: List of PDG IDs of long lived particles so that one can keep their children.

Definition at line 133 of file MenuTruthThinning.h.

◆ m_particlesKey

SG::ThinningHandleKey<xAOD::TruthParticleContainer> DerivationFramework::MenuTruthThinning::m_particlesKey { this, "ParticlesKey", "TruthParticles", "TruthParticle container name" }
private

Definition at line 55 of file MenuTruthThinning.h.

56{ this, "ParticlesKey", "TruthParticles", "TruthParticle container name" };

◆ m_partonPtThresh

float DerivationFramework::MenuTruthThinning::m_partonPtThresh
private

Parameter: Write partons with Pt above this threshold.

Definition at line 90 of file MenuTruthThinning.h.

◆ m_pdgIdsToKeep

std::vector<int> DerivationFramework::MenuTruthThinning::m_pdgIdsToKeep
private

Parameter: List of PDG IDs to always keep.

Definition at line 130 of file MenuTruthThinning.h.

◆ m_photonPtCut

float DerivationFramework::MenuTruthThinning::m_photonPtCut
private

Definition at line 87 of file MenuTruthThinning.h.

◆ m_preserveAncestors

bool DerivationFramework::MenuTruthThinning::m_preserveAncestors
private

Definition at line 125 of file MenuTruthThinning.h.

◆ m_preserveDescendants

bool DerivationFramework::MenuTruthThinning::m_preserveDescendants
private

Parameter: preserve descendant/ancestor graph completeness.

Definition at line 123 of file MenuTruthThinning.h.

◆ m_preserveGeneratorDescendants

bool DerivationFramework::MenuTruthThinning::m_preserveGeneratorDescendants
private

Definition at line 124 of file MenuTruthThinning.h.

◆ m_preserveHadVtx

bool DerivationFramework::MenuTruthThinning::m_preserveHadVtx
private

Definition at line 127 of file MenuTruthThinning.h.

◆ m_preserveImmediate

bool DerivationFramework::MenuTruthThinning::m_preserveImmediate
private

Definition at line 126 of file MenuTruthThinning.h.

◆ m_removedpart

std::atomic<unsigned int> DerivationFramework::MenuTruthThinning::m_removedpart
mutableprivate

Definition at line 137 of file MenuTruthThinning.h.

◆ m_streamName

StringProperty DerivationFramework::MenuTruthThinning::m_streamName { this, "StreamName", "", "Name of the stream being thinned" }
private

Definition at line 53 of file MenuTruthThinning.h.

54{ this, "StreamName", "", "Name of the stream being thinned" };

◆ m_totpart

std::atomic<unsigned int> DerivationFramework::MenuTruthThinning::m_totpart
mutableprivate

Definition at line 136 of file MenuTruthThinning.h.

◆ m_verticesKey

SG::ThinningHandleKey<xAOD::TruthVertexContainer> DerivationFramework::MenuTruthThinning::m_verticesKey { this, "VerticesKey", "TruthVertices", "TruthVertex container name" }
private

Definition at line 57 of file MenuTruthThinning.h.

58{ this, "VerticesKey", "TruthVertices", "TruthVertex container name" };

◆ m_writeAllLeptons

bool DerivationFramework::MenuTruthThinning::m_writeAllLeptons
private

Parameter: Write all leptons.

Definition at line 105 of file MenuTruthThinning.h.

◆ m_writeAllStable

bool DerivationFramework::MenuTruthThinning::m_writeAllStable
private

Definition at line 111 of file MenuTruthThinning.h.

◆ m_writeBHadrons

bool DerivationFramework::MenuTruthThinning::m_writeBHadrons
private

Parameter: Keep b-hadrons?

Definition at line 67 of file MenuTruthThinning.h.

◆ m_writeBosonProducts

bool DerivationFramework::MenuTruthThinning::m_writeBosonProducts
private

Parameter: Write boson decay products.

Definition at line 96 of file MenuTruthThinning.h.

◆ m_writeBosons

bool DerivationFramework::MenuTruthThinning::m_writeBosons
private

Parameter: Keep bosons?

Definition at line 86 of file MenuTruthThinning.h.

◆ m_writeBSM

bool DerivationFramework::MenuTruthThinning::m_writeBSM
private

Parameter: Keep BSM particles?

Definition at line 83 of file MenuTruthThinning.h.

◆ m_writeBSMProducts

bool DerivationFramework::MenuTruthThinning::m_writeBSMProducts
private

Parameter: Write BSM decay products.

Definition at line 99 of file MenuTruthThinning.h.

◆ m_writeCHadrons

bool DerivationFramework::MenuTruthThinning::m_writeCHadrons
private

Parameter: Keep c-hadrons?

Definition at line 70 of file MenuTruthThinning.h.

◆ m_writeEverything

bool DerivationFramework::MenuTruthThinning::m_writeEverything
private

Parameter: Write absolutely everything.

Definition at line 93 of file MenuTruthThinning.h.

◆ m_writeFirstN

int DerivationFramework::MenuTruthThinning::m_writeFirstN
private

Parameter: First N particles to write.

Definition at line 120 of file MenuTruthThinning.h.

◆ m_writeGeant

bool DerivationFramework::MenuTruthThinning::m_writeGeant
private

Parameter: Keep geant particles?

Definition at line 73 of file MenuTruthThinning.h.

◆ m_writeHadrons

bool DerivationFramework::MenuTruthThinning::m_writeHadrons
private

Parameter: Keep hadrons?

Definition at line 64 of file MenuTruthThinning.h.

◆ m_writeLeptonsNotFromHadrons

bool DerivationFramework::MenuTruthThinning::m_writeLeptonsNotFromHadrons
private

Parameter: Write all leptons.

Definition at line 108 of file MenuTruthThinning.h.

◆ m_writeNotPhysical

bool DerivationFramework::MenuTruthThinning::m_writeNotPhysical
private

Parameter: Write particles with status code 3.

Definition at line 114 of file MenuTruthThinning.h.

◆ m_writePartons

bool DerivationFramework::MenuTruthThinning::m_writePartons
private

Parameter: Keep partons?

Definition at line 61 of file MenuTruthThinning.h.

◆ m_writeTauHad

bool DerivationFramework::MenuTruthThinning::m_writeTauHad
private

Parameter: Keep hadronic tau decays?

Definition at line 80 of file MenuTruthThinning.h.

◆ m_writeTopAndDecays

bool DerivationFramework::MenuTruthThinning::m_writeTopAndDecays
private

Parameter: Write top and decay products.

Definition at line 102 of file MenuTruthThinning.h.

◆ m_writettHFHadrons

bool DerivationFramework::MenuTruthThinning::m_writettHFHadrons
private

Parameter: Write particles for tt+HF classification.

Definition at line 117 of file MenuTruthThinning.h.


The documentation for this class was generated from the following files: