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