CGEM BOSS 6.6.5.h
BESIII Offline Software System
Loading...
Searching...
No Matches
QtCalib Class Reference

#include <QtCalib.h>

+ Inheritance diagram for QtCalib:

Public Member Functions

 QtCalib ()
 
 ~QtCalib ()
 
void init (TObjArray *hlist, MdcCosGeom *pGeom)
 
void mergeHist (TFile *fhist)
 
void calib (MdcCalibConst *calconst, TObjArray *newXtList, TObjArray *r2tList)
 
- Public Member Functions inherited from CalibBase
 CalibBase ()
 
virtual ~CalibBase ()
 

Static Public Member Functions

static Double_t qtFun (Double_t *x, Double_t *par)
 

Detailed Description

Definition at line 12 of file QtCalib.h.

Constructor & Destructor Documentation

◆ QtCalib()

QtCalib::QtCalib ( )

Definition at line 6 of file QtCalib.cpp.

6 {
7 cout << "Calibration type: QtCalib" << endl;
8}

◆ ~QtCalib()

QtCalib::~QtCalib ( )

Definition at line 10 of file QtCalib.cpp.

10 {
11}

Member Function Documentation

◆ calib()

void QtCalib::calib ( MdcCalibConst * calconst,
TObjArray * newXtList,
TObjArray * r2tList )
virtual

Implements CalibBase.

Definition at line 77 of file QtCalib.cpp.

77 {
78 CalibBase::calib(calconst, newXtList, r2tList);
79
80 double vdr = 0.03;
81 Stat_t entry;
82 double qtpar;
83 double qbcen;
84 double tw;
85 double deltw;
86 double qterr;
87 TF1* funQt = new TF1("funQt", qtFun, 200, 2000, 2);
88
89 ofstream fqtlog("qtlog");
90 for(int lay=0; lay<NLAYER; lay++){
91 if(0 == gFgCalib[lay]) continue;
92
93 fqtlog << "Layer" << lay << endl;
94 double qtini[2];
95 for(int ord=0; ord<QtOrd; ord++) qtini[ord] = calconst->getQtpar(lay, ord);
96 for(int bin=0; bin<NQBin; bin++){
97 entry = m_hqt[lay][bin]->GetEntries();
98 if(entry > 300){
99 deltw = m_hqt[lay][bin] -> GetMean();
100 qterr = ( m_hqt[lay][bin]->GetRMS() ) / sqrt((double)entry);
101 deltw /= vdr;
102 qterr /= vdr;
103 } else{
104 continue;
105 }
106
107 qbcen = ( (double)bin + 0.5 ) * m_qbinw[lay] + m_qmin[lay];
108// tw = qtFun(qbcen, m_qtpar[lay]) + deltw;
109// tw = (m_mdcFunSvc->getTimeWalk(lay, qbcen)) + deltw;
110 tw = qtini[1] / sqrt(qbcen) + qtini[0] + deltw;
111
112 m_grqt[lay]->SetPoint(bin, qbcen, tw);
113 m_grqt[lay]->SetPointError(bin, 0, qterr);
114
115 m_grqdt[lay]->SetPoint(bin, qbcen, deltw);
116 m_grqdt[lay]->SetPointError(bin, 0, qterr);
117
118 fqtlog << setw(3) << bin << setw(12) << deltw << setw(12) << tw
119 << setw(12) << qbcen << setw(12) << qterr << endl;
120 }
121
122 m_grqt[lay]->Fit("funQt", "Q+", "", m_qmin[lay], m_qmax[lay]);
123
124 fqtlog << "Qtpar: ";
125 for(int ord=0; ord<QtOrd; ord++){
126 qtpar = funQt->GetParameter(ord);
127 qterr = funQt->GetParError(ord);
128 calconst -> resetQtpar(lay, ord, qtpar);
129
130 fqtlog << setw(12) << qtpar << setw(12) << qterr << endl;
131 }
132 } // end of layer loop
133 fqtlog.close();
134 renameHist();
135 delete funQt;
136}
*******INTEGER m_nBinMax INTEGER m_NdiMax !No of bins in histogram for cell exploration division $ !Last vertex $ !Last active cell $ !Last cell in buffer $ !No of sampling when dividing cell $ !No of function total $ !Flag for random ceel for $ !Flag for type of for WtMax $ !Flag which decides whether vertices are included in the sampling $ entire domain is hyp !Maximum effective eevents per bin
Definition FoamA.h:85
virtual void calib(MdcCalibConst *calconst, TObjArray *newXtList, TObjArray *r2tList)=0
double getQtpar(int lay, int order) const
static Double_t qtFun(Double_t *x, Double_t *par)
Definition QtCalib.cpp:138

◆ init()

void QtCalib::init ( TObjArray * hlist,
MdcCosGeom * pGeom )
virtual

Implements CalibBase.

Definition at line 13 of file QtCalib.cpp.

13 {
14 CalibBase::init(hlist, pGeom);
15 m_pGeom = pGeom;
16
17 char hname[200];
18 for(int lay=0; lay<NLAYER; lay++){
19 m_qmin[lay] = gQmin[lay];
20 m_qmax[lay] = gQmax[lay];
21 m_qbinw[lay] = (m_qmax[lay] - m_qmin[lay]) / (double)NQBin;
22 }
23
24 m_fdQt = new TFolder("mfdQt", "fdQt");
25 m_fdQ_T = new TFolder("mQtPlot", "QtPlot");
26 hlist -> Add(m_fdQt);
27 hlist -> Add(m_fdQ_T);
28
29 for(int lay=0; lay<NLAYER; lay++){
30 sprintf(hname, "mHQ_Layer%02d", lay);
31 m_hqhit[lay] = new TH1F(hname, "", 1500, 0, 3000);
32 m_fdQt -> Add(m_hqhit[lay]);
33
34 sprintf(hname, "mHQT_Plot_lay%02d", lay);
35 m_grqt[lay] = new TGraphErrors();
36 m_grqt[lay]->SetName(hname);
37 m_grqt[lay]->SetMarkerStyle(20);
38 m_grqt[lay]->SetMarkerColor(1);
39 m_grqt[lay]->SetLineColor(10);
40 m_fdQ_T->Add(m_grqt[lay]);
41
42 sprintf(hname, "mHQdelT_Plot_lay%02d", lay);
43 m_grqdt[lay] = new TGraphErrors();
44 m_grqdt[lay]->SetName(hname);
45 m_grqdt[lay]->SetMarkerStyle(10);
46 m_grqdt[lay]->SetMarkerColor(1);
47 m_grqdt[lay]->SetLineColor(10);
48 m_fdQ_T->Add(m_grqdt[lay]);
49
50 for(int bin=0; bin<NQBin; bin++){
51 sprintf(hname, "mHQT_Lay%02d_Bin%02d", lay, bin);
52 m_hqt[lay][bin] = new TH1F(hname, "", 200, -1, 1);
53 m_fdQt -> Add(m_hqt[lay][bin]);
54 }
55 }
56}
virtual void init(TObjArray *hlist, MdcCosGeom *pGeom)=0
Definition CalibBase.cpp:14

◆ mergeHist()

void QtCalib::mergeHist ( TFile * fhist)
virtual

Implements CalibBase.

Definition at line 58 of file QtCalib.cpp.

58 {
60
61 char hname[200];
62 TH1F* hist;
63 TFolder* fdQt = (TFolder*)fhist->Get("fdQt");
64 for(int lay=0; lay<NLAYER; lay++){
65 sprintf(hname, "HQ_Layer%02d", lay);
66 hist = (TH1F*)fdQt->FindObjectAny(hname);
67 m_hqhit[lay]->Add(hist);
68
69 for(int bin=0; bin<NQBin; bin++){
70 sprintf(hname, "HQT_Lay%02d_Bin%02d", lay, bin);
71 hist = (TH1F*)fdQt->FindObjectAny(hname);
72 m_hqt[lay][bin]->Add(hist);
73 }
74 }
75}
virtual void mergeHist(TFile *fhist)=0

◆ qtFun()

Double_t QtCalib::qtFun ( Double_t * x,
Double_t * par )
static

Definition at line 138 of file QtCalib.cpp.

138 {
139 Double_t tw = par[1] / sqrt(x[0]) + par[0];
140 return tw;
141}

Referenced by calib().


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