Belle II Software release-09-00-04
CDCDedxValidationAlgorithm.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/CDCDedxValidationAlgorithm.h>
10#include <cdc/calibration/CDCdEdx/CDCDedxWireGainAlgorithm.h>
11
12#include <cdc/dbobjects/CDCDedxWireGain.h>
13#include <cdc/dbobjects/CDCDedxCosineCor.h>
14#include <cdc/dbobjects/CDCDedx1DCell.h>
15#include <cdc/dbobjects/CDCDedxRunGain.h>
16#include <cdc/dbobjects/CDCDedxBadWires.h>
17#include <framework/database/IntervalOfValidity.h>
18
19#include <framework/database/Database.h>
20#include <framework/database/DBStore.h>
21#include <framework/database/Configuration.h>
22
23#include <cmath>
24#include <TTree.h>
25#include <TMap.h>
26#include <TLegend.h>
27#include <TF1.h>
28
29
30using namespace Belle2;
31using namespace CDC;
32
33//-----------------------------------------------------------------
34// Implementation
35//-----------------------------------------------------------------
36
38 CalibrationAlgorithm("ElectronValCollector"),
39 m_sigmaR(2.0),
40 m_dedxBins(600),
41 m_dedxMin(0.0),
42 m_dedxMax(5.0),
43 m_cosBins(100),
44 m_cosMin(-1.0),
45 m_cosMax(1.0),
46 m_momBins(80),
47 m_momMin(0.0),
48 m_momMax(8.0),
49 m_eaBin(316),
50 m_eaMin(-TMath::Pi() / 2),
51 m_eaMax(+TMath::Pi() / 2),
52 m_suffix("")
53{
54 // Set module properties
55 setDescription("A validation algorithm for CDC dE/dx electron");
56}
57
58//-----------------------------------------------------------------
59// Run the calibration
60//-----------------------------------------------------------------
61
63{
64
66
67 std::vector<std::string> subdirs = {"run", "costh", "mom", "wire", "injection", "oneD"};
68 for (const auto& dir : subdirs) {
69 gSystem->Exec(Form("mkdir -p plots/%s", dir.c_str()));
70 }
71
72 // Get data objects
73 auto tBhabha = getObjectPtr<TTree>("tBhabha");
74
75 // require at least 100 tracks
76 if (tBhabha->GetEntries() < 100) return c_NotEnoughData;
77
78 // Get data objects
79 auto tRadee = getObjectPtr<TTree>("tRadee");
80
81 // require at least 100 tracks
82 if (tRadee->GetEntries() < 100) return c_NotEnoughData;
83
87 return c_OK;
88}
89
90//------------------------------------
92{
93
94 int cruns = 0;
95 for (auto expRun : getRunList()) {
96 if (cruns == 0) B2INFO("start exp " << expRun.first << " and run " << expRun.second << "");
97 cruns++;
98 }
99
100 const auto erStart = getRunList()[0];
101 int estart = erStart.first;
102 int rstart = erStart.second;
103
104 updateDBObjPtrs(1, rstart, estart);
105
106 if (m_suffix.length() > 0) m_suffix = Form("%s_e%d_r%d", m_suffix.data(), estart, rstart);
107 else m_suffix = Form("e%d_r%d", estart, rstart);
108}
109
110
112{
113 auto ttree = getObjectPtr<TTree>("tRadee");
114
115 double dedx, costh, p, injtime = 0.0, injring = 1.0;
116 int charge;
117
118 std::vector<double>* dedxhit = 0, *enta = 0;
119 std::vector<int>* layer = 0;
120
121 ttree->SetBranchAddress("dedx", &dedx);
122 ttree->SetBranchAddress("p", &p);
123 ttree->SetBranchAddress("costh", &costh);
124 ttree->SetBranchAddress("charge", &charge);
125 ttree->SetBranchAddress("injtime", &injtime);
126 ttree->SetBranchAddress("injring", &injring);
127 ttree->SetBranchAddress("dedxhit", &dedxhit);
128 ttree->SetBranchAddress("entaRS", &enta);
129 ttree->SetBranchAddress("layer", &layer);
130
131 std::vector<double> vtlocaledges;
132 defineTimeBins(vtlocaledges);
133 m_tbins = vtlocaledges.size() - 1;
134 m_tedges = &vtlocaledges[0];
135
136 std::array<std::array<std::vector<TH1D*>, 2>, 13> hdedx_mom;
137 std::array<std::vector<TH1D*>, 2> hdedx_mom_peaks, hdedx_inj, hdedx_oned;
138 TH1D* htimes = new TH1D("htimes", "", m_tbins, m_tedges);
139
140 const double momBinWidth = (m_momMax - m_momMin) / m_momBins;
141 const double momBinW = (4.0 - m_momMin) / 4;
142
143 std::string scos[13] = {"acos", "cos#theta > 0.0", "cos#theta < 0.0", "cos#theta <= -0.8",
144 "cos#theta > -0.8 and cos#theta <= -0.6",
145 "cos#theta > -0.6 and cos#theta <= -0.4", "cos#theta > -0.4 and cos#theta <= -0.2",
146 "cos#theta > -0.2 and cos#theta <= 0", "cos#theta > 0 and cos#theta <= 0.2",
147 "cos#theta > 0.2 and cos#theta <= 0.4", "cos#theta > 0.4 and cos#theta <= 0.6",
148 "cos#theta > 0.6 and cos#theta <= 0.8", "cos#theta > 0.8"
149 };
150 std::string stype[2] = {"posi", "elec"};
151 std::string sLayer[2] = {"IL", "OL"};
152
153 // Define histograms for momentum bins and charge types
154 for (int ic = 0; ic < 13; ic++) {
155 for (int it = 0; it < 2; ++it) {
156 hdedx_mom[ic][it].resize(m_momBins);
157 defineHisto(hdedx_mom[ic][it], "mom", Form("%d_%s", ic, stype[it].data()));
158 }
159 }
160
161 // Define histograms for injection time bins and rings
162 for (unsigned int ir = 0; ir < 2; ir++) {
163 hdedx_inj[ir].resize(m_tbins);
164 hdedx_mom_peaks[ir].resize(4);
165 hdedx_oned[ir].resize(m_eaBin);
166
167 defineHisto(hdedx_inj[ir], "inj", m_sring[ir].data());
168 defineHisto(hdedx_mom_peaks[ir], "mom_peaks", Form("%s", stype[ir].data()));
169 defineHisto(hdedx_oned[ir], "oned", Form("%s", sLayer[ir].data()));
170 }
171
172 double eaBW = (m_eaMax - m_eaMin) / m_eaBin;
173 double icos[3] = {0, -1, -1};
174 double chgtype;
175
176 // Loop over all the entries in the tree
177 for (int i = 0; i < ttree->GetEntries(); ++i) {
178 ttree->GetEvent(i);
179
180 // Skip invalid events
181 if (dedx <= 0 || injtime < 0 || injring < 0) continue;
182
183 // Calculate momentum bin index for hdedx_mom
184 int binIndex = static_cast<int>((abs(p) - m_momMin) / momBinWidth);
185
186 // Determine cos(theta) category
187
188 icos[1] = (costh > 0) ? 1 : 2;
189 icos[2] = int((costh + 1.0) / 0.2) + 3;
190 if (icos[2] < 3) icos[2] = 3;
191 if (icos[2] > 12) icos[2] = 12;
192
193 // Determine charge type
194 chgtype = (charge > 0) ? 0 : 1;
195
196 // Fill momentum histograms (only if binIndex is valid)
197 if (binIndex >= 0 && binIndex < m_momBins) {
198 hdedx_mom[icos[0]][chgtype][binIndex]->Fill(dedx);
199 hdedx_mom[icos[1]][chgtype][binIndex]->Fill(dedx);
200 hdedx_mom[icos[2]][chgtype][binIndex]->Fill(dedx);
201 }
202
203 // Add larger times to the last bin
204 if (injtime > m_tedges[m_tbins]) injtime = m_tedges[m_tbins] - 10.0;
205
206 // Injection ring type
207 int wr = (injring > 0.5) ? 1 : 0;
208
209 // Injection time bin
210 unsigned int tb = htimes->GetXaxis()->FindBin(injtime);
211 tb = std::min(tb, static_cast<unsigned int>(m_tbins)) - 1;
212
213 // Fill injection time and dE/dx histograms
214 htimes->Fill(injtime);
215 hdedx_inj[wr][tb]->Fill(dedx);
216
217 // Fill hdedx_mom_peaks with its own binning
218 int binI = static_cast<int>((abs(p) - m_momMin) / momBinW);
219 if (binI >= 0 && binI < 4) {
220 hdedx_mom_peaks[chgtype][binI]->Fill(dedx);
221 }
222
223 // Fill dE/dx in enta bins from hits
224 for (unsigned int j = 0; j < dedxhit->size(); ++j) {
225 if (dedxhit->at(j) == 0) continue;
226
227 double entaval = enta->at(j);
228 int ibin = std::floor((entaval - m_eaMin) / eaBW);
229 if (ibin < 0 || ibin >= m_eaBin) continue;
230
231 int mL = (layer->at(j) < 8) ? 0 : 1;
232 hdedx_oned[mL][ibin]->Fill(dedxhit->at(j));
233 }
234 }
235
236
237 for (int ic = 0; ic < 13; ic++) {
238 for (int it = 0; it < 2; ++it) {
239 printCanvas(hdedx_mom[ic][it], Form("plots/mom/dedx_vs_mom_%d_%s_%s", ic, stype[it].data(), m_suffix.data()), "mom");
240 }
241 }
242 for (int it = 0; it < 2; ++it) {
243 printCanvas(hdedx_inj[it], Form("plots/injection/dedx_vs_inj_%s_%s", m_sring[it].data(), m_suffix.data()), "inj");
244 printCanvas(hdedx_oned[it], Form("plots/oneD/dedx_vs_1D_%s_%s", sLayer[it].data(), m_suffix.data()), "oned");
245 }
246
247 printCanvasdEdx(hdedx_mom_peaks, Form("plots/mom/dedxpeaks_vs_mom_%s", m_suffix.data()), "mom");
248
249}
250
252{
253 auto ttree = getObjectPtr<TTree>("tBhabha");
254
255 double dedx, costh;
256 int run, charge;
257
258 std::vector<int>* wire = 0;
259 ttree->SetBranchAddress("wire", &wire);
260
261 std::vector<double>* dedxhit = 0;
262 ttree->SetBranchAddress("dedxhit", &dedxhit);
263
264 ttree->SetBranchAddress("dedx", &dedx);
265 ttree->SetBranchAddress("run", &run);
266 ttree->SetBranchAddress("charge", &charge);
267 ttree->SetBranchAddress("costh", &costh);
268
269 std::map<int, TH1D*> hdedx_run;
270 std::array<std::vector<TH1D*>, 3> hdedx_cos;
271 std::array<std::vector<TH1D*>, 2> hdedx_cos_peaks;
272 std::vector<TH1D*> hdedxhit(c_nSenseWires);
273
274 const double cosBinWidth = (m_cosMax - m_cosMin) / m_cosBins;
275 const double cosBinW = (m_cosMax - m_cosMin) / 4;
276
277 std::string stype[3] = {"all", "posi", "elec"};
278
279 for (int it = 0; it < 3; ++it) {
280 hdedx_cos[it].resize(m_cosBins);
281 defineHisto(hdedx_cos[it], "costh", stype[it]);
282 }
283
284 for (int ir = 0; ir < 2; ir++) {
285 hdedx_cos_peaks[ir].resize(4);
286 defineHisto(hdedx_cos_peaks[ir], "cos_peaks", Form("%s", stype[ir + 1].data()));
287 }
288
289 defineHisto(hdedxhit, "wire", "wire");
290
291 // Loop over all the entries in the tree
292 for (int i = 0; i < ttree->GetEntries(); ++i) {
293 ttree->GetEvent(i);
294 if (dedx <= 0) continue;
295
296 // Check if a dE/dx histogram for this run number already exists
297 if (hdedx_run.find(run) == hdedx_run.end()) {
298 std::string histName = Form("hist_dedx_run_%d", run);
299 std::string histTitle = Form("dE/dx Histogram for Run %d", run);
300 hdedx_run[run] = new TH1D(histName.data(), histTitle.data(), m_dedxBins, m_dedxMin, m_dedxMax);
301 }
302
303 // Fill run-specific histogram
304 hdedx_run[run]->Fill(dedx);
305
306 // Fill cos(theta) histograms (all charge + by charge sign)
307 int binIndex = static_cast<int>((costh - m_cosMin) / cosBinWidth);
308 if (binIndex >= 0 && binIndex < m_cosBins) {
309 hdedx_cos[0][binIndex]->Fill(dedx); // All charge
310
311 if (charge > 0)
312 hdedx_cos[1][binIndex]->Fill(dedx);
313 else if (charge < 0)
314 hdedx_cos[2][binIndex]->Fill(dedx);
315 }
316
317 // Fill dE/dx for each wire hit
318 for (unsigned int j = 0; j < wire->size(); ++j) {
319 int jwire = wire->at(j);
320 double jhitdedx = dedxhit->at(j);
321 hdedxhit[jwire]->Fill(jhitdedx);
322 }
323
324 // Fill cos(theta) peaks histograms
325 int binI = static_cast<int>((costh - m_cosMin) / cosBinW);
326 if (binI >= 0 && binI < 4) {
327 if (charge > 0)
328 hdedx_cos_peaks[0][binI]->Fill(dedx);
329 else if (charge < 0)
330 hdedx_cos_peaks[1][binI]->Fill(dedx);
331 }
332 }
333
334 printCanvasRun(hdedx_run, Form("plots/run/dedx_vs_run_%s", m_suffix.data()));
335 printCanvas(hdedx_cos[0], Form("plots/costh/dedx_vs_cos_all_%s", m_suffix.data()), "costh");
336 printCanvas(hdedx_cos[1], Form("plots/costh/dedx_vs_cos_positrons_%s", m_suffix.data()), "costh");
337 printCanvas(hdedx_cos[2], Form("plots/costh/dedx_vs_cos_electrons_%s", m_suffix.data()), "costh");
338 printCanvasdEdx(hdedx_cos_peaks, Form("plots/costh/dedxpeaks_vs_cos_%s", m_suffix.data()), "costh");
339 wireGain(hdedxhit);
340}
341
342//------------------------------------
343void CDCDedxValidationAlgorithm::defineHisto(std::vector<TH1D*>& htemp, std::string var, std::string stype)
344{
345 int xbins;
346 double xmin, xmax;
347 double binWidth = 0.0;
348
349 if (var == "mom") {
350 xbins = m_momBins; xmin = m_momMin; xmax = m_momMax;
351 } else if (var == "oned") {
352 xbins = m_eaBin; xmin = m_eaMin; xmax = m_eaMax; m_dedxBins = 250;
353 } else if (var == "costh") {
354 xbins = m_cosBins; xmin = m_cosMin; xmax = m_cosMax;
355 } else if (var == "inj") {
356 xbins = m_tbins;
357 } else if (var == "mom_peaks") {
358 xbins = 4; xmin = m_momMin; xmax = 4.0;
359 } else if (var == "cos_peaks") {
360 xbins = 4; xmin = m_cosMin; xmax = m_cosMax;
361 } else {
362 xbins = c_nSenseWires; m_dedxBins = 250;
363 }
364
365 if (var == "costh" || var == "mom" || var == "mom_peaks" || var == "cos_peaks" || var == "oned") {
366 binWidth = (xmax - xmin) / xbins;
367 }
368
369 for (int ic = 0; ic < xbins; ic++) {
370 std::string title = Form("dedxhit-dist, wire:%d", ic);
371 std::string name = Form("hdedx_%s_%d", var.data(), ic);
372
373 if (var == "costh" || var == "mom" || var == "mom_peaks" || var == "cos_peaks" || var == "oned") {
374 double min = ic * binWidth + xmin;
375 double max = min + binWidth;
376 title = Form("%s: (%0.02f, %0.02f) %s", var.data(), min, max, stype.data());
377 name = Form("hdedx_%s_%s_%d", var.data(), stype.data(), ic);
378 } else if (var == "inj") {
379 std::string label = getTimeBinLabel(m_tedges[ic], ic);
380 title = Form("%s, time(%s)", stype.data(), label.data());
381 name = Form("h%s_%s_t%d", var.data(), stype.data(), ic);
382 }
383 htemp[ic] = new TH1D(name.data(), "", m_dedxBins, m_dedxMin, m_dedxMax);
384 htemp[ic]->SetTitle(Form("%s;dedx;entries", title.data()));
385 }
386}
387
388void CDCDedxValidationAlgorithm::printCanvasdEdx(std::array<std::vector<TH1D*>, 2>& htemp, std::string namesfx, std::string svar)
389{
390 int xbins = 4;
391 double xmin, xmax;
392
393 if (svar == "mom") {
394 xmin = m_momMin; xmax = 4.0;
395 } else if (svar == "costh") {
396 xmin = m_cosMin; xmax = m_cosMax;
397 } else {
398 B2FATAL("wrong input");
399 }
400 double binWidth = (xmax - xmin) / xbins;
401
402 // Set up the TCanvas with 2x2 grid
403 TCanvas* ctmp = new TCanvas("tmp", "tmp", 1200, 1200);
404 ctmp->Divide(2, 2); // Divide into 2x2 grid
405
406 // Prepare the PDF output
407 std::stringstream psname;
408 psname << Form("%s.pdf[", namesfx.data());
409 ctmp->Print(psname.str().c_str());
410 psname.str("");
411 psname << Form("%s.pdf", namesfx.data());
412
413 // Iterate through the histograms and plot them in the canvas
414 for (int i = 0; i < xbins; ++i) {
415
416 ctmp->cd(i % 4 + 1);
417
418 double emean, emeanErr, esigma, esigmaErr;
419 double pmean, pmeanErr, psigma, psigmaErr;
420
421 fit(htemp[0][i], emean, emeanErr, esigma, esigmaErr);
422 fit(htemp[1][i], pmean, pmeanErr, psigma, psigmaErr);
423
424 double min = i * binWidth + xmin;
425 double max = min + binWidth;
426
427 TPaveText pt(0.6, 0.63, 0.85, 0.89, "NBNDC");
428 setTextCosmetics(pt, kBlack);
429 pt.AddText("e+");
430 pt.AddText(Form("#mu_{fit}: %0.03f#pm%0.03f", emean, emeanErr));
431 pt.AddText(Form("#sigma_{fit}: %0.03f#pm%0.03f", esigma, esigmaErr));
432
433 pt.AddText("e-");
434 pt.AddText(Form("#mu_{fit}: %0.03f#pm%0.03f", pmean, pmeanErr));
435 pt.AddText(Form("#sigma_{fit}: %0.03f#pm%0.03f", psigma, psigmaErr));
436
437 htemp[0][i]->SetStats(0);
438 htemp[1][i]->SetStats(0);
439 htemp[0][i]->SetFillColor(0);
440 htemp[1][i]->SetFillColor(0);
441 htemp[0][i]->SetLineColor(8);
442 htemp[1][i]->SetLineColor(9);
443 htemp[0][i]->SetTitle(Form("%s: (%0.02f, %0.02f)", svar.data(), min, max));
444 if (htemp[0][i]->GetEntries() > 0)
445 htemp[0][i]->Scale(1.0 / htemp[0][i]->GetEntries());
446 if (htemp[1][i]->GetEntries() > 0)
447 htemp[1][i]->Scale(1.0 / htemp[1][i]->GetEntries());
448
449 if (htemp[1][i]->GetMaximum() > htemp[0][i]->GetMaximum())
450 htemp[0][i]->SetMaximum(htemp[1][i]->GetMaximum());
451
452 htemp[0][i]->DrawCopy("HIST");
453 htemp[1][i]->DrawCopy("same HIST");
454 pt.DrawClone("same");
455
456 TLegend* lego = new TLegend(0.15, 0.67, 0.3, 0.8);
457 lego->AddEntry(htemp[0][i], "e+", "l");
458 lego->AddEntry(htemp[1][i], "e-", "l");
459 lego->Draw("same");
460
461 if ((i + 1) % 4 == 0 || i == xbins - 1) {
462 ctmp->SetBatch(kTRUE);
463 ctmp->Print(psname.str().c_str());
464 if ((i + 1) % 4 == 0) ctmp->Clear("D");
465 }
466 }
467
468 psname.str("");
469 psname << Form("%s.pdf]", namesfx.data());
470 ctmp->Print(psname.str().c_str());
471
472 delete ctmp;
473}
474
475void CDCDedxValidationAlgorithm::printCanvas(std::vector<TH1D*>& htemp, std::string namesfx, std::string svar)
476{
477 int xbins;
478 double xmin, xmax;
479
480 if (svar == "mom") {
481 xbins = m_momBins; xmin = m_momMin; xmax = m_momMax;
482 } else if (svar == "oned") {
483 xbins = m_eaBin; xmin = m_eaMin; xmax = m_eaMax;
484 } else if (svar == "costh") {
485 xbins = m_cosBins; xmin = m_cosMin; xmax = m_cosMax;
486 } else if (svar == "inj") {
487 xbins = m_tbins;
488 } else if (svar == "mom_peaks") {
489 xbins = 4; xmin = m_momMin; xmax = 4.0;
490 } else {
491 B2FATAL("wrong input");
492 }
493
494 // Set up the TCanvas with 4x4 grid
495 TCanvas* ctmp = new TCanvas("tmp", "tmp", 1200, 1200);
496 ctmp->Divide(4, 4); // Divide into 4x4 grid
497
498 // Prepare the PDF output
499 std::stringstream psname;
500 psname << Form("%s.pdf[", namesfx.data());
501 ctmp->Print(psname.str().c_str());
502 psname.str("");
503 psname << Form("%s.pdf", namesfx.data());
504
505 std::ofstream outFile;
506 outFile.open(Form("%s.txt", namesfx.data()));
508
509 // Iterate through the histograms and plot them in the canvas
510 for (int i = 0; i < xbins; ++i) {
511
512 ctmp->cd(i % 16 + 1);
513 TPaveText pt(0.6, 0.73, 0.85, 0.89, "NBNDC");
514 setTextCosmetics(pt, kBlack);
515
516 if (svar == "oned") {
517 unsigned int minbin, maxbin;
518 wireg.getTruncatedBins(htemp[i], minbin, maxbin);
519 htemp[i]->SetTitle(Form("dedxhit-dist, entabin: %d ;%d;%d", i, minbin, maxbin));
520
521 double dedxmean = wireg.getTruncationMean(htemp[i], minbin, maxbin);
522
523 const double binWidth = (xmax - xmin) / xbins;
524 double binCenter = xmin + (i + 0.5) * binWidth; // Calculate bin center for cos(theta) or mom
525
526 outFile << binCenter << " " << dedxmean << std::endl;
527 } else {
528 double mean, meanErr, sigma, sigmaErr;
529 fit(htemp[i], mean, meanErr, sigma, sigmaErr);
530
531 if (svar == "mom" || svar == "costh" || svar == "mom_peaks") {
532 const double binWidth = (xmax - xmin) / xbins;
533 double binCenter = xmin + (i + 0.5) * binWidth; // Calculate bin center for cos(theta) or mom
534
535 outFile << binCenter << " " << mean << " " << meanErr << " " << sigma << " " << sigmaErr << std::endl;
536 } else {
537 std::string label = getTimeBinLabel(m_tedges[i], i);
538 outFile << i << " " << label << " " << mean << " " << meanErr << " " << sigma << " " << sigmaErr << std::endl;
539 }
540
541 pt.AddText(Form("#mu_{fit}: %0.03f#pm%0.03f", mean, meanErr));
542 pt.AddText(Form("#sigma_{fit}: %0.03f#pm%0.03f", sigma, sigmaErr));
543 }
544 htemp[i]->SetStats(0);
545 htemp[i]->DrawCopy("");
546 pt.DrawClone("same");
547
548 if ((i + 1) % 16 == 0 || ((i + 1) == xbins)) {
549 ctmp->SetBatch(kTRUE);
550 ctmp->Print(psname.str().c_str());
551 ctmp->Clear("D");
552 }
553 }
554
555 // ctmp->Print(psname.str().c_str());
556 psname.str("");
557 psname << Form("%s.pdf]", namesfx.data());
558 ctmp->Print(psname.str().c_str());
559
560 outFile.close();
561
562 delete ctmp;
563}
564
565void CDCDedxValidationAlgorithm::fit(TH1D*& hist, double& mean, double& meanErr, double& sigma, double& sigmaErr)
566{
567
568 std::string status = "";
569
570 if (hist->Integral() > 100)
571 fitGaussianWRange(hist, status);
572
573 if (status != "fitOK") {
574 hist->SetFillColor(kOrange);
575 mean = 0.0, meanErr = 0.0, sigma = 0.0, sigmaErr = 0.0;
576 } else {
577 mean = hist->GetFunction("gaus")->GetParameter(1);
578 meanErr = hist->GetFunction("gaus")->GetParError(1);
579 sigma = hist->GetFunction("gaus")->GetParameter(2);
580 sigmaErr = hist->GetFunction("gaus")->GetParError(2);
581 hist->SetFillColor(kYellow);
582 }
583}
584
585void CDCDedxValidationAlgorithm::printCanvasRun(std::map<int, TH1D*>& htemp, std::string namesfx)
586{
587 // Set up the TCanvas with 4x4 grid
588 TCanvas* ctmp = new TCanvas("tmp", "tmp", 1200, 1200);
589 ctmp->Divide(4, 4); // Divide into 4x4 grid
590
591 // Prepare the PDF output
592 std::stringstream psname;
593 psname << Form("%s.pdf[", namesfx.data());
594 ctmp->Print(psname.str().c_str());
595 psname.str("");
596 psname << Form("%s.pdf", namesfx.data());
597
598 std::ofstream outFile;
599 outFile.open(Form("%s.txt", namesfx.data()));
600
601 // Iterate through the histograms and plot them in the canvas
602 int irun = 0;
603 for (const auto& entry : htemp) {
604 int run = entry.first;
605 TH1D* hist = entry.second;
606
607 ctmp->cd(irun % 16 + 1);
608
609 TPaveText pt(0.6, 0.73, 0.85, 0.89, "NBNDC");
610 setTextCosmetics(pt, kBlack);
611
612 double mean, meanErr, sigma, sigmaErr;
613 fit(hist, mean, meanErr, sigma, sigmaErr);
614
615 outFile << run << " " << mean << " " << meanErr << " " << sigma << " " << sigmaErr << std::endl;
616
617 pt.AddText(Form("#mu_{fit}: %0.03f#pm%0.03f", mean, meanErr));
618 pt.AddText(Form("#sigma_{fit}: %0.03f#pm%0.03f", sigma, sigmaErr));
619
620 hist->SetStats(0);
621 hist->DrawCopy("");
622 pt.DrawClone("same");
623
624 if ((irun + 1) % 16 == 0 || irun == int(htemp.size() - 1)) {
625 ctmp->SetBatch(kTRUE);
626 ctmp->Print(psname.str().c_str());
627 ctmp->Clear("D");
628 }
629 irun++;
630 }
631
632 ctmp->Print(psname.str().c_str());
633 psname.str("");
634 psname << Form("%s.pdf]", namesfx.data());
635 ctmp->Print(psname.str().c_str());
636
637 outFile.close();
638
639 delete ctmp;
640}
641
642//----------------------------------------
643void CDCDedxValidationAlgorithm::fitGaussianWRange(TH1D*& temphist, std::string& status)
644{
645 double histmean = temphist->GetMean();
646 double histrms = temphist->GetRMS();
647 temphist->GetXaxis()->SetRangeUser(histmean - 5.0 * histrms, histmean + 5.0 * histrms);
648
649 int fs = temphist->Fit("gaus", "Q0");
650 if (fs != 0) {
651 B2INFO(Form("\tFit (round 1) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
652 status = "fitFailed";
653 return;
654 } else {
655 double mean = temphist->GetFunction("gaus")->GetParameter(1);
656 double width = temphist->GetFunction("gaus")->GetParameter(2);
657 temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
658 fs = temphist->Fit("gaus", "QR", "", mean - m_sigmaR * width, mean + m_sigmaR * width);
659 if (fs != 0) {
660 B2INFO(Form("\tFit (round 2) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
661 status = "fitFailed";
662 return;
663 } else {
664 temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
665 B2INFO(Form("\tFit for hist (%s) successful (status = %d)", temphist->GetName(), fs));
666 status = "fitOK";
667 }
668 }
669}
670
671void CDCDedxValidationAlgorithm::wireGain(std::vector<TH1D*>& hdedxhit)
672{
673
674 std::vector<double> vdedx_means;
675 std::vector<double> layermean(c_maxNSenseLayers);
676 std::vector<double> lgmean(c_maxNSenseLayers);
677
678 std::ofstream outFile, outFileLayer, outFileAvg, outFilebdwire;
679 outFile.open(Form("plots/wire/dedx_mean_gwire_%s.txt", m_suffix.data()));
680 outFilebdwire.open(Form("plots/wire/dedx_mean_badwire_%s.txt", m_suffix.data()));
681 outFileLayer.open(Form("plots/wire/dedx_mean_layer_%s.txt", m_suffix.data()));
682 outFileAvg.open(Form("plots/wire/dedx_mean_layer_avg_%s.txt", m_suffix.data()));
683
684 int activelayers = 0;
685 double layeravg = 0.0;
686
689
691
692 int jwire = -1;
693 for (unsigned int il = 0; il < c_maxNSenseLayers; ++il) {
694
695 int activewires = 0, goodwires = 0;
696 layermean[il] = 0.0;
697 lgmean[il] = 0.0;
698
699 for (unsigned int iw = 0; iw < cdcgeo.nWiresInLayer(il); ++iw) {
700 jwire++;
701
702 unsigned int minbin, maxbin;
703 wireg.getTruncatedBins(hdedxhit[jwire], minbin, maxbin);
704 hdedxhit[jwire]->SetTitle(Form("dedxhit-dist, wire: %d ;%d;%d", jwire, minbin, maxbin));
705
706 double dedxmean = wireg.getTruncationMean(hdedxhit[jwire], minbin, maxbin);
707 vdedx_means.push_back(dedxmean);
708
709 if (Badwire->getBadWireStatus(jwire) == kTRUE)
710 outFilebdwire << jwire << " " << dedxmean << std::endl;
711 else
712 outFile << jwire << " " << dedxmean << std::endl;
713
714 if (vdedx_means.at(jwire) > 0) {
715 layermean[il] += vdedx_means.at(jwire);
716 activewires++;
717 if (Badwire->getBadWireStatus(jwire) != kTRUE) {
718 lgmean[il] += vdedx_means.at(jwire);
719 goodwires++;
720 }
721 }
722 }
723
724 if (activewires > 0) layermean[il] /= activewires;
725 else layermean[il] = 1.0;
726
727 if (goodwires > 0) lgmean[il] /= goodwires;
728 else lgmean[il] = 1.0;
729
730 outFileLayer << il << " " << layermean[il] << " " << lgmean[il] << std::endl;
731
732 //calculate outer layer average for active layer
733 if (il >= 8 && layermean[il] > 0) {
734 layeravg += layermean[il];
735 activelayers++;
736 }
737 }
738
739 //normalisation of wiregains to get outergain ~1.0
740 if (activelayers > 0) layeravg /= activelayers;
741 outFileAvg << layeravg << std::endl;
742
743 outFile.close();
744 outFilebdwire.close();
745 outFileLayer.close();
746 outFileAvg.close();
747 printCanvasWire(hdedxhit, Form("plots/wire/dedx_vs_wire_%s", m_suffix.data()), vdedx_means);
748}
749
750void CDCDedxValidationAlgorithm::printCanvasWire(std::vector<TH1D*> temp, std::string namesfx,
751 const std::vector<double>& vdedx_mean)
752{
753 TCanvas* ctmp = new TCanvas("tmp", "tmp", 900, 900);
754 ctmp->Divide(4, 4);
755
756 std::stringstream psname;
757 psname << Form("%s.pdf[", namesfx.data());
758 ctmp->Print(psname.str().c_str());
759 psname.str("");
760 psname << Form("%s.pdf", namesfx.data());
761
762 for (unsigned int ip = 0; ip < c_nwireCDC; ip++) {
763 int minbin = std::stoi(temp[ip]->GetXaxis()->GetTitle());
764 int maxbin = std::stoi(temp[ip]->GetYaxis()->GetTitle());
765 temp[ip]->SetFillColor(kYellow - 9);
766 temp[ip]->SetTitle(Form("%s, #mu_{trunc} %0.03f;dedxhit;entries", temp[ip]->GetTitle(), vdedx_mean.at(ip)));
767
768 ctmp->cd(ip % 16 + 1);
769 gPad->cd();
770 temp[ip]->DrawCopy("hist");
771 TH1D* hdedxhitC = (TH1D*)temp[ip]->Clone(Form("%sC", temp[ip]->GetName()));
772 hdedxhitC->GetXaxis()->SetRange(minbin, maxbin);
773 hdedxhitC->SetFillColor(kAzure + 1);
774 hdedxhitC->DrawCopy("same histo");
775
776 if ((ip + 1) % 16 == 0) {
777 ctmp->SetBatch(kTRUE);
778 ctmp->Print(psname.str().c_str());
779 ctmp->Clear("D");
780 }
781 delete temp[ip];
782 delete hdedxhitC;
783
784 }
785
786 psname.str("");
787 psname << Form("%s.pdf]", namesfx.data());
788 ctmp->Print(psname.str().c_str());
789 delete ctmp;
790}
791
792void CDCDedxValidationAlgorithm::defineTimeBins(std::vector<double>& vtlocaledges)
793{
794 double fixedges[69];
795 for (int ib = 0; ib < 69; ib++) {
796 fixedges[ib] = ib * 0.5 * 1e3;
797 if (ib > 40 && ib <= 60) fixedges[ib] = fixedges[ib - 1] + 1.0 * 1e3;
798 else if (ib > 60 && ib <= 64) fixedges[ib] = fixedges[ib - 1] + 10.0 * 1e3;
799 else if (ib > 64 && ib <= 65) fixedges[ib] = fixedges[ib - 1] + 420.0 * 1e3;
800 else if (ib > 65 && ib <= 66) fixedges[ib] = fixedges[ib - 1] + 500.0 * 1e3;
801 else if (ib > 66) fixedges[ib] = fixedges[ib - 1] + 2e6;
802 vtlocaledges.push_back(fixedges[ib]);
803 }
804}
805
806//--------------------------
808{
809
810 TCanvas cstats("cstats", "cstats", 800, 400);
811 cstats.SetBatch(kTRUE);
812 cstats.Divide(2, 1);
813
814 cstats.cd(1);
815 auto hestats = getObjectPtr<TH1I>("hestats");
816 if (hestats) {
817 hestats->SetName(Form("hestats_%s", m_suffix.data()));
818 hestats->SetStats(0);
819 hestats->DrawCopy("");
820 }
821
822 cstats.cd(2);
823 auto htstats = getObjectPtr<TH1I>("htstats");
824 if (htstats) {
825 htstats->SetName(Form("htstats_%s", m_suffix.data()));
826 htstats->SetStats(0);
827 htstats->DrawCopy("");
828 }
829
830 cstats.Print(Form("cdcdedx_stats_%s.pdf", m_suffix.data()));
831}
832
833void CDCDedxValidationAlgorithm::DatabaseIN(int experiment, int run)
834{
835 if (m_EventMetaData.isValid()) {
836 m_EventMetaData->setExperiment(experiment);
837 m_EventMetaData->setRun(run);
838 }
839
840 auto& dbConfiguration = Conditions::Configuration::getInstance();
841 dbConfiguration.setGlobalTags({"online"});
842 if (!m_testingPayloadName.empty() && m_GlobalTagName.empty()) {
843 dbConfiguration.prependTestingPayloadLocation(m_testingPayloadName);
844 } else if (m_testingPayloadName.empty() && !m_GlobalTagName.empty()) {
845 dbConfiguration.prependGlobalTag(m_GlobalTagName);
846 } else
847 B2FATAL("Setting both testing payload and Global Tag or setting no one of them.");
848
849 /* Mimic a module initialization. */
851 m_EventMetaData.registerInDataStore();
853 if (!m_EventMetaData.isValid())
854 m_EventMetaData.construct(1, run, experiment);
855
856 /* Database instance and configuration. */
857 DBStore& dbStore = DBStore::Instance();
858 dbStore.update();
859 dbStore.updateEvent();
860}
861
863{
864
865 DatabaseIN(experiment, run);
866
867 std::vector<double> wiregain;
868 std::vector<double> layermean(c_maxNSenseLayers);
869
870 DBObjPtr<CDCDedxWireGain> DBWireGains;
871 if (!DBWireGains.isValid()) B2FATAL("Wire gain data are not valid.");
872
874
875 int jwire = -1;
876 for (unsigned int il = 0; il < c_maxNSenseLayers; ++il) {
877
878 int activewires = 0;
879 layermean[il] = 0.0;
880
881 for (unsigned int iw = 0; iw < cdcgeo.nWiresInLayer(il); ++iw) {
882 jwire++;
883
884 wiregain.push_back(DBWireGains->getWireGain(jwire));
885
886 if (wiregain.at(jwire) > 0) {
887 layermean[il] += wiregain.at(jwire);
888 activewires++;
889 }
890 }
891
892 if (activewires > 0) layermean[il] /= activewires;
893 else layermean[il] = 1.0;
894 }
895
897 return { wiregain, layermean };
898}
899
901{
902
903 DatabaseIN(experiment, run);
904
905 std::vector<double> cosgain, cos;
906
907 DBObjPtr<CDCDedxCosineCor> DBCosineCor;
908 if (!DBCosineCor.isValid()) B2FATAL("Cosine gain data are not valid.");
909
910 unsigned int nCosBins = DBCosineCor->getSize();
911
912 for (unsigned int il = 0; il < nCosBins; ++il) {
913 double costh = -1.0 + (il + 0.5) * 2.0 / nCosBins;
914 costh += .000001;
915 cosgain.push_back(DBCosineCor->getMean(il));
916 cos.push_back(costh);
917 }
918
920 return {cosgain, cos};
921}
922
924{
925
926 DatabaseIN(experiment, run);
927
928 std::vector<double> inner1D, outer1D, Enta;
929
930 DBObjPtr<CDCDedx1DCell> DBOneDCell;
931 if (!DBOneDCell.isValid()) B2FATAL("OneD cell gain data are not valid.");
932
933 for (int i = 0; i < 2; i++) {
934
935 unsigned int nBins = DBOneDCell->getNBins(i);
936 double binSize = TMath::Pi() / nBins;
937
938 for (unsigned int nbin = 0; nbin < nBins; nbin++) {
939
940 double enta = (-1 * TMath::Pi() / 2.0) + binSize * nbin;
941 if (i == 0) {
942 Enta.push_back(enta);
943 inner1D.push_back(DBOneDCell->getMean(0, nbin));
944 } else
945 outer1D.push_back(DBOneDCell->getMean(17, nbin));
946 }
947 }
949 return {inner1D, outer1D, Enta};
950
951}
952
953double CDCDedxValidationAlgorithm::getrungain(int experiment, int run)
954{
955
956 DatabaseIN(experiment, run);
957
959 if (!RunGain.isValid()) B2FATAL("Run gain data are not valid.");
960 double gain = RunGain->getRunGain();
961 return gain;
962}
963
965{
966 /* Reset both DataStore and Database. */
968 Database::Instance().reset(false);
969 DBStore::Instance().reset(false);
970}
void bhabhaValidation()
Validate dE/dx using bhabha sample (vs run, cosine)
double m_eaMax
upper edge of entrance angle
void resetDatabase()
Clear current DB pointers and state.
void wireGain(std::vector< TH1D * > &hdedxhit)
Validate wire gain data using dE/dx histograms.
void printCanvasdEdx(std::array< std::vector< TH1D * >, 2 > &htemp, std::string namesfx, std::string svar)
Draw dE/dx histograms for momentum and cosine bins.
void printCanvasRun(std::map< int, TH1D * > &htemp, std::string namesfx)
Draw dE/dx per run histogram canvas.
void radeeValidation()
Validate dE/dx using radee sample (vs momentum, injection time)
double m_sigmaR
fit dedx dist in sigma range
void setTextCosmetics(TPaveText pt, Color_t color)
Set text cosmetics for TPaveText.
WireGainData getwiregain(int experiment, int run)
Retrieve wire gain data from DB.
void getExpRunInfo()
function to get extract calibration run/exp
OnedData getonedgain(int experiment, int run)
Retrieve 1D gain data from DB.
double * m_tedges
internal time array (copy of vtlocaledges)
void DatabaseIN(int experiment, int run)
Load database payload for given run.
void defineTimeBins(std::vector< double > &vtlocaledges)
Set bin edges for injection time.
std::array< std::string, 2 > m_sring
injection ring name
CosGainData getcosgain(int experiment, int run)
Retrieve cosine gain data from DB.
void printCanvasWire(std::vector< TH1D * > temp, std::string namesfx, const std::vector< double > &vdedx_mean)
Plot dE/dx vs wire number.
std::string m_testingPayloadName
Testing payload location.
std::string m_suffix
suffix string to separate plots
DBObjPtr< CDCGeometry > m_cdcGeo
Geometry of CDC.
CDCDedxValidationAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
double getrungain(int experiment, int run)
Retrieve run gain data from DB.
void fitGaussianWRange(TH1D *&temphist, std::string &status)
Perform Gaussian fit with range on a histogram.
void plotEventStats()
Plot summary statistics of selected events.
virtual EResult calibrate() override
Main calibration method.
void printCanvas(std::vector< TH1D * > &htemp, std::string namesfx, std::string svar)
Draw dE/dx histograms across bins.
StoreObjPtr< EventMetaData > m_EventMetaData
Event metadata.
unsigned int m_tbins
internal time bins
void defineHisto(std::vector< TH1D * > &htemp, std::string var, std::string stype)
Define dE/dx histograms for plotting.
std::string getTimeBinLabel(const double &tedges, const int &it)
Get time bin label string.
void fit(TH1D *&hist, double &mean, double &meanErr, double &sigma, double &sigmaErr)
Perform full Gaussian fit and extract parameters.
double m_eaMin
lower edge of entrance angle
A calibration algorithm for CDC dE/dx wire gains.
void getTruncatedBins(TH1D *hdedxhit, unsigned int &binlow, unsigned int &binhigh)
function to get bins of trunction from histogram
double getTruncationMean(TH1D *hdedxhit, int binlow, int binhigh)
function to get mean of trunction from histogram
The Class for CDC Geometry Parameters.
unsigned nWiresInLayer(int layerId) const
Returns wire numbers in a layer.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
Base class for calibration algorithms.
void updateDBObjPtrs(const unsigned int event, const int run, const int experiment)
Updates any DBObjPtrs by calling update(event) for DBStore.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
static Configuration & getInstance()
Get a reference to the instance which will be used when the Database is initialized.
bool isValid() const
Check whether a valid object was obtained from the database.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Singleton class to cache database objects.
Definition: DBStore.h:31
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:54
void setInitializeActive(bool active)
Setter for m_initializeActive.
Definition: DataStore.cc:94
void reset(EDurability durability)
Frees memory occupied by data store items and removes all objects from the map.
Definition: DataStore.cc:86
void reset(bool keepEntries=false)
Invalidate all payloads.
Definition: DBStore.cc:177
static Database & Instance()
Instance of a singleton Database.
Definition: Database.cc:42
static DBStore & Instance()
Instance of a singleton DBStore.
Definition: DBStore.cc:28
void updateEvent()
Updates all intra-run dependent objects.
Definition: DBStore.cc:142
void update()
Updates all objects that are outside their interval of validity.
Definition: DBStore.cc:79
static void reset(bool keepConfig=false)
Reset the database instance.
Definition: Database.cc:50
Abstract base class for different kinds of events.
Container for cosine gain data.
Container for 1D gain data.
Container for wire gain data.