CMS 3D CMS Logo

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

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

Static Public Member Functions

static double flat ()
 

Private Member Functions

bool addToHepMC (HepMC::GenParticle *partHep, const EvtId &idEvt, HepMC::GenEvent *theEvent, bool del_daug)
 
bool findLastinChain (HepMC::GenParticle *&p)
 
void go_through_daughters (EvtParticle *part)
 
bool hasnoDaughter (HepMC::GenParticle *p)
 
void SetDefault_m_PDGs ()
 
void update_particles (HepMC::GenParticle *partHep, HepMC::GenEvent *theEvent, HepMC::GenParticle *p)
 

Private Attributes

int BmixingOption = 1
 
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< std::string > fSpecialSettings
 
std::vector< int > m_PDGs
 

Detailed Description

Definition at line 36 of file EvtGenInterface.h.

Constructor & Destructor Documentation

◆ EvtGenInterface()

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

Definition at line 55 of file EvtGenInterface.cc.

References muonDTDigis_cfi::pset.

55  {
56  fPSet = new ParameterSet(pset);
57  the_engine = new myEvtRandomEngine(nullptr);
58 }
edm::ParameterSet * fPSet
myEvtRandomEngine * the_engine

◆ ~EvtGenInterface()

EvtGenInterface::~EvtGenInterface ( )
override

Definition at line 287 of file EvtGenInterface.cc.

287 {}

Member Function Documentation

◆ addToHepMC()

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

Definition at line 536 of file EvtGenInterface.cc.

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

539  {
540  // Set up the parent particle from the HepMC GenEvent tree.
541  //EvtVector4R pInit(EvtPDL::getMass(idEvt),partHep->momentum().px(),partHep->momentum().py(),partHep->momentum().pz());
542  EvtVector4R pInit(
543  partHep->momentum().e(), partHep->momentum().px(), partHep->momentum().py(), partHep->momentum().pz());
544  EvtParticle* parent = EvtParticleFactory::particleFactory(idEvt, pInit);
545  HepMC::FourVector posHep = (partHep->production_vertex())->position();
546  EvtVector4R vInit(posHep.t(), posHep.x(), posHep.y(), posHep.z());
547  // Reset polarization if requested....
548  if (EvtPDL::getSpinType(idEvt) == EvtSpinType::DIRAC && polarizations.count(partHep->pdg_id()) > 0) {
549  const HepMC::FourVector& momHep = partHep->momentum();
550  EvtVector4R momEvt;
551  momEvt.set(momHep.t(), momHep.x(), momHep.y(), momHep.z());
552  // Particle is spin 1/2, so we can polarize it. Check polarizations map for particle, grab its polarization if it exists
553  // and make the spin density matrix
554  float pol = polarizations.find(partHep->pdg_id())->second;
555  GlobalVector pPart(momHep.x(), momHep.y(), momHep.z());
556  GlobalVector zHat(0., 0., 1.);
557  GlobalVector zCrossP = zHat.cross(pPart);
558  GlobalVector polVec = pol * zCrossP.unit();
559  EvtSpinDensity theSpinDensity;
560  theSpinDensity.setDim(2);
561  theSpinDensity.set(0, 0, EvtComplex(1. / 2. + polVec.z() / 2., 0.));
562  theSpinDensity.set(0, 1, EvtComplex(polVec.x() / 2., -polVec.y() / 2.));
563  theSpinDensity.set(1, 0, EvtComplex(polVec.x() / 2., polVec.y() / 2.));
564  theSpinDensity.set(1, 1, EvtComplex(1. / 2. - polVec.z() / 2., 0.));
565  parent->setSpinDensityForwardHelicityBasis(theSpinDensity);
566  }
567  if (parent) {
568  // Generate the event
569  m_EvtGen->generateDecay(parent);
570  if (del_daug)
571  go_through_daughters(parent); // will delete all daugthers which are listed as forced, to allow decay them later
572 
573  // create HepMCTree
574  EvtHepMCEvent evtHepMCEvent;
575  evtHepMCEvent.constructEvent(parent, vInit);
576  HepMC::GenEvent* evtGenHepMCTree = evtHepMCEvent.getEvent();
577  parent->deleteTree();
578 
579  // update the event using a recursive function
580  if (!evtGenHepMCTree->particles_empty())
581  update_particles(partHep, theEvent, (*evtGenHepMCTree->particles_begin()));
582 
583  } else {
584  // this should never hapend, in such case, speak out.
585  return false;
586  }
587  return true;
588 }
Vector3DBase< typename PreciseFloatType< T, U >::Type, FrameTag > cross(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:110
T z() const
Definition: PV3DBase.h:61
void go_through_daughters(EvtParticle *part)
U second(std::pair< T, U > const &p)
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
std::map< int, float > polarizations
void update_particles(HepMC::GenParticle *partHep, HepMC::GenEvent *theEvent, HepMC::GenParticle *p)
static int position[264][3]
Definition: ReadPGInfo.cc:289
Vector3DBase unit() const
Definition: Vector3DBase.h:54

◆ decay()

HepMC::GenEvent * EvtGenInterface::decay ( HepMC::GenEvent evt)
overridevirtual

Reimplemented from gen::EvtGenInterfaceBase.

Definition at line 451 of file EvtGenInterface.cc.

References funct::abs(), Exception, mps_fire::i, heavyIonCSV_trainingSettings::idx, createfilelist::int, dqmiolumiharvest::j, edm::errors::LogicError, gen::p, and eostools::which().

451  {
452  if (the_engine->engine() == nullptr) {
454  << "The EvtGen code attempted to use a random number engine while\n"
455  << "the engine pointer was null in EvtGenInterface::decay. This might\n"
456  << "mean that the code was modified to generate a random number outside\n"
457  << "the event and beginLuminosityBlock methods, which is not allowed.\n";
458  }
459  CLHEP::RandFlat m_flat(*the_engine->engine(), 0., 1.);
460 
461  // decay all request unforced particles and store the forced decays to later decay one per event
462  unsigned int nisforced = 0;
463  std::vector<std::vector<HepMC::GenParticle*> > forcedparticles;
464  for (unsigned int i = 0; i < forced_pdgids.size(); i++)
465  forcedparticles.push_back(std::vector<HepMC::GenParticle*>());
466 
467  // notice this is a dynamic loop
468  for (HepMC::GenEvent::particle_const_iterator p = evt->particles_begin(); p != evt->particles_end(); ++p) {
469  if ((*p)->status() == 1) { // all particles to be decays are set to status 1 by generator.hadronizer
470  int idHep = (*p)->pdg_id();
471  bool isignore = false;
472  for (unsigned int i = 0; i < ignore_pdgids.size(); i++) {
473  if (idHep == ignore_pdgids[i])
474  isignore = true;
475  }
476  if (!isignore) {
477  bool isforced = false;
478  for (unsigned int i = 0; i < forced_pdgids.size(); i++) {
479  if (idHep == forced_pdgids[i]) {
480  forcedparticles.at(i).push_back(*p); // storing and counting forced decays.
481  nisforced++;
482  isforced = true;
483  break;
484  }
485  }
486  bool isDefaultEvtGen = false;
487  for (unsigned int i = 0; i < m_PDGs.size(); i++) {
488  if (abs(idHep) == abs(m_PDGs[i])) {
489  isDefaultEvtGen = true;
490  break;
491  }
492  }
493  EvtId idEvt = EvtPDL::evtIdFromStdHep(idHep);
494  int ipart = idEvt.getId();
495  EvtDecayTable* evtDecayTable = EvtDecayTable::getInstance();
496  if (!isforced && isDefaultEvtGen && ipart != -1 && evtDecayTable->getNMode(ipart) != 0) {
497  addToHepMC(*p, idEvt, evt, true); // generate decay and remove daugther if they are forced
498  }
499  }
500  }
501  }
502 
503  // decay all forced particles (only 1/event is forced)... with no mixing allowed
504  unsigned int which = (unsigned int)(nisforced * flat());
505  if (which == nisforced && nisforced > 0)
506  which = nisforced - 1;
507 
508  unsigned int idx = 0;
509  for (unsigned int i = 0; i < forcedparticles.size(); i++) {
510  for (unsigned int j = 0; j < forcedparticles.at(i).size(); j++) {
511  EvtId idEvt = EvtPDL::evtIdFromStdHep(forcedparticles.at(i).at(j)->pdg_id()); // "standard" decay Id
512  if (idx == which) {
513  idEvt = forced_id[i]; // force decay Id
514  edm::LogInfo("EvtGenInterface::decay ")
515  << EvtPDL::getStdHep(idEvt) << " will force to decay " << idx + 1 << " out of " << nisforced << std::endl;
516  }
517  bool decayed = false;
518  while (!decayed) {
519  decayed = addToHepMC(forcedparticles.at(i).at(j), idEvt, evt, false);
520  }
521  idx++;
522  }
523  }
524 
525  // add code to ensure all particles have an end vertex and if they are undecayed with no end vertes set to status 1
526  for (HepMC::GenEvent::particle_const_iterator p = evt->particles_begin(); p != evt->particles_end(); ++p) {
527  if ((*p)->end_vertex() && (*p)->status() == 1)
528  edm::LogWarning("EvtGenInterface::decay error: incorrect status!"); //(*p)->set_status(2);
529  if ((*p)->end_vertex() && (*p)->end_vertex()->particles_out_size() == 0)
530  edm::LogWarning("EvtGenInterface::decay error: empty end vertex!");
531  }
532  return evt;
533 }
CLHEP::HepRandomEngine * engine() const
static double flat()
std::vector< int > forced_pdgids
bool addToHepMC(HepMC::GenParticle *partHep, const EvtId &idEvt, HepMC::GenEvent *theEvent, bool del_daug)
double p[5][pyjets_maxn]
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< EvtId > forced_id
Log< level::Info, false > LogInfo
myEvtRandomEngine * the_engine
std::vector< int > ignore_pdgids
def which(cmd)
Definition: eostools.py:336

◆ findLastinChain()

bool EvtGenInterface::findLastinChain ( HepMC::GenParticle *&  p)
private

Definition at line 628 of file EvtGenInterface.cc.

References funct::abs(), ztail::d, and gen::p.

628  {
629  if (p->end_vertex()) {
630  if (p->end_vertex()->particles_out_size() != 0) {
631  for (HepMC::GenVertex::particles_out_const_iterator d = p->end_vertex()->particles_out_const_begin();
632  d != p->end_vertex()->particles_out_const_end();
633  d++) {
634  if (abs((*d)->pdg_id()) == abs(p->pdg_id())) {
635  p = *d;
637  return false;
638  }
639  }
640  }
641  }
642  return true;
643 }
double p[5][pyjets_maxn]
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
d
Definition: ztail.py:151
bool findLastinChain(HepMC::GenParticle *&p)

◆ flat()

double EvtGenInterface::flat ( )
static

Definition at line 618 of file EvtGenInterface.cc.

References Exception.

618  {
619  if (!fRandomEngine) {
620  throw cms::Exception("LogicError")
621  << "EvtGenInterface::flat: Attempt to generate random number when engine pointer is null\n"
622  << "This might mean that the code was modified to generate a random number outside the\n"
623  << "event and beginLuminosityBlock methods, which is not allowed.\n";
624  }
625  return fRandomEngine->flat();
626 }
static CLHEP::HepRandomEngine * fRandomEngine

◆ go_through_daughters()

void EvtGenInterface::go_through_daughters ( EvtParticle *  part)
private

Definition at line 654 of file EvtGenInterface.cc.

References newFWLiteAna::found, mps_fire::i, and dqmiolumiharvest::j.

654  {
655  int NDaug = part->getNDaug();
656  if (NDaug) {
657  EvtParticle* Daughter;
658  for (int i = 0; i < NDaug; i++) {
659  Daughter = part->getDaug(i);
660  int idHep = Daughter->getPDGId();
661  int found = 0;
662  for (unsigned int j = 0; j < forced_pdgids.size(); j++) {
663  if (idHep == forced_pdgids[j]) {
664  found = 1;
665  Daughter->deleteDaughters();
666  break;
667  }
668  }
669  if (!found)
670  go_through_daughters(Daughter);
671  }
672  }
673 }
void go_through_daughters(EvtParticle *part)
std::vector< int > forced_pdgids
part
Definition: HCALResponse.h:20

◆ hasnoDaughter()

bool EvtGenInterface::hasnoDaughter ( HepMC::GenParticle *  p)
private

Definition at line 645 of file EvtGenInterface.cc.

References gen::p.

645  {
646  if (p->end_vertex()) {
647  if (p->end_vertex()->particles_out_size() != 0) {
648  return false;
649  }
650  }
651  return true;
652 }
double p[5][pyjets_maxn]

◆ init()

void EvtGenInterface::init ( void  )
overridevirtual

Reimplemented from gen::EvtGenInterfaceBase.

Definition at line 289 of file EvtGenInterface.cc.

References Exception, beamvalidation::exit(), newFWLiteAna::found, edm::FileInPath::fullPath(), getId(), EvtModelUserReg::getUserModels(), mps_fire::i, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, dqmiolumiharvest::j, mergeVDriftHistosByStation::name, objects.autophobj::photonType, and AlCaHLTBitMon_QueryRunRegistry::string.

289  {
290  // flags for pythia8
291  fSpecialSettings.push_back("Pythia8:ParticleDecays:mixB = off");
292  //
293 
294  edm::FileInPath decay_table(fPSet->getParameter<std::string>("decay_table"));
295  edm::FileInPath pdt(fPSet->getParameter<edm::FileInPath>("particle_property_file"));
296 
297  bool usePythia = fPSet->getUntrackedParameter<bool>("use_internal_pythia", true);
298  bool useTauola = fPSet->getUntrackedParameter<bool>("use_internal_tauola", true);
299  bool usePhotos = fPSet->getUntrackedParameter<bool>("use_internal_photos", true);
300 
301  //Setup evtGen following instructions on http://evtgen.warwick.ac.uk/docs/external/
302  bool convertPythiaCodes = fPSet->getUntrackedParameter<bool>(
303  "convertPythiaCodes", true); // Specify if we want to use Pythia 6 physics codes for decays
304  std::string pythiaDir =
305  std::getenv("PYTHIA8DATA"); // Specify the pythia xml data directory to use the default PYTHIA8DATA location
306  if (pythiaDir.empty()) {
307  edm::LogError("EvtGenInterface::~EvtGenInterface")
308  << "EvtGenInterface::init() PYTHIA8DATA not defined. Terminating program ";
309  exit(0);
310  }
311  std::string photonType("gamma"); // Specify the photon type for Photos
312  bool useEvtGenRandom(true); // Specify if we want to use the EvtGen random number engine for these generators
313 
314  // Set up the default external generator list: Photos, Pythia and/or Tauola
315  EvtExternalGenList genList(convertPythiaCodes, pythiaDir, photonType, useEvtGenRandom);
316  EvtAbsRadCorr* radCorrEngine = nullptr;
317  if (usePhotos)
318  radCorrEngine = genList.getPhotosModel(); // Get interface to radiative correction engine
319  std::list<EvtDecayBase*> extraModels = genList.getListOfModels(); // get interface to Pythia and Tauola
320  std::list<EvtDecayBase*> myExtraModels;
321  for (unsigned int i = 0; i < extraModels.size(); i++) {
322  std::list<EvtDecayBase*>::iterator it = extraModels.begin();
323  std::advance(it, i);
324  TString name = (*it)->getName();
325  if (name.Contains("PYTHIA") && usePythia)
326  myExtraModels.push_back(*it);
327  if (name.Contains("TAUOLA") && useTauola)
328  myExtraModels.push_back(*it);
329  }
330 
331  //Set up user evtgen models
332 
333  EvtModelUserReg userList;
334  std::list<EvtDecayBase*> userModels = userList.getUserModels(); // get interface to user models
335  for (unsigned int i = 0; i < userModels.size(); i++) {
336  std::list<EvtDecayBase*>::iterator it = userModels.begin();
337  std::advance(it, i);
338  TString name = (*it)->getName();
339  edm::LogInfo("EvtGenInterface::~EvtGenInterface") << "Adding user model: " << name;
340  myExtraModels.push_back(*it);
341  }
342 
343  // Set up the incoherent (1) or coherent (0) B mixing option
344  BmixingOption = fPSet->getUntrackedParameter<int>("B_Mixing", 1);
345  if (BmixingOption != 0 && BmixingOption != 1) {
346  throw cms::Exception("Configuration") << "EvtGenProducer requires B_Mixing to be 0 (coherent) or 1 (incoherent) \n"
347  "Please fix this in your configuration.";
348  }
349 
351  // Create the EvtGen generator object, passing the external generators
352  m_EvtGen = new EvtGen(
353  decay_table.fullPath().c_str(), pdt.fullPath().c_str(), the_engine, radCorrEngine, &myExtraModels, BmixingOption);
354 
355  // Add additional user information
356  if (fPSet->exists("user_decay_file")) {
357  std::vector<std::string> user_decays = fPSet->getParameter<std::vector<std::string> >("user_decay_file");
358  for (unsigned int i = 0; i < user_decays.size(); i++) {
359  edm::FileInPath user_decay(user_decays.at(i));
360  m_EvtGen->readUDecay(user_decay.fullPath().c_str());
361  }
362  }
363 
364  if (fPSet->exists("user_decay_embedded")) {
365  std::vector<std::string> user_decay_lines = fPSet->getParameter<std::vector<std::string> >("user_decay_embedded");
366  char user_decay_tmp[] = "user_decay_tmpfileXXXXXX";
367  int tmp_creation = mkstemp(user_decay_tmp);
368  FILE* tmpf = std::fopen(user_decay_tmp, "w");
369  if (!tmpf || (tmp_creation == -1)) {
370  edm::LogError("EvtGenInterface::~EvtGenInterface")
371  << "EvtGenInterface::init() fails when trying to open a temporary file for embedded user.dec. Terminating "
372  "program ";
373  exit(0);
374  }
375  for (unsigned int i = 0; i < user_decay_lines.size(); i++) {
376  user_decay_lines.at(i) += "\n";
377  std::fputs(user_decay_lines.at(i).c_str(), tmpf);
378  }
379  std::fclose(tmpf);
380  m_EvtGen->readUDecay(user_decay_tmp);
381  }
382 
383  // setup pdgid which the generator/hadronizer should not decay
384  if (fPSet->exists("operates_on_particles")) {
385  std::vector<int> tmpPIDs = fPSet->getParameter<std::vector<int> >("operates_on_particles");
386  m_PDGs.clear();
387  bool goodinput = false;
388  if (!tmpPIDs.empty()) {
389  if (tmpPIDs.size() == 1 && tmpPIDs[0] == 0)
390  goodinput = false;
391  else
392  goodinput = true;
393  } else {
394  goodinput = false;
395  }
396  if (goodinput)
397  m_PDGs = tmpPIDs;
398  else
400  } else
402 
403  for (unsigned int i = 0; i < m_PDGs.size(); i++) {
404  edm::LogInfo("EvtGenInterface::~EvtGenInterface")
405  << "EvtGenInterface::init() Particles to Operate on: " << m_PDGs[i];
406  }
407 
408  // Obtain information to set polarization of particles
409  polarize_ids = fPSet->getUntrackedParameter<std::vector<int> >("particles_to_polarize", std::vector<int>());
410  polarize_pol = fPSet->getUntrackedParameter<std::vector<double> >("particle_polarizations", std::vector<double>());
411  if (polarize_ids.size() != polarize_pol.size()) {
412  throw cms::Exception("Configuration")
413  << "EvtGenProducer requires that the particles_to_polarize and particle_polarization\n"
414  "vectors be the same size. Please fix this in your configuration.";
415  }
416  for (unsigned int ndx = 0; ndx < polarize_ids.size(); ndx++) {
417  if (polarize_pol[ndx] < -1. || polarize_pol[ndx] > 1.) {
418  throw cms::Exception("Configuration")
419  << "EvtGenProducer error: particle polarizations must be in the range -1 < P < 1";
420  }
421  polarizations.insert(std::pair<int, float>(polarize_ids[ndx], polarize_pol[ndx]));
422  }
423 
424  // Forced decays are particles that are aliased and forced to be decayed by EvtGen
425  if (fPSet->exists("list_forced_decays")) {
426  std::vector<std::string> forced_names = fPSet->getParameter<std::vector<std::string> >("list_forced_decays");
427  for (unsigned int i = 0; i < forced_names.size(); i++) {
428  EvtId found = EvtPDL::getId(forced_names[i]);
429  if (found.getId() == -1)
430  throw cms::Exception("Configuration") << "name in part list for ignored decays not found: " << forced_names[i];
431  if (found.getId() == found.getAlias())
432  throw cms::Exception("Configuration") << "name of ignored decays is not an alias: " << forced_names[i];
433  forced_id.push_back(found);
434  forced_pdgids.push_back(EvtPDL::getStdHep(found)); // force_pdgids is the list of stdhep codes
435  }
436  }
437  edm::LogInfo("EvtGenInterface::~EvtGenInterface")
438  << "Number of Forced Paricles is: " << forced_pdgids.size() << std::endl;
439  for (unsigned int j = 0; j < forced_id.size(); j++) {
440  edm::LogInfo("EvtGenInterface::~EvtGenInterface")
441  << "Forced Paricles are: " << forced_pdgids.at(j) << " " << forced_id.at(j) << std::endl;
442  }
443  // Ignore decays are particles that are not to be decayed by EvtGen
444  if (fPSet->exists("list_ignored_pdgids")) {
445  ignore_pdgids = fPSet->getUntrackedParameter<std::vector<int> >("list_ignored_pdgids");
446  }
447 
448  return;
449 }
edm::ParameterSet * fPSet
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< int > polarize_ids
Log< level::Error, false > LogError
photonType
PHOTON
Definition: autophobj.py:157
std::vector< int > forced_pdgids
T getUntrackedParameter(std::string const &, T const &) const
std::map< int, float > polarizations
std::vector< EvtId > forced_id
std::vector< double > polarize_pol
Log< level::Info, false > LogInfo
static unsigned int getId()
myEvtRandomEngine * the_engine
std::vector< std::string > fSpecialSettings
std::list< EvtDecayBase * > getUserModels()
std::vector< int > ignore_pdgids
def exit(msg="")

◆ operatesOnParticles()

const std::vector<int>& gen::EvtGenInterface::operatesOnParticles ( )
inlineoverridevirtual

Reimplemented from gen::EvtGenInterfaceBase.

Definition at line 43 of file EvtGenInterface.h.

References gen::EvtGenInterfaceBase::m_PDGs.

43 { return m_PDGs; }

◆ SetDefault_m_PDGs()

void EvtGenInterface::SetDefault_m_PDGs ( )
private

Definition at line 60 of file EvtGenInterface.cc.

References mps_fire::i.

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

◆ setRandomEngine()

void EvtGenInterface::setRandomEngine ( CLHEP::HepRandomEngine *  v)
overridevirtual

Implements gen::EvtGenInterfaceBase.

Definition at line 613 of file EvtGenInterface.cc.

References gen::v.

613  {
615  fRandomEngine = v;
616 }
double v[5][pyjets_maxn]
static CLHEP::HepRandomEngine * fRandomEngine
myEvtRandomEngine * the_engine
void setRandomEngine(CLHEP::HepRandomEngine *v)

◆ update_particles()

void EvtGenInterface::update_particles ( HepMC::GenParticle *  partHep,
HepMC::GenEvent theEvent,
HepMC::GenParticle *  p 
)
private

Definition at line 591 of file EvtGenInterface.cc.

References ztail::d, GenParticle::GenParticle, gen::p, and L1BJetProducer_cff::vtx.

591  {
592  if (p->end_vertex()) {
593  if (!partHep->end_vertex()) {
594  HepMC::GenVertex* vtx = new HepMC::GenVertex(p->end_vertex()->position());
595  theEvent->add_vertex(vtx);
596  vtx->add_particle_in(partHep);
597  }
598  if (p->end_vertex()->particles_out_size() != 0) {
599  for (HepMC::GenVertex::particles_out_const_iterator d = p->end_vertex()->particles_out_const_begin();
600  d != p->end_vertex()->particles_out_const_end();
601  d++) {
602  HepMC::GenParticle* daughter = new HepMC::GenParticle((*d)->momentum(), (*d)->pdg_id(), (*d)->status());
603  daughter->suggest_barcode(theEvent->particles_size() + 1);
604  partHep->end_vertex()->add_particle_out(daughter);
605  if ((*d)->end_vertex())
606  update_particles(daughter, theEvent, (*d)); // if daugthers add them as well
607  }
608  partHep->set_status(p->status());
609  }
610  }
611 }
double p[5][pyjets_maxn]
d
Definition: ztail.py:151
void update_particles(HepMC::GenParticle *partHep, HepMC::GenEvent *theEvent, HepMC::GenParticle *p)

Member Data Documentation

◆ BmixingOption

int gen::EvtGenInterface::BmixingOption = 1
private

Definition at line 67 of file EvtGenInterface.h.

◆ forced_id

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

Definition at line 58 of file EvtGenInterface.h.

◆ forced_pdgids

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

Definition at line 59 of file EvtGenInterface.h.

◆ fPSet

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

Definition at line 68 of file EvtGenInterface.h.

◆ fRandomEngine

CLHEP::HepRandomEngine * EvtGenInterface::fRandomEngine
staticprivate

Definition at line 70 of file EvtGenInterface.h.

◆ ignore_pdgids

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

Definition at line 61 of file EvtGenInterface.h.

◆ m_EvtGen

EvtGen* gen::EvtGenInterface::m_EvtGen
private

Definition at line 56 of file EvtGenInterface.h.

◆ polarizations

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

Definition at line 66 of file EvtGenInterface.h.

◆ polarize_ids

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

Definition at line 64 of file EvtGenInterface.h.

◆ polarize_pol

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

Definition at line 65 of file EvtGenInterface.h.

◆ the_engine

myEvtRandomEngine* gen::EvtGenInterface::the_engine
private

Definition at line 71 of file EvtGenInterface.h.