CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
PythiaDecays Class Reference

#include <PythiaDecays.h>

Public Member Functions

const DaughterParticleListparticleDaughters (ParticlePropagator &particle, CLHEP::HepRandomEngine *)
 
 PythiaDecays ()
 
 ~PythiaDecays ()
 

Private Attributes

std::auto_ptr< Pythia8::Pythia > decayer
 
std::unique_ptr< gen::P8RndmEnginep8RndmEngine
 
DaughterParticleList theList
 

Detailed Description

Definition at line 28 of file PythiaDecays.h.

Constructor & Destructor Documentation

PythiaDecays::PythiaDecays ( )

Definition at line 8 of file PythiaDecays.cc.

References decayer, p8RndmEngine, and sysUtil::pid.

9 {
10  // inspired by method Pythia8Hadronizer::residualDecay() in GeneratorInterface/Pythia8Interface/src/Py8GunBase.cc
11  decayer.reset(new Pythia8::Pythia);
13  decayer->setRndmEnginePtr(p8RndmEngine.get());
14  decayer->settings.flag("ProcessLevel:all",false);
15  decayer->settings.flag("PartonLevel:FSRinResonances",false);
16  decayer->settings.flag("ProcessLevel:resonanceDecays",false);
17  decayer->init();
18 
19  // forbid all decays
20  // (decays are allowed selectively in the particleDaughters function)
21  Pythia8::ParticleData & pdt = decayer->particleData;
22  int pid = 1;
23  while(pdt.nextId(pid) > pid){
24  pid = pdt.nextId(pid);
25  pdt.mayDecay(pid,false);
26  }
27 }
std::auto_ptr< Pythia8::Pythia > decayer
Definition: PythiaDecays.h:39
HepPDT::ParticleData ParticleData
std::unique_ptr< gen::P8RndmEngine > p8RndmEngine
Definition: PythiaDecays.h:40
PythiaDecays::~PythiaDecays ( )

Definition at line 30 of file PythiaDecays.cc.

30  {
31 }

Member Function Documentation

const DaughterParticleList & PythiaDecays::particleDaughters ( ParticlePropagator particle,
CLHEP::HepRandomEngine *  engine 
)

Definition at line 34 of file PythiaDecays.cc.

References decayer, RawParticle::mass(), RawParticle::momentum(), p8RndmEngine, sysUtil::pid, RawParticle::pid(), RawParticle::T(), theList, RawParticle::X(), RawParticle::Y(), and RawParticle::Z().

Referenced by TrajectoryManager::updateWithDaughters().

35 {
37 
38  theList.clear();
39 
40  // inspired by method Pythia8Hadronizer::residualDecay() in GeneratorInterface/Pythia8Interface/src/Py8GunBase.cc
41  int pid = particle.pid();
42  decayer->event.reset();
43  Pythia8::Particle py8part( pid , 93, 0, 0, 0, 0, 0, 0,
44  particle.momentum().x(), // note: momentum().x() and Px() are the same
45  particle.momentum().y(),
46  particle.momentum().z(),
47  particle.momentum().t(),
48  particle.mass() );
49  py8part.vProd( particle.X(), particle.Y(),
50  particle.Z(), particle.T() );
51  decayer->event.append( py8part );
52 
53  int nentries_before = decayer->event.size();
54  decayer->particleData.mayDecay(pid,true); // switch on the decay of this and only this particle (avoid double decays)
55  decayer->next(); // do the decay
56  decayer->particleData.mayDecay(pid,false); // switch it off again
57  int nentries_after = decayer->event.size();
58  if ( nentries_after <= nentries_before ) return theList;
59 
60  theList.resize(nentries_after - nentries_before,RawParticle());
61 
62 
63  for ( int ipart=nentries_before; ipart<nentries_after; ipart++ )
64  {
65  Pythia8::Particle& py8daughter = decayer->event[ipart];
66  theList[ipart-nentries_before].SetXYZT( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
67  theList[ipart-nentries_before].setVertex( py8daughter.xProd(),
68  py8daughter.yProd(),
69  py8daughter.zProd(),
70  py8daughter.tProd() );
71  theList[ipart-nentries_before].setID( py8daughter.id() );
72  theList[ipart-nentries_before].setMass( py8daughter.m() );
73  }
74 
75  return theList;
76 }
int pid() const
get the HEP particle ID number
Definition: RawParticle.h:265
double mass() const
get the MEASURED mass
Definition: RawParticle.h:283
DaughterParticleList theList
Definition: PythiaDecays.h:38
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:286
double Y() const
y of vertex
Definition: RawParticle.h:275
double Z() const
z of vertex
Definition: RawParticle.h:276
std::auto_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:274
double T() const
vertex time
Definition: RawParticle.h:277

Member Data Documentation

std::auto_ptr<Pythia8::Pythia> PythiaDecays::decayer
private

Definition at line 39 of file PythiaDecays.h.

Referenced by particleDaughters(), and PythiaDecays().

std::unique_ptr<gen::P8RndmEngine> PythiaDecays::p8RndmEngine
private

Definition at line 40 of file PythiaDecays.h.

Referenced by particleDaughters(), and PythiaDecays().

DaughterParticleList PythiaDecays::theList
private

Definition at line 38 of file PythiaDecays.h.

Referenced by particleDaughters().