CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EvtGenInterface.cc
Go to the documentation of this file.
1 // Class Based on EvtGenInterface(LHC).
2 //
3 // Created March 2014
4 //
5 // This class is a modification of the original EvtGenInterface which was developed for EvtGenLHC 9.1.
6 // The modifications for EvtGen 1.3.0 are implemented by Ian M. Nugent
7 // I would like to thank the EvtGen developers, in particular John Black, and Mikhail Kirsanov for their assistance.
8 //
9 // January 2015: Setting of coherent or incoherent B mixing included by Eduard Burelo
10 // January 2015: Adding new feature to allow users to provide new evtgen models
11 
14 
17 
29 
30 #include "CLHEP/Random/RandFlat.h"
32 
33 #include "EvtGen/EvtGen.hh"
34 #include "EvtGenExternal/EvtExternalGenList.hh"
35 #include "EvtGenBase/EvtAbsRadCorr.hh"
36 #include "EvtGenBase/EvtDecayBase.hh"
37 #include "EvtGenBase/EvtId.hh"
38 #include "EvtGenBase/EvtDecayTable.hh"
39 #include "EvtGenBase/EvtParticle.hh"
40 #include "EvtGenBase/EvtParticleFactory.hh"
41 #include "EvtGenBase/EvtHepMCEvent.hh"
42 
43 #include "HepPID/ParticleIDTranslations.hh"
44 
45 #include "TString.h"
46 #include <string>
47 #include <stdlib.h>
48 
49 using namespace gen;
50 using namespace edm;
51 
52 CLHEP::HepRandomEngine* EvtGenInterface::fRandomEngine;
53 
55  fPSet=new ParameterSet(pset);
56  the_engine = new myEvtRandomEngine(nullptr);
57 }
58 
60  // fill up default list of particles to be declared stable in the "master generator"
61  // these are assumed to be PDG ID's
62  //
63  // Note: Pythia6's kc=43, 44, and 84 commented out because they're obsolete (per S.Mrenna)
64  //
65  m_PDGs.push_back( 300553 ) ;
66  m_PDGs.push_back( 511 ) ;
67  m_PDGs.push_back( 521 ) ;
68  m_PDGs.push_back( 523 ) ;
69  m_PDGs.push_back( 513 ) ;
70  m_PDGs.push_back( 533 ) ;
71  m_PDGs.push_back( 531 ) ;
72 
73  m_PDGs.push_back( 15 ) ;
74 
75  m_PDGs.push_back( 413 ) ;
76  m_PDGs.push_back( 423 ) ;
77  m_PDGs.push_back( 433 ) ;
78  m_PDGs.push_back( 411 ) ;
79  m_PDGs.push_back( 421 ) ;
80  m_PDGs.push_back( 431 ) ;
81  m_PDGs.push_back( 10411 );
82  m_PDGs.push_back( 10421 );
83  m_PDGs.push_back( 10413 );
84  m_PDGs.push_back( 10423 );
85  m_PDGs.push_back( 20413 );
86  m_PDGs.push_back( 20423 );
87 
88  m_PDGs.push_back( 415 );
89  m_PDGs.push_back( 425 );
90  m_PDGs.push_back( 10431 );
91  m_PDGs.push_back( 20433 );
92  m_PDGs.push_back( 10433 );
93  m_PDGs.push_back( 435 );
94 
95  m_PDGs.push_back( 310 );
96  m_PDGs.push_back( 311 );
97  m_PDGs.push_back( 313 );
98  m_PDGs.push_back( 323 );
99  m_PDGs.push_back( 10321 );
100  m_PDGs.push_back( 10311 );
101  m_PDGs.push_back( 10313 );
102  m_PDGs.push_back( 10323 );
103  m_PDGs.push_back( 20323 );
104  m_PDGs.push_back( 20313 );
105  m_PDGs.push_back( 325 );
106  m_PDGs.push_back( 315 );
107 
108  m_PDGs.push_back( 100313 );
109  m_PDGs.push_back( 100323 );
110  m_PDGs.push_back( 30313 );
111  m_PDGs.push_back( 30323 );
112  m_PDGs.push_back( 30343 );
113  m_PDGs.push_back( 30353 );
114  m_PDGs.push_back( 30363 );
115 
116  m_PDGs.push_back( 111 );
117  m_PDGs.push_back( 221 );
118  m_PDGs.push_back( 113 );
119  m_PDGs.push_back( 213 );
120  m_PDGs.push_back( 223 );
121  m_PDGs.push_back( 331 );
122  m_PDGs.push_back( 333 );
123  m_PDGs.push_back( 20213 );
124  m_PDGs.push_back( 20113 );
125  m_PDGs.push_back( 215 );
126  m_PDGs.push_back( 115 );
127  m_PDGs.push_back( 10213 );
128  m_PDGs.push_back( 10113 );
129  m_PDGs.push_back( 9000111 ); // PDG ID = 9000111, Pythia6 ID = 10111
130  m_PDGs.push_back( 9000211 ); // PDG ID = 9000211, Pythia6 ID = 10211
131  m_PDGs.push_back( 9010221 ); // PDG ID = 9010211, Pythia6 ID = ???
132  m_PDGs.push_back( 10221 );
133  m_PDGs.push_back( 20223 );
134  m_PDGs.push_back( 20333 );
135  m_PDGs.push_back( 225 );
136  m_PDGs.push_back( 9020221 ); // PDG ID = 9020211, Pythia6 ID = ???
137  m_PDGs.push_back( 335 );
138  m_PDGs.push_back( 10223 );
139  m_PDGs.push_back( 10333 );
140  m_PDGs.push_back( 100213 );
141  m_PDGs.push_back( 100113 );
142 
143  m_PDGs.push_back( 441 );
144  m_PDGs.push_back( 100441 );
145  m_PDGs.push_back( 443 );
146  m_PDGs.push_back( 100443 );
147  m_PDGs.push_back( 9000443 );
148  m_PDGs.push_back( 9010443 );
149  m_PDGs.push_back( 9020443 );
150  m_PDGs.push_back( 10441 );
151  m_PDGs.push_back( 20443 );
152  m_PDGs.push_back( 445 );
153 
154  m_PDGs.push_back( 30443 );
155  m_PDGs.push_back( 551 );
156  m_PDGs.push_back( 553 );
157  m_PDGs.push_back( 100553 );
158  m_PDGs.push_back( 200553 );
159  m_PDGs.push_back( 10551 );
160  m_PDGs.push_back( 20553 );
161  m_PDGs.push_back( 555 );
162  m_PDGs.push_back( 10553 );
163 
164  m_PDGs.push_back( 110551 );
165  m_PDGs.push_back( 120553 );
166  m_PDGs.push_back( 100555 );
167  m_PDGs.push_back( 210551 );
168  m_PDGs.push_back( 220553 );
169  m_PDGs.push_back( 200555 );
170  m_PDGs.push_back( 30553 );
171  m_PDGs.push_back( 20555 );
172 
173  m_PDGs.push_back( 557 );
174  m_PDGs.push_back( 130553 );
175  m_PDGs.push_back( 120555 );
176  m_PDGs.push_back( 100557 );
177  m_PDGs.push_back( 110553 );
178  m_PDGs.push_back( 210553 );
179  m_PDGs.push_back( 10555 );
180  m_PDGs.push_back( 110555 );
181 
182  m_PDGs.push_back( 4122 );
183  m_PDGs.push_back( 4132 );
184  // m_PDGs.push_back( 84 ); // obsolete
185  m_PDGs.push_back( 4112 );
186  m_PDGs.push_back( 4212 );
187  m_PDGs.push_back( 4232 );
188  m_PDGs.push_back( 4222 );
189  m_PDGs.push_back( 4322 );
190  m_PDGs.push_back( 4312 );
191 
192  m_PDGs.push_back( 13122 );
193  m_PDGs.push_back( 13124 );
194  m_PDGs.push_back( 23122 );
195  m_PDGs.push_back( 33122 );
196  m_PDGs.push_back( 43122 );
197  m_PDGs.push_back( 53122 );
198  m_PDGs.push_back( 13126 );
199  m_PDGs.push_back( 13212 );
200  //m_PDGs.push_back( 13241 ); unknown particle -typo?
201 
202  m_PDGs.push_back( 3126 );
203  m_PDGs.push_back( 3124 );
204  m_PDGs.push_back( 3122 );
205  m_PDGs.push_back( 3222 );
206  m_PDGs.push_back( 2214 );
207  m_PDGs.push_back( 2224 );
208  m_PDGs.push_back( 3324 );
209  m_PDGs.push_back( 2114 );
210  m_PDGs.push_back( 1114 );
211  m_PDGs.push_back( 3112 );
212  m_PDGs.push_back( 3212 );
213  m_PDGs.push_back( 3114 );
214  m_PDGs.push_back( 3224 );
215  m_PDGs.push_back( 3214 );
216  m_PDGs.push_back( 3216 );
217  m_PDGs.push_back( 3322 );
218  m_PDGs.push_back( 3312 );
219  m_PDGs.push_back( 3314 );
220  m_PDGs.push_back( 3334 );
221 
222  m_PDGs.push_back( 4114 );
223  m_PDGs.push_back( 4214 );
224  m_PDGs.push_back( 4224 );
225  m_PDGs.push_back( 4314 );
226  m_PDGs.push_back( 4324 );
227  m_PDGs.push_back( 4332 );
228  m_PDGs.push_back( 4334 );
229  //m_PDGs.push_back( 43 ); // obsolete (?)
230  //m_PDGs.push_back( 44 ); // obsolete (?)
231  m_PDGs.push_back( 10443 );
232 
233  m_PDGs.push_back( 5122 );
234  m_PDGs.push_back( 5132 );
235  m_PDGs.push_back( 5232 );
236  m_PDGs.push_back( 5332 );
237  m_PDGs.push_back( 5222 );
238  m_PDGs.push_back( 5112 );
239  m_PDGs.push_back( 5212 );
240  m_PDGs.push_back( 541 );
241  m_PDGs.push_back( 14122 );
242  m_PDGs.push_back( 14124 );
243  m_PDGs.push_back( 5312 );
244  m_PDGs.push_back( 5322 );
245  m_PDGs.push_back( 10521 );
246  m_PDGs.push_back( 20523 );
247  m_PDGs.push_back( 10523 );
248 
249  m_PDGs.push_back( 525 );
250  m_PDGs.push_back( 10511 );
251  m_PDGs.push_back( 20513 );
252  m_PDGs.push_back( 10513 );
253  m_PDGs.push_back( 515 );
254  m_PDGs.push_back( 10531 );
255  m_PDGs.push_back( 20533 );
256  m_PDGs.push_back( 10533 );
257  m_PDGs.push_back( 535 );
258  m_PDGs.push_back( 543 );
259  m_PDGs.push_back( 545 );
260  m_PDGs.push_back( 5114 );
261  m_PDGs.push_back( 5224 );
262  m_PDGs.push_back( 5214 );
263  m_PDGs.push_back( 5314 );
264  m_PDGs.push_back( 5324 );
265  m_PDGs.push_back( 5334 );
266  m_PDGs.push_back( 10541 );
267  m_PDGs.push_back( 10543 );
268  m_PDGs.push_back( 20543 );
269 
270  m_PDGs.push_back( 4424 );
271  m_PDGs.push_back( 4422 );
272  m_PDGs.push_back( 4414 );
273  m_PDGs.push_back( 4412 );
274  m_PDGs.push_back( 4432 );
275  m_PDGs.push_back( 4434 );
276 
277  m_PDGs.push_back( 130 );
278 
279 
280  for(unsigned int i=0; i<m_PDGs.size();i++){
281  int pdt=HepPID::translatePythiatoPDT (m_PDGs.at(i));
282  if(pdt!=m_PDGs.at(i))m_PDGs.at(i)=pdt;
283  }
284 }
285 
287 }
288 
290  edm::FileInPath decay_table(fPSet->getParameter<std::string>("decay_table"));
291  edm::FileInPath pdt(fPSet->getParameter<edm::FileInPath>("particle_property_file"));
292 
293  bool usePythia = fPSet->getUntrackedParameter<bool>("use_internal_pythia",true);
294  bool useTauola = fPSet->getUntrackedParameter<bool>("use_internal_tauola",true);
295  bool usePhotos = fPSet->getUntrackedParameter<bool>("use_internal_photos",true);
296 
297  //Setup evtGen following instructions on http://evtgen.warwick.ac.uk/docs/external/
298  bool convertPythiaCodes=fPSet->getUntrackedParameter<bool>("convertPythiaCodes",true); // Specify if we want to use Pythia 6 physics codes for decays
299  std::string pythiaDir = getenv ("PYTHIA8DATA"); // Specify the pythia xml data directory to use the default PYTHIA8DATA location
300  if(pythiaDir==NULL){
301  std::cout << "EvtGenInterface::init() PYTHIA8DATA not defined. Terminating program " << std::endl; exit(0);
302  }
303  std::string photonType("gamma"); // Specify the photon type for Photos
304  bool useEvtGenRandom(true); // Specify if we want to use the EvtGen random number engine for these generators
305 
306  // Set up the default external generator list: Photos, Pythia and/or Tauola
307  EvtExternalGenList genList(convertPythiaCodes, pythiaDir, photonType, useEvtGenRandom);
308  EvtAbsRadCorr* radCorrEngine=0;
309  if(usePhotos) radCorrEngine = genList.getPhotosModel(); // Get interface to radiative correction engine
310  std::list<EvtDecayBase*> extraModels = genList.getListOfModels(); // get interface to Pythia and Tauola
311  std::list<EvtDecayBase*> myExtraModels;
312  for(unsigned int i=0; i<extraModels.size();i++){
313  std::list<EvtDecayBase*>::iterator it = extraModels.begin();
314  std::advance(it,i);
315  TString name=(*it)->getName();
316  if(name.Contains("PYTHIA") && usePythia) myExtraModels.push_back(*it);
317  if(name.Contains("TAUOLA") && useTauola) myExtraModels.push_back(*it);
318  }
319 
320  //Set up user evtgen models
321 
322  EvtModelUserReg userList;
323  std::list<EvtDecayBase*> userModels = userList.getUserModels(); // get interface to user models
324  for(unsigned int i=0; i<userModels.size();i++){
325  std::list<EvtDecayBase*>::iterator it = userModels.begin();
326  std::advance(it,i);
327  TString name=(*it)->getName();
328  std::cout<<" Adding user model: "<<name<<std::endl;
329  myExtraModels.push_back(*it);
330  }
331 
332 
333 
334  // Set up the incoherent (1) or coherent (0) B mixing option
335  BmixingOption = fPSet->getUntrackedParameter<int>("B_Mixing",1);
336  if(BmixingOption!=0 && BmixingOption!=1){
337  throw cms::Exception("Configuration") << "EvtGenProducer requires B_Mixing to be 0 (coherent) or 1 (incoherent) \n"
338  "Please fix this in your configuration.";
339  }
340 
342  // Create the EvtGen generator object, passing the external generators
343  m_EvtGen = new EvtGen(decay_table.fullPath().c_str(),pdt.fullPath().c_str(),the_engine,radCorrEngine,&myExtraModels,BmixingOption);
344 
345  // Add additional user information
346  if (fPSet->exists("user_decay_file")){
347  std::vector<std::string> user_decays = fPSet->getParameter<std::vector<std::string> >("user_decay_file");
348  for(unsigned int i=0;i<user_decays.size();i++){
349  edm::FileInPath user_decay(user_decays.at(i));
350  m_EvtGen->readUDecay(user_decay.fullPath().c_str());
351  }
352  }
353 
354  // setup pdgid which the generator/hadronizer should not decay
355  if (fPSet->exists("operates_on_particles")){
356  std::vector<int> tmpPIDs = fPSet->getParameter< std::vector<int> >("operates_on_particles");
357  m_PDGs.clear();
358  bool goodinput=false;
359  if(tmpPIDs.size()>0){ if(tmpPIDs.size()==1 && tmpPIDs[0]==0) goodinput=false;}
360  else{goodinput=false;}
361  if(goodinput) m_PDGs = tmpPIDs;
362  else SetDefault_m_PDGs();
363  }
364  else SetDefault_m_PDGs();
365 
366 
367 
368  for(unsigned int i=0;i<m_PDGs.size();i++){
369  std::cout << "EvtGenInterface::init() Particles to Operate on: " << m_PDGs[i] << std::endl;
370  }
371 
372  // Obtain information to set polarization of particles
373  polarize_ids = fPSet->getUntrackedParameter<std::vector<int> >("particles_to_polarize",std::vector<int>());
374  polarize_pol = fPSet->getUntrackedParameter<std::vector<double> >("particle_polarizations",std::vector<double>());
375  if (polarize_ids.size() != polarize_pol.size()) {
376  throw cms::Exception("Configuration") << "EvtGenProducer requires that the particles_to_polarize and particle_polarization\n"
377  "vectors be the same size. Please fix this in your configuration.";
378  }
379  for (unsigned int ndx = 0; ndx < polarize_ids.size(); ndx++) {
380  if (polarize_pol[ndx] < -1. || polarize_pol[ndx] > 1.) {
381  throw cms::Exception("Configuration") << "EvtGenProducer error: particle polarizations must be in the range -1 < P < 1";
382  }
383  polarizations.insert(std::pair<int, float>(polarize_ids[ndx], polarize_pol[ndx]));
384  }
385 
386  // Forced decays are particles that are aliased and forced to be decayed by EvtGen
387  if (fPSet->exists("list_forced_decays")){
388  std::vector<std::string> forced_names = fPSet->getParameter< std::vector<std::string> >("list_forced_decays");
389  for(unsigned int i=0;i<forced_names.size();i++){
390  EvtId found = EvtPDL::getId(forced_names[i]);
391  if(found.getId() == -1) throw cms::Exception("Configuration") << "name in part list for ignored decays not found: " << forced_names[i];
392  if(found.getId() == found.getAlias()) throw cms::Exception("Configuration") << "name of ignored decays is not an alias: " << forced_names[i];
393  forced_id.push_back(found);
394  forced_pdgids.push_back(EvtPDL::getStdHep(found)); // force_pdgids is the list of stdhep codes
395  }
396  }
397 
398  // Ignore decays are particles that are not to be decayed by EvtGen
399  if (fPSet->exists("list_ignored_pdgids")){
400  ignore_pdgids = fPSet->getUntrackedParameter< std::vector<int> >("list_ignored_pdgids");
401  }
402 
403  return;
404 }
405 
406 HepMC::GenEvent* EvtGenInterface::decay( HepMC::GenEvent* evt ){
407  if(the_engine->engine() == nullptr){
409  << "The EvtGen code attempted to use a random number engine while\n"
410  << "the engine pointer was null in EvtGenInterface::decay. This might\n"
411  << "mean that the code was modified to generate a random number outside\n"
412  << "the event and beginLuminosityBlock methods, which is not allowed.\n";
413  }
414  CLHEP::RandFlat m_flat(*the_engine->engine(), 0., 1.);
415 
416  for (HepMC::GenEvent::particle_const_iterator p= evt->particles_begin(); p != evt->particles_end(); ++p){
417  if((*p)->status()==1){ // all particles to be decays are set to status 1 by generator.hadronizer
418  int idHep = (*p)->pdg_id();
419  EvtId idEvt = EvtPDL::evtIdFromStdHep(idHep);
420  for(unsigned int i=0;i<forced_pdgids.size();i++){ // First check if part with forced decay in that case do not decay immediately
421  if(idHep == forced_pdgids[i]){
422  idEvt = forced_id[i];
423  }
424  }
425  int ipart = idEvt.getId();
426  if (ipart==-1) continue; // particle not known to EvtGen
427  EvtDecayTable *evtDecayTable=EvtDecayTable::getInstance();
428  if (evtDecayTable->getNMode(ipart)==0) continue; // particles stable for EvtGen
429  addToHepMC(*p,idEvt,evt); // generate decay
430  }
431  }
432 
433  // add code to ensure all particles have an end vertex and if they are undecayed with no end vertes set to status 1
434  for (HepMC::GenEvent::particle_const_iterator p= evt->particles_begin(); p != evt->particles_end(); ++p){
435  if((*p)->end_vertex() && (*p)->status() == 1)(*p)->set_status(2);
436  if((*p)->end_vertex() && (*p)->end_vertex()->particles_out_size()==0) edm::LogWarning("EvtGenInterface::decay error: empty end vertex!");
437  }
438  return evt;
439 }
440 
441 // Add particle to MC
442 void EvtGenInterface::addToHepMC(HepMC::GenParticle* partHep,const EvtId &idEvt, HepMC::GenEvent* theEvent){
443  // Set up the parent particle from the HepMC GenEvent tree.
444  //EvtVector4R pInit(EvtPDL::getMass(idEvt),partHep->momentum().px(),partHep->momentum().py(),partHep->momentum().pz());
445  EvtVector4R pInit(partHep->momentum().e(),partHep->momentum().px(),partHep->momentum().py(),partHep->momentum().pz());
446  EvtParticle* parent = EvtParticleFactory::particleFactory(idEvt, pInit);
447  // Reset polarization if requested....
448  if(EvtPDL::getSpinType(idEvt) == EvtSpinType::DIRAC && polarizations.count(partHep->pdg_id())>0){
449  HepMC::FourVector momHep = partHep->momentum();
450  EvtVector4R momEvt;
451  momEvt.set(momHep.t(),momHep.x(),momHep.y(),momHep.z());
452  // Particle is spin 1/2, so we can polarize it. Check polarizations map for particle, grab its polarization if it exists
453  // and make the spin density matrix
454  float pol = polarizations.find(partHep->pdg_id())->second;
455  GlobalVector pPart(momHep.x(), momHep.y(), momHep.z());
456  GlobalVector zHat(0., 0., 1.);
457  GlobalVector zCrossP = zHat.cross(pPart);
458  GlobalVector polVec = pol * zCrossP.unit();
459  EvtSpinDensity theSpinDensity;
460  theSpinDensity.setDim(2);
461  theSpinDensity.set(0, 0, EvtComplex(1./2. + polVec.z()/2., 0.));
462  theSpinDensity.set(0, 1, EvtComplex(polVec.x()/2., -polVec.y()/2.));
463  theSpinDensity.set(1, 0, EvtComplex(polVec.x()/2., polVec.y()/2.));
464  theSpinDensity.set(1, 1, EvtComplex(1./2. - polVec.z()/2., 0.));
465  parent->setSpinDensityForwardHelicityBasis(theSpinDensity);
466  }
467  if(parent){
468  // Generate the event
469  m_EvtGen->generateDecay(parent);
470 
471  // Write out the results
472  EvtHepMCEvent evtHepMCEvent;
473  evtHepMCEvent.constructEvent(parent);
474  HepMC::GenEvent* evtGenHepMCTree = evtHepMCEvent.getEvent();
475 
476  // update the event using a recursive function
477  if(!evtGenHepMCTree->particles_empty()) update_particles(partHep,theEvent,(*evtGenHepMCTree->particles_begin()));
478 
479  //clean up
480  parent->deleteTree();
481  }
482 }
483 
484 //Recursivley add EvtGen decay to to Event Decy tree
486  if(p->end_vertex()){
487  if(!partHep->end_vertex()){
488  HepMC::GenVertex* vtx = new HepMC::GenVertex(p->end_vertex()->position());
489  theEvent->add_vertex(vtx);
490  vtx->add_particle_in(partHep);
491  }
492  if(p->end_vertex()->particles_out_size()!=0){
493  for(HepMC::GenVertex::particles_out_const_iterator d=p->end_vertex()->particles_out_const_begin(); d!=p->end_vertex()->particles_out_const_end();d++){
494 
495  // Create daughter and add to event
496  HepMC::GenParticle *daughter = new HepMC::GenParticle((*d)->momentum(),(*d)->pdg_id(),(*d)->status());
497  daughter->suggest_barcode(theEvent->particles_size()+1);
498  partHep->end_vertex()->add_particle_out(daughter);
499 
500  // Ensure forced decays are done with the alias
501  bool isforced=false;
502  for(unsigned int i=0;i<forced_pdgids.size();i++){
503  if(daughter->pdg_id()==forced_pdgids[i]){
504  addToHepMC(daughter,forced_id[i],theEvent);
505  break;
506  }
507  }
508  if(isforced)continue;
509 
510  // check to see if particle is on list of patricles not be decayed by EvtGen
511  bool isignore=false;
512  for(unsigned int i=0;i<ignore_pdgids.size();i++){
513  if(daughter->pdg_id()==ignore_pdgids[i])isignore=true;
514  }
515 
516  // Recursively add daughters
517  if((*d)->end_vertex() && !isignore) update_particles(daughter,theEvent,(*d));
518  }
519  }
520  }
521 }
522 
523 void EvtGenInterface::setRandomEngine(CLHEP::HepRandomEngine* v) {
524  the_engine->setRandomEngine(v);
525  fRandomEngine=v;
526 }
527 
529  if ( !fRandomEngine ) {
530  throw cms::Exception("LogicError")
531  << "EvtGenInterface::flat: Attempt to generate random number when engine pointer is null\n"
532  << "This might mean that the code was modified to generate a random number outside the\n"
533  << "event and beginLuminosityBlock methods, which is not allowed.\n";
534  }
535  return fRandomEngine->flat();
536 }
537 
int i
Definition: DBlmapReader.cc:9
virtual HepMC::GenEvent * decay(HepMC::GenEvent *)
list parent
Definition: dbtoconf.py:74
static unsigned int getId(void)
T y() const
Definition: PV3DBase.h:63
#define NULL
Definition: scimark2.h:8
static double flat()
double v[5][pyjets_maxn]
virtual void setRandomEngine(CLHEP::HepRandomEngine *v)
U second(std::pair< T, U > const &p)
double p[5][pyjets_maxn]
static CLHEP::HepRandomEngine * fRandomEngine
Vector3DBase< typename PreciseFloatType< T, U >::Type, FrameTag > cross(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:119
T z() const
Definition: PV3DBase.h:64
const reco::GenParticle * daughter(const reco::GenParticle &p, unsigned int idau)
Vector3DBase unit() const
Definition: Vector3DBase.h:57
void update_particles(HepMC::GenParticle *partHep, HepMC::GenEvent *theEvent, HepMC::GenParticle *p)
EvtGenInterface(const edm::ParameterSet &)
tuple photonType
PHOTON
Definition: autophobj.py:139
std::list< EvtDecayBase * > getUserModels()
tuple cout
Definition: gather_cfg.py:121
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::string fullPath() const
Definition: FileInPath.cc:165
T x() const
Definition: PV3DBase.h:62
void addToHepMC(HepMC::GenParticle *partHep, const EvtId &idEvt, HepMC::GenEvent *theEvent)