1 |
schiefer |
1.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/plot.h"
|
12 |
|
|
|
13 |
|
|
#include <TROOT.h>
|
14 |
|
|
#include <TRint.h>
|
15 |
|
|
#include <TFile.h>
|
16 |
|
|
#include <TTree.h>
|
17 |
|
|
#include <TEventList.h>
|
18 |
|
|
#include <TH1F.h>
|
19 |
|
|
#include <TF1.h>
|
20 |
|
|
#include <TGraphErrors.h>
|
21 |
|
|
#include <TMultiGraph.h>
|
22 |
|
|
#include <TCanvas.h>
|
23 |
|
|
#include <TLegend.h>
|
24 |
|
|
|
25 |
|
|
#include <iostream>
|
26 |
|
|
#include <iomanip>
|
27 |
|
|
#include <fstream>
|
28 |
|
|
#include <cmath>
|
29 |
|
|
#include <string>
|
30 |
|
|
#include <vector>
|
31 |
|
|
|
32 |
|
|
|
33 |
|
|
using namespace std;
|
34 |
|
|
|
35 |
|
|
|
36 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
37 |
|
|
// defaults
|
38 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
39 |
|
|
string dflt_ptbins ="7,9,12,15,20,25,30,45,65,90,120";
|
40 |
|
|
//string dflt_etabins=".2,.4,.7,1.,1.3,1.5,1.75,2.,2.3,2.5,2.7,2.9,3.2";
|
41 |
|
|
string dflt_etabins=".3,.6,.9,1.2,1.5,1.8,2.1,2.4,2.7";
|
42 |
|
|
string dflt_phibins="-2.75,-2.25,-1.75.,-1.25,-.75,-.25,.25,.75,1.25,1.75,2.25,2.75";
|
43 |
|
|
|
44 |
|
|
|
45 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
46 |
|
|
// declare global functions
|
47 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
48 |
|
|
|
49 |
|
|
/*
|
50 |
|
|
TH1F** initRspHistos(const string& jetalg,const string& varname,
|
51 |
|
|
const vector<float>& bins,const string& xtitle);
|
52 |
|
|
TH1F*** initRspHistos(const string& jetalg,
|
53 |
|
|
const string& varname1,const string& varname2,
|
54 |
|
|
const vector<float>& bins1,const vector<float>& bins2,
|
55 |
|
|
const string& xtitle);
|
56 |
|
|
TH1F** initVarHistos(const string& jetalg,const string& varname,
|
57 |
|
|
const vector<float>& bins,const string& xtitle);
|
58 |
|
|
TH1F*** initVarHistos(const string& jetalg,
|
59 |
|
|
const string& varname1,const string& varname2,
|
60 |
|
|
const vector<float>& bins1,const vector<float>& bins2,
|
61 |
|
|
const string& xtitle);
|
62 |
|
|
TH1F** initPtHistos (const string& jetalg,
|
63 |
|
|
const string& varname,const vector<float>& bins);
|
64 |
|
|
TH1F*** initPtHistos (const string& jetalg,
|
65 |
|
|
const string& varname1,const string& varname2,
|
66 |
|
|
const vector<float>& bins1,const vector<float>& bins2);
|
67 |
|
|
*/
|
68 |
|
|
|
69 |
|
|
void fillRspAndSgm(TGraphErrors* gRsp,TGraphErrors* gSgm,float nsigma,
|
70 |
|
|
const vector<float>& bins,TH1F** hVar,TH1F** hPt,TH1F** hRsp);
|
71 |
|
|
TCanvas* plotVarHistos(const vector<float>& bins,TH1F** hVar,bool logy=false);
|
72 |
|
|
TCanvas* plotRspHistos(const vector<float>& bins,TH1F** hRsp);
|
73 |
|
|
TCanvas* plotResponse (const vector<string>& jetalgs,short applyjes,
|
74 |
|
|
const vector<TGraphErrors*> gRsp,
|
75 |
|
|
const string& varname,const string& xtitle);
|
76 |
|
|
TCanvas* plotResponse (const string& jetalg,short applyjes,bool fitcorr,
|
77 |
|
|
const vector<float>& bins1,TGraphErrors** gRsp,
|
78 |
|
|
const string& varname1,const string vartitle1,
|
79 |
|
|
const string& varname2,const string& xtitle);
|
80 |
|
|
TCanvas* plotResolution(const vector<string>& jetalgs,
|
81 |
|
|
const vector<TGraphErrors*> gSgm,
|
82 |
|
|
const string& varname,const string& xtitle);
|
83 |
|
|
TCanvas* plotResolution(const string& jetalg,const vector<float>& bins1,
|
84 |
|
|
TGraphErrors** gSgm,
|
85 |
|
|
const string& varname1,const string vartitle1,
|
86 |
|
|
const string& varname2,const string& xtitle);
|
87 |
|
|
|
88 |
|
|
/*
|
89 |
|
|
int get_ibin(const vector<float>& bins,float value);
|
90 |
|
|
string get_bin_title(const string& varname,int ibin,const vector<float>& bins);
|
91 |
|
|
string get_bin_name(const string& varname,int ibin,const vector<float>& bins);
|
92 |
|
|
void get_bin_limits(float& min,float& max,int ibin,const vector<float>& bins);
|
93 |
|
|
*/
|
94 |
|
|
float get_xmax(TGraph* g);
|
95 |
|
|
|
96 |
|
|
|
97 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
98 |
|
|
// main
|
99 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
100 |
|
|
|
101 |
|
|
//______________________________________________________________________________
|
102 |
|
|
int main(int argc,char**argv)
|
103 |
|
|
{
|
104 |
|
|
cmdline cl;
|
105 |
|
|
cl.parse(argc,argv);
|
106 |
|
|
|
107 |
|
|
vector<string> input = cl.get_vector<string>("input");
|
108 |
|
|
string selection = cl.get_value<string> ("selection", "njt>0");
|
109 |
|
|
float drmax = cl.get_value<float> ("drmax", 0.3);
|
110 |
|
|
string treename = cl.get_value<string> ("treename", "t");
|
111 |
|
|
float nsigma = cl.get_value<float> ("nsigma", 2.0);
|
112 |
|
|
short applyjes = cl.get_value<short> ("applyjes", 0);
|
113 |
|
|
bool fitcorr = cl.get_value<bool> ("fitcorr", false);
|
114 |
|
|
vector<float> ptbins = cl.get_vector<float> ("ptbins", dflt_ptbins);
|
115 |
|
|
vector<float> etabins = cl.get_vector<float> ("etabins",dflt_etabins);
|
116 |
|
|
vector<float> phibins = cl.get_vector<float> ("phibins",dflt_phibins);
|
117 |
|
|
bool abseta = cl.get_value<bool> ("abseta", true);
|
118 |
|
|
|
119 |
|
|
if (!cl.check()) return 0;
|
120 |
|
|
cl.print();
|
121 |
|
|
|
122 |
|
|
|
123 |
|
|
// etabins
|
124 |
|
|
if (!abseta) {
|
125 |
|
|
int neta=(int)etabins.size();
|
126 |
|
|
std::reverse(etabins.begin(),etabins.end());
|
127 |
|
|
for (int ieta=neta-1;ieta>=0;ieta--) etabins.push_back(etabins[ieta]);
|
128 |
|
|
for (int ieta=0;ieta<neta;ieta++) etabins[ieta]*=-1;
|
129 |
|
|
}
|
130 |
|
|
|
131 |
|
|
|
132 |
|
|
// jetalgs
|
133 |
|
|
vector<string> jetalgs;
|
134 |
|
|
for (vector<string>::const_iterator it=input.begin();it!=input.end();++it) {
|
135 |
|
|
string jetalg = *it;
|
136 |
|
|
string::size_type pos;
|
137 |
|
|
while ((pos=jetalg.find('/'))!=string::npos) jetalg = jetalg.substr(pos+1);
|
138 |
|
|
jetalg = jetalg.substr(0,jetalg.find(".root"));
|
139 |
|
|
jetalgs.push_back(jetalg);
|
140 |
|
|
}
|
141 |
|
|
|
142 |
|
|
// prepare branch variables and histograms / graphs
|
143 |
|
|
float weight;
|
144 |
|
|
char njt;
|
145 |
|
|
float jtpt[100];
|
146 |
|
|
float jteta[100];
|
147 |
|
|
float jtphi[100];
|
148 |
|
|
float jtgenpt[100];
|
149 |
|
|
float jtgendr[100];
|
150 |
|
|
float jtjes[100][3];
|
151 |
|
|
|
152 |
|
|
int nptbins = ptbins.size()+1;
|
153 |
|
|
int netabins = etabins.size()+1;
|
154 |
|
|
int nphibins = phibins.size()+1;
|
155 |
|
|
|
156 |
|
|
vector<TH1F**> hPt;
|
157 |
|
|
vector<TH1F**> hPtCorr;
|
158 |
|
|
vector<TH1F**> hRspVsPt;
|
159 |
|
|
vector<TH1F**> hEta;
|
160 |
|
|
vector<TH1F**> hEtaPtCorr;
|
161 |
|
|
vector<TH1F**> hRspVsEta;
|
162 |
|
|
vector<TH1F**> hPhi;
|
163 |
|
|
vector<TH1F**> hPhiPtCorr;
|
164 |
|
|
vector<TH1F**> hRspVsPhi;
|
165 |
|
|
|
166 |
|
|
vector<TH1F***> hEtaPtPt;
|
167 |
|
|
vector<TH1F***> hEtaPtPtCorr;
|
168 |
|
|
vector<TH1F***> hRspVsEtaPt;
|
169 |
|
|
|
170 |
|
|
vector<TGraphErrors*> gRspVsPt;
|
171 |
|
|
vector<TGraphErrors*> gSgmVsPt;
|
172 |
|
|
vector<TGraphErrors*> gRspVsEta;
|
173 |
|
|
vector<TGraphErrors*> gSgmVsEta;
|
174 |
|
|
vector<TGraphErrors*> gRspVsPhi;
|
175 |
|
|
vector<TGraphErrors*> gSgmVsPhi;
|
176 |
|
|
|
177 |
|
|
vector<TGraphErrors**> gRspVsEtaPt;
|
178 |
|
|
vector<TGraphErrors**> gSgmVsEtaPt;
|
179 |
|
|
|
180 |
|
|
argc=1; //argv[1]="-b";
|
181 |
|
|
TRint* app = new TRint(argv[0],&argc,argv);
|
182 |
|
|
|
183 |
|
|
|
184 |
|
|
//
|
185 |
|
|
// loop over all input files
|
186 |
|
|
//
|
187 |
|
|
for (unsigned int i=0;i<input.size();++i) {
|
188 |
|
|
TFile* f = new TFile(input[i].c_str(),"READ"); if (!f->IsOpen()) return 0;
|
189 |
|
|
TTree* t = (TTree*)f->Get(treename.c_str());
|
190 |
|
|
|
191 |
|
|
TEventList* el = new TEventList("sel","sel");
|
192 |
|
|
t->Draw(">>sel",selection.c_str());
|
193 |
|
|
|
194 |
|
|
t->SetBranchAddress("weight", &weight);
|
195 |
|
|
t->SetBranchAddress("njt", &njt);
|
196 |
|
|
t->SetBranchAddress("jtpt", jtpt);
|
197 |
|
|
t->SetBranchAddress("jteta", jteta);
|
198 |
|
|
t->SetBranchAddress("jtphi", jtphi);
|
199 |
|
|
t->SetBranchAddress("jtgenpt",jtgenpt);
|
200 |
|
|
t->SetBranchAddress("jtgendr",jtgendr);
|
201 |
|
|
|
202 |
|
|
if (t->FindBranch("jtjes")) t->SetBranchAddress("jtjes", jtjes);
|
203 |
|
|
else
|
204 |
|
|
for (int i1=0;i1<100;i1++)
|
205 |
|
|
for (int i2=0;i2<3;i2++) jtjes[i1][i2]=1.0;
|
206 |
|
|
|
207 |
|
|
|
208 |
|
|
string rsp_xtitle = "E_{T}^{gen}-E_{T}^{rec} [GeV]";
|
209 |
|
|
|
210 |
|
|
/*
|
211 |
|
|
hPt .push_back(initVarHistos(jetalgs[i],"Pt", ptbins,"jet p_{T} [GeV]"));
|
212 |
|
|
hPtCorr .push_back(initPtHistos (jetalgs[i],"Pt", ptbins));
|
213 |
|
|
hRspVsPt .push_back(initRspHistos(jetalgs[i],"Pt", ptbins, rsp_xtitle));
|
214 |
|
|
hEta .push_back(initVarHistos(jetalgs[i],"Eta",etabins,"jet #eta"));
|
215 |
|
|
hEtaPtCorr.push_back(initPtHistos (jetalgs[i],"Eta",etabins));
|
216 |
|
|
hRspVsEta .push_back(initRspHistos(jetalgs[i],"Eta",etabins,rsp_xtitle));
|
217 |
|
|
hPhi .push_back(initVarHistos(jetalgs[i],"Phi",phibins,"jet #phi"));
|
218 |
|
|
hPhiPtCorr.push_back(initPtHistos (jetalgs[i],"Phi",phibins));
|
219 |
|
|
hRspVsPhi .push_back(initRspHistos(jetalgs[i],"Phi",phibins,rsp_xtitle));
|
220 |
|
|
|
221 |
|
|
hEtaPtPt.push_back(initVarHistos(jetalgs[i],
|
222 |
|
|
"Eta","Pt",etabins,ptbins,"jet p_{T} [GeV]"));
|
223 |
|
|
hEtaPtPtCorr.push_back(initPtHistos(jetalgs[i],
|
224 |
|
|
"Eta","Pt",etabins,ptbins));
|
225 |
|
|
hRspVsEtaPt.push_back(initRspHistos(jetalgs[i],
|
226 |
|
|
"Eta","Pt",etabins,ptbins,rsp_xtitle));
|
227 |
|
|
*/
|
228 |
|
|
|
229 |
|
|
hPt .push_back(initHistos("JetPt",jetalgs[i],
|
230 |
|
|
1000,0,4000,true,"jet p_{T} [GeV]",
|
231 |
|
|
"Pt",ptbins));
|
232 |
|
|
hPtCorr .push_back(initHistos("JetPtPtCorr",jetalgs[i],
|
233 |
|
|
1000,0,4000,true,"jet p_{T} [GeV]",
|
234 |
|
|
"Pt",ptbins));
|
235 |
|
|
hRspVsPt .push_back(initHistos("RspVsPt",jetalgs[i],
|
236 |
|
|
100,-100,100,true,rsp_xtitle,
|
237 |
|
|
"Pt",ptbins));
|
238 |
|
|
hEta .push_back(initHistos("JetEta",jetalgs[i],
|
239 |
|
|
100,0,5,false,"jet #eta",
|
240 |
|
|
"Eta",etabins));
|
241 |
|
|
hEtaPtCorr.push_back(initHistos("JetEtaPtCorr",jetalgs[i],
|
242 |
|
|
100,0,5,true,"jet p_{T} [GeV]",
|
243 |
|
|
"Eta",etabins));
|
244 |
|
|
hRspVsEta .push_back(initHistos("RspVsEta",jetalgs[i],
|
245 |
|
|
100,-100,100,true,rsp_xtitle,
|
246 |
|
|
"Eta",etabins));
|
247 |
|
|
hPhi .push_back(initHistos(jetalgs[i],"Phi",phibins,"jet #phi"));
|
248 |
|
|
hPhiPtCorr.push_back(initHistos (jetalgs[i],"Phi",phibins));
|
249 |
|
|
hRspVsPhi .push_back(initHistos(jetalgs[i],"Phi",phibins,rsp_xtitle));
|
250 |
|
|
|
251 |
|
|
hEtaPtPt.push_back(initVarHistos(jetalgs[i],
|
252 |
|
|
"Eta","Pt",etabins,ptbins,"jet p_{T} [GeV]"));
|
253 |
|
|
hEtaPtPtCorr.push_back(initPtHistos(jetalgs[i],
|
254 |
|
|
"Eta","Pt",etabins,ptbins));
|
255 |
|
|
hRspVsEtaPt.push_back(initRspHistos(jetalgs[i],
|
256 |
|
|
"Eta","Pt",etabins,ptbins,rsp_xtitle));
|
257 |
|
|
|
258 |
|
|
|
259 |
|
|
int nevts = el->GetN();
|
260 |
|
|
cout<<jetalgs[i]<<": "<<nevts<<" events selected."<<endl;
|
261 |
|
|
|
262 |
|
|
|
263 |
|
|
// loop over all events
|
264 |
|
|
for (int ievt=0;ievt<nevts;ievt++) {
|
265 |
|
|
|
266 |
|
|
t->GetEntry(el->GetEntry(ievt));
|
267 |
|
|
|
268 |
|
|
// loop over all jets in the event
|
269 |
|
|
for (int ijt=0;ijt<njt;ijt++) {
|
270 |
|
|
|
271 |
|
|
if (jtgendr[ijt]<0.0||jtgendr[ijt]>drmax) continue;
|
272 |
|
|
|
273 |
|
|
if (ijt>0) continue;
|
274 |
|
|
|
275 |
|
|
float pt = jtpt[ijt];
|
276 |
|
|
float genpt = jtgenpt[ijt];
|
277 |
|
|
float eta = (abseta) ? std::abs(jteta[ijt]) : jteta[ijt];
|
278 |
|
|
float phi = jtphi[ijt];
|
279 |
|
|
|
280 |
|
|
float ptcorr = (applyjes==0) ? pt : pt*jtjes[ijt][applyjes-1];
|
281 |
|
|
float deltapt = genpt - ptcorr;
|
282 |
|
|
|
283 |
|
|
// pt
|
284 |
|
|
int ipt = get_ibin(ptbins,pt);
|
285 |
|
|
hPt[i][ipt] ->Fill(pt,weight);
|
286 |
|
|
hPtCorr[i][ipt] ->Fill(ptcorr,weight);
|
287 |
|
|
hRspVsPt[i][ipt]->Fill(deltapt,weight);
|
288 |
|
|
|
289 |
|
|
// eta
|
290 |
|
|
int ieta = get_ibin(etabins,eta);
|
291 |
|
|
hEta[i][ieta] ->Fill(eta, weight);
|
292 |
|
|
hEtaPtCorr[i][ieta]->Fill(ptcorr, weight);
|
293 |
|
|
hRspVsEta[i][ieta] ->Fill(deltapt,weight);
|
294 |
|
|
|
295 |
|
|
// phi
|
296 |
|
|
int iphi = get_ibin(phibins,phi);
|
297 |
|
|
hPhi[i][iphi] ->Fill(phi, weight);
|
298 |
|
|
hPhiPtCorr[i][iphi]->Fill(ptcorr, weight);
|
299 |
|
|
hRspVsPhi[i][iphi] ->Fill(deltapt,weight);
|
300 |
|
|
|
301 |
|
|
// pt-eta
|
302 |
|
|
hEtaPtPt[i][ieta][ipt] ->Fill(pt, weight);
|
303 |
|
|
hEtaPtPtCorr[i][ieta][ipt]->Fill(ptcorr, weight);
|
304 |
|
|
hRspVsEtaPt[i][ieta][ipt] ->Fill(deltapt,weight);
|
305 |
|
|
|
306 |
|
|
|
307 |
|
|
} // jets
|
308 |
|
|
} // evts
|
309 |
|
|
|
310 |
|
|
cout<<jetalgs[i]<<": all events processed / histograms filled."<<endl;
|
311 |
|
|
|
312 |
|
|
|
313 |
|
|
// response / resolution vs *pt*
|
314 |
|
|
gRspVsPt.push_back(new TGraphErrors(nptbins));
|
315 |
|
|
gSgmVsPt.push_back(new TGraphErrors(nptbins));
|
316 |
|
|
gRspVsPt.back()->SetName(("RspVsPt_"+jetalgs[i]).c_str());
|
317 |
|
|
gSgmVsPt.back()->SetName(("SgmVsPt_"+jetalgs[i]).c_str());
|
318 |
|
|
fillRspAndSgm(gRspVsPt.back(),gSgmVsPt.back(),nsigma,
|
319 |
|
|
ptbins,hPt[i],hPtCorr[i],hRspVsPt[i]);
|
320 |
|
|
|
321 |
|
|
// response / resolution vs *eta*
|
322 |
|
|
gRspVsEta.push_back(new TGraphErrors(netabins));
|
323 |
|
|
gSgmVsEta.push_back(new TGraphErrors(netabins));
|
324 |
|
|
gRspVsEta.back()->SetName(("RspVsEta_"+jetalgs[i]).c_str());
|
325 |
|
|
gSgmVsEta.back()->SetName(("SgmVsEta_"+jetalgs[i]).c_str());
|
326 |
|
|
fillRspAndSgm(gRspVsEta.back(),gSgmVsEta.back(),nsigma,
|
327 |
|
|
etabins,hEta[i],hEtaPtCorr[i],hRspVsEta[i]);
|
328 |
|
|
|
329 |
|
|
// response / resolution vs *phi*
|
330 |
|
|
gRspVsPhi.push_back(new TGraphErrors(nphibins));
|
331 |
|
|
gSgmVsPhi.push_back(new TGraphErrors(nphibins));
|
332 |
|
|
gRspVsPhi.back()->SetName(("RspVsPhi_"+jetalgs[i]).c_str());
|
333 |
|
|
gSgmVsPhi.back()->SetName(("SgmVsPhi_"+jetalgs[i]).c_str());
|
334 |
|
|
fillRspAndSgm(gRspVsPhi.back(),gSgmVsPhi.back(),nsigma,
|
335 |
|
|
phibins,hPhi[i],hPhiPtCorr[i],hRspVsPhi[i]);
|
336 |
|
|
|
337 |
|
|
// response / resolution vs *pt-eta*
|
338 |
|
|
TGraphErrors** g;
|
339 |
|
|
g=new TGraphErrors*[netabins];
|
340 |
|
|
for (int ieta=0;ieta<netabins;ieta++) g[ieta] = new TGraphErrors(nptbins);
|
341 |
|
|
gRspVsEtaPt.push_back(g);
|
342 |
|
|
g=new TGraphErrors*[netabins];
|
343 |
|
|
for (int ieta=0;ieta<netabins;ieta++) g[ieta] = new TGraphErrors(nptbins);
|
344 |
|
|
gSgmVsEtaPt.push_back(g);
|
345 |
|
|
for (int ieta=0;ieta<netabins;ieta++) {
|
346 |
|
|
gRspVsEtaPt[i][ieta]->SetName(("RspVsPhi_"+get_bin_name("Eta",ieta,etabins)+
|
347 |
|
|
"_"+jetalgs[i]).c_str());
|
348 |
|
|
gSgmVsEtaPt[i][ieta]->SetName(("SgmVsPhi_"+get_bin_name("Eta",ieta,etabins)+
|
349 |
|
|
"_"+jetalgs[i]).c_str());
|
350 |
|
|
fillRspAndSgm(gRspVsEtaPt[i][ieta],gSgmVsEtaPt[i][ieta],nsigma,ptbins,
|
351 |
|
|
hEtaPtPt[i][ieta],hEtaPtPtCorr[i][ieta],hRspVsEtaPt[i][ieta]);
|
352 |
|
|
}
|
353 |
|
|
}
|
354 |
|
|
|
355 |
|
|
|
356 |
|
|
//
|
357 |
|
|
// make plots
|
358 |
|
|
//
|
359 |
|
|
TFile* f = new TFile("mcresponse.root","RECREATE");
|
360 |
|
|
|
361 |
|
|
for (unsigned int ialg=0;ialg<jetalgs.size();ialg++) {
|
362 |
|
|
|
363 |
|
|
// pt
|
364 |
|
|
TCanvas* cPt = plotVarHistos(ptbins,hPt[ialg],true);
|
365 |
|
|
TCanvas* cPtCorr = plotVarHistos(ptbins,hPtCorr[ialg],true);
|
366 |
|
|
TCanvas* cRspPt = plotRspHistos(ptbins,hRspVsPt[ialg]);
|
367 |
|
|
|
368 |
|
|
delete cPt;
|
369 |
|
|
delete cPtCorr;
|
370 |
|
|
delete cRspPt;
|
371 |
|
|
|
372 |
|
|
// eta
|
373 |
|
|
TCanvas* cEta = plotVarHistos(etabins,hEta[ialg]);
|
374 |
|
|
TCanvas* cEtaPtCorr = plotVarHistos(etabins,hEtaPtCorr[ialg],true);
|
375 |
|
|
TCanvas* cRspEta = plotRspHistos(etabins,hRspVsEta[ialg]);
|
376 |
|
|
|
377 |
|
|
delete cEta;
|
378 |
|
|
delete cEtaPtCorr;
|
379 |
|
|
delete cRspEta;
|
380 |
|
|
|
381 |
|
|
// phi
|
382 |
|
|
TCanvas* cPhi = plotVarHistos(phibins,hPhi[ialg]);
|
383 |
|
|
TCanvas* cPhiPtCorr = plotVarHistos(phibins,hPhiPtCorr[ialg],true);
|
384 |
|
|
TCanvas* cRspPhi = plotRspHistos(phibins,hRspVsPhi[ialg]);
|
385 |
|
|
|
386 |
|
|
delete cPhi;
|
387 |
|
|
delete cPhiPtCorr;
|
388 |
|
|
delete cRspPhi;
|
389 |
|
|
|
390 |
|
|
// eta-pt
|
391 |
|
|
for (int ieta=0;ieta<netabins;ieta++) {
|
392 |
|
|
TCanvas* cEtaPtPt = plotVarHistos(ptbins,hEtaPtPt[ialg][ieta],true);
|
393 |
|
|
TCanvas* cEtaPtPtCorr = plotVarHistos(ptbins,hEtaPtPtCorr[ialg][ieta],true);
|
394 |
|
|
TCanvas* cRspEtaPt = plotRspHistos(ptbins,hRspVsEtaPt[ialg][ieta]);
|
395 |
|
|
|
396 |
|
|
delete cEtaPtPt;
|
397 |
|
|
delete cEtaPtPtCorr;
|
398 |
|
|
delete cRspEtaPt;
|
399 |
|
|
}
|
400 |
|
|
|
401 |
|
|
}
|
402 |
|
|
|
403 |
|
|
// response & resolution vs pT
|
404 |
|
|
TCanvas* cRspVsPt = plotResponse(jetalgs,applyjes,gRspVsPt,"Pt","jet p_{T} [GeV]");
|
405 |
|
|
TCanvas* cSgmVsPt = plotResolution(jetalgs,gSgmVsPt,"Pt","jet p_{T} [GeV]");
|
406 |
|
|
|
407 |
|
|
// response & resolution vs eta
|
408 |
|
|
TCanvas* cRspVsEta = plotResponse(jetalgs,applyjes,gRspVsEta,"Eta","jet #eta");
|
409 |
|
|
TCanvas* cSgmVsEta = plotResolution(jetalgs,gSgmVsEta,"Eta","jet #eta");
|
410 |
|
|
|
411 |
|
|
// response & resolution vs phi
|
412 |
|
|
TCanvas* cRspVsPhi = plotResponse(jetalgs,applyjes,gRspVsPhi,"Phi","jet #phi");
|
413 |
|
|
TCanvas* cSgmVsPhi = plotResolution(jetalgs,gSgmVsPhi,"Phi","jet #phi");
|
414 |
|
|
|
415 |
|
|
// response & resolution vs eta/pT
|
416 |
|
|
for (unsigned int ialg=0;ialg<jetalgs.size();ialg++) {
|
417 |
|
|
|
418 |
|
|
TCanvas* cRspVsEtaPt = plotResponse(jetalgs[ialg],applyjes,fitcorr,
|
419 |
|
|
etabins,gRspVsEtaPt[ialg],
|
420 |
|
|
"Eta","#eta","Pt","jet p_{T} [GeV]");
|
421 |
|
|
TCanvas* cSgmVsEtaPt = plotResolution(jetalgs[ialg],etabins,gSgmVsEtaPt[ialg],
|
422 |
|
|
"Eta","#eta","Pt","jet p_{T} [GeV]");
|
423 |
|
|
}
|
424 |
|
|
|
425 |
|
|
f->Write();
|
426 |
|
|
f->Close();
|
427 |
|
|
delete f;
|
428 |
|
|
|
429 |
|
|
app->Run();
|
430 |
|
|
|
431 |
|
|
return 0;
|
432 |
|
|
}
|
433 |
|
|
|
434 |
|
|
|
435 |
|
|
|
436 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
437 |
|
|
// implementation of global functions
|
438 |
|
|
////////////////////////////////////////////////////////////////////////////////
|
439 |
|
|
|
440 |
|
|
//______________________________________________________________________________
|
441 |
|
|
TH1F** initRspHistos(const string& jetalg,const string& varname,
|
442 |
|
|
const vector<float>& bins,const string& xtitle)
|
443 |
|
|
{
|
444 |
|
|
unsigned int nhistos = bins.size()+1;
|
445 |
|
|
TH1F** result = new TH1F*[nhistos];
|
446 |
|
|
for (unsigned int i=0;i<nhistos;i++) {
|
447 |
|
|
string hname = "Rsp_"+get_bin_name(varname,i,bins)+"_"+jetalg;
|
448 |
|
|
result[i] = new TH1F(hname.c_str(),hname.c_str(),51,-100.,100.);
|
449 |
|
|
result[i]->SetXTitle(xtitle.c_str());
|
450 |
|
|
result[i]->SetYTitle("number of jets");
|
451 |
|
|
result[i]->Sumw2();
|
452 |
|
|
}
|
453 |
|
|
|
454 |
|
|
return result;
|
455 |
|
|
}
|
456 |
|
|
|
457 |
|
|
|
458 |
|
|
//______________________________________________________________________________
|
459 |
|
|
TH1F*** initRspHistos(const string& jetalg,
|
460 |
|
|
const string& varname1,const string& varname2,
|
461 |
|
|
const vector<float>& bins1,const vector<float>& bins2,
|
462 |
|
|
const string& xtitle)
|
463 |
|
|
{
|
464 |
|
|
unsigned int n1 = bins1.size()+1;
|
465 |
|
|
TH1F*** result = new TH1F**[n1];
|
466 |
|
|
for (unsigned int i1=0;i1<n1;i1++) {
|
467 |
|
|
unsigned int n2 = bins2.size()+1;
|
468 |
|
|
result[i1] = new TH1F*[n2];
|
469 |
|
|
for (unsigned int i2=0;i2<n2;i2++) {
|
470 |
|
|
string hname =
|
471 |
|
|
"Rsp_"+
|
472 |
|
|
get_bin_name(varname1,i1,bins1)+"_"+
|
473 |
|
|
get_bin_name(varname2,i2,bins2)+"_"+
|
474 |
|
|
jetalg;
|
475 |
|
|
result[i1][i2] = new TH1F(hname.c_str(),hname.c_str(),51,-100.,100.);
|
476 |
|
|
result[i1][i2]->SetXTitle(xtitle.c_str());
|
477 |
|
|
result[i1][i2]->SetYTitle("number of jets");
|
478 |
|
|
result[i1][i2]->Sumw2();
|
479 |
|
|
}
|
480 |
|
|
}
|
481 |
|
|
return result;
|
482 |
|
|
}
|
483 |
|
|
|
484 |
|
|
|
485 |
|
|
//______________________________________________________________________________
|
486 |
|
|
TH1F** initVarHistos(const string& jetalg,const string& varname,
|
487 |
|
|
const vector<float>& bins,const string& xtitle)
|
488 |
|
|
{
|
489 |
|
|
unsigned int nhistos = bins.size()+1;
|
490 |
|
|
TH1F** result = new TH1F*[nhistos];
|
491 |
|
|
for (unsigned int i=0;i<nhistos;i++) {
|
492 |
|
|
string hname = get_bin_name(varname,i,bins)+"_"+jetalg;
|
493 |
|
|
float xmin,xmax;
|
494 |
|
|
get_bin_limits(xmin,xmax,i,bins);
|
495 |
|
|
result[i] = new TH1F(hname.c_str(),hname.c_str(),21,xmin,xmax);
|
496 |
|
|
result[i]->SetXTitle(xtitle.c_str());
|
497 |
|
|
result[i]->SetYTitle("number of jets");
|
498 |
|
|
result[i]->Sumw2();
|
499 |
|
|
}
|
500 |
|
|
|
501 |
|
|
return result;
|
502 |
|
|
}
|
503 |
|
|
|
504 |
|
|
|
505 |
|
|
//______________________________________________________________________________
|
506 |
|
|
TH1F** initPtHistos(const string& jetalg,
|
507 |
|
|
const string& varname,const vector<float>& bins)
|
508 |
|
|
{
|
509 |
|
|
unsigned int nhistos = bins.size()+1;
|
510 |
|
|
TH1F** result = new TH1F*[nhistos];
|
511 |
|
|
for (unsigned int i=0;i<nhistos;i++) {
|
512 |
|
|
string hname = "PtCorr_"+get_bin_name(varname,i,bins)+"_"+jetalg;
|
513 |
|
|
result[i] = new TH1F(hname.c_str(),hname.c_str(),100,0.,300.);
|
514 |
|
|
result[i]->SetXTitle("jet p_{T} [GeV]");
|
515 |
|
|
result[i]->SetYTitle("number of jets");
|
516 |
|
|
result[i]->Sumw2();
|
517 |
|
|
}
|
518 |
|
|
|
519 |
|
|
return result;
|
520 |
|
|
}
|
521 |
|
|
|
522 |
|
|
|
523 |
|
|
//______________________________________________________________________________
|
524 |
|
|
TH1F*** initVarHistos(const string& jetalg,
|
525 |
|
|
const string& varname1,const string& varname2,
|
526 |
|
|
const vector<float>& bins1,const vector<float>& bins2,
|
527 |
|
|
const string& xtitle)
|
528 |
|
|
{
|
529 |
|
|
unsigned int n1 = bins1.size()+1;
|
530 |
|
|
unsigned int n2 = bins2.size()+1;
|
531 |
|
|
TH1F*** result = new TH1F**[n1];
|
532 |
|
|
for (unsigned int i1=0;i1<n1;i1++) {
|
533 |
|
|
result[i1] = new TH1F*[n2];
|
534 |
|
|
for (unsigned int i2=0;i2<n2;i2++) {
|
535 |
|
|
string hname =
|
536 |
|
|
get_bin_name(varname1,i1,bins1)+"_"+
|
537 |
|
|
get_bin_name(varname2,i2,bins2)+"_"+
|
538 |
|
|
jetalg;
|
539 |
|
|
float xmin,xmax;
|
540 |
|
|
get_bin_limits(xmin,xmax,i2,bins2);
|
541 |
|
|
result[i1][i2] = new TH1F(hname.c_str(),hname.c_str(),21,xmin,xmax);
|
542 |
|
|
result[i1][i2]->SetXTitle(xtitle.c_str());
|
543 |
|
|
result[i1][i2]->SetYTitle("number of jets");
|
544 |
|
|
result[i1][i2]->Sumw2();
|
545 |
|
|
}
|
546 |
|
|
}
|
547 |
|
|
|
548 |
|
|
return result;
|
549 |
|
|
}
|
550 |
|
|
|
551 |
|
|
|
552 |
|
|
//______________________________________________________________________________
|
553 |
|
|
TH1F*** initPtHistos(const string& jetalg,
|
554 |
|
|
const string& varname1,const string& varname2,
|
555 |
|
|
const vector<float>& bins1,const vector<float>& bins2)
|
556 |
|
|
|
557 |
|
|
{
|
558 |
|
|
unsigned int n1 = bins1.size()+1;
|
559 |
|
|
unsigned int n2 = bins2.size()+1;
|
560 |
|
|
TH1F*** result = new TH1F**[n1];
|
561 |
|
|
for (unsigned int i1=0;i1<n1;i1++) {
|
562 |
|
|
result[i1] = new TH1F*[n2];
|
563 |
|
|
for (unsigned int i2=0;i2<n2;i2++) {
|
564 |
|
|
string hname =
|
565 |
|
|
"PtCorr_"+
|
566 |
|
|
get_bin_name(varname1,i1,bins1)+"_"+
|
567 |
|
|
get_bin_name(varname2,i2,bins2)+"_"+
|
568 |
|
|
jetalg;
|
569 |
|
|
result[i1][i2] = new TH1F(hname.c_str(),hname.c_str(),100,0,300.);
|
570 |
|
|
result[i1][i2]->SetXTitle("jet p_{T} [GeV]");
|
571 |
|
|
result[i1][i2]->SetYTitle("number of jets");
|
572 |
|
|
result[i1][i2]->Sumw2();
|
573 |
|
|
}
|
574 |
|
|
}
|
575 |
|
|
|
576 |
|
|
return result;
|
577 |
|
|
}
|
578 |
|
|
|
579 |
|
|
|
580 |
|
|
//______________________________________________________________________________
|
581 |
|
|
void fillRspAndSgm(TGraphErrors* gRsp,TGraphErrors* gSgm,float nsigma,
|
582 |
|
|
const vector<float>& bins,TH1F** hVar,TH1F** hPt,TH1F** hRsp)
|
583 |
|
|
{
|
584 |
|
|
int nbins = bins.size()+1;
|
585 |
|
|
for (int ibin=0;ibin<nbins;ibin++) {
|
586 |
|
|
|
587 |
|
|
float var = hVar[ibin]->GetMean();
|
588 |
|
|
float pt = hPt[ibin]->GetMean();
|
589 |
|
|
TH1F* h = hRsp[ibin];
|
590 |
|
|
|
591 |
|
|
float de = h->GetMean();
|
592 |
|
|
float errde = h->GetMeanError();
|
593 |
|
|
float sgmde = h->GetRMS();
|
594 |
|
|
float errsgmde = h->GetRMSError();
|
595 |
|
|
|
596 |
|
|
if (h->GetEntries()>50) {
|
597 |
|
|
string fname = "fitfnc" + (string)h->GetName();
|
598 |
|
|
TF1* f = new TF1(fname.c_str(),"gaus",de-nsigma*sgmde,de+nsigma*sgmde);
|
599 |
|
|
f->SetLineColor(kRed);
|
600 |
|
|
f->SetLineWidth(1);
|
601 |
|
|
f->SetParameter(1,de);
|
602 |
|
|
f->SetParameter(2,sgmde);
|
603 |
|
|
h->Fit(f,"QR");
|
604 |
|
|
|
605 |
|
|
de = f->GetParameter(1);
|
606 |
|
|
errde = f->GetParError(1);
|
607 |
|
|
sgmde = f->GetParameter(2);
|
608 |
|
|
errsgmde = f->GetParError(2);
|
609 |
|
|
}
|
610 |
|
|
|
611 |
|
|
float rsp = pt/(pt+de);
|
612 |
|
|
float errrsp = pt/(pt+de)/(pt+de)*errde;
|
613 |
|
|
float sgm = pt/(pt+de)/(pt+de)*sgmde;
|
614 |
|
|
float errsgm = 2*pt*sgmde/(pt+de)/(pt+de)/(pt+de)*errde;
|
615 |
|
|
|
616 |
|
|
sgm /= rsp;
|
617 |
|
|
errsgm /= rsp;
|
618 |
|
|
|
619 |
|
|
gRsp->SetPoint(ibin,var,rsp);
|
620 |
|
|
gRsp->SetPointError(ibin,0.0,errrsp);
|
621 |
|
|
gSgm->SetPoint(ibin,var,sgm);
|
622 |
|
|
gSgm->SetPointError(ibin,0.0,errsgm);
|
623 |
|
|
}
|
624 |
|
|
|
625 |
|
|
}
|
626 |
|
|
|
627 |
|
|
|
628 |
|
|
//______________________________________________________________________________
|
629 |
|
|
TCanvas* plotVarHistos(const vector<float>& bins,TH1F** hVar,bool logy)
|
630 |
|
|
{
|
631 |
|
|
string::size_type pos;
|
632 |
|
|
string cname,jetalg;
|
633 |
|
|
|
634 |
|
|
cname = hVar[0]->GetName();
|
635 |
|
|
pos = cname.find_last_of('_');
|
636 |
|
|
jetalg = cname.substr(pos+1);
|
637 |
|
|
cname = cname.substr(0,pos);
|
638 |
|
|
pos = cname.find_last_of(':');
|
639 |
|
|
cname = cname.substr(0,pos);
|
640 |
|
|
cname += "_"+jetalg;
|
641 |
|
|
|
642 |
|
|
TCanvas* c = new TCanvas(cname.c_str(),cname.c_str());
|
643 |
|
|
|
644 |
|
|
int nbins = bins.size()+1;
|
645 |
|
|
int ncx= int(std::sqrt((float)nbins));
|
646 |
|
|
int ncy= ncx;
|
647 |
|
|
if (ncx*ncy<nbins) ncx++;
|
648 |
|
|
if (ncx*ncy<nbins) ncy++;
|
649 |
|
|
|
650 |
|
|
c->Divide(ncx,ncy);
|
651 |
|
|
|
652 |
|
|
for (int ibin=0;ibin<nbins;ibin++) {
|
653 |
|
|
c->cd(ibin+1);
|
654 |
|
|
gPad->SetLogy(logy);
|
655 |
|
|
hVar[ibin]->Draw("EHIST");
|
656 |
|
|
}
|
657 |
|
|
|
658 |
|
|
c->Write();
|
659 |
|
|
|
660 |
|
|
return c;
|
661 |
|
|
}
|
662 |
|
|
|
663 |
|
|
|
664 |
|
|
//______________________________________________________________________________
|
665 |
|
|
TCanvas* plotRspHistos(const vector<float>& bins,TH1F** hRsp)
|
666 |
|
|
{
|
667 |
|
|
string::size_type pos;
|
668 |
|
|
string cname,jetalg;
|
669 |
|
|
|
670 |
|
|
cname = hRsp[0]->GetName();
|
671 |
|
|
pos = cname.find_last_of('_');
|
672 |
|
|
jetalg = cname.substr(pos+1);
|
673 |
|
|
cname = cname.substr(0,pos);
|
674 |
|
|
pos = cname.find_last_of(':');
|
675 |
|
|
cname = cname.substr(0,pos);
|
676 |
|
|
cname += "_"+jetalg;
|
677 |
|
|
|
678 |
|
|
TCanvas* c = new TCanvas(cname.c_str(),cname.c_str());
|
679 |
|
|
|
680 |
|
|
int nbins = bins.size()+1;
|
681 |
|
|
int ncx= int(std::sqrt((float)nbins));
|
682 |
|
|
int ncy= ncx;
|
683 |
|
|
if (ncx*ncy<nbins) ncx++;
|
684 |
|
|
if (ncx*ncy<nbins) ncy++;
|
685 |
|
|
|
686 |
|
|
c->Divide(ncx,ncy);
|
687 |
|
|
|
688 |
|
|
for (int ibin=0;ibin<nbins;ibin++) {
|
689 |
|
|
c->cd(ibin+1);
|
690 |
|
|
if (hRsp[ibin]->GetEntries()>0) gPad->SetLogy();
|
691 |
|
|
hRsp[ibin]->Draw("EHIST");
|
692 |
|
|
string fname = "fitfnc" + (string)hRsp[ibin]->GetName();
|
693 |
|
|
TF1* fitfnc=hRsp[ibin]->GetFunction(fname.c_str());
|
694 |
|
|
if (0!=fitfnc) fitfnc->Draw("SAME");
|
695 |
|
|
}
|
696 |
|
|
|
697 |
|
|
c->Write();
|
698 |
|
|
|
699 |
|
|
return c;
|
700 |
|
|
}
|
701 |
|
|
|
702 |
|
|
|
703 |
|
|
//______________________________________________________________________________
|
704 |
|
|
TCanvas* plotResponse(const vector<string>& jetalgs,short applyjes,
|
705 |
|
|
const vector<TGraphErrors*> gRsp,
|
706 |
|
|
const string& varname,const string& xtitle)
|
707 |
|
|
{
|
708 |
|
|
assert(jetalgs.size()==gRsp.size());
|
709 |
|
|
|
710 |
|
|
string cname = "RspVs"+varname;
|
711 |
|
|
TCanvas* c = new TCanvas(cname.c_str(),cname.c_str());
|
712 |
|
|
TLegend* l = new TLegend(0.7,0.25+jetalgs.size()*0.045,0.85,0.25);
|
713 |
|
|
|
714 |
|
|
l->SetFillColor(10);
|
715 |
|
|
Color_t color=kBlack;
|
716 |
|
|
TMultiGraph* mg = new TMultiGraph();
|
717 |
|
|
|
718 |
|
|
for (unsigned int i=0;i<gRsp.size();i++) {
|
719 |
|
|
TGraphErrors* g = gRsp[i];
|
720 |
|
|
l->AddEntry(g,jetalgs[i].c_str(),"l");
|
721 |
|
|
while (color==kYellow||color==kWhite||color==10) color++;
|
722 |
|
|
g->SetMarkerColor(color);
|
723 |
|
|
g->SetLineColor(color);
|
724 |
|
|
g->SetMarkerStyle(kFullCircle);
|
725 |
|
|
g->SetMarkerSize(0.5);
|
726 |
|
|
mg->Add(g,"P");
|
727 |
|
|
|
728 |
|
|
if (varname=="Pt") {
|
729 |
|
|
string fname = "fitfnc"+cname;
|
730 |
|
|
string fitfncAsString =
|
731 |
|
|
//(applyjes==0) ? "[0]+[1]*log(x)+[2]*log(x)*log(x)+[3]*x+[4]*x*x" : "pol0";
|
732 |
|
|
(applyjes==0) ? "[0]*(pow(log(x),[1])-1/x)" : "pol0";
|
733 |
|
|
|
734 |
|
|
TF1* fitfnc = new TF1(fname.c_str(),fitfncAsString.c_str(),1,get_xmax(g));
|
735 |
|
|
|
736 |
|
|
if (applyjes==0) {
|
737 |
|
|
fitfnc->SetParameter(0,0.4);
|
738 |
|
|
fitfnc->SetParameter(1,0.4);
|
739 |
|
|
}
|
740 |
|
|
else {
|
741 |
|
|
fitfnc->SetParameter(0,1.0);
|
742 |
|
|
}
|
743 |
|
|
|
744 |
|
|
fitfnc->SetLineWidth(1);
|
745 |
|
|
fitfnc->SetLineColor(color);
|
746 |
|
|
g->Fit(fitfnc,"QR");
|
747 |
|
|
}
|
748 |
|
|
|
749 |
|
|
g->Write();
|
750 |
|
|
|
751 |
|
|
color++;
|
752 |
|
|
}
|
753 |
|
|
|
754 |
|
|
string mgtitle = "Jet Response Vs "+varname;
|
755 |
|
|
mg->Draw("a");
|
756 |
|
|
mg->SetTitle(mgtitle.c_str());
|
757 |
|
|
mg->GetXaxis()->SetTitle(xtitle.c_str());
|
758 |
|
|
mg->GetYaxis()->SetTitle("E_{T}^{rec}/E_{T}^{gen}");
|
759 |
|
|
if (varname=="Pt") mg->GetXaxis()->SetLimits(0.0,mg->GetXaxis()->GetXmax());
|
760 |
|
|
mg->SetMinimum(0.0);
|
761 |
|
|
mg->SetMaximum(1.1);
|
762 |
|
|
l->Draw();
|
763 |
|
|
|
764 |
|
|
c->Write();
|
765 |
|
|
|
766 |
|
|
return c;
|
767 |
|
|
}
|
768 |
|
|
|
769 |
|
|
|
770 |
|
|
//______________________________________________________________________________
|
771 |
|
|
TCanvas* plotResponse(const string& jetalg,short applyjes,bool fitcorr,
|
772 |
|
|
const vector<float>& bins1,TGraphErrors** gRsp,
|
773 |
|
|
const string& varname1,const string vartitle1,
|
774 |
|
|
const string& varname2,const string& xtitle)
|
775 |
|
|
|
776 |
|
|
{
|
777 |
|
|
int nbins = bins1.size()+1;
|
778 |
|
|
|
779 |
|
|
string fitfncAsString = "pol0";
|
780 |
|
|
int fitfncNbParams = 1;
|
781 |
|
|
|
782 |
|
|
ofstream fout;
|
783 |
|
|
|
784 |
|
|
if (applyjes==0) {
|
785 |
|
|
//fitfncAsString = "[0]+[1]*log(x)+[2]*log(x)*log(x)+[3]*x+[4]*x*x";
|
786 |
|
|
//fitfncNbParams = 5;
|
787 |
|
|
fitfncAsString = "[0]*(pow(log(x),[1])-1/x)";
|
788 |
|
|
fitfncNbParams = 2;
|
789 |
|
|
}
|
790 |
|
|
if (fitcorr) {
|
791 |
|
|
fout.open((jetalg+".txt").c_str()); assert(fout.is_open());
|
792 |
|
|
fout<<"#\n# jet correction parameters for "<<jetalg<<" algorithm\n#\n\n"
|
793 |
|
|
<<"# number of eta bins\n"<<nbins
|
794 |
|
|
<<"\n\n# fit function\n"<<fitfncAsString<<"\n"<<fitfncNbParams
|
795 |
|
|
<<endl<<endl;
|
796 |
|
|
}
|
797 |
|
|
|
798 |
|
|
string cname = "RspVs"+varname1+varname2+"_"+jetalg;
|
799 |
|
|
TCanvas* c = new TCanvas(cname.c_str(),cname.c_str());
|
800 |
|
|
TLegend* l = new TLegend(0.7,0.25+nbins*0.045,0.85,0.25);
|
801 |
|
|
|
802 |
|
|
l->SetFillColor(10);
|
803 |
|
|
Color_t color=kBlack;
|
804 |
|
|
TMultiGraph* mg = new TMultiGraph();
|
805 |
|
|
|
806 |
|
|
for (int ibin=0;ibin<nbins;ibin++) {
|
807 |
|
|
|
808 |
|
|
TGraphErrors* g = gRsp[ibin];
|
809 |
|
|
l->AddEntry(g,get_bin_title(vartitle1,ibin,bins1).c_str(),"l");
|
810 |
|
|
while (color==kYellow||color==kWhite||color==10) color++;
|
811 |
|
|
g->SetMarkerColor(color);
|
812 |
|
|
g->SetLineColor(color);
|
813 |
|
|
g->SetMarkerStyle(kFullCircle);
|
814 |
|
|
g->SetMarkerSize(0.5);
|
815 |
|
|
mg->Add(g,"P");
|
816 |
|
|
|
817 |
|
|
if (varname2=="Pt") {
|
818 |
|
|
string fname = "fitfnc"+cname;
|
819 |
|
|
TF1* fitfnc = new TF1(fname.c_str(),fitfncAsString.c_str(),1,get_xmax(g));
|
820 |
|
|
|
821 |
|
|
if (applyjes==0) {
|
822 |
|
|
fitfnc->SetParameter(0,0.4);
|
823 |
|
|
fitfnc->SetParameter(1,0.4);
|
824 |
|
|
}
|
825 |
|
|
else {
|
826 |
|
|
fitfnc->SetParameter(0,1.0);
|
827 |
|
|
}
|
828 |
|
|
|
829 |
|
|
fitfnc->SetLineWidth(1);
|
830 |
|
|
fitfnc->SetLineColor(color);
|
831 |
|
|
g->Fit(fitfnc,"QR");
|
832 |
|
|
|
833 |
|
|
if (fitcorr==0) {
|
834 |
|
|
if (ibin==nbins-1)
|
835 |
|
|
fout<<"# eta > "<<bins1.back()<<"\n999.9"<<endl;
|
836 |
|
|
else
|
837 |
|
|
fout<<"# eta < "<<bins1[ibin]<<"\n"<<bins1[ibin]<<endl;
|
838 |
|
|
for (int ipar=0;ipar<fitfnc->GetNpar();ipar++)
|
839 |
|
|
fout<<setw(12)<<fitfnc->GetParameter(ipar)<<" "
|
840 |
|
|
<<setw(12)<<fitfnc->GetParError(ipar)<<endl;
|
841 |
|
|
fout<<endl;
|
842 |
|
|
}
|
843 |
|
|
|
844 |
|
|
}
|
845 |
|
|
|
846 |
|
|
g->Write();
|
847 |
|
|
|
848 |
|
|
color++;
|
849 |
|
|
}
|
850 |
|
|
|
851 |
|
|
if (fout.is_open()) fout.close();
|
852 |
|
|
|
853 |
|
|
string mgtitle = jetalg+": Jet Response vs "+varname1+" and "+varname2;
|
854 |
|
|
mg->Draw("a");
|
855 |
|
|
mg->SetTitle(mgtitle.c_str());
|
856 |
|
|
mg->GetXaxis()->SetTitle(xtitle.c_str());
|
857 |
|
|
mg->GetYaxis()->SetTitle("E_{T}^{rec}/E_{T}^{gen}");
|
858 |
|
|
if (varname2=="Pt") mg->GetXaxis()->SetLimits(0.0,mg->GetXaxis()->GetXmax());
|
859 |
|
|
mg->SetMinimum(0.0);
|
860 |
|
|
mg->SetMaximum(1.1);
|
861 |
|
|
l->Draw();
|
862 |
|
|
|
863 |
|
|
c->Write();
|
864 |
|
|
|
865 |
|
|
return c;
|
866 |
|
|
}
|
867 |
|
|
|
868 |
|
|
|
869 |
|
|
//______________________________________________________________________________
|
870 |
|
|
TCanvas* plotResolution(const vector<string>& jetalgs,
|
871 |
|
|
const vector<TGraphErrors*> gSgm,
|
872 |
|
|
const string& varname,const string& xtitle)
|
873 |
|
|
{
|
874 |
|
|
string cname = "SgmVs"+varname;
|
875 |
|
|
TCanvas* c = new TCanvas(cname.c_str(),cname.c_str());
|
876 |
|
|
TLegend* l = new TLegend(0.7,0.85-jetalgs.size()*0.045,0.85,0.85);
|
877 |
|
|
|
878 |
|
|
l->SetFillColor(10);
|
879 |
|
|
Color_t color=kBlack;
|
880 |
|
|
TMultiGraph* mg = new TMultiGraph();
|
881 |
|
|
|
882 |
|
|
for (unsigned int i=0;i<gSgm.size();i++) {
|
883 |
|
|
TGraphErrors* g = gSgm[i];
|
884 |
|
|
l->AddEntry(g,jetalgs[i].c_str(),"l");
|
885 |
|
|
while (color==kYellow||color==kWhite||color==10) color++;
|
886 |
|
|
g->SetMarkerColor(color);
|
887 |
|
|
g->SetLineColor(color);
|
888 |
|
|
g->SetMarkerStyle(kFullCircle);
|
889 |
|
|
g->SetMarkerSize(0.5);
|
890 |
|
|
mg->Add(g,"P");
|
891 |
|
|
|
892 |
|
|
if (varname=="Pt") {
|
893 |
|
|
string fname = "fitfnc"+cname;
|
894 |
|
|
TF1* fitfnc = new TF1(fname.c_str(),"sqrt([0]*[0]/x*x+[1]*[1]/x+[2]*[2])");
|
895 |
|
|
fitfnc->SetParameter(0,0.0);
|
896 |
|
|
fitfnc->SetParameter(1,0.5);
|
897 |
|
|
fitfnc->SetParameter(2,0.1);
|
898 |
|
|
fitfnc->SetLineWidth(1);
|
899 |
|
|
fitfnc->SetLineColor(color);
|
900 |
|
|
fitfnc->SetRange(2.0,get_xmax(g));
|
901 |
|
|
g->Fit(fitfnc,"QR");
|
902 |
|
|
}
|
903 |
|
|
|
904 |
|
|
g->Write();
|
905 |
|
|
|
906 |
|
|
color++;
|
907 |
|
|
}
|
908 |
|
|
|
909 |
|
|
string mgtitle = "Jet Resolution Vs "+varname;
|
910 |
|
|
mg->Draw("a");
|
911 |
|
|
mg->SetTitle(mgtitle.c_str());
|
912 |
|
|
mg->GetXaxis()->SetTitle(xtitle.c_str());
|
913 |
|
|
mg->GetYaxis()
|
914 |
|
|
->SetTitle("#sigma(E_{T}^{rec}/E_{T}^{gen})/<E_{T}^{rec}/E_{T}^{gen}>");
|
915 |
|
|
if (varname=="Pt") mg->GetXaxis()->SetLimits(0.0,mg->GetXaxis()->GetXmax());
|
916 |
|
|
mg->SetMinimum(0.0);
|
917 |
|
|
mg->SetMaximum(0.5);
|
918 |
|
|
l->Draw();
|
919 |
|
|
|
920 |
|
|
c->Write();
|
921 |
|
|
|
922 |
|
|
return c;
|
923 |
|
|
}
|
924 |
|
|
|
925 |
|
|
|
926 |
|
|
//______________________________________________________________________________
|
927 |
|
|
TCanvas* plotResolution(const string& jetalg,const vector<float>& bins1,
|
928 |
|
|
TGraphErrors** gSgm,
|
929 |
|
|
const string& varname1,const string vartitle1,
|
930 |
|
|
const string& varname2,const string& xtitle)
|
931 |
|
|
{
|
932 |
|
|
int nbins = bins1.size()+1;
|
933 |
|
|
|
934 |
|
|
string cname = "SgmVs"+varname1+varname2+"_"+jetalg;
|
935 |
|
|
TCanvas* c = new TCanvas(cname.c_str(),cname.c_str());
|
936 |
|
|
TLegend* l = new TLegend(0.7,0.85-nbins*0.045,0.85,0.85);
|
937 |
|
|
|
938 |
|
|
l->SetFillColor(10);
|
939 |
|
|
Color_t color=kBlack;
|
940 |
|
|
TMultiGraph* mg = new TMultiGraph();
|
941 |
|
|
|
942 |
|
|
for (int ibin=0;ibin<nbins;ibin++) {
|
943 |
|
|
TGraphErrors* g = gSgm[ibin];
|
944 |
|
|
l->AddEntry(g,get_bin_title(vartitle1,ibin,bins1).c_str(),"l");
|
945 |
|
|
while (color==kYellow||color==kWhite||color==10) color++;
|
946 |
|
|
g->SetMarkerColor(color);
|
947 |
|
|
g->SetLineColor(color);
|
948 |
|
|
g->SetMarkerStyle(kFullCircle);
|
949 |
|
|
g->SetMarkerSize(0.5);
|
950 |
|
|
mg->Add(g,"P");
|
951 |
|
|
|
952 |
|
|
if (varname2=="Pt") {
|
953 |
|
|
string fname = "fitfnc"+cname;
|
954 |
|
|
TF1* fitfnc = new TF1(fname.c_str(),"sqrt([0]*[0]/x*x+[1]*[1]/x+[2]*[2])");
|
955 |
|
|
fitfnc->SetParameter(0,0.0);
|
956 |
|
|
fitfnc->SetParameter(1,0.5);
|
957 |
|
|
fitfnc->SetParameter(2,0.1);
|
958 |
|
|
fitfnc->SetLineWidth(1);
|
959 |
|
|
fitfnc->SetLineColor(color);
|
960 |
|
|
fitfnc->SetRange(2.0,get_xmax(g));
|
961 |
|
|
g->Fit(fitfnc,"QR");
|
962 |
|
|
}
|
963 |
|
|
|
964 |
|
|
g->Write();
|
965 |
|
|
|
966 |
|
|
color++;
|
967 |
|
|
}
|
968 |
|
|
|
969 |
|
|
string mgtitle=jetalg+": Jet Resolution Vs "+varname1+" and "+varname2;
|
970 |
|
|
mg->Draw("a");
|
971 |
|
|
mg->SetTitle(mgtitle.c_str());
|
972 |
|
|
mg->GetXaxis()->SetTitle(xtitle.c_str());
|
973 |
|
|
mg->GetYaxis()
|
974 |
|
|
->SetTitle("#sigma(E_{T}^{rec}/E_{T}^{gen})/<E_{T}^{rec}/E_{T}^{gen}>");
|
975 |
|
|
if (varname2=="Pt") mg->GetXaxis()->SetLimits(0.0,mg->GetXaxis()->GetXmax());
|
976 |
|
|
mg->SetMinimum(0.0);
|
977 |
|
|
mg->SetMaximum(0.5);
|
978 |
|
|
l->Draw();
|
979 |
|
|
|
980 |
|
|
c->Write();
|
981 |
|
|
|
982 |
|
|
return c;
|
983 |
|
|
}
|
984 |
|
|
|
985 |
|
|
|
986 |
|
|
//______________________________________________________________________________
|
987 |
|
|
int get_ibin(const vector<float>& bins,float value)
|
988 |
|
|
{
|
989 |
|
|
int result(-1);
|
990 |
|
|
if (value< bins.front()) result=0;
|
991 |
|
|
else if (value>=bins.back()) result=(int)bins.size();
|
992 |
|
|
else {
|
993 |
|
|
for (int i=0;i<(int)bins.size()-1;i++)
|
994 |
|
|
if (value>=bins[i]&&value<bins[i+1]) { result=i+1; break; }
|
995 |
|
|
if (result<0) cout<<"get_ibin: FATAL ERROR!"<<endl;
|
996 |
|
|
}
|
997 |
|
|
return result;
|
998 |
|
|
}
|
999 |
|
|
|
1000 |
|
|
|
1001 |
|
|
//______________________________________________________________________________
|
1002 |
|
|
string get_bin_name(const string& varname,int ibin,const vector<float>& bins)
|
1003 |
|
|
{
|
1004 |
|
|
stringstream ss;
|
1005 |
|
|
ss<<varname<<":";
|
1006 |
|
|
if (ibin==0) ss<<"lt"<<bins.front();
|
1007 |
|
|
else if (ibin==(int)bins.size()) ss<<"gt"<<bins.back();
|
1008 |
|
|
else ss<<bins[ibin-1]<<"to"<<bins[ibin];
|
1009 |
|
|
return ss.str();
|
1010 |
|
|
}
|
1011 |
|
|
|
1012 |
|
|
|
1013 |
|
|
//______________________________________________________________________________
|
1014 |
|
|
string get_bin_title(const string& varname,int ibin,const vector<float>& bins)
|
1015 |
|
|
{
|
1016 |
|
|
stringstream ss;
|
1017 |
|
|
if (ibin==0) ss<<varname<<" < "<<setw(4)<<bins.front();
|
1018 |
|
|
else if (ibin==(int)bins.size()) ss<<varname<<" > "<<setw(4)<<bins.back();
|
1019 |
|
|
else ss<<setw(4)<<bins[ibin-1]<<" #leq "+varname+" < "<<setw(4)<<bins[ibin];
|
1020 |
|
|
return ss.str();
|
1021 |
|
|
}
|
1022 |
|
|
|
1023 |
|
|
|
1024 |
|
|
//______________________________________________________________________________
|
1025 |
|
|
void get_bin_limits(float& min,float& max,int ibin,const vector<float>& bins)
|
1026 |
|
|
{
|
1027 |
|
|
if (ibin==0) {
|
1028 |
|
|
min=0.0;
|
1029 |
|
|
max=bins.front();
|
1030 |
|
|
}
|
1031 |
|
|
else if (ibin==(int)bins.size()) {
|
1032 |
|
|
min=bins.back();
|
1033 |
|
|
max=min+2.*(bins[ibin-1]-bins[ibin-2]);
|
1034 |
|
|
}
|
1035 |
|
|
else {
|
1036 |
|
|
min = bins[ibin-1];
|
1037 |
|
|
max = bins[ibin];
|
1038 |
|
|
}
|
1039 |
|
|
}
|
1040 |
|
|
|
1041 |
|
|
|
1042 |
|
|
//______________________________________________________________________________
|
1043 |
|
|
float get_xmax(TGraph* g)
|
1044 |
|
|
{
|
1045 |
|
|
float result(-1e200);
|
1046 |
|
|
for (int i=0;i<g->GetN();++i) if (g->GetX()[i]>result) result=g->GetX()[i];
|
1047 |
|
|
return result;
|
1048 |
|
|
}
|