StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StarPythia8Decayer.cxx
1 #include "StarPythia8Decayer.h"
2 #include "TLorentzVector.h"
3 #include <assert.h>
4 #include "TParticle.h"
5 #include "StarGenerator/UTIL/StarRandom.h"
6 #include "TSystem.h"
7 #include "StMessMgr.h"
8 #include "ParticleData.h"
9 
10 #ifndef Pythia8_version
11 #error "Pythia8_version is not defined"
12 #endif
13 
14 // ----------------------------------------------------------------------------
15 // Remap pythia's random number generator to StarRandom
16 class PyRand : public Pythia8::RndmEngine {
17 public:
18  double flat() { return StarRandom::Instance().flat(); }
19 };
20 
21 //..................................................................................................
22 StarPythia8Decayer::StarPythia8Decayer( Pythia8::Pythia *pythia ) : mPythia(pythia), mOwner(false), mDebug(0)
23 {
24 
25  if ( pythia ) return; // If provided, we expect pythia to be correctly initialized
26 
27 
28  // Create the new instance of pythia, specifying the path to the
29  // auxilliary files required by pythia.
30  // NOTE: When adding new versions of Pythia8, we need to specify
31  // the version of the code in the source code
32  TString path = "StRoot/StarGenerator/"; path+= Pythia8_version; path+="/xmldoc/";
33  {
34  ifstream in(path);
35  if (!in.good()) { path = "$(STAR)/"+path; }
36  path = gSystem->ExpandPathName(path.Data());
37  }
38 
39  Info(GetName(),Form("MC version is %s data at %s",Pythia8_version,path.Data()));
40  Info(GetName(),Form("Configuration files at %s",path.Data()));
41 
42  mPythia = new Pythia8::Pythia( path.Data() );
43  mPythia -> setRndmEnginePtr( new PyRand() );
44  mPythia->readString("SoftQCD:elastic = on");
45 
46  // Initialize for pp @ 510.0 GeV
47  mPythia->init( 2212, 2212, 510.0 );
48 
49 }
50 //..................................................................................................
52 {
53 
54 }
55 //..................................................................................................
56 void StarPythia8Decayer::Decay( int pdgid, TLorentzVector *_p )
57 {
58 
59  LOG_INFO << "Decay pdgid=" << pdgid << endm;
60 
61  // Clear the event from the last run
62  ClearEvent();
63  // Add the particle to the pythia stack
64  AppendParticle( pdgid, _p );
65  // Get the particle ID and allow it to decay
66  int id = mPythia -> event[0].id();
67  mPythia->particleData.mayDecay( id, true );
68  mPythia->moreDecays();
69  // Print list of pythia lines on debug
70  if ( mDebug ) mPythia->event.list();
71 }
72 //..................................................................................................
73 int StarPythia8Decayer::ImportParticles( TClonesArray *_array )
74 {
75  // Save the decay products
76  assert(_array);
77  TClonesArray &array = *_array;
78  array.Clear();
79 
80  int nparts = 0;
81 
82  for ( int i=0; i<mPythia->event.size();i++ )
83  {
84  // if ( mPythia->event[i].id() == 90 ) continue; //??
85  new(array[nparts++]) TParticle (
86  mPythia->event[i].id(),
87  mPythia->event[i].status(),
88  mPythia->event[i].mother1(),
89  mPythia->event[i].mother2(),
90  mPythia->event[i].daughter1(),
91  mPythia->event[i].daughter2(),
92  mPythia->event[i].px(), // [GeV/c]
93  mPythia->event[i].py(), // [GeV/c]
94  mPythia->event[i].pz(), // [GeV/c]
95  mPythia->event[i].e(), // [GeV]
96  mPythia->event[i].xProd(), // [mm]
97  mPythia->event[i].yProd(), // [mm]
98  mPythia->event[i].zProd(), // [mm]
99  mPythia->event[i].tProd()); // [mm/c]
100  };
101  return nparts;
102 }
103 //..................................................................................................
104 // not implemented. complain about it. loudly.
105 void StarPythia8Decayer::SetForceDecay( int type ){ assert(0); }
106 void StarPythia8Decayer::ForceDecay(){ assert(0); }
107 float StarPythia8Decayer::GetPartialBranchingRatio( int ipdg ){ assert(0); }
108 void StarPythia8Decayer::ReadDecayTable(){ assert(0); }
109 //..................................................................................................
111 {
112  // return lifetime in seconds of teh particle with PDG number pdg
113  return (mPythia->particleData.tau0(pdg) * 3.3333e-12) ;
114 }
115 //..................................................................................................
116 void StarPythia8Decayer::AppendParticle(int pdg, TLorentzVector* p)
117 {
118  // Append a particle to the stack
119  const int status = 11;
120 
121  double pt2 = p->Perp2(); // pt squared
122  double px = p->Px();
123  double py = p->Py();
124  double pz = p->Pz(); // z-component
125  double M2 = p->M2(); // mass squared
126 
127  // Lookup particle entry in pythia8
128  const auto* pde = mPythia->particleData.particleDataEntryPtr( pdg );
129  if ( pde ) {
130  M2 = pde->m0() * pde->m0();
131  }
132 
133  double E = TMath::Sqrt( pt2 + M2 + pz*pz );
134 
135  mPythia->event.append(pdg, status, 0, 0, px, py, pz, E, TMath::Sqrt(M2) );
136 }
137 //..................................................................................................
139 {
140  mPythia->event.clear();
141 }
142 //..................................................................................................
144 {
145  if (mPythia) delete mPythia; // delete
146 }
static StarRandom & Instance()
Obtain the single instance of the random number generator.
Definition: StarRandom.cxx:87
float GetPartialBranchingRatio(int pdgid)
Return the branching ratio for the spdcified PDG ID.
~StarPythia8Decayer()
Class destructor.
int ImportParticles(TClonesArray *array=0)
Returns the decay products in a TClonesArray of TParticle.
Double_t flat() const
Return a random number uniformly distributed between 0 and 1.
Definition: StarRandom.cxx:68
void Init()
Initializes the decayer.
StarPythia8Decayer()
Class constructor.
void AppendParticle(int pdgid, TLorentzVector *p=0)
Add a particle with specified PDG ID to the stack.
void ClearEvent()
Clear the event.
float GetLifetime(int pdgid)
Return teh lifetime in seconds for the specified particle.
void Decay(int pdg, TLorentzVector *p=0)
Decays the particle specified by PDG id and lorentz vector.