Belle II Software  release-05-02-19
SVD3SampleCoGTimeCalibrationAlgorithm.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2017 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Luigi Corona, Giulia Casarosa *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 #include <svd/calibration/SVD3SampleCoGTimeCalibrationAlgorithm.h>
11 
12 #include <svd/dbobjects/SVDCoGCalibrationFunction.h>
13 #include <svd/calibration/SVD3SampleCoGTimeCalibrations.h>
14 
15 #include <TF1.h>
16 #include <TProfile.h>
17 #include <TH2F.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 
36  gROOT->SetBatch(true);
37 
38  int ladderOfLayer[4] = {7, 10, 12, 16};
39  int sensorOnLayer[4] = {2, 3, 4, 5};
40 
41  auto timeCal = new Belle2::SVDCoGCalibrationFunction();
42  auto payload = new Belle2::SVD3SampleCoGTimeCalibrations::t_payload(*timeCal, m_id);
43 
44  std::unique_ptr<TF1> pol3(new TF1("pol3", "[0] + [1]*x + [2]*x*x + [3]*x*x*x", -10,
45  80)); // 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.) However, original value (-10,80) also seems working.
46 
47  FileStat_t info;
48  int cal_rev = 1;
49  while (gSystem->GetPathInfo(Form("algorithm_3SampleCoG_output_rev_%d.root", cal_rev), info) == 0)
50  cal_rev++;
51  std::unique_ptr<TFile> f(new TFile(Form("algorithm_3SampleCoG_output_rev_%d.root", cal_rev), "RECREATE"));
52 
53  auto m_tree = new TTree(Form("rev_%d", cal_rev), "RECREATE");
54  int layer_num, ladder_num, sensor_num, view, ndf;
55  float a, b, c, d, a_err, b_err, c_err, d_err, chi2, p;
56  m_tree->Branch("layer", &layer_num, "layer/I");
57  m_tree->Branch("ladder", &ladder_num, "ladder/I");
58  m_tree->Branch("sensor", &sensor_num, "sensor/I");
59  m_tree->Branch("isU", &view, "isU/I");
60  m_tree->Branch("a", &a, "a/F");
61  m_tree->Branch("b", &b, "b/F");
62  m_tree->Branch("c", &c, "c/F");
63  m_tree->Branch("d", &d, "d/F");
64  m_tree->Branch("a_err", &a_err, "a_err/F");
65  m_tree->Branch("b_err", &b_err, "b_err/F");
66  m_tree->Branch("c_err", &c_err, "c_err/F");
67  m_tree->Branch("d_err", &d_err, "d_err/F");
68  m_tree->Branch("chi2", &chi2, "chi2/F");
69  m_tree->Branch("ndf", &ndf, "ndf/I");
70  m_tree->Branch("p", &p, "p/F");
71 
72  for (int layer = 0; layer < 4; layer++) {
73  layer_num = layer + 3;
74  for (int ladder = 0; ladder < (int)ladderOfLayer[layer]; ladder++) {
75  ladder_num = ladder + 1;
76  for (int sensor = 0; sensor < (int)sensorOnLayer[layer]; sensor++) {
77  sensor_num = sensor + 1;
78  for (view = 1; view > -1; view--) {
79  char side = 'U';
80  if (view == 0)
81  side = 'V';
82  auto hEventT0vsCoG = getObjectPtr<TH2F>(Form("eventT0vsCoG__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
83  auto hEventT0 = getObjectPtr<TH1F>(Form("eventT0__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
84  auto hEventT0nosync = getObjectPtr<TH1F>(Form("eventT0nosync__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
85  B2INFO("Histogram: " << hEventT0vsCoG->GetName() <<
86  " Entries (n. clusters): " << hEventT0vsCoG->GetEntries());
87  if (layer_num == 3 && hEventT0vsCoG->GetEntries() < m_minEntries) {
88  B2INFO("Histogram: " << hEventT0vsCoG->GetName() <<
89  " Entries (n. clusters): " << hEventT0vsCoG->GetEntries() <<
90  " Entries required: " << m_minEntries);
91  B2WARNING("Not enough data, adding one run to the collector");
92  f->Close();
93  gSystem->Unlink(Form("algorithm_3SampleCoG_output_rev_%d.root", cal_rev));
94  return c_NotEnoughData;
95  }
96  for (int i = 1; i <= hEventT0vsCoG->GetNbinsX(); i++) {
97  for (int j = 1; j <= hEventT0vsCoG->GetNbinsY(); j++) {
98  if (hEventT0vsCoG->GetBinContent(i, j) < max(2, int(hEventT0vsCoG->GetEntries() * 0.001))) {
99  hEventT0vsCoG->SetBinContent(i, j, 0);
100  }
101  }
102  }
103  TProfile* pfx = hEventT0vsCoG->ProfileX();
104  std::string name = "pfx_" + std::string(hEventT0vsCoG->GetName());
105  pfx->SetName(name.c_str());
106  TFitResultPtr tfr = pfx->Fit("pol3", "QMRS");
107  double par[4];
108  pol3->GetParameters(par);
109  // double meanT0 = hEventT0->GetMean();
110  // double meanT0NoSync = hEventT0nosync->GetMean();
111  timeCal->set_current(1);
112  // timeCal->set_current(2);
113  timeCal->set_pol3parameters(par[0], par[1], par[2], par[3]);
114  payload->set(layer_num, ladder_num, sensor_num, bool(view), 1, *timeCal);
115  f->cd();
116  hEventT0->Write();
117  hEventT0vsCoG->Write();
118  hEventT0nosync->Write();
119  pfx->Write();
120 
121  a = par[0]; b = par[1]; c = par[2]; d = par[3];
122  a_err = tfr->ParError(0); b_err = tfr->ParError(1); c_err = tfr->ParError(2); d_err = tfr->ParError(3);
123  chi2 = tfr->Chi2();
124  ndf = tfr->Ndf();
125  p = tfr->Prob();
126  m_tree->Fill();
127 
128  }
129  }
130  }
131  }
132  m_tree->Write();
133  f->Close();
134  saveCalibration(payload, "SVD3SampleCoGTimeCalibrations");
135 
136  //delete f;
137 
138  // probably not needed - would trigger re-doing the collection
139  // if ( ... too large corrections ... ) return c_Iterate;
140  return c_OK;
141 }
142 
143 bool SVD3SampleCoGTimeCalibrationAlgorithm::isBoundaryRequired(const Calibration::ExpRun& currentRun)
144 {
145  float meanRawTimeL3V = 0;
146  // auto eventT0Hist = getObjectPtr<TH1F>("hEventT0FromCDC");
147  auto rawTimeL3V = getObjectPtr<TH1F>("hRawTimeL3V");
148  // float meanEventT0 = eventT0Hist->GetMean();
149  if (!rawTimeL3V) {
151  meanRawTimeL3V = m_previousRawTimeMeanL3V.value();
152  } else {
153  if (rawTimeL3V->GetEntries() > m_minEntries)
154  meanRawTimeL3V = rawTimeL3V->GetMean();
155  else {
157  meanRawTimeL3V = m_previousRawTimeMeanL3V.value();
158  }
159  }
161  B2INFO("Setting start payload boundary to be the first run ("
162  << currentRun.first << "," << currentRun.second << ")");
163  m_previousRawTimeMeanL3V.emplace(meanRawTimeL3V);
164 
165  return true;
166  } else if (abs(meanRawTimeL3V - m_previousRawTimeMeanL3V.value()) > m_allowedTimeShift) {
167  B2INFO("Histogram mean has shifted from " << m_previousRawTimeMeanL3V.value()
168  << " to " << meanRawTimeL3V << ". We are requesting a new payload boundary for ("
169  << currentRun.first << "," << currentRun.second << ")");
170  m_previousRawTimeMeanL3V.emplace(meanRawTimeL3V);
171  return true;
172  } else {
173  return false;
174  }
175 }
Belle2::CalibrationAlgorithm::saveCalibration
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
Definition: CalibrationAlgorithm.cc:290
Belle2::SVD3SampleCoGTimeCalibrationAlgorithm::m_previousRawTimeMeanL3V
std::optional< float > m_previousRawTimeMeanL3V
Raw CoG of the previous run.
Definition: SVD3SampleCoGTimeCalibrationAlgorithm.h:73
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::CalibrationAlgorithm::setDescription
void setDescription(const std::string &description)
Set algorithm description (in constructor)
Definition: CalibrationAlgorithm.h:331
Belle2::SVD3SampleCoGTimeCalibrationAlgorithm::m_allowedTimeShift
float m_allowedTimeShift
Allowed Raw CoGshift.
Definition: SVD3SampleCoGTimeCalibrationAlgorithm.h:74
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::SVD3SampleCoGTimeCalibrationAlgorithm::calibrate
virtual EResult calibrate() override
Run algo on data.
Definition: SVD3SampleCoGTimeCalibrationAlgorithm.cc:33
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::SVD3SampleCoGTimeCalibrations::t_payload
SVDCalibrationsBase< SVDCalibrationsScalar< SVDCoGCalibrationFunction > > t_payload
typedef for the SVDCoGCalibrationFunction payload of all SVD sensors
Definition: SVD3SampleCoGTimeCalibrations.h:45
Belle2::CalibrationAlgorithm::c_NotEnoughData
@ c_NotEnoughData
Needs more data =2 in Python.
Definition: CalibrationAlgorithm.h:53
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47
Belle2::SVDCoGCalibrationFunction
class to contain the CoG Time calibrations
Definition: SVDCoGCalibrationFunction.h:33
Belle2::SVD3SampleCoGTimeCalibrationAlgorithm::m_minEntries
float m_minEntries
Set the minimun number of entries required in the histograms of layer 3.
Definition: SVD3SampleCoGTimeCalibrationAlgorithm.h:75
Belle2::SVD3SampleCoGTimeCalibrationAlgorithm::isBoundaryRequired
virtual bool isBoundaryRequired(const Calibration::ExpRun &currentRun) override
If the event T0 changes significantly return true.
Definition: SVD3SampleCoGTimeCalibrationAlgorithm.cc:143
Belle2::SVD3SampleCoGTimeCalibrationAlgorithm::m_id
std::string m_id
Parameter given to set the UniqueID of the payload.
Definition: SVD3SampleCoGTimeCalibrationAlgorithm.h:72