CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
PythiaDecays Class Reference

#include <PythiaDecays.h>

Public Member Functions

const DaughterParticleListparticleDaughtersPy6 (ParticlePropagator &particle)
 
const DaughterParticleListparticleDaughtersPy8 (ParticlePropagator &particle)
 
 PythiaDecays (std::string program)
 
 ~PythiaDecays ()
 

Private Attributes

std::auto_ptr< Pythia8::Pythia > decayer
 
std::string program_
 
Pythia6jetspyjets
 
gen::Pythia6Servicepyservice
 
std::auto_ptr< Pythia8::Pythia > pythia
 
DaughterParticleList theList
 

Detailed Description

Definition at line 27 of file PythiaDecays.h.

Constructor & Destructor Documentation

PythiaDecays::PythiaDecays ( std::string  program)

Definition at line 21 of file PythiaDecays.cc.

References gather_cfg::cout, decayer, program_, pyjets, pyservice, pythia, and Pythia6jets::Pythia6jets().

22 {
23  program_=program;
24  if (program_ == "pythia6") {
26  pyjets = new Pythia6jets();
28  // The PYTHIA decay tables will be initialized later
29  } else if (program_ == "pythia8") {
31  pythia.reset(new Pythia8::Pythia);
32  decayer.reset(new Pythia8::Pythia);
33 
34  RandomP8* RP8 = new RandomP8();
35  pythia->setRndmEnginePtr(RP8);
36  decayer->setRndmEnginePtr(RP8);
37  } else {
38  std::cout << "WARNING: you are requesting an option which is not available in PythiaDecays::PythiaDecays " << std::endl;
39  }
40 
41 }
gen::Pythia6Service * pyservice
Definition: PythiaDecays.h:42
std::auto_ptr< Pythia8::Pythia > decayer
Definition: PythiaDecays.h:46
std::string program_
Definition: PythiaDecays.h:40
Pythia6jets * pyjets
Definition: PythiaDecays.h:43
tuple cout
Definition: gather_cfg.py:121
std::auto_ptr< Pythia8::Pythia > pythia
Definition: PythiaDecays.h:45
PythiaDecays::~PythiaDecays ( )

Definition at line 43 of file PythiaDecays.cc.

References program_, pyjets, and pyservice.

43  {
44  if (program_ == "pythia6") {
45  delete pyjets;
46  delete pyservice;
47  }
48 }
gen::Pythia6Service * pyservice
Definition: PythiaDecays.h:42
std::string program_
Definition: PythiaDecays.h:40
Pythia6jets * pyjets
Definition: PythiaDecays.h:43

Member Function Documentation

const DaughterParticleList & PythiaDecays::particleDaughtersPy6 ( ParticlePropagator particle)

Definition at line 108 of file PythiaDecays.cc.

References i, Pythia6jets::k(), RawParticle::mass(), max(), Pythia6jets::n(), Pythia6jets::p(), RawParticle::pid(), pyjets, pyservice, PYTHIA6PYDECY, RawParticle::T(), theList, Pythia6jets::v(), RawParticle::X(), RawParticle::Y(), and RawParticle::Z().

Referenced by TrajectoryManager::updateWithDaughters().

109 {
110  gen::Pythia6Service::InstanceWrapper guard(pyservice); // grab Py6 context
111 
112  // Pythia6jets pyjets;
113  int ip;
114 
115  pyjets->k(1,1) = 1;
116  pyjets->k(1,2) = particle.pid();
117  pyjets->p(1,1) = particle.Px();
118  pyjets->p(1,2) = particle.Py();
119  pyjets->p(1,3) = particle.Pz();
120  pyjets->p(1,4) = std::max(particle.mass(),particle.e());
121  pyjets->p(1,5) = particle.mass();
122  pyjets->v(1,1) = particle.X();
123  pyjets->v(1,2) = particle.Y();
124  pyjets->v(1,3) = particle.Z();
125  pyjets->v(1,4) = particle.T();
126  pyjets->n() = 1;
127 
128  ip = 1;
129  PYTHIA6PYDECY(&ip);
130 
131  // Fill the list of daughters
132  theList.clear();
133  if ( pyjets->n()==1 ) return theList;
134 
135  theList.resize(pyjets->n()-1,RawParticle());
136 
137  for (int i=2;i<=pyjets->n();++i) {
138 
139  theList[i-2].SetXYZT(pyjets->p(i,1),pyjets->p(i,2),
140  pyjets->p(i,3),pyjets->p(i,4));
141  theList[i-2].setVertex(pyjets->v(i,1),pyjets->v(i,2),
142  pyjets->v(i,3),pyjets->v(i,4));
143  theList[i-2].setID(pyjets->k(i,2));
144  theList[i-2].setMass(pyjets->p(i,5));
145 
146  }
147 
148  return theList;
149 
150 }
int i
Definition: DBlmapReader.cc:9
int & n(void)
Definition: Pythia6jets.cc:37
gen::Pythia6Service * pyservice
Definition: PythiaDecays.h:42
int pid() const
get the HEP particle ID number
Definition: RawParticle.h:264
double mass() const
get the MEASURED mass
Definition: RawParticle.h:282
double & v(int i, int j)
Definition: Pythia6jets.cc:73
DaughterParticleList theList
Definition: PythiaDecays.h:39
const T & max(const T &a, const T &b)
double Y() const
y of vertex
Definition: RawParticle.h:274
double Z() const
z of vertex
Definition: RawParticle.h:275
int & k(int i, int j)
Definition: Pythia6jets.cc:49
double X() const
x of vertex
Definition: RawParticle.h:273
#define PYTHIA6PYDECY
Definition: PythiaDecays.cc:15
Pythia6jets * pyjets
Definition: PythiaDecays.h:43
double T() const
vertex time
Definition: RawParticle.h:276
double & p(int i, int j)
Definition: Pythia6jets.cc:61
const DaughterParticleList & PythiaDecays::particleDaughtersPy8 ( ParticlePropagator particle)

Definition at line 51 of file PythiaDecays.cc.

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

Referenced by TrajectoryManager::updateWithDaughters().

52 {
53 
54  theList.clear();
55 
56  /*
57  std::cout << "###" << std::endl;
58  std::cout << "particle.PDGname() == " << particle.PDGname() << std::endl;
59  std::cout << "particle.status() == " << particle.status() << std::endl;
60  std::cout << "particle.pid() == " << particle.pid() << std::endl;
61  std::cout << "particle.momentum().x() == " << particle.momentum().x() << std::endl;
62  std::cout << "particle.Px() == " << particle.Px() << std::endl;
63  std::cout << "particle.X() == " << particle.X() << std::endl;
64  std::cout << "particle.mass() == " << particle.mass() << std::endl;
65  std::cout << "particle.PDGmass() == " << particle.PDGmass() << std::endl;
66  std::cout << "particle.PDGcTau() == " << particle.PDGcTau() << std::endl;
67  std::cout << "(decayer->particleData).tau0( particle.pid() ) == " << (decayer->particleData).tau0( particle.pid() ) << std::endl;
68  */
69 
70  // inspired by method Pythia8Hadronizer::residualDecay() in GeneratorInterface/Pythia8Interface/src/Pythia8Hadronizer.cc
71  decayer->event.reset();
72  Pythia8::Particle py8part( particle.pid(), 93, 0, 0, 0, 0, 0, 0,
73  particle.momentum().x(), // note: momentum().x() and Px() are the same
74  particle.momentum().y(),
75  particle.momentum().z(),
76  particle.momentum().t(),
77  particle.mass() );
78  py8part.vProd( particle.X(), particle.Y(),
79  particle.Z(), particle.T() );
80  py8part.tau( (decayer->particleData).tau0( particle.pid() ) );
81  decayer->event.append( py8part );
82 
83  int nentries = decayer->event.size();
84  if ( !decayer->event[nentries-1].mayDecay() ) return theList;
85  decayer->next();
86  int nentries1 = decayer->event.size();
87  if ( nentries1 <= nentries ) return theList; //same number of particles, no decays...
88  Pythia8::Particle& py8daughter = decayer->event[nentries]; // the 1st daughter // DO I NEED THIS LINE?
89 
90  theList.resize(nentries,RawParticle());
91 
92  for ( int ipart=nentries+1; ipart<nentries1; ipart++ )
93  {
94  py8daughter = decayer->event[ipart];
95  theList[ipart-nentries-1].SetXYZT( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
96  theList[ipart-nentries-1].setVertex( py8daughter.xProd(),
97  py8daughter.yProd(),
98  py8daughter.zProd(),
99  py8daughter.tProd() );
100  theList[ipart-nentries-1].setID( py8daughter.id() );
101  theList[ipart-nentries-1].setMass( py8daughter.m() );
102  }
103 
104  return theList;
105 }
int pid() const
get the HEP particle ID number
Definition: RawParticle.h:264
double mass() const
get the MEASURED mass
Definition: RawParticle.h:282
DaughterParticleList theList
Definition: PythiaDecays.h:39
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:285
double Y() const
y of vertex
Definition: RawParticle.h:274
double Z() const
z of vertex
Definition: RawParticle.h:275
std::auto_ptr< Pythia8::Pythia > decayer
Definition: PythiaDecays.h:46
double X() const
x of vertex
Definition: RawParticle.h:273
double T() const
vertex time
Definition: RawParticle.h:276

Member Data Documentation

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

Definition at line 46 of file PythiaDecays.h.

Referenced by particleDaughtersPy8(), and PythiaDecays().

std::string PythiaDecays::program_
private

Definition at line 40 of file PythiaDecays.h.

Referenced by PythiaDecays(), and ~PythiaDecays().

Pythia6jets* PythiaDecays::pyjets
private

Definition at line 43 of file PythiaDecays.h.

Referenced by particleDaughtersPy6(), PythiaDecays(), and ~PythiaDecays().

gen::Pythia6Service* PythiaDecays::pyservice
private

Definition at line 42 of file PythiaDecays.h.

Referenced by particleDaughtersPy6(), PythiaDecays(), and ~PythiaDecays().

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

Definition at line 45 of file PythiaDecays.h.

Referenced by PythiaDecays().

DaughterParticleList PythiaDecays::theList
private

Definition at line 39 of file PythiaDecays.h.

Referenced by particleDaughtersPy6(), and particleDaughtersPy8().