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