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  declareProperty("forced_momentum", m_forced_momentum = "", "Forced momentum unit");
30  declareProperty("forced_length", m_forced_length = "", "Forced length unit");
31 }
32 #ifndef HEPMC3
33 //---->//This is copied from MCUtils
35 
36 
38 static inline void reduce(HepMC::GenEvent* ge, HepMC::GenParticle* gp) {
39  // Do nothing if for some reason this particle is not actually in this event
40  if (gp->parent_event() != ge) return;
41 
42  // Get start and end vertices
43  HepMC::GenVertex* vstart = gp->production_vertex();
44  HepMC::GenVertex* vend = gp->end_vertex();
45 
46  // Disconnect the unwanted particle from its vertices and delete it
47  if (vstart != nullptr) vstart->remove_particle(gp);
48  if (vend != nullptr) vend->remove_particle(gp);
49  delete gp;
50 
51  // If start/end vertices are valid and distinct, and this was the only particle that
52  // connected them, then reassign the end vertex decay products to the start vertex
53  // and rewrite the vertex position as most appropriate.
56  if (vstart != nullptr && vend != nullptr && vend != vstart) {
57  bool is_only_link = true;
58  for (auto pchild=vstart->particles_out_const_begin();pchild!=vstart->particles_out_const_end();++pchild) {
59  if ((*pchild)->end_vertex() == vend) is_only_link = false;
60  }
61  if (is_only_link) {
62  if (vend->position() != HepMC::FourVector())
63  vstart->set_position(vend->position()); //< @todo Always use end position if defined... ok?
64  while (vend->particles_out_size() > 0) {
65  vstart->add_particle_out(*vend->particles_out_const_begin());
66  }
67  while (vend->particles_in_size() > 0) {
68  vstart->add_particle_in(*vend->particles_in_const_begin());
69  }
70  }
71  }
72 
73  // Sweep up any vertices orphaned by the particle removal
78  std::vector<HepMC::GenVertex*> orphaned_vtxs;
79  for (HepMC::GenEvent::vertex_const_iterator vi = ge->vertices_begin(); vi != ge->vertices_end(); ++vi) {
80  if ((*vi)->particles_in_size() == 0 && (*vi)->particles_out_size() == 0) orphaned_vtxs.push_back(*vi);
81  }
82  for (HepMC::GenVertex* gv : orphaned_vtxs) delete gv;
83  }
84 
86  inline void reduce(HepMC::GenEvent* ge, std::vector<HepMC::GenParticlePtr> toremove) {
87  while (toremove.size()) {
88  auto gp = toremove.back();
89  toremove.pop_back();
90  reduce(ge, gp);
91  }
92  }
93 //<----//This is copied from MCUtils
94 #endif
95 
97  for (McEventCollection::const_iterator ievt = events()->begin(); ievt != events()->end(); ++ievt) {
98  // FIXME: const_cast
99  HepMC::GenEvent* evt = const_cast<HepMC::GenEvent*>(*ievt);
100  m_looper.findLoops(evt,true);
101  /* AV: To understand why some algorithms fail one would have to request debug explicitly, but that is the reviwers consensus.*/
102  if (!m_looper.loop_particles().empty() || !m_looper.loop_vertices().empty()) {
103  ATH_MSG_DEBUG("Found " << m_looper.loop_vertices().size() << " vertices in loops");
104  ATH_MSG_DEBUG("Found " << m_looper.loop_particles().size() << " particles in loops");
105  ATH_MSG_DEBUG("Please use MC::Loops::findLoops for this event to obtain all particles and vertices in the loops");
106  }
107 #ifdef HEPMC3
108  auto old_momentum = evt->momentum_unit();
109  auto old_length = evt->length_unit();
110  evt->set_units(
111  m_forced_momentum != "" ? HepMC3::Units::momentum_unit(m_forced_momentum): old_momentum,
113  if ( m_forced_momentum != "" && m_forced_momentum != HepMC3::Units::name(old_momentum) ) ATH_MSG_WARNING("Updated momentum units " << HepMC3::Units::name(old_momentum) << "->" << m_forced_momentum);
114  if ( m_forced_length != "" && m_forced_length != HepMC3::Units::name(old_length) ) ATH_MSG_WARNING("Updated length units " << HepMC3::Units::name(old_length) << "->" << m_forced_length);
115 #else
116  std::string old_momentum = (evt->momentum_unit() == HepMC::Units::MomentumUnit::MEV ? "MEV" : "GEV" );
117  std::string old_length = (evt->length_unit() == HepMC::Units::LengthUnit::MM ? "MM" : "CM" );
118  evt->define_units(
119  m_forced_momentum != "" ? m_forced_momentum: old_momentum,
120  m_forced_length != "" ? m_forced_length: old_length);
121  if ( m_forced_momentum != "" && m_forced_momentum != old_momentum ) ATH_MSG_WARNING("Updated momentum units " << old_momentum << "->" << m_forced_momentum);
122  if ( m_forced_length != "" && m_forced_length != old_length ) ATH_MSG_WARNING("Updated length units " << old_length << "->" << m_forced_length);
123 #endif
124  if (!m_pidmap.empty()) {
125  for (auto ip: *evt) {
126  // Skip this particle if (somehow) its pointer is null
127  if (!ip) continue;
128  auto newpid = m_pidmap.find(ip->pdg_id());
129  if (newpid == m_pidmap.end()) continue;
130  ip->set_pdg_id(newpid->second);
131  // Increment the counter for replaced PDG IDs
132  ++m_replacedPIDs;
133  m_replacedpid_counts[ip->pdg_id()]++;
134  }
135  }
136 #ifdef HEPMC3
137  // Add a unit entry to the event weight vector if it's currently empty
138  if (evt->weights().empty()) {
139  ATH_MSG_DEBUG("Adding a unit weight to empty event weight vector");
140  evt->weights().push_back(1);
141  }
142  // Set a (0,0,0) vertex to be the signal vertex if not already set
144  const HepMC::FourVector nullpos;
145  for (auto iv: evt->vertices()) {
146  if (iv->position() == nullpos) {
147  ATH_MSG_DEBUG("Setting representative event position vertex");
149  break;
150  }
151  }
152  }
153 
154  // Catch cases with more than 2 beam particles (16.11.2021)
155  std::vector<HepMC::GenParticlePtr> beams_t;
156  for (const HepMC::GenParticlePtr& p : evt->beams()) {
157  if (p->status() == 4) beams_t.push_back(p);
158  }
159  if (beams_t.size() > 2) {
160  ATH_MSG_INFO("Invalid number of beam particles " << beams_t.size() << ". Will try to fix.");
161  std::vector<HepMC::GenParticlePtr> bparttoremove;
162  for (const auto& bpart : beams_t) {
163  if (bpart->id() == 0 && bpart->production_vertex()) bparttoremove.push_back(bpart);
164  }
165  for (auto bpart: bparttoremove) {
166  bpart->production_vertex()->remove_particle_out(bpart);
167  }
168  }
169 
170  // Some heuristics to catch problematic cases
171  std::vector<HepMC::GenParticlePtr> semi_disconnected, decay_loop_particles;
172  for (auto ip : evt->particles()) {
173  // Skip this particle if (somehow) its pointer is null
174  if (!ip) continue;
175  bool particle_to_fix = false;
176  int abspid = std::abs(ip->pdg_id());
177  auto vProd = ip->production_vertex();
178  auto vEnd = ip->end_vertex();
180  if ( (!vProd || vProd->id() == 0) && vEnd && ip->status() != 4) {
181  particle_to_fix = true;
182  ATH_MSG_DEBUG("Found particle " << ip->pdg_id() << " without production vertex! HepMC status = " << ip->status());
183  }
185  if (vProd && !vEnd && ip->status() != 1) {
186  particle_to_fix = true;
187  ATH_MSG_DEBUG("Found particle " << ip->pdg_id() << " without decay vertex! HepMC status = " << ip->status());
188  }
189  if (particle_to_fix) semi_disconnected.push_back(ip);
190  // Case 3: keep track of loop particles inside decay chains (seen in H7+EvtGen)
191  if (abspid == 43 || abspid == 44 || abspid == 30353 || abspid == 30343) {
192  decay_loop_particles.push_back(std::move(ip));
193  }
194  }
195 
199 
204  if ( !m_ignoreSemiDisconnected && (semi_disconnected.size() == 4 || semi_disconnected.size() == 3 || semi_disconnected.size() == 2)) {
205  size_t no_endv = 0;
206  size_t no_prov = 0;
207  HepMC::FourVector sum(0,0,0,0);
208  std::set<HepMC::GenVertexPtr> standalone;
209  for (const auto& part : semi_disconnected) {
210  if (!part->production_vertex() || !part->production_vertex()->id()) {
211  no_prov++; sum += part->momentum();
212  }
213  if (!part->end_vertex()) {
214  no_endv++; sum -= part->momentum();
215  }
216  if (part->production_vertex()) standalone.insert(part->production_vertex());
217  if (part->end_vertex()) standalone.insert(part->end_vertex());
218  }
219  ATH_MSG_INFO("Heuristics: found " << semi_disconnected.size() << " semi-disconnected particles. Momentum sum is " << sum << " Standalone vertices " << standalone.size());
220  bool standalonevertex = (standalone.size() == 1 && (*standalone.begin())->particles_in().size() + (*standalone.begin())->particles_out().size() == semi_disconnected.size());
222  if (! standalonevertex && no_endv && no_prov && ( no_endv + no_prov == semi_disconnected.size() )) {
223  if (std::abs(sum.px()) < 1e-2 && std::abs(sum.py()) < 1e-2 && std::abs(sum.pz()) < 1e-2 ) {
224  ATH_MSG_INFO("Try " << no_endv << "->" << no_prov << " splitting/merging.");
225  auto v = HepMC::newGenVertexPtr();
226  for (auto part : semi_disconnected) {
227  if (!part->production_vertex() || part->production_vertex()->id() == 0) v->add_particle_out(std::move(part));
228  }
229  for (auto part : semi_disconnected) {
230  if (!part->end_vertex()) v->add_particle_in(std::move(part));
231  }
232  evt->add_vertex(std::move(v));
233  }
234  }
235  }
236 
241  for (auto part: decay_loop_particles) {
243  auto vend = part->end_vertex();
244  auto vprod = part->production_vertex();
245  if (!vprod || !vend) continue;
246  bool loop_in_decay = true;
249  auto sisters = vend->particles_in();
250  for (auto sister: sisters) {
251  if (vprod != sister->production_vertex()) loop_in_decay = false;
252  }
253  if (!loop_in_decay) continue;
254 
256  auto daughters = vend->particles_out();
257  for (auto p : daughters) vprod->add_particle_out(std::move(p));
258  for (auto sister : sisters) {
259  vprod->remove_particle_out(sister);
260  vend->remove_particle_in(sister);
261  evt->remove_particle(std::move(sister));
262  }
263  evt->remove_vertex(std::move(vend));
264 
265  }
266 
267  // Event particle content cleaning -- remove "bad" structures
268  std::vector<HepMC::GenParticlePtr> toremove;
269  for (auto ip: evt->particles()) {
270  // Skip this particle if (somehow) its pointer is null
271  if (!ip) continue;
272  m_totalSeen += 1;
273  // Flag to declare if a particle should be removed
274  bool bad_particle = false;
275  // Check for loops
276  if ( m_killLoops && isSimpleLoop(ip) ) {
277  bad_particle = true;
278  m_loopKilled += 1;
279  ATH_MSG_DEBUG( "Found a looper : " );
281  }
282  // Check on PDG ID 0
283  if ( m_killPDG0 && isPID0(ip) ) {
284  bad_particle = true;
285  m_pdg0Killed += 1;
286  ATH_MSG_DEBUG( "Found PDG ID 0 : " );
288  }
289  // Clean decays
290  int abs_pdg_id = std::abs(ip->pdg_id());
291  bool is_decayed_weak_boson = ( abs_pdg_id == 23 || abs_pdg_id == 24 || abs_pdg_id == 25 ) && ip->end_vertex();
292  if ( m_cleanDecays && isNonTransportableInDecayChain(ip) && !is_decayed_weak_boson ) {
293  bad_particle = true;
294  m_decayCleaned += 1;
295  ATH_MSG_DEBUG( "Found a bad particle in a decay chain : " );
297  }
298  // Only add to the toremove vector once, even if multiple tests match
299  if (bad_particle) toremove.push_back(ip);
300  }
301 
302  // Properties before cleaning
303  const int num_particles_orig = evt->particles().size();
304 
305  // Do the cleaning
306  if (!toremove.empty()) {
307  ATH_MSG_DEBUG("Cleaning event record of " << toremove.size() << " bad particles");
308  for (auto part: toremove) evt->remove_particle(part);
309  }
310 
312  int purged=0;
313  do {
314  purged=0;
315  const std::vector <HepMC::GenParticlePtr> allParticles=evt->particles();
316  for(auto p : allParticles) {
317  HepMC::ConstGenVertexPtr end_v=p->end_vertex();
318  if(p->status() == 2 && !end_v) {
319  evt->remove_particle(std::move(p));
320  ++purged;
322  }
323  }
324  }
325  while (purged>0);
326  }
327 
328  const int num_particles_filt = evt->particles().size();
329 
330  if(num_particles_orig!=num_particles_filt) {
331  // Write out the change in the number of particles
332  ATH_MSG_INFO("Particles filtered: " << num_particles_orig << " -> " << num_particles_filt);
333  }
334  #else
335 
336  // Add a unit entry to the event weight vector if it's currently empty
337  if (evt->weights().empty()) {
338  ATH_MSG_DEBUG("Adding a unit weight to empty event weight vector");
339  evt->weights().push_back(1);
340  }
341 
342  // Set a (0,0,0) vertex to be the signal vertex if not already set
343  if (evt->signal_process_vertex() == NULL) {
344  const HepMC::FourVector nullpos;
345  for (HepMC::GenEvent::vertex_const_iterator iv = evt->vertices_begin(); iv != evt->vertices_end(); ++iv) {
346  if ((*iv)->position() == nullpos) {
347  ATH_MSG_DEBUG("Setting representative event position vertex");
348  evt->set_signal_process_vertex(const_cast<HepMC::GenVertex*>(*iv));
349  break;
350  }
351  }
352  }
353 
354 
355  // Some heuristics to catch problematic cases
356  std::vector<HepMC::GenParticlePtr> semi_disconnected, decay_loop_particles;
357  for (auto ip : *evt) {
358  // Skip this particle if (somehow) its pointer is null
359  if (!ip) continue;
360  bool particle_to_fix = false;
361  int abspid = std::abs(ip->pdg_id());
362  auto vProd = ip->production_vertex();
363  auto vEnd = ip->end_vertex();
365  if ( (!vProd || vProd->id() == 0) && vEnd && ip->status() != 4) {
366  particle_to_fix = true;
367  ATH_MSG_DEBUG("Found particle " << ip->pdg_id() << " without production vertex! HepMC status = " << ip->status());
368  }
370  if (vProd && !vEnd && ip->status() != 1) {
371  particle_to_fix = true;
372  ATH_MSG_DEBUG("Found particle " << ip->pdg_id() << " without decay vertex! HepMC status = " << ip->status());
373  }
374  if (particle_to_fix) semi_disconnected.push_back(ip);
375  // Case 3: keep track of loop particles inside decay chains (seen in H7+EvtGen)
376  if (abspid == 43 || abspid == 44 || abspid == 30353 || abspid == 30343) {
377  decay_loop_particles.push_back(ip);
378  }
379  }
380 
383  if (semi_disconnected.size() == 4 ||semi_disconnected.size() == 3 || semi_disconnected.size() == 2) {
384  size_t no_endv = 0;
385  size_t no_prov = 0;
386  double vsum[4] = {0,0,0,0};
387  for (auto part : semi_disconnected) {
388  if (!part->production_vertex() ) {
389  no_prov++;
390  vsum[0] += part->momentum().px();
391  vsum[1] += part->momentum().py();
392  vsum[2] += part->momentum().pz();
393  vsum[3] += part->momentum().e();
394  }
395  if (!part->end_vertex()) {
396  no_endv++;
397  vsum[0] -= part->momentum().px();
398  vsum[1] -= part->momentum().py();
399  vsum[2] -= part->momentum().pz();
400  vsum[3] -= part->momentum().e();
401  }
402  }
403  HepMC::FourVector sum(vsum[0],vsum[1],vsum[2],vsum[3]);
404  ATH_MSG_INFO("Heuristics: found " << semi_disconnected.size() << " semi-disconnected particles. Momentum sum is " << vsum[0] << " " << vsum[1] << " " << vsum[2] << " " << vsum[3]);
406  if (no_endv && no_prov && ( no_endv + no_prov == semi_disconnected.size() )) {
407  if (std::abs(sum.px()) < 1e-2 && std::abs(sum.py()) < 1e-2 && std::abs(sum.pz()) < 1e-2 ) {
408  ATH_MSG_INFO("Try " << no_endv << "->" << no_prov << " splitting/merging.");
409  auto v = HepMC::newGenVertexPtr();
410  for (auto part : semi_disconnected) {
411  if (!part->production_vertex()) v->add_particle_out(part);
412  }
413  for (auto part : semi_disconnected) {
414  if (!part->end_vertex()) v->add_particle_in(part);
415  }
416  evt->add_vertex(v);
417  }
418  }
419  }
420 
425  for (auto part: decay_loop_particles) {
427  auto vend = part->end_vertex();
428  auto vprod = part->production_vertex();
429  if (!vend || !vprod) continue;
430  bool loop_in_decay = true;
432  std::vector<HepMC::GenParticlePtr> sisters;
433  for (auto p = vend->particles_begin(HepMC::parents); p!= vend->particles_end(HepMC::parents); ++p) sisters.push_back(*p);
434  for (auto sister : sisters) {
435  if (vprod != sister->production_vertex()) loop_in_decay = false;
436  }
437  if (!loop_in_decay) continue;
438 
439  std::vector<HepMC::GenParticlePtr> daughters;
440  for (auto p = vend->particles_begin(HepMC::children); p!= vend->particles_end(HepMC::children); ++p) daughters.push_back(*p);
441  for (auto p : daughters) vprod->add_particle_out(p);
442  for (auto sister : sisters) {
443  vprod->remove_particle(sister);
444  vend->remove_particle(sister);
445  }
446  evt->remove_vertex(vend);
447 
448  }
449 
450 
451  // Event particle content cleaning -- remove "bad" structures
452  std::vector<HepMC::GenParticlePtr> toremove; toremove.reserve(10);
453  for (HepMC::GenEvent::particle_const_iterator ip = evt->particles_begin(); ip != evt->particles_end(); ++ip) {
454  // Skip this particle if (somehow) its pointer is null
455  if (*ip == NULL) continue;
456  m_totalSeen += 1;
457 
458  // Flag to declare if a particle should be removed
459  bool bad_particle = false;
460 
461  // Check for loops
462  if ( m_killLoops && isSimpleLoop(*ip) ) {
463  bad_particle = true;
464  m_loopKilled += 1;
465  ATH_MSG_DEBUG( "Found a looper : " );
466  if ( msgLvl( MSG::DEBUG ) ) (*ip)->print();
467  }
468 
469  // Check on PDG ID 0
470  if ( m_killPDG0 && isPID0(*ip) ) {
471  bad_particle = true;
472  m_pdg0Killed += 1;
473  ATH_MSG_DEBUG( "Found PDG ID 0 : " );
474  if ( msgLvl( MSG::DEBUG ) ) (*ip)->print();
475  }
476 
477  // Clean decays
478  int abs_pdg_id = std::abs((*ip)->pdg_id());
479  bool is_decayed_weak_boson = ( abs_pdg_id == 23 || abs_pdg_id == 24 || abs_pdg_id == 25 ) && (*ip)->end_vertex();
480  if ( m_cleanDecays && isNonTransportableInDecayChain(*ip) && !is_decayed_weak_boson ) {
481  bad_particle = true;
482  m_decayCleaned += 1;
483  ATH_MSG_DEBUG( "Found a bad particle in a decay chain : " );
484  if ( msgLvl( MSG::DEBUG ) ) (*ip)->print();
485  }
486 
487  // Only add to the toremove vector once, even if multiple tests match
488  if (bad_particle) toremove.push_back(*ip);
489  }
490 
491  // Properties before cleaning
492  const int num_particles_orig = evt->particles_size();
493  int num_orphan_vtxs_orig = 0;
494  int num_noparent_vtxs_orig = 0;
495  int num_nochild_vtxs_orig = 0;
496  for (auto v = evt->vertices_begin(); v != evt->vertices_end(); ++v) {
497  if ((*v)->particles_in_size()==0&&(*v)->particles_out_size()==0) num_orphan_vtxs_orig++;
498  if ((*v)->particles_in_size()==0) num_noparent_vtxs_orig++;
499  if ((*v)->particles_out_size()==0) num_nochild_vtxs_orig++;
500  }
501 
502  // Do the cleaning
503  if (!toremove.empty()) {
504  ATH_MSG_DEBUG("Cleaning event record of " << toremove.size() << " bad particles");
505  // Clean!
506  int signal_vertex_bc = evt->signal_process_vertex() ? evt->signal_process_vertex()->barcode() : 0;
507  //This is the only place where reduce is used.
508  reduce(evt , toremove);
509  if (evt->barcode_to_vertex (signal_vertex_bc) == nullptr) {
510  evt->set_signal_process_vertex (nullptr);
511  }
512  }
513 
515  int purged=0;
516  do {
517  for (HepMC::GenParticle* p : *evt) {
518  HepMC::ConstGenVertexPtr end_v = p->end_vertex();
519  if (p->status() == 2 && !end_v) {
520  delete p->production_vertex()->remove_particle(p);
521  ++purged;
523  }
524  }
525  }
526  while (purged>0);
527  }
528 
529  // Properties after cleaning
530  const int num_particles_filt = evt->particles_size();
531  if(num_particles_orig!=num_particles_filt) {
532  int num_orphan_vtxs_filt = 0;
533  int num_noparent_vtxs_filt = 0;
534  int num_nochild_vtxs_filt = 0;
535  for (auto v = evt->vertices_begin(); v != evt->vertices_end(); ++v) {
536  if ((*v)->particles_in_size()==0&&(*v)->particles_out_size()==0) num_orphan_vtxs_filt++;
537  if ((*v)->particles_in_size()==0) num_noparent_vtxs_filt++;
538  if ((*v)->particles_out_size()==0) num_nochild_vtxs_filt++;
539  }
540 
541  // Write out the change in the number of particles
542  ATH_MSG_INFO("Particles filtered: " << num_particles_orig << " -> " << num_particles_filt);
543  // Warn if the numbers of "strange" vertices have changed
544  if (num_orphan_vtxs_filt != num_orphan_vtxs_orig)
545  ATH_MSG_WARNING("Change in orphaned vertices: " << num_orphan_vtxs_orig << " -> " << num_orphan_vtxs_filt);
546  if (num_noparent_vtxs_filt != num_noparent_vtxs_orig)
547  ATH_MSG_WARNING("Change in no-parent vertices: " << num_noparent_vtxs_orig << " -> " << num_noparent_vtxs_filt);
548  if (num_nochild_vtxs_filt != num_nochild_vtxs_orig)
549  ATH_MSG_WARNING("Change in no-parent vertices: " << num_nochild_vtxs_orig << " -> " << num_nochild_vtxs_filt);
550  }
551 #endif
552  }
553  return StatusCode::SUCCESS;
554 }
555 
556 
558  if (m_killLoops ) ATH_MSG_INFO( "Removed " << m_loopKilled << " of " << m_totalSeen << " particles because of loops." );
559  if (m_killPDG0 ) ATH_MSG_INFO( "Removed " << m_pdg0Killed << " of " << m_totalSeen << " particles because of PDG ID 0." );
560  if (m_cleanDecays) ATH_MSG_INFO( "Removed " << m_decayCleaned << " of " << m_totalSeen << " particles while cleaning decay chains." );
561  if(m_purgeUnstableWithoutEndVtx) ATH_MSG_INFO( "Removed " << m_unstablePurged << " of " << m_totalSeen << " unstable particles because they had no decay vertex." );
562  if (!m_pidmap.empty()) {
563  ATH_MSG_INFO( "Replaced " << m_replacedPIDs << " PIDs of particles." );
564  for (const auto& p: m_replacedpid_counts) ATH_MSG_INFO( "Replaced " << p.first << " PIDs " << p.second << "times." );
565  }
566  return StatusCode::SUCCESS;
567 }
568 
569 
571 
572 
573 // Identify PDG ID = 0 particles, usually from HEPEVT padding
575  return p->pdg_id() == 0;
576 }
577 
578 // Identify the particles from
579 bool FixHepMC::fromDecay(const HepMC::ConstGenParticlePtr& p, std::shared_ptr<std::set<int> >& storage) const {
580  if (!p) return false;
581  auto v=p->production_vertex();
582  if (!v) return false;
583 #ifdef HEPMC3
584  for ( const auto& anc: v->particles_in()) {
585  if (MC::isDecayed(anc) && (MC::isTau(anc->pdg_id()) || MC::isHadron(anc->pdg_id()))) return true;
586  }
587  if (storage->find(p->id()) != storage->end()) return false;
588  storage->insert(p->id());
589  for ( const auto& anc: v->particles_in()) {
590  if (fromDecay(anc, storage)) return true;
591  }
592 #else
593  for (auto anc=v->particles_in_const_begin(); anc != v->particles_in_const_end(); ++anc) {
594  if (MC::isDecayed((*anc)) && (MC::isTau((*anc)->pdg_id()) || MC::isHadron((*anc)->pdg_id()))) return true;
595  }
596  if (storage->find(p->barcode()) != storage->end()) return false;
597  storage->insert(p->barcode());
598  for (auto anc=v->particles_in_const_begin(); anc != v->particles_in_const_end(); ++anc){
599  if (fromDecay(*anc, storage)) return true;
600  }
601 #endif
602  return false;
603 }
604 
605 // Identify non-transportable stuff _after_ hadronisation
607  auto storage = std::make_shared<std::set<int>>();
608  return !MC::isTransportable(p->pdg_id()) && fromDecay(p, storage);
609 }
610 
611 // Identify internal "loop" particles
613  return (p->production_vertex() == p->end_vertex() && p->end_vertex() != nullptr);
614 }
615 
617 
618 #endif
619 
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
AllowedVariables::e
e
Definition: AsgElectronSelectorTool.cxx:37
FixHepMC::m_replacedPIDs
long m_replacedPIDs
Definition: FixHepMC.h:67
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:557
FixHepMC::m_pdg0Killed
long m_pdg0Killed
Definition: FixHepMC.h:63
python.DecayParser.parents
parents
print ("==> buf:",buf)
Definition: DecayParser.py:30
GenVertex.h
HepMC::GenParticlePtr
GenParticle * GenParticlePtr
Definition: GenParticle.h:37
FixHepMC::m_cleanDecays
bool m_cleanDecays
Definition: FixHepMC.h:53
ReadBchFromCrest.begin
begin
Definition: ReadBchFromCrest.py:80
MM
@ MM
Definition: RegSelEnums.h:38
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:96
HepMC::Print::line
void line(std::ostream &os, const GenEvent &e)
Definition: GenEvent.h:677
length_unit
constexpr double length_unit
Definition: ActsDetectorElement.cxx:51
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:86
FixHepMC::m_replacedpid_counts
std::map< int, int > m_replacedpid_counts
map of counters of replacements.
Definition: FixHepMC.h:71
HepMC::set_signal_process_vertex
void set_signal_process_vertex(GenEvent *e, T v)
Definition: GenEvent.h:651
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:64
FixHepMC::fromDecay
bool fromDecay(const HepMC::ConstGenParticlePtr &p, std::shared_ptr< std::set< int > > &storage) const
Definition: FixHepMC.cxx:579
FixHepMC::m_forced_momentum
std::string m_forced_momentum
Definition: FixHepMC.h:56
isTau
bool isTau(const T &p)
Definition: AtlasPID.h:208
FixHepMC::m_totalSeen
long m_totalSeen
Definition: FixHepMC.h:66
FixHepMC::m_looper
MC::Loops< HepMC::GenEvent, HepMC::ConstGenParticlePtr, HepMC::ConstGenVertexPtr > m_looper
member to detect loops
Definition: FixHepMC.h:74
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:612
MC::Loops::loop_particles
const std::vector< Prt > & loop_particles() const
Accessor: return the full list of particles in loops.
Definition: Loops.h:36
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:240
isHadron
bool isHadron(const T &p)
Definition: AtlasPID.h:351
FixHepMC::isNonTransportableInDecayChain
bool isNonTransportableInDecayChain(const HepMC::ConstGenParticlePtr &p) const
Definition: FixHepMC.cxx:606
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
Accessor: return the full list of vertices in loops.
Definition: Loops.h:41
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
isTransportable
bool isTransportable(const T &p)
Definition: AtlasPID.h:873
FixHepMC::isPID0
bool isPID0(const HepMC::ConstGenParticlePtr &p) const
Definition: FixHepMC.cxx:574
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:31
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)
Function that does the work to identify loops Pass the event of interest.
Definition: Loops.h:47
FixHepMC::m_loopKilled
long m_loopKilled
Definition: FixHepMC.h:62
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:70
FixHepMC::m_forced_length
std::string m_forced_length
Definition: FixHepMC.h:57
FixHepMC::m_unstablePurged
long m_unstablePurged
Definition: FixHepMC.h:65
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:626