BOSS 7.0.8
BESIII Offline Software System
Loading...
Searching...
No Matches
QtMdcCalib.cxx
Go to the documentation of this file.
2
3#include "GaudiKernel/MsgStream.h"
4#include "GaudiKernel/IMessageSvc.h"
5#include "GaudiKernel/StatusCode.h"
6#include "GaudiKernel/ISvcLocator.h"
7#include "GaudiKernel/Bootstrap.h"
8
9#include <iostream>
10#include <iomanip>
11#include <fstream>
12#include <string>
13#include <cstring>
14#include <math.h>
15
16#include "TF1.h"
17#include "TMinuit.h"
18
19using namespace std;
20
22 m_vdr = 0.03;
23
24 m_nlayer = MdcCalNLayer;
25 m_qtorder = MdcCalQtOrd;
26 m_nbin = MdcCalNQBin;
27 m_innNLay = MdcCalInnNLay;
28}
29
31}
32
34 int bin;
35 for(int lay=0; lay<MdcCalNLayer; lay++){
36 delete m_hqhit[lay];
37 delete m_grqt[lay];
38 delete m_grqdt[lay];
39 for(bin=0; bin<MdcCalQtOrd; bin++){
40 delete m_hqt[lay][bin];
41 }
42 }
43 delete m_fdQt;
44 delete m_fdQ_T;
45
47}
48
49void QtMdcCalib::initialize(TObjArray* hlist, IMdcGeomSvc* mdcGeomSvc,
50 IMdcCalibFunSvc* mdcFunSvc, IMdcUtilitySvc* mdcUtilitySvc) {
51 IMessageSvc* msgSvc;
52 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
53 MsgStream log(msgSvc, "QtMdcCalib");
54 log << MSG::INFO << "QtMdcCalib::initialize()" << endreq;
55
56 m_hlist = hlist;
57 m_mdcGeomSvc = mdcGeomSvc;
58 m_mdcFunSvc = mdcFunSvc;
59 m_mdcUtilitySvc = mdcUtilitySvc;
60
61 MdcCalib::initialize(m_hlist, m_mdcGeomSvc, m_mdcFunSvc, m_mdcUtilitySvc);
62
63 int bin;
64 int lay;
65 double qbinw;
66 char hname[200];
67
68 for(lay=0; lay<m_nlayer; lay++){
69 m_qmin[lay] = m_param.qmin[lay];
70 m_qmax[lay] = m_param.qmax[lay];
71 m_qbinw[lay] = (m_qmax[lay] - m_qmin[lay]) / (double)m_nbin;
72 }
73
74 m_fdQt = new TFolder("fdQt", "fdQt");
75 m_fdQ_T = new TFolder("QtPlot", "QtPlot");
76 m_hlist -> Add(m_fdQt);
77 m_hlist -> Add(m_fdQ_T);
78
79 for(lay=0; lay<m_nlayer; lay++){
80 sprintf(hname, "HQ_Layer%02d", lay);
81 m_hqhit[lay] = new TH1F(hname, "", 1500, 0, 3000);
82 m_fdQt -> Add(m_hqhit[lay]);
83
84 sprintf(hname, "HQT_Plot_lay%02d", lay);
85 m_grqt[lay] = new TGraphErrors();
86 m_grqt[lay]->SetName(hname);
87 m_grqt[lay]->SetMarkerStyle(20);
88 m_grqt[lay]->SetMarkerColor(1);
89 m_fdQ_T->Add(m_grqt[lay]);
90
91 sprintf(hname, "HQdelT_Plot_lay%02d", lay);
92 m_grqdt[lay] = new TGraphErrors();
93 m_grqdt[lay]->SetName(hname);
94 m_grqdt[lay]->SetMarkerStyle(10);
95 m_grqdt[lay]->SetMarkerColor(1);
96 m_fdQ_T->Add(m_grqdt[lay]);
97
98 for(bin=0; bin<m_nbin; bin++){
99 sprintf(hname, "HQT_Lay%02d_Bin%02d", lay, bin);
100 m_hqt[lay][bin] = new TH1F(hname, "", 200, -1, 1);
101 m_fdQt -> Add(m_hqt[lay][bin]);
102 }
103 }
104}
105
107 IMessageSvc* msgSvc;
108 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
109 MsgStream log(msgSvc, "QtMdcCalib");
110 log << MSG::DEBUG << "QtMdcCalib::fillHist()" << endreq;
111
112 MdcCalib::fillHist(event);
113
114 // get EsTimeCol
115 bool esCutFg = event->getEsCutFlag();
116 if( ! esCutFg ) return -1;
117
118 int i;
119 int k;
120 int bin;
121 int lay;
122
123 double doca;
124 double dmeas;
125
126 int ntrk;
127 int nhit;
128
129 bool fgHitLay[MdcCalNLayer];
130
131 MdcCalRecTrk* rectrk;
132 MdcCalRecHit* rechit;
133
134 ntrk = event -> getNTrk();
135 for(i=0; i<ntrk; i++){
136 rectrk = event -> getRecTrk(i);
137 nhit = rectrk -> getNHits();
138
139 // dr cut
140 double dr = rectrk->getDr();
141 if(fabs(dr) > m_param.drCut) continue;
142
143 // dz cut
144 double dz = rectrk->getDz();
145 if(fabs(dz) > m_param.dzCut) continue;
146
147 for(lay=0; lay<MdcCalNLayer; lay++) fgHitLay[lay] = false;
148 for(k=0; k<nhit; k++){
149 rechit = rectrk -> getRecHit(k);
150 lay = rechit -> getLayid();
151 fgHitLay[lay] = true;
152 }
153
154 int nhitlay = 0;
155 for(lay=0; lay<MdcCalNLayer; lay++) if(fgHitLay[lay]) nhitlay++;
156 if(nhitlay < m_param.nHitLayCut) continue;
157
158 bool fgNoise = rectrk->getFgNoiseRatio();
159 if(m_param.noiseCut && (!fgNoise)) continue;
160
161 for(k=0; k<nhit; k++){
162 rechit = rectrk -> getRecHit(k);
163 lay = rechit -> getLayid();
164 doca = rechit -> getDocaInc();
165 dmeas = rechit -> getDmeas();
166 m_resi = rechit -> getResiInc();
167 m_qhit = rechit -> getQhit();
168
169 m_hqhit[lay] -> Fill(m_qhit);
170
171 bin = (int)((m_qhit - m_qmin[lay]) / m_qbinw[lay]);
172 if( (bin >= 0) && (bin < m_nbin) ){
173 m_hqt[lay][bin] -> Fill( m_resi );
174 }
175 }
176 }
177 return 1;
178}
179
181 IMessageSvc* msgSvc;
182 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
183 MsgStream log(msgSvc, "QtMdcCalib");
184 log << MSG::INFO << "QtMdcCalib::updateConst()" << endreq;
185
186 MdcCalib::updateConst(calconst);
187
188 int lay;
189 int bin;
190 int ord;
191
192 Stat_t entry;
193 double qtpar;
194 double qbcen;
195 double tw;
196 double deltw;
197 double qterr;
198
199 TF1* funQt = new TF1("funQt", qtFun, 200, 2000, 2);
200
201 ofstream fqtlog("qtlog");
202 for(lay=0; lay<m_nlayer; lay++){
203 if(0 == m_param.fgCalib[lay]) continue;
204
205 fqtlog << "Layer" << lay << endl;
206
207 for(ord=0; ord<m_qtorder; ord++){
208 m_qtpar[lay][ord] = calconst -> getQtpar(lay, ord);
209 }
210
211 for(bin=0; bin<m_nbin; bin++){
212 entry = m_hqt[lay][bin] -> GetEntries();
213
214 if(entry > 300){
215 deltw = m_hqt[lay][bin] -> GetMean();
216 qterr = ( m_hqt[lay][bin]->GetRMS() ) / sqrt((double)entry);
217 deltw /= m_vdr;
218 qterr /= m_vdr;
219 } else{
220 continue;
221 }
222
223 qbcen = ( (double)bin + 0.5 ) * m_qbinw[lay] + m_qmin[lay];
224// tw = qtFun(qbcen, m_qtpar[lay]) + deltw;
225 tw = (m_mdcFunSvc->getTimeWalk(lay, qbcen)) + deltw;
226
227 m_grqt[lay]->SetPoint(bin, qbcen, tw);
228 m_grqt[lay]->SetPointError(bin, 0, qterr);
229
230 m_grqdt[lay]->SetPoint(bin, qbcen, deltw);
231 m_grqdt[lay]->SetPointError(bin, 0, qterr);
232
233 fqtlog << setw(3) << bin
234 << setw(12) << deltw
235 << setw(12) << tw
236 << setw(12) << qbcen
237 << setw(12) << qterr
238 << endl;
239 }
240
241 m_grqt[lay]->Fit("funQt", "Q+", "", m_qmin[lay], m_qmax[lay]);
242
243 fqtlog << "Qtpar: ";
244 for(ord=0; ord<m_qtorder; ord++){
245 qtpar = funQt->GetParameter(ord);
246 qterr = funQt->GetParError(ord);
247 calconst -> resetQtpar(lay, ord, qtpar);
248
249 fqtlog << setw(12) << qtpar
250 << setw(12) << qterr
251 << endl;
252 }
253
254// if( (0 == ierflg) && (3 == istat) ){
255// for(ord=0; ord<m_qtorder; ord++){
256// gmqt -> GetParameter(ord, qtpar, qterr);
257// calconst -> resetQtpar(lay, ord, qtpar);
258
259// fqtlog << setw(12) << qtpar
260// << setw(12) << qterr
261// << endl;
262// }
263// } else{
264// fqtlog << setw(12) << m_qtpar[lay][0]
265// << setw(12) << "0"
266// << setw(12) << m_qtpar[lay][1]
267// << setw(12) << "0"
268// << endl;
269// }
270
271 } // end of layer loop
272
273 fqtlog.close();
274 delete funQt;
275 return 1;
276}
277
278Double_t QtMdcCalib::qtFun(Double_t *x, Double_t *par){
279 double tw = par[1] / sqrt(x[0]) + par[0];
280 return tw;
281}
282
curve Fill()
sprintf(cut,"kal_costheta0_em>-0.93&&kal_costheta0_em<0.93&&kal_pxy0_em>=0.05+%d*0.1&&kal_pxy0_em<0.15+%d*0.1&&NGch>=2", j, j)
mg Add(gr3)
Double_t x[10]
*******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
const int MdcCalNLayer
Definition: MdcCalParams.h:6
const int MdcCalQtOrd
Definition: MdcCalParams.h:16
const int MdcCalNQBin
Definition: MdcCalParams.h:17
const int MdcCalInnNLay
Definition: MdcCalParams.h:7
IMessageSvc * msgSvc()
virtual double getTimeWalk(int layid, double Q) const =0
int fgCalib[MdcCalNLayer]
Definition: MdcCalParams.h:75
double qmax[MdcCalNLayer]
Definition: MdcCalParams.h:81
double dzCut
Definition: MdcCalParams.h:71
double qmin[MdcCalNLayer]
Definition: MdcCalParams.h:80
double drCut
Definition: MdcCalParams.h:70
double getDz() const
Definition: MdcCalRecTrk.h:33
bool getFgNoiseRatio() const
Definition: MdcCalRecTrk.h:46
double getDr() const
Definition: MdcCalRecTrk.h:30
virtual void clear()=0
Definition: MdcCalib.cxx:76
virtual void initialize(TObjArray *hlist, IMdcGeomSvc *mdcGeomSvc, IMdcCalibFunSvc *mdcFunSvc, IMdcUtilitySvc *mdcUtilitySvc)=0
Definition: MdcCalib.cxx:202
virtual int updateConst(MdcCalibConst *calconst)=0
Definition: MdcCalib.cxx:1322
virtual int fillHist(MdcCalEvent *event)=0
Definition: MdcCalib.cxx:692
static Double_t qtFun(Double_t *x, Double_t *par)
Definition: QtMdcCalib.cxx:278
void clear()
Definition: QtMdcCalib.cxx:33
void initialize(TObjArray *hlist, IMdcGeomSvc *mdcGeomSvc, IMdcCalibFunSvc *mdcFunSvc, IMdcUtilitySvc *mdcUtilitySvc)
Definition: QtMdcCalib.cxx:49
int fillHist(MdcCalEvent *event)
Definition: QtMdcCalib.cxx:106
int updateConst(MdcCalibConst *calconst)
Definition: QtMdcCalib.cxx:180