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

#include <XtMdcCalib.h>

+ Inheritance diagram for XtMdcCalib:

Public Member Functions

 XtMdcCalib ()
 
 ~XtMdcCalib ()
 
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 ()
 
int getHxtKey (int lay, int entr, int lr, int bin) const
 
- Public Member Functions inherited from MdcCalib
 MdcCalib ()
 
virtual ~MdcCalib ()
 

Static Public Member Functions

static void fcnXT (Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag)
 
static void fcnXtEdge (Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag)
 
static double xtFun (double t, double xtpar[])
 

Static Public Attributes

static std::vector< double > XMEAS
 
static std::vector< double > TBINCEN
 
static std::vector< double > ERR
 
static double Tmax
 
static double Dmax
 
static std::vector< double > XMEASED
 
static std::vector< double > TBINCENED
 
static std::vector< double > ERRED
 

Detailed Description

Definition at line 10 of file XtMdcCalib.h.

Constructor & Destructor Documentation

◆ XtMdcCalib()

XtMdcCalib::XtMdcCalib ( )

Definition at line 34 of file XtMdcCalib.cxx.

34 {
35 m_nbin = 50; // m_nbin=50 before 2011-11-16
36 m_nLR = 3;
37 m_nxtpar = 8;
38
39 for(int lay=0; lay<43; lay++){
40 if(lay < 15) m_nEntr[lay] = 1;
41 else m_nEntr[lay] = 2; // for entr>0 and entr<0
42 }
43
44 m_tbinw = 10.0;
45 m_fgIni = false;
46}

◆ ~XtMdcCalib()

XtMdcCalib::~XtMdcCalib ( )

Definition at line 48 of file XtMdcCalib.cxx.

48 {
49}

Member Function Documentation

◆ clear()

void XtMdcCalib::clear ( )
virtual

Implements MdcCalib.

Definition at line 51 of file XtMdcCalib.cxx.

51 {
52 cout << "~XtMdcCalib" << endl;
53
54 unsigned int i;
55 for(i=0; i<m_hxt.size(); i++){
56 delete m_hxt[i];
57 if(0 == (i % 1000)) cout << "delete hxt[" << i << "]" << endl;
58 }
59 delete m_fdXt;
60 m_hxt.clear();
61 m_hxtmap.clear();
62
64}
virtual void clear()=0
Definition MdcCalib.cxx:80

◆ fcnXT()

void XtMdcCalib::fcnXT ( Int_t & npar,
Double_t * gin,
Double_t & f,
Double_t * par,
Int_t iflag )
static

Definition at line 531 of file XtMdcCalib.cxx.

532 {
533 unsigned int i;
534 int ord;
535 Double_t deta;
536 Double_t chisq = 0.0;
537 Double_t dfit;
538
539 for(i=0; i<TBINCEN.size(); i++){
540 dfit = 0;
541 for(ord=0; ord<=5; ord++){
542 dfit += par[ord] * pow(TBINCEN[i], ord);
543 }
544 deta = (dfit - XMEAS[i]) / ERR[i];
545 chisq += deta * deta;
546 }
547
548 f = chisq;
549}
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")
static std::vector< double > ERR
Definition XtMdcCalib.h:33
static std::vector< double > XMEAS
Definition XtMdcCalib.h:31
static std::vector< double > TBINCEN
Definition XtMdcCalib.h:32

Referenced by updateConst().

◆ fcnXtEdge()

void XtMdcCalib::fcnXtEdge ( Int_t & npar,
Double_t * gin,
Double_t & f,
Double_t * par,
Int_t iflag )
static

Definition at line 569 of file XtMdcCalib.cxx.

570 {
571 unsigned int i;
572 Double_t deta;
573 Double_t chisq = 0.0;
574 Double_t dfit;
575
576 for(i=0; i<TBINCENED.size(); i++){
577 dfit = par[0] * (TBINCENED[i] - Tmax) + Dmax;
578 deta = (dfit - XMEASED[i]) / ERRED[i];
579 chisq += deta * deta;
580 }
581
582 f = chisq;
583}
static double Dmax
Definition XtMdcCalib.h:36
static double Tmax
Definition XtMdcCalib.h:35
static std::vector< double > ERRED
Definition XtMdcCalib.h:39
static std::vector< double > XMEASED
Definition XtMdcCalib.h:37
static std::vector< double > TBINCENED
Definition XtMdcCalib.h:38

Referenced by updateConst().

◆ fillHist()

int XtMdcCalib::fillHist ( MdcCalEvent * event)
virtual

Implements MdcCalib.

Definition at line 114 of file XtMdcCalib.cxx.

114 {
115 IMessageSvc* msgSvc;
116 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
117 MsgStream log(msgSvc, "XtMdcCalib");
118 log << MSG::DEBUG << "XtMdcCalib::fillHist()" << endreq;
119
120 MdcCalib::fillHist(event);
121
122 // get EsTimeCol
123 bool esCutFg = event->getEsCutFlag();
124 if( ! esCutFg ) return -1;
125
126 int i;
127 int k;
128 int lay;
129 int iLR;
130 int iEntr;
131 int bin;
132 int key;
133 int hid;
134
135 double dr;
136 double dz;
137 double tdr;
138 double doca;
139 double resi;
140 double entrance;
141 int stat;
142
143 double xtpar[MdcCalXtNPars];
144 int nhitlay;
145 bool fgHitLay[MdcCalNLayer];
146
147 if(! m_fgIni){
148 for(lay=0; lay<MdcCalNLayer; lay++){
149 if(lay < 8) m_docaMax[lay] = m_param.maxDocaInner;
150 else m_docaMax[lay] = m_param.maxDocaOuter;
151 }
152
153 for(lay=0; lay<MdcCalNLayer; lay++){
154 for(iEntr=0; iEntr<m_nEntr[lay]; iEntr++){
155 for(iLR=0; iLR<MdcCalLR; iLR++){
156// m_mdcFunSvc -> getXtpar(lay, iEntr, iLR, xtpar);
157 if(0 == iEntr) m_mdcFunSvc -> getXtpar(lay, 8, iLR, xtpar);
158 else if(1 == iEntr) m_mdcFunSvc -> getXtpar(lay, 9, iLR, xtpar);
159 m_tm[lay][iEntr][iLR] = xtpar[6];
160 }
161 }
162 }
163 m_fgIni = true;
164 }
165
166 MdcCalRecTrk* rectrk;
167 MdcCalRecHit* rechit;
168
169 int nhit;
170 int ntrk = event -> getNTrk();
171 if((ntrk < m_param.nTrkCut[0]) || (ntrk > m_param.nTrkCut[1])) return -1;
172 for(i=0; i<ntrk; i++){
173 rectrk = event->getRecTrk(i);
174 nhit = rectrk -> getNHits();
175
176 // dr cut
177 dr = rectrk->getDr();
178 if(fabs(dr) > m_param.drCut) continue;
179
180 // dz cut
181 dz = rectrk->getDz();
182 if(fabs(dz) > m_param.dzCut) continue;
183
184 // momentum cut
185 double p = rectrk -> getP();
186 if((fabs(p) < m_param.pCut[0]) || (fabs(p) > m_param.pCut[1])) continue;
187
188 // select cosmic track with y<0
189 double phi0 = rectrk -> getPhi0();
190 double phiTrk = phi0 + CLHEP::halfpi;
191 if(phiTrk > CLHEP::twopi) phiTrk -= CLHEP::twopi;
192 if(m_param.cosmicDwTrk && (phiTrk<CLHEP::pi)) continue;
193
194 for(lay=0; lay<MdcCalNLayer; lay++){
195 fgHitLay[lay] = false;
196 }
197
198 for(k=0; k<nhit; k++){
199 rechit = rectrk -> getRecHit(k);
200 lay = rechit -> getLayid();
201 fgHitLay[lay] = true;
202 }
203
204 nhitlay = 0;
205 for(lay=0; lay<MdcCalNLayer; lay++){
206 if(fgHitLay[lay]) nhitlay++;
207 }
208 if(nhitlay < m_param.nHitLayCut) continue;
209
210// bool fgNoise = rectrk->getFgNoiseRatio();
211// if(m_param.noiseCut && (!fgNoise)) continue;
212
213 for(k=0; k<nhit; k++){
214 rechit = rectrk -> getRecHit(k);
215 lay = rechit -> getLayid();
216 doca = rechit -> getDocaInc();
217 iLR = rechit -> getLR();
218 entrance = rechit -> getEntra();
219 resi = rechit -> getResiInc();
220 tdr = rechit -> getTdrift();
221 stat = rechit -> getStat();
222
223 if(1 != stat) continue;
224
225 if( (fabs(doca) > m_docaMax[lay]) ||
226 (fabs(resi) > m_param.resiCut[lay]) ){
227 continue;
228 }
229
230 if(0 == lay){
231 if( ! fgHitLay[1] ) continue;
232 } else if(42 == lay){
233 if( ! fgHitLay[41] ) continue;
234 } else{
235 if( (!fgHitLay[lay-1]) && (!fgHitLay[lay+1]) ) continue;
236 }
237
238 iEntr = m_mdcFunSvc -> getXtEntrIndex(entrance);
239 if(1 == m_nEntr[lay]){
240 iEntr = 0;
241 } else if(2 == m_nEntr[lay]){
242 if(entrance > 0.0) iEntr = 1;
243 else iEntr = 0;
244 }
245
246 bin = (int)(tdr / m_tbinw);
247
248 // left-right separately
249 key = getHxtKey(lay, iEntr, iLR, bin);
250 if((bin < m_nbin) && (1 == m_hxtmap.count(key))){
251 hid = m_hxtmap[key];
252 m_hxt[hid] -> Fill( resi );
253 }
254
255 // left-right in common
256 key = getHxtKey(lay, iEntr, 2, bin);
257 if((bin < m_nbin) && (1 == m_hxtmap.count(key))){
258 hid = m_hxtmap[key];
259 m_hxt[hid] -> Fill( resi );
260 }
261
262 if(fabs(tdr - m_tm[lay][iEntr][iLR]) < 5){
263 key = getHxtKey(lay, iEntr, iLR, m_nbin);
264 if( 1 == m_hxtmap.count(key) ){
265 hid = m_hxtmap[key];
266 m_hxt[hid] -> Fill( resi );
267 }
268
269 key = getHxtKey(lay, iEntr, 2, m_nbin);
270 if( 1 == m_hxtmap.count(key) ){
271 hid = m_hxtmap[key];
272 m_hxt[hid] -> Fill( resi );
273 }
274 }
275
276 } // hit loop
277 } // track loop
278 return 1;
279}
curve Fill()
*******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 MdcCalXtNPars
const int MdcCalLR
IMessageSvc * msgSvc()
*************DOUBLE PRECISION m_pi *DOUBLE PRECISION m_HvecTau2 DOUBLE PRECISION m_HvClone2 DOUBLE PRECISION m_gamma1 DOUBLE PRECISION m_gamma2 DOUBLE PRECISION m_thet1 DOUBLE PRECISION m_thet2 INTEGER m_IFPHOT *COMMON c_Taupair $ !Spin Polarimeter vector first Tau $ !Spin Polarimeter vector second Tau $ !Clone Spin Polarimeter vector first Tau $ !Clone Spin Polarimeter vector second Tau $ !Random Euler angle for cloning st tau $ !Random Euler angle for cloning st tau $ !Random Euler angle for cloning st tau $ !Random Euler angle for cloning nd tau $ !Random Euler angle for cloning nd tau $ !Random Euler angle for cloning nd tau $ !phi of HvecTau1 $ !theta of HvecTau1 $ !phi of HvecTau2 $ !theta of HvecTau2 $ !super key
Definition Taupair.h:42
double maxDocaOuter
double pCut[2]
double maxDocaInner
double resiCut[MdcCalNLayer]
double getDz() const
double getDr() const
virtual int fillHist(MdcCalEvent *event)=0
Definition MdcCalib.cxx:730
int getHxtKey(int lay, int entr, int lr, int bin) const

◆ getHxtKey()

int XtMdcCalib::getHxtKey ( int lay,
int entr,
int lr,
int bin ) const

Definition at line 520 of file XtMdcCalib.cxx.

520 {
521 int key;
522
523 key = ( (lay << HXTLAYER_INDEX) & HXTLAYER_MASK ) |
524 ( (entr << HXTENTRA_INDEX) & HXTENTRA_MASK ) |
525 ( (lr << HXTLR_INDEX) & HXTLR_MASK ) |
526 ( (bin << HXTBIN_INDEX) & HXTBIN_MASK );
527
528 return key;
529}

Referenced by fillHist(), initialize(), and updateConst().

◆ initialize()

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

Implements MdcCalib.

Definition at line 67 of file XtMdcCalib.cxx.

68 {
69 IMessageSvc* msgSvc;
70 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
71 MsgStream log(msgSvc, "XtMdcCalib");
72 log << MSG::INFO << "XtMdcCalib::initialize()" << endreq;
73
74 m_hlist = hlist;
75 m_mdcGeomSvc = mdcGeomSvc;
76 m_mdcFunSvc = mdcFunSvc;
77 m_mdcUtilitySvc = mdcUtilitySvc;
78
79 MdcCalib::initialize(m_hlist, m_mdcGeomSvc, m_mdcFunSvc, m_mdcUtilitySvc);
80
81 int lay;
82 int iLR;
83 int iEntr;
84 int bin;
85 int hid;
86 int key;
87 char hname[200];
88 TH1D* hist;
89
90 m_fdXt = new TFolder("fdXt", "fdXt");
91 m_hlist -> Add(m_fdXt);
92
93 m_nlayer = m_mdcGeomSvc -> getLayerSize();
94
95 hid = 0;
96 for(lay=0; lay<m_nlayer; lay++){
97 for(iEntr=0; iEntr<m_nEntr[lay]; iEntr++){
98 for(iLR=0; iLR<m_nLR; iLR++){
99 for(bin=0; bin<=m_nbin; bin++){
100 sprintf(hname, "Hxt%02d_E%02d_LR%01d_B%02d", lay, iEntr, iLR, bin);
101 hist = new TH1D(hname, "", 300, -1.5, 1.5);
102 m_hxt.push_back( hist );
103 m_fdXt -> Add( hist );
104
105 key = getHxtKey(lay, iEntr, iLR, bin);
106 m_hxtmap.insert( valType2( key, hid ) );
107 hid++;
108 }
109 }
110 }
111 }
112}
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)
map< int, int >::value_type valType2
virtual void initialize(TObjArray *hlist, IMdcGeomSvc *mdcGeomSvc, IMdcCalibFunSvc *mdcFunSvc, IMdcUtilitySvc *mdcUtilitySvc)=0
Definition MdcCalib.cxx:214

◆ printCut()

void XtMdcCalib::printCut ( ) const
virtual

Implements MdcCalib.

Definition at line 281 of file XtMdcCalib.cxx.

281 {
283}
virtual void printCut() const =0

◆ setParam()

void XtMdcCalib::setParam ( MdcCalParams & param)
inlinevirtual

Implements MdcCalib.

Definition at line 80 of file XtMdcCalib.h.

80 {
81 MdcCalib::setParam(param);
82 m_param = param;
83}
virtual void setParam(MdcCalParams &param)=0
Definition MdcCalib.h:306

◆ updateConst()

int XtMdcCalib::updateConst ( MdcCalibConst * calconst)
virtual

Implements MdcCalib.

Definition at line 285 of file XtMdcCalib.cxx.

285 {
286 IMessageSvc* msgSvc;
287 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
288 MsgStream log(msgSvc, "XtMdcCalib");
289 log << MSG::INFO << "XtMdcCalib::updateConst()" << endreq;
290
291 MdcCalib::updateConst(calconst);
292
293 int i;
294 int hid;
295 int lay;
296 int iLR;
297 int iEntr;
298 int bin;
299 int ord;
300 int key;
301
302 Stat_t entry;
303 double xtpar;
304 double xterr;
305 double tbcen;
306 double deltx;
307 double xcor;
308 double xerr;
309 double xtini[8];
310 double xtfit[8];
311
312 Int_t ierflg;
313 Int_t istat;
314 Int_t nvpar;
315 Int_t nparx;
316 Double_t fmin;
317 Double_t edm;
318 Double_t errdef;
319 Double_t arglist[10];
320
321 TMinuit* gmxt = new TMinuit(6);
322 gmxt -> SetPrintLevel(-1);
323 gmxt -> SetFCN(fcnXT);
324 gmxt -> SetErrorDef(1.0);
325 gmxt -> mnparm(0, "xtpar0", 0, 0.1, 0, 0, ierflg);
326 gmxt -> mnparm(1, "xtpar1", 0, 0.1, 0, 0, ierflg);
327 gmxt -> mnparm(2, "xtpar2", 0, 0.1, 0, 0, ierflg);
328 gmxt -> mnparm(3, "xtpar3", 0, 0.1, 0, 0, ierflg);
329 gmxt -> mnparm(4, "xtpar4", 0, 0.1, 0, 0, ierflg);
330 gmxt -> mnparm(5, "xtpar5", 0, 0.1, 0, 0, ierflg);
331 arglist[0] = 0;
332 gmxt -> mnexcm("SET NOW", arglist, 0, ierflg);
333
334 TMinuit* gmxtEd = new TMinuit(1);
335 gmxtEd -> SetPrintLevel(-1);
336 gmxtEd -> SetFCN(fcnXtEdge);
337 gmxtEd -> SetErrorDef(1.0);
338 gmxtEd -> mnparm(0, "xtpar0", 0, 0.1, 0, 0, ierflg);
339 arglist[0] = 0;
340 gmxtEd -> mnexcm("SET NOW", arglist, 0, ierflg);
341
342 ofstream fxtlog("xtlog");
343 for(lay=0; lay<m_nlayer; lay++){
344 if(0 == m_param.fgCalib[lay]) continue;
345
346 for(iEntr=0; iEntr<m_nEntr[lay]; iEntr++){
347 for(iLR=0; iLR<m_nLR; iLR++){
348
349 fxtlog << "Layer " << setw(3) << lay << setw(3) << iEntr
350 << setw(3) << iLR << endl;
351
352 for(ord=0; ord<m_nxtpar; ord++){
353// xtpar = calconst -> getXtpar(lay, iEntr, iLR, ord);
354 if(0 == iEntr) xtpar = calconst -> getXtpar(lay, 8, iLR, ord);
355 else if(1 == iEntr) xtpar = calconst -> getXtpar(lay, 9, iLR, ord);
356
357 xtini[ord] = xtpar;
358 xtfit[ord] = xtpar;
359 }
360
361 Tmax = xtini[6];
362
363 for(bin=0; bin<=m_nbin; bin++){
364 key = getHxtKey(lay, iEntr, iLR, bin);
365 hid = m_hxtmap[key];
366
367 entry = m_hxt[hid] -> GetEntries();
368 if(entry > 100){
369 deltx = m_hxt[hid] -> GetMean();
370 xerr = m_hxt[hid]->GetRMS();
371 } else{
372 continue;
373 }
374
375 if(bin < m_nbin)
376 tbcen = ( (double)bin + 0.5 ) * m_tbinw;
377 else tbcen = m_tm[lay][iEntr][iLR];
378 xcor = xtFun(tbcen, xtini) - deltx;
379
380 if((tbcen <= Tmax) || (bin == m_nbin)){
381 TBINCEN.push_back( tbcen );
382 XMEAS.push_back( xcor );
383 ERR.push_back( xerr );
384 } else{
385 TBINCENED.push_back( tbcen );
386 XMEASED.push_back( xcor );
387 ERRED.push_back( xerr );
388 }
389 fxtlog << setw(3) << bin
390 << setw(15) << deltx
391 << setw(15) << xcor
392 << setw(15) << tbcen
393 << setw(15) << xerr
394 << endl;
395 } // end of bin loop
396
397 if( XMEAS.size() < 12 ){
398 TBINCEN.clear();
399 XMEAS.clear();
400 ERR.clear();
401
402 TBINCENED.clear();
403 XMEASED.clear();
404 ERRED.clear();
405
406 continue;
407 }
408
409 for(ord=0; ord<=5; ord++){
410 arglist[0] = ord + 1;
411 arglist[1] = xtini[ord];
412 gmxt -> mnexcm("SET PARameter", arglist, 2, ierflg);
413 }
414
415 // fix the xtpar[0] at 0
416 if(1 == m_param.fixXtC0){
417 arglist[0] = 1;
418 arglist[1] = 0.0;
419 gmxt -> mnexcm("SET PARameter", arglist, 2, ierflg);
420 gmxt -> mnexcm("FIX", arglist, 1, ierflg);
421 }
422
423 arglist[0] = 1000;
424 arglist[1] = 0.1;
425 gmxt -> mnexcm("MIGRAD", arglist, 2, ierflg);
426 gmxt -> mnstat(fmin, edm, errdef, nvpar, nparx, istat);
427
428 fxtlog << "Xtpar: " << endl;
429 if( (0 == ierflg) && (istat >= 2) ){
430 for(ord=0; ord<=5; ord++){
431 gmxt -> GetParameter(ord, xtpar, xterr);
432// calconst -> resetXtpar(lay, iEntr, iLR, ord, xtpar);
433 xtfit[ord] = xtpar;
434
435 if(1 == m_nEntr[lay]){
436 for(i=0; i<18; i++)
437 calconst -> resetXtpar(lay, i, iLR, ord, xtpar);
438 } else if(2 == m_nEntr[lay]){
439 if(0 == iEntr){
440 for(i=0; i<9; i++) // entr<0
441 calconst->resetXtpar(lay, i, iLR, ord, xtpar);
442 } else{
443 for(i=9; i<18; i++) // entr>0
444 calconst->resetXtpar(lay, i, iLR, ord, xtpar);
445 }
446 }
447 fxtlog << setw(15) << xtpar << setw(15) << xterr << endl;
448 }
449 } else{
450 for(ord=0; ord<=5; ord++){
451 fxtlog << setw(15) << xtini[ord] << setw(15) << "0" << endl;
452 }
453 }
454 fxtlog << setw(15) << Tmax << setw(15) << "0" << endl;
455
456 // release the first parameter
457 if(1 == m_param.fixXtC0){
458 arglist[0] = 1;
459 gmxt -> mnexcm("REL", arglist, 1, ierflg);
460 }
461
462 Dmax = xtFun(Tmax, xtfit);
463
464 if( XMEASED.size() >= 3 ){
465 // fit xt in the edge area
466 arglist[0] = 1;
467 arglist[1] = xtini[7];
468 gmxtEd -> mnexcm("SET PARameter", arglist, 2, ierflg);
469
470 arglist[0] = 1000;
471 arglist[1] = 0.1;
472 gmxtEd -> mnexcm("MIGRAD", arglist, 2, ierflg);
473 gmxtEd -> mnstat(fmin, edm, errdef, nvpar, nparx, istat);
474
475 if( (0 == ierflg) && (istat >=2) ){
476 gmxtEd -> GetParameter(0, xtpar, xterr);
477 if(xtpar < 0.0) xtpar = 0.0;
478 if(m_param.fixXtEdge) xtpar = xtini[7];
479// calconst -> resetXtpar(lay, iEntr, iLR, 7, xtpar);
480
481 if(1 == m_nEntr[lay]){
482 for(i=0; i<18; i++)
483 calconst -> resetXtpar(lay, i, iLR, 7, xtpar);
484 } else if(2 == m_nEntr[lay]){
485 if(0 == iEntr){
486 for(i=0; i<9; i++)
487 calconst->resetXtpar(lay, i, iLR, 7, xtpar);
488 } else{
489 for(i=9; i<18; i++)
490 calconst->resetXtpar(lay, i, iLR, 7, xtpar);
491 }
492 }
493 fxtlog << setw(15) << xtpar << setw(15) << xterr << endl;
494 } else {
495 fxtlog << setw(15) << xtini[7] << setw(15) << "0" << endl;
496 }
497 } else {
498 fxtlog << setw(15) << xtini[7] << setw(15) << "0" << endl;
499 }
500 fxtlog << "Tm " << setw(15) << Tmax
501 << " Dmax " << setw(15) << Dmax << endl;
502
503 TBINCEN.clear();
504 XMEAS.clear();
505 ERR.clear();
506 TBINCENED.clear();
507 XMEASED.clear();
508 ERRED.clear();
509 } // end of l-r loop
510 } // end of entrance angle loop
511 } // end of layer loop
512
513 fxtlog.close();
514 delete gmxt;
515 delete gmxtEd;
516 cout << "Xt update finished" << endl;
517 return 1;
518}
int fgCalib[MdcCalNLayer]
void resetXtpar(int lay, int entr, int lr, int order, double val)
virtual int updateConst(MdcCalibConst *calconst)=0
static double xtFun(double t, double xtpar[])
static void fcnXT(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag)
static void fcnXtEdge(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag)
double xerr[1000]

◆ xtFun()

double XtMdcCalib::xtFun ( double t,
double xtpar[] )
static

Definition at line 551 of file XtMdcCalib.cxx.

551 {
552 int ord;
553 double dist = 0.0;
554 double tmax = xtpar[6];
555
556 if(t < tmax ){
557 for(ord=0; ord<=5; ord++){
558 dist += xtpar[ord] * pow(t, ord);
559 }
560 }else{
561 for(ord=0; ord<=5; ord++){
562 dist += xtpar[ord] * pow(tmax, ord);
563 }
564 dist += xtpar[7] * (t - tmax);
565 }
566
567 return dist;
568}
TTree * t
Definition binning.cxx:23

Referenced by updateConst().

Member Data Documentation

◆ Dmax

double XtMdcCalib::Dmax
static

Definition at line 36 of file XtMdcCalib.h.

Referenced by fcnXtEdge(), and updateConst().

◆ ERR

vector< double > XtMdcCalib::ERR
static

Definition at line 33 of file XtMdcCalib.h.

Referenced by fcnXT(), and updateConst().

◆ ERRED

vector< double > XtMdcCalib::ERRED
static

Definition at line 39 of file XtMdcCalib.h.

Referenced by fcnXtEdge(), and updateConst().

◆ TBINCEN

vector< double > XtMdcCalib::TBINCEN
static

Definition at line 32 of file XtMdcCalib.h.

Referenced by fcnXT(), and updateConst().

◆ TBINCENED

vector< double > XtMdcCalib::TBINCENED
static

Definition at line 38 of file XtMdcCalib.h.

Referenced by fcnXtEdge(), and updateConst().

◆ Tmax

double XtMdcCalib::Tmax
static

Definition at line 35 of file XtMdcCalib.h.

Referenced by fcnXtEdge(), and updateConst().

◆ XMEAS

vector< double > XtMdcCalib::XMEAS
static

Definition at line 31 of file XtMdcCalib.h.

Referenced by fcnXT(), and updateConst().

◆ XMEASED

vector< double > XtMdcCalib::XMEASED
static

Definition at line 37 of file XtMdcCalib.h.

Referenced by fcnXtEdge(), and updateConst().


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