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