StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StarPythia6.cxx
1 #include "StarPythia6.h"
2 ClassImp(StarPythia6);
3 
4 #include "StarCallf77.h"
5 #include "StarGenerator/EVENT/StarGenPPEvent.h"
6 #include "StarGenerator/EVENT/StarGenEPEvent.h"
7 #include "StarGenerator/EVENT/StarGenParticle.h"
8 
9 #include "StarGenerator/UTIL/StarRandom.h"
10 #include <map>
11 #include <iostream>
12 
13 // ----------------------------------------------------------------------------
14 // Remap pythia's random number generator to StarRandom
15 extern "C" {
16  Double_t pyrstar_( Int_t *idummy ){ return StarRandom::Instance().flat(); };
17 };
18 // ----------------------------------------------------------------------------
19 // ----------------------------------------------------------------------------
20 // ----------------------------------------------------------------------------
21 StarPythia6::StarPythia6( const Char_t *name ) : StarGenerator(name)
22 {
23 
24  // Setup a map between pythia6 status codes and the HepMC status codes
25  // used in StarEvent
26  mStatusCode[0] = StarGenParticle::kNull;
27  mStatusCode[1] = StarGenParticle::kFinal;
28  const Int_t decays[] = { 2, 3, 4, 5, 11, 12, 13, 14, 15 };
29  for ( UInt_t i=0;i<sizeof(decays)/sizeof(Int_t); i++ )
30  {
31  mStatusCode[decays[i]]=StarGenParticle::kDecayed;
32  }
33  const Int_t docum[] = { 21, 31, 32, 41, 42, 51, 52 };
34  for ( UInt_t i=0;i<sizeof(docum)/sizeof(Int_t); i++ )
35  {
36  mStatusCode[docum[i]]=StarGenParticle::kDocumentation;
37  }
38 
39 
40 }
41 // ----------------------------------------------------------------------------
42 // ----------------------------------------------------------------------------
43 // ----------------------------------------------------------------------------
44 void StarPythia6::PyTune( Int_t tune ){ ::PyTune(tune); }
45 void StarPythia6::PyStat( Int_t stat ){ ::PyStat(stat); }
46 void StarPythia6::PyList( Int_t list ){ ::PyList(list); }
47 
48 // ----------------------------------------------------------------------------
50 {
51 
52  //
53  // Create a new event record for either pp or ep events
54  //
55  if ( mBlue == "proton" ) mEvent = new StarGenPPEvent();
56  else mEvent = new StarGenEPEvent();
57 
63  pydat3().mdcy(102,1)=0; // PI0 111
64  pydat3().mdcy(106,1)=0; // PI+ 211
65  pydat3().mdcy(109,1)=0; // ETA 221
66  pydat3().mdcy(116,1)=0; // K+ 321
67  pydat3().mdcy(112,1)=0; // K_SHORT 310
68  pydat3().mdcy(105,1)=0; // K_LONG 130
69  pydat3().mdcy(164,1)=0; // LAMBDA0 3122
70  pydat3().mdcy(167,1)=0; // SIGMA0 3212
71  pydat3().mdcy(162,1)=0; // SIGMA- 3112
72  pydat3().mdcy(169,1)=0; // SIGMA+ 3222
73  pydat3().mdcy(172,1)=0; // Xi- 3312
74  pydat3().mdcy(174,1)=0; // Xi0 3322
75  pydat3().mdcy(176,1)=0; // OMEGA- 3334
76 
77 
78  // Remap the ROOT names to pythia6 names
79  std::map< TString, string > particle;
80  particle["proton"]="p";
81  particle["e-"]="e-";
82 
83  // For frames other than CMS, we need to set the beam momenta in pyjets
84  if ( mFrame=="3MOM" || mFrame=="4MOM" || mFrame=="5MOM" )
85  {
86  for ( Int_t i=0;i<3;i++ ) pyjets().p(1,i+1) = mBlueMomentum[i];
87  for ( Int_t i=0;i<3;i++ ) pyjets().p(2,i+1) = mYellMomentum[i];
88  }
89  if ( mFrame=="4MOM" || mFrame=="5MOM" )
90  { Int_t i=3;
91  pyjets().p(1,i+1) = mBlueMomentum[i];
92  pyjets().p(2,i+1) = mYellMomentum[i];
93  }
94 
95 
96  // Initialize pythia
97  PyInit( mFrame.Data(), particle[mBlue], particle[mYell], mRootS );
98 
99  return StMaker::Init();
100 }
101 // ----------------------------------------------------------------------------
102 //
103 // ----------------------------------------------------------------------------
105 {
106 
107  // Generate the event
108  PyEvnt();
109 
110  // Blue beam is a proton, running PP
111  if ( mBlue == "proton" ) FillPP( mEvent );
112  // Otherwise, runnin EP
113  else FillEP( mEvent );
114 
115  // Loop over all particles in the event
116  mNumberOfParticles = pyjets().n;
117  for ( Int_t idx=1; idx<=mNumberOfParticles; idx++ )
118  {
119 
120  Int_t id = pyjets().k(idx,2);
121  Int_t stat = mStatusCode[ pyjets().k(idx,1) ];
122  if ( !stat ) {
124  }
125  Int_t m1 = pyjets().k(idx,3);
126  Int_t m2 = 0;
127  Int_t d1 = pyjets().k(idx,4);
128  Int_t d2 = pyjets().k(idx,5);
129  Double_t px = pyjets().p(idx,1);
130  Double_t py = pyjets().p(idx,2);
131  Double_t pz = pyjets().p(idx,3);
132  Double_t E = pyjets().p(idx,4);
133  Double_t M = pyjets().p(idx,5);
134  Double_t vx = pyjets().v(idx,1);
135  Double_t vy = pyjets().v(idx,2);
136  Double_t vz = pyjets().v(idx,3);
137  Double_t vt = pyjets().v(idx,4);
138 
139  mEvent -> AddParticle( stat, id, m1, m2, d1, d2, px, py, pz, E, M, vx, vy, vz, vt );
140 
141  }
142 
143  return kStOK;
144 }
145 // ----------------------------------------------------------------------------
146 //
147 // ----------------------------------------------------------------------------
149 {
150 
151  // Fill the event-wise information
152  StarGenPPEvent *myevent = (StarGenPPEvent *)event;
153 
154  myevent -> idBlue = pypars().msti(11);
155  myevent -> idYell = pypars().msti(12);
156  myevent -> process = pysubs().msel;
157  myevent -> subprocess = pypars().msti(1);
158 
159  myevent -> idParton1 = pypars().msti(15);
160  myevent -> idParton2 = pypars().msti(16);
161  myevent -> xParton1 = pypars().pari(31);
162  myevent -> xParton2 = pypars().pari(32);
163  myevent -> xPdf1 = 0;
164  myevent -> xPdf2 = 0;
165  myevent -> Q2fac = pypars().pari(22);
166  myevent -> Q2ren = 0.;
167  myevent -> valence1 = -1;
168  myevent -> valence2 = -1;
169 
170  myevent -> sHat = pypars().pari(14);
171  myevent -> tHat = pypars().pari(15);
172  myevent -> uHat = pypars().pari(16);
173  myevent -> ptHat = pypars().pari(17);
174  myevent -> thetaHat = TMath::ACos( pypars().pari(41) );
175  myevent -> phiHat = -999;
176 
177  myevent -> weight = pypars().pari(7);
178 
179  myevent -> mstu72 = pydat1().mstu(72);
180  myevent -> mstu73 = pydat1().mstu(73);
181 
182 }
183 // ----------------------------------------------------------------------------
184 //
185 // ----------------------------------------------------------------------------
187 {
188 
189  // Fill the event-wise information
190  StarGenEPEvent *myevent = (StarGenEPEvent *)event;
191 
192  myevent -> idBlue = pypars().msti(11);
193  myevent -> idYell = pypars().msti(12);
194  myevent -> process = pysubs().msel;
195  myevent -> subprocess = pypars().msti(1);
196 
197  // ID of the colliding parton. If MSTI(15) isn't a parton, grab MSTI(16)
198  Int_t id = pypars().msti(15);
199  Double_t x = pypars().pari(31);
200  if ( TMath::Abs(id)>6 ) { id = pypars().msti(16); x = pypars().pari(32); }
201  myevent -> idParton = id;
202  myevent -> xParton = x;
203 
204  myevent -> xPdf = 0;
205 
206  myevent -> Q2 = pypars().pari(22);
207  myevent -> weight = pypars().pari(7);
208 
209 }
210 
212 {
213  StarGenStats stats("Pythia6Stats","Pythia 6 Run Statistics");
214  PyStat(1);
215 
216  stats.nTried = pyint5().ngen(0,1);
217  stats.nSelected = pyint5().ngen(0,2);
218  stats.nAccepted = pyint5().ngen(0,3);
219  stats.sigmaGen = pyint5().xsec(0,3); //xsection
220  // ... we may want to extend StarGenStats to add these, if they are useful ...
221  // stats.xxxxxxxx = pyint5().xsec(0,1); //estimated maximum differential cross section
222  // stats.xxxxxxxxxxxx = pyint5().xsec(0,2); //gives the sum of differential cross sections phase-space points evaluated so far
223 
224  stats.nFilterSeen = stats.nAccepted;
225  stats.nFilterAccept = stats.nAccepted;
226 
227  //stats.Dump();
228 
229  // Return a copy of the class we just created
230  return stats;
231 }
void FillEP(StarGenEvent *event)
(Optional) Method to fill a DIS event
void PyTune(Int_t tune)
Calls the pytune function.
Definition: StarPythia6.cxx:44
Double_t mRootS
CMS energy or incident beam momentum for fixed target collisions.
static StarRandom & Instance()
Obtain the single instance of the random number generator.
Definition: StarRandom.cxx:87
TString mYell
Name of the yellow beam particle (-z)
PyDat3_t & pydat3()
Returns a reference to the /PYDAT3/ common block.
Definition: StarPythia6.h:65
Event record class tailored to PP kinematics.
TLorentzVector mYellMomentum
4-Momentum of the yellow beam
ABC for defining event generator interfaces.
Definition: StarGenerator.h:34
StarGenEvent * mEvent
Generated event.
void FillPP(StarGenEvent *event)
(Optional) Method to fill a PP event
Int_t Init()
Definition: StarPythia6.cxx:49
Int_t Generate()
PyInt5_t & pyint5()
Returns a reference to the /PYINT5/ common block.
Definition: StarPythia6.h:69
void PyStat(Int_t stat)
Calls the pystat function.
Definition: StarPythia6.cxx:45
TLorentzVector mBlueMomentum
4-Momentum of the blue beam
PyDat1_t & pydat1()
Returns a reference to the /PYDAT1/ common block.
Definition: StarPythia6.h:63
PySubs_t & pysubs()
Returns a reference to the /PYSUBS/ common block.
Definition: StarPythia6.h:61
Double_t flat() const
Return a random number uniformly distributed between 0 and 1.
Definition: StarRandom.cxx:68
End of run statistics for event generators.
Definition: StarGenStats.h:21
StarGenStats Stats()
Create and retrieve end-of-run statistics.
PyJets_t & pyjets()
Returns a reference to the /PYJETS/ common block.
Definition: StarPythia6.h:59
Base class for event records.
Definition: StarGenEvent.h:81
Interface to pythia 6.
Definition: StarPythia6.h:48
Event record class tailored to DIS kinemaics.
Definition: Stypes.h:40
PyPars_t & pypars()
Returns a reference to the /PYPARS/ common block.
Definition: StarPythia6.h:67
TString mFrame
Frame of the collision, i.e. CMS, FIXT, 3MOM, 4MOM, 5MOM.
TString mBlue
Name of the blue beam particle (+z)
void PyList(Int_t list)
Calls the pylist function.
Definition: StarPythia6.cxx:46