Belle II Software  release-05-02-19
XTCalibrationAlgorithm.cc
1 #include <cdc/calibration/XTCalibrationAlgorithm.h>
2 #include <cdc/calibration/XTFunction.h>
3 #include <cdc/geometry/CDCGeometryPar.h>
4 #include <cdc/dbobjects/CDCXtRelations.h>
5 
6 #include <TError.h>
7 #include <TROOT.h>
8 #include <TProfile.h>
9 #include <TF1.h>
10 #include <TFile.h>
11 #include <TTree.h>
12 #include <iostream>
13 #include <memory>
14 
15 #include <framework/database/DBObjPtr.h>
16 
17 using namespace std;
18 using namespace Belle2;
19 using namespace CDC;
20 
21 typedef std::array<float, 3> array3;
22 XTCalibrationAlgorithm::XTCalibrationAlgorithm() : CalibrationAlgorithm("CDCCalibrationCollector")
23 {
25  " -------------------------- XT Calibration Algorithm -------------------------\n"
26  );
27 }
28 
30 {
31 
32  B2INFO("create and fill histo");
33 
34  /*Create histogram*/
35  for (int i = 0; i < 56; ++i) {
36  for (int lr = 0; lr < 2; ++lr) {
37  for (int al = 0; al < m_nAlphaBins; ++al) {
38  for (int th = 0; th < m_nThetaBins; ++th) {
39  m_hProf[i][lr][al][th] = new TProfile(Form("m_hProf%d_%d_%d_%d", i, lr, al, th),
40  Form("(L=%d)-(lr=%d)-(#alpha=%3.0f)-(#theta=%3.0f); Drift time (ns);Drift Length (cm)",
41  i, lr, m_iAlpha[al], m_iTheta[th]), 210, -20, 600, 0, 1.2, "i");
42  m_hist2d[i][lr][al][th] = new TH2F(Form("h%d_%d_%d_%d", i, lr, al, th),
43  Form("(L=%d)-(lr=%d)-(#alpha=%3.0f)-(#theta=%3.0f); Drift time (ns);Drift Length (cm)",
44  i, lr, m_iAlpha[al], m_iTheta[th]), 210, -20, 600, 110, 0, 1.2);
45  if (lr == 1)
46  m_hist2dDraw[i][al][th] = new TH2F(Form("h_draw%d_%d_%d", i, al, th),
47  Form("(L=%d)-(#alpha=%3.0f)-(#theta=%3.0f); Drift time (ns);Drift Length (cm)",
48  i, m_iAlpha[al], m_iTheta[th]), 210, -20, 600, 220, -1.2, 1.2);
49  }
50  }
51  }
52  }
53 
54  /* Read data and make histo*/
55 
56  auto tree = getObjectPtr<TTree>("tree");
57 
58  UChar_t lay;
59  Float_t dt;
60  Float_t dx;
61  Float_t Pval, alpha, theta;
62  Float_t ndf;
63 
64  tree->SetBranchAddress("lay", &lay);
65  tree->SetBranchAddress("t", &dt);
66  tree->SetBranchAddress("x_u", &dx);
67  tree->SetBranchAddress("alpha", &alpha);
68  tree->SetBranchAddress("theta", &theta);
69  tree->SetBranchAddress("Pval", &Pval);
70  tree->SetBranchAddress("ndf", &ndf);
71 
72  /* Disable unused branch */
73  std::vector<TString> list_vars = {"lay", "t", "x_u", "alpha", "theta", "Pval", "ndf"};
74  tree->SetBranchStatus("*", 0);
75 
76  for (TString brname : list_vars) {
77  tree->SetBranchStatus(brname, 1);
78  }
79 
80 
81  int al = 0;
82  int th = 0;
83  const Long64_t nEntries = tree->GetEntries();
84  B2INFO("Number of entries " << nEntries);
85  for (Long64_t i = 0; i < nEntries; ++i) {
86  tree->GetEntry(i);
87  /* protect in case |alpha|>90 */
88  if (fabs(alpha > 90)) {
89  if (alpha < 0) alpha += 180;
90  if (alpha > 0) alpha -= 180;
91  }
92 
93  if (Pval < m_minPval || ndf < m_minNdf) continue;
94 
95  for (int k = 0; k < m_nAlphaBins; ++k) {
96  if (alpha < m_upperAlpha[k]) {
97  al = k;
98  break;
99  }
100  }
101  for (int j = 0; j < m_nThetaBins; ++j) {
102  if (theta < m_upperTheta[j]) {
103  th = j;
104  break;
105  }
106  }
107  int lr = dx > 0 ? c_Right : c_Left;
108  if (m_LRseparate) {
109  m_hProf[lay][lr][al][th]->Fill(dt, abs(dx));
110  m_hist2d[lay][lr][al][th]->Fill(dt, abs(dx));
111  } else {
112  m_hProf[lay][0][al][th]->Fill(dt, abs(dx));
113  m_hist2d[lay][0][al][th]->Fill(dt, abs(dx));
114  m_hProf[lay][1][al][th]->Fill(dt, abs(dx));
115  m_hist2d[lay][1][al][th]->Fill(dt, abs(dx));
116  }
117  m_hist2dDraw[lay][al][th]->Fill(dt, dx);
118  }
119 }
120 
122 {
123  gROOT->SetBatch(1);
124  gPrintViaErrorHandler = true; // Suppress huge log output from TMinuit
125  gErrorIgnoreLevel = 3001;
126  B2INFO("Start calibration");
127 
128 
129  const auto exprun = getRunList()[0];
130  B2INFO("ExpRun used for DB Geometry : " << exprun.first << " " << exprun.second);
131  updateDBObjPtrs(1, exprun.second, exprun.first);
132  B2INFO("Creating CDCGeometryPar object");
134 
135  prepare();
136  createHisto();
137 
138  B2INFO("Start Fitting");
139  std::unique_ptr<XTFunction> xt;
140  for (int l = 0; l < 56; ++l) {
141  for (int lr = 0; lr < 2; ++lr) {
142  for (int al = 0; al < m_nAlphaBins; ++al) {
143  for (int th = 0; th < m_nThetaBins; ++th) {
144  if (m_hist2d[l][lr][al][th]->GetEntries() < m_minEntriesRequired) {
145  m_fitStatus[l][lr][al][th] = FitStatus::c_lowStat;
146  continue;
147  }
148  double p0, p1, tmin;
149  TF1* fpol1;
150  if (m_useSliceFit) {
151  m_hist2d[l][lr][al][th]->FitSlicesY(0, 0, -1, 5);
152  m_hist2d_1[l][lr][al][th] = (TH1F*)gDirectory->Get(Form("h%d_%d_%d_%d_1", l, lr, al, th));
153  if (!m_hist2d_1[l][lr][al][th]) {
154  m_fitStatus[l][lr][al][th] = FitStatus::c_lowStat;
155  B2WARNING("Error, not found results of slices fit");
156  continue;
157  }
158  m_hist2d_1[l][lr][al][th]->Fit("pol1", "Q", "", 30, 60);
159  fpol1 = (TF1*)m_hProf[l][lr][al][th]->GetFunction("pol1");
160  } else {
161  /*Set Error for low statistic bin*/
162  for (int n = 0; n < m_hProf[l][lr][al][th]->GetNbinsX(); ++n) {
163  if (m_hProf[l][lr][al][th]->GetBinEntries(n) < 5 && m_hProf[l][lr][al][th]->GetBinEntries(n) > 1) {
164  m_hProf[l][lr][al][th]->SetBinError(n, 0.3 / m_hProf[l][lr][al][th]->GetBinEntries(n));
165  }
166  }
167  m_hProf[l][lr][al][th]->Fit("pol1", "Q", "", 30, 60);
168  fpol1 = (TF1*)m_hProf[l][lr][al][th]->GetFunction("pol1");
169  }
170 
171  if (fpol1) {
172  //determine tmin in fitting
173  p0 = fpol1->GetParameter(0);
174  p1 = fpol1->GetParameter(1);
175  tmin = -1 * p0 / p1 + 15;
176  } else {
177  p0 = 0;
178  p1 = 0.005;
179  tmin = 12;
180  }
181 
182  // B2INFO("layer " << l << ", lr " << lr << ", alpha " << m_iAlpha[al] << ", theta " << m_iTheta[th]);
183  if (m_useSliceFit) { // if slice fit results exist.
184  xt.reset(new XTFunction(m_hist2d_1[l][lr][al][th], m_xtMode));
185  } else { // from TProfile.
186  xt.reset(new XTFunction(m_hProf[l][lr][al][th], m_xtMode));
187  }
188 
189  if (m_bField) {
190  int ial_old = 0;
191  int ith_old = 0;
192  for (int k = 0; k < m_nAlphaBins; ++k) {
193  if (m_iAlpha[al] < m_upperAlpha[k]) {
194  ial_old = k;
195  break;
196  }
197  }
198  for (int j = 0; j < m_nThetaBins; ++j) {
199  if (m_iTheta[th] < m_upperTheta[j]) {
200  ith_old = j;
201  break;
202  }
203  }
204 
205  double p6 = m_xtPrior[l][lr][ial_old][ith_old][6];
206  if (p6 > 400) {
207  p6 = 400;
208  }
209 
210  if (m_xtMode == m_xtModePrior) {
211  xt->setXTParams(m_xtPrior[l][lr][ial_old][ith_old]);
212  xt->setP6(p6);
213  } else {
214  xt->setXTParams(p0, p1, 0., 0., 0., 0., p6, m_xtPrior[l][lr][ial_old][ith_old][7]);
215  }
216  xt->setFitRange(tmin, p6 + 100);
217  } else {
218  xt->setXTParams(p0, p1, 0., 0., 0., 0., m_par6[l], 0.0001);
219  xt->setFitRange(tmin, m_par6[l] + 100);
220  }
221  xt->setDebug(m_debug);
222  xt->setBField(m_bField);
223  xt->fitXT();
224  if (xt->isValid() == false) {
225  B2WARNING("Empty xt");
226  m_fitStatus[l][lr][al][th] = c_fitFailure;
227  continue;
228  }
229  if (xt->getFitStatus() != 1) {
230  B2WARNING("Fit failed");
231  m_fitStatus[l][lr][al][th] = c_fitFailure;
232  continue;
233  }
234  if (xt->validate() == true) {
235  m_fitStatus[l][lr][al][th] = xt->getFitStatus();
236  m_xtFunc[l][lr][al][th] = (TF1*)xt->getXTFunction();
237 
238  if (m_useSliceFit) {
239  m_hist2d_1[l][lr][al][th] = (TH1F*)xt->getFittedHisto();
240  } else {
241  m_hProf[l][lr][al][th] = (TProfile*)xt->getFittedHisto();
242  }
243  } else {
244  m_fitStatus[l][lr][al][th] = c_fitFailure;
245  }
246  }
247  }
248  }
249  }
250  sanitaryCheck();
251  write();
252  storeHisto();
253  return checkConvergence();
254 }
255 
257 {
258  const double tMax = 500; // max drift time (nsec)
259  for (int l = 0; l < 56; ++l) {
260  for (int lr = 0; lr < 2; ++lr) {
261  for (int al = 0; al < m_nAlphaBins; ++al) {
262  for (int th = 0; th < m_nThetaBins; ++th) {
263  if (m_fitStatus[l][lr][al][th] == FitStatus::c_OK) {
264  TF1* fun = m_xtFunc[l][lr][al][th];
265  double y = fun->Eval(tMax);
266  if (y < 0) {
267  B2INFO("Strange XT function l " << l << " lr " << lr << " alpha " << al << " theta " << th
268  << ", replaced by initial one");
269  fun->SetParameters(m_xtPrior[l][lr][al][th]);
270  }
271  }
272  }
273  }
274  }
275  }
276 }
278 {
279 
280  const int nTotal = 56 * 2 * m_nAlphaBins * m_nThetaBins;
281  int nFitCompleted = 0;
282  for (int l = 0; l < 56; ++l) {
283  for (int lr = 0; lr < 2; ++lr) {
284  for (int al = 0; al < m_nAlphaBins; ++al) {
285  for (int th = 0; th < m_nThetaBins; ++th) {
286  if (m_fitStatus[l][lr][al][th] == FitStatus::c_OK) {
287  nFitCompleted++;
288  }
289  }
290  }
291  }
292  }
293 
294  if (static_cast<double>(nFitCompleted) / nTotal < m_threshold) {
295  B2WARNING("Less than " << m_threshold * 100 << " % of XTs were fitted.");
296  return c_NotEnoughData;
297  }
298  return c_OK;
299 }
300 
302 {
303  B2INFO("Prepare calibration of XT");
304  const double rad2deg = 180 / M_PI;
305 
307  m_nAlphaBins = dbXT->getNoOfAlphaBins();
308  m_nThetaBins = dbXT->getNoOfThetaBins();
309  for (unsigned short i = 0; i < m_nAlphaBins; ++i) {
310  array3 alpha = dbXT->getAlphaBin(i);
311  m_lowerAlpha[i] = alpha[0] * rad2deg;
312  m_upperAlpha[i] = alpha[1] * rad2deg;
313  m_iAlpha[i] = alpha[2] * rad2deg;
314  }
315 
316  for (unsigned short i = 0; i < m_nThetaBins; ++i) {
317  array3 theta = dbXT->getThetaBin(i);
318  m_lowerTheta[i] = theta[0] * rad2deg;
319  m_upperTheta[i] = theta[1] * rad2deg;
320  m_iTheta[i] = theta[2] * rad2deg;
321  }
322 
323  m_xtModePrior = dbXT->getXtParamMode();
324  if (!(m_xtModePrior == c_Chebyshev || m_xtModePrior == c_Polynomial)) {
325  B2FATAL("Function type before calibration is wrong " << m_xtModePrior);
326  }
327 
328  B2INFO("Number of alpha bins " << m_nAlphaBins);
329  B2INFO("Number of theta bins " << m_nThetaBins);
330  B2INFO("Function type " << m_xtMode);
331 
332  for (unsigned short iCL = 0; iCL < 56; ++iCL) {
333  for (unsigned short iLR = 0; iLR < 2; ++iLR) {
334  for (unsigned short iA = 0; iA < m_nAlphaBins; ++iA) {
335  for (unsigned short iT = 0; iT < m_nThetaBins; ++iT) {
336  const std::vector<float> params = dbXT->getXtParams(iCL, iLR, iA, iT);
337  unsigned short np = params.size();
338  for (unsigned short i = 0; i < np; ++i) {
339  m_xtPrior[iCL][iLR][iA][iT][i] = params[i];
340  }
341  }
342  }
343  }
344  }
345 }
346 
348 {
349  B2INFO("write calibrated XT");
350  double par[8];
351 
352 
353  int nfitted = 0;
354  int nfailure = 0;
355 
356  //
357  // Save to the localDB
358  //
359 
360  CDCXtRelations* xtRel = new CDCXtRelations();
361  const float deg2rad = static_cast<float>(Unit::deg);
362 
363  for (int i = 0; i < m_nAlphaBins; ++i) {
364  std::array<float, 3> alpha3 = {m_lowerAlpha[i]* deg2rad,
365  m_upperAlpha[i]* deg2rad,
366  m_iAlpha[i]* deg2rad
367  };
368  xtRel->setAlphaBin(alpha3);
369  }
370 
371  for (int i = 0; i < m_nThetaBins; ++i) {
372  std::array<float, 3> theta3 = {m_lowerTheta[i]* deg2rad,
373  m_upperTheta[i]* deg2rad,
374  m_iTheta[i]* deg2rad
375  };
376  xtRel->setThetaBin(theta3);
377  }
378 
379  xtRel->setXtParamMode(m_xtMode);
380 
381  for (int th = 0; th < m_nThetaBins; ++th) {
382  for (int al = 0; al < m_nAlphaBins; ++al) {
383  for (int l = 0; l < 56; ++l) {
384  for (int lr = 0; lr < 2; ++lr) {
385  if (m_fitStatus[l][lr][al][th] != FitStatus::c_OK) {
386  nfailure += 1;
387  B2DEBUG(21, "fit failure status = " << m_fitStatus[l][lr][al][th]);
388  B2DEBUG(21, "layer " << l << ", r " << lr << ", alpha " << m_iAlpha[al] << ", theta " << m_iTheta[th]);
389  B2DEBUG(21, "number of event: " << m_hProf[l][lr][al][th]->GetEntries());
390  if (m_fitStatus[l][lr][al][th] != FitStatus::c_lowStat) {
391  if (m_xtFunc[l][lr][al][th]) {
392  B2DEBUG(21, "Probability of fit: " << m_xtFunc[l][lr][al][th]->GetProb());
393  }
394  }
395  // If fit is failed
396  // and mode of input xt (prior) is same as output, previous xt is used.
397  if (m_xtMode == m_xtModePrior) {
398  for (int i = 0; i < 8; ++i) {
399  par[i] = m_xtPrior[l][lr][al][th][i];
400  }
401  } else {
402  B2FATAL("XT mode before/after calibration is different!");
403  }
404 
405  } else {
406  if (par[1] < 0) { // if negative c1, privious xt is kept.
407  for (int i = 0; i < 8; ++i) {
408  par[i] = m_xtPrior[l][lr][al][th][i];
409  }
410  } else {
411  m_xtFunc[l][lr][al][th]->GetParameters(par);
412  nfitted += 1;
413  }
414  }
415  std::vector<float> xtbuff;
416  for (int i = 0; i < 8; ++i) {
417  xtbuff.push_back(par[i]);
418  }
419  xtRel->setXtParams(l, lr, al, th, xtbuff);
420  }//lr
421  }//layer
422  }//alpha
423  }//theta
424 
425  if (m_textOutput == true) {
427  }
428 
429  saveCalibration(xtRel, "CDCXtRelations");
430 
431  B2RESULT("Total number of xt fit: " << m_nAlphaBins * m_nThetaBins * 2 * 56);
432  B2RESULT("Successfully Fitted: " << nfitted);
433  B2RESULT("Failure Fit: " << nfailure);
434 
435 }
436 
438 {
439 
440  auto hNDF = getObjectPtr<TH1F>("hNDF");
441  auto hPval = getObjectPtr<TH1F>("hPval");
442  auto hEvtT0 = getObjectPtr<TH1F>("hEventT0");
443  B2INFO("saving histograms");
444  TFile* fout = new TFile(m_histName.c_str(), "RECREATE");
445  TDirectory* top = gDirectory;
446  //store NDF, P-val. EventT0 histogram for monitoring during calibration
447  if (hNDF && hPval && hEvtT0) {
448  hEvtT0->Write();
449  hPval->Write();
450  hNDF->Write();
451  }
452  // for each layer
453 
454  TDirectory* Direct[56];
455  int nhisto = 0;
456  for (int l = 0; l < 56; ++l) {
457  top->cd();
458  Direct[l] = gDirectory->mkdir(Form("lay_%d", l));
459  Direct[l]->cd();
460  for (int th = 0; th < m_nThetaBins; ++th) {
461  for (int al = 0; al < m_nAlphaBins; ++al) {
462  m_hist2dDraw[l][al][th]->Write();
463  for (int lr = 0; lr < 2; ++lr) {
464  m_hist2d[l][lr][al][th]->Write();
465  if (m_fitStatus[l][lr][al][th] != 1) continue;
466  if (m_useSliceFit) {
467  if (m_hist2d_1[l][lr][al][th]) {
468  m_hist2d_1[l][lr][al][th]->Write();
469  nhisto += 1;
470  }
471  } else {
472  m_hProf[l][lr][al][th]->Write();
473  nhisto += 1;
474  }
475  }
476  }
477  }
478  }
479  top->cd();
480 
481  fout->Close();
482  B2RESULT(" " << nhisto << " histograms was stored.");
483 }
484 
Belle2::CDC::XTCalibrationAlgorithm::m_LRseparate
bool m_LRseparate
Separate LR in calibration or mix.
Definition: XTCalibrationAlgorithm.h:130
Belle2::CDC::XTFunction::setXTParams
void setXTParams(double p[8])
Set Paramerters for fit.
Definition: XTFunction.h:211
Belle2::CDCXtRelations
Database object for xt-relations.
Definition: CDCXtRelations.h:38
Belle2::CDC::XTFunction::isValid
bool isValid()
Is valid.
Definition: XTFunction.h:185
Belle2::CDC::XTFunction::getFitStatus
int getFitStatus()
get fitted flag.
Definition: XTFunction.h:253
Belle2::CDC::XTFunction
Class to perform fitting for each xt function.
Definition: XTFunction.h:79
Belle2::CDCXtRelations::setThetaBin
void setThetaBin(const array3 &theta)
Set theta-angle bin (rad)
Definition: CDCXtRelations.h:76
Belle2::CDC::XTCalibrationAlgorithm::m_lowerTheta
float m_lowerTheta[7]
Lower boundays of theta bins.
Definition: XTCalibrationAlgorithm.h:152
Belle2::CDC::XTCalibrationAlgorithm::calibrate
EResult calibrate() override
Run algo on data.
Definition: XTCalibrationAlgorithm.cc:121
Belle2::CDCXtRelations::setAlphaBin
void setAlphaBin(const array3 &alpha)
Set alpha-angle bin (rad)
Definition: CDCXtRelations.h:62
Belle2::CDC::XTCalibrationAlgorithm::prepare
void prepare()
Prepare the calibration of XT.
Definition: XTCalibrationAlgorithm.cc:301
Belle2::CalibrationAlgorithm::saveCalibration
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
Definition: CalibrationAlgorithm.cc:290
Belle2::CDC::XTCalibrationAlgorithm::m_xtFunc
TF1 * m_xtFunc[56][2][20][10]
XTFunction.
Definition: XTCalibrationAlgorithm.h:138
Belle2::CDC::XTCalibrationAlgorithm::m_xtMode
int m_xtMode
Mode of xt; 0 is polynomial;1 is Chebyshev.
Definition: XTCalibrationAlgorithm.h:147
Belle2::CDC::XTCalibrationAlgorithm::m_cdcGeo
DBObjPtr< CDCGeometry > m_cdcGeo
Geometry of CDC.
Definition: XTCalibrationAlgorithm.h:170
Belle2::CDC::XTCalibrationAlgorithm::m_hist2d
TH2F * m_hist2d[56][2][20][10]
2D histo of xt
Definition: XTCalibrationAlgorithm.h:135
Belle2::CDC::XTFunction::setFitRange
void setFitRange(double tmin, double tmax)
Set Fit range.
Definition: XTFunction.h:231
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::CDC::XTFunction::getXTFunction
TF1 * getXTFunction()
Get XT function.
Definition: XTFunction.h:275
Belle2::CDC::XTCalibrationAlgorithm::m_minNdf
double m_minNdf
minimum ndf required
Definition: XTCalibrationAlgorithm.h:126
Belle2::CalibrationAlgorithm::setDescription
void setDescription(const std::string &description)
Set algorithm description (in constructor)
Definition: CalibrationAlgorithm.h:331
Belle2::CalibrationAlgorithm::updateDBObjPtrs
void updateDBObjPtrs(const unsigned int event, const int run, const int experiment)
Updates any DBObjPtrs by calling update(event) for DBStore.
Definition: CalibrationAlgorithm.cc:398
Belle2::CDC::XTCalibrationAlgorithm::m_textOutput
bool m_textOutput
output text file if true
Definition: XTCalibrationAlgorithm.h:167
Belle2::CDC::XTCalibrationAlgorithm::m_nThetaBins
int m_nThetaBins
number of theta bins
Definition: XTCalibrationAlgorithm.h:146
Belle2::CDC::XTCalibrationAlgorithm::createHisto
void createHisto()
Create histogram for calibration.
Definition: XTCalibrationAlgorithm.cc:29
Belle2::CDC::XTFunction::setDebug
void setDebug(bool debug)
Set Debug.
Definition: XTFunction.h:246
Belle2::CDC::XTCalibrationAlgorithm::m_hist2d_1
TH1F * m_hist2d_1[56][2][20][10]
1D xt histo, results of slice fit
Definition: XTCalibrationAlgorithm.h:137
Belle2::CDC::XTCalibrationAlgorithm::m_par6
double m_par6[56]
boundary parameter for fitting, semi-experiment number
Definition: XTCalibrationAlgorithm.h:156
Belle2::DBObjPtr
Class for accessing objects in the database.
Definition: DBObjPtr.h:31
Belle2::CDC::XTCalibrationAlgorithm::m_bField
bool m_bField
with b field or none
Definition: XTCalibrationAlgorithm.h:131
Belle2::CDC::XTCalibrationAlgorithm::m_threshold
double m_threshold
minimal requirement for the fraction of fitted results
Definition: XTCalibrationAlgorithm.h:132
Belle2::CDC::XTCalibrationAlgorithm::m_hist2dDraw
TH2F * m_hist2dDraw[56][20][10]
2d histo for draw
Definition: XTCalibrationAlgorithm.h:136
Belle2::CDC::XTCalibrationAlgorithm::m_minPval
double m_minPval
minimum pvalue required
Definition: XTCalibrationAlgorithm.h:127
Belle2::CDC::XTCalibrationAlgorithm::m_hProf
TProfile * m_hProf[56][2][20][10]
Profile xt histo.
Definition: XTCalibrationAlgorithm.h:134
Belle2::CDC::CDCGeometryPar::Instance
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
Definition: CDCGeometryPar.cc:41
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::CDC::XTCalibrationAlgorithm::m_iAlpha
float m_iAlpha[18]
Represented alpha in alpha bins.
Definition: XTCalibrationAlgorithm.h:151
Belle2::CDC::XTCalibrationAlgorithm::m_lowerAlpha
float m_lowerAlpha[18]
Lower boundays of alpha bins.
Definition: XTCalibrationAlgorithm.h:149
Belle2::CDC::XTFunction::validate
bool validate()
Validate the xt has proper shape.
Definition: XTFunction.h:427
Belle2::CDCXtRelations::setXtParamMode
void setXtParamMode(unsigned short mode)
Set xt parameterization mode.
Definition: CDCXtRelations.h:97
Belle2::CDC::XTCalibrationAlgorithm::m_fitStatus
int m_fitStatus[56][2][20][10]
Fit flag.
Definition: XTCalibrationAlgorithm.h:142
Belle2::CDC::XTCalibrationAlgorithm::m_useSliceFit
bool m_useSliceFit
Use slice fit or profile.
Definition: XTCalibrationAlgorithm.h:143
Belle2::CDC::XTCalibrationAlgorithm::storeHisto
void storeHisto()
Store histogram to file.
Definition: XTCalibrationAlgorithm.cc:437
Belle2::CDC::XTFunction::setBField
void setBField(bool bfield)
set to use BField
Definition: XTFunction.h:206
Belle2::CDC::XTCalibrationAlgorithm::m_minEntriesRequired
int m_minEntriesRequired
minimum number of hit per hitosgram.
Definition: XTCalibrationAlgorithm.h:144
Belle2::Unit::deg
static const double deg
degree to radians
Definition: Unit.h:119
Belle2::CDCXtRelations::setXtParams
void setXtParams(const XtID xtID, const std::vector< float > &params)
Set xt parameters for the specified id.
Definition: CDCXtRelations.h:105
Belle2::CDC::XTCalibrationAlgorithm::m_upperAlpha
float m_upperAlpha[18]
Upper boundays of alpha bins.
Definition: XTCalibrationAlgorithm.h:150
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::CDC::XTCalibrationAlgorithm::m_xtModePrior
int m_xtModePrior
Mode of xt before calibration; 0 is polynomial;1 is Chebyshev.
Definition: XTCalibrationAlgorithm.h:148
Belle2::CDC::XTCalibrationAlgorithm::m_iTheta
float m_iTheta[7]
Represented theta in theta bins.
Definition: XTCalibrationAlgorithm.h:154
Belle2::CalibrationAlgorithm::getRunList
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
Definition: CalibrationAlgorithm.h:276
Belle2::CDC::XTFunction::fitXT
void fitXT()
Do fitting.
Definition: XTFunction.h:287
Belle2::CDC::XTCalibrationAlgorithm::checkConvergence
EResult checkConvergence()
Check the convergence of XT fit.
Definition: XTCalibrationAlgorithm.cc:277
Belle2::CalibrationAlgorithm::c_NotEnoughData
@ c_NotEnoughData
Needs more data =2 in Python.
Definition: CalibrationAlgorithm.h:53
Belle2::CDC::XTCalibrationAlgorithm::m_upperTheta
float m_upperTheta[7]
Upper boundays of theta bins.
Definition: XTCalibrationAlgorithm.h:153
Belle2::CDC::XTCalibrationAlgorithm::m_nAlphaBins
int m_nAlphaBins
number of alpha bins
Definition: XTCalibrationAlgorithm.h:145
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47
Belle2::CDC::XTCalibrationAlgorithm::m_debug
bool m_debug
run in debug or silent
Definition: XTCalibrationAlgorithm.h:128
Belle2::CDC::XTCalibrationAlgorithm::sanitaryCheck
void sanitaryCheck()
Check if there are any wrong xt functions.
Definition: XTCalibrationAlgorithm.cc:256
Belle2::CDC::XTCalibrationAlgorithm::m_histName
std::string m_histName
root file name
Definition: XTCalibrationAlgorithm.h:169
Belle2::CDC::XTCalibrationAlgorithm::m_xtPrior
double m_xtPrior[56][2][18][7][8]
paremeters of XT before calibration
Definition: XTCalibrationAlgorithm.h:140
Belle2::CDC::XTFunction::setP6
void setP6(double p6)
Set Parameter 6 for polynomia fit.
Definition: XTFunction.h:177
Belle2::CDC::XTCalibrationAlgorithm::m_outputFileName
std::string m_outputFileName
Output xt filename.
Definition: XTCalibrationAlgorithm.h:168
Belle2::CDC::XTFunction::getFittedHisto
TProfile * getFittedHisto()
Get histogram.
Definition: XTFunction.h:282
Belle2::CDCXtRelations::outputToFile
void outputToFile(std::string fileName) const
Output the contents in test file format.
Definition: CDCXtRelations.h:373
Belle2::CDC::XTCalibrationAlgorithm::write
void write()
Store calibrated constand.
Definition: XTCalibrationAlgorithm.cc:347