Belle II Software  release-06-02-00
SVDCoGTimeCalibrationAlgorithm.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 #include <svd/calibration/SVDCoGTimeCalibrationAlgorithm.h>
9 
10 #include <svd/dbobjects/SVDCoGCalibrationFunction.h>
11 #include <svd/calibration/SVDCoGTimeCalibrations.h>
12 
13 #include <TF1.h>
14 #include <TProfile.h>
15 #include <TH2F.h>
16 #include <framework/logging/Logger.h>
17 #include <iostream>
18 #include <TString.h>
19 #include <TFitResult.h>
20 
21 using namespace std;
22 using namespace Belle2;
23 
24 SVDCoGTimeCalibrationAlgorithm::SVDCoGTimeCalibrationAlgorithm(const std::string& str) :
25  CalibrationAlgorithm("SVDTimeCalibrationCollector")
26 {
27  setDescription("SVDCoGTimeCalibration calibration algorithm");
28  m_id = str;
29 }
30 
32 {
33 
34  gROOT->SetBatch(true);
35 
36  int ladderOfLayer[4] = {7, 10, 12, 16};
37  int sensorOnLayer[4] = {2, 3, 4, 5};
38 
39  auto timeCal = new Belle2::SVDCoGCalibrationFunction();
40  auto payload = new Belle2::SVDCoGTimeCalibrations::t_payload(*timeCal, m_id);
41 
42  std::unique_ptr<TF1> pol1(new TF1("pol1", "[0] + [1]*x", -10, 80));
43  pol1->SetParameters(-40, 0.9);
44  std::unique_ptr<TF1> pol3(new TF1("pol3", "[0] + [1]*x + [2]*x*x + [3]*x*x*x", -10, 80));
45  pol3->SetParLimits(0, -200, 0);
46  pol3->SetParLimits(1, 0, 10);
47  pol3->SetParLimits(2, -1, 0);
48  pol3->SetParLimits(3, 0, 1);
49  std::unique_ptr<TF1> pol5(new TF1("pol5", "[0] + [1]*x + [2]*x*x + [3]*x*x*x + [4]*x*x*x*x + [5]*x*x*x*x*x", -100, 100));
50  pol5->SetParameters(-50, 1.5, 0.01, 0.0001, 0.00001, 0.000001);
51 
52  FileStat_t info;
53  int cal_rev = 1;
54  while (gSystem->GetPathInfo(Form("algorithm_6SampleCoG_output_rev_%d.root", cal_rev), info) == 0)
55  cal_rev++;
56  std::unique_ptr<TFile> f(new TFile(Form("algorithm_6SampleCoG_output_rev_%d.root", cal_rev), "RECREATE"));
57 
58  auto m_tree = new TTree(Form("rev_%d", cal_rev), "RECREATE");
59  int layer_num, ladder_num, sensor_num, view, ndf;
60  float a, b, c, d, a_err, b_err, c_err, d_err, chi2, p;
61  m_tree->Branch("layer", &layer_num, "layer/I");
62  m_tree->Branch("ladder", &ladder_num, "ladder/I");
63  m_tree->Branch("sensor", &sensor_num, "sensor/I");
64  m_tree->Branch("isU", &view, "isU/I");
65  m_tree->Branch("a", &a, "a/F");
66  m_tree->Branch("b", &b, "b/F");
67  m_tree->Branch("c", &c, "c/F");
68  m_tree->Branch("d", &d, "d/F");
69  m_tree->Branch("a_err", &a_err, "a_err/F");
70  m_tree->Branch("b_err", &b_err, "b_err/F");
71  m_tree->Branch("c_err", &c_err, "c_err/F");
72  m_tree->Branch("d_err", &d_err, "d_err/F");
73  m_tree->Branch("chi2", &chi2, "chi2/F");
74  m_tree->Branch("ndf", &ndf, "ndf/I");
75  m_tree->Branch("p", &p, "p/F");
76 
77  for (int layer = 0; layer < 4; layer++) {
78  layer_num = layer + 3;
79  for (int ladder = 0; ladder < (int)ladderOfLayer[layer]; ladder++) {
80  ladder_num = ladder + 1;
81  for (int sensor = 0; sensor < (int)sensorOnLayer[layer]; sensor++) {
82  sensor_num = sensor + 1;
83  for (view = 1; view > -1; view--) {
84  char side = 'U';
85  if (view == 0)
86  side = 'V';
87  auto hEventT0vsCoG = getObjectPtr<TH2F>(Form("eventT0vsCoG__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
88  auto hEventT0 = getObjectPtr<TH1F>(Form("eventT0__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
89  auto hEventT0nosync = getObjectPtr<TH1F>(Form("eventT0nosync__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
90  B2INFO("Histogram: " << hEventT0vsCoG->GetName() <<
91  " Entries (n. clusters): " << hEventT0vsCoG->GetEntries());
92  if (layer_num == 3 && hEventT0vsCoG->GetEntries() < m_minEntries) {
93  B2INFO("Histogram: " << hEventT0vsCoG->GetName() <<
94  " Entries (n. clusters): " << hEventT0vsCoG->GetEntries() <<
95  " Entries required: " << m_minEntries);
96  B2WARNING("Not enough data, adding one run to the collector");
97  f->Close();
98  gSystem->Unlink(Form("algorithm_6SampleCoG_output_rev_%d.root", cal_rev));
99  return c_NotEnoughData;
100  }
101  for (int i = 1; i <= hEventT0vsCoG->GetNbinsX(); i++) {
102  for (int j = 1; j <= hEventT0vsCoG->GetNbinsY(); j++) {
103  if (hEventT0vsCoG->GetBinContent(i, j) < max(2, int(hEventT0vsCoG->GetEntries() * 0.001))) {
104  hEventT0vsCoG->SetBinContent(i, j, 0);
105  }
106  }
107  }
108  TProfile* pfx = hEventT0vsCoG->ProfileX();
109  std::string name = "pfx_" + std::string(hEventT0vsCoG->GetName());
110  pfx->SetName(name.c_str());
111  TFitResultPtr tfr = pfx->Fit("pol3", "RQSM");
112  double par[4];
113  pol3->GetParameters(par);
115  /*
116  pfx->Fit("pol1", "RQ");
117  double par[4];
118  pol1->GetParameters(par);
119  par[2] = 0;
120  par[3] = 0;
121  */
122  // double meanT0 = hEventT0->GetMean();
123  // double meanT0NoSync = hEventT0nosync->GetMean();
124  timeCal->set_current(1);
125  // timeCal->set_current(2);
126  timeCal->set_pol3parameters(par[0], par[1], par[2], par[3]);
127  payload->set(layer_num, ladder_num, sensor_num, bool(view), 1, *timeCal);
128  f->cd();
129  hEventT0->Write();
130  hEventT0vsCoG->Write();
131  hEventT0nosync->Write();
132  pfx->Write();
133 
134  if (!tfr) {
135  f->Close();
136  B2FATAL("Fit to the histogram failed in SVDCoGTimeCalibrationAlgorithm. "
137  << "Check the 2-D histogram to clarify the reason.");
138  } else {
139  a = par[0]; b = par[1]; c = par[2]; d = par[3];
140  a_err = tfr->ParError(0); b_err = tfr->ParError(1); c_err = tfr->ParError(2); d_err = tfr->ParError(3);
141  chi2 = tfr->Chi2();
142  ndf = tfr->Ndf();
143  p = tfr->Prob();
144  m_tree->Fill();
145  }
146 
147  }
148  }
149  }
150  }
151  m_tree->Write();
152  f->Close();
153  saveCalibration(payload, "SVDCoGTimeCalibrations");
154 
155  //delete f;
156 
157  // probably not needed - would trigger re-doing the collection
158  // if ( ... too large corrections ... ) return c_Iterate;
159  return c_OK;
160 }
161 
162 bool SVDCoGTimeCalibrationAlgorithm::isBoundaryRequired(const Calibration::ExpRun& currentRun)
163 {
164  float meanRawTimeL3V = 0;
165  // auto eventT0Hist = getObjectPtr<TH1F>("hEventT0FromCDC");
166  auto rawTimeL3V = getObjectPtr<TH1F>("hRawTimeL3V");
167  // float meanEventT0 = eventT0Hist->GetMean();
168  if (!rawTimeL3V) {
170  meanRawTimeL3V = m_previousRawTimeMeanL3V.value();
171  } else {
172  if (rawTimeL3V->GetEntries() > m_minEntries)
173  meanRawTimeL3V = rawTimeL3V->GetMean();
174  else {
176  meanRawTimeL3V = m_previousRawTimeMeanL3V.value();
177  }
178  }
180  B2INFO("Setting start payload boundary to be the first run ("
181  << currentRun.first << "," << currentRun.second << ")");
182  m_previousRawTimeMeanL3V.emplace(meanRawTimeL3V);
183 
184  return true;
185  } else if (abs(meanRawTimeL3V - m_previousRawTimeMeanL3V.value()) > m_allowedTimeShift) {
186  B2INFO("Histogram mean has shifted from " << m_previousRawTimeMeanL3V.value()
187  << " to " << meanRawTimeL3V << ". We are requesting a new payload boundary for ("
188  << currentRun.first << "," << currentRun.second << ")");
189  m_previousRawTimeMeanL3V.emplace(meanRawTimeL3V);
190  return true;
191  } else {
192  return false;
193  }
194 }
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.
@ c_NotEnoughData
Needs more data =2 in Python.
class to contain the CoG Time calibrations
std::string m_id
Parameter given to set the UniqueID of the payload.
std::optional< float > m_previousRawTimeMeanL3V
CoG time mean of the previous run for V side of layer 3.
float m_minEntries
Set the minimun number of entries required in the histograms of layer 3.
virtual EResult calibrate() override
Run algo on data.
virtual bool isBoundaryRequired(const Calibration::ExpRun &currentRun) override
If the event T0 changes significantly return true.
SVDCalibrationsBase< SVDCalibrationsScalar< SVDCoGCalibrationFunction > > t_payload
typedef for the SVDCoGCalibrationFunction payload of all SVD sensors
Abstract base class for different kinds of events.