CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 <assert.h>
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  //work around missing initialization inside pythia8
46  infoPtr->eventAttributes = new std::map<std::string, std::string >;
47 
48 
49  return true;
50 }
51 
52 
53 bool LHAupLesHouches::setEvent(int inProcId)
54 {
55  if (!event) return false;
56 
57  if ( event->getReadAttempts() > 0 ) return false; // record already tried
58 
59  const lhef::HEPEUP &hepeup = *event->getHEPEUP();
60 
61  if ( !hepeup.NUP ) return false;
62 
63  setProcess(hepeup.IDPRUP, hepeup.XWGTUP, hepeup.SCALUP,
64  hepeup.AQEDUP, hepeup.AQCDUP);
65 
66  const std::vector<float> &scales = event->scales();
67 
68  unsigned int iscale = 0;
69  for(int i = 0; i < hepeup.NUP; i++) {
70  //retrieve scale corresponding to each particle
71  double scalein = -1.;
72 
73  //handle clustering scales if present,
74  //applies to outgoing partons only
75  if (setScalesFromLHEF_ && scales.size()>0 && hepeup.ISTUP[i]==1) {
76  if (iscale>=scales.size()) {
77  edm::LogError("InvalidLHEInput") << "Pythia8 requires"
78  << "cluster scales for all outgoing partons or for none"
79  << std::endl;
80  }
81  scalein = scales[iscale];
82  ++iscale;
83  }
84 
85  addParticle(hepeup.IDUP[i], hepeup.ISTUP[i],
86  hepeup.MOTHUP[i].first, hepeup.MOTHUP[i].second,
87  hepeup.ICOLUP[i].first, hepeup.ICOLUP[i].second,
88  hepeup.PUP[i][0], hepeup.PUP[i][1],
89  hepeup.PUP[i][2], hepeup.PUP[i][3],
90  hepeup.PUP[i][4], hepeup.VTIMUP[i],
91  hepeup.SPINUP[i],scalein);
92  }
93 
94  infoPtr->eventAttributes->clear();
95 
96  //fill parton multiplicities if available
97  int npLO = event->npLO();
98  int npNLO = event->npNLO();
99 
100  //default value of -99 indicates tags were not present in the original LHE file
101  //don't pass to pythia in this case to emulate pythia internal lhe reader behaviour
102  if (npLO!=-99) {
103  char buffer [100];
104  snprintf( buffer, 100, "%i",npLO);
105  (*infoPtr->eventAttributes)["npLO"] = buffer;
106  }
107  if (npNLO!=-99) {
108  char buffer [100];
109  snprintf( buffer, 100, "%i",npNLO);
110  (*infoPtr->eventAttributes)["npNLO"] = buffer;
111  }
112 
113  const lhef::LHEEvent::PDF *pdf = event->getPDF();
114  if (pdf) {
115  this->setPdf(pdf->id.first, pdf->id.second,
116  pdf->x.first, pdf->x.second,
117  pdf->scalePDF,
118  pdf->xPDF.first, pdf->xPDF.second, true);
119  }
120  else {
121  this->setPdf(hepeup.IDUP[0], hepeup.IDUP[1],
122  hepeup.PUP[0][3] / runInfo->getHEPRUP()->EBMUP.first,
123  hepeup.PUP[1][3] / runInfo->getHEPRUP()->EBMUP.second,
124  0., 0., 0., false);
125  }
126 
127  event->attempted();
128 
129  return true;
130 }
int i
Definition: DBlmapReader.cc:9
std::pair< double, double > EBMUP
Definition: LesHouches.h:78
std::pair< double, double > x
Definition: PdfInfo.h:11
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:12
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
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
std::vector< int > IDUP
Definition: LesHouches.h:225
std::vector< double > XERRUP
Definition: LesHouches.h:114
std::pair< int, int > id
Definition: PdfInfo.h:10
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)
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:13
std::vector< std::pair< int, int > > ICOLUP
Definition: LesHouches.h:242
std::vector< int > LPRUP
Definition: LesHouches.h:124
double SCALUP
Definition: LesHouches.h:210