ATLAS Offline Software
Loading...
Searching...
No Matches
GenEvent.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4/* Author: Andrii Verbytskyi andrii.verbytskyi@mpp.mpg.de */
5
6#ifndef ATLASHEPMC_GENEVENT_H
7#define ATLASHEPMC_GENEVENT_H
8#ifdef HEPMC3
9#undef private
10#undef protected
11#include "HepMC3/GenEvent.h"
12#include "HepMC3/GenHeavyIon.h"
13#include "HepMC3/GenPdfInfo.h"
14#include "HepMC3/PrintStreams.h"
15#include "AtlasHepMC/Barcode.h"
19
20#include <unordered_map>
21
22namespace HepMC3 {
23inline std::vector<HepMC3::GenParticlePtr>::const_iterator begin(HepMC3::GenEvent& e) { return e.particles().begin(); }
24inline std::vector<HepMC3::GenParticlePtr>::const_iterator end(HepMC3::GenEvent& e) { return e.particles().end(); }
25inline std::vector<HepMC3::ConstGenParticlePtr>::const_iterator begin(const HepMC3::GenEvent& e) { return e.particles().begin(); }
26inline std::vector<HepMC3::ConstGenParticlePtr>::const_iterator end(const HepMC3::GenEvent& e) { return e.particles().end(); }
27}
28
29namespace HepMC {
30using Print=HepMC3::Print;
31using GenHeavyIon=HepMC3::GenHeavyIon;
32using GenEvent=HepMC3::GenEvent;
33
34class ShortEventAttribute : public HepMC3::Attribute {
35public:
36 ShortEventAttribute():HepMC3::Attribute(){}
37 ShortEventAttribute(const HepMC3::GenEvent* e):HepMC3::Attribute(){ from_event(e); }
38
39 bool from_event(const HepMC3::GenEvent* e){
40 NUP=e->particles().size();
41 resize();
42 XWGTUP = e->weights().size() ? e->weights()[0] : 1.0;
43 auto A_signal_process_id=e->attribute<HepMC3::IntAttribute>("signal_process_id");
44 IDPRUP = A_signal_process_id?A_signal_process_id->value() : 0;
45 auto A_event_scale=e->attribute<HepMC3::DoubleAttribute>("event_scale");
46 SCALUP = A_event_scale? A_event_scale->value():0;
47 auto A_alphaQCD=e->attribute<HepMC3::DoubleAttribute>("alphaQCD");
48 AQCDUP = A_alphaQCD? A_alphaQCD->value():0;
49 auto A_alphaQED=e->attribute<HepMC3::DoubleAttribute>("alphaQED");
50 AQEDUP = A_alphaQED? A_alphaQED->value():0;
51
52 for ( int i = 0; i < NUP; ++i ){
53 PUP[i][0] = e->particles().at(i)->momentum().px();
54 PUP[i][1] = e->particles().at(i)->momentum().py();
55 PUP[i][2] = e->particles().at(i)->momentum().pz();
56 PUP[i][3] = e->particles().at(i)->momentum().e();
57 PUP[i][4] = e->particles().at(i)->momentum().m();
58 IDUP[i] = e->particles().at(i)->pdg_id();
59 auto pv = e->particles().at(i)->production_vertex();
60 auto ev = e->particles().at(i)->end_vertex();
61 if (pv && ev ) ISTUP[i] = 2;
62 if (pv && !ev ) ISTUP[i] = 1;
63 if (e->particles().at(i)->status() == 4 || !pv ) ISTUP[i] = -1;
64 auto flow1 = e->particles().at(i)->attribute<HepMC3::IntAttribute>("flow1");
65 auto flow2 = e->particles().at(i)->attribute<HepMC3::IntAttribute>("flow2");
66 ICOLUP[i].first = flow1 ? flow1->value() : 0;
67 ICOLUP[i].second = flow2 ? flow2->value() : 0;
68 int l = 0,h = 0;
69 if (pv) for (auto p: pv->particles_in()) { l= std::min(l,p->id()); h=std::max(h,p->id());}
70 MOTHUP[i].first = h;
71 MOTHUP[i].second = l;
72 }
73 return true;
74 }
75
76 bool from_string(const std::string &att) override {
77 std::istringstream iss(att);
78 iss >> NUP;
79 iss >> IDPRUP;
80 iss >> XWGTUP;
81 iss >> SCALUP;
82 iss >> AQEDUP;
83 iss >> AQCDUP;
84 //assume input is a trusted source
85 //coverity[tainted_data]
86 resize();
87 for ( int i = 0; i < NUP; ++i ){
88 iss >> IDUP[i];
89 iss >> ISTUP[i];
90 iss >> MOTHUP[i].first;
91 iss >> MOTHUP[i].second;
92 iss >> ICOLUP[i].first;
93 iss >> ICOLUP[i].second;
94 iss >> PUP[i][0];
95 iss >> PUP[i][1];
96 iss >> PUP[i][2];
97 iss >> PUP[i][3];
98 iss >> PUP[i][4];
99 iss >> VTIMUP[i];
100 iss >> SPINUP[i];
101 }
102 set_is_parsed(true);
103 return true;
104 }
105
106 bool to_string(std::string &fl) const override {
107 std::ostringstream file;
108 file << " " << std::setw(4) << NUP
109 << " " << std::setw(6) << IDPRUP
110 << " " << std::setw(14) << XWGTUP
111 << " " << std::setw(14) << SCALUP
112 << " " << std::setw(14) << AQEDUP
113 << " " << std::setw(14) << AQCDUP << "\n";
114 for ( int i = 0; i < NUP; ++i )
115 file << " " << std::setw(8) << IDUP[i]
116 << " " << std::setw(2) << ISTUP[i]
117 << " " << std::setw(4) << MOTHUP[i].first
118 << " " << std::setw(4) << MOTHUP[i].second
119 << " " << std::setw(4) << ICOLUP[i].first
120 << " " << std::setw(4) << ICOLUP[i].second
121 << " " << std::setw(14) << PUP[i][0]
122 << " " << std::setw(14) << PUP[i][1]
123 << " " << std::setw(14) << PUP[i][2]
124 << " " << std::setw(14) << PUP[i][3]
125 << " " << std::setw(14) << PUP[i][4]
126 << " " << std::setw(1) << VTIMUP[i]
127 << " " << std::setw(1) << SPINUP[i] << std::endl;
128 fl+=file.str();
129 return true;
130 }
131
132 void resize() {
133 IDUP.resize(NUP);
134 ISTUP.resize(NUP);
135 MOTHUP.resize(NUP);
136 ICOLUP.resize(NUP);
137 PUP.resize(NUP, std::vector<double>(5));
138 VTIMUP.resize(NUP);
139 SPINUP.resize(NUP);
140 }
141
142 int NUP=0;
143 int IDPRUP=0;
144 double XWGTUP=0;
145 double SCALUP=0;
146 double AQEDUP=0;
147 double AQCDUP=0;
148 std::vector<long> IDUP{};
149 std::vector<int> ISTUP{};
150 std::vector< std::pair<int,int> > MOTHUP{};
151 std::vector< std::pair<int,int> > ICOLUP{};
152 std::vector< std::vector<double> > PUP{};
153 std::vector<double> VTIMUP{};
154 std::vector<double> SPINUP{};
155};
156
157
158class GenEventBarcodes : public HepMC3::Attribute
159{
160public:
161 virtual bool from_string(const std::string & /*att*/) override { return false; }
162 virtual bool to_string(std::string &att) const override
163 {
164 att = "GenEventBarcodes";
165 return true;
166 }
167 ConstGenVertexPtr barcode_to_vertex (int id) const {
168 auto it = m_vertexBC.find (id);
169 if (it != m_vertexBC.end()) return it->second;
170 return nullptr;
171 }
173 auto it = m_vertexBC.find (id);
174 if (it != m_vertexBC.end()) return it->second;
175 return nullptr;
176 }
178 auto it = m_particleBC.find (id);
179 if (it != m_particleBC.end()) return it->second;
180 return nullptr;
181 }
183 auto it = m_particleBC.find (id);
184 if (it != m_particleBC.end()) return it->second;
185 return nullptr;
186 }
187
188 void add (GenVertexPtr p) {
189 if (!p) return;
190 auto barcode = p->attribute<HepMC3::IntAttribute>("barcode");
191 if (barcode) {
192 m_vertexBC[barcode->value()] = std::move(p);
193 }
194 }
195
196 void remove (GenVertexPtr p) {
197 if (!p) return;
198 auto barcode = p->attribute<HepMC3::IntAttribute>("barcode");
199 if (barcode) {
200 auto it = m_vertexBC.find (barcode->value());
201 if (it != m_vertexBC.end()) {
202 m_vertexBC.erase (it);
203 }
204 }
205 }
206
207 void add (GenParticlePtr p) {
208 if (!p) return;
209 auto barcode = p->attribute<HepMC3::IntAttribute>("barcode");
210 if (barcode) {
211 m_particleBC[barcode->value()] = std::move(p);
212 }
213 }
214
215 void remove (GenParticlePtr p) {
216 if (!p) return;
217 auto barcode = p->attribute<HepMC3::IntAttribute>("barcode");
218 if (barcode) {
219 auto it = m_particleBC.find (barcode->value());
220 if (it != m_particleBC.end()) {
221 m_particleBC.erase (it);
222 }
223 }
224 }
225
226 std::map<int, ConstGenVertexPtr> barcode_to_vertex_map() const {
227 std::map<int, ConstGenVertexPtr> ret;
228 for (const auto &bcvertpair: m_vertexBC)
229 ret.insert({bcvertpair.first,std::const_pointer_cast<const HepMC3::GenVertex>(bcvertpair.second)});
230 return ret;
231 }
232 std::map<int, ConstGenParticlePtr> barcode_to_particle_map() const {
233 std::map<int, ConstGenParticlePtr> ret;
234 for (const auto &bcpartpair: m_particleBC)
235 ret.insert({bcpartpair.first,std::const_pointer_cast<const HepMC3::GenParticle>(bcpartpair.second)});
236 return ret;
237 }
238 std::map<int,int> id_to_barcode_map() const {
239 std::map<int, int> ret;
240 for (const auto &bcvertpair: m_vertexBC) ret.insert({bcvertpair.second->id(),bcvertpair.first});
241 for (const auto &bcpartpair: m_particleBC) ret.insert({bcpartpair.second->id(),bcpartpair.first});
242 return ret;
243 }
244
245
246 void fillAttribute(GenEvent* e) {
247 const auto eventAttributes = e->attributes(); // this makes a copy
248 const auto barcodeAttributeIt = eventAttributes.find("barcode");
249 const bool hasBarcodeAttribute = barcodeAttributeIt != eventAttributes.end();
250
251 const auto &particles = e->particles();
252 for (size_t i = 1; i <= particles.size(); i++) {
253 if (hasBarcodeAttribute && barcodeAttributeIt->second.count(i) > 0) {
254 const auto &ptr = barcodeAttributeIt->second.at(i);
255 if (ptr->is_parsed()) {
256 m_particleBC[static_cast<HepMC3::IntAttribute*>(ptr.get())->value()] = ptr->particle();
257 }
258 else {
259 m_particleBC[std::atoi(ptr->unparsed_string().c_str())] = ptr->particle();
260 }
261 } else {
262 m_particleBC[i] = particles[i-1];
263 }
264 }
265 const auto &vertices = e->vertices();
266 for (size_t i = 1; i <= vertices.size(); i++) {
267 if (hasBarcodeAttribute && barcodeAttributeIt->second.count(-i) > 0) {
268 const auto &ptr = barcodeAttributeIt->second.at(-i);
269 if (ptr->is_parsed()) {
270 m_vertexBC[static_cast<HepMC3::IntAttribute*>(ptr.get())->value()] = ptr->vertex();
271 }
272 else {
273 m_vertexBC[std::atoi(ptr->unparsed_string().c_str())] = ptr->vertex();
274 }
275 } else {
276 m_vertexBC[-i] = vertices[i-1];
277 }
278 }
279 set_is_parsed(true);
280 }
281
282private:
283 std::unordered_map<int, GenVertexPtr> m_vertexBC;
284 std::unordered_map<int, GenParticlePtr> m_particleBC;
285};
286
292class ShadowParticle : public HepMC3::Attribute {
293public:
294
296 ShadowParticle() {}
297
299 ShadowParticle(ConstGenParticlePtr p)
300 : Attribute(), m_shadow(p) {}
301
303 virtual bool from_string(const std::string &) override { return false; }
304
306 virtual bool to_string(std::string &att) const override
307 {
308 att = "ShadowParticle";
309 return true;
310 }
312 ConstGenParticlePtr value() const {
313 return m_shadow;
314 }
315
316private:
317
318 ConstGenParticlePtr m_shadow;
319};
320
321inline bool set_ll_event_number(HepMC3::GenEvent* e, long long int num){
322 e->add_attribute("long_long_event_number", std::make_shared<HepMC3::LongLongAttribute>(num));
323 return true;
324}
325inline long long int get_ll_event_number(const HepMC3::GenEvent* e){
326 auto at = e->attribute<HepMC3::LongLongAttribute>("long_long_event_number");
327 return at?at->value():e->event_number();
328}
329
330inline std::map<std::string, std::size_t> weights_map(const HepMC3::GenEvent* e) {
331 std::map<std::string, std::size_t> ret;
332 auto run = e->run_info();
333 if (!run) return ret;
334 std::vector<std::string> names = run->weight_names();
335 for (const auto& name: names) ret[name] = run->weight_index(name);
336 return ret;
337}
338
339inline std::vector<HepMC3::GenParticlePtr>::const_iterator begin(HepMC3::GenEvent& e) { return e.particles().begin(); }
340inline std::vector<HepMC3::GenParticlePtr>::const_iterator end(HepMC3::GenEvent& e) { return e.particles().end(); }
341inline std::vector<HepMC3::ConstGenParticlePtr>::const_iterator begin(const HepMC3::GenEvent& e) { return e.particles().begin(); }
342inline std::vector<HepMC3::ConstGenParticlePtr>::const_iterator end(const HepMC3::GenEvent& e) { return e.particles().end(); }
343
344inline GenEvent* newGenEvent(const int signal_process_id, const int event_number ) { // TODO Update event_number to long long int?
345 GenEvent* e= new GenEvent();
346 std::shared_ptr<HepMC3::IntAttribute> signal_process_id_A = std::make_shared<HepMC3::IntAttribute>(signal_process_id);
347 e->add_attribute("signal_process_id",signal_process_id_A);
348 e->add_attribute("barcodes", std::make_shared<GenEventBarcodes>());
349 e->set_event_number(event_number);
350 return e;
351}
352
353inline GenEvent* copyemptyGenEvent(const GenEvent* inEvt) {
354 GenEvent* e= new GenEvent();
355 e->set_event_number(inEvt->event_number());
356 e->weights()=inEvt->weights();
357 auto a_mpi = inEvt->attribute<HepMC3::IntAttribute>("mpi");
358 if (a_mpi) e->add_attribute("mpi",std::make_shared<HepMC3::IntAttribute>(*a_mpi));
359 auto a_signal_process_id = inEvt->attribute<HepMC3::IntAttribute>("signal_process_id");
360 if (a_signal_process_id) e->add_attribute("signal_process_id",std::make_shared<HepMC3::IntAttribute>(*a_signal_process_id));
361 auto a_event_scale = inEvt->attribute<HepMC3::DoubleAttribute>("event_scale");
362 if (a_event_scale) e->add_attribute("event_scale",std::make_shared<HepMC3::DoubleAttribute>(*a_event_scale));
363 auto a_alphaQCD = inEvt->attribute<HepMC3::DoubleAttribute>("alphaQCD");
364 if (a_alphaQCD) e->add_attribute("alphaQCD",std::make_shared<HepMC3::DoubleAttribute>(*a_alphaQCD));
365 auto a_alphaQED = inEvt->attribute<HepMC3::DoubleAttribute>("alphaQED");
366 if (a_alphaQED) e->add_attribute("alphaQED",std::make_shared<HepMC3::DoubleAttribute>(*a_alphaQED));
367 auto a_pi = inEvt->pdf_info();
368 if (a_pi) e->set_pdf_info(std::make_shared<HepMC3::GenPdfInfo>(*a_pi));
369 auto a_hi = inEvt->heavy_ion();
370 if (a_hi) e->set_heavy_ion(std::make_shared<HepMC3::GenHeavyIon>(*a_hi));
371 auto a_random_states = inEvt->attribute<HepMC3::VectorLongIntAttribute>("random_states");
372 if (a_random_states) e->add_attribute("random_states",std::make_shared<HepMC3::VectorLongIntAttribute>(*a_random_states));
373 e->add_attribute("barcodes", std::make_shared<GenEventBarcodes>());
374 return e;
375}
376
377inline void fillBarcodesAttribute(GenEvent* e) {
378 auto barcodes = e->attribute<GenEventBarcodes> ("barcodes");
379 if (!barcodes) {
380 barcodes = std::make_shared<GenEventBarcodes>();
381 e->add_attribute("barcodes", barcodes);
382 }
383 // force re-parsing as calling barcodes->is_parsed() returns true here
384 barcodes->fillAttribute(e);
385}
386
387inline ConstGenVertexPtr barcode_to_vertex(const GenEvent* e, int id ) {
388 // Prefer to use optimized GenEvent barcodes attribute
389 const auto &barcodes = e->attribute<GenEventBarcodes> ("barcodes");
390 if (barcodes) {
391 ConstGenVertexPtr ptr = barcodes->barcode_to_vertex (id);
392 if (ptr) return ptr;
393 }
394 // Fallback to unoptimized GenVertex barcode attribute
395 const auto eventAttributes = e->attributes(); // this makes a copy
396 const auto barcodeAttributeIt = eventAttributes.find("barcode");
397 const bool hasBarcodeAttribute = barcodeAttributeIt != eventAttributes.end();
398
399 const auto &vertices = e->vertices();
400 if (hasBarcodeAttribute) {
401 for (size_t i = 1; i <= vertices.size(); i++) {
402 const auto &ptrIt = barcodeAttributeIt->second.find(-i);
403 if (ptrIt != barcodeAttributeIt->second.end()) {
404 const auto &ptr = ptrIt->second;
405 if (ptr->is_parsed()) {
406 if (id == static_cast<HepMC3::IntAttribute*>(ptr.get())->value()) {
407 return ptr->vertex();
408 }
409 }
410 else {
411 if (id == std::atoi(ptr->unparsed_string().c_str())) {
412 return ptr->vertex();
413 }
414 }
415 }
416 }
417 }
418 // No barcodes attribute, so assume that we are passing the id member variable instead of a barcode
419 if (-id > 0 && -id <= static_cast<int>(vertices.size())) {
420 if (!vertices[-id-1]->attribute<HepMC3::IntAttribute>("barcode")) {
421 return vertices[-id-1];
422 }
423 }
424 return HepMC3::ConstGenVertexPtr();
425}
426
427inline ConstGenParticlePtr barcode_to_particle(const GenEvent* e, int id ) {
428 // Prefer to use optimized GenEvent barcodes attribute
429 const auto &barcodes = e->attribute<GenEventBarcodes> ("barcodes");
430 if (barcodes) {
431 ConstGenParticlePtr ptr = barcodes->barcode_to_particle (id);
432 if (ptr) return ptr;
433 }
434 // Fallback to unoptimized GenParticle barcode attribute
435 const auto eventAttributes = e->attributes(); // this makes a copy
436 const auto barcodeAttributeIt = eventAttributes.find("barcode");
437 const bool hasBarcodeAttribute = barcodeAttributeIt != eventAttributes.end();
438
439 const auto &particles = e->particles();
440 if (hasBarcodeAttribute) {
441 for (size_t i = 1; i <= particles.size(); i++) {
442 const auto &ptrIt = barcodeAttributeIt->second.find(i);
443 if (ptrIt != barcodeAttributeIt->second.end()) {
444 const auto &ptr = ptrIt->second;
445 if (ptr->is_parsed()) {
446 if (id == static_cast<HepMC3::IntAttribute*>(ptr.get())->value()) {
447 return ptr->particle();
448 }
449 }
450 else {
451 if (id == std::atoi(ptr->unparsed_string().c_str())) {
452 return ptr->particle();
453 }
454 }
455 }
456 }
457 }
458 // No barcodes attribute, so assume that we are passing the id member variable instead of a barcode
459 if (id > 0 && id <= static_cast<int>(particles.size())) {
460 if (!particles[id-1]->attribute<HepMC3::IntAttribute>("barcode")) {
461 return particles[id-1];
462 }
463 }
464 return HepMC3::ConstGenParticlePtr();
465}
466
467inline GenVertexPtr barcode_to_vertex(GenEvent* e, int id ) {
468 // Prefer to use optimized GenEvent barcodes attribute
469 const auto &barcodes = e->attribute<GenEventBarcodes> ("barcodes");
470 if (barcodes) {
471 GenVertexPtr ptr = barcodes->barcode_to_vertex (id);
472 if (ptr) return ptr;
473 }
474 // Fallback to unoptimized GenVertex barcode attribute
475 const auto eventAttributes = e->attributes(); // this makes a copy
476 const auto barcodeAttributeIt = eventAttributes.find("barcode");
477 const bool hasBarcodeAttribute = barcodeAttributeIt != eventAttributes.end();
478
479 const auto &vertices = e->vertices();
480 if (hasBarcodeAttribute) {
481 for (size_t i = 1; i <= vertices.size(); i++) {
482 const auto &ptrIt = barcodeAttributeIt->second.find(-i);
483 if (ptrIt != barcodeAttributeIt->second.end()) {
484 const auto &ptr = ptrIt->second;
485 if (ptr->is_parsed()) {
486 if (id == static_cast<HepMC3::IntAttribute*>(ptr.get())->value()) {
487 return ptr->vertex();
488 }
489 }
490 else {
491 if (id == std::atoi(ptr->unparsed_string().c_str())) {
492 return ptr->vertex();
493 }
494 }
495 }
496 }
497 }
498 // No barcodes attribute, so assume that we are passing the id member variable instead of a barcode
499 if (-id > 0 && -id <= static_cast<int>(vertices.size())) {
500 if (!vertices[-id-1]->attribute<HepMC3::IntAttribute>("barcode")) {
501 return vertices[-id-1];
502 }
503 }
504 return HepMC3::GenVertexPtr();
505}
506
507inline GenParticlePtr barcode_to_particle(GenEvent* e, int id ) {
508 // Prefer to use optimized GenEvent barcodes attribute
509 const auto &barcodes = e->attribute<GenEventBarcodes> ("barcodes");
510 if (barcodes) {
511 GenParticlePtr ptr = barcodes->barcode_to_particle (id);
512 if (ptr) return ptr;
513 }
514 // Fallback to unoptimized GenParticle barcode attribute
515 const auto eventAttributes = e->attributes(); // this makes a copy
516 const auto barcodeAttributeIt = eventAttributes.find("barcode");
517 const bool hasBarcodeAttribute = barcodeAttributeIt != eventAttributes.end();
518
519 const auto &particles = e->particles();
520 if (hasBarcodeAttribute) {
521 for (size_t i = 1; i <= particles.size(); i++) {
522 const auto &ptrIt = barcodeAttributeIt->second.find(i);
523 if (ptrIt != barcodeAttributeIt->second.end()) {
524 const auto &ptr = ptrIt->second;
525 if (ptr->is_parsed()) {
526 if (id == static_cast<HepMC3::IntAttribute*>(ptr.get())->value()) {
527 return ptr->particle();
528 }
529 }
530 else {
531 if (id == std::atoi(ptr->unparsed_string().c_str())) {
532 return ptr->particle();
533 }
534 }
535 }
536 }
537 }
538 // No barcodes attribute, so assume that we are passing the id member variable instead of a barcode
539 if (id > 0 && id <= static_cast<int>(particles.size())) {
540 if (!particles[id-1]->attribute<HepMC3::IntAttribute>("barcode")) {
541 return particles[id-1];
542 }
543 }
544 return HepMC3::GenParticlePtr();
545}
546
547inline int mpi(const GenEvent& evt) {
548 std::shared_ptr<HepMC3::IntAttribute> A_mpi=evt.attribute<HepMC3::IntAttribute>("mpi");
549 return A_mpi?(A_mpi->value()):0;
550}
551inline int mpi(const GenEvent* evt) {
552 std::shared_ptr<HepMC3::IntAttribute> A_mpi=evt->attribute<HepMC3::IntAttribute>("mpi");
553 return A_mpi?(A_mpi->value()):0;
554}
555
556inline int signal_process_id(const GenEvent& evt) {
557 std::shared_ptr<HepMC3::IntAttribute> A_signal_process_id=evt.attribute<HepMC3::IntAttribute>("signal_process_id");
558 return A_signal_process_id?(A_signal_process_id->value()):0;
559}
560inline int signal_process_id(const GenEvent* evt) {
561 std::shared_ptr<HepMC3::IntAttribute> A_signal_process_id=evt->attribute<HepMC3::IntAttribute>("signal_process_id");
562 return A_signal_process_id?(A_signal_process_id->value()):0;
563}
564inline void set_signal_process_id(GenEvent* e, const int i=0) {
565 std::shared_ptr<HepMC3::IntAttribute> signal_process_id = std::make_shared<HepMC3::IntAttribute>(i);
566 e->add_attribute("signal_process_id",std::move(signal_process_id));
567}
568inline void set_mpi(GenEvent* e, const int i=0) {
569 std::shared_ptr<HepMC3::IntAttribute> mpi = std::make_shared<HepMC3::IntAttribute>(i);
570 e->add_attribute("mpi",std::move(mpi));
571}
572inline void set_random_states(GenEvent* e, std::vector<long int>& a) {
573 e->add_attribute("random_states",std::make_shared<HepMC3::VectorLongIntAttribute>(a));
574}
575template <class T> void set_signal_process_vertex(GenEvent* e, T& v) {
576 if (!v || !e) return;
577/* AV: HepMC2 adds the vertex to event */
578 e->add_vertex(v);
579 v->add_attribute("signal_process_vertex",std::make_shared<HepMC3::IntAttribute>(1));
580}
581inline ConstGenVertexPtr signal_process_vertex(const GenEvent* e) { for (auto v: e->vertices()) if (v->attribute<HepMC3::IntAttribute>("signal_process_vertex")) return v; return nullptr; }
582inline GenVertexPtr signal_process_vertex(GenEvent* e) { for (auto v: e->vertices()) if (v->attribute<HepMC3::IntAttribute>("signal_process_vertex")) return v; return nullptr; }
583inline bool valid_beam_particles(const GenEvent* e) {
584 if (!e) return false;
585 size_t nBeams = 0;
586 for (const auto& p : e->beams()) { if (p->status() == 4) ++nBeams; }
587 if (nBeams != 2) return false;
588 return true;
589}
590
591template <class T> bool suggest_barcode(T& p, int i) {
592 if (!p->parent_event()) return false;
593 auto barcodes = p->parent_event()->template attribute<GenEventBarcodes> ("barcodes");
594 if (!barcodes) {
595 barcodes = std::make_shared<GenEventBarcodes>();
596 p->parent_event()->add_attribute("barcodes", barcodes);
597 }
598 barcodes->remove(p);
599 bool ret = p->add_attribute("barcode",std::make_shared<HepMC3::IntAttribute>(i));
600 if (ret) barcodes->add(p);
601 return ret;
602}
603
604}
605
606#else
607
608#include "HepMC/GenEvent.h"
609#include "HepMC/GenVertex.h"
610#include "AtlasHepMC/GenVertex.h"
611#include "AtlasHepMC/Barcode.h"
612#include <memory>
613namespace HepMC {
614inline bool set_ll_event_number(HepMC::GenEvent* e, long long int num){
615 if (num > std::numeric_limits<int>::max()) return false;
616 e->set_event_number((int)num);
617 return true;
618}
619inline long long int get_ll_event_number(const HepMC::GenEvent* e){
620 return e->event_number();
621}
622inline GenEvent::particle_iterator begin(HepMC::GenEvent& e) { return e.particles_begin(); }
623inline GenEvent::particle_iterator end(HepMC::GenEvent& e) { return e.particles_end(); }
624inline GenEvent::particle_const_iterator begin(const HepMC::GenEvent& e) { return e.particles_begin(); }
625inline GenEvent::particle_const_iterator end(const HepMC::GenEvent& e) { return e.particles_end(); }
626inline GenEvent* newGenEvent(const int a, const int b ) { return new GenEvent(a,b); }
627inline GenVertex* signal_process_vertex(const GenEvent* e) { return e->signal_process_vertex(); }
628inline void fillBarcodesAttribute(GenEvent* ) { }
629inline GenVertex* barcode_to_vertex(const GenEvent* e, int id ) {return e->barcode_to_vertex(id);}
630inline GenParticle* barcode_to_particle(const GenEvent* e, int id ) {return e->barcode_to_particle(id);}
631inline int mpi(const GenEvent& e) {
632 return e.mpi();
633}
634inline int mpi(const GenEvent* e) {
635 return e->mpi();
636}
637inline int signal_process_id(const GenEvent& e) {
638 return e.signal_process_id();
639}
640inline int signal_process_id(const GenEvent* e) {
641 return e->signal_process_id();
642}
643inline void set_signal_process_id(GenEvent* e, const int i) {
644 e->set_signal_process_id(i);
645}
646inline void set_mpi(GenEvent* e, const int i) {
647 e->set_mpi(i);
648}
649template <class T> void set_random_states(GenEvent* e, std::vector<T> a) {
650 e->set_random_states(a);
651}
652template <class T> void set_signal_process_vertex(GenEvent* e, T v) {
653 e->set_signal_process_vertex(v);
654}
655inline GenEvent* copyemptyGenEvent(const GenEvent* inEvt) {
656 HepMC::GenEvent* outEvt = new HepMC::GenEvent( inEvt->signal_process_id(), inEvt->event_number() );
657 outEvt->set_mpi ( inEvt->mpi() );
658 outEvt->set_event_scale ( inEvt->event_scale() );
659 outEvt->set_alphaQCD ( inEvt->alphaQCD() );
660 outEvt->set_alphaQED ( inEvt->alphaQED() );
661 outEvt->weights() = inEvt->weights();
662 outEvt->set_random_states( inEvt->random_states() );
663 if ( nullptr != inEvt->heavy_ion() ) {
664 outEvt->set_heavy_ion ( *inEvt->heavy_ion() );
665 }
666 if ( nullptr != inEvt->pdf_info() ) {
667 outEvt->set_pdf_info ( *inEvt->pdf_info() );
668 }
669 return outEvt;
670}
671
672template <class T> bool suggest_barcode(T& p, int i) {return p.suggest_barcode(i);}
673template <class T> bool suggest_barcode(T* p, int i) {return p->suggest_barcode(i);}
674//Smart pointers should not be used with HepMC2. But it happens.
675template <> inline bool suggest_barcode<std::unique_ptr<HepMC::GenParticle> >(std::unique_ptr<HepMC::GenParticle>& p, int i) {return p->suggest_barcode(i);}
676
677namespace Print {
678inline void line(std::ostream& os,const GenEvent& e) {e.print(os);}
679inline void line(std::ostream& os,const GenEvent* e) {e->print(os);}
680inline void content(std::ostream& os,const GenEvent& e) {e.print(os);}
681inline void content(std::ostream& os,const GenEvent* e) {e->print(os);}
682}
683inline bool valid_beam_particles(const GenEvent* e) {return e->valid_beam_particles();}
684}
686#endif
687#endif
static std::string to_string(const std::vector< T > &v)
static Double_t a
@ GenParticle
bool add(const std::string &hname, TKey *tobj)
Definition fastadd.cxx:55
int ev
Definition globals.cxx:25
void line(std::ostream &os, const GenEvent &e)
Definition GenEvent.h:678
void content(std::ostream &os, const GenEvent &e)
Definition GenEvent.h:680
int signal_process_id(const GenEvent &e)
Definition GenEvent.h:637
void set_mpi(GenEvent *e, const int i)
Definition GenEvent.h:646
bool set_ll_event_number(HepMC::GenEvent *e, long long int num)
Definition GenEvent.h:614
void set_signal_process_vertex(GenEvent *e, T v)
Definition GenEvent.h:652
GenParticle * barcode_to_particle(const GenEvent *e, int id)
Definition GenEvent.h:630
int barcode(const T *p)
Definition Barcode.h:16
GenEvent::particle_iterator begin(HepMC::GenEvent &e)
Definition GenEvent.h:622
GenEvent::particle_iterator end(HepMC::GenEvent &e)
Definition GenEvent.h:623
void set_random_states(GenEvent *e, std::vector< T > a)
Definition GenEvent.h:649
GenVertex * barcode_to_vertex(const GenEvent *e, int id)
Definition GenEvent.h:629
long long int get_ll_event_number(const HepMC::GenEvent *e)
Definition GenEvent.h:619
HepMC::GenVertex * GenVertexPtr
Definition GenVertex.h:59
bool suggest_barcode(T &p, int i)
Definition GenEvent.h:672
void set_signal_process_id(GenEvent *e, const int i)
Definition GenEvent.h:643
GenEvent * copyemptyGenEvent(const GenEvent *inEvt)
Definition GenEvent.h:655
void fillBarcodesAttribute(GenEvent *)
Definition GenEvent.h:628
GenParticle * GenParticlePtr
Definition GenParticle.h:37
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38
GenVertex * signal_process_vertex(const GenEvent *e)
Definition GenEvent.h:627
int mpi(const GenEvent &e)
Definition GenEvent.h:631
bool valid_beam_particles(const GenEvent *e)
Definition GenEvent.h:683
const HepMC::GenVertex * ConstGenVertexPtr
Definition GenVertex.h:60
GenEvent * newGenEvent(const int a, const int b)
Definition GenEvent.h:626
l
Printing final latex table to .tex output file.
void * ptr(T *p)
Definition SGImplSvc.cxx:74
barcodes(beg, end, sz)
Definition Dumpers.py:2831
TFile * file