eic-smear  1.0.3
A collection of ROOT classes for Monte Carlo events and a fast-smearing code simulating detector effects for the Electron-Ion Collider task force
read.cxx
Go to the documentation of this file.
1 // read.cxx
2 //
3 // Created by TB on 6/13/11.
4 // Copyright 2011 BNL. All rights reserved.
5 //
6 // Example of how to read a file produced by BuildTree for a simple analysis.
7 // To run, in ROOT do:
8 // root [0] .L /path/to/read.cxx
9 // root [1] read("myInputFile.root", 10000 )
10 
11 void read(TString inFileNames, int nEvents ) {
12 
13  // If the analysis solely uses TTree::Draw statements,
14  // you don't need to load
15  // the shared library. You will receive warnings such as
16  // Warning in <TClass::TClass>: no dictionary for class Particle
17  // is available
18  // but these can be ignored. However if you wish to work with the event
19  // objects themselves, the shared library must be loaded:
20  // Load the shared library, if not done automaticlly:
21  // gSystem->Load("/path/to/libeicsmear.so" );
22 
23  // The TTrees are named EICTree.
24  // Create a TChain for trees with this name.
25  TChain tree("EICTree");
26 
27  // Add the file(s) we want to analyse to the chain.
28  // We could add multiple files if we wanted.
29  tree.Add(inFileNames); // Wild cards are allowed e.g. tree.Add("*.root" );
30 // tree.Add(/path/to/otherFileNames ); // etc...
31 
32  // Create an object to store the current event from the tree.
33  // This is how we access the values in the tree.
34  // If you want to use generator-specific values, then
35  // the event type should match the type in the TTree. Valid types are
36  // EventPythia, EventPepsi, EventRapgap, EventDjangoh, EventMilou.
37  // If you only need shared quantities like x, Q2 and the particle list
38  // you can use EventBase and the macro will be general for any Monte Carlo.
39  erhic::EventPythia* event(NULL);// = new EventPythia;
40 // EventBase* event(NULL);
41 
42  // Now associate the contents of the branch with the buffer.
43  // The events are stored in a branch named event:
44  tree.SetBranchAddress("event", &event ); // Note &event, not event.
45 
46  // Now we can do some analysis...
47 
48  // We record the largest particle pT we find here:
49  double highestPt(-1. );
50 
51  // Histograms for our analysis.
52  TH2D Q2VsX("Q2VsX",
53  "Q^{2} vs. Bjorken x;log_{10}(x);log_{10}(Q^{2})",
54  100, -5., 0., 100, -2., 3. );
55  TH1D ptHist("ptHist",
56  "pT of charged pions",
57  50, 0.0, 2.0 );
58  TH1D deltaPhi("deltaPhi",
59  "Delta-phi of hadrons",
60  40, 0.0, 3.1415 );
61 
62  // Loop over events:
63  for(int i(0); i < nEvents; ++i ) {
64 
65  // Read the next entry from the tree.
66  tree.GetEntry(i);
67 
68  // Fill the Q2 vs. x histogram:
69  Q2VsX.Fill(TMath::Log10(event->GetX()),
70  TMath::Log10(event->GetQ2()));
71 
72  // The event contains a vector (array) of particles.
73  int nParticles = event->GetNTracks();
74 
75  // We now know the number of particles in the event, so loop over
76  // the particles:
77  for(int j(0); j < nParticles; ++j ) {
78  const Particle* particle = event->GetTrack(j);
79  // Let's just select charged pions for this example:
80  int pdg = particle->GetPdgCode();
81  if(abs(pdg) != 211 ) continue;
82 
83  ptHist.Fill(particle->GetPt());
84 
85  // Update the highest pT:
86  if(particle->GetPt() > highestPt ) {
87  highestPt = particle->GetPt();
88  } // if
89  } // for
90  } // for
91 
92  std::cout << "The highest pT was " << highestPt << " GeV/c" << std::endl;
93 
94  TCanvas canvas;
95  Q2VsX.Draw("colz" );
96  canvas.Print("Q2VsX.png" );
97  ptHist.Draw();
98  canvas.Print("pt.png" );
99 }
erhic::ParticleMCbase::GetPdgCode
virtual Pid GetPdgCode() const
Returns the ID of the particle.
Definition: erhic/ParticleMC.h:252
erhic::ParticleMCbase::GetPt
virtual Double_t GetPt() const
Returns momentum transverse to the beam direction.
Definition: erhic/ParticleMC.h:496
tree
Definition: tree.py:1
erhic::EventPythia
Pythia 6 DIS event.
Definition: erhic/EventPythia.h:28
read
void read(TString inFileNames, int nEvents)
Definition: read.cxx:11
erhic::EventDis::GetX
virtual Double_t GetX() const
Returns Bjorken-x of the event.
Definition: EventDis.h:198
erhic::ParticleMC
Definition: erhic/ParticleMC.h:403
erhic::EventDis::GetQ2
virtual Double_t GetQ2() const
Returns the four-momentum transfer (exchange boson mass) Q2.
Definition: EventDis.h:206