Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files 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 /// \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 #if PYTHIA_VERSION_INTEGER == 8235
39  // Pythia 8.235 is known to cause crashes:
40  printf("ABORTING PYTHIA8 TUTORIAL!\n");
41  printf("The version of Pythia you use is known to case crashes due to memory errors.\n");
42  printf("They have been reported to the authors; the Pythia versions 8.1... are known to work.\n");
43  return;
44 #endif
45 
46 // Configure
47  pythia8->ReadString("HardQCD:all = on");
48  pythia8->ReadString("Random:setSeed = on");
49  // use a reproducible seed: always the same results for the tutorial.
50  pythia8->ReadString("Random:seed = 42");
51 
52 
53 // Initialize
54 
55  pythia8->Initialize(2212 /* p */, 2212 /* p */, 14000. /* TeV */);
56 
57 // Event loop
58  for (Int_t iev = 0; iev < nev; iev++) {
59  pythia8->GenerateEvent();
60  if (iev < ndeb) pythia8->EventListing();
61  pythia8->ImportParticles(particles,"All");
62  Int_t np = particles->GetEntriesFast();
63 // Particle loop
64  for (Int_t ip = 0; ip < np; ip++) {
65  TParticle* part = (TParticle*) particles->At(ip);
66  Int_t ist = part->GetStatusCode();
67  // Positive codes are final particles.
68  if (ist <= 0) continue;
69  Int_t pdg = part->GetPdgCode();
70  Float_t charge = TDatabasePDG::Instance()->GetParticle(pdg)->Charge();
71  if (charge == 0.) continue;
72  Float_t eta = part->Eta();
73  Float_t pt = part->Pt();
74 
75  etaH->Fill(eta);
76  if (pt > 0.) ptH->Fill(pt, 1./(2. * pt));
77  }
78  }
79 
80  pythia8->PrintStatistics();
81 
82  TCanvas* c1 = new TCanvas("c1","Pythia8 test example",800,800);
83  c1->Divide(1, 2);
84  c1->cd(1);
85  etaH->Scale(5./Float_t(nev));
86  etaH->Draw();
87  etaH->SetXTitle("#eta");
88  etaH->SetYTitle("dN/d#eta");
89 
90  c1->cd(2);
91  gPad->SetLogy();
92  ptH->Scale(5./Float_t(nev));
93  ptH->Draw();
94  ptH->SetXTitle("p_{t} [GeV/c]");
95  ptH->SetYTitle("dN/dp_{t}^{2} [GeV/c]^{-2}");
96  }