ATLAS Offline Software
FixHepMC.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 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_totalSeen(0)
20  , m_replacedPIDs(0)
21 {
22  declareProperty("KillLoops", m_killLoops = true, "Remove particles in loops?");
23  declareProperty("KillPDG0", m_killPDG0 = true, "Remove particles with PDG ID 0?");
24  declareProperty("CleanDecays", m_cleanDecays = true, "Clean decay chains from non-propagating particles?");
25  declareProperty("LoopsByBarcode", m_loopByBC = false, "Detect loops based on barcodes as well as vertices?");
26  declareProperty("PIDmap", m_pidmap = std::map<int,int>(), "Map of PDG IDs to replace");
27 }
28 #ifndef HEPMC3
29 //---->//This is copied from MCUtils
31 
32 
34 static inline void reduce(HepMC::GenEvent* ge, HepMC::GenParticle* gp) {
35  // Do nothing if for some reason this particle is not actually in this event
36  if (gp->parent_event() != ge) return;
37 
38  // Get start and end vertices
39  HepMC::GenVertex* vstart = gp->production_vertex();
40  HepMC::GenVertex* vend = gp->end_vertex();
41 
42  // Disconnect the unwanted particle from its vertices and delete it
43  if (vstart != nullptr) vstart->remove_particle(gp);
44  if (vend != nullptr) vend->remove_particle(gp);
45  delete gp;
46 
47  // If start/end vertices are valid and distinct, and this was the only particle that
48  // connected them, then reassign the end vertex decay products to the start vertex
49  // and rewrite the vertex position as most appropriate.
52  if (vstart != nullptr && vend != nullptr && vend != vstart) {
53  bool is_only_link = true;
54  for (auto pchild=vstart->particles_out_const_begin();pchild!=vstart->particles_out_const_end();++pchild) {
55  if ((*pchild)->end_vertex() == vend) is_only_link = false;
56  }
57  if (is_only_link) {
58  if (vend->position() != HepMC::FourVector())
59  vstart->set_position(vend->position()); //< @todo Always use end position if defined... ok?
60  while (vend->particles_out_size() > 0) {
61  vstart->add_particle_out(*vend->particles_out_const_begin());
62  }
63  while (vend->particles_in_size() > 0) {
64  vstart->add_particle_in(*vend->particles_in_const_begin());
65  }
66  }
67  }
68 
69  // Sweep up any vertices orphaned by the particle removal
74  std::vector<HepMC::GenVertex*> orphaned_vtxs;
75  for (HepMC::GenEvent::vertex_const_iterator vi = ge->vertices_begin(); vi != ge->vertices_end(); ++vi) {
76  if ((*vi)->particles_in_size() == 0 && (*vi)->particles_out_size() == 0) orphaned_vtxs.push_back(*vi);
77  }
78  for (HepMC::GenVertex* gv : orphaned_vtxs) delete gv;
79  }
80 
82  inline void reduce(HepMC::GenEvent* ge, std::vector<HepMC::GenParticlePtr> toremove) {
83  while (toremove.size()) {
84  auto gp = toremove.back();
85  toremove.pop_back();
86  reduce(ge, gp);
87  }
88  }
89 //<----//This is copied from MCUtils
90 #endif
91 
93  for (McEventCollection::const_iterator ievt = events()->begin(); ievt != events()->end(); ++ievt) {
94  // FIXME: const_cast
95  HepMC::GenEvent* evt = const_cast<HepMC::GenEvent*>(*ievt);
96  if (!m_pidmap.empty()) {
97  for (auto ip: *evt) {
98  // Skip this particle if (somehow) its pointer is null
99  if (!ip) continue;
100  auto newpid = m_pidmap.find(ip->pdg_id());
101  if (newpid == m_pidmap.end()) continue;
102  ip->set_pdg_id(newpid->second);
103  }
104  }
105 #ifdef HEPMC3
106  // Add a unit entry to the event weight vector if it's currently empty
107  if (evt->weights().empty()) {
108  ATH_MSG_DEBUG("Adding a unit weight to empty event weight vector");
109  evt->weights().push_back(1);
110  }
111  // Set a (0,0,0) vertex to be the signal vertex if not already set
113  const HepMC::FourVector nullpos;
114  for (auto iv: evt->vertices()) {
115  if (iv->position() == nullpos) {
116  ATH_MSG_DEBUG("Setting representative event position vertex");
118  break;
119  }
120  }
121  }
122 
123  // Catch cases with more than 2 beam particles (16.11.2021)
124  std::vector<HepMC::GenParticlePtr> beams_t;
125  for (const HepMC::GenParticlePtr& p : evt->beams()) {
126  if (p->status() == 4) beams_t.push_back(p);
127  }
128  if (beams_t.size() > 2) {
129  ATH_MSG_INFO("Invalid number of beam particles " << beams_t.size() << ". Will try to fix.");
130  std::vector<HepMC::GenParticlePtr> bparttoremove;
131  for (const auto& bpart : beams_t) {
132  if (bpart->id() == 0 && bpart->production_vertex()) bparttoremove.push_back(bpart);
133  }
134  for (auto bpart: bparttoremove) {
135  bpart->production_vertex()->remove_particle_out(bpart);
136  }
137  }
138 
139  // Some heuristics to catch problematic cases
140  std::vector<HepMC::GenParticlePtr> semi_disconnected, decay_loop_particles;
141  for (auto ip : evt->particles()) {
142  // Skip this particle if (somehow) its pointer is null
143  if (!ip) continue;
144  bool particle_to_fix = false;
145  int abspid = std::abs(ip->pdg_id());
146  auto vProd = ip->production_vertex();
147  auto vEnd = ip->end_vertex();
149  if ( (!vProd || vProd->id() == 0) && vEnd && ip->status() != 4) {
150  particle_to_fix = true;
151  ATH_MSG_DEBUG("Found particle " << ip->pdg_id() << " without production vertex! HepMC status = " << ip->status());
152  }
154  if (vProd && !vEnd && ip->status() != 1) {
155  particle_to_fix = true;
156  ATH_MSG_DEBUG("Found particle " << ip->pdg_id() << " without decay vertex! HepMC status = " << ip->status());
157  }
158  if (particle_to_fix) semi_disconnected.push_back(ip);
159  // Case 3: keep track of loop particles inside decay chains (seen in H7+EvtGen)
160  if (abspid == 43 || abspid == 44 || abspid == 30353 || abspid == 30343) {
161  decay_loop_particles.push_back(ip);
162  }
163  }
164 
168  if ( semi_disconnected.size() == 4 || semi_disconnected.size() == 3 || semi_disconnected.size() == 2) {
169  size_t no_endv = 0;
170  size_t no_prov = 0;
171  HepMC::FourVector sum(0,0,0,0);
172  for (const auto& part : semi_disconnected) {
173  if (!part->production_vertex() || !part->production_vertex()->id()) {
174  no_prov++; sum += part->momentum();
175  }
176  if (!part->end_vertex()) {
177  no_endv++; sum -= part->momentum();
178  }
179  }
180  ATH_MSG_INFO("Heuristics: found " << semi_disconnected.size() << " semi-disconnected particles. Momentum sum is " << sum);
182  if (no_endv && no_prov && ( no_endv + no_prov == semi_disconnected.size() )) {
183  if (std::abs(sum.px()) < 1e-2 && std::abs(sum.py()) < 1e-2 && std::abs(sum.pz()) < 1e-2 ) {
184  ATH_MSG_INFO("Try " << no_endv << "->" << no_prov << " splitting/merging.");
185  auto v = HepMC::newGenVertexPtr();
186  for (auto part : semi_disconnected) {
187  if (!part->production_vertex() || part->production_vertex()->id() == 0) v->add_particle_out(part);
188  }
189  for (auto part : semi_disconnected) {
190  if (!part->end_vertex()) v->add_particle_in(part);
191  }
192  evt->add_vertex(v);
193  }
194  }
195  }
196 
201  for (auto part: decay_loop_particles) {
203  auto vend = part->end_vertex();
204  auto vprod = part->production_vertex();
205  if (!vprod || !vend) continue;
206  bool loop_in_decay = true;
209  auto sisters = vend->particles_in();
210  for (auto sister: sisters) {
211  if (vprod != sister->production_vertex()) loop_in_decay = false;
212  }
213  if (!loop_in_decay) continue;
214 
216  auto daughters = vend->particles_out();
217  for (auto p : daughters) vprod->add_particle_out(p);
218  for (auto sister : sisters) {
219  vprod->remove_particle_out(sister);
220  vend->remove_particle_in(sister);
221  evt->remove_particle(sister);
222  }
223  evt->remove_vertex(vend);
224 
225  }
226 
227  // Event particle content cleaning -- remove "bad" structures
228  std::vector<HepMC::GenParticlePtr> toremove;
229  for (auto ip: evt->particles()) {
230  // Skip this particle if (somehow) its pointer is null
231  if (!ip) continue;
232  m_totalSeen += 1;
233  // Flag to declare if a particle should be removed
234  bool bad_particle = false;
235  // Check for loops
236  if ( m_killLoops && isLoop(ip) ) {
237  bad_particle = true;
238  m_loopKilled += 1;
239  ATH_MSG_DEBUG( "Found a looper : " );
241  }
242  // Check on PDG ID 0
243  if ( m_killPDG0 && isPID0(ip) ) {
244  bad_particle = true;
245  m_pdg0Killed += 1;
246  ATH_MSG_DEBUG( "Found PDG ID 0 : " );
248  }
249  // Clean decays
250  int abs_pdg_id = std::abs(ip->pdg_id());
251  bool is_decayed_weak_boson = ( abs_pdg_id == 23 || abs_pdg_id == 24 || abs_pdg_id == 25 ) && ip->end_vertex();
252  if ( m_cleanDecays && isNonTransportableInDecayChain(ip) && !is_decayed_weak_boson ) {
253  bad_particle = true;
254  m_decayCleaned += 1;
255  ATH_MSG_DEBUG( "Found a bad particle in a decay chain : " );
257  }
258  // Only add to the toremove vector once, even if multiple tests match
259  if (bad_particle) toremove.push_back(ip);
260  }
261 
262  // Escape here if there's nothing more to do, otherwise do the cleaning
263  if (toremove.empty()) continue;
264  ATH_MSG_DEBUG("Cleaning event record of " << toremove.size() << " bad particles");
265  // Properties before cleaning
266  const int num_particles_orig = evt->particles().size();
267  for (auto part: toremove) evt->remove_particle(part);
268  const int num_particles_filt = evt->particles().size();
269  // Write out the change in the number of particles
270  ATH_MSG_INFO("Particles filtered: " << num_particles_orig << " -> " << num_particles_filt);
271  #else
272 
273  // Add a unit entry to the event weight vector if it's currently empty
274  if (evt->weights().empty()) {
275  ATH_MSG_DEBUG("Adding a unit weight to empty event weight vector");
276  evt->weights().push_back(1);
277  }
278 
279  // Set a (0,0,0) vertex to be the signal vertex if not already set
280  if (evt->signal_process_vertex() == NULL) {
281  const HepMC::FourVector nullpos;
282  for (HepMC::GenEvent::vertex_const_iterator iv = evt->vertices_begin(); iv != evt->vertices_end(); ++iv) {
283  if ((*iv)->position() == nullpos) {
284  ATH_MSG_DEBUG("Setting representative event position vertex");
285  evt->set_signal_process_vertex(const_cast<HepMC::GenVertex*>(*iv));
286  break;
287  }
288  }
289  }
290 
291 
292  // Some heuristics to catch problematic cases
293  std::vector<HepMC::GenParticlePtr> semi_disconnected, decay_loop_particles;
294  for (auto ip : *evt) {
295  // Skip this particle if (somehow) its pointer is null
296  if (!ip) continue;
297  bool particle_to_fix = false;
298  int abspid = std::abs(ip->pdg_id());
299  auto vProd = ip->production_vertex();
300  auto vEnd = ip->end_vertex();
302  if ( (!vProd || vProd->id() == 0) && vEnd && ip->status() != 4) {
303  particle_to_fix = true;
304  ATH_MSG_DEBUG("Found particle " << ip->pdg_id() << " without production vertex! HepMC status = " << ip->status());
305  }
307  if (vProd && !vEnd && ip->status() != 1) {
308  particle_to_fix = true;
309  ATH_MSG_DEBUG("Found particle " << ip->pdg_id() << " without decay vertex! HepMC status = " << ip->status());
310  }
311  if (particle_to_fix) semi_disconnected.push_back(ip);
312  // Case 3: keep track of loop particles inside decay chains (seen in H7+EvtGen)
313  if (abspid == 43 || abspid == 44 || abspid == 30353 || abspid == 30343) {
314  decay_loop_particles.push_back(ip);
315  }
316  }
317 
320  if (semi_disconnected.size() == 4 ||semi_disconnected.size() == 3 || semi_disconnected.size() == 2) {
321  size_t no_endv = 0;
322  size_t no_prov = 0;
323  double vsum[4] = {0,0,0,0};
324  for (auto part : semi_disconnected) {
325  if (!part->production_vertex() ) {
326  no_prov++;
327  vsum[0] += part->momentum().px();
328  vsum[1] += part->momentum().py();
329  vsum[2] += part->momentum().pz();
330  vsum[3] += part->momentum().e();
331  }
332  if (!part->end_vertex()) {
333  no_endv++;
334  vsum[0] -= part->momentum().px();
335  vsum[1] -= part->momentum().py();
336  vsum[2] -= part->momentum().pz();
337  vsum[3] -= part->momentum().e();
338  }
339  }
340  HepMC::FourVector sum(vsum[0],vsum[1],vsum[2],vsum[3]);
341  ATH_MSG_INFO("Heuristics: found " << semi_disconnected.size() << " semi-disconnected particles. Momentum sum is " << vsum[0] << " " << vsum[1] << " " << vsum[2] << " " << vsum[3]);
343  if (no_endv && no_prov && ( no_endv + no_prov == semi_disconnected.size() )) {
344  if (std::abs(sum.px()) < 1e-2 && std::abs(sum.py()) < 1e-2 && std::abs(sum.pz()) < 1e-2 ) {
345  ATH_MSG_INFO("Try " << no_endv << "->" << no_prov << " splitting/merging.");
346  auto v = HepMC::newGenVertexPtr();
347  for (auto part : semi_disconnected) {
348  if (!part->production_vertex()) v->add_particle_out(part);
349  }
350  for (auto part : semi_disconnected) {
351  if (!part->end_vertex()) v->add_particle_in(part);
352  }
353  evt->add_vertex(v);
354  }
355  }
356  }
357 
362  for (auto part: decay_loop_particles) {
364  auto vend = part->end_vertex();
365  auto vprod = part->production_vertex();
366  if (!vend || !vprod) continue;
367  bool loop_in_decay = true;
369  std::vector<HepMC::GenParticlePtr> sisters;
370  for (auto p = vend->particles_begin(HepMC::parents); p!= vend->particles_end(HepMC::parents); ++p) sisters.push_back(*p);
371  for (auto sister : sisters) {
372  if (vprod != sister->production_vertex()) loop_in_decay = false;
373  }
374  if (!loop_in_decay) continue;
375 
376  std::vector<HepMC::GenParticlePtr> daughters;
377  for (auto p = vend->particles_begin(HepMC::children); p!= vend->particles_end(HepMC::children); ++p) daughters.push_back(*p);
378  for (auto p : daughters) vprod->add_particle_out(p);
379  for (auto sister : sisters) {
380  vprod->remove_particle(sister);
381  vend->remove_particle(sister);
382  }
383  evt->remove_vertex(vend);
384 
385  }
386 
387 
388  // Event particle content cleaning -- remove "bad" structures
389  std::vector<HepMC::GenParticlePtr> toremove; toremove.reserve(10);
390  for (HepMC::GenEvent::particle_const_iterator ip = evt->particles_begin(); ip != evt->particles_end(); ++ip) {
391  // Skip this particle if (somehow) its pointer is null
392  if (*ip == NULL) continue;
393  m_totalSeen += 1;
394 
395  // Flag to declare if a particle should be removed
396  bool bad_particle = false;
397 
398  // Check for loops
399  if ( m_killLoops && isLoop(*ip) ) {
400  bad_particle = true;
401  m_loopKilled += 1;
402  ATH_MSG_DEBUG( "Found a looper : " );
403  if ( msgLvl( MSG::DEBUG ) ) (*ip)->print();
404  }
405 
406  // Check on PDG ID 0
407  if ( m_killPDG0 && isPID0(*ip) ) {
408  bad_particle = true;
409  m_pdg0Killed += 1;
410  ATH_MSG_DEBUG( "Found PDG ID 0 : " );
411  if ( msgLvl( MSG::DEBUG ) ) (*ip)->print();
412  }
413 
414  // Clean decays
415  int abs_pdg_id = std::abs((*ip)->pdg_id());
416  bool is_decayed_weak_boson = ( abs_pdg_id == 23 || abs_pdg_id == 24 || abs_pdg_id == 25 ) && (*ip)->end_vertex();
417  if ( m_cleanDecays && isNonTransportableInDecayChain(*ip) && !is_decayed_weak_boson ) {
418  bad_particle = true;
419  m_decayCleaned += 1;
420  ATH_MSG_DEBUG( "Found a bad particle in a decay chain : " );
421  if ( msgLvl( MSG::DEBUG ) ) (*ip)->print();
422  }
423 
424  // Only add to the toremove vector once, even if multiple tests match
425  if (bad_particle) toremove.push_back(*ip);
426  }
427 
428  // Escape here if there's nothing more to do, otherwise do the cleaning
429  if (toremove.empty()) continue;
430  ATH_MSG_DEBUG("Cleaning event record of " << toremove.size() << " bad particles");
431 
432  // Properties before cleaning
433  const int num_particles_orig = evt->particles_size();
434  int num_orphan_vtxs_orig = 0;
435  int num_noparent_vtxs_orig = 0;
436  int num_nochild_vtxs_orig = 0;
437  for (auto v = evt->vertices_begin(); v != evt->vertices_end(); ++v) {
438  if ((*v)->particles_in_size()==0&&(*v)->particles_out_size()==0) num_orphan_vtxs_orig++;
439  if ((*v)->particles_in_size()==0) num_noparent_vtxs_orig++;
440  if ((*v)->particles_out_size()==0) num_nochild_vtxs_orig++;
441  }
442  // Clean!
443  int signal_vertex_bc = evt->signal_process_vertex() ? evt->signal_process_vertex()->barcode() : 0;
444  //This is the only place where reduce is used.
445  reduce(evt , toremove);
446  if (evt->barcode_to_vertex (signal_vertex_bc) == nullptr) {
447  evt->set_signal_process_vertex (nullptr);
448  }
449 
450  // Properties after cleaning
451  const int num_particles_filt = evt->particles_size();
452  int num_orphan_vtxs_filt = 0;
453  int num_noparent_vtxs_filt = 0;
454  int num_nochild_vtxs_filt = 0;
455  for (auto v = evt->vertices_begin(); v != evt->vertices_end(); ++v) {
456  if ((*v)->particles_in_size()==0&&(*v)->particles_out_size()==0) num_orphan_vtxs_filt++;
457  if ((*v)->particles_in_size()==0) num_noparent_vtxs_filt++;
458  if ((*v)->particles_out_size()==0) num_nochild_vtxs_filt++;
459  }
460 
461  // Write out the change in the number of particles
462  ATH_MSG_INFO("Particles filtered: " << num_particles_orig << " -> " << num_particles_filt);
463  // Warn if the numbers of "strange" vertices have changed
464  if (num_orphan_vtxs_filt != num_orphan_vtxs_orig)
465  ATH_MSG_WARNING("Change in orphaned vertices: " << num_orphan_vtxs_orig << " -> " << num_orphan_vtxs_filt);
466  if (num_noparent_vtxs_filt != num_noparent_vtxs_orig)
467  ATH_MSG_WARNING("Change in no-parent vertices: " << num_noparent_vtxs_orig << " -> " << num_noparent_vtxs_filt);
468  if (num_nochild_vtxs_filt != num_nochild_vtxs_orig)
469  ATH_MSG_WARNING("Change in no-parent vertices: " << num_nochild_vtxs_orig << " -> " << num_nochild_vtxs_filt);
470 
471 #endif
472  }
473  return StatusCode::SUCCESS;
474 }
475 
476 
478  if (m_killLoops ) ATH_MSG_INFO( "Removed " << m_loopKilled << " of " << m_totalSeen << " particles because of loops." );
479  if (m_killPDG0 ) ATH_MSG_INFO( "Removed " << m_pdg0Killed << " of " << m_totalSeen << " particles because of PDG ID 0." );
480  if (m_cleanDecays) ATH_MSG_INFO( "Removed " << m_decayCleaned << " of " << m_totalSeen << " particles while cleaning decay chains." );
481  if (!m_pidmap.empty()) ATH_MSG_INFO( "Replaced " << m_replacedPIDs << "PIDs of particles." );
482  return StatusCode::SUCCESS;
483 }
484 
485 
487 
488 
489 // Identify PDG ID = 0 particles, usually from HEPEVT padding
491  return p->pdg_id() == 0;
492 }
493 
494 // Identify the particles from
495 bool FixHepMC::fromDecay(const HepMC::ConstGenParticlePtr& p, std::shared_ptr<std::set<int> >& storage) const {
496  if (!p) return false;
497  auto v=p->production_vertex();
498  if (!v) return false;
499 #ifdef HEPMC3
500  for ( const auto& anc: v->particles_in()) {
501  if (MC::isDecayed(anc) && (MC::isTau(anc->pdg_id()) || MC::isHadron(anc->pdg_id()))) return true;
502  }
503  if (storage->find(p->id()) != storage->end()) return false;
504  storage->insert(p->id());
505  for ( const auto& anc: v->particles_in()) {
506  if (fromDecay(anc, storage)) return true;
507  }
508 #else
509  for (auto anc=v->particles_in_const_begin(); anc != v->particles_in_const_end(); ++anc) {
510  if (MC::isDecayed((*anc)) && (MC::isTau((*anc)->pdg_id()) || MC::isHadron((*anc)->pdg_id()))) return true;
511  }
512  if (storage->find(p->barcode()) != storage->end()) return false;
513  storage->insert(p->barcode());
514  for (auto anc=v->particles_in_const_begin(); anc != v->particles_in_const_end(); ++anc){
515  if (fromDecay(*anc, storage)) return true;
516  }
517 #endif
518  return false;
519 }
520 
521 // Identify non-transportable stuff _after_ hadronisation
523  auto storage = std::make_shared<std::set<int>>();
524  return !MC::isTransportable(p->pdg_id()) && fromDecay(p, storage);
525 }
526 
527 // Identify internal "loop" particles
529  if (p->production_vertex() == p->end_vertex() && p->end_vertex() != NULL) return true;
530  if (m_loopByBC && p->production_vertex()) {
532  int barcodep = HepMC::barcode(p);
533  for (auto itrParent: *(p->production_vertex())) {
534  if ( HepMC::barcode(itrParent) > barcodep) {
535  ATH_MSG_VERBOSE("Found a loop (a la Sherpa sample) via barcode.");
536  return true; // Cannot vectorize, but this is a pretty short loop
537  } // Check on barcodes
538  } // Loop over parent particles
539  } // Has a production vertex
540  return false;
541 }
542 
544 
545 #endif
546 
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
FixHepMC::m_replacedPIDs
long m_replacedPIDs
Definition: FixHepMC.h:62
DataModel_detail::const_iterator
Const iterator class for DataVector/DataList.
Definition: DVLIterator.h:82
GenEvent.h
python.PerfMonSerializer.p
def p
Definition: PerfMonSerializer.py:743
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
FixHepMC::finalize
StatusCode finalize()
Definition: FixHepMC.cxx:477
FixHepMC::m_pdg0Killed
long m_pdg0Killed
Definition: FixHepMC.h:59
AthCommonDataStore< AthCommonMsg< Algorithm > >::declareProperty
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T > &t)
Definition: AthCommonDataStore.h:145
python.DecayParser.parents
parents
print ("==> buf:",buf)
Definition: DecayParser.py:31
FixHepMC::isLoop
bool isLoop(const HepMC::ConstGenParticlePtr &p) const
Definition: FixHepMC.cxx:528
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:52
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:92
ATH_MSG_VERBOSE
#define ATH_MSG_VERBOSE(x)
Definition: AthMsgStreamMacros.h:28
HepMC::Print::line
void line(std::ostream &os, const GenEvent &e)
Definition: GenEvent.h:554
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:82
FixHepMC::m_loopByBC
bool m_loopByBC
Definition: FixHepMC.h:53
HepMC::set_signal_process_vertex
void set_signal_process_vertex(GenEvent *e, T v)
Definition: GenEvent.h:528
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
HepMC::barcode
int barcode(const T *p)
Definition: Barcode.h:16
find_tgc_unfilled_channelids.ip
ip
Definition: find_tgc_unfilled_channelids.py:3
FixHepMC::m_decayCleaned
long m_decayCleaned
Definition: FixHepMC.h:60
FixHepMC::fromDecay
bool fromDecay(const HepMC::ConstGenParticlePtr &p, std::shared_ptr< std::set< int > > &storage) const
Definition: FixHepMC.cxx:495
isTau
bool isTau(const T &p)
Definition: AtlasPID.h:148
FixHepMC::m_totalSeen
long m_totalSeen
Definition: FixHepMC.h:61
HepMC::ConstGenParticlePtr
const GenParticle * ConstGenParticlePtr
Definition: GenParticle.h:38
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:192
isHadron
bool isHadron(const T &p)
Definition: AtlasPID.h:207
FixHepMC::isNonTransportableInDecayChain
bool isNonTransportableInDecayChain(const HepMC::ConstGenParticlePtr &p) const
Definition: FixHepMC.cxx:522
python.PyAthena.v
v
Definition: PyAthena.py:157
DiTauMassTools::MaxHistStrategyV2::e
e
Definition: PhysicsAnalysis/TauID/DiTauMassTools/DiTauMassTools/HelperFunctions.h:26
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
isTransportable
bool isTransportable(const T &p)
Definition: AtlasPID.h:216
FixHepMC::isPID0
bool isPID0(const HepMC::ConstGenParticlePtr &p) const
Definition: FixHepMC.cxx:490
MC::isDecayed
bool isDecayed(const T &p)
Definition: HepMCHelpers.h:29
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:50
FixHepMC::m_loopKilled
long m_loopKilled
Definition: FixHepMC.h:58
FixHepMC::m_killPDG0
bool m_killPDG0
Definition: FixHepMC.h:51
FixHepMC::m_pidmap
std::map< int, int > m_pidmap
map of pids to change.
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:503