Loading [MathJax]/extensions/tex2jax.js
ROOT  6.06/09
Reference Guide
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
TPythia8.cxx
Go to the documentation of this file.
1 // @(#)root/pythia8:$Name$:$Id$
2 // Author: Andreas Morsch 27/10/2007
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 
13 ////////////////////////////////////////////////////////////////////////////////
14 // //
15 // TPythia8 //
16 // //
17 // TPythia is an interface class to C++ version of Pythia 8.1 //
18 // event generators, written by T.Sjostrand. //
19 // //
20 // The user is assumed to be familiar with the Pythia package. //
21 // This class includes only a basic interface to Pythia8. Because Pythia8 is //
22 // also written in C++, its functions/classes can be called directly from a //
23 // compiled C++ script. //
24 // To call Pythia functions not available in this interface a dictionary must //
25 // be generated. //
26 // see $ROOTSYS/tutorials/pythia/pythia8.C for an example of use from CINT. //
27 // //
28 ////////////////////////////////////////////////////////////////////////////////
29 /*
30 *------------------------------------------------------------------------------------*
31  | |
32  | *------------------------------------------------------------------------------* |
33  | | | |
34  | | | |
35  | | PPP Y Y TTTTT H H III A Welcome to the Lund Monte Carlo! | |
36  | | P P Y Y T H H I A A This is PYTHIA version 8.100 | |
37  | | PPP Y T HHHHH I AAAAA Last date of change: 20 Oct 2007 | |
38  | | P Y T H H I A A | |
39  | | P Y T H H III A A Now is 27 Oct 2007 at 18:26:53 | |
40  | | | |
41  | | Main author: Torbjorn Sjostrand; CERN/PH, CH-1211 Geneva, Switzerland, | |
42  | | and Department of Theoretical Physics, Lund University, Lund, Sweden; | |
43  | | phone: + 41 - 22 - 767 82 27; e-mail: torbjorn@thep.lu.se | |
44  | | Author: Stephen Mrenna; Computing Division, Simulations Group, | |
45  | | Fermi National Accelerator Laboratory, MS 234, Batavia, IL 60510, USA; | |
46  | | phone: + 1 - 630 - 840 - 2556; e-mail: mrenna@fnal.gov | |
47  | | Author: Peter Skands; CERN/PH, CH-1211 Geneva, Switzerland, | |
48  | | and Theoretical Physics Department, | |
49  | | Fermi National Accelerator Laboratory, MS 106, Batavia, IL 60510, USA; | |
50  | | phone: + 41 - 22 - 767 24 59; e-mail: skands@fnal.gov | |
51  | | | |
52  | | The main program reference is the 'Brief Introduction to PYTHIA 8.1', | |
53  | | T. Sjostrand, S. Mrenna and P. Skands, arXiv:0710.3820 | |
54  | | | |
55  | | The main physics reference is the 'PYTHIA 6.4 Physics and Manual', | |
56  | | T. Sjostrand, S. Mrenna and P. Skands, JHEP05 (2006) 026 [hep-ph/0603175]. | |
57  | | | |
58  | | An archive of program versions and documentation is found on the web: | |
59  | | http://www.thep.lu.se/~torbjorn/Pythia.html | |
60  | | | |
61  | | This program is released under the GNU General Public Licence version 2. | |
62  | | Please respect the MCnet Guidelines for Event Generator Authors and Users. | |
63  | | | |
64  | | Disclaimer: this program comes without any guarantees. | |
65  | | Beware of errors and use common sense when interpreting results. | |
66  | | | |
67  | | Copyright (C) 2007 Torbjorn Sjostrand | |
68  | | | |
69  | | | |
70  | *------------------------------------------------------------------------------* |
71  | |
72  *------------------------------------------------------------------------------------*
73 */
74 
75 #include "TPythia8.h"
76 
77 #include "TClonesArray.h"
78 #include "TParticle.h"
79 #include "TDatabasePDG.h"
80 #include "TLorentzVector.h"
81 
82 ClassImp(TPythia8)
83 
84 TPythia8* TPythia8::fgInstance = 0;
85 
86 ////////////////////////////////////////////////////////////////////////////////
87 /// Constructor
88 
89 TPythia8::TPythia8():
90  TGenerator("TPythia8", "TPythia8"),
91  fPythia(0),
92  fNumberOfParticles(0)
93 {
94  if (fgInstance)
95  Fatal("TPythia8", "There's already an instance of TPythia8");
96 
97  delete fParticles; // was allocated as TObjArray in TGenerator
98 
99  fParticles = new TClonesArray("TParticle",50);
100  fPythia = new Pythia8::Pythia();
101 }
102 
103 ////////////////////////////////////////////////////////////////////////////////
104 /// Constructor with an xmlDir (eg "../xmldoc"
105 
106 TPythia8::TPythia8(const char *xmlDir):
107  TGenerator("TPythia8", "TPythia8"),
108  fPythia(0),
109  fNumberOfParticles(0)
110 {
111  if (fgInstance)
112  Fatal("TPythia8", "There's already an instance of TPythia8");
113 
114  delete fParticles; // was allocated as TObjArray in TGenerator
115 
116  fParticles = new TClonesArray("TParticle",50);
117  fPythia = new Pythia8::Pythia(xmlDir);
118 }
119 
120 ////////////////////////////////////////////////////////////////////////////////
121 /// Destructor
122 
124 {
125  if (fParticles) {
126  fParticles->Delete();
127  delete fParticles;
128  fParticles = 0;
129  }
130  delete fPythia;
131 }
132 
133 ////////////////////////////////////////////////////////////////////////////////
134 /// Return an instance of TPythia8
135 
137 {
138  return fgInstance ? fgInstance : (fgInstance = new TPythia8()) ;
139 }
140 
141 ////////////////////////////////////////////////////////////////////////////////
142 /// Initialization
143 
145 {
147 
148  // Set arguments in Settings database.
149  fPythia->settings.mode("Beams:idA", idAin);
150  fPythia->settings.mode("Beams:idB", idBin);
151  fPythia->settings.mode("Beams:frameType", 1);
152  fPythia->settings.parm("Beams:eCM", ecms);
153 
154  return fPythia->init();
155 
156  //return fPythia->init(idAin, idBin, ecms);
157 }
158 
159 ////////////////////////////////////////////////////////////////////////////////
160 /// Initialization
161 
163 {
165 
166  // Set arguments in Settings database.
167  fPythia->settings.mode("Beams:idA", idAin);
168  fPythia->settings.mode("Beams:idB", idBin);
169  fPythia->settings.mode("Beams:frameType", 2);
170  fPythia->settings.parm("Beams:eA", eAin);
171  fPythia->settings.parm("Beams:eB", eBin);
172 
173  // Send on to common initialization.
174  return fPythia->init();
175 
176  //return fPythia->init(idAin, idBin, eAin, eBin);
177 }
178 
179 ////////////////////////////////////////////////////////////////////////////////
180 /// Generate the next event
181 
183 {
184  fPythia->next();
185  fNumberOfParticles = fPythia->event.size() - 1;
186  ImportParticles();
187 }
188 ////////////////////////////////////////////////////////////////////////////////
189 /// Import particles from Pythia stack
190 
192 {
193  if (particles == 0) return 0;
194  TClonesArray &clonesParticles = *particles;
195  clonesParticles.Clear();
196  Int_t nparts=0;
197  Int_t i;
198  Int_t ioff = 0;
199  fNumberOfParticles = fPythia->event.size();
200  if (fPythia->event[0].id() == 90) {
201  ioff = -1;
202  }
203 
204  if (!strcmp(option,"") || !strcmp(option,"Final")) {
205  for (i = 0; i < fNumberOfParticles; i++) {
206  if (fPythia->event[i].id() == 90) continue;
207  if (fPythia->event[i].isFinal()) {
208  new(clonesParticles[nparts]) TParticle(
209  fPythia->event[i].id(),
210  fPythia->event[i].isFinal(),
211  fPythia->event[i].mother1() + ioff,
212  fPythia->event[i].mother2() + ioff,
213  fPythia->event[i].daughter1() + ioff,
214  fPythia->event[i].daughter2() + ioff,
215  fPythia->event[i].px(), // [GeV/c]
216  fPythia->event[i].py(), // [GeV/c]
217  fPythia->event[i].pz(), // [GeV/c]
218  fPythia->event[i].e(), // [GeV]
219  fPythia->event[i].xProd(), // [mm]
220  fPythia->event[i].yProd(), // [mm]
221  fPythia->event[i].zProd(), // [mm]
222  fPythia->event[i].tProd()); // [mm/c]
223  nparts++;
224  } // final state partice
225  } // particle loop
226  } else if (!strcmp(option,"All")) {
227  for (i = 0; i < fNumberOfParticles; i++) {
228  if (fPythia->event[i].id() == 90) continue;
229  new(clonesParticles[nparts]) TParticle(
230  fPythia->event[i].id(),
231  fPythia->event[i].isFinal(),
232  fPythia->event[i].mother1() + ioff,
233  fPythia->event[i].mother2() + ioff,
234  fPythia->event[i].daughter1() + ioff,
235  fPythia->event[i].daughter2() + ioff,
236  fPythia->event[i].px(), // [GeV/c]
237  fPythia->event[i].py(), // [GeV/c]
238  fPythia->event[i].pz(), // [GeV/c]
239  fPythia->event[i].e(), // [GeV]
240  fPythia->event[i].xProd(), // [mm]
241  fPythia->event[i].yProd(), // [mm]
242  fPythia->event[i].zProd(), // [mm]
243  fPythia->event[i].tProd()); // [mm/c]
244  nparts++;
245  } // particle loop
246  }
247  if(ioff==-1) fNumberOfParticles--;
248  return nparts;
249 }
250 
251 ////////////////////////////////////////////////////////////////////////////////
252 /// Import particles from Pythia stack
253 
255 {
256  fParticles->Clear();
257  Int_t ioff = 0;
258  Int_t numpart = fPythia->event.size();
259  if (fPythia->event[0].id() == 90) {
260  numpart--;
261  ioff = -1;
262  }
263 
264 
266  for (Int_t i = 1; i <= numpart; i++) {
267  new(a[i]) TParticle(
268  fPythia->event[i].id(),
269  fPythia->event[i].isFinal(),
270  fPythia->event[i].mother1() + ioff,
271  fPythia->event[i].mother2() + ioff,
272  fPythia->event[i].daughter1() + ioff,
273  fPythia->event[i].daughter2() + ioff,
274  fPythia->event[i].px(), // [GeV/c]
275  fPythia->event[i].py(), // [GeV/c]
276  fPythia->event[i].pz(), // [GeV/c]
277  fPythia->event[i].e(), // [GeV]
278  fPythia->event[i].xProd(), // [mm]
279  fPythia->event[i].yProd(), // [mm]
280  fPythia->event[i].zProd(), // [mm]
281  fPythia->event[i].tProd()); // [mm/c]
282  }
283  return fParticles;
284 }
285 
286 ////////////////////////////////////////////////////////////////////////////////
287 /// Initialization
288 
290 {
291  return (fPythia->event.size() - 1);
292 }
293 
294 ////////////////////////////////////////////////////////////////////////////////
295 /// Configuration
296 
297 void TPythia8::ReadString(const char* string) const
298 {
299  fPythia->readString(string);
300 }
301 
302 ////////////////////////////////////////////////////////////////////////////////
303 /// Configuration
304 
305 void TPythia8::ReadConfigFile(const char* string) const
306 {
307  fPythia->readFile(string);
308 }
309 
310 ////////////////////////////////////////////////////////////////////////////////
311 /// Event listing
312 
313 void TPythia8::ListAll() const
314 {
315  fPythia->settings.listAll();
316 }
317 
318 ////////////////////////////////////////////////////////////////////////////////
319 /// Event listing
320 
322 {
323  fPythia->settings.listChanged();
324 }
325 
326 ////////////////////////////////////////////////////////////////////////////////
327 /// Event listing
328 
329 void TPythia8::Plist(Int_t id) const
330 {
331  fPythia->particleData.list(id);
332 }
333 
334 ////////////////////////////////////////////////////////////////////////////////
335 /// Event listing
336 
337 void TPythia8::PlistAll() const
338 {
339  fPythia->particleData.listAll();
340 }
341 
342 ////////////////////////////////////////////////////////////////////////////////
343 /// Event listing
344 
346 {
347  fPythia->particleData.listChanged();
348 }
349 
350 ////////////////////////////////////////////////////////////////////////////////
351 /// Print end of run statistics
352 
354 {
355  fPythia->stat();
356 }
357 
358 ////////////////////////////////////////////////////////////////////////////////
359 /// Event listing
360 
362 {
363  fPythia->event.list();
364 }
365 
366 ////////////////////////////////////////////////////////////////////////////////
367 /// Add some pythia specific particle code to the data base
368 
370 {
372  pdgDB->AddParticle("string","string", 0, kTRUE,
373  0, 0, "QCD string", 90);
374  pdgDB->AddParticle("rho_diff0", "rho_diff0", 0, kTRUE,
375  0, 0, "QCD diffr. state", 9900110);
376  pdgDB->AddParticle("pi_diffr+", "pi_diffr+", 0, kTRUE,
377  0, 1, "QCD diffr. state", 9900210);
378  pdgDB->AddParticle("omega_di", "omega_di", 0, kTRUE,
379  0, 0, "QCD diffr. state", 9900220);
380  pdgDB->AddParticle("phi_diff","phi_diff", 0, kTRUE,
381  0, 0, "QCD diffr. state", 9900330);
382  pdgDB->AddParticle("J/psi_di", "J/psi_di", 0, kTRUE,
383  0, 0, "QCD diffr. state", 9900440);
384  pdgDB->AddParticle("n_diffr0","n_diffr0",0,kTRUE,
385  0, 0, "QCD diffr. state", 9902110);
386  pdgDB->AddParticle("p_diffr+","p_diffr+", 0, kTRUE,
387  0, 1, "QCD diffr. state", 9902210);
388 }
389 
void PlistAll() const
Event listing.
Definition: TPythia8.cxx:337
void PrintStatistics() const
Print end of run statistics.
Definition: TPythia8.cxx:353
An array of TObjects.
Definition: TObjArray.h:39
virtual Int_t ImportParticles(TClonesArray *particles, Option_t *option="")
Import particles from Pythia stack.
Definition: TPythia8.cxx:191
Pythia8::Pythia * fPythia
singleton instance
Definition: TPythia8.h:82
TPythia8()
Number of particles.
static TPythia8 * Instance()
Return an instance of TPythia8.
Definition: TPythia8.cxx:136
void PlistChanged() const
Event listing.
Definition: TPythia8.cxx:345
virtual ~TPythia8()
Destructor.
Definition: TPythia8.cxx:123
void Fatal(const char *location, const char *msgfmt,...)
virtual void Clear(Option_t *option="")
Remove all objects from the array.
Definition: TObjArray.cxx:297
const char Option_t
Definition: RtypesCore.h:62
virtual void Delete(Option_t *option="")
Remove all objects from the array AND delete all heap based objects.
Definition: TObjArray.cxx:328
Bool_t Initialize(Int_t idAin, Int_t idBin, Double_t ecms)
Initialization.
Definition: TPythia8.cxx:144
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
void ListChanged() const
Event listing.
Definition: TPythia8.cxx:321
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:946
virtual TParticlePDG * AddParticle(const char *Name, const char *Title, Double_t Mass, Bool_t Stable, Double_t DecayWidth, Double_t Charge, const char *ParticleClass, Int_t PdgCode, Int_t Anti=-1, Int_t TrackingCode=0)
Particle definition normal constructor.
friend class TClonesArray
Definition: TObject.h:214
static TDatabasePDG * Instance()
static function
void Plist(Int_t id) const
Event listing.
Definition: TPythia8.cxx:329
void ListAll() const
Event listing.
Definition: TPythia8.cxx:313
void AddParticlesToPdgDataBase()
Add some pythia specific particle code to the data base.
Definition: TPythia8.cxx:369
virtual void Clear(Option_t *option="")
Clear the clones array.
void EventListing() const
Event listing.
Definition: TPythia8.cxx:361
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
TObjArray * fParticles
display neutrons if true
Definition: TGenerator.h:151
Int_t GetN() const
Initialization.
Definition: TPythia8.cxx:289
Int_t fNumberOfParticles
The pythia8 instance.
Definition: TPythia8.h:83
An array of clone (identical) objects.
Definition: TClonesArray.h:32
virtual void GenerateEvent()
Generate the next event.
Definition: TPythia8.cxx:182
const Bool_t kTRUE
Definition: Rtypes.h:91
void ReadConfigFile(const char *string) const
Configuration.
Definition: TPythia8.cxx:305
void ReadString(const char *string) const
Configuration.
Definition: TPythia8.cxx:297