15 #include <TGraphErrors.h>
16 #include <TGraphAsymmErrors.h>
41 TFile*
f = TFile::Open(inFilename.c_str(),
"UPDATE");
43 if (
f == 0 || !
f->IsOpen()) {
49 if (
f->GetSize() < 1000.) {
56 bool dirExists =
false;
58 TIter next_run(
f->GetListOfKeys());
60 while ((key_run =
dynamic_cast<TKey*
>(next_run())) != 0) {
61 TObject* obj_run = key_run->ReadObj();
62 TDirectory* tdir_run =
dynamic_cast<TDirectory*
>(obj_run);
64 std::string tdir_run_name(tdir_run->GetName());
65 if (tdir_run_name.find(
"run") != std::string::npos) {
66 run_dir = tdir_run_name;
68 dirExists =
f->GetDirectory(run_dir +
"/InDetGlobal/PrimaryVertex");
71 f->cd(run_dir +
"/InDetGlobal/PrimaryVertex");
82 dirExists =
f->GetDirectory(
"InDetGlobal/PrimaryVertex");
83 if (dirExists)
f->cd(
"InDetGlobal/PrimaryVertex");
110 TH2F* h_Vrt_pullVsSomething_split(0);
111 TH2F* h_Vrt_err_vs_Something(0);
113 TString xAxisLabel(
"");
115 if (versus ==
"Ntrk") {
119 xAxisLabel =
"Number of fitted tracks";
120 }
else if (versus ==
"SumPt2") {
124 xAxisLabel =
"#sqrt{#sum p_{T}^{2}} [GeV]";
129 if (h_Vrt_pullVsSomething_split == 0 or h_Vrt_err_vs_Something == 0)
return;
131 int n_bins = h_Vrt_pullVsSomething_split->GetNbinsX();
132 std::vector<float> rms_z;
133 std::vector<float> rms_z_er;
134 std::vector<float> sigma_z;
135 std::vector<float> sigma_z_er;
136 std::vector<float> bins_z_nt;
137 std::vector<float> bins_z_nt_er;
147 const Int_t minEntriesForKFactorBin = 1000;
148 for (
int bin_count = 1; bin_count < n_bins + 1; bin_count++) {
155 TH1D* profileZTmp = h_Vrt_pullVsSomething_split->ProjectionY(
"projectionPulls", bin_count, bin_count,
"e");
158 startBin = bin_count;
159 profileZ = (TH1D*) profileZTmp->Clone(
"projectionPulls_Integrated");
162 profileZ->Add(profileZTmp);
166 if ((profileZ->GetEntries() < minEntriesForKFactorBin) && (bin_count < n_bins))
169 Double_t lowEdge = h_Vrt_pullVsSomething_split->GetXaxis()->GetBinLowEdge(startBin);
170 Double_t highEdge = h_Vrt_pullVsSomething_split->GetXaxis()->GetBinLowEdge(bin_count) +
171 h_Vrt_pullVsSomething_split->GetXaxis()->GetBinWidth(bin_count);
172 Double_t binCenter = (lowEdge + highEdge) / 2;
173 Double_t
binWidth = (highEdge - lowEdge) / 2;
177 bins_z_nt.push_back(binCenter);
180 rms_z.push_back(profileZ->GetRMS());
181 rms_z_er.push_back(profileZ->GetRMSError());
184 if (profileZ->GetEntries() > 100.) {
186 sigma_z.push_back(fit_res[0]);
187 sigma_z_er.push_back(fit_res[1]);
189 sigma_z.push_back(0.);
190 sigma_z_er.push_back(0.);
197 TGraphErrors* krms_z_vs_ntrk =
new TGraphErrors(
198 bins_z_nt.size(), &(bins_z_nt[0]), &(rms_z[0]), &(bins_z_nt_er[0]), &(rms_z_er[0]));
199 krms_z_vs_ntrk->GetYaxis()->SetTitle(
coordinate +
" scale factor from RMS");
200 krms_z_vs_ntrk->GetXaxis()->SetTitle(xAxisLabel);
201 krms_z_vs_ntrk->SetTitle(
"scaleFactor" +
coordinate +
"_RMS");
202 krms_z_vs_ntrk->SetName(
"scaleFactor" +
coordinate +
"_" + versus +
"_RMS");
204 TGraphErrors* kgs_z_vs_ntrk =
new TGraphErrors(
205 bins_z_nt.size(), &(bins_z_nt[0]), &(sigma_z[0]), &(bins_z_nt_er[0]), &(sigma_z_er[0]));
206 kgs_z_vs_ntrk->GetYaxis()->SetTitle(
coordinate +
" scale factor from gauss fit");
207 kgs_z_vs_ntrk->GetXaxis()->SetTitle(xAxisLabel);
208 kgs_z_vs_ntrk->SetTitle(
"scaleFactor" +
coordinate +
"_Fit");
209 kgs_z_vs_ntrk->SetName(
"scaleFactor_" +
coordinate +
"_" + versus +
"_Fit");
212 float maxFitRange(100.);
213 float minFitRange(2.);
214 if (versus ==
"SumPt2") {
219 const Double_t* kgs_z_ntrk_fit_er;
220 int fitResKFactorMethod = 2;
221 if (fitResKFactorMethod == 1) {
223 kgs_z_vs_ntrk->Fit(
"pol2",
"Q",
"", minFitRange, maxFitRange);
224 kgs_z_vs_ntrk->GetFunction(
"pol2")->SetLineColor(kRed);
225 kgs_z_ntrk_fit = kgs_z_vs_ntrk->GetFunction(
"pol2");
226 kgs_z_ntrk_fit_er = kgs_z_ntrk_fit->GetParErrors();
227 }
else if (fitResKFactorMethod == 2) {
229 kgs_z_vs_ntrk->Fit(
"pol1",
"Q",
"", minFitRange, maxFitRange);
230 kgs_z_vs_ntrk->GetFunction(
"pol1")->SetLineColor(kRed);
231 kgs_z_ntrk_fit = kgs_z_vs_ntrk->GetFunction(
"pol1");
232 kgs_z_ntrk_fit_er = kgs_z_ntrk_fit->GetParErrors();
233 }
else if (fitResKFactorMethod == 3) {
234 TF1* kgsFitFcn =
new TF1(
"kgsFitFcn",
scaleFactorFitFcn, minFitRange, maxFitRange, 3);
235 kgsFitFcn->SetParameter(0, minFitRange);
236 kgsFitFcn->SetParameter(1, 1.0);
237 kgsFitFcn->SetParameter(2, 1.0);
238 for (
int ifit = 0; ifit < 1; ifit++)
239 kgs_z_vs_ntrk->Fit(kgsFitFcn,
"Q");
240 kgs_z_vs_ntrk->Fit(kgsFitFcn,
"Q");
241 kgs_z_ntrk_fit = kgsFitFcn;
242 kgs_z_ntrk_fit_er = kgsFitFcn->GetParErrors();
247 }
else if (fitResKFactorMethod == 4) {
249 kgs_z_vs_ntrk->Fit(
"pol0",
"Q",
"", minFitRange, maxFitRange);
250 kgs_z_vs_ntrk->GetFunction(
"pol0")->SetLineColor(kRed);
251 kgs_z_ntrk_fit = kgs_z_vs_ntrk->GetFunction(
"pol0");
252 kgs_z_ntrk_fit_er = kgs_z_ntrk_fit->GetParErrors();
258 int nbins_z_err_ntrk = h_Vrt_err_vs_Something->GetNbinsX();
260 std::vector<float> av_err_z;
261 std::vector<float> av_err_z_er;
264 std::vector<float> err_bins_z_nt;
265 std::vector<float> err_bins_z_nt_er;
266 std::vector<float> res_z;
267 std::vector<float> res_z_er;
271 for (
int bin_count = 1; bin_count <= nbins_z_err_ntrk; ++bin_count) {
272 err_bins_z_nt.push_back(h_Vrt_err_vs_Something->GetXaxis()->GetBinCenter(bin_count));
273 err_bins_z_nt_er.push_back(h_Vrt_err_vs_Something->GetXaxis()->GetBinWidth(bin_count) / 2.);
275 TH1D* profileY = h_Vrt_err_vs_Something->ProjectionY(
"projectionErrors", bin_count, bin_count,
"e");
280 float mean = profileY->GetMean();
281 float mean_error = profileY->GetMeanError();
291 av_err_z.push_back(
mean);
292 av_err_z_er.push_back(mean_error);
299 if (fitResKFactorMethod == 1) {
300 pr_er =
error_func(bin_count, kgs_z_ntrk_fit_er);
301 }
else if (fitResKFactorMethod == 2) {
302 val = h_Vrt_err_vs_Something->GetXaxis()->GetBinCenter(bin_count);
303 pr_er = TMath::Power(kgs_z_ntrk_fit_er[1] *
val, 2) + TMath::Power(kgs_z_ntrk_fit_er[0], 2);
304 pr_er = TMath::Sqrt(pr_er);
307 }
else if (fitResKFactorMethod == 3) {
308 val = h_Vrt_err_vs_Something->GetXaxis()->GetBinCenter(bin_count);
310 pr_er = kgs_z_ntrk_fit_er[2];
311 }
else if (fitResKFactorMethod == 4) {
312 pr_er = kgs_z_ntrk_fit_er[0];
315 res_z.push_back(
mean * kgs_z_ntrk_fit->Eval(h_Vrt_err_vs_Something->GetXaxis()->GetBinCenter(bin_count)));
316 res_z_er.push_back(TMath::Sqrt(TMath::Power(mean_error *
317 kgs_z_ntrk_fit->Eval(h_Vrt_err_vs_Something->GetXaxis()->GetBinCenter(
319 2) + TMath::Power(pr_er *
mean, 2)));
325 TGraphErrors* res_z_vs_ntrk =
326 new TGraphErrors(err_bins_z_nt.size(), &(err_bins_z_nt[0]), &(res_z[0]), &(err_bins_z_nt_er[0]), &(res_z_er[0]));
327 res_z_vs_ntrk->GetYaxis()->SetTitle(
coordinate +
" Vertex Resolution [mm]");
328 res_z_vs_ntrk->GetXaxis()->SetTitle(xAxisLabel);
329 res_z_vs_ntrk->SetTitle(
coordinate +
" Vertex Resolution");
330 res_z_vs_ntrk->SetName(
"resolution_" +
coordinate +
"_" + versus);
349 if (versus ==
"Ntrk") res_z_vs_ntrk->GetXaxis()->SetRangeUser(0., 100.);
350 else res_z_vs_ntrk->GetXaxis()->SetRangeUser(0., 20.);
351 res_z_vs_ntrk->GetYaxis()->SetRangeUser(0.0025, 1.);
352 res_z_vs_ntrk->Write(
"", TObject::kOverwrite);
353 delete res_z_vs_ntrk;
355 if (versus ==
"Ntrk") krms_z_vs_ntrk->GetXaxis()->SetRangeUser(0., 100.);
356 else krms_z_vs_ntrk->GetXaxis()->SetRangeUser(0., 20.);
357 krms_z_vs_ntrk->GetYaxis()->SetRangeUser(0.5, 1.3);
358 krms_z_vs_ntrk->Write(
"", TObject::kOverwrite);
359 delete krms_z_vs_ntrk;
361 if (versus ==
"Ntrk") kgs_z_vs_ntrk->GetXaxis()->SetRangeUser(0., 100.);
362 else kgs_z_vs_ntrk->GetXaxis()->SetRangeUser(0., 20.);
363 kgs_z_vs_ntrk->GetYaxis()->SetRangeUser(0.5, 1.3);
364 kgs_z_vs_ntrk->Write(
"", TObject::kOverwrite);
365 delete kgs_z_vs_ntrk;
381 if (h_Vrt_split_tag_ntrk == 0)
return;
384 if (h_Vrt_split_probe_ntrk == 0)
return;
386 TH1F* h_Vrt_split_matched_tag_ntrk = (
TH1F*)
gDirectory->Get(
"Vrt_split_matched_tag_ntrk");
387 if (h_Vrt_split_matched_tag_ntrk == 0)
return;
389 TH1F* h_Vrt_split_matched_probe_ntrk = (
TH1F*)
gDirectory->Get(
"Vrt_split_matched_probe_ntrk");
390 if (h_Vrt_split_matched_probe_ntrk == 0)
return;
393 if (h_Vrt_split_dist_tag == 0)
return;
396 if (h_Vrt_split_dist_probe == 0)
return;
399 TGraphAsymmErrors* g_Vrt_rec_eff_m1_split_vs_ntrk =
new TGraphAsymmErrors();
401 g_Vrt_rec_eff_m1_split_vs_ntrk->BayesDivide(h_Vrt_split_probe_ntrk, h_Vrt_split_tag_ntrk);
402 g_Vrt_rec_eff_m1_split_vs_ntrk->SetName(
"g_RecEff_M1");
404 TGraphAsymmErrors* g_Vrt_sel_eff_m1_split_vs_ntrk =
new TGraphAsymmErrors();
405 g_Vrt_sel_eff_m1_split_vs_ntrk->BayesDivide(h_Vrt_split_matched_probe_ntrk, h_Vrt_split_matched_tag_ntrk);
406 g_Vrt_sel_eff_m1_split_vs_ntrk->SetName(
"g_SelEff_M1");
409 g_Vrt_rec_eff_m1_split_vs_ntrk->GetHistogram()->GetXaxis()->SetTitle(
"Number of tracks");
410 g_Vrt_rec_eff_m1_split_vs_ntrk->GetHistogram()->GetYaxis()->SetTitle(
"Reconstruction efficiency");
411 g_Vrt_rec_eff_m1_split_vs_ntrk->SetMarkerStyle(20);
412 g_Vrt_rec_eff_m1_split_vs_ntrk->Write(
"", TObject::kOverwrite);
413 delete g_Vrt_rec_eff_m1_split_vs_ntrk;
415 g_Vrt_sel_eff_m1_split_vs_ntrk->GetHistogram()->GetXaxis()->SetTitle(
"Number of tracks");
416 g_Vrt_sel_eff_m1_split_vs_ntrk->GetHistogram()->GetYaxis()->SetTitle(
"Selection Efficiency");
417 g_Vrt_sel_eff_m1_split_vs_ntrk->SetMarkerStyle(20);
418 g_Vrt_sel_eff_m1_split_vs_ntrk->Write(
"", TObject::kOverwrite);
419 delete g_Vrt_sel_eff_m1_split_vs_ntrk;
429 histo->Fit(
"gaus",
"Q0",
"", -2, 2);
430 TF1* func =
histo->GetFunction(
"gaus");
431 double actualSigma = func->GetParameter(2);
432 double actualSigmaErr = func->GetParError(2);
434 for (
int u = 0;
u < 10;
u++) {
435 histo->Fit(
"gaus",
"Q0",
"", -sigmas * actualSigma, sigmas * actualSigma);
436 func =
histo->GetFunction(
"gaus");
437 actualSigma = func->GetParameter(2);
438 actualSigmaErr = func->GetParError(2);
441 results.push_back(actualSigma);
442 results.push_back(actualSigmaErr);
453 return(TMath::Power(
par[0], 2) +
x * TMath::Power(
par[1], 2) + TMath::Power(
x *
par[2], 2));
461 return (
x[0] >=
par[0]) * TMath::Erf(
par[1] *
x[0] -
par[0]) *
par[2];