Logo ROOT   6.14/05
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 
28 void 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 // Configure
44  pythia8->ReadString("HardQCD:all = on");
45 
46 
47 // Initialize
48 
49  pythia8->Initialize(2212 /* p */, 2212 /* p */, 14000. /* TeV */);
50 
51 // Event loop
52  for (Int_t iev = 0; iev < nev; iev++) {
53  pythia8->GenerateEvent();
54  if (iev < ndeb) pythia8->EventListing();
55  pythia8->ImportParticles(particles,"All");
56  Int_t np = particles->GetEntriesFast();
57 // Particle loop
58  for (Int_t ip = 0; ip < np; ip++) {
59  TParticle* part = (TParticle*) particles->At(ip);
60  Int_t ist = part->GetStatusCode();
61  // Positive codes are final particles.
62  if (ist <= 0) continue;
63  Int_t pdg = part->GetPdgCode();
65  if (charge == 0.) continue;
66  Float_t eta = part->Eta();
67  Float_t pt = part->Pt();
68 
69  etaH->Fill(eta);
70  if (pt > 0.) ptH->Fill(pt, 1./(2. * pt));
71  }
72  }
73 
74  pythia8->PrintStatistics();
75 
76  TCanvas* c1 = new TCanvas("c1","Pythia8 test example",800,800);
77  c1->Divide(1, 2);
78  c1->cd(1);
79  etaH->Scale(5./Float_t(nev));
80  etaH->Draw();
81  etaH->SetXTitle("#eta");
82  etaH->SetYTitle("dN/d#eta");
83 
84  c1->cd(2);
85  gPad->SetLogy();
86  ptH->Scale(5./Float_t(nev));
87  ptH->Draw();
88  ptH->SetXTitle("p_{t} [GeV/c]");
89  ptH->SetYTitle("dN/dp_{t}^{2} [GeV/c]^{-2}");
90  }
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:6101
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3251
virtual Int_t ImportParticles(TClonesArray *particles, Option_t *option="")
Import particles from Pythia stack.
Definition: TPythia8.cxx:193
float Float_t
Definition: RtypesCore.h:53
return c1
Definition: legend1.C:41
THist< 1, float, THistStatContent, THistStatUncertainty > TH1F
Definition: THist.hxx:285
Description of the dynamic properties of a particle.
Definition: TParticle.h:26
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:688
Bool_t Initialize(Int_t idAin, Int_t idBin, Double_t ecms)
Initialization.
Definition: TPythia8.cxx:146
void PrintStatistics() const
Print end of run statistics.
Definition: TPythia8.cxx:355
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
Load a shared library.
Definition: TSystem.cxx:1829
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:567
int Int_t
Definition: RtypesCore.h:41
virtual void SetYTitle(const char *title)
Definition: TH1.h:406
TObject * At(Int_t idx) const
Definition: TObjArray.h:165
Int_t GetStatusCode() const
Definition: TParticle.h:82
static TDatabasePDG * Instance()
static function
Double_t Eta() const
Definition: TParticle.h:137
TPaveText * pt
R__EXTERN TSystem * gSystem
Definition: TSystem.h:540
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2974
Int_t GetEntriesFast() const
Definition: TObjArray.h:64
Double_t Charge() const
Definition: TParticlePDG.h:68
void EventListing() const
Event listing.
Definition: TPythia8.cxx:363
void ReadString(const char *string) const
Configuration.
Definition: TPythia8.cxx:299
The Canvas class.
Definition: TCanvas.h:31
TPythia8 is an interface class to C++ version of Pythia 8.1 event generators, written by T...
Definition: TPythia8.h:76
TParticlePDG * GetParticle(Int_t pdgCode) const
Get a pointer to the particle object according to the MC code number.
Double_t Pt() const
Definition: TParticle.h:135
An array of clone (identical) objects.
Definition: TClonesArray.h:32
virtual void SetXTitle(const char *title)
Definition: TH1.h:405
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1162
#define gPad
Definition: TVirtualPad.h:285
virtual void GenerateEvent()
Generate the next event.
Definition: TPythia8.cxx:184
Int_t GetPdgCode() const
Definition: TParticle.h:83