]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFcalibHisto.cxx
Minor bugfixes and extensions for the onlineDiplay Interface
[u/mrichter/AliRoot.git] / TOF / AliTOFcalibHisto.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*************************************************************************
17  *
18  * AliTOFcalibHisto - class to handle TOF calibration histograms,
19  *                    map histograms and more
20  *
21  * autors:   Roberto Preghenella (R+)
22  * concacts: preghenella@bo.infn.it
23  *
24  *************************************************************************/
25
26 #include "AliTOFcalibHisto.h"
27 #include "AliLog.h"
28 #include "TH1F.h"
29 #include "TFile.h"
30 #include "AliTOFRawStream.h"
31 #include "AliTOFCableLengthMap.h"
32 #include "AliESDtrack.h"
33
34 ClassImp(AliTOFcalibHisto)
35
36 //__________________________________________________________________________
37
38 TFile *AliTOFcalibHisto::fgCalibHistoFile = NULL;
39 TFile *AliTOFcalibHisto::fgCalibParFile = NULL;
40
41 //__________________________________________________________________________
42
43 TString AliTOFcalibHisto::fgCalibHistoFileName = "$ALICE_ROOT/TOF/data/AliTOFcalibHisto.root";
44 TString AliTOFcalibHisto::fgCalibParFileName = "$ALICE_ROOT/TOF/data/AliTOFcalibPar.root";
45
46 //__________________________________________________________________________
47
48 const TString AliTOFcalibHisto::fgkCalibConstName[kNcalibConsts] = {
49   "LHCperiod",
50   "AmphenolCableDelay",
51   "FlatCableDelay",
52   "InterfaceCardDelay"
53 };
54
55 //__________________________________________________________________________
56
57 const TString AliTOFcalibHisto::fgkCalibMapName[kNcalibMaps] = {
58   /* main index */
59   "Index",
60   /* EO indices */
61   "DDL",
62   "TRM", 
63   "Chain", 
64   "TDC", 
65   "Channel", 
66   /* DO indices */
67   "Sector", 
68   "Plate", 
69   "Strip", 
70   "SectorStrip", 
71   "PadZ", 
72   "PadX", 
73   "Pad",
74   "InterfaceCardIndex",
75   /* calib constants */
76   "DDLBCshift",
77   "FlatCableLength",
78   "InterfaceCardLength",
79   "AmphenolCableLength"
80 };
81
82 //__________________________________________________________________________
83
84 const TString AliTOFcalibHisto::fgkCalibParName[kNcalibPars] = {
85   "hDDLDelay",
86   "hHPTDCDelay",
87   "hLeftFEAchDelay",
88   "hRightFEAchDelay",
89   "hFEADelay",
90   "hTRMDelay",
91   "hSlew"
92 };
93
94 //__________________________________________________________________________
95
96 /* LHC clock period [ns] */
97 const Float_t AliTOFcalibHisto::fgkLHCperiod = 25.; /* SET THE CORRECT VALUE !!! */
98
99 //__________________________________________________________________________
100
101 /* Amphenol cable delay [ns/cm] */
102 const Float_t AliTOFcalibHisto::fgkAmphenolCableDelay = 5.13e-2; /* from measurement */
103
104 //__________________________________________________________________________
105
106 /* flat cable delay [ns/cm] */
107 //const Float_t AliTOFcalibHisto::fgkFlatCableDelay = 5.3e-2; /* from Amphenol 132-2829 series data-sheet */
108 const Float_t AliTOFcalibHisto::fgkFlatCableDelay = 5.124e-2; /* from LHC08d calibration */
109
110 //__________________________________________________________________________
111
112 /* interface card delay [ns/cm] */
113 //const Float_t AliTOFcalibHisto::fgkInterfaceCardDelay = 6.9e-2; /* from HyperLinx simulation */
114 const Float_t AliTOFcalibHisto::fgkInterfaceCardDelay = 5.7898e-2; /* from LHC08d calibration */
115
116 //__________________________________________________________________________
117
118 /* number of readout channels (DO/EO) */
119 const Int_t AliTOFcalibHisto::fgkNchannels = 157248;
120 const Int_t AliTOFcalibHisto::fgkNchannelsEO = 172800;
121
122 //__________________________________________________________________________
123
124 /* DDL BC shifts due to TTC fibers [LHCperiod] */
125 const Int_t AliTOFcalibHisto::fgkDDLBCshift[72] = {
126   2, 2, -1, -1,
127   2, 2, 0, 0,
128   2, 2, 0, 0,
129   2, 2, 0, 0,
130   2, 2, 0, 0,
131   2, 2, 0, 0,
132   2, 2, 0, 0,
133   2, 2, 0, 0,
134   2, 2, 0, 0,
135   2, 2, 0, 0,
136   2, 2, -1, -1,
137   2, 2, -1, -1,
138   2, 2, -2, -2,
139   2, 2, -2, -2,
140   2, 2, -2, -2,
141   2, 2, -1, -1,
142   2, 2, -1, -1,
143   2, 2, -1, -1,
144 };
145
146 //__________________________________________________________________________
147
148 /* strip flat-cable length (preliminary) [cm] */
149 const Float_t AliTOFcalibHisto::fgkFlatCableLength[91] = {
150   18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 17.,
151   21., 21., 21., 21., 21., 17., 17., 21., 21., 17., 21., 21., 21., 17., 21., 21., 17., 21., 23.,
152   17., 19., 17., 19., 17., 19., 17., 19., 17., 19., 17., 19., 17., 19., 17.,
153   23., 21., 17., 21., 21., 17., 21., 21., 21., 17., 21., 21., 17., 17., 21., 21., 21., 21., 21.,
154   17., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18.
155 };
156
157 //__________________________________________________________________________
158
159 /* interface card lenght (preliminary) [cm] */
160 const Float_t AliTOFcalibHisto::fgkInterfaceCardLength[48] = {
161   13.97, 12.57, 14.52, 13.10, 15.44, 13.60, 10.58, 9.14, 
162   11.21, 9.76, 12.11, 10.76, 8.67, 7.58, 9.32, 8.09,
163   10.24, 8.4, 5.51, 4.31, 6.54, 5.23, 7.48, 6.28,
164   10.43, 8.76, 11.05, 9.43, 11.72, 10.14, 7.2, 5.69,
165   7.71, 6.26, 8.36, 7.19, 4.85, 4.09, 5.57, 4.35, 
166   6.59, 5.12, 2.49, 2.96, 2.70, 2.76, 2.91, 2.55
167 };
168
169 //__________________________________________________________________________
170
171 AliTOFcalibHisto::AliTOFcalibHisto() :
172   TObject(),
173   fCalibConst(),
174   fCalibMap(),
175   fCalibPar()
176 {
177   /* default constructor */
178 }
179
180 //__________________________________________________________________________
181
182 AliTOFcalibHisto::~AliTOFcalibHisto()
183 {
184   /* default destructor */
185 }
186
187 //__________________________________________________________________________
188
189 void 
190 AliTOFcalibHisto::LoadHisto(TFile* file, TH1F **histo, const Char_t *name) 
191 {
192   /* load histo */
193   *histo = (TH1F *)file->Get(name);
194   if (!*histo)
195     AliWarning(Form("error while getting %s histo", name));
196 }
197
198 //__________________________________________________________________________
199
200 void 
201 AliTOFcalibHisto::CreateHisto(TH1F **histo, const Char_t *name, Int_t size) 
202 {
203   /* create histo */
204   *histo = new TH1F(name, Form(";index;%s", name), size, 0, size);
205   if (!*histo)
206     AliWarning(Form("error while creating %s histo", name));
207 }
208
209 //__________________________________________________________________________
210
211 void 
212 AliTOFcalibHisto::WriteHisto(TFile *file, TH1F *histo) 
213 {
214   /* write histo */
215   if (!file || !file->IsOpen() || !histo)
216     return;
217   file->cd(); 
218   histo->Write();
219 }
220
221 //__________________________________________________________________________
222
223 void
224 AliTOFcalibHisto::SetHisto(TH1F *histo, Int_t index, Float_t value)
225 {
226   /* set histo */
227   if (!histo)
228     return;
229   histo->SetBinContent(index + 1, value);
230 }
231
232 //__________________________________________________________________________
233
234 Float_t
235 AliTOFcalibHisto::GetHisto(TH1F *histo, Int_t index)
236 {
237   /* get histo */
238   if (!histo) {
239     AliWarning("cannot get histo");
240     return 0.;
241   }
242   return histo->GetBinContent(index + 1);
243 }
244
245 //__________________________________________________________________________
246
247 void
248 AliTOFcalibHisto::LoadCalibHisto()
249 {
250   /* load calib histo */
251
252   if (fgCalibHistoFile && fgCalibHistoFile->IsOpen())
253     AliWarning("calib histo file already open: reloading"); 
254
255   /* open input file */
256   TFile *fileIn = TFile::Open(GetCalibHistoFileName());
257   if (!fileIn || !fileIn->IsOpen())
258     AliFatal(Form("cannot open input file %s", GetCalibHistoFileName()));
259
260   /* set calib histo file */
261   fgCalibHistoFile = fileIn;
262
263   /* load consts */
264   for (Int_t iConst = 0; iConst < kNcalibConsts; iConst++)
265     LoadHisto(fileIn, &fCalibConst[iConst], fgkCalibConstName[iConst].Data());
266   /* load maps */
267   for (Int_t iMap = 0; iMap < kNcalibMaps; iMap++)
268     LoadHisto(fileIn, &fCalibMap[iMap], fgkCalibMapName[iMap].Data());
269 }
270
271 //__________________________________________________________________________
272
273 void
274 AliTOFcalibHisto::LoadCalibPar()
275 {
276   /* load calib par */
277
278   if (fgCalibParFile && fgCalibParFile->IsOpen())
279     AliWarning("calib par file already open: reloading"); 
280
281   /* load calib histo */
282   LoadCalibHisto();
283
284   /* open input file */
285   TFile *fileIn = TFile::Open(GetCalibParFileName());
286   if (!fileIn || !fileIn->IsOpen())
287     AliError(Form("cannot open input file %s", GetCalibParFileName()));
288
289   /* set calib par file */
290   fgCalibParFile = fileIn;
291
292   /* load pars */
293   for (Int_t i = 0; i < kNcalibPars; i++)
294     LoadHisto(fileIn, &fCalibPar[i], fgkCalibParName[i].Data());
295
296 }
297
298 //__________________________________________________________________________
299
300 void
301 AliTOFcalibHisto::WriteCalibHisto()
302 {
303   /* write calib histo */
304
305   /* open output file */
306   TFile *fileOut = TFile::Open(GetCalibHistoFileName(), "RECREATE");
307   if (!fileOut || !fileOut->IsOpen())
308     AliFatal(Form("cannot open output file %s", GetCalibHistoFileName()));
309
310   /* create consts */
311   for (Int_t iConst = 0; iConst < kNcalibConsts; iConst++)
312     CreateHisto(&fCalibConst[iConst], fgkCalibConstName[iConst].Data(), 1);
313   /* create maps */
314   for (Int_t iMap = 0; iMap < kNcalibMaps; iMap++)
315     if (iMap == kIndex)
316       CreateHisto(&fCalibMap[iMap], fgkCalibMapName[iMap].Data(), fgkNchannelsEO);
317     else
318       CreateHisto(&fCalibMap[iMap], fgkCalibMapName[iMap].Data(), fgkNchannels);
319
320   /*** SETUP CONSTS ***/
321
322   SetHisto(fCalibConst[kLHCperiod], 0, fgkLHCperiod);
323   SetHisto(fCalibConst[kAmphenolCableDelay], 0, fgkAmphenolCableDelay);
324   SetHisto(fCalibConst[kFlatCableDelay], 0, fgkFlatCableDelay);
325   SetHisto(fCalibConst[kInterfaceCardDelay], 0, fgkInterfaceCardDelay);
326   
327   /***  SETUP MAPS  ***/
328
329   AliTOFRawStream rawStream;
330   Int_t indexEO, det[5], dummy, index, sector, plate, strip, sectorStrip, padz, padx, pad, icIndex;
331
332   /* temporarly disable warnings */
333   AliLog::EType_t logLevel = (AliLog::EType_t)AliLog::GetGlobalLogLevel();
334   AliLog::SetGlobalLogLevel(AliLog::kError);
335
336   /* loop over electronics oriented (EO) indices */
337   for (Int_t ddl = 0; ddl < 72; ddl++)
338     for (Int_t trm = 0; trm < 10; trm++)
339       for (Int_t chain = 0; chain < 2; chain++)
340         for (Int_t tdc = 0; tdc < 15; tdc++)
341           for (Int_t channel = 0; channel < 8; channel++) {
342             
343             /* compute index EO */
344             indexEO = GetIndexEO(ddl, trm, chain, tdc, channel);
345
346             /* convert EO indices into detector oriented (DO) indices
347                (this call causes some warnings because the loop includes
348                EO indices which are not connected to physical channels) */
349             rawStream.EquipmentId2VolumeId(ddl, trm + 3, chain, tdc, channel, det);
350             
351             /* swap det[3] and det[4] */
352             dummy = det[3]; det[3] = det[4]; det[4] = dummy;
353             
354             /* check detector indices */
355             if (det[0] < 0 || det[0] > 71 ||
356                 det[1] < 0 || det[1] > 4 ||
357                 det[2] < 0 || det[2] > 18 ||
358                 det[3] < 0 || det[3] > 1 ||
359                 det[4] < 0 || det[4] > 47) {
360               SetHisto(fCalibMap[kIndex], indexEO, -1);
361               continue;
362             }
363             
364             /* setup information */
365             index = AliTOFGeometry::GetIndex(det);
366             sector = det[0];
367             plate = det[1];
368             strip = det[2];
369             sectorStrip = plate < 3 ? plate * 19 + strip : plate * 19 - 4 + strip;
370             padz = det[3];
371             padx = det[4];
372             pad = padz + 2 * padx;
373             icIndex = pad < 48 ? pad : 95 - pad;
374
375             /* set maps */
376
377             /* main index */
378             SetHisto(fCalibMap[kIndex], indexEO, index);
379             /* EO indices */
380             SetHisto(fCalibMap[kDDL], index, ddl);
381             SetHisto(fCalibMap[kTRM], index, trm);
382             SetHisto(fCalibMap[kChain], index, chain);
383             SetHisto(fCalibMap[kTDC], index, tdc);
384             SetHisto(fCalibMap[kChannel], index, channel);
385             /* DO indices */
386             SetHisto(fCalibMap[kSector], index, sector);
387             SetHisto(fCalibMap[kPlate], index, plate);
388             SetHisto(fCalibMap[kStrip], index, strip);
389             SetHisto(fCalibMap[kSectorStrip], index, sectorStrip);
390             SetHisto(fCalibMap[kPadZ], index, padz);
391             SetHisto(fCalibMap[kPadX], index, padx);
392             SetHisto(fCalibMap[kPad], index, pad);
393             SetHisto(fCalibMap[kInterfaceCardIndex], index, icIndex);
394             /* calib constants */
395             SetHisto(fCalibMap[kDDLBCshift], index, fgkDDLBCshift[ddl]);
396             SetHisto(fCalibMap[kFlatCableLength], index, fgkFlatCableLength[sectorStrip]);
397             SetHisto(fCalibMap[kInterfaceCardLength], index, fgkInterfaceCardLength[icIndex]);
398             SetHisto(fCalibMap[kAmphenolCableLength], index, AliTOFCableLengthMap::GetCableLength(ddl, trm + 3, chain, tdc));
399             
400           } /* loop over electronics oriented (EO) indices */
401
402   /* re-enable warnings */
403   AliLog::SetGlobalLogLevel(logLevel);
404
405   /* write consts */
406   for (Int_t iConst = 0; iConst < kNcalibConsts; iConst++)
407     WriteHisto(fileOut, fCalibConst[iConst]);
408   /* write maps */
409   for (Int_t iMap = 0; iMap < kNcalibMaps; iMap++)
410     WriteHisto(fileOut, fCalibMap[iMap]);
411
412   /* close output file */
413   fileOut->Close();
414 }
415
416 //__________________________________________________________________________
417
418 Float_t
419 AliTOFcalibHisto::GetCorrection(Int_t corr, Int_t index, Float_t tot)
420 {
421   /* apply correction */
422
423   Int_t ddl, chain, tdc, channel, hptdc, pbCh, feaIndex, sector, plate, strip, padx, trm;
424   Float_t slewing;
425   
426   switch (corr) {
427   case kDDLBCcorr:
428     return -GetCalibConst(kLHCperiod) * GetCalibMap(kDDLBCshift, index);
429   case kAmphenolCableCorr:
430     return GetCalibConst(kAmphenolCableDelay) * GetCalibMap(kAmphenolCableLength, index);
431   case kFlatCableCorr:
432     return GetCalibConst(kFlatCableDelay) * GetCalibMap(kFlatCableLength, index);
433   case kInterfaceCardCorr:
434     return GetCalibConst(kInterfaceCardDelay) * GetCalibMap(kInterfaceCardLength, index);
435   case kDDLdelayCorr:
436     ddl = (Int_t)GetCalibMap(kDDL, index);
437     return GetCalibPar(kDDLdelayPar, ddl);
438   case kHPTDCdelayCorr:
439     chain = (Int_t)GetCalibMap(kChain, index);
440     tdc = (Int_t)GetCalibMap(kTDC, index);
441     hptdc = tdc + 15 * chain;
442     return GetCalibPar(kHPTDCdelayPar, hptdc);
443   case kFEAchDelayCorr:
444     ddl = (Int_t)GetCalibMap(kDDL, index);
445     tdc = (Int_t)GetCalibMap(kTDC, index);
446     channel = (Int_t)GetCalibMap(kChannel, index);
447     pbCh = channel + 8 * (tdc % 3);
448     if (ddl % 2 == 0)
449       return GetCalibPar(kRightFEAchDelayPar, pbCh);
450     else
451       return GetCalibPar(kLeftFEAchDelayPar, pbCh);
452   case kFEAdelayCorr:
453     sector = (Int_t)GetCalibMap(kSector, index);
454     plate = (Int_t)GetCalibMap(kPlate, index);
455     strip = (Int_t)GetCalibMap(kStrip, index);
456     padx = (Int_t)GetCalibMap(kPadX, index);
457     feaIndex = padx / 12 + 4 * strip + 4 * 19 * plate + 4 * 19 * 5 * sector;      
458     return GetCalibPar(kFEAdelayPar, feaIndex);
459   case kTRMdelayCorr:
460     trm = (Int_t)GetCalibMap(kTRM, index);
461     return GetCalibPar(kTRMdelayPar, trm);
462   case kTimeSlewingCorr:
463     slewing = 0.;
464     for (Int_t i = 0; i < fCalibPar[kTimeSlewingPar]->GetNbinsX(); i++)
465       slewing += GetCalibPar(kTimeSlewingPar, i) * TMath::Power(tot, i);
466     return slewing;
467   default:
468     AliWarning(Form("unknown correction flag (%d)", corr));
469     return 0.;
470   }
471 }
472
473 //__________________________________________________________________________
474
475 Float_t
476 AliTOFcalibHisto::GetNominalCorrection(Int_t index)
477 {
478   /* get nominal correction */
479   Float_t corr = 0;
480   for (Int_t iCorr = 0; iCorr < kNcorrections; iCorr++)
481     corr += GetCorrection(iCorr, index);
482   return corr;
483 }
484
485 //__________________________________________________________________________
486
487 void
488 AliTOFcalibHisto::ApplyNominalCorrection(AliESDtrack *track)
489 {
490   /* apply nominal correction */
491   
492   Double_t rawTime = track->GetTOFsignalRaw();
493   Int_t index = track->GetTOFCalChannel();
494   Double_t time = rawTime - 1.e3 * GetNominalCorrection(index);
495   track->SetTOFsignal(time);
496 }
497