Belle II Software development
LaserCalibratorFit.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/modules/TOPLaserCalibrator/LaserCalibratorFit.h>
10
11#include <framework/logging/Logger.h>
12
13#include <Math/PdfFuncMathCore.h>
14#include <Math/MinimizerOptions.h>
15#include <TFile.h>
16#include <TTree.h>
17
18namespace Belle2 {
23 namespace TOP {
24 //Crystal Ball function
25 // cppcheck-suppress constParameter
26 double fcnCB(double* x, double* par)
27 {
28 double m1 = par[1];
29 double s1 = par[2];
30 double a1 = par[3];
31 double n1 = par[4];
32 double f1 = ROOT::Math::crystalball_pdf(x[0], a1, n1, s1, m1);
33 return par[0] * f1;
34 }
35
36 //double Crystal Ball function
37 // cppcheck-suppress constParameter
38 double fcnCB2(double* x, double* par)
39 {
40 double m1 = par[1];
41 double s1 = par[2];
42 double a1 = par[3];
43 double n1 = par[4];
44 double m2 = m1 - par[6];
45 double s2 = par[7];
46 double a2 = par[8];
47 double n2 = par[9];
48
49 double f1 = ROOT::Math::crystalball_pdf(x[0], a1, n1, s1, m1);
50 double f2 = ROOT::Math::crystalball_pdf(x[0], a2, n2, s2, m2);
51 return par[0] * f1 + par[5] * f2;
52 }
53
55 m_moduleID(moduleID)
56 {
57 }
58
60 {
61 for (auto f : m_func) {
62 delete f;
63 }
64 m_func.clear();
65
66 for (auto h : m_hist) {
67 delete h;
68 }
69 m_hist.clear();
70 }
71
72//initial TF1 of all channels in one slot
73 void LaserCalibratorFit::setHist(const std::vector<TH1F*>& hist)
74 {
75 m_hist = hist;
76 for (unsigned i = 0; i < hist.size(); i++) {
77 m_func.push_back(0);
78 m_maxpos.push_back(0);
79 m_maxpos_error.push_back(0);
80 }
81 }
82
83 double LaserCalibratorFit::getFitChisq(unsigned channel)
84 {
85 if (m_func[channel]) {
86 return m_func[channel]->GetChisquare();
87 } else {
88 return -9;
89 };
90 }
91
92 int LaserCalibratorFit::fitChannel(unsigned channel)
93 {
94 if (channel > m_hist.size()) {
95 B2WARNING("Wrong channel!");
96 return 0;
97 }
98 if (!m_hist[channel]) {
99 B2WARNING("No Evts. in channel " << channel);
100 return 0;
101 }
102 //require at least 10 evts for each fit
103 if (m_hist[channel]->GetEntries() < 10) {
104 B2WARNING("Too little Evts. for fitting (Evts. < 10)!");
105 return 0;
106 }
107 m_hist[channel]->SetAxisRange(m_xmin, m_xmax);
108 m_maxpos[channel] = m_hist[channel]->GetXaxis()->GetBinCenter(m_hist[channel]->GetMaximumBin());
109 m_hist[channel]->SetAxisRange(m_maxpos[channel] - 1.5, m_maxpos[channel] + 2);
110 m_maxpos_error[channel] = m_hist[channel]->GetXaxis()->GetBinWidth(m_hist[channel]->GetMaximumBin()); // Error = bin width ?
111 if (m_fitMethod == "gauss" || m_fitMethod == "gaus") {
112 m_func[channel] = makeGFit(channel);
113 } else if (m_fitMethod == "cb") {
114 m_func[channel] = makeCBFit(channel);
115 } else if (m_fitMethod == "cb2") {
116 m_func[channel] = makeCB2Fit(channel, 0);
117 } else {
118 B2WARNING("No matched fit method!");
119 return 0;
120 }
121
122 return 1;
123 }
124
125 void LaserCalibratorFit::writeFile(const std::string& outfile)
126 {
127 auto file = new TFile(outfile.c_str(), "RECREATE");
128 auto otree = new TTree("fits", "fitted times");
129
130 unsigned channel = 0;
131 double maxpos = 0;
132
133 otree->Branch("maxpos", &maxpos, "maxpos/D");
134
135 if (m_fitMethod == "gauss") {
136 double parms[3];
137 double parmErrs[3];
138 otree->Branch("channel", &channel, "channel/I");
139 otree->Branch("norm", &(parms[0]), "norm/D");
140 otree->Branch("time", &(parms[1]), "time/D");
141 otree->Branch("reso", &(parms[2]), "reso/D");
142 otree->Branch("normErr", &(parmErrs[0]), "normErr/D");
143 otree->Branch("timeErr", &(parmErrs[1]), "timeErr/D");
144 otree->Branch("resoErr", &(parmErrs[2]), "resoErr/D");
145
146 for (auto& f : m_func) {
147 maxpos = m_maxpos[channel];
148 if (f) {
149 f->GetParameters(parms);
150 for (int iPar = 0; iPar < 3; iPar++)
151 parmErrs[iPar] = f->GetParError(iPar); // cppcheck-suppress unreadVariable
152 otree->Fill();
153 channel++;
154 }
155 }
156 } else if (m_fitMethod == "cb2") {
157 double parms[8];
158 double parmErrs[8];
159
160 double time2 = 0;
161 otree->Branch("channel", &channel, "channel/I");
162 otree->Branch("norm1", &(parms[0]), "norm1/D");
163 otree->Branch("time1", &(parms[1]), "time1/D");
164 otree->Branch("reso1", &(parms[2]), "reso1/D");
165 otree->Branch("alpha_CB", &(parms[3]), "alpha_CB/D");
166 otree->Branch("n_CB", &(parms[4]), "n_CB/D");
167 otree->Branch("norm2", &(parms[5]), "norm2/D");
168 otree->Branch("dt12", &(parms[6]), "dt12/D");
169 otree->Branch("reso2", &(parms[7]), "reso2/D");
170 otree->Branch("time2", &time2, "time2/D");
171
172 otree->Branch("norm1Err", &(parmErrs[0]), "norm1Err/D");
173 otree->Branch("time1Err", &(parmErrs[1]), "time1Err/D");
174 otree->Branch("reso1Err", &(parmErrs[2]), "reso1Err/D");
175 otree->Branch("alpha_CBErr", &(parmErrs[3]), "alpha_CBErr/D");
176 otree->Branch("n_CBErr", &(parmErrs[4]), "n_CBErr/D");
177 otree->Branch("norm2Err", &(parmErrs[5]), "norm2Err/D");
178 otree->Branch("dt12Err", &(parmErrs[6]), "dt12Err/D");
179 otree->Branch("reso2Err", &(parmErrs[7]), "reso2Err/D");
180
181
182 for (auto& f : m_func) {
183 maxpos = m_maxpos[channel];
184 if (f) {
185 f->GetParameters(parms);
186 for (int iPar = 0; iPar < 8; iPar++)
187 parmErrs[iPar] = f->GetParError(iPar); // cppcheck-suppress unreadVariable
188 otree->Fill();
189 channel++;
190 }
191 }
192 } else if (m_fitMethod == "cb") {
193 double parms[5];
194 double parmErrs[5];
195
196 otree->Branch("channel", &channel, "channel/I");
197 otree->Branch("norm", &(parms[0]), "norm/D");
198 otree->Branch("time", &(parms[1]), "time/D");
199 otree->Branch("reso", &(parms[2]), "reso/D");
200 otree->Branch("alpha_CB", &(parms[3]), "alpha_CB/D");
201 otree->Branch("n_CB", &(parms[4]), "n_CB/D");
202
203 otree->Branch("normErr", &(parmErrs[0]), "normErr/D");
204 otree->Branch("timeErr", &(parmErrs[1]), "timeErr/D");
205 otree->Branch("resoErr", &(parmErrs[2]), "resoErr/D");
206 otree->Branch("alpha_CBErr", &(parmErrs[3]), "alpha_CBErr/D");
207 otree->Branch("n_CBErr", &(parmErrs[4]), "n_CBErr/D");
208
209 for (auto& f : m_func) {
210 maxpos = m_maxpos[channel];
211 if (f) {
212 f->GetParameters(parms);
213 for (int iPar = 0; iPar < 5; iPar++)
214 parmErrs[iPar] = f->GetParError(iPar); // cppcheck-suppress unreadVariable
215 otree->Fill();
216 channel++;
217 }
218 }
219 } else {
220 B2WARNING("Nothing be written to files!");
221 return;
222 }
223 for (auto& h : m_hist) {
224 if (h && h->GetEntries() >= 10) {
225 h->Write();
226 }
227 }
228 otree->Write();
229 delete otree;
230 file->Close();
231 delete file;
232 }
233
234// single Gaussian fit
235 TF1* LaserCalibratorFit::makeGFit(unsigned channel)
236 {
237 TH1F* h = m_hist[channel];
238 double m = h->GetMean();
239 double w = h->GetRMS();
240 h->Fit("gaus", "Q", "", m - 2.*w, m + 4.*w);
241 double parms[3];
242 TF1* func = h->GetFunction("gaus");
243 func->GetParameters(parms);
244 func->SetNpx(1000);
245 h->Fit("gaus", "Q", "", m - w, m + w);
246 func = h->GetFunction("gaus");
247 m_fitT = parms[1];
248 m_fitTErr = func->GetParError(1);
249 return func;
250 }
251
252// single Crystal Ball fit
253 TF1* LaserCalibratorFit::makeCBFit(unsigned channel)
254 {
255 TH1F* h = m_hist[channel];
256 auto func = new TF1("fcnCB", fcnCB, m_xmin, m_xmax, 5);
257
258 double parms[5];
259 parms[0] = 1;
260 parms[1] = h->GetMean();
261 parms[2] = h->GetRMS();
262 parms[3] = -0.5;
263 parms[4] = 4;
264
265 func->SetParameters(parms);
266 func->SetNpx(1000);
267
268 func->SetParName(0, "norm");
269 func->SetParName(1, "time");
270 func->SetParName(2, "sigma");
271 func->SetParName(3, "alpha_CB");
272 func->SetParName(4, "n_CB");
273
274 //func->SetParLimits(1,par[1]-0.2,par[1]+0.2);
275 //func->SetParLimits(2,par[2]-0.02,par[2]+0.08);
276 func->SetParLimits(3, -5, 0);
277 func->SetParLimits(4, 0, 10);
278
279 h->Fit(func, "Q", "");
280 if (fabs(m_maxpos[channel] - parms[1]) < parms[2]) {
281 m_fitT = parms[1];
282 m_fitTErr = func->GetParError(1);
283 } else {
284 m_fitT = m_maxpos[channel];
285 m_fitTErr = m_maxpos_error[channel];
286 }
287 return func;
288 }
289
290// double Crystal Ball fit
291 TF1* LaserCalibratorFit::makeCB2Fit(unsigned channel, bool minOut)
292 {
293 TH1F* h = m_hist[channel];
294 auto func = new TF1("fcnCB2", fcnCB2, m_xmin, m_xmax, 8);
295
296 const double vdt[8] = {0.272, 0.242, 0.208, 0.178, 0.113, 0.082, 0.0485, 0.017}; //input para. from MC study, need further studies
297 double parms[8];
298 //the input para. need further studies
299 parms[0] = 1;
300 parms[1] = 0.6;
301 parms[2] = 0.1;
302 parms[3] = -0.5;
303 parms[4] = 4;
304 parms[5] = 1;
305 parms[6] = vdt[(channel) / 64]; //refers to a typical time separation of two main peaks
306 //parms[6] = 0.3;
307 parms[7] = 0.1;
308
309 func->SetParameters(parms);
310 func->SetNpx(1000);
311
312 func->SetParName(0, "norm1");
313 func->SetParName(1, "time1");
314 func->SetParName(2, "sigma1");
315 func->SetParName(3, "alpha1_CB");
316 func->SetParName(4, "n1_CB");
317 func->SetParName(5, "norm2");
318 func->SetParName(6, "dt12");
319 func->SetParName(7, "sigma2");
320 func->SetParName(8, "alpha2_CB");
321 func->SetParName(9, "n2_CB");
322
323 func->SetParLimits(1, parms[1] - 0.25, parms[1] + 0.15);
324 //func->SetParLimits(2,par[2]-0.02,par[2]+0.08);
325 func->SetParLimits(3, -5, 0);
326 func->SetParLimits(4, 0, 10);
327 func->SetParLimits(6, parms[6] * 0.7, parms[6] * 1.5);
328 func->SetParLimits(8, -5, 0);
329 func->SetParLimits(9, 0, 10);
330
331 ROOT::Math::MinimizerOptions::SetDefaultMaxFunctionCalls(50000); //set to a large number for a test
332 if (minOut) {
333 h->Fit(func, "Q", "");
334 } else {
335 h->Fit(func);
336 }
337 m_fitT = parms[1];
338 m_fitTErr = func->GetParError(1);
339 return func;
340 }
341 } // TOP namespace
343}// Belle2 namespace
TF1 * makeCB2Fit(unsigned channel, bool minOut)
Fit process using double Crystal Ball fuction.
void setHist(const std::vector< TH1F * > &hist)
set time hist of all channels in one moduleID
std::vector< TH1F * > m_hist
time hist of 512 channels
TF1 * makeCBFit(unsigned channel)
Fit process using single Crystal Ball fuction.
std::vector< TF1 * > m_func
fitting function
double m_fitTErr
error on the mean position estimated by the fit
LaserCalibratorFit(unsigned moduleID)
Constructor.
double m_xmax
fitting upper-edge
std::vector< double > m_maxpos_error
error on the center positon of hist max bin
int fitChannel(unsigned channel)
fit for a specific channel
void writeFile(const std::string &outfile)
write fit result to a root file
double m_fitT
mean position after fit
std::string m_fitMethod
fitting method
TF1 * makeGFit(unsigned channel)
Fit process using single gaussian function.
std::vector< double > m_maxpos
center positon of hist max bin
double getFitChisq(unsigned channel)
get chi^2 in the fit
Abstract base class for different kinds of events.