19#include "HepMC3/ReaderFactory.h"
20#include "HepMC3/GenEvent.h"
21#include "HepMC3/GenParticle.h"
22#include "HepMC3/GenVertex.h"
23#include "HepMC3/FourVector.h"
24#include "HepMC3/Version.h"
27#include <fairlogger/Logger.h>
28#include "FairPrimaryGenerator.h"
52 mEvent =
new HepMC3::GenEvent();
62 LOG(info) <<
"Destructing GeneratorHepMC";
69 if (not
mCmd.empty()) {
83 if (not
param.fileName.empty()) {
84 LOG(warn) <<
"The use of the key \"HepMC.fileName\" is "
85 <<
"deprecated, use \"GeneratorFileOrCmd.fileNames\" instead";
89 if (not
param.fileName.empty()) {
102 if (
param.version != 0 and
mCmd.empty()) {
103 LOG(warn) <<
"The key \"HepMC.version\" is no longer needed when "
104 <<
"reading from files. The format version of the input files "
105 <<
"are automatically deduced. However, it is mandatory when reading "
106 <<
"from a pipe containing HepMC2 data.";
115 if (not
param.fileName.empty()) {
116 LOG(warn) <<
"The use of the key \"HepMC.fileName\" is "
117 <<
"deprecated, use \"GeneratorFileOrCmd.fileNames\" instead";
121 if (not
param.fileName.empty()) {
134 if (
param.version != 0 and
mCmd.empty()) {
135 LOG(warn) <<
"The key \"HepMC.version\" is no longer needed when "
136 <<
"reading from files. The format version of the input files "
137 <<
"are automatically deduced. However, it is mandatory when reading "
138 <<
"from a pipe containing HepMC2 data.";
145 LOG(
debug) <<
"Generating an event";
148 constexpr int max_tries = 3;
160 mEvent->set_units(HepMC3::Units::GEV, HepMC3::Units::MM);
164 LOG(error) <<
"Event reading from HepMC failed ...";
167 }
while (tries < max_tries);
169 LOG(error) <<
"HepMC event gen failed (Does the file/stream have enough events)?";
178 HepMC3::GenEvent&
event = *
mEvent;
180 auto particles =
event.particles();
181 auto vertices =
event.vertices();
182 std::list<HepMC3::GenParticlePtr> toRemove;
184 LOG(
debug) <<
"HepMC events has " << particles.size()
185 <<
" particles and " << vertices.size()
186 <<
" vertices" << std::endl;
189 for (
size_t i = 0;
i < particles.size(); ++
i) {
190 auto particle = particles[
i];
197 toRemove.push_back(particle);
198 LOG(
debug) <<
" Remove " << std::setw(3) << particle->id();
200 auto endVtx = particle->end_vertex();
201 auto prdVtx = particle->production_vertex();
204 endVtx->remove_particle_in(particle);
205 LOG(
debug) <<
" end " << std::setw(3) << endVtx->id();
207 if (prdVtx and prdVtx->id() != endVtx->id()) {
208 auto outbound = endVtx->particles_out();
209 auto inbound = endVtx->particles_in();
210 LOG(
debug) <<
" prd " << std::setw(3) << prdVtx->id() <<
" "
211 << std::setw(3) << outbound.size() <<
" out "
213 << std::setw(3) << inbound.size() <<
" in ";
217 for (
auto outgoing : outbound) {
221 auto ee = outgoing->end_vertex();
222 if (not ee or ee->id() != prdVtx->id()) {
223 prdVtx->add_particle_out(outgoing);
225 LOG(
debug) <<
" " << std::setw(3) << outgoing->id();
232 for (
auto incoming : inbound) {
234 auto pp = incoming->production_vertex();
235 if (not pp or pp->id() != prdVtx->id()) {
236 prdVtx->add_particle_in(incoming);
239 LOG(
debug) <<
" " << std::setw(3) << incoming->id();
245 prdVtx->remove_particle_out(particle);
249 LOG(
debug) <<
"Selected " << nSelect <<
" particles\n"
250 <<
"Removing " << toRemove.size() <<
" particles";
251 size_t oldSize = particles.size();
252 for (
auto particle : toRemove) {
253 event.remove_particle(particle);
256 std::list<HepMC3::GenVertexPtr> remVtx;
257 for (
auto vtx :
event.vertices()) {
259 (vtx->particles_out().empty() and
260 vtx->particles_in().empty())) {
261 remVtx.push_back(vtx);
264 LOG(
debug) <<
"Removing " << remVtx.size() <<
" vertexes";
265 for (
auto vtx : remVtx) {
266 event.remove_vertex(vtx);
269 LOG(
debug) <<
"HepMC events was pruned from " << oldSize
270 <<
" particles to " <<
event.particles().size()
271 <<
" particles and " <<
event.vertices().size()
281 auto select = [](HepMC3::ConstGenParticlePtr particle) {
282 switch (particle->status()) {
297 auto particles =
mEvent->particles();
298 for (
int i = 0;
i < particles.size(); ++
i) {
301 auto particle = particles.at(
i);
302 auto momentum = particle->momentum();
303 auto vertex = particle->production_vertex()->position();
304 auto parents = particle->parents();
305 auto children = particle->children();
308 auto m1 = parents.empty() ? -1 : parents.front()->id() - 1;
309 auto m2 = parents.empty() ? -1 : parents.back()->id() - 1;
312 auto d1 = children.empty() ? -1 : children.front()->id() - 1;
313 auto d2 = children.empty() ? -1 : children.back()->id() - 1;
316 mParticles.push_back(TParticle(particle->pid(),
332 particle->status() == 1);
342template <
typename AttributeType,
typename TargetType>
344 const std::string&
name,
345 const std::shared_ptr<HepMC3::Attribute>&
a)
347 if (
auto* p =
dynamic_cast<AttributeType*
>(
a.get())) {
348 eventHeader->
putInfo<TargetType>(
name, p->value());
355 const std::string&
name,
356 const std::shared_ptr<HepMC3::Attribute>&
a)
358 using IntAttribute = HepMC3::IntAttribute;
359 using LongAttribute = HepMC3::LongAttribute;
360 using FloatAttribute = HepMC3::FloatAttribute;
361 using DoubleAttribute = HepMC3::DoubleAttribute;
362 using StringAttribute = HepMC3::StringAttribute;
363 using CharAttribute = HepMC3::CharAttribute;
364 using LongLongAttribute = HepMC3::LongLongAttribute;
365 using LongDoubleAttribute = HepMC3::LongDoubleAttribute;
366 using UIntAttribute = HepMC3::UIntAttribute;
367 using ULongAttribute = HepMC3::ULongAttribute;
368 using ULongLongAttribute = HepMC3::ULongLongAttribute;
369 using BoolAttribute = HepMC3::BoolAttribute;
371 if (putAttributeInfoImpl<IntAttribute, int>(eventHeader,
name,
a)) {
374 if (putAttributeInfoImpl<LongAttribute, int>(eventHeader,
name,
a)) {
377 if (putAttributeInfoImpl<FloatAttribute, float>(eventHeader,
name,
a)) {
380 if (putAttributeInfoImpl<DoubleAttribute, float>(eventHeader,
name,
a)) {
383 if (putAttributeInfoImpl<StringAttribute, std::string>(eventHeader,
name,
a)) {
386 if (putAttributeInfoImpl<CharAttribute, char>(eventHeader,
name,
a)) {
389 if (putAttributeInfoImpl<LongLongAttribute, int>(eventHeader,
name,
a)) {
392 if (putAttributeInfoImpl<LongDoubleAttribute, float>(eventHeader,
name,
a)) {
395 if (putAttributeInfoImpl<UIntAttribute, int>(eventHeader,
name,
a)) {
398 if (putAttributeInfoImpl<ULongAttribute, int>(eventHeader,
name,
a)) {
401 if (putAttributeInfoImpl<ULongLongAttribute, int>(eventHeader,
name,
a)) {
404 if (putAttributeInfoImpl<BoolAttribute, bool>(eventHeader,
name,
a)) {
417 eventHeader->
putInfo<std::string>(Key::generator,
"hepmc");
418 eventHeader->
putInfo<
int>(Key::generatorVersion, HEPMC3_VERSION_CODE);
420 auto xSection =
mEvent->cross_section();
421 auto pdfInfo =
mEvent->pdf_info();
422 auto hiInfo =
mEvent->heavy_ion();
430 auto attStr =
mEvent->attribute_as_string(
"GenHeavyIon");
431 if (!attStr.empty() && attStr[0] ==
'v') {
432 std::istringstream is(attStr);
436 auto hi = std::make_shared<HepMC3::GenHeavyIon>();
437 double spectNeutrons, spectProtons, eccentricity, userCentEst;
438 is >> hi->Ncoll_hard >> hi->Npart_proj >> hi->Npart_targ >> hi->Ncoll >> spectNeutrons >> spectProtons
439 >> hi->N_Nwounded_collisions >> hi->Nwounded_N_collisions >> hi->Nwounded_Nwounded_collisions >> hi->impact_parameter >> hi->event_plane_angle >> eccentricity
440 >> hi->sigma_inel_NN >> hi->centrality >> userCentEst
441 >> hi->Nspec_proj_n >> hi->Nspec_targ_n >> hi->Nspec_proj_p >> hi->Nspec_targ_p;
443 LOG(
debug) <<
"GenHeavyIon: using manual v0 parser (workaround for HepMC3 from_string bug)";
446 LOG(warn) <<
"GenHeavyIon: manual v0 parser also failed on: [" << attStr <<
"]";
454 eventHeader->
putInfo<
float>(Key::xSection, xSection->xsec());
455 eventHeader->
putInfo<
float>(Key::xSectionError, xSection->xsec_err());
456 eventHeader->
putInfo<
int>(Key::acceptedEvents,
457 xSection->get_accepted_events());
458 eventHeader->
putInfo<
int>(Key::attemptedEvents,
459 xSection->get_attempted_events());
464 for (
auto w :
mEvent->weights()) {
466 eventHeader->
putInfo<
float>(Key::weight + post,
w);
468 eventHeader->
putInfo<
float>(Key::xSection, xSection->xsec(iw));
469 eventHeader->
putInfo<
float>(Key::xSectionError, xSection->xsec_err(iw));
476 eventHeader->
putInfo<
int>(Key::pdfParton1Id, pdfInfo->parton_id[0]);
477 eventHeader->
putInfo<
int>(Key::pdfParton2Id, pdfInfo->parton_id[1]);
478 eventHeader->
putInfo<
float>(Key::pdfX1, pdfInfo->x[0]);
479 eventHeader->
putInfo<
float>(Key::pdfX2, pdfInfo->x[1]);
480 eventHeader->
putInfo<
float>(Key::pdfScale, pdfInfo->scale);
481 eventHeader->
putInfo<
float>(Key::pdfXF1, pdfInfo->xf[0]);
482 eventHeader->
putInfo<
float>(Key::pdfXF2, pdfInfo->xf[1]);
483 eventHeader->
putInfo<
int>(Key::pdfCode1, pdfInfo->pdf_id[0]);
484 eventHeader->
putInfo<
int>(Key::pdfCode2, pdfInfo->pdf_id[1]);
489 eventHeader->SetB(hiInfo->impact_parameter);
490 eventHeader->
putInfo<
float>(Key::impactParameter,
491 hiInfo->impact_parameter);
492 eventHeader->
putInfo<
int>(Key::nPart,
493 hiInfo->Npart_proj + hiInfo->Npart_targ);
494 eventHeader->
putInfo<
int>(Key::nPartProjectile, hiInfo->Npart_proj);
495 eventHeader->
putInfo<
int>(Key::nPartTarget, hiInfo->Npart_targ);
496 eventHeader->
putInfo<
int>(Key::nColl, hiInfo->Ncoll);
497 eventHeader->
putInfo<
int>(Key::nCollHard, hiInfo->Ncoll_hard);
498 eventHeader->
putInfo<
int>(Key::nCollNNWounded,
499 hiInfo->N_Nwounded_collisions);
500 eventHeader->
putInfo<
int>(Key::nCollNWoundedN,
501 hiInfo->Nwounded_N_collisions);
502 eventHeader->
putInfo<
int>(Key::nCollNWoundedNwounded,
503 hiInfo->Nwounded_Nwounded_collisions);
504 eventHeader->
putInfo<
double>(Key::planeAngle, hiInfo->event_plane_angle);
505 eventHeader->
putInfo<
float>(Key::sigmaInelNN, hiInfo->sigma_inel_NN);
506 eventHeader->
putInfo<
float>(Key::centrality, hiInfo->centrality);
507 eventHeader->
putInfo<
int>(Key::nSpecProjectileProton, hiInfo->Nspec_proj_p);
508 eventHeader->
putInfo<
int>(Key::nSpecProjectileNeutron, hiInfo->Nspec_proj_n);
509 eventHeader->
putInfo<
int>(Key::nSpecTargetProton, hiInfo->Nspec_targ_p);
510 eventHeader->
putInfo<
int>(Key::nSpecTargetNeutron, hiInfo->Nspec_targ_n);
513 for (
auto na :
mEvent->attributes()) {
514 std::string
name = na.first;
515 if (
name ==
"GenPdfInfo" ||
516 name ==
"GenCrossSection" ||
517 name ==
"GenHeavyIon") {
521 for (
auto ia : na.second) {
526 putAttributeInfo(eventHeader,
name + post, at);
536 LOG(
debug) <<
"Reseting the reader";
541 LOG(
debug) <<
"No more files to read, return false";
552 if (not
mCmd.empty()) {
558 LOG(info) <<
"Creating ASCII reader of " <<
filename;
570 LOG(info) <<
"Reader is " <<
mReader.get() <<
" " << ret;
611 if (not
mCmd.empty()) {
631 LOG(
debug) <<
"EG command line is \"" << cmd <<
"\"";
635 LOG(fatal) <<
"Failed to spawn \"" << cmd <<
"\"";
o2::monitoring::tags::Key Key
Utility functions for MC particles.
void updateHeader(o2::dataformats::MCEventHeader *eventHeader) override
void setEventsToSkip(uint64_t val)
Bool_t importParticles() override
~GeneratorHepMC() override
void pruneEvent(Select select)
Bool_t generateEvent() override
HepMC3::GenEvent * mEvent
std::shared_ptr< HepMC3::Reader > mReader
void setup(const GeneratorFileOrCmdParam ¶m0, const GeneratorHepMCParam ¶m, const conf::SimConfig &config)
std::string mInterfaceName
std::vector< TParticle > mParticles
static void encodeParticleStatusAndTracking(TParticle &particle, bool wanttracking=true)
GLuint const GLchar * name
GLboolean GLboolean GLboolean GLboolean a
GLubyte GLubyte GLubyte GLubyte w
std::vector< InputSpec > select(char const *matcher="")
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::string to_string(gsl::span< T, Size > span)
virtual bool ensureFiles()
virtual bool terminateCmd()
void setFileNames(const std::string &filenames)
void setup(const GeneratorFileOrCmdParam ¶m, const conf::SimConfig &config)
std::list< std::string > mFileNames
virtual std::string makeCmdLine() const
virtual bool removeTemp() const
virtual bool executeCmdLine(const std::string &cmd)
virtual bool makeTemp(const bool &)
virtual bool makeFifo() const
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"