StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StarPythia8.cxx
1 #include "StarPythia8.h"
2 ClassImp(StarPythia8);
3 
4 #include "StarGenerator/UTIL/StarParticleData.h"
5 #include "TParticlePDG.h"
6 
7 #include "StarGenerator/UTIL/StarRandom.h"
8 #include "StarGenerator/EVENT/StarGenPPEvent.h"
9 #include "StarGenerator/EVENT/StarGenEPEvent.h"
10 #include "TString.h"
11 #include "TSystem.h"
12 
13 #ifndef Pythia8_version
14 #error "Pythia8_version is not defined"
15 #endif
16 
17 // ----------------------------------------------------------------------------
18 // Remap pythia's random number generator to StarRandom
19 class PyRand : public Pythia8::RndmEngine {
20 public:
21  double flat() { return StarRandom::Instance().flat(); }
22 };
23 // ----------------------------------------------------------------------------
24 // ----------------------------------------------------------------------------
25 // ----------------------------------------------------------------------------
26 StarPythia8::StarPythia8(const char *name) : StarGenerator(name)
27 {
28 
29  // Create the new instance of pythia, specifying the path to the
30  // auxilliary files required by pythia.
31  // NOTE: When adding new versions of Pythia8, we need to specify
32  // the version of the code in the source code
33  TString path = "StRoot/StarGenerator/"; path+= Pythia8_version; path+="/share/Pythia8/xmldoc/";
34  {
35  ifstream in(path);
36  if (!in.good()) { path = "$(STAR)/"+path; }
37  path = gSystem->ExpandPathName(path.Data());
38  }
39 
40 
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()));
43 
44  mPythia = new Pythia8::Pythia( path.Data() );
45  mPythia -> setRndmEnginePtr( new PyRand() );
46 
47 }
48 // ----------------------------------------------------------------------------
49 // ----------------------------------------------------------------------------
50 // ----------------------------------------------------------------------------
52 {
53 
54  //
55  // Remap pythia8 to pdg particles
56  //
57  map<TString,TString> particles;
58  particles["electron"]="e-";
59  particles["proton"]="proton";
60 
61  TString myBlue = particles[mBlue]; if ( myBlue == "" ) myBlue=mBlue;
62  TString myYell = particles[mYell]; if ( myYell == "" ) myYell=mYell;
63 
64  //
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");
71 
72  //
73  // Initialize pythia based on the frame and registered beam momenta
74  //
76 
77  TParticlePDG *blue = pdg.GetParticle(myBlue); assert(blue);
78  TParticlePDG *yell = pdg.GetParticle(myYell); assert(yell);
79  //
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() );
85  //
86  // Setup event record based upon the beam species
87  //
88  mEvent = new StarGenPPEvent();
89  if ( mBlue == "proton" ) {
90  LOG_INFO << "pp or pA mode detected" << endm;
91  }
92  else {
93  LOG_INFO << "AA (or eA) mode detected... event record will still be a pp event record." << endm;
94  }
95  //
96  // Make several particles stable which may cross the beam pipe,
97  // and so the simulation package must be allowed to decide to
98  // decay them or not.
99  //
100  Set("111:onMode=0"); // pi0 stable to permit mother/daughter in star record
101  Set("211:onMode=0"); // pi+/- stable
102  Set("221:onMode=0"); // eta stable
103  Set("321:onMode=0"); // K+/- stable
104  Set("310:onMode=0"); // K short
105  Set("130:onMode=0"); // K long
106  Set("3122:onMode=0"); // Lambda 0
107  Set("3112:onMode=0"); // Sigma -
108  Set("3222:onMode=0"); // Sigma +
109  Set("3212:onMode=0"); // Sigma 0
110  Set("3312:onMode=0"); // Xi -
111  Set("3322:onMode=0"); // Xi 0
112  Set("3334:onMode=0"); // Omega -
113 
114 
115  return StMaker::Init();
116  //
117 }
118 // ----------------------------------------------------------------------------
120 {
121 
122  //
123  // Generate the event. This is where the action happens. The rest of this
124  // method is just copying data from pythia into the event record.
125  //
126  mPythia -> next();
127 
128  //
129  // Get the pythis event record
130  //
131  Pythia8::Event &event = mPythia->event;
132 
133  //
134  // Get the number of particles in the event. Pythia 8 include a "0th" entry,
135  // which represents the system in the event record. We will skip over this.
136  //
137  mNumberOfParticles = event.size() - 1;
138 
139 
140  if ( mBlue == "proton" ) FillPP( mEvent );
141  else FillEP( mEvent );
142 
143  // Loop over all particles in the event
144  for ( int idx=1; idx <= mNumberOfParticles; idx++ )
145  {
146 
147  int id = event[idx].id();
148 // int stat = event.statusHepMC(idx);
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(); // mm
160  double vy = event[idx].yProd(); // mm
161  double vz = event[idx].zProd(); // mm
162  double vt = event[idx].tProd(); // mm/c
163 
164  mEvent -> AddParticle( stat, id, mother1, mother2, daughter1, daughter2, px, py, pz, energy, mass, vx, vy, vz, vt );
165 
166  }
167 
168  return kStOK;
169 }
170 // ----------------------------------------------------------------------------
171 // ----------------------------------------------------------------------------
172 // ----------------------------------------------------------------------------
173 void StarPythia8::FillEP( StarGenEvent *myevent )
174 {
175  //
176  // Fill event-wise information
177  //
178  StarGenEPEvent *event = (StarGenEPEvent *)myevent;
179  const Pythia8::Info &info = mPythia->info;
180 
181  event -> idBlue = info.idA(); // Beam particle
182  event -> idYell = info.idB(); // Beam particle
183  event -> process = info.code();
184  event -> subprocess = (info.hasSub())? 0 : info.codeSub();
185 
186  int id = info.id1();
187  double x = info.x1();
188  double xPdf = info.pdf1();
189  int valence = info.isValence1();
190  if ( TMath::Abs(id)>6 )
191  {
192  id = info.id2();
193  x = info.x2();
194  xPdf = info.pdf2();
195  valence = info.isValence2();
196  }
197 
198  event -> idParton = id;
199  event -> xParton = x;
200  event -> xPdf = xPdf;
201 
202  event -> Q2 = info.Q2Fac(); // factorization scale
203  event -> valence = valence;
204 
205 // event -> sHat = info.sHat();
206 // event -> tHat = info.tHat();
207 // event -> uHat = info.uHat();
208 // event -> ptHat = info.pTHat();
209 // event -> thetaHat = info.thetaHat();
210 // event -> phiHat = info.phiHat();
211 
212  event -> weight = info.weight();
213 
214 }
215 // ----------------------------------------------------------------------------
216 // ----------------------------------------------------------------------------
217 // ----------------------------------------------------------------------------
218 void StarPythia8::FillPP( StarGenEvent *myevent )
219 {
220  //
221  // Fill event-wise information
222  //
223  StarGenPPEvent *event = (StarGenPPEvent *)myevent;
224  const Pythia8::Info &info = mPythia->info;
225 
226  event -> idBlue = info.idA(); // Beam particle
227  event -> idYell = info.idB(); // Beam particle
228  event -> process = info.code();
229  event -> subprocess = (info.hasSub())? 0 : info.codeSub();
230 
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(); // factorization scale
238  event -> Q2ren = info.Q2Ren(); // renormalization scale
239  event -> valence1 = info.isValence1();
240  event -> valence2 = info.isValence2();
241 
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();
248 
249  event -> weight = info.weight();
250 
251 }
252 // ----------------------------------------------------------------------------
253 // ----------------------------------------------------------------------------
254 // ----------------------------------------------------------------------------
256 {
257 
258  StarGenStats stats( GetName(), "Pythia 8 Run Statistics" );
259  const Pythia8::Info &info = mPythia->info;
260 
261  // Print pythia statistics
262  mPythia->stat();
263 
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();
270 
271  stats.nFilterSeen = stats.nAccepted;
272  stats.nFilterAccept = stats.nAccepted;
273 
274  stats.Dump();
275 
276  // Return a copy of the class we just created
277  return stats;
278 
279 }
280 
281 
282 // ----------------------------------------------------------------------------
283 int StarPythia8::InitCMS( int blue, int yell )
284 {
285  Set( Form("Beams:idA=%i",blue) ); // +z --> towards EEMC
286  Set( Form("Beams:idB=%i",yell) ); // -z --> towards ETOF
287  Set( "Beams:frameType=1"); // 1=CMS
288  Set( Form("Beams:eCM=%f",mRootS) );
289 //mPythia->init( blue, yell, mRootS );
290  mPythia->init();
291  return 0;
292 }
293 // ----------------------------------------------------------------------------
294 int StarPythia8::InitFIXT( int blue, int yell )
295 {
296  Set( Form("Beams:idA=%i",blue) ); // +z --> towards EEMC
297  Set( Form("Beams:idB=%i",yell) ); // -z --> towards ETOF
298  Set( "Beams:frameType=2"); // 2=fixed target mode
299  double eA = (mRootS>0)? mRootS : 0.0;
300  double eB = (mRootS<0)? mRootS : 0.0;
301  Set( Form("Beams:eA=%f", eA) );
302  Set( Form("Beams:eB=%f", eB) );
303  mPythia->init();
304  return 0;
305 }
306 // ----------------------------------------------------------------------------
307 int StarPythia8::Init3MOM( int blue, int yell )
308 {
309  LOG_INFO << "Init3/4/5MOM not implemented yet" << endm;
310  assert(0);
311  // cout << "Initializing 3MOM: " << endl;
312  // mBlueMomentum.Print();
313  // mYellMomentum.Print();
314 
315 
316  // mPythia -> init( blue, yell,
317  // mBlueMomentum.X(), mBlueMomentum.Y(), mBlueMomentum.Z(),
318  // mYellMomentum.X(), mYellMomentum.Y(), mYellMomentum.Z() );
319 
320  return 0;
321 }
322 // ----------------------------------------------------------------------------
323 int StarPythia8::Init4MOM( int blue, int yell )
324 {
325  return Init3MOM( blue, yell );
326 }
327 // ----------------------------------------------------------------------------
328 int StarPythia8::Init5MOM( int blue, int yell )
329 {
330  return Init3MOM( blue, yell );
331 }
332 // ----------------------------------------------------------------------------
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.
Definition: StarRandom.cxx:87
TString mYell
Name of the yellow beam particle (-z)
Event record class tailored to PP kinematics.
ABC for defining event generator interfaces.
Definition: StarGenerator.h:34
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.
Definition: StarRandom.cxx:68
End of run statistics for event generators.
Definition: StarGenStats.h:21
Base class for event records.
Definition: StarGenEvent.h:81
virtual const char * GetName() const
special overload
Definition: StMaker.cxx:237
Event record class tailored to DIS kinemaics.
Definition: Stypes.h:40
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.
Definition: StarPythia8.cxx:51
StarGenStats Stats()
Return end-of-run statistics.
void Set(const char *s)
Pass a string to Pythia8::Pythia::readString(), for user configuration.
Definition: StarPythia8.h:91