Belle II Software  release-08-01-10
TOPLocalCalFitter.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 <top/calibration/TOPLocalCalFitter.h>
10 
11 #include "TH2F.h"
12 #include "TFile.h"
13 #include "TMath.h"
14 #include "TF1.h"
15 #include "TROOT.h"
16 #include "TTree.h"
17 #include "math.h"
18 #include <limits>
19 #include <framework/logging/Logger.h>
20 #include <top/dbobjects/TOPCalChannelT0.h>
21 
22 
23 //using namespace std;
24 using namespace Belle2;
25 using namespace TOP;
26 
27 double crystalball_function(double x, double alpha, double n, double sigma, double mean)
28 {
29  // evaluate the crystal ball function
30  if (sigma < 0.) return 0.;
31  double z = (x - mean) / sigma;
32  if (alpha < 0) z = -z;
33  double abs_alpha = std::abs(alpha);
34  if (z > - abs_alpha)
35  return std::exp(- 0.5 * z * z);
36  else {
37  double nDivAlpha = n / abs_alpha;
38  double AA = std::exp(-0.5 * abs_alpha * abs_alpha);
39  double B = nDivAlpha - abs_alpha;
40  double arg = nDivAlpha / (B - z);
41  return AA * std::pow(arg, n);
42  }
43 }
44 
45 double crystalball_pdf(double x, double alpha, double n, double sigma, double mean)
46 {
47  // evaluation of the PDF ( is defined only for n >1)
48  if (sigma < 0.) return 0.;
49  if (n <= 1) return std::numeric_limits<double>::quiet_NaN(); // pdf is not normalized for n <=1
50  double abs_alpha = std::abs(alpha);
51  double C = n / abs_alpha * 1. / (n - 1.) * std::exp(-alpha * alpha / 2.);
52  double D = std::sqrt(M_PI / 2.) * (1. + erf(abs_alpha / std::sqrt(2.)));
53  double N = 1. / (sigma * (C + D));
54  return N * crystalball_function(x, alpha, n, sigma, mean);
55 }
56 
57 
58 
59 // Two gaussians to model the TTS of the MCP-PMTs.
60 double TTSPDF(double x, double time, double deltaT, double sigma1, double sigma2, double f1)
61 {
62  return f1 * TMath::Gaus(x, time, sigma1, kTRUE) + (1 - f1) * TMath::Gaus(x, time + deltaT, sigma2, kTRUE) ;
63 }
64 
65 
66 // Full PDF to fit the laser, made of:
67 // 2 TTSPDF
68 // 1 crystal ball PDF for the extra peak at +1 ns we don't understand
69 // 1 gaussian to help modelling the tail
70 // cppcheck-suppress constParameter
71 double laserPDF(double* x, double* p)
72 {
73  double time = p[0];
74  double sigma = p[1];
75  double fraction = p[2];
76  double deltaTLaser = p[3];
77  double sigmaRatio = p[4];
78  double deltaTTS = p[5];
79  double f1 = p[6];
80  double norm = p[7];
81 
82  double deltaTExtra = p[8];
83  double sigmaExtra = p[9];
84  double normExtra = p[10];
85 
86  double deltaTBkg = p[11];
87  double sigmaBkg = p[12];
88  double bkg = p[13];
89 
90  double alpha = p[14];
91  double n = p[15];
92 
93  double mainPeak = fraction * TTSPDF(x[0], time, deltaTTS, sigma, TMath::Sqrt(sigma * sigma + sigmaRatio * sigmaRatio), f1);
94  double secondaryPeak = (1. - fraction) * TTSPDF(x[0], time + deltaTLaser, deltaTTS, sigma,
95  TMath::Sqrt(sigma * sigma + sigmaRatio * sigmaRatio), f1);
96  double extraPeak = crystalball_pdf(x[0], alpha, n, sigmaExtra, time + deltaTExtra);
97  double background = TMath::Gaus(x[0], time + deltaTBkg + deltaTTS, sigmaBkg, kTRUE);
98 
99  return norm * (mainPeak + secondaryPeak) + normExtra * extraPeak + bkg * background;
100 }
101 
102 
104 {
106  "Perform the fit of the laser and pulser runs"
107  );
108 
109 }
110 
111 
112 
114 {
115  m_inputTTS = TFile::Open(m_TTSData.c_str());
116  m_inputConstraints = TFile::Open(m_fitConstraints.c_str());
117 
118  B2INFO("Getting the TTS parameters from " << m_TTSData);
119  m_inputTTS->cd();
120  m_inputTTS->GetObject("tree", m_treeTTS);
121  m_treeTTS->SetBranchAddress("mean2", &m_mean2);
122  m_treeTTS->SetBranchAddress("sigma1", &m_sigma1);
123  m_treeTTS->SetBranchAddress("sigma2", &m_sigma2);
124  m_treeTTS->SetBranchAddress("fraction1", &m_f1);
125  m_treeTTS->SetBranchAddress("fraction2", &m_f2);
126  m_treeTTS->SetBranchAddress("pixelRow", &m_pixelRow);
127  m_treeTTS->SetBranchAddress("pixelCol", &m_pixelCol);
128 
129  if (m_fitterMode == "MC")
130  std::cout << "Running in MC mode, not constraints will be set" << std::endl;
131  else {
132  B2INFO("Getting the laser fit parameters from " << m_fitConstraints);
133  m_inputConstraints->cd();
134  m_inputConstraints->GetObject("fitTree", m_treeConstraints);
135  m_treeConstraints->SetBranchAddress("peakTime", &m_peakTimeConstraints);
136  m_treeConstraints->SetBranchAddress("deltaT", &m_deltaTConstraints);
137  m_treeConstraints->SetBranchAddress("fraction", &m_fractionConstraints);
138  if (m_fitterMode == "monitoring") {
139  m_treeConstraints->SetBranchAddress("timeExtra", &m_timeExtraConstraints);
140  m_treeConstraints->SetBranchAddress("sigmaExtra", &m_sigmaExtraConstraints);
141  m_treeConstraints->SetBranchAddress("alphaExtra", &m_alphaExtraConstraints);
142  m_treeConstraints->SetBranchAddress("nExtra", &m_nExtraConstraints);
143  m_treeConstraints->SetBranchAddress("timeBackground", &m_timeBackgroundConstraints);
144  m_treeConstraints->SetBranchAddress("sigmaBackground", &m_sigmaBackgroundConstraints);
145  }
146  }
147  return;
148 }
149 
150 
152 {
153  m_histFile = new TFile(m_output.c_str(), "recreate");
154  m_histFile->cd();
155  m_fitTree = new TTree("fitTree", "fitTree");
156  m_fitTree->Branch<short>("channel", &m_channel);
157  m_fitTree->Branch<short>("slot", &m_slot);
158  m_fitTree->Branch<short>("row", &m_row);
159  m_fitTree->Branch<short>("col", &m_col);
160  m_fitTree->Branch<float>("peakTime", &m_peakTime);
161  m_fitTree->Branch<float>("peakTimeErr", &m_peakTimeErr);
162  m_fitTree->Branch<float>("deltaT", &m_deltaT);
163  m_fitTree->Branch<float>("deltaTErr", &m_deltaTErr);
164  m_fitTree->Branch<float>("sigma", &m_sigma);
165  m_fitTree->Branch<float>("sigmaErr", &m_sigmaErr);
166  m_fitTree->Branch<float>("fraction", &m_fraction);
167  m_fitTree->Branch<float>("fractionErr", &m_fractionErr);
168  m_fitTree->Branch<float>("yieldLaser", &m_yieldLaser);
169  m_fitTree->Branch<float>("yieldLaserErr", &m_yieldLaserErr);
170  m_fitTree->Branch<float>("timeExtra", &m_timeExtra);
171  m_fitTree->Branch<float>("sigmaExtra", &m_sigmaExtra);
172  m_fitTree->Branch<float>("nExtra", &m_nExtra);
173  m_fitTree->Branch<float>("alphaExtra", &m_alphaExtra);
174  m_fitTree->Branch<float>("yieldLaserExtra", &m_yieldLaserExtra);
175  m_fitTree->Branch<float>("timeBackground", &m_timeBackground);
176  m_fitTree->Branch<float>("sigmaBackground", &m_sigmaBackground);
177  m_fitTree->Branch<float>("yieldLaserBackground", &m_yieldLaserBackground);
178 
179  m_fitTree->Branch<float>("fractionMC", &m_fractionMC);
180  m_fitTree->Branch<float>("deltaTMC", &m_deltaTMC);
181  m_fitTree->Branch<float>("peakTimeMC", &m_peakTimeMC);
182  m_fitTree->Branch<float>("firstPulserTime", &m_firstPulserTime);
183  m_fitTree->Branch<float>("firstPulserSigma", &m_firstPulserSigma);
184  m_fitTree->Branch<float>("secondPulserTime", &m_secondPulserTime);
185  m_fitTree->Branch<float>("secondPulserSigma", &m_secondPulserSigma);
186  m_fitTree->Branch<short>("fitStatus", &m_fitStatus);
187 
188  m_fitTree->Branch<float>("chi2", &m_chi2);
189  m_fitTree->Branch<float>("rms", &m_rms);
190 
191 
193  m_timewalkTree = new TTree("timewalkTree", "timewalkTree");
194 
195  m_timewalkTree->Branch<float>("binLowerEdge", &m_binLowerEdge);
196  m_timewalkTree->Branch<float>("binUpperEdge", &m_binUpperEdge);
197  m_timewalkTree->Branch<short>("channel", &m_channel);
198  m_timewalkTree->Branch<short>("slot", &m_slot);
199  m_timewalkTree->Branch<short>("row", &m_row);
200  m_timewalkTree->Branch<short>("col", &m_col);
201  m_timewalkTree->Branch<float>("histoIntegral", &m_histoIntegral);
202  m_timewalkTree->Branch<float>("peakTime", &m_peakTime);
203  m_timewalkTree->Branch<float>("peakTimeErr", &m_peakTimeErr);
204  m_timewalkTree->Branch<float>("deltaT", &m_deltaT);
205  m_timewalkTree->Branch<float>("deltaTErr", &m_deltaTErr);
206  m_timewalkTree->Branch<float>("sigma", &m_sigma);
207  m_timewalkTree->Branch<float>("sigmaErr", &m_sigmaErr);
208  m_timewalkTree->Branch<float>("fraction", &m_fraction);
209  m_timewalkTree->Branch<float>("fractionErr", &m_fractionErr);
210  m_timewalkTree->Branch<float>("yieldLaser", &m_yieldLaser);
211  m_timewalkTree->Branch<float>("yieldLaserErr", &m_yieldLaserErr);
212  m_timewalkTree->Branch<float>("timeExtra", &m_timeExtra);
213  m_timewalkTree->Branch<float>("sigmaExtra", &m_sigmaExtra);
214  m_timewalkTree->Branch<float>("nExtra", &m_nExtra);
215  m_timewalkTree->Branch<float>("alphaExtra", &m_alphaExtra);
216  m_timewalkTree->Branch<float>("yieldLaserExtra", &m_yieldLaserExtra);
217  m_timewalkTree->Branch<float>("timeBackground", &m_timeBackground);
218  m_timewalkTree->Branch<float>("sigmaBackground", &m_sigmaBackground);
219  m_timewalkTree->Branch<float>("yieldLaserBackground", &m_yieldLaserBackground);
220 
221  m_timewalkTree->Branch<float>("fractionMC", &m_fractionMC);
222  m_timewalkTree->Branch<float>("deltaTMC", &m_deltaTMC);
223  m_timewalkTree->Branch<float>("peakTimeMC", &m_peakTimeMC);
224  m_timewalkTree->Branch<float>("firstPulserTime", &m_firstPulserTime);
225  m_timewalkTree->Branch<float>("firstPulserSigma", &m_firstPulserSigma);
226  m_timewalkTree->Branch<float>("secondPulserTime", &m_secondPulserTime);
227  m_timewalkTree->Branch<float>("secondPulserSigma", &m_secondPulserSigma);
228  m_timewalkTree->Branch<short>("fitStatus", &m_fitStatus);
229 
230  m_timewalkTree->Branch<float>("chi2", &m_chi2);
231  m_timewalkTree->Branch<float>("rms", &m_rms);
232 
233  }
234 
235  return;
236 }
237 
238 
239 
240 
241 void TOPLocalCalFitter::fitChannel(short iSlot, short iChannel, TH1* h_profile)
242 {
243 
244  // loads the TTS infos and the fit constraint for the given channel and slot
245  if (m_fitterMode == "monitoring")
246  m_treeConstraints->GetEntry(iChannel + 512 * iSlot);
247  else if (m_fitterMode == "calibration") // The MC-based constraint file has only slot 1 at the moment
248  m_treeConstraints->GetEntry(iChannel);
249 
250  m_treeTTS->GetEntry(iChannel + 512 * iSlot);
251  // finds the maximum of the hit timing histogram and adjust the histogram range around it (3 ns window)
252  double maxpos = h_profile->GetBinCenter(h_profile->GetMaximumBin());
253  h_profile->GetXaxis()->SetRangeUser(maxpos - 1, maxpos + 2.);
254 
255  // gets the histogram integral to give a starting value to the fitter
256  double integral = h_profile->Integral();
257 
258  // creates the fit function
259  TF1 laser = TF1("laser", laserPDF, maxpos - 1, maxpos + 2., 16);
260 
261  // par[0] = peakTime
262  laser.SetParameter(0, maxpos);
263  laser.SetParLimits(0, maxpos - 0.06, maxpos + 0.06);
264 
265  // par[1] = sigma
266  laser.SetParameter(1, 0.1);
267  laser.SetParLimits(1, 0.05, 0.25);
268  if (m_fitterMode == "MC") {
269  laser.SetParameter(1, 0.02);
270  laser.SetParLimits(1, 0., 0.04);
271  }
272  // par[2] = fraction of the main peak respect to the total
273  laser.SetParameter(2, m_fractionConstraints);
274  laser.SetParLimits(2, 0.5, 1.);
275 
276  // par[3]= time difference between the main and secondary path. fixed to the MC value
277  laser.FixParameter(3, m_deltaTConstraints);
278 
279  // This is an hack: in some channels the MC sees one peak only, while in the data there are clearly
280  // two well distinguished peaks. This will disappear if we'll ever get a better laser simulation.
281  if (m_deltaTConstraints > -0.001) {
282  laser.SetParameter(3, -0.3);
283  laser.SetParLimits(3, -0.4, -0.2);
284  }
285 
286  // par[4] is the quadratic difference of the sigmas of the two TTS gaussians (tail - core)
287  laser.FixParameter(4, TMath::Sqrt(m_sigma2 * m_sigma2 - m_sigma1 * m_sigma1));
288  // par[5] is the position of the second TTS gaussian w/ respect to the first one
289  laser.FixParameter(5, m_mean2);
290  // par[6] is the relative contribution of the second TTS gaussian
291  laser.FixParameter(6, m_f1);
292  if (m_fitterMode == "MC")
293  laser.FixParameter(6, 0);
294 
295  // par[7] is the PDF normalization, = integral*bin width
296  laser.SetParameter(7, integral * 0.005);
297  laser.SetParLimits(7, 0.2 * integral * 0.005, 2.*integral * 0.005);
298 
299  // par[8-10] are the relative position, the sigma and the integral of the extra peak
300  laser.SetParameter(8, 1.);
301  laser.SetParLimits(8, 0.3, 2.);
302  laser.SetParameter(9, 0.2);
303  laser.SetParLimits(9, 0.08, 1.);
304  laser.SetParameter(10, 0.1 * integral * 0.005);
305  laser.SetParLimits(10, 0., 0.2 * integral * 0.005);
306  // par[14-15] are the tail parameters of the crystal ball function used to describe the extra peak
307  laser.SetParameter(14, -2.);
308  laser.SetParameter(15, 2.);
309  laser.SetParLimits(15, 1.01, 20.);
310 
311  // par[11-13] are relative position, sigma and integral of the broad gaussian added to better describe the tail at high times
312  laser.SetParameter(11, 1.);
313  laser.SetParLimits(11, 0.1, 5.);
314  laser.SetParameter(12, 0.8);
315  laser.SetParLimits(12, 0., 5.);
316  laser.SetParameter(13, 0.01 * integral * 0.005);
317  laser.SetParLimits(13, 0., 0.2 * integral * 0.005);
318 
319  // if it's a monitoring fit, fix a buch more parameters.
320  if (m_fitterMode == "monitoring") {
321  laser.FixParameter(2, m_fractionConstraints);
322  laser.FixParameter(3, m_deltaTConstraints);
323  laser.FixParameter(8, m_timeExtraConstraints);
324  laser.FixParameter(9, m_sigmaExtraConstraints);
325  laser.FixParameter(14, m_alphaExtraConstraints);
326  laser.FixParameter(15, m_nExtraConstraints);
327  laser.FixParameter(11, m_timeBackgroundConstraints);
328  laser.FixParameter(12, m_sigmaBackgroundConstraints);
329  }
330 
331 
332  // if it's a MC fit, fix a buch more parameters.
333  if (m_fitterMode == "MC") {
334  laser.SetParameter(2, 0.8);
335  laser.SetParLimits(2, 0., 1.);
336  laser.SetParameter(3, -0.1);
337  laser.SetParLimits(3, -0.4, -0.);
338  // The following are just random reasonable number, only to pin-point the tail components to some value and remove them form the fit
339  laser.FixParameter(8, 0);
340  laser.FixParameter(9, 0.1);
341  laser.FixParameter(14, -2.);
342  laser.FixParameter(15, 2);
343  laser.FixParameter(11, 1.);
344  laser.FixParameter(12, 0.1);
345  laser.FixParameter(13, 0.);
346  laser.FixParameter(10, 0.);
347  }
348 
349 
350 
351  // make the plot of the fit function nice setting 2000 sampling points
352  laser.SetNpx(2000);
353 
354 
355  // do the fit!
356  h_profile->Fit("laser", "R L Q");
357 
358 
359  // Add by hand the different fit components to the histogram, mostly for debugging/presentation purposes
360  TF1* peak1 = new TF1("peak1", laserPDF, maxpos - 1, maxpos + 2., 16);
361  TF1* peak2 = new TF1("peak2", laserPDF, maxpos - 1, maxpos + 2., 16);
362  TF1* extra = new TF1("extra", laserPDF, maxpos - 1, maxpos + 2., 16);
363  TF1* background = new TF1("background", laserPDF, maxpos - 1, maxpos + 2., 16);
364  for (int iPar = 0; iPar < 16; iPar++) {
365  peak1->FixParameter(iPar, laser.GetParameter(iPar));
366  peak2->FixParameter(iPar, laser.GetParameter(iPar));
367  extra->FixParameter(iPar, laser.GetParameter(iPar));
368  background->FixParameter(iPar, laser.GetParameter(iPar));
369  }
370  peak1->FixParameter(2, 0.);
371  peak1->FixParameter(7, (1 - laser.GetParameter(2))*laser.GetParameter(7));
372  peak1->FixParameter(10, 0.);
373  peak1->FixParameter(13, 0.);
374 
375  peak2->FixParameter(2, 1.);
376  peak2->FixParameter(7, laser.GetParameter(2)*laser.GetParameter(7));
377  peak2->FixParameter(10, 0.);
378  peak2->FixParameter(13, 0.);
379 
380  extra->FixParameter(7, 0.);
381  extra->FixParameter(13, 0.);
382 
383  background->FixParameter(7, 0.);
384  background->FixParameter(10, 0.);
385 
386  h_profile->GetListOfFunctions()->Add(peak1);
387  h_profile->GetListOfFunctions()->Add(peak2);
388  h_profile->GetListOfFunctions()->Add(extra);
389  h_profile->GetListOfFunctions()->Add(background);
390 
391 
392  // save the results in the variables linked to the tree branches
393  m_channel = iChannel;
394  m_row = m_pixelRow;
395  m_col = m_pixelCol;
396  m_slot = iSlot + 1;
397  m_peakTime = laser.GetParameter(0);
398  m_peakTimeErr = laser.GetParError(0);
399  m_deltaT = laser.GetParameter(3);
400  m_deltaTErr = laser.GetParError(3);
401  m_sigma = laser.GetParameter(1);
402  m_sigmaErr = laser.GetParError(1);
403  m_fraction = laser.GetParameter(2);
404  m_fractionErr = laser.GetParError(2);
405  m_yieldLaser = laser.GetParameter(7) / 0.005;
406  m_yieldLaserErr = laser.GetParError(7) / 0.005;
407  m_timeExtra = laser.GetParameter(8);
408  m_sigmaExtra = laser.GetParameter(9);
409  m_yieldLaserExtra = laser.GetParameter(10) / 0.005;
410  m_alphaExtra = laser.GetParameter(14);
411  m_nExtra = laser.GetParameter(15);
412  m_timeBackground = laser.GetParameter(11);
413  m_sigmaBackground = laser.GetParameter(12);
414  m_yieldLaserBackground = laser.GetParameter(13) / 0.005;
415  m_chi2 = laser.GetChisquare() / laser.GetNDF();
416 
417  // copy some MC information to the output tree
421 
422  return;
423 }
424 
425 
426 void TOPLocalCalFitter::fitPulser(TH1* h_profileFirstPulser, TH1* h_profileSecondPulser)
427 {
428  float maxpos = h_profileFirstPulser->GetBinCenter(h_profileFirstPulser->GetMaximumBin());
429  h_profileFirstPulser->GetXaxis()->SetRangeUser(maxpos - 1, maxpos + 1.);
430  if (h_profileFirstPulser->Integral() > 1000) {
431  TF1 pulser1 = TF1("pulser1", "[0]*TMath::Gaus(x, [1], [2], kTRUE)", maxpos - 1, maxpos + 1.);
432  pulser1.SetParameter(0, 1.);
433  pulser1.SetParameter(1, maxpos);
434  pulser1.SetParameter(2, 0.05);
435  h_profileFirstPulser->Fit("pulser1", "R Q");
436  m_firstPulserTime = pulser1.GetParameter(1);
437  m_firstPulserSigma = pulser1.GetParameter(2);
438  h_profileFirstPulser->Write();
439  } else {
440  m_firstPulserTime = -999;
441  m_firstPulserSigma = -999;
442  }
443 
444  maxpos = h_profileSecondPulser->GetBinCenter(h_profileSecondPulser->GetMaximumBin());
445  h_profileSecondPulser->GetXaxis()->SetRangeUser(maxpos - 1, maxpos + 1.);
446  if (h_profileSecondPulser->Integral() > 1000) {
447  TF1 pulser2 = TF1("pulser2", "[0]*TMath::Gaus(x, [1], [2], kTRUE)", maxpos - 1, maxpos + 1.);
448  pulser2.SetParameter(0, 1.);
449  pulser2.SetParameter(1, maxpos);
450  pulser2.SetParameter(2, 0.05);
451  h_profileSecondPulser->Fit("pulser2", "R Q");
452  m_secondPulserTime = pulser2.GetParameter(1);
453  m_secondPulserSigma = pulser2.GetParameter(2);
454  h_profileSecondPulser->Write();
455  } else {
456  m_secondPulserTime = -999;
457  m_secondPulserSigma = -999;
458  }
459  return;
460 }
461 
462 
463 
465 {
466  if (m_chi2 < 4 && m_sigma < 0.2 && m_yieldLaser > 1000) {
467  m_fitStatus = 0;
468  } else {
469  m_fitStatus = 1;
470  }
471  return;
472 }
473 
474 
476 {
477  int nEntries = m_fitTree->GetEntries();
478  if (nEntries != 8192) {
479  B2ERROR("fitTree does not contain an entry wit a fit result for each channel. Found " << nEntries <<
480  " instead of 8192. Perhaps you tried to run the commonT0 calculation before finishing the fitting?");
481  return;
482  }
483 
484  // Create and fill the TOPCalChannelT0 object.
485  // This part is mostly copy-pasted from the DB importer used up to Jan 2020
486  auto* channelT0 = new TOPCalChannelT0();
487  short nCal[16] = {0};
488  for (int i = 0; i < nEntries; i++) {
489  m_fitTree->GetEntry(i);
490  channelT0->setT0(m_slot, m_channel, m_peakTime - m_peakTimeMC, m_peakTimeErr);
491  if (m_fitStatus == 0) {
492  nCal[m_slot - 1]++;
493  } else {
494  channelT0->setUnusable(m_slot, m_channel);
495  }
496  }
497 
498  // Normalize the constants
499  channelT0->suppressAverage();
500 
501  // create the localDB
502  saveCalibration(channelT0);
503 
504  short nCalTot = 0;
505  B2INFO("Summary: ");
506  for (int iSlot = 1; iSlot < 17; iSlot++) {
507  B2INFO("--> Number of calibrated channels on Slot " << iSlot << " : " << nCal[iSlot - 1] << "/512");
508  B2INFO("--> Cal on ch 1, 256 and 511: " << channelT0->getT0(iSlot, 0) << ", " << channelT0->getT0(iSlot,
509  257) << ", " << channelT0->getT0(iSlot, 511));
510  nCalTot += nCal[iSlot - 1];
511  }
512 
513  B2RESULT("Channel T0 calibration constants imported to database, calibrated channels: " << nCalTot << "/ 8192");
514 
515 
516  // Loop again on the output tree to save the constants there too, adding two more branches.
517  TBranch* channelT0Branch = m_fitTree->Branch<float>("channelT0", &m_channelT0);
518  TBranch* channelT0ErrBranch = m_fitTree->Branch<float>("channelT0Err", &m_channelT0Err);
519 
520  for (int i = 0; i < nEntries; i++) {
521  m_fitTree->GetEntry(i);
522  m_channelT0 = channelT0->getT0(m_slot, m_channel);
523  m_channelT0Err = channelT0->getT0Error(m_slot, m_channel);
524  channelT0Branch->Fill();
525  channelT0ErrBranch->Fill();
526  }
527  return ;
528 
529 }
530 
531 
532 
534 {
535  gROOT->SetBatch();
536 
537  loadMCInfoTrees();
538 
539 
541 
542  // Loads the tree with the hits
543  auto hitTree = getObjectPtr<TTree>("hitTree");
544  TH2F* h_hitTime = new TH2F("h_hitTime", " ", 512 * 16, 0., 512 * 16, 22000, -70, 40.); // 5 ps bins
545 
546  float amplitude, hitTime;
547  short channel, slot;
548  bool refTimeValid;
549  hitTree->SetBranchAddress("amplitude", &amplitude);
550  hitTree->SetBranchAddress("hitTime", &hitTime);
551  hitTree->SetBranchAddress("channel", &channel);
552  hitTree->SetBranchAddress("slot", &slot);
553  hitTree->SetBranchAddress("refTimeValid", &refTimeValid);
554 
555  // An attempt to speed things up looping over the tree only once.
556  std::vector<TH2F*> h_hitTimeLaserHistos = {};
557  for (int iLowerEdge = 0; iLowerEdge < (int)m_binEdges.size() - 1; iLowerEdge++) {
558  TH2F* h_hitTimeLaser = new TH2F(("h_hitTimeLaser_" + std::to_string(iLowerEdge + 1)).c_str(), " ", 512 * 16, 0., 512 * 16, 14000,
559  -70, 0.); // 5 ps bins
560  h_hitTimeLaserHistos.push_back(h_hitTimeLaser);
561  }
562 
563  for (unsigned int i = 0; i < hitTree->GetEntries(); i++) {
564  auto onepc = (unsigned int)(hitTree->GetEntries() / 100);
565  if (i % onepc == 0)
566  std::cout << "processing hit " << i << " of " << hitTree->GetEntries() << " (" << i / (onepc * 10) << " %)" << std::endl;
567  hitTree->GetEntry(i);
568  auto it = std::lower_bound(m_binEdges.cbegin(), m_binEdges.cend(), amplitude);
569  int iLowerEdge = std::distance(m_binEdges.cbegin(), it) - 1;
570 
571  if (iLowerEdge >= 0 && iLowerEdge < static_cast<int>(m_binEdges.size()) - 1 && refTimeValid)
572  h_hitTimeLaserHistos[iLowerEdge]->Fill(channel + (slot - 1) * 512, hitTime);
573  if (amplitude > 80. && refTimeValid)
574  h_hitTime->Fill(channel + (slot - 1) * 512, hitTime);
575  }
576 
577  m_histFile->cd();
578  h_hitTime->Write();
579 
580  for (short iSlot = 0; iSlot < 16; iSlot++) {
581  std::cout << "fitting slot " << iSlot + 1 << std::endl;
582  for (short iChannel = 0; iChannel < 512; iChannel++) {
583  TH1D* h_profile = h_hitTime->ProjectionY(("profile_" + std::to_string(iSlot + 1) + "_" + std::to_string(iChannel)).c_str(),
584  iSlot * 512 + iChannel + 1, iSlot * 512 + iChannel + 1);
585  if (m_fitterMode == "MC")
586  h_profile->GetXaxis()->SetRangeUser(-10, -10);
587  else
588  h_profile->GetXaxis()->SetRangeUser(-65,
589  -5); // if you will even change it, make sure not to include the h_hitTime overflow bins in this range
590  fitChannel(iSlot, iChannel, h_profile);
591 
593 
594  // Now let's fit the pulser
595  TH1D* h_profileFirstPulser = h_hitTime->ProjectionY(("profileFirstPulser_" + std::to_string(iSlot + 1) + "_" + std::to_string(
596  iChannel)).c_str(), iSlot * 512 + iChannel + 1, iSlot * 512 + iChannel + 1);
597  TH1D* h_profileSecondPulser = h_hitTime->ProjectionY(("profileSecondPulser_" + std::to_string(iSlot + 1) + "_" + std::to_string(
598  iChannel)).c_str(), iSlot * 512 + iChannel + 1, iSlot * 512 + iChannel + 1);
599  h_profileFirstPulser->GetXaxis()->SetRangeUser(-10, 10);
600  h_profileSecondPulser->GetXaxis()->SetRangeUser(10, 40);
601 
602  fitPulser(h_profileFirstPulser, h_profileSecondPulser);
603 
604  m_fitTree->Fill();
605 
606  h_profile->Write();
607  h_profileFirstPulser->Write();
608  h_profileSecondPulser->Write();
609  }
610 
611  h_hitTime->Write();
612  }
613 
615 
616  m_fitTree->Write();
617 
618 
620  std::cout << "Fitting in bins" << std::endl;
621  for (int iLowerEdge = 0; iLowerEdge < (int)m_binEdges.size() - 1; iLowerEdge++) {
622  m_binLowerEdge = m_binEdges[iLowerEdge];
623  m_binUpperEdge = m_binEdges[iLowerEdge + 1];
624  std::cout << "Fitting the amplitude interval (" << m_binLowerEdge << ", " << m_binUpperEdge << " )" << std::endl;
625 
626  for (short iSlot = 0; iSlot < 16; iSlot++) {
627  std::cout << " Fitting slot " << iSlot + 1 << std::endl;
628  for (short iChannel = 0; iChannel < 512; iChannel++) {
629  TH1D* h_profile = h_hitTimeLaserHistos[iLowerEdge]->ProjectionY(("profile_" + std::to_string(iSlot + 1) + "_" + std::to_string(
630  iChannel) + "_" + std::to_string(iLowerEdge)).c_str(),
631  iSlot * 512 + iChannel + 1, iSlot * 512 + iChannel + 1);
632  if (m_fitterMode == "MC")
633  h_profile->GetXaxis()->SetRangeUser(-10, -10);
634  else
635  h_profile->GetXaxis()->SetRangeUser(-65,
636  -5); // if you will even change it, make sure not to include the h_hitTime overflow bins in this range
637  fitChannel(iSlot, iChannel, h_profile);
638  m_histoIntegral = h_profile->Integral();
640 
641  m_timewalkTree->Fill();
642  h_profile->Write();
643  }
644  }
645  }
646 
647  m_timewalkTree->Write();
648  }
649 
650  m_histFile->Close();
651 
652  return c_OK;
653 }
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
Channel T0 calibration constants for all 512 channels of 16 modules.
short m_fitStatus
Fit quality flag, propagated to the constants.
float m_yieldLaserErr
Statistical error on yield.
float m_binLowerEdge
Lower edge of the amplitude bin in which this fit is performed.
void fitChannel(short, short, TH1 *)
Fits the laser light on one channel.
float m_nExtraConstraints
parameter n of the tail of the extra peak
float m_chi2
Reduced chi2 of the fit.
float m_timeExtraConstraints
Position of the guassian used to describe the extra peak on the timing distribution tail.
float m_f1
Fraction of the first gaussian on the TTS parametrization.
float m_fraction
Fraction of events in the secondary peak.
float m_sigmaExtra
Gaussian sigma of the extra peak in the timing tail
float m_secondPulserSigma
Time resolution from the fit of the first electronic pulse, from a Gaussian fit.
float m_yieldLaserBackground
Integral of the background gaussian.
float m_sigma
Gaussian time resolution, fitted.
float m_peakTimeMC
Time of the main peak in teh MC simulation, i.e.
std::string m_output
Name of the output file.
float m_timeBackground
Position of the gaussian used to describe the background, w/ respect to peakTime.
void loadMCInfoTrees()
loads the TTS parameters and the MC truth info
float m_deltaTMC
Time difference between the main peak and the secondary peak in the MC simulation.
float m_peakTimeErr
Statistical error on peakTime.
short m_channel
Channel number (0-511)
float m_channelT0Err
Statistical error on channelT0.
std::string m_TTSData
File with the Fit constraints and MC info.
std::vector< float > m_binEdges
Amplitude bins.
float m_sigmaBackground
Sigma of the gaussian used to describe the background.
float m_peakTime
Fitted time of the main (i.e.
float m_fractionErr
Statistical error on fraction.
bool m_isFitInAmplitudeBins
Enables the fit in amplitude bins.
void calculateChennelT0()
Calculates the commonT0 calibration after the fits have been done.
float m_alphaExtra
alpha parameter of the tail of the extra peak.
float m_rms
RMS of the histogram used for the fit.
float m_mean2
Position of the second gaussian of the TTS parametrization with respect to the first one.
std::string m_fitConstraints
File with the TTS parametrization.
TTree * m_timewalkTree
Output of the fitter.
float m_yieldLaser
Total number of laser hits from the fitting function integral.
TFile * m_inputTTS
File containing m_treeTTS.
float m_nExtra
parameter n of the tail of the extra peak
float m_alphaExtraConstraints
alpha parameter of the tail of the extra peak.
float m_channelT0
Raw, channelT0 calibration, defined as peakTime-peakTimeMC.
float m_sigmaErr
Statistical error on sigma.
std::string m_fitterMode
Fit mode.
float m_binUpperEdge
Upper edge of the amplitude bin in which this fit is performed.
TFile * m_histFile
Output of the fitter.
TTree * m_treeConstraints
Input to the fitter.
float m_sigmaExtraConstraints
Width of the guassian used to describe the extra peak on the timing distribution tail.
float m_firstPulserSigma
Time resolution from the fit of the first electronic pulse, from a Gaussian fit.
float m_deltaT
Time difference between the main peak and the secondary peak.
TTree * m_fitTree
Output of the fitter.
float m_f2
Fraction of the second gaussian on the TTS parametrization.
void determineFitStatus()
determines if the constant obtained by the fit are good or not
float m_secondPulserTime
Average time of the second electronic pulse respect to the reference pulse, from a gaussian fit.
float m_peakTimeConstraints
Time of the main laser peak in the MC simulation (aka MC correction)
float m_sigma1
Width of the first gaussian on the TTS parametrization.
float m_fractionMC
Fraction of events in the secondary peak form the MC simulation.
float m_sigmaBackgroundConstraints
Sigma of the gaussian used to describe the background.
virtual EResult calibrate() override
Runs the algorithm on events.
TTree * m_treeTTS
Input to the fitter.
void setupOutputTreeAndFile()
prepares the output tree
float m_timeBackgroundConstraints
Position of the gaussian used to describe the background, w/ respect to peakTime.
TFile * m_inputConstraints
File containing m_treeConstraints.
float m_deltaTErr
Statistical error on deltaT.
float m_firstPulserTime
Average time of the first electronic pulse respect to the reference pulse, from a Gaussian fit.
float m_histoIntegral
Integral of the fitted histogram.
float m_yieldLaserExtra
Integral of the extra peak.
float m_sigma2
Width of the second gaussian on the TTS parametrization.
void fitPulser(TH1 *, TH1 *)
Fits the two pulsers.
float m_timeExtra
Position of the extra peak seen in the timing tail, w/ respect to peakTime.
float m_fractionConstraints
Fraction of the main peak.
float m_deltaTConstraints
Distance between the main and the secondary laser peak.
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.