BOSS 7.1.2
BESIII Offline Software System
Loading...
Searching...
No Matches
QtMdcCalib Class Reference

#include <QtMdcCalib.h>

+ Inheritance diagram for QtMdcCalib:

Public Member Functions

 QtMdcCalib ()
 
 ~QtMdcCalib ()
 
void initialize (TObjArray *hlist, IMdcGeomSvc *mdcGeomSvc, IMdcCalibFunSvc *mdcFunSvc, IMdcUtilitySvc *mdcUtilitySvc)
 
void setParam (MdcCalParams &param)
 
int fillHist (MdcCalEvent *event)
 
int updateConst (MdcCalibConst *calconst)
 
void printCut () const
 
void clear ()
 
- Public Member Functions inherited from MdcCalib
 MdcCalib ()
 
virtual ~MdcCalib ()
 

Static Public Member Functions

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

Detailed Description

Definition at line 8 of file QtMdcCalib.h.

Constructor & Destructor Documentation

◆ QtMdcCalib()

QtMdcCalib::QtMdcCalib ( )

Definition at line 21 of file QtMdcCalib.cxx.

21 {
22 m_vdr = 0.03;
23
24 m_nlayer = MdcCalNLayer;
25 m_qtorder = MdcCalQtOrd;
26 m_nbin = MdcCalNQBin;
27 m_innNLay = MdcCalInnNLay;
28}
const int MdcCalNLayer
Definition MdcCalParams.h:6
const int MdcCalQtOrd
const int MdcCalNQBin
const int MdcCalInnNLay
Definition MdcCalParams.h:7

◆ ~QtMdcCalib()

QtMdcCalib::~QtMdcCalib ( )

Definition at line 30 of file QtMdcCalib.cxx.

30 {
31}

Member Function Documentation

◆ clear()

void QtMdcCalib::clear ( )
virtual

Implements MdcCalib.

Definition at line 33 of file QtMdcCalib.cxx.

33 {
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}
*******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 clear()=0
Definition MdcCalib.cxx:80

◆ fillHist()

int QtMdcCalib::fillHist ( MdcCalEvent * event)
virtual

Implements MdcCalib.

Definition at line 106 of file QtMdcCalib.cxx.

106 {
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 // momentum cut
148 double p = rectrk -> getP();
149 if((fabs(p) < m_param.pCut[0]) || (fabs(p) > m_param.pCut[1])) continue;
150
151 for(lay=0; lay<MdcCalNLayer; lay++) fgHitLay[lay] = false;
152 for(k=0; k<nhit; k++){
153 rechit = rectrk -> getRecHit(k);
154 lay = rechit -> getLayid();
155 fgHitLay[lay] = true;
156 }
157
158 int nhitlay = 0;
159 for(lay=0; lay<MdcCalNLayer; lay++) if(fgHitLay[lay]) nhitlay++;
160 if(nhitlay < m_param.nHitLayCut) continue;
161
162 bool fgNoise = rectrk->getFgNoiseRatio();
163 if(m_param.noiseCut && (!fgNoise)) continue;
164
165 for(k=0; k<nhit; k++){
166 rechit = rectrk -> getRecHit(k);
167 lay = rechit -> getLayid();
168 doca = rechit -> getDocaInc();
169 dmeas = rechit -> getDmeas();
170 m_resi = rechit -> getResiInc();
171 m_qhit = rechit -> getQhit();
172
173 m_hqhit[lay] -> Fill(m_qhit);
174
175 bin = (int)((m_qhit - m_qmin[lay]) / m_qbinw[lay]);
176 if( (bin >= 0) && (bin < m_nbin) ){
177 m_hqt[lay][bin] -> Fill( m_resi );
178 }
179 }
180 }
181 return 1;
182}
curve Fill()
IMessageSvc * msgSvc()
double pCut[2]
double getDz() const
bool getFgNoiseRatio() const
double getDr() const
virtual int fillHist(MdcCalEvent *event)=0
Definition MdcCalib.cxx:730

◆ initialize()

void QtMdcCalib::initialize ( TObjArray * hlist,
IMdcGeomSvc * mdcGeomSvc,
IMdcCalibFunSvc * mdcFunSvc,
IMdcUtilitySvc * mdcUtilitySvc )
virtual

Implements MdcCalib.

Definition at line 49 of file QtMdcCalib.cxx.

50 {
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}
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 qmax[MdcCalNLayer]
double qmin[MdcCalNLayer]
virtual void initialize(TObjArray *hlist, IMdcGeomSvc *mdcGeomSvc, IMdcCalibFunSvc *mdcFunSvc, IMdcUtilitySvc *mdcUtilitySvc)=0
Definition MdcCalib.cxx:214

◆ printCut()

void QtMdcCalib::printCut ( ) const
virtual

Implements MdcCalib.

Definition at line 184 of file QtMdcCalib.cxx.

184 {
186}
virtual void printCut() const =0

◆ qtFun()

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

Definition at line 286 of file QtMdcCalib.cxx.

286 {
287 double tw = par[1] / sqrt(x[0]) + par[0];
288 return tw;
289}

Referenced by updateConst().

◆ setParam()

void QtMdcCalib::setParam ( MdcCalParams & param)
inlinevirtual

Implements MdcCalib.

Definition at line 53 of file QtMdcCalib.h.

53 {
54 MdcCalib::setParam(param);
55 m_param = param;
56}
virtual void setParam(MdcCalParams &param)=0
Definition MdcCalib.h:306

◆ updateConst()

int QtMdcCalib::updateConst ( MdcCalibConst * calconst)
virtual

Implements MdcCalib.

Definition at line 188 of file QtMdcCalib.cxx.

188 {
189 IMessageSvc* msgSvc;
190 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
191 MsgStream log(msgSvc, "QtMdcCalib");
192 log << MSG::INFO << "QtMdcCalib::updateConst()" << endreq;
193
194 MdcCalib::updateConst(calconst);
195
196 int lay;
197 int bin;
198 int ord;
199
200 Stat_t entry;
201 double qtpar;
202 double qbcen;
203 double tw;
204 double deltw;
205 double qterr;
206
207 TF1* funQt = new TF1("funQt", qtFun, 200, 2000, 2);
208
209 ofstream fqtlog("qtlog");
210 for(lay=0; lay<m_nlayer; lay++){
211 if(0 == m_param.fgCalib[lay]) continue;
212
213 fqtlog << "Layer" << lay << endl;
214
215 for(ord=0; ord<m_qtorder; ord++){
216 m_qtpar[lay][ord] = calconst -> getQtpar(lay, ord);
217 }
218
219 for(bin=0; bin<m_nbin; bin++){
220 entry = m_hqt[lay][bin] -> GetEntries();
221
222 if(entry > 300){
223 deltw = m_hqt[lay][bin] -> GetMean();
224 qterr = ( m_hqt[lay][bin]->GetRMS() ) / sqrt((double)entry);
225 deltw /= m_vdr;
226 qterr /= m_vdr;
227 } else{
228 continue;
229 }
230
231 qbcen = ( (double)bin + 0.5 ) * m_qbinw[lay] + m_qmin[lay];
232// tw = qtFun(qbcen, m_qtpar[lay]) + deltw;
233 tw = (m_mdcFunSvc->getTimeWalk(lay, qbcen)) + deltw;
234
235 m_grqt[lay]->SetPoint(bin, qbcen, tw);
236 m_grqt[lay]->SetPointError(bin, 0, qterr);
237
238 m_grqdt[lay]->SetPoint(bin, qbcen, deltw);
239 m_grqdt[lay]->SetPointError(bin, 0, qterr);
240
241 fqtlog << setw(3) << bin
242 << setw(12) << deltw
243 << setw(12) << tw
244 << setw(12) << qbcen
245 << setw(12) << qterr
246 << endl;
247 }
248
249 m_grqt[lay]->Fit("funQt", "Q+", "", m_qmin[lay], m_qmax[lay]);
250
251 fqtlog << "Qtpar: ";
252 for(ord=0; ord<m_qtorder; ord++){
253 qtpar = funQt->GetParameter(ord);
254 qterr = funQt->GetParError(ord);
255 calconst -> resetQtpar(lay, ord, qtpar);
256
257 fqtlog << setw(12) << qtpar
258 << setw(12) << qterr
259 << endl;
260 }
261
262// if( (0 == ierflg) && (3 == istat) ){
263// for(ord=0; ord<m_qtorder; ord++){
264// gmqt -> GetParameter(ord, qtpar, qterr);
265// calconst -> resetQtpar(lay, ord, qtpar);
266
267// fqtlog << setw(12) << qtpar
268// << setw(12) << qterr
269// << endl;
270// }
271// } else{
272// fqtlog << setw(12) << m_qtpar[lay][0]
273// << setw(12) << "0"
274// << setw(12) << m_qtpar[lay][1]
275// << setw(12) << "0"
276// << endl;
277// }
278
279 } // end of layer loop
280
281 fqtlog.close();
282 delete funQt;
283 return 1;
284}
virtual double getTimeWalk(int layid, double Q) const =0
int fgCalib[MdcCalNLayer]
virtual int updateConst(MdcCalibConst *calconst)=0
static Double_t qtFun(Double_t *x, Double_t *par)

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