CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PythiaDecays.cc
Go to the documentation of this file.
1 // HepMC Headers
2 #include "HepMC/PythiaWrapper6_4.h"
3 
4 // Pythia6 framework integration service Headers
6 
7 // FAMOS Headers
11 
13 
15 
16 // Needed for Pythia6
17 #define PYTHIA6PYDECY pythia6pydecy_
18 
19 extern "C" {
20  void PYTHIA6PYDECY(int *ip);
21 }
22 
24 {
25  program_=program;
26  if (program_ == "pythia6") {
28  pyjets = new Pythia6jets();
30  // The PYTHIA decay tables will be initialized later
31  } else if (program_ == "pythia8") {
33  pythia.reset(new Pythia8::Pythia);
34  decayer.reset(new Pythia8::Pythia);
35 
37  pythia->setRndmEnginePtr(p8RndmEngine.get());
38  decayer->setRndmEnginePtr(p8RndmEngine.get());
39  } else {
40  std::cout << "WARNING: you are requesting an option which is not available in PythiaDecays::PythiaDecays " << std::endl;
41  }
42 
43 }
44 
46  if (program_ == "pythia6") {
47  delete pyjets;
48  delete pyservice;
49  }
50 }
51 
53 PythiaDecays::particleDaughtersPy8(ParticlePropagator& particle, CLHEP::HepRandomEngine* engine)
54 {
56 
57  theList.clear();
58 
59  /*
60  std::cout << "###" << std::endl;
61  std::cout << "particle.PDGname() == " << particle.PDGname() << std::endl;
62  std::cout << "particle.status() == " << particle.status() << std::endl;
63  std::cout << "particle.pid() == " << particle.pid() << std::endl;
64  std::cout << "particle.momentum().x() == " << particle.momentum().x() << std::endl;
65  std::cout << "particle.Px() == " << particle.Px() << std::endl;
66  std::cout << "particle.X() == " << particle.X() << std::endl;
67  std::cout << "particle.mass() == " << particle.mass() << std::endl;
68  std::cout << "particle.PDGmass() == " << particle.PDGmass() << std::endl;
69  std::cout << "particle.PDGcTau() == " << particle.PDGcTau() << std::endl;
70  std::cout << "(decayer->particleData).tau0( particle.pid() ) == " << (decayer->particleData).tau0( particle.pid() ) << std::endl;
71  */
72 
73  // inspired by method Pythia8Hadronizer::residualDecay() in GeneratorInterface/Pythia8Interface/src/Pythia8Hadronizer.cc
74  decayer->event.reset();
75  Pythia8::Particle py8part( particle.pid(), 93, 0, 0, 0, 0, 0, 0,
76  particle.momentum().x(), // note: momentum().x() and Px() are the same
77  particle.momentum().y(),
78  particle.momentum().z(),
79  particle.momentum().t(),
80  particle.mass() );
81  py8part.vProd( particle.X(), particle.Y(),
82  particle.Z(), particle.T() );
83  py8part.tau( (decayer->particleData).tau0( particle.pid() ) );
84  decayer->event.append( py8part );
85 
86  int nentries = decayer->event.size();
87  if ( !decayer->event[nentries-1].mayDecay() ) return theList;
88  decayer->next();
89  int nentries1 = decayer->event.size();
90  if ( nentries1 <= nentries ) return theList; //same number of particles, no decays...
91  Pythia8::Particle& py8daughter = decayer->event[nentries]; // the 1st daughter // DO I NEED THIS LINE?
92 
93  theList.resize(nentries,RawParticle());
94 
95  for ( int ipart=nentries+1; ipart<nentries1; ipart++ )
96  {
97  py8daughter = decayer->event[ipart];
98  theList[ipart-nentries-1].SetXYZT( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
99  theList[ipart-nentries-1].setVertex( py8daughter.xProd(),
100  py8daughter.yProd(),
101  py8daughter.zProd(),
102  py8daughter.tProd() );
103  theList[ipart-nentries-1].setID( py8daughter.id() );
104  theList[ipart-nentries-1].setMass( py8daughter.m() );
105  }
106 
107  return theList;
108 }
109 
111 PythiaDecays::particleDaughtersPy6(ParticlePropagator& particle, CLHEP::HepRandomEngine* engine)
112 {
114 
115  gen::Pythia6Service::InstanceWrapper guard(pyservice); // grab Py6 context
116 
117  // Pythia6jets pyjets;
118  int ip;
119 
120  pyjets->k(1,1) = 1;
121  pyjets->k(1,2) = particle.pid();
122  pyjets->p(1,1) = particle.Px();
123  pyjets->p(1,2) = particle.Py();
124  pyjets->p(1,3) = particle.Pz();
125  pyjets->p(1,4) = std::max(particle.mass(),particle.e());
126  pyjets->p(1,5) = particle.mass();
127  pyjets->v(1,1) = particle.X();
128  pyjets->v(1,2) = particle.Y();
129  pyjets->v(1,3) = particle.Z();
130  pyjets->v(1,4) = particle.T();
131  pyjets->n() = 1;
132 
133  ip = 1;
134  PYTHIA6PYDECY(&ip);
135 
136  // Fill the list of daughters
137  theList.clear();
138  if ( pyjets->n()==1 ) return theList;
139 
140  theList.resize(pyjets->n()-1,RawParticle());
141 
142  for (int i=2;i<=pyjets->n();++i) {
143 
144  theList[i-2].SetXYZT(pyjets->p(i,1),pyjets->p(i,2),
145  pyjets->p(i,3),pyjets->p(i,4));
146  theList[i-2].setVertex(pyjets->v(i,1),pyjets->v(i,2),
147  pyjets->v(i,3),pyjets->v(i,4));
148  theList[i-2].setID(pyjets->k(i,2));
149  theList[i-2].setMass(pyjets->p(i,5));
150 
151  }
152 
153  return theList;
154 
155 }
int i
Definition: DBlmapReader.cc:9
const DaughterParticleList & particleDaughtersPy6(ParticlePropagator &particle, CLHEP::HepRandomEngine *)
int & n(void)
Definition: Pythia6jets.cc:37
Pythia6jets(void)
Definition: Pythia6jets.cc:23
gen::Pythia6Service * pyservice
Definition: PythiaDecays.h:49
int pid() const
get the HEP particle ID number
Definition: RawParticle.h:265
double mass() const
get the MEASURED mass
Definition: RawParticle.h:283
double & v(int i, int j)
Definition: Pythia6jets.cc:73
DaughterParticleList theList
Definition: PythiaDecays.h:46
const DaughterParticleList & particleDaughtersPy8(ParticlePropagator &particle, CLHEP::HepRandomEngine *)
Definition: PythiaDecays.cc:53
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:286
const T & max(const T &a, const T &b)
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:53
std::string program_
Definition: PythiaDecays.h:47
int & k(int i, int j)
Definition: Pythia6jets.cc:49
std::unique_ptr< gen::P8RndmEngine > p8RndmEngine
Definition: PythiaDecays.h:54
double X() const
x of vertex
Definition: RawParticle.h:274
#define PYTHIA6PYDECY
Definition: PythiaDecays.cc:17
Pythia6jets * pyjets
Definition: PythiaDecays.h:50
tuple cout
Definition: gather_cfg.py:121
double T() const
vertex time
Definition: RawParticle.h:277
std::auto_ptr< Pythia8::Pythia > pythia
Definition: PythiaDecays.h:52
double & p(int i, int j)
Definition: Pythia6jets.cc:61
PythiaDecays(std::string program)
Definition: PythiaDecays.cc:23
std::vector< RawParticle > DaughterParticleList
Definition: PythiaDecays.h:28