Belle II Software development
TOPModuleT0DeltaTAlgorithm.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
9#include <top/calibration/TOPModuleT0DeltaTAlgorithm.h>
10#include <top/dbobjects/TOPCalModuleT0.h>
11#include <math.h>
12#include <TROOT.h>
13#include <TFile.h>
14#include <TF1.h>
15#include <TH2F.h>
16#include <TMatrixDSym.h>
17#include <TDecompChol.h>
18#include <string>
19#include <vector>
20#include <sstream>
21
22using namespace std;
23
24namespace Belle2 {
29 namespace TOP {
30
32 CalibrationAlgorithm("TOPModuleT0DeltaTCollector")
33 {
34 setDescription("Calibration algorithm for method DeltaT");
35 }
36
38 {
39 gROOT->SetBatch();
40
41 // get basic histogram
42
43 auto slotPairs = getObjectPtr<TH2F>("slots");
44 if (not slotPairs) {
45 B2ERROR("TOPModuleT0DeltaTAlgorithm: histogram 'slots' not found");
46 return c_Failure;
47 }
48
49 // construct file name and open output root file
50
51 const auto& expRun = getRunList();
52 string expNo = to_string(expRun[0].first);
53 while (expNo.length() < 4) expNo.insert(0, "0");
54 string runNo = to_string(expRun[0].second);
55 while (runNo.length() < 5) runNo.insert(0, "0");
56 string outputFileName = "moduleT0_rough-e" + expNo + "-r" + runNo + ".root";
57 auto* file = TFile::Open(outputFileName.c_str(), "recreate");
58
59 slotPairs->Write();
60
61 // create vectors and a histogram to store the results of fits
62
63 std::vector<double> deltaT0;
64 std::vector<double> sigma;
65 std::vector<std::vector<double> > A;
66 auto* chi2Fits = new TH2F("chi2_of_fits", "normalized chi2 of succesfull fits",
67 16, 0.5, 16.5, 16, 0.5, 16.5);
68 chi2Fits->SetXTitle("first slot number");
69 chi2Fits->SetYTitle("second slot number");
70
71 // fit histograms of time differences and store the results
72
73 for (int i = 0; i < slotPairs->GetNbinsX(); i++) {
74 for (int k = 0; k < slotPairs->GetNbinsY(); k++) {
75 if (slotPairs->GetBinContent(i + 1, k + 1) < m_minEntries) continue;
76 int slot1 = slotPairs->GetXaxis()->GetBinCenter(i + 1);
77 int slot2 = slotPairs->GetYaxis()->GetBinCenter(k + 1);
78 string name = "deltaT0_" + to_string(slot1) + "-" + to_string(slot2);
79 auto h = getObjectPtr<TH1F>(name);
80 if (not h) continue;
81 int status = fitHistogram(h);
82 if (status != 0) continue;
83 h->Write();
84 deltaT0.push_back(m_delT0);
85 sigma.push_back(m_error);
86 std::vector<double> a(16, 0);
87 a[slot1 - 1] = 1;
88 a[slot2 - 1] = -1;
89 A.push_back(a);
90 chi2Fits->SetBinContent(i + 1, k + 1, m_chi2 / m_ndf);
91 }
92 }
93
94 // append the bound (sum of all calibration constants equals to 0)
95
96 A.push_back(std::vector<double>(16, 1));
97 deltaT0.push_back(0);
98 sigma.push_back(0.001); // arbitrary but small compared to calibration precision
99
100 // check degrees of freedom and return if not positive
101
102 int m = deltaT0.size();
103 int ndf = m - 16;
104 if (ndf <= 0) {
105 file->Write();
106 file->Close();
107 B2INFO("TOPModuleT0DeltaTAlgorithm: NDF < 0");
108 return c_NotEnoughData;
109 }
110
111 // construct the matrix of a linear system of equations
112
113 TMatrixDSym B(16);
114 for (int i = 0; i < 16; i++) {
115 for (int j = 0; j < 16; j++) {
116 for (int k = 0; k < m; k++) {
117 B[i][j] += A[k][i] * A[k][j] / (sigma[k] * sigma[k]);
118 }
119 }
120 }
121
122 // invert the matrix and return if not positive definite
123
124 TDecompChol C(B); // Choletski decomposition to check also for positive definitness
125 bool ok = C.Invert(B);
126 if (not ok) {
127 file->Write();
128 file->Close();
129 B2INFO("TOPModuleT0DeltaTAlgorithm: matrix is not positive definite");
130 return c_NotEnoughData;
131 }
132
133 // construct the right side of a linear system of equations
134
135 std::vector<double> b(16, 0);
136 for (int i = 0; i < 16; i++) {
137 for (int k = 0; k < m; k++) {
138 b[i] += A[k][i] * deltaT0[k] / (sigma[k] * sigma[k]);
139 }
140 }
141
142 // solve for unknown module T0's
143
144 std::vector<double> x(16, 0); // module T0's
145 std::vector<double> e(16, 0); // uncertainties on module T0
146 for (int i = 0; i < 16; i++) {
147 for (int j = 0; j < 16; j++) {
148 x[i] += B[i][j] * b[j];
149 }
150 e[i] = sqrt(B[i][i]);
151 }
152
153 // calculate chi^2
154
155 double chi2 = 0;
156 for (int k = 0; k < m; k++) {
157 double s = 0;
158 for (int i = 0; i < 16; i++) {
159 s += A[k][i] * x[i];
160 }
161 chi2 += pow((s - deltaT0[k]) / sigma[k], 2);
162 }
163
164 // store the results in a histogram
165
166 stringstream ss;
167 ss.precision(3);
168 ss << "Module T0, chi^2/NDF = " << chi2 << "/" << ndf;
169 string title = ss.str();
170
171 auto* h_moduleT0 = new TH1F("moduleT0", title.c_str(), 16, 0.5, 16.5);
172 h_moduleT0->SetXTitle("slot number");
173 h_moduleT0->SetYTitle("module T0 [ns]");
174 for (int i = 0; i < 16; i++) {
175 h_moduleT0->SetBinContent(i + 1, x[i]);
176 h_moduleT0->SetBinError(i + 1, e[i]);
177 }
178
179 // write the results and close the file
180
181 file->Write();
182 file->Close();
183
184 // check the results and return if not good enough
185
186 for (auto err : e) {
187 if (err > m_minError) return c_NotEnoughData;
188 }
189
190 // otherwise create and import payload to DB
191
192 auto* moduleT0 = new TOPCalModuleT0();
193 for (int i = 0; i < 16; i++) {
194 moduleT0->setT0(i + 1, x[i], e[i]);
195 }
196 moduleT0->suppressAverage();
197 saveCalibration(moduleT0);
198
199 return c_OK;
200 }
201
202
203 int TOPModuleT0DeltaTAlgorithm::fitHistogram(std::shared_ptr<TH1F> h)
204 {
205 int numEntries = h->GetSumOfWeights();
206 int status = 0;
207 if (numEntries > m_cutoffEntries and m_tailFractInit > 0) {
208 status = fitDoubleGaus(h);
209 if (status != 0) status = fitSingleGaus(h);
210 } else {
211 status = fitSingleGaus(h);
212 }
213 return status;
214 }
215
216
217 int TOPModuleT0DeltaTAlgorithm::fitSingleGaus(std::shared_ptr<TH1F> h)
218 {
219 double sum = h->GetSumOfWeights();
220 if (sum < 5) return 5;
221 double maxPosition = h->GetBinCenter(h->GetMaximumBin());
222 double binWidth = h->GetBinWidth(1);
223 double xmin = h->GetXaxis()->GetXmin();
224 double xmax = h->GetXaxis()->GetXmax();
225
226 auto* func = new TF1("func1g",
227 "[0] + [1]/sqrt(2*pi)/[3]*exp(-0.5*((x-[2])/[3])**2)",
228 xmin, xmax);
229 func->SetParameter(0, 0);
230 func->SetParameter(1, sum * binWidth);
231 func->SetParameter(2, maxPosition);
232 func->SetParameter(3, m_sigmaCoreInit);
233
234 int status = h->Fit(func, "LERSQ");
235
236 m_delT0 = func->GetParameter(2);
237 m_error = func->GetParError(2);
238 m_chi2 = func->GetChisquare();
239 m_ndf = func->GetNDF();
240
241 return status;
242 }
243
244
245 int TOPModuleT0DeltaTAlgorithm::fitDoubleGaus(std::shared_ptr<TH1F> h)
246 {
247 double sum = h->GetSumOfWeights();
248 if (sum < 7) return 7;
249 double maxPosition = h->GetBinCenter(h->GetMaximumBin());
250 double binWidth = h->GetBinWidth(1);
251 double xmin = h->GetXaxis()->GetXmin();
252 double xmax = h->GetXaxis()->GetXmax();
253
254 auto* func = new TF1("func2g",
255 "[0] + [1]*((1-[4])/sqrt(2*pi)/[3]*exp(-0.5*((x-[2])/[3])**2)"
256 "+ [4]/sqrt(2*pi)/[5]*exp(-0.5*((x-[2])/[5])**2))",
257 xmin, xmax);
258 func->SetParameter(0, 0);
259 func->SetParameter(1, sum * binWidth);
260 func->SetParameter(2, maxPosition);
261 func->SetParameter(3, m_sigmaCoreInit);
262 func->SetParameter(4, m_tailFractInit);
263 func->SetParameter(5, m_sigmaTailInit);
264
265 int status = h->Fit(func, "LERSQ");
266
267 m_delT0 = func->GetParameter(2);
268 m_error = func->GetParError(2);
269 m_chi2 = func->GetChisquare();
270 m_ndf = func->GetNDF();
271
272 return status;
273 }
274
275 } // end namespace TOP
277} // end namespace Belle2
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)
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
@ c_Failure
Failed =3 in Python.
Module T0 calibration constants for all 16 modules.
int fitDoubleGaus(std::shared_ptr< TH1F > h)
Fit double gaus w/ same mean + constant.
int m_minEntries
minimal number of histogram entries to perform fit
virtual EResult calibrate() final
algorithm implementation
int fitHistogram(std::shared_ptr< TH1F > h)
Fit histogram.
int m_cutoffEntries
cutoff entries for single/double gaussian fit
double m_minError
minimal moduleT0 uncertainty [ns] to declare c_OK
double m_tailFractInit
fraction of tail gaussian
double m_sigmaCoreInit
core gaussian sigma [ns]
int fitSingleGaus(std::shared_ptr< TH1F > h)
Fit single gaus + constant.
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.
STL namespace.