1 #include "StarPythia8.h"
4 #include "TDatabasePDG.h"
5 #include "TParticlePDG.h"
7 #include "StarGenerator/UTIL/StarRandom.h"
8 #include "StarGenerator/EVENT/StarGenPPEvent.h"
9 #include "StarGenerator/EVENT/StarGenEPEvent.h"
13 #ifndef Pythia8_version
14 #error "Pythia8_version is not defined"
26 StarPythia8::StarPythia8(
const char *name) :
StarGenerator(name)
33 TString path =
"StRoot/StarGenerator/"; path+= Pythia8_version; path+=
"/xmldoc/";
36 if (!in.good()) { path =
"$(STAR)/"+path; }
37 path = gSystem->ExpandPathName(path.Data());
41 Info(GetName(),Form(
"MC version is %s data at %s",Pythia8_version,path.Data()));
42 Info(GetName(),Form(
"Configuration files at %s",path.Data()));
45 mPythia -> setRndmEnginePtr(
new PyRand() );
57 map<TString,TString> particles;
58 particles[
"electron"]=
"e-";
59 particles[
"proton"]=
"proton";
61 TString myBlue = particles[
mBlue];
if ( myBlue ==
"" ) myBlue=
mBlue;
62 TString myYell = particles[
mYell];
if ( myYell ==
"" ) myYell=
mYell;
68 TDatabasePDG &pdg = (*TDatabasePDG::Instance());
69 TParticlePDG *blue = pdg.GetParticle(myBlue); assert(blue);
70 TParticlePDG *yell = pdg.GetParticle(myYell); assert(yell);
72 if (
mFrame ==
"CMS" ) InitCMS ( blue->PdgCode(), yell->PdgCode() );
73 if (
mFrame ==
"FIXT" ) InitFIXT( blue->PdgCode(), yell->PdgCode() );
74 if (
mFrame ==
"3MOM" ) Init3MOM( blue->PdgCode(), yell->PdgCode() );
75 if (
mFrame ==
"4MOM" ) Init4MOM( blue->PdgCode(), yell->PdgCode() );
76 if (
mFrame ==
"5MOM" ) Init5MOM( blue->PdgCode(), yell->PdgCode() );
101 return StMaker::Init();
123 mNumberOfParticles =
event.size() - 1;
130 for (
int idx=1; idx <= mNumberOfParticles; idx++ )
133 int id =
event[idx].id();
134 int stat =
event.statusHepMC(idx);
135 int mother1 =
event[idx].mother1();
136 int mother2 =
event[idx].mother2();
137 int daughter1 =
event[idx].daughter1();
138 int daughter2 =
event[idx].daughter2();
139 double px =
event[idx].px();
140 double py =
event[idx].py();
141 double pz =
event[idx].pz();
142 double energy =
event[idx].e();
143 double mass =
event[idx].m();
144 double vx =
event[idx].xProd();
145 double vy =
event[idx].yProd();
146 double vz =
event[idx].zProd();
147 double vt =
event[idx].tProd();
149 mEvent -> AddParticle( stat,
id, mother1, mother2, daughter1, daughter2, px, py, pz, energy, mass, vx, vy, vz, vt );
166 event -> idBlue = info.idA();
167 event -> idYell = info.idB();
168 event -> process = info.code();
169 event -> subprocess = (info.hasSub())? 0 : info.codeSub();
172 double x = info.x1();
173 double xPdf = info.pdf1();
174 int valence = info.isValence1();
175 if ( TMath::Abs(
id)>6 )
180 valence = info.isValence2();
183 event -> idParton = id;
184 event -> xParton = x;
185 event -> xPdf = xPdf;
187 event -> Q2 = info.Q2Fac();
188 event -> valence = valence;
197 event -> weight = info.weight();
211 event -> idBlue = info.idA();
212 event -> idYell = info.idB();
213 event -> process = info.code();
214 event -> subprocess = (info.hasSub())? 0 : info.codeSub();
216 event -> idParton1 = info.id1();
217 event -> idParton2 = info.id2();
218 event -> xParton1 = info.x1();
219 event -> xParton2 = info.x2();
220 event -> xPdf1 = info.pdf1();
221 event -> xPdf2 = info.pdf2();
222 event -> Q2fac = info.Q2Fac();
223 event -> Q2ren = info.Q2Ren();
224 event -> valence1 = info.isValence1();
225 event -> valence2 = info.isValence2();
227 event -> sHat = info.sHat();
228 event -> tHat = info.tHat();
229 event -> uHat = info.uHat();
230 event -> ptHat = info.pTHat();
231 event -> thetaHat = info.thetaHat();
232 event -> phiHat = info.phiHat();
234 event -> weight = info.weight();
249 stats.nTried = info.nTried();
250 stats.nSelected = info.nSelected();
251 stats.nAccepted = info.nAccepted();
252 stats.sigmaGen = info.sigmaGen();
253 stats.sigmaErr = info.sigmaErr();
254 stats.sumWeightGen = info.weightSum();
256 stats.nFilterSeen = stats.nAccepted;
257 stats.nFilterAccept = stats.nAccepted;
268 int StarPythia8::InitCMS(
int blue,
int yell )
270 mPythia->init( blue, yell,
mRootS );
274 int StarPythia8::InitFIXT(
int blue,
int yell )
277 mPythia->init( blue, yell,
mRootS, 0.0 );
279 mPythia->init( blue, yell, 0,
mRootS );
284 int StarPythia8::Init3MOM(
int blue,
int yell )
287 cout <<
"Initializing 3MOM: " << endl;
292 mPythia -> init( blue, yell,
299 int StarPythia8::Init4MOM(
int blue,
int yell )
301 return Init3MOM( blue, yell );
304 int StarPythia8::Init5MOM(
int blue,
int yell )
306 return Init3MOM( blue, yell );
Double_t mRootS
CMS energy or incident beam momentum for fixed target collisions.
int Generate()
Generate one event.
static StarRandom & Instance()
Obtain the single instance of the random number generator.
TString mYell
Name of the yellow beam particle (-z)
Event record class tailored to PP kinematics.
TLorentzVector mYellMomentum
4-Momentum of the yellow beam
ABC for defining event generator interfaces.
StarGenEvent * mEvent
Generated event.
void FillEP(StarGenEvent *event)
(Optional) Method to fill a DIS event
TLorentzVector mBlueMomentum
4-Momentum of the blue beam
Double_t flat() const
Return a random number uniformly distributed between 0 and 1.
End of run statistics for event generators.
Base class for event records.
virtual const char * GetName() const
special overload
Event record class tailored to DIS kinemaics.
TString mFrame
Frame of the collision, i.e. CMS, FIXT, 3MOM, 4MOM, 5MOM.
TString mBlue
Name of the blue beam particle (+z)
void FillPP(StarGenEvent *event)
(Optional) Method to fill a PP event
int Init()
Initialize the event generator.
StarGenStats Stats()
Return end-of-run statistics.
void Set(const char *s)
Pass a string to Pythia8::Pythia::readString(), for user configuration.