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

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

References pythia_, and pythiaRandomEngine_.

◆ ~Decayer()

fastsim::Decayer::~Decayer ( )

Default destructor.

Definition at line 10 of file Decayer.cc.

10 { ; }

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

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

References funct::abs(), fastsim::Particle::charge(), fastsim::Particle::getMotherDeltaR(), fastsim::Particle::getMotherPdgId(), isExotic(), 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
isExotic
bool isExotic(int pdgid_)
Definition: ParticleManager.h:144
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