Logo ROOT   6.16/01
Reference Guide
pythia8.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_pythia
3/// pythia8 basic example
4///
5/// to run, do:
6///
7/// ~~~{.cpp}
8/// root > .x pythia8.C
9/// ~~~
10///
11/// Note that before executing this script, for some Pythia8 builds:
12///
13/// - the env variable PYTHIA8 must point to the pythia8100 (or newer) directory
14/// - the env variable PYTHIA8DATA must be defined and it must point to $PYTHIA8/xmldoc
15///
16/// \macro_code
17///
18/// \author Andreas Morsch
19
20#include "TSystem.h"
21#include "TH1F.h"
22#include "TClonesArray.h"
23#include "TPythia8.h"
24#include "TParticle.h"
25#include "TDatabasePDG.h"
26#include "TCanvas.h"
27
28void pythia8(Int_t nev = 100, Int_t ndeb = 1)
29{
30// Load libraries
31 gSystem->Load("libEG");
32 gSystem->Load("libEGPythia8");
33// Histograms
34 TH1F* etaH = new TH1F("etaH", "Pseudorapidity", 120, -12., 12.);
35 TH1F* ptH = new TH1F("ptH", "pt", 50, 0., 10.);
36
37
38// Array of particles
39 TClonesArray* particles = new TClonesArray("TParticle", 1000);
40// Create pythia8 object
41 TPythia8* pythia8 = new TPythia8();
42
43#if PYTHIA_VERSION_INTEGER == 8235
44 // Pythia 8.235 is known to cause crashes:
45 printf("ABORTING PYTHIA8 TUTORIAL!\n");
46 printf("The version of Pythia you use is known to case crashes due to memory errors.\n");
47 printf("They have been reported to the authors; the Pythia versions 8.1... are known to work.\n");
48 return;
49#endif
50
51// Configure
52 pythia8->ReadString("HardQCD:all = on");
53 pythia8->ReadString("Random:setSeed = on");
54 // use a reproducible seed: always the same results for the tutorial.
55 pythia8->ReadString("Random:seed = 42");
56
57
58// Initialize
59
60 pythia8->Initialize(2212 /* p */, 2212 /* p */, 14000. /* TeV */);
61
62// Event loop
63 for (Int_t iev = 0; iev < nev; iev++) {
64 pythia8->GenerateEvent();
65 if (iev < ndeb) pythia8->EventListing();
66 pythia8->ImportParticles(particles,"All");
67 Int_t np = particles->GetEntriesFast();
68// Particle loop
69 for (Int_t ip = 0; ip < np; ip++) {
70 TParticle* part = (TParticle*) particles->At(ip);
71 Int_t ist = part->GetStatusCode();
72 // Positive codes are final particles.
73 if (ist <= 0) continue;
74 Int_t pdg = part->GetPdgCode();
76 if (charge == 0.) continue;
77 Float_t eta = part->Eta();
78 Float_t pt = part->Pt();
79
80 etaH->Fill(eta);
81 if (pt > 0.) ptH->Fill(pt, 1./(2. * pt));
82 }
83 }
84
85 pythia8->PrintStatistics();
86
87 TCanvas* c1 = new TCanvas("c1","Pythia8 test example",800,800);
88 c1->Divide(1, 2);
89 c1->cd(1);
90 etaH->Scale(5./Float_t(nev));
91 etaH->Draw();
92 etaH->SetXTitle("#eta");
93 etaH->SetYTitle("dN/d#eta");
94
95 c1->cd(2);
96 gPad->SetLogy();
97 ptH->Scale(5./Float_t(nev));
98 ptH->Draw();
99 ptH->SetXTitle("p_{t} [GeV/c]");
100 ptH->SetYTitle("dN/dp_{t}^{2} [GeV/c]^{-2}");
101 }
int Int_t
Definition: RtypesCore.h:41
float Float_t
Definition: RtypesCore.h:53
R__EXTERN TSystem * gSystem
Definition: TSystem.h:540
#define gPad
Definition: TVirtualPad.h:286
The Canvas class.
Definition: TCanvas.h:31
An array of clone (identical) objects.
Definition: TClonesArray.h:32
static TDatabasePDG * Instance()
static function
TParticlePDG * GetParticle(Int_t pdgCode) const
Get a pointer to the particle object according to the MC code number.
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
virtual void SetXTitle(const char *title)
Definition: TH1.h:409
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3251
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2974
virtual void SetYTitle(const char *title)
Definition: TH1.h:410
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:6126
Int_t GetEntriesFast() const
Definition: TObjArray.h:64
TObject * At(Int_t idx) const
Definition: TObjArray.h:165
Double_t Charge() const
Definition: TParticlePDG.h:68
Description of the dynamic properties of a particle.
Definition: TParticle.h:26
Int_t GetStatusCode() const
Definition: TParticle.h:82
Int_t GetPdgCode() const
Definition: TParticle.h:83
Double_t Pt() const
Definition: TParticle.h:135
Double_t Eta() const
Definition: TParticle.h:137
TPythia8 is an interface class to C++ version of Pythia 8.1 event generators, written by T....
Definition: TPythia8.h:77
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 PrintStatistics() const
Print end of run statistics.
Definition: TPythia8.cxx:355
virtual void GenerateEvent()
Generate the next event.
Definition: TPythia8.cxx:184
virtual Int_t ImportParticles(TClonesArray *particles, Option_t *option="")
Import particles from Pythia stack.
Definition: TPythia8.cxx:193
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
Load a shared library.
Definition: TSystem.cxx:1843
TPaveText * pt
return c1
Definition: legend1.C:41