Belle II Software  release-08-01-10
SpaceResolutionCalibration Class Reference

Class for Space resolution calibration. More...

#include <SpaceResolutionCalibration.h>

Collaboration diagram for SpaceResolutionCalibration:

Public Member Functions

 SpaceResolutionCalibration ()
 Constructor.
 
virtual ~SpaceResolutionCalibration ()
 Destructor.
 
virtual void setDebug (bool debug=false)
 Debug or not.
 
virtual void setUseDB (bool useDB=false)
 Use database or text mode.
 
virtual void setMinimumNDF (double minndf)
 minimum NDF required for track
 
virtual void setMinimumPval (double minPval)
 Minimum Pval required.
 
virtual void setBinWidth (double bw)
 Bin width of each slide.
 
virtual void BField (bool bfield)
 Work with B field or not;.
 
virtual void inputFileNames (std::string inputname)
 Input root file names, results of collector module.
 
virtual void setStoreHisto (bool storeHist=false)
 Store histograms durring the calibration or not.
 
virtual void ProfileFileNames (std::string profileFileName)
 File name describe theta/alpha bin, if don't want to use default from input sigma.
 
virtual void useProfileFromInputSigma (bool useProfileFromInputSigma)
 use sigma bin profile form input sigma or new one from input file
 
virtual void setSigmaFileName (std::string name)
 Output sigma file name, for text mode.
 
void execute ()
 execute all, make the interface the same as CAF
 

Protected Member Functions

virtual bool calibrate ()
 Run algo on data. More...
 
virtual void createHisto ()
 create histogram
 
virtual void readProfile ()
 read sigma bining (alpha, theta bining)
 
virtual void readSigma ()
 read sigma from previous calibration, (input sigma)
 
virtual void readSigmaFromDB ()
 read sigma from DB More...
 
virtual void readSigmaFromText ()
 read sigma from text file
 
virtual void storeHisto ()
 store histogram
 
virtual void write ()
 save calibration, in text file or db
 

Private Member Functions

double getUpperBoundaryForFit (TGraphErrors *graph)
 search max point at boundary region
 

Private Attributes

double m_ndfmin = 5
 Minimum NDF

 
double m_Pvalmin = 0.
 Minimum Prob(chi2) of track.
 
double m_binWidth = 0.05
 width of each bin, unit cm
 
bool m_debug = false
 Debug or not.
 
bool m_storeHisto = false
 Store histogram or not.
 
bool m_useDB = false
 use db or text mode
 
bool m_useProfileFromInputSigma = true
 Use binning from old sigma or new one form input.
 
bool m_BField = true
 Work with BField, fit range and initial parameters is different incase B and noB.
 
double sigma_old [56][2][18][7][8]
 old sigma prameters.
 
double sigma_new [56][2][18][7][8]
 new sigma prameters.
 
TGraphErrors * gfit [56][2][18][7]
 sigma*sigma graph for fit
 
TGraphErrors * gr [56][2][18][7]
 sigma graph.
 
TH2F * hist_b [56][2][Max_nalpha][Max_ntheta]
 2D histogram of biased residual
 
TH2F * hist_u [56][2][Max_nalpha][Max_ntheta]
 2D histogram of unbiased residual
 
TH1F * hu_m [56][2][Max_nalpha][Max_ntheta]
 mean histogram biased residual
 
TH1F * hu_s [56][2][Max_nalpha][Max_ntheta]
 sigma histogram of biased residual
 
TH1F * hb_m [56][2][Max_nalpha][Max_ntheta]
 mean histogram of unbiased residual
 
TH1F * hb_s [56][2][Max_nalpha][Max_ntheta]
 sigma histogram of ubiased residual
 
int m_fitflag [56][2][Max_nalpha][Max_ntheta] = {{{{0}}}}
 Fit flag; 1:OK ; 0:error.
 
std::string m_outputSigmaFileName = "sigma_new.dat"
 Output sigma file name.
 
std::string m_inputRootFileNames = "rootfile/output*"
 Input root file names.
 
std::string m_ProfileFileName = "sigma_profile"
 Profile file name.
 
DBObjPtr< CDCSpaceResols > * m_sResolFromDB
 Database for sigma.
 
std::string m_sigmafile = "cdc/data/sigma.dat"
 Sigma file name, for text mode.
 
int m_nalpha
 number of alpha bins
 
int m_ntheta
 number of theta bins
 
double l_alpha [18]
 Lower boundays of alpha bins.
 
double u_alpha [18]
 Upper boundays of alpha bins.
 
double ialpha [18]
 represented alphas of alpha bins.
 
double l_theta [7]
 Lower boundays of theta bins.
 
double u_theta [7]
 Upper boundays of theta bins.
 
double itheta [7]
 represented alphas of theta bins.
 
int nalpha_old
 number of alpha bins from input
 
int ntheta_old
 number of theta bins from input
 
double l_alpha_old [18]
 Lower boundays of alpha bins from input.
 
double u_alpha_old [18]
 Upper boundays of alpha bins from input.
 
double ialpha_old [18]
 represented alphas of alpha bins from input.
 
double l_theta_old [7]
 Lower boundays of theta bins from input.
 
double u_theta_old [7]
 Upper boundays of theta bins from input.
 
double itheta_old [7]
 represented alphas of theta bins from input.
 
unsigned short m_sigmaParamMode_old
 sigma mode from input.
 

Static Private Attributes

static const int Max_nalpha = 18
 Maximum alpha bin.
 
static const int Max_ntheta = 7
 maximum theta bin

 
static const unsigned short Max_np = 40
 Maximum number of point =1/binwidth.
 

Detailed Description

Class for Space resolution calibration.

Definition at line 25 of file SpaceResolutionCalibration.h.

Member Function Documentation

◆ calibrate()

bool calibrate ( )
protectedvirtual

Run algo on data.

Upper limit of fitting.

Definition at line 286 of file SpaceResolutionCalibration.cc.

287 {
288  gROOT->SetBatch(1);
289  gErrorIgnoreLevel = 3001;
290  createHisto();
291  B2INFO("Start to calibrate");
292  gSystem->Exec("mkdir -p Sigma_Fit_Err"); //create a folder to store error histo
293 
294  TF1* func = new TF1("func", "[0]/(x*x + [1])+[2]* x+[3]+[4]*exp([5]*(x-[6])*(x-[6]))", 0, 1.);
295  TH1F* hprob = new TH1F("h1", "", 20, 0, 1);
296  double upFit;
297  double intp6;
298 
299  for (int i = 0; i < 56; ++i) {
300  for (int lr = 0; lr < 2; ++lr) {
301  for (int al = 0; al < m_nalpha; ++al) {
302  for (int th = 0; th < m_ntheta; ++th) {
303  if (m_fitflag[i][lr][al][th] != -1) { /*if graph exist, do fitting*/
304 
305 
306  upFit = getUpperBoundaryForFit(gfit[i][lr][al][th]);
307  intp6 = upFit + 0.2;
308  B2DEBUG(199, "xmax for fitting: " << upFit);
309  // func->SetLineColor(1 + lr + al * 2 + th * 3);
310  func->SetParameters(5E-6, 0.007, 1E-4, 1E-5, 0.00008, -30, intp6);
311  func->SetParLimits(0, 1E-7, 1E-4);
312  func->SetParLimits(1, 0.00001, 0.02);
313  func->SetParLimits(2, 1E-6, 0.0005);
314  func->SetParLimits(3, 1E-8, 0.0005);
315  func->SetParLimits(4, 0., 0.001);
316  func->SetParLimits(5, -40, 0.);
317  func->SetParLimits(6, intp6 - 0.5, intp6 + 0.3);
318  B2DEBUG(199, "FITTING for layer: " << i << "lr: " << lr << " ial" << al << " ith:" << th);
319  B2DEBUG(199, "Fit flag before fit:" << m_fitflag[i][lr][al][th]);
320  // if(!(gfit[i][lr][al][th]->isValid())) continue;
321  // m_fitflag[i][lr][al][th] = 0;
322  for (int j = 0; j < 10; j++) {
323 
324  B2DEBUG(199, "loop: " << j);
325  B2DEBUG(199, "Int p6: " << intp6);
326  B2DEBUG(199, "Number of Point: " << gfit[i][lr][al][th]->GetN());
327  Int_t stat = gfit[i][lr][al][th]->Fit("func", "MQE", "", 0.05, upFit);
328  B2DEBUG(199, "stat of fit" << stat);
329  std::string Fit_status = gMinuit->fCstatu.Data();
330  B2DEBUG(199, "FIT STATUS: " << Fit_status);
331  // stat=gfit[i]->Fit(Form("ffit[%d]",i),"M "+Q,"",0.0,cellsize(i)+0.05+j*0.005);
332  if (Fit_status == "OK" ||
333  Fit_status == "SUCCESSFUL" ||
334  Fit_status == "CALL LIMIT" ||
335  Fit_status == "PROBLEMS") {
336  if (fabs(func->Eval(0.3)) > 0.00035 || func->Eval(0.3) < 0) {
337  func->SetParameters(5E-6, 0.007, 1E-4, 1E-7, 0.0007, -30, intp6 + 0.05 * j);
338  // func->SetParameters(defaultparsmall);
339  m_fitflag[i][lr][al][th] = 0;
340  } else {
341  B2DEBUG(199, "Prob of fit: " << func->GetProb());
342  m_fitflag[i][lr][al][th] = 1;
343  break;
344  }
345  } else {
346  m_fitflag[i][lr][al][th] = 0;
347  func->SetParameters(5E-6, 0.007, 1E-4, 1E-7, 0.0007, -30, intp6 + 0.05 * j);
348  // upFit += 0.025;
349  if (j == 9) {
350  TCanvas* c1 = new TCanvas("c1", "", 600, 600);
351  gfit[i][lr][al][th]->Draw();
352  c1->SaveAs(Form("Sigma_Fit_Err/%d_%d_al%d_th%d_%s.png", i, lr, al, th, Fit_status.c_str()));
353  B2WARNING("Fit error: " << i << " " << lr << " " << al << " " << th);
354  }
355  }
356  }
357  if (m_fitflag[i][lr][al][th] == 1) {
358  B2DEBUG(199, "ProbFit: Lay_lr_al_th: " << i << " " << lr << " " << al << " " << th << func->GetProb());
359  hprob->Fill(func->GetProb());
360  func->GetParameters(sigma_new[i][lr][al][th]);
361  }
362  }
363  }
364  }
365  }
366  }
367  storeHisto();
368  write();
369 
370  return true;
371 }
R E
internal precision of FFTW codelets
int m_fitflag[56][2][Max_nalpha][Max_ntheta]
Fit flag; 1:OK ; 0:error.
TGraphErrors * gfit[56][2][18][7]
sigma*sigma graph for fit
double getUpperBoundaryForFit(TGraphErrors *graph)
search max point at boundary region
virtual void write()
save calibration, in text file or db
double sigma_new[56][2][18][7][8]
new sigma prameters.

◆ readSigmaFromDB()

void readSigmaFromDB ( )
protectedvirtual

read sigma from DB

< angle bin info.

Definition at line 570 of file SpaceResolutionCalibration.cc.


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