ATLAS Offline Software
FixHepMC.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 #ifndef XAOD_ANALYSIS
6 
9 #include "AtlasHepMC/GenVertex.h"
10 #include "AtlasHepMC/GenEvent.h"
11 #include <vector>
12 
13 
14 FixHepMC::FixHepMC(const std::string& name, ISvcLocator* pSvcLocator)
15  : GenBase(name, pSvcLocator)
16  , m_loopKilled(0)
17  , m_pdg0Killed(0)
18  , m_decayCleaned(0)
19  , m_unstablePurged(0)
20  , m_totalSeen(0)
21  , m_replacedPIDs(0)
22 {
23  declareProperty("KillLoops", m_killLoops = true, "Remove particles in loops?");
24  declareProperty("KillPDG0", m_killPDG0 = true, "Remove particles with PDG ID 0?");
25  declareProperty("CleanDecays", m_cleanDecays = true, "Clean decay chains from non-propagating particles?");
26  declareProperty("PurgeUnstableWithoutEndVtx", m_purgeUnstableWithoutEndVtx = false, "Remove unstable particles without decay vertex?");
27  declareProperty("IgnoreSemiDisconnected", m_ignoreSemiDisconnected = false, "Ignore semi-disconnected particles (normal in Sherpa)");
28  declareProperty("PIDmap", m_pidmap = std::map<int,int>(), "Map of PDG IDs to replace");
29 }
30 #ifndef HEPMC3
31 //---->//This is copied from MCUtils
33 
34 
36 static inline void reduce(HepMC::GenEvent* ge, HepMC::GenParticle* gp) {
37  // Do nothing if for some reason this particle is not actually in this event
38  if (gp->parent_event() != ge) return;
39 
40  // Get start and end vertices
41  HepMC::GenVertex* vstart = gp->production_vertex();
42  HepMC::GenVertex* vend = gp->end_vertex();
43 
44  // Disconnect the unwanted particle from its vertices and delete it
45  if (vstart != nullptr) vstart->remove_particle(gp);
46  if (vend != nullptr) vend->remove_particle(gp);
47  delete gp;
48 
49  // If start/end vertices are valid and distinct, and this was the only particle that
50  // connected them, then reassign the end vertex decay products to the start vertex
51  // and rewrite the vertex position as most appropriate.
54  if (vstart != nullptr && vend != nullptr && vend != vstart) {
55  bool is_only_link = true;
56  for (auto pchild=vstart->particles_out_const_begin();pchild!=vstart->particles_out_const_end();++pchild) {
57  if ((*pchild)->end_vertex() == vend) is_only_link = false;
58  }
59  if (is_only_link) {
60  if (vend->position() != HepMC::FourVector())
61  vstart->set_position(vend->position()); //< @todo Always use end position if defined... ok?
62  while (vend->particles_out_size() > 0) {
63  vstart->add_particle_out(*vend->particles_out_const_begin());
64  }
65  while (vend->particles_in_size() > 0) {
66  vstart->add_particle_in(*vend->particles_in_const_begin());
67  }
68  }
69  }
70 
71  // Sweep up any vertices orphaned by the particle removal
76  std::vector<HepMC::GenVertex*> orphaned_vtxs;
77  for (HepMC::GenEvent::vertex_const_iterator vi = ge->vertices_begin(); vi != ge->vertices_end(); ++vi) {
78  if ((*vi)->particles_in_size() == 0 && (*vi)->particles_out_size() == 0) orphaned_vtxs.push_back(*vi);
79  }
80  for (HepMC::GenVertex* gv : orphaned_vtxs) delete gv;
81  }
82 
84  inline void reduce(HepMC::GenEvent* ge, std::vector<HepMC::GenParticlePtr> toremove) {
85  while (toremove.size()) {
86  auto gp = toremove.back();
87  toremove.pop_back();
88  reduce(ge, gp);
89  }
90  }
91 //<----//This is copied from MCUtils
92 #endif
93 
95  for (McEventCollection::const_iterator ievt = events()->begin(); ievt != events()->end(); ++ievt) {
96  // FIXME: const_cast
97  HepMC::GenEvent* evt = const_cast<HepMC::GenEvent*>(*ievt);
98  m_looper.findLoops(evt,true);
99  /* AV: To understand why some algorithms fail one would have to request debug explicitly, but that is the reviwers consensus.*/
100  if (!m_looper.loop_particles().empty() || !m_looper.loop_vertices().empty()) {
101  ATH_MSG_DEBUG("Found " << m_looper.loop_vertices().size() << " vertices in loops");
102  ATH_MSG_DEBUG("Found " << m_looper.loop_particles().size() << " particles in loops");
103  ATH_MSG_DEBUG("Please use MC::Loops::findLoops for this event to obtain all particles and vertices in the loops");
104  }
105  if (!m_pidmap.empty()) {
106  for (auto ip: *evt) {
107  // Skip this particle if (somehow) its pointer is null
108  if (!ip) continue;
109  auto newpid = m_pidmap.find(ip->pdg_id());
110  if (newpid == m_pidmap.end()) continue;
111  ip->set_pdg_id(newpid->second);
112  }
113  }
114 #ifdef HEPMC3
115  // Add a unit entry to the event weight vector if it's currently empty
116  if (evt->weights().empty()) {
117  ATH_MSG_DEBUG("Adding a unit weight to empty event weight vector");
118  evt->weights().push_back(1);
119  }
120  // Set a (0,0,0) vertex to be the signal vertex if not already set
122  const HepMC::FourVector nullpos;
123  for (auto iv: evt->vertices()) {
124  if (iv->position() == nullpos) {
125  ATH_MSG_DEBUG("Setting representative event position vertex");
127  break;
128  }
129  }
130  }
131 
132  // Catch cases with more than 2 beam particles (16.11.2021)
133  std::vector<HepMC::GenParticlePtr> beams_t;
134  for (const HepMC::GenParticlePtr& p : evt->beams()) {
135  if (p->status() == 4) beams_t.push_back(p);
136  }
137  if (beams_t.size() > 2) {
138  ATH_MSG_INFO("Invalid number of beam particles " << beams_t.size() << ". Will try to fix.");
139  std::vector<HepMC::GenParticlePtr> bparttoremove;
140  for (const auto& bpart : beams_t) {
141  if (bpart->id() == 0 && bpart->production_vertex()) bparttoremove.push_back(bpart);
142  }
143  for (auto bpart: bparttoremove) {
144  bpart->production_vertex()->remove_particle_out(bpart);
145  }
146  }
147 
148  // Some heuristics to catch problematic cases
149  std::vector<HepMC::GenParticlePtr> semi_disconnected, decay_loop_particles;
150  for (auto ip : evt->particles()) {
151  // Skip this particle if (somehow) its pointer is null
152  if (!ip) continue;
153  bool particle_to_fix = false;
154  int abspid = std::abs(ip->pdg_id());
155  auto vProd = ip->production_vertex();
156  auto vEnd = ip->end_vertex();
158  if ( (!vProd || vProd->id() == 0) && vEnd && ip->status() != 4) {
159  particle_to_fix = true;
160  ATH_MSG_DEBUG("Found particle " << ip->pdg_id() << " without production vertex! HepMC status = " << ip->status());
161  }
163  if (vProd && !vEnd && ip->status() != 1) {
164  particle_to_fix = true;
165  ATH_MSG_DEBUG("Found particle " << ip->pdg_id() << " without decay vertex! HepMC status = " << ip->status());
166  }
167  if (particle_to_fix) semi_disconnected.push_back(ip);
168  // Case 3: keep track of loop particles inside decay chains (seen in H7+EvtGen)
169  if (abspid == 43 || abspid == 44 || abspid == 30353 || abspid == 30343) {
170  decay_loop_particles.push_back(ip);
171  }
172  }
173 
177 
182  if ( !m_ignoreSemiDisconnected && (semi_disconnected.size() == 4 || semi_disconnected.size() == 3 || semi_disconnected.size() == 2)) {
183  size_t no_endv = 0;
184  size_t no_prov = 0;
185  HepMC::FourVector sum(0,0,0,0);
186  std::set<HepMC::GenVertexPtr> standalone;
187  for (const auto& part : semi_disconnected) {
188  if (!part->production_vertex() || !part->production_vertex()->id()) {
189  no_prov++; sum += part->momentum();
190  }
191  if (!part->end_vertex()) {
192  no_endv++; sum -= part->momentum();
193  }
194  if (part->production_vertex()) standalone.insert(part->production_vertex());
195  if (part->end_vertex()) standalone.insert(part->end_vertex());
196  }
197  ATH_MSG_INFO("Heuristics: found " << semi_disconnected.size() << " semi-disconnected particles. Momentum sum is " << sum << " Standalone vertices " << standalone.size());
198  bool standalonevertex = (standalone.size() == 1 && (*standalone.begin())->particles_in().size() + (*standalone.begin())->particles_out().size() == semi_disconnected.size());
200  if (! standalonevertex && no_endv && no_prov && ( no_endv + no_prov == semi_disconnected.size() )) {
201  if (std::abs(sum.px()) < 1e-2 && std::abs(sum.py()) < 1e-2 && std::abs(sum.pz()) < 1e-2 ) {
202  ATH_MSG_INFO("Try " << no_endv << "->" << no_prov << " splitting/merging.");
203  auto v = HepMC::newGenVertexPtr();
204  for (auto part : semi_disconnected) {
205  if (!part->production_vertex() || part->production_vertex()->id() == 0) v->add_particle_out(part);
206  }
207  for (auto part : semi_disconnected) {
208  if (!part->end_vertex()) v->add_particle_in(std::move(part));
209  }
210  evt->add_vertex(std::move(v));
211  }
212  }
213  }
214 
219  for (auto part: decay_loop_particles) {
221  auto vend = part->end_vertex();
222  auto vprod = part->production_vertex();
223  if (!vprod || !vend) continue;
224  bool loop_in_decay = true;
227  auto sisters = vend->particles_in();
228  for (auto sister: sisters) {
229  if (vprod != sister->production_vertex()) loop_in_decay = false;
230  }
231  if (!loop_in_decay) continue;
232 
234  auto daughters = vend->particles_out();
235  for (auto p : daughters) vprod->add_particle_out(p);
236  for (auto sister : sisters) {
237  vprod->remove_particle_out(sister);
238  vend->remove_particle_in(sister);
239  evt->remove_particle(std::move(sister));
240  }
241  evt->remove_vertex(std::move(vend));
242 
243  }
244 
245  // Event particle content cleaning -- remove "bad" structures
246  std::vector<HepMC::GenParticlePtr> toremove;
247  for (auto ip: evt->particles()) {
248  // Skip this particle if (somehow) its pointer is null
249  if (!ip) continue;
250  m_totalSeen += 1;
251  // Flag to declare if a particle should be removed
252  bool bad_particle = false;
253  // Check for loops
254  if ( m_killLoops && isSimpleLoop(ip) ) {
255  bad_particle = true;
256  m_loopKilled += 1;
257  ATH_MSG_DEBUG( "Found a looper : " );
259  }
260  // Check on PDG ID 0
261  if ( m_killPDG0 && isPID0(ip) ) {
262  bad_particle = true;
263  m_pdg0Killed += 1;
264  ATH_MSG_DEBUG( "Found PDG ID 0 : " );
266  }
267  // Clean decays
268  int abs_pdg_id = std::abs(ip->pdg_id());
269  bool is_decayed_weak_boson = ( abs_pdg_id == 23 || abs_pdg_id == 24 || abs_pdg_id == 25 ) && ip->end_vertex();
270  if ( m_cleanDecays && isNonTransportableInDecayChain(ip) && !is_decayed_weak_boson ) {
271  bad_particle = true;
272  m_decayCleaned += 1;
273  ATH_MSG_DEBUG( "Found a bad particle in a decay chain : " );
275  }
276  // Only add to the toremove vector once, even if multiple tests match
277  if (bad_particle) toremove.push_back(ip);
278  }
279 
280  // Properties before cleaning
281  const int num_particles_orig = evt->particles().size();
282 
283  // Do the cleaning
284  if (!toremove.empty()) {
285  ATH_MSG_DEBUG("Cleaning event record of " << toremove.size() << " bad particles");
286  for (auto part: toremove) evt->remove_particle(part);
287  }
288 
290  int purged=0;
291  do {
292  purged=0;
293  const std::vector <HepMC::GenParticlePtr> allParticles=evt->particles();
294  for(auto p : allParticles) {
295  HepMC::ConstGenVertexPtr end_v=p->end_vertex();
296  if(p->status() == 2 && !end_v) {
297  evt->remove_particle(std::move(p));
298  ++purged;
300  }
301  }
302  }
303  while (purged>0);
304  }
305 
306  const int num_particles_filt = evt->particles().size();
307 
308  if(num_particles_orig!=num_particles_filt) {
309  // Write out the change in the number of particles
310  ATH_MSG_INFO("Particles filtered: " << num_particles_orig << " -> " << num_particles_filt);
311  }
312  #else
313 
314  // Add a unit entry to the event weight vector if it's currently empty
315  if (evt->weights().empty()) {
316  ATH_MSG_DEBUG("Adding a unit weight to empty event weight vector");
317  evt->weights().push_back(1);
318  }
319 
320  // Set a (0,0,0) vertex to be the signal vertex if not already set
321  if (evt->signal_process_vertex() == NULL) {
322  const HepMC::FourVector nullpos;
323  for (HepMC::GenEvent::vertex_const_iterator iv = evt->vertices_begin(); iv != evt->vertices_end(); ++iv) {
324  if ((*iv)->position() == nullpos) {
325  ATH_MSG_DEBUG("Setting representative event position vertex");
326  evt->set_signal_process_vertex(const_cast<HepMC::GenVertex*>(*iv));
327  break;
328  }
329  }
330  }
331 
332 
333  // Some heuristics to catch problematic cases
334  std::vector<HepMC::GenParticlePtr> semi_disconnected, decay_loop_particles;
335  for (auto ip : *evt) {
336  // Skip this particle if (somehow) its pointer is null
337  if (!ip) continue;
338  bool particle_to_fix = false;
339  int abspid = std::abs(ip->pdg_id());
340  auto vProd = ip->production_vertex();
341  auto vEnd = ip->end_vertex();
343  if ( (!vProd || vProd->id() == 0) && vEnd && ip->status() != 4) {
344  particle_to_fix = true;
345  ATH_MSG_DEBUG("Found particle " << ip->pdg_id() << " without production vertex! HepMC status = " << ip->status());
346  }
348  if (vProd && !vEnd && ip->status() != 1) {
349  particle_to_fix = true;
350  ATH_MSG_DEBUG("Found particle " << ip->pdg_id() << " without decay vertex! HepMC status = " << ip->status());
351  }
352  if (particle_to_fix) semi_disconnected.push_back(ip);
353  // Case 3: keep track of loop particles inside decay chains (seen in H7+EvtGen)
354  if (abspid == 43 || abspid == 44 || abspid == 30353 || abspid == 30343) {
355  decay_loop_particles.push_back(ip);
356  }
357  }
358 
361  if (semi_disconnected.size() == 4 ||semi_disconnected.size() == 3 || semi_disconnected.size() == 2) {
362  size_t no_endv = 0;
363  size_t no_prov = 0;
364  double vsum[4] = {0,0,0,0};
365  for (auto part : semi_disconnected) {
366  if (!part->production_vertex() ) {
367  no_prov++;
368  vsum[0] += part->momentum().px();
369  vsum[1] += part->momentum().py();
370  vsum[2] += part->momentum().pz();
371  vsum[3] += part->momentum().e();
372  }
373  if (!part->end_vertex()) {
374  no_endv++;
375  vsum[0] -= part->momentum().px();
376  vsum[1] -= part->momentum().py();
377  vsum[2] -= part->momentum().pz();
378  vsum[3] -= part->momentum().e();
379  }
380  }
381  HepMC::FourVector sum(vsum[0],vsum[1],vsum[2],vsum[3]);
382  ATH_MSG_INFO("Heuristics: found " << semi_disconnected.size() << " semi-disconnected particles. Momentum sum is " << vsum[0] << " " << vsum[1] << " " << vsum[2] << " " << vsum[3]);
384  if (no_endv && no_prov && ( no_endv + no_prov == semi_disconnected.size() )) {
385  if (std::abs(sum.px()) < 1e-2 && std::abs(sum.py()) < 1e-2 && std::abs(sum.pz()) < 1e-2 ) {
386  ATH_MSG_INFO("Try " << no_endv << "->" << no_prov << " splitting/merging.");
387  auto v = HepMC::newGenVertexPtr();
388  for (auto part : semi_disconnected) {
389  if (!part->production_vertex()) v->add_particle_out(part);
390  }
391  for (auto part : semi_disconnected) {
392  if (!part->end_vertex()) v->add_particle_in(part);
393  }
394  evt->add_vertex(v);
395  }
396  }
397  }
398 
403  for (auto part: decay_loop_particles) {
405  auto vend = part->end_vertex();
406  auto vprod = part->production_vertex();
407  if (!vend || !vprod) continue;
408  bool loop_in_decay = true;
410  std::vector<HepMC::GenParticlePtr> sisters;
411  for (auto p = vend->particles_begin(HepMC::parents); p!= vend->particles_end(HepMC::parents); ++p) sisters.push_back(*p);
412  for (auto sister : sisters) {
413  if (vprod != sister->production_vertex()) loop_in_decay = false;
414  }
415  if (!loop_in_decay) continue;
416 
417  std::vector<HepMC::GenParticlePtr> daughters;
418  for (auto p = vend->particles_begin(HepMC::children); p!= vend->particles_end(HepMC::children); ++p) daughters.push_back(*p);
419  for (auto p : daughters) vprod->add_particle_out(p);
420  for (auto sister : sisters) {
421  vprod->remove_particle(sister);
422  vend->remove_particle(sister);
423  }
424  evt->remove_vertex(vend);
425 
426  }
427 
428 
429  // Event particle content cleaning -- remove "bad" structures
430  std::vector<HepMC::GenParticlePtr> toremove; toremove.reserve(10);
431  for (HepMC::GenEvent::particle_const_iterator ip = evt->particles_begin(); ip != evt->particles_end(); ++ip) {
432  // Skip this particle if (somehow) its pointer is null
433  if (*ip == NULL) continue;
434  m_totalSeen += 1;
435 
436  // Flag to declare if a particle should be removed
437  bool bad_particle = false;
438 
439  // Check for loops
440  if ( m_killLoops && isSimpleLoop(*ip) ) {
441  bad_particle = true;
442  m_loopKilled += 1;
443  ATH_MSG_DEBUG( "Found a looper : " );
444  if ( msgLvl( MSG::DEBUG ) ) (*ip)->print();
445  }
446 
447  // Check on PDG ID 0
448  if ( m_killPDG0 && isPID0(*ip) ) {
449  bad_particle = true;
450  m_pdg0Killed += 1;
451  ATH_MSG_DEBUG( "Found PDG ID 0 : " );
452  if ( msgLvl( MSG::DEBUG ) ) (*ip)->print();
453  }
454 
455  // Clean decays
456  int abs_pdg_id = std::abs((*ip)->pdg_id());
457  bool is_decayed_weak_boson = ( abs_pdg_id == 23 || abs_pdg_id == 24 || abs_pdg_id == 25 ) && (*ip)->end_vertex();
458  if ( m_cleanDecays && isNonTransportableInDecayChain(*ip) && !is_decayed_weak_boson ) {
459  bad_particle = true;
460  m_decayCleaned += 1;
461  ATH_MSG_DEBUG( "Found a bad particle in a decay chain : " );
462  if ( msgLvl( MSG::DEBUG ) ) (*ip)->print();
463  }
464 
465  // Only add to the toremove vector once, even if multiple tests match
466  if (bad_particle) toremove.push_back(*ip);
467  }
468 
469  // Properties before cleaning
470  const int num_particles_orig = evt->particles_size();
471  int num_orphan_vtxs_orig = 0;
472  int num_noparent_vtxs_orig = 0;
473  int num_nochild_vtxs_orig = 0;
474  for (auto v = evt->vertices_begin(); v != evt->vertices_end(); ++v) {
475  if ((*v)->particles_in_size()==0&&(*v)->particles_out_size()==0) num_orphan_vtxs_orig++;
476  if ((*v)->particles_in_size()==0) num_noparent_vtxs_orig++;
477  if ((*v)->particles_out_size()==0) num_nochild_vtxs_orig++;
478  }
479 
480  // Do the cleaning
481  if (!toremove.empty()) {
482  ATH_MSG_DEBUG("Cleaning event record of " << toremove.size() << " bad particles");
483  // Clean!
484  int signal_vertex_bc = evt->signal_process_vertex() ? evt->signal_process_vertex()->barcode() : 0;
485  //This is the only place where reduce is used.
486  reduce(evt , toremove);
487  if (evt->barcode_to_vertex (signal_vertex_bc) == nullptr) {
488  evt->set_signal_process_vertex (nullptr);
489  }
490  }
491 
493  int purged=0;
494  do {
495  for (HepMC::GenParticle* p : *evt) {
496  HepMC::ConstGenVertexPtr end_v = p->end_vertex();
497  if (p->status() == 2 && !end_v) {
498  delete p->production_vertex()->remove_particle(p);
499  ++purged;
501  }
502  }
503  }
504  while (purged>0);
505  }
506 
507  // Properties after cleaning
508  const int num_particles_filt = evt->particles_size();
509  if(num_particles_orig!=num_particles_filt) {
510  int num_orphan_vtxs_filt = 0;
511  int num_noparent_vtxs_filt = 0;
512  int num_nochild_vtxs_filt = 0;
513  for (auto v = evt->vertices_begin(); v != evt->vertices_end(); ++v) {
514  if ((*v)->particles_in_size()==0&&(*v)->particles_out_size()==0) num_orphan_vtxs_filt++;
515  if ((*v)->particles_in_size()==0) num_noparent_vtxs_filt++;
516  if ((*v)->particles_out_size()==0) num_nochild_vtxs_filt++;
517  }
518 
519  // Write out the change in the number of particles
520  ATH_MSG_INFO("Particles filtered: " << num_particles_orig << " -> " << num_particles_filt);
521  // Warn if the numbers of "strange" vertices have changed
522  if (num_orphan_vtxs_filt != num_orphan_vtxs_orig)
523  ATH_MSG_WARNING("Change in orphaned vertices: " << num_orphan_vtxs_orig << " -> " << num_orphan_vtxs_filt);
524  if (num_noparent_vtxs_filt != num_noparent_vtxs_orig)
525  ATH_MSG_WARNING("Change in no-parent vertices: " << num_noparent_vtxs_orig << " -> " << num_noparent_vtxs_filt);
526  if (num_nochild_vtxs_filt != num_nochild_vtxs_orig)
527  ATH_MSG_WARNING("Change in no-parent vertices: " << num_nochild_vtxs_orig << " -> " << num_nochild_vtxs_filt);
528  }
529 #endif
530  }
531  return StatusCode::SUCCESS;
532 }
533 
534 
536  if (m_killLoops ) ATH_MSG_INFO( "Removed " << m_loopKilled << " of " << m_totalSeen << " particles because of loops." );
537  if (m_killPDG0 ) ATH_MSG_INFO( "Removed " << m_pdg0Killed << " of " << m_totalSeen << " particles because of PDG ID 0." );
538  if (m_cleanDecays) ATH_MSG_INFO( "Removed " << m_decayCleaned << " of " << m_totalSeen << " particles while cleaning decay chains." );
539  if(m_purgeUnstableWithoutEndVtx) ATH_MSG_INFO( "Removed " << m_unstablePurged << " of " << m_totalSeen << " unstable particles because they had no decay vertex." );
540  if (!m_pidmap.empty()) ATH_MSG_INFO( "Replaced " << m_replacedPIDs << "PIDs of particles." );
541  return StatusCode::SUCCESS;
542 }
543 
544 
546 
547 
548 // Identify PDG ID = 0 particles, usually from HEPEVT padding
550  return p->pdg_id() == 0;
551 }
552 
553 // Identify the particles from
554 bool FixHepMC::fromDecay(const HepMC::ConstGenParticlePtr& p, std::shared_ptr<std::set<int> >& storage) const {
555  if (!p) return false;
556  auto v=p->production_vertex();
557  if (!v) return false;
558 #ifdef HEPMC3
559  for ( const auto& anc: v->particles_in()) {
560  if (MC::isDecayed(anc) && (MC::isTau(anc->pdg_id()) || MC::isHadron(anc->pdg_id()))) return true;
561  }
562  if (storage->find(p->id()) != storage->end()) return false;
563  storage->insert(p->id());
564  for ( const auto& anc: v->particles_in()) {
565  if (fromDecay(anc, storage)) return true;
566  }
567 #else
568  for (auto anc=v->particles_in_const_begin(); anc != v->particles_in_const_end(); ++anc) {
569  if (MC::isDecayed((*anc)) && (MC::isTau((*anc)->pdg_id()) || MC::isHadron((*anc)->pdg_id()))) return true;
570  }
571  if (storage->find(p->barcode()) != storage->end()) return false;
572  storage->insert(p->barcode());
573  for (auto anc=v->particles_in_const_begin(); anc != v->particles_in_const_end(); ++anc){
574  if (fromDecay(*anc, storage)) return true;
575  }
576 #endif
577  return false;
578 }
579 
580 // Identify non-transportable stuff _after_ hadronisation
582  auto storage = std::make_shared<std::set<int>>();
583  return !MC::isTransportable(p->pdg_id()) && fromDecay(p, storage);
584 }
585 
586 // Identify internal "loop" particles
588  return (p->production_vertex() == p->end_vertex() && p->end_vertex() != nullptr);
589 }
590 
592 
593 #endif
594 
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
AllowedVariables::e
e
Definition: AsgElectronSelectorTool.cxx:37
FixHepMC::m_replacedPIDs
long m_replacedPIDs
Definition: FixHepMC.h:65
DataModel_detail::const_iterator
Const iterator class for DataVector/DataList.
Definition: DVLIterator.h:82
GenEvent.h
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
L1CaloPhase1Monitoring.standalone
standalone
Definition: L1CaloPhase1Monitoring.py:139
FixHepMC::finalize
StatusCode finalize()
Definition: FixHepMC.cxx:535
FixHepMC::m_pdg0Killed
long m_pdg0Killed
Definition: FixHepMC.h:61
python.DecayParser.parents
parents
print ("==> buf:",buf)
Definition: DecayParser.py:31
GenVertex.h
HepMC::GenParticlePtr
GenParticle * GenParticlePtr
Definition: GenParticle.h:37
PlotCalibFromCool.begin
begin
Definition: PlotCalibFromCool.py:94
FixHepMC::m_cleanDecays
bool m_cleanDecays
Definition: FixHepMC.h:53
AthCommonMsg< Algorithm >::msgLvl
bool msgLvl(const MSG::Level lvl) const
Definition: AthCommonMsg.h:30
LArG4FSStartPointFilter.evt
evt
Definition: LArG4FSStartPointFilter.py:42
FixHepMC::execute
StatusCode execute()
Definition: FixHepMC.cxx:94
HepMC::Print::line
void line(std::ostream &os, const GenEvent &e)
Definition: GenEvent.h:676
python.DataFormatRates.events
events
Definition: DataFormatRates.py:105
reduce
void reduce(HepMC::GenEvent *ge, std::vector< HepMC::GenParticlePtr > toremove)
Remove unwanted particles from the event, collapsing the graph structure consistently.
Definition: FixHepMC.cxx:84
HepMC::set_signal_process_vertex
void set_signal_process_vertex(GenEvent *e, T v)
Definition: GenEvent.h:650
python.utils.AtlRunQueryDQUtils.p
p
Definition: AtlRunQueryDQUtils.py:209
GenBase
Base class for common behaviour of MC truth algorithms.
Definition: GenBase.h:47
convertTimingResiduals.sum
sum
Definition: convertTimingResiduals.py:55
HepMC::newGenVertexPtr
GenVertexPtr newGenVertexPtr(const HepMC::FourVector &pos=HepMC::FourVector(0.0, 0.0, 0.0, 0.0), const int i=0)
Definition: GenVertex.h:64
FixHepMC.h
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
AthCommonDataStore< AthCommonMsg< Algorithm > >::declareProperty
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
Definition: AthCommonDataStore.h:145
find_tgc_unfilled_channelids.ip
ip
Definition: find_tgc_unfilled_channelids.py:3
FixHepMC::m_decayCleaned
long m_decayCleaned
Definition: FixHepMC.h:62
FixHepMC::fromDecay
bool fromDecay(const HepMC::ConstGenParticlePtr &p, std::shared_ptr< std::set< int > > &storage) const
Definition: FixHepMC.cxx:554
isTau
bool isTau(const T &p)
Definition: AtlasPID.h:205
FixHepMC::m_totalSeen
long m_totalSeen
Definition: FixHepMC.h:64
FixHepMC::m_looper
MC::Loops< HepMC::GenEvent, HepMC::ConstGenParticlePtr, HepMC::ConstGenVertexPtr > m_looper
member to detect loops
Definition: FixHepMC.h:71
HepMC::ConstGenParticlePtr
const GenParticle * ConstGenParticlePtr
Definition: GenParticle.h:38
FixHepMC::m_purgeUnstableWithoutEndVtx
bool m_purgeUnstableWithoutEndVtx
Definition: FixHepMC.h:54
FixHepMC::isSimpleLoop
bool isSimpleLoop(const HepMC::ConstGenParticlePtr &p) const
Definition: FixHepMC.cxx:587
MC::Loops::loop_particles
const std::vector< Prt > & loop_particles() const
Definition: Loops.h:25
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:240
isHadron
bool isHadron(const T &p)
Definition: AtlasPID.h:348
FixHepMC::isNonTransportableInDecayChain
bool isNonTransportableInDecayChain(const HepMC::ConstGenParticlePtr &p) const
Definition: FixHepMC.cxx:581
python.PyAthena.v
v
Definition: PyAthena.py:154
FixHepMC::m_ignoreSemiDisconnected
bool m_ignoreSemiDisconnected
Definition: FixHepMC.h:55
MC::Loops::loop_vertices
const std::vector< Vtx > & loop_vertices() const
Definition: Loops.h:28
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
isTransportable
bool isTransportable(const T &p)
Definition: AtlasPID.h:862
FixHepMC::isPID0
bool isPID0(const HepMC::ConstGenParticlePtr &p) const
Definition: FixHepMC.cxx:549
MC::isDecayed
bool isDecayed(const T &p)
Identify if the particle decayed.
Definition: HepMCHelpers.h:42
DEBUG
#define DEBUG
Definition: page_access.h:11
python.DecayParser.children
children
Definition: DecayParser.py:32
FixHepMC::m_killLoops
bool m_killLoops
Definition: FixHepMC.h:51
HepMC::ConstGenVertexPtr
const HepMC::GenVertex * ConstGenVertexPtr
Definition: GenVertex.h:60
MC::Loops::findLoops
int findLoops(const Evt *evt, bool force)
Definition: Loops.h:31
FixHepMC::m_loopKilled
long m_loopKilled
Definition: FixHepMC.h:60
FixHepMC::m_killPDG0
bool m_killPDG0
Definition: FixHepMC.h:52
FixHepMC::m_pidmap
std::map< int, int > m_pidmap
map of pids to change.
Definition: FixHepMC.h:68
FixHepMC::m_unstablePurged
long m_unstablePurged
Definition: FixHepMC.h:63
HepMCHelpers.h
GenParticle
@ GenParticle
Definition: TruthClasses.h:30
FixHepMC::FixHepMC
FixHepMC(const std::string &name, ISvcLocator *pSvcLocator)
Definition: FixHepMC.cxx:14
HepMC::signal_process_vertex
GenVertex * signal_process_vertex(const GenEvent *e)
Definition: GenEvent.h:625