Belle II Software  release-05-02-19
XT.h
1 #include "TProfile.h"
2 #include "TCanvas.h"
3 #include "TF1.h"
4 #include "TH1D.h"
5 #include "Math/ChebyshevPol.h"
6 #include "iostream"
7 
12 Double_t pol5pol1(Double_t* x, Double_t* par)
13 {
14  Double_t xx = x[0];
15  Double_t x6, x2;
16  Double_t f, ctp;
17  if (xx < par[6]) {
18  f = par[0] + par[1] * xx + par[2] * xx * xx + par[3] * xx * xx * xx + par[4] * xx * xx * xx * xx + par[5] * xx * xx * xx * xx * xx;
19  } else {
20  x6 = par[6];
21  x2 = xx - x6;
22  ctp = par[0] + par[1] * x6 + par[2] * x6 * x6 + par[3] * x6 * x6 * x6 + par[4] * x6 * x6 * x6 * x6 + par[5] * x6 * x6 * x6 * x6 *
23  x6;
24  f = par[7] * x2 + ctp;
25  }
26  return f;
27 }
28 
33 Double_t Cheb5pol1(Double_t* x, Double_t* par)
34 {
35  Double_t xx = x[0];
36  Double_t x6, x2;
37  Double_t f, ctp;
38 
39  if (xx < par[6]) {
40  f = ROOT::Math::Chebyshev5(xx, par[0], par[1], par[2], par[3], par[4], par[5]);
41  } else {
42  x6 = par[6];
43  x2 = xx - x6;
44  ctp = ROOT::Math::Chebyshev5(x6, par[0], par[1], par[2], par[3], par[4], par[5]);
45  f = par[7] * x2 + ctp;
46  }
47  return f;
48 }
53 class XT {
54 public:
58  explicit XT(TProfile* h1)
59  {
60  m_h1 = (TProfile*)h1->Clone();
61  m_h1->SetDirectory(0);
62  }
66  XT(TProfile* h1, int mode)
67  {
68  m_h1 = (TProfile*)h1->Clone();
69  m_h1->SetDirectory(0);
70  m_mode = mode;
71  }
72 
76  XT(TH1D* h1, int mode)
77  {
78  m_h1 = (TProfile*)h1->Clone();
79  m_h1->SetDirectory(0);
80  m_mode = mode;
81  }
82 
86  void setP6(double p6)
87  {
88  m_XTParam[6] = p6;
89  }
90 
96  void setMode(int mode)
97  {
98  m_mode = mode;
99  }
103  void BField(bool bfield) {m_BField = bfield;}
104 
108  void setXTParams(double p[8])
109  {
110  for (int i = 0; i < 8; ++i) {m_XTParam[i] = p[i];}
111  m_tmax = p[6] + 250;
112  }
116  void setXTParams(double p0, double p1, double p2, double p3,
117  double p4, double p5, double p6, double p7)
118  {
119  m_XTParam[0] = p0; m_XTParam[1] = p1; m_XTParam[2] = p2;
120  m_XTParam[3] = p3; m_XTParam[4] = p4; m_XTParam[5] = p5;
121  m_XTParam[6] = p6; m_XTParam[7] = p7;
122  m_tmax = p6 + 250;
123  }
124 
128  void setFitRange(double tmin, double tmax)
129  {
130  m_tmin = tmin;
131  m_tmax = tmax;
132  }
137  {
138  m_minRequiredEntry = min;
139  }
143  void setDebug(bool debug)
144  {
145  m_debug = debug;
146  }
151  {
152  return m_fitflag;
153  }
154 
158  double getProb()
159  {
160  return m_Prob;
161  }
165  void getFittedXTParams(double pa[8])
166  {
167  for (int i = 0; i < 8; ++i) {pa[i] = m_FittedXTParams[i];}
168  }
172  TF1* getXTFunction(int mode)
173  {
174  if (mode == 0) return xtpol5;
175  else return xtCheb5;
176  }
181  {
182  if (m_mode == 0) return xtpol5;
183  else return xtCheb5;
184  };
188  TProfile* getFittedHisto() {return m_h1;}
189 
193  void FitXT()
194  {
195  if (m_mode == 0) FitPol5();
196  else FitChebyshev();
197  }
198 
202  void FitXT(int mode)
203  {
204  if (mode == 0) FitPol5();
205  else FitChebyshev();
206  }
210  void FitPol5();
214  void FitChebyshev();
215 private:
216 
217  TProfile* m_h1;
218  TF1* xtpol5 = new TF1("xtpol5", pol5pol1, 0.0, 700, 8);
219  TF1* xtCheb5 = new TF1("xtCheb5", Cheb5pol1, 0.0, 700, 8);
221  int m_mode = 1;
222  bool m_debug = true;
223  bool m_draw = false;
224  bool m_BField = true;
225  int m_minRequiredEntry = 800;
226  double m_XTParam[8] = {};
227  double m_FittedXTParams[8] = {};
237  int m_fitflag = 0;
238  double m_Prob = 0;
239  double m_tmin = 20;
240  double m_tmax = m_XTParam[6] + 250;
241 };
242 
244 {
245  double max_dif = 0.12;
246  double max_dif2 = 0.05;
247  double par[8];
248  m_h1->Fit("pol1", "MQ", "", m_tmin, 50);
249  TF1* f1 = (TF1*)m_h1->GetFunction("pol1");
250  double p0 = f1->GetParameter(0);
251  double p1 = f1->GetParameter(1);
252  double f10 = f1->Eval(10);
253  /****************************/
254  int out = 0; /*how many time outer part change fit limit*/
255  xtpol5->SetParameters(p0, p1, 0, 0, 0, 0, m_XTParam[6], 0);
256  double p6default = m_XTParam[6];
257  if (m_BField) {
258  xtpol5->SetParLimits(7, 0.0, 0.001);
259  xtpol5->SetParLimits(1, 0.0, 0.01);
260  } else {
261  xtpol5->SetParLimits(7, 0.0, 0.01);
262  xtpol5->SetParLimits(1, 0.0, 0.01);
263  }
264  // xtpol5->SetParLimits(6, p6default - 30, p6default + 30);
265  for (int i = 0; i < 10; ++i) {
266  m_fitflag = 1;
267  std::cout << "Fitting" << std::endl;
268  double stat = m_h1->Fit(xtpol5, "M", "0", m_tmin, m_tmax);
269  xtpol5->GetParameters(par);
270 
271  /*Eval outer region,*/
272  double fp6 = xtpol5->Eval(par[6]);
273  double fbehindp6 = xtpol5->Eval(par[6] - 12) - 0.01;
274  if (fp6 < fbehindp6 || fp6 > 1) { /*may be change to good value*/
275  m_fitflag = 2;
276  out += 1;
277  xtpol5->SetParameters(p0, p1, 0, 0, 0, 0, p6default, 0);
278  xtpol5->SetParLimits(6, p6default - 10, p6default + 10 - out / 2);
279  m_tmax += 2;
280  if (m_tmax < p6default + 30) {
281  m_tmax = p6default + 30;
282  }
283  }
284  /* EVal inner region*/
285  /* compare p0 of linear fit vs p0 of xtfun fit and eval point t=10)*/
286  if (fabs(par[0] - p0) > max_dif || fabs(f10 - xtpol5->Eval(10)) > max_dif2) {
287  m_fitflag = 3;
288  if (i == 9) std::cout << "ERROR XT FIT inner part" << std::endl;
289  xtpol5->SetParameters(p0, p1, 0, 0, 0, 0, p6default, 0);
290  xtpol5->SetParLimits(1, 0, 0.08);
291  m_tmin -= 0.5;
292 
293  if (m_tmin < 14) {
294  m_tmin = 14; std::cout << "ERROR: tmin so small, fit iter: " << std::endl;
295  }
296  }
297 
298  if (stat == 0) {
299  m_fitflag = 0;
300  if (m_debug) {std::cout << "FIT Failure; Layer: " << std::endl;}
301  }
302 
303  if (m_debug) {printf("P6default= %3.2f, p6 fit = %3.2f", p6default, par[6]);}
304 
305  if (m_fitflag == 1) {
306  if (m_debug) std::cout << "Fit success" << std::endl;
307  xtpol5->GetParameters(m_FittedXTParams);
308  m_Prob = xtpol5->GetProb();
309  break;
310  }
311 
312  } //end loop of fitting
313  if (m_draw) {
314  TString hname = m_h1->GetName();
315  TString name = hname + ".pdf";
316  TCanvas* c1 = new TCanvas("c1", "", 800, 600);
317  m_h1->Draw();
318  xtpol5->SetLineColor(kBlack);
319  xtpol5->DrawF1(0, 400, "same");
320  c1->SaveAs(name);
321 
322  }
323 }
324 
326 {
327  if (m_h1->GetEntries() < m_minRequiredEntry) {
328  m_fitflag = -1;
329  return;
330  }
331  // m_tmax = m_XTParam[6] + 100;
332  //xtCheb5->SetParameters(0.0, 0.005, 0., 0., 0., 0., m_XTParam[6], 0.001);
333  double p[6];
334  double par[8];
335  xtCheb5->SetParLimits(7, 0., 0.001);
336  m_h1->Fit("chebyshev5", "QME", "", m_tmin, m_XTParam[6]);
337  m_h1->GetFunction("chebyshev5")->GetParameters(p);
338  xtCheb5->SetParameters(p[0], p[1], p[2], p[3], p[4], p[5], m_XTParam[6], 0.000);
339 
340  double stat;
341  for (int i = 0; i < 10; ++i) {
342  stat = m_h1->Fit("xtCheb5", "MQ", "0", m_tmin, m_tmax);
343  if (stat == 0) {
344  xtCheb5->SetParameters(p[0], p[1], p[2], p[3], p[4], p[5], m_XTParam[6] - 20, 0.000);
345  m_tmax -= 10;
346  continue;
347  }
348  xtCheb5->GetParameters(par);
349  /*Eval outer region,*/
350  double fp6 = xtCheb5->Eval(par[6]);
351  double fbehindp6 = xtCheb5->Eval(par[6] - 10) - 0.005;
352  if (fp6 < fbehindp6 || fp6 > 1) { /*may be change to good value*/
353  m_fitflag = 2;
354  // out += 1;
355  xtCheb5->SetParameters(p[0], p[1], p[2], p[3], p[4], p[5], par[6] - 20, 0.000);
356  xtCheb5->SetParLimits(6, par[6] - 50, par[6] - 10);
357  m_tmax -= 10;
358  // if (m_tmax < p6default + 30) {
359  // m_tmax = p6default + 30;
360  // }
361  } else {
362  break;
363  }
364  // m_tmax +=10;
365  //if (stat != 0) break;
366  }
367  xtCheb5->GetParameters(m_FittedXTParams);
368  m_Prob = xtCheb5->GetProb();
369  if (stat == 0)
370  m_fitflag = 0;
371  else
372  m_fitflag = 1;
373  // xtCheb5->DrawF1(0,400,"same");
374  if (m_draw) {
375  TString hname = m_h1->GetName();
376  TString name = hname + ".pdf";
377  TCanvas* c1 = new TCanvas("c1", "", 800, 600);
378  m_h1->Draw();
379  xtCheb5->SetLineColor(kBlack);
380  xtCheb5->DrawF1(0, 400, "same");
381  c1->SaveAs(name);
382  }
383 }
XT::getFittedHisto
TProfile * getFittedHisto()
Get histogram.
Definition: XT.h:188
XT::getProb
double getProb()
Get the chi2 probability.
Definition: XT.h:158
XT::xtCheb5
TF1 * xtCheb5
5th order Cheb.
Definition: XT.h:219
XT::getXTFunction
TF1 * getXTFunction()
Get XT function.
Definition: XT.h:180
XT::setP6
void setP6(double p6)
Set Parameter 6 for polynomia fit.
Definition: XT.h:86
XT::m_Prob
double m_Prob
Chi2 prob of fitting.
Definition: XT.h:238
XT::m_BField
bool m_BField
With magnetic field or not.
Definition: XT.h:224
XT::setXTParams
void setXTParams(double p[8])
Set Paramerters for fit.
Definition: XT.h:108
XT::m_fitflag
int m_fitflag
Fit Flag =-1: low statitic =1: good =0: Fit failure =2: Error Outer =3: Error Inner part;.
Definition: XT.h:237
XT::m_debug
bool m_debug
Print debug durring fitting or not.
Definition: XT.h:222
XT::xtpol5
TF1 * xtpol5
5th order polynomial function
Definition: XT.h:218
XT::FitXT
void FitXT()
Do fitting.
Definition: XT.h:193
XT::setXTParams
void setXTParams(double p0, double p1, double p2, double p3, double p4, double p5, double p6, double p7)
Set Initial parameters for fitting.
Definition: XT.h:116
XT::m_tmax
double m_tmax
upper boundary of fit range
Definition: XT.h:240
XT::BField
void BField(bool bfield)
set to use BField
Definition: XT.h:103
XT
Class to perform fitting for each xt function.
Definition: XT.h:53
XT::setMode
void setMode(int mode)
Set XT mode.
Definition: XT.h:96
XT::XT
XT(TProfile *h1, int mode)
Initialized with TProfile histogram and mode.
Definition: XT.h:66
XT::m_XTParam
double m_XTParam[8]
Parameter fo xt.
Definition: XT.h:226
XT::XT
XT(TH1D *h1, int mode)
Initialized with TH1D histogram and mode.
Definition: XT.h:76
XT::m_tmin
double m_tmin
lower boundary of fit range
Definition: XT.h:239
XT::m_FittedXTParams
double m_FittedXTParams[8]
Fitted parameters.
Definition: XT.h:227
XT::getFitStatus
int getFitStatus()
get fitted flag.
Definition: XT.h:150
XT::m_minRequiredEntry
int m_minRequiredEntry
Minimum entry required for each histo.
Definition: XT.h:225
XT::setDebug
void setDebug(bool debug)
Set Debug.
Definition: XT.h:143
XT::FitXT
void FitXT(int mode)
Do fitting.
Definition: XT.h:202
XT::m_mode
int m_mode
XT mode, 0 is for 5th order polynomial, 1 is Chebshev polynomial.
Definition: XT.h:221
XT::m_h1
TProfile * m_h1
Histogram of xt relation.
Definition: XT.h:217
XT::FitChebyshev
void FitChebyshev()
Fit xt histogram incase 5th order Chebeshev polynomial is used.
Definition: XT.h:325
XT::FitPol5
void FitPol5()
Fit xt histogram incase 5th order polynomial is used.
Definition: XT.h:243
XT::getXTFunction
TF1 * getXTFunction(int mode)
Get XT function.
Definition: XT.h:172
XT::setSmallestEntryRequired
void setSmallestEntryRequired(int min)
Set minimum number of entry required for fit.
Definition: XT.h:136
XT::getFittedXTParams
void getFittedXTParams(double pa[8])
get fit parameters.
Definition: XT.h:165
XT::XT
XT(TProfile *h1)
Initialized with TProfile histogram.
Definition: XT.h:58
XT::setFitRange
void setFitRange(double tmin, double tmax)
Set Fit range.
Definition: XT.h:128
XT::m_draw
bool m_draw
Draw and store png plot of each histo or not.
Definition: XT.h:223