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 #include "TGenericTable.h"
14 
15 StMaker *_maker = 0;
16 
17 TGenericTable *regtable( const Char_t *type, const Char_t *name, void *address )
18 {
19  TGenericTable *table = new TGenericTable(type,name);
20  table->Adopt( 1, address );
21  _maker -> AddData( table, ".const" );
22  return table;
23 };
24 
25 // ----------------------------------------------------------------------------
26 // Remap pythia's random number generator to StarRandom
27 extern "C" {
28  Double_t pyrstar_( Int_t *idummy ){ return StarRandom::Instance().flat(); };
29 };
30 // ----------------------------------------------------------------------------
31 // ----------------------------------------------------------------------------
32 // ----------------------------------------------------------------------------
33 StarPythia6::StarPythia6( const Char_t *name ) : StarGenerator(name)
34 {
35 
36  // Setup a map between pythia6 status codes and the HepMC status codes
37  // used in StarEvent
38  mStatusCode[0] = StarGenParticle::kNull;
39  mStatusCode[1] = StarGenParticle::kFinal;
40  const Int_t decays[] = { 2, 3, 4, 5, 11, 12, 13, 14, 15 };
41  for ( UInt_t i=0;i<sizeof(decays)/sizeof(Int_t); i++ )
42  {
43  mStatusCode[decays[i]]=StarGenParticle::kDecayed;
44  }
45  const Int_t docum[] = { 21, 31, 32, 41, 42, 51, 52 };
46  for ( UInt_t i=0;i<sizeof(docum)/sizeof(Int_t); i++ )
47  {
48  mStatusCode[docum[i]]=StarGenParticle::kDocumentation;
49  }
50 
51  _maker = this;
52 
53  regtable("PyJets_t", "pyjets", address_of_pyjets() );
54  regtable("PySubs_t", "pysubs", address_of_pysubs() );
55  regtable("PyDat1_t", "pydat1", address_of_pydat1() );
56  regtable("PyDat3_t", "pydat3", address_of_pydat3() );
57  regtable("PyPars_t", "pypars", address_of_pypars() );
58  regtable("PyInt5_t", "pyint5", address_of_pyint5() );
59 
60 }
61 // ----------------------------------------------------------------------------
62 // ----------------------------------------------------------------------------
63 // ----------------------------------------------------------------------------
64 void StarPythia6::PyTune( Int_t tune ){ ::PyTune(tune); }
65 void StarPythia6::PyStat( Int_t stat ){ ::PyStat(stat); }
66 void StarPythia6::PyList( Int_t list ){ ::PyList(list); }
67 
68 // ----------------------------------------------------------------------------
69 Int_t StarPythia6::Init()
70 {
71 
72  //
73  // Create a new event record for either pp or ep events
74  //
75  if ( mBlue == "proton" ) mEvent = new StarGenPPEvent();
76  else mEvent = new StarGenEPEvent();
77 
83  pydat3().mdcy(102,1)=0; // PI0 111
84  pydat3().mdcy(106,1)=0; // PI+ 211
85  pydat3().mdcy(109,1)=0; // ETA 221
86  pydat3().mdcy(116,1)=0; // K+ 321
87  pydat3().mdcy(112,1)=0; // K_SHORT 310
88  pydat3().mdcy(105,1)=0; // K_LONG 130
89  pydat3().mdcy(164,1)=0; // LAMBDA0 3122
90  pydat3().mdcy(167,1)=0; // SIGMA0 3212
91  pydat3().mdcy(162,1)=0; // SIGMA- 3112
92  pydat3().mdcy(169,1)=0; // SIGMA+ 3222
93  pydat3().mdcy(172,1)=0; // Xi- 3312
94  pydat3().mdcy(174,1)=0; // Xi0 3322
95  pydat3().mdcy(176,1)=0; // OMEGA- 3334
96 
97 
98  // Remap the ROOT names to pythia6 names
99  std::map< TString, string > particle;
100  particle["proton"]="p";
101  particle["e-"]="e-";
102 
103  // For frames other than CMS, we need to set the beam momenta in pyjets
104  if ( mFrame=="3MOM" || mFrame=="4MOM" || mFrame=="5MOM" )
105  {
106  for ( Int_t i=0;i<3;i++ ) pyjets().p(1,i+1) = mBlueMomentum[i];
107  for ( Int_t i=0;i<3;i++ ) pyjets().p(2,i+1) = mYellMomentum[i];
108  }
109  if ( mFrame=="4MOM" || mFrame=="5MOM" )
110  { Int_t i=3;
111  pyjets().p(1,i+1) = mBlueMomentum[i];
112  pyjets().p(2,i+1) = mYellMomentum[i];
113  }
114 
115 
116  // Initialize pythia
117  PyInit( mFrame.Data(), particle[mBlue], particle[mYell], mRootS );
118 
119  return StMaker::Init();
120 }
121 // ----------------------------------------------------------------------------
122 //
123 // ----------------------------------------------------------------------------
124 Int_t StarPythia6::Generate()
125 {
126 
127  // Generate the event
128  PyEvnt();
129 
130  // Blue beam is a proton, running PP
131  if ( mBlue == "proton" ) FillPP( mEvent );
132  // Otherwise, runnin EP
133  else FillEP( mEvent );
134 
135  // Loop over all particles in the event
136  mNumberOfParticles = pyjets().n;
137  for ( Int_t idx=1; idx<=mNumberOfParticles; idx++ )
138  {
139 
140  Int_t id = pyjets().k(idx,2);
141  Int_t stat = mStatusCode[ pyjets().k(idx,1) ];
142  if ( !stat ) {
144  }
145  Int_t m1 = pyjets().k(idx,3);
146  Int_t m2 = 0;
147  Int_t d1 = pyjets().k(idx,4);
148  Int_t d2 = pyjets().k(idx,5);
149  Double_t px = pyjets().p(idx,1);
150  Double_t py = pyjets().p(idx,2);
151  Double_t pz = pyjets().p(idx,3);
152  Double_t E = pyjets().p(idx,4);
153  Double_t M = pyjets().p(idx,5);
154  Double_t vx = pyjets().v(idx,1);
155  Double_t vy = pyjets().v(idx,2);
156  Double_t vz = pyjets().v(idx,3);
157  Double_t vt = pyjets().v(idx,4);
158 
159  mEvent -> AddParticle( stat, id, m1, m2, d1, d2, px, py, pz, E, M, vx, vy, vz, vt );
160 
161  }
162 
163  return kStOK;
164 }
165 // ----------------------------------------------------------------------------
166 //
167 // ----------------------------------------------------------------------------
168 void StarPythia6::FillPP( StarGenEvent *event )
169 {
170 
171  // Fill the event-wise information
172  StarGenPPEvent *myevent = (StarGenPPEvent *)event;
173 
174  myevent -> idBlue = pypars().msti(11);
175  myevent -> idYell = pypars().msti(12);
176  myevent -> process = pysubs().msel;
177  myevent -> subprocess = pypars().msti(1);
178 
179  myevent -> idParton1 = pypars().msti(15);
180  myevent -> idParton2 = pypars().msti(16);
181  myevent -> xParton1 = pypars().pari(31);
182  myevent -> xParton2 = pypars().pari(32);
183  myevent -> xPdf1 = 0;
184  myevent -> xPdf2 = 0;
185  myevent -> Q2fac = pypars().pari(22);
186  myevent -> Q2ren = 0.;
187  myevent -> valence1 = -1;
188  myevent -> valence2 = -1;
189 
190  myevent -> sHat = pypars().pari(14);
191  myevent -> tHat = pypars().pari(15);
192  myevent -> uHat = pypars().pari(16);
193  myevent -> ptHat = pypars().pari(17);
194  myevent -> thetaHat = TMath::ACos( pypars().pari(41) );
195  myevent -> phiHat = -999;
196 
197  myevent -> weight = pypars().pari(7);
198 
199  myevent -> mstu72 = pydat1().mstu(72);
200  myevent -> mstu73 = pydat1().mstu(73);
201 
202 }
203 // ----------------------------------------------------------------------------
204 //
205 // ----------------------------------------------------------------------------
206 void StarPythia6::FillEP( StarGenEvent *event )
207 {
208 
209  // Fill the event-wise information
210  StarGenEPEvent *myevent = (StarGenEPEvent *)event;
211 
212  myevent -> idBlue = pypars().msti(11);
213  myevent -> idYell = pypars().msti(12);
214  myevent -> process = pysubs().msel;
215  myevent -> subprocess = pypars().msti(1);
216 
217  // ID of the colliding parton. If MSTI(15) isn't a parton, grab MSTI(16)
218  Int_t id = pypars().msti(15);
219  Double_t x = pypars().pari(31);
220  if ( TMath::Abs(id)>6 ) { id = pypars().msti(16); x = pypars().pari(32); }
221  myevent -> idParton = id;
222  myevent -> xParton = x;
223 
224  myevent -> xPdf = 0;
225 
226  myevent -> Q2 = pypars().pari(22);
227  myevent -> weight = pypars().pari(7);
228 
229 }
230 
232 {
233  StarGenStats stats("Pythia6Stats","Pythia 6 Run Statistics");
234  PyStat(1);
235 
236  stats.nTried = pyint5().ngen(0,1);
237  stats.nSelected = pyint5().ngen(0,2);
238  stats.nAccepted = pyint5().ngen(0,3);
239  stats.sigmaGen = pyint5().xsec(0,3); //xsection
240  // ... we may want to extend StarGenStats to add these, if they are useful ...
241  // stats.xxxxxxxx = pyint5().xsec(0,1); //estimated maximum differential cross section
242  // stats.xxxxxxxxxxxx = pyint5().xsec(0,2); //gives the sum of differential cross sections phase-space points evaluated so far
243 
244  stats.nFilterSeen = stats.nAccepted;
245  stats.nFilterAccept = stats.nAccepted;
246 
247  //stats.Dump();
248 
249  // Return a copy of the class we just created
250  return stats;
251 }
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
virtual void Adopt(Int_t n, void *array)
Definition: TTable.cxx:1107