Logo ROOT   6.18/05
Reference Guide
ProofPythia.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_proofpythia
3///
4/// Selector to generate Monte Carlo events with Pythia8
5///
6/// \macro_code
7///
8/// \author Gerardo Ganis (gerardo.ganis@cern.ch)
9
10#define ProofPythia_cxx
11
12#include <TCanvas.h>
13#include <TFrame.h>
14#include <TPaveText.h>
15#include <TFormula.h>
16#include <TF1.h>
17#include <TH1F.h>
18#include <TMath.h>
19#include <TString.h>
20#include <TStyle.h>
21#include <TSystem.h>
22#include <TParameter.h>
23#include "TClonesArray.h"
24#include "TParticle.h"
25#include "TDatabasePDG.h"
26
27#include "ProofPythia.h"
28#include "TPythia8.h"
29
30//_____________________________________________________________________________
31ProofPythia::ProofPythia()
32{
33 // Constructor
34
35 fHist = 0;
36 fPt = 0;
37 fEta = 0;
38 fPythia = 0;
39 fP = 0;
40}
41
42//_____________________________________________________________________________
43ProofPythia::~ProofPythia()
44{
45 // Destructor
46
47 SafeDelete(fPythia);
48 SafeDelete(fP);
49}
50
51//_____________________________________________________________________________
52void ProofPythia::Begin(TTree * /*tree*/)
53{
54 // The Begin() function is called at the start of the query.
55 // When running with PROOF Begin() is only called on the client.
56 // The tree argument is deprecated (on PROOF 0 is passed).
57
58 TString option = GetOption();
59 Info("Begin", "starting a simple exercise with process option: %s", option.Data());
60}
61
62//_____________________________________________________________________________
63void ProofPythia::SlaveBegin(TTree * /*tree*/)
64{
65 // The SlaveBegin() function is called after the Begin() function.
66 // When running with PROOF SlaveBegin() is called on each slave server.
67 // The tree argument is deprecated (on PROOF 0 is passed).
68
69 TString option = GetOption();
70
71 // Histograms
72 fTot = new TH1F("histo1", "total multiplicity", 25, 0.5, 2500.5);
73 fHist = new TH1F("histo2", "charged multiplicity", 20, 0.5, 500.5);
74 fPt = new TH1F("histo3", "particles pT", 100, 0., 10);
75 fEta = new TH1F("histo4", "particles Eta", 100, -10., 10);
76 fTot->SetFillColor(kBlue);
77 fHist->SetFillColor(kRed);
78 fOutput->Add(fTot);
79 fOutput->Add(fHist);
80 fOutput->Add(fPt);
81 fOutput->Add(fEta);
82
83 fPythia = new TPythia8();
84 // Configure
85 fPythia->SetName("pythia8");
86 fPythia->ReadConfigFile("pythia8/main03.cmnd");
87
88 // Initialize
89 fPythia->Initialize( 2212, 2212, 14000.);
90 fP = new TClonesArray("TParticle", 1000);
91
92}
93
94//_____________________________________________________________________________
95Bool_t ProofPythia::Process(Long64_t entry)
96{
97 // Main event loop
98
99 fPythia->GenerateEvent();
100 if (entry < 2)
101 fPythia->EventListing();
102 fPythia->ImportParticles(fP, "All");
103 Int_t nTot = fPythia->GetN();
104 fPythia->ImportParticles(fP, "All");
105 Int_t np = fP->GetEntriesFast();
106 // Particle loop
107 Int_t nCharged = 0;
108 for (Int_t ip = 0; ip < np; ip++) {
109 TParticle* part = (TParticle*) fP->At(ip);
110 Int_t ist = part->GetStatusCode();
111 Int_t pdg = part->GetPdgCode();
112 if (ist != 1) continue;
114 if (charge == 0.) continue;
115 nCharged++;
116 Float_t eta = part->Eta();
117 Float_t pt = part->Pt();
118 if (pt > 0.) fPt->Fill(pt);
119 if ((eta > -10) && (eta < 10)) fEta->Fill(eta);
120 }
121 fHist->Fill(nCharged);
122 fTot->Fill(nTot);
123
124 return kTRUE;
125}
126
127//_____________________________________________________________________________
128void ProofPythia::SlaveTerminate()
129{
130 // The SlaveTerminate() function is called after all entries or objects
131 // have been processed. When running with PROOF SlaveTerminate() is called
132 // on each slave server.
133}
134
135//_____________________________________________________________________________
136void ProofPythia::Terminate()
137{
138 // The Terminate() function is the last function to be called during
139 // a query. It always runs on the client, it can be used to present
140 // the results graphically or save the results to file.
141
142 //
143 // Create canvas
144 //
145 TCanvas *c1 = new TCanvas("c1","Proof ProofPythia canvas",200,10,700,700);
146 c1->Divide(2, 2);
147
148 if ((fTot = dynamic_cast<TH1F *>(fOutput->FindObject("histo1")))) {
149 c1->cd(1);
150 fTot->Draw("h");
151 }
152
153 if ((fHist = dynamic_cast<TH1F *>(fOutput->FindObject("histo2")))) {
154 c1->cd(2);
155 fHist->Draw("h");
156 }
157
158 if ((fPt = dynamic_cast<TH1F *>(fOutput->FindObject("histo3")))) {
159 c1->cd(3);
160 fPt->Draw("h");
161 }
162
163 if ((fEta = dynamic_cast<TH1F *>(fOutput->FindObject("histo4")))) {
164 c1->cd(4);
165 fEta->Draw("h");
166 }
167
168 // Final update
169 c1->cd();
170 c1->Update();
171}
Selector to generate Monte Carlo events with Pythia8.
#define SafeDelete(p)
Definition: RConfig.hxx:543
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
long long Long64_t
Definition: RtypesCore.h:69
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
@ kRed
Definition: Rtypes.h:64
@ kBlue
Definition: Rtypes.h:64
void Info(const char *location, const char *msgfmt,...)
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
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
Basic string class.
Definition: TString.h:131
const char * Data() const
Definition: TString.h:364
A TTree represents a columnar dataset.
Definition: TTree.h:71
TPaveText * pt
return c1
Definition: legend1.C:41
void Begin(Int_t type)