StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StTofCalibMaker.cxx
1 /*******************************************************************
2  *
3  * $Id: StTofCalibMaker.cxx,v 1.24 2018/02/26 23:26:50 smirnovd Exp $
4  *
5  * Author: Xin Dong
6  *****************************************************************
7  *
8  * Description: - Tof Calibration Maker to do the calibration for pVPD
9  * (start timing) , TOFp and TOFr
10  * - store into StTofHit
11  * - if no valid calibration parameters, store matched hits
12  *
13  *
14  *******************************************************************/
15 #include <iostream>
16 #include "TFile.h"
17 #include "TProfile.h"
18 #include "TNtuple.h"
19 #include "TH1.h"
20 #include "TH2.h"
21 #include "TF1.h"
22 
23 #include "StEvent.h"
24 #include "StTofHit.h"
25 #include "StTofCell.h"
26 #include "StTofSlat.h"
27 #include "StTofData.h"
28 #include "StTofPidTraits.h"
29 #include "StEventTypes.h"
30 #include "Stypes.h"
31 #include "StMessMgr.h"
32 #include "StThreeVectorD.hh"
33 #include "StHelix.hh"
34 #include "StTrackGeometry.h"
35 #include "StEventUtilities/StuRefMult.hh"
36 #include "PhysicalConstants.h"
37 #include "phys_constants.h"
38 #include "StPhysicalHelixD.hh"
39 #include "tables/St_tofTzero_Table.h"
40 #include "tables/St_tofTACorr_Table.h"
41 #include "tables/St_tofCorrection_Table.h"
42 #include "tables/St_tofAdcRange_Table.h"
43 #include "tables/St_tofResolution_Table.h"
44 
45 #include "tables/St_tofr5INLtable_Table.h"
46 #include "tables/St_tofTotCorr_Table.h"
47 #include "tables/St_tofZCorr_Table.h"
48 
49 #include "tables/St_vertexSeed_Table.h"
50 
51 #include "tables/St_tofTOffset_Table.h"
52 #include "tables/St_tofPhaseOffset_Table.h"
53 
54 #include "StTofUtil/tofPathLength.hh"
55 #include "StTofUtil/StTofDataCollection.h"
56 #include "StTofUtil/StTofSlatCollection.h"
57 #include "StTofUtil/StTofCellCollection.h"
58 #include "StTofUtil/StTofHitCollection.h"
59 #include "StTofUtil/StTofGeometry.h"
60 #include "StTofCalibMaker.h"
61 
62 #include "StMessMgr.h"
63 #include "StMemoryInfo.hh"
64 #include "StTimer.hh"
65 
66 
67 //_____________________________________________________________________________
68 StTofCalibMaker::StTofCalibMaker(const char *name) : StMaker(name)
69 {
73  setTDCLimits(20, 1500); // TDC range
74  setADCCut(300); // TA correction adc cut
75  setTDCWidth(0.05); // 50 ps
76 
77  setPVPDADCLimits(30,1100);
78  setPVPDTDCLimits(1,2000);
79  setPVPDHitsCut(1,1);
80  setOuterGeometry(true);
81 
82  mValidStartTime = kTRUE;
83  mEastPVPDValid = kTRUE;
84  mSlewingCorr = kTRUE;
85 
86  StThreeVectorD MomFstPt(0.,0.,9999.);
87  StThreeVectorD origin(0.,0.,0.);
88  mBeamHelix = new StPhysicalHelixD(MomFstPt,origin,0.5*tesla,1.);
89 
90  resetPars();
91 }
92 
93 //_____________________________________________________________________________
95 {
96  resetPars();
97 }
98 
99 //_____________________________________________________________________________
101 {
102  for(int i=0;i<mNPar;i++) {
103  mTofrZPar[i] = 0.;
104  mTofpZPar[i] = 0.;
105  }
106  for(int i=0;i<mNTOFr;i++) {
107  mTofrT0[i] = 0.;
108  }
109  for(int i=0;i<mNTOFr*mNPar;i++) {
110  mTofrTAPar[i] = 0.0;
111  }
112  for(int i=0;i<mNTOFp;i++) {
113  mTofpT0[i] = 0.;
114  }
115  for(int i=0;i<mNTOFp*mNPar;i++) {
116  mTofpTAPar[i] = 0.0;
117  }
118  for(int i=0;i<mNPVPD*mNPar;i++) {
119  mPVPDTAPar[i] = 0.0;
120  }
121 
122  mValidCalibPar = kFALSE;
123 
125  for(int i=0;i<mNTOFr;i++) {
126  mTofrADCMin[i] = 0.;
127  mTofrADCMax[i] = 1024.;
128  mTofrRes[i] = 999.;
129  }
130  for(int i=0;i<mNTOFp;i++) {
131  mTofpADCMin[i] = 0.;
132  mTofpADCMax[i] = 1024.;
133  mTofpRes[i] = 999.;
134  }
135  for(int i=0;i<mNPVPD;i++) {
136  mPVPDRes[i] = 999.;
137  }
138 
139  // run5
140  for(int i=0;i<mTdigBoard;i++) {
141  for(int j=0;j<mTdcOnBoard;j++) {
142  for(int k=0;k<mTdcChannel;k++) {
143  mINLtable[i][j][k] = 0.0;
144  }
145  }
146  }
147 
148  for(int i=0;i<mNTOFr5;i++) {
149  for(int j=0;j<mNBinMax;j++) {
150  mTofr5TotEdge[i][j] = 0.0;
151  mTofr5TotCorr[i][j] = 0.0;
152  mTofr5ZEdge[i][j] = 0.0;
153  mTofr5ZCorr[i][j] = 0.0;
154  }
155  }
156 
157  for(int i=0;i<mNPVPD;i++) {
158  for(int j=0;j<mNBinMax;j++) {
159  mPVPDTotEdge[i][j] = 0.0;
160  mPVPDTotCorr[i][j] = 0.0;
161  }
162  }
163 
164  // run8
165  for(int i=0;i<mNTray;i++) {
166  for(int j=0;j<mNTDIG;j++) {
167  for(int k=0;k<mNBinMax;k++) {
168  mTofTotEdge[i][j][k] = 0.0;
169  mTofTotCorr[i][j][k] = 0.0;
170  mTofZEdge[i][j][k] = 0.0;
171  mTofZCorr[i][j][k] = 0.0;
172  }
173  }
174  for(int j=0;j<mNModule;j++) {
175  for(int k=0;k<mNCell;k++) {
176  mTofTZero[i][j][k] = 0.0;
177  }
178  }
179  }
180  for(int i=0;i<2*mNVPD;i++) {
181  for(int j=0;j<mNBinMax;j++) {
182  mVPDTotEdge[i][j] = 0.0;
183  mVPDTotCorr[i][j] = 0.0;
184  }
185  mVPDTZero[i] = 0.0;
186  mVPDLeTime[i] = 0.0;
187  mVPDTot[i] = 0.0;
188  }
189 
190  mTStart = -9999.;
191  mTDiff = -9999.;
192  mVPDVtxZ = -9999.;
193  mProjVtxZ = -9999.;
194  mVPDHitPatternEast = 0;
195  mVPDHitPatternWest = 0;
196 }
197 
198 //____________________________________________________________________________
199 Int_t StTofCalibMaker::Init()
200 {
201  initFormulas();
202  return kStOK;
203 }
204 
205 //_____________________________________________________________________________
207 {
209  mTofrSlewing = new TF1("TofrSlewing", "[0]+[1]/sqrt(x)+[2]/x+[3]/sqrt(x)/x+[4]/x/x");
210  // changed back
211  mTofpSlewing = new TF1("TofpSlewing", "[0]+[1]/sqrt(x)+[2]/x+[3]/sqrt(x)/x+[4]/x/x");
212  // Run 4, AuAu200GeV, Jiansong's calibration function -- removed later
213  // mTofpSlewing = new TF1("TofpSlewing","[0]+[1]*sqrt(x)+[2]*x+[3]*x*sqrt(x)");
214  mTofrZCorr = new TF1("TofrZCorr", "pol7");
215  // mTofpZCorr = new TF1("TofpZCorr", "pol7");
216  // Run 4, AuAu200GeV, Jiansong's calibration function
217  mTofpZCorr = new TF1("TofpZCorr", "[0]+[1]*sqrt(x)+[2]*x+[3]*x*sqrt(x)");
218  mPVPDSlewing = new TF1("pVPDSlewing","[0]+[1]/sqrt(x)+[2]/x+[3]*x");
219 }
220 
221 //____________________________________________________________________________
222 Int_t StTofCalibMaker::InitRun(int runnumber)
223 {
224  // tof run configurations
225 
226  mTofpGeom = new StTofGeometry();
227  mTofpGeom->init(this);
228 
229  Int_t val = kStOK;
230  val = initParameters(runnumber);
231  if(val==kStOK) {
232  mValidCalibPar = kTRUE;
233  } else {
234  mValidCalibPar = kFALSE;
235  }
236 
237  if(mValidCalibPar) {
238  gMessMgr->Info(" ==> Good! Valid cali parameters! ","OS");
239  } else {
240  gMessMgr->Info(" ==> No valid cali parameters! ","OS");
241  }
242 
243  // RUN IV 023-035, east PVPD dead
244  if(runnumber>5023000&&runnumber<5035000) {
245  mEastPVPDValid = kFALSE;
246  } else {
247  mEastPVPDValid = kTRUE;
248  }
249 
250  return kStOK;
251 
252 }
253 
254 //_____________________________________________________________________________
256 {
257  mYear2 = (runnumber<4000000);
258  mYear3 = (runnumber>4000000&&runnumber<5000000);
259  mYear4 = (runnumber>5000000&&runnumber<6000000);
260  mYear5 = (runnumber>6000000&&runnumber<7000000);
261  mYear8 = (runnumber>9000000&&runnumber<10000000);
264  gMessMgr->Info("","OS") << " -- retrieving run parameters from Calibrations_tof" << endm;
265  TDataSet *mDbDataSet = GetDataBase("Calibrations/tof/tofTzero");
266  if (!mDbDataSet){
267  gMessMgr->Error("unable to get TOF run parameters","OS");
268  return kStErr;
269  }
270 
274  if(mYear2||mYear3||mYear4) {
275  gMessMgr->Info(" loading parameters for Run II/III/IV", "OS");
276  // -- T0 parameters --
277  St_tofTzero* tofT0 = static_cast<St_tofTzero*>(mDbDataSet->Find("tofTzero"));
278  if(!tofT0) {
279  gMessMgr->Error("unable to get Tzero table","OS");
280  return kStErr;
281  } else {
282  tofTzero_st* t0 = static_cast<tofTzero_st*>(tofT0->GetArray());
283  Int_t nRows = t0[0].entries;
284  if(nRows<0||nRows>mNMax) {
285  gMessMgr->Error("# of Tzero out of range","OS");
286  return kStErr;
287  }
288  for(Int_t i=0;i<nRows;i++) {
289  int daqId = t0[0].daqChannel[i];
290  int tdcChan = t0[0].tdcChan[i];
291  // check the daqId
292  if(daqId<0) continue;
293  if(tdcChan<42) {
294  mTofpT0[daqId] = (Double_t)(t0[0].Tzero[i]);
295  if(Debug()) {
296  gMessMgr->Info("","OS") << " -TOFp- daqId=" << daqId << " tdcChan=" << tdcChan << " T0=" << mTofpT0[daqId] << endm;
297  }
298  } else {
299  mTofrT0[daqId] = (Double_t)(t0[0].Tzero[i]);
300  if(Debug()) {
301  gMessMgr->Info("","OS") << " -TOFr- daqId=" << daqId << " tdcChan=" << tdcChan << " T0=" << mTofrT0[daqId] << endm;
302  }
303  }
304  }
305  }
306  // -- TA slewing parameters --
307  St_tofTACorr *tofTA = static_cast<St_tofTACorr*>(mDbDataSet->Find("tofTACorr"));
308  if(!tofTA) {
309  gMessMgr->Error("unable to find TA slewing parameters","OS");
310  return kStErr;
311  }
312  tofTACorr_st *tofta = static_cast<tofTACorr_st*>(tofTA->GetArray());
313  for(Int_t i=0;i<mNMax;i++) {
314  int daqId = tofta[0].daqChannel[i];
315  int tdcChan = tofta[0].tdcChan[i];
316  // check the daqId
317  if(daqId<0) continue;
318  for(int j=0;j<mNPar;j++) {
319  int ijdaq = daqId*mNPar+j;
320  int ij = i*mNPar+j;
321  if(tdcChan<42) {
322  if(daqId>=mNTOFp) {
323  gMessMgr->Warning("More than expected TOFp channels read in","OS");
324  } else {
325  mTofpTAPar[ijdaq] = (Double_t)(tofta[0].a[ij]);
326  if(Debug()) {
327  gMessMgr->Info("","OS") << " -TOFp- daqId=" << daqId << " tdcChan=" << tdcChan << " TA Corr[" << j << "]=" << mTofpTAPar[ijdaq] << endm;
328  }
329  }
330  } else if(tdcChan<48) {
331  if(daqId>=mNPVPD) {
332  gMessMgr->Warning("More than expected pVPD channels read in","OS");
333  } else {
334  mPVPDTAPar[ijdaq] = (Double_t)(tofta[0].a[ij]);
335  if(Debug()) {
336  gMessMgr->Info("","OS") << " -pVPD- daqId=" << daqId << " tdcChan=" << tdcChan << " TA Corr[" << j << "]=" << mPVPDTAPar[ijdaq] << endm;
337  }
338  }
339  } else {
340  if(daqId>=mNTOFr) {
341  gMessMgr->Warning("More than expected TOFr channels read in","OS");
342  } else {
343  mTofrTAPar[ijdaq] = (Double_t)(tofta[0].a[ij]);
344  if(Debug()) {
345  gMessMgr->Info("","OS") << " -TOFr- daqId=" << daqId << " tdcChan=" << tdcChan << " TA Corr[" << j << "]=" << mTofrTAPar[ij] << endm;
346  }
347  }
348  }
349  }
350  }
351 
352  // -- Z corr parameters --
353  St_tofCorrection *tofrZCorr = static_cast<St_tofCorrection*>(mDbDataSet->Find("tofrZhitCorr"));
354  St_tofCorrection *tofpZCorr = static_cast<St_tofCorrection*>(mDbDataSet->Find("tofpZhitCorr"));
355  if(!tofrZCorr && !tofpZCorr) {
356  gMessMgr->Error("unable to find Zhit corr parameters","OS");
357  return kStErr;
358  }
359  if(tofrZCorr) {
360  tofCorrection_st *tofrZ = static_cast<tofCorrection_st*>(tofrZCorr->GetArray());
361  for(int i=0;i<mNPar;i++) {
362  mTofrZPar[i] = tofrZ[0].a[i];
363  if(Debug()) {
364  gMessMgr->Info("","OS") << " -TOFr- Zcorr[" << i << "]=" << mTofrZPar[i] << endm;
365  }
366  }
367  }
368  if(tofpZCorr) {
369  tofCorrection_st *tofpZ = static_cast<tofCorrection_st*>(tofpZCorr->GetArray());
370  for(int i=0;i<mNPar;i++) {
371  mTofpZPar[i] = tofpZ[0].a[i];
372  if(Debug()) {
373  gMessMgr->Info("","OS") << " -TOFp- Zcorr[" << i << "]=" << mTofpZPar[i] << endm;
374  }
375  }
376  }
377 
378  // ADC range parameters ( specially for TOFp in Run IV )
379  St_tofAdcRange *tofAdc = static_cast<St_tofAdcRange*>(mDbDataSet->Find("tofAdcRange"));
380  if(!tofAdc) {
381  gMessMgr->Warning("unable to find ADC range parameters, use default values!","OS");
382  }
383  tofAdcRange_st *tofadc = static_cast<tofAdcRange_st*>(tofAdc->GetArray());
384  for(Int_t i=0;i<mNMax;i++) {
385  int daqId = tofadc[0].daqChannel[i];
386  int adcChan = tofadc[0].adcChan[i];
387  if(daqId<0) continue;
388  if(adcChan<42) {
389  if(daqId>=mNTOFp) {
390  gMessMgr->Warning("More than expected TOFp channels read in","OS");
391  } else {
392  mTofpADCMin[daqId] = tofadc[0].adcMin[i];
393  mTofpADCMax[daqId] = tofadc[0].adcMax[i];
394  if(Debug()) {
395  gMessMgr->Info("","OS") << " -TOFp- daqId=" << daqId << " adcChan=" << adcChan << " min=" << mTofpADCMin[daqId] << " max=" << mTofpADCMax[daqId] << endm;
396  }
397  }
398  }
399  if(adcChan>=60) {
400  if(daqId>=mNTOFr) {
401  gMessMgr->Warning("More than expected TOFr channels read in","OS");
402  } else {
403  mTofrADCMin[daqId] = tofadc[0].adcMin[i];
404  mTofrADCMax[daqId] = tofadc[0].adcMax[i];
405  if(Debug()) {
406  gMessMgr->Info("","OS") << " -TOFr- daqId=" << daqId << " adcChan=" << adcChan << " min=" << mTofrADCMin[daqId] << " max=" << mTofrADCMax[daqId] << endm;
407  }
408  }
409  }
410  }
411 
412  // res parameters ( specially for TOFp in Run IV )
413  St_tofResolution *tofRes = static_cast<St_tofResolution*>(mDbDataSet->Find("tofResolution"));
414  if(!tofRes) {
415  gMessMgr->Warning("unable to find resolution parameters, nSimgaTof UNAVAILABLE!","OS");
416  }
417  tofResolution_st *tofres = static_cast<tofResolution_st*>(tofRes->GetArray());
418  for(Int_t i=0;i<mNMax;i++) {
419  int daqId = tofres[0].daqChannel[i];
420  int tdcChan = tofres[0].tdcChan[i];
421  if(daqId<0) continue;
422  if(tdcChan<42) {
423  if(daqId>=mNTOFp) {
424  gMessMgr->Warning("More than expected TOFp channels read in","OS");
425  } else {
426  mTofpRes[daqId] = tofres[0].resolution[i];
427  if(Debug()) {
428  gMessMgr->Info("","OS") << " -TOFp- daqId=" << daqId << " tdcChan=" << tdcChan << " resolution=" << mTofpRes[daqId] << endm;
429  }
430  }
431  } else if(tdcChan<48) {
432  if(daqId>=mNPVPD) {
433  gMessMgr->Warning("More than expected PVPD channels read in","OS");
434  } else {
435  mPVPDRes[daqId] = tofres[0].resolution[i];
436  if(Debug()) {
437  gMessMgr->Info("","OS") << " -PVPD- daqId=" << daqId << " tdcChan=" << tdcChan << " resolution=" << mPVPDRes[daqId] << endm;
438  }
439  }
440  } else {
441  if(daqId>=mNTOFr) {
442  gMessMgr->Warning("More than expected TOFr channels read in","OS");
443  } else {
444  mTofrRes[daqId] = tofres[0].resolution[i];
445  if(Debug()) {
446  gMessMgr->Info("","OS") << " -TOFr- daqId=" << daqId << " tdcChan=" << tdcChan << " resolution=" << mTofrRes[daqId] << endm;
447  }
448  }
449  }
450  }
454  } else if(mYear5) {
455  gMessMgr->Info(" loading parameters for Run V","OS");
456 
457  // read INL table
458  St_tofr5INLtable* tofr5INLtable = static_cast<St_tofr5INLtable*>(mDbDataSet->Find("tofr5INLtable"));
459  if(!tofr5INLtable) {
460  gMessMgr->Error("unable to get tofr5 INL table parameters","OS");
461  // assert(tofr5INLtable);
462  return kStErr;
463  }
464  tofr5INLtable_st* inltable = static_cast<tofr5INLtable_st*>(tofr5INLtable->GetArray());
465  Int_t numRows = tofr5INLtable->GetNRows();
466 
467  if(Debug()) gMessMgr->Info("","OS") << " Number of rows read in: " << numRows << " for INL tables" << endm;
468 
469  const Char_t *boardName;
470  Short_t boardId;
471  Short_t tdcId;
472  Float_t INLcorr[1024];
473 
474  for (Int_t i=0;i<numRows;i++) {
475  boardName = "";
476  boardId = -1;
477  tdcId = -1;
478  for(int j=0;j<mTdcChannel;j++) {
479  INLcorr[j] = 0.0;
480  }
481 
482  boardName = (Char_t *)(inltable[i].boardID);
483  boardId = inltable[i].boardNumber;
484  tdcId = inltable[i].TDCID;
485  if(Debug())
486  gMessMgr->Info("","OS") << " name = " << boardName << " bId = " << boardId << " tdcId = " << tdcId << endm;
487  for(int j=0;j<mTdcChannel;j++) {
488  INLcorr[j] = inltable[i].INLcorrection[j];
489  if(Debug()&&j%100==0) gMessMgr->Info("","OS") << " j=" << j << " inlcorr=" << INLcorr[j] << endm;
490  mINLtable[boardId][tdcId][j] = INLcorr[j];
491  }
492 
493  }
494 
495  // read tofTotCorr table
496  St_tofTotCorr* tofTotCorr = static_cast<St_tofTotCorr*>(mDbDataSet->Find("tofTotCorr"));
497  if(!tofTotCorr) {
498  gMessMgr->Error("unable to get tofr5 TotCorr table parameters","OS");
499  // assert(tofTotCorr);
500  return kStErr;
501  }
502  tofTotCorr_st* totCorr = static_cast<tofTotCorr_st*>(tofTotCorr->GetArray());
503  numRows = tofTotCorr->GetNRows();
504 
505  if(Debug()) gMessMgr->Info("","OS") << " Number of rows read in: " << numRows << " for ToT correction" << endm;
506 
507  if(numRows!=mNTOFr5+mNPVPD) {
508  gMessMgr->Warning("","OS") << " Mis-matched number of rows in tofTotCorr table! Return! " << endm;
509  // return kStErr;
510  }
511 
512  // for (Int_t i=0;i<numRows;i++) {
513  for(Int_t i=0;i<mNTOFr5+mNPVPD;i++) {
514  short trayId = totCorr[i].trayId;
515  short moduleId = totCorr[i].moduleId;
516  short cellId = totCorr[i].cellId;
517  short tdcId = totCorr[i].tdcId;
518 
519  if(Debug()) gMessMgr->Info("","OS") << " module " << moduleId << " cell " << cellId << " tdcId " << tdcId << endm;
520  for(Int_t j=0;j<mNBinMax;j++) {
521  if(trayId==-1||trayId==-2) { // pVPD east west
522  mPVPDTotEdge[cellId-1][j] = totCorr[i].tot[j];
523  mPVPDTotCorr[cellId-1][j] = totCorr[i].corr[j];
524  } else if(trayId==93) { // TOFr5 tray
525  mTofr5TotEdge[(moduleId-1)*6+cellId-1][j] = totCorr[i].tot[j];
526  mTofr5TotCorr[(moduleId-1)*6+cellId-1][j] = totCorr[i].corr[j];
527  if(Debug()&&j%10==0) gMessMgr->Info("","OS") << " j=" << j << " tot " << mTofr5TotEdge[(moduleId-1)*6+cellId-1][j] << " corr " << mTofr5TotCorr[(moduleId-1)*6+cellId-1][j] << endm;
528  } else { // wrong trayId
529  }
530  } // end j 0->mNBinMax
531  } // end i 0->numRows
532 
533 
534  // read tofZCorr table
535  St_tofZCorr* tofZCorr = static_cast<St_tofZCorr*>(mDbDataSet->Find("tofZCorr"));
536  if(!tofZCorr) {
537  gMessMgr->Error("unable to get tofr5 ZCorr table parameters","OS");
538  // assert(tofZCorr);
539  return kStErr;
540  }
541  tofZCorr_st* zCorr = static_cast<tofZCorr_st*>(tofZCorr->GetArray());
542  numRows = tofZCorr->GetNRows();
543 
544  if(Debug()) gMessMgr->Info("","OS") << " Number of rows read in: " << numRows << " for Z correction" << endm;
545 
546  if(numRows!=mNTOFr5) { // only for TOFr5 tray
547  gMessMgr->Warning("","OS") << " Mis-matched number of rows in tofZCorr table! Return! " << endm;
548  // return kStErr;
549  }
550 
551  // for (Int_t i=0;i<numRows;i++) {
552  for (Int_t i=0;i<mNTOFr5;i++) {
553  short trayId = zCorr[i].trayId;
554  short moduleId = zCorr[i].moduleId;
555  short cellId = zCorr[i].cellId;
556 // short tdcId = zCorr[i].tdcId;
557 
558  if(Debug()) gMessMgr->Info("","OS") << " module " << moduleId << " cell " << cellId << endm;
559  for(Int_t j=0;j<mNBinMax;j++) {
560  if(trayId==93) { // TOFr5 tray
561  mTofr5ZEdge[(moduleId-1)*6+cellId-1][j] = zCorr[i].z[j];
562  mTofr5ZCorr[(moduleId-1)*6+cellId-1][j] = zCorr[i].corr[j];
563  if(Debug()&&j%10==0) gMessMgr->Info("","OS") << " j=" << j << " z " << mTofr5ZEdge[(moduleId-1)*6+cellId-1][j] << " corr " << mTofr5ZCorr[(moduleId-1)*6+cellId-1][j] << endm;
564  } else { // wrong trayId
565  }
566  } // end j 0->mNBinMax
567  } // end i 0->numRows
568 
572  } else if(mYear8) {
573  gMessMgr->Info("","OS") << " loading parameters for Run VIII" << endm;
574 
575  // read tofTotCorr table
576  St_tofTotCorr* tofTotCorr = static_cast<St_tofTotCorr*>(mDbDataSet->Find("tofTotCorr"));
577  if(!tofTotCorr) {
578  gMessMgr->Error("unable to get tof TotCorr table parameters","OS");
579  // assert(tofTotCorr);
580  return kStErr;
581  }
582  tofTotCorr_st* totCorr = static_cast<tofTotCorr_st*>(tofTotCorr->GetArray());
583  Int_t numRows = tofTotCorr->GetNRows();
584 
585  if(numRows!=mNTray8*mNTDIG+mNVPD*2) {
586  gMessMgr->Warning("","OS") << " Mis-matched number of rows in tofTotCorr table! " << endm;
587  // return kStErr;
588  }
589 
590  if(Debug()) gMessMgr->Info("","OS") << " Number of rows read in: " << numRows << " for ToT correction" << endm;
591 
592  for (Int_t i=0;i<mNTray8*mNTDIG+mNVPD*2;i++) {
593  short trayId = totCorr[i].trayId;
594  short moduleId = totCorr[i].moduleId;
595  short boardId = (moduleId-1)/4+1; // used for trays
596  short cellId = totCorr[i].cellId; // used for vpds
597 
598  int index = (trayId-mNTray-1)*mNVPD+(cellId-1); // used for vpd index
599 
600  if(Debug()) gMessMgr->Info("","OS") << " tray " << trayId << " board " << boardId << " cell " << cellId << endm;
601  for(Int_t j=0;j<mNBinMax;j++) {
602  if(trayId==mWestVpdTrayId||trayId==mEastVpdTrayId) { // upVPD east west
603  mVPDTotEdge[index][j] = totCorr[i].tot[j];
604  mVPDTotCorr[index][j] = totCorr[i].corr[j];
605  } else if(trayId>0&&trayId<=mNTray){ // trays
606  mTofTotEdge[trayId-1][boardId-1][j] = totCorr[i].tot[j];
607  mTofTotCorr[trayId-1][boardId-1][j] = totCorr[i].corr[j];
608  if(Debug()&&j%10==0) gMessMgr->Info("","OS") << " j=" << j << " tot " << mTofTotEdge[trayId-1][boardId-1][j] << " corr " << mTofTotCorr[trayId-1][boardId-1][j] << endm;
609  }
610  } // end j 0->mNBinMax
611  } // end i 0->numRows
612 
613 
614  // read tofZCorr table
615  St_tofZCorr* tofZCorr = static_cast<St_tofZCorr*>(mDbDataSet->Find("tofZCorr"));
616  if(!tofZCorr) {
617  gMessMgr->Error("unable to get tof ZCorr table parameters","OS");
618  // assert(tofZCorr);
619  return kStErr;
620  }
621  tofZCorr_st* zCorr = static_cast<tofZCorr_st*>(tofZCorr->GetArray());
622  numRows = tofZCorr->GetNRows();
623 
624  if(numRows!=mNTray8*mNTDIG) {
625  gMessMgr->Warning("","OS") << " Mis-matched number of rows in tofZCorr table! " << endm;
626  // return kStErr;
627  }
628  if(Debug()) gMessMgr->Info("","OS") << " Number of rows read in: " << numRows << " for Z correction" << endm;
629 
630  for (Int_t i=0;i<mNTray8*mNTDIG;i++) {
631  short trayId = totCorr[i].trayId;
632  short moduleId = totCorr[i].moduleId;
633  short boardId = (moduleId-1)/4+1; // used for trays
634  short cellId = totCorr[i].cellId; // used for vpds
635 
636  if(Debug()) gMessMgr->Info("","OS") << " tray " << trayId << " board " << boardId << " cell " << cellId << endm;
637  for(Int_t j=0;j<mNBinMax;j++) {
638  if(trayId>0&&trayId<=120) { // trays
639  mTofZEdge[trayId-1][boardId-1][j] = zCorr[i].z[j];
640  mTofZCorr[trayId-1][boardId-1][j] = zCorr[i].corr[j];
641  if(Debug()&&j%10==0) gMessMgr->Info("","OS") << " j=" << j << " tot " << mTofZEdge[trayId-1][boardId-1][j] << " corr " << mTofZCorr[trayId-1][boardId-1][j] << endm;
642  }
643  } // end j 0->mNBinMax
644  } // end i 0->numRows
645 
646  // read tofTOffset table
647  St_tofTOffset* tofTOffset = static_cast<St_tofTOffset*>(mDbDataSet->Find("tofTOffset"));
648  if(!tofTOffset) {
649  gMessMgr->Error("unable to get tof TOffset table parameters","OS");
650  // assert(tofTOffset);
651  return kStErr;
652  }
653  tofTOffset_st* tZero = static_cast<tofTOffset_st*>(tofTOffset->GetArray());
654  numRows = tofTOffset->GetNRows();
655 
656  if(Debug()) gMessMgr->Info("","OS") << " Number of rows read in: " << numRows << " for TOffset correction" << endm;
657 
658  if(numRows!=mNTray8) {
659  gMessMgr->Warning("","OS") << " Mis-matched number of rows in tofTOffset table! " << endm;
660  // return kStErr;
661  }
662  for (Int_t i=0;i<mNTray8;i++) {
663  short trayId = tZero[i].trayId;
664  if(Debug()) gMessMgr->Info("","OS") << " tray " << trayId << endm;
665 
666  if(trayId>0&&trayId<=mNTray) {
667  for(int j=0;j<mNTOF;j++) {
668  mTofTZero[trayId-1][j/6][j%6] = tZero[i].T0[j];
669  if(Debug()&&j%10==0) gMessMgr->Info("","OS") << " j=" << j << " T0 " << mTofTZero[trayId-1][j/6][j%6] << endm;
670  }
671  }
672  }
673 
674  // read tofPhaseOffset table
675  St_tofPhaseOffset* tofPhaseOffset = static_cast<St_tofPhaseOffset*>(mDbDataSet->Find("tofPhaseOffset"));
676  if(!tofPhaseOffset) {
677  gMessMgr->Error("unable to get tof PhaseOffset table parameters","OS");
678  // assert(tofPhaseOffset);
679  return kStErr;
680  }
681  tofPhaseOffset_st* tPhaseDiff = static_cast<tofPhaseOffset_st*>(tofPhaseOffset->GetArray());
682 
683  mPhaseOffset8 = tPhaseDiff[0].T0[0] + tPhaseDiff[0].T0[1]; // run 8, only e/w difference, to double precision
684  if(Debug()) gMessMgr->Info("","OS") << " PhaseOffset = " << mPhaseOffset8 << endm;
685 
686  /*
687  // test -- read in from data files
688  ifstream inData;
689 
690  inData.open("/star/u/dongx/TOFr/y8update/dbase/dat/pvpdCali_9050074.dat");
691  for (Int_t i=0;i<mNVPD*2;i++) {
692  int tubeId, nbin;
693  inData>>tubeId;
694  inData>>nbin; // nbin intervals, nbin+1 edges
695 
696  short trayId = (tubeId-1)/mNVPD + mNTray + 1;
697  short cellId = (tubeId-1)%mNVPD + 1; // used for vpds
698 
699  if(Debug()) gMessMgr->Info("","OS") << " tray " << trayId << " cell " << cellId << endm;
700  int index = (trayId-mNTray-1)*mNVPD+(cellId-1);
701  for(Int_t j=0;j<=nbin;j++) {
702  double xedge;
703  inData>>xedge;
704  if(trayId==mWestVpdTrayId||trayId==mEastVpdTrayId) { // upVPD east west
705  mVPDTotEdge[index][j] = xedge;
706  }
707  }
708  for(Int_t j=0;j<=nbin;j++) {
709  double ycorr;
710  inData>>ycorr;
711  if(trayId==mWestVpdTrayId||trayId==mEastVpdTrayId) { // upVPD east west
712  mVPDTotCorr[index][j] = ycorr;
713  }
714  } // end j 0->mNBinMax
715  for(Int_t j=0;j<=nbin;j++) {
716  if(Debug()&&j%10==0) gMessMgr->Info("","OS") << " j=" << j << " tot " << mVPDTotEdge[index][j] << " corr " << mVPDTotCorr[index][j] << endm;
717  } // end j 0->mNBinMax
718  } // end i 0->numRows
719  inData.close();
720 
721  inData.open("/star/u/dongx/TOFr/y8update/dbase/dat/totCali_9050074.dat");
722  for (Int_t i=0;i<mNTray8*mNTDIG;i++) {
723  int trayId, boardId, nbin;
724  inData >> trayId >> boardId;
725  inData >> nbin;
726 
727  if(Debug()) gMessMgr->Info("","OS") << " tray " << trayId << " board " << boardId << endm;
728  for(Int_t j=0;j<=nbin;j++) {
729  if (trayId>0&&trayId<=mNTray){ // trays
730  double xedge;
731  inData>>xedge;
732  mTofTotEdge[trayId-1][boardId-1][j] = xedge;
733  }
734  }
735  for(Int_t j=0;j<=nbin;j++) {
736  if (trayId>0&&trayId<=mNTray){ // trays
737  double ycorr;
738  inData>>ycorr;
739  mTofTotCorr[trayId-1][boardId-1][j] = ycorr;
740  }
741  }
742  for(Int_t j=0;j<=nbin;j++) {
743  if(Debug()&&j%10==0) gMessMgr->Info("","OS") << " j=" << j << " tot " << mTofTotEdge[trayId-1][boardId-1][j] << " corr " << mTofTotCorr[trayId-1][boardId-1][j] << endm;
744  } // end j 0->mNBinMax
745  } // end i 0->numRows
746  inData.close();
747 
748  inData.open("/star/u/dongx/TOFr/y8update/dbase/dat/zCali_9050074.dat");
749  for (Int_t i=0;i<mNTray8*mNTDIG;i++) {
750  int trayId, boardId, nbin;
751  inData >> trayId >> boardId;
752  inData >> nbin;
753 
754  if(Debug()) gMessMgr->Info("","OS") << " tray " << trayId << " board " << boardId << endm;
755  for(Int_t j=0;j<=nbin;j++) {
756  if (trayId>0&&trayId<=mNTray){ // trays
757  double xedge;
758  inData>>xedge;
759  mTofZEdge[trayId-1][boardId-1][j] = xedge;
760  }
761  }
762  for(Int_t j=0;j<=nbin;j++) {
763  if (trayId>0&&trayId<=mNTray){ // trays
764  double ycorr;
765  inData>>ycorr;
766  mTofZCorr[trayId-1][boardId-1][j] = ycorr;
767  }
768  }
769  for(Int_t j=0;j<=nbin;j++) {
770  if(Debug()&&j%10==0) gMessMgr->Info("","OS") << " j=" << j << " z " << mTofZEdge[trayId-1][boardId-1][j] << " corr " << mTofZCorr[trayId-1][boardId-1][j] << endm;
771  } // end j 0->mNBinMax
772  } // end i 0->numRows
773  inData.close();
774 
775  inData.open("/star/u/dongx/TOFr/y8update/dbase/dat/t0_9050074.dat");
776  for (Int_t i=0;i<mNTray8;i++) {
777  for (Int_t j=0;j<mNModule;j++) {
778  for (Int_t k=0;k<mNCell;k++) {
779  int trayId, moduleId, cellId;
780  inData>>trayId>>moduleId>>cellId;
781 
782  double t0;
783  inData >> t0;
784  mTofTZero[trayId-1][moduleId-1][cellId-1] = t0;
785  if(Debug()&&k==0) gMessMgr->Info("","OS") << " trayId=" << trayId << " moduleId=" << moduleId << " cellId=" << cellId << " T0 " << mTofTZero[trayId-1][moduleId-1][cellId-1] << endm;
786  }
787  }
788  }
789  inData.close();
790 
791  //mPhaseOffset8 = -33796.; // run 8, only e/w difference
792  mPhaseOffset8 = - 27945.8; //9058055
793  */
794 
795  // ========== Set Beam Line =====================
796  double x0 = 0.;
797  double y0 = 0.;
798  double dxdz = 0.;
799  double dydz = 0.;
800 
801  // Get Current Beam Line Constraint from database
802  TDataSet* dbDataSet = this->GetDataBase("Calibrations/rhic/vertexSeed");
803 
804  if (dbDataSet) {
805  vertexSeed_st* vSeed = ((St_vertexSeed*) (dbDataSet->FindObject("vertexSeed")))->GetTable();
806 
807  x0 = vSeed->x0;
808  y0 = vSeed->y0;
809  dxdz = vSeed->dxdz;
810  dydz = vSeed->dydz;
811  }
812  else {
813  LOG_INFO << "StTofCalibMaker -- No Database for beamline" << endm;
814  }
815 
816  LOG_INFO << "BeamLine Constraint: " << endm;
817  LOG_INFO << "x(z) = " << x0 << " + " << dxdz << " * z" << endm;
818  LOG_INFO << "y(z) = " << y0 << " + " << dydz << " * z" << endm;
819 
820  //**********
821  //beam line not be calibrated yet
822  //x0 shift by 0.5
823  //x0 = 0.5;
824  //**********
825  StThreeVectorD origin(x0,y0,0.0);
826  double pt = 88889999;
827  double nxy=::sqrt(dxdz*dxdz + dydz*dydz);
828  if(nxy<1.e-5){ // beam line _MUST_ be tilted
829  LOG_WARN << "StTofCalibMaker:: Beam line must be tilted!" << endm;
830  nxy=dxdz=1.e-5;
831  }
832  double p0=pt/nxy;
833  double px = p0*dxdz;
834  double py = p0*dydz;
835  double pz = p0; // approximation: nx,ny<<0
836  StThreeVectorD MomFstPt(px*GeV, py*GeV, pz*GeV);
837  if(mBeamHelix) delete mBeamHelix;
838  mBeamHelix = new StPhysicalHelixD(MomFstPt,origin,0.5*tesla,1.);
839 
840  } // end if (mYear8)
841 
842  return kStOK;
843 }
844 
845 //____________________________________________________________________________
846 Int_t StTofCalibMaker::FinishRun(int runnumber)
847 {
848 
849  if(mTofpGeom) delete mTofpGeom;
850  mTofpGeom = 0;
851 
852  if(mBeamHelix) delete mBeamHelix;
853  mBeamHelix = 0;
854 
855  return kStOK;
856 }
857 
858 //_____________________________________________________________________________
860 {
861  clearFormulars();
862  return kStOK;
863 }
864 
865 //_____________________________________________________________________________
867 {
868  if (mTofrSlewing) delete mTofrSlewing;
869  if (mTofpSlewing) delete mTofpSlewing;
870  if (mTofrZCorr) delete mTofrZCorr;
871  if (mTofpZCorr) delete mTofpZCorr;
872  if (mPVPDSlewing) delete mPVPDSlewing;
873 }
874 
875 //_____________________________________________________________________________
877 {
878  gMessMgr->Info(" StTofCalibMaker::Maker: starting ...","OS");
879 
880  mTStart = -9999.;
881  mTDiff = -9999.;
882  mVPDVtxZ = -9999.;
883  mProjVtxZ = -9999.;
884  mVPDHitPatternEast = 0;
885  mVPDHitPatternWest = 0;
886 
887  Int_t iret = kStOK;
888  if(mYear2||mYear3||mYear4){
889  iret = processEventYear2to4();
890  } else if(mYear5) {
891  iret = processEventYear5();
892  } else if(mYear8) {
893  iret = processEventYear8();
894  }
895  return kStOK;
896 }
897 
898 //____________________________________________________________________________
899 Int_t StTofCalibMaker::processEventYear2to4(){
900  mEvent = (StEvent *) GetInputDS("StEvent");
901 
902  // event selection
903  if( !mEvent || !mEvent->primaryVertex() ||
904  !mEvent->tofCollection() ||
905  !mEvent->tofCollection()->dataPresent() ||
906  (!mEvent->tofCollection()->cellsPresent() &&
907  !mEvent->tofCollection()->slatsPresent()) ) {
908  gMessMgr->Info("","OS") << "StTofCalibMaker -- nothing to do ... bye-bye" << endm;
909  return kStOK;
910  }
911 
912  StThreeVectorD vtx = mEvent->primaryVertex()->position();
913  Double_t vz = vtx.z();
914 
915  //-------------------------------------------------
916  // pVPD calibration and calculate the start timing
917  //-------------------------------------------------
918  StTofCollection *theTof = mEvent->tofCollection();
919  StSPtrVecTofData &tofData = theTof->tofData();
920  for(int i=0;i<mNPVPD;i++) {
921  mPVPDAdc[i] = tofData[42+i]->adc();
922  mPVPDTdc[i] = tofData[42+i]->tdc();
923  if(mYear3||mYear4) {
924  mPVPDAdc[i] = tofData[54+i]->adc();
925  }
926  }
927 
928  mValidStartTime = kTRUE;
929  Double_t T0 = tstart(mPVPDAdc, mPVPDTdc, vz);
930  if(T0>1000.) {
931  gMessMgr->Info(" Not a good PVPD-required event!","OS");
932  mValidStartTime = kFALSE;
933  // return kStOK;
934  }
935  gMessMgr->Info("","OS") << " pVPD start timing T0 = " << T0 << endm;
936 
937  StTofHitCollection *tofHit = new StTofHitCollection();
938 
939  //---------------------------------------
940  // Tofr calibrations cells -> hits
941  //---------------------------------------
942  StSPtrVecTofCell &tofCell = theTof->tofCells();
943  Int_t ncells = tofCell.size();
944  gMessMgr->Info("","OS") << " TOFr matched cells : " << ncells << endm;
945  for(int i=0;i<ncells;i++) {
946  StTofHit *aHit = new StTofHit;
947  StTofCell *aCell = dynamic_cast<StTofCell*>(tofCell[i]);
948  aHit->setTrayIndex(aCell->trayIndex());
949  aHit->setModuleIndex(aCell->moduleIndex());
950  aHit->setCellIndex(aCell->cellIndex());
951  aHit->setCellCollIndex(i);
952 
953  int daqId = aCell->daqIndex();
954  aHit->setDaqIndex(daqId);
955  Double_t tof = (Double_t)aCell->tdc()*mTDCWidth;
956  Double_t adc = (Double_t)aCell->adc();
957  Double_t zhit = (Double_t)aCell->zHit();
958 
959  StTrack *thisTrack = aCell->associatedTrack();
960  aHit->setAssociatedTrack(thisTrack);
961  StTrackGeometry *theTrackGeometry =
962  (mOuterGeometry)?thisTrack->outerGeometry():thisTrack->geometry();
963  Double_t L = tofPathLength(&vtx, &aCell->position(), theTrackGeometry->helix().curvature());
964  aHit->setPathLength((Float_t)L);
965 
966  if(adc<mTofrADCMin[daqId]||adc>=mTofrADCMax[daqId]) {
967  delete aHit;
968  continue;
969  }
970  if(mValidCalibPar&&mValidStartTime) {
971  Double_t tofcorr = tofrAllCorr(tof, T0, adc, zhit, daqId);
972  aHit->setTimeOfFlight((Float_t)tofcorr);
973 
974  Double_t beta = L/(tofcorr*(C_C_LIGHT/1.e9));
975  aHit->setBeta((Float_t)beta);
976 
977  const StThreeVectorF momentum = thisTrack->geometry()->momentum();
978  Double_t ptot = momentum.mag();
979  Double_t tofe = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_ELECTRON*M_ELECTRON)/ptot;
980  Double_t tofpi = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_PION_PLUS*M_PION_PLUS)/ptot;
981  Double_t tofk = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_KAON_PLUS*M_KAON_PLUS)/ptot;
982  Double_t tofp = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_PROTON*M_PROTON)/ptot;
983 
984  aHit->setTofExpectedAsElectron((Float_t)tofe);
985  aHit->setTofExpectedAsPion((Float_t)tofpi);
986  aHit->setTofExpectedAsKaon((Float_t)tofk);
987  aHit->setTofExpectedAsProton((Float_t)tofp);
988 
989  float sigmae = 999.;
990  float sigmapi = 999.;
991  float sigmak = 999.;
992  float sigmap = 999.;
993  float res = mTofrRes[daqId];
994  if(fabs(res)>1.e-5) {
995  sigmae = (Float_t)((tofcorr-tofe)/res);
996  sigmapi = (Float_t)((tofcorr-tofpi)/res);
997  sigmak = (Float_t)((tofcorr-tofk)/res);
998  sigmap = (Float_t)((tofcorr-tofp)/res);
999  }
1000  aHit->setSigmaElectron(sigmae);
1001  aHit->setSigmaPion(sigmapi);
1002  aHit->setSigmaKaon(sigmak);
1003  aHit->setSigmaProton(sigmap);
1004 
1005  } else {
1006  aHit->setTimeOfFlight(9999.);
1007  aHit->setBeta(9999.);
1008  }
1009 
1010  tofHit->push_back(aHit);
1011  }
1012 
1013  //------------------------------------
1014  // Tofp calibrations slats -> hits
1015  // -----------------------------------
1016  StSPtrVecTofSlat &tofSlat = theTof->tofSlats();
1017  Int_t nslats = tofSlat.size();
1018  for(int i=0;i<nslats;i++) {
1019  StTofHit *aHit = new StTofHit;
1020  StTofSlat *aSlat = tofSlat[i];
1021  aHit->setTrayIndex(0);
1022  aHit->setModuleIndex(0);
1023  //
1024  // slatIndex() fake slatId, = daqId+1 actually -- in TofpMatchMaker
1025  //
1026  int daqId = aSlat->slatIndex()-1;
1027  int slatId = mTofpGeom->daqToSlatId(daqId);
1028 
1029  aHit->setCellIndex(slatId);
1030  aHit->setCellCollIndex(i);
1031  aHit->setDaqIndex(daqId);
1032 
1033  Double_t tof = (Double_t)aSlat->tdc()*mTDCWidth;
1034  Double_t adc = (Double_t)aSlat->adc();
1035  Double_t zhit = (Double_t)aSlat->zHit();
1036 
1037  StTrack *thisTrack = aSlat->associatedTrack();
1038  aHit->setAssociatedTrack(thisTrack);
1039  StTrackGeometry *theTrackGeometry =
1040  (mOuterGeometry)?thisTrack->outerGeometry():thisTrack->geometry();
1041  Double_t L = tofPathLength(&vtx, &aSlat->position(), theTrackGeometry->helix().curvature());
1042  aHit->setPathLength((Float_t)L);
1043 
1044  // if(adc<mTofpADCMin[daqId]||adc>=mTofpADCMax[daqId]) {
1045  if(adc<mTofpADCMin[daqId]) {
1046  delete aHit;
1047  continue;
1048  }
1049  if(mValidCalibPar&&mValidStartTime) {
1050  Double_t tofcorr = tofpAllCorr(tof, T0, adc, zhit, daqId);
1051  aHit->setTimeOfFlight((Float_t)tofcorr);
1052 
1053  Double_t beta = L/(tofcorr*(C_C_LIGHT/1.e9));
1054  aHit->setBeta((Float_t)beta);
1055 
1056  const StThreeVectorF momentum = thisTrack->geometry()->momentum();
1057  Double_t ptot = momentum.mag();
1058  Double_t tofe = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_ELECTRON*M_ELECTRON)/ptot;
1059  Double_t tofpi = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_PION_PLUS*M_PION_PLUS)/ptot;
1060  Double_t tofk = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_KAON_PLUS*M_KAON_PLUS)/ptot;
1061  Double_t tofp = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_PROTON*M_PROTON)/ptot;
1062  aHit->setTofExpectedAsElectron((Float_t)tofe);
1063  aHit->setTofExpectedAsPion((Float_t)tofpi);
1064  aHit->setTofExpectedAsKaon((Float_t)tofk);
1065  aHit->setTofExpectedAsProton((Float_t)tofp);
1066 
1067  float sigmae = 999.;
1068  float sigmapi = 999.;
1069  float sigmak = 999.;
1070  float sigmap = 999.;
1071  float res = mTofpRes[daqId];
1072  if(fabs(res)>1.e-5) {
1073  sigmae = (Float_t)((tofcorr-tofe)/res);
1074  sigmapi = (Float_t)((tofcorr-tofpi)/res);
1075  sigmak = (Float_t)((tofcorr-tofk)/res);
1076  sigmap = (Float_t)((tofcorr-tofp)/res);
1077  }
1078  aHit->setSigmaElectron(sigmae);
1079  aHit->setSigmaPion(sigmapi);
1080  aHit->setSigmaKaon(sigmak);
1081  aHit->setSigmaProton(sigmap);
1082 
1083  } else {
1084  aHit->setTimeOfFlight(9999.);
1085  aHit->setBeta(9999.);
1086  }
1087  tofHit->push_back(aHit);
1088  }
1089 
1090 
1091  // store the hit collection in StEvent
1092  for (size_t j=0;j<tofHit->size();j++){
1093  theTof->addHit(tofHit->getHit(j));
1094  if (Debug())
1095  gMessMgr->Info("","OS") << "storing " << j << " " << " tray:"
1096  << tofHit->getHit(j)->trayIndex() << " module:"
1097  << tofHit->getHit(j)->moduleIndex() << " cell "
1098  << tofHit->getHit(j)->cellIndex() << endm;
1099  }
1100 
1101  // add the tof pidtraits
1102  StSPtrVecTofHit& tofHitVec = theTof->tofHits();
1103  StSPtrVecTrackNode& nodes = mEvent->trackNodes();
1104  for (unsigned int iNode=0; iNode<nodes.size(); iNode++){
1105  StTrack *theTrack = nodes[iNode]->track(primary);
1106  if(!theTrack) continue;
1107  Int_t trkId = theTrack->key();
1108  for (size_t j=0;j<tofHitVec.size();j++){
1109  StTrack *aTrack = tofHitVec[j]->associatedTrack();
1110  if(!aTrack) continue;
1111  if(aTrack->key()!=trkId) continue;
1112  StTofPidTraits* pidTof = new StTofPidTraits(tofHitVec[j]->trayIndex(), tofHitVec[j]->moduleIndex(), tofHitVec[j]->cellIndex(), tofHitVec[j]->timeOfFlight(), tofHitVec[j]->pathLength(), tofHitVec[j]->beta());
1113  pidTof->setSigmaElectron(tofHitVec[j]->sigmaElectron());
1114  pidTof->setSigmaPion(tofHitVec[j]->sigmaPion());
1115  pidTof->setSigmaKaon(tofHitVec[j]->sigmaKaon());
1116  pidTof->setSigmaProton(tofHitVec[j]->sigmaProton());
1117 
1118  theTrack->addPidTraits(pidTof);
1119  }
1120  }
1121 
1122  delete tofHit;
1123 
1124  return kStOK;
1125 }
1126 
1127 //____________________________________________________________________________
1128 Int_t StTofCalibMaker::processEventYear5(){
1129 
1130  mEvent = (StEvent *) GetInputDS("StEvent");
1131 
1132  // event selection
1133  if( !mEvent || !mEvent->primaryVertex() ||
1134  !mEvent->tofCollection() ||
1135  (!mEvent->tofCollection()->dataPresent() &&
1136  !mEvent->tofCollection()->rawdataPresent()) ||
1137  (!mEvent->tofCollection()->cellsPresent()) ) {
1138  gMessMgr->Info("","OS") << "StTofCalibMaker -- nothing to do ... bye-bye" << endm;
1139  return kStOK;
1140  }
1141 
1142  StThreeVectorD vtx = mEvent->primaryVertex()->position();
1143  Double_t vz = vtx.z();
1144 
1145  //-------------------------------------------------
1146  // pVPD calibration and calculate the start timing
1147  //-------------------------------------------------
1148  StTofCollection *theTof = mEvent->tofCollection();
1149  mSortTofRawData = new StSortTofRawData(theTof);
1150  IntVec validchannel = mSortTofRawData->GetValidChannel();
1151 
1152  int used[mNPVPD] = {0,0,0,0,0,0};
1153  int channum[mNPVPD] ={-1,-1,-1,-1,-1,-1};
1154 
1155  for(unsigned int ich=0;ich<validchannel.size();ich++){
1156  int chan = validchannel[ich];
1157  if(chan<mNTOFr5) continue; // need only pvpd
1158  int ichan = chan - mNTOFr5;
1159  if(ichan<0||ichan>=mNPVPD) continue;
1160  if(used[ichan]>0) continue; // skip multi hits
1161  used[ichan]++;
1162  // leading edge
1163  int tmptdc = (mSortTofRawData->GetLeadingTdc(chan))[0];
1164  // do inl correction
1165  int bin = int(tmptdc)&0x03ff;
1166  float pvpdletdc = tmptdc + GetINLcorr(4,chan,bin);
1167  mPVPDLeTime[ichan] = pvpdletdc * VHRBIN2PS/1000.;
1168  // Trailing edge
1169  tmptdc = (mSortTofRawData->GetTrailingTdc(chan))[0];
1170  // do inl correction
1171  bin = int(tmptdc)&0x0ff;
1172  float pvpdtetdc = tmptdc + GetINLcorr(5,chan,bin);
1173  float tetime = pvpdtetdc * HRBIN2PS/1000.;
1174  mPVPDTot[ichan] = tetime - mPVPDLeTime[ichan];
1175  channum[ichan]=ichan;
1176  }
1177 
1178  // sort out the piled-up hit - remove those hits not from this event
1179  int nPVPDFired = 0;
1180  for(int i=0;i<mNPVPD;i++) {
1181  if(channum[i]!=i) {
1182  mPVPDLeTime[i] = 0.;
1183  mPVPDTot[i] = 0.;
1184  } else {
1185  nPVPDFired++;
1186  }
1187  }
1188  for(int i=0;i<mNPVPD;i++) {
1189  if(channum[i]!=i) continue;
1190  int n0 = 0;
1191  for(int j=0;j<mNPVPD;j++) {
1192  if(channum[j]!=j) continue;
1193  float dt = mPVPDLeTime[j] - mPVPDLeTime[i];
1194  if(fabs(dt)<200.) n0++;
1195  }
1196 
1197  if(n0>=nPVPDFired/2) { // OK
1198  } else { // not from this event
1199  mPVPDLeTime[i] = 0.;
1200  mPVPDTot[i] = 0.;
1201  }
1202  }
1203 
1204  mValidStartTime = kTRUE;
1205  Double_t T0 = tstart5(mPVPDTot, mPVPDLeTime, vz);
1206  if(T0<0.) {
1207  gMessMgr->Info(" Not a good PVPD-required event!","OS");
1208  mValidStartTime = kFALSE;
1209  // return kStOK;
1210  }
1211  gMessMgr->Info("","OS") << " pVPD start timing T0 = " << T0 << endm;
1212 
1213  StTofHitCollection *tofHit = new StTofHitCollection();
1214 
1215  //---------------------------------------
1216  // Tofr calibrations cells -> hits
1217  //---------------------------------------
1218  StSPtrVecTofCell &tofCell = theTof->tofCells();
1219  Int_t ncells = tofCell.size();
1220  gMessMgr->Info("","OS") << " TOFr matched cells : " << ncells << endm;
1221  for(int i=0;i<ncells;i++) {
1222  StTofHit *aHit = new StTofHit;
1223  StTofCell *aCell = dynamic_cast<StTofCell*>(tofCell[i]);
1224  aHit->setTrayIndex(aCell->trayIndex());
1225  aHit->setModuleIndex(aCell->moduleIndex());
1226  aHit->setCellIndex(aCell->cellIndex());
1227  aHit->setCellCollIndex(i);
1228 
1229  int daqId = aCell->daqIndex();
1230  aHit->setDaqIndex(daqId);
1231 
1232 
1233 
1234  Int_t letdc = (Int_t)aCell->tdc();
1235  int bin = int(letdc)&0x03ff;
1236  double tmptdc_f = letdc + GetINLcorr(4, daqId, bin);
1237  double letime = tmptdc_f * VHRBIN2PS / 1000.; // ns
1238 
1239  Int_t tetdc = (Int_t)aCell->adc();
1240  bin = int(tetdc)&0x0ff;
1241  tmptdc_f = tetdc + GetINLcorr(5, daqId, bin);
1242  double tetime = tmptdc_f * HRBIN2PS / 1000.; // ns
1243 
1244  double tot = tetime - letime; // ns
1245  double tof = letime - T0;
1246  Double_t zhit = (Double_t)aCell->zHit();
1247 
1248  StTrack *thisTrack = aCell->associatedTrack();
1249  aHit->setAssociatedTrack(thisTrack);
1250  StTrackGeometry *theTrackGeometry =
1251  (mOuterGeometry)?thisTrack->outerGeometry():thisTrack->geometry();
1252  Double_t L = tofPathLength(&vtx, &aCell->position(), theTrackGeometry->helix().curvature());
1253  aHit->setPathLength((Float_t)L);
1254 
1255  if(fabs(tof)>200.) { // +/-200 ns window
1256  gMessMgr->Info("","OS") << " the hit is not from this event!" << endm;
1257  delete aHit;
1258  continue;
1259  }
1260  if(mValidCalibPar&&mValidStartTime) {
1261  int moduleChan = (aCell->moduleIndex()-1)*6 + (aCell->cellIndex()-1);
1262  Double_t tofcorr = tofr5AllCorr(tof, tot, zhit, moduleChan);
1263  aHit->setTimeOfFlight((Float_t)tofcorr);
1264 
1265  Double_t beta = L/(tofcorr*(C_C_LIGHT/1.e9));
1266  aHit->setBeta((Float_t)beta);
1267 
1268  const StThreeVectorF momentum = thisTrack->geometry()->momentum();
1269  Double_t ptot = momentum.mag();
1270  Double_t tofe = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_ELECTRON*M_ELECTRON)/ptot;
1271  Double_t tofpi = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_PION_PLUS*M_PION_PLUS)/ptot;
1272  Double_t tofk = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_KAON_PLUS*M_KAON_PLUS)/ptot;
1273  Double_t tofp = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_PROTON*M_PROTON)/ptot;
1274 
1275  aHit->setTofExpectedAsElectron((Float_t)tofe);
1276  aHit->setTofExpectedAsPion((Float_t)tofpi);
1277  aHit->setTofExpectedAsKaon((Float_t)tofk);
1278  aHit->setTofExpectedAsProton((Float_t)tofp);
1279 
1280  float sigmae = 999.;
1281  float sigmapi = 999.;
1282  float sigmak = 999.;
1283  float sigmap = 999.;
1284  float res = 0.095; // 95 ps by default
1285  if(fabs(res)>1.e-5) {
1286  sigmae = (Float_t)((tofcorr-tofe)/res);
1287  sigmapi = (Float_t)((tofcorr-tofpi)/res);
1288  sigmak = (Float_t)((tofcorr-tofk)/res);
1289  sigmap = (Float_t)((tofcorr-tofp)/res);
1290  }
1291  aHit->setSigmaElectron(sigmae);
1292  aHit->setSigmaPion(sigmapi);
1293  aHit->setSigmaKaon(sigmak);
1294  aHit->setSigmaProton(sigmap);
1295 
1296  } else {
1297  aHit->setTimeOfFlight(9999.);
1298  aHit->setBeta(9999.);
1299  }
1300 
1301  tofHit->push_back(aHit);
1302  }
1303 
1304 
1305  // store the hit collection in StEvent
1306  for (size_t j=0;j<tofHit->size();j++){
1307  theTof->addHit(tofHit->getHit(j));
1308  if (Debug())
1309  gMessMgr->Info("","OS") << "storing " << j << " " << " tray:"
1310  << tofHit->getHit(j)->trayIndex() << " module:"
1311  << tofHit->getHit(j)->moduleIndex() << " cell:"
1312  << tofHit->getHit(j)->cellIndex() << endm;
1313  }
1314 
1315  // add the tof pidtraits
1316  StSPtrVecTofHit& tofHitVec = theTof->tofHits();
1317  StSPtrVecTrackNode& nodes = mEvent->trackNodes();
1318  for (unsigned int iNode=0; iNode<nodes.size(); iNode++){
1319  StTrack *theTrack = nodes[iNode]->track(primary);
1320  if(!theTrack) continue;
1321  Int_t trkId = theTrack->key();
1322  for (size_t j=0;j<tofHitVec.size();j++){
1323  StTrack *aTrack = tofHitVec[j]->associatedTrack();
1324  if(!aTrack) continue;
1325  if(aTrack->key()!=trkId) continue;
1326  StTofPidTraits* pidTof = new StTofPidTraits(tofHitVec[j]->trayIndex(), tofHitVec[j]->moduleIndex(), tofHitVec[j]->cellIndex(), tofHitVec[j]->timeOfFlight(), tofHitVec[j]->pathLength(), tofHitVec[j]->beta());
1327  pidTof->setSigmaElectron(tofHitVec[j]->sigmaElectron());
1328  pidTof->setSigmaPion(tofHitVec[j]->sigmaPion());
1329  pidTof->setSigmaKaon(tofHitVec[j]->sigmaKaon());
1330  pidTof->setSigmaProton(tofHitVec[j]->sigmaProton());
1331 
1332  theTrack->addPidTraits(pidTof);
1333  }
1334  }
1335 
1336  delete tofHit;
1337 
1338  return kStOK;
1339 }
1340 
1341 //____________________________________________________________________________
1343 
1344  mEvent = (StEvent *) GetInputDS("StEvent");
1345 
1346  // event selection // no primary vertex in run8
1347  if( !mEvent ||
1348  !mEvent->tofCollection() ||
1349  (!mEvent->tofCollection()->dataPresent() &&
1350  !mEvent->tofCollection()->rawdataPresent()) ||
1351  (!mEvent->tofCollection()->cellsPresent()) ) {
1352  gMessMgr->Info("","OS") << "StTofCalibMaker -- nothing to do ... bye-bye" << endm;
1353  return kStOK;
1354  }
1355 
1356  StTofCollection *theTof = mEvent->tofCollection();
1357  StSPtrVecTofCell &tofCell = theTof->tofCells();
1358  Int_t ncells = tofCell.size();
1359  gMessMgr->Info("","OS") << " TOFr matched cells + upVPD fired tubes : " << ncells << endm;
1360 
1361 
1362  //added by Zebo, projVtxZ from track cloest to beam line
1363  float dcaRmin = 9999;
1364  for(int i=0;i<ncells;i++) {
1365  StTofHit *aHit = new StTofHit;
1366  StTofCell *aCell = dynamic_cast<StTofCell*>(tofCell[i]);
1367  int trayId = aCell->trayIndex();
1368  if(trayId<=0 || trayId>120) continue;
1369  StTrack *thisTrack = aCell->associatedTrack();
1370  aHit->setAssociatedTrack(thisTrack);
1371  StTrackGeometry *theTrackGeometry = thisTrack->geometry();
1372 
1373  StThreeVectorD tofPos = theTrackGeometry->helix().at(theTrackGeometry->helix().pathLengths(*mBeamHelix).first);
1374  LOG_INFO<<"tofPos(x,y,z)= "<<tofPos.x()<<" "<<tofPos.y()<<" "<<tofPos.z()<<endm;
1375  StThreeVectorD dcatof = tofPos - mBeamHelix->at(theTrackGeometry->helix().pathLengths(*mBeamHelix).second);
1376  LOG_INFO<<"dcatof(x,y)= "<<dcatof.x()<<" "<<dcatof.y()<<endm;
1377 
1378  if(dcaRmin>dcatof.perp()) {
1379  mProjVtxZ = tofPos.z();
1380  dcaRmin = dcatof.perp();
1381  }
1382  delete aHit;
1383  }
1384  //end
1385 
1386  //-------------------------------------------------
1387  // pVPD calibration and calculate the start timing
1388  //-------------------------------------------------
1389  for(int i=0;i<2*mNVPD;i++) {
1390  mVPDLeTime[i] = -999.;
1391  mVPDTot[i] = -999.;
1392  }
1393  for(int i=0;i<ncells;i++) {
1394  StTofCell *aCell = dynamic_cast<StTofCell*>(tofCell[i]);
1395  if(!aCell) continue;
1396  int trayId = aCell->trayIndex();
1397 
1398  if(trayId!=mWestVpdTrayId && trayId!=mEastVpdTrayId) continue;
1399 
1400  int tubeId = aCell->cellIndex();
1401  mVPDLeTime[(trayId-121)*mNVPD+(tubeId-1)] = aCell->leadingEdgeTime();
1402  mVPDTot[(trayId-121)*mNVPD+(tubeId-1)] = aCell->tot();
1403 
1404  // T0 is calculated one-by-one according to different vz.
1405 // Double_t T0 = tstart5(mVPDTot, mVPDLeTime, vz);
1406 // if(T0<0.) {
1407 // gMessMgr->Info(" Not a good PVPD-required event!","OS");
1408 // mValidStartTime = kFALSE;
1409 // // return kStOK;
1410 // }
1411  }
1412  tsum8(mVPDTot, mVPDLeTime);
1413  mTStart = tstart8(mProjVtxZ);
1414 
1415  gMessMgr->Info("","OS") << " NWest = " << mNWest << " NEast = " << mNEast << " TdcSum West = " << mTSumWest << " East = " << mTSumEast << endm;
1416  if(mTStart<-1000.) {
1417  gMessMgr->Info("","OS") << " mTStart not available!" << endm;
1418  mValidStartTime = kFALSE;
1419  } else {
1420  mValidStartTime = kTRUE;
1421  }
1422  gMessMgr->Info("","OS") << " mValidCalibPar = " << mValidCalibPar << " mValidStartTime = " << mValidStartTime << endm;
1423 
1425  theTof->setVpdEast(mVPDHitPatternEast);
1426  theTof->setVpdWest(mVPDHitPatternWest);
1427  theTof->setTdiff(mTDiff);
1428  theTof->setVzVpd(mVPDVtxZ);
1429  theTof->setTstart(mTStart);
1430 
1431  gMessMgr->Info("","OS") << " TofCollection: NWest = " << theTof->numberOfVpdWest() << " NEast = " << theTof->numberOfVpdEast() << endm;
1432  gMessMgr->Info("","OS") << "Tdiff = " << mTDiff <<" vpd vz = " << mVPDVtxZ << " proj vz = " << mProjVtxZ<<endm;
1433 
1434  StTofHitCollection *tofHit = new StTofHitCollection();
1435 
1436  //---------------------------------------
1437  // Tofr calibrations cells -> hits
1438  //---------------------------------------
1439 
1440  for(int i=0;i<ncells;i++) {
1441  StTofHit *aHit = new StTofHit;
1442  StTofCell *aCell = dynamic_cast<StTofCell*>(tofCell[i]);
1443  int trayId = aCell->trayIndex();
1444  if(trayId<=0 || trayId>120) continue;
1445 
1446  aHit->setTrayIndex(aCell->trayIndex());
1447  aHit->setModuleIndex(aCell->moduleIndex());
1448  aHit->setCellIndex(aCell->cellIndex());
1449  aHit->setCellCollIndex(i);
1450 
1451  int daqId = aCell->daqIndex();
1452  aHit->setDaqIndex(daqId);
1453 
1454  StTrack *thisTrack = aCell->associatedTrack();
1455  aHit->setAssociatedTrack(thisTrack);
1456  StTrackGeometry *theTrackGeometry = thisTrack->geometry();
1457 
1458  StThreeVectorD tofPos = theTrackGeometry->helix().at(theTrackGeometry->helix().pathLengths(*mBeamHelix).first);
1459  StThreeVectorD dcatof = tofPos - mBeamHelix->at(theTrackGeometry->helix().pathLengths(*mBeamHelix).second);
1460 
1461  Double_t L = tofPathLength(&tofPos, &aCell->position(), theTrackGeometry->helix().curvature());
1462  aHit->setPathLength((Float_t)L);
1463 
1464  //double mTStart = tstart8(dcatof.z());
1465  gMessMgr->Info("","OS") << " p(x,y,z) = "<<theTrackGeometry->momentum().x()<<" "<<theTrackGeometry->momentum().y()<<" "<<theTrackGeometry->momentum().z()<<endm;
1466  gMessMgr->Info("","OS") << " Project Z = " << tofPos.z() << " mTStart = " << mTStart << endm;
1467 
1468  double tot = aCell->tot(); // ns
1469  double tdc = aCell->leadingEdgeTime();
1470  tdc -= mPhaseOffset8;
1471  while(tdc>TMAX) tdc -= TMAX;
1472  double tof = tdc - mTStart;
1473  Double_t zhit = (Double_t)aCell->zHit();
1474 
1475 
1476  if(mValidCalibPar&&mValidStartTime) {
1477  int moduleChan = (aCell->moduleIndex()-1)*6 + (aCell->cellIndex()-1);
1478  Double_t tofcorr = tofr8AllCorr(tof, tot, zhit, trayId, moduleChan);
1479  aHit->setTimeOfFlight((Float_t)tofcorr);
1480 
1481  Double_t beta = L/(tofcorr*(C_C_LIGHT/1.e9));
1482  aHit->setBeta((Float_t)beta);
1483 
1484  const StThreeVectorF momentum = thisTrack->geometry()->momentum();
1485  Double_t ptot = momentum.mag();
1486  Double_t tofe = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_ELECTRON*M_ELECTRON)/ptot;
1487  Double_t tofpi = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_PION_PLUS*M_PION_PLUS)/ptot;
1488  Double_t tofk = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_KAON_PLUS*M_KAON_PLUS)/ptot;
1489  Double_t tofp = L/(C_C_LIGHT/1.e9)*sqrt(ptot*ptot+M_PROTON*M_PROTON)/ptot;
1490 
1491  aHit->setTofExpectedAsElectron((Float_t)tofe);
1492  aHit->setTofExpectedAsPion((Float_t)tofpi);
1493  aHit->setTofExpectedAsKaon((Float_t)tofk);
1494  aHit->setTofExpectedAsProton((Float_t)tofp);
1495 
1496  float sigmae = 999.;
1497  float sigmapi = 999.;
1498  float sigmak = 999.;
1499  float sigmap = 999.;
1500  float res = 0.085; // 85 ps by default
1501  if(fabs(res)>1.e-5) {
1502  sigmae = (Float_t)((tofcorr-tofe)/res);
1503  sigmapi = (Float_t)((tofcorr-tofpi)/res);
1504  sigmak = (Float_t)((tofcorr-tofk)/res);
1505  sigmap = (Float_t)((tofcorr-tofp)/res);
1506  }
1507  aHit->setSigmaElectron(sigmae);
1508  aHit->setSigmaPion(sigmapi);
1509  aHit->setSigmaKaon(sigmak);
1510  aHit->setSigmaProton(sigmap);
1511 
1512  } else {
1513  aHit->setTimeOfFlight(9999.);
1514  aHit->setBeta(9999.);
1515  }
1516 
1517  tofHit->push_back(aHit);
1518  }
1519 
1520 
1521  // store the hit collection in StEvent
1522  for (size_t j=0;j<tofHit->size();j++){
1523  theTof->addHit(tofHit->getHit(j));
1524  if (Debug())
1525  gMessMgr->Info("","OS") << "storing " << j << " " << " tray:"
1526  << tofHit->getHit(j)->trayIndex() << " module:"
1527  << tofHit->getHit(j)->moduleIndex() << " cell:"
1528  << tofHit->getHit(j)->cellIndex() << " TOF:"
1529  << tofHit->getHit(j)->timeOfFlight() << " beta:"
1530  << tofHit->getHit(j)->beta() << endm;
1531  }
1532 
1533  // add the tof pidtraits
1534  StSPtrVecTofHit& tofHitVec = theTof->tofHits();
1535  StSPtrVecTrackNode& nodes = mEvent->trackNodes();
1536  for (unsigned int iNode=0; iNode<nodes.size(); iNode++){
1537  StTrack *theTrack = nodes[iNode]->track(primary);
1538  if(!theTrack) continue;
1539  Int_t trkId = theTrack->key();
1540  for (size_t j=0;j<tofHitVec.size();j++){
1541  StTrack *aTrack = tofHitVec[j]->associatedTrack();
1542  if(!aTrack) continue;
1543  if(aTrack->key()!=trkId) continue;
1544  StTofPidTraits* pidTof = new StTofPidTraits(tofHitVec[j]->trayIndex(), tofHitVec[j]->moduleIndex(), tofHitVec[j]->cellIndex(), tofHitVec[j]->timeOfFlight(), tofHitVec[j]->pathLength(), tofHitVec[j]->beta());
1545  pidTof->setSigmaElectron(tofHitVec[j]->sigmaElectron());
1546  pidTof->setSigmaPion(tofHitVec[j]->sigmaPion());
1547  pidTof->setSigmaKaon(tofHitVec[j]->sigmaKaon());
1548  pidTof->setSigmaProton(tofHitVec[j]->sigmaProton());
1549 
1550  theTrack->addPidTraits(pidTof);
1551  }
1552  }
1553 
1554  delete tofHit;
1555 
1556  return kStOK;
1557 }
1558 
1559 //_____________________________________________________________________________
1560 Double_t StTofCalibMaker::tofrT0Corr(const Double_t tof, const Double_t Tstart, const Int_t iDaqChan)
1561 {
1562  return tof - Tstart - mTofrT0[iDaqChan];
1563 }
1564 
1565 //_____________________________________________________________________________
1566 Double_t StTofCalibMaker::tofrSlewingCorr(const Double_t tof, const Double_t adc, const Int_t iDaqChan)
1567 {
1568  Double_t par[mNPar];
1569  for(int i=0;i<mNPar;i++) {
1570  par[i] = mTofrTAPar[iDaqChan*mNPar+i];
1571  }
1572 
1573  mTofrSlewing->SetParameters(par);
1574  return tof - mTofrSlewing->Eval(adc);
1575 
1576 }
1577 
1578 //_____________________________________________________________________________
1579 Double_t StTofCalibMaker::tofrZCorr(const Double_t tof, const Double_t zhit)
1580 {
1581  mTofrZCorr->SetParameters(mTofrZPar);
1582  return tof - mTofrZCorr->Eval(zhit);
1583 }
1584 
1585 //_____________________________________________________________________________
1586 Double_t StTofCalibMaker::tofrAllCorr(const Double_t tof, const Double_t T0, const Double_t adc, const Double_t z, const Int_t iDaqChan)
1587 {
1588  gMessMgr->Info("","OS") << "\nStTofCalibMaker::tofrAllCorr: Tofr calibrating...\n"
1589  << "\tDoing Calibration in TOFr Channel " << iDaqChan
1590  << "\n\tinput tof = " << tof << " Tstart = "
1591  << T0 << " ADC = " << adc << " Zlocal = " << z
1592  << "\n" << endm;
1593 
1594  Double_t tofcorr = tofrT0Corr( tof, T0, iDaqChan );
1595  gMessMgr->Info("","OS") << " T0 corr: tofcorr = " << tofcorr << endm;
1596  Double_t tofcorr2 = tofrSlewingCorr( tofcorr, adc, iDaqChan );
1597  gMessMgr->Info("","OS") << " Slewing corr: tofcorr2 = " << tofcorr2 << endm;
1598  Double_t tofcorr3 = tofrZCorr( tofcorr2, z);
1599  gMessMgr->Info("","OS") << " Z position corr: tofcorr3 = " << tofcorr3 << endm;
1600  return tofcorr3;
1601 }
1602 
1603 //_____________________________________________________________________________
1604 Double_t StTofCalibMaker::tofpT0Corr(const Double_t tof, const Double_t Tstart, const Int_t iDaqChan)
1605 {
1606  return tof - Tstart - mTofpT0[iDaqChan];
1607 }
1608 
1609 //_____________________________________________________________________________
1610 Double_t StTofCalibMaker::tofpSlewingCorr(const Double_t tof, const Double_t adc, const Int_t iDaqChan)
1611 {
1612  Double_t par[mNPar];
1613  for(int i=0;i<mNPar;i++) {
1614  par[i] = mTofpTAPar[iDaqChan*mNPar+i];
1615  }
1616 
1617  mTofpSlewing->SetParameters(par);
1618  return tof - mTofpSlewing->Eval(adc);
1619 
1620 }
1621 
1622 //_____________________________________________________________________________
1623 Double_t StTofCalibMaker::tofpZCorr(const Double_t tof, const Double_t zhit)
1624 {
1625  mTofpZCorr->SetParameters(mTofpZPar);
1626  return tof - mTofpZCorr->Eval(zhit);
1627 }
1628 
1629 //_____________________________________________________________________________
1630 Double_t StTofCalibMaker::tofpAllCorr(const Double_t tof, const Double_t T0, const Double_t adc, const Double_t z, const Int_t iDaqChan)
1631 {
1632  gMessMgr->Info("","OS") << "\nStTofCalibMaker::tofpAllCorr: Tofp calibrating...\n"
1633  << "\tDoing Calibration in TOFp Channel " << iDaqChan
1634  << "\n\tinput tof = " << tof << " Tstart = "
1635  << T0 << " ADC = " << adc << " Zlocal = " << z
1636  << "\n" << endm;
1637 
1638  Double_t tofcorr = tofpT0Corr( tof, T0, iDaqChan );
1639  gMessMgr->Info("","OS") << " T0 corr: tofcorr = " << tofcorr << endm;
1640  Double_t tofcorr2 = tofpSlewingCorr( tofcorr, adc, iDaqChan );
1641  gMessMgr->Info("","OS") << " Slewing corr: tofcorr2 = " << tofcorr2 << endm;
1642  Double_t tofcorr3 = tofpZCorr( tofcorr2, z);
1643  gMessMgr->Info("","OS") << " Z position corr: tofcorr3 = " << tofcorr3 << endm;
1644  return tofcorr3;
1645 }
1646 
1647 //_____________________________________________________________________________
1648 Double_t StTofCalibMaker::tofr5AllCorr(const Double_t tof, const Double_t tot, const Double_t z, const Int_t iModuleChan)
1649 {
1650  int module = iModuleChan/6 + 1;
1651  int cell = iModuleChan%6 + 1;
1652  gMessMgr->Info("","OS") << "\nStTofCalibMaker::tofr5AllCorr: Tofr5 calibrating...\n"
1653  << "\tDoing Calibration in TOFr5 Module " << module << " Cell " << cell
1654  << "\n\tinput tof = " << tof
1655  << " TOT = " << tot << " Zlocal = " << z
1656  << "\n" << endm;
1657 
1658  Double_t tofcorr = tof;
1659 
1660  int iTotBin = -1;
1661  for(int i=0;i<mNBinMax-1;i++) {
1662  if(tot>=mTofr5TotEdge[iModuleChan][i] && tot<mTofr5TotEdge[iModuleChan][i+1]) {
1663  iTotBin = i;
1664  break;
1665  }
1666  }
1667  if(iTotBin>=0&&iTotBin<mNBinMax) {
1668  tofcorr -= mTofr5TotCorr[iModuleChan][iTotBin];
1669  } else {
1670  tofcorr = -9999.;
1671  }
1672 
1673  int iZBin = -1;
1674  for(int i=0;i<mNBinMax-1;i++) {
1675  if(z>=mTofr5ZEdge[iModuleChan][i] && z<mTofr5ZEdge[iModuleChan][i+1]) {
1676  iZBin = i;
1677  break;
1678  }
1679  }
1680  if(iZBin>=0&&iZBin<mNBinMax) {
1681  tofcorr -= mTofr5ZCorr[iModuleChan][iZBin];
1682  } else {
1683  tofcorr = -9999.;
1684  }
1685 
1686  gMessMgr->Info("","OS") << " Corrected tof: tofcorr = " << tofcorr << endm;
1687  return tofcorr;
1688 }
1689 
1690 //_____________________________________________________________________________
1691 Double_t StTofCalibMaker::tofr8AllCorr(const Double_t tof, const Double_t tot, const Double_t z, const Int_t iTray, const Int_t iModuleChan)
1692 {
1693  int tray = iTray;
1694  int module = iModuleChan/6 + 1;
1695  int cell = iModuleChan%6 + 1;
1696  int board = iModuleChan/24 + 1;
1697  gMessMgr->Info("","OS") << "\nStTofCalibMaker::tofr8AllCorr: Tofr8 calibrating...\n"
1698  << "\tDoing Calibration in TOFr8 Tray " << tray << " Module " << module << " Cell " << cell
1699  << "\n\tinput tof = " << tof
1700  << " TOT = " << tot << " Zlocal = " << z
1701  << "\n" << endm;
1702 
1703  Double_t tofcorr = tof;
1704 
1705  tofcorr -= mTofTZero[tray-1][module-1][cell-1];
1706 
1707  if(Debug()) gMessMgr->Info("","OS") << "T0 correction: "<<mTofTZero[tray-1][module-1][cell-1]<<endm;
1708 
1709  if(mSlewingCorr) {
1710  int iTotBin = -1;
1711  for(int i=0;i<mNBinMax-1;i++) {
1712  if(tot>=mTofTotEdge[tray-1][board-1][i] && tot<mTofTotEdge[tray-1][board-1][i+1]) {
1713  iTotBin = i;
1714  break;
1715  }
1716  }
1717  if(iTotBin>=0&&iTotBin<mNBinMax) {
1718  double x1 = mTofTotEdge[tray-1][board-1][iTotBin];
1719  double x2 = mTofTotEdge[tray-1][board-1][iTotBin+1];
1720  double y1 = mTofTotCorr[tray-1][board-1][iTotBin];
1721  double y2 = mTofTotCorr[tray-1][board-1][iTotBin+1];
1722  double dcorr = y1 + (tot-x1)*(y2-y1)/(x2-x1);
1723  if(Debug()) gMessMgr->Info("","OS") << "TOT correction: "<<dcorr<<endm;
1724 
1725  tofcorr -= dcorr;
1726  } else {
1727  gMessMgr->Info("","OS") << " TOT out of range! EXIT! " << endm;
1728  return -9999.;
1729  }
1730 
1731  int iZBin = -1;
1732  for(int i=0;i<mNBinMax-1;i++) {
1733  if(z>=mTofZEdge[tray-1][board-1][i] && z<mTofZEdge[tray-1][board-1][i+1]) {
1734  iZBin = i;
1735  break;
1736  }
1737  }
1738  if(iZBin>=0&&iZBin<mNBinMax) {
1739  double x1 = mTofZEdge[tray-1][board-1][iZBin];
1740  double x2 = mTofZEdge[tray-1][board-1][iZBin+1];
1741  double y1 = mTofZCorr[tray-1][board-1][iZBin];
1742  double y2 = mTofZCorr[tray-1][board-1][iZBin+1];
1743 // double dcorr = y1 + (tot-x1)*(y2-y1)/(x2-x1);
1744  double dcorr = y1 + (z-x1)*(y2-y1)/(x2-x1);
1745 
1746  tofcorr -= dcorr;
1747  if(Debug()) gMessMgr->Info("","OS") << "zHit correction: "<<dcorr<<endm;
1748 
1749  } else {
1750  gMessMgr->Info("","OS") << " Z our of range! EXIT! " << endm;
1751  return -9999.;
1752  }
1753 
1754  }
1755 
1756  gMessMgr->Info("","OS") << " Corrected tof: tofcorr = " << tofcorr << endm;
1757  return tofcorr;
1758 }
1759 
1760 //_____________________________________________________________________________
1761 Double_t StTofCalibMaker::tstart(const Double_t *adc, const Double_t *tdc, const Double_t vz)
1762 {
1765  Double_t Tstart = 9999.;
1766 
1767  Int_t Ieast = 0, Iwest = 0;
1768  Double_t TdcSumEast = 0., TdcSumWest = 0.;
1769  Double_t ped = 0.0;
1770  for(int i=0;i<3;i++) {
1771  if( validPVPDADC(adc[i]) && validPVPDTDC(tdc[i]) ) {
1772  Ieast++;
1773  Double_t par[mNPar];
1774  for(int j=0;j<mNPar;j++) {
1775  par[j] = mPVPDTAPar[i*mNPar+j];
1776  }
1777  mPVPDSlewing->SetParameters(par);
1778  mPVPDTdc[i] = tdc[i] - mPVPDSlewing->Eval(adc[i]-ped);
1779  TdcSumEast += mPVPDTdc[i];
1780  }
1781 
1782  if( validPVPDADC(adc[i+3]) && validPVPDTDC(tdc[i+3]) ) {
1783  Iwest++;
1784  Double_t par[mNPar];
1785  for(int j=0;j<mNPar;j++) {
1786  par[j] = mPVPDTAPar[(i+3)*mNPar+j];
1787  }
1788  mPVPDSlewing->SetParameters(par);
1789  mPVPDTdc[i+3] = tdc[i+3] - mPVPDSlewing->Eval(adc[i+3]-ped);
1790  TdcSumWest += mPVPDTdc[i+3];
1791  }
1792  }
1793 
1794  Double_t TdcSum = TdcSumEast + TdcSumWest;
1795 
1796  // if ( Ieast>=mPVPDEastHitsCut && Iwest>=mPVPDWestHitsCut ) {
1798  if ( (Ieast>=mPVPDEastHitsCut || !mEastPVPDValid) && Iwest>=mPVPDWestHitsCut ) {
1799  Tstart = (TdcSum*mTDCWidth-(Ieast-Iwest)*vz/(C_C_LIGHT/1.e9))/(Ieast+Iwest);
1800  }
1801 
1802  return Tstart;
1803 }
1804 //_____________________________________________________________________________
1805 Double_t StTofCalibMaker::tstart5(const Double_t *tot, const Double_t *time, const Double_t vz)
1806 {
1809  Double_t Tstart = -999999.;
1810 
1811  Int_t Ieast = 0, Iwest = 0;
1812  Double_t TSumEast = 0., TSumWest = 0.;
1813  for(int i=0;i<3;i++) {
1814  if(Debug()) gMessMgr->Info("","OS") << " East pVPD tot = " << tot[i] << " time = " << time[i] << endm;
1815  if(Debug()) gMessMgr->Info("","OS") << " West pVPD tot = " << tot[i+3] << " time = " << time[i+3] << endm;
1816 
1817  if( time[i]>0. ) {
1818  int ibin = -1;
1819  for(int j=0;j<mNBinMax-1;j++) {
1820  if(tot[i]>=mPVPDTotEdge[i][j] && tot[i]<mPVPDTotEdge[i][j+1]) {
1821  ibin = j;
1822  break;
1823  }
1824  }
1825  if(ibin>=0&&ibin<mNBinMax) {
1826  Ieast++;
1827  mPVPDLeTime[i] = time[i] - mPVPDTotCorr[i][ibin];
1828  TSumEast += mPVPDLeTime[i];
1829  }
1830  }
1831 
1832  if( time[i+3]>0. ) {
1833  int ibin = -1;
1834  for(int j=0;j<mNBinMax-1;j++) {
1835  if(tot[i+3]>=mPVPDTotEdge[i+3][j] && tot[i+3]<mPVPDTotEdge[i+3][j+1]) {
1836  ibin = j;
1837  break;
1838  }
1839  }
1840  if(ibin>=0&&ibin<mNBinMax) {
1841  Iwest++;
1842  mPVPDLeTime[i+3] = time[i+3] - mPVPDTotCorr[i+3][ibin];
1843  TSumWest += mPVPDLeTime[i+3];
1844  }
1845  }
1846  }
1847 
1848  Double_t TSum = TSumEast + TSumWest;
1849 
1850  if ( Ieast>=mPVPDEastHitsCut && Iwest>=mPVPDWestHitsCut ) {
1851  Tstart = (TSum-(Ieast-Iwest)*vz/(C_C_LIGHT/1.e9))/(Ieast+Iwest);
1852  }
1853 
1854  return Tstart;
1855 }
1856 
1857 //_____________________________________________________________________________
1858 void StTofCalibMaker::tsum8(const Double_t *tot, const Double_t *time)
1859 {
1861  mTSumEast = 0.;
1862  mTSumWest = 0.;
1863  mNEast = 0;
1864  mNWest = 0;
1865  mVPDHitPatternWest = 0;
1866  mVPDHitPatternEast = 0;
1867 
1868  // West
1869  for(int i=0;i<mNVPD;i++) {
1870  if( time[i]>0. && tot[i]>0. ) {
1871  // west no offset
1872  if(mSlewingCorr) { // a switch
1873  int ibin = -1;
1874  for(int j=0;j<mNBinMax-1;j++) {
1875  if(tot[i]>=mVPDTotEdge[i][j] && tot[i]<mVPDTotEdge[i][j+1]) {
1876  ibin = j;
1877  break;
1878  }
1879  }
1880  if(ibin>=0&&ibin<mNBinMax) {
1881  mNWest++;
1882 
1883  double x1 = mVPDTotEdge[i][ibin];
1884  double x2 = mVPDTotEdge[i][ibin+1];
1885  double y1 = mVPDTotCorr[i][ibin];
1886  double y2 = mVPDTotCorr[i][ibin+1];
1887  double dcorr = y1 + (tot[i]-x1)*(y2-y1)/(x2-x1);
1888 
1889  mVPDLeTime[i] = time[i] - dcorr;
1890  mTSumWest += mVPDLeTime[i];
1891 
1892  mVPDHitPatternWest |= 1<<i;
1893  }
1894  } else {
1895  gMessMgr->Info("","OS") << " Vpd West tube " << i+1 << " TOT out of range!" << endm;
1896  // out of range, remove this hit
1897  }
1898  }
1899 
1900  // East
1901  if( time[i+mNVPD]>0. && tot[i+mNVPD]>0. ) {
1902  double tmp = time[i+mNVPD] - mPhaseOffset8;
1903  while(tmp>TMAX) tmp -= TMAX;
1904  if(mSlewingCorr) { // a switch
1905  int ibin = -1;
1906  for(int j=0;j<mNBinMax-1;j++) {
1907  if(tot[i+mNVPD]>=mVPDTotEdge[i+mNVPD][j] && tot[i+mNVPD]<mVPDTotEdge[i+mNVPD][j+1]) {
1908  ibin = j;
1909  break;
1910  }
1911  }
1912  if(ibin>=0&&ibin<mNBinMax) {
1913  mNEast++;
1914 
1915  double x1 = mVPDTotEdge[i+mNVPD][ibin];
1916  double x2 = mVPDTotEdge[i+mNVPD][ibin+1];
1917  double y1 = mVPDTotCorr[i+mNVPD][ibin];
1918  double y2 = mVPDTotCorr[i+mNVPD][ibin+1];
1919  double dcorr = y1 + (tot[i+mNVPD]-x1)*(y2-y1)/(x2-x1);
1920 
1921  mVPDLeTime[i+mNVPD] = tmp - dcorr;
1922  mTSumEast += mVPDLeTime[i+mNVPD];
1923 
1924  mVPDHitPatternEast |= 1<<i;
1925 
1926  }
1927  } else {
1928  gMessMgr->Info("","OS") << " Vpd East tube " << i+1 << " TOT out of range!" << endm;
1929  // out of range, remove this hit
1930  }
1931  }
1932  }
1933 
1934  if ( mNEast>=mPVPDEastHitsCut && mNWest>=mPVPDWestHitsCut ) {
1935  mTDiff = (mTSumEast/mNEast - mTSumWest/mNWest)/2. - mProjVtxZ/(C_C_LIGHT/1.e9);
1936  mVPDVtxZ = (mTSumEast/mNEast - mTSumWest/mNWest)/2.*(C_C_LIGHT/1.e9);
1937  }
1938 
1939  return;
1940 }
1941 //_____________________________________________________________________________
1942 Double_t StTofCalibMaker::tstart8(const Double_t vz)
1943 {
1944  Double_t Tstart = -9999.;
1945 
1946  Double_t TSum = mTSumEast + mTSumWest;
1947 
1948  if ( mNEast>=mPVPDEastHitsCut && mNWest>=mPVPDWestHitsCut ) {
1949  Tstart = (TSum-(mNEast-mNWest)*vz/(C_C_LIGHT/1.e9))/(mNEast+mNWest);
1950  }
1951 
1952  return Tstart;
1953 }
1954 
1955 //_____________________________________________________________________________
1956 Double_t StTofCalibMaker::GetINLcorr(const int edgeflag,const int tdcchan,const int bin)
1957 {
1958  int iboard=tdcchan/24;
1959  int chan = (tdcchan%24);
1960  int itdc=0;
1961  if(edgeflag==4) itdc = chan/8; // leading edge
1962  if(edgeflag==5) itdc = 3; // trailing edge
1963  if(tdcchan==192&&edgeflag==4) itdc=0;
1964  if(tdcchan==193&&edgeflag==4) itdc=1;
1965  if(tdcchan==194&&edgeflag==4) itdc=1;
1966  if(tdcchan==195&&edgeflag==4) itdc=0;
1967  if(tdcchan==196&&edgeflag==4) itdc=1;
1968  if(tdcchan==197&&edgeflag==4) itdc=1;
1969 
1970  return mINLtable[iboard][itdc][bin];
1971 }
1972 
1973 /*****************************************************************
1974  *
1975  * $Log: StTofCalibMaker.cxx,v $
1976  * Revision 1.24 2018/02/26 23:26:50 smirnovd
1977  * StTof: Remove outdated ClassImp macro
1978  *
1979  * Revision 1.23 2018/02/26 23:13:19 smirnovd
1980  * Move embedded CVS log messages to the end of file
1981  *
1982  * Revision 1.22 2012/12/14 06:35:48 geurts
1983  * Changed global database calls to direct table access and/or removed deprecated database access code.
1984  *
1985  * Revision 1.21 2011/05/27 18:25:32 genevb
1986  * Propagate StTrack::key => Int_t to other codes
1987  *
1988  * Revision 1.20 2008/09/03 22:30:43 dongx
1989  * mTStart added, applicable for Run8
1990  *
1991  * Revision 1.19 2008/07/30 20:03:03 dongx
1992  * mBeamLine initialized in constructor to avoid null pointer when no calibration tables
1993  *
1994  * Revision 1.18 2008/07/22 00:16:47 dongx
1995  * initialization added for global variables in StTofCollection in each event
1996  *
1997  * Revision 1.17 2008/06/25 21:45:29 dongx
1998  * Reset time values when tot/z is out of range
1999  *
2000  * Revision 1.16 2008/06/24 21:33:37 dongx
2001  * Added the filling of vzvpd into StTofCollection
2002  *
2003  * Revision 1.15 2008/06/17 17:49:19 dongx
2004  * Update for Run 8 - first release
2005  *
2006  * Revision 1.14 2007/11/30 17:11:23 dongx
2007  * removed tdcId in tofZCorr for tofZCorr.tdcId is not defined in db, and removed from idl definition
2008  *
2009  * Revision 1.13 2007/03/13 15:09:10 dongx
2010  * Remove breaking of failure on number of return rows during db I/O for tofTotCorr and tofZCorr
2011  *
2012  * Revision 1.12 2007/03/05 18:51:02 dongx
2013  * updated for Run V CuCu calibration
2014  * - INL correction moved in this maker
2015  * - Tot Corr and Z Corr use new tables in data base
2016  * - pVPD calibrated information cannot be fully stored within current infrastructure, need update on TofCollection. Configurations better than (1,1) are all selected.
2017  *
2018  * Revision 1.11 2005/04/12 17:33:47 dongx
2019  * update for year 5 data. not completed, leave as empty now.
2020  *
2021  * Revision 1.10 2004/09/20 16:07:35 dongx
2022  * correct the nsigma in StTofPidTraits
2023  *
2024  * Revision 1.9 2004/08/13 00:15:03 dongx
2025  * correct the nSigmaXXX calculation
2026  *
2027  * Revision 1.8 2004/08/11 19:35:40 dongx
2028  * loose the ADC cut of tofp
2029  *
2030  * Revision 1.7 2004/08/11 18:58:40 dongx
2031  * missing nSigmaXX in tofHit implemented
2032  *
2033  * Revision 1.6 2004/07/24 03:33:56 dongx
2034  * Tofp slewing function changed back
2035  *
2036  * Revision 1.5 2004/07/16 18:28:17 dongx
2037  * -Tofp Slewing function changed in AuAu200 GeV Run IV
2038  * -Include those runs with eastern PVPD dead
2039  *
2040  * Revision 1.4 2004/07/16 15:06:08 dongx
2041  * Z correction function separated for TOFp and TOFr.
2042  * Use a new one for RunIV AuAu 200GeV runs
2043  *
2044  * Revision 1.3 2004/07/15 18:11:22 dongx
2045  * -introduce two new tables in dbase: tofAdcRange & tofResolution
2046  * -continue update on writing StTofPidTraits
2047  *
2048  * Revision 1.2 2004/07/08 18:26:09 dongx
2049  * filling StTofPidTraits added (not completed, null nsigmaXXX now)
2050  *
2051  * Revision 1.1 2004/07/01 17:23:48 dongx
2052  * first release
2053  */
Double_t tofpSlewingCorr(const Double_t tof, const Double_t adc, const Int_t iDaqChan)
tofp TA slewing correction function
void resetPars()
Reset the calibration parameters.
Int_t initParameters(int)
Initialize the calibration parameters from dbase.
void setOuterGeometry(const bool)
switch to select the Helix Geometry
pair< double, double > pathLengths(const StHelix &, double minStepSize=10 *micrometer, double minRange=10 *centimeter) const
path lengths at dca between two helices
Definition: StHelix.cc:511
Double_t tofrAllCorr(const Double_t tof, const Double_t Tstart, const Double_t adc, const Double_t zlocal, const Int_t iDaqChan)
tofr overall correction function
Double_t GetINLcorr(const int edgeid, const int tempchan, const int bin)
Run 5 -&gt;
Double_t tstart(const Double_t *adc, const Double_t *tdc, const Double_t Vz)
Start timing calculation function.
Definition: tof.h:15
void tsum8(const Double_t *tot, const Double_t *time)
Run 8 -&gt;
void initFormulas()
Initialize the formulas of calibration functions.
Double_t tofrSlewingCorr(const Double_t tof, const Double_t adc, const Int_t iDaqChan)
tofr TA slewing correction function
StTofCalibMaker(const char *name="tofCalib")
Default constructor.
Time-of-Flight Geometry Utilities.
Double_t tofpZCorr(const Double_t tof, const Double_t zlocal)
tofp loca Zhit correction function
Double_t tofpT0Corr(const Double_t tof, const Double_t Tstart, const Int_t iDaqChan)
tofp T0 calibration function
Definition: Stypes.h:40
void clearFormulars()
Delete the calibration functions.
Double_t tofrZCorr(const Double_t tof, const Double_t zlocal)
tofr local Zhit correction function
Double_t tstart5(const Double_t *tot, const Double_t *time, const Double_t Vz)
Definition: Stypes.h:44
~StTofCalibMaker()
Destructor.
Double_t tofrT0Corr(const Double_t tof, const Double_t Tstart, const Int_t iDaqChan)
tofr T0 calibration function
Double_t tofr8AllCorr(const Double_t tof, const Double_t tot, const Double_t zlocal, const Int_t iTray, const Int_t iModuleChan)
tstart calculation splitted into 2 steps in Run 8 since no primary vertex in the event ...
Double_t tofpAllCorr(const Double_t tof, const Double_t Tstart, const Double_t adc, const Double_t zlocal, const Int_t iDaqChan)
tofp overall correction function
virtual TDataSet * Find(const char *path) const
Definition: TDataSet.cxx:362
void init()
initialize geometry class from XDF file and set-up DAQ/Slat mappings