CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Pythia8Hadronizer.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <sstream>
3 #include <string>
4 #include <memory>
5 #include <stdint.h>
6 
7 #include <HepMC/GenEvent.h>
8 #include <HepMC/GenParticle.h>
9 
10 #include <Pythia.h>
11 #include <HepMCInterface.h>
12 
14 
16 
20 
23 
29 
31 
32 #include "HepPID/ParticleIDTranslations.hh"
33 
35 
36 using namespace gen;
37 using namespace Pythia8;
38 
40  public:
41  Pythia8Hadronizer(const edm::ParameterSet &params);
43 
44  bool initializeForInternalPartons();
45  bool initializeForExternalPartons();
46 
47  bool declareStableParticles(const std::vector<int> &pdgIds);
48  bool declareSpecialSettings( const std::vector<std::string> );
49 
50  void statistics();
51 
52  bool generatePartonsAndHadronize();
53  bool hadronize();
54  bool decay();
55  bool residualDecay();
56  void finalizeEvent();
57 
58  const char *classname() const { return "Pythia8Hadronizer"; }
59 
60  private:
62 
64  double comEnergy;
66  unsigned int pythiaPylistVerbosity;
70  unsigned int maxEventsToPrint;
71 
73 
76 
77  std::auto_ptr<LHAupLesHouches> lhaUP;
78 
79  std::auto_ptr<Pythia> pythia;
80  Event *pythiaEvent;
81  HepMC::I_Pythia8 toHepMC;
82 
83  UserHooks* UserHookPointer;
84 };
85 
86 
88  BaseHadronizer(params),
89  parameters(params.getParameter<edm::ParameterSet>("PythiaParameters")),
90  comEnergy(params.getParameter<double>("comEnergy")),
91  pythiaPylistVerbosity(params.getUntrackedParameter<int>("pythiaPylistVerbosity", 0)),
92  pythiaHepMCVerbosity(params.getUntrackedParameter<bool>("pythiaHepMCVerbosity", false)),
93  maxEventsToPrint(params.getUntrackedParameter<int>("maxEventsToPrint", 0)),
94  LHEInputFileName(params.getUntrackedParameter<string>("LHEInputFileName","")),
95  useUserHook(false),
96  UserHookPointer(0)
97 {
98  if( params.exists( "useUserHook" ) )
99  useUserHook = params.getParameter<bool>("useUserHook");
101 }
102 
104 {
105 }
106 
108 {
109  //Old code that used Pythia8 own random engine
110  //edm::Service<edm::RandomNumberGenerator> rng;
111  //uint32_t seed = rng->mySeed();
112  //Pythia8::Rndm::init(seed);
113 
114  RandomP8* RP8 = new RandomP8();
115 
116  pythia.reset(new Pythia);
117 
118  pythia->setRndmEnginePtr(RP8);
119  if(useUserHook){
121  pythia->setUserHooksPtr(UserHookPointer);
122  }
123 
124  pythiaEvent = &pythia->event;
125 
127  line != parameters.end(); ++line) {
128  if (line->find("Random:") != std::string::npos)
129  throw cms::Exception("PythiaError")
130  << "Attempted to set random number "
131  "using Pythia commands. Please use "
132  "the RandomNumberGeneratorService."
133  << std::endl;
134 
135  if (!pythia->readString(*line))
136  throw cms::Exception("PythiaError")
137  << "Pythia 8 did not accept \""
138  << *line << "\"." << std::endl;
139  }
140 
141  if(pythiaPylistVerbosity > 10) {
143  pythia->settings.listAll();
145  pythia->particleData.listAll();
146  }
147 
148  pythia->init(2212, 2212, comEnergy);
149 
150  pythia->settings.listChanged();
151 
152  return true;
153 }
154 
155 
157 {
158 
159  std::cout << "Initializing for external partons" << std::endl;
160 
161  RandomP8* RP8 = new RandomP8();
162 
163  pythia.reset(new Pythia);
164 
165  pythia->setRndmEnginePtr(RP8);
166 
167  pythiaEvent = &pythia->event;
168 
170  line != parameters.end(); ++line) {
171  if (line->find("Random:") != std::string::npos)
172  throw cms::Exception("PythiaError")
173  << "Attempted to set random number "
174  "using Pythia commands. Please use "
175  "the RandomNumberGeneratorService."
176  << std::endl;
177 
178  if (!pythia->readString(*line))
179  throw cms::Exception("PythiaError")
180  << "Pythia 8 did not accept \""
181  << *line << "\"." << std::endl;
182  }
183 
184  if(pythiaPylistVerbosity > 10) {
186  pythia->settings.listAll();
188  pythia->particleData.listAll();
189  }
190 
191  if(LHEInputFileName != string()) {
192 
193  cout << endl;
194  cout << "LHE Input File Name = " << LHEInputFileName << endl;
195  cout << endl;
196  pythia->init(LHEInputFileName);
197 
198  } else {
199 
200  lhaUP.reset(new LHAupLesHouches());
201  lhaUP->loadRunInfo(lheRunInfo());
202  pythia->init(lhaUP.get());
203 
204  }
205 
206  return true;
207 }
208 
209 
210 #if 0
211 // naive Pythia8 HepMC status fixup
212 static int getStatus(const HepMC::GenParticle *p)
213 {
214  int status = p->status();
215  if (status > 0)
216  return status;
217  else if (status > -30 && status < 0)
218  return 3;
219  else
220  return 2;
221 }
222 #endif
223 
224 bool Pythia8Hadronizer::declareStableParticles(const std::vector<int> &pdgIds)
225 {
226 
227  for ( size_t i=0; i<pdgIds.size(); i++ )
228  {
229  // FIXME: need to double check if PID's are the same in Py6 & Py8,
230  // because the HepPDT translation tool is actually for **Py6**
231  //
232  int PyID = HepPID::translatePDTtoPythia( pdgIds[i] );
233  std::ostringstream pyCard ;
234  pyCard << PyID <<":mayDecay=false";
235  pythia->readString( pyCard.str() );
236  // alternative:
237  // set the 2nd input argument warn=false
238  // - this way Py8 will NOT print warnings about unknown particle code(s)
239  // pythia->readString( pyCard.str(), false )
240  }
241 
242  return true;
243 
244 }
245 bool Pythia8Hadronizer::declareSpecialSettings( const std::vector<std::string> )
246 {
247  return true;
248 }
249 
250 
252 {
253  pythia->statistics();
254 
255  double xsec = pythia->info.sigmaGen(); // cross section in mb
256  xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
257  runInfo().setInternalXSec(xsec);
258 }
259 
261 {
262  if (!pythia->next())
263  return false;
264 
265  event().reset(new HepMC::GenEvent);
266  toHepMC.fill_next_event(*pythiaEvent, event().get());
267 
268  return true;
269 }
270 
272 {
273  if(LHEInputFileName == string()) {
274  //cout << "start loading event" << endl;
275  lhaUP->loadEvent(lheEvent());
276  //cout << "finish loading event" << endl;
277  }
278 
279  if (!pythia->next())
280  return false;
281 
282  // update LHE matching statistics
283  //
285 
286  event().reset(new HepMC::GenEvent);
287  toHepMC.fill_next_event(*pythiaEvent, event().get());
288 
289  return true;
290 }
291 
293 {
294  return true;
295 }
296 
298 {
299  return true;
300 }
301 
303 {
304  bool lhe = lheEvent() != 0;
305 
306  event()->set_signal_process_id(pythia->info.code());
307  event()->set_event_scale(pythia->info.pTHat()); //FIXME
308 
309  // Putting pdf info into the HepMC record
310  // There is the overloaded pythia8 HepMCInterface method fill_next_event
311  // that does this, but CMSSW GeneratorInterface does not fill HepMC
312  // record according to the HepMC convention (stores f(x) instead of x*f(x)
313  // and converts gluon PDG ID to zero). For this reason we use the
314  // method fill_next_event (above) that does NOT this and fill pdf info here
315  //
316  int id1 = pythia->info.id1();
317  int id2 = pythia->info.id2();
318  if (id1 == 21) id1 = 0;
319  if (id2 == 21) id2 = 0;
320  double x1 = pythia->info.x1();
321  double x2 = pythia->info.x2();
322  double Q = pythia->info.QRen();
323  double pdf1 = pythia->info.pdf1() / pythia->info.x1();
324  double pdf2 = pythia->info.pdf2() / pythia->info.x2();
325  event()->set_pdf_info(HepMC::PdfInfo(id1,id2,x1,x2,Q,pdf1,pdf2));
326 
327  // Storing weights. Will be moved to pythia8 HepMCInterface
328  //
329  if (lhe && std::abs(lheRunInfo()->getHEPRUP()->IDWTUP) == 4)
330  // translate mb to pb (CMS/Gen "convention" as of May 2009)
331  event()->weights().push_back( pythia->info.weight() * 1.0e9 );
332  else
333  event()->weights().push_back( pythia->info.weight() );
334  //
335  // Below is the weight to be used in the case of reweighting by UserHook
336  // This is a temporary solution. It requires the explicit conversion of
337  // the pointer. In future versions (>150) the additional factor introduced
338  // in the UserHook will be stored by pythia. Also, putting the weights
339  // into the HepMC record will be made in the pythia8 HepMCInterface
340  //
341  if(useUserHook)
342  event()->weights().push_back(1./((PtHatReweightUserHook*)UserHookPointer)->getFactor());
343 
344  // now create the GenEventInfo product from the GenEvent and fill
345  // the missing pieces
346  eventInfo().reset( new GenEventInfoProduct( event().get() ) );
347 
348  // in pythia pthat is used to subdivide samples into different bins
349  // in LHE mode the binning is done by the external ME generator
350  // which is likely not pthat, so only filling it for Py6 internal mode
351  if (!lhe) {
352  eventInfo()->setBinningValues(std::vector<double>(1, pythia->info.pTHat()));
353  }
354 
355  //******** Verbosity ********
356 
357  if (maxEventsToPrint > 0 &&
360  if (pythiaPylistVerbosity) {
361  pythia->info.list(std::cout);
362  pythia->event.list(std::cout);
363  }
364 
365  if (pythiaHepMCVerbosity) {
366  std::cout << "Event process = "
367  << pythia->info.code() << "\n"
368  << "----------------------" << std::endl;
369  event()->print();
370  }
371  }
372 }
373 
376 
bool pythiaHepMCVerbosity
HepMC verbosity flag.
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
dictionary parameters
Definition: Parameters.py:2
double comEnergy
Center-of-Mass energy.
edm::GeneratorFilter< Pythia8Hadronizer, ExternalDecayDriver > Pythia8GeneratorFilter
bool declareStableParticles(const std::vector< int > &pdgIds)
std::auto_ptr< Pythia > pythia
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
#define abs(x)
Definition: mlp_lapack.h:159
void count(LHERunInfo::CountMode count, double weight=1.0, double matchWeight=1.0)
Definition: LHEEvent.cc:198
bool declareSpecialSettings(const std::vector< std::string >)
static int getStatus(const HepMC::GenParticle *p)
void setInternalXSec(const XSec &xsec)
bool useUserHook
Switch User Hook flag.
std::auto_ptr< HepMC::GenEvent > & event()
CLHEP::HepRandomEngine & getEngineReference()
GenRunInfoProduct & runInfo()
std::auto_ptr< LHAupLesHouches > lhaUP
lhef::LHEEvent * lheEvent()
CLHEP::HepRandomEngine * randomEngine
Definition: PYR.cc:4
const char * classname() const
unsigned int pythiaPylistVerbosity
Pythia PYLIST Verbosity flag.
std::auto_ptr< GenEventInfoProduct > & eventInfo()
lhef::LHERunInfo * lheRunInfo()
ParameterCollector parameters
unsigned int maxEventsToPrint
Events to print if verbosity.
Pythia8Hadronizer(const edm::ParameterSet &params)
const_iterator end() const
const_iterator begin() const
tuple cout
Definition: gather_cfg.py:41
HepMC::I_Pythia8 toHepMC
UserHooks * UserHookPointer
tuple status
Definition: ntuplemaker.py:245
edm::HadronizerFilter< Pythia8Hadronizer, ExternalDecayDriver > Pythia8HadronizerFilter