Belle II Software development
HadronBgPrep.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <cdc/calibration/CDCdEdx/HadronBgPrep.h>
10using namespace Belle2;
11
13{
14 m_bgBins = 15;
15 m_bgMin = 0.1;
16 m_bgMax = 40;
17 m_cosBins = 18;
18 m_cosMin = -1.0;
19 m_cosMax = 1.0;
20 m_injBins = 20;
21 m_injMin = 0;
22 m_injMax = 80000;
23 m_nhitBins = 10;
24 m_nhitMin = 7;
25 m_nhitMax = 39;
26 m_cut = 0.5;
27}
28HadronBgPrep::HadronBgPrep(int bgbins, double lowerbg, double upperbg, int cosbins, double lowercos, double uppercos, int injbins,
29 double lowerinj, double upperinj, int nhitbins, double lowernhit, double uppernhit, double cut)
30{
31 m_bgBins = bgbins;
32 m_bgMin = lowerbg;
33 m_bgMax = upperbg;
34 m_cosBins = cosbins;
35 m_cosMin = lowercos;
36 m_cosMax = uppercos;
37 m_injBins = injbins;
38 m_injMin = lowerinj;
39 m_injMax = upperinj;
40 m_nhitBins = nhitbins;
41 m_nhitMin = lowernhit;
42 m_nhitMax = uppernhit;
43 m_cut = cut;
44}
45
46void HadronBgPrep::prepareSample(std::shared_ptr<TTree> hadron, TFile*& outfile, const std::string& suffix,
47 const std::string& bgcurvefile,
48 const std::string& bgsigmafile, const std::string& pdg, bool ismakePlots)
49{
50
51 std::vector<TH1F*> hdedx_bg, hchi_bg, hionzsigma_bg;
52 std::map<int, std::vector<TH1F*>> hchi_inj, hchicos_allbg, hchicos_1by3bg, hchicos_2by3bg, hchicos_3by3bg;
53
54 //define histograms
55 defineHisto(hdedx_bg, "dedx", "bg", pdg.data());
56 defineHisto(hchi_bg, "chi", "bg", pdg.data());
57 defineHisto(hionzsigma_bg, "ionzsigma", "bg", pdg.data());
58 for (int i = 0; i < 2; ++i) {
59
60 defineHisto(hchi_inj[i], "chi", Form("inj_%d", i), pdg.data());
61
62 std::string charge = "pos";
63 if (i == 1) charge = "neg";
64
65 defineHisto(hchicos_allbg[i], "chi", Form("%s_allbg_cos", charge.data()), pdg.data());
66 defineHisto(hchicos_1by3bg[i], "chi", Form("%s_1b3bg_cos", charge.data()), pdg.data());
67 defineHisto(hchicos_2by3bg[i], "chi", Form("%s_2b3bg_cos", charge.data()), pdg.data());
68 defineHisto(hchicos_3by3bg[i], "chi", Form("%s_3b3bg_cos", charge.data()), pdg.data());
69 }
70
71
72 // Create some containers to calculate averages
73 for (int i = 0; i < m_bgBins; ++i) {
74 m_sumcos.push_back(0.);
75 m_sumbg.push_back(0.);
76 m_sumres_square.push_back(0.);
77 m_sumsize.push_back(0.);
78 m_means.push_back(0.);
79 m_errors.push_back(0.);
80 }
81
82 for (int i = 0; i < m_injBins; ++i) {
83 m_injsize.push_back(0.);
84 m_suminj.push_back(0.);
85 }
86
87 // --------------------------------------------------
88 // LOOP OVER EVENTS AND FILL CONTAINERS
89 // --------------------------------------------------
90 // Fill the histograms to be fitted
91
92 double dedxnosat; // dE/dx w/o HS correction (use to get new parameters)
93 double p; // track momentum
94 double costh; // cosine of track polar angle
95 double timereso;
96 int nhit; // number of hits on this track
97 double injtime; //injection time
98 double isher;
99 int charge;
100
101 hadron->SetBranchAddress("dedxnosat", &dedxnosat); //must be without any HS correction
102 hadron->SetBranchAddress("p", &p);
103 hadron->SetBranchAddress("costh", &costh);
104 hadron->SetBranchAddress("timereso", &timereso);
105 hadron->SetBranchAddress("nhits", &nhit);
106 hadron->SetBranchAddress("injtime", &injtime);
107 hadron->SetBranchAddress("injring", &isher);
108 hadron->SetBranchAddress("charge", &charge);
109
110 double mass = getParticleMass(pdg);
111 if (mass == 0.0) B2FATAL("Mass of particle " << pdg.data() << " is zero");
112
113 const double cosstep = (m_cosMax - m_cosMin) / m_cosBins;
114 const double tstep = (m_injMax - m_injMin) / m_injBins;
115
116 CDCDedxMeanPred mbg;
117 mbg.setParameters(bgcurvefile);
119 sbg.setParameters(bgsigmafile);
120
121 CDCDedxHadSat had;
122 had.setParameters("sat-pars.fit.txt");
123
124 unsigned int entries = hadron->GetEntries();
125
126 for (unsigned int index = 0; index < entries; ++index) {
127
128 hadron->GetEvent(index);
129
130 int chg = (charge < 0) ? 1 : 0;
131 double bg = fabs(p) / mass;
132
133 // clean up bad events and restrict the momentum range
134 if (nhit < 0 || nhit > 100)continue;
135 if (injtime <= 0 || isher < 0)continue;
136
137 if (dedxnosat <= 0)continue;
138 if (costh != costh)continue;
139 if (bg < m_bgMin || bg > m_bgMax)continue;
140
141 if (fabs(p) > 7.0)continue;
142
143 if (pdg == "proton") if ((dedxnosat - 0.45) * abs(p) * abs(p) < m_cut)continue;
144
145 int bgBin = (int)((bg - m_bgMin) / (m_bgMax - m_bgMin) * m_bgBins);
146
147 double dedx_new = had.D2I(costh, had.I2D(costh, 1.0) * dedxnosat);
148
149 hdedx_bg[bgBin]->Fill(dedx_new);
150
151 //get predicted value of dEdx mean
152 double dedx_cur = mbg.getMean(bg);
153 double dedx_res = sbg.getSigma(dedx_new, nhit, costh, timereso);
154
155 if (dedx_res == 0) {
156 B2INFO("RESOLUTION IS ZERO!!!");
157 continue;
158 }
159
160 //Modified chi_new values
161 double chi_new = (dedx_new - dedx_cur) / dedx_res;
162
163 hchi_bg[bgBin]->Fill(chi_new);
164
165 double ionz_res = sbg.cosPrediction(costh) * sbg.nhitPrediction(nhit) * timereso;
166
167 if (ionz_res == 0) {
168 B2INFO("RESOLUTION costh * nhit * timereso IS ZERO!!!");
169 continue;
170 }
171
172 hionzsigma_bg[bgBin]->Fill((dedx_new - dedx_cur) / ionz_res);
173
174 m_sumcos[bgBin] += costh;
175 m_sumbg[bgBin] += bg;
176 m_sumres_square[bgBin] += pow(dedx_res, 2);
177 m_sumsize[bgBin] += 1;
178
179 // make histograms of dE/dx vs. cos(theta) for validation
180 int icos = (int)((costh + 1) / cosstep);
181 hchicos_allbg[chg][icos]->Fill(chi_new);
182
183 if (bgBin <= int(m_bgBins / 3)) hchicos_1by3bg[chg][icos]->Fill(chi_new);
184 else if (bgBin <= int(2 * m_bgBins / 3)) hchicos_2by3bg[chg][icos]->Fill(chi_new);
185 else hchicos_3by3bg[chg][icos]->Fill(chi_new);
186
187 if (injtime > m_injMax) injtime = m_injMax - 10.0;
188 int wr = 0;
189 if (isher > 0.5) wr = 1;
190 int injBin = (int)((injtime - m_injMin) / tstep);
191 hchi_inj[wr][injBin]->Fill(chi_new);
192
193 m_suminj[injBin] += injtime;
194 m_injsize[injBin] += 1;
195
196 } // end of event loop
197
198 // --------------------------------------------------
199 // FIT IN BINS OF BETA-GAMMA AND COS(THETA)
200 // --------------------------------------------------
201 // fit the histograms with Gaussian functions
202 // and extract the means and errors
203
204 setPars(outfile, pdg, hdedx_bg, hchi_bg, hionzsigma_bg, hchi_inj);
205
206 // Plot the histograms
207 if (ismakePlots) {
208
209 plotDist(hdedx_bg, Form("fits_dedx_inbg_%s_%s", suffix.data(), pdg.data()), m_bgBins);
210 plotDist(hchi_bg, Form("fits_chi_inbg_%s_%s", suffix.data(), pdg.data()), m_bgBins);
211 plotDist(hionzsigma_bg, Form("fits_ionzreso_inbg_%s_%s", suffix.data(), pdg.data()), m_bgBins);
212 plotDist(hchi_inj, Form("fits_chi_ininj_%s_%s", suffix.data(), pdg.data()), m_injBins);
213 printCanvasCos(hchicos_allbg, hchicos_1by3bg, hchicos_2by3bg, hchicos_3by3bg, pdg, suffix);
214 }
215
216 // delete histograms
217 clearVars();
218 deleteHistos(hdedx_bg);
219 deleteHistos(hchi_bg);
220 deleteHistos(hionzsigma_bg);
221 for (int i = 0; i < 2; ++i) {
222 deleteHistos(hchi_inj[i]);
223 deleteHistos(hchicos_allbg[i]);
224 deleteHistos(hchicos_1by3bg[i]);
225 deleteHistos(hchicos_2by3bg[i]);
226 deleteHistos(hchicos_3by3bg[i]);
227 }
228
229}
230
231//------------------------------------
232void HadronBgPrep::defineHisto(std::vector<TH1F*>& htemp, const std::string& svar, const std::string& stype, const std::string& pdg)
233{
234 int nbdEdx = 200, nbchi = 200;
235 double dedxMax = 20.0;
236
237 if (pdg == "pion") {
238 nbdEdx = 400, dedxMax = 4.0;
239 } else if (pdg == "kaon") {
240 nbdEdx = 500, dedxMax = 5.0;
241 } else if (pdg == "proton") {
242 nbdEdx = 1500, dedxMax = 30.0;
243 nbchi = 350;
244 } else if (pdg == "muon") {
245 nbdEdx = 300, dedxMax = 3.0;
246 } else if (pdg == "electron") {
247 nbdEdx = 200, dedxMax = 2.0;
248 }
249
250 int bins;
251 double min, max;
252
253 if (stype == "bg") bins = m_bgBins, min = m_bgMin, max = m_bgMax ;
254 else if (stype == "inj_0" || stype == "inj_1") bins = m_injBins, min = m_injMin, max = m_injMax ;
255 else if (stype == "nhit") bins = m_nhitBins, min = m_nhitMin, max = m_nhitMax ;
256 else bins = m_cosBins, min = m_cosMin, max = m_cosMax ;
257
258 double step = (max - min) / bins;
259 if (stype == "nhit") step = (max - min + 1) / bins;
260
261 for (int j = 0; j < bins; ++j) {
262
263 double start = min + j * step;
264 double end = start + step;
265 if (stype == "nhit") end = int((start + step) * 0.99999);
266
267 // initialize the histograms
268 std::string histname = Form("%s_%s_%s_%d", pdg.data(), svar.data(), stype.data(), j);
269 std::string title = Form("%s_%s_%s (%.02f, %.02f)", pdg.data(), svar.data(), stype.data(), start, end);
270
271 if (svar == "dedx")
272 htemp.push_back(new TH1F(histname.data(), title.data(), nbdEdx, 0, dedxMax));
273 else if (svar == "chi")
274 htemp.push_back(new TH1F(histname.data(), title.data(), nbchi, -10.0, 10.0));
275 else
276 htemp.push_back(new TH1F(histname.data(), title.data(), 300, -3.0, 3.0));
277
278 htemp[j]->GetXaxis()->SetTitle(Form("%s", svar.data()));
279 htemp[j]->GetYaxis()->SetTitle("Entries");
280 }
281}
282
283//------------------------------------
284void HadronBgPrep::plotDist(std::map<int, std::vector<TH1F*>>& hist, const std::string& sname, int bins)
285{
286
287 TCanvas ctmp(Form("cdcdedx_%s", sname.data()), "", 1200, 600);
288 ctmp.Divide(2, 1);
289 ctmp.SetBatch(kTRUE);
290
291 std::stringstream psname;
292 psname << Form("plots/HadronPrep/%s.pdf[", sname.data());
293 ctmp.Print(psname.str().c_str());
294 psname.str("");
295 psname << Form("plots/HadronPrep/%s.pdf", sname.data());
296 for (int j = 0; j < bins; ++j) {
297
298 for (int i = 0 ; i < 2; ++i) {
299 ctmp.cd(i + 1);
300 hist[i][j]->SetFillColorAlpha(i + 5, 0.25);
301 hist[i][j]->Draw();
302 }
303 ctmp.Print(psname.str().c_str());
304 }
305 psname.str("");
306 psname << Form("plots/HadronPrep/%s.pdf]", sname.data());
307 ctmp.Print(psname.str().c_str());
308}
309
310//------------------------------------
311void HadronBgPrep::plotDist(std::vector<TH1F*>& hist, const std::string& sname, int nbins)
312{
313
314 TCanvas ctmp(Form("cdcdedx_%s", sname.data()), "", 1200, 1200);
315 ctmp.Divide(2, 2);
316 ctmp.SetBatch(kTRUE);
317
318 std::stringstream psname;
319 psname << Form("plots/HadronPrep/%s.pdf[", sname.data());
320 ctmp.Print(psname.str().c_str());
321 psname.str("");
322 psname << Form("plots/HadronPrep/%s.pdf", sname.data());
323
324 for (int i = 0 ; i < nbins; ++i) {
325 ctmp.cd(i % 4 + 1);
326 hist[i]->SetFillColor(kYellow - 9);
327 hist[i]->Draw();
328
329 if ((i + 1) % 4 == 0 || (i + 1) == nbins) {
330 ctmp.Print(psname.str().c_str());
331 ctmp.Clear("D");
332 }
333 }
334 psname.str("");
335 psname << Form("plots/HadronPrep/%s.pdf]", sname.data());
336 ctmp.Print(psname.str().c_str());
337}
338
339//----------------------------------------
340void HadronBgPrep::setPars(TFile*& outfile, std::string pdg, std::vector<TH1F*>& hdedx_bg, std::vector<TH1F*>& hchi_bg,
341 std::vector<TH1F*>& hionzsigma_bg, std::map<int, std::vector<TH1F*>>& hchi_inj)
342{
343 outfile->cd();
344
345 TTree* satTree = new TTree(Form("%s", pdg.data()), "dE/dx m_means and m_errors");
346 double satbg; // beta-gamma value for this bin
347 double satcosth; // cos(theta) value for this bin
348 double satdedx; // mean dE/dx value for this bin
349 double satdedxerr; // error on ^
350 double satdedxwidth;// width of ^ distribution
351
352 double satchi; // mean chi value for this bin
353 double satchierr; // error on ^
354 double satchiwidth; // width of ^ distribution
355 double satchiwidth_err;
356 double sationzres; // width of dedx reso
357
358 double satbg_avg; // average beta-gamma value for this sample
359 double satcosth_avg; // average cos(theta) value for this sample
360 double satdedxres_avg; // average dE/dx error squared for this sample
361
362 satTree->Branch("bg", &satbg, "bg/D");
363 satTree->Branch("costh", &satcosth, "costh/D");
364
365 //dEdx related variables
366 satTree->Branch("dedx", &satdedx, "dedx/D");
367 satTree->Branch("dedxerr", &satdedxerr, "dedxerr/D");
368 satTree->Branch("dedxwidth", &satdedxwidth, "dedxwidth/D");
369
370 // Chi related variables
371 satTree->Branch("chimean", &satchi, "chimean/D");
372 satTree->Branch("chimean_err", &satchierr, "chimean_err/D");
373 satTree->Branch("chisigma", &satchiwidth, "chisigma/D");
374 satTree->Branch("chisigma_err", &satchiwidth_err, "chisigma_err/D");
375
376 satTree->Branch("ionzres", &sationzres, "ionzres/D");
377
378 // Other variables
379 satTree->Branch("bg_avg", &satbg_avg, "bg_avg/D");
380 satTree->Branch("costh_avg", &satcosth_avg, "costh_avg/D");
381 satTree->Branch("dedxres_avg", &satdedxres_avg, "dedxres_avg/D");
382
383 double bgstep = (m_bgMax - m_bgMin) / m_bgBins;
384
385 for (int i = 0; i < m_bgBins; ++i) {
386
387 // fill some details for this bin
388 satbg = m_bgMin + 0.5 * bgstep + i * bgstep;
389 satcosth = 0.0;
390
391 satbg_avg = m_sumbg[i] / m_sumsize[i];
392 satcosth_avg = m_sumcos[i] / m_sumsize[i];
393 satdedxres_avg = m_sumres_square[i] / m_sumsize[i];
394
395 //1. -------------------------
396 // fit the dE/dx distribution in bins of beta-gamma
397 fit(hdedx_bg[i], pdg.data());
398 satdedx = m_means[i] = hdedx_bg[i]->GetFunction("gaus")->GetParameter(1);
399 satdedxerr = m_errors[i] = hdedx_bg[i]->GetFunction("gaus")->GetParError(1);
400 satdedxwidth = hdedx_bg[i]->GetFunction("gaus")->GetParameter(2);
401
402 //2. -------------------------
403 // fit the chi distribution in bins of beta-gamma
404 fit(hchi_bg[i], pdg.data());
405 satchi = hchi_bg[i]->GetFunction("gaus")->GetParameter(1);
406 satchierr = hchi_bg[i]->GetFunction("gaus")->GetParError(1);
407 satchiwidth = hchi_bg[i]->GetFunction("gaus")->GetParameter(2);
408 satchiwidth_err = hchi_bg[i]->GetFunction("gaus")->GetParError(2);
409
410 //3. -------------------------
411 // fit the chi distribution in bins of beta-gamma
412 fit(hionzsigma_bg[i], pdg.data());
413 sationzres = hionzsigma_bg[i]->GetFunction("gaus")->GetParameter(2);
414
415 // fill the tree for this bin
416 satTree->Fill();
417 }
418 // write out the data to file
419 satTree->Write();
420
421 // --------------------------------------------------
422 // FIT IN BINS OF INJECTION TIME FOR VALIDATION
423 // --------------------------------------------------
424
425 std::string svar = "ler";
426 for (int ir = 0; ir < 2; ir++) {
427
428 if (ir == 1) svar = "her";
429
430 TTree* injTree = new TTree(Form("%s_%s", pdg.data(), svar.data()), "chi m_means and m_errors");
431 double inj_avg;
432 double mean;
433 double mean_err;
434 double sigma;
435 double sigma_err;
436
437 injTree->Branch("inj_avg", &inj_avg, "inj_avg/D");
438 injTree->Branch("chimean", &mean, "chimean/D");
439 injTree->Branch("chimean_err", &mean_err, "chimean_err/D");
440 injTree->Branch("chisigma", &sigma, "chisigma/D");
441 injTree->Branch("chisigma_err", &sigma_err, "chisigma_err/D");
442
443 for (int i = 0; i < m_injBins; ++i) {
444
445 inj_avg = m_suminj[i] / m_injsize[i];
446
447 // fit the dE/dx distribution in bins of injection time'
448 fit(hchi_inj[ir][i], pdg.data());
449
450 mean = hchi_inj[ir][i]->GetFunction("gaus")->GetParameter(1);
451 mean_err = hchi_inj[ir][i]->GetFunction("gaus")->GetParError(1);
452 sigma = hchi_inj[ir][i]->GetFunction("gaus")->GetParameter(2);
453 sigma_err = hchi_inj[ir][i]->GetFunction("gaus")->GetParError(2);
454
455 injTree->Fill();
456 }
457
458 injTree->Write();
459 }
460 outfile->Write();
461}
462
463//----------------------------------------
464void HadronBgPrep::fit(TH1F*& hist, const std::string& pdg)
465{
466 gstatus status;
467 if (pdg == "pion") fitGaussianWRange(hist, status, 1.0);
468 else fitGaussianWRange(hist, status, 2.0);
469
470 hist->SetFillColorAlpha(kAzure + 1, 0.30);
471
472 if (status == OK) {
473 double mean = hist->GetFunction("gaus")->GetParameter(1);
474 double meanerr = hist->GetFunction("gaus")->GetParError(1);
475 double width = hist->GetFunction("gaus")->GetParameter(2);
476
477 std::string title = Form("#mu_{fit}: %0.03f #pm %0.03f, #sigma_{fit}: %0.03f", mean, meanerr, width);
478 hist->SetTitle(Form("%s, %s", hist->GetTitle(), title.data()));
479 }
480
481}
482
483//----------------------------------------
484void HadronBgPrep::fitGaussianWRange(TH1F*& temphist, gstatus& status, double sigmaR)
485{
486 double histmean = temphist->GetMean();
487 double histrms = temphist->GetRMS();
488 temphist->GetXaxis()->SetRangeUser(histmean - 5.0 * histrms, histmean + 5.0 * histrms);
489
490 int fs = temphist->Fit("gaus", "Q");
491 if (fs != 0) {
492 B2INFO(Form("\tFit (round 1) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
493 status = Failed;
494 return;
495 } else {
496 double mean = temphist->GetFunction("gaus")->GetParameter(1);
497 double width = temphist->GetFunction("gaus")->GetParameter(2);
498 temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
499 fs = temphist->Fit("gaus", "QR", "", mean - sigmaR * width, mean + sigmaR * width);
500 if (fs != 0) {
501 B2INFO(Form("\tFit (round 2) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
502 status = Failed;
503 return;
504 } else {
505 temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
506 B2INFO(Form("\tFit for hist (%s) successful (status = %d)", temphist->GetName(), fs));
507 status = OK;
508 }
509 }
510}
511
512
513// --------------------------------------------------
514// FIT IN BINS OF COS(THETA) FOR VALIDATION
515// --------------------------------------------------
516void HadronBgPrep::printCanvasCos(std::map<int, std::vector<TH1F*>>& hchicos_allbg,
517 std::map<int, std::vector<TH1F*>>& hchicos_1by3bg, std::map<int, std::vector<TH1F*>>& hchicos_2by3bg,
518 std::map<int, std::vector<TH1F*>>& hchicos_3by3bg, const std::string& pdg, const std::string& suffix)
519{
520
521 std::vector<std::vector<double>> chicos(2, std::vector<double>(m_cosBins));
522 std::vector<std::vector<double>> sigmacos(2, std::vector<double>(m_cosBins));
523 std::vector<std::vector<double>> chicoserr(2, std::vector<double>(m_cosBins));
524 std::vector<std::vector<double>> sigmacoserr(2, std::vector<double>(m_cosBins));
525
526 std::vector<std::vector<double>> chicos_1b3bg(2, std::vector<double>(m_cosBins));
527 std::vector<std::vector<double>> sigmacos_1b3bg(2, std::vector<double>(m_cosBins));
528 std::vector<std::vector<double>> chicos_1b3bgerr(2, std::vector<double>(m_cosBins));
529 std::vector<std::vector<double>> sigmacos_1b3bgerr(2, std::vector<double>(m_cosBins));
530
531 std::vector<std::vector<double>> chicos_2b3bg(2, std::vector<double>(m_cosBins));
532 std::vector<std::vector<double>> sigmacos_2b3bg(2, std::vector<double>(m_cosBins));
533 std::vector<std::vector<double>> chicos_2b3bgerr(2, std::vector<double>(m_cosBins));
534 std::vector<std::vector<double>> sigmacos_2b3bgerr(2, std::vector<double>(m_cosBins));
535
536 std::vector<std::vector<double>> chicos2(2, std::vector<double>(m_cosBins));
537 std::vector<std::vector<double>> sigmacos2(2, std::vector<double>(m_cosBins));
538 std::vector<std::vector<double>> chicos2err(2, std::vector<double>(m_cosBins));
539 std::vector<std::vector<double>> sigmacos2err(2, std::vector<double>(m_cosBins));
540
541 for (int c = 0; c < 2; ++c) {
542 for (int i = 0; i < m_cosBins; ++i) {
543 if (hchicos_allbg[c][i]->Integral() > 100) {
544 fit(hchicos_allbg[c][i], pdg.data());
545 chicos[c][i] = hchicos_allbg[c][i]->GetFunction("gaus")->GetParameter(1);
546 chicoserr[c][i] = hchicos_allbg[c][i]->GetFunction("gaus")->GetParError(1);
547 sigmacos[c][i] = hchicos_allbg[c][i]->GetFunction("gaus")->GetParameter(2);
548 sigmacoserr[c][i] = hchicos_allbg[c][i]->GetFunction("gaus")->GetParError(2);
549 }
550
551 if (hchicos_1by3bg[c][i]->Integral() > 100) {
552 fit(hchicos_1by3bg[c][i], pdg.data());
553 chicos_1b3bg[c][i] = hchicos_1by3bg[c][i]->GetFunction("gaus")->GetParameter(1);
554 chicos_1b3bgerr[c][i] = hchicos_1by3bg[c][i]->GetFunction("gaus")->GetParError(1);
555 sigmacos_1b3bg[c][i] = hchicos_1by3bg[c][i]->GetFunction("gaus")->GetParameter(2);
556 sigmacos_1b3bgerr[c][i] = hchicos_1by3bg[c][i]->GetFunction("gaus")->GetParError(2);
557 }
558
559 if (hchicos_2by3bg[c][i]->Integral() > 100) {
560 fit(hchicos_2by3bg[c][i], pdg.data());
561 chicos_2b3bg[c][i] = hchicos_2by3bg[c][i]->GetFunction("gaus")->GetParameter(1);
562 chicos_2b3bgerr[c][i] = hchicos_2by3bg[c][i]->GetFunction("gaus")->GetParError(1);
563 sigmacos_2b3bg[c][i] = hchicos_2by3bg[c][i]->GetFunction("gaus")->GetParameter(2);
564 sigmacos_2b3bgerr[c][i] = hchicos_2by3bg[c][i]->GetFunction("gaus")->GetParError(2);
565 }
566
567 if (hchicos_3by3bg[c][i]->Integral() > 100) {
568 fit(hchicos_3by3bg[c][i], pdg.data());
569 chicos2[c][i] = hchicos_3by3bg[c][i]->GetFunction("gaus")->GetParameter(1);
570 chicos2err[c][i] = hchicos_3by3bg[c][i]->GetFunction("gaus")->GetParError(1);
571 sigmacos2[c][i] = hchicos_3by3bg[c][i]->GetFunction("gaus")->GetParameter(2);
572 sigmacos2err[c][i] = hchicos_3by3bg[c][i]->GetFunction("gaus")->GetParError(2);
573 }
574 }
575 }
576
577 plotDist(hchicos_allbg, Form("fits_chi_vscos_allbg_%s_%s", pdg.data(), suffix.data()), m_cosBins);
578 plotDist(hchicos_1by3bg, Form("fits_chi_vscos_1by3bg_%s_%s", pdg.data(), suffix.data()), m_cosBins);
579 plotDist(hchicos_2by3bg, Form("fits_chi_vscos_2by3bg_%s_%s", pdg.data(), suffix.data()), m_cosBins);
580 plotDist(hchicos_3by3bg, Form("fits_chi_vscos_3by3bg_%s_%s", pdg.data(), suffix.data()), m_cosBins);
581
582 const double cosstep = 2.0 / m_cosBins;
583 std::vector<double> cosArray(m_cosBins), cosArrayErr(m_cosBins);
584 for (int i = 0; i < m_cosBins; ++i) {
585 cosArray[i] = -1.0 + (i * cosstep + cosstep / 2.0); //finding bin centre
586 cosArrayErr[i] = 0.0;
587 }
588
589 TGraphErrors grchicos(m_cosBins, cosArray.data(), chicos[0].data(), cosArrayErr.data(), chicoserr[0].data());
590 TGraphErrors grchicos_1b3bg(m_cosBins, cosArray.data(), chicos_1b3bg[0].data(), cosArrayErr.data(), chicos_1b3bgerr[0].data());
591 TGraphErrors grchicos_2b3bg(m_cosBins, cosArray.data(), chicos_2b3bg[0].data(), cosArrayErr.data(), chicos_2b3bgerr[0].data());
592 TGraphErrors grchicos2(m_cosBins, cosArray.data(), chicos2[0].data(), cosArrayErr.data(), chicos2err[0].data());
593
594 TGraphErrors grchicosn(m_cosBins, cosArray.data(), chicos[1].data(), cosArrayErr.data(), chicoserr[1].data());
595 TGraphErrors grchicos_1b3bgn(m_cosBins, cosArray.data(), chicos_1b3bg[1].data(), cosArrayErr.data(), chicos_1b3bgerr[1].data());
596 TGraphErrors grchicos_2b3bgn(m_cosBins, cosArray.data(), chicos_2b3bg[1].data(), cosArrayErr.data(), chicos_2b3bgerr[1].data());
597 TGraphErrors grchicos2n(m_cosBins, cosArray.data(), chicos2[1].data(), cosArrayErr.data(), chicos2err[1].data());
598
599 TGraphErrors grsigmacos(m_cosBins, cosArray.data(), sigmacos[0].data(), cosArrayErr.data(), sigmacoserr[0].data());
600 TGraphErrors grsigmacos_1b3bg(m_cosBins, cosArray.data(), sigmacos_1b3bg[0].data(), cosArrayErr.data(),
601 sigmacos_1b3bgerr[0].data());
602 TGraphErrors grsigmacos_2b3bg(m_cosBins, cosArray.data(), sigmacos_2b3bg[0].data(), cosArrayErr.data(),
603 sigmacos_2b3bgerr[0].data());
604 TGraphErrors grsigmacos2(m_cosBins, cosArray.data(), sigmacos2[0].data(), cosArrayErr.data(), sigmacos2err[0].data());
605
606 TGraphErrors grsigmacosn(m_cosBins, cosArray.data(), sigmacos[1].data(), cosArrayErr.data(), sigmacoserr[1].data());
607 TGraphErrors grsigmacos_1b3bgn(m_cosBins, cosArray.data(), sigmacos_1b3bg[1].data(), cosArrayErr.data(),
608 sigmacos_1b3bgerr[1].data());
609 TGraphErrors grsigmacos_2b3bgn(m_cosBins, cosArray.data(), sigmacos_2b3bg[1].data(), cosArrayErr.data(),
610 sigmacos_2b3bgerr[1].data());
611 TGraphErrors grsigmacos2n(m_cosBins, cosArray.data(), sigmacos2[1].data(), cosArrayErr.data(), sigmacos2err[1].data());
612
613 TLine line0(-1, 0, 1, 0);
614 line0.SetLineStyle(kDashed);
615 line0.SetLineColor(kRed);
616
617 TLine line1(-1, 1, 1, 1);
618 line1.SetLineStyle(kDashed);
619 line1.SetLineColor(kRed);
620
621 TString ptype("");
622 double mass = 0.0;
623 if (pdg == "pion") {ptype += "#pi"; mass = Const::pion.getMass();}
624 else if (pdg == "kaon") {ptype += "K"; mass = Const::kaon.getMass();}
625 else if (pdg == "proton") { ptype += "p"; mass = Const::proton.getMass();}
626 else if (pdg == "muon") {ptype += "#mu"; mass = Const::muon.getMass();}
627 else if (pdg == "electron") {ptype += "e"; mass = Const::electron.getMass();}
628 if (mass == 0.0) B2FATAL("Mass of particle " << pdg.data() << " is zero");
629
630 int bglow = 0, bghigh = 0;
631 double bgstep = (m_bgMax - m_bgMin) / m_bgBins;
632
633 TLegend lchi(0.7, 0.75, 0.8, 0.85);
634 lchi.SetBorderSize(0);
635 lchi.SetFillColor(kYellow);
636 lchi.AddEntry(&grchicos, ptype + "^{+}", "pc");
637 lchi.AddEntry(&grchicosn, ptype + "^{-}", "pc");
638
639 TCanvas cchi(Form("cchi_%s", suffix.data()), "cchi", 700, 600);
640 cchi.Divide(2, 2);
641 cchi.cd(1);
642 FormatGraph(grchicos, 0, Form("all (%s) bg bins, p =(%0.02f, %0.02f)", pdg.data(), m_bgMin * mass, m_bgMax * mass));
643 FormatGraph(grchicosn, 1);
644 grchicos.Draw("AP");
645 grchicosn.Draw("P,same");
646 line0.Draw("same");
647 lchi.Draw("same");
648
649 cchi.cd(2);
650 bghigh = int(m_bgBins / 3);
651 FormatGraph(grchicos_1b3bg, 0, Form("first 1/3 bg bins, p =(%0.02f, %0.02f)", m_bgMin * mass, bghigh * bgstep * mass));
652 FormatGraph(grchicos_1b3bgn, 1);
653 grchicos_1b3bg.Draw("AP");
654 grchicos_1b3bgn.Draw("P,same");
655 line0.Draw("same");
656 lchi.Draw("same");
657
658 cchi.cd(3);
659 bglow = int(m_bgBins / 3);
660 bghigh = int(2 * m_bgBins / 3);
661 FormatGraph(grchicos_2b3bg, 0, Form("second 1/3 bg bins, p =(%0.02f, %0.02f)", bglow * bgstep * mass, bghigh * bgstep * mass));
662 FormatGraph(grchicos_2b3bgn, 1);
663 grchicos_2b3bg.Draw("AP");
664 grchicos_2b3bgn.Draw("P,same");
665 line0.Draw("same");
666 lchi.Draw("same");
667
668 cchi.cd(4);
669 bglow = int(2 * m_bgBins / 3);
670 FormatGraph(grchicos2, 0, Form("third 1/3 bg bins, p =(%0.02f, %0.02f)", bglow * bgstep * mass, m_bgMax * mass));
671 FormatGraph(grchicos2n, 1);
672 grchicos2.Draw("AP");
673 grchicos2n.Draw("P,same");
674 line0.Draw("same");
675 lchi.Draw("same");
676 cchi.SaveAs(Form("plots/HadronCal/Monitoring/mean_chi_vscos_inbg_%s_%s.pdf", pdg.data(), suffix.data()));
677
678 cchi.cd(1);
679 FormatGraph(grsigmacos, 2, Form("all (%s) bg bins, p =(%0.02f, %0.02f)", pdg.data(), m_bgMin * mass, m_bgMax * mass));
680 FormatGraph(grsigmacosn, 1);
681 grsigmacos.Draw("AP");
682 grsigmacosn.Draw("P,same");
683 line1.Draw("same");
684 lchi.Draw("same");
685
686 cchi.cd(2);
687 bghigh = int(m_bgBins / 3);
688 FormatGraph(grsigmacos_1b3bg, 2, Form("first 1/3 bg bins, p =(%0.02f, %0.02f)", m_bgMin * mass, bghigh * bgstep * mass));
689 FormatGraph(grsigmacos_1b3bgn, 1);
690 grsigmacos_1b3bg.Draw("AP");
691 grsigmacos_1b3bgn.Draw("P,same");
692 line1.Draw("same");
693 lchi.Draw("same");
694
695 cchi.cd(3);
696 bglow = int(m_bgBins / 3);
697 bghigh = int(2 * m_bgBins / 3);
698 FormatGraph(grsigmacos_2b3bg, 2, Form("second 1/3 bg bins, p =(%0.02f, %0.02f)", bglow * bgstep * mass, bghigh * bgstep * mass));
699 FormatGraph(grsigmacos_2b3bgn, 1);
700 grsigmacos_2b3bg.Draw("AP");
701 grsigmacos_2b3bgn.Draw("P,same");
702 line1.Draw("same");
703 lchi.Draw("same");
704
705 cchi.cd(4);
706 bglow = int(2 * m_bgBins / 3);
707 FormatGraph(grsigmacos2, 2, Form("third 1/3 bg bins, p =(%0.02f, %0.02f)", bglow * bgstep * mass, m_bgMax * mass));
708 FormatGraph(grsigmacos2n, 1);
709 grsigmacos2.Draw("AP");
710 grsigmacos2n.Draw("P,same");
711 line1.Draw("same");
712 lchi.Draw("same");
713 cchi.SaveAs(Form("plots/HadronCal/Monitoring/sigma_chi_vscos_inbg_%s_%s.pdf", pdg.data(), suffix.data()));
714
715}
Class to hold the hadron saturation functions.
Definition: CDCDedxHadSat.h:31
double I2D(double cosTheta, double I) const
hadron saturation parameterization part 2
double D2I(double cosTheta, double D) const
hadron saturation parameterization part 1
void setParameters()
set the parameters
Class to hold the prediction of mean as a function of beta-gamma (bg)
double getMean(double bg)
Return the predicted mean value as a function of beta-gamma (bg)
void setParameters(std::string infile)
set the parameters from file
Class to hold the prediction of resolution depending dE/dx, nhit, and cos(theta)
double getSigma(double dedx, double nhit, double cos, double timereso)
Return the predicted resolution depending on dE/dx, nhit, and cos(theta)
double cosPrediction(double cos)
Return sigma from the cos parameterization.
double nhitPrediction(double nhit)
Return sigma from the nhit parameterization.
void setParameters(std::string infile)
set the parameters from file
double getMass() const
Particle mass.
Definition: UnitConst.cc:353
static const ChargedStable muon
muon particle
Definition: Const.h:660
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
static const ChargedStable proton
proton particle
Definition: Const.h:663
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:662
static const ChargedStable electron
electron particle
Definition: Const.h:659
void clearVars()
function to clear the variables
Definition: HadronBgPrep.h:138
void defineHisto(std::vector< TH1F * > &htemp, const std::string &svar, const std::string &stype, const std::string &pdg)
function to define histograms
std::vector< double > m_sumcos
variables to add cos values
Definition: HadronBgPrep.h:174
std::vector< double > m_suminj
variables to add injection time
Definition: HadronBgPrep.h:177
int m_nhitBins
bins for nhits
Definition: HadronBgPrep.h:196
int m_injBins
bins for injection time
Definition: HadronBgPrep.h:188
int m_cosBins
bins for cosine
Definition: HadronBgPrep.h:192
double getParticleMass(const std::string &particle)
function to get the particle mass
Definition: HadronBgPrep.h:160
void deleteHistos(std::vector< TH1F * > &htemp)
function to delete the histograms
Definition: HadronBgPrep.h:152
double m_cosMax
max range of cosine
Definition: HadronBgPrep.h:194
void printCanvasCos(std::map< int, std::vector< TH1F * > > &hchicos_allbg, std::map< int, std::vector< TH1F * > > &hchicos_1by3bg, std::map< int, std::vector< TH1F * > > &hchicos_2by3bg, std::map< int, std::vector< TH1F * > > &hchicos_3by3bg, const std::string &particle, const std::string &suffix)
function to draw the dedx vs costh histograms
int m_bgBins
bins for bg
Definition: HadronBgPrep.h:184
void fit(TH1F *&hist, const std::string &pdg)
function to fit the histograms
std::vector< double > m_means
mean variable
Definition: HadronBgPrep.h:178
std::vector< double > m_errors
error variable
Definition: HadronBgPrep.h:179
std::vector< int > m_sumsize
size of the bg bins
Definition: HadronBgPrep.h:181
void plotDist(std::map< int, std::vector< TH1F * > > &hist, const std::string &suffix, int bins)
function to plot the map of histograms
HadronBgPrep()
Constructor: Sets the description, the properties and the parameters of the algorithm.
Definition: HadronBgPrep.cc:12
void FormatGraph(TGraphErrors &gr, int flag, const std::string &name="")
function to set graph cosmetics
Definition: HadronBgPrep.h:110
double m_injMax
max range of injection time
Definition: HadronBgPrep.h:190
double m_nhitMax
max range of nhits
Definition: HadronBgPrep.h:198
std::vector< double > m_sumres_square
variables to add square of resolution
Definition: HadronBgPrep.h:176
double m_injMin
min range of injection time
Definition: HadronBgPrep.h:189
void setPars(TFile *&outfile, std::string pdg, std::vector< TH1F * > &hdedx_bg, std::vector< TH1F * > &hchi_bg, std::vector< TH1F * > &hionzsigma_bg, std::map< int, std::vector< TH1F * > > &hchi_inj)
function to fill the parameters like mean and reso in the tree
double m_bgMax
max range of bg
Definition: HadronBgPrep.h:186
std::vector< int > m_injsize
size of the injection bins
Definition: HadronBgPrep.h:182
double m_bgMin
min range of bg
Definition: HadronBgPrep.h:185
double m_cosMin
min range of cosine
Definition: HadronBgPrep.h:193
double m_cut
cut to clean protons
Definition: HadronBgPrep.h:200
std::vector< double > m_sumbg
variables to add bg values
Definition: HadronBgPrep.h:175
void fitGaussianWRange(TH1F *&temphist, gstatus &status, double sigmaR)
function to perform gauss fit for input histogram
double m_nhitMin
min range of nhits
Definition: HadronBgPrep.h:197
void prepareSample(std::shared_ptr< TTree > hadron, TFile *&outfile, const std::string &suffix, const std::string &bgcurvefile, const std::string &bgsigmafile, const std::string &pdg, bool ismakePlots)
function to prepare sample for monitoring plots, bg curve fitting and sigma vs ionz fitting
Definition: HadronBgPrep.cc:46
Abstract base class for different kinds of events.