CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
gen::EvtGenInterface Class Reference

#include <EvtGenInterface.h>

Inheritance diagram for gen::EvtGenInterface:
gen::EvtGenInterfaceBase

Public Member Functions

virtual HepMC::GenEvent * decay (HepMC::GenEvent *)
 
 EvtGenInterface (const edm::ParameterSet &)
 
virtual void init ()
 
virtual const std::vector< int > & operatesOnParticles ()
 
virtual void setRandomEngine (CLHEP::HepRandomEngine *v)
 
 ~EvtGenInterface ()
 
- Public Member Functions inherited from gen::EvtGenInterfaceBase
 EvtGenInterfaceBase ()
 
virtual void SetPhotosDecayRandomEngine (CLHEP::HepRandomEngine *decayRandomEngine)
 
virtual ~EvtGenInterfaceBase ()
 

Static Public Member Functions

static double flat ()
 

Private Member Functions

void addToHepMC (HepMC::GenParticle *partHep, const EvtId &idEvt, HepMC::GenEvent *theEvent)
 
void SetDefault_m_PDGs ()
 
void update_particles (HepMC::GenParticle *partHep, HepMC::GenEvent *theEvent, HepMC::GenParticle *p)
 

Private Attributes

std::vector< EvtId > forced_id
 
std::vector< int > forced_pdgids
 
edm::ParameterSetfPSet
 
std::vector< int > ignore_pdgids
 
EvtGen * m_EvtGen
 
std::map< int, float > polarizations
 
std::vector< int > polarize_ids
 
std::vector< double > polarize_pol
 
myEvtRandomEnginethe_engine
 

Static Private Attributes

static CLHEP::HepRandomEngine * fRandomEngine
 

Additional Inherited Members

- Protected Attributes inherited from gen::EvtGenInterfaceBase
std::vector< int > m_PDGs
 

Detailed Description

Definition at line 37 of file EvtGenInterface.h.

Constructor & Destructor Documentation

EvtGenInterface::EvtGenInterface ( const edm::ParameterSet pset)

Definition at line 51 of file EvtGenInterface.cc.

51  {
52  fPSet=new ParameterSet(pset);
53  the_engine = new myEvtRandomEngine(nullptr);
54 }
edm::ParameterSet * fPSet
myEvtRandomEngine * the_engine
EvtGenInterface::~EvtGenInterface ( )

Definition at line 283 of file EvtGenInterface.cc.

283  {
284 }

Member Function Documentation

void EvtGenInterface::addToHepMC ( HepMC::GenParticle *  partHep,
const EvtId &  idEvt,
HepMC::GenEvent *  theEvent 
)
private

Definition at line 412 of file EvtGenInterface.cc.

References Vector3DBase< T, FrameTag >::cross(), dbtoconf::parent, edm::second(), Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

412  {
413  // Set up the parent particle from the HepMC GenEvent tree.
414  //EvtVector4R pInit(EvtPDL::getMass(idEvt),partHep->momentum().px(),partHep->momentum().py(),partHep->momentum().pz());
415  EvtVector4R pInit(partHep->momentum().e(),partHep->momentum().px(),partHep->momentum().py(),partHep->momentum().pz());
416  EvtParticle* parent = EvtParticleFactory::particleFactory(idEvt, pInit);
417  // Reset polarization if requested....
418  if(EvtPDL::getSpinType(idEvt) == EvtSpinType::DIRAC && polarizations.count(partHep->pdg_id())>0){
419  HepMC::FourVector momHep = partHep->momentum();
420  EvtVector4R momEvt;
421  momEvt.set(momHep.t(),momHep.x(),momHep.y(),momHep.z());
422  // Particle is spin 1/2, so we can polarize it. Check polarizations map for particle, grab its polarization if it exists
423  // and make the spin density matrix
424  float pol = polarizations.find(partHep->pdg_id())->second;
425  GlobalVector pPart(momHep.x(), momHep.y(), momHep.z());
426  GlobalVector zHat(0., 0., 1.);
427  GlobalVector zCrossP = zHat.cross(pPart);
428  GlobalVector polVec = pol * zCrossP.unit();
429  EvtSpinDensity theSpinDensity;
430  theSpinDensity.setDim(2);
431  theSpinDensity.set(0, 0, EvtComplex(1./2. + polVec.z()/2., 0.));
432  theSpinDensity.set(0, 1, EvtComplex(polVec.x()/2., -polVec.y()/2.));
433  theSpinDensity.set(1, 0, EvtComplex(polVec.x()/2., polVec.y()/2.));
434  theSpinDensity.set(1, 1, EvtComplex(1./2. - polVec.z()/2., 0.));
435  parent->setSpinDensityForwardHelicityBasis(theSpinDensity);
436  }
437  if(parent){
438  // Generate the event
439  m_EvtGen->generateDecay(parent);
440 
441  // Write out the results
442  EvtHepMCEvent evtHepMCEvent;
443  evtHepMCEvent.constructEvent(parent);
444  HepMC::GenEvent* evtGenHepMCTree = evtHepMCEvent.getEvent();
445 
446  // update the event using a recursive function
447  if(!evtGenHepMCTree->particles_empty()) update_particles(partHep,theEvent,(*evtGenHepMCTree->particles_begin()));
448 
449  //clean up
450  parent->deleteTree();
451  }
452 }
list parent
Definition: dbtoconf.py:74
T y() const
Definition: PV3DBase.h:63
U second(std::pair< T, U > const &p)
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
std::map< int, float > polarizations
Vector3DBase unit() const
Definition: Vector3DBase.h:57
void update_particles(HepMC::GenParticle *partHep, HepMC::GenEvent *theEvent, HepMC::GenParticle *p)
T x() const
Definition: PV3DBase.h:62
HepMC::GenEvent * EvtGenInterface::decay ( HepMC::GenEvent *  evt)
virtual

Reimplemented from gen::EvtGenInterfaceBase.

Definition at line 382 of file EvtGenInterface.cc.

References edm::hlt::Exception, i, edm::errors::LogicError, and gen::p.

382  {
383  if(the_engine->engine() == nullptr){
385  << "The EvtGen code attempted to use a random number engine while\n"
386  << "the engine pointer was null in EvtGenInterface::decay. This might\n"
387  << "mean that the code was modified to generate a random number outside\n"
388  << "the event and beginLuminosityBlock methods, which is not allowed.\n";
389  }
390  CLHEP::RandFlat m_flat(*the_engine->engine(), 0., 1.);
391 
392  for (HepMC::GenEvent::particle_const_iterator p= evt->particles_begin(); p != evt->particles_end(); ++p){
393  if((*p)->status()==1){ // all particles to be decays are set to status 1 by generator.hadronizer
394  int idHep = (*p)->pdg_id();
395  EvtId idEvt = EvtPDL::evtIdFromStdHep(idHep);
396  for(unsigned int i=0;i<forced_pdgids.size();i++){ // First check if part with forced decay in that case do not decay immediately
397  if(idHep == forced_pdgids[i]){
398  idEvt = forced_id[i];
399  }
400  }
401  int ipart = idEvt.getId();
402  if (ipart==-1) continue; // particle not known to EvtGen
403  EvtDecayTable *evtDecayTable=EvtDecayTable::getInstance();
404  if (evtDecayTable->getNMode(ipart)==0) continue; // particles stable for EvtGen
405  addToHepMC(*p,idEvt,evt); // generate decay
406  }
407  }
408  return evt;
409 }
int i
Definition: DBlmapReader.cc:9
std::vector< int > forced_pdgids
double p[5][pyjets_maxn]
CLHEP::HepRandomEngine * engine() const
std::vector< EvtId > forced_id
myEvtRandomEngine * the_engine
void addToHepMC(HepMC::GenParticle *partHep, const EvtId &idEvt, HepMC::GenEvent *theEvent)
double EvtGenInterface::flat ( )
static

Definition at line 498 of file EvtGenInterface.cc.

References edm::hlt::Exception.

498  {
499  if ( !fRandomEngine ) {
500  throw cms::Exception("LogicError")
501  << "EvtGenInterface::flat: Attempt to generate random number when engine pointer is null\n"
502  << "This might mean that the code was modified to generate a random number outside the\n"
503  << "event and beginLuminosityBlock methods, which is not allowed.\n";
504  }
505  return fRandomEngine->flat();
506 }
static CLHEP::HepRandomEngine * fRandomEngine
void EvtGenInterface::init ( void  )
virtual

Reimplemented from gen::EvtGenInterfaceBase.

Definition at line 286 of file EvtGenInterface.cc.

References gather_cfg::cout, edm::hlt::Exception, cmsRelvalreport::exit, newFWLiteAna::found, edm::FileInPath::fullPath(), getId(), i, mergeVDriftHistosByStation::name, NULL, and AlCaHLTBitMon_QueryRunRegistry::string.

286  {
287  edm::FileInPath decay_table(fPSet->getParameter<std::string>("decay_table"));
288  edm::FileInPath pdt(fPSet->getParameter<edm::FileInPath>("particle_property_file"));
289 
290  bool usePythia = fPSet->getUntrackedParameter<bool>("use_internal_pythia",true);
291  bool useTauola = fPSet->getUntrackedParameter<bool>("use_internal_tauola",true);
292  bool usePhotos = fPSet->getUntrackedParameter<bool>("use_internal_photos",true);
293 
294  //Setup evtGen following instructions on http://evtgen.warwick.ac.uk/docs/external/
295  bool convertPythiaCodes=fPSet->getUntrackedParameter<bool>("convertPythiaCodes",true); // Specify if we want to use Pythia 6 physics codes for decays
296  std::string pythiaDir = getenv ("PYTHIA8DATA"); // Specify the pythia xml data directory to use the default PYTHIA8DATA location
297  if(pythiaDir==NULL){
298  std::cout << "EvtGenInterface::init() PYTHIA8DATA not defined. Terminating program " << std::endl; exit(0);
299  }
300  std::string photonType("gamma"); // Specify the photon type for Photos
301  bool useEvtGenRandom(true); // Specify if we want to use the EvtGen random number engine for these generators
302 
303  // Set up the default external generator list: Photos, Pythia and/or Tauola
304  EvtExternalGenList genList(convertPythiaCodes, pythiaDir, photonType, useEvtGenRandom);
305  EvtAbsRadCorr* radCorrEngine=0;
306  if(usePhotos) radCorrEngine = genList.getPhotosModel(); // Get interface to radiative correction engine
307  std::list<EvtDecayBase*> extraModels = genList.getListOfModels(); // get interface to Pythia and Tauola
308  std::list<EvtDecayBase*> myExtraModels;
309  for(unsigned int i=0; i<extraModels.size();i++){
310  std::list<EvtDecayBase*>::iterator it = extraModels.begin();
311  std::advance(it,i);
312  TString name=(*it)->getName();
313  if(name.Contains("PYTHIA") && usePythia) myExtraModels.push_back(*it);
314  if(name.Contains("TAUOLA") && useTauola) myExtraModels.push_back(*it);
315  }
316 
318  // Create the EvtGen generator object, passing the external generators
319  m_EvtGen = new EvtGen(decay_table.fullPath().c_str(),pdt.fullPath().c_str(),the_engine,radCorrEngine,&myExtraModels);
320 
321  // Add additional user information
322  if (fPSet->exists("user_decay_file")){
323  std::vector<std::string> user_decays = fPSet->getParameter<std::vector<std::string> >("user_decay_files");
324  for(unsigned int i=0;i<user_decays.size();i++){
325  edm::FileInPath user_decay(user_decays.at(i));
326  m_EvtGen->readUDecay(user_decay.fullPath().c_str());
327  }
328  }
329 
330  // setup pdgid which the generator/hadronizer should not decay
331  if (fPSet->exists("operates_on_particles")){
332  std::vector<int> tmpPIDs = fPSet->getParameter< std::vector<int> >("operates_on_particles");
333  m_PDGs.clear();
334  bool goodinput=false;
335  if(tmpPIDs.size()>0){ if(tmpPIDs.size()==1 && tmpPIDs[0]==0) goodinput=false;}
336  else{goodinput=false;}
337  if(goodinput) m_PDGs = tmpPIDs;
338  else SetDefault_m_PDGs();
339  }
340  else SetDefault_m_PDGs();
341 
342 
343 
344  for(unsigned int i=0;i<m_PDGs.size();i++){
345  std::cout << "EvtGenInterface::init() Particles to Operate on: " << m_PDGs[i] << std::endl;
346  }
347 
348  // Obtain information to set polarization of particles
349  polarize_ids = fPSet->getUntrackedParameter<std::vector<int> >("particles_to_polarize",std::vector<int>());
350  polarize_pol = fPSet->getUntrackedParameter<std::vector<double> >("particle_polarizations",std::vector<double>());
351  if (polarize_ids.size() != polarize_pol.size()) {
352  throw cms::Exception("Configuration") << "EvtGenProducer requires that the particles_to_polarize and particle_polarization\n"
353  "vectors be the same size. Please fix this in your configuration.";
354  }
355  for (unsigned int ndx = 0; ndx < polarize_ids.size(); ndx++) {
356  if (polarize_pol[ndx] < -1. || polarize_pol[ndx] > 1.) {
357  throw cms::Exception("Configuration") << "EvtGenProducer error: particle polarizations must be in the range -1 < P < 1";
358  }
359  polarizations.insert(std::pair<int, float>(polarize_ids[ndx], polarize_pol[ndx]));
360  }
361 
362  // Forced decays are particles that are aliased and forced to be decayed by EvtGen
363  if (fPSet->exists("list_forced_decays")){
364  std::vector<std::string> forced_names = fPSet->getParameter< std::vector<std::string> >("list_forced_decays");
365  for(unsigned int i=0;i<forced_names.size();i++){
366  EvtId found = EvtPDL::getId(forced_names[i]);
367  if(found.getId() == -1) throw cms::Exception("Configuration") << "name in part list for ignored decays not found: " << forced_names[i];
368  if(found.getId() == found.getAlias()) throw cms::Exception("Configuration") << "name of ignored decays is not an alias: " << forced_names[i];
369  forced_id.push_back(found);
370  forced_pdgids.push_back(EvtPDL::getStdHep(found)); // force_pdgids is the list of stdhep codes
371  }
372  }
373 
374  // Ignore decays are particles that are not to be decayed by EvtGen
375  if (fPSet->exists("list_ignored_pdgids")){
376  ignore_pdgids = fPSet->getUntrackedParameter< std::vector<int> >("list_ignored_pdgids");
377  }
378 
379  return;
380 }
T getParameter(std::string const &) const
edm::ParameterSet * fPSet
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
static unsigned int getId(void)
bool exists(std::string const &parameterName) const
checks if a parameter exists
#define NULL
Definition: scimark2.h:8
std::vector< int > polarize_ids
std::vector< int > forced_pdgids
std::map< int, float > polarizations
std::vector< EvtId > forced_id
std::vector< double > polarize_pol
myEvtRandomEngine * the_engine
std::vector< int > ignore_pdgids
tuple cout
Definition: gather_cfg.py:121
virtual const std::vector<int>& gen::EvtGenInterface::operatesOnParticles ( )
inlinevirtual

Reimplemented from gen::EvtGenInterfaceBase.

Definition at line 44 of file EvtGenInterface.h.

References gen::EvtGenInterfaceBase::m_PDGs.

44 { return m_PDGs; }
void EvtGenInterface::SetDefault_m_PDGs ( )
private

Definition at line 56 of file EvtGenInterface.cc.

References i.

56  {
57  // fill up default list of particles to be declared stable in the "master generator"
58  // these are assumed to be PDG ID's
59  //
60  // Note: Pythia6's kc=43, 44, and 84 commented out because they're obsolete (per S.Mrenna)
61  //
62  m_PDGs.push_back( 300553 ) ;
63  m_PDGs.push_back( 511 ) ;
64  m_PDGs.push_back( 521 ) ;
65  m_PDGs.push_back( 523 ) ;
66  m_PDGs.push_back( 513 ) ;
67  m_PDGs.push_back( 533 ) ;
68  m_PDGs.push_back( 531 ) ;
69 
70  m_PDGs.push_back( 15 ) ;
71 
72  m_PDGs.push_back( 413 ) ;
73  m_PDGs.push_back( 423 ) ;
74  m_PDGs.push_back( 433 ) ;
75  m_PDGs.push_back( 411 ) ;
76  m_PDGs.push_back( 421 ) ;
77  m_PDGs.push_back( 431 ) ;
78  m_PDGs.push_back( 10411 );
79  m_PDGs.push_back( 10421 );
80  m_PDGs.push_back( 10413 );
81  m_PDGs.push_back( 10423 );
82  m_PDGs.push_back( 20413 );
83  m_PDGs.push_back( 20423 );
84 
85  m_PDGs.push_back( 415 );
86  m_PDGs.push_back( 425 );
87  m_PDGs.push_back( 10431 );
88  m_PDGs.push_back( 20433 );
89  m_PDGs.push_back( 10433 );
90  m_PDGs.push_back( 435 );
91 
92  m_PDGs.push_back( 310 );
93  m_PDGs.push_back( 311 );
94  m_PDGs.push_back( 313 );
95  m_PDGs.push_back( 323 );
96  m_PDGs.push_back( 10321 );
97  m_PDGs.push_back( 10311 );
98  m_PDGs.push_back( 10313 );
99  m_PDGs.push_back( 10323 );
100  m_PDGs.push_back( 20323 );
101  m_PDGs.push_back( 20313 );
102  m_PDGs.push_back( 325 );
103  m_PDGs.push_back( 315 );
104 
105  m_PDGs.push_back( 100313 );
106  m_PDGs.push_back( 100323 );
107  m_PDGs.push_back( 30313 );
108  m_PDGs.push_back( 30323 );
109  m_PDGs.push_back( 30343 );
110  m_PDGs.push_back( 30353 );
111  m_PDGs.push_back( 30363 );
112 
113  m_PDGs.push_back( 111 );
114  m_PDGs.push_back( 221 );
115  m_PDGs.push_back( 113 );
116  m_PDGs.push_back( 213 );
117  m_PDGs.push_back( 223 );
118  m_PDGs.push_back( 331 );
119  m_PDGs.push_back( 333 );
120  m_PDGs.push_back( 20213 );
121  m_PDGs.push_back( 20113 );
122  m_PDGs.push_back( 215 );
123  m_PDGs.push_back( 115 );
124  m_PDGs.push_back( 10213 );
125  m_PDGs.push_back( 10113 );
126  m_PDGs.push_back( 9000111 ); // PDG ID = 9000111, Pythia6 ID = 10111
127  m_PDGs.push_back( 9000211 ); // PDG ID = 9000211, Pythia6 ID = 10211
128  m_PDGs.push_back( 9010221 ); // PDG ID = 9010211, Pythia6 ID = ???
129  m_PDGs.push_back( 10221 );
130  m_PDGs.push_back( 20223 );
131  m_PDGs.push_back( 20333 );
132  m_PDGs.push_back( 225 );
133  m_PDGs.push_back( 9020221 ); // PDG ID = 9020211, Pythia6 ID = ???
134  m_PDGs.push_back( 335 );
135  m_PDGs.push_back( 10223 );
136  m_PDGs.push_back( 10333 );
137  m_PDGs.push_back( 100213 );
138  m_PDGs.push_back( 100113 );
139 
140  m_PDGs.push_back( 441 );
141  m_PDGs.push_back( 100441 );
142  m_PDGs.push_back( 443 );
143  m_PDGs.push_back( 100443 );
144  m_PDGs.push_back( 9000443 );
145  m_PDGs.push_back( 9010443 );
146  m_PDGs.push_back( 9020443 );
147  m_PDGs.push_back( 10441 );
148  m_PDGs.push_back( 20443 );
149  m_PDGs.push_back( 445 );
150 
151  m_PDGs.push_back( 30443 );
152  m_PDGs.push_back( 551 );
153  m_PDGs.push_back( 553 );
154  m_PDGs.push_back( 100553 );
155  m_PDGs.push_back( 200553 );
156  m_PDGs.push_back( 10551 );
157  m_PDGs.push_back( 20553 );
158  m_PDGs.push_back( 555 );
159  m_PDGs.push_back( 10553 );
160 
161  m_PDGs.push_back( 110551 );
162  m_PDGs.push_back( 120553 );
163  m_PDGs.push_back( 100555 );
164  m_PDGs.push_back( 210551 );
165  m_PDGs.push_back( 220553 );
166  m_PDGs.push_back( 200555 );
167  m_PDGs.push_back( 30553 );
168  m_PDGs.push_back( 20555 );
169 
170  m_PDGs.push_back( 557 );
171  m_PDGs.push_back( 130553 );
172  m_PDGs.push_back( 120555 );
173  m_PDGs.push_back( 100557 );
174  m_PDGs.push_back( 110553 );
175  m_PDGs.push_back( 210553 );
176  m_PDGs.push_back( 10555 );
177  m_PDGs.push_back( 110555 );
178 
179  m_PDGs.push_back( 4122 );
180  m_PDGs.push_back( 4132 );
181  // m_PDGs.push_back( 84 ); // obsolete
182  m_PDGs.push_back( 4112 );
183  m_PDGs.push_back( 4212 );
184  m_PDGs.push_back( 4232 );
185  m_PDGs.push_back( 4222 );
186  m_PDGs.push_back( 4322 );
187  m_PDGs.push_back( 4312 );
188 
189  m_PDGs.push_back( 13122 );
190  m_PDGs.push_back( 13124 );
191  m_PDGs.push_back( 23122 );
192  m_PDGs.push_back( 33122 );
193  m_PDGs.push_back( 43122 );
194  m_PDGs.push_back( 53122 );
195  m_PDGs.push_back( 13126 );
196  m_PDGs.push_back( 13212 );
197  //m_PDGs.push_back( 13241 ); unknown particle -typo?
198 
199  m_PDGs.push_back( 3126 );
200  m_PDGs.push_back( 3124 );
201  m_PDGs.push_back( 3122 );
202  m_PDGs.push_back( 3222 );
203  m_PDGs.push_back( 2214 );
204  m_PDGs.push_back( 2224 );
205  m_PDGs.push_back( 3324 );
206  m_PDGs.push_back( 2114 );
207  m_PDGs.push_back( 1114 );
208  m_PDGs.push_back( 3112 );
209  m_PDGs.push_back( 3212 );
210  m_PDGs.push_back( 3114 );
211  m_PDGs.push_back( 3224 );
212  m_PDGs.push_back( 3214 );
213  m_PDGs.push_back( 3216 );
214  m_PDGs.push_back( 3322 );
215  m_PDGs.push_back( 3312 );
216  m_PDGs.push_back( 3314 );
217  m_PDGs.push_back( 3334 );
218 
219  m_PDGs.push_back( 4114 );
220  m_PDGs.push_back( 4214 );
221  m_PDGs.push_back( 4224 );
222  m_PDGs.push_back( 4314 );
223  m_PDGs.push_back( 4324 );
224  m_PDGs.push_back( 4332 );
225  m_PDGs.push_back( 4334 );
226  //m_PDGs.push_back( 43 ); // obsolete (?)
227  //m_PDGs.push_back( 44 ); // obsolete (?)
228  m_PDGs.push_back( 10443 );
229 
230  m_PDGs.push_back( 5122 );
231  m_PDGs.push_back( 5132 );
232  m_PDGs.push_back( 5232 );
233  m_PDGs.push_back( 5332 );
234  m_PDGs.push_back( 5222 );
235  m_PDGs.push_back( 5112 );
236  m_PDGs.push_back( 5212 );
237  m_PDGs.push_back( 541 );
238  m_PDGs.push_back( 14122 );
239  m_PDGs.push_back( 14124 );
240  m_PDGs.push_back( 5312 );
241  m_PDGs.push_back( 5322 );
242  m_PDGs.push_back( 10521 );
243  m_PDGs.push_back( 20523 );
244  m_PDGs.push_back( 10523 );
245 
246  m_PDGs.push_back( 525 );
247  m_PDGs.push_back( 10511 );
248  m_PDGs.push_back( 20513 );
249  m_PDGs.push_back( 10513 );
250  m_PDGs.push_back( 515 );
251  m_PDGs.push_back( 10531 );
252  m_PDGs.push_back( 20533 );
253  m_PDGs.push_back( 10533 );
254  m_PDGs.push_back( 535 );
255  m_PDGs.push_back( 543 );
256  m_PDGs.push_back( 545 );
257  m_PDGs.push_back( 5114 );
258  m_PDGs.push_back( 5224 );
259  m_PDGs.push_back( 5214 );
260  m_PDGs.push_back( 5314 );
261  m_PDGs.push_back( 5324 );
262  m_PDGs.push_back( 5334 );
263  m_PDGs.push_back( 10541 );
264  m_PDGs.push_back( 10543 );
265  m_PDGs.push_back( 20543 );
266 
267  m_PDGs.push_back( 4424 );
268  m_PDGs.push_back( 4422 );
269  m_PDGs.push_back( 4414 );
270  m_PDGs.push_back( 4412 );
271  m_PDGs.push_back( 4432 );
272  m_PDGs.push_back( 4434 );
273 
274  m_PDGs.push_back( 130 );
275 
276 
277  for(unsigned int i=0; i<m_PDGs.size();i++){
278  int pdt=HepPID::translatePythiatoPDT (m_PDGs.at(i));
279  if(pdt!=m_PDGs.at(i))m_PDGs.at(i)=pdt;
280  }
281 }
int i
Definition: DBlmapReader.cc:9
void EvtGenInterface::setRandomEngine ( CLHEP::HepRandomEngine *  v)
virtual

Implements gen::EvtGenInterfaceBase.

Definition at line 493 of file EvtGenInterface.cc.

References gen::v.

493  {
496 }
double v[5][pyjets_maxn]
static CLHEP::HepRandomEngine * fRandomEngine
myEvtRandomEngine * the_engine
void setRandomEngine(CLHEP::HepRandomEngine *v)
void EvtGenInterface::update_particles ( HepMC::GenParticle *  partHep,
HepMC::GenEvent *  theEvent,
HepMC::GenParticle *  p 
)
private

Definition at line 455 of file EvtGenInterface.cc.

References configurableAnalysis::GenParticle, and i.

455  {
456  if(p->end_vertex()){
457  if(!partHep->end_vertex()){
458  HepMC::GenVertex* vtx = new HepMC::GenVertex(p->end_vertex()->position());
459  theEvent->add_vertex(vtx);
460  vtx->add_particle_in(partHep);
461  }
462  if(p->end_vertex()->particles_out_size()!=0){
463  for(HepMC::GenVertex::particles_out_const_iterator d=p->end_vertex()->particles_out_const_begin(); d!=p->end_vertex()->particles_out_const_end();d++){
464 
465  // Create daughter and add to event
466  HepMC::GenParticle *daughter = new HepMC::GenParticle((*d)->momentum(),(*d)->pdg_id(),(*d)->status());
467  daughter->suggest_barcode(theEvent->particles_size()+1);
468  partHep->end_vertex()->add_particle_out(daughter);
469 
470  // Ensure forced decays are done with the alias
471  bool isforced=false;
472  for(unsigned int i=0;i<forced_pdgids.size();i++){
473  if(daughter->pdg_id()==forced_pdgids[i]){
474  addToHepMC(daughter,forced_id[i],theEvent);
475  break;
476  }
477  }
478  if(isforced)continue;
479 
480  // check to see if particle is on list of patricles not be decayed by EvtGen
481  bool isignore=false;
482  for(unsigned int i=0;i<ignore_pdgids.size();i++){
483  if(daughter->pdg_id()==ignore_pdgids[i])isignore=true;
484  }
485 
486  // Recursively add daughters
487  if((*d)->end_vertex() && !isignore) update_particles(daughter,theEvent,(*d));
488  }
489  }
490  }
491 }
int i
Definition: DBlmapReader.cc:9
std::vector< int > forced_pdgids
double p[5][pyjets_maxn]
std::vector< EvtId > forced_id
void update_particles(HepMC::GenParticle *partHep, HepMC::GenEvent *theEvent, HepMC::GenParticle *p)
std::vector< int > ignore_pdgids
void addToHepMC(HepMC::GenParticle *partHep, const EvtId &idEvt, HepMC::GenEvent *theEvent)

Member Data Documentation

std::vector<EvtId> gen::EvtGenInterface::forced_id
private

Definition at line 56 of file EvtGenInterface.h.

std::vector<int> gen::EvtGenInterface::forced_pdgids
private

Definition at line 57 of file EvtGenInterface.h.

edm::ParameterSet* gen::EvtGenInterface::fPSet
private

Definition at line 66 of file EvtGenInterface.h.

CLHEP::HepRandomEngine * EvtGenInterface::fRandomEngine
staticprivate

Definition at line 68 of file EvtGenInterface.h.

std::vector<int> gen::EvtGenInterface::ignore_pdgids
private

Definition at line 59 of file EvtGenInterface.h.

EvtGen* gen::EvtGenInterface::m_EvtGen
private

Definition at line 54 of file EvtGenInterface.h.

std::map<int, float> gen::EvtGenInterface::polarizations
private

Definition at line 64 of file EvtGenInterface.h.

std::vector<int> gen::EvtGenInterface::polarize_ids
private

Definition at line 62 of file EvtGenInterface.h.

std::vector<double> gen::EvtGenInterface::polarize_pol
private

Definition at line 63 of file EvtGenInterface.h.

myEvtRandomEngine* gen::EvtGenInterface::the_engine
private

Definition at line 69 of file EvtGenInterface.h.