Belle II Software  release-08-01-10
SpaceResolutionCalibration.h
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 <TH1F.h>
9 #include <TH2F.h>
10 #include <TGraphErrors.h>
11 #include <TF1.h>
12 #include <framework/database/DBObjPtr.h>
13 #include <cdc/dbobjects/CDCSpaceResols.h>
14 
15 namespace Belle2 {
20  namespace CDC {
21 
26  // typedef std::array<float, 3> array3; /**< angle bin info. */
27  public:
33  virtual void setDebug(bool debug = false) {m_debug = debug; }
35  virtual void setUseDB(bool useDB = false) {m_useDB = useDB; }
37  virtual void setMinimumNDF(double minndf) {m_ndfmin = minndf;}
39  virtual void setMinimumPval(double minPval) {m_Pvalmin = minPval;}
41  virtual void setBinWidth(double bw) {m_binWidth = bw;}
43  virtual void BField(bool bfield) {m_BField = bfield;}
45  virtual void inputFileNames(std::string inputname)
46  {
47  m_inputRootFileNames.assign(inputname);
48  }
50  virtual void setStoreHisto(bool storeHist = false) {m_storeHisto = storeHist;}
52  virtual void ProfileFileNames(std::string profileFileName)
53  {
54  m_ProfileFileName.assign(profileFileName);
55  }
58  {
60  }
62  virtual void setSigmaFileName(std::string name) {m_sigmafile.assign(name);}
63  // virtual void setMode(unsigned short mode=1){m_SigmaParamMode = mode;}
65  void execute()
66  {
67  calibrate();
68  }
69 
70  protected:
72  virtual bool calibrate();
74  virtual void createHisto();
76  virtual void readProfile();
78  virtual void readSigma();
80  virtual void readSigmaFromDB();
82  virtual void readSigmaFromText();
84  virtual void storeHisto();
86  virtual void write();
87 
88  private:
89  static const int Max_nalpha = 18;
90  static const int Max_ntheta = 7;
91  static const unsigned short Max_np = 40;
93  double m_ndfmin = 5;
94  double m_Pvalmin = 0.;
95  double m_binWidth = 0.05;
96  bool m_debug = false;
97  //bool m_draw = false; /**< print out histogram in pdf file or not*/
98  bool m_storeHisto = false;
99  bool m_useDB = false;
101  bool m_BField = true;
102  // bool m_LRseparate = true;
103  double sigma_old[56][2][18][7][8];
104  double sigma_new[56][2][18][7][8];
105  //TF1* ffit[56][2][18][7]; /**< fitting function*/
106  TGraphErrors* gfit[56][2][18][7];
107  TGraphErrors* gr[56][2][18][7];
108  TH2F* hist_b[56][2][Max_nalpha][Max_ntheta];
109  TH2F* hist_u[56][2][Max_nalpha][Max_ntheta];
110  TH1F* hu_m[56][2][Max_nalpha][Max_ntheta];
111  TH1F* hu_s[56][2][Max_nalpha][Max_ntheta];
112  TH1F* hb_m[56][2][Max_nalpha][Max_ntheta];
113  TH1F* hb_s[56][2][Max_nalpha][Max_ntheta];
114  int m_fitflag[56][2][Max_nalpha][Max_ntheta] = {{{{0}}}} ;
116  std::string m_outputSigmaFileName = "sigma_new.dat";
117  std::string m_inputRootFileNames = "rootfile/output*";
118  std::string m_ProfileFileName = "sigma_profile";
120  std::string m_sigmafile = "cdc/data/sigma.dat";
121  //int m_firstExperiment; /**< First experiment. */
122  //int m_firstRun; /**< First run. */
123  //int m_lastExperiment; /**< Last experiment */
124  //int m_lastRun; /**< Last run. */
125 
126  int m_nalpha;
127  int m_ntheta;
128  double l_alpha[18];
129  double u_alpha[18];
130  double ialpha[18];
131  double l_theta[7];
132  double u_theta[7];
133  double itheta[7];
138  double l_alpha_old[18];
139  double u_alpha_old[18];
140  double ialpha_old[18];
141  double l_theta_old[7];
142  double u_theta_old[7];
143  double itheta_old[7];
145  unsigned short m_sigmaParamMode_old;
146  //unsigned short m_sigmaParamMode = 1; /**< sigma mode for this calibration.*/
147 
151  double getUpperBoundaryForFit(TGraphErrors* graph)
152  {
153  double ymax = 0;
154  double xmax = 0;
155  int imax = 0;
156  double x, y;
157  int unCount = floor(0.05 / m_binWidth);
158  int N = graph->GetN();
159  int Nstart = floor(0.5 * (N - unCount));
160  int Nend = N - unCount;
161  for (int i = Nstart; i < Nend; ++i) {
162  graph->GetPoint(i, x, y);
163  if (graph->GetErrorY(i) > 0.06E-3) continue;
164  if (y > ymax) {
165  xmax = x; ymax = y;
166  imax = i;
167  }
168  }
169  if (imax <= Nstart) {
170  graph->GetPoint(Nend, x, y);
171  xmax = x;
172  }
173  return xmax;
174  }
175 
176  };
177  }
179 }
Class for Space resolution calibration.
bool m_BField
Work with BField, fit range and initial parameters is different incase B and noB.
double ialpha[18]
represented alphas of alpha bins.
virtual void setStoreHisto(bool storeHist=false)
Store histograms durring the calibration or not.
std::string m_inputRootFileNames
Input root file names.
virtual void setMinimumPval(double minPval)
Minimum Pval required.
TH1F * hb_s[56][2][Max_nalpha][Max_ntheta]
sigma histogram of ubiased residual
TGraphErrors * gr[56][2][18][7]
sigma graph.
int m_fitflag[56][2][Max_nalpha][Max_ntheta]
Fit flag; 1:OK ; 0:error.
virtual void setSigmaFileName(std::string name)
Output sigma file name, for text mode.
static const int Max_ntheta
maximum theta bin
TH2F * hist_u[56][2][Max_nalpha][Max_ntheta]
2D histogram of unbiased residual
bool m_useProfileFromInputSigma
Use binning from old sigma or new one form input.
TH1F * hu_m[56][2][Max_nalpha][Max_ntheta]
mean histogram biased residual
double u_alpha[18]
Upper boundays of alpha bins.
virtual void setMinimumNDF(double minndf)
minimum NDF required for track
virtual void BField(bool bfield)
Work with B field or not;.
virtual void readProfile()
read sigma bining (alpha, theta bining)
double sigma_old[56][2][18][7][8]
old sigma prameters.
void execute()
execute all, make the interface the same as CAF
virtual void inputFileNames(std::string inputname)
Input root file names, results of collector module.
TGraphErrors * gfit[56][2][18][7]
sigma*sigma graph for fit
DBObjPtr< CDCSpaceResols > * m_sResolFromDB
Database for sigma.
virtual void useProfileFromInputSigma(bool useProfileFromInputSigma)
use sigma bin profile form input sigma or new one from input file
double l_theta_old[7]
Lower boundays of theta bins from input.
double l_alpha[18]
Lower boundays of alpha bins.
virtual void readSigmaFromDB()
read sigma from DB
double m_binWidth
width of each bin, unit cm
double getUpperBoundaryForFit(TGraphErrors *graph)
search max point at boundary region
virtual void setBinWidth(double bw)
Bin width of each slide.
TH2F * hist_b[56][2][Max_nalpha][Max_ntheta]
2D histogram of biased residual
int nalpha_old
number of alpha bins from input
std::string m_outputSigmaFileName
Output sigma file name.
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
double itheta_old[7]
represented alphas of theta bins from input.
virtual void setUseDB(bool useDB=false)
Use database or text mode.
static const unsigned short Max_np
Maximum number of point =1/binwidth.
std::string m_sigmafile
Sigma file name, for text mode.
double l_alpha_old[18]
Lower boundays of alpha bins from input.
virtual void ProfileFileNames(std::string profileFileName)
File name describe theta/alpha bin, if don't want to use default from input sigma.
double u_theta_old[7]
Upper boundays of theta bins from input.
virtual void write()
save calibration, in text file or db
double itheta[7]
represented alphas of theta bins.
double ialpha_old[18]
represented alphas of alpha bins from input.
double l_theta[7]
Lower boundays of theta bins.
double sigma_new[56][2][18][7][8]
new sigma prameters.
virtual void readSigmaFromText()
read sigma from text file
static const int Max_nalpha
Maximum alpha bin.
double m_Pvalmin
Minimum Prob(chi2) of track.
double u_alpha_old[18]
Upper boundays of alpha bins from input.
virtual void setDebug(bool debug=false)
Debug or not.
virtual void readSigma()
read sigma from previous calibration, (input sigma)
double u_theta[7]
Upper boundays of theta bins.
unsigned short m_sigmaParamMode_old
sigma mode from input.
int ntheta_old
number of theta bins from input
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Abstract base class for different kinds of events.