Logo ROOT  
Reference Guide
jets.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_tree
3///
4/// Usage of a Tree using the JetEvent class.
5///
6/// The JetEvent class has several collections (TClonesArray)
7/// and other collections (TRefArray) referencing objects
8/// in the TClonesArrays.
9/// The JetEvent class is in $ROOTSYS/tutorials/tree/JetEvent.h,cxx
10/// to execute the script, do
11/// ~~~
12/// .x jets.C
13/// ~~~
14///
15/// \macro_code
16///
17/// \author Rene Brun
18
19#ifdef JETS_SECOND_RUN
20
21#include "TFile.h"
22#include "TTree.h"
23#include "TRandom.h"
24#include "TROOT.h"
25#include "TSystem.h"
26#include "JetEvent.h"
27#include "Riostream.h"
28
29
30void write(Int_t nev=100) {
31 //write nev Jet events
32 TFile f("JetEvent.root","recreate");
33 TTree *T = new TTree("T","Event example with Jets");
34 JetEvent *event = new JetEvent;
35 T->Branch("event","JetEvent",&event,8000,2);
36
37 for (Int_t ev=0;ev<nev;ev++) {
38 event->Build();
39 T->Fill();
40 }
41
42 T->Print();
43 T->Write();
44}
45
46void read() {
47 //read the JetEvent file
48 TFile f("JetEvent.root");
49 TTree *T = (TTree*)f.Get("T");
50 JetEvent *event = 0;
51 T->SetBranchAddress("event", &event);
52 Long64_t nentries = T->GetEntries();
53
54 for (Long64_t ev=0;ev<nentries;ev++) {
55 T->GetEntry(ev);
56 if (ev) continue; //dump first event only
57 cout << " Event: "<< ev
58 << " Jets: " << event->GetNjet()
59 << " Tracks: " << event->GetNtrack()
60 << " Hits A: " << event->GetNhitA()
61 << " Hits B: " << event->GetNhitB() << endl;
62 }
63}
64
65void pileup(Int_t nev=200) {
66 //make nev pileup events, each build with LOOPMAX events selected
67 //randomly among the nentries
68 TFile f("JetEvent.root");
69 TTree *T = (TTree*)f.Get("T");
70 // Long64_t nentries = T->GetEntries();
71
72 const Int_t LOOPMAX=10;
73 JetEvent *events[LOOPMAX];
74 Int_t loop;
75 for (loop=0;loop<LOOPMAX;loop++) events[loop] = 0;
76 for (Long64_t ev=0;ev<nev;ev++) {
77 if (ev%10 == 0) printf("building pileup: %lld\n",ev);
78 for (loop=0;loop<LOOPMAX;loop++) {
79 Int_t rev = gRandom->Uniform(LOOPMAX);
80 T->SetBranchAddress("event", &events[loop]);
81 T->GetEntry(rev);
82 }
83 }
84}
85
86void jets(Int_t nev=100, Int_t npileup=200, Bool_t secondrun = true) {
87 // Embedding these loads inside the first run of the script is not yet
88 // supported in v6
89 // gROOT->ProcessLine(".L $ROOTSYS/tutorials/tree/JetEvent.cxx+");
90 write(nev);
91 read();
92 pileup(npileup);
93}
94
95#else
96
97//void jets(Int_t nev=100, Int_t npileup=200, Bool_t secondrun);
98void jets(Int_t nev=100, Int_t npileup=200) {
99 TString tutdir = gROOT->GetTutorialDir();
100 gROOT->ProcessLine(".L " + tutdir + "/tree/JetEvent.cxx+");
101 gROOT->ProcessLine("#define JETS_SECOND_RUN yes");
102 gROOT->ProcessLine("#include \"" __FILE__ "\"");
103 gROOT->ProcessLine("jets(100,200,true)");
104}
105
106#endif
#define f(i)
Definition: RSha256.hxx:104
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
long long Long64_t
Definition: RtypesCore.h:69
int nentries
Definition: THbookFile.cxx:89
#define gROOT
Definition: TROOT.h:415
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:48
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition: TRandom.cxx:635
Basic string class.
Definition: TString.h:131
A TTree represents a columnar dataset.
Definition: TTree.h:72
void jets()
Definition: jets.C:38
double T(double x)
Definition: ChebyshevPol.h:34