9#include <top/calibration/TOPModuleT0DeltaTAlgorithm.h>
10#include <top/dbobjects/TOPCalModuleT0.h>
16#include <TMatrixDSym.h>
17#include <TDecompChol.h>
43 auto slotPairs = getObjectPtr<TH2F>(
"slots");
45 B2ERROR(
"TOPModuleT0DeltaTAlgorithm: histogram 'slots' not found");
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");
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");
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);
82 if (status != 0)
continue;
86 std::vector<double> a(16, 0);
90 chi2Fits->SetBinContent(i + 1, k + 1,
m_chi2 /
m_ndf);
96 A.push_back(std::vector<double>(16, 1));
98 sigma.push_back(0.001);
102 int m = deltaT0.size();
107 B2INFO(
"TOPModuleT0DeltaTAlgorithm: NDF < 0");
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]);
125 bool ok = C.Invert(B);
129 B2INFO(
"TOPModuleT0DeltaTAlgorithm: matrix is not positive definite");
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]);
144 std::vector<double> x(16, 0);
145 std::vector<double> e(16, 0);
146 for (
int i = 0; i < 16; i++) {
147 for (
int j = 0; j < 16; j++) {
148 x[i] += B[i][j] * b[j];
150 e[i] =
sqrt(B[i][i]);
156 for (
int k = 0; k < m; k++) {
158 for (
int i = 0; i < 16; i++) {
161 chi2 += pow((s - deltaT0[k]) / sigma[k], 2);
168 ss <<
"Module T0, chi^2/NDF = " << chi2 <<
"/" << ndf;
169 string title = ss.str();
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]);
193 for (
int i = 0; i < 16; i++) {
194 moduleT0->setT0(i + 1, x[i], e[i]);
196 moduleT0->suppressAverage();
205 int numEntries = h->GetSumOfWeights();
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();
226 auto* func =
new TF1(
"func1g",
227 "[0] + [1]/sqrt(2*pi)/[3]*exp(-0.5*((x-[2])/[3])**2)",
229 func->SetParameter(0, 0);
230 func->SetParameter(1, sum * binWidth);
231 func->SetParameter(2, maxPosition);
234 int status = h->Fit(func,
"LERSQ");
236 m_delT0 = func->GetParameter(2);
237 m_error = func->GetParError(2);
238 m_chi2 = func->GetChisquare();
239 m_ndf = func->GetNDF();
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();
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))",
258 func->SetParameter(0, 0);
259 func->SetParameter(1, sum * binWidth);
260 func->SetParameter(2, maxPosition);
265 int status = h->Fit(func,
"LERSQ");
267 m_delT0 = func->GetParameter(2);
268 m_error = func->GetParError(2);
269 m_chi2 = func->GetChisquare();
270 m_ndf = func->GetNDF();
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 successfully =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
double m_delT0
fitted delta T0
virtual EResult calibrate() final
algorithm implementation
int fitHistogram(std::shared_ptr< TH1F > h)
Fit histogram.
double m_chi2
chi2 of the fit
int m_cutoffEntries
cutoff entries for single/double gaussian fit
double m_sigmaTailInit
tail gaussian sigma [ns]
double m_minError
minimal moduleT0 uncertainty [ns] to declare c_OK
double m_ndf
NDF of the fit.
double m_tailFractInit
fraction of tail gaussian
double m_error
error on fitted delta T0
double m_sigmaCoreInit
core gaussian sigma [ns]
TOPModuleT0DeltaTAlgorithm()
Constructor.
int fitSingleGaus(std::shared_ptr< TH1F > h)
Fit single gaus + constant.
double sqrt(double a)
sqrt for double
Abstract base class for different kinds of events.