CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
fastsim::Decayer Class Reference

Implementation of non-stable particle decays. More...

#include <Decayer.h>

Public Member Functions

void decay (const Particle &particle, std::vector< std::unique_ptr< Particle > > &secondaries, CLHEP::HepRandomEngine &engine) const
 Decay particle using pythia. More...
 
 Decayer ()
 Default Constructor. More...
 
 ~Decayer ()
 Default destructor. More...
 

Private Attributes

std::unique_ptr< Pythia8::Pythia > pythia_
 Instance of pythia. More...
 
std::unique_ptr< gen::P8RndmEnginepythiaRandomEngine_
 Instance of pythia Random Engine. More...
 

Detailed Description

Implementation of non-stable particle decays.

Inspired by method Pythia8Hadronizer::residualDecay() in GeneratorInterface/Pythia8Interface/src/Py8GunBase.cc

Definition at line 37 of file Decayer.h.

Constructor & Destructor Documentation

fastsim::Decayer::Decayer ( )

Default Constructor.

Definition at line 11 of file Decayer.cc.

References sysUtil::pid, pythia_, and pythiaRandomEngine_.

12  : pythia_(new Pythia8::Pythia())
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 }
std::unique_ptr< gen::P8RndmEngine > pythiaRandomEngine_
Instance of pythia Random Engine.
Definition: Decayer.h:56
HepPDT::ParticleData ParticleData
std::unique_ptr< Pythia8::Pythia > pythia_
Instance of pythia.
Definition: Decayer.h:55
fastsim::Decayer::~Decayer ( )

Default destructor.

Definition at line 9 of file Decayer.cc.

9 {;}

Member Function Documentation

void fastsim::Decayer::decay ( const Particle particle,
std::vector< std::unique_ptr< Particle > > &  secondaries,
CLHEP::HepRandomEngine &  engine 
) const

Decay particle using pythia.

Parameters
particleThe particle that should be decayed.
secondariesThe decay products.
engineThe Random Engine.

Definition at line 33 of file Decayer.cc.

References funct::abs(), fastsim::Particle::charge(), fastsim::Particle::getMotherDeltaR(), fastsim::Particle::getMotherPdgId(), fastsim::Particle::momentum(), fastsim::Particle::pdgId(), sysUtil::pid, fastsim::Particle::position(), pythia_, pythiaRandomEngine_, and fastsim::Particle::simTrackIndex().

Referenced by FastSimProducer::createFSimTrack(), and FastSimProducer::produce().

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
double T() const
Definition: Particle.h:55
double X() const
Definition: Particle.h:49
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double Z() const
Definition: Particle.h:53
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
Definition: Particle.h:19
double E() const
Definition: Particle.h:63
std::unique_ptr< Pythia8::Pythia > pythia_
Instance of pythia.
Definition: Decayer.h:55
double Y() const
Definition: Particle.h:51

Member Data Documentation

std::unique_ptr<Pythia8::Pythia> fastsim::Decayer::pythia_
private

Instance of pythia.

Definition at line 55 of file Decayer.h.

Referenced by decay(), and Decayer().

std::unique_ptr<gen::P8RndmEngine> fastsim::Decayer::pythiaRandomEngine_
private

Instance of pythia Random Engine.

Definition at line 56 of file Decayer.h.

Referenced by decay(), and Decayer().