1 |
#define JetResolutions_cxx
|
2 |
// The class definition in JetResolutions.h has been generated automatically
|
3 |
// by the ROOT utility TTree::MakeSelector(). This class is derived
|
4 |
// from the ROOT class TSelector. For more information on the TSelector
|
5 |
// framework see $ROOTSYS/README/README.SELECTOR or the ROOT User Manual.
|
6 |
|
7 |
// The following methods are defined in this file:
|
8 |
// Begin(): called every time a loop on the tree starts,
|
9 |
// a convenient place to create your histograms.
|
10 |
// SlaveBegin(): called after Begin(), when on PROOF called only on the
|
11 |
// slave servers.
|
12 |
// Process(): called for each event, in this function you decide what
|
13 |
// to read and fill your histograms.
|
14 |
// SlaveTerminate: called at the end of the loop on the tree, when on PROOF
|
15 |
// called only on the slave servers.
|
16 |
// Terminate(): called at the end of the loop on the tree,
|
17 |
// a convenient place to draw/fit your histograms.
|
18 |
//
|
19 |
// To use this file, try the following session on your Tree T:
|
20 |
//
|
21 |
// Root > T->Process("JetResolutions.C")
|
22 |
// Root > T->Process("JetResolutions.C","some options")
|
23 |
// Root > T->Process("JetResolutions.C+")
|
24 |
//
|
25 |
|
26 |
#include "JetResolutions.h"
|
27 |
|
28 |
void JetResolutions::Begin(TTree * /*tree*/) {
|
29 |
// The Begin() function is called at the start of the query.
|
30 |
// When running with PROOF Begin() is only called on the client.
|
31 |
// The tree argument is deprecated (on PROOF 0 is passed).
|
32 |
|
33 |
TString option = GetOption();
|
34 |
EBinEdges.push_back(0);
|
35 |
EBinEdges.push_back(20);
|
36 |
EBinEdges.push_back(30);
|
37 |
EBinEdges.push_back(50);
|
38 |
EBinEdges.push_back(80);
|
39 |
EBinEdges.push_back(120);
|
40 |
EBinEdges.push_back(170);
|
41 |
EBinEdges.push_back(230);
|
42 |
EBinEdges.push_back(300);
|
43 |
EBinEdges.push_back(380);
|
44 |
EBinEdges.push_back(470);
|
45 |
EBinEdges.push_back(570);
|
46 |
EBinEdges.push_back(680);
|
47 |
EBinEdges.push_back(800);
|
48 |
EBinEdges.push_back(1000);
|
49 |
EBinEdges.push_back(1300);
|
50 |
EBinEdges.push_back(1700);
|
51 |
EBinEdges.push_back(2200);
|
52 |
|
53 |
PtBinEdges.push_back(0);
|
54 |
PtBinEdges.push_back(20);
|
55 |
PtBinEdges.push_back(30);
|
56 |
PtBinEdges.push_back(50);
|
57 |
PtBinEdges.push_back(80);
|
58 |
PtBinEdges.push_back(120);
|
59 |
PtBinEdges.push_back(170);
|
60 |
PtBinEdges.push_back(230);
|
61 |
PtBinEdges.push_back(300);
|
62 |
PtBinEdges.push_back(380);
|
63 |
PtBinEdges.push_back(470);
|
64 |
PtBinEdges.push_back(570);
|
65 |
PtBinEdges.push_back(680);
|
66 |
PtBinEdges.push_back(800);
|
67 |
PtBinEdges.push_back(1000);
|
68 |
PtBinEdges.push_back(1300);
|
69 |
PtBinEdges.push_back(1700);
|
70 |
PtBinEdges.push_back(2200);
|
71 |
|
72 |
EtaBinEdges.push_back(0.0);
|
73 |
EtaBinEdges.push_back(1.4);
|
74 |
EtaBinEdges.push_back(2.6);
|
75 |
EtaBinEdges.push_back(3.2);
|
76 |
EtaBinEdges.push_back(4.5);
|
77 |
EtaBinEdges.push_back(5.0);
|
78 |
|
79 |
h_JetResPt_Pt.resize(EtaBinEdges.size() - 1);
|
80 |
for (vector<vector<TH1F*> >::iterator it = h_JetResPt_Pt.begin(); it != h_JetResPt_Pt.end(); ++it) {
|
81 |
it->resize(PtBinEdges.size() - 1);
|
82 |
}
|
83 |
h_JetResEta_Pt.resize(EtaBinEdges.size() - 1);
|
84 |
for (vector<vector<TH1F*> >::iterator it = h_JetResEta_Pt.begin(); it != h_JetResEta_Pt.end(); ++it) {
|
85 |
it->resize(PtBinEdges.size() - 1);
|
86 |
}
|
87 |
h_JetResPhi_Pt.resize(EtaBinEdges.size() - 1);
|
88 |
for (vector<vector<TH1F*> >::iterator it = h_JetResPhi_Pt.begin(); it != h_JetResPhi_Pt.end(); ++it) {
|
89 |
it->resize(PtBinEdges.size() - 1);
|
90 |
}
|
91 |
|
92 |
h_JetResPt_E.resize(EtaBinEdges.size() - 1);
|
93 |
for (vector<vector<TH1F*> >::iterator it = h_JetResPt_E.begin(); it != h_JetResPt_E.end(); ++it) {
|
94 |
it->resize(EBinEdges.size() - 1);
|
95 |
}
|
96 |
h_JetResEta_E.resize(EtaBinEdges.size() - 1);
|
97 |
for (vector<vector<TH1F*> >::iterator it = h_JetResEta_E.begin(); it != h_JetResEta_E.end(); ++it) {
|
98 |
it->resize(EBinEdges.size() - 1);
|
99 |
}
|
100 |
h_JetResPhi_E.resize(EtaBinEdges.size() - 1);
|
101 |
for (vector<vector<TH1F*> >::iterator it = h_JetResPhi_E.begin(); it != h_JetResPhi_E.end(); ++it) {
|
102 |
it->resize(EBinEdges.size() - 1);
|
103 |
}
|
104 |
|
105 |
for (unsigned int i_pt = 0; i_pt < PtBinEdges.size() - 1; ++i_pt) {
|
106 |
for (unsigned int i_eta = 0; i_eta < EtaBinEdges.size() - 1; ++i_eta) {
|
107 |
char hname[100];
|
108 |
sprintf(hname, "hResponsePt_Pt%i_Eta%i", i_pt, i_eta);
|
109 |
h_JetResPt_Pt.at(i_eta).at(i_pt) = new TH1F(hname, hname, 300, 0., 3.);
|
110 |
h_JetResPt_Pt.at(i_eta).at(i_pt)->Sumw2();
|
111 |
sprintf(hname, "hResponseEta_Pt%i_Eta%i", i_pt, i_eta);
|
112 |
h_JetResEta_Pt.at(i_eta).at(i_pt) = new TH1F(hname, hname, 200, -0.5, 0.5);
|
113 |
h_JetResEta_Pt.at(i_eta).at(i_pt)->Sumw2();
|
114 |
sprintf(hname, "hResponsePhi_Pt%i_Eta%i", i_pt, i_eta);
|
115 |
h_JetResPhi_Pt.at(i_eta).at(i_pt) = new TH1F(hname, hname, 200, -0.5, 0.5);
|
116 |
h_JetResPhi_Pt.at(i_eta).at(i_pt)->Sumw2();
|
117 |
}
|
118 |
}
|
119 |
|
120 |
for (unsigned int i_e = 0; i_e < EBinEdges.size() - 1; ++i_e) {
|
121 |
for (unsigned int i_eta = 0; i_eta < EtaBinEdges.size() - 1; ++i_eta) {
|
122 |
char hname[100];
|
123 |
sprintf(hname, "hResponsePt_E%i_Eta%i", i_e, i_eta);
|
124 |
h_JetResPt_E.at(i_eta).at(i_e) = new TH1F(hname, hname, 300, 0., 3.);
|
125 |
h_JetResPt_E.at(i_eta).at(i_e)->Sumw2();
|
126 |
sprintf(hname, "hResponseEta_E%i_Eta%i", i_e, i_eta);
|
127 |
h_JetResEta_E.at(i_eta).at(i_e) = new TH1F(hname, hname, 200, -0.5, 0.5);
|
128 |
h_JetResEta_E.at(i_eta).at(i_e)->Sumw2();
|
129 |
sprintf(hname, "hResponsePhi_E%i_Eta%i", i_e, i_eta);
|
130 |
h_JetResPhi_E.at(i_eta).at(i_e) = new TH1F(hname, hname, 200, -0.5, 0.5);
|
131 |
h_JetResPhi_E.at(i_eta).at(i_e)->Sumw2();
|
132 |
}
|
133 |
}
|
134 |
|
135 |
}
|
136 |
|
137 |
void JetResolutions::SlaveBegin(TTree * /*tree*/) {
|
138 |
// The SlaveBegin() function is called after the Begin() function.
|
139 |
// When running with PROOF SlaveBegin() is called on each slave server.
|
140 |
// The tree argument is deprecated (on PROOF 0 is passed).
|
141 |
|
142 |
TString option = GetOption();
|
143 |
|
144 |
}
|
145 |
|
146 |
Bool_t JetResolutions::Process(Long64_t entry) {
|
147 |
// The Process() function is called for each entry in the tree (or possibly
|
148 |
// keyed object in the case of PROOF) to be processed. The entry argument
|
149 |
// specifies which entry in the currently loaded tree is to be processed.
|
150 |
// It can be passed to either JetResolutions::GetEntry() or TBranch::GetEntry()
|
151 |
// to read either all or the required parts of the data. When processing
|
152 |
// keyed objects with PROOF, the object is already loaded and is available
|
153 |
// via the fObject pointer.
|
154 |
//
|
155 |
// This function should contain the "body" of the analysis. It can contain
|
156 |
// simple or elaborate selection criteria, run algorithms on the data
|
157 |
// of the event and typically fill histograms.
|
158 |
//
|
159 |
// The processing can be stopped by calling Abort().
|
160 |
//
|
161 |
// Use fStatus to set the return value of TTree::Process().
|
162 |
//
|
163 |
// The return value is currently not used.
|
164 |
|
165 |
bool SumJets = false;
|
166 |
|
167 |
GetEntry(entry);
|
168 |
|
169 |
if (!SumJets) {
|
170 |
//// Exact and unique 1:1 matching
|
171 |
//// Loop over GenJets (up to 4 leading GenJets)
|
172 |
for (int i = 0; i < NobjGenJet && i < 50; ++i) {
|
173 |
double pt = GenJetColPt[i];
|
174 |
double eta = GenJetColEta[i];
|
175 |
double phi = GenJetColPhi[i];
|
176 |
double e = GenJetColE[i];
|
177 |
if (pt > 10. && fabs(eta) < 5.) {
|
178 |
TLorentzVector gj;
|
179 |
gj.SetPtEtaPhiE(pt, eta, phi, e);
|
180 |
bool badJet = false;
|
181 |
//// Check that no further GenJet is within DeltaR = 1
|
182 |
for (int j = 0; j < NobjGenJet && j < 50; ++j) {
|
183 |
if (i == j)
|
184 |
continue;
|
185 |
double pttmp = GenJetColPt[j];
|
186 |
double etatmp = GenJetColEta[j];
|
187 |
if ((pttmp > 0.05 || pt && pttmp > 20.) && fabs(etatmp) < 5.) {
|
188 |
double phitmp = GenJetColPhi[j];
|
189 |
double etmp = GenJetColE[j];
|
190 |
TLorentzVector gjtmp;
|
191 |
gjtmp.SetPtEtaPhiE(pttmp, etatmp, phitmp, etmp);
|
192 |
double dR = gj.DeltaR(gjtmp);
|
193 |
if (dR < 1.0) {
|
194 |
badJet = true;
|
195 |
break;
|
196 |
}
|
197 |
}
|
198 |
}
|
199 |
//// Matching and check that no further CaloJet is within DeltaR = 1
|
200 |
double dRmin = 999;
|
201 |
double dRnext = 999;
|
202 |
double pt_min = 0;
|
203 |
double eta_min = 0;
|
204 |
double phi_min = 0;
|
205 |
double dPhi_min = 0;
|
206 |
double e_min = 0;
|
207 |
for (int j = 0; j < NobjJet && j < 50; ++j) {
|
208 |
double pttmp = JetCorrL2L3[j] * JetPt[j];
|
209 |
double etatmp = JetEta[j];
|
210 |
double phitmp = JetPhi[j];
|
211 |
double etmp = JetCorrL2L3[j] * JetE[j];
|
212 |
if (pttmp > 0.05 || pt && pttmp > 20.) {
|
213 |
TLorentzVector jtmp;
|
214 |
jtmp.SetPtEtaPhiE(pttmp, etatmp, phitmp, etmp);
|
215 |
double dR = gj.DeltaR(jtmp);
|
216 |
if (dR > dRmin && dR < dRnext) {
|
217 |
dRnext = dR;
|
218 |
}
|
219 |
if (dR < dRmin) {
|
220 |
dRnext = dRmin;
|
221 |
pt_min = pttmp;
|
222 |
eta_min = etatmp;
|
223 |
phi_min = phitmp;
|
224 |
e_min = etmp;
|
225 |
dRmin = dR;
|
226 |
dPhi_min = gj.DeltaPhi(jtmp);
|
227 |
}
|
228 |
}
|
229 |
}
|
230 |
if (!badJet && dRmin < 0.2 && dRnext > 1.0) {
|
231 |
//if (!badJet && dRmin < 0.2) {
|
232 |
h_JetResPt_Pt.at(EtaBin(eta)).at(PtBin(pt))->Fill(pt_min / pt);
|
233 |
h_JetResPt_E.at(EtaBin(eta)).at(EBin(e))->Fill(pt_min / pt);
|
234 |
h_JetResEta_Pt.at(EtaBin(eta)).at(PtBin(pt))->Fill(eta_min - eta);
|
235 |
h_JetResEta_E.at(EtaBin(eta)).at(EBin(e))->Fill(eta_min - eta);
|
236 |
h_JetResPhi_Pt.at(EtaBin(eta)).at(PtBin(pt))->Fill(dPhi_min);
|
237 |
h_JetResPhi_E.at(EtaBin(eta)).at(EBin(e))->Fill(dPhi_min);
|
238 |
}
|
239 |
}
|
240 |
}
|
241 |
} else {
|
242 |
//// Possible matching of one GenJet to many CaloJets
|
243 |
//// Loop over GenJets (up to 4 leading GenJets)
|
244 |
for (int i = 0; i < NobjGenJet && i < 50; ++i) {
|
245 |
double pt = GenJetColPt[i];
|
246 |
double eta = GenJetColEta[i];
|
247 |
double phi = GenJetColPhi[i];
|
248 |
double e = GenJetColE[i];
|
249 |
|
250 |
TLorentzVector MonsterJet(0., 0., 0., 0.);
|
251 |
if (pt > 10. && fabs(eta) < 5.) {
|
252 |
TLorentzVector gj;
|
253 |
gj.SetPtEtaPhiE(pt, eta, phi, e);
|
254 |
bool badJet = false;
|
255 |
//// Check that no further GenJet is within DeltaR = 1
|
256 |
for (int j = 0; j < NobjGenJet; ++j) {
|
257 |
if (i == j)
|
258 |
continue;
|
259 |
double pttmp = GenJetColPt[j];
|
260 |
double etatmp = GenJetColEta[j];
|
261 |
if ((pttmp > 0.05 || pt && pttmp > 20.) && fabs(etatmp) < 5.) {
|
262 |
double phitmp = GenJetColPhi[j];
|
263 |
double etmp = GenJetColE[j];
|
264 |
TLorentzVector gjtmp;
|
265 |
gjtmp.SetPtEtaPhiE(pttmp, etatmp, phitmp, etmp);
|
266 |
double dR = gj.DeltaR(gjtmp);
|
267 |
if (dR < 1.5) {
|
268 |
badJet = true;
|
269 |
break;
|
270 |
}
|
271 |
}
|
272 |
}
|
273 |
//// Matching and check that no further CaloJet is within DeltaR = 1
|
274 |
double dRmin = 999;
|
275 |
double dRnext = 999;
|
276 |
double pt_min = 0;
|
277 |
double eta_min = 0;
|
278 |
double phi_min = 0;
|
279 |
double e_min = 0;
|
280 |
for (int j = 0; j < NobjJet; ++j) {
|
281 |
double pttmp = JetCorrL2L3[j] * JetPt[j];
|
282 |
double etatmp = JetEta[j];
|
283 |
double phitmp = JetPhi[j];
|
284 |
double etmp = JetCorrL2L3[j] * JetE[j];
|
285 |
TLorentzVector jtmp;
|
286 |
jtmp.SetPtEtaPhiE(pttmp, etatmp, phitmp, etmp);
|
287 |
double dR = gj.DeltaR(jtmp);
|
288 |
if (dR > 0.7 && dR < 1.5) {
|
289 |
badJet = true;
|
290 |
}
|
291 |
if (dR < 0.7) {
|
292 |
MonsterJet += jtmp;
|
293 |
}
|
294 |
}
|
295 |
if (!badJet) {
|
296 |
h_JetResPt_Pt.at(EtaBin(eta)).at(PtBin(pt))->Fill(MonsterJet.Pt() / pt);
|
297 |
h_JetResPt_E.at(EtaBin(eta)).at(EBin(e))->Fill(MonsterJet.Pt() / pt);
|
298 |
h_JetResEta_Pt.at(EtaBin(eta)).at(PtBin(pt))->Fill(MonsterJet.Eta() - eta);
|
299 |
h_JetResEta_E.at(EtaBin(eta)).at(EBin(e))->Fill(MonsterJet.Eta() - eta);
|
300 |
h_JetResPhi_Pt.at(EtaBin(eta)).at(PtBin(pt))->Fill(gj.DeltaPhi(MonsterJet));
|
301 |
h_JetResPhi_E.at(EtaBin(eta)).at(EBin(e))->Fill(gj.DeltaPhi(MonsterJet));
|
302 |
}
|
303 |
}
|
304 |
}
|
305 |
}
|
306 |
|
307 |
return kTRUE;
|
308 |
}
|
309 |
|
310 |
void JetResolutions::SlaveTerminate() {
|
311 |
// The SlaveTerminate() function is called after all entries or objects
|
312 |
// have been processed. When running with PROOF SlaveTerminate() is called
|
313 |
// on each slave server.
|
314 |
|
315 |
}
|
316 |
|
317 |
void JetResolutions::Terminate() {
|
318 |
// The Terminate() function is the last function to be called during
|
319 |
// a query. It always runs on the client, it can be used to present
|
320 |
// the results graphically or save the results to file.
|
321 |
TFile hfile("output.root", "RECREATE", "Jet response in pT and eta bins");
|
322 |
|
323 |
// Save all objects in this file
|
324 |
for (unsigned int i_pt = 0; i_pt < PtBinEdges.size() - 1; ++i_pt) {
|
325 |
for (unsigned int i_eta = 0; i_eta < EtaBinEdges.size() - 1; ++i_eta) {
|
326 |
|
327 |
char cname[100];
|
328 |
sprintf(cname, "c_Pt%i_Eta%i.eps", i_pt, i_eta);
|
329 |
TCanvas* c = new TCanvas(cname, cname, 800, 800);
|
330 |
c->Divide(2,2);
|
331 |
c->cd(1)->SetLogy();
|
332 |
char fname[100];
|
333 |
sprintf(fname, "f_Pt%i_Eta%i", i_pt, i_eta);
|
334 |
//// Gauss
|
335 |
TF1* fitfunction = new TF1(fname, "gaus(0)", 0.5, 2.);
|
336 |
TF1* fitfunction_eta = new TF1(fname, "gaus(0)", -0.5, 0.5);
|
337 |
TF1* fitfunction_phi = new TF1(fname, "gaus(0)", -0.5, 0.5);
|
338 |
//// Single sided crystal ball
|
339 |
// TF1* fitfunction = new TF1(fname, CrystalBall, 0., 3., 5);
|
340 |
// fitfunction->SetParNames("#alpha", "n", "Mean", "#sigma", "N");
|
341 |
//// Double sided crystal ball
|
342 |
// TF1* fitfunction = new TF1(fname, DSCrystalBall, 0., 3., 7);
|
343 |
// fitfunction->SetParNames("#alpha_l", "n_l", "#alpha_r", "n_r", "Mean", "#sigma", "N");
|
344 |
double integral = h_JetResPt_Pt.at(i_eta).at(i_pt)->Integral();
|
345 |
if (integral > 0) {
|
346 |
h_JetResPt_Pt.at(i_eta).at(i_pt)->Scale(1. / integral);
|
347 |
h_JetResEta_Pt.at(i_eta).at(i_pt)->Scale(1. / integral);
|
348 |
h_JetResPhi_Pt.at(i_eta).at(i_pt)->Scale(1. / integral);
|
349 |
}
|
350 |
h_JetResPt_Pt.at(i_eta).at(i_pt)->Draw();
|
351 |
c->cd(2)->SetLogy();
|
352 |
h_JetResEta_Pt.at(i_eta).at(i_pt)->Draw();
|
353 |
c->cd(3)->SetLogy();
|
354 |
h_JetResPhi_Pt.at(i_eta).at(i_pt)->Draw();
|
355 |
double ptmean = 0.5 * (PtBinEdges.at(i_pt) + PtBinEdges.at(i_pt + 1));
|
356 |
double etamean = 0.5 * (EtaBinEdges.at(i_eta) + EtaBinEdges.at(i_eta + 1));
|
357 |
double emean = ptmean * TMath::CosH(etamean);
|
358 |
double sigma = TMath::Sqrt(pow(1.2 / sqrt(emean), 2) + pow(0.03, 2));
|
359 |
double norm = 1. / sigma / sqrt(2* TMath ::Pi()) / 100.;
|
360 |
//// For Gauss
|
361 |
fitfunction->SetParameters(norm, 1., sigma);
|
362 |
fitfunction_eta->SetParameters(norm, 0., sigma);
|
363 |
fitfunction_phi->SetParameters(norm, 0., sigma);
|
364 |
//// For single sided CB
|
365 |
// fitfunction->SetParameters(2., 2., 1., sigma, norm);
|
366 |
//// For double sided CB
|
367 |
// fitfunction->SetParameters(3., 3., 3., 3., 1., sigma, norm);
|
368 |
if (integral > 10) {
|
369 |
h_JetResPt_Pt.at(i_eta).at(i_pt)->Fit(fitfunction, "LLRQN");
|
370 |
h_JetResEta_Pt.at(i_eta).at(i_pt)->Fit(fitfunction_eta, "LLRQN");
|
371 |
h_JetResPhi_Pt.at(i_eta).at(i_pt)->Fit(fitfunction_phi, "LLRQN");
|
372 |
}
|
373 |
// for (int i = 1; i <= h_JetRes_Pt.at(i_eta).at(i_pt)->GetNbinsX(); ++i) {
|
374 |
// h_JetRes_Pt.at(i_eta).at(i_pt)->SetBinContent(i, fitfunction->Eval(
|
375 |
// h_JetRes_Pt.at(i_eta).at(i_pt)->GetBinCenter(i)));
|
376 |
// h_JetRes_Pt.at(i_eta).at(i_pt)->SetBinError(i, 0.);
|
377 |
// }
|
378 |
c->cd(1);
|
379 |
fitfunction->Draw("same");
|
380 |
c->cd(2);
|
381 |
fitfunction_eta->Draw("same");
|
382 |
c->cd(3);
|
383 |
fitfunction_phi->Draw("same");
|
384 |
c->SaveAs(cname);
|
385 |
h_JetResPt_Pt.at(i_eta).at(i_pt)->Write();
|
386 |
std::cout << "Pt = " << PtBinEdges.at(i_pt) << " ..." << PtBinEdges.at(i_pt + 1) << ", Eta = "
|
387 |
<< EtaBinEdges.at(i_eta) << "..." << EtaBinEdges.at(i_eta + 1) << ", mean response = "
|
388 |
<< fitfunction->GetParameter(1) << ", sigma pt = " << fitfunction->GetParameter(2)
|
389 |
<< ", delta Eta mean = " << fitfunction_eta->GetParameter(1) << ", delta Eta sigma = "
|
390 |
<< fitfunction_eta->GetParameter(2) << ", delta Phi mean = " << fitfunction_phi->GetParameter(1)
|
391 |
<< ", delta Phi sigma = " << fitfunction_phi->GetParameter(2) << std::endl;
|
392 |
}
|
393 |
}
|
394 |
for (unsigned int i_e = 0; i_e < EBinEdges.size() - 1; ++i_e) {
|
395 |
for (unsigned int i_eta = 0; i_eta < EtaBinEdges.size() - 1; ++i_eta) {
|
396 |
char cname[100];
|
397 |
sprintf(cname, "c_E%i_Eta%i.eps", i_e, i_eta);
|
398 |
TCanvas* c = new TCanvas(cname, cname, 800, 800);
|
399 |
c->Divide(2,2);
|
400 |
c->cd(1)->SetLogy();
|
401 |
char fname[100];
|
402 |
sprintf(fname, "f_E%i_Eta%i", i_e, i_eta);
|
403 |
//// Gauss
|
404 |
TF1* fitfunction = new TF1(fname, "gaus(0)", 0.5, 2.);
|
405 |
TF1* fitfunction_eta = new TF1(fname, "gaus(0)", -0.5, 0.5);
|
406 |
TF1* fitfunction_phi = new TF1(fname, "gaus(0)", -0.5, 0.5);
|
407 |
//// Single sided crystal ball
|
408 |
// TF1* fitfunction = new TF1(fname, CrystalBall, 0., 3., 5);
|
409 |
// fitfunction->SetParNames("#alpha", "n", "Mean", "#sigma", "N");
|
410 |
//// Double sided crystal ball
|
411 |
// TF1* fitfunction = new TF1(fname, DSCrystalBall, 0., 3., 7);
|
412 |
// fitfunction->SetParNames("#alpha_l", "n_l", "#alpha_r", "n_r", "Mean", "#sigma", "N");
|
413 |
double integral = h_JetResPt_E.at(i_eta).at(i_e)->Integral();
|
414 |
if (integral > 0) {
|
415 |
h_JetResPt_E.at(i_eta).at(i_e)->Scale(1. / integral);
|
416 |
}
|
417 |
h_JetResPt_E.at(i_eta).at(i_e)->Draw();
|
418 |
c->cd(2)->SetLogy();
|
419 |
h_JetResEta_E.at(i_eta).at(i_e)->Draw();
|
420 |
c->cd(3)->SetLogy();
|
421 |
h_JetResPhi_E.at(i_eta).at(i_e)->Draw();
|
422 |
double emean = (EBinEdges.at(i_e) + EBinEdges.at(i_e + 1)) / 2.;
|
423 |
double sigma = TMath::Sqrt(pow(1.2 / sqrt(emean), 2) + pow(0.03, 2));
|
424 |
double norm = 1. / sigma / sqrt(2* TMath ::Pi()) / 100.;
|
425 |
//// For Gauss
|
426 |
fitfunction->SetParameters(norm, 1., sigma);
|
427 |
fitfunction_eta->SetParameters(norm, 0., sigma);
|
428 |
fitfunction_phi->SetParameters(norm, 0., sigma);
|
429 |
//// For single sided CB
|
430 |
// fitfunction->SetParameters(2., 2., 1., sigma, norm);
|
431 |
//// For double sided CB
|
432 |
// fitfunction->SetParameters(3., 3., 3., 3., 1., sigma, norm);
|
433 |
if (integral > 10) {
|
434 |
h_JetResPt_E.at(i_eta).at(i_e)->Fit(fitfunction, "LLRQN");
|
435 |
h_JetResEta_E.at(i_eta).at(i_e)->Fit(fitfunction_eta, "LLRQN");
|
436 |
h_JetResPhi_E.at(i_eta).at(i_e)->Fit(fitfunction_phi, "LLRQN");
|
437 |
}
|
438 |
// for (int i = 1; i <= h_JetRes_E.at(i_eta).at(i_e)->GetNbinsX(); ++i) {
|
439 |
// h_JetRes_E.at(i_eta).at(i_e)->SetBinContent(i, fitfunction->Eval(
|
440 |
// h_JetRes_E.at(i_eta).at(i_e)->GetBinCenter(i)));
|
441 |
// h_JetRes_E.at(i_eta).at(i_e)->SetBinError(i, 0.);
|
442 |
// }
|
443 |
c->cd(1);
|
444 |
fitfunction->Draw("same");
|
445 |
c->cd(2);
|
446 |
fitfunction_eta->Draw("same");
|
447 |
c->cd(3);
|
448 |
fitfunction_phi->Draw("same");
|
449 |
c->SaveAs(cname);
|
450 |
h_JetResPt_E.at(i_eta).at(i_e)->Write();
|
451 |
std::cout << "E = " << EBinEdges.at(i_e) << " ..." << EBinEdges.at(i_e + 1) << ", Eta = " << EtaBinEdges.at(
|
452 |
i_eta) << "..." << EtaBinEdges.at(i_eta + 1) << ", mean response = " << fitfunction->GetParameter(1)
|
453 |
<< ", sigma pt = " << fitfunction->GetParameter(2) << ", delta Eta mean = "
|
454 |
<< fitfunction_eta->GetParameter(1) << ", delta Eta sigma = " << fitfunction_eta->GetParameter(2)
|
455 |
<< ", delta Phi mean = " << fitfunction_phi->GetParameter(1) << ", delta Phi sigma = "
|
456 |
<< fitfunction_phi->GetParameter(2) << std::endl;
|
457 |
}
|
458 |
}
|
459 |
hfile.WriteObject(&EBinEdges, "EBinEdges");
|
460 |
hfile.WriteObject(&PtBinEdges, "PtBinEdges");
|
461 |
hfile.WriteObject(&EtaBinEdges, "EtaBinEdges");
|
462 |
|
463 |
// Close the file.
|
464 |
hfile.Close();
|
465 |
|
466 |
}
|
467 |
|
468 |
int JetResolutions::EBin(const double& e) {
|
469 |
int i_e = -1;
|
470 |
for (vector<double>::const_iterator it = EBinEdges.begin(); it != EBinEdges.end(); ++it) {
|
471 |
if ((*it) > e)
|
472 |
break;
|
473 |
++i_e;
|
474 |
}
|
475 |
if (i_e < 0)
|
476 |
i_e = 0;
|
477 |
if (i_e > EBinEdges.size() - 2)
|
478 |
i_e = EBinEdges.size() - 2;
|
479 |
|
480 |
return i_e;
|
481 |
}
|
482 |
|
483 |
int JetResolutions::PtBin(const double& pt) {
|
484 |
int i_pt = -1;
|
485 |
for (vector<double>::const_iterator it = PtBinEdges.begin(); it != PtBinEdges.end(); ++it) {
|
486 |
if ((*it) > pt)
|
487 |
break;
|
488 |
++i_pt;
|
489 |
}
|
490 |
if (i_pt < 0)
|
491 |
i_pt = 0;
|
492 |
if (i_pt > PtBinEdges.size() - 2)
|
493 |
i_pt = PtBinEdges.size() - 2;
|
494 |
|
495 |
return i_pt;
|
496 |
}
|
497 |
|
498 |
int JetResolutions::EtaBin(const double& eta) {
|
499 |
int i_eta = -1;
|
500 |
for (vector<double>::const_iterator it = EtaBinEdges.begin(); it != EtaBinEdges.end(); ++it) {
|
501 |
if ((*it) > fabs(eta))
|
502 |
break;
|
503 |
++i_eta;
|
504 |
}
|
505 |
if (i_eta < 0)
|
506 |
i_eta = 0;
|
507 |
if (i_eta > EtaBinEdges.size() - 2)
|
508 |
i_eta = EtaBinEdges.size() - 2;
|
509 |
return i_eta;
|
510 |
}
|
511 |
|
512 |
int main() {
|
513 |
TChain* ch = new TChain("DiJetTree");
|
514 |
JetResolutions *S = new JetResolutions();
|
515 |
for (unsigned int i = 1; i <= 92; ++i) {
|
516 |
char fname[100];
|
517 |
sprintf(fname, "Summer09QCDFlat-MC_31X_V9_7TeV-v1_E/Summer09QCDFlat_Pt15to3000MC_31X_V9_7TeV-v1_%i.root", i);
|
518 |
ch->Add(fname);
|
519 |
}
|
520 |
ch->Process(S);
|
521 |
return 0;
|
522 |
}
|