CMS 3D CMS Logo

PythiaDecays.cc
Go to the documentation of this file.
5 
6 #include <Pythia8/Pythia.h>
7 #include "Pythia8Plugins/HepMC2.h"
8 
10 {
11  // inspired by method Pythia8Hadronizer::residualDecay() in GeneratorInterface/Pythia8Interface/src/Py8GunBase.cc
12  decayer.reset(new Pythia8::Pythia);
14  decayer->setRndmEnginePtr(p8RndmEngine.get());
15  decayer->settings.flag("ProcessLevel:all",false);
16  decayer->settings.flag("PartonLevel:FSRinResonances",false);
17  decayer->settings.flag("ProcessLevel:resonanceDecays",false);
18  decayer->init();
19 
20  // forbid all decays
21  // (decays are allowed selectively in the particleDaughters function)
22  Pythia8::ParticleData & pdt = decayer->particleData;
23  int pid = 1;
24  while(pdt.nextId(pid) > pid){
25  pid = pdt.nextId(pid);
26  pdt.mayDecay(pid,false);
27  }
28 }
29 
30 
32 }
33 
35 PythiaDecays::particleDaughters(ParticlePropagator& particle, CLHEP::HepRandomEngine* engine)
36 {
38 
39  theList.clear();
40 
41  // inspired by method Pythia8Hadronizer::residualDecay() in GeneratorInterface/Pythia8Interface/src/Py8GunBase.cc
42  int pid = particle.particle().pid();
43  decayer->event.reset();
44  Pythia8::Particle py8part( pid , 93, 0, 0, 0, 0, 0, 0,
45  particle.particle().momentum().x(), // note: momentum().x() and Px() are the same
46  particle.particle().momentum().y(),
47  particle.particle().momentum().z(),
48  particle.particle().momentum().t(),
49  particle.particle().mass() );
50  py8part.vProd( particle.particle().X(), particle.particle().Y(),
51  particle.particle().Z(), particle.particle().T() );
52  decayer->event.append( py8part );
53 
54  int nentries_before = decayer->event.size();
55  decayer->particleData.mayDecay(pid,true); // switch on the decay of this and only this particle (avoid double decays)
56  decayer->next(); // do the decay
57  decayer->particleData.mayDecay(pid,false); // switch it off again
58  int nentries_after = decayer->event.size();
59  if ( nentries_after <= nentries_before ) return theList;
60 
61  theList.reserve(nentries_after - nentries_before);
62 
63 
64  for ( int ipart=nentries_before; ipart<nentries_after; ipart++ )
65  {
66  Pythia8::Particle& py8daughter = decayer->event[ipart];
67  theList.emplace_back(
69  py8daughter.id(),
70  XYZTLorentzVector( py8daughter.px(),
71  py8daughter.py(),
72  py8daughter.pz(),
73  py8daughter.e() ),
74  XYZTLorentzVector( py8daughter.xProd(),
75  py8daughter.yProd(),
76  py8daughter.zProd(),
77  py8daughter.tProd() ) ) ).setMass( py8daughter.m() );
78  }
79 
80  return theList;
81 }
const HepPDT::ParticleDataTable * particleDataTable() const
const DaughterParticleList & particleDaughters(ParticlePropagator &particle, CLHEP::HepRandomEngine *)
Definition: PythiaDecays.cc:35
int pid() const
get the HEP particle ID number
Definition: RawParticle.h:296
RawParticle const & particle() const
The particle being propagated.
double mass() const
get the MEASURED mass
Definition: RawParticle.h:314
DaughterParticleList theList
Definition: PythiaDecays.h:38
RawParticle makeParticle(HepPDT::ParticleDataTable const *, int id, const math::XYZTLorentzVector &p)
Definition: makeParticle.cc:29
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:340
double Y() const
y of vertex
Definition: RawParticle.h:306
double Z() const
z of vertex
Definition: RawParticle.h:307
HepPDT::ParticleData ParticleData
std::unique_ptr< Pythia8::Pythia > decayer
Definition: PythiaDecays.h:39
std::unique_ptr< gen::P8RndmEngine > p8RndmEngine
Definition: PythiaDecays.h:40
double X() const
x of vertex
Definition: RawParticle.h:305
double T() const
vertex time
Definition: RawParticle.h:308
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:27
std::vector< RawParticle > DaughterParticleList
Definition: PythiaDecays.h:25