Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
91TPythia8::TPythia8(bool printBanner /*= true*/):
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 // If only we could skip stating the default of the first parameter...:
103 fPythia = new Pythia8::Pythia("../share/Pythia8/xmldoc", printBanner);
104}
105
106////////////////////////////////////////////////////////////////////////////////
107/// Constructor with an xmlDir (eg "../xmldoc"
108
109TPythia8::TPythia8(const char *xmlDir, bool printBanner /*= true*/):
110 TGenerator("TPythia8", "TPythia8"),
111 fPythia(0),
112 fNumberOfParticles(0)
113{
114 if (fgInstance)
115 Fatal("TPythia8", "There's already an instance of TPythia8");
116
117 delete fParticles; // was allocated as TObjArray in TGenerator
118
119 fParticles = new TClonesArray("TParticle",50);
120 fPythia = new Pythia8::Pythia(xmlDir, printBanner);
121}
122
123////////////////////////////////////////////////////////////////////////////////
124/// Destructor
125
127{
128 if (fParticles) {
130 delete fParticles;
131 fParticles = 0;
132 }
133 delete fPythia;
134}
135
136////////////////////////////////////////////////////////////////////////////////
137/// Return an instance of TPythia8
138
140{
141 return fgInstance ? fgInstance : (fgInstance = new TPythia8()) ;
142}
143
144////////////////////////////////////////////////////////////////////////////////
145/// Initialization
146
148{
150
151 // Set arguments in Settings database.
152 fPythia->settings.mode("Beams:idA", idAin);
153 fPythia->settings.mode("Beams:idB", idBin);
154 fPythia->settings.mode("Beams:frameType", 1);
155 fPythia->settings.parm("Beams:eCM", ecms);
156
157 return fPythia->init();
158
159 //return fPythia->init(idAin, idBin, ecms);
160}
161
162////////////////////////////////////////////////////////////////////////////////
163/// Initialization
164
166{
168
169 // Set arguments in Settings database.
170 fPythia->settings.mode("Beams:idA", idAin);
171 fPythia->settings.mode("Beams:idB", idBin);
172 fPythia->settings.mode("Beams:frameType", 2);
173 fPythia->settings.parm("Beams:eA", eAin);
174 fPythia->settings.parm("Beams:eB", eBin);
175
176 // Send on to common initialization.
177 return fPythia->init();
178
179 //return fPythia->init(idAin, idBin, eAin, eBin);
180}
181
182////////////////////////////////////////////////////////////////////////////////
183/// Generate the next event
184
186{
187 fPythia->next();
188 fNumberOfParticles = fPythia->event.size() - 1;
190}
191////////////////////////////////////////////////////////////////////////////////
192/// Import particles from Pythia stack
193
195{
196 if (particles == 0) return 0;
197 TClonesArray &clonesParticles = *particles;
198 clonesParticles.Clear();
199 Int_t nparts=0;
200 Int_t i;
201 Int_t ioff = 0;
202 fNumberOfParticles = fPythia->event.size();
203 if (fPythia->event[0].id() == 90) {
204 ioff = -1;
205 }
206
207 if (!strcmp(option,"") || !strcmp(option,"Final")) {
208 for (i = 0; i < fNumberOfParticles; i++) {
209 if (fPythia->event[i].id() == 90) continue;
210 if (fPythia->event[i].isFinal()) {
211 new(clonesParticles[nparts]) TParticle(
212 fPythia->event[i].id(),
213 fPythia->event[i].isFinal(),
214 fPythia->event[i].mother1() + ioff,
215 fPythia->event[i].mother2() + ioff,
216 fPythia->event[i].daughter1() + ioff,
217 fPythia->event[i].daughter2() + ioff,
218 fPythia->event[i].px(), // [GeV/c]
219 fPythia->event[i].py(), // [GeV/c]
220 fPythia->event[i].pz(), // [GeV/c]
221 fPythia->event[i].e(), // [GeV]
222 fPythia->event[i].xProd(), // [mm]
223 fPythia->event[i].yProd(), // [mm]
224 fPythia->event[i].zProd(), // [mm]
225 fPythia->event[i].tProd()); // [mm/c]
226 nparts++;
227 } // final state partice
228 } // particle loop
229 } else if (!strcmp(option,"All")) {
230 for (i = 0; i < fNumberOfParticles; i++) {
231 if (fPythia->event[i].id() == 90) continue;
232 new(clonesParticles[nparts]) TParticle(
233 fPythia->event[i].id(),
234 fPythia->event[i].isFinal(),
235 fPythia->event[i].mother1() + ioff,
236 fPythia->event[i].mother2() + ioff,
237 fPythia->event[i].daughter1() + ioff,
238 fPythia->event[i].daughter2() + ioff,
239 fPythia->event[i].px(), // [GeV/c]
240 fPythia->event[i].py(), // [GeV/c]
241 fPythia->event[i].pz(), // [GeV/c]
242 fPythia->event[i].e(), // [GeV]
243 fPythia->event[i].xProd(), // [mm]
244 fPythia->event[i].yProd(), // [mm]
245 fPythia->event[i].zProd(), // [mm]
246 fPythia->event[i].tProd()); // [mm/c]
247 nparts++;
248 } // particle loop
249 }
250 if(ioff==-1) fNumberOfParticles--;
251 return nparts;
252}
253
254////////////////////////////////////////////////////////////////////////////////
255/// Import particles from Pythia stack
256
258{
259 fParticles->Clear();
260 Int_t ioff = 0;
261 Int_t numpart = fPythia->event.size();
262 if (fPythia->event[0].id() == 90) {
263 numpart--;
264 ioff = -1;
265 }
266
267
269 for (Int_t i = 1; i <= numpart; i++) {
270 new(a[i]) TParticle(
271 fPythia->event[i].id(),
272 fPythia->event[i].isFinal(),
273 fPythia->event[i].mother1() + ioff,
274 fPythia->event[i].mother2() + ioff,
275 fPythia->event[i].daughter1() + ioff,
276 fPythia->event[i].daughter2() + ioff,
277 fPythia->event[i].px(), // [GeV/c]
278 fPythia->event[i].py(), // [GeV/c]
279 fPythia->event[i].pz(), // [GeV/c]
280 fPythia->event[i].e(), // [GeV]
281 fPythia->event[i].xProd(), // [mm]
282 fPythia->event[i].yProd(), // [mm]
283 fPythia->event[i].zProd(), // [mm]
284 fPythia->event[i].tProd()); // [mm/c]
285 }
286 return fParticles;
287}
288
289////////////////////////////////////////////////////////////////////////////////
290/// Initialization
291
293{
294 return (fPythia->event.size() - 1);
295}
296
297////////////////////////////////////////////////////////////////////////////////
298/// Configuration
299
300void TPythia8::ReadString(const char* string) const
301{
302 fPythia->readString(string);
303}
304
305////////////////////////////////////////////////////////////////////////////////
306/// Configuration
307
308void TPythia8::ReadConfigFile(const char* string) const
309{
310 fPythia->readFile(string);
311}
312
313////////////////////////////////////////////////////////////////////////////////
314/// Event listing
315
317{
318 fPythia->settings.listAll();
319}
320
321////////////////////////////////////////////////////////////////////////////////
322/// Event listing
323
325{
326 fPythia->settings.listChanged();
327}
328
329////////////////////////////////////////////////////////////////////////////////
330/// Event listing
331
332void TPythia8::Plist(Int_t id) const
333{
334 fPythia->particleData.list(id);
335}
336
337////////////////////////////////////////////////////////////////////////////////
338/// Event listing
339
341{
342 fPythia->particleData.listAll();
343}
344
345////////////////////////////////////////////////////////////////////////////////
346/// Event listing
347
349{
350 fPythia->particleData.listChanged();
351}
352
353////////////////////////////////////////////////////////////////////////////////
354/// Print end of run statistics
355
357{
358 fPythia->stat();
359}
360
361////////////////////////////////////////////////////////////////////////////////
362/// Event listing
363
365{
366 fPythia->event.list();
367}
368
369////////////////////////////////////////////////////////////////////////////////
370/// Add some pythia specific particle code to the data base
371
373{
375 pdgDB->AddParticle("string","string", 0, kTRUE,
376 0, 0, "QCD string", 90);
377 pdgDB->AddParticle("rho_diff0", "rho_diff0", 0, kTRUE,
378 0, 0, "QCD diffr. state", 9900110);
379 pdgDB->AddParticle("pi_diffr+", "pi_diffr+", 0, kTRUE,
380 0, 1, "QCD diffr. state", 9900210);
381 pdgDB->AddParticle("omega_di", "omega_di", 0, kTRUE,
382 0, 0, "QCD diffr. state", 9900220);
383 pdgDB->AddParticle("phi_diff","phi_diff", 0, kTRUE,
384 0, 0, "QCD diffr. state", 9900330);
385 pdgDB->AddParticle("J/psi_di", "J/psi_di", 0, kTRUE,
386 0, 0, "QCD diffr. state", 9900440);
387 pdgDB->AddParticle("n_diffr0","n_diffr0",0,kTRUE,
388 0, 0, "QCD diffr. state", 9902110);
389 pdgDB->AddParticle("p_diffr+","p_diffr+", 0, kTRUE,
390 0, 1, "QCD diffr. state", 9902210);
391}
392
#define a(i)
Definition RSha256.hxx:99
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t option
An array of clone (identical) objects.
void Clear(Option_t *option="") override
Clear the clones array.
Particle database manager class.
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:31
void Clear(Option_t *option="") override
Remove all objects from the array.
void Delete(Option_t *option="") override
Remove all objects from the array AND delete all heap based objects.
friend class TClonesArray
Definition TObject.h:242
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:1004
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:332
void PlistChanged() const
Event listing.
Definition TPythia8.cxx:348
Int_t ImportParticles(TClonesArray *particles, Option_t *option="") override
Import particles from Pythia stack.
Definition TPythia8.cxx:194
void ListAll() const
Event listing.
Definition TPythia8.cxx:316
void ReadString(const char *string) const
Configuration.
Definition TPythia8.cxx:300
Bool_t Initialize(Int_t idAin, Int_t idBin, Double_t ecms)
Initialization.
Definition TPythia8.cxx:147
void EventListing() const
Event listing.
Definition TPythia8.cxx:364
void GenerateEvent() override
Generate the next event.
Definition TPythia8.cxx:185
void ListChanged() const
Event listing.
Definition TPythia8.cxx:324
Int_t fNumberOfParticles
The pythia8 instance.
Definition TPythia8.h:83
void PrintStatistics() const
Print end of run statistics.
Definition TPythia8.cxx:356
Int_t GetN() const
Initialization.
Definition TPythia8.cxx:292
void PlistAll() const
Event listing.
Definition TPythia8.cxx:340
static TPythia8 * Instance()
Return an instance of TPythia8.
Definition TPythia8.cxx:139
TPythia8(bool printBanner=true)
Number of particles.
Definition TPythia8.cxx:91
static TPythia8 * fgInstance
Definition TPythia8.h:81
void ReadConfigFile(const char *string) const
Configuration.
Definition TPythia8.cxx:308
void AddParticlesToPdgDataBase()
Add some pythia specific particle code to the data base.
Definition TPythia8.cxx:372
Pythia8::Pythia * fPythia
singleton instance
Definition TPythia8.h:82
~TPythia8() override
Destructor.
Definition TPythia8.cxx:126