1 |
////////////////////////////////////////////////////////////////////////////////
|
2 |
//
|
3 |
// mcresponse_x
|
4 |
// ------------
|
5 |
//
|
6 |
// 06/16/2007 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
|
7 |
////////////////////////////////////////////////////////////////////////////////
|
8 |
|
9 |
|
10 |
#include "utils/cmdline.h"
|
11 |
#include "utils/hist.h"
|
12 |
|
13 |
#include <TSystem.h>
|
14 |
#include <TROOT.h>
|
15 |
#include <TRint.h>
|
16 |
#include <TFile.h>
|
17 |
#include <TTree.h>
|
18 |
#include <TEventList.h>
|
19 |
#include <TH1F.h>
|
20 |
#include <TF1.h>
|
21 |
#include <TMath.h>
|
22 |
|
23 |
#include <iostream>
|
24 |
#include <iomanip>
|
25 |
#include <cmath>
|
26 |
#include <string>
|
27 |
#include <vector>
|
28 |
|
29 |
|
30 |
using namespace std;
|
31 |
|
32 |
|
33 |
////////////////////////////////////////////////////////////////////////////////
|
34 |
// declare global functions
|
35 |
////////////////////////////////////////////////////////////////////////////////
|
36 |
void replaceHistos(int nbins,
|
37 |
TH1F**& hVar,TH1F**& hEtCorr,TH1F**& hAbsRsp,TH1F**& hRelRsp,
|
38 |
TTree**& tVar);
|
39 |
void replaceHistos(int nbins1,int nbins2,
|
40 |
TH1F***&hVar,TH1F***&hEtCorr,TH1F***&hAbsRsp,TH1F***&hRelRsp,
|
41 |
TTree***& tVar);
|
42 |
void replaceHistos(int nbins1,int nbins2,int nbins3,
|
43 |
TH1F****&hVar,TH1F****&hEtCorr,
|
44 |
TH1F****&hAbsRsp,TH1F****&hRelRsp,
|
45 |
TTree****& tVar);
|
46 |
void fitHistos (TH1F** histos,int n,double nsigma);
|
47 |
void fitHistos (TH1F*** histos,int n1,int n2,double nsigma);
|
48 |
void fitHistos (TH1F**** histos,int n1,int n2,int n3,double nsigma);
|
49 |
void saveHistos (TH1F**& histos,int n);
|
50 |
void saveHistos (TH1F***& histos,int n1,int n2);
|
51 |
void saveHistos (TH1F****& histos,int n1,int n2,int n3);
|
52 |
|
53 |
|
54 |
|
55 |
////////////////////////////////////////////////////////////////////////////////
|
56 |
// main
|
57 |
////////////////////////////////////////////////////////////////////////////////
|
58 |
|
59 |
//______________________________________________________________________________
|
60 |
int main(int argc,char**argv)
|
61 |
{
|
62 |
cmdline cl;
|
63 |
cl.parse(argc,argv);
|
64 |
|
65 |
vector<string> input = cl.get_vector<string>("input");
|
66 |
string output = cl.get_value<string> ("output", "mcrsp.root");
|
67 |
string selection = cl.get_value<string> ("selection", "njt>0");
|
68 |
float drmax = cl.get_value<float> ("drmax", 0.3);
|
69 |
float jtetmin = cl.get_value<float> ("jtetmin", 0.0);
|
70 |
float jtetmax = cl.get_value<float> ("jtetmax", 1e06);
|
71 |
float jtetamin = cl.get_value<float> ("jtetamin", -9.9);
|
72 |
float jtetamax = cl.get_value<float> ("jtetamax", 9.9);
|
73 |
short applyjes = cl.get_value<short> ("applyjes", 0);
|
74 |
string treename = cl.get_value<string> ("treename", "t");
|
75 |
float nsigma = cl.get_value<float> ("nsigma", 3.0);
|
76 |
vector<float> etbins = cl.get_vector<float> ("etbins", "");
|
77 |
vector<float> etabins = cl.get_vector<float> ("etabins", "");
|
78 |
vector<float> phibins = cl.get_vector<float> ("phibins", "");
|
79 |
vector<float> emfbins = cl.get_vector<float> ("emfbins", "");
|
80 |
bool abseta = cl.get_value<bool> ("abseta", true);
|
81 |
|
82 |
if (!cl.check()) return 0;
|
83 |
cl.print();
|
84 |
|
85 |
int netbins = (etbins.size() >0) ? etbins.size() +1 : 0;
|
86 |
int netabins = (etabins.size()>0) ? etabins.size()+1 : 0;
|
87 |
int nphibins = (phibins.size()>0) ? phibins.size()+1 : 0;
|
88 |
int nemfbins = (emfbins.size()>0) ? emfbins.size()+1 : 0;
|
89 |
|
90 |
if (netbins+netabins+nphibins+nemfbins==0) {
|
91 |
cout<<"Must specify bins: etbins, etabins, phibins, AND/OR emfbins!"<<endl;
|
92 |
return 0;
|
93 |
}
|
94 |
|
95 |
// etabins
|
96 |
if (netabins>0&&!abseta) {
|
97 |
int neta=(int)etabins.size();
|
98 |
std::reverse(etabins.begin(),etabins.end());
|
99 |
for (int ieta=neta-1;ieta>=0;ieta--) etabins.push_back(etabins[ieta]);
|
100 |
for (int ieta=0;ieta<neta;ieta++) etabins[ieta]*=-1;
|
101 |
}
|
102 |
|
103 |
|
104 |
// jetalgs
|
105 |
vector<string> jetalgs;
|
106 |
for (vector<string>::const_iterator it=input.begin();it!=input.end();++it) {
|
107 |
string jetalg = *it;
|
108 |
string::size_type pos;
|
109 |
while ((pos=jetalg.find('/'))!=string::npos) jetalg = jetalg.substr(pos+1);
|
110 |
jetalg = jetalg.substr(0,jetalg.find(".root"));
|
111 |
jetalgs.push_back(jetalg);
|
112 |
}
|
113 |
|
114 |
// prepare branch variables and histograms / graphs
|
115 |
float weight;
|
116 |
char njt;
|
117 |
float jtet[100];
|
118 |
float jteta[100];
|
119 |
float jtphi[100];
|
120 |
float jtemf[100];
|
121 |
float jtgenet[100];
|
122 |
float jtgendr[100];
|
123 |
short jtgenflv[100];
|
124 |
float jtjes[100][3];
|
125 |
|
126 |
argc=3; argv[1]="-l"; argv[2]="-b";
|
127 |
TRint* app = new TRint(argv[0],&argc,argv); app->Argc();
|
128 |
|
129 |
|
130 |
//
|
131 |
// loop over all input files
|
132 |
//
|
133 |
TFile* plotfile = new TFile(output.c_str(),"UPDATE");
|
134 |
|
135 |
for (unsigned int i=0;i<input.size();++i) {
|
136 |
TFile* f = new TFile(input[i].c_str(),"READ"); if (!f->IsOpen()) return 0;
|
137 |
TTree* t = (TTree*)f->Get(treename.c_str());
|
138 |
|
139 |
TEventList* el = new TEventList("sel","sel");
|
140 |
t->Draw(">>sel",selection.c_str());
|
141 |
|
142 |
t->SetBranchAddress("weight", &weight);
|
143 |
t->SetBranchAddress("njt", &njt);
|
144 |
t->SetBranchAddress("jtet", jtet);
|
145 |
t->SetBranchAddress("jteta", jteta);
|
146 |
t->SetBranchAddress("jtphi", jtphi);
|
147 |
t->SetBranchAddress("jtemf", jtemf);
|
148 |
t->SetBranchAddress("jtgenet", jtgenet);
|
149 |
t->SetBranchAddress("jtgendr", jtgendr);
|
150 |
t->SetBranchAddress("jtgenflv",jtgenflv);
|
151 |
|
152 |
vector<TBranch*> branches;
|
153 |
branches.push_back(t->GetBranch("weight"));
|
154 |
branches.push_back(t->GetBranch("njt"));
|
155 |
branches.push_back(t->GetBranch("jtet"));
|
156 |
branches.push_back(t->GetBranch("jteta"));
|
157 |
branches.push_back(t->GetBranch("jtphi"));
|
158 |
branches.push_back(t->GetBranch("jtemf"));
|
159 |
branches.push_back(t->GetBranch("jtgenet"));
|
160 |
branches.push_back(t->GetBranch("jtgendr"));
|
161 |
branches.push_back(t->GetBranch("jtgenflv"));
|
162 |
|
163 |
if (t->FindBranch("jtjes")) {
|
164 |
t->SetBranchAddress("jtjes", jtjes);
|
165 |
branches.push_back(t->GetBranch("jtjes"));
|
166 |
}
|
167 |
else {
|
168 |
for (int i1=0;i1<100;i1++)
|
169 |
for (int i2=0;i2<3;i2++) jtjes[i1][i2]=1.0;
|
170 |
}
|
171 |
string absrsp_xtitle = "E_{T}^{gen}-E_{T} [GeV]";
|
172 |
string relrsp_xtitle = "E_{T}/E_{T}^{gen}";
|
173 |
|
174 |
|
175 |
TH1F** hGenEt(0);
|
176 |
TH1F** hEtCorrGenEt(0);
|
177 |
TH1F** hAbsRspGenEt(0);
|
178 |
TH1F** hRelRspGenEt(0);
|
179 |
|
180 |
if (netbins>0) {
|
181 |
hGenEt = hist::initHistos("GenEt",jetalgs[i],100,
|
182 |
"E_{T}^{gen} [GeV]",
|
183 |
"GenEt",etbins);
|
184 |
hEtCorrGenEt = hist::initHistos("EtCorrGenEt",jetalgs[i],100,
|
185 |
"E_{T}^{gen} [GeV]",
|
186 |
"GenEt",etbins);
|
187 |
hAbsRspGenEt = hist::initHistos("AbsRspGenEt",jetalgs[i],
|
188 |
100,-50,150,
|
189 |
absrsp_xtitle,
|
190 |
"GenEt",etbins);
|
191 |
hRelRspGenEt = hist::initHistos("RelRspGenEt",jetalgs[i],
|
192 |
100,0,2,
|
193 |
relrsp_xtitle,
|
194 |
"GenEt",etbins);
|
195 |
}
|
196 |
|
197 |
|
198 |
TH1F** hEt(0);
|
199 |
TH1F** hEtCorrEt(0);
|
200 |
TH1F** hAbsRspEt(0);
|
201 |
TH1F** hRelRspEt(0);
|
202 |
|
203 |
if (netbins>0) {
|
204 |
hEt = hist::initHistos("Et",jetalgs[i],100,
|
205 |
"E_{T} [GeV]","Et",etbins);
|
206 |
hEtCorrEt = hist::initHistos("EtCorrEt",jetalgs[i],100,
|
207 |
"E_{T} [GeV]","Et",etbins);
|
208 |
hAbsRspEt = hist::initHistos("AbsRspEt",jetalgs[i],100,-50,150,
|
209 |
absrsp_xtitle,"Et",etbins);
|
210 |
hRelRspEt = hist::initHistos("RelRspEt",jetalgs[i],100,0,2,
|
211 |
relrsp_xtitle,"Et",etbins);
|
212 |
}
|
213 |
|
214 |
|
215 |
TH1F** hEta(0);
|
216 |
TH1F** hEtCorrEta(0);
|
217 |
TH1F** hAbsRspEta(0);
|
218 |
TH1F** hRelRspEta(0);
|
219 |
|
220 |
if (netabins>0) {
|
221 |
hEta = hist::initHistos("Eta",jetalgs[i],100,
|
222 |
"#eta","Eta",etabins);
|
223 |
hEtCorrEta = hist::initHistos("EtCorrEta",jetalgs[i],100,
|
224 |
"E_{T} [GeV]","Eta",etabins);
|
225 |
hAbsRspEta = hist::initHistos("AbsRspEta",jetalgs[i],100,-50,150,
|
226 |
absrsp_xtitle,"Eta",etabins);
|
227 |
hRelRspEta = hist::initHistos("RelRspEta",jetalgs[i],100,0,2,
|
228 |
relrsp_xtitle,"Eta",etabins);
|
229 |
}
|
230 |
|
231 |
|
232 |
TH1F** hPhi(0);
|
233 |
TH1F** hEtCorrPhi(0);
|
234 |
TH1F** hAbsRspPhi(0);
|
235 |
TH1F** hRelRspPhi(0);
|
236 |
|
237 |
if (nphibins>0) {
|
238 |
hPhi = hist::initHistos("Phi",jetalgs[i],100,
|
239 |
"#phi","Phi",phibins);
|
240 |
hEtCorrPhi = hist::initHistos("EtCorrPhi",jetalgs[i],100,
|
241 |
"E_{T} [Gev]","Phi",phibins);
|
242 |
hAbsRspPhi = hist::initHistos("AbsRspPhi",jetalgs[i],100,-50,150,
|
243 |
absrsp_xtitle,"Phi",phibins);
|
244 |
hRelRspPhi = hist::initHistos("RelRspPhi",jetalgs[i],100,0,2,
|
245 |
relrsp_xtitle,"Phi",phibins);
|
246 |
}
|
247 |
|
248 |
|
249 |
TH1F** hEmf(0);
|
250 |
TH1F** hEtCorrEmf(0);
|
251 |
TH1F** hAbsRspEmf(0);
|
252 |
TH1F** hRelRspEmf(0);
|
253 |
|
254 |
if (nemfbins>0) {
|
255 |
hEmf = hist::initHistos("Emf",jetalgs[i],100,
|
256 |
"emf","Emf",emfbins);
|
257 |
hEtCorrEmf = hist::initHistos("EtCorrEmf",jetalgs[i],100,
|
258 |
"E_{T} [Gev]","Emf",emfbins);
|
259 |
hAbsRspEmf = hist::initHistos("AbsRspEmf",jetalgs[i],100,-50,150,
|
260 |
absrsp_xtitle,"Emf",emfbins);
|
261 |
hRelRspEmf = hist::initHistos("RelRspEmf",jetalgs[i],100,0,2,
|
262 |
relrsp_xtitle,"Emf",emfbins);
|
263 |
}
|
264 |
|
265 |
|
266 |
TH1F*** hEtEtaEt(0);
|
267 |
TH1F*** hEtCorrEtaEt(0);
|
268 |
TH1F*** hAbsRspEtaEt(0);
|
269 |
TH1F*** hRelRspEtaEt(0);
|
270 |
|
271 |
if (netbins>0&&netabins>0) {
|
272 |
hEtEtaEt = hist::initHistos("Et",jetalgs[i],100,"E_{T} [GeV]",
|
273 |
"Eta",etabins,"Et",etbins);
|
274 |
hEtCorrEtaEt = hist::initHistos("EtCorrEtaEt",jetalgs[i],100,
|
275 |
"E_{T} [GeV]",
|
276 |
"Eta",etabins,"Et",etbins);
|
277 |
hAbsRspEtaEt = hist::initHistos("AbsRspEtaEt",jetalgs[i],100,-50,150,
|
278 |
absrsp_xtitle,
|
279 |
"Eta",etabins,"Et",etbins);
|
280 |
hRelRspEtaEt = hist::initHistos("RelRspEtaEt",jetalgs[i],100,0,2,
|
281 |
relrsp_xtitle,
|
282 |
"Eta",etabins,"Et",etbins);
|
283 |
}
|
284 |
|
285 |
|
286 |
TH1F*** hEtEmfEt(0);
|
287 |
TH1F*** hEtCorrEmfEt(0);
|
288 |
TH1F*** hAbsRspEmfEt(0);
|
289 |
TH1F*** hRelRspEmfEt(0);
|
290 |
|
291 |
if (netbins>0&&nemfbins>0) {
|
292 |
hEtEmfEt = hist::initHistos("Et",jetalgs[i],100,"E_{T} [GeV]",
|
293 |
"Emf",emfbins,"Et",etbins);
|
294 |
hEtCorrEmfEt = hist::initHistos("EtCorrEmfEt",jetalgs[i],100,
|
295 |
"E_{T} [GeV]",
|
296 |
"Emf",emfbins,"Et",etbins);
|
297 |
hAbsRspEmfEt = hist::initHistos("AbsRspEmfEt",jetalgs[i],100,-50,150,
|
298 |
absrsp_xtitle,
|
299 |
"Emf",emfbins,"Et",etbins);
|
300 |
hRelRspEmfEt = hist::initHistos("RelRspEmfEt",jetalgs[i],100,0,2,
|
301 |
relrsp_xtitle,
|
302 |
"Emf",emfbins,"Et",etbins);
|
303 |
}
|
304 |
|
305 |
|
306 |
TH1F**** hEtEmfEtaEt(0);
|
307 |
TH1F**** hEtCorrEmfEtaEt(0);
|
308 |
TH1F**** hAbsRspEmfEtaEt(0);
|
309 |
TH1F**** hRelRspEmfEtaEt(0);
|
310 |
|
311 |
if (netbins>0&&netabins>0&&nemfbins>0) {
|
312 |
hEtEmfEtaEt = hist::initHistos("Et",jetalgs[i],100,"E_{T} [GeV]",
|
313 |
"Emf",emfbins,
|
314 |
"Eta",etabins,
|
315 |
"Et", etbins);
|
316 |
hEtCorrEmfEtaEt = hist::initHistos("EtCorrEmfEtaEt",jetalgs[i],100,
|
317 |
"E_{T} [GeV]",
|
318 |
"Emf",emfbins,
|
319 |
"Eta",etabins,
|
320 |
"Et", etbins);
|
321 |
hAbsRspEmfEtaEt = hist::initHistos("AbsRspEmfEtaEt",jetalgs[i],
|
322 |
100,-50,150,absrsp_xtitle,
|
323 |
"Emf",emfbins,
|
324 |
"Eta",etabins,
|
325 |
"Et", etbins);
|
326 |
hRelRspEmfEtaEt = hist::initHistos("RelRspEmfEtaEt",jetalgs[i],
|
327 |
100,0,2,relrsp_xtitle,
|
328 |
"Emf",emfbins,
|
329 |
"Eta",etabins,
|
330 |
"Et", etbins);
|
331 |
}
|
332 |
|
333 |
int nevts = el->GetN();
|
334 |
cout<<jetalgs[i]<<": "<<nevts<<" events selected."<<endl;
|
335 |
|
336 |
// prepare intermediate trees
|
337 |
float dr,genet,et,eta,phi,emf,etcorr,absrsp,relrsp;
|
338 |
short flv;
|
339 |
|
340 |
TTree** tGenEt(0);
|
341 |
TTree** tEt(0);
|
342 |
TTree** tEta(0);
|
343 |
TTree** tPhi(0);
|
344 |
TTree** tEmf(0);
|
345 |
TTree*** tEtEtaEt(0);
|
346 |
TTree*** tEtEmfEt(0);
|
347 |
TTree**** tEtEmfEtaEt(0);
|
348 |
|
349 |
gROOT->cd();
|
350 |
|
351 |
if (netbins>0) tGenEt = new TTree*[netbins];
|
352 |
if (netbins>0) tEt = new TTree*[netbins];
|
353 |
if (netabins>0) tEta = new TTree*[netabins];
|
354 |
if (nphibins>0) tPhi = new TTree*[nphibins];
|
355 |
if (nemfbins>0) tEmf = new TTree*[nemfbins];
|
356 |
if (netbins>0&&netabins>0) tEtEtaEt = new TTree**[netabins];
|
357 |
if (netbins>0&&nemfbins>0) tEtEmfEt = new TTree**[nemfbins];
|
358 |
if (netbins>0&&netabins>0&&nemfbins>0) tEtEmfEtaEt = new TTree***[nemfbins];
|
359 |
|
360 |
for (int iet=0;iet<netbins;iet++) {
|
361 |
|
362 |
tGenEt[iet]=new TTree("tGenEt","tGenEt");
|
363 |
tGenEt[iet]->Branch("val", &genet, "val/F");
|
364 |
tGenEt[iet]->Branch("etcorr",&genet, "etcorr/F");
|
365 |
tGenEt[iet]->Branch("absrsp",&absrsp,"absrsp/F");
|
366 |
tGenEt[iet]->Branch("relrsp",&relrsp,"relrsp/F");
|
367 |
tGenEt[iet]->Branch("weight",&weight,"weight/F");
|
368 |
|
369 |
tEt[iet]=new TTree("tEt","tEt");
|
370 |
tEt[iet]->Branch("val", &et, "val/F");
|
371 |
tEt[iet]->Branch("etcorr",&etcorr,"etcorr/F");
|
372 |
tEt[iet]->Branch("absrsp",&absrsp,"absrsp/F");
|
373 |
tEt[iet]->Branch("relrsp",&relrsp,"relrsp/F");
|
374 |
tEt[iet]->Branch("weight",&weight,"weight/F");
|
375 |
}
|
376 |
for (int ieta=0;ieta<netabins;ieta++) {
|
377 |
tEta[ieta]=new TTree("tEta","tEta");
|
378 |
tEta[ieta]->Branch("val", &eta, "val/F");
|
379 |
tEta[ieta]->Branch("etcorr",&etcorr,"etcorr/F");
|
380 |
tEta[ieta]->Branch("absrsp",&absrsp,"absrsp/F");
|
381 |
tEta[ieta]->Branch("relrsp",&relrsp,"relrsp/F");
|
382 |
tEta[ieta]->Branch("weight",&weight,"weight/F");
|
383 |
}
|
384 |
for (int iphi=0;iphi<nphibins;iphi++) {
|
385 |
tPhi[iphi]=new TTree("tPhi","tPhi");
|
386 |
tPhi[iphi]->Branch("val", &phi, "val/F");
|
387 |
tPhi[iphi]->Branch("etcorr",&etcorr,"etcorr/F");
|
388 |
tPhi[iphi]->Branch("absrsp",&absrsp,"absrsp/F");
|
389 |
tPhi[iphi]->Branch("relrsp",&relrsp,"relrsp/F");
|
390 |
tPhi[iphi]->Branch("weight",&weight,"weight/F");
|
391 |
}
|
392 |
for (int iemf=0;iemf<nemfbins;iemf++) {
|
393 |
tEmf[iemf]=new TTree("tEmf","tEmf");
|
394 |
tEmf[iemf]->Branch("val", &emf, "val/F");
|
395 |
tEmf[iemf]->Branch("etcorr",&etcorr,"etcorr/F");
|
396 |
tEmf[iemf]->Branch("absrsp",&absrsp,"absrsp/F");
|
397 |
tEmf[iemf]->Branch("relrsp",&relrsp,"relrsp/F");
|
398 |
tEmf[iemf]->Branch("weight",&weight,"weight/F");
|
399 |
}
|
400 |
|
401 |
if (netbins>0) {
|
402 |
for (int ieta=0;ieta<netabins;ieta++) {
|
403 |
tEtEtaEt[ieta] = new TTree*[netbins];
|
404 |
for (int iet=0;iet<netbins;iet++) {
|
405 |
tEtEtaEt[ieta][iet] = new TTree("tEtEtaEt","tEtEtaEt");
|
406 |
tEtEtaEt[ieta][iet]->Branch("val", &et, "val/F");
|
407 |
tEtEtaEt[ieta][iet]->Branch("etcorr",&etcorr,"etcorr/F");
|
408 |
tEtEtaEt[ieta][iet]->Branch("absrsp",&absrsp,"absrsp/F");
|
409 |
tEtEtaEt[ieta][iet]->Branch("relrsp",&relrsp,"relrsp/F");
|
410 |
tEtEtaEt[ieta][iet]->Branch("weight",&weight,"weight/F");
|
411 |
}
|
412 |
}
|
413 |
for (int iemf=0;iemf<nemfbins;iemf++) {
|
414 |
tEtEmfEt[iemf] = new TTree*[netbins];
|
415 |
for (int iet=0;iet<netbins;iet++) {
|
416 |
tEtEmfEt[iemf][iet] = new TTree("tEtEmfEt","tEtEmfEt");
|
417 |
tEtEmfEt[iemf][iet]->Branch("val", &et, "val/F");
|
418 |
tEtEmfEt[iemf][iet]->Branch("etcorr",&etcorr,"etcorr/F");
|
419 |
tEtEmfEt[iemf][iet]->Branch("absrsp",&absrsp,"absrsp/F");
|
420 |
tEtEmfEt[iemf][iet]->Branch("relrsp",&relrsp,"relrsp/F");
|
421 |
tEtEmfEt[iemf][iet]->Branch("weight",&weight,"weight/F");
|
422 |
}
|
423 |
}
|
424 |
if (netabins>0) {
|
425 |
for (int iemf=0;iemf<nemfbins;iemf++) {
|
426 |
tEtEmfEtaEt[iemf] = new TTree**[netabins];
|
427 |
for (int ieta=0;ieta<netabins;ieta++) {
|
428 |
tEtEmfEtaEt[iemf][ieta] = new TTree*[netbins];
|
429 |
for (int iet=0;iet<netbins;iet++) {
|
430 |
tEtEmfEtaEt[iemf][ieta][iet] = new TTree("tEtEmfEtaEt",
|
431 |
"tEtEmfEtaEt");
|
432 |
tEtEmfEtaEt[iemf][ieta][iet]->Branch("val", &et, "val/F");
|
433 |
tEtEmfEtaEt[iemf][ieta][iet]->Branch("etcorr",&etcorr,"etcorr/F");
|
434 |
tEtEmfEtaEt[iemf][ieta][iet]->Branch("absrsp",&absrsp,"absrsp/F");
|
435 |
tEtEmfEtaEt[iemf][ieta][iet]->Branch("relrsp",&relrsp,"relrsp/F");
|
436 |
tEtEmfEtaEt[iemf][ieta][iet]->Branch("weight",&weight,"weight/F");
|
437 |
}
|
438 |
}
|
439 |
}
|
440 |
}
|
441 |
}
|
442 |
|
443 |
cout<<jetalgs[i]<<": trees created."<<endl;
|
444 |
|
445 |
// loop over all events
|
446 |
for (int ievt=0;ievt<nevts;ievt++) {
|
447 |
|
448 |
//t->GetEntry(el->GetEntry(ievt));
|
449 |
for (int ibrnch=0;ibrnch<(int)branches.size();ibrnch++)
|
450 |
branches[ibrnch]->GetEntry(el->GetEntry(ievt));
|
451 |
|
452 |
// loop over all jets in the event
|
453 |
for (int ijt=0;ijt<njt;ijt++) {
|
454 |
|
455 |
et = jtet[ijt];
|
456 |
genet = jtgenet[ijt];
|
457 |
eta = (abseta) ? std::abs(jteta[ijt]) : jteta[ijt];
|
458 |
phi = jtphi[ijt];
|
459 |
emf = jtemf[ijt];
|
460 |
dr = jtgendr[ijt];
|
461 |
flv = std::abs(jtgenflv[ijt]);
|
462 |
|
463 |
if (dr <0.0 ||dr >drmax) continue;
|
464 |
if (et <jtetmin ||et >jtetmax) continue;
|
465 |
if (eta<jtetamin||eta>jtetamax) continue;
|
466 |
// TEST
|
467 |
//if (flv!=21) continue;
|
468 |
|
469 |
etcorr = (applyjes==0) ? et : et*jtjes[ijt][applyjes-1];
|
470 |
|
471 |
absrsp = genet - etcorr;
|
472 |
relrsp = etcorr/genet;
|
473 |
|
474 |
if (TMath::IsNaN(emf)) emf = 0.0;
|
475 |
|
476 |
int igenet = (netbins>0) ? hist::get_ibin(genet,etbins) : -1;
|
477 |
int iet = (netbins>0) ? hist::get_ibin(et,etbins) : -1;
|
478 |
int ieta = (netabins>0) ? hist::get_ibin(eta,etabins) : -1;
|
479 |
int iphi = (nphibins>0) ? hist::get_ibin(phi,phibins) : -1;
|
480 |
int iemf = (nemfbins>0) ? hist::get_ibin(emf,emfbins) : -1;
|
481 |
|
482 |
if (netbins>0) tGenEt[igenet]->Fill();
|
483 |
if (netbins>0) tEt[iet]->Fill();
|
484 |
if (netabins>0) tEta[ieta]->Fill();
|
485 |
if (nphibins>0) tPhi[iphi]->Fill();
|
486 |
if (nemfbins>0) tEmf[iemf]->Fill();
|
487 |
|
488 |
if (netbins>0&&netabins>0) tEtEtaEt[ieta][iet]->Fill();
|
489 |
if (netbins>0&&nemfbins>0) tEtEmfEt[iemf][iet]->Fill();
|
490 |
|
491 |
if (netbins>0&&netabins>0&&nemfbins>0)
|
492 |
tEtEmfEtaEt[iemf][ieta][iet]->Fill();
|
493 |
|
494 |
} // jets
|
495 |
} // evts
|
496 |
|
497 |
cout<<jetalgs[i]<<": all events processed / histograms filled."<<endl;
|
498 |
|
499 |
|
500 |
// replace histograms
|
501 |
|
502 |
replaceHistos(netbins,
|
503 |
hGenEt,hEtCorrGenEt,hAbsRspGenEt,hRelRspGenEt,tGenEt);
|
504 |
replaceHistos(netbins,hEt,hEtCorrEt,hAbsRspEt,hRelRspEt,tEt);
|
505 |
replaceHistos(netabins,hEta,hEtCorrEta,hAbsRspEta,hRelRspEta,tEta);
|
506 |
replaceHistos(nphibins,hPhi,hEtCorrPhi,hAbsRspPhi,hRelRspPhi,tPhi);
|
507 |
replaceHistos(nemfbins,hEmf,hEtCorrEmf,hAbsRspEmf,hRelRspEmf,tEmf);
|
508 |
replaceHistos(netabins,netbins,
|
509 |
hEtEtaEt,hEtCorrEtaEt,hAbsRspEtaEt,hRelRspEtaEt,
|
510 |
tEtEtaEt);
|
511 |
replaceHistos(nemfbins,netbins,
|
512 |
hEtEmfEt,hEtCorrEmfEt,hAbsRspEmfEt,hRelRspEmfEt,
|
513 |
tEtEmfEt);
|
514 |
replaceHistos(nemfbins,netabins,netbins,
|
515 |
hEtEmfEtaEt,hEtCorrEmfEtaEt,
|
516 |
hAbsRspEmfEtaEt,hRelRspEmfEtaEt,
|
517 |
tEtEmfEtaEt);
|
518 |
|
519 |
cout<<jetalgs[i]<<": all histograms replaced, and tree deleted."<<endl;
|
520 |
|
521 |
|
522 |
plotfile->cd();
|
523 |
TDirectory* d = plotfile->GetDirectory(jetalgs[i].c_str());
|
524 |
if (0!=d) plotfile->rmdir(jetalgs[i].c_str());
|
525 |
d = plotfile->mkdir(jetalgs[i].c_str());
|
526 |
d->cd();
|
527 |
|
528 |
// fit response histos
|
529 |
fitHistos(hAbsRspGenEt,netbins,nsigma);
|
530 |
fitHistos(hRelRspGenEt,netbins,nsigma);
|
531 |
|
532 |
fitHistos(hAbsRspEt, netbins,nsigma);
|
533 |
fitHistos(hRelRspEt, netbins,nsigma);
|
534 |
|
535 |
fitHistos(hAbsRspEta, netabins,nsigma);
|
536 |
fitHistos(hRelRspEta, netabins,nsigma);
|
537 |
|
538 |
fitHistos(hAbsRspPhi, nphibins,nsigma);
|
539 |
fitHistos(hRelRspPhi, nphibins,nsigma);
|
540 |
|
541 |
fitHistos(hAbsRspEmf, nemfbins,nsigma);
|
542 |
fitHistos(hRelRspEmf, nemfbins,nsigma);
|
543 |
|
544 |
fitHistos(hAbsRspEtaEt,netabins,netbins,nsigma);
|
545 |
fitHistos(hRelRspEtaEt,netabins,netbins,nsigma);
|
546 |
|
547 |
fitHistos(hAbsRspEmfEt,nemfbins,netbins,nsigma);
|
548 |
fitHistos(hRelRspEmfEt,nemfbins,netbins,nsigma);
|
549 |
|
550 |
fitHistos(hAbsRspEmfEtaEt,nemfbins,netabins,netbins,nsigma);
|
551 |
fitHistos(hRelRspEmfEtaEt,nemfbins,netabins,netbins,nsigma);
|
552 |
|
553 |
cout<<jetalgs[i]<<": fits done."<<endl;
|
554 |
|
555 |
|
556 |
// save histograms
|
557 |
saveHistos(hGenEt, netbins);
|
558 |
saveHistos(hEtCorrGenEt,netbins);
|
559 |
saveHistos(hAbsRspGenEt,netbins);
|
560 |
saveHistos(hRelRspGenEt,netbins);
|
561 |
|
562 |
saveHistos(hEt, netbins);
|
563 |
saveHistos(hEtCorrEt, netbins);
|
564 |
saveHistos(hAbsRspEt, netbins);
|
565 |
saveHistos(hRelRspEt, netbins);
|
566 |
|
567 |
saveHistos(hEta, netabins);
|
568 |
saveHistos(hEtCorrEta, netabins);
|
569 |
saveHistos(hAbsRspEta, netabins);
|
570 |
saveHistos(hRelRspEta, netabins);
|
571 |
|
572 |
saveHistos(hPhi, nphibins);
|
573 |
saveHistos(hEtCorrPhi, nphibins);
|
574 |
saveHistos(hAbsRspPhi, nphibins);
|
575 |
saveHistos(hRelRspPhi, nphibins);
|
576 |
|
577 |
saveHistos(hEmf, nemfbins);
|
578 |
saveHistos(hEtCorrEmf, nemfbins);
|
579 |
saveHistos(hAbsRspEmf, nemfbins);
|
580 |
saveHistos(hRelRspEmf, nemfbins);
|
581 |
|
582 |
saveHistos(hEtEtaEt, netabins,netbins);
|
583 |
saveHistos(hEtCorrEtaEt,netabins,netbins);
|
584 |
saveHistos(hAbsRspEtaEt,netabins,netbins);
|
585 |
saveHistos(hRelRspEtaEt,netabins,netbins);
|
586 |
|
587 |
saveHistos(hEtEmfEt, nemfbins,netbins);
|
588 |
saveHistos(hEtCorrEmfEt,nemfbins,netbins);
|
589 |
saveHistos(hAbsRspEmfEt,nemfbins,netbins);
|
590 |
saveHistos(hRelRspEmfEt,nemfbins,netbins);
|
591 |
|
592 |
saveHistos(hEtEmfEtaEt, nemfbins,netabins,netbins);
|
593 |
saveHistos(hEtCorrEmfEtaEt,nemfbins,netabins,netbins);
|
594 |
saveHistos(hAbsRspEmfEtaEt,nemfbins,netabins,netbins);
|
595 |
saveHistos(hRelRspEmfEtaEt,nemfbins,netabins,netbins);
|
596 |
|
597 |
cout<<jetalgs[i]<<": histograms saved."<<endl;
|
598 |
}
|
599 |
|
600 |
plotfile->Write();
|
601 |
plotfile->Close();
|
602 |
delete plotfile;
|
603 |
|
604 |
//app->Run();
|
605 |
|
606 |
return 0;
|
607 |
}
|
608 |
|
609 |
|
610 |
|
611 |
////////////////////////////////////////////////////////////////////////////////
|
612 |
// implementation of global functions
|
613 |
////////////////////////////////////////////////////////////////////////////////
|
614 |
|
615 |
//______________________________________________________________________________
|
616 |
void replaceHistos(int nbins,
|
617 |
TH1F**& hVar,TH1F**& hEtCorr,TH1F**& hAbsRsp,TH1F**& hRelRsp,
|
618 |
TTree**& tVar)
|
619 |
{
|
620 |
if (nbins==0) return;
|
621 |
|
622 |
TH1::SetDefaultSumw2();
|
623 |
for (int ibin=0;ibin<nbins;ibin++) {
|
624 |
|
625 |
TH1F* h(0);
|
626 |
|
627 |
if (tVar[ibin]->GetEntries()==0) continue;
|
628 |
|
629 |
tVar[ibin]->Project("h(50)","val","weight*(1)");
|
630 |
h = (TH1F*)gROOT->FindObject("h");
|
631 |
h->SetNameTitle(hVar[ibin]->GetName(),hVar[ibin]->GetTitle());
|
632 |
h->SetXTitle(hVar[ibin]->GetXaxis()->GetTitle());
|
633 |
h->SetYTitle(hVar[ibin]->GetYaxis()->GetTitle());
|
634 |
delete hVar[ibin];
|
635 |
hVar[ibin] = h;
|
636 |
|
637 |
tVar[ibin]->Project("h(50)","etcorr","weight*(1)");
|
638 |
h = (TH1F*)gROOT->FindObject("h");
|
639 |
h->SetNameTitle(hEtCorr[ibin]->GetName(),hEtCorr[ibin]->GetTitle());
|
640 |
h->SetXTitle(hEtCorr[ibin]->GetXaxis()->GetTitle());
|
641 |
h->SetYTitle(hEtCorr[ibin]->GetYaxis()->GetTitle());
|
642 |
delete hEtCorr[ibin];
|
643 |
hEtCorr[ibin] = h;
|
644 |
|
645 |
tVar[ibin]->Project("h(50)","absrsp","weight*(1)");
|
646 |
h = (TH1F*)gROOT->FindObject("h");
|
647 |
h->SetNameTitle(hAbsRsp[ibin]->GetName(),hAbsRsp[ibin]->GetTitle());
|
648 |
h->SetXTitle(hAbsRsp[ibin]->GetXaxis()->GetTitle());
|
649 |
h->SetYTitle(hAbsRsp[ibin]->GetYaxis()->GetTitle());
|
650 |
delete hAbsRsp[ibin];
|
651 |
hAbsRsp[ibin] = h;
|
652 |
|
653 |
tVar[ibin]->Project("h(50)","relrsp","weight*(1)");
|
654 |
h = (TH1F*)gROOT->FindObject("h");
|
655 |
h->SetNameTitle(hRelRsp[ibin]->GetName(),hRelRsp[ibin]->GetTitle());
|
656 |
h->SetXTitle(hRelRsp[ibin]->GetXaxis()->GetTitle());
|
657 |
h->SetYTitle(hRelRsp[ibin]->GetYaxis()->GetTitle());
|
658 |
delete hRelRsp[ibin];
|
659 |
hRelRsp[ibin] = h;
|
660 |
|
661 |
delete tVar[ibin];
|
662 |
}
|
663 |
|
664 |
delete [] tVar;
|
665 |
}
|
666 |
|
667 |
|
668 |
//______________________________________________________________________________
|
669 |
void replaceHistos(int nbins1,int nbins2,
|
670 |
TH1F***& hVar,TH1F***& hEtCorr,TH1F***& hAbsRsp,TH1F***& hRelRsp,
|
671 |
TTree***& tVar)
|
672 |
{
|
673 |
if (nbins1==0) return;
|
674 |
|
675 |
for (int i1=0;i1<nbins1;i1++)
|
676 |
replaceHistos(nbins2,hVar[i1],hEtCorr[i1],hAbsRsp[i1],hRelRsp[i1],tVar[i1]);
|
677 |
delete [] tVar;
|
678 |
}
|
679 |
|
680 |
|
681 |
//______________________________________________________________________________
|
682 |
void replaceHistos(int nbins1,int nbins2,int nbins3,
|
683 |
TH1F****& hVar,TH1F****& hEtCorr,TH1F****&hAbsRsp,TH1F****&hRelRsp,
|
684 |
TTree****& tVar)
|
685 |
{
|
686 |
if (nbins1==0) return;
|
687 |
for (int i1=0;i1<nbins1;i1++)
|
688 |
replaceHistos(nbins2,nbins3,
|
689 |
hVar[i1],hEtCorr[i1],hAbsRsp[i1],hRelRsp[i1],tVar[i1]);
|
690 |
delete [] tVar;
|
691 |
}
|
692 |
|
693 |
|
694 |
//______________________________________________________________________________
|
695 |
void fitHistos(TH1F** histos,int n,double nsigma)
|
696 |
{
|
697 |
for (int i=0;i<n;i++) {
|
698 |
float nrm = histos[i]->Integral();
|
699 |
float mean = histos[i]->GetMean();
|
700 |
float rms = histos[i]->GetRMS();
|
701 |
TF1* f = new TF1("gaus","gaus",mean-nsigma*rms,mean+nsigma*rms);
|
702 |
f->SetParameter(0,nrm);
|
703 |
f->SetParameter(1,mean);
|
704 |
f->SetParameter(2,rms);
|
705 |
f->SetLineWidth(1);
|
706 |
f->SetLineColor(kRed);
|
707 |
histos[i]->Fit(f,"QR"); // Avoid TCanvas?
|
708 |
}
|
709 |
}
|
710 |
|
711 |
|
712 |
//______________________________________________________________________________
|
713 |
void fitHistos(TH1F*** histos,int n1,int n2,double nsigma)
|
714 |
{
|
715 |
if (n1==0||n2==0) return;
|
716 |
for (int i1=0;i1<n1;i1++) fitHistos(histos[i1],n2,nsigma);
|
717 |
}
|
718 |
|
719 |
|
720 |
//______________________________________________________________________________
|
721 |
void fitHistos(TH1F**** histos,int n1,int n2,int n3,double nsigma)
|
722 |
{
|
723 |
if (n1==0||n2==0||n3==0) return;
|
724 |
for (int i1=0;i1<n1;i1++) fitHistos(histos[i1],n2,n3,nsigma);
|
725 |
}
|
726 |
|
727 |
|
728 |
//______________________________________________________________________________
|
729 |
void saveHistos(TH1F**& histos,int n)
|
730 |
{
|
731 |
if (n==0) return;
|
732 |
for (int i=0;i<n;i++) {
|
733 |
histos[i]->Write();
|
734 |
delete histos[i];
|
735 |
}
|
736 |
delete [] histos;
|
737 |
}
|
738 |
|
739 |
|
740 |
//______________________________________________________________________________
|
741 |
void saveHistos(TH1F***& histos,int n1,int n2)
|
742 |
{
|
743 |
if (n1==0) return;
|
744 |
for (int i1=0;i1<n1;i1++) saveHistos(histos[i1],n2);
|
745 |
delete [] histos;
|
746 |
}
|
747 |
|
748 |
|
749 |
//______________________________________________________________________________
|
750 |
void saveHistos(TH1F****& histos,int n1,int n2,int n3)
|
751 |
{
|
752 |
if (n1==0) return;
|
753 |
for (int i1=0;i1<n1;i1++) saveHistos(histos[i1],n2,n3);
|
754 |
delete [] histos;
|
755 |
}
|