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 34 of file Decayer.h.

Constructor & Destructor Documentation

◆ Decayer()

fastsim::Decayer::Decayer ( )

Default Constructor.

Definition at line 11 of file Decayer.cc.

11  : 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 }

References pythia_, and pythiaRandomEngine_.

◆ ~Decayer()

fastsim::Decayer::~Decayer ( )

Default destructor.

Definition at line 9 of file Decayer.cc.

9 { ; }

Member Function Documentation

◆ decay()

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 28 of file Decayer.cc.

30  {
31  // make sure pythia takes random numbers from the engine past through via the function arguments
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 }

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

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

Member Data Documentation

◆ pythia_

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

Instance of pythia.

Definition at line 53 of file Decayer.h.

Referenced by Decayer().

◆ pythiaRandomEngine_

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

Instance of pythia Random Engine.

Definition at line 54 of file Decayer.h.

Referenced by Decayer().

ParticleData
HepPDT::ParticleData ParticleData
Definition: ParticleDataTable.h:9
Particle::E
double E() const
Definition: Particle.h:94
gen::P8RndmEngine
Definition: P8RndmEngine.h:27
Particle::Z
double Z() const
Definition: Particle.h:69
fastsim::Decayer::pythia_
std::unique_ptr< Pythia8::Pythia > pythia_
Instance of pythia.
Definition: Decayer.h:53
Particle::Y
double Y() const
Definition: Particle.h:64
Particle::X
double X() const
Definition: Particle.h:59
Particle::T
double T() const
Definition: Particle.h:74
fastsim::Particle
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
Definition: Particle.h:16
fastsim::Decayer::pythiaRandomEngine_
std::unique_ptr< gen::P8RndmEngine > pythiaRandomEngine_
Instance of pythia Random Engine.
Definition: Decayer.h:54
edm::RandomEngineSentry
Definition: RandomEngineSentry.h:28
math::XYZTLorentzVector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22