Belle II Software development
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
23using namespace Belle2;
24
26 CalibrationAlgorithm("SVDTimeCalibrationCollector")
27{
28 setDescription("SVD3SampleCoGTimeCalibration calibration algorithm");
29 m_id = str;
30}
31
33{
34 gROOT->SetBatch(true);
35
36 auto timeCal = new Belle2::SVDCoGCalibrationFunction();
37 auto payload = new Belle2::SVD3SampleCoGTimeCalibrations::t_payload(*timeCal, m_id);
38
39 // 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.)
40 // However, original value (-10,80) also seems working.
41 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,
42 80));
43 // apply parameter limits to be monotonic: [1] > 0 and [2] > 0. upper limits are loose on [2].
44 // lower limit on [1] is the minimum tilt (= tilt at the pole)
45 pol3->SetParLimits(1, 0, 100);
46 pol3->SetParLimits(2, 0, 0.1);
47 // the pole position [3] is different in the data (~40) and the current MC (~45)
48 // apply parameter limits on [3] to help non-linear minimization
49 pol3->SetParLimits(3, 30, 60);
50
51 FileStat_t info;
52 int cal_rev = 1;
53 while (gSystem->GetPathInfo(Form("algorithm_3SampleCoG_output_rev_%d.root", cal_rev), info) == 0)
54 cal_rev++;
55 std::unique_ptr<TFile> f(new TFile(Form("algorithm_3SampleCoG_output_rev_%d.root", cal_rev), "RECREATE"));
56
57 auto m_tree = new TTree(Form("rev_%d", cal_rev), "RECREATE");
58 int layer_num, ladder_num, sensor_num, view, ndf;
59 float a, b, c, d, a_err, b_err, c_err, d_err, chi2, p;
60 m_tree->Branch("layer", &layer_num, "layer/I");
61 m_tree->Branch("ladder", &ladder_num, "ladder/I");
62 m_tree->Branch("sensor", &sensor_num, "sensor/I");
63 m_tree->Branch("isU", &view, "isU/I");
64 m_tree->Branch("a", &a, "a/F");
65 m_tree->Branch("b", &b, "b/F");
66 m_tree->Branch("c", &c, "c/F");
67 m_tree->Branch("d", &d, "d/F");
68 m_tree->Branch("a_err", &a_err, "a_err/F");
69 m_tree->Branch("b_err", &b_err, "b_err/F");
70 m_tree->Branch("c_err", &c_err, "c_err/F");
71 m_tree->Branch("d_err", &d_err, "d_err/F");
72 m_tree->Branch("chi2", &chi2, "chi2/F");
73 m_tree->Branch("ndf", &ndf, "ndf/I");
74 m_tree->Branch("p", &p, "p/F");
75
77 B2INFO("--------- Applyingselection, 2D-region selection parameters: ");
78 B2INFO("Upper Line (q, m): " << m_interceptUpperLine << ", " << m_angularCoefficientUpperLine);
79 B2INFO("Lower Line (q, m): " << m_interceptLowerLine << ", " << m_angularCoefficientLowerLine);
80 } //B2INFO("Selecton applied : " << m_applyLinearCutsToRemoveBkg);
81
82 auto __hEventT0vsCoG__ = getObjectPtr<TH3F>("__hEventT0vsCoG__");
83 auto __hEventT0__ = getObjectPtr<TH2F>("__hEventT0__");
84 auto __hEventT0NoSync__ = getObjectPtr<TH2F>("__hEventT0NoSync__");
85 auto __hBinToSensorMap__ = getObjectPtr<TH1F>("__hBinToSensorMap__");
86
87 for (int ij = 0; ij < (__hBinToSensorMap__->GetNbinsX()); ij++) {
88 {
89 {
90 {
91
92 auto binLabel = __hBinToSensorMap__->GetXaxis()->GetBinLabel(ij + 1);
93 char side;
94 std::sscanf(binLabel, "L%dL%dS%d%c", &layer_num, &ladder_num, &sensor_num, &side);
95 view = 0;
96 if (side == 'U')
97 view = 1;
98
99 B2INFO("Projecting for Sensor: " << binLabel << " with Bin Number: " << ij + 1);
100
101 __hEventT0vsCoG__->GetZaxis()->SetRange(ij + 1, ij + 1);
102 auto hEventT0vsCoG = (TH2D*)__hEventT0vsCoG__->Project3D("yxe");
103 auto hEventT0 = (TH1D*)__hEventT0__->ProjectionX("hEventT0_tmp", ij + 1, ij + 1);
104 auto hEventT0nosync = (TH1D*)__hEventT0NoSync__->ProjectionX("hEventT0NoSync_tmp", ij + 1, ij + 1);
105
106 hEventT0vsCoG->SetName(Form("eventT0vsCoG__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
107 hEventT0->SetName(Form("eventT0__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
108 hEventT0nosync->SetName(Form("eventT0nosync__L%dL%dS%d%c", layer_num, ladder_num, sensor_num, side));
109
110 char sidePN = (side == 'U' ? 'P' : 'N');
111 hEventT0vsCoG->SetTitle(Form("EventT0Sync vs rawTime in %d.%d.%d %c/%c", layer_num, ladder_num, sensor_num, side, sidePN));
112 hEventT0->SetTitle(Form("EventT0Sync in %d.%d.%d %c/%c", layer_num, ladder_num, sensor_num, side, sidePN));
113 hEventT0nosync->SetTitle(Form("EventT0NoSync in %d.%d.%d %c/%c", layer_num, ladder_num, sensor_num, side, sidePN));
114
115 hEventT0vsCoG->SetDirectory(0);
116 hEventT0->SetDirectory(0);
117 hEventT0nosync->SetDirectory(0);
118
119 int nEntriesForFilter = hEventT0vsCoG->GetEntries();
120 B2INFO("Histogram: " << hEventT0vsCoG->GetName() <<
121 " Entries (n. clusters): " << hEventT0vsCoG->GetEntries());
122 if (layer_num == 3 && hEventT0vsCoG->GetEntries() < m_minEntries) {
123 B2INFO("Histogram: " << hEventT0vsCoG->GetName() <<
124 " Entries (n. clusters): " << hEventT0vsCoG->GetEntries() <<
125 " Entries required: " << m_minEntries);
126 B2WARNING("Not enough data, adding one run to the collector");
127 f->Close();
128 gSystem->Unlink(Form("algorithm_3SampleCoG_output_rev_%d.root", cal_rev));
129 return c_NotEnoughData;
130 }
131 if (layer_num != 3 && hEventT0vsCoG->GetEntries() < m_minEntries / 10) {
132 B2INFO("Histogram: " << hEventT0vsCoG->GetName() <<
133 " Entries (n. clusters): " << hEventT0vsCoG->GetEntries() <<
134 " Entries required: " << m_minEntries / 10);
135 B2WARNING("Not enough data, adding one run to the collector");
136 f->Close();
137 gSystem->Unlink(Form("algorithm_3SampleCoG_output_rev_%d.root", cal_rev));
138 return c_NotEnoughData;
139 }
140 TF1* f1 = new TF1("f1", "[0]+[1]*x", -100, 100);
142 TF1* f2 = new TF1("f1", "[0]+[1]*x", -100, 100);
144 for (int i = 1; i <= hEventT0vsCoG->GetNbinsX(); i++) {
145 for (int j = 1; j <= hEventT0vsCoG->GetNbinsY(); j++) {
146 double bcx = ((TAxis*)hEventT0vsCoG->GetXaxis())->GetBinCenter(i);
147 double bcy = ((TAxis*)hEventT0vsCoG->GetYaxis())->GetBinCenter(j);
148 if (m_applyLinearCutsToRemoveBkg && (hEventT0vsCoG->GetBinContent(i, j) > 0 && (bcy > f1->Eval(bcx) || bcy < f2->Eval(bcx)))) {
149 hEventT0vsCoG->SetBinContent(i, j, 0);
150 } else if (hEventT0vsCoG->GetBinContent(i, j) > 0
151 && (hEventT0vsCoG->GetBinContent(i, j) < std::max(2, int(nEntriesForFilter * 0.001)))) {
152 hEventT0vsCoG->SetBinContent(i, j, 0);
153 }
154 }
155 }
156
157 TProfile* pfx = hEventT0vsCoG->ProfileX();
158 std::string name = "pfx_" + std::string(hEventT0vsCoG->GetName());
159 pfx->SetName(name.c_str());
160 // set initial value for d, which is quadratic in the formula.
161 // also for b and c, with the nominal value in data.
162 pol3->SetParameter(1, 1.75);
163 pol3->SetParameter(2, 0.005);
164 pol3->SetParameter(3, 40);
165 TFitResultPtr tfr = pfx->Fit("pol3", "QMRS");
166 double par[4];
167 pol3->GetParameters(par);
168 // double meanT0 = hEventT0->GetMean();
169 // double meanT0NoSync = hEventT0nosync->GetMean();
170 timeCal->set_current(1);
171 timeCal->set_pol3parameters(par[0], par[1] + par[2]*par[3]*par[3], -par[2]*par[3], par[2] / 3);
172 payload->set(layer_num, ladder_num, sensor_num, bool(view), 1, *timeCal);
173 f->cd();
174 hEventT0->Write();
175 hEventT0vsCoG->Write();
176 hEventT0nosync->Write();
177 pfx->Write();
178
179 delete pfx;
180 delete hEventT0vsCoG;
181 delete hEventT0;
182 delete hEventT0nosync;
183
184
185 if (tfr.Get() == nullptr || (tfr->Status() != 0 && tfr->Status() != 4 && tfr->Status() != 4000)) {
186 f->Close();
187 B2FATAL("Fit to the histogram failed in SVD3SampleCoGTimeCalibrationAlgorithm. "
188 << "Check the 2-D histogram to clarify the reason.");
189 } else {
190 a = par[0]; b = par[1]; c = par[2]; d = par[3];
191 a_err = tfr->ParError(0); b_err = tfr->ParError(1); c_err = tfr->ParError(2); d_err = tfr->ParError(3);
192 chi2 = tfr->Chi2();
193 ndf = tfr->Ndf();
194 p = tfr->Prob();
195 m_tree->Fill();
196 }
197 }
198 }
199 }
200 }
201 m_tree->Write();
202 f->Close();
203 saveCalibration(payload, "SVD3SampleCoGTimeCalibrations");
204
205 //delete f;
206
207 // probably not needed - would trigger re-doing the collection
208 // if ( ... too large corrections ... ) return c_Iterate;
209 return c_OK;
210}
211
212bool SVD3SampleCoGTimeCalibrationAlgorithm::isBoundaryRequired(const Calibration::ExpRun& currentRun)
213{
214 float meanRawTimeL3V = 0;
215 // auto eventT0Hist = getObjectPtr<TH1F>("hEventT0FromCDC");
216 auto rawTimeL3V = getObjectPtr<TH1F>("hRawTimeL3V");
217 // float meanEventT0 = eventT0Hist->GetMean();
218 if (!rawTimeL3V) {
220 meanRawTimeL3V = m_previousRawTimeMeanL3V.value();
221 } else {
222 if (rawTimeL3V->GetEntries() > m_minEntries)
223 meanRawTimeL3V = rawTimeL3V->GetMean();
224 else {
226 meanRawTimeL3V = m_previousRawTimeMeanL3V.value();
227 }
228 }
230 B2INFO("Setting start payload boundary to be the first run ("
231 << currentRun.first << "," << currentRun.second << ")");
232 m_previousRawTimeMeanL3V.emplace(meanRawTimeL3V);
233
234 return true;
235 } else if (abs(meanRawTimeL3V - m_previousRawTimeMeanL3V.value()) > m_allowedTimeShift) {
236 B2INFO("Histogram mean has shifted from " << m_previousRawTimeMeanL3V.value()
237 << " to " << meanRawTimeL3V << ". We are requesting a new payload boundary for ("
238 << currentRun.first << "," << currentRun.second << ")");
239 m_previousRawTimeMeanL3V.emplace(meanRawTimeL3V);
240 return true;
241 } else {
242 return false;
243 }
244}
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 successfully =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.
SVD3SampleCoGTimeCalibrationAlgorithm(const std::string &str)
Constructor set the prefix to SVDTimeCalibrationCollector.
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.