Loading [MathJax]/extensions/tex2jax.js
ATLAS Offline Software
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Pythia8_i.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4 #include "Pythia8_i/Pythia8_i.h"
8 
10 #include <boost/algorithm/string.hpp>
11 #include <boost/lexical_cast.hpp>
12 
13 // calls to fortran routines
16 
17 #include <sstream>
18 // For limits
19 #include <limits>
20 
21 
22 namespace {
23  // Name of random number stream
24  std::string s_pythia_stream{"PYTHIA8_INIT"};
25 }
26 
27 
28 // fix Pythia8 shower weights change in conventions
29 #define PYTHIA8_NWEIGHTS nWeights
30 #define PYTHIA8_WEIGHT weight
31 #define PYTHIA8_WLABEL weightLabel
32 #define PYTHIA8_CONVERSION 1.0
33 
34 #ifdef PYTHIA_VERSION_INTEGER
35  #undef PYTHIA8_NWEIGHTS
36  #undef PYTHIA8_WEIGHT
37  #undef PYTHIA8_WLABEL
38  #if PYTHIA_VERSION_INTEGER > 8303
39  #define PYTHIA8_NWEIGHTS nWeightGroups
40  #else
41  #define PYTHIA8_NWEIGHTS nVariationGroups
42  #endif
43  #define PYTHIA8_WEIGHT getGroupWeight
44  #define PYTHIA8_WLABEL getGroupName
45 #endif
46 
51 namespace {
52 
53 std::string py8version()
54 {
55  static const std::string incdir (PY8INCLUDE_DIR);
56  std::string::size_type pos = incdir.find ("/pythia8/");
57  if (pos == std::string::npos) return "";
58  pos += 9;
59  std::string::size_type pos2 = incdir.find ("/", pos);
60  if (pos2 == std::string::npos) pos2 = incdir.size();
61  return incdir.substr (pos, pos2-pos);
62 }
63 
64 
65 }
66 
68 Pythia8_i::Pythia8_i(const std::string &name, ISvcLocator *pSvcLocator)
69 : GenModule(name, pSvcLocator)
70 {
71  m_particleIDs["PROTON"] = PROTON;
72  m_particleIDs["ANTIPROTON"] = ANTIPROTON;
73  m_particleIDs["ELECTRON"] = ELECTRON;
74  m_particleIDs["POSITRON"] = POSITRON;
75  m_particleIDs["NEUTRON"] = NEUTRON;
76  m_particleIDs["ANTINEUTRON"] = ANTINEUTRON;
77  m_particleIDs["MUON"] = MUON;
78  m_particleIDs["ANTIMUON"] = ANTIMUON;
79  m_particleIDs["LEAD"] = LEAD;
80  m_particleIDs["OXYGEN"] = OXYGEN;
81 
82 }
83 
86 
87  #ifndef PYTHIA8_3SERIES
88 // if(m_userHookPtr != 0) delete m_userHookPtr;
89 
90 // for(UserHooksPtrType ptr: m_userHooksPtrs){
91 // delete ptr;
92 // }
93  #endif
94 }
95 
98 
99  ATH_MSG_DEBUG("Pythia8_i from genInitialize()");
100 
101  ATH_MSG_INFO("XML Path is " + xmlpath());
102  m_pythia = std::make_unique<Pythia8::Pythia> (xmlpath());
103 #ifdef HEPMC3
104  m_runinfo = std::make_shared<HepMC3::GenRunInfo>();
106  struct HepMC3::GenRunInfo::ToolInfo generator={std::string("Pythia8"),py8version(),std::string("Used generator")};
107  m_runinfo->tools().push_back(generator);
108 #endif
109 
110  bool canInit = true;
111 
112  m_version = m_pythia->settings.parm("Pythia:versionNumber");
113 
114  s_pythia_stream = "PYTHIA8_INIT";
115 
116  // By default add "nominal" to the list of shower weight names
118  m_showerWeightNames.insert(m_showerWeightNames.begin(), "nominal");
119 
120  // We do explicitly set tune 4C, since it is the starting point for many other tunes
121  // Tune 4C for pp collisions
122  m_pythia->readString("Tune:pp = 5");
123 
124  // also use CTEQ6L1 from LHAPDF 6 by default
125  // can be over-written using JO
126  m_pythia->readString("PDF:pSet= LHAPDF6:cteq6l1");
127 
128  // switch off verbose event print out
129  m_pythia->readString("Next:numberShowEvent = 0");
130 
131  // Add flag to switch off from JO the Pythia8ToHepMC::print_inconsistency internal variable
132  m_pythia->settings.addFlag("AthenaPythia8ToHepMC:print_inconsistency",true);
133 
134  // Add UserHooks first because these potentially add new settings that must exist prior to parsing commands
135 
136  for(const auto &hook: m_userHooks){
137  ATH_MSG_INFO("Adding user hook " + hook + ".");
138  bool canSetHook = true;
139  if (hook == "SuppressSmallPT") {
141  canSetHook=m_pythia->setUserHooksPtr(PYTHIA8_PTRWRAP(m_SuppressSmallPT));
142 }
143  else {
145  canSetHook = m_pythia->addUserHooksPtr(m_userHooksPtrs.back());
146 }
147  if(!canSetHook){
148  ATH_MSG_ERROR("Unable to set requested user hook.");
149  ATH_MSG_ERROR("Pythia 8 initialisation will FAIL!");
150  canInit=false;
151  }
152  }
153 
154  for(const std::pair<const std::string, double> &param : Pythia8_UserHooks::UserHooksFactory::userSettings<double>()){
155  m_pythia->settings.addParm(param.first, param.second, false, false, 0., 0.);
156  }
157 
158  for(const std::pair<const std::string, int> &param : Pythia8_UserHooks::UserHooksFactory::userSettings<int>()){
159  m_pythia->settings.addMode(param.first, param.second, false, false, 0., 0.);
160  }
161 
162  for(const std::pair<const std::string, bool> &param : Pythia8_UserHooks::UserHooksFactory::userSettings<bool>()){
163  m_pythia->settings.addFlag(param.first, param.second);
164  }
165 
166  for(const std::pair<const std::string, std::string> &param : Pythia8_UserHooks::UserHooksFactory::userSettings<std::string>()){
167  m_pythia->settings.addWord(param.first, param.second);
168  }
169 
170  if( ! m_athenaTool.empty() ){
171  if(m_athenaTool.retrieve().isFailure()){
172  ATH_MSG_ERROR("Unable to retrieve Athena Tool for custom Pythia processing");
173  return StatusCode::FAILURE;
174  }
175  else {
176  StatusCode status = m_athenaTool->InitializePythiaInfo(*m_pythia);
177  if(status != StatusCode::SUCCESS) return status;
178  }
179  }
180 
181  // Now apply the settings from the JO
182  for(const std::string &cmd : m_commands){
183 
184  if(cmd.compare("")==0) continue;
185  else if (cmd.find("Beams:id") != std::string::npos ) {
186  ATH_MSG_ERROR("With command '" << cmd << "' you are trying to set a beam different from p: please use the Beam1/Beam2 properties instead:");
187  ATH_MSG_ERROR(" example: genSeq.Pythia8.Beam1 = 'LEAD'");
188  ATH_MSG_ERROR(" genSeq.Pythia8.Beam2 = 'ANTINEUTRON'");
189  return StatusCode::FAILURE;
190  }
191  else if (cmd.find("Beams:frameType") != std::string::npos ) {
193  ATH_MSG_WARNING(" Found an explicit 'Beams:frameType' command: this will override transform beams momenta/energy parameters, regardless of its requested value. ");
194  }
195 
196  bool read = m_pythia->readString(cmd);
197 
198  if(!read){
199  ATH_MSG_ERROR("Pythia could not understand the command '"<< cmd<<"'");
200  return StatusCode::FAILURE;
201  }
202  }
203 
206 
207  ATH_MSG_INFO("Beam1 = "<<beam1);
208  ATH_MSG_INFO("Beam2 = "<<beam2);
209 
210  if(beam1 == 0 || beam2 == 0){
211  ATH_MSG_ERROR("Invalid beam particle!");
212  return StatusCode::FAILURE;
213  }
214 
215 
216 
217  if(m_useRndmGenSvc){
218 
219  ATH_MSG_INFO(" !!!!!!!!!!!! WARNING ON PYTHIA RANDOM NUMBERS !!!!!!!!!! ");
220  ATH_MSG_INFO(" THE ATHENA SERVICE AthRNGenSvc IS USED.");
221  ATH_MSG_INFO(" !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! ");
222 
223  m_atlasRndmEngine = std::make_shared<customRndm>();
224  CLHEP::HepRandomEngine* rndmEngine = getRandomEngineDuringInitialize(s_pythia_stream, m_randomSeed, m_dsid); // NOT THREAD-SAFE
225  m_atlasRndmEngine->init(rndmEngine);
226 #if PYTHIA_VERSION_INTEGER >= 8310
227  m_pythia->setRndmEnginePtr(m_atlasRndmEngine);
228 #else
229  m_pythia->setRndmEnginePtr(m_atlasRndmEngine.get());
230 #endif
231  s_pythia_stream = "PYTHIA8";
232  }else{
233  ATH_MSG_INFO(" !!!!!!!!!!!! WARNING ON PYTHIA RANDOM NUMBERS !!!!!!!!!! ");
234  ATH_MSG_INFO(" THE STANDARD PYTHIA8 RANDOM NUMBER SERVICE IS USED.");
235  ATH_MSG_INFO(" THE ATHENA SERVICE AthRNGSvc IS ***NOT*** USED.");
236  ATH_MSG_INFO(" !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! ");
237  }
238 
239  if(m_userProcess != ""){
241  }
242 
243 
244 
245  if(m_userResonances.value() != ""){
246 
247  std::vector<std::string> resonanceArgs;
248 
249  boost::split(resonanceArgs, m_userResonances.value(), boost::is_any_of(":"));
250  if(resonanceArgs.size() != 2){
251  ATH_MSG_ERROR("Cannot Understand UserResonance job option!");
252  ATH_MSG_ERROR("You should specify it as a 'name:id1,id2,id3...'");
253  ATH_MSG_ERROR("Where name is the name of your UserResonance, and id1,id2,id3 are a comma separated list of PDG IDs to which it is applied");
254  canInit = false;
255  }
256  std::vector<std::string> resonanceIds;
257  boost::split(resonanceIds, resonanceArgs.back(), boost::is_any_of(","));
258  if(resonanceIds.size()==0){
259  ATH_MSG_ERROR("You did not specifiy any PDG ids to which your user resonance width should be applied!");
260  ATH_MSG_ERROR("You should specify a list as 'name:id1,id2,id3...'");
261  ATH_MSG_ERROR("Where name is the name of your UserResonance, and id1,id2,id3 are a comma separated list of PDG IDs to which it is applied");
262  canInit=false;
263  }
264 
265  for(std::vector<std::string>::const_iterator sId = resonanceIds.begin();
266  sId != resonanceIds.end(); ++sId){
267  int idResIn = boost::lexical_cast<int>(*sId);
268  m_userResonancePtrs.push_back(Pythia8_UserResonance::UserResonanceFactory::create(resonanceArgs.front(), idResIn));
269  }
270 
271  for(std::shared_ptr<Pythia8::ResonanceWidths>& resonance : m_userResonancePtrs) {
272 #if PYTHIA_VERSION_INTEGER >= 8310
273  m_pythia->setResonancePtr(resonance);
274 #else
275  m_pythia->setResonancePtr(resonance.get());
276 #endif
277  }
278 
279  }
280 
281  if(m_particleDataFile != "") {
282  if(!m_pythia->particleData.reInit(m_particleDataFile, true)){
283  ATH_MSG_ERROR("Unable to read requested particle data table: " + m_particleDataFile + " !!");
284  ATH_MSG_ERROR("Pythia 8 initialisation will FAIL!");
285  canInit = false;
286  }
287  }
288 
289  if(m_lheFile != ""){
290  if(m_procPtr){
291  ATH_MSG_ERROR("Both LHE file and user process have been specified");
292  ATH_MSG_ERROR("LHE input does not make sense with a user process!");
293  canInit = false;
294  }
295 
296  canInit = canInit && m_pythia->readString("Beams:frameType = 4");
297  canInit = canInit && m_pythia->readString("Beams:LHEF = " + m_lheFile);
298  if(!canInit){
299  ATH_MSG_ERROR("Unable to read requested LHE file: " + m_lheFile + " !");
300  ATH_MSG_ERROR("Pythia 8 initialisation will FAIL!");
301  }
302  }else{
304  canInit = canInit && m_pythia->readString("Beams:frameType = 1");
305  canInit = canInit && m_pythia->readString("Beams:eCM = " + std::to_string(m_collisionEnergy));
306  }
307  canInit = canInit && m_pythia->readString("Beams:idA = " + std::to_string(beam1));
308  canInit = canInit && m_pythia->readString("Beams:idB = " + std::to_string(beam2));
309  }
310 
311  if(m_procPtr){
312 #if PYTHIA_VERSION_INTEGER >= 8310
313  if(!m_pythia->setSigmaPtr(m_procPtr))
314 #else
315  if(!m_pythia->setSigmaPtr(m_procPtr.get()))
316 #endif
317  {
318  ATH_MSG_ERROR("Unable to set requested user process: " + m_userProcess + " !!");
319  ATH_MSG_ERROR("Pythia 8 initialisation will FAIL!");
320  canInit = false;
321  }
322  }
323 
324  StatusCode returnCode = StatusCode::SUCCESS;
325 
326  m_pythia->particleData.listXML(m_outputParticleDataFile.value().substr(0,m_outputParticleDataFile.value().find("xml"))+"orig.xml");
327  m_pythia->settings.writeFile("Settings_before.log",true);
328 
329  if(canInit){
330  canInit = m_pythia->init();
331  }
332 
333  if(!canInit){
334  returnCode = StatusCode::FAILURE;
335  ATH_MSG_ERROR(" *** Unable to initialise Pythia !! ***");
336  }
337 
338  m_pythia->particleData.listXML(m_outputParticleDataFile);
339  m_pythia->settings.writeFile("Settings_after.log",true);
340 
341  //counter for event failures;
342  m_failureCount = 0;
343 
345 
346  // Set set_print_inconsistency to Athena corresponding flag (allowing to change it from JO)
347  m_pythiaToHepMC.set_print_inconsistency( m_pythia->settings.flag("AthenaPythia8ToHepMC:print_inconsistency") );
348 
349  m_pythiaToHepMC.set_store_pdf(true);
350  m_pythiaToHepMC.set_store_weights(false);
351 
352  return returnCode;
353 }
354 
357 
358  ATH_MSG_DEBUG(">>> Pythia8_i from callGenerator");
360  //Re-seed the random number stream
361  long seeds[7];
362  const EventContext& ctx = Gaudi::Hive::currentContext();
363  ATHRNG::calculateSeedsMC21(seeds, s_pythia_stream, ctx.eventID().event_number(), m_dsid, m_randomSeed);
364  m_atlasRndmEngine->getEngine()->setSeeds(seeds, 0); // NOT THREAD-SAFE
365  m_seeds.clear();
366  m_seeds.push_back(seeds[0]);
367  m_seeds.push_back(seeds[1]);
368  }
369 
370  bool status = m_pythia->next();
371 
372  StatusCode returnCode = StatusCode::SUCCESS;
373 
374  if(!status){
375  ++m_failureCount;
377  ATH_MSG_ERROR("Exceeded the max number of consecutive event failures.");
378  returnCode = StatusCode::FAILURE;
379  }else{
380  ATH_MSG_INFO("Event generation failed - re-trying.");
381  returnCode = this->callGenerator();
382  }
383  }
384 
385  ATH_MSG_DEBUG("Now checking with Tool.empty() second time");
386  if( ! m_athenaTool.empty() ){
387  StatusCode stat = m_athenaTool->ModifyPythiaEvent(*m_pythia);
388  if(stat != StatusCode::SUCCESS) returnCode = stat;
389  }
390 
391  m_failureCount = 0;
392 
393  m_nAccepted += 1.;
394  // some CKKWL merged events have zero weight (or unfilled event).
395  // start again with such events
396 
397  double eventWeight = m_pythia->info.mergingWeight()*m_pythia->info.weight();
398 
399 
400  if(returnCode != StatusCode::FAILURE &&
401  // Here this is a double, but it may be converted into float at some point downstream
402  (std::abs(eventWeight) < std::numeric_limits<float>::min() ||
403  m_pythia->event.size() < 2)){
404 
405  returnCode = this->callGenerator();
406  } else if ( std::abs(eventWeight) < std::numeric_limits<float>::min() &&
407  std::abs(eventWeight) > std::numeric_limits<double>::min() ){
408  ATH_MSG_WARNING("Found event weight " << eventWeight << " between the float and double precision limits. Rejecting event.");
409  } else {
410  m_nMerged += eventWeight;
412 
413  // For FxFx cross section:
414  m_sigmaTotal+=m_pythia->info.weight();
415  }
416 
417  return returnCode;
418 }
419 
422 
423  ATH_MSG_DEBUG(">>> Pythia8_i from fillEvt");
424 
425  evt->set_event_number(m_internal_event_number);
426 
427  // if using "getGroupWeight" and | lhastrategy | = 4, then need to convert mb to pb ( done otherwise when calling info.weight(), [...] )
428  if( m_internal_event_number == 1 && std::abs(m_pythia->info.lhaStrategy()) == 4 ) {
430  ATH_MSG_DEBUG(" LHA strategy needs a conversion to fix Pythia8 shower weights bug(s) equal to " << m_conversion);
431  }
432 
433 
434  if(m_pythia->event.size() < 2){
435  ATH_MSG_ERROR("Something wrong with this event - it contains fewer than 2 particles!");
436  ATH_MSG_ERROR("internal event number is "<<m_internal_event_number);
437  return StatusCode::FAILURE;
438  }
439 
441 
442  // In case we're asked, save the LHE event as an attribute. This only works in HEPMC3.
443  if(m_saveLHE){
444 #ifdef HEPMC3
445  auto evtlhe = std::make_shared<HepMC::GenEvent>();
446  m_pythiaToHepMC.fill_next_event(m_pythia->process, evtlhe.get(), m_internal_event_number, &m_pythia->info, &m_pythia->settings);
447  auto extra = std::make_shared<HepMC::ShortEventAttribute>(evtlhe.get());
448  evt->add_attribute("LHERecord", extra);
449 #else
450  ATH_MSG_WARNING("LHE record saving requested, but not implemented with HEPMC2 functionality.");
451 #endif
452  }
453 
454  // in debug mode you can check whether the pdf information is stored
455  if(evt->pdf_info()){
456 #ifdef HEPMC3
457  ATH_MSG_DEBUG("PDFinfo id1:" << evt->pdf_info()->parton_id[0]);
458  ATH_MSG_DEBUG("PDFinfo id2:" << evt->pdf_info()->parton_id[1]);
459  ATH_MSG_DEBUG("PDFinfo x1:" << evt->pdf_info()->x[0]);
460  ATH_MSG_DEBUG("PDFinfo x2:" << evt->pdf_info()->x[1]);
461  ATH_MSG_DEBUG("PDFinfo scalePDF:" << evt->pdf_info()->scale);
462  ATH_MSG_DEBUG("PDFinfo pdf1:" << evt->pdf_info()->pdf_id[0]);
463  ATH_MSG_DEBUG("PDFinfo pdf2:" << evt->pdf_info()->pdf_id[1]);
464 #else
465  ATH_MSG_DEBUG("PDFinfo id1:" << evt->pdf_info()->id1());
466  ATH_MSG_DEBUG("PDFinfo id2:" << evt->pdf_info()->id2());
467  ATH_MSG_DEBUG("PDFinfo x1:" << evt->pdf_info()->x1());
468  ATH_MSG_DEBUG("PDFinfo x2:" << evt->pdf_info()->x2());
469  ATH_MSG_DEBUG("PDFinfo scalePDF:" << evt->pdf_info()->scalePDF());
470  ATH_MSG_DEBUG("PDFinfo pdf1:" << evt->pdf_info()->pdf1());
471  ATH_MSG_DEBUG("PDFinfo pdf2:" << evt->pdf_info()->pdf2());
472 #endif
473  }
474  else
475  ATH_MSG_DEBUG("No PDF information available in HepMC::GenEvent!");
476 
477  // set the randomseeds
479 
480  // fill weights
482 
483  return StatusCode::SUCCESS;
484 }
485 
488 
489  // reset
490  m_mergingWeight = m_pythia->info.mergingWeightNLO(); // first initialisation, good for the nominal weight in any case (see Merging:includeWeightInXsection).
491  m_enhanceWeight = 1.0; // better to keep the enhancement from UserHooks in a separate variable, to keep the code clear
492 
493 #ifndef PYTHIA8_304SERIES
494  // include Enhance userhook weight
495  for(const auto &hook: m_userHooksPtrs) {
496  if (hook->canEnhanceEmission()) {
497  m_enhanceWeight *= hook->getEnhancedEventWeight();
498  }
499  }
500 #endif // not PYTHIA8_304SERIES
501 
502  // DO NOT try to distinguish between phase space and merging contributions: only their product contains both, so always multiply them in order to be robust
503  double eventWeight = m_pythia->info.weight() * m_mergingWeight * m_enhanceWeight;
504 
505  std::map<std::string,double> fWeights; // temporary weight vector, needed while managing hepmc2 and 3 simultaneously
506 
507  // save as Default the usual pythia8 weight, including enhancements
508  // NB: info.weight() is built upon info.eventWeightLHEF (but includes shower "decorations"), not from the rwgt LHE vector.
509  size_t atlas_specific_weights = 1;
510  fWeights["Default"]=eventWeight;
511  if(m_internal_event_number == 1) {
512  ATH_MSG_INFO("found LHEF version "<< m_pythia->info.LHEFversion() );
513  if (m_pythia->info.LHEFversion() > 1) m_doLHE3Weights = true;
514  m_weightIDs.clear();
515  m_weightIDs.push_back("Default"); // notice that we are assigning a weightname = weightID
516  }
517 
518  // Now systematic weights: a) generator and b) shower variations
519  // merging implies reading from lhe files (while the opposite is not true, ok, but still works)
520  if (m_lheFile!="") {
521  // remember that the true merging component can be disentangled only when there are LHE input events
522 
523  // make sure to get the merging weight correctly: this is needed -only- for systematic weights (generator and parton shower variations)
524  // NB: most robust way, regardless of the specific value of Merging:includeWeightInXsection (if on, info.mergingWeight() is always 1.0 ! )
525  m_mergingWeight *= ( m_pythia->info.weight() / m_pythia->info.eventWeightLHEF ); // this is the actual merging component, regardless of conventions,
526  // needed by systematic variation weights
527  // this includes rare compensation or selection-biased weights
528  }
529 
530 
531  ATH_MSG_DEBUG("Event weights: enhancing weight, merging weight, total weight = "<<m_enhanceWeight<<", "<<m_mergingWeight<<", "<<eventWeight);
532 
533  // we already filled the "Default" weight
534  std::vector<std::string>::const_iterator id = m_weightIDs.begin()+atlas_specific_weights;
535 
536  // a) fill generator weights
537  if(m_pythia->info.getWeightsDetailedSize() != 0){
538  for(std::map<std::string, Pythia8::LHAwgt>::const_iterator wgt = m_pythia->info.rwgt->wgts.begin();
539  wgt != m_pythia->info.rwgt->wgts.end(); ++wgt){
540 
541  if(m_internal_event_number == 1){
542  m_weightIDs.push_back(wgt->first);
543  }else{
544  if(*id != wgt->first){
545  ATH_MSG_ERROR("Mismatch in LHE3 weight id. Found "<<wgt->first<<", expected "<<*id);
546  return StatusCode::FAILURE;
547  }
548  ++id;
549  }
550 
551  std::map<std::string, Pythia8::LHAweight>::const_iterator weightName = m_pythia->info.init_weights->find(wgt->first);
552  if(weightName != m_pythia->info.init_weights->end()){
553  fWeights[weightName->second.contents] = m_mergingWeight * m_enhanceWeight * wgt->second.contents;
554  }else{
555  // if no name was assigned/found, assign a weightname = weightid
556  fWeights[wgt->first] = m_mergingWeight * m_enhanceWeight * wgt->second.contents;
557  }
558 
559  }
560  }
561 
562  // b) shower weight variations
563  // always start from second shower weight: the first was saved already as "Default"
564 
565  for(int iw = 1; iw < m_pythia->info.PYTHIA8_NWEIGHTS(); ++iw){
566 
567  std::string wtName = ((int)m_showerWeightNames.size() == m_pythia->info.PYTHIA8_NWEIGHTS())? m_showerWeightNames[iw]: "ShowerWt_" +std::to_string(iw);
568 
569  fWeights[wtName] = m_mergingWeight * m_enhanceWeight * m_pythia->info.PYTHIA8_WEIGHT(iw)*m_conversion;
570 
571  if(m_internal_event_number == 1) {
572  m_weightIDs.push_back(wtName);
573  ATH_MSG_DEBUG("Shower weight name, value, conversion: "<<wtName<<", "<< fWeights[wtName] <<","<<m_conversion );
574  }
575  }
576 
577  // Now save, as a backup, also the LHEF weight without enhancements (useful for possible, rare, non-default cases)
578  // to make clear that it should not be used in analyses, save it always with a -10 factor (so that its goal is clear even if not checking its name)
579  if (m_lheFile!="")
580  {
581  fWeights["EXTRA_bare_LHE_weight"]=(-10.0)*m_pythia->info.eventWeightLHEF;
582  if(m_internal_event_number == 1) m_weightIDs.push_back("EXTRA_bare_LHE_weight");
583  }
584 
585  // Sad, but needed: create a string vector with acceptable order of weight names
586  if(m_internal_event_number == 1){
587  m_weightNames.clear();
588  for(const std::string &id : m_weightIDs){
589  if(m_doLHE3Weights) {
590  std::map<std::string, Pythia8::LHAweight>::const_iterator weight = m_pythia->info.init_weights->find(id);
591  if(weight != m_pythia->info.init_weights->end()) m_weightNames.push_back(weight->second.contents);
592  else m_weightNames.push_back(id);
593  }
594  else {
595  m_weightNames.push_back(id);
596  }
597  }
598  if (m_weightNames.size() != fWeights.size() ) {
599  ATH_MSG_ERROR("Something wrong when building list of weight names: " << m_weightNames.size() << " vs "<< fWeights.size() << ", exiting ...");
600  return StatusCode::FAILURE;
601  }
602  }
603 
604  // The following depends on the specific hepmc2/3 implementation
605 #ifdef HEPMC3
606  if (!evt->run_info()) {
607  evt->set_run_info(m_runinfo);
608  }
609  evt->run_info()->set_weight_names(m_weightNames);
610 
611  // for the first event, weight EXTRA_bare_LHE_weight is not present in evt->weights(), so we need to book a place for it by hand
612  if (m_internal_event_number == 1 && evt->run_info()->weight_names().size() == evt->weights().size()+1 ) {
613  evt->weights().push_back(1.0);
614  }
615 
616  // added conversion GeV -> MeV to ensure correct units
617  evt->set_units(HepMC3::Units::MEV, HepMC3::Units::MM);
618 
619  evt->weights().resize(fWeights.size(), 1.0);
620  for (auto w: fWeights) {
621  evt->weight(w.first)=w.second;
622  }
623 
624  auto beams=evt->beams();
625  ATH_MSG_DEBUG( " Energy of the beams " << beams[0]->momentum().e() );
626 
627 //uncomment to list HepMC events
628 // std::cout << " print::listing Pythia8 " << std::endl;
629 // HepMC3::Print::listing(std::cout, *evt);
630 
631 #else
632  evt->weights().clear();
633  for (const auto& w: m_weightNames ) {evt->weights()[w]=fWeights[w];}
634  auto beams=evt->beam_particles();
635  ATH_MSG_DEBUG( " Energy of the beams " << beams.first->momentum().e() );
636 
637 //uncomment to list HepMC events
638 // std::cout << " print::listing Pythia8 " << std::endl;
639 // evt->print();
640 #endif
641 
642  return StatusCode::SUCCESS;
643 }
644 
647 
648  ATH_MSG_INFO(">>> Pythia8_i from genFinalize");
649  m_pythia->stat();
650 
651  Pythia8::Info info = m_pythia->info;
652  double xs = info.sigmaGen(); // in mb
653 
655  const double accfactor = m_nMerged / info.nAccepted();
656  ATH_MSG_DEBUG("Multiplying cross-section by CKKWL merging acceptance of "<<m_nMerged <<"/" <<info.nAccepted() << " = " << accfactor
657  << ": " << xs << " -> " << xs*accfactor);
658  xs *= accfactor;
659  }
660 
661  if(m_doFxFxXS){
662  ATH_MSG_DEBUG("Using FxFx cross section recipe: xs = "<< m_sigmaTotal << " / " << 1e9*info.nTried());
663  xs = m_sigmaTotal / (1e9*info.nTried());
664  std::cout << "Using FxFx cross section recipe: xs = "<< m_sigmaTotal << " / " << 1e9*info.nTried() << std::endl;
665  }
666 
667  if( ! m_athenaTool.empty()){
668  double xsmod = m_athenaTool->CrossSectionScaleFactor();
669  ATH_MSG_DEBUG("Multiplying cross-section by Pythia Modifier tool factor " << xsmod );
670  xs *= xsmod;
671  }
672 
673  xs *= 1000. * 1000.;//convert to nb
674 
675  std::cout << "MetaData: cross-section (nb)= " << xs <<std::endl;
676  std::cout << "MetaData: generator= Pythia 8." << std::string(py8version()) <<std::endl;
677 
678  if(m_doLHE3Weights || m_weightIDs.size()>1 ){
679  std::cout<<"MetaData: weights = ";
680  for (const auto& w : m_weightNames ) std::cout<< w <<" | ";
681  std::cout<<std::endl;
682  }
683 
684  if (m_computeEfficiency){
685  if (info.nTried()>0) ATH_MSG_INFO("Pythia8 efficiency (nAccepted/nTried %) = " << info.nAccepted()*100./info.nTried());
686  else ATH_MSG_INFO("Pythia8 efficiency cannot be computed, nTried <=0");
687  }
688 
689 if(m_useRndmGenSvc){
690  ATH_MSG_INFO("Number of random numbers produced " << m_atlasRndmEngine->getRNCalls());
691  if (m_useReseed) ATH_MSG_INFO("Each event was reseeded ");
692  else ATH_MSG_INFO("Events were not reseeded ");
693 }
694 
695  return StatusCode::SUCCESS;
696 }
697 
700  return m_version;
701 }
702 
704 const std::string& Pythia8_i::pythia_stream() {
705  return s_pythia_stream;
706 }
707 
709 std::string Pythia8_i::xmlpath(){
710 
711 
712  std::string foundpath = "";
713 
714 // Try to find the xmldoc directory using PathResolver:
715  foundpath = PathResolverFindCalibDirectory( "Pythia8/xmldoc" );
716 
717  return foundpath;
718 }
719 
721 
722 // fix Pythia8 shower weights change in conventions
723 #ifdef PYTHIA8_NWEIGHTS
724  #undef PYTHIA8_NWEIGHTS
725  #undef PYTHIA8_WEIGHT
726  #undef PYTHIA8_WLABEL
727  #undef PYTHIA8_CONVERSION
728 #endif
grepfile.info
info
Definition: grepfile.py:38
Pythia8_i::NEUTRON
@ NEUTRON
Definition: Pythia8_i.h:120
read
IovVectorMap_t read(const Folder &theFolder, const SelectionCriterion &choice, const unsigned int limit=10)
Definition: openCoraCool.cxx:569
AllowedVariables::e
e
Definition: AsgElectronSelectorTool.cxx:37
Pythia8_i::m_maxFailures
UnsignedIntegerProperty m_maxFailures
Definition: Pythia8_i.h:94
Pythia8_i::m_enhanceWeight
double m_enhanceWeight
Definition: Pythia8_i.h:162
Pythia8_i::m_nMerged
double m_nMerged
Definition: Pythia8_i.h:135
Pythia8_i::m_saveLHE
BooleanProperty m_saveLHE
Definition: Pythia8_i.h:172
ForwardTracker::beam1
@ beam1
Definition: ForwardTrackerConstants.h:13
Pythia8_i::fillWeights
virtual StatusCode fillWeights(HepMC::GenEvent *evt)
Definition: Pythia8_i.cxx:487
Pythia8_i::MUON
@ MUON
Definition: Pythia8_i.h:120
Pythia8_i::m_version
double m_version
Definition: Pythia8_i.h:114
Pythia8_i::m_procPtr
std::shared_ptr< Pythia8::Sigma2Process > m_procPtr
Definition: Pythia8_i.h:148
ForwardTracker::beam2
@ beam2
Definition: ForwardTrackerConstants.h:13
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
Pythia8_i::xmlpath
static std::string xmlpath()
Definition: Pythia8_i.cxx:709
Pythia8_i::m_particleDataFile
StringProperty m_particleDataFile
Definition: Pythia8_i.h:158
Pythia8_i::m_weightNames
std::vector< std::string > m_weightNames
Definition: Pythia8_i.h:164
Pythia8_i::m_useReseed
BooleanProperty m_useReseed
Definition: Pythia8_i.h:99
rerun_display.cmd
string cmd
Definition: rerun_display.py:67
min
constexpr double min()
Definition: ap_fixedTest.cxx:26
Pythia8_i::m_failureCount
unsigned int m_failureCount
Definition: Pythia8_i.h:139
Pythia8_i::ELECTRON
@ ELECTRON
Definition: Pythia8_i.h:120
Pythia8_i::m_athenaTool
PublicToolHandle< IPythia8Custom > m_athenaTool
Definition: Pythia8_i.h:170
MM
@ MM
Definition: RegSelEnums.h:38
Pythia8_i::pythiaVersion
double pythiaVersion() const
Definition: Pythia8_i.cxx:699
Pythia8_i::m_userHooks
StringArrayProperty m_userHooks
Definition: Pythia8_i.h:102
LArG4FSStartPointFilter.evt
evt
Definition: LArG4FSStartPointFilter.py:42
Pythia8_i::m_pt0timesMPI
DoubleProperty m_pt0timesMPI
Definition: Pythia8_i.h:104
Pythia8_i.h
Pythia8_i::m_showerWeightNamesProp
StringArrayProperty m_showerWeightNamesProp
Definition: Pythia8_i.h:168
Pythia8_i::ANTIPROTON
@ ANTIPROTON
Definition: Pythia8_i.h:120
Pythia8_i::ANTINEUTRON
@ ANTINEUTRON
Definition: Pythia8_i.h:120
Pythia8_i::m_userHooksPtrs
std::vector< UserHooksPtrType > m_userHooksPtrs
Definition: Pythia8_i.h:150
PathResolverFindCalibDirectory
std::string PathResolverFindCalibDirectory(const std::string &logical_file_name)
Definition: PathResolver.cxx:432
Pythia8_i::m_doFxFxXS
BooleanProperty m_doFxFxXS
Definition: Pythia8_i.h:132
dqt_zlumi_pandas.weight
int weight
Definition: dqt_zlumi_pandas.py:189
Pythia8_UserResonance::UserResonanceFactory::create
static std::shared_ptr< ResonanceWidths > create(const std::string &name, int pdgid)
Call this with the name of the ResonanceWidth and PDG ID to which it will be applied e....
Definition: UserResonanceFactory.cxx:13
Pythia8_i::m_doLHE3Weights
bool m_doLHE3Weights
Definition: Pythia8_i.h:165
GenModule
Base class for common behaviour of generator interfaces.
Definition: GenModule.h:39
Pythia8_i::m_userResonancePtrs
std::vector< std::shared_ptr< Pythia8::ResonanceWidths > > m_userResonancePtrs
Definition: Pythia8_i.h:154
Pythia8_i::fillEvt
virtual StatusCode fillEvt(HepMC::GenEvent *evt)
For filling the HepMC event object.
Definition: Pythia8_i.cxx:421
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
ParticleGun_EoverP_Config.momentum
momentum
Definition: ParticleGun_EoverP_Config.py:63
Pythia8_i::m_numberAlphaS
DoubleProperty m_numberAlphaS
Definition: Pythia8_i.h:105
McEventCollection.h
Pythia8_i::m_showerWeightNames
std::vector< std::string > m_showerWeightNames
Definition: Pythia8_i.h:167
Pythia8_UserProcess::UserProcessFactory::create
static std::shared_ptr< Sigma2Process > create(const std::string &procName)
Definition: UserProcessFactory.cxx:11
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
Pythia8_i::OXYGEN
@ OXYGEN
Definition: Pythia8_i.h:120
Pythia8_i::m_atlasRndmEngine
std::shared_ptr< customRndm > m_atlasRndmEngine
Definition: Pythia8_i.h:97
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
python.LArMinBiasAlgConfig.int
int
Definition: LArMinBiasAlgConfig.py:59
PYTHIA8_CONVERSION
#define PYTHIA8_CONVERSION
Definition: Pythia8_i.cxx:32
Pythia8_i::callGenerator
virtual StatusCode callGenerator()
For calling the generator on each iteration of the event loop.
Definition: Pythia8_i.cxx:356
UserResonanceFactory.h
UserProcessFactory.h
Pythia8_i::genFinalize
virtual StatusCode genFinalize()
For finalising the generator, if required.
Definition: Pythia8_i.cxx:646
Pythia8_i::m_nAccepted
double m_nAccepted
Definition: Pythia8_i.h:134
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
xAOD::double
double
Definition: CompositeParticle_v1.cxx:159
ppToJpsimu0mu0.SuppressSmallPT
SuppressSmallPT
Definition: ppToJpsimu0mu0.py:26
beamspotman.stat
stat
Definition: beamspotman.py:266
Pythia8_i::m_beam2
StringProperty m_beam2
Definition: Pythia8_i.h:126
Pythia8_i::m_lheFile
StringProperty m_lheFile
Definition: Pythia8_i.h:129
Pythia8_i::m_sigmaTotal
double m_sigmaTotal
Definition: Pythia8_i.h:136
Pythia8_i::m_SuppressSmallPT
Pythia8::SuppressSmallPT * m_SuppressSmallPT
Definition: Pythia8_i.h:176
python.handimod.extra
int extra
Definition: handimod.py:522
PathResolver.h
id
SG::auxid_t id
Definition: Control/AthContainers/Root/debug.cxx:227
Pythia8_i::m_weightIDs
std::vector< std::string > m_weightIDs
Definition: Pythia8_i.h:163
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:228
MagicNumbers.h
ActsTrk::to_string
std::string to_string(const DetectorType &type)
Definition: GeometryDefs.h:34
Pythia8_i::~Pythia8_i
~Pythia8_i()
Definition: Pythia8_i.cxx:85
Pythia8_i::LEAD
@ LEAD
Definition: Pythia8_i.h:120
PYTHIA8_PTRWRAP
#define PYTHIA8_PTRWRAP(A)
Definition: UserHooksFactory.h:14
RNGWrapper.h
Pythia8_i::PDGID
PDGID
Definition: Pythia8_i.h:120
python.LumiBlobConversion.pos
pos
Definition: LumiBlobConversion.py:18
Pythia8_i::Pythia8_i
Pythia8_i(const std::string &name, ISvcLocator *pSvcLocator)
Definition: Pythia8_i.cxx:68
Pythia8_i::m_userProcess
StringProperty m_userProcess
Definition: Pythia8_i.h:145
Pythia8_i::m_seeds
std::vector< long int > m_seeds
Definition: Pythia8_i.h:143
Pythia8_i::m_particleIDs
std::map< std::string, PDGID > m_particleIDs
Definition: Pythia8_i.h:141
Pythia8_i::POSITRON
@ POSITRON
Definition: Pythia8_i.h:120
Pythia8_i::m_mergingWeight
double m_mergingWeight
Definition: Pythia8_i.h:161
mc.generator
generator
Configure Herwig7 These are the commands corresponding to what would go into the regular Herwig infil...
Definition: mc.MGH7_FxFx_H71-DEFAULT_test.py:18
Pythia8_i::ANTIMUON
@ ANTIMUON
Definition: Pythia8_i.h:120
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
GenModule::m_randomSeed
IntegerProperty m_randomSeed
Seed for random number engine.
Definition: GenModule.h:84
Pythia8_i::m_outputParticleDataFile
StringProperty m_outputParticleDataFile
Definition: Pythia8_i.h:159
Pythia8_i::m_conversion
double m_conversion
Definition: Pythia8_i.h:137
Pythia8_i::m_dsid
IntegerProperty m_dsid
Definition: Pythia8_i.h:101
Pythia8_i::m_override_transform_beamenergy
bool m_override_transform_beamenergy
Definition: Pythia8_i.h:127
ATHRNG::calculateSeedsMC21
void calculateSeedsMC21(long *seeds, const std::string &algName, uint64_t ev, uint64_t run, uint64_t offset=0)
Set the random seed using a string (e.g.
Definition: RNGWrapper.cxx:37
Pythia8_i::m_sameAlphaSAsMPI
BooleanProperty m_sameAlphaSAsMPI
Definition: Pythia8_i.h:106
Pythia8_i::m_userResonances
StringProperty m_userResonances
Definition: Pythia8_i.h:152
Pythia8_i::m_doCKKWLAcceptance
BooleanProperty m_doCKKWLAcceptance
Definition: Pythia8_i.h:131
Pythia8_i::m_commands
StringArrayProperty m_commands
Definition: Pythia8_i.h:116
merge.status
status
Definition: merge.py:17
Pythia8_i::m_collisionEnergy
DoubleProperty m_collisionEnergy
Definition: Pythia8_i.h:122
Pythia8_i::m_useRndmGenSvc
BooleanProperty m_useRndmGenSvc
Definition: Pythia8_i.h:96
python.IoTestsLib.w
def w
Definition: IoTestsLib.py:200
Pythia8_i::PROTON
@ PROTON
Definition: Pythia8_i.h:120
Pythia8_i::m_beam1
StringProperty m_beam1
Definition: Pythia8_i.h:125
Pythia8_UserHooks::UserHooksFactory::create
static UserHooks * create(const std::string &hookName)
Definition: UserHooksFactory.cxx:10
Pythia8_i::genInitialize
virtual StatusCode genInitialize()
For initializing the generator, if required.
Definition: Pythia8_i.cxx:97
Trk::split
@ split
Definition: LayerMaterialProperties.h:38
Pythia8_i::m_internal_event_number
int m_internal_event_number
Definition: Pythia8_i.h:112
Pythia8_i::pythia_stream
static const std::string & pythia_stream()
Definition: Pythia8_i.cxx:704
Pythia8_i::m_pythia
std::unique_ptr< Pythia8::Pythia > m_pythia
Definition: Pythia8_i.h:92
Pythia8_i::m_computeEfficiency
BooleanProperty m_computeEfficiency
Definition: Pythia8_i.h:133
GenModule::getRandomEngineDuringInitialize
CLHEP::HepRandomEngine * getRandomEngineDuringInitialize(const std::string &streamName, unsigned long int randomSeedOffset, unsigned int conditionsRun=1, unsigned int lbn=1) const
Definition: GenModule.cxx:53
Pythia8_i::m_pythiaToHepMC
HepMC::Pythia8ToHepMC m_pythiaToHepMC
Definition: Pythia8_i.h:93
HepMC::set_random_states
void set_random_states(GenEvent *e, std::vector< T > a)
Definition: GenEvent.h:647