Belle II Software development
XTFunction Class Reference

Class to perform fitting for each xt function. More...

#include <XTFunction.h>

Public Member Functions

 XTFunction (TProfile *h1)
 Initialized with TProfile histogram.
 
 XTFunction (TProfile *h1, int mode)
 Initialized with TProfile histogram and mode.
 
 XTFunction (TH1F *h1, int mode)
 Initialized with TH1D histogram and mode.
 
 XTFunction (const XTFunction &x)
 Copy constructor.
 
XTFunctionoperator= (const XTFunction &x)
 Assignment operator.
 
void setP6 (double p6)
 Set Parameter 6 for polynomia fit.
 
bool isValid ()
 Is valid.
 
void setMode (int mode)
 Set XT mode.
 
void setBField (bool bfield)
 set to use BField
 
void setXTParams (const double p[8])
 Set Parameters for fit.
 
void setXTParams (double p0, double p1, double p2, double p3, double p4, double p5, double p6, double p7)
 Set Initial parameters for fitting.
 
void setFitRange (double tmin, double tmax)
 Set Fit range.
 
void setSmallestEntryRequired (int min)
 Set minimum number of entry required for fit.
 
void setDebug (bool debug)
 Set Debug.
 
int getFitStatus ()
 get fitted flag.
 
double getProb ()
 Get the chi2 probability.
 
void getFittedXTParams (double pa[8])
 get fit parameters.
 
TF1 * getXTFunction ()
 Get XT function.
 
TProfile * getFittedHisto ()
 Get histogram.
 
void fitXT ()
 Do fitting.
 
void FitPol5 ()
 Fit xt histogram incase 5th order polynomial is used.
 
void FitChebyshev ()
 Fit xt histogram incase 5th order Chebeshev polynomial is used.
 
bool validate ()
 Validate the xt has proper shape.
 

Private Attributes

TProfile * m_h1
 Histogram of xt relation.
 
TF1 * m_fitFunc
 Fit function.
 
int m_mode = c_Chebyshev
 XT mode, 0 is for 5th order polynomial, 1 is Chebshev polynomial.
 
bool m_debug = true
 Print debug durring fitting or not.
 
bool m_draw = false
 Draw and store png plot of each histo or not.
 
bool m_bField = true
 With magnetic field or not.
 
int m_minRequiredEntry = 800
 Minimum entry required for each histo.
 
double m_XTParam [8] = {}
 Parameter fo xt.
 
double m_FittedXTParams [8] = {}
 Fitted parameters.
 
int m_fitflag = 0
 Fit Flag =-1: low statitic =1: good =0: Fit failure =2: Error Outer =3: Error Inner part;.
 
double m_Prob = 0
 Chi2 prob of fitting.
 
double m_tmin = 20
 lower boundary of fit range
 
double m_tmax = m_XTParam[6] + 50
 upper boundary of fit range
 

Detailed Description

Class to perform fitting for each xt function.

Definition at line 69 of file XTFunction.h.

Constructor & Destructor Documentation

◆ XTFunction() [1/4]

XTFunction ( TProfile *  h1)
inlineexplicit

Initialized with TProfile histogram.

Definition at line 74 of file XTFunction.h.

75 {
76 m_h1 = (TProfile*)h1->Clone();
77 m_h1->SetDirectory(0);
78 if (m_mode == c_Chebyshev) {
79 m_fitFunc = new TF1("xtCheb5", cheby5pol1, 0.0, 700, 8);
80 } else {
81 m_fitFunc = new TF1("xtpol5", pol5pol1, 0.0, 700, 8);
82 }
83 }
TF1 * m_fitFunc
Fit function.
Definition: XTFunction.h:304
int m_mode
XT mode, 0 is for 5th order polynomial, 1 is Chebshev polynomial.
Definition: XTFunction.h:309
TProfile * m_h1
Histogram of xt relation.
Definition: XTFunction.h:303
Double_t cheby5pol1(Double_t *x, Double_t *par)
helper function to initialize xt function with 5th order Chebshev Polynomial + linear.
Definition: XTFunction.h:48
Double_t pol5pol1(Double_t *x, Double_t *par)
helper function to initialize xt function with 5th order polynomial + linear.
Definition: XTFunction.h:27

◆ XTFunction() [2/4]

XTFunction ( TProfile *  h1,
int  mode 
)
inline

Initialized with TProfile histogram and mode.

Definition at line 87 of file XTFunction.h.

88 {
89 m_h1 = (TProfile*)h1->Clone();
90 m_h1->SetDirectory(0);
91 m_mode = mode;
92 if (m_mode == c_Chebyshev) {
93 m_fitFunc = new TF1("xtCheb5", cheby5pol1, 0.0, 700, 8);
94 } else {
95 m_fitFunc = new TF1("xtpol5", pol5pol1, 0.0, 700, 8);
96 }
97
98 }

◆ XTFunction() [3/4]

XTFunction ( TH1F *  h1,
int  mode 
)
inline

Initialized with TH1D histogram and mode.

Definition at line 103 of file XTFunction.h.

104 {
105 m_h1 = (TProfile*)h1->Clone();
106 m_h1->SetDirectory(0);
107 m_mode = mode;
108 if (m_mode == c_Chebyshev) {
109 m_fitFunc = new TF1("xtCheb5", cheby5pol1, 0.0, 700, 8);
110 } else {
111 m_fitFunc = new TF1("xtpol5", pol5pol1, 0.0, 700, 8);
112 }
113 }

◆ XTFunction() [4/4]

XTFunction ( const XTFunction x)
inline

Copy constructor.

Definition at line 118 of file XTFunction.h.

118 :
119 m_h1(x.m_h1),
120 m_mode(x.m_mode),
121 m_debug(x.m_debug),
122 m_draw(x.m_draw),
123 m_bField(x.m_bField),
124 m_minRequiredEntry(x.m_minRequiredEntry),
125 m_fitflag(x.m_fitflag),
126 m_Prob(x.m_Prob),
127 m_tmin(x.m_tmin),
128 m_tmax(x.m_tmax)
129 {
130 m_fitFunc = (TF1*) x.m_fitFunc->Clone();
131 for (int i = 0; i < 8; ++i) {
132 m_XTParam[i] = x.m_XTParam[i];
133 m_FittedXTParams[i] = x.m_XTParam[i];
134 }
135 }
int m_minRequiredEntry
Minimum entry required for each histo.
Definition: XTFunction.h:313
double m_FittedXTParams[8]
Fitted parameters.
Definition: XTFunction.h:315
int m_fitflag
Fit Flag =-1: low statitic =1: good =0: Fit failure =2: Error Outer =3: Error Inner part;.
Definition: XTFunction.h:325
bool m_debug
Print debug durring fitting or not.
Definition: XTFunction.h:310
double m_Prob
Chi2 prob of fitting.
Definition: XTFunction.h:326
bool m_bField
With magnetic field or not.
Definition: XTFunction.h:312
double m_tmin
lower boundary of fit range
Definition: XTFunction.h:327
double m_tmax
upper boundary of fit range
Definition: XTFunction.h:328
bool m_draw
Draw and store png plot of each histo or not.
Definition: XTFunction.h:311
double m_XTParam[8]
Parameter fo xt.
Definition: XTFunction.h:314

Member Function Documentation

◆ FitChebyshev()

void FitChebyshev ( )

Fit xt histogram incase 5th order Chebeshev polynomial is used.

Definition at line 436 of file XTFunction.h.

437 {
438 if (m_mode != c_Chebyshev) {
439 B2ERROR("Fitting function is wrong");
440 }
441
442 if (m_h1->GetEntries() < m_minRequiredEntry) {
443 m_fitflag = -1;
444 return;
445 }
446 // m_tmax = m_XTParam[6] + 100;
447 //xtCheb5->SetParameters(0.0, 0.005, 0., 0., 0., 0., m_XTParam[6], 0.001);
448 double par[8] = {0.0};
449 m_fitFunc->SetParLimits(7, 0., 0.001);
450 int fitresult = m_h1->Fit("chebyshev5", "QME", "", m_tmin, m_XTParam[6]);
451 if (fitresult >= 0) {
452 m_h1->GetFunction("chebyshev5")->GetParameters(par);
453 m_fitFunc->SetParameters(par[0], par[1], par[2], par[3], par[4], par[5], m_XTParam[6], 0.000);
454 }
455 double stat;
456 for (int i = 0; i < 10; ++i) {
457 stat = m_h1->Fit(m_fitFunc, "MQ", "0", m_tmin, m_tmax);
458 if (stat == 0) {
459 m_fitFunc->SetParameters(par[0], par[1], par[2], par[3], par[4], par[5], m_XTParam[6] - 20, 0.000);
460 m_tmax -= 10;
461 continue;
462 }
463 m_fitFunc->GetParameters(par);
464 if (par[1] < 0) { // negative c1
465 // std::cout << " neg c1 converted" << std::endl;
466 par[1] *= -1.0;
467 m_fitFunc->SetParLimits(1, 0., 0.01);
468 m_tmin += 10.0;
469 continue;
470 }
471
472 /*Eval outer region,*/
473 double fp6 = m_fitFunc->Eval(par[6]);
474 double fbehindp6 = m_fitFunc->Eval(par[6] - 10) - 0.005;
475 if (fp6 < fbehindp6 || fp6 > 1) { /*may be change to good value*/
476 m_fitflag = 2;
477 // out += 1;
478 m_fitFunc->SetParameters(par[0], par[1], par[2], par[3], par[4], par[5], par[6] - 20, 0.000);
479 m_fitFunc->SetParLimits(6, par[6] - 50, par[6] - 10);
480 m_tmax -= 10;
481 // if (m_tmax < p6default + 30) {
482 // m_tmax = p6default + 30;
483 // }
484 } else {
485 break;
486 }
487 // m_tmax +=10;
488 //if (stat != 0) break;
489 }
490 if (par[1] > 0) {
491 m_fitFunc->GetParameters(m_FittedXTParams);
492 m_Prob = m_fitFunc->GetProb();
493 }
494 if (stat == 0)
495 m_fitflag = 0;
496 else
497 m_fitflag = 1;
498 // m_fitFunc->DrawF1(0,400,"same");
499 if (m_draw) {
500 TString hname = m_h1->GetName();
501 TString name = hname + ".pdf";
502 TCanvas* c1 = new TCanvas("c1", "", 800, 600);
503 m_h1->Draw();
504 m_fitFunc->SetLineColor(kBlack);
505 m_fitFunc->DrawF1(0, 400, "same");
506 c1->SaveAs(name);
507 }
508 }

◆ FitPol5()

void FitPol5 ( )

Fit xt histogram incase 5th order polynomial is used.

Definition at line 331 of file XTFunction.h.

332 {
333 if (m_mode != c_Polynomial) {
334 B2ERROR("Fitting function is wrong");
335 }
336 double max_dif = 0.12;
337 double max_dif2 = 0.05;
338 double par[8];
339 m_h1->Fit("pol1", "MQ", "", m_tmin, 50);
340 TF1* f1 = (TF1*)m_h1->GetFunction("pol1");
341 double p0 = f1->GetParameter(0);
342 double p1 = f1->GetParameter(1);
343 double f10 = f1->Eval(10);
344 /****************************/
345 int in = 0; /*how many time inner part change fit limit*/
346 int out = 0; /*how many time outer part change fit limit*/
347 m_fitFunc->SetParameters(p0, p1, 0, 0, 0, 0, m_XTParam[6], 0);
348 double p6default = m_XTParam[6];
349 if (m_bField) {
350 m_fitFunc->SetParLimits(7, 0.0, 0.001);
351 m_fitFunc->SetParLimits(1, 0.0, 0.01);
352 } else {
353 m_fitFunc->SetParLimits(7, 0.0, 0.01);
354 m_fitFunc->SetParLimits(1, 0.0, 0.01);
355 }
356
357 for (int i = 0; i < 10; ++i) {
358 m_fitflag = 1;
359 std::cout << "Fitting" << std::endl;
360 double stat = m_h1->Fit(m_fitFunc, "M", "0", m_tmin, m_tmax);
361 m_fitFunc->GetParameters(par);
362
363 /*Eval outer region,*/
364 double fp6 = m_fitFunc->Eval(par[6]);
365 double fbehindp6 = m_fitFunc->Eval(par[6] - 12) - 0.01;
366 if (fp6 < fbehindp6 || fp6 > 1) { /*may be change to good value*/
367 m_fitflag = 2;
368 out += 1;
369 m_fitFunc->SetParameters(p0, p1, 0, 0, 0, 0, p6default, 0);
370 m_fitFunc->SetParLimits(6, p6default - 10, p6default + 10 - out / 2);
371 m_tmax += 2;
372 if (m_tmax < p6default + 30) {
373 m_tmax = p6default + 30;
374 }
375 }
376 /* EVal inner region*/
377 /* compare p0 of linear fit vs p0 of xtfun fit and eval point t=10)*/
378 if (fabs(par[0] - p0) > max_dif || fabs(f10 - m_fitFunc->Eval(10)) > max_dif2) {
379 m_fitflag = 3;
380 if (i == 9) std::cout << "ERROR XT FIT inner part" << std::endl;
381 in += 1;
382 m_fitFunc->SetParameters(p0, p1, 0, 0, 0, 0, p6default, 0);
383 m_fitFunc->SetParLimits(1, 0, 0.08);
384 m_tmin -= 0.5;
385
386 if (m_tmin < 14) {
387 m_tmin = 14; std::cout << "ERROR: tmin so small, fit iter: " << std::endl;
388 }
389 }
390
391 if (stat == 0) {
392 m_fitflag = 0;
393 if (m_debug) {std::cout << "FIT Failure; Layer: " << std::endl;}
394 }
395
396 if (m_debug) {printf("P6default= %3.2f, p6 fit = %3.2f", p6default, par[6]);}
397
398 if (m_fitflag == 1) {
399 if (m_debug) std::cout << "Fit success" << std::endl;
400 m_fitFunc->GetParameters(m_FittedXTParams);
401 m_Prob = m_fitFunc->GetProb();
402 break;
403 }
404
405 } //end loop of fitting
406 if (m_debug) B2INFO("Number of failures due to inner (outer) regions " << in << "(" << out << ")");
407
408 if (m_draw) {
409 TString hname = m_h1->GetName();
410 TString name = hname + ".pdf";
411 TCanvas* c1 = new TCanvas("c1", "", 800, 600);
412 m_h1->Draw();
413 m_fitFunc->SetLineColor(kBlack);
414 m_fitFunc->DrawF1(0, 400, "same");
415 c1->SaveAs(name);
416 }
417 }

◆ fitXT()

void fitXT ( )
inline

Do fitting.

Definition at line 277 of file XTFunction.h.

278 {
279 if (m_mode == c_Polynomial) {
280 FitPol5();
281 } else if (m_mode == c_Chebyshev) {
282 FitChebyshev();
283 } else {
284 B2ERROR("Undefined fitting function");
285 }
286 }
void FitChebyshev()
Fit xt histogram incase 5th order Chebeshev polynomial is used.
Definition: XTFunction.h:436
void FitPol5()
Fit xt histogram incase 5th order polynomial is used.
Definition: XTFunction.h:331

◆ getFitStatus()

int getFitStatus ( )
inline

get fitted flag.

Definition at line 243 of file XTFunction.h.

244 {
245 return m_fitflag;
246 }

◆ getFittedHisto()

TProfile * getFittedHisto ( )
inline

Get histogram.

Definition at line 272 of file XTFunction.h.

272{return m_h1;}

◆ getFittedXTParams()

void getFittedXTParams ( double  pa[8])
inline

get fit parameters.

Definition at line 258 of file XTFunction.h.

259 {
260 for (int i = 0; i < 8; ++i) {pa[i] = m_FittedXTParams[i];}
261 }

◆ getProb()

double getProb ( )
inline

Get the chi2 probability.

Definition at line 251 of file XTFunction.h.

252 {
253 return m_Prob;
254 }

◆ getXTFunction()

TF1 * getXTFunction ( )
inline

Get XT function.

Definition at line 265 of file XTFunction.h.

266 {
267 return m_fitFunc;
268 }

◆ isValid()

bool isValid ( )
inline

Is valid.

Definition at line 175 of file XTFunction.h.

176 {
177 if (m_fitFunc->IsValid() == true) {
178 return true;
179 } else {
180 return false;
181 }
182 }

◆ operator=()

XTFunction & operator= ( const XTFunction x)
inline

Assignment operator.

Definition at line 140 of file XTFunction.h.

141 {
142 if (this != &x) {
143 m_h1 = x.m_h1;
144 m_fitFunc = x.m_fitFunc;
145 m_mode = x.m_mode;
146 m_debug = x.m_debug;
147 m_draw = x.m_draw;
148 m_bField = x.m_bField;
149 m_minRequiredEntry = x.m_minRequiredEntry;
150 m_fitflag = x.m_fitflag;
151 m_Prob = x.m_Prob;
152 m_tmin = x.m_tmin;
153 m_tmax = x.m_tmax;
154
155 for (int i = 0; i < 8; ++i) {
156 m_XTParam[i] = x.m_XTParam[i];
157 m_FittedXTParams[i] = x.m_XTParam[i];
158 }
159 }
160 return *this;
161 }

◆ setBField()

void setBField ( bool  bfield)
inline

set to use BField

Definition at line 196 of file XTFunction.h.

196{m_bField = bfield;}

◆ setDebug()

void setDebug ( bool  debug)
inline

Set Debug.

Definition at line 236 of file XTFunction.h.

237 {
238 m_debug = debug;
239 }

◆ setFitRange()

void setFitRange ( double  tmin,
double  tmax 
)
inline

Set Fit range.

Definition at line 221 of file XTFunction.h.

222 {
223 m_tmin = tmin;
224 m_tmax = tmax;
225 }

◆ setMode()

void setMode ( int  mode)
inline

Set XT mode.

1 is 5th order Chebshev polynomial. 0 is 5th order polynomial.

Definition at line 189 of file XTFunction.h.

190 {
191 m_mode = mode;
192 }

◆ setP6()

void setP6 ( double  p6)
inline

Set Parameter 6 for polynomia fit.

Definition at line 167 of file XTFunction.h.

168 {
169 m_XTParam[6] = p6;
170 }

◆ setSmallestEntryRequired()

void setSmallestEntryRequired ( int  min)
inline

Set minimum number of entry required for fit.

Definition at line 229 of file XTFunction.h.

230 {
231 m_minRequiredEntry = min;
232 }

◆ setXTParams() [1/2]

void setXTParams ( const double  p[8])
inline

Set Parameters for fit.

Definition at line 201 of file XTFunction.h.

202 {
203 for (int i = 0; i < 8; ++i) {m_XTParam[i] = p[i];}
204 m_tmax = p[6] + 50;
205 }

◆ setXTParams() [2/2]

void setXTParams ( double  p0,
double  p1,
double  p2,
double  p3,
double  p4,
double  p5,
double  p6,
double  p7 
)
inline

Set Initial parameters for fitting.

Definition at line 209 of file XTFunction.h.

211 {
212 m_XTParam[0] = p0; m_XTParam[1] = p1; m_XTParam[2] = p2;
213 m_XTParam[3] = p3; m_XTParam[4] = p4; m_XTParam[5] = p5;
214 m_XTParam[6] = p6; m_XTParam[7] = p7;
215 m_tmax = p6 + 50;
216 }

◆ validate()

bool validate ( )

Validate the xt has proper shape.

Suppose to be bad xt if |xt(0)| > 0.2.

Definition at line 419 of file XTFunction.h.

420 {
421
422 const double p6 = m_fitFunc->GetParameter(6);
423 if (fabs(m_fitFunc->Eval(0)) > 0.3) {
424 B2WARNING("Bad xt function");
425 m_fitflag = 0;
426 return false;
427 } else if (p6 < 50.0) {
428 B2WARNING("Unrealistic p6");
429 m_fitflag = 0;
430 return false;
431 } else {
432 return true;
433 }
434 }

Member Data Documentation

◆ m_bField

bool m_bField = true
private

With magnetic field or not.

Definition at line 312 of file XTFunction.h.

◆ m_debug

bool m_debug = true
private

Print debug durring fitting or not.

Definition at line 310 of file XTFunction.h.

◆ m_draw

bool m_draw = false
private

Draw and store png plot of each histo or not.

Definition at line 311 of file XTFunction.h.

◆ m_fitflag

int m_fitflag = 0
private

Fit Flag =-1: low statitic =1: good =0: Fit failure =2: Error Outer =3: Error Inner part;.

Definition at line 325 of file XTFunction.h.

◆ m_fitFunc

TF1* m_fitFunc
private

Fit function.

Definition at line 304 of file XTFunction.h.

◆ m_FittedXTParams

double m_FittedXTParams[8] = {}
private

Fitted parameters.

Definition at line 315 of file XTFunction.h.

◆ m_h1

TProfile* m_h1
private

Histogram of xt relation.

Definition at line 303 of file XTFunction.h.

◆ m_minRequiredEntry

int m_minRequiredEntry = 800
private

Minimum entry required for each histo.

Definition at line 313 of file XTFunction.h.

◆ m_mode

int m_mode = c_Chebyshev
private

XT mode, 0 is for 5th order polynomial, 1 is Chebshev polynomial.

Definition at line 309 of file XTFunction.h.

◆ m_Prob

double m_Prob = 0
private

Chi2 prob of fitting.

Definition at line 326 of file XTFunction.h.

◆ m_tmax

double m_tmax = m_XTParam[6] + 50
private

upper boundary of fit range

Definition at line 328 of file XTFunction.h.

◆ m_tmin

double m_tmin = 20
private

lower boundary of fit range

Definition at line 327 of file XTFunction.h.

◆ m_XTParam

double m_XTParam[8] = {}
private

Parameter fo xt.

Definition at line 314 of file XTFunction.h.


The documentation for this class was generated from the following file: