1 | /* |
---|
2 | root |
---|
3 | gSystem->Load("libDelphes"); |
---|
4 | .X examples/Example2.C("delphes_output.root"); |
---|
5 | */ |
---|
6 | |
---|
7 | //------------------------------------------------------------------------------ |
---|
8 | |
---|
9 | struct MyPlots |
---|
10 | { |
---|
11 | TH1 *fTrackDeltaPT; |
---|
12 | TH1 *fJetPT[2]; |
---|
13 | TH1 *fMissingET; |
---|
14 | TH1 *fElectronPT; |
---|
15 | }; |
---|
16 | |
---|
17 | //------------------------------------------------------------------------------ |
---|
18 | |
---|
19 | void BookHistograms(ExRootResult *result, MyPlots *plots) |
---|
20 | { |
---|
21 | THStack *stack; |
---|
22 | TLegend *legend; |
---|
23 | TPaveText *comment; |
---|
24 | |
---|
25 | // book 2 histograms for PT of 1st and 2nd leading jets |
---|
26 | |
---|
27 | plots->fJetPT[0] = result->AddHist1D("jet_pt_0", "leading jet P_{T}", |
---|
28 | "jet P_{T}, GeV/c", "number of jets", |
---|
29 | 50, 0.0, 100.0); |
---|
30 | |
---|
31 | plots->fJetPT[1] = result->AddHist1D("jet_pt_1", "2nd leading jet P_{T}", |
---|
32 | "jet P_{T}, GeV/c", "number of jets", |
---|
33 | 50, 0.0, 100.0); |
---|
34 | |
---|
35 | plots->fJetPT[0]->SetLineColor(kRed); |
---|
36 | plots->fJetPT[1]->SetLineColor(kBlue); |
---|
37 | |
---|
38 | // book 1 stack of 2 histograms |
---|
39 | |
---|
40 | stack = result->AddHistStack("jet_pt_all", "1st and 2nd jets P_{T}"); |
---|
41 | stack->Add(plots->fJetPT[0]); |
---|
42 | stack->Add(plots->fJetPT[1]); |
---|
43 | |
---|
44 | // book legend for stack of 2 histograms |
---|
45 | |
---|
46 | legend = result->AddLegend(0.72, 0.86, 0.98, 0.98); |
---|
47 | legend->AddEntry(plots->fJetPT[0], "leading jet", "l"); |
---|
48 | legend->AddEntry(plots->fJetPT[1], "second jet", "l"); |
---|
49 | |
---|
50 | // attach legend to stack (legend will be printed over stack in .eps file) |
---|
51 | |
---|
52 | result->Attach(stack, legend); |
---|
53 | |
---|
54 | // book more histograms |
---|
55 | |
---|
56 | plots->fMissingET = result->AddHist1D("missing_et", "Missing E_{T}", |
---|
57 | "Missing E_{T}, GeV", |
---|
58 | "number of events", |
---|
59 | 60, 0.0, 30.0); |
---|
60 | |
---|
61 | plots->fElectronPT = result->AddHist1D("electron_pt", "electron P_{T}", |
---|
62 | "electron P_{T}, GeV/c", |
---|
63 | "number of electrons", |
---|
64 | 50, 0.0, 100.0); |
---|
65 | |
---|
66 | // book general comment |
---|
67 | |
---|
68 | comment = result->AddComment(0.54, 0.72, 0.98, 0.98); |
---|
69 | comment->AddText("demonstration plot"); |
---|
70 | comment->AddText("produced by Example.C"); |
---|
71 | |
---|
72 | // attach comment to single histograms |
---|
73 | |
---|
74 | result->Attach(plots->fJetPT[0], comment); |
---|
75 | result->Attach(plots->fJetPT[1], comment); |
---|
76 | result->Attach(plots->fMissingET, comment); |
---|
77 | result->Attach(plots->fElectronPT, comment); |
---|
78 | } |
---|
79 | |
---|
80 | //------------------------------------------------------------------------------ |
---|
81 | |
---|
82 | void AnalyseEvents(ExRootTreeReader *treeReader, MyPlots *plots) |
---|
83 | { |
---|
84 | TClonesArray *branchJet = treeReader->UseBranch("Jet"); |
---|
85 | TClonesArray *branchElectron = treeReader->UseBranch("Electron"); |
---|
86 | TClonesArray *branchMissingET = treeReader->UseBranch("MissingET"); |
---|
87 | |
---|
88 | Long64_t allEntries = treeReader->GetEntries(); |
---|
89 | |
---|
90 | cout << "** Chain contains " << allEntries << " events" << endl; |
---|
91 | |
---|
92 | Jet *jet[2]; |
---|
93 | MissingET *met; |
---|
94 | Electron *electron; |
---|
95 | |
---|
96 | TIter itElectron(branchElectron); |
---|
97 | |
---|
98 | Long64_t entry; |
---|
99 | |
---|
100 | // Loop over all events |
---|
101 | for(entry = 0; entry < allEntries; ++entry) |
---|
102 | { |
---|
103 | // Load selected branches with data from specified event |
---|
104 | treeReader->ReadEntry(entry); |
---|
105 | |
---|
106 | // Analyse two leading jets |
---|
107 | if(branchJet->GetEntriesFast() >= 2) |
---|
108 | { |
---|
109 | jet[0] = (Jet*) branchJet->At(0); |
---|
110 | jet[1] = (Jet*) branchJet->At(1); |
---|
111 | |
---|
112 | plots->fJetPT[0]->Fill(jet[0]->PT); |
---|
113 | plots->fJetPT[1]->Fill(jet[1]->PT); |
---|
114 | } |
---|
115 | |
---|
116 | // Analyse missing ET |
---|
117 | if(branchMissingET->GetEntriesFast() > 0) |
---|
118 | { |
---|
119 | met = (MissingET*) branchMissingET->At(0); |
---|
120 | plots->fMissingET->Fill(met->MET); |
---|
121 | } |
---|
122 | |
---|
123 | // Loop over all electrons in event |
---|
124 | itElectron.Reset(); |
---|
125 | while((electron = (Electron*) itElectron.Next())) |
---|
126 | { |
---|
127 | plots->fElectronPT->Fill(electron->PT); |
---|
128 | } |
---|
129 | } |
---|
130 | } |
---|
131 | |
---|
132 | //------------------------------------------------------------------------------ |
---|
133 | |
---|
134 | void PrintHistograms(ExRootResult *result, MyPlots *plots) |
---|
135 | { |
---|
136 | result->Print(); |
---|
137 | } |
---|
138 | |
---|
139 | //------------------------------------------------------------------------------ |
---|
140 | |
---|
141 | void Example2(const char *inputFile) |
---|
142 | { |
---|
143 | TChain *chain = new TChain("Delphes"); |
---|
144 | chain->Add(inputFile); |
---|
145 | |
---|
146 | ExRootTreeReader *treeReader = new ExRootTreeReader(chain); |
---|
147 | ExRootResult *result = new ExRootResult(); |
---|
148 | |
---|
149 | MyPlots *plots = new MyPlots; |
---|
150 | |
---|
151 | BookHistograms(result, plots); |
---|
152 | |
---|
153 | AnalyseEvents(treeReader, plots); |
---|
154 | |
---|
155 | PrintHistograms(result, plots); |
---|
156 | |
---|
157 | result->Write("results.root"); |
---|
158 | |
---|
159 | cout << "** Exiting..." << endl; |
---|
160 | |
---|
161 | delete plots; |
---|
162 | delete result; |
---|
163 | delete treeReader; |
---|
164 | delete chain; |
---|
165 | } |
---|
166 | |
---|
167 | //------------------------------------------------------------------------------ |
---|
168 | |
---|