1 |
#if !defined(__CINT__) && !defined(__MAKECINT__)
|
2 |
#include "DataFormats/FWLite/interface/Handle.h"
|
3 |
#include "DataFormats/FWLite/interface/Event.h"
|
4 |
#include "DataFormats/FWLite/interface/ChainEvent.h"
|
5 |
//Headers for the data items
|
6 |
#include "DataFormats/VertexReco/interface/Vertex.h"
|
7 |
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
|
8 |
|
9 |
#endif
|
10 |
|
11 |
#include "TROOT.h"
|
12 |
#include "TH2F.h"
|
13 |
#include "TFile.h"
|
14 |
#include "TF1.h"
|
15 |
#include "TCanvas.h"
|
16 |
#include "TMath.h"
|
17 |
#include "TBox.h"
|
18 |
#include <vector>
|
19 |
#include <map>
|
20 |
#include <string>
|
21 |
|
22 |
using namespace std;
|
23 |
typedef map<unsigned int, unsigned int> RunNumberMap;
|
24 |
|
25 |
//
|
26 |
// iterative fit of +- 2 sigma Gaussian core on a 1d histogram
|
27 |
//
|
28 |
bool fitHistogram (TH1* slice, TF1& fg, double* pars, double* epars, bool dbg = false) {
|
29 |
// initialization of parameters
|
30 |
pars[0] = slice->GetMaximum();
|
31 |
pars[1] = slice->GetMean();
|
32 |
pars[2] = slice->GetRMS();
|
33 |
epars[0] = 0;
|
34 |
epars[1] = 0;
|
35 |
epars[2] = 0;
|
36 |
// some iterations
|
37 |
bool converged(false);
|
38 |
for ( unsigned int it=0; it<10; ++it) {
|
39 |
const double nSigRange = 2.;
|
40 |
fg.SetParameters(pars);
|
41 |
double oldMean = pars[1];
|
42 |
double oldWidth = pars[2];
|
43 |
if ( it>0 ) {
|
44 |
fg.SetRange(pars[1]-nSigRange*pars[2],
|
45 |
pars[1]+nSigRange*pars[2]);
|
46 |
// if ( dbg ) {
|
47 |
// cout << "Setting range to " << pars[1]-nSigRange*pars[2]
|
48 |
// << " " << pars[1]+nSigRange*pars[2] << endl;
|
49 |
// }
|
50 |
}
|
51 |
else {
|
52 |
fg.SetRange(slice->GetXaxis()->GetXmin(),
|
53 |
slice->GetXaxis()->GetXmax());
|
54 |
// if ( dbg ) {
|
55 |
// cout << "Setting range to " << slice->GetXaxis()->GetXmin()
|
56 |
// << " " << slice->GetXaxis()->GetXmax() << endl;
|
57 |
// }
|
58 |
}
|
59 |
if ( dbg )
|
60 |
slice->Fit(&fg,"QLR+");
|
61 |
else
|
62 |
slice->Fit(&fg,"QLR");
|
63 |
fg.GetParameters(pars);
|
64 |
epars[0] = fg.GetParError(0);
|
65 |
epars[1] = fg.GetParError(1);
|
66 |
epars[2] = fg.GetParError(2);
|
67 |
if ( converged ) {
|
68 |
cout << "terminating at iteration " << it << " "
|
69 |
<< pars[1] << " " << oldMean << " " << epars[1] << " "
|
70 |
<< pars[2] << " " << oldWidth << " " << epars[2] << endl;
|
71 |
break;
|
72 |
}
|
73 |
else if ( dbg ) {
|
74 |
cout << " at iteration " << it << " "
|
75 |
<< pars[1] << " " << oldMean << " " << epars[1] << " "
|
76 |
<< pars[2] << " " << oldWidth << " " << epars[2] << endl;
|
77 |
}
|
78 |
// convergence criterion
|
79 |
if ( it>2 &&
|
80 |
fabs((pars[1]-oldMean))<max(fabs(0.01*epars[1]),0.0010) &&
|
81 |
fabs((pars[2]-oldWidth))<max(fabs(0.01*epars[2]),0.0010) ) {
|
82 |
converged = true;
|
83 |
}
|
84 |
}
|
85 |
|
86 |
return converged;
|
87 |
}
|
88 |
|
89 |
//
|
90 |
// fit Gaussian to slices of the 2d "run" histogram
|
91 |
//
|
92 |
void fitRunHistogram (TH2* histo, TH1*& hmean, TH1*& hwidth) {
|
93 |
//
|
94 |
// create histograms for mean and width vs. run
|
95 |
//
|
96 |
TDirectory* saveDir = gDirectory;
|
97 |
gROOT->cd();
|
98 |
int nbin = histo->GetNbinsX();
|
99 |
string title = histo->GetName();
|
100 |
TAxis* xhisto = histo->GetXaxis();
|
101 |
hmean = new TH1F((title+"Mean").c_str(),title.c_str(),nbin,-0.5,nbin-0.5);
|
102 |
TAxis* xhmean = hmean->GetXaxis();
|
103 |
hwidth = new TH1F((title+"Width").c_str(),title.c_str(),nbin,-0.5,nbin-0.5);
|
104 |
TAxis* xhwidth = hwidth->GetXaxis();
|
105 |
//
|
106 |
// fit function
|
107 |
//
|
108 |
TF1 fg("fg","gaus");
|
109 |
//
|
110 |
// loop over runs
|
111 |
//
|
112 |
bool converged;
|
113 |
double pars[3];
|
114 |
double epars[3];
|
115 |
for ( unsigned int ib=1; ib<=nbin; ++ib ) {
|
116 |
//
|
117 |
// initialize mean and width to 0 error
|
118 |
//
|
119 |
hmean->SetBinError(ib,0.);
|
120 |
hwidth->SetBinError(ib,0.);
|
121 |
//
|
122 |
// get slice
|
123 |
//
|
124 |
TH1* slice =
|
125 |
histo->ProjectionY((title+xhisto->GetBinLabel(ib)).c_str(),ib,ib);
|
126 |
// drop histograms with too few entries
|
127 |
if ( slice->GetEntries()<100 ) {
|
128 |
delete slice;
|
129 |
continue;
|
130 |
}
|
131 |
//
|
132 |
// fit slice
|
133 |
//
|
134 |
converged = fitHistogram(slice,fg,pars,epars);
|
135 |
if ( converged ) {
|
136 |
// copy parameters and errors
|
137 |
hmean->SetBinContent(ib,pars[1]);
|
138 |
hmean->SetBinError(ib,epars[1]);
|
139 |
hwidth->SetBinContent(ib,pars[2]);
|
140 |
hwidth->SetBinError(ib,epars[2]);
|
141 |
}
|
142 |
else {
|
143 |
cout << "No convergence for " << histo->GetName() << " / run "
|
144 |
<< histo->GetXaxis()->GetBinLabel(ib) << endl;
|
145 |
// refit with debug flag
|
146 |
// fitHistogram(slice,fg,pars,epars,true);
|
147 |
hmean->SetBinContent(ib,pars[1]);
|
148 |
hmean->SetBinError(ib,epars[1]);
|
149 |
hwidth->SetBinContent(ib,pars[2]);
|
150 |
hwidth->SetBinError(ib,epars[2]);
|
151 |
}
|
152 |
// keep failing fits for debugging
|
153 |
// if ( converged ) delete slice;
|
154 |
xhmean->SetBinLabel(ib,xhisto->GetBinLabel(ib));
|
155 |
xhwidth->SetBinLabel(ib,xhisto->GetBinLabel(ib));
|
156 |
}
|
157 |
|
158 |
saveDir->cd();
|
159 |
}
|
160 |
|
161 |
//
|
162 |
// set x-axis labels to run number
|
163 |
//
|
164 |
void setRunNumberLabels (TH1* histo, RunNumberMap& runNumbers) {
|
165 |
char label[32];
|
166 |
TAxis* axis = histo->GetXaxis();
|
167 |
for ( RunNumberMap::const_iterator i=runNumbers.begin();
|
168 |
i!=runNumbers.end(); ++i ) {
|
169 |
sprintf(label,"%d",(*i).first);
|
170 |
axis->SetBinLabel((*i).second+1,label);
|
171 |
}
|
172 |
}
|
173 |
|
174 |
|
175 |
//
|
176 |
// plot PV distribution for a single run with indication of BS position
|
177 |
//
|
178 |
void plotRun (char coord, int run)
|
179 |
{
|
180 |
char title[64];
|
181 |
|
182 |
sprintf(title,"pv%cVsRun%d",coord,run);
|
183 |
TH1* hpv = (TH1*)gROOT->FindObject(title);
|
184 |
if ( hpv==0 ) return;
|
185 |
|
186 |
sprintf(title,"pv%cRun%d",coord,run);
|
187 |
TCanvas* c = new TCanvas(title);
|
188 |
|
189 |
hpv->Draw();
|
190 |
|
191 |
sprintf(title,"bs%cVsRunMean",coord);
|
192 |
TH1* hbs = (TH1*)gROOT->FindObject(title);
|
193 |
sprintf(title,"%d",run);
|
194 |
double vbs(0.);
|
195 |
double ebs(0.);
|
196 |
for ( unsigned int ib=1; ib<=hbs->GetNbinsX(); ++ib ) {
|
197 |
const char* label = hbs->GetXaxis()->GetBinLabel(ib);
|
198 |
if ( strcmp(title,label)==0 ) {
|
199 |
vbs = hbs->GetBinContent(ib);
|
200 |
ebs = hbs->GetBinError(ib);
|
201 |
break;
|
202 |
}
|
203 |
}
|
204 |
TBox* box = new TBox(vbs-ebs,0,vbs+ebs,hpv->GetMaximum());
|
205 |
box->SetFillStyle(1001);
|
206 |
box->SetFillColor(7);
|
207 |
box->Draw();
|
208 |
|
209 |
hpv->Draw("same");
|
210 |
|
211 |
sprintf(title,"pv%cRun%d.png",coord,run);
|
212 |
c->SaveAs(title);
|
213 |
delete c;
|
214 |
}
|
215 |
//
|
216 |
// summary plots (PV + BS) vs run
|
217 |
//
|
218 |
void plotPvVsBs () {
|
219 |
char coord[] = { "xyz" };
|
220 |
for ( int ic=0; ic<3; ++ic ) {
|
221 |
string canvasTitle("pvVsBs");
|
222 |
canvasTitle += coord[ic];
|
223 |
TCanvas* c = new TCanvas(canvasTitle.c_str(),canvasTitle.c_str(),
|
224 |
700,500);
|
225 |
string bsTitle("bs.VsRunMean");
|
226 |
bsTitle.replace(2,1,1,coord[ic]);
|
227 |
TH1* bsVsRun = (TH1*)gROOT->FindObject(bsTitle.c_str());
|
228 |
string pvTitle("pv.VsRunMean");
|
229 |
pvTitle.replace(2,1,1,coord[ic]);
|
230 |
TH1* pvVsRun = (TH1*)gROOT->FindObject(pvTitle.c_str());
|
231 |
pvVsRun->SetMarkerStyle(20);
|
232 |
pvVsRun->Draw();
|
233 |
bsVsRun->SetFillStyle(1001);
|
234 |
bsVsRun->SetFillColor(7);
|
235 |
bsVsRun->Draw("E2 same");
|
236 |
pvVsRun->Draw("same");
|
237 |
}
|
238 |
}
|
239 |
|
240 |
//
|
241 |
// main entry point
|
242 |
//
|
243 |
void pvHistory (const vector<string>& fileNames) {
|
244 |
|
245 |
fwlite::ChainEvent ev(fileNames);
|
246 |
|
247 |
fwlite::Handle< reco::BeamSpot > bsHandle;
|
248 |
fwlite::Handle< vector<reco::Vertex> > pvHandle;
|
249 |
//
|
250 |
// create list of runs with > 100 PVs
|
251 |
//
|
252 |
RunNumberMap runNumberCounts;
|
253 |
int lastRun(0);
|
254 |
for( ev.toBegin(); !ev.atEnd(); ++ev) {
|
255 |
unsigned int run = ev.id().run();
|
256 |
pvHandle.getByLabel(ev,"offlinePrimaryVertices");
|
257 |
|
258 |
const reco::Vertex& pv = (*pvHandle)[0];
|
259 |
if ( pv.tracksSize()>0 && pv.ndof()>2 &&
|
260 |
(pv.ndof()+3.)/2./pv.tracksSize()>0.5 ) {
|
261 |
++runNumberCounts[run];
|
262 |
}
|
263 |
}
|
264 |
// create list of bin indices
|
265 |
unsigned int nRun(0);
|
266 |
RunNumberMap runNumberIndices;
|
267 |
for ( RunNumberMap::iterator i=runNumberCounts.begin();
|
268 |
i!=runNumberCounts.end(); ++i ) {
|
269 |
if ( (*i).second>100 ) runNumberIndices[(*i).first] = nRun++;
|
270 |
}
|
271 |
cout << "Found " << nRun << " run numbers" << endl;
|
272 |
//
|
273 |
// create histograms
|
274 |
//
|
275 |
TDirectory* saveDir = gDirectory;
|
276 |
gROOT->cd();
|
277 |
TH2* pvxVsRun = new TH2F("pvxVsRun","pvxVsRun",nRun,-0.5,nRun-0.5,80,-0.1,0.3);
|
278 |
TH2* pvyVsRun = new TH2F("pvyVsRun","pvyVsRun",nRun,-0.5,nRun-0.5,80,-0.2,0.2);
|
279 |
TH2* pvzVsRun = new TH2F("pvzVsRun","pvzVsRun",nRun,-0.5,nRun-0.5,80,-20,20);
|
280 |
TH2* bsxVsRun = new TH2F("bsxVsRun","bsxVsRun",nRun,-0.5,nRun-0.5,80,0,0.4);
|
281 |
TH2* bsyVsRun = new TH2F("bsyVsRun","bsyVsRun",nRun,-0.5,nRun-0.5,80,0,0.4);
|
282 |
TH2* bszVsRun = new TH2F("bszVsRun","bszVsRun",nRun,-0.5,nRun-0.5,80,-20,20);
|
283 |
TH1* bsxVsRunMean = new TH1F("bsxVsRunMean","bsxVsRunMean",nRun,-0.5,nRun-0.5);
|
284 |
TH1* bsyVsRunMean = new TH1F("bsyVsRunMean","bsyVsRunMean",nRun,-0.5,nRun-0.5);
|
285 |
TH1* bszVsRunMean = new TH1F("bszVsRunMean","bszVsRunMean",nRun,-0.5,nRun-0.5);
|
286 |
saveDir->cd();
|
287 |
|
288 |
//
|
289 |
// fill histograms
|
290 |
//
|
291 |
lastRun = -1;
|
292 |
int idx(-1);
|
293 |
int nev(0);
|
294 |
for( ev.toBegin(); !ev.atEnd(); ++ev) {
|
295 |
int run = ev.id().run();
|
296 |
// if ( (++nev)%1000 == 0 ) cout << "processing record " << nev << " ( run "
|
297 |
// << ev.id().run() << " " << " evt "
|
298 |
// << ev.id().event() << " )" << endl;
|
299 |
//
|
300 |
// update bin index if run changed
|
301 |
//
|
302 |
if ( run!=lastRun ) {
|
303 |
RunNumberMap::const_iterator irun = runNumberIndices.find(run);
|
304 |
if ( irun == runNumberIndices.end() ) {
|
305 |
idx = -1;
|
306 |
}
|
307 |
else {
|
308 |
idx = (*irun).second;
|
309 |
}
|
310 |
lastRun = run;
|
311 |
// cout << "Run " << run << " idx " << idx << endl;
|
312 |
}
|
313 |
if ( idx<0 ) continue;
|
314 |
|
315 |
//
|
316 |
// retrieve beam spot and PV
|
317 |
//
|
318 |
bsHandle.getByLabel(ev,"offlineBeamSpot");
|
319 |
pvHandle.getByLabel(ev,"offlinePrimaryVertices");
|
320 |
//
|
321 |
// fill PV histograms (after quality cut)
|
322 |
//
|
323 |
const reco::Vertex& pv = (*pvHandle)[0];
|
324 |
if ( pv.tracksSize()>0 && pv.ndof()>2 &&
|
325 |
(pv.ndof()+3.)/2./pv.tracksSize()>0.5 ) {
|
326 |
pvxVsRun->Fill(idx,pv.x());
|
327 |
pvyVsRun->Fill(idx,pv.y());
|
328 |
pvzVsRun->Fill(idx,pv.z());
|
329 |
}
|
330 |
//
|
331 |
// fill BS histograms
|
332 |
//
|
333 |
bsxVsRun->Fill(idx,bsHandle->x0());
|
334 |
bsyVsRun->Fill(idx,bsHandle->y0());
|
335 |
bszVsRun->Fill(idx,bsHandle->z0());
|
336 |
if ( idx>=0 && idx<bsxVsRun->GetNbinsX() ) {
|
337 |
// assume BS is constant / run: fill current value into 1d histogram
|
338 |
bsxVsRunMean->SetBinContent(idx+1,bsHandle->x0());
|
339 |
bsxVsRunMean->SetBinError(idx+1,bsHandle->x0Error());
|
340 |
bsyVsRunMean->SetBinContent(idx+1,bsHandle->y0());
|
341 |
bsyVsRunMean->SetBinError(idx+1,bsHandle->y0Error());
|
342 |
bszVsRunMean->SetBinContent(idx+1,bsHandle->z0());
|
343 |
bszVsRunMean->SetBinError(idx+1,bsHandle->z0Error());
|
344 |
}
|
345 |
|
346 |
}
|
347 |
saveDir = gDirectory;
|
348 |
//
|
349 |
// number of entries / run
|
350 |
//
|
351 |
gROOT->cd();
|
352 |
TH1* npvVsRun = pvzVsRun->ProjectionX("npvVsRun");
|
353 |
npvVsRun->SetTitle("npvVsRun");
|
354 |
TH1* nbsVsRun = bszVsRun->ProjectionX("nbsVsRun");
|
355 |
npvVsRun->SetTitle("nbsVsRun");
|
356 |
saveDir->cd();
|
357 |
//
|
358 |
// set run numbers on x-axis labels
|
359 |
//
|
360 |
setRunNumberLabels(npvVsRun,runNumberIndices);
|
361 |
setRunNumberLabels(pvxVsRun,runNumberIndices);
|
362 |
setRunNumberLabels(pvyVsRun,runNumberIndices);
|
363 |
setRunNumberLabels(pvzVsRun,runNumberIndices);
|
364 |
setRunNumberLabels(nbsVsRun,runNumberIndices);
|
365 |
setRunNumberLabels(bsxVsRun,runNumberIndices);
|
366 |
setRunNumberLabels(bsyVsRun,runNumberIndices);
|
367 |
setRunNumberLabels(bszVsRun,runNumberIndices);
|
368 |
setRunNumberLabels(bsxVsRunMean,runNumberIndices);
|
369 |
setRunNumberLabels(bsyVsRunMean,runNumberIndices);
|
370 |
setRunNumberLabels(bszVsRunMean,runNumberIndices);
|
371 |
//
|
372 |
// fit histograms
|
373 |
//
|
374 |
TH1* pvxMean; TH1* pvxWidth;
|
375 |
fitRunHistogram(pvxVsRun,pvxMean,pvxWidth);
|
376 |
TH1* pvyMean; TH1* pvyWidth;
|
377 |
fitRunHistogram(pvyVsRun,pvyMean,pvyWidth);
|
378 |
TH1* pvzMean; TH1* pvzWidth;
|
379 |
fitRunHistogram(pvzVsRun,pvzMean,pvzWidth);
|
380 |
|
381 |
// cout << pvxMean << endl;
|
382 |
// cout << pvxMean->GetNbinsX() << " " << pvxMean->GetXaxis() << endl;
|
383 |
//
|
384 |
// summary table
|
385 |
//
|
386 |
for ( unsigned int i=1; i<=pvxMean->GetNbinsX(); ++i ) {
|
387 |
if ( pvxMean->GetBinError(i)>1.e-10 ||
|
388 |
pvxMean->GetBinError(i)>1.e-10 ||
|
389 |
pvxMean->GetBinError(i)>1.e-10 ) {
|
390 |
printf(" %-6s : %8.4f +- %6.4f %6.4f +- %6.4f %8.4f +- %6.4f %6.4f +- %6.4f %9.4f +- %7.4f %6.4f +- %7.4f\n",
|
391 |
pvxMean->GetXaxis()->GetBinLabel(i),
|
392 |
pvxMean->GetBinContent(i),pvxMean->GetBinError(i),
|
393 |
pvxWidth->GetBinContent(i),pvxWidth->GetBinError(i),
|
394 |
pvyMean->GetBinContent(i),pvyMean->GetBinError(i),
|
395 |
pvyWidth->GetBinContent(i),pvyWidth->GetBinError(i),
|
396 |
pvzMean->GetBinContent(i),pvzMean->GetBinError(i),
|
397 |
pvzWidth->GetBinContent(i),pvzWidth->GetBinError(i));
|
398 |
}
|
399 |
}
|
400 |
|
401 |
//
|
402 |
// save histograms
|
403 |
//
|
404 |
TFile* fo = new TFile("../skim_BSC0_AND_40OR41_MinBias_Dec19th_v1.root","recreate");
|
405 |
pvxVsRun->Write();
|
406 |
pvyVsRun->Write();
|
407 |
pvzVsRun->Write();
|
408 |
pvxMean->Write();
|
409 |
pvxWidth->Write();
|
410 |
pvyMean->Write();
|
411 |
pvyWidth->Write();
|
412 |
pvzMean->Write();
|
413 |
pvzWidth->Write();
|
414 |
fo->Write();
|
415 |
delete fo;
|
416 |
|
417 |
plotPvVsBs();
|
418 |
//
|
419 |
// plot inidividual histograms / run
|
420 |
//
|
421 |
// for ( RunNumberMap::iterator i=runNumberIndices.begin();
|
422 |
// i!=runNumberIndices.end(); ++i ) {
|
423 |
// plotRun('x',(*i).first);
|
424 |
// plotRun('y',(*i).first);
|
425 |
// plotRun('z',(*i).first);
|
426 |
// }
|
427 |
|
428 |
}
|
429 |
//
|
430 |
// call with single file name
|
431 |
//
|
432 |
void pvHistory (const char* fileName) {
|
433 |
vector<string> fileNames(1,fileName);
|
434 |
pvHistory(fileNames);
|
435 |
}
|