1 #include "StarPythia8.h"
4 #include "StarGenerator/UTIL/StarParticleData.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+=
"/share/Pythia8/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;
65 LOG_INFO <<
"Default configuration with the Detroit Tune" << endm;
66 Set(
"MultipartonInteractions:pT0Ref=1.40");
67 Set(
"MultipartonInteractions:ecmPow=0.135");
68 Set(
"MultipartonInteractions:coreRadius=0.56");
69 Set(
"MultipartonInteractions:coreFraction=0.78");
70 Set(
"ColourReconnection:range=5.4");
77 TParticlePDG *blue = pdg.
GetParticle(myBlue); assert(blue);
78 TParticlePDG *yell = pdg.
GetParticle(myYell); assert(yell);
80 if (
mFrame ==
"CMS" ) InitCMS ( blue->PdgCode(), yell->PdgCode() );
81 if (
mFrame ==
"FIXT" ) InitFIXT( blue->PdgCode(), yell->PdgCode() );
82 if (
mFrame ==
"3MOM" ) Init3MOM( blue->PdgCode(), yell->PdgCode() );
83 if (
mFrame ==
"4MOM" ) Init4MOM( blue->PdgCode(), yell->PdgCode() );
84 if (
mFrame ==
"5MOM" ) Init5MOM( blue->PdgCode(), yell->PdgCode() );
89 if (
mBlue ==
"proton" ) {
90 LOG_INFO <<
"pp or pA mode detected" << endm;
93 LOG_INFO <<
"AA (or eA) mode detected... event record will still be a pp event record." << endm;
106 Set(
"3122:onMode=0");
107 Set(
"3112:onMode=0");
108 Set(
"3222:onMode=0");
109 Set(
"3212:onMode=0");
110 Set(
"3312:onMode=0");
111 Set(
"3322:onMode=0");
112 Set(
"3334:onMode=0");
115 return StMaker::Init();
137 mNumberOfParticles =
event.size() - 1;
144 for (
int idx=1; idx <= mNumberOfParticles; idx++ )
147 int id =
event[idx].id();
149 int stat =
event[idx].statusHepMC();
150 int mother1 =
event[idx].mother1();
151 int mother2 =
event[idx].mother2();
152 int daughter1 =
event[idx].daughter1();
153 int daughter2 =
event[idx].daughter2();
154 double px =
event[idx].px();
155 double py =
event[idx].py();
156 double pz =
event[idx].pz();
157 double energy =
event[idx].e();
158 double mass =
event[idx].m();
159 double vx =
event[idx].xProd();
160 double vy =
event[idx].yProd();
161 double vz =
event[idx].zProd();
162 double vt =
event[idx].tProd();
164 mEvent -> AddParticle( stat,
id, mother1, mother2, daughter1, daughter2, px, py, pz, energy, mass, vx, vy, vz, vt );
181 event -> idBlue = info.idA();
182 event -> idYell = info.idB();
183 event -> process = info.code();
184 event -> subprocess = (info.hasSub())? 0 : info.codeSub();
187 double x = info.x1();
188 double xPdf = info.pdf1();
189 int valence = info.isValence1();
190 if ( TMath::Abs(
id)>6 )
195 valence = info.isValence2();
198 event -> idParton = id;
199 event -> xParton = x;
200 event -> xPdf = xPdf;
202 event -> Q2 = info.Q2Fac();
203 event -> valence = valence;
212 event -> weight = info.weight();
226 event -> idBlue = info.idA();
227 event -> idYell = info.idB();
228 event -> process = info.code();
229 event -> subprocess = (info.hasSub())? 0 : info.codeSub();
231 event -> idParton1 = info.id1();
232 event -> idParton2 = info.id2();
233 event -> xParton1 = info.x1();
234 event -> xParton2 = info.x2();
235 event -> xPdf1 = info.pdf1();
236 event -> xPdf2 = info.pdf2();
237 event -> Q2fac = info.Q2Fac();
238 event -> Q2ren = info.Q2Ren();
239 event -> valence1 = info.isValence1();
240 event -> valence2 = info.isValence2();
242 event -> sHat = info.sHat();
243 event -> tHat = info.tHat();
244 event -> uHat = info.uHat();
245 event -> ptHat = info.pTHat();
246 event -> thetaHat = info.thetaHat();
247 event -> phiHat = info.phiHat();
249 event -> weight = info.weight();
264 stats.nTried = info.nTried();
265 stats.nSelected = info.nSelected();
266 stats.nAccepted = info.nAccepted();
267 stats.sigmaGen = info.sigmaGen();
268 stats.sigmaErr = info.sigmaErr();
269 stats.sumWeightGen = info.weightSum();
271 stats.nFilterSeen = stats.nAccepted;
272 stats.nFilterAccept = stats.nAccepted;
283 int StarPythia8::InitCMS(
int blue,
int yell )
285 Set( Form(
"Beams:idA=%i",blue) );
286 Set( Form(
"Beams:idB=%i",yell) );
287 Set(
"Beams:frameType=1");
294 int StarPythia8::InitFIXT(
int blue,
int yell )
296 Set( Form(
"Beams:idA=%i",blue) );
297 Set( Form(
"Beams:idB=%i",yell) );
298 Set(
"Beams:frameType=2");
301 Set( Form(
"Beams:eA=%f", eA) );
302 Set( Form(
"Beams:eB=%f", eB) );
307 int StarPythia8::Init3MOM(
int blue,
int yell )
309 LOG_INFO <<
"Init3/4/5MOM not implemented yet" << endm;
323 int StarPythia8::Init4MOM(
int blue,
int yell )
325 return Init3MOM( blue, yell );
328 int StarPythia8::Init5MOM(
int blue,
int yell )
330 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.
ABC for defining event generator interfaces.
StarGenEvent * mEvent
Generated event.
static StarParticleData & instance()
Returns a reference to the single instance of this class.
TParticlePDG * GetParticle(const Char_t *name)
Get a particle by name.
Interface to PDG information.
void FillEP(StarGenEvent *event)
(Optional) Method to fill a DIS event
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.