BOSS 7.1.1
BESIII Offline Software System
Loading...
Searching...
No Matches
MdcCalibFunSvc Class Reference

#include <MdcCalibFunSvc.h>

+ Inheritance diagram for MdcCalibFunSvc:

Public Member Functions

 MdcCalibFunSvc (const std::string &name, ISvcLocator *svcloc)
 
 ~MdcCalibFunSvc ()
 
virtual StatusCode initialize ()
 
virtual StatusCode finalize ()
 
void handle (const Incident &)
 
double getVprop (int lay) const
 
double getTprop (int lay, double z) const
 
double driftTimeToDist (double drifttime, int layid, int cellid, int lr, double entrance=0.0) const
 
double distToDriftTime (double dist, int layid, int cellid, int lr, double entrance=0.0) const
 
double getSigma (int layid, int lr, double dist, double entrance=0.0, double tanlam=0.0, double z=0.0, double Q=1000.0) const
 
double getSigmaLR (int layid, int lr, double dist, double entrance=0.0, double tanlam=0.0, double z=0.0, double Q=1000.0) const
 
double getSigma1 (int layid, int lr, double dist, double entrance=0.0, double tanlam=0.0, double z=0.0, double Q=1000.0) const
 
double getSigma2 (int layid, int lr, double dist, double entrance=0.0, double tanlam=0.0, double z=0.0, double Q=1000.0) const
 
double getF (int layid, int lr, double dist, double entrance=0.0, double tanlam=0.0, double z=0.0, double Q=1000.0) const
 
double getSigmaToT (int layid, int lr, double tdr, double entrance=0.0, double tanlam=0.0, double z=0.0, double Q=1000.0) const
 
double getSigmaToTLR (int layid, int lr, double tdr, double entrance=0.0, double tanlam=0.0, double z=0.0, double Q=1000.0) const
 
void setXtBegin ()
 
int getNextXtpar (int &key, double &par)
 
void getXtpar (int layid, int entr, int lr, double par[]) const
 
bool getNewXtpar ()
 
TTree * getNewXtparTree (int layid, int entr, int lr) const
 
double getT0 (int layid, int cellid) const
 
double getT0 (int wireid) const
 
double getDelT0 (int wireid) const
 
double getTimeWalk (int layid, double Q) const
 
double getQtpar (int layid, int ord) const
 
double getWireEff (int layid, int cellid) const
 
double getWireEff (int wireid) const
 
void setSdBegin ()
 
int getNextSdpar (int &key, double &par)
 
void getSdpar (int layid, int entr, int lr, double par[]) const
 
bool getR2tpar ()
 
TTree * getR2tTree (int layid) const
 
int getXtEntrIndex (double entrance) const
 
int getSdEntrIndex (double entrance) const
 

Public Attributes

int m_run
 

Detailed Description

Definition at line 21 of file MdcCalibFunSvc.h.

Constructor & Destructor Documentation

◆ MdcCalibFunSvc()

MdcCalibFunSvc::MdcCalibFunSvc ( const std::string & name,
ISvcLocator * svcloc )

Definition at line 40 of file MdcCalibFunSvc.cxx.

40 :
41 base_class (name, svcloc), m_layInfSig(-1) {
42
43 // declare properties
44 declareProperty("CheckConst", m_checkConst = false);
45 declareProperty("LayerInfSig", m_layInfSig);
46 declareProperty("XtMode", m_xtMode = 1);
47 declareProperty("NewXtFile", m_xtfile);
48 declareProperty("ReadWireEffDb", m_readWireEffDb = true);
49 declareProperty("WireEffFile", m_wireEffFile);
50 declareProperty("LinearXT", m_linearXT = false);
51 declareProperty("FixSigma", m_fixSigma = false);
52 declareProperty("FixSigmaValue", m_fixSigmaValue = 130.0); // micron
53 m_outputXtMode = true;
54}

◆ ~MdcCalibFunSvc()

MdcCalibFunSvc::~MdcCalibFunSvc ( )

Definition at line 56 of file MdcCalibFunSvc.cxx.

56 {
57}

Member Function Documentation

◆ distToDriftTime()

double MdcCalibFunSvc::distToDriftTime ( double dist,
int layid,
int cellid,
int lr,
double entrance = 0.0 ) const

Definition at line 216 of file MdcCalibFunSvc.cxx.

217 {
218 int i = 0;
219 double time;
220 int ord;
221 double xtpar[8];
222 double dxdtpar[5];
223 double x;
224 double dxdt;
225 double deltax;
226
227 int entr = getXtEntrIndex(entrance);
228 getXtpar(layid, entr, lr, xtpar);
229
230 double tm1 = xtpar[6];
231 double tm2 = 2000.0;
232 double dm1 = driftTimeToDist(tm1, layid, cellid, lr, entrance);
233 double dm2 = driftTimeToDist(tm2, layid, cellid, lr, entrance);
234
235 if(dist < 0){
236 cout << "Warning in MdcCalibFunSvc: driftDist < 0" << endl;
237 time = 0.0;
238 } else if(dist < xtpar[0]){
239 time = 0.0;
240 } else if(dist < dm1){
241 for(ord=0; ord<5; ord++){
242 dxdtpar[ord] = (double)(ord+1) * xtpar[ord+1];
243 }
244 time = dist / 0.03;
245 while(1){
246 if( i > 50 ){
247 cout << "Warning in MdcCalibFunSvc: "
248 << "can not get the exact value in the dist-to-time conversion."
249 << endl;
250 time = dist / 0.03;
251 break;
252 }
253
254 x = 0.0;
255 for(ord=0; ord<6; ord++)
256 x += xtpar[ord] * pow(time, ord);
257
258 deltax = x - dist;
259 if( fabs(deltax) < 0.001 ){
260 break;
261 }
262
263 dxdt = 0.0;
264 for(ord=0; ord<5; ord++)
265 dxdt += dxdtpar[ord] * pow(time, ord);
266
267 time = time - (deltax / dxdt);
268 i++;
269 }
270 } else if(dist < dm2){
271 time = (dist - dm1) * (tm2 - tm1) / (dm2 - dm1) + tm1;
272 } else{
273 time = tm2;
274 }
275
276 if(m_linearXT) time = dist / 0.03;
277 return time;
278}
Double_t x[10]
Double_t time
double driftTimeToDist(double drifttime, int layid, int cellid, int lr, double entrance=0.0) const
void getXtpar(int layid, int entr, int lr, double par[]) const
int getXtEntrIndex(double entrance) const

Referenced by EsTimeAlg::execute().

◆ driftTimeToDist()

double MdcCalibFunSvc::driftTimeToDist ( double drifttime,
int layid,
int cellid,
int lr,
double entrance = 0.0 ) const

Definition at line 151 of file MdcCalibFunSvc.cxx.

152 {
153 double dist;
154 if(0 == m_xtMode){
155 dist = t2dPoly(drifttime, layid, cellid, lr, entrance);
156 } else{
157 if((0==lr) || (1==lr)) dist = t2dInter(drifttime, layid, cellid, lr, entrance);
158 else{
159 double dl = t2dInter(drifttime, layid, cellid, lr, entrance);
160 double dr = t2dInter(drifttime, layid, cellid, lr, entrance);
161 dist = (dl + dr) * 0.5;
162 }
163 }
164// cout << setw(15) << drifttime << setw(15) << dist << endl;
165 if(m_linearXT) dist = 0.03 * drifttime;
166 return dist;
167}

Referenced by HoughHit::calDriftDist(), MdcxHit::d(), distToDriftTime(), MdcHit::driftDist(), TrkReco::execute(), and KalFitTrack::getDriftDist().

◆ finalize()

StatusCode MdcCalibFunSvc::finalize ( )
virtual

Definition at line 108 of file MdcCalibFunSvc.cxx.

108 {
109 MsgStream log(messageService(), name());
110 log << MSG::INFO << "MdcCalibFunSvc::finalize()" << endreq;
111
112 m_xtmap.clear();
113 m_t0.clear();
114 m_delt0.clear();
115 m_qtpar0.clear();
116 m_qtpar1.clear();
117 m_sdmap.clear();
118
119 return StatusCode::SUCCESS;
120}

◆ getDelT0()

double MdcCalibFunSvc::getDelT0 ( int wireid) const
inline

Definition at line 87 of file MdcCalibFunSvc.h.

87{ return m_delt0[wireid]; }

◆ getF()

double MdcCalibFunSvc::getF ( int layid,
int lr,
double dist,
double entrance = 0.0,
double tanlam = 0.0,
double z = 0.0,
double Q = 1000.0 ) const

Definition at line 354 of file MdcCalibFunSvc.cxx.

356 {
357
358 return 1.0;
359}

◆ getNewXtpar()

bool MdcCalibFunSvc::getNewXtpar ( )

Definition at line 418 of file MdcCalibFunSvc.cxx.

418 {
419 MsgStream log(messageService(), name());
420 log << MSG::INFO << "read calib const from TCDS" << endreq;
421
422 for(int layid=0; layid<NLAYER; layid++){
423 for(int entr=0; entr<NXTENTR; entr++){
424 for(int lr=0; lr<2; lr++){
425 double br_t,br_d;
426 TTree* newXtTree = getNewXtparTree(layid,entr,lr);
427 if(!newXtTree) return false;
428 newXtTree -> SetBranchAddress("t", &br_t);
429 newXtTree -> SetBranchAddress("d", &br_d);
430 int nEntries = newXtTree -> GetEntries();
431 if((nEntries<10) || (nEntries>=200)){
432 log << MSG::ERROR << "wrong X-T constants: layer " << layid
433 << ", iEntr " << entr << ", lr " << lr << endreq;
434 return false;
435 }
436 m_nNewXt[layid][entr][lr] = nEntries;
437 for(int i=0; i<nEntries; i++){
438 newXtTree->GetEntry(i);
439 m_vt[layid][entr][lr][i] = br_t;
440 m_vd[layid][entr][lr][i] = br_d;
441 }//end loop entries
442 }//end lr
443 }//end entr
444 }//end layid
445
446 return true;
447}
data SetBranchAddress("time",&time)
TTree * getNewXtparTree(int layid, int entr, int lr) const

◆ getNewXtparTree()

TTree * MdcCalibFunSvc::getNewXtparTree ( int layid,
int entr,
int lr ) const

Definition at line 449 of file MdcCalibFunSvc.cxx.

449 {
450 MsgStream log(messageService(), name());
451 string fullPath = "/Calib/MdcCal";
452 SmartDataPtr<CalibData::MdcCalibData> calConst(m_pCalDataSvc, fullPath);
453 if( ! calConst ){
454 log << MSG::ERROR << "can not get MdcCalibConst via SmartPtr" << endreq;
455 return NULL;
456 }
457
458 TTree* newXtTree = calConst->getNewXtpar(layid,entr,lr);
459 return newXtTree;
460}
#define NULL

Referenced by getNewXtpar().

◆ getNextSdpar()

int MdcCalibFunSvc::getNextSdpar ( int & key,
double & par )

Definition at line 558 of file MdcCalibFunSvc.cxx.

558 {
559 if( m_sditer != m_sdmap.end() ){
560 key = (*m_sditer).first;
561 par = (*m_sditer).second;
562 m_sditer++;
563 return 1;
564 }
565 else return 0;
566}
*************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

◆ getNextXtpar()

int MdcCalibFunSvc::getNextXtpar ( int & key,
double & par )

Definition at line 400 of file MdcCalibFunSvc.cxx.

400 {
401 if( m_xtiter != m_xtmap.end() ){
402 key = (*m_xtiter).first;
403 par = (*m_xtiter).second;
404 m_xtiter++;
405 return 1;
406 }
407 else return 0;
408}

◆ getQtpar()

double MdcCalibFunSvc::getQtpar ( int layid,
int ord ) const

Definition at line 538 of file MdcCalibFunSvc.cxx.

538 {
539 if(0 == ord) return m_qtpar0[layid];
540 else if(1 == ord) return m_qtpar1[layid];
541 else {
542 cout << "wrong order number" << endl;
543 return 0.0;
544 }
545}

Referenced by getTimeWalk().

◆ getR2tpar()

bool MdcCalibFunSvc::getR2tpar ( )

Definition at line 462 of file MdcCalibFunSvc.cxx.

462 {
463 for(int layid=0; layid<NLAYER; layid++){
464 int br_iEntr,br_lr;
465 double br_s,br_t;
466 TTree* r2tTree = getR2tTree(layid);
467 if(!r2tTree) return false;
468 r2tTree -> SetBranchAddress("iEntr", &br_iEntr);
469 r2tTree -> SetBranchAddress("lr", &br_lr);
470 r2tTree -> SetBranchAddress("s", &br_s);
471 r2tTree -> SetBranchAddress("t", &br_t);
472 int nEntries = r2tTree -> GetEntries();
473 for(int i=0; i<nEntries; i++){
474 r2tTree->GetEntry(i);
475 int bin = m_nR2t[layid][br_iEntr][br_lr];
476 if(bin>=200){
477 cout << "Error: number of sigma-time bins overflow" << endl;
478 return false;
479 }
480 m_tR2t[layid][br_iEntr][br_lr][bin] = br_t;
481 m_sR2t[layid][br_iEntr][br_lr][bin] = br_s;
482 m_nR2t[layid][br_iEntr][br_lr]++;
483 }
484 }
485 for(int layid=0; layid<NLAYER; layid++){
486 for(int iEntr=0; iEntr<NXTENTR; iEntr++){
487 for(int lr=0; lr<2; lr++){
488 if((m_nR2t[layid][iEntr][lr]<10) || (m_nR2t[layid][iEntr][lr]>=200)) return false;
489 }
490 }
491 }
492 return true;
493}
*******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
TTree * getR2tTree(int layid) const

◆ getR2tTree()

TTree * MdcCalibFunSvc::getR2tTree ( int layid) const

Definition at line 495 of file MdcCalibFunSvc.cxx.

495 {
496 MsgStream log(messageService(), name());
497 string fullPath = "/Calib/MdcCal";
498 SmartDataPtr<CalibData::MdcCalibData> calConst(m_pCalDataSvc, fullPath);
499 if( ! calConst ){
500 log << MSG::ERROR << "can not get MdcCalibConst via SmartPtr" << endreq;
501 return NULL;
502 }
503
504 TTree* r2tTree = calConst->getR2tpar(layid);
505 return r2tTree;
506}

Referenced by getR2tpar().

◆ getSdEntrIndex()

int MdcCalibFunSvc::getSdEntrIndex ( double entrance) const

Definition at line 599 of file MdcCalibFunSvc.cxx.

599 {
600 int i;
601 int index;
602 int idmax = 5;
603 double aglpi = 3.141592653;
604 double aglmin = -1.570796327; // -90 degree
605 double aglmax = 1.570796327; // 90 degree
606 double delAngle = 0.523598776; // 30 degree
607
608 MsgStream log(messageService(), name());
609 if(entrance < aglmin){
610 log << MSG::WARNING << "entrance angle < -pi/2" << endreq;
611 while(1){
612 entrance += aglpi;
613 if(entrance >= aglmin) break;
614 }
615 } else if(entrance > aglmax){
616 log << MSG::WARNING << "entrance angle > pi/2" << endreq;
617 while(1){
618 entrance -= aglpi;
619 if(entrance <= aglmax) break;
620 }
621 }
622
623 index = (int)((entrance-aglmin) / delAngle);
624 if(index < 0) index = 0;
625 else if(index > idmax) index = idmax;
626
627 return index;
628}

Referenced by getSigmaLR().

◆ getSdpar()

void MdcCalibFunSvc::getSdpar ( int layid,
int entr,
int lr,
double par[] ) const

Definition at line 547 of file MdcCalibFunSvc.cxx.

547 {
548 int parId;
549 if( (entr < 0) || (entr >= 18) ){
550 entr = 17;
551 }
552 for(int bin=0; bin<NSDBIN; bin++){
553 parId = getSdparId(layid, entr, lr, bin);
554 par[bin] = m_sdpar[parId];
555 }
556}

Referenced by getSigmaLR().

◆ getSigma()

double MdcCalibFunSvc::getSigma ( int layid,
int lr,
double dist,
double entrance = 0.0,
double tanlam = 0.0,
double z = 0.0,
double Q = 1000.0 ) const

Definition at line 280 of file MdcCalibFunSvc.cxx.

282 {
283 double sigma;
284 if( (0 == lr) || (1 == lr) ){
285 sigma = getSigmaLR(layid, lr, dist, entrance, tanlam, z, Q);
286 } else{
287 double sl = getSigmaLR(layid, 0, dist, entrance, tanlam, z, Q);
288 double sr = getSigmaLR(layid, 1, dist, entrance, tanlam, z, Q);
289 sigma = (sl + sr) * 0.5;
290 }
291
292 if(m_fixSigma) sigma = 0.001 * m_fixSigmaValue;
293 if(layid == m_layInfSig) sigma = 9999.0;
294 return sigma;
295}
TTree * sigma
double getSigmaLR(int layid, int lr, double dist, double entrance=0.0, double tanlam=0.0, double z=0.0, double Q=1000.0) const

Referenced by TrkReco::execute(), KalFitTrack::getSigma(), getSigma1(), and MdcHit::sigma().

◆ getSigma1()

double MdcCalibFunSvc::getSigma1 ( int layid,
int lr,
double dist,
double entrance = 0.0,
double tanlam = 0.0,
double z = 0.0,
double Q = 1000.0 ) const

Definition at line 340 of file MdcCalibFunSvc.cxx.

342 {
343 double sigma1 = getSigma(layid, lr, dist, entrance, tanlam, z, Q);
344 return sigma1;
345}
double getSigma(int layid, int lr, double dist, double entrance=0.0, double tanlam=0.0, double z=0.0, double Q=1000.0) const

◆ getSigma2()

double MdcCalibFunSvc::getSigma2 ( int layid,
int lr,
double dist,
double entrance = 0.0,
double tanlam = 0.0,
double z = 0.0,
double Q = 1000.0 ) const

Definition at line 347 of file MdcCalibFunSvc.cxx.

349 {
350
351 return 0.0;
352}

◆ getSigmaLR()

double MdcCalibFunSvc::getSigmaLR ( int layid,
int lr,
double dist,
double entrance = 0.0,
double tanlam = 0.0,
double z = 0.0,
double Q = 1000.0 ) const

Definition at line 297 of file MdcCalibFunSvc.cxx.

299 {
300
301 double sigma = 9999.0;
302 double par[NSDBIN];
303
304 int entr = getSdEntrIndex(entrance);
305 getSdpar(layid, entr, lr, par);
306
307 int nmaxBin;
308 double dw = 0.5; // width of each distance bin
309 double dmin = 0.25; // mm
310 if(layid < 8){
311 nmaxBin = 20; // 11->20 2011-12-10
312 } else{
313 nmaxBin = 20; // 15->20 2011-12-10
314 }
315
316 double dref[2];
317 double distAbs = fabs(dist);
318 if(distAbs < dmin){
319 sigma = par[0];
320 } else{
321 int bin = (int)((distAbs - dmin) / dw);
322 if(bin >= nmaxBin){
323 sigma = par[nmaxBin];
324 } else if(bin < 0){
325 sigma = par[0];
326 } else{
327 dref[0] = (double)bin * dw + 0.25;
328 dref[1] = (double)(bin+1) * dw + 0.25;
329 if((dref[1] - dref[0]) <= 0){
330 sigma = 9999.0;
331 } else{
332 sigma = (par[bin+1] - par[bin]) * (distAbs - dref[0]) /
333 (dref[1] - dref[0]) + par[bin];
334 }
335 }
336 }
337 return sigma;
338}
int getSdEntrIndex(double entrance) const
void getSdpar(int layid, int entr, int lr, double par[]) const

Referenced by getSigma().

◆ getSigmaToT()

double MdcCalibFunSvc::getSigmaToT ( int layid,
int lr,
double tdr,
double entrance = 0.0,
double tanlam = 0.0,
double z = 0.0,
double Q = 1000.0 ) const

Definition at line 361 of file MdcCalibFunSvc.cxx.

362 {
363 if(!m_fgR2t){
364 cout << "ERROR: can not get sigma-time functions" << endl;
365 return -999.0;
366 } else if( (0 == lr) || (1 == lr) ){
367 return getSigmaToTLR(layid, lr, tdr, entrance, tanlam, z, Q);
368 } else{
369 double sl = getSigmaToTLR(layid, 0, tdr, entrance, tanlam, z, Q);
370 double sr = getSigmaToTLR(layid, 1, tdr, entrance, tanlam, z, Q);
371 double sigma = (sl + sr) * 0.5;
372 return sigma;
373 }
374}
double getSigmaToTLR(int layid, int lr, double tdr, double entrance=0.0, double tanlam=0.0, double z=0.0, double Q=1000.0) const

◆ getSigmaToTLR()

double MdcCalibFunSvc::getSigmaToTLR ( int layid,
int lr,
double tdr,
double entrance = 0.0,
double tanlam = 0.0,
double z = 0.0,
double Q = 1000.0 ) const

Definition at line 376 of file MdcCalibFunSvc.cxx.

377 {
378 double sigma;
379 int iEntr = getXtEntrIndex(entrance);
380 int nBin = m_nR2t[layid][iEntr][lr];
381 if(tdr < m_tR2t[layid][iEntr][lr][0]){
382 sigma = m_sR2t[layid][iEntr][lr][0];
383 } else if(tdr < m_tR2t[layid][iEntr][lr][nBin-1]){
384 for(int i=0; i<(nBin-1); i++){
385 if((tdr>=m_tR2t[layid][iEntr][lr][i]) && (tdr<m_tR2t[layid][iEntr][lr][i+1])){
386 double t1 = m_tR2t[layid][iEntr][lr][i];
387 double t2 = m_tR2t[layid][iEntr][lr][i+1];
388 double s1 = m_sR2t[layid][iEntr][lr][i];
389 double s2 = m_sR2t[layid][iEntr][lr][i+1];
390 sigma = (tdr-t1) * (s2-s1) / (t2-t1) + s1;
391 break;
392 }
393 }
394 } else{
395 sigma = m_sR2t[layid][iEntr][lr][nBin-1];
396 }
397 return sigma;
398}

Referenced by getSigmaToT().

◆ getT0() [1/2]

double MdcCalibFunSvc::getT0 ( int layid,
int cellid ) const

Definition at line 509 of file MdcCalibFunSvc.cxx.

509 {
510 int wireid = m_pMdcGeomSvc->Wire(layid, cellid)->Id();
511 double t0 = getT0(wireid);
512
513 return t0;
514}
virtual const MdcGeoWire *const Wire(unsigned id)=0
double getT0(int layid, int cellid) const
int Id(void) const
Definition MdcGeoWire.h:127

Referenced by HoughHit::driftTime(), HoughHit::driftTime(), EsTimeAlg::execute(), TrkReco::execute(), Hough2D::fit(), Hough3D::fit(), KalFitTrack::getDriftTime(), getT0(), MdcxHit::process(), and MdcHit::setCalibSvc().

◆ getT0() [2/2]

double MdcCalibFunSvc::getT0 ( int wireid) const
inline

Definition at line 86 of file MdcCalibFunSvc.h.

86{ return m_t0[wireid]; }

◆ getTimeWalk()

double MdcCalibFunSvc::getTimeWalk ( int layid,
double Q ) const

Definition at line 516 of file MdcCalibFunSvc.cxx.

516 {
517 double tw = 0.0;
518 double qtpar[2];
519 int ord;
520
521 if(Q < 0.0001) Q = 0.0001;
522
523 for(ord=0; ord<2; ord++){
524 qtpar[ord] = getQtpar(layid, ord);
525 }
526
527 tw = qtpar[0] + qtpar[1] / sqrt( Q );
528 if(m_run < 0) tw = 0.0; // for MC
529
530 return tw;
531}
double getQtpar(int layid, int ord) const

Referenced by HoughHit::driftTime(), HoughHit::driftTime(), EsTimeAlg::execute(), TrkReco::execute(), Hough2D::fit(), Hough3D::fit(), KalFitTrack::getDriftTime(), MdcxHit::process(), and MdcHit::setCalibSvc().

◆ getTprop()

double MdcCalibFunSvc::getTprop ( int lay,
double z ) const

Definition at line 145 of file MdcCalibFunSvc.cxx.

145 {
146 double vp = getVprop(lay);
147 double tp = fabs(z - m_zst[lay]) / vp;
148 return tp;
149}
double getVprop(int lay) const

Referenced by HoughHit::driftTime(), HoughHit::driftTime(), MdcHit::driftTime(), and MdcxHit::tcor().

◆ getVprop()

double MdcCalibFunSvc::getVprop ( int lay) const
inline

Definition at line 215 of file MdcCalibFunSvc.h.

215 {
216 if(lay<8) return 220.0;
217 else return 240.0;
218}

Referenced by getTprop().

◆ getWireEff() [1/2]

double MdcCalibFunSvc::getWireEff ( int layid,
int cellid ) const

Definition at line 533 of file MdcCalibFunSvc.cxx.

533 {
534 int wireid = m_pMdcGeomSvc->Wire(layid, cellid)->Id();
535 return m_wireEff[wireid];
536}

◆ getWireEff() [2/2]

double MdcCalibFunSvc::getWireEff ( int wireid) const
inline

Definition at line 93 of file MdcCalibFunSvc.h.

93{ return m_wireEff[wireid]; }

◆ getXtEntrIndex()

int MdcCalibFunSvc::getXtEntrIndex ( double entrance) const

Definition at line 568 of file MdcCalibFunSvc.cxx.

568 {
569 int i;
570 int index;
571 int idmax = 17;
572 double aglpi = 3.141592653;
573 double aglmin = -1.570796327; // -90 degree
574 double aglmax = 1.570796327; // 90 degree
575 double delAngle = 0.174532925; // 10 degree
576
577 MsgStream log(messageService(), name());
578 if(entrance < aglmin){
579 log << MSG::WARNING << "entrance angle < -pi/2" << endreq;
580 while(1){
581 entrance += aglpi;
582 if(entrance >= aglmin) break;
583 }
584 } else if(entrance > aglmax){
585 log << MSG::WARNING << "entrance angle > pi/2" << endreq;
586 while(1){
587 entrance -= aglpi;
588 if(entrance <= aglmax) break;
589 }
590 }
591
592 index = (int)((entrance-aglmin) / delAngle);
593 if(index < 0) index = 0;
594 else if(index > idmax) index = idmax;
595
596 return index;
597}

Referenced by distToDriftTime(), and getSigmaToTLR().

◆ getXtpar()

void MdcCalibFunSvc::getXtpar ( int layid,
int entr,
int lr,
double par[] ) const

Definition at line 410 of file MdcCalibFunSvc.cxx.

410 {
411 int parId;
412 for(int ord=0; ord<8; ord++){
413 parId = getXtparId(layid, entr, lr, ord);
414 par[ord] = m_xtpar[parId];
415 }
416}

Referenced by distToDriftTime().

◆ handle()

void MdcCalibFunSvc::handle ( const Incident & inc)

Definition at line 122 of file MdcCalibFunSvc.cxx.

122 {
123 MsgStream log( messageService(), name() );
124 log << MSG::DEBUG << "handle: " << inc.type() << endreq;
125
126 if ( inc.type() == "NewRun" ){
127 log << MSG::DEBUG << "NewRun" << endreq;
128
129 if( ! initCalibConst() ){
130 log << MSG::ERROR
131 << "can not initilize Mdc Calib Constants" << endl
132 << " Please insert the following statement "
133 << "in your \"jobOption.txt\" "
134 << "before the include file of Mdc Reconstruction: "
135 << endl << " "
136 << "#include \"$CALIBSVCROOT/share/job-CalibData.txt\""
137 << endl
138 << " If still error, please contact with Wu Linghui "
139 << "([email protected])."
140 << endreq;
141 }
142 }
143}

◆ initialize()

StatusCode MdcCalibFunSvc::initialize ( )
virtual

Definition at line 68 of file MdcCalibFunSvc.cxx.

68 {
69 MsgStream log(messageService(), name());
70 log << MSG::INFO << "MdcCalibFunSvc::initialize()" << endreq;
71
72 StatusCode sc = Service::initialize();
73 if( sc.isFailure() ) return sc;
74
75 IIncidentSvc* incsvc;
76 sc = service("IncidentSvc", incsvc);
77 int priority = 100;
78 if( sc.isSuccess() ){
79 incsvc -> addListener(this, "NewRun", priority);
80 }
81
82 sc = service("CalibDataSvc", m_pCalDataSvc, true);
83 if( sc == StatusCode::SUCCESS ){
84 log << MSG::INFO << "Retrieve IDataProviderSvc" << endreq;
85 }else{
86 log << MSG::FATAL << "can not get IDataProviderSvc" << endreq;
87 }
88
89 sc = service("MdcGeomSvc", m_pMdcGeomSvc);
90 if( sc != StatusCode::SUCCESS ){
91 log << MSG::ERROR << "can not use MdcGeomSvc" << endreq;
92 return StatusCode::FAILURE;
93 }
94
95 if(m_fixSigma) cout << "Fix MDC sigma to " << m_fixSigmaValue << " micron." << endl;
96
97 m_updateNum = 0;
98 for(int wir=0; wir<6796; wir++) m_wireEff[wir] = 1.0;
99 for(int lay=0; lay<NLAYER; lay++){
100 for(int iEntr=0; iEntr<NXTENTR; iEntr++){
101 for(int lr=0; lr<2; lr++) m_nR2t[lay][iEntr][lr] = 0;
102 }
103 }
104
105 return StatusCode::SUCCESS;
106}

◆ setSdBegin()

void MdcCalibFunSvc::setSdBegin ( )
inline

Definition at line 95 of file MdcCalibFunSvc.h.

95{ m_sditer = m_sdmap.begin(); }

◆ setXtBegin()

void MdcCalibFunSvc::setXtBegin ( )
inline

Definition at line 79 of file MdcCalibFunSvc.h.

79{ m_xtiter = m_xtmap.begin(); }

Member Data Documentation

◆ m_run

int MdcCalibFunSvc::m_run

Definition at line 34 of file MdcCalibFunSvc.h.

Referenced by getTimeWalk().


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