Belle II Software  release-05-02-19
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_draw = false
 print out histogram in pdf file 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.
 
TF1 * ffit [56][2][18][7]
 fitting function
 
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_firstExperiment
 First experiment.
 
int m_firstRun
 First run.
 
int m_lastExperiment
 Last experiment.
 
int m_lastRun
 Last run.
 
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.
 
unsigned short m_sigmaParamMode = 1
 sigma mode for this calibration.
 

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 18 of file SpaceResolutionCalibration.h.

Member Function Documentation

◆ calibrate()

bool calibrate ( )
protectedvirtual

Run algo on data.

Upper limit of fitting.

Definition at line 280 of file SpaceResolutionCalibration.cc.

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

◆ readSigmaFromDB()

void readSigmaFromDB ( )
protectedvirtual

read sigma from DB

< angle bin info.

Definition at line 566 of file SpaceResolutionCalibration.cc.


The documentation for this class was generated from the following files:
Belle2::CDC::SpaceResolutionCalibration::createHisto
virtual void createHisto()
create histogram
Definition: SpaceResolutionCalibration.cc:31
Belle2::CDC::SpaceResolutionCalibration::write
virtual void write()
save calibration, in text file or db
Definition: SpaceResolutionCalibration.cc:399
Belle2::CDC::SpaceResolutionCalibration::sigma_new
double sigma_new[56][2][18][7][8]
new sigma prameters.
Definition: SpaceResolutionCalibration.h:97
Belle2::CDC::SpaceResolutionCalibration::m_fitflag
int m_fitflag[56][2][Max_nalpha][Max_ntheta]
Fit flag; 1:OK ; 0:error.
Definition: SpaceResolutionCalibration.h:107
Belle2::CDC::SpaceResolutionCalibration::storeHisto
virtual void storeHisto()
store histogram
Definition: SpaceResolutionCalibration.cc:366
Belle2::CDC::SpaceResolutionCalibration::gfit
TGraphErrors * gfit[56][2][18][7]
sigma*sigma graph for fit
Definition: SpaceResolutionCalibration.h:99
Belle2::CDC::SpaceResolutionCalibration::m_ntheta
int m_ntheta
number of theta bins
Definition: SpaceResolutionCalibration.h:120
Belle2::CDC::SpaceResolutionCalibration::m_nalpha
int m_nalpha
number of alpha bins
Definition: SpaceResolutionCalibration.h:119
Belle2::CDC::SpaceResolutionCalibration::getUpperBoundaryForFit
double getUpperBoundaryForFit(TGraphErrors *graph)
search max point at boundary region
Definition: SpaceResolutionCalibration.h:144