CMS 3D CMS Logo

Decayer.cc
Go to the documentation of this file.
6 
7 #include <Pythia8/Pythia.h>
8 #include "Pythia8Plugins/HepMC2.h"
9 
11 
12 fastsim::Decayer::Decayer() : pythia_(new Pythia8::Pythia()), pythiaRandomEngine_(new gen::P8RndmEngine()) {
13  pythia_->setRndmEnginePtr(pythiaRandomEngine_.get());
14  pythia_->settings.flag("ProcessLevel:all", false);
15  pythia_->settings.flag("PartonLevel:FSRinResonances", false);
16  pythia_->settings.flag("ProcessLevel:resonanceDecays", false);
17  pythia_->init();
18 
19  // forbid all decays
20  // (decays are allowed selectively in the decay function)
21  Pythia8::ParticleData& pdt = pythia_->particleData;
22  int pid = 0;
23  while (pdt.nextId(pid) > pid) {
24  pid = pdt.nextId(pid);
25  pdt.mayDecay(pid, false);
26  }
27 }
28 
29 void fastsim::Decayer::decay(const Particle& particle,
30  std::vector<std::unique_ptr<fastsim::Particle> >& secondaries,
31  CLHEP::HepRandomEngine& engine) const {
32  // make sure pythia takes random numbers from the engine past through via the function arguments
33  edm::RandomEngineSentry<gen::P8RndmEngine> sentry(pythiaRandomEngine_.get(), &engine);
34 
35  // inspired by method Pythia8Hadronizer::residualDecay() in GeneratorInterface/Pythia8Interface/src/Py8GunBase.cc
36  int pid = particle.pdgId();
37  // snip decay products of exotic particles or their children. These decay products are preserved from the event record.
38  // limitation: if exotic incurs heavy energy loss during propagation, the saved decay products could be too hard.
39 
40  if (isExotic(pid) || isExotic(particle.getMotherPdgId())) {
41  return;
42  }
43 
44  pythia_->event.reset();
45 
46  // create a pythia particle which has the same properties as the FastSim particle
47  Pythia8::Particle pythiaParticle(pid,
48  93,
49  0,
50  0,
51  0,
52  0,
53  0,
54  0,
55  particle.momentum().X(),
56  particle.momentum().Y(),
57  particle.momentum().Z(),
58  particle.momentum().E(),
59  particle.momentum().M());
60  pythiaParticle.vProd(
61  particle.position().X(), particle.position().Y(), particle.position().Z(), particle.position().T());
62  pythia_->event.append(pythiaParticle);
63 
64  int nentries_before = pythia_->event.size();
65  // switch on the decay of this and only this particle (avoid double decays)
66  pythia_->particleData.mayDecay(pid, true);
67  // do the decay
68  pythia_->next();
69  // switch it off again
70  pythia_->particleData.mayDecay(pid, false);
71  int nentries_after = pythia_->event.size();
72 
73  if (nentries_after <= nentries_before)
74  return;
75 
76  // add decay products back to the event
77  for (int ipart = nentries_before; ipart < nentries_after; ipart++) {
78  Pythia8::Particle& daughter = pythia_->event[ipart];
79 
80  secondaries.emplace_back(new fastsim::Particle(
81  daughter.id(),
82  math::XYZTLorentzVector(daughter.xProd(), daughter.yProd(), daughter.zProd(), daughter.tProd()),
83  math::XYZTLorentzVector(daughter.px(), daughter.py(), daughter.pz(), daughter.e())));
84 
85  // daughter can inherit the SimTrackIndex of mother (if both charged): necessary for FastSim (cheat) tracking
86  if (particle.charge() != 0 && std::abs(particle.charge() - daughter.charge()) < 1E-10) {
87  secondaries.back()->setMotherDeltaR(particle.momentum());
88  secondaries.back()->setMotherPdgId(particle.getMotherDeltaR() == -1 ? particle.pdgId()
89  : particle.getMotherPdgId());
90  secondaries.back()->setMotherSimTrackIndex(particle.simTrackIndex());
91  }
92  }
93 
94  return;
95 }
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
isExotic
bool isExotic(int pdgid_)
Definition: ParticleManager.h:144
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:12
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:29
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
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
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
ParticleManager.h
fastsim::Decayer::~Decayer
~Decayer()
Default destructor.
Definition: Decayer.cc:10
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