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