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