Belle II Software  release-06-01-15
XTCalibration.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 <cdc/calibration/XT.h>
9 #include <cdc/calibration/XTCalibration.h>
10 
11 #include <cdc/dbobjects/CDCXtRelations.h>
12 
13 #include <TError.h>
14 #include <TROOT.h>
15 #include <TH1D.h>
16 #include <TProfile.h>
17 #include <TF1.h>
18 #include <TFile.h>
19 #include <TChain.h>
20 #include <TSystem.h>
21 #include <iostream>
22 #include <iomanip>
23 
24 #include <framework/database/DBObjPtr.h>
25 #include <framework/logging/Logger.h>
26 #include <framework/utilities/FileSystem.h>
27 #include <cdc/calibration/CDCDatabaseImporter.h>
28 
29 #include <boost/iostreams/filtering_stream.hpp>
30 #include <boost/iostreams/device/file.hpp>
31 #include <boost/iostreams/filter/gzip.hpp>
32 
33 using namespace std;
34 using namespace Belle2;
35 using namespace CDC;
36 typedef std::array<float, 3> array3;
37 XTCalibration::XTCalibration()
38 // : m_firstExperiment(0), m_firstRun(0), m_lastExperiment(-1), m_lastRun(-1)
39 {
40  /*
41  setDescription(
42  " -------------------------- Test Calibration Algoritm -------------------------\n"
43  " \n"
44  " Testing algorithm which just gets mean of a test histogram collected by \n"
45  0 " CaTest module and provides a DB object with another histogram with one \n"
46  " entry at calibrated value. \n"
47  " ------------------------------------------------------------------------------\n"
48  );
49  */
50 }
51 
52 void XTCalibration::CreateHisto()
53 {
54  readXT();
55  readProfile();
56  /* read data from tree, make histo for fit*/
57  TChain* tree = new TChain("tree");
58  tree->Add(m_inputRootFileNames.c_str());
59  B2INFO("Open Files: " << m_inputRootFileNames);
60  if (!tree->GetBranch("ndf")) {
61  cout << "input data do not exits, please check!" << endl;
62  B2FATAL("echo rootfile do not exits or something wrong");
63  gSystem->Exec("echo rootfile do not exits or something wrong >> error");
64  return;
65  }
66  int lay;
67  double dt;
68  double dx;
69  double Pval, alpha, theta;
70  double ndf;
71 
72  tree->SetBranchAddress("lay", &lay);
73  tree->SetBranchAddress("t", &dt);
74  tree->SetBranchAddress("x_u", &dx);
75  tree->SetBranchAddress("alpha", &alpha);
76  tree->SetBranchAddress("theta", &theta);
77  tree->SetBranchAddress("Pval", &Pval);
78  tree->SetBranchAddress("ndf", &ndf);
79 
80  /* Disable unused branch */
81  std::vector<TString> list_vars = {"lay", "t", "x_u", "alpha", "theta", "Pval", "ndf"};
82  tree->SetBranchStatus("*", 0);
83 
84  for (TString brname : list_vars) {
85  tree->SetBranchStatus(brname, 1);
86  }
87 
88 
89  /*Create histogram*/
90  for (int i = 0; i < 56; ++i) {
91  for (int lr = 0; lr < 2; ++lr) {
92  for (int al = 0; al < m_nalpha; ++al) {
93  for (int th = 0; th < m_ntheta; ++th) {
94  hprof[i][lr][al][th] = new TProfile(Form("hprof%d_%d_%d_%d", i, lr, al, th),
95  Form("(L=%d)-(lr=%d)-(#alpha=%3.0f)-(#theta=%3.0f); Drift time (ns);Drift Length (cm)",
96  i, lr, ialpha[al], itheta[th]), 210, -20, 600, 0, 1.2, "i");
97  hist2d[i][lr][al][th] = new TH2D(Form("h%d_%d_%d_%d", i, lr, al, th),
98  Form("(L=%d)-(lr=%d)-(#alpha=%3.0f)-(#theta=%3.0f); Drift time (ns);Drift Length (cm)",
99  i, lr, ialpha[al], itheta[th]), 210, -20, 600, 110, 0, 1.2);
100  if (lr == 1)
101  hist2d_draw[i][al][th] = new TH2D(Form("h_draw%d_%d_%d", i, al, th),
102  Form("(L=%d)-(#alpha=%3.0f)-(#theta=%3.0f); Drift time (ns);Drift Length (cm)",
103  i, ialpha[al], itheta[th]), 210, -20, 600, 2200, -1.2, 1.2);
104  }
105  }
106  }
107  }
108 
109  /*Now read data and make histo*/
110  int al = 0;
111  int th = 0;
112  int lr = 0;
113  Long64_t nbytes = 0;
114  int nEntries = tree->GetEntries();
115  cout << "Number of Hit: " << nEntries << endl;
116 
117  for (int i = 0; i < nEntries; ++i) {
118  nbytes += tree->GetEntry(i);
119  /* protect in case |alpha|>90*/
120  if (fabs(alpha) > 90) {
121  if (alpha < 0) alpha += 180;
122  if (alpha > 0) alpha -= 180;
123  }
124 
125  if (Pval < m_Pvalmin) continue;
126  if (ndf < m_ndfmin) continue;
127 
128  for (int k = 0; k < m_nalpha; ++k) {
129  if (alpha < u_alpha[k]) {al = k; break;}
130  }
131  for (int j = 0; j < m_ntheta; ++j) {
132  if (theta < u_theta[j]) {th = j; break;}
133  }
134  if (dx > 0)
135  lr = 1;
136  else
137  lr = 0;
138  if (m_LRseparate) {
139  hprof[lay][lr][al][th]->Fill(dt, abs(dx));
140  hist2d[lay][lr][al][th]->Fill(dt, abs(dx));
141  } else {
142  hprof[lay][0][al][th]->Fill(dt, abs(dx));
143  hist2d[lay][0][al][th]->Fill(dt, abs(dx));
144  hprof[lay][1][al][th]->Fill(dt, abs(dx));
145  hist2d[lay][1][al][th]->Fill(dt, abs(dx));
146  }
147  hist2d_draw[lay][al][th]->Fill(dt, dx);
148  }
149 }
150 void XTCalibration::readProfile()
151 {
152  /*Read profile for xt*/
153  if (m_useProfileXTFromInputXT) {
154  B2INFO("use XT bining from input XT");
155  m_nalpha = nalpha_old;
156  m_ntheta = ntheta_old;
157  B2INFO("Nalpha: " << m_nalpha << "\n Ntheta: " << m_ntheta);
158  for (int i = 0; i < m_nalpha; ++i) {
159  l_alpha[i] = l_alpha_old[i];
160  u_alpha[i] = u_alpha_old[i];
161  ialpha[i] = ialpha_old[i];
162  B2INFO("-" << l_alpha[i] << " " << u_alpha[i] << " " << ialpha[i]);
163  }
164  for (int i = 0; i < m_ntheta; ++i) {
165  l_theta[i] = l_theta_old[i];
166  u_theta[i] = u_theta_old[i];
167  itheta[i] = itheta_old[i];
168  B2INFO("-" << l_theta[i] << " " << u_theta[i] << " " << itheta[i]);
169  }
170  } else {
171  B2INFO("use XT bining from profile file");
172  ifstream proxt(m_profileFileName.c_str());
173  if (!proxt) {
174  B2FATAL("file not found: " << m_profileFileName);
175  }
176  double dumy1, dumy2, dumy3;
177  proxt >> m_nalpha;
178  B2DEBUG(99, "Number of alpha bin" << m_nalpha);
179  if (m_nalpha > m_MAXalpha) {
180  B2FATAL("number of alpha bin excess limit; please increse uplimit: " << m_nalpha << " > " << m_MAXalpha);
181  }
182  for (int i = 0; i < m_nalpha; ++i) {
183  proxt >> dumy1 >> dumy2 >> dumy3;
184  l_alpha[i] = dumy1;
185  u_alpha[i] = dumy2;
186  ialpha[i] = dumy3;
187  }
188  proxt >> m_ntheta;
189  B2DEBUG(99, "Number of theta bin" << m_nalpha);
190  if (m_ntheta > m_MAXtheta) {B2FATAL("number of theta bin excess limit; please increse uplimit: " << m_ntheta << " > " << m_MAXtheta);}
191  for (int i = 0; i < m_ntheta; ++i) {
192  proxt >> dumy1 >> dumy2 >> dumy3;
193  l_theta[i] = dumy1;
194  u_theta[i] = dumy2;
195  itheta[i] = dumy3;
196  }
197  }
198  B2INFO("Finish asssign XT bining");
199 }
200 
201 bool XTCalibration::calibrate()
202 {
203  gROOT->SetBatch(1);
204  gErrorIgnoreLevel = 3001;
205 
206  CreateHisto();
207  B2INFO("Start Fitting");
208  for (int l = 0; l < 56; ++l) {
209  for (int lr = 0; lr < 2; ++lr) {
210  for (int al = 0; al < m_nalpha; ++al) {
211  for (int th = 0; th < m_ntheta; ++th) {
212 
213  if (hist2d[l][lr][al][th]->GetEntries() < m_smallestEntryRequire) {
214  fitflag[l][lr][al][th] = -1;
215  continue;
216  }
217  double p0, p1, tmin;
218  TF1* fpol1;
219  if (m_useSliceFit) {
220  hist2d[l][lr][al][th]->FitSlicesY(0, 0, -1, 5);
221  hist2d_1[l][lr][al][th] = (TH1D*)gDirectory->Get(Form("h%d_%d_%d_%d_1", l, lr, al, th));
222  if (!hist2d_1[l][lr][al][th]) {
223  fitflag[l][lr][al][th] = -1;
224  B2WARNING("error, not found results of slices fit");
225  continue;
226  }
227  hist2d_1[l][lr][al][th]->Fit("pol1", "Q", "", 30, 60);
228  fpol1 = (TF1*)hprof[l][lr][al][th]->GetFunction("pol1");
229  } else {
230  /*Set Error for low statistic bin*/
231  for (int n = 0; n < hprof[l][lr][al][th]->GetNbinsX(); ++n) {
232  if (hprof[l][lr][al][th]->GetBinEntries(n) < 5 && hprof[l][lr][al][th]->GetBinEntries(n) > 1) {
233  hprof[l][lr][al][th]->SetBinError(n, 0.3 / hprof[l][lr][al][th]->GetBinEntries(n));
234  }
235  }
236  hprof[l][lr][al][th]->Fit("pol1", "Q", "", 30, 60);
237  fpol1 = (TF1*)hprof[l][lr][al][th]->GetFunction("pol1");
238  }
239  if (fpol1) {
240  //determine tmin in fitting
241  p0 = fpol1->GetParameter(0);
242  p1 = fpol1->GetParameter(1);
243  tmin = -1 * p0 / p1 + 15;
244  } else {
245  p0 = 0; p1 = 0.005;
246  tmin = 12;
247  }
248  XT* xt;
249  if (m_useSliceFit) {
250  xt = new XT(hist2d_1[l][lr][al][th], m_xtmode);
251  } else {
252  xt = new XT(hprof[l][lr][al][th], m_xtmode);
253  }
254  // xt->setSmallestEntryRequired(m_smallestEntryRequire);
255  if (m_BField) {
256  int ial_old = 0;
257  int ith_old = 0;
258  for (int k = 0; k < nalpha_old; ++k) {
259  if (ialpha[al] < u_alpha_old[k]) {ial_old = k; break;}
260  }
261  for (int j = 0; j < ntheta_old; ++j) {
262  if (itheta[th] < u_theta_old[j]) {ith_old = j; break;}
263  }
264  double p6 = xtold[l][lr][ial_old][ith_old][6];
265  if (p6 > 400)
266  p6 = 400;
267  if (m_xtmode == xtmode_old) {
268  xt->setXTParams(xtold[l][lr][ial_old][ith_old]);
269  xt->setP6(p6);
270  } else {
271  xt->setXTParams(p0, p1, 0., 0., 0., 0., p6, xtold[l][lr][ial_old][ith_old][7]);
272  }
273  xt->setFitRange(tmin, p6 + 100);
274  } else {
275  xt->setXTParams(p0, p1, 0., 0., 0., 0., m_par6[l], 0.0001);
276  xt->setFitRange(tmin, m_par6[l] + 100);
277  }
278  xt->setDebug(m_debug);
279  xt->BField(m_BField);
280  xt->FitXT(m_xtmode);
281  /*get result*/
282  fitflag[l][lr][al][th] = xt->getFitStatus();
283  xtf5r[l][lr][al][th] = (TF1*)xt->getXTFunction();
284  if (m_useSliceFit) {
285  hist2d_1[l][lr][al][th] = (TH1D*)xt->getFittedHisto();
286  } else {
287  hprof[l][lr][al][th] = (TProfile*)xt->getFittedHisto();
288  }
289  }
290  }
291  }
292  }
293  Write();
294  storeHisto();
295  return true;
296 }
297 void XTCalibration::Write()
298 {
299  /*Set parameter for layer that fit do not success*/
300  /* and then write output file*/
301  double par[8];
302  ofstream xtout(m_OutputXTFileName.c_str());
303  xtout << m_nalpha << endl;
304  for (int i = 0; i < m_nalpha; ++i) {
305  xtout << std::setprecision(3) << l_alpha[i] << " "
306  << std::setprecision(3) << u_alpha[i] << " "
307  << std::setprecision(3) << ialpha[i] << endl;
308  }
309  xtout << m_ntheta << endl;
310  for (int i = 0; i < m_ntheta; ++i) {
311  xtout << std::setprecision(3) << l_theta[i] << " "
312  << std::setprecision(3) << u_theta[i] << " "
313  << std::setprecision(3) << itheta[i] << endl;
314  }
315  xtout << m_xtmode << " " << 8 << endl;
316 
317  int nfitted = 0;
318  int nfailure = 0;
319 
320  for (int th = 0; th < m_ntheta; ++th) {
321  for (int al = 0; al < m_nalpha; ++al) {
322  for (int l = 0; l < 56; ++l) {
323  for (int lr = 0; lr < 2; ++lr) {
324  /*Set Parameter for bad fit*/
325  if (fitflag[l][lr][al][th] != 1) {
326  nfailure += 1;
327  printf("fit failure status = %d \n", fitflag[l][lr][al][th]);
328  printf("layer %d, r %d, alpha %3.1f, theta %3.1f \n", l, lr, ialpha[al], itheta[th]);
329  printf("number of event: %3.2f", hprof[l][lr][al][th]->GetEntries());
330  if (fitflag[l][lr][al][th] != -1) {
331  printf("Probability of fit: %3.4f", xtf5r[l][lr][al][th]->GetProb());
332  }
333  if (m_xtmode == xtmode_old) {
334  int ial_old = 0;
335  int ith_old = 0;
336  for (int k = 0; k < nalpha_old; ++k) {
337  if (ialpha[al] < u_alpha_old[k]) {ial_old = k; break;}
338  }
339  for (int j = 0; j < ntheta_old; ++j) {
340  if (itheta[th] < u_theta_old[j]) {ith_old = j; break;}
341  }
342  for (int p = 0; p < 8; ++p) {
343  par[p] = xtold[l][lr][ial_old][ith_old][p];
344  }
345  } else {
346  //if mode of input xt is different from output, simple xt is used.
347  par[0] = 0; par[1] = 0.004; par[2] = 0; par[3] = 0; par[4] = 0; par[5] = 0; par[6] = m_par6[l]; par[7] = 0.00001;
348  }
349  } else {
350  xtf5r[l][lr][al][th]->GetParameters(par);
351  nfitted += 1;
352  }
353  /*Write params*/
354  xtout << l << std::setw(5) << itheta[th] << std::setw(5) << ialpha[al] << std::setw(5) << "0.0" << std::setw(4) << lr << std::setw(
355  15);
356  for (int p = 0; p < 8; ++p) {
357  if (p != 7) { xtout << std::setprecision(7) << par[p] << std::setw(15);}
358  if (p == 7) { xtout << std::setprecision(7) << par[p] << std::endl;}
359  }
360  }//lr
361  }//th
362  }//al
363  }//lay
364  xtout.close();
365  B2RESULT(" Total number of xt fit: " << m_nalpha * m_ntheta * 2 * 56);
366  B2RESULT(" Successfully Fitted: " << nfitted);
367  B2RESULT(" Failure Fit: " << nfailure);
368  B2RESULT("Finish export xt to text file");
369  if (m_useDB) {
370  CDCDatabaseImporter import(0, 0, -1, -1);
371  import.importXT(m_OutputXTFileName.c_str());
372  }
373 }
374 
375 void XTCalibration::storeHisto()
376 {
377  B2INFO("start store histogram");
378  TFile* fout = new TFile("XTFIT.root", "RECREATE");
379  TDirectory* top = gDirectory;
380  TDirectory* Direct[56];
381  int nhisto = 0;
382  for (int l = 0; l < 56; ++l) {
383  top->cd();
384  Direct[l] = gDirectory->mkdir(Form("lay_%d", l));
385  Direct[l]->cd();
386  for (int th = 0; th < m_ntheta; ++th) {
387  for (int al = 0; al < m_nalpha; ++al) {
388  hist2d_draw[l][al][th]->Write();
389  for (int lr = 0; lr < 2; ++lr) {
390  hist2d[l][lr][al][th]->Write();
391  if (fitflag[l][lr][al][th] != 1) continue;
392  if (m_useSliceFit) {
393  if (hist2d_1[l][lr][al][th]) {
394  hist2d_1[l][lr][al][th]->Write();
395  nhisto += 1;
396  }
397  } else {
398  hprof[l][lr][al][th]->Write();
399  nhisto += 1;
400  }
401  }
402  }
403  }
404  }
405  fout->Close();
406  B2RESULT(" " << nhisto << " histograms was stored.");
407 }
408 
409 void XTCalibration::readXT()
410 {
411  if (m_useDB) {
412  B2INFO("reading xt from DB");
413 
414  /*
415  ReadXT:readXTFromDB(&xtold,dbXT_old,
416  &nalpha_old,l_alpha_old,u_alpha_old,ialpha_old,
417  &ntheta_old,l_theta_old,u_theta_old, itheta_old);
418  */
419  readXTFromDB();
420  B2INFO("Number of theta bin from xt: " << ntheta_old);
421  B2INFO("Theta 0: " << itheta_old[0]);
422  // if(!a){
424  } else {
425  B2INFO("Read Xt from text");
426  /*
427  ReadXT::readXTFromText(xtold,m_xtfile,
428  nalpha_old,l_alpha_old,u_alpha_old,ialpha_old,
429  ntheta_old,l_theta_old,u_theta_old, itheta_old);
430  */
431  readXTFromText();
432  B2INFO("nalpha: " << nalpha_old);
433  // if(!a)
434  // {B2FATAL("Error reading xt from text");return;}
435  }
436 }
437 
438 void XTCalibration::readXTFromText()
439 {
440  std::string fileName1 = "/data/cdc" + m_xtfile;
441  std::string fileName = FileSystem::findFile(fileName1);
442  boost::iostreams::filtering_istream ifs;
443  if (fileName == "") {
444  fileName = FileSystem::findFile(m_xtfile);
445  }
446  if (fileName == "") {
447  B2FATAL("CDCGeometryPar: " << fileName1 << " not exist!");
448  } else {
449  B2INFO("CDCGeometryPar: open " << fileName1);
450  if ((fileName.rfind(".gz") != string::npos) && (fileName.length() - fileName.rfind(".gz") == 3)) {
451  ifs.push(boost::iostreams::gzip_decompressor());
452  }
453  ifs.push(boost::iostreams::file_source(fileName));
454  if (!ifs) B2FATAL("CDCGeometryPar: cannot open " << fileName1 << " !");
455 
456 
457  }
458  int npar = 8;
459  //read alpha bin info.
460  // unsigned short nAlphaBins = 0;
461  ifs >> nalpha_old;
462 
463  for (unsigned short i = 0; i < nalpha_old; ++i) {
464  ifs >> l_alpha_old[i] >> u_alpha_old[i] >> ialpha_old[i];
465  // ifs >> alpha0 >> alpha1 >> alpha2;
466  //m_alphaPoints[i] = alpha2;
467  }
468 
469  //read theta bin info.
470  // unsigned short nThetaBins = 0;
471  ifs >> ntheta_old;
472 
473  for (unsigned short i = 0; i < ntheta_old; ++i) {
474  ifs >> l_theta_old[i] >> u_theta_old[i] >> itheta_old[i];
475  //ifs >> theta0 >> theta1 >> theta2;
476  //m_thetaPoints[i] = theta2;
477  }
478 
479  B2INFO("number of alpha - theta bin" << nalpha_old << " - " << ntheta_old);
480  short np = 0;
481  unsigned short iCL, iLR;
482  // const unsigned short npx = nXTParams - 1;
483  double xtc[npar];
484  double theta, alpha, dummy1;
485  unsigned nRead = 0;
486  // unsigned m_xtParamMode_old;
487  ifs >> xtmode_old >> np;
488 
489  const double epsi = 0.1;
490 
491  while (ifs >> iCL) {
492  ifs >> theta >> alpha >> dummy1 >> iLR;
493  for (int i = 0; i < np; ++i) {
494  ifs >> xtc[i];
495  }
496  ++nRead;
497 
498  int ith = -99;
499  for (unsigned short i = 0; i < ntheta_old; ++i) {
500  if (fabs(theta - itheta_old[i]) < epsi) {
501  ith = i;
502  break;
503  }
504  }
505  if (ith < 0) {
506  gSystem->Exec("echo xt_theta error binning>> error");
507  return;
508  }
509 
510  int ial = -99;
511  for (unsigned short i = 0; i < nalpha_old; ++i) {
512  if (fabs(alpha - ialpha_old[i]) < epsi) {
513  ial = i;
514  break;
515  }
516  }
517  if (ial < 0) {
518  gSystem->Exec("echo xt_alpha error binning>> error");
519  return;
520  }
521 
522  for (int i = 0; i < np; ++i) {
523  xtold[iCL][iLR][ial][ith][i] = xtc[i];
524  }
525 
526  } //end of while loop
527 
528  //convert unit
529  /*
530  const double degrad = M_PI / 180.;
531  for (unsigned i = 0; i < nAlphaBins; ++i) {
532  m_alphaPoints[i] *= degrad;
533  }
534  for (unsigned i = 0; i < nThetaBins; ++i) {
535  m_thetaPoints[i] *= degrad;
536  }
537  */
538  // return true;
539 }
540 void XTCalibration::readXTFromDB()
541 {
542  DBObjPtr<CDCXtRelations> dbXT_old;
543  nalpha_old = dbXT_old->getNoOfAlphaBins();
544  B2INFO("Number of alpha" << nalpha_old);
545  double rad2deg = 180 / M_PI;
546  for (unsigned short i = 0; i < nalpha_old; ++i) {
547  array3 alpha = dbXT_old->getAlphaBin(i);
548  l_alpha_old[i] = alpha[0] * rad2deg;
549  u_alpha_old[i] = alpha[1] * rad2deg;
550  ialpha_old[i] = alpha[2] * rad2deg;
551  // std::cout << m_alphaPoints[i]*180./M_PI << std::endl;
552  }
553 
554  ntheta_old = dbXT_old->getNoOfThetaBins();
555  B2INFO("Ntheta: " << ntheta_old);
556  for (unsigned short i = 0; i < ntheta_old; ++i) {
557  // m_thetaPoints[i] = (*dbXT_old).getThetaPoint(i);
558  array3 theta = dbXT_old->getThetaBin(i);
559  l_theta_old[i] = theta[0] * rad2deg;
560  u_theta_old[i] = theta[1] * rad2deg;
561  itheta_old[i] = theta[2] * rad2deg;
562 
563 
564  // std::cout << m_thetaPoints[i]*180./M_PI << std::endl;
565  }
566 
567  xtmode_old = dbXT_old->getXtParamMode();
568 
569  for (unsigned short iCL = 0; iCL < 56; ++iCL) {
570  for (unsigned short iLR = 0; iLR < 2; ++iLR) {
571  for (unsigned short iA = 0; iA < nalpha_old; ++iA) {
572  for (unsigned short iT = 0; iT < ntheta_old; ++iT) {
573  const std::vector<float> params = dbXT_old->getXtParams(iCL, iLR, iA, iT);
574  unsigned short np = params.size();
575  // std::cout <<"np4xt= " << np << std::endl;
576  for (unsigned short i = 0; i < np; ++i) {
577  xtold[iCL][iLR][iA][iT][i] = params[i];
578  }
579  }
580  }
581  }
582  }
583  // return true;
584 }
CDC database importer.
void importXT(std::string fileName)
Import xt table to the database.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Class to perform fitting for each xt function.
Definition: XT.h:60
void setXTParams(double p[8])
Set Paramerters for fit.
Definition: XT.h:115
void setP6(double p6)
Set Parameter 6 for polynomia fit.
Definition: XT.h:93
void setFitRange(double tmin, double tmax)
Set Fit range.
Definition: XT.h:135
TF1 * getXTFunction(int mode)
Get XT function.
Definition: XT.h:179
void FitXT()
Do fitting.
Definition: XT.h:200
void setDebug(bool debug)
Set Debug.
Definition: XT.h:150
TProfile * getFittedHisto()
Get histogram.
Definition: XT.h:195
int getFitStatus()
get fitted flag.
Definition: XT.h:157
void BField(bool bfield)
set to use BField
Definition: XT.h:110
Abstract base class for different kinds of events.