ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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,
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 /// \author Andreas Morsch
17 
18 void pythia8(Int_t nev = 100, Int_t ndeb = 1)
19 {
20  const char *p8dataenv = gSystem->Getenv("PYTHIA8DATA");
21  if (!p8dataenv) {
22  const char *p8env = gSystem->Getenv("PYTHIA8");
23  if (!p8env) {
24  Error("pythia8.C",
25  "Environment variable PYTHIA8 must contain path to pythia directory!");
26  return;
27  }
28  TString p8d = p8env;
29  p8d += "/xmldoc";
30  gSystem->Setenv("PYTHIA8DATA", p8d);
31  }
32 
33  const char* path = gSystem->ExpandPathName("$PYTHIA8DATA");
34  if (gSystem->AccessPathName(path)) {
35  Error("pythia8.C",
36  "Environment variable PYTHIA8DATA must contain path to $PYTHIA8/xmldoc directory !");
37  return;
38  }
39 
40 // Load libraries
41 #ifndef G__WIN32 // Pythia8 is a static library on Windows
42  gSystem->Load("$PYTHIA8/lib/libpythia8");
43 #endif
44  gSystem->Load("libEG");
45  gSystem->Load("libEGPythia8");
46 // Histograms
47  TH1F* etaH = new TH1F("etaH", "Pseudorapidity", 120, -12., 12.);
48  TH1F* ptH = new TH1F("ptH", "pt", 50, 0., 10.);
49 
50 
51 // Array of particles
52  TClonesArray* particles = new TClonesArray("TParticle", 1000);
53 // Create pythia8 object
54  TPythia8* pythia8 = new TPythia8();
55 
56 // Configure
57  pythia8->ReadString("HardQCD:all = on");
58 
59 
60 // Initialize
61 
62  pythia8->Initialize(2212 /* p */, 2212 /* p */, 14000. /* TeV */);
63 
64 // Event loop
65  for (Int_t iev = 0; iev < nev; iev++) {
66  pythia8->GenerateEvent();
67  if (iev < ndeb) pythia8->EventListing();
68  pythia8->ImportParticles(particles,"All");
69  Int_t np = particles->GetEntriesFast();
70 // Particle loop
71  for (Int_t ip = 0; ip < np; ip++) {
72  TParticle* part = (TParticle*) particles->At(ip);
73  Int_t ist = part->GetStatusCode();
74  // Positive codes are final particles.
75  if (ist <= 0) continue;
76  Int_t pdg = part->GetPdgCode();
78  if (charge == 0.) continue;
79  Float_t eta = part->Eta();
80  Float_t pt = part->Pt();
81 
82  etaH->Fill(eta);
83  if (pt > 0.) ptH->Fill(pt, 1./(2. * pt));
84  }
85  }
86 
87  pythia8->PrintStatistics();
88 
89  TCanvas* c1 = new TCanvas("c1","Pythia8 test example",800,800);
90  c1->Divide(1, 2);
91  c1->cd(1);
92  etaH->Scale(5./Float_t(nev));
93  etaH->Draw();
94  etaH->SetXTitle("#eta");
95  etaH->SetYTitle("dN/d#eta");
96 
97  c1->cd(2);
98  gPad->SetLogy();
99  ptH->Scale(5./Float_t(nev));
100  ptH->Draw();
101  ptH->SetXTitle("p_{t} [GeV/c]");
102  ptH->SetYTitle("dN/dp_{t}^{2} [GeV/c]^{-2}");
103  }
void PrintStatistics() const
Print end of run statistics.
Definition: TPythia8.cxx:355
Double_t Charge() const
Definition: TParticlePDG.h:72
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1213
Double_t Pt() const
Definition: TParticle.h:122
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition: TH1.cxx:6174
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
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
Double_t Eta() const
Definition: TParticle.h:124
Description of the dynamic properties of a particle.
Definition: TParticle.h:34
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
TParticlePDG * GetParticle(Int_t pdgCode) const
Get a pointer to the particle object according to the MC code number.
Bool_t Initialize(Int_t idAin, Int_t idBin, Double_t ecms)
Initialization.
Definition: TPythia8.cxx:146
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
Load a shared library.
Definition: TSystem.cxx:1766
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
int Int_t
Definition: RtypesCore.h:41
virtual void SetYTitle(const char *title)
Definition: TH1.h:409
Int_t GetEntriesFast() const
Definition: TObjArray.h:66
static TDatabasePDG * Instance()
static function
virtual const char * Getenv(const char *env)
Get environment variable.
Definition: TSystem.cxx:1575
tuple np
Definition: multifit.py:30
void Error(const char *location, const char *msgfmt,...)
virtual void Setenv(const char *name, const char *value)
Set environment variable.
Definition: TSystem.cxx:1559
R__EXTERN TSystem * gSystem
Definition: TSystem.h:545
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
Int_t GetStatusCode() const
Definition: TParticle.h:89
Int_t GetPdgCode() const
Definition: TParticle.h:90
The Canvas class.
Definition: TCanvas.h:48
void EventListing() const
Event listing.
Definition: TPythia8.cxx:363
TPythia8 is an interface class to C++ version of Pythia 8.1 event generators, written by T...
Definition: TPythia8.h:76
An array of clone (identical) objects.
Definition: TClonesArray.h:32
virtual void SetXTitle(const char *title)
Definition: TH1.h:408
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:1073
#define gPad
Definition: TVirtualPad.h:288
virtual void GenerateEvent()
Generate the next event.
Definition: TPythia8.cxx:184
virtual Bool_t ExpandPathName(TString &path)
Expand a pathname getting rid of special shell characters like ~.
Definition: TSystem.cxx:1191
TObject * At(Int_t idx) const
Definition: TObjArray.h:167
void ReadString(const char *string) const
Configuration.
Definition: TPythia8.cxx:299