CMS 3D CMS Logo

Decayer.cc
Go to the documentation of this file.
5 
6 #include <Pythia8/Pythia.h>
7 #include "Pythia8Plugins/HepMC2.h"
8 
10 
11 fastsim::Decayer::Decayer() : pythia_(new Pythia8::Pythia()), pythiaRandomEngine_(new gen::P8RndmEngine()) {
12  pythia_->setRndmEnginePtr(pythiaRandomEngine_.get());
13  pythia_->settings.flag("ProcessLevel:all", false);
14  pythia_->settings.flag("PartonLevel:FSRinResonances", false);
15  pythia_->settings.flag("ProcessLevel:resonanceDecays", false);
16  pythia_->init();
17 
18  // forbid all decays
19  // (decays are allowed selectively in the decay function)
20  Pythia8::ParticleData& pdt = pythia_->particleData;
21  int pid = 0;
22  while (pdt.nextId(pid) > pid) {
23  pid = pdt.nextId(pid);
24  pdt.mayDecay(pid, false);
25  }
26 }
27 
28 void fastsim::Decayer::decay(const Particle& particle,
29  std::vector<std::unique_ptr<fastsim::Particle> >& secondaries,
30  CLHEP::HepRandomEngine& engine) const {
31  // make sure pythia takes random numbers from the engine past through via the function arguments
32  edm::RandomEngineSentry<gen::P8RndmEngine> sentry(pythiaRandomEngine_.get(), &engine);
33 
34  // inspired by method Pythia8Hadronizer::residualDecay() in GeneratorInterface/Pythia8Interface/src/Py8GunBase.cc
35  int pid = particle.pdgId();
36  pythia_->event.reset();
37 
38  // create a pythia particle which has the same properties as the FastSim particle
39  Pythia8::Particle pythiaParticle(pid,
40  93,
41  0,
42  0,
43  0,
44  0,
45  0,
46  0,
47  particle.momentum().X(),
48  particle.momentum().Y(),
49  particle.momentum().Z(),
50  particle.momentum().E(),
51  particle.momentum().M());
52  pythiaParticle.vProd(
53  particle.position().X(), particle.position().Y(), particle.position().Z(), particle.position().T());
54  pythia_->event.append(pythiaParticle);
55 
56  int nentries_before = pythia_->event.size();
57  // switch on the decay of this and only this particle (avoid double decays)
58  pythia_->particleData.mayDecay(pid, true);
59  // do the decay
60  pythia_->next();
61  // switch it off again
62  pythia_->particleData.mayDecay(pid, false);
63  int nentries_after = pythia_->event.size();
64 
65  if (nentries_after <= nentries_before)
66  return;
67 
68  // add decay products back to the event
69  for (int ipart = nentries_before; ipart < nentries_after; ipart++) {
70  Pythia8::Particle& daughter = pythia_->event[ipart];
71 
72  secondaries.emplace_back(new fastsim::Particle(
73  daughter.id(),
74  math::XYZTLorentzVector(daughter.xProd(), daughter.yProd(), daughter.zProd(), daughter.tProd()),
75  math::XYZTLorentzVector(daughter.px(), daughter.py(), daughter.pz(), daughter.e())));
76 
77  // daughter can inherit the SimTrackIndex of mother (if both charged): necessary for FastSim (cheat) tracking
78  if (particle.charge() != 0 && std::abs(particle.charge() - daughter.charge()) < 1E-10) {
79  secondaries.back()->setMotherDeltaR(particle.momentum());
80  secondaries.back()->setMotherPdgId(particle.getMotherDeltaR() == -1 ? particle.pdgId()
81  : particle.getMotherPdgId());
82  secondaries.back()->setMotherSimTrackIndex(particle.simTrackIndex());
83  }
84  }
85 
86  return;
87 }
fastsim::Particle::momentum
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:143
P8RndmEngine.h
fastsim::Particle::simTrackIndex
int simTrackIndex() const
Return index of the SimTrack.
Definition: Particle.h:153
ParticleData
HepPDT::ParticleData ParticleData
Definition: ParticleDataTable.h:9
Particle.h
Pythia8
Definition: PythiaDecays.h:21
fastsim::Decayer::pythia_
std::unique_ptr< Pythia8::Pythia > pythia_
Instance of pythia.
Definition: Decayer.h:53
gen
Definition: PythiaDecays.h:13
fastsim::Particle::position
const math::XYZTLorentzVector & position() const
Return position of the particle.
Definition: Particle.h:140
fastsim::Particle::getMotherPdgId
int getMotherPdgId() const
Get pdgIdto mother particle (only makes sense if mother and daughter charged).
Definition: Particle.h:209
fastsim::Decayer::Decayer
Decayer()
Default Constructor.
Definition: Decayer.cc:11
fastsim::Particle::pdgId
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:134
fastsim::Decayer::decay
void decay(const Particle &particle, std::vector< std::unique_ptr< Particle > > &secondaries, CLHEP::HepRandomEngine &engine) const
Decay particle using pythia.
Definition: Decayer.cc:28
fastsim::Particle::charge
double charge() const
Return charge of the particle.
Definition: Particle.h:137
fastsim::Particle
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
Definition: Particle.h:16
Decayer.h
fastsim::Decayer::pythiaRandomEngine_
std::unique_ptr< gen::P8RndmEngine > pythiaRandomEngine_
Instance of pythia Random Engine.
Definition: Decayer.h:54
RandomEngineSentry.h
edm::RandomEngineSentry
Definition: RandomEngineSentry.h:28
fastsim::Decayer::~Decayer
~Decayer()
Default destructor.
Definition: Decayer.cc:9
math::XYZTLorentzVector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
Particle
Definition: Particle.py:1
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
fastsim::Particle::getMotherDeltaR
double getMotherDeltaR() const
Get delta R to mother particle (only makes sense if mother and daughter charged).
Definition: Particle.h:202