Belle II Software  release-08-01-10
SVD3SampleCoGTimeCalibrationAlgorithm.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/SVD3SampleCoGTimeCalibrationAlgorithm.h>
9 
10 #include <svd/dbobjects/SVDCoGCalibrationFunction.h>
11 #include <svd/calibration/SVD3SampleCoGTimeCalibrations.h>
12 
13 #include <TF1.h>
14 #include <TProfile.h>
15 #include <TH1F.h>
16 #include <TH2F.h>
17 #include <TH3F.h>
18 #include <framework/logging/Logger.h>
19 #include <iostream>
20 #include <TString.h>
21 #include <TFitResult.h>
22 
23 using namespace std;
24 using namespace Belle2;
25 
26 SVD3SampleCoGTimeCalibrationAlgorithm::SVD3SampleCoGTimeCalibrationAlgorithm(const std::string& str) :
27  CalibrationAlgorithm("SVDTimeCalibrationCollector")
28 {
29  setDescription("SVD3SampleCoGTimeCalibration calibration algorithm");
30  m_id = str;
31 }
32 
34 {
35  gROOT->SetBatch(true);
36 
37  auto timeCal = new Belle2::SVDCoGCalibrationFunction();
38  auto payload = new Belle2::SVD3SampleCoGTimeCalibrations::t_payload(*timeCal, m_id);
39 
40  // In the local study, Y. Uematsu tuned the range to (31.5,48), because the correlation is not exactly like pol3. (The deviation appears at around 48, very naive.)
41  // However, original value (-10,80) also seems working.
42  std::unique_ptr<TF1> pol3(new TF1("pol3", "[0] + [1]*x + [2]*[3]*[3]*x - [2]*[3]*x*x + [2]*x*x*x/3", -10,
43  80));
44  // apply parameter limits to be monotonic: [1] > 0 and [2] > 0. upper limits are loose on [2].
45  // lower limit on [1] is the minimum tilt (= tilt at the pole)
46  pol3->SetParLimits(1, 0, 100);
47  pol3->SetParLimits(2, 0, 0.1);
48  // the pole position [3] is different in the data (~40) and the current MC (~45)
49  // apply parameter limits on [3] to help non-linear minimization
50  pol3->SetParLimits(3, 30, 60);
51 
52  FileStat_t info;
53  int cal_rev = 1;
54  while (gSystem->GetPathInfo(Form("algorithm_3SampleCoG_output_rev_%d.root", cal_rev), info) == 0)
55  cal_rev++;
56  std::unique_ptr<TFile> f(new TFile(Form("algorithm_3SampleCoG_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 
78  B2INFO("--------- Applyingselection, 2D-region selection parameters: ");
79  B2INFO("Upper Line (q, m): " << m_interceptUpperLine << ", " << m_angularCoefficientUpperLine);
80  B2INFO("Lower Line (q, m): " << m_interceptLowerLine << ", " << m_angularCoefficientLowerLine);
81  } //B2INFO("Selecton applied : " << m_applyLinearCutsToRemoveBkg);
82 
83  auto __hEventT0vsCoG__ = getObjectPtr<TH3F>("__hEventT0vsCoG__");
84  auto __hEventT0__ = getObjectPtr<TH2F>("__hEventT0__");
85  auto __hEventT0NoSync__ = getObjectPtr<TH2F>("__hEventT0NoSync__");
86  auto __hBinToSensorMap__ = getObjectPtr<TH1F>("__hBinToSensorMap__");
87 
88  for (int ij = 0; ij < (__hBinToSensorMap__->GetNbinsX()); ij++) {
89  {
90  {
91  {
92 
93  auto binLabel = __hBinToSensorMap__->GetXaxis()->GetBinLabel(ij + 1);
94  char side;
95  std::sscanf(binLabel, "L%dL%dS%d%c", &layer_num, &ladder_num, &sensor_num, &side);
96  view = 0;
97  if (side == 'U')
98  view = 1;
99 
100  B2INFO("Projecting for Sensor: " << binLabel << " with Bin Number: " << ij + 1);
101 
102  __hEventT0vsCoG__->GetZaxis()->SetRange(ij + 1, ij + 1);
103  auto hEventT0vsCoG = (TH2D*)__hEventT0vsCoG__->Project3D("yxe");
104  auto hEventT0 = (TH1D*)__hEventT0__->ProjectionX("hEventT0_tmp", ij + 1, ij + 1);
105  auto hEventT0nosync = (TH1D*)__hEventT0NoSync__->ProjectionX("hEventT0NoSync_tmp", ij + 1, ij + 1);
106 
107  hEventT0vsCoG->SetName(Form("eventT0vsCoG__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
108  hEventT0->SetName(Form("eventT0__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
109  hEventT0nosync->SetName(Form("eventT0nosync__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
110 
111  char sidePN = (side == 'U' ? 'P' : 'N');
112  hEventT0vsCoG->SetTitle(Form("EventT0Sync vs rawTime in %d.%d.%d %c/%c", layer_num, ladder_num, sensor_num, side, sidePN));
113  hEventT0->SetTitle(Form("EventT0Sync in %d.%d.%d %c/%c", layer_num, ladder_num, sensor_num, side, sidePN));
114  hEventT0nosync->SetTitle(Form("EventT0NoSync in %d.%d.%d %c/%c", layer_num, ladder_num, sensor_num, side, sidePN));
115 
116  hEventT0vsCoG->SetDirectory(0);
117  hEventT0->SetDirectory(0);
118  hEventT0nosync->SetDirectory(0);
119 
120  int nEntriesForFilter = hEventT0vsCoG->GetEntries();
121  B2INFO("Histogram: " << hEventT0vsCoG->GetName() <<
122  " Entries (n. clusters): " << hEventT0vsCoG->GetEntries());
123  if (layer_num == 3 && hEventT0vsCoG->GetEntries() < m_minEntries) {
124  B2INFO("Histogram: " << hEventT0vsCoG->GetName() <<
125  " Entries (n. clusters): " << hEventT0vsCoG->GetEntries() <<
126  " Entries required: " << m_minEntries);
127  B2WARNING("Not enough data, adding one run to the collector");
128  f->Close();
129  gSystem->Unlink(Form("algorithm_3SampleCoG_output_rev_%d.root", cal_rev));
130  return c_NotEnoughData;
131  }
132  if (layer_num != 3 && hEventT0vsCoG->GetEntries() < m_minEntries / 10) {
133  B2INFO("Histogram: " << hEventT0vsCoG->GetName() <<
134  " Entries (n. clusters): " << hEventT0vsCoG->GetEntries() <<
135  " Entries required: " << m_minEntries / 10);
136  B2WARNING("Not enough data, adding one run to the collector");
137  f->Close();
138  gSystem->Unlink(Form("algorithm_3SampleCoG_output_rev_%d.root", cal_rev));
139  return c_NotEnoughData;
140  }
141  TF1* f1 = new TF1("f1", "[0]+[1]*x", -100, 100);
143  TF1* f2 = new TF1("f1", "[0]+[1]*x", -100, 100);
145  for (int i = 1; i <= hEventT0vsCoG->GetNbinsX(); i++) {
146  for (int j = 1; j <= hEventT0vsCoG->GetNbinsY(); j++) {
147  double bcx = ((TAxis*)hEventT0vsCoG->GetXaxis())->GetBinCenter(i);
148  double bcy = ((TAxis*)hEventT0vsCoG->GetYaxis())->GetBinCenter(j);
149  if (m_applyLinearCutsToRemoveBkg && (hEventT0vsCoG->GetBinContent(i, j) > 0 && (bcy > f1->Eval(bcx) || bcy < f2->Eval(bcx)))) {
150  hEventT0vsCoG->SetBinContent(i, j, 0);
151  } else if (hEventT0vsCoG->GetBinContent(i, j) > 0
152  && (hEventT0vsCoG->GetBinContent(i, j) < max(2, int(nEntriesForFilter * 0.001)))) {
153  hEventT0vsCoG->SetBinContent(i, j, 0);
154  }
155  }
156  }
157 
158  TProfile* pfx = hEventT0vsCoG->ProfileX();
159  std::string name = "pfx_" + std::string(hEventT0vsCoG->GetName());
160  pfx->SetName(name.c_str());
161  // set initial value for d, which is quadratic in the formula.
162  // also for b and c, with the nominal value in data.
163  pol3->SetParameter(1, 1.75);
164  pol3->SetParameter(2, 0.005);
165  pol3->SetParameter(3, 40);
166  TFitResultPtr tfr = pfx->Fit("pol3", "QMRS");
167  double par[4];
168  pol3->GetParameters(par);
169  // double meanT0 = hEventT0->GetMean();
170  // double meanT0NoSync = hEventT0nosync->GetMean();
171  timeCal->set_current(1);
172  timeCal->set_pol3parameters(par[0], par[1] + par[2]*par[3]*par[3], -par[2]*par[3], par[2] / 3);
173  payload->set(layer_num, ladder_num, sensor_num, bool(view), 1, *timeCal);
174  f->cd();
175  hEventT0->Write();
176  hEventT0vsCoG->Write();
177  hEventT0nosync->Write();
178  pfx->Write();
179 
180  delete pfx;
181  delete hEventT0vsCoG;
182  delete hEventT0;
183  delete hEventT0nosync;
184 
185 
186  if (tfr.Get() == nullptr || (tfr->Status() != 0 && tfr->Status() != 4 && tfr->Status() != 4000)) {
187  f->Close();
188  B2FATAL("Fit to the histogram failed in SVD3SampleCoGTimeCalibrationAlgorithm. "
189  << "Check the 2-D histogram to clarify the reason.");
190  } else {
191  a = par[0]; b = par[1]; c = par[2]; d = par[3];
192  a_err = tfr->ParError(0); b_err = tfr->ParError(1); c_err = tfr->ParError(2); d_err = tfr->ParError(3);
193  chi2 = tfr->Chi2();
194  ndf = tfr->Ndf();
195  p = tfr->Prob();
196  m_tree->Fill();
197  }
198  }
199  }
200  }
201  }
202  m_tree->Write();
203  f->Close();
204  saveCalibration(payload, "SVD3SampleCoGTimeCalibrations");
205 
206  //delete f;
207 
208  // probably not needed - would trigger re-doing the collection
209  // if ( ... too large corrections ... ) return c_Iterate;
210  return c_OK;
211 }
212 
213 bool SVD3SampleCoGTimeCalibrationAlgorithm::isBoundaryRequired(const Calibration::ExpRun& currentRun)
214 {
215  float meanRawTimeL3V = 0;
216  // auto eventT0Hist = getObjectPtr<TH1F>("hEventT0FromCDC");
217  auto rawTimeL3V = getObjectPtr<TH1F>("hRawTimeL3V");
218  // float meanEventT0 = eventT0Hist->GetMean();
219  if (!rawTimeL3V) {
221  meanRawTimeL3V = m_previousRawTimeMeanL3V.value();
222  } else {
223  if (rawTimeL3V->GetEntries() > m_minEntries)
224  meanRawTimeL3V = rawTimeL3V->GetMean();
225  else {
227  meanRawTimeL3V = m_previousRawTimeMeanL3V.value();
228  }
229  }
231  B2INFO("Setting start payload boundary to be the first run ("
232  << currentRun.first << "," << currentRun.second << ")");
233  m_previousRawTimeMeanL3V.emplace(meanRawTimeL3V);
234 
235  return true;
236  } else if (abs(meanRawTimeL3V - m_previousRawTimeMeanL3V.value()) > m_allowedTimeShift) {
237  B2INFO("Histogram mean has shifted from " << m_previousRawTimeMeanL3V.value()
238  << " to " << meanRawTimeL3V << ". We are requesting a new payload boundary for ("
239  << currentRun.first << "," << currentRun.second << ")");
240  m_previousRawTimeMeanL3V.emplace(meanRawTimeL3V);
241  return true;
242  } else {
243  return false;
244  }
245 }
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.
bool m_applyLinearCutsToRemoveBkg
if true turns on the selection to remove background clusters from CoG3 calibration
std::string m_id
Parameter given to set the UniqueID of the payload.
float m_angularCoefficientUpperLine
Angular coefficient of one of the two lines that define the signal region used in the CoG3 calibratio...
std::optional< float > m_previousRawTimeMeanL3V
Raw CoG of the previous run.
float m_interceptUpperLine
Intercept of one of the two lines that define the signal region used in the CoG3 calibration (to remo...
float m_minEntries
Set the minimun number of entries required in the histograms of layer 3.
float m_angularCoefficientLowerLine
Angulat coefficienct of one of the two lines that define the signal region used in the CoG3 calibrati...
virtual bool isBoundaryRequired(const Calibration::ExpRun &currentRun) override
If the event T0 changes significantly return true.
float m_interceptLowerLine
Intercept of one of the two lines that define the signal region used in the CoG3 calibration (to remo...
SVDCalibrationsBase< SVDCalibrationsScalar< SVDCoGCalibrationFunction > > t_payload
typedef for the SVDCoGCalibrationFunction payload of all SVD sensors
class to contain the CoG Time calibrations
Abstract base class for different kinds of events.