Belle II Software release-09-00-00
SVD3SampleELSTimeCalibrationAlgorithm.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/SVD3SampleELSTimeCalibrationAlgorithm.h>
9
10#include <svd/dbobjects/SVDCoGCalibrationFunction.h>
11#include <svd/calibration/SVD3SampleELSTimeCalibrations.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
23using namespace std;
24using namespace Belle2;
25
27 CalibrationAlgorithm("SVDTimeCalibrationCollector")
28{
29 setDescription("SVD3SampleELSTimeCalibration calibration algorithm");
30 m_id = str;
31}
32
34{
35
36 gROOT->SetBatch(true);
37
38 auto timeCal = new Belle2::SVDCoGCalibrationFunction();
39 auto payload = new Belle2::SVD3SampleELSTimeCalibrations::t_payload(*timeCal, m_id);
40
41 std::unique_ptr<TF1> pol1pole(new TF1("pol1pole", "[0] + [1]*x + [2]/(x - [3])", -60,
42 0)); // In the local study, Y. Uematsu tuned the range to (-21.5,0). Original value is (-10,80).
43 pol1pole->SetParameter(1, 0.6);
44 pol1pole->SetParameter(2, -150);
45 pol1pole->SetParLimits(1, 0, 1);
46 pol1pole->SetParLimits(2, -300, 0);
47
48 FileStat_t info;
49 int cal_rev = 1;
50 while (gSystem->GetPathInfo(Form("algorithm_3SampleELS_output_rev_%d.root", cal_rev), info) == 0)
51 cal_rev++;
52 std::unique_ptr<TFile> f(new TFile(Form("algorithm_3SampleELS_output_rev_%d.root", cal_rev), "RECREATE"));
53
54 auto m_tree = new TTree(Form("rev_%d", cal_rev), "RECREATE");
55 int layer_num, ladder_num, sensor_num, view, ndf;
56 float a, b, c, d, a_err, b_err, c_err, d_err, chi2, p;
57 m_tree->Branch("layer", &layer_num, "layer/I");
58 m_tree->Branch("ladder", &ladder_num, "ladder/I");
59 m_tree->Branch("sensor", &sensor_num, "sensor/I");
60 m_tree->Branch("isU", &view, "isU/I");
61 m_tree->Branch("a", &a, "a/F");
62 m_tree->Branch("b", &b, "b/F");
63 m_tree->Branch("c", &c, "c/F");
64 m_tree->Branch("d", &d, "d/F");
65 m_tree->Branch("a_err", &a_err, "a_err/F");
66 m_tree->Branch("b_err", &b_err, "b_err/F");
67 m_tree->Branch("c_err", &c_err, "c_err/F");
68 m_tree->Branch("d_err", &d_err, "d_err/F");
69 m_tree->Branch("chi2", &chi2, "chi2/F");
70 m_tree->Branch("ndf", &ndf, "ndf/I");
71 m_tree->Branch("p", &p, "p/F");
72
73 auto __hEventT0vsCoG__ = getObjectPtr<TH3F>("__hEventT0vsCoG__");
74 auto __hEventT0__ = getObjectPtr<TH2F>("__hEventT0__");
75 auto __hEventT0NoSync__ = getObjectPtr<TH2F>("__hEventT0NoSync__");
76 auto __hBinToSensorMap__ = getObjectPtr<TH1F>("__hBinToSensorMap__");
77
78 for (int ij = 0; ij < (__hBinToSensorMap__->GetNbinsX()); ij++) {
79 {
80 {
81 {
82
83 auto binLabel = __hBinToSensorMap__->GetXaxis()->GetBinLabel(ij + 1);
84 char side;
85 std::sscanf(binLabel, "L%dL%dS%d%c", &layer_num, &ladder_num, &sensor_num, &side);
86 view = 0;
87 if (side == 'U')
88 view = 1;
89
90 B2INFO("Projecting for Sensor: " << binLabel << " with Bin Number: " << ij + 1);
91
92 __hEventT0vsCoG__->GetZaxis()->SetRange(ij + 1, ij + 1);
93 auto hEventT0vsELS = (TH2D*)__hEventT0vsCoG__->Project3D("yxe");
94 auto hEventT0 = (TH1D*)__hEventT0__->ProjectionX("hEventT0_tmp", ij + 1, ij + 1);
95 auto hEventT0nosync = (TH1D*)__hEventT0NoSync__->ProjectionX("hEventT0NoSync_tmp", ij + 1, ij + 1);
96
97 hEventT0vsELS->SetName(Form("eventT0vsCoG__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
98 hEventT0->SetName(Form("eventT0__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
99 hEventT0nosync->SetName(Form("eventT0nosync__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
100
101 char sidePN = (side == 'U' ? 'P' : 'N');
102 hEventT0vsELS->SetTitle(Form("EventT0Sync vs rawTime in %d.%d.%d %c/%c", layer_num, ladder_num, sensor_num, side, sidePN));
103 hEventT0->SetTitle(Form("EventT0Sync in %d.%d.%d %c/%c", layer_num, ladder_num, sensor_num, side, sidePN));
104 hEventT0nosync->SetTitle(Form("EventT0NoSync in %d.%d.%d %c/%c", layer_num, ladder_num, sensor_num, side, sidePN));
105
106 hEventT0vsELS->SetDirectory(0);
107 hEventT0->SetDirectory(0);
108 hEventT0nosync->SetDirectory(0);
109
110 B2INFO("Histogram: " << hEventT0vsELS->GetName() <<
111 " Entries (n. clusters): " << hEventT0vsELS->GetEntries());
112 if (layer_num == 3 && hEventT0vsELS->GetEntries() < m_minEntries) {
113 B2INFO("Histogram: " << hEventT0vsELS->GetName() <<
114 " Entries (n. clusters): " << hEventT0vsELS->GetEntries() <<
115 " Entries required: " << m_minEntries);
116 B2WARNING("Not enough data, adding one run to the collector");
117 f->Close();
118 gSystem->Unlink(Form("algorithm_3SampleELS_output_rev_%d.root", cal_rev));
119 return c_NotEnoughData;
120 }
121 if (layer_num != 3 && hEventT0vsELS->GetEntries() < m_minEntries / 10) {
122 B2INFO("Histogram: " << hEventT0vsELS->GetName() <<
123 " Entries (n. clusters): " << hEventT0vsELS->GetEntries() <<
124 " Entries required: " << m_minEntries / 10);
125 B2WARNING("Not enough data, adding one run to the collector");
126 f->Close();
127 gSystem->Unlink(Form("algorithm_3SampleELS_output_rev_%d.root", cal_rev));
128 return c_NotEnoughData;
129 }
130 for (int i = 1; i <= hEventT0vsELS->GetNbinsX(); i++) {
131 for (int j = 1; j <= hEventT0vsELS->GetNbinsY(); j++) {
132 if (hEventT0vsELS->GetBinContent(i, j) < max(2, int(hEventT0vsELS->GetEntries() * 0.001))) {
133 hEventT0vsELS->SetBinContent(i, j, 0);
134 }
135 }
136 }
137 TProfile* pfx = hEventT0vsELS->ProfileX();
138 std::string name = "pfx_" + std::string(hEventT0vsELS->GetName());
139 pfx->SetName(name.c_str());
140 // Remove non-liniarity from the fit by fixing the pole position
141 pol1pole->FixParameter(3, 7);
142 pfx->Fit("pol1pole", "QMRS");
143 pol1pole->SetParLimits(3, 0, 15);
144 TFitResultPtr tfr = pfx->Fit("pol1pole", "QMRS");
145 // There is small possibility that the chi2 minimization ends at the parameter boundary.
146 // In such case, we may get a completely wrong fit function, thus we refit w/o boundaries.
147 // if (tfr->Chi2() > 500) {
148 // B2WARNING("Correlation fit Chi2 is too large, refit with very wide parameter boundaries.");
149 // pol1pole->SetParLimits(1, 0, 0);
150 // pol1pole->SetParLimits(2, 0, 0);
151 // pol1pole->SetParLimits(3, 0, 0);
152 // tfr = pfx->Fit("pol1pole", "QMRS");
153 // pol1pole->SetParLimits(1, 0, 1);
154 // pol1pole->SetParLimits(2, -300, 0);
155 // }
156 double par[4];
157 pol1pole->GetParameters(par);
158 // double meanT0 = hEventT0->GetMean();
159 // double meanT0NoSync = hEventT0nosync->GetMean();
160 timeCal->set_current(3);
161 timeCal->set_elsparameters(par[0], par[1], par[2], par[3]);
162 payload->set(layer_num, ladder_num, sensor_num, bool(view), 1, *timeCal);
163 f->cd();
164 hEventT0->Write();
165 hEventT0vsELS->Write();
166 hEventT0nosync->Write();
167 pfx->Write();
168
169 delete pfx;
170 delete hEventT0vsELS;
171 delete hEventT0;
172 delete hEventT0nosync;
173
174 if (tfr.Get() == nullptr || (tfr->Status() != 0 && tfr->Status() != 4 && tfr->Status() != 4000)) {
175 f->Close();
176 B2FATAL("Fit to the histogram failed in SVD3SampleELSTimeCalibrationAlgorithm. "
177 << "Check the 2-D histogram to clarify the reason.");
178 } else {
179 a = par[0]; b = par[1]; c = par[2]; d = par[3];
180 a_err = tfr->ParError(0); b_err = tfr->ParError(1); c_err = tfr->ParError(2); d_err = tfr->ParError(3);
181 chi2 = tfr->Chi2();
182 ndf = tfr->Ndf();
183 p = tfr->Prob();
184 m_tree->Fill();
185 }
186 }
187 }
188 }
189 }
190 m_tree->Write();
191 f->Close();
192 saveCalibration(payload, "SVD3SampleELSTimeCalibrations");
193
194 //delete f;
195
196 // probably not needed - would trigger re-doing the collection
197 // if ( ... too large corrections ... ) return c_Iterate;
198 return c_OK;
199}
200
201bool SVD3SampleELSTimeCalibrationAlgorithm::isBoundaryRequired(const Calibration::ExpRun& currentRun)
202{
203 float meanRawTimeL3V = 0;
204 // auto eventT0Hist = getObjectPtr<TH1F>("hEventT0FromCDC");
205 auto rawTimeL3V = getObjectPtr<TH1F>("hRawTimeL3V");
206 // float meanEventT0 = eventT0Hist->GetMean();
207 if (!rawTimeL3V) {
209 meanRawTimeL3V = m_previousRawTimeMeanL3V.value();
210 } else {
211 if (rawTimeL3V->GetEntries() > m_minEntries)
212 meanRawTimeL3V = rawTimeL3V->GetMean();
213 else {
215 meanRawTimeL3V = m_previousRawTimeMeanL3V.value();
216 }
217 }
219 B2INFO("Setting start payload boundary to be the first run ("
220 << currentRun.first << "," << currentRun.second << ")");
221 m_previousRawTimeMeanL3V.emplace(meanRawTimeL3V);
222
223 return true;
224 } else if (abs(meanRawTimeL3V - m_previousRawTimeMeanL3V.value()) > m_allowedTimeShift) {
225 B2INFO("Histogram mean has shifted from " << m_previousRawTimeMeanL3V.value()
226 << " to " << meanRawTimeL3V << ". We are requesting a new payload boundary for ("
227 << currentRun.first << "," << currentRun.second << ")");
228 m_previousRawTimeMeanL3V.emplace(meanRawTimeL3V);
229 return true;
230 } else {
231 return false;
232 }
233}
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.
std::string m_id
Parameter given to set the UniqueID of the payload.
std::optional< float > m_previousRawTimeMeanL3V
Raw CoG of the previous run.
float m_minEntries
Set the minimun number of entries required in the histograms of layer 3.
virtual bool isBoundaryRequired(const Calibration::ExpRun &currentRun) override
If the event T0 changes significantly return true.
SVD3SampleELSTimeCalibrationAlgorithm(const std::string &str)
Constructor set the prefix to SVDTimeCalibrationCollector.
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.
STL namespace.