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