Belle II Software development
XTCalibration Class Reference

Class to perform xt calibration for drift chamber. More...

#include <XTCalibration.h>

Public Member Functions

 XTCalibration ()
 Constructor.
 
virtual ~XTCalibration ()
 Destructor.
 
virtual void BField (bool bfield)
 set to use BField
 
virtual void setDebug (bool debug=false)
 Run in debug or silent.
 
virtual void setUseDB (bool useDB=false)
 Set to run with database mode or text mode.
 
virtual void setMinimumNDF (double minndf)
 set minimum number of degree of freedom requirement
 
virtual void setMinimumPval (double minPval)
 set minimum Prob(Chi2) requirement
 
virtual void inputFileNames (std::string inputname)
 Input root file name, output of collector.
 
virtual void profileFileNames (std::string profileFileName)
 Profile file name incase you want to change the xt binning.
 
virtual void useProfileFromInputXT (bool useProfileXTFromInputXT)
 if you want to change xt bining, you have to set this to true
 
virtual void setXTFileName (std::string name)
 input xt file name incase text mode is used.
 
virtual void setMode (unsigned short mode=1)
 set xt mode, 0 is polynimial, 1 is Chebshev polynomial
 
virtual void setStoreHisto (bool storeHist=false)
 set to store histogram or not.
 
void setLRSeparate (bool lr=true)
 Set LR separate mode (default is true).
 
void execute ()
 Run calibration.
 

Protected Member Functions

virtual bool calibrate ()
 Run algo on data.
 
virtual void readXTFromDB ()
 Read old xt parameter from database.
 
virtual void readXTFromText ()
 Read old xt parameter from text file, incase text mode is used.
 
virtual void CreateHisto ()
 Create histogram for calibration.
 
virtual void readProfile ()
 Read profile xt file.
 
virtual void readXT ()
 read xt paramter (wrap text mode and database mode)
 
virtual void Write ()
 Store calibrated constand.
 
virtual void storeHisto ()
 Store histogram to file.
 

Private Attributes

double m_ndfmin = 5
 minimum ndf required
 
double m_Pvalmin = 0.
 minimum pvalue required
 
bool m_debug = false
 run in debug or silent
 
bool m_storeHisto = false
 Store histogram or not.
 
bool m_useDB = false
 Use Database or text mode.
 
bool m_useProfileXTFromInputXT = true
 use profile from text file or default in input xt
 
bool m_LRseparate = true
 Separate LR in calibration or mix.
 
bool m_useSliceFit = false
 Use slice fit or profile.
 
bool m_BField = true
 with b field or none
 
double xtold [56][2][18][7][8]
 Old paremeter.
 
int fitflag [56][2][20][10]
 Fit flag.
 
TF1 * xtf5r [56][2][20][10]
 XTFunction.
 
TProfile * hprof [56][2][20][10]
 Profile xt histo.
 
TH2D * hist2d [56][2][20][10]
 2D histo of xt
 
TH2D * hist2d_draw [56][20][10]
 2d histo for draw
 
TH1D * hist2d_1 [56][2][20][10]
 1D xt histo, results of slice fit
 
std::string m_OutputXTFileName = "xt_new.dat"
 Out put xt filename.
 
std::string m_inputRootFileNames = "rootfile/output*"
 input root filename
 
std::string m_profileFileName = "xt_profile"
 profile file name
 
std::string m_xtfile = "cdc/data/xt.dat"
 Input xt file name, incase 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 xtmode_old
 XT mode old, 0-polynomial, 1 Cheb.
 
int m_MAXalpha = 18
 max alpha bin
 
int m_MAXtheta = 7
 max theta bin
 
unsigned short m_xtmode = 1
 Mode of xt; 0 is polynomial;1 is Chebyshev.
 
int m_smallestEntryRequire = 1000
 minimum number of hit per hitosgram.
 
double m_par6 [56]
 boundary parameter for fitting, semi-experiment number
 

Detailed Description

Class to perform xt calibration for drift chamber.

Definition at line 21 of file XTCalibration.h.

Constructor & Destructor Documentation

◆ XTCalibration()

Constructor.

Definition at line 37 of file XTCalibration.cc.

39{
40 /*
41 setDescription(
42 " -------------------------- Test Calibration Algoritm -------------------------\n"
43 " \n"
44 " Testing algorithm which just gets mean of a test histogram collected by \n"
45 0 " CaTest module and provides a DB object with another histogram with one \n"
46 " entry at calibrated value. \n"
47 " ------------------------------------------------------------------------------\n"
48 );
49 */
50}

◆ ~XTCalibration()

virtual ~XTCalibration ( )
inlinevirtual

Destructor.

Definition at line 26 of file XTCalibration.h.

26{}

Member Function Documentation

◆ BField()

virtual void BField ( bool  bfield)
inlinevirtual

set to use BField

Definition at line 28 of file XTCalibration.h.

28{m_BField = bfield;}
bool m_BField
with b field or none
Definition: XTCalibration.h:86

◆ calibrate()

bool calibrate ( )
protectedvirtual

Run algo on data.

Definition at line 200 of file XTCalibration.cc.

201{
202 gROOT->SetBatch(1);
203 gErrorIgnoreLevel = 3001;
204
205 CreateHisto();
206 B2INFO("Start Fitting");
207 for (int l = 0; l < 56; ++l) {
208 for (int lr = 0; lr < 2; ++lr) {
209 for (int al = 0; al < m_nalpha; ++al) {
210 for (int th = 0; th < m_ntheta; ++th) {
211
212 if (hist2d[l][lr][al][th]->GetEntries() < m_smallestEntryRequire) {
213 fitflag[l][lr][al][th] = -1;
214 continue;
215 }
216 double p0, p1, tmin;
217 TF1* fpol1;
218 if (m_useSliceFit) {
219 hist2d[l][lr][al][th]->FitSlicesY(0, 0, -1, 5);
220 hist2d_1[l][lr][al][th] = (TH1D*)gDirectory->Get(Form("h%d_%d_%d_%d_1", l, lr, al, th));
221 if (!hist2d_1[l][lr][al][th]) {
222 fitflag[l][lr][al][th] = -1;
223 B2WARNING("error, not found results of slices fit");
224 continue;
225 }
226 hist2d_1[l][lr][al][th]->Fit("pol1", "Q", "", 30, 60);
227 fpol1 = (TF1*)hprof[l][lr][al][th]->GetFunction("pol1");
228 } else {
229 /*Set Error for low statistic bin*/
230 for (int n = 0; n < hprof[l][lr][al][th]->GetNbinsX(); ++n) {
231 if (hprof[l][lr][al][th]->GetBinEntries(n) < 5 && hprof[l][lr][al][th]->GetBinEntries(n) > 1) {
232 hprof[l][lr][al][th]->SetBinError(n, 0.3 / hprof[l][lr][al][th]->GetBinEntries(n));
233 }
234 }
235 hprof[l][lr][al][th]->Fit("pol1", "Q", "", 30, 60);
236 fpol1 = (TF1*)hprof[l][lr][al][th]->GetFunction("pol1");
237 }
238 if (fpol1) {
239 //determine tmin in fitting
240 p0 = fpol1->GetParameter(0);
241 p1 = fpol1->GetParameter(1);
242 tmin = -1 * p0 / p1 + 15;
243 } else {
244 p0 = 0; p1 = 0.005;
245 tmin = 12;
246 }
247 XT* xt(nullptr);
248 if (m_useSliceFit) {
249 xt = new XT(hist2d_1[l][lr][al][th], m_xtmode);
250 } else {
251 xt = new XT(hprof[l][lr][al][th], m_xtmode);
252 }
253 // xt->setSmallestEntryRequired(m_smallestEntryRequire);
254 if (m_BField) {
255 int ial_old = 0;
256 int ith_old = 0;
257 for (int k = 0; k < nalpha_old; ++k) {
258 if (ialpha[al] < u_alpha_old[k]) {ial_old = k; break;}
259 }
260 for (int j = 0; j < ntheta_old; ++j) {
261 if (itheta[th] < u_theta_old[j]) {ith_old = j; break;}
262 }
263 double p6 = xtold[l][lr][ial_old][ith_old][6];
264 if (p6 > 400)
265 p6 = 400;
266 if (m_xtmode == xtmode_old) {
267 xt->setXTParams(xtold[l][lr][ial_old][ith_old]);
268 xt->setP6(p6);
269 } else {
270 xt->setXTParams(p0, p1, 0., 0., 0., 0., p6, xtold[l][lr][ial_old][ith_old][7]);
271 }
272 xt->setFitRange(tmin, p6 + 100);
273 } else {
274 xt->setXTParams(p0, p1, 0., 0., 0., 0., m_par6[l], 0.0001);
275 xt->setFitRange(tmin, m_par6[l] + 100);
276 }
277 xt->setDebug(m_debug);
278 xt->BField(m_BField);
279 xt->FitXT(m_xtmode);
280 /*get result*/
281 fitflag[l][lr][al][th] = xt->getFitStatus();
282 xtf5r[l][lr][al][th] = (TF1*)xt->getXTFunction();
283 if (m_useSliceFit) {
284 hist2d_1[l][lr][al][th] = (TH1D*)xt->getFittedHisto();
285 } else {
286 hprof[l][lr][al][th] = (TProfile*)xt->getFittedHisto();
287 }
288 delete xt;
289 }
290 }
291 }
292 }
293 Write();
294 storeHisto();
295 return true;
296}
unsigned short xtmode_old
XT mode old, 0-polynomial, 1 Cheb.
virtual void storeHisto()
Store histogram to file.
double ialpha[18]
represented alphas of alpha bins.
double m_par6[56]
boundary parameter for fitting, semi-experiment number
double xtold[56][2][18][7][8]
Old paremeter.
Definition: XTCalibration.h:89
int m_nalpha
number of alpha bins
TProfile * hprof[56][2][20][10]
Profile xt histo.
Definition: XTCalibration.h:93
bool m_debug
run in debug or silent
Definition: XTCalibration.h:80
int m_ntheta
number of theta bins
virtual void Write()
Store calibrated constand.
int nalpha_old
number of alpha bins from input
bool m_useSliceFit
Use slice fit or profile.
Definition: XTCalibration.h:85
TH1D * hist2d_1[56][2][20][10]
1D xt histo, results of slice fit
Definition: XTCalibration.h:96
double u_theta_old[7]
Upper boundays of theta bins from input.
virtual void CreateHisto()
Create histogram for calibration.
double itheta[7]
represented alphas of theta bins.
unsigned short m_xtmode
Mode of xt; 0 is polynomial;1 is Chebyshev.
int fitflag[56][2][20][10]
Fit flag.
Definition: XTCalibration.h:90
TF1 * xtf5r[56][2][20][10]
XTFunction.
Definition: XTCalibration.h:91
TH2D * hist2d[56][2][20][10]
2D histo of xt
Definition: XTCalibration.h:94
double u_alpha_old[18]
Upper boundays of alpha bins from input.
int ntheta_old
number of theta bins from input
int m_smallestEntryRequire
minimum number of hit per hitosgram.
Class to perform fitting for each xt function.
Definition: XT.h:60

◆ CreateHisto()

void CreateHisto ( )
protectedvirtual

Create histogram for calibration.

Definition at line 52 of file XTCalibration.cc.

53{
54 readXT();
56 /* read data from tree, make histo for fit*/
57 TChain* tree = new TChain("tree");
58 tree->Add(m_inputRootFileNames.c_str());
59 B2INFO("Open Files: " << m_inputRootFileNames);
60 if (!tree->GetBranch("ndf")) {
61 cout << "input data do not exits, please check!" << endl;
62 B2FATAL("echo rootfile do not exits or something wrong");
63 gSystem->Exec("echo rootfile do not exits or something wrong >> error");
64 return;
65 }
66 int lay;
67 double dt;
68 double dx;
69 double Pval, alpha, theta;
70 double ndf;
71
72 tree->SetBranchAddress("lay", &lay);
73 tree->SetBranchAddress("t", &dt);
74 tree->SetBranchAddress("x_u", &dx);
75 tree->SetBranchAddress("alpha", &alpha);
76 tree->SetBranchAddress("theta", &theta);
77 tree->SetBranchAddress("Pval", &Pval);
78 tree->SetBranchAddress("ndf", &ndf);
79
80 /* Disable unused branch */
81 std::vector<TString> list_vars = {"lay", "t", "x_u", "alpha", "theta", "Pval", "ndf"};
82 tree->SetBranchStatus("*", 0);
83
84 for (TString brname : list_vars) {
85 tree->SetBranchStatus(brname, 1);
86 }
87
88
89 /*Create histogram*/
90 for (int i = 0; i < 56; ++i) {
91 for (int lr = 0; lr < 2; ++lr) {
92 for (int al = 0; al < m_nalpha; ++al) {
93 for (int th = 0; th < m_ntheta; ++th) {
94 hprof[i][lr][al][th] = new TProfile(Form("hprof%d_%d_%d_%d", i, lr, al, th),
95 Form("(L=%d)-(lr=%d)-(#alpha=%3.0f)-(#theta=%3.0f); Drift time (ns);Drift Length (cm)",
96 i, lr, ialpha[al], itheta[th]), 210, -20, 600, 0, 1.2, "i");
97 hist2d[i][lr][al][th] = new TH2D(Form("h%d_%d_%d_%d", i, lr, al, th),
98 Form("(L=%d)-(lr=%d)-(#alpha=%3.0f)-(#theta=%3.0f); Drift time (ns);Drift Length (cm)",
99 i, lr, ialpha[al], itheta[th]), 210, -20, 600, 110, 0, 1.2);
100 if (lr == 1)
101 hist2d_draw[i][al][th] = new TH2D(Form("h_draw%d_%d_%d", i, al, th),
102 Form("(L=%d)-(#alpha=%3.0f)-(#theta=%3.0f); Drift time (ns);Drift Length (cm)",
103 i, ialpha[al], itheta[th]), 210, -20, 600, 2200, -1.2, 1.2);
104 }
105 }
106 }
107 }
108
109 /*Now read data and make histo*/
110 int al = 0;
111 int th = 0;
112 int lr = 0;
113 int nEntries = tree->GetEntries();
114 cout << "Number of Hit: " << nEntries << endl;
115
116 for (int i = 0; i < nEntries; ++i) {
117 tree->GetEntry(i);
118 /* protect in case |alpha|>90*/
119 if (fabs(alpha) > 90) {
120 if (alpha < 0) alpha += 180;
121 if (alpha > 0) alpha -= 180;
122 }
123
124 if (Pval < m_Pvalmin) continue;
125 if (ndf < m_ndfmin) continue;
126
127 for (int k = 0; k < m_nalpha; ++k) {
128 if (alpha < u_alpha[k]) {al = k; break;}
129 }
130 for (int j = 0; j < m_ntheta; ++j) {
131 if (theta < u_theta[j]) {th = j; break;}
132 }
133 if (dx > 0)
134 lr = 1;
135 else
136 lr = 0;
137 if (m_LRseparate) {
138 hprof[lay][lr][al][th]->Fill(dt, abs(dx));
139 hist2d[lay][lr][al][th]->Fill(dt, abs(dx));
140 } else {
141 hprof[lay][0][al][th]->Fill(dt, abs(dx));
142 hist2d[lay][0][al][th]->Fill(dt, abs(dx));
143 hprof[lay][1][al][th]->Fill(dt, abs(dx));
144 hist2d[lay][1][al][th]->Fill(dt, abs(dx));
145 }
146 hist2d_draw[lay][al][th]->Fill(dt, dx);
147 }
148}
std::string m_inputRootFileNames
input root filename
double u_alpha[18]
Upper boundays of alpha bins.
virtual void readXT()
read xt paramter (wrap text mode and database mode)
bool m_LRseparate
Separate LR in calibration or mix.
Definition: XTCalibration.h:84
virtual void readProfile()
Read profile xt file.
double m_ndfmin
minimum ndf required
Definition: XTCalibration.h:78
TH2D * hist2d_draw[56][20][10]
2d histo for draw
Definition: XTCalibration.h:95
double m_Pvalmin
minimum pvalue required
Definition: XTCalibration.h:79
double u_theta[7]
Upper boundays of theta bins.

◆ execute()

void execute ( )
inline

Run calibration.

Definition at line 54 of file XTCalibration.h.

55 {
56 calibrate();
57 }
virtual bool calibrate()
Run algo on data.

◆ inputFileNames()

virtual void inputFileNames ( std::string  inputname)
inlinevirtual

Input root file name, output of collector.

Definition at line 38 of file XTCalibration.h.

38{m_inputRootFileNames.assign(inputname);}

◆ profileFileNames()

virtual void profileFileNames ( std::string  profileFileName)
inlinevirtual

Profile file name incase you want to change the xt binning.

Definition at line 40 of file XTCalibration.h.

40{m_profileFileName.assign(profileFileName);}
std::string m_profileFileName
profile file name

◆ readProfile()

void readProfile ( )
protectedvirtual

Read profile xt file.

Definition at line 149 of file XTCalibration.cc.

150{
151 /*Read profile for xt*/
153 B2INFO("use XT bining from input XT");
156 B2INFO("Nalpha: " << m_nalpha << "\n Ntheta: " << m_ntheta);
157 for (int i = 0; i < m_nalpha; ++i) {
158 l_alpha[i] = l_alpha_old[i];
159 u_alpha[i] = u_alpha_old[i];
160 ialpha[i] = ialpha_old[i];
161 B2INFO("-" << l_alpha[i] << " " << u_alpha[i] << " " << ialpha[i]);
162 }
163 for (int i = 0; i < m_ntheta; ++i) {
164 l_theta[i] = l_theta_old[i];
165 u_theta[i] = u_theta_old[i];
166 itheta[i] = itheta_old[i];
167 B2INFO("-" << l_theta[i] << " " << u_theta[i] << " " << itheta[i]);
168 }
169 } else {
170 B2INFO("use XT bining from profile file");
171 ifstream proxt(m_profileFileName.c_str());
172 if (!proxt) {
173 B2FATAL("file not found: " << m_profileFileName);
174 }
175 double dumy1, dumy2, dumy3;
176 proxt >> m_nalpha;
177 B2DEBUG(99, "Number of alpha bin" << m_nalpha);
178 if (m_nalpha > m_MAXalpha) {
179 B2FATAL("number of alpha bin excess limit; please increse uplimit: " << m_nalpha << " > " << m_MAXalpha);
180 }
181 for (int i = 0; i < m_nalpha; ++i) {
182 proxt >> dumy1 >> dumy2 >> dumy3;
183 l_alpha[i] = dumy1;
184 u_alpha[i] = dumy2;
185 ialpha[i] = dumy3;
186 }
187 proxt >> m_ntheta;
188 B2DEBUG(99, "Number of theta bin" << m_nalpha);
189 if (m_ntheta > m_MAXtheta) {B2FATAL("number of theta bin excess limit; please increse uplimit: " << m_ntheta << " > " << m_MAXtheta);}
190 for (int i = 0; i < m_ntheta; ++i) {
191 proxt >> dumy1 >> dumy2 >> dumy3;
192 l_theta[i] = dumy1;
193 u_theta[i] = dumy2;
194 itheta[i] = dumy3;
195 }
196 }
197 B2INFO("Finish asssign XT bining");
198}
bool m_useProfileXTFromInputXT
use profile from text file or default in input xt
Definition: XTCalibration.h:83
int m_MAXalpha
max alpha bin
double l_theta_old[7]
Lower boundays of theta bins from input.
double l_alpha[18]
Lower boundays of alpha bins.
double itheta_old[7]
represented alphas of theta bins from input.
double l_alpha_old[18]
Lower boundays of alpha bins from input.
double ialpha_old[18]
represented alphas of alpha bins from input.
double l_theta[7]
Lower boundays of theta bins.
int m_MAXtheta
max theta bin

◆ readXT()

void readXT ( )
protectedvirtual

read xt paramter (wrap text mode and database mode)

B2FATAL("Error reading xt from DB");return;}

Definition at line 409 of file XTCalibration.cc.

410{
411 if (m_useDB) {
412 B2INFO("reading xt from DB");
413
414 /*
415 ReadXT:readXTFromDB(&xtold,dbXT_old,
416 &nalpha_old,l_alpha_old,u_alpha_old,ialpha_old,
417 &ntheta_old,l_theta_old,u_theta_old, itheta_old);
418 */
419 readXTFromDB();
420 B2INFO("Number of theta bin from xt: " << ntheta_old);
421 B2INFO("Theta 0: " << itheta_old[0]);
422 // if(!a){
424 } else {
425 B2INFO("Read Xt from text");
426 /*
427 ReadXT::readXTFromText(xtold,m_xtfile,
428 nalpha_old,l_alpha_old,u_alpha_old,ialpha_old,
429 ntheta_old,l_theta_old,u_theta_old, itheta_old);
430 */
432 B2INFO("nalpha: " << nalpha_old);
433 // if(!a)
434 // {B2FATAL("Error reading xt from text");return;}
435 }
436}
virtual void readXTFromText()
Read old xt parameter from text file, incase text mode is used.
virtual void readXTFromDB()
Read old xt parameter from database.
bool m_useDB
Use Database or text mode.
Definition: XTCalibration.h:82

◆ readXTFromDB()

void readXTFromDB ( )
protectedvirtual

Read old xt parameter from database.

Definition at line 538 of file XTCalibration.cc.

539{
541 nalpha_old = dbXT_old->getNoOfAlphaBins();
542 B2INFO("Number of alpha" << nalpha_old);
543 double rad2deg = 180 / M_PI;
544 for (unsigned short i = 0; i < nalpha_old; ++i) {
545 array3 alpha = dbXT_old->getAlphaBin(i);
546 l_alpha_old[i] = alpha[0] * rad2deg;
547 u_alpha_old[i] = alpha[1] * rad2deg;
548 ialpha_old[i] = alpha[2] * rad2deg;
549 // std::cout << m_alphaPoints[i]*180./M_PI << std::endl;
550 }
551
552 ntheta_old = dbXT_old->getNoOfThetaBins();
553 B2INFO("Ntheta: " << ntheta_old);
554 for (unsigned short i = 0; i < ntheta_old; ++i) {
555 // m_thetaPoints[i] = (*dbXT_old).getThetaPoint(i);
556 array3 theta = dbXT_old->getThetaBin(i);
557 l_theta_old[i] = theta[0] * rad2deg;
558 u_theta_old[i] = theta[1] * rad2deg;
559 itheta_old[i] = theta[2] * rad2deg;
560
561
562 // std::cout << m_thetaPoints[i]*180./M_PI << std::endl;
563 }
564
565 xtmode_old = dbXT_old->getXtParamMode();
566
567 for (unsigned short iCL = 0; iCL < 56; ++iCL) {
568 for (unsigned short iLR = 0; iLR < 2; ++iLR) {
569 for (unsigned short iA = 0; iA < nalpha_old; ++iA) {
570 for (unsigned short iT = 0; iT < ntheta_old; ++iT) {
571 const std::vector<float> params = dbXT_old->getXtParams(iCL, iLR, iA, iT);
572 unsigned short np = params.size();
573 // std::cout <<"np4xt= " << np << std::endl;
574 for (unsigned short i = 0; i < np; ++i) {
575 xtold[iCL][iLR][iA][iT][i] = params[i];
576 }
577 }
578 }
579 }
580 }
581 // return true;
582}
Class for accessing objects in the database.
Definition: DBObjPtr.h:21

◆ readXTFromText()

void readXTFromText ( )
protectedvirtual

Read old xt parameter from text file, incase text mode is used.

Definition at line 438 of file XTCalibration.cc.

439{
440 std::string fileName1 = "/data/cdc" + m_xtfile;
441 std::string fileName = FileSystem::findFile(fileName1);
442 boost::iostreams::filtering_istream ifs;
443 if (fileName == "") {
444 fileName = FileSystem::findFile(m_xtfile);
445 }
446 if (fileName == "") {
447 B2FATAL("CDCGeometryPar: " << fileName1 << " not exist!");
448 } else {
449 B2INFO("CDCGeometryPar: open " << fileName1);
450 if ((fileName.rfind(".gz") != string::npos) && (fileName.length() - fileName.rfind(".gz") == 3)) {
451 ifs.push(boost::iostreams::gzip_decompressor());
452 }
453 ifs.push(boost::iostreams::file_source(fileName));
454 if (!ifs) B2FATAL("CDCGeometryPar: cannot open " << fileName1 << " !");
455
456
457 }
458 const int npar = 8;
459 //read alpha bin info.
460 // unsigned short nAlphaBins = 0;
461 ifs >> nalpha_old;
462
463 for (unsigned short i = 0; i < nalpha_old; ++i) {
464 ifs >> l_alpha_old[i] >> u_alpha_old[i] >> ialpha_old[i];
465 // ifs >> alpha0 >> alpha1 >> alpha2;
466 //m_alphaPoints[i] = alpha2;
467 }
468
469 //read theta bin info.
470 // unsigned short nThetaBins = 0;
471 ifs >> ntheta_old;
472
473 for (unsigned short i = 0; i < ntheta_old; ++i) {
474 ifs >> l_theta_old[i] >> u_theta_old[i] >> itheta_old[i];
475 //ifs >> theta0 >> theta1 >> theta2;
476 //m_thetaPoints[i] = theta2;
477 }
478
479 B2INFO("number of alpha - theta bin" << nalpha_old << " - " << ntheta_old);
480 short np = 0;
481 unsigned short iCL, iLR;
482 // const unsigned short npx = nXTParams - 1;
483 double xtc[npar]; // cppcheck-suppress constVariable
484 double theta, alpha, dummy1;
485 // unsigned m_xtParamMode_old;
486 ifs >> xtmode_old >> np;
487
488 const double epsi = 0.1;
489
490 while (ifs >> iCL) {
491 ifs >> theta >> alpha >> dummy1 >> iLR;
492 for (int i = 0; i < np; ++i) {
493 ifs >> xtc[i];
494 }
495
496 int ith = -99;
497 for (unsigned short i = 0; i < ntheta_old; ++i) {
498 if (fabs(theta - itheta_old[i]) < epsi) {
499 ith = i;
500 break;
501 }
502 }
503 if (ith < 0) {
504 gSystem->Exec("echo xt_theta error binning>> error");
505 return;
506 }
507
508 int ial = -99;
509 for (unsigned short i = 0; i < nalpha_old; ++i) {
510 if (fabs(alpha - ialpha_old[i]) < epsi) {
511 ial = i;
512 break;
513 }
514 }
515 if (ial < 0) {
516 gSystem->Exec("echo xt_alpha error binning>> error");
517 return;
518 }
519
520 for (int i = 0; i < np; ++i) {
521 xtold[iCL][iLR][ial][ith][i] = xtc[i];
522 }
523
524 } //end of while loop
525
526 //convert unit
527 /*
528 const double degrad = M_PI / 180.;
529 for (unsigned i = 0; i < nAlphaBins; ++i) {
530 m_alphaPoints[i] *= degrad;
531 }
532 for (unsigned i = 0; i < nThetaBins; ++i) {
533 m_thetaPoints[i] *= degrad;
534 }
535 */
536 // return true;
537}
std::string m_xtfile
Input xt file name, incase text mode.
static std::string findFile(const std::string &path, bool silent=false)
Search for given file or directory in local or central release directory, and return absolute path if...
Definition: FileSystem.cc:151

◆ setDebug()

virtual void setDebug ( bool  debug = false)
inlinevirtual

Run in debug or silent.

Definition at line 30 of file XTCalibration.h.

30{m_debug = debug; }

◆ setLRSeparate()

void setLRSeparate ( bool  lr = true)
inline

Set LR separate mode (default is true).

Definition at line 50 of file XTCalibration.h.

50{m_LRseparate = lr;}

◆ setMinimumNDF()

virtual void setMinimumNDF ( double  minndf)
inlinevirtual

set minimum number of degree of freedom requirement

Definition at line 34 of file XTCalibration.h.

34{m_ndfmin = minndf;}

◆ setMinimumPval()

virtual void setMinimumPval ( double  minPval)
inlinevirtual

set minimum Prob(Chi2) requirement

Definition at line 36 of file XTCalibration.h.

36{m_Pvalmin = minPval;}

◆ setMode()

virtual void setMode ( unsigned short  mode = 1)
inlinevirtual

set xt mode, 0 is polynimial, 1 is Chebshev polynomial

Definition at line 46 of file XTCalibration.h.

46{m_xtmode = mode;}

◆ setStoreHisto()

virtual void setStoreHisto ( bool  storeHist = false)
inlinevirtual

set to store histogram or not.

Definition at line 48 of file XTCalibration.h.

48{m_storeHisto = storeHist;}
bool m_storeHisto
Store histogram or not.
Definition: XTCalibration.h:81

◆ setUseDB()

virtual void setUseDB ( bool  useDB = false)
inlinevirtual

Set to run with database mode or text mode.

Definition at line 32 of file XTCalibration.h.

32{m_useDB = useDB; }

◆ setXTFileName()

virtual void setXTFileName ( std::string  name)
inlinevirtual

input xt file name incase text mode is used.

Definition at line 44 of file XTCalibration.h.

44{m_xtfile.assign(name);}

◆ storeHisto()

void storeHisto ( )
protectedvirtual

Store histogram to file.

Definition at line 375 of file XTCalibration.cc.

376{
377 B2INFO("start store histogram");
378 TFile* fout = new TFile("XTFIT.root", "RECREATE");
379 TDirectory* top = gDirectory;
380 TDirectory* Direct[56];
381 int nhisto = 0;
382 for (int l = 0; l < 56; ++l) {
383 top->cd();
384 Direct[l] = gDirectory->mkdir(Form("lay_%d", l));
385 Direct[l]->cd();
386 for (int th = 0; th < m_ntheta; ++th) {
387 for (int al = 0; al < m_nalpha; ++al) {
388 hist2d_draw[l][al][th]->Write();
389 for (int lr = 0; lr < 2; ++lr) {
390 hist2d[l][lr][al][th]->Write();
391 if (fitflag[l][lr][al][th] != 1) continue;
392 if (m_useSliceFit) {
393 if (hist2d_1[l][lr][al][th]) {
394 hist2d_1[l][lr][al][th]->Write();
395 nhisto += 1;
396 }
397 } else {
398 hprof[l][lr][al][th]->Write();
399 nhisto += 1;
400 }
401 }
402 }
403 }
404 }
405 fout->Close();
406 B2RESULT(" " << nhisto << " histograms was stored.");
407}

◆ useProfileFromInputXT()

virtual void useProfileFromInputXT ( bool  useProfileXTFromInputXT)
inlinevirtual

if you want to change xt bining, you have to set this to true

Definition at line 42 of file XTCalibration.h.

42{m_useProfileXTFromInputXT = useProfileXTFromInputXT;}

◆ Write()

void Write ( )
protectedvirtual

Store calibrated constand.

Definition at line 297 of file XTCalibration.cc.

298{
299 /*Set parameter for layer that fit do not success*/
300 /* and then write output file*/
301 double par[8];
302 ofstream xtout(m_OutputXTFileName.c_str());
303 xtout << m_nalpha << endl;
304 for (int i = 0; i < m_nalpha; ++i) {
305 xtout << std::setprecision(3) << l_alpha[i] << " "
306 << std::setprecision(3) << u_alpha[i] << " "
307 << std::setprecision(3) << ialpha[i] << endl;
308 }
309 xtout << m_ntheta << endl;
310 for (int i = 0; i < m_ntheta; ++i) {
311 xtout << std::setprecision(3) << l_theta[i] << " "
312 << std::setprecision(3) << u_theta[i] << " "
313 << std::setprecision(3) << itheta[i] << endl;
314 }
315 xtout << m_xtmode << " " << 8 << endl;
316
317 int nfitted = 0;
318 int nfailure = 0;
319
320 for (int th = 0; th < m_ntheta; ++th) {
321 for (int al = 0; al < m_nalpha; ++al) {
322 for (int l = 0; l < 56; ++l) {
323 for (int lr = 0; lr < 2; ++lr) {
324 /*Set Parameter for bad fit*/
325 if (fitflag[l][lr][al][th] != 1) {
326 nfailure += 1;
327 printf("fit failure status = %d \n", fitflag[l][lr][al][th]);
328 printf("layer %d, r %d, alpha %3.1f, theta %3.1f \n", l, lr, ialpha[al], itheta[th]);
329 printf("number of event: %3.2f", hprof[l][lr][al][th]->GetEntries());
330 if (fitflag[l][lr][al][th] != -1) {
331 printf("Probability of fit: %3.4f", xtf5r[l][lr][al][th]->GetProb());
332 }
333 if (m_xtmode == xtmode_old) {
334 int ial_old = 0;
335 int ith_old = 0;
336 for (int k = 0; k < nalpha_old; ++k) {
337 if (ialpha[al] < u_alpha_old[k]) {ial_old = k; break;}
338 }
339 for (int j = 0; j < ntheta_old; ++j) {
340 if (itheta[th] < u_theta_old[j]) {ith_old = j; break;}
341 }
342 for (int p = 0; p < 8; ++p) {
343 par[p] = xtold[l][lr][ial_old][ith_old][p];
344 }
345 } else {
346 //if mode of input xt is different from output, simple xt is used.
347 par[0] = 0; par[1] = 0.004; par[2] = 0; par[3] = 0; par[4] = 0; par[5] = 0; par[6] = m_par6[l]; par[7] = 0.00001;
348 }
349 } else {
350 xtf5r[l][lr][al][th]->GetParameters(par);
351 nfitted += 1;
352 }
353 /*Write params*/
354 xtout << l << std::setw(5) << itheta[th] << std::setw(5) << ialpha[al] << std::setw(5) << "0.0" << std::setw(4) << lr << std::setw(
355 15);
356 for (int p = 0; p < 8; ++p) {
357 if (p != 7) { xtout << std::setprecision(7) << par[p] << std::setw(15);}
358 if (p == 7) { xtout << std::setprecision(7) << par[p] << std::endl;}
359 }
360 }//lr
361 }//th
362 }//al
363 }//lay
364 xtout.close();
365 B2RESULT(" Total number of xt fit: " << m_nalpha * m_ntheta * 2 * 56);
366 B2RESULT(" Successfully Fitted: " << nfitted);
367 B2RESULT(" Failure Fit: " << nfailure);
368 B2RESULT("Finish export xt to text file");
369 if (m_useDB) {
370 CDCDatabaseImporter import(0, 0, -1, -1);
371 import.importXT(m_OutputXTFileName.c_str());
372 }
373}
CDC database importer.
void importXT(std::string fileName)
Import xt table to the database.
std::string m_OutputXTFileName
Out put xt filename.

Member Data Documentation

◆ fitflag

int fitflag[56][2][20][10]
private

Fit flag.

Definition at line 90 of file XTCalibration.h.

◆ hist2d

TH2D* hist2d[56][2][20][10]
private

2D histo of xt

Definition at line 94 of file XTCalibration.h.

◆ hist2d_1

TH1D* hist2d_1[56][2][20][10]
private

1D xt histo, results of slice fit

Definition at line 96 of file XTCalibration.h.

◆ hist2d_draw

TH2D* hist2d_draw[56][20][10]
private

2d histo for draw

Definition at line 95 of file XTCalibration.h.

◆ hprof

TProfile* hprof[56][2][20][10]
private

Profile xt histo.

Definition at line 93 of file XTCalibration.h.

◆ ialpha

double ialpha[18]
private

represented alphas of alpha bins.

Definition at line 124 of file XTCalibration.h.

◆ ialpha_old

double ialpha_old[18]
private

represented alphas of alpha bins from input.

Definition at line 133 of file XTCalibration.h.

◆ itheta

double itheta[7]
private

represented alphas of theta bins.

Definition at line 127 of file XTCalibration.h.

◆ itheta_old

double itheta_old[7]
private

represented alphas of theta bins from input.

Definition at line 136 of file XTCalibration.h.

◆ l_alpha

double l_alpha[18]
private

Lower boundays of alpha bins.

Definition at line 122 of file XTCalibration.h.

◆ l_alpha_old

double l_alpha_old[18]
private

Lower boundays of alpha bins from input.

Definition at line 131 of file XTCalibration.h.

◆ l_theta

double l_theta[7]
private

Lower boundays of theta bins.

Definition at line 125 of file XTCalibration.h.

◆ l_theta_old

double l_theta_old[7]
private

Lower boundays of theta bins from input.

Definition at line 134 of file XTCalibration.h.

◆ m_BField

bool m_BField = true
private

with b field or none

Definition at line 86 of file XTCalibration.h.

◆ m_debug

bool m_debug = false
private

run in debug or silent

Definition at line 80 of file XTCalibration.h.

◆ m_inputRootFileNames

std::string m_inputRootFileNames = "rootfile/output*"
private

input root filename

Definition at line 107 of file XTCalibration.h.

◆ m_LRseparate

bool m_LRseparate = true
private

Separate LR in calibration or mix.

Definition at line 84 of file XTCalibration.h.

◆ m_MAXalpha

int m_MAXalpha = 18
private

max alpha bin

Definition at line 139 of file XTCalibration.h.

◆ m_MAXtheta

int m_MAXtheta = 7
private

max theta bin

Definition at line 140 of file XTCalibration.h.

◆ m_nalpha

int m_nalpha
private

number of alpha bins

Definition at line 120 of file XTCalibration.h.

◆ m_ndfmin

double m_ndfmin = 5
private

minimum ndf required

Definition at line 78 of file XTCalibration.h.

◆ m_ntheta

int m_ntheta
private

number of theta bins

Definition at line 121 of file XTCalibration.h.

◆ m_OutputXTFileName

std::string m_OutputXTFileName = "xt_new.dat"
private

Out put xt filename.

Definition at line 106 of file XTCalibration.h.

◆ m_par6

double m_par6[56]
private
Initial value:
= {89, 91, 94, 99, 104, 107, 110, 117,
126, 144, 150, 157, 170, 180,
160, 167, 183, 205, 200, 194,
177, 189, 192, 206, 224, 234,
193, 206, 209, 215, 222, 239,
204, 212, 217, 227, 235, 240,
215, 222, 230, 239, 246, 253,
227, 232, 239, 243, 253, 258,
231, 243, 246, 256, 263, 300
}

boundary parameter for fitting, semi-experiment number

Definition at line 144 of file XTCalibration.h.

◆ m_profileFileName

std::string m_profileFileName = "xt_profile"
private

profile file name

Definition at line 108 of file XTCalibration.h.

◆ m_Pvalmin

double m_Pvalmin = 0.
private

minimum pvalue required

Definition at line 79 of file XTCalibration.h.

◆ m_smallestEntryRequire

int m_smallestEntryRequire = 1000
private

minimum number of hit per hitosgram.

Definition at line 142 of file XTCalibration.h.

◆ m_storeHisto

bool m_storeHisto = false
private

Store histogram or not.

Definition at line 81 of file XTCalibration.h.

◆ m_useDB

bool m_useDB = false
private

Use Database or text mode.

Definition at line 82 of file XTCalibration.h.

◆ m_useProfileXTFromInputXT

bool m_useProfileXTFromInputXT = true
private

use profile from text file or default in input xt

Definition at line 83 of file XTCalibration.h.

◆ m_useSliceFit

bool m_useSliceFit = false
private

Use slice fit or profile.

Definition at line 85 of file XTCalibration.h.

◆ m_xtfile

std::string m_xtfile = "cdc/data/xt.dat"
private

Input xt file name, incase text mode.

Definition at line 109 of file XTCalibration.h.

◆ m_xtmode

unsigned short m_xtmode = 1
private

Mode of xt; 0 is polynomial;1 is Chebyshev.

Definition at line 141 of file XTCalibration.h.

◆ nalpha_old

int nalpha_old
private

number of alpha bins from input

Definition at line 129 of file XTCalibration.h.

◆ ntheta_old

int ntheta_old
private

number of theta bins from input

Definition at line 130 of file XTCalibration.h.

◆ u_alpha

double u_alpha[18]
private

Upper boundays of alpha bins.

Definition at line 123 of file XTCalibration.h.

◆ u_alpha_old

double u_alpha_old[18]
private

Upper boundays of alpha bins from input.

Definition at line 132 of file XTCalibration.h.

◆ u_theta

double u_theta[7]
private

Upper boundays of theta bins.

Definition at line 126 of file XTCalibration.h.

◆ u_theta_old

double u_theta_old[7]
private

Upper boundays of theta bins from input.

Definition at line 135 of file XTCalibration.h.

◆ xtf5r

TF1* xtf5r[56][2][20][10]
private

XTFunction.

Definition at line 91 of file XTCalibration.h.

◆ xtmode_old

unsigned short xtmode_old
private

XT mode old, 0-polynomial, 1 Cheb.

Definition at line 138 of file XTCalibration.h.

◆ xtold

double xtold[56][2][18][7][8]
private

Old paremeter.

Definition at line 89 of file XTCalibration.h.


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