Logo ROOT   6.12/06
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 /// \macro_code
12 ///
13 /// \author Andreas Morsch
14 
15 #include "TSystem.h"
16 #include "TH1F.h"
17 #include "TClonesArray.h"
18 #include "TPythia8.h"
19 #include "TParticle.h"
20 #include "TDatabasePDG.h"
21 #include "TCanvas.h"
22 
23 void pythia8(Int_t nev = 100, Int_t ndeb = 1)
24 {
25 // Load libraries
26  gSystem->Load("libEG");
27  gSystem->Load("libEGPythia8");
28 // Histograms
29  TH1F* etaH = new TH1F("etaH", "Pseudorapidity", 120, -12., 12.);
30  TH1F* ptH = new TH1F("ptH", "pt", 50, 0., 10.);
31 
32 
33 // Array of particles
34  TClonesArray* particles = new TClonesArray("TParticle", 1000);
35 // Create pythia8 object
36  TPythia8* pythia8 = new TPythia8();
37 
38 // Configure
39  pythia8->ReadString("HardQCD:all = on");
40 
41 
42 // Initialize
43 
44  pythia8->Initialize(2212 /* p */, 2212 /* p */, 14000. /* TeV */);
45 
46 // Event loop
47  for (Int_t iev = 0; iev < nev; iev++) {
48  pythia8->GenerateEvent();
49  if (iev < ndeb) pythia8->EventListing();
50  pythia8->ImportParticles(particles,"All");
51  Int_t np = particles->GetEntriesFast();
52 // Particle loop
53  for (Int_t ip = 0; ip < np; ip++) {
54  TParticle* part = (TParticle*) particles->At(ip);
55  Int_t ist = part->GetStatusCode();
56  // Positive codes are final particles.
57  if (ist <= 0) continue;
58  Int_t pdg = part->GetPdgCode();
60  if (charge == 0.) continue;
61  Float_t eta = part->Eta();
62  Float_t pt = part->Pt();
63 
64  etaH->Fill(eta);
65  if (pt > 0.) ptH->Fill(pt, 1./(2. * pt));
66  }
67  }
68 
69  pythia8->PrintStatistics();
70 
71  TCanvas* c1 = new TCanvas("c1","Pythia8 test example",800,800);
72  c1->Divide(1, 2);
73  c1->cd(1);
74  etaH->Scale(5./Float_t(nev));
75  etaH->Draw();
76  etaH->SetXTitle("#eta");
77  etaH->SetYTitle("dN/d#eta");
78 
79  c1->cd(2);
80  gPad->SetLogy();
81  ptH->Scale(5./Float_t(nev));
82  ptH->Draw();
83  ptH->SetXTitle("p_{t} [GeV/c]");
84  ptH->SetYTitle("dN/dp_{t}^{2} [GeV/c]^{-2}");
85  }
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:6061
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3242
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
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:556
int Int_t
Definition: RtypesCore.h:41
virtual void SetYTitle(const char *title)
Definition: TH1.h:395
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:2967
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:394
#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