]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFcalibHisto.cxx
correct mask for V0 charge decoding in STU payload
[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 // *
22 // * autors:   Roberto Preghenella (R+)
23 // * concacts: preghenella@bo.infn.it
24 // *
25 // *
26
27 #include "AliTOFcalibHisto.h"
28 #include "AliLog.h"
29 #include "TH1D.h"
30 #include "TFile.h"
31 #include "AliTOFRawStream.h"
32 #include "AliTOFCableLengthMap.h"
33 #include "AliESDtrack.h"
34
35 #define SLEW_TOTMIN 10.
36 #define SLEW_TOTMAX 16.
37
38 ClassImp(AliTOFcalibHisto)
39
40 //__________________________________________________________________________
41
42 TFile *AliTOFcalibHisto::fgCalibHistoFile = NULL;
43 TFile *AliTOFcalibHisto::fgCalibParFile = NULL;
44 TFile *AliTOFcalibHisto::fgCalibStatFile = NULL;
45
46 //__________________________________________________________________________
47
48 TString AliTOFcalibHisto::fgCalibHistoFileName = "$ALICE_ROOT/TOF/data/AliTOFcalibHisto.root";
49 TString AliTOFcalibHisto::fgCalibParFileName = "$ALICE_ROOT/TOF/data/AliTOFcalibPar.root";
50 TString AliTOFcalibHisto::fgCalibStatFileName = "$ALICE_ROOT/TOF/data/AliTOFcalibStat.root";
51
52 //__________________________________________________________________________
53
54 const TString AliTOFcalibHisto::fgkCalibConstName[kNcalibConsts] = {
55   "LHCperiod",
56   "AmphenolCableDelay",
57   "FlatCableDelay",
58   "InterfaceCardDelay"
59 };
60
61 //__________________________________________________________________________
62
63 const TString AliTOFcalibHisto::fgkCalibMapName[kNcalibMaps] = {
64   /* main index */
65   "Index",
66   /* EO indices */
67   "DDL",
68   "TRM", 
69   "Chain", 
70   "TDC", 
71   "Channel", 
72   /* DO indices */
73   "Sector", 
74   "Plate", 
75   "Strip", 
76   "SectorStrip", 
77   "PadZ", 
78   "PadX", 
79   "Pad",
80   "InterfaceCardIndex",
81   /* calib constants */
82   "DDLBCshift",
83   "FlatCableLength",
84   "InterfaceCardLength",
85   "AmphenolCableLength"
86 };
87
88 //__________________________________________________________________________
89
90 const TString AliTOFcalibHisto::fgkCalibParName[kNcalibPars] = {
91   "hDDLDelay",
92   "hHPTDCDelay",
93   "hLeftFEAchDelay",
94   "hRightFEAchDelay",
95   "hFEADelay",
96   "hICDelay",
97   "hTRMDelay",
98   "hStripDelay",
99   "hIndexDelay",
100   "hSlewing"
101 };
102
103 //__________________________________________________________________________
104
105 const TString AliTOFcalibHisto::fgkCalibStatName[kNcalibStats] = {
106   "hStripStat"
107 };
108
109 //__________________________________________________________________________
110
111 /* LHC clock period [ns] */
112 const Double_t AliTOFcalibHisto::fgkLHCperiod = (24.4e-3 * 1024); /* ns */
113
114 //__________________________________________________________________________
115
116 /* Amphenol cable delay [ns/cm] */
117 const Double_t AliTOFcalibHisto::fgkAmphenolCableDelay = 5.13e-2; /* from measurement */
118
119 //__________________________________________________________________________
120
121 /* flat cable delay [ns/cm] */
122 //const Double_t AliTOFcalibHisto::fgkFlatCableDelay = 5.3e-2; /* from Amphenol 132-2829 series data-sheet */
123 const Double_t AliTOFcalibHisto::fgkFlatCableDelay = 5.124e-2; /* from LHC08d calibration */
124
125 //__________________________________________________________________________
126
127 /* interface card delay [ns/cm] */
128 //const Double_t AliTOFcalibHisto::fgkInterfaceCardDelay = 6.9e-2; /* from HyperLinx simulation */
129 //const Double_t AliTOFcalibHisto::fgkInterfaceCardDelay = 5.7898e-2; /* from LHC08d calibration */
130 const Double_t AliTOFcalibHisto::fgkInterfaceCardDelay = 6.31360207815420404e-02; /* from LHC09c calibration */
131
132 //__________________________________________________________________________
133
134 /* number of readout channels (DO/EO) */
135 const Int_t AliTOFcalibHisto::fgkNchannels = 157248;
136 const Int_t AliTOFcalibHisto::fgkNchannelsEO = 172800;
137
138 //__________________________________________________________________________
139
140 /* DDL BC shifts due to TTC fibers [LHCperiod] */
141 const Int_t AliTOFcalibHisto::fgkDDLBCshift[72] = {
142   2, 2, -1, -1,
143   2, 2, 0, 0,
144   2, 2, 0, 0,
145   2, 2, 0, 0,
146   2, 2, 0, 0,
147   2, 2, 0, 0,
148   2, 2, 0, 0,
149   2, 2, 0, 0,
150   2, 2, 0, 0,
151   2, 2, 0, 0,
152   2, 2, -1, -1,
153   2, 2, -1, -1,
154   2, 2, -2, -2,
155   2, 2, -2, -2,
156   2, 2, -2, -2,
157   2, 2, -1, -1,
158   2, 2, -1, -1,
159   2, 2, -1, -1
160 };
161
162 //__________________________________________________________________________
163
164 /* strip flat-cable length (preliminary) [cm] */
165 const Double_t AliTOFcalibHisto::fgkFlatCableLength[91] = {
166   18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 17.,
167   21., 21., 21., 21., 21., 17., 17., 21., 21., 17., 21., 21., 21., 17., 21., 21., 17., 21., 23.,
168   17., 19., 17., 19., 17., 19., 17., 19., 17., 19., 17., 19., 17., 19., 17.,
169   23., 21., 17., 21., 21., 17., 21., 21., 21., 17., 21., 21., 17., 17., 21., 21., 21., 21., 21.,
170   17., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18., 18.
171 };
172
173 //__________________________________________________________________________
174
175 /* interface card lenght (preliminary) [cm] */
176 const Double_t AliTOFcalibHisto::fgkInterfaceCardLength[48] = {
177   13.97, 12.57, 14.52, 13.10, 15.44, 13.60, 10.58, 9.14, 
178   11.21, 9.76, 12.11, 10.76, 8.67, 7.58, 9.32, 8.09,
179   10.24, 8.4, 5.51, 4.31, 6.54, 5.23, 7.48, 6.28,
180   10.43, 8.76, 11.05, 9.43, 11.72, 10.14, 7.2, 5.69,
181   7.71, 6.26, 8.36, 7.19, 4.85, 4.09, 5.57, 4.35, 
182   6.59, 5.12, 2.49, 2.96, 2.70, 2.76, 2.91, 2.55
183 };
184
185 //__________________________________________________________________________
186
187 Bool_t AliTOFcalibHisto::fgCableCorrectionFlag[kNcorrections] = {
188   kFALSE, // kDDLBCcorr
189   kTRUE, // kAmphenolCableCorr
190   kTRUE, // kFlatCableCorr
191   kTRUE, // kInterfaceCardCorr
192   kFALSE, // kDDLdelayCorr
193   kFALSE, // kHPTDCdelayCorr
194   kFALSE, // kFEAchDelayCorr
195   kFALSE, // kFEAdelayCorr
196   kFALSE, // kTRMdelayCorr
197   kFALSE, // kICdelayCorr
198   kFALSE, // kStripDelayCorr
199   kFALSE, // kIndexDelayCorr
200   kFALSE, // kTimeSlewingCorr
201 };
202
203 //__________________________________________________________________________
204
205 Bool_t AliTOFcalibHisto::fgFullCorrectionFlag[kNcorrections] = {
206   kFALSE, // kDDLBCcorr
207   kTRUE, // kAmphenolCableCorr
208   kTRUE, // kFlatCableCorr
209   kTRUE, // kInterfaceCardCorr
210   kTRUE, // kDDLdelayCorr
211   kTRUE, // kHPTDCdelayCorr
212   kTRUE, // kFEAchDelayCorr
213   kTRUE, // kFEAdelayCorr
214   kTRUE, // kTRMdelayCorr
215   kFALSE, // kICdelayCorr
216   kTRUE, // kStripDelayCorr
217   kTRUE, // kIndexDelayCorr
218   kTRUE, // kTimeSlewingCorr
219 };
220
221 //__________________________________________________________________________
222
223 const Int_t AliTOFcalibHisto::fgkStripStat[18][91] = {
224
225   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
226    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
227    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
228    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
229    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S00 */
230   
231   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
232    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
233    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
234    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
235    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S01 */
236
237   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
238    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
239    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
240    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
241    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S02 */
242   
243   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
244    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
245    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 0, 0, 0, 
246    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
247    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S03 */
248   
249   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
250    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
251    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
252    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
253    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S04 */
254
255   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
256    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
257    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
258    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
259    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S05 */
260
261   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
262    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
263    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
264    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
265    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S06 */
266   
267   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
268    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
269    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
270    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
271    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S07 */
272   
273   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
274    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
275    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
276    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
277    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S08 */
278   
279   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
280    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
281    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
282    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
283    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S09 */
284   
285   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
286    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
287    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
288    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
289    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S10 */
290   
291   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
292    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
293    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
294    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
295    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S11 */
296   
297   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
298    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
299    1, 1, 1, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 1, 
300    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
301    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S12 */
302
303   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
304    1, 1, 1, 1, 1, 1, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
305    0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 
306    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
307    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S13 */
308
309   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
310    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
311    0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 
312    0, 0, 0, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
313    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S14 */
314   
315   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
316    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
317    0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 
318    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
319    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S15 */
320
321   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
322    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
323    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
324    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
325    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}, /* S16 */
326   
327   {1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
328    1, 1, 1, 1, 1, 1, 1, 1, 0, 0, 0, 0, 1, 1, 1, 1, 1, 1, 1,
329    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
330    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 
331    1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}  /* S17 */
332
333 };
334
335 //__________________________________________________________________________
336
337 AliTOFcalibHisto::AliTOFcalibHisto() :
338   TObject(),
339   fCalibConst(),
340   fCalibMap(),
341   fCalibPar(),
342   fCalibStat()
343 {
344   /* default constructor */
345 }
346
347 //__________________________________________________________________________
348
349 AliTOFcalibHisto::~AliTOFcalibHisto()
350 {
351   /* default destructor */
352 }
353
354 //__________________________________________________________________________
355
356 void 
357 AliTOFcalibHisto::LoadHisto(TFile * const file, TH1D **histo, const Char_t *name) 
358 {
359   /* load histo */
360   *histo = (TH1D *)file->Get(name);
361   if (!*histo)
362     AliWarning(Form("error while getting %s histo", name));
363 }
364
365 //__________________________________________________________________________
366
367 void 
368 AliTOFcalibHisto::CreateHisto(TH1D **histo, const Char_t *name, Int_t size) 
369 {
370   /* create histo */
371   *histo = new TH1D(name, Form(";index;%s", name), size, 0, size);
372   if (!*histo)
373     AliWarning(Form("error while creating %s histo", name));
374 }
375
376 //__________________________________________________________________________
377
378 void 
379 AliTOFcalibHisto::WriteHisto(TFile *file, TH1D *histo) 
380 {
381   /* write histo */
382   if (!file || !file->IsOpen() || !histo)
383     return;
384   file->cd(); 
385   histo->Write();
386 }
387
388 //__________________________________________________________________________
389
390 void
391 AliTOFcalibHisto::SetHisto(TH1D *histo, Int_t index, Double_t value)
392 {
393   /* set histo */
394   if (!histo)
395     return;
396   histo->SetBinContent(index + 1, value);
397 }
398
399 //__________________________________________________________________________
400
401 Double_t
402 AliTOFcalibHisto::GetHisto(TH1D *histo, Int_t index)
403 {
404   /* get histo */
405   if (!histo) {
406     AliWarning("cannot get histo");
407     return 0.;
408   }
409   return histo->GetBinContent(index + 1);
410 }
411
412 //__________________________________________________________________________
413
414 void
415 AliTOFcalibHisto::LoadCalibHisto()
416 {
417   /* load calib histo */
418
419   if (fgCalibHistoFile && fgCalibHistoFile->IsOpen())
420     AliWarning("calib histo file already open: reloading"); 
421
422   /* open input file */
423   TFile *fileIn = TFile::Open(GetCalibHistoFileName());
424   if (!fileIn || !fileIn->IsOpen()) {
425     AliFatal(Form("cannot open input file %s", GetCalibHistoFileName()));
426     return;
427   }
428
429   /* set calib histo file */
430   fgCalibHistoFile = fileIn;
431
432   /* load consts */
433   for (Int_t iConst = 0; iConst < kNcalibConsts; iConst++)
434     LoadHisto(fileIn, &fCalibConst[iConst], fgkCalibConstName[iConst].Data());
435   /* load maps */
436   for (Int_t iMap = 0; iMap < kNcalibMaps; iMap++)
437     LoadHisto(fileIn, &fCalibMap[iMap], fgkCalibMapName[iMap].Data());
438 }
439
440 //__________________________________________________________________________
441
442 void
443 AliTOFcalibHisto::LoadCalibPar()
444 {
445   /* load calib par */
446
447   if (fgCalibParFile && fgCalibParFile->IsOpen())
448     AliWarning("calib par file already open: reloading"); 
449
450   /* load calib histo */
451   LoadCalibHisto();
452
453   /* open input file */
454   TFile *fileIn = TFile::Open(GetCalibParFileName());
455   if (!fileIn || !fileIn->IsOpen()) {
456     AliError(Form("cannot open input file %s", GetCalibParFileName()));
457     return;
458   }
459
460   /* set calib par file */
461   fgCalibParFile = fileIn;
462
463   /* load pars */
464   for (Int_t i = 0; i < kNcalibPars; i++)
465     LoadHisto(fileIn, &fCalibPar[i], fgkCalibParName[i].Data());
466
467 }
468
469 //__________________________________________________________________________
470
471 void
472 AliTOFcalibHisto::LoadCalibStat()
473 {
474   /* load calib stat */
475
476   if (fgCalibStatFile && fgCalibStatFile->IsOpen())
477     AliWarning("calib par file already open: reloading"); 
478
479   /* load calib histo */
480   LoadCalibHisto();
481
482   /* open input file */
483   TFile *fileIn = TFile::Open(GetCalibStatFileName());
484   if (!fileIn || !fileIn->IsOpen()) {
485     AliError(Form("cannot open input file %s", GetCalibStatFileName()));
486     return;
487   }
488
489   /* set calib par file */
490   fgCalibStatFile = fileIn;
491
492   /* load pars */
493   for (Int_t i = 0; i < kNcalibStats; i++)
494     LoadHisto(fileIn, &fCalibStat[i], fgkCalibStatName[i].Data());
495
496 }
497
498 //__________________________________________________________________________
499
500 void
501 AliTOFcalibHisto::WriteCalibHisto()
502 {
503   /* write calib histo */
504
505   /* open output file */
506   TFile *fileOut = TFile::Open(GetCalibHistoFileName(), "RECREATE");
507   if (!fileOut || !fileOut->IsOpen()) {
508     AliFatal(Form("cannot open output file %s", GetCalibHistoFileName()));
509     return;
510   }
511
512   /* create consts */
513   for (Int_t iConst = 0; iConst < kNcalibConsts; iConst++)
514     CreateHisto(&fCalibConst[iConst], fgkCalibConstName[iConst].Data(), 1);
515   /* create maps */
516   for (Int_t iMap = 0; iMap < kNcalibMaps; iMap++)
517     if (iMap == kIndex)
518       CreateHisto(&fCalibMap[iMap], fgkCalibMapName[iMap].Data(), fgkNchannelsEO);
519     else
520       CreateHisto(&fCalibMap[iMap], fgkCalibMapName[iMap].Data(), fgkNchannels);
521
522   /*** SETUP CONSTS ***/
523
524   SetHisto(fCalibConst[kLHCperiod], 0, fgkLHCperiod);
525   SetHisto(fCalibConst[kAmphenolCableDelay], 0, fgkAmphenolCableDelay);
526   SetHisto(fCalibConst[kFlatCableDelay], 0, fgkFlatCableDelay);
527   SetHisto(fCalibConst[kInterfaceCardDelay], 0, fgkInterfaceCardDelay);
528   
529   /***  SETUP MAPS  ***/
530
531   AliTOFRawStream rawStream;
532   Int_t indexEO, det[5], dummy, index, sector, plate, strip, sectorStrip, padz, padx, pad, icIndex;
533
534   /* temporarly disable warnings */
535   AliLog::EType_t logLevel = (AliLog::EType_t)AliLog::GetGlobalLogLevel();
536   AliLog::SetGlobalLogLevel(AliLog::kError);
537
538   /* loop over electronics oriented (EO) indices */
539   for (Int_t ddl = 0; ddl < 72; ddl++)
540     for (Int_t trm = 0; trm < 10; trm++)
541       for (Int_t chain = 0; chain < 2; chain++)
542         for (Int_t tdc = 0; tdc < 15; tdc++)
543           for (Int_t channel = 0; channel < 8; channel++) {
544             
545             /* compute index EO */
546             indexEO = GetIndexEO(ddl, trm, chain, tdc, channel);
547
548             /* convert EO indices into detector oriented (DO) indices
549                (this call causes some warnings because the loop includes
550                EO indices which are not connected to physical channels) */
551             rawStream.EquipmentId2VolumeId(ddl, trm + 3, chain, tdc, channel, det);
552             
553             /* swap det[3] and det[4] */
554             dummy = det[3]; det[3] = det[4]; det[4] = dummy;
555             
556             /* check detector indices */
557             if (det[0] < 0 || det[0] > 17 ||
558                 det[1] < 0 || det[1] > 4 ||
559                 det[2] < 0 || det[2] > 18 ||
560                 det[3] < 0 || det[3] > 1 ||
561                 det[4] < 0 || det[4] > 47) {
562               SetHisto(fCalibMap[kIndex], indexEO, -1);
563               continue;
564             }
565             
566             /* setup information */
567             index = AliTOFGeometry::GetIndex(det);
568             sector = det[0];
569             plate = det[1];
570             strip = det[2];
571             sectorStrip = plate < 3 ? plate * 19 + strip : plate * 19 - 4 + strip;
572             padz = det[3];
573             padx = det[4];
574             pad = padz + 2 * padx;
575             icIndex = pad < 48 ? pad : 95 - pad;
576
577             /* set maps */
578
579             /* main index */
580             SetHisto(fCalibMap[kIndex], indexEO, index);
581             /* EO indices */
582             SetHisto(fCalibMap[kDDL], index, ddl);
583             SetHisto(fCalibMap[kTRM], index, trm);
584             SetHisto(fCalibMap[kChain], index, chain);
585             SetHisto(fCalibMap[kTDC], index, tdc);
586             SetHisto(fCalibMap[kChannel], index, channel);
587             /* DO indices */
588             SetHisto(fCalibMap[kSector], index, sector);
589             SetHisto(fCalibMap[kPlate], index, plate);
590             SetHisto(fCalibMap[kStrip], index, strip);
591             SetHisto(fCalibMap[kSectorStrip], index, sectorStrip);
592             SetHisto(fCalibMap[kPadZ], index, padz);
593             SetHisto(fCalibMap[kPadX], index, padx);
594             SetHisto(fCalibMap[kPad], index, pad);
595             SetHisto(fCalibMap[kInterfaceCardIndex], index, icIndex);
596             /* calib constants */
597             SetHisto(fCalibMap[kDDLBCshift], index, fgkDDLBCshift[ddl]);
598             SetHisto(fCalibMap[kFlatCableLength], index, fgkFlatCableLength[sectorStrip]);
599             SetHisto(fCalibMap[kInterfaceCardLength], index, fgkInterfaceCardLength[icIndex]);
600             SetHisto(fCalibMap[kAmphenolCableLength], index, AliTOFCableLengthMap::GetCableLength(ddl, trm + 3, chain, tdc));
601             
602           } /* loop over electronics oriented (EO) indices */
603
604   /* re-enable warnings */
605   AliLog::SetGlobalLogLevel(logLevel);
606
607   /* write consts */
608   for (Int_t iConst = 0; iConst < kNcalibConsts; iConst++)
609     WriteHisto(fileOut, fCalibConst[iConst]);
610   /* write maps */
611   for (Int_t iMap = 0; iMap < kNcalibMaps; iMap++)
612     WriteHisto(fileOut, fCalibMap[iMap]);
613
614   /* close output file */
615   fileOut->Close();
616 }
617
618 //__________________________________________________________________________
619
620 void
621 AliTOFcalibHisto::WriteCalibStat()
622 {
623   /* write calib stat */
624
625   /* open output file */
626   TFile *fileOut = TFile::Open(GetCalibStatFileName(), "RECREATE");
627   if (!fileOut || !fileOut->IsOpen()) {
628     AliFatal(Form("cannot open output file %s", GetCalibStatFileName()));
629     return;
630   }
631
632   /* create stats */
633   for (Int_t iStat = 0; iStat < kNcalibStats; iStat++)
634     CreateHisto(&fCalibStat[iStat], fgkCalibStatName[iStat].Data(), fgkNchannels);
635
636   /***  SETUP STATS  ***/
637
638   Int_t sector, sectorStrip;
639
640   /* load calib histo */
641   LoadCalibHisto();
642
643   /* loop over channels */
644   for (Int_t index = 0; index < fgkNchannels; index++) {
645     sector = (Int_t)GetCalibMap(kSector, index);
646     sectorStrip = (Int_t)GetCalibMap(kSectorStrip, index);
647     /* strip stat */
648     SetHisto(fCalibStat[kStripStat], index, (Double_t)fgkStripStat[sector][sectorStrip]);
649   } /* loop over channels */
650
651   /* write stats */
652   for (Int_t iStat = 0; iStat < kNcalibStats; iStat++)
653     WriteHisto(fileOut, fCalibStat[iStat]);
654
655   /* close output file */
656   fileOut->Close();
657 }
658
659 //__________________________________________________________________________
660
661 Double_t
662 AliTOFcalibHisto::GetCorrection(Int_t corr, Int_t index, Double_t tot)
663 {
664   /* apply correction */
665
666   Int_t ddl, chain, tdc, channel, hptdc, pbCh, feaIndex, sector, plate, padx, trm, icIndex, sectorStrip;
667   Double_t slewing;
668   
669   switch (corr) {
670   case kDDLBCcorr:
671     return -GetCalibConst(kLHCperiod) * GetCalibMap(kDDLBCshift, index);
672   case kAmphenolCableCorr:
673     return GetCalibConst(kAmphenolCableDelay) * GetCalibMap(kAmphenolCableLength, index);
674   case kFlatCableCorr:
675     return GetCalibConst(kFlatCableDelay) * GetCalibMap(kFlatCableLength, index);
676   case kInterfaceCardCorr:
677     return GetCalibConst(kInterfaceCardDelay) * GetCalibMap(kInterfaceCardLength, index);
678   case kDDLdelayCorr:
679     ddl = (Int_t)GetCalibMap(kDDL, index);
680     return GetCalibPar(kDDLdelayPar, ddl);
681   case kHPTDCdelayCorr:
682     chain = (Int_t)GetCalibMap(kChain, index);
683     tdc = (Int_t)GetCalibMap(kTDC, index);
684     hptdc = tdc + 15 * chain;
685     return GetCalibPar(kHPTDCdelayPar, hptdc);
686   case kFEAchDelayCorr:
687     ddl = (Int_t)GetCalibMap(kDDL, index);
688     tdc = (Int_t)GetCalibMap(kTDC, index);
689     channel = (Int_t)GetCalibMap(kChannel, index);
690     pbCh = channel + 8 * (tdc % 3);
691     if (ddl % 2 == 0)
692       return GetCalibPar(kRightFEAchDelayPar, pbCh);
693     else
694       return GetCalibPar(kLeftFEAchDelayPar, pbCh);
695   case kFEAdelayCorr:
696     sector = (Int_t)GetCalibMap(kSector, index);
697     plate = (Int_t)GetCalibMap(kPlate, index);
698     sectorStrip = (Int_t)GetCalibMap(kSectorStrip, index);
699     padx = (Int_t)GetCalibMap(kPadX, index);
700     feaIndex = padx / 12 + 4 * sectorStrip + 364 * sector;      
701     return GetCalibPar(kFEAdelayPar, feaIndex);
702   case kTRMdelayCorr:
703     ddl = (Int_t)GetCalibMap(kDDL, index);
704     trm = (Int_t)GetCalibMap(kTRM, index);
705     return GetCalibPar(kTRMdelayPar, trm + 10 * ddl);
706   case kICdelayCorr:
707     icIndex = (Int_t)GetCalibMap(kInterfaceCardIndex, index);
708     return GetCalibPar(kICdelayPar, icIndex);
709   case kStripDelayCorr:
710     sector = (Int_t)GetCalibMap(kSector, index);
711     sectorStrip = (Int_t)GetCalibMap(kSectorStrip, index);
712     return GetCalibPar(kStripDelayPar, sectorStrip + 91 * sector);
713   case kIndexDelayCorr:
714     return GetCalibPar(kIndexDelayPar, index);
715   case kTimeSlewingCorr:
716     tot = tot < SLEW_TOTMIN ? SLEW_TOTMIN : tot;
717     tot = tot > SLEW_TOTMAX ? SLEW_TOTMAX : tot;
718     slewing = 0.;
719     for (Int_t i = 0; i < fCalibPar[kTimeSlewingPar]->GetNbinsX(); i++)
720       slewing += GetCalibPar(kTimeSlewingPar, i) * TMath::Power(tot, i);
721     return slewing;
722   default:
723     AliWarning(Form("unknown correction flag (%d)", corr));
724     return 0.;
725   }
726 }
727
728 //__________________________________________________________________________
729
730 Double_t
731 AliTOFcalibHisto::GetNominalCorrection(Int_t index)
732 {
733   /* get nominal correction */
734   Double_t corr = 0;
735   for (Int_t iCorr = 0; iCorr < kNcorrections; iCorr++)
736     corr += GetCorrection(iCorr, index);
737   return corr;
738 }
739
740 //__________________________________________________________________________
741
742 void
743 AliTOFcalibHisto::ApplyNominalCorrection(AliESDtrack *track)
744 {
745   /* apply nominal correction */
746   
747   Double_t rawTime = track->GetTOFsignalRaw();
748   Int_t index = track->GetTOFCalChannel();
749   Double_t time = rawTime - 1.e3 * GetNominalCorrection(index);
750   track->SetTOFsignal(time);
751 }
752
753 //__________________________________________________________________________
754
755 Double_t
756 AliTOFcalibHisto::GetCableCorrection(Int_t index)
757 {
758   /* get cable correction */
759   Double_t corr = 0;
760   for (Int_t iCorr = 0; iCorr < kNcorrections; iCorr++)
761     if (fgCableCorrectionFlag[iCorr])
762       corr += GetCorrection(iCorr, index);
763   return corr;
764 }
765
766 //__________________________________________________________________________
767
768 Double_t
769 AliTOFcalibHisto::GetFullCorrection(Int_t index, Double_t tot)
770 {
771   /* get full correction */
772   Double_t corr = 0;
773   for (Int_t iCorr = 0; iCorr < kNcorrections; iCorr++)
774     if (fgFullCorrectionFlag[iCorr]) {
775       corr += GetCorrection(iCorr, index, tot);
776     }
777   return corr;
778 }
779
780 //__________________________________________________________________________
781
782 Bool_t
783 AliTOFcalibHisto::GetStatus(Int_t index)
784 {
785   /* get status */
786
787   Bool_t status = kTRUE;
788
789   for (Int_t iStat = 0; iStat < kNcalibStats; iStat++)
790     status &= GetCalibStat(iStat, index);
791
792   return status;
793 }