CMS 3D CMS Logo

LHAupLesHouches.cc
Go to the documentation of this file.
1 #include <algorithm>
2 #include <iostream>
3 #include <iterator>
4 #include <sstream>
5 #include <string>
6 #include <memory>
7 #include <cassert>
8 
11 
12 
13 using namespace Pythia8;
14 
15 
17 {
18  if (!runInfo) return false;
19  const lhef::HEPRUP &heprup = *runInfo->getHEPRUP();
20 
21  setBeamA(heprup.IDBMUP.first, heprup.EBMUP.first,
22  heprup.PDFGUP.first, heprup.PDFSUP.first);
23  setBeamB(heprup.IDBMUP.second, heprup.EBMUP.second,
24  heprup.PDFGUP.second, heprup.PDFSUP.second);
25  setStrategy(heprup.IDWTUP);
26 
27  for(int i = 0; i < heprup.NPRUP; i++)
28  addProcess(heprup.LPRUP[i], heprup.XSECUP[i],
29  heprup.XERRUP[i], heprup.XMAXUP[i]);
30 
31  //hadronisation->onInit().emit();
32 
33  //runInfo.reset();
34 
35  //fill SLHA header information if available
36  std::vector<std::string> slha = runInfo->findHeader("slha");
37  if (!slha.empty()) {
38  std::string slhaheader;
39  for(std::vector<std::string>::const_iterator iter = slha.begin(); iter != slha.end(); ++iter) {
40  slhaheader.append(*iter);
41  }
42  infoPtr->setHeader("slha",slhaheader);
43  }
44 
45  //will be used to work around missing initialization inside pythia8
46  if(!fEvAttributes) {
47  fEvAttributes = new std::map<std::string, std::string >;
48  } else {
49  fEvAttributes->clear();
50  }
51 
52  return true;
53 }
54 
55 
56 bool LHAupLesHouches::setEvent(int inProcId)
57 {
58  if (!event) return false;
59 
60  if ( event->getReadAttempts() > 0 ) return false; // record already tried
61 
62  const lhef::HEPEUP &hepeup = *event->getHEPEUP();
63 
64  if ( !hepeup.NUP ) return false;
65 
66  setProcess(hepeup.IDPRUP, hepeup.XWGTUP, hepeup.SCALUP,
67  hepeup.AQEDUP, hepeup.AQCDUP);
68 
69  const std::vector<float> &scales = event->scales();
70 
71  unsigned int iscale = 0;
72  for(int i = 0; i < hepeup.NUP; i++) {
73  //retrieve scale corresponding to each particle
74  double scalein = -1.;
75 
76  //handle clustering scales if present,
77  //applies to outgoing partons only
78  if (setScalesFromLHEF_ && !scales.empty() && hepeup.ISTUP[i]==1) {
79  if (iscale>=scales.size()) {
80  edm::LogError("InvalidLHEInput") << "Pythia8 requires"
81  << "cluster scales for all outgoing partons or for none"
82  << std::endl;
83  }
84  scalein = scales[iscale];
85  ++iscale;
86  }
87 
88  addParticle(hepeup.IDUP[i], hepeup.ISTUP[i],
89  hepeup.MOTHUP[i].first, hepeup.MOTHUP[i].second,
90  hepeup.ICOLUP[i].first, hepeup.ICOLUP[i].second,
91  hepeup.PUP[i][0], hepeup.PUP[i][1],
92  hepeup.PUP[i][2], hepeup.PUP[i][3],
93  hepeup.PUP[i][4], hepeup.VTIMUP[i],
94  hepeup.SPINUP[i],scalein);
95  }
96 
97  if(!infoPtr->eventAttributes) {
98  fEvAttributes->clear();
99  infoPtr->eventAttributes = fEvAttributes;
100  } else {
101  infoPtr->eventAttributes = fEvAttributes; // make sure still there
102  infoPtr->eventAttributes->clear();
103  }
104 
105  //fill parton multiplicities if available
106  int npLO = event->npLO();
107  int npNLO = event->npNLO();
108 
109  //default value of -99 indicates tags were not present in the original LHE file
110  //don't pass to pythia in this case to emulate pythia internal lhe reader behaviour
111  if (npLO!=-99) {
112  char buffer [100];
113  snprintf( buffer, 100, "%i",npLO);
114  (*infoPtr->eventAttributes)["npLO"] = buffer;
115  }
116  if (npNLO!=-99) {
117  char buffer [100];
118  snprintf( buffer, 100, "%i",npNLO);
119  (*infoPtr->eventAttributes)["npNLO"] = buffer;
120  }
121 
122  //add #rwgt info from comments
123  const std::vector<std::string> &comments = event->getComments();
124  for (unsigned i=0; i<comments.size(); i++){
125  if (comments[i].rfind("#rwgt", 0)==0)
126  (*infoPtr->eventAttributes)["#rwgt"] = comments[i];
127  }
128 
129  const lhef::LHEEvent::PDF *pdf = event->getPDF();
130  if (pdf) {
131  this->setPdf(pdf->id.first, pdf->id.second,
132  pdf->x.first, pdf->x.second,
133  pdf->scalePDF,
134  pdf->xPDF.first, pdf->xPDF.second, true);
135  }
136  else {
137  this->setPdf(hepeup.IDUP[0], hepeup.IDUP[1],
138  hepeup.PUP[0][3] / runInfo->getHEPRUP()->EBMUP.first,
139  hepeup.PUP[1][3] / runInfo->getHEPRUP()->EBMUP.second,
140  0., 0., 0., false);
141  }
142 
143  event->attempted();
144 
145  return true;
146 }
std::pair< double, double > EBMUP
Definition: LesHouches.h:78
std::pair< double, double > x
Definition: PdfInfo.h:13
std::vector< double > VTIMUP
Definition: LesHouches.h:254
std::pair< int, int > IDBMUP
Definition: LesHouches.h:73
std::pair< double, double > xPDF
Definition: PdfInfo.h:14
std::pair< int, int > PDFGUP
Definition: LesHouches.h:84
std::vector< std::pair< int, int > > MOTHUP
Definition: LesHouches.h:236
std::vector< FiveVector > PUP
Definition: LesHouches.h:248
std::vector< double > SPINUP
Definition: LesHouches.h:261
std::vector< int > ISTUP
Definition: LesHouches.h:230
std::vector< int > IDUP
Definition: LesHouches.h:225
std::vector< double > XERRUP
Definition: LesHouches.h:114
std::pair< int, int > id
Definition: PdfInfo.h:12
std::vector< double > XMAXUP
Definition: LesHouches.h:119
double AQCDUP
Definition: LesHouches.h:220
std::pair< int, int > PDFSUP
Definition: LesHouches.h:90
bool setEvent(int idProcIn) override
bool setInit() override
double AQEDUP
Definition: LesHouches.h:215
std::vector< double > XSECUP
Definition: LesHouches.h:108
double XWGTUP
Definition: LesHouches.h:196
double scalePDF
Definition: PdfInfo.h:15
std::vector< std::pair< int, int > > ICOLUP
Definition: LesHouches.h:242
Definition: event.py:1
std::vector< int > LPRUP
Definition: LesHouches.h:124
double SCALUP
Definition: LesHouches.h:210