]>
Commit | Line | Data |
---|---|---|
b3fcfd96 | 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 | /* | |
4c865c34 | 19 | Responsible: Raphaelle Bailhache (rbailhache@ikf.uni-frankfurt.de) |
20 | Code to analyze the TRD calibration and to produce OCDB entries | |
b3fcfd96 | 21 | |
22 | ||
23 | .x ~/rootlogon.C | |
24 | gSystem->Load("libANALYSIS"); | |
4c865c34 | 25 | gSystem->Load("libTRDcalib"); |
b3fcfd96 | 26 | |
27 | AliTRDPreprocessorOffline proces; | |
28 | TString ocdbPath="local:////" | |
29 | ocdbPath+=gSystem->GetFromPipe("pwd"); | |
30 | ||
31 | proces.CalibTimeGain("CalibObjects.root",run0,run1,ocdbPath); | |
32 | proces.CalibTimeVdrift("CalibObjects.root",run0,run1,ocdbPath); | |
33 | // take the raw calibration data from the file CalibObjects.root | |
34 | // and make a OCDB entry with run validity run0-run1 | |
35 | // results are stored at the ocdbPath - local or alien ... | |
36 | // default storage ""- data stored at current working directory | |
37 | ||
38 | */ | |
00d203b6 | 39 | #include "AliLog.h" |
b3fcfd96 | 40 | #include "Riostream.h" |
41 | #include <fstream> | |
42 | #include "TFile.h" | |
43 | #include "TCanvas.h" | |
44 | #include "TLegend.h" | |
45 | #include "TH2I.h" | |
46 | #include "TH1I.h" | |
47 | #include "TH2F.h" | |
48 | #include "TH1F.h" | |
82b413fd | 49 | #include "TMath.h" |
83d0cc79 | 50 | #include "THnSparse.h" |
b3fcfd96 | 51 | #include "TProfile2D.h" |
52 | #include "AliTRDCalDet.h" | |
53 | #include "AliTRDCalPad.h" | |
54 | #include "AliCDBMetaData.h" | |
55 | #include "AliCDBId.h" | |
56 | #include "AliCDBManager.h" | |
57 | #include "AliCDBStorage.h" | |
58 | #include "AliTRDCalibraMode.h" | |
59 | #include "AliTRDCalibraFit.h" | |
60 | #include "AliTRDCalibraVdriftLinearFit.h" | |
a0bb5615 | 61 | #include "AliTRDCalibraExbAltFit.h" |
b3fcfd96 | 62 | #include "AliTRDPreprocessorOffline.h" |
81a5aeca | 63 | #include "AliTRDCalChamberStatus.h" |
83d0cc79 | 64 | #include "AliTRDCalibChamberStatus.h" |
840ec79d | 65 | #include "AliTRDCommonParam.h" |
83d0cc79 | 66 | #include "AliCDBManager.h" |
67 | #include "AliCDBEntry.h" | |
b3fcfd96 | 68 | |
69 | ||
70 | ClassImp(AliTRDPreprocessorOffline) | |
71 | ||
54f2ff1c | 72 | AliTRDPreprocessorOffline::AliTRDPreprocessorOffline(): |
01239968 | 73 | TNamed("TPCPreprocessorOffline","TPCPreprocessorOffline"), |
b3fcfd96 | 74 | fMethodSecond(kTRUE), |
01239968 | 75 | fNameList("TRDCalib"), |
76 | fCalDetGainUsed(0x0), | |
4c865c34 | 77 | fCalDetVdriftUsed(0x0), |
840ec79d | 78 | fCalDetExBUsed(0x0), |
b3fcfd96 | 79 | fCH2d(0x0), |
80 | fPH2d(0x0), | |
81 | fPRF2d(0x0), | |
83d0cc79 | 82 | fSparse(0x0), |
b3fcfd96 | 83 | fAliTRDCalibraVdriftLinearFit(0x0), |
a0bb5615 | 84 | fAliTRDCalibraExbAltFit(0x0), |
b3fcfd96 | 85 | fNEvents(0x0), |
86 | fAbsoluteGain(0x0), | |
a0bb5615 | 87 | fPlots(new TObjArray(kNumCalibObjs)), |
88 | fCalibObjects(new TObjArray(kNumCalibObjs)), | |
83d0cc79 | 89 | fFirstRunGainUsed(0), |
4c865c34 | 90 | fVersionGainUsed(0), |
91 | fSubVersionGainUsed(0), | |
ca7e6e64 | 92 | fFirstRunVdriftUsed(0), |
4c865c34 | 93 | fVersionVdriftUsed(0), |
00d203b6 | 94 | fSubVersionVdriftUsed(0), |
83d0cc79 | 95 | fFirstRunExBUsed(0), |
96 | fVersionExBUsed(0), | |
97 | fSubVersionExBUsed(0), | |
98 | fNoExBUsedInReco(kFALSE), | |
00d203b6 | 99 | fSwitchOnValidation(kTRUE), |
100 | fVdriftValidated(kFALSE), | |
840ec79d | 101 | fExBValidated(kFALSE), |
a2a4ec8e | 102 | fT0Validated(kFALSE), |
103 | fMinStatsVdriftT0PH(800*20), | |
104 | fMinStatsVdriftLinear(800), | |
105 | fMinStatsGain(800), | |
54f2ff1c | 106 | fMinStatsPRF(600), |
83d0cc79 | 107 | fMinStatsChamberStatus(20), |
54f2ff1c | 108 | fBackCorrectGain(kFALSE), |
109 | fBackCorrectVdrift(kTRUE), | |
110 | fNotEnoughStatisticsForTheGain(kFALSE), | |
111 | fNotEnoughStatisticsForTheVdriftLinear(kFALSE), | |
82b413fd | 112 | fStatusNeg(0), |
83d0cc79 | 113 | fStatusPos(0), |
114 | fBadCalibValidate(20), | |
115 | fNoDataValidate(20), | |
116 | fRMSBadCalibratedGain(20.0), | |
117 | fRMSBadCalibratedVdrift(20.0), | |
118 | fRMSBadCalibratedExB(50.0) | |
b3fcfd96 | 119 | { |
120 | // | |
121 | // default constructor | |
122 | // | |
83d0cc79 | 123 | |
124 | memset(fBadCalib, 0, sizeof(Int_t) * 18); | |
125 | memset(fNoData, 0, sizeof(Int_t) * 18); | |
b3fcfd96 | 126 | } |
00d203b6 | 127 | //_________________________________________________________________________________________________________________ |
b3fcfd96 | 128 | AliTRDPreprocessorOffline::~AliTRDPreprocessorOffline() { |
129 | // | |
130 | // Destructor | |
131 | // | |
132 | ||
01239968 | 133 | if(fCalDetGainUsed) delete fCalDetGainUsed; |
4c865c34 | 134 | if(fCalDetVdriftUsed) delete fCalDetVdriftUsed; |
840ec79d | 135 | if(fCalDetExBUsed) delete fCalDetExBUsed; |
b3fcfd96 | 136 | if(fCH2d) delete fCH2d; |
137 | if(fPH2d) delete fPH2d; | |
138 | if(fPRF2d) delete fPRF2d; | |
83d0cc79 | 139 | if(fSparse) delete fSparse; |
b3fcfd96 | 140 | if(fAliTRDCalibraVdriftLinearFit) delete fAliTRDCalibraVdriftLinearFit; |
a0bb5615 | 141 | if(fAliTRDCalibraExbAltFit) delete fAliTRDCalibraExbAltFit; |
b3fcfd96 | 142 | if(fNEvents) delete fNEvents; |
143 | if(fAbsoluteGain) delete fAbsoluteGain; | |
144 | if(fPlots) delete fPlots; | |
145 | if(fCalibObjects) delete fCalibObjects; | |
146 | ||
83d0cc79 | 147 | } |
148 | //___________________________________________________________________________________ | |
149 | void AliTRDPreprocessorOffline::Process(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage) | |
150 | { | |
151 | // | |
152 | // Process to the gain, vdrift, timeoffset, exb and chamber status calibration | |
153 | // | |
154 | ||
155 | if(SetCalDetGain(startRunNumber,fVersionGainUsed,fSubVersionGainUsed) && SetCalDetVdriftExB(startRunNumber,fVersionVdriftUsed,fSubVersionVdriftUsed,fVersionExBUsed,fSubVersionExBUsed)) { | |
156 | ||
157 | CalibVdriftT0(file,startRunNumber,endRunNumber,ocdbStorage); | |
158 | CalibGain(file,startRunNumber,endRunNumber,ocdbStorage); | |
159 | CalibChamberStatus(file,startRunNumber,endRunNumber,ocdbStorage); | |
160 | ||
161 | } | |
162 | ||
163 | PrintStatus(); | |
164 | ||
b3fcfd96 | 165 | } |
00d203b6 | 166 | //___________________________________________________________________________________________________________________ |
b3fcfd96 | 167 | |
168 | void AliTRDPreprocessorOffline::CalibVdriftT0(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){ | |
169 | // | |
170 | // make calibration of the drift velocity | |
171 | // Input parameters: | |
172 | // file - the location of input file | |
173 | // startRunNumber, endRunNumber - run validity period | |
174 | // ocdbStorage - path to the OCDB storage | |
175 | // - if empty - local storage 'pwd' uesed | |
176 | if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB"; | |
177 | // | |
178 | // 1. Initialization | |
179 | // | |
008817a3 | 180 | fVdriftValidated = kTRUE; |
181 | fT0Validated = kTRUE; | |
840ec79d | 182 | fExBValidated = kTRUE; |
54f2ff1c | 183 | fNotEnoughStatisticsForTheVdriftLinear = kFALSE; |
b3fcfd96 | 184 | // |
185 | // 2. extraction of the information | |
186 | // | |
840ec79d | 187 | if(ReadVdriftLinearFitGlobal(file) && fCalDetVdriftUsed && fCalDetExBUsed) AnalyzeVdriftLinearFit(); |
54f2ff1c | 188 | if(ReadVdriftT0Global(file)) AnalyzeVdriftT0(); |
b3fcfd96 | 189 | // |
190 | // 3. Append QA plots | |
191 | // | |
192 | //MakeDefaultPlots(fVdriftArray,fVdriftArray); | |
193 | // | |
194 | // | |
00d203b6 | 195 | // 4. validate OCDB entries |
b3fcfd96 | 196 | // |
00d203b6 | 197 | if(fSwitchOnValidation==kTRUE && ValidateVdrift()==kFALSE) { |
82b413fd | 198 | //AliError("TRD vdrift OCDB parameters out of range!"); |
00d203b6 | 199 | fVdriftValidated = kFALSE; |
200 | } | |
201 | if(fSwitchOnValidation==kTRUE && ValidateT0()==kFALSE) { | |
82b413fd | 202 | //AliError("TRD t0 OCDB parameters out of range!"); |
00d203b6 | 203 | fT0Validated = kFALSE; |
204 | } | |
840ec79d | 205 | if(fSwitchOnValidation==kTRUE && ValidateExB()==kFALSE) { |
206 | //AliError("TRD t0 OCDB parameters out of range!"); | |
207 | fExBValidated = kFALSE; | |
208 | } | |
b3fcfd96 | 209 | // |
00d203b6 | 210 | // 5. update of OCDB |
211 | // | |
212 | // | |
213 | if(fVdriftValidated) UpdateOCDBVdrift(startRunNumber,endRunNumber,ocdbStorage); | |
214 | if(fT0Validated) UpdateOCDBT0(startRunNumber,endRunNumber,ocdbStorage); | |
840ec79d | 215 | if(fExBValidated) UpdateOCDBExB(startRunNumber,endRunNumber,ocdbStorage); |
00d203b6 | 216 | |
b3fcfd96 | 217 | } |
a0bb5615 | 218 | //___________________________________________________________________________________________________________________ |
219 | ||
220 | void AliTRDPreprocessorOffline::CalibExbAlt(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){ | |
221 | // | |
222 | // make calibration of the drift velocity | |
223 | // Input parameters: | |
224 | // file - the location of input file | |
225 | // startRunNumber, endRunNumber - run validity period | |
226 | // ocdbStorage - path to the OCDB storage | |
227 | // - if empty - local storage 'pwd' uesed | |
228 | if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB"; | |
229 | // | |
230 | // 1. Initialization | |
231 | // | |
232 | ||
233 | // | |
234 | // 2. extraction of the information | |
235 | // | |
236 | if(ReadExbAltFitGlobal(file)) AnalyzeExbAltFit(); | |
237 | // | |
238 | // 3. Append QA plots | |
239 | // | |
240 | //MakeDefaultPlots(fVdriftArray,fVdriftArray); | |
241 | // | |
242 | // | |
243 | // 4. validate OCDB entries | |
244 | // | |
245 | // | |
246 | // 5. update of OCDB | |
247 | // | |
248 | // | |
249 | UpdateOCDBExBAlt(startRunNumber,endRunNumber,ocdbStorage); | |
250 | ||
251 | } | |
252 | ||
00d203b6 | 253 | //_________________________________________________________________________________________________________________ |
b3fcfd96 | 254 | |
255 | void AliTRDPreprocessorOffline::CalibGain(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){ | |
256 | // | |
257 | // make calibration of the drift velocity | |
258 | // Input parameters: | |
259 | // file - the location of input file | |
260 | // startRunNumber, endRunNumber - run validity period | |
261 | // ocdbStorage - path to the OCDB storage | |
262 | // - if empty - local storage 'pwd' uesed | |
263 | if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB"; | |
264 | // | |
54f2ff1c | 265 | fNotEnoughStatisticsForTheGain = kFALSE; |
266 | // | |
b3fcfd96 | 267 | // 1. Initialization |
268 | if(!ReadGainGlobal(file)) return; | |
269 | // | |
270 | // | |
271 | // 2. extraction of the information | |
272 | // | |
273 | AnalyzeGain(); | |
54f2ff1c | 274 | if(fBackCorrectGain) CorrectFromDetGainUsed(); |
a0bb5615 | 275 | //if(fBackCorrectVdrift) CorrectFromDetVdriftUsed(); |
b3fcfd96 | 276 | // |
277 | // 3. Append QA plots | |
278 | // | |
279 | //MakeDefaultPlots(fVdriftArray,fVdriftArray); | |
280 | // | |
281 | // | |
00d203b6 | 282 | // 4. validate OCDB entries |
283 | // | |
284 | if(fSwitchOnValidation==kTRUE && ValidateGain()==kFALSE) { | |
82b413fd | 285 | //AliError("TRD gain OCDB parameters out of range!"); |
00d203b6 | 286 | return; |
287 | } | |
b3fcfd96 | 288 | // |
00d203b6 | 289 | // 5. update of OCDB |
b3fcfd96 | 290 | // |
00d203b6 | 291 | // |
292 | if((!fCalDetVdriftUsed) || (fCalDetVdriftUsed && fVdriftValidated)) UpdateOCDBGain(startRunNumber,endRunNumber,ocdbStorage); | |
293 | ||
b3fcfd96 | 294 | |
295 | } | |
00d203b6 | 296 | //________________________________________________________________________________________________________________ |
b3fcfd96 | 297 | |
298 | void AliTRDPreprocessorOffline::CalibPRF(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){ | |
299 | // | |
300 | // make calibration of the drift velocity | |
301 | // Input parameters: | |
302 | // file - the location of input file | |
303 | // startRunNumber, endRunNumber - run validity period | |
304 | // ocdbStorage - path to the OCDB storage | |
305 | // - if empty - local storage 'pwd' uesed | |
306 | if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB"; | |
307 | // | |
308 | // 1. Initialization | |
309 | if(!ReadPRFGlobal(file)) return; | |
310 | // | |
311 | // | |
312 | // 2. extraction of the information | |
313 | // | |
314 | AnalyzePRF(); | |
315 | // | |
316 | // 3. Append QA plots | |
317 | // | |
318 | //MakeDefaultPlots(fVdriftArray,fVdriftArray); | |
319 | // | |
320 | // | |
00d203b6 | 321 | // |
322 | // 4. validate OCDB entries | |
323 | // | |
324 | if(fSwitchOnValidation==kTRUE && ValidatePRF()==kFALSE) { | |
82b413fd | 325 | //AliError("TRD prf OCDB parameters out of range!"); |
00d203b6 | 326 | return; |
327 | } | |
328 | // | |
329 | // 5. update of OCDB | |
b3fcfd96 | 330 | // |
331 | // | |
332 | UpdateOCDBPRF(startRunNumber,endRunNumber,ocdbStorage); | |
333 | ||
a5dcf618 | 334 | } |
335 | //________________________________________________________________________________________________________________ | |
336 | ||
83d0cc79 | 337 | void AliTRDPreprocessorOffline::CalibChamberStatus(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){ |
a5dcf618 | 338 | // |
339 | // make calibration of the chamber status | |
340 | // Input parameters: | |
341 | // startRunNumber, endRunNumber - run validity period | |
342 | // ocdbStorage - path to the OCDB storage | |
343 | // - if empty - local storage 'pwd' uesed | |
344 | if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB"; | |
345 | // | |
346 | // | |
83d0cc79 | 347 | // 1. Initialization |
348 | if(!ReadStatusGlobal(file)) return; | |
349 | // | |
350 | // | |
351 | // | |
a5dcf618 | 352 | // 2. extraction of the information |
353 | // | |
354 | if(!AnalyzeChamberStatus()) return; | |
355 | // | |
356 | // 3. Append QA plots | |
357 | // | |
358 | //MakeDefaultPlots(fVdriftArray,fVdriftArray); | |
359 | // | |
360 | // | |
361 | // | |
362 | // 4. validate OCDB entries | |
363 | // | |
364 | if(fSwitchOnValidation==kTRUE && ValidateChamberStatus()==kFALSE) { | |
82b413fd | 365 | //AliError("TRD Chamber status OCDB parameters not ok!"); |
a5dcf618 | 366 | return; |
367 | } | |
368 | // | |
369 | // 5. update of OCDB | |
370 | // | |
371 | // | |
83d0cc79 | 372 | if((!fNotEnoughStatisticsForTheVdriftLinear) && (!fNotEnoughStatisticsForTheGain)) UpdateOCDBChamberStatus(startRunNumber,endRunNumber,ocdbStorage); |
373 | //UpdateOCDBChamberStatus(startRunNumber,endRunNumber,ocdbStorage); | |
a5dcf618 | 374 | |
b3fcfd96 | 375 | } |
00d203b6 | 376 | //______________________________________________________________________________________________________ |
4c865c34 | 377 | Bool_t AliTRDPreprocessorOffline::Init(const Char_t* fileName){ |
378 | // | |
379 | // read the calibration used during the reconstruction | |
380 | // | |
381 | ||
382 | if(ReadVdriftT0Global(fileName)) { | |
383 | ||
384 | TString nameph = fPH2d->GetTitle(); | |
ca7e6e64 | 385 | fFirstRunVdriftUsed = GetFirstRun(nameph); |
4c865c34 | 386 | fVersionVdriftUsed = GetVersion(nameph); |
387 | fSubVersionVdriftUsed = GetSubVersion(nameph); | |
388 | ||
389 | //printf("Found Version %d, Subversion %d for vdrift\n",fVersionVdriftUsed,fSubVersionVdriftUsed); | |
390 | ||
391 | } | |
392 | ||
393 | if(ReadGainGlobal(fileName)) { | |
394 | ||
395 | TString namech = fCH2d->GetTitle(); | |
83d0cc79 | 396 | fFirstRunGainUsed = GetFirstRun(namech); |
4c865c34 | 397 | fVersionGainUsed = GetVersion(namech); |
398 | fSubVersionGainUsed = GetSubVersion(namech); | |
399 | ||
400 | //printf("Found Version %d, Subversion %d for gain\n",fVersionGainUsed,fSubVersionGainUsed); | |
401 | ||
402 | } | |
83d0cc79 | 403 | |
404 | if(ReadVdriftLinearFitGlobal(fileName)) { | |
405 | ||
406 | TString namelinear = fAliTRDCalibraVdriftLinearFit->GetNameCalibUsed(); | |
407 | fFirstRunExBUsed = GetFirstRun(namelinear); | |
408 | fVersionExBUsed = GetVersion(namelinear); | |
409 | fSubVersionExBUsed = GetSubVersion(namelinear); | |
410 | ||
411 | //printf("Found Version %d, Subversion %d, run %d for ExB\n",fVersionExBUsed,fSubVersionExBUsed,fFirstRunExBUsed); | |
412 | ||
413 | } | |
4c865c34 | 414 | |
82b413fd | 415 | if(fVersionVdriftUsed == 0) fStatusPos = fStatusPos |kVdriftErrorOld; |
416 | if(fVersionGainUsed == 0) fStatusPos = fStatusPos | kGainErrorOld; | |
840ec79d | 417 | |
4c865c34 | 418 | return kTRUE; |
419 | ||
420 | } | |
00d203b6 | 421 | //___________________________________________________________________________________________________________________ |
b3fcfd96 | 422 | |
83d0cc79 | 423 | Bool_t AliTRDPreprocessorOffline::ReadStatusGlobal(const Char_t* fileName){ |
424 | // | |
425 | // read calibration entries from file | |
426 | // | |
427 | if(fSparse) return kTRUE; | |
428 | TFile fcalib(fileName); | |
429 | TObjArray * array = (TObjArray*)fcalib.Get(fNameList); | |
430 | if (array){ | |
431 | fSparse = (THnSparseI *) array->FindObject("NumberOfEntries"); | |
432 | if(!fSparse) return kFALSE; | |
433 | } | |
434 | else | |
435 | return kFALSE; | |
436 | ||
437 | return kTRUE; | |
438 | ||
439 | } | |
440 | //___________________________________________________________________________________________________________________ | |
441 | ||
b3fcfd96 | 442 | Bool_t AliTRDPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){ |
443 | // | |
444 | // read calibration entries from file | |
445 | // | |
4c865c34 | 446 | if(fCH2d) return kTRUE; |
b3fcfd96 | 447 | TFile fcalib(fileName); |
01239968 | 448 | TObjArray * array = (TObjArray*)fcalib.Get(fNameList); |
b3fcfd96 | 449 | if (array){ |
450 | TH2I *ch2d = (TH2I *) array->FindObject("CH2d"); | |
451 | if(!ch2d) return kFALSE; | |
452 | fCH2d = (TH2I*)ch2d->Clone(); | |
453 | //fNEvents = (TH1I *) array->FindObject("NEvents"); | |
454 | //fAbsoluteGain = (TH2F *) array->FindObject("AbsoluteGain"); | |
455 | }else{ | |
456 | TH2I *ch2d = (TH2I *) fcalib.Get("CH2d"); | |
457 | if(!ch2d) return kFALSE; | |
458 | fCH2d = (TH2I*)ch2d->Clone(); | |
459 | //fNEvents = (TH1I *) fcalib.Get("NEvents"); | |
460 | //fAbsoluteGain = (TH2F *) fcalib.Get("AbsoluteGain"); | |
461 | } | |
462 | fCH2d->SetDirectory(0); | |
463 | //printf("title of CH2d %s\n",fCH2d->GetTitle()); | |
464 | ||
465 | return kTRUE; | |
466 | ||
467 | } | |
00d203b6 | 468 | //_________________________________________________________________________________________________________________ |
b3fcfd96 | 469 | |
470 | Bool_t AliTRDPreprocessorOffline::ReadVdriftT0Global(const Char_t* fileName){ | |
471 | // | |
472 | // read calibration entries from file | |
473 | // | |
4c865c34 | 474 | if(fPH2d) return kTRUE; |
b3fcfd96 | 475 | TFile fcalib(fileName); |
01239968 | 476 | TObjArray * array = (TObjArray*)fcalib.Get(fNameList); |
b3fcfd96 | 477 | if (array){ |
478 | TProfile2D *ph2d = (TProfile2D *) array->FindObject("PH2d"); | |
479 | if(!ph2d) return kFALSE; | |
480 | fPH2d = (TProfile2D*)ph2d->Clone(); | |
481 | //fNEvents = (TH1I *) array->FindObject("NEvents"); | |
482 | }else{ | |
483 | TProfile2D *ph2d = (TProfile2D *) fcalib.Get("PH2d"); | |
484 | if(!ph2d) return kFALSE; | |
485 | fPH2d = (TProfile2D*)ph2d->Clone(); | |
486 | //fNEvents = (TH1I *) fcalib.Get("NEvents"); | |
487 | } | |
488 | fPH2d->SetDirectory(0); | |
489 | //printf("title of PH2d %s\n",fPH2d->GetTitle()); | |
490 | ||
491 | return kTRUE; | |
492 | ||
493 | } | |
00d203b6 | 494 | //___________________________________________________________________________________________________________________ |
b3fcfd96 | 495 | |
496 | Bool_t AliTRDPreprocessorOffline::ReadVdriftLinearFitGlobal(const Char_t* fileName){ | |
497 | // | |
498 | // read calibration entries from file | |
499 | // | |
4c865c34 | 500 | if(fAliTRDCalibraVdriftLinearFit) return kTRUE; |
b3fcfd96 | 501 | TFile fcalib(fileName); |
01239968 | 502 | TObjArray * array = (TObjArray*)fcalib.Get(fNameList); |
b3fcfd96 | 503 | if (array){ |
504 | fAliTRDCalibraVdriftLinearFit = (AliTRDCalibraVdriftLinearFit *) array->FindObject("AliTRDCalibraVdriftLinearFit"); | |
505 | //fNEvents = (TH1I *) array->FindObject("NEvents"); | |
506 | }else{ | |
507 | fAliTRDCalibraVdriftLinearFit = (AliTRDCalibraVdriftLinearFit *) fcalib.Get("AliTRDCalibraVdriftLinearFit"); | |
508 | //fNEvents = (TH1I *) fcalib.Get("NEvents"); | |
509 | } | |
510 | if(!fAliTRDCalibraVdriftLinearFit) { | |
511 | //printf("No AliTRDCalibraVdriftLinearFit\n"); | |
512 | return kFALSE; | |
513 | } | |
514 | return kTRUE; | |
515 | ||
a0bb5615 | 516 | } |
517 | //_____________________________________________________________________________________________________________ | |
518 | Bool_t AliTRDPreprocessorOffline::ReadExbAltFitGlobal(const Char_t* fileName){ | |
519 | // | |
520 | // read calibration entries from file | |
521 | // | |
522 | if(fAliTRDCalibraExbAltFit) return kTRUE; | |
523 | TFile fcalib(fileName); | |
524 | TObjArray * array = (TObjArray*)fcalib.Get(fNameList); | |
525 | if (array){ | |
526 | fAliTRDCalibraExbAltFit = (AliTRDCalibraExbAltFit *) array->FindObject("AliTRDCalibraExbAltFit"); | |
527 | //fNEvents = (TH1I *) array->FindObject("NEvents"); | |
528 | }else{ | |
529 | fAliTRDCalibraExbAltFit = (AliTRDCalibraExbAltFit *) fcalib.Get("AliTRDCalibraExbAltFit"); | |
530 | //fNEvents = (TH1I *) fcalib.Get("NEvents"); | |
531 | } | |
532 | if(!fAliTRDCalibraExbAltFit) { | |
533 | //printf("No AliTRDCalibraExbAltFit\n"); | |
534 | return kFALSE; | |
535 | } | |
536 | return kTRUE; | |
537 | ||
b3fcfd96 | 538 | } |
00d203b6 | 539 | //_____________________________________________________________________________________________________________ |
b3fcfd96 | 540 | |
541 | Bool_t AliTRDPreprocessorOffline::ReadPRFGlobal(const Char_t* fileName){ | |
542 | // | |
543 | // read calibration entries from file | |
544 | // | |
4c865c34 | 545 | if(fPRF2d) return kTRUE; |
b3fcfd96 | 546 | TFile fcalib(fileName); |
01239968 | 547 | TObjArray * array = (TObjArray*)fcalib.Get(fNameList); |
b3fcfd96 | 548 | if (array){ |
549 | TProfile2D *prf2d = (TProfile2D *) array->FindObject("PRF2d"); | |
550 | if(!prf2d) return kFALSE; | |
551 | fPRF2d = (TProfile2D*)prf2d->Clone(); | |
552 | //fNEvents = (TH1I *) array->FindObject("NEvents"); | |
553 | }else{ | |
554 | TProfile2D *prf2d = (TProfile2D *) fcalib.Get("PRF2d"); | |
555 | if(!prf2d) return kFALSE; | |
556 | fPRF2d = (TProfile2D*)prf2d->Clone(); | |
557 | //fNEvents = (TH1I *) fcalib.Get("NEvents"); | |
558 | } | |
559 | fPRF2d->SetDirectory(0); | |
560 | //printf("title of PRF2d %s\n",fPRF2d->GetTitle()); | |
561 | ||
562 | return kTRUE; | |
563 | ||
564 | } | |
00d203b6 | 565 | //__________________________________________________________________________________________________________ |
b3fcfd96 | 566 | |
567 | Bool_t AliTRDPreprocessorOffline::AnalyzeGain(){ | |
568 | // | |
569 | // Analyze gain - produce the calibration objects | |
570 | // | |
571 | ||
572 | AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance(); | |
a2a4ec8e | 573 | calibra->SetMinEntries(fMinStatsGain); // If there is less than 1000 entries in the histo: no fit |
b3fcfd96 | 574 | calibra->AnalyseCH(fCH2d); |
575 | ||
576 | Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(0)) | |
577 | + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(0)); | |
578 | Int_t nbfit = calibra->GetNumberFit(); | |
579 | Int_t nbE = calibra->GetNumberEnt(); | |
580 | ||
581 | ||
582 | Bool_t ok = kFALSE; | |
4c865c34 | 583 | Bool_t meanother = kFALSE; |
b3fcfd96 | 584 | // enough statistics |
585 | if ((nbtg > 0) && | |
586 | (nbfit >= 0.5*nbE) && (nbE > 30)) { | |
587 | // create the cal objects | |
54f2ff1c | 588 | if(!fBackCorrectGain) { |
4c865c34 | 589 | calibra->PutMeanValueOtherVectorFit(1,kTRUE); |
590 | meanother = kTRUE; | |
591 | } | |
b3fcfd96 | 592 | TObjArray object = calibra->GetVectorFit(); |
4c865c34 | 593 | AliTRDCalDet *calDetGain = calibra->CreateDetObjectGain(&object,meanother); |
b3fcfd96 | 594 | TH1F *coefGain = calDetGain->MakeHisto1DAsFunctionOfDet(); |
595 | // Put them in the array | |
596 | fCalibObjects->AddAt(calDetGain,kGain); | |
597 | fPlots->AddAt(coefGain,kGain); | |
598 | // | |
599 | ok = kTRUE; | |
600 | } | |
54f2ff1c | 601 | else { |
602 | fNotEnoughStatisticsForTheGain = kTRUE; | |
dee5f636 | 603 | Int_t minStatsGain = fMinStatsGain*30; |
604 | calibra->SetMinEntries(minStatsGain); // Because we do it for all, we increase this | |
54f2ff1c | 605 | Double_t gainoverallnotnormalized = calibra->AnalyseCHAllTogether(fCH2d); |
606 | if(fCalDetGainUsed && (gainoverallnotnormalized > 0.0)) { | |
607 | AliTRDCalDet *calDetGain = new AliTRDCalDet(*fCalDetGainUsed); | |
608 | Double_t oldmean = fCalDetGainUsed->CalcMean(kFALSE); | |
609 | //printf("oldmean %f\n",oldmean); | |
610 | if(oldmean > 0.0) { | |
611 | Double_t scalefactor = calibra->GetScaleFactorGain(); | |
612 | //printf("Correction factor %f\n",gainoverallnotnormalized*scalefactor); | |
613 | calDetGain->Multiply(gainoverallnotnormalized*scalefactor/oldmean); | |
614 | //printf("newmean %f\n",calDetGain->CalcMean(kFALSE)); | |
615 | TH1F *coefGain = calDetGain->MakeHisto1DAsFunctionOfDet(); | |
616 | fCalibObjects->AddAt(calDetGain,kGain); | |
617 | fPlots->AddAt(coefGain,kGain); | |
618 | // | |
619 | ok = kTRUE; | |
82b413fd | 620 | fStatusNeg = fStatusNeg | kGainNotEnoughStatsButFill; |
54f2ff1c | 621 | } |
622 | else { | |
82b413fd | 623 | fStatusPos = fStatusPos | kGainErrorOld; |
54f2ff1c | 624 | } |
625 | } | |
626 | else { | |
82b413fd | 627 | if(gainoverallnotnormalized <= 0.0) fStatusNeg = fStatusNeg | kGainNotEnoughStatsNotFill; |
628 | if(!fCalDetGainUsed) fStatusPos = fStatusPos | kGainErrorOld; | |
54f2ff1c | 629 | } |
630 | } | |
b3fcfd96 | 631 | |
632 | calibra->ResetVectorFit(); | |
633 | ||
634 | return ok; | |
635 | ||
636 | } | |
00d203b6 | 637 | //_____________________________________________________________________________________________________ |
b3fcfd96 | 638 | Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftT0(){ |
639 | // | |
640 | // Analyze VdriftT0 - produce the calibration objects | |
641 | // | |
642 | ||
643 | AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance(); | |
a2a4ec8e | 644 | calibra->SetMinEntries(fMinStatsVdriftT0PH); // If there is less than 1000 entries in the histo: no fit |
b3fcfd96 | 645 | calibra->AnalysePH(fPH2d); |
646 | ||
647 | Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1)) | |
648 | + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1)); | |
649 | Int_t nbfit = calibra->GetNumberFit(); | |
650 | Int_t nbfitSuccess = calibra->GetNumberFitSuccess(); | |
651 | Int_t nbE = calibra->GetNumberEnt(); | |
652 | ||
653 | //printf("nbtg %d, nbfit %d, nbE %d, nbfitSuccess %d\n",nbtg,nbfit,nbE,nbfitSuccess); | |
654 | ||
655 | Bool_t ok = kFALSE; | |
656 | if ((nbtg > 0) && | |
657 | (nbfit >= 0.5*nbE) && (nbE > 30) && (nbfitSuccess > 30)) { | |
658 | //printf("Pass the cut for VdriftT0\n"); | |
659 | // create the cal objects | |
01239968 | 660 | calibra->RemoveOutliers(1,kFALSE); |
661 | calibra->PutMeanValueOtherVectorFit(1,kFALSE); | |
662 | calibra->RemoveOutliers2(kFALSE); | |
663 | calibra->PutMeanValueOtherVectorFit2(1,kFALSE); | |
664 | // | |
b3fcfd96 | 665 | TObjArray object = calibra->GetVectorFit(); |
666 | AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object); | |
667 | TH1F *coefVdriftPH = calDetVdrift->MakeHisto1DAsFunctionOfDet(); | |
668 | AliTRDCalPad *calPadVdrift = (AliTRDCalPad *)calibra->CreatePadObjectVdrift(&object,calDetVdrift); | |
669 | TH1F *coefPadVdrift = calPadVdrift->MakeHisto1D(); | |
670 | object = calibra->GetVectorFit2(); | |
671 | AliTRDCalDet *calDetT0 = calibra->CreateDetObjectT0(&object); | |
672 | TH1F *coefT0 = calDetT0->MakeHisto1DAsFunctionOfDet(); | |
673 | AliTRDCalPad *calPadT0 = (AliTRDCalPad *)calibra->CreatePadObjectT0(&object,calDetT0); | |
674 | TH1F *coefPadT0 = calPadT0->MakeHisto1D(); | |
675 | // Put them in the array | |
676 | fCalibObjects->AddAt(calDetT0,kT0PHDet); | |
677 | fCalibObjects->AddAt(calDetVdrift,kVdriftPHDet); | |
678 | fCalibObjects->AddAt(calPadT0,kT0PHPad); | |
679 | fCalibObjects->AddAt(calPadVdrift,kVdriftPHPad); | |
680 | fPlots->AddAt(coefVdriftPH,kVdriftPHDet); | |
681 | fPlots->AddAt(coefT0,kT0PHDet); | |
682 | fPlots->AddAt(coefPadVdrift,kVdriftPHPad); | |
683 | fPlots->AddAt(coefPadT0,kT0PHPad); | |
684 | // | |
685 | ok = kTRUE; | |
686 | } | |
54f2ff1c | 687 | else { |
840ec79d | 688 | //printf("Not enough stats timeoffset\n"); |
82b413fd | 689 | fStatusNeg = fStatusNeg | kTimeOffsetNotEnoughStatsNotFill; |
54f2ff1c | 690 | } |
b3fcfd96 | 691 | calibra->ResetVectorFit(); |
692 | ||
693 | return ok; | |
694 | ||
695 | } | |
00d203b6 | 696 | //____________________________________________________________________________________________________________________ |
b3fcfd96 | 697 | Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftLinearFit(){ |
698 | // | |
699 | // Analyze vdrift linear fit - produce the calibration objects | |
700 | // | |
701 | ||
a5dcf618 | 702 | //printf("Analyse linear fit\n"); |
703 | ||
f558cb62 | 704 | |
b3fcfd96 | 705 | AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance(); |
840ec79d | 706 | calibra->SetCalDetVdriftExB(fCalDetVdriftUsed,fCalDetExBUsed); |
a2a4ec8e | 707 | calibra->SetMinEntries(fMinStatsVdriftLinear); // If there is less than 1000 entries in the histo: no fit |
a5dcf618 | 708 | //printf("Fill PE Array\n"); |
01239968 | 709 | fAliTRDCalibraVdriftLinearFit->FillPEArray(); |
a5dcf618 | 710 | //printf("AliTRDCalibraFit\n"); |
b3fcfd96 | 711 | calibra->AnalyseLinearFitters(fAliTRDCalibraVdriftLinearFit); |
a5dcf618 | 712 | //printf("After\n"); |
b3fcfd96 | 713 | |
714 | //Int_t nbtg = 540; | |
715 | Int_t nbfit = calibra->GetNumberFit(); | |
716 | Int_t nbE = calibra->GetNumberEnt(); | |
717 | ||
718 | ||
719 | Bool_t ok = kFALSE; | |
720 | // enough statistics | |
721 | if ((nbfit >= 0.5*nbE) && (nbE > 30)) { | |
722 | // create the cal objects | |
01239968 | 723 | //calibra->RemoveOutliers(1,kTRUE); |
724 | calibra->PutMeanValueOtherVectorFit(1,kTRUE); | |
725 | //calibra->RemoveOutliers2(kTRUE); | |
726 | calibra->PutMeanValueOtherVectorFit2(1,kTRUE); | |
727 | // | |
b3fcfd96 | 728 | TObjArray object = calibra->GetVectorFit(); |
729 | AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object,kTRUE); | |
730 | TH1F *coefDriftLinear = calDetVdrift->MakeHisto1DAsFunctionOfDet(); | |
731 | object = calibra->GetVectorFit2(); | |
732 | AliTRDCalDet *calDetLorentz = calibra->CreateDetObjectLorentzAngle(&object); | |
733 | TH1F *coefLorentzAngle = calDetLorentz->MakeHisto1DAsFunctionOfDet(); | |
eec29d56 | 734 | //if(!calDetLorentz) printf("No lorentz created\n"); |
b3fcfd96 | 735 | // Put them in the array |
736 | fCalibObjects->AddAt(calDetVdrift,kVdriftLinear); | |
737 | fCalibObjects->AddAt(calDetLorentz,kLorentzLinear); | |
738 | fPlots->AddAt(coefDriftLinear,kVdriftLinear); | |
739 | fPlots->AddAt(coefLorentzAngle,kLorentzLinear); | |
740 | // | |
741 | ok = kTRUE; | |
742 | } | |
54f2ff1c | 743 | else { |
744 | fNotEnoughStatisticsForTheVdriftLinear = kTRUE; | |
dee5f636 | 745 | Int_t minNumberOfEntriesForAll = fMinStatsVdriftLinear*30; |
746 | calibra->SetMinEntries(minNumberOfEntriesForAll); // Because we do it for all, we increase this | |
840ec79d | 747 | Double_t vdriftoverall = -100.0; |
748 | Double_t exboverall = 100.0; | |
749 | calibra->AnalyseLinearFittersAllTogether(fAliTRDCalibraVdriftLinearFit,vdriftoverall,exboverall); | |
750 | //printf("Found mean vdrift %f and exb %f\n",vdriftoverall,exboverall); | |
751 | if(fCalDetVdriftUsed && (vdriftoverall > 0.0) && (exboverall < 70.0)) { | |
54f2ff1c | 752 | AliTRDCalDet *calDetVdrift = new AliTRDCalDet(*fCalDetVdriftUsed); |
840ec79d | 753 | AliTRDCalDet *calDetLorentz = new AliTRDCalDet(*fCalDetExBUsed); |
754 | Double_t oldmeanvdrift = fCalDetVdriftUsed->CalcMean(kFALSE); | |
755 | Double_t oldmeanexb = fCalDetExBUsed->CalcMean(kFALSE); | |
54f2ff1c | 756 | //printf("oldmean %f\n",oldmean); |
840ec79d | 757 | if((oldmeanvdrift > 0.0) && (oldmeanexb < 70.0)) { |
54f2ff1c | 758 | //printf("Correction factor %f\n",vdriftoverall); |
840ec79d | 759 | calDetVdrift->Multiply(vdriftoverall/oldmeanvdrift); |
760 | calDetLorentz->Multiply(exboverall/oldmeanexb); | |
54f2ff1c | 761 | //printf("newmean %f\n",calDetVdrift->CalcMean(kFALSE)); |
762 | TH1F *coefDriftLinear = calDetVdrift->MakeHisto1DAsFunctionOfDet(); | |
840ec79d | 763 | TH1F *coefLorentzAngle = calDetLorentz->MakeHisto1DAsFunctionOfDet(); |
54f2ff1c | 764 | // Put them in the array |
765 | fCalibObjects->AddAt(calDetVdrift,kVdriftLinear); | |
840ec79d | 766 | fCalibObjects->AddAt(calDetLorentz,kLorentzLinear); |
54f2ff1c | 767 | fPlots->AddAt(coefDriftLinear,kVdriftLinear); |
840ec79d | 768 | fPlots->AddAt(coefLorentzAngle,kLorentzLinear); |
54f2ff1c | 769 | // |
770 | ok = kTRUE; | |
82b413fd | 771 | fStatusNeg = fStatusNeg | kVdriftNotEnoughStatsButFill; |
54f2ff1c | 772 | } |
82b413fd | 773 | else { |
840ec79d | 774 | if(oldmeanvdrift) fStatusPos = fStatusPos | kVdriftErrorOld; |
775 | if(oldmeanexb) fStatusPos = fStatusPos | kExBErrorOld; | |
82b413fd | 776 | } |
777 | } | |
778 | else { | |
840ec79d | 779 | if((vdriftoverall <= 0.0) && (exboverall > 70.0)) fStatusNeg = fStatusNeg | kVdriftNotEnoughStatsNotFill; |
82b413fd | 780 | if(!fCalDetVdriftUsed) fStatusPos = fStatusPos | kVdriftErrorOld; |
840ec79d | 781 | if(!fCalDetExBUsed) fStatusPos = fStatusPos | kExBErrorOld; |
54f2ff1c | 782 | } |
54f2ff1c | 783 | } |
b3fcfd96 | 784 | |
785 | calibra->ResetVectorFit(); | |
786 | ||
787 | return ok; | |
788 | ||
789 | } | |
00d203b6 | 790 | //________________________________________________________________________________________________________________ |
b3fcfd96 | 791 | |
a0bb5615 | 792 | Bool_t AliTRDPreprocessorOffline::AnalyzeExbAltFit(){ |
793 | // | |
794 | // Analyze vdrift linear fit - produce the calibration objects | |
795 | // | |
796 | ||
797 | //printf("Analyse linear fit\n"); | |
798 | ||
799 | ||
800 | AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance(); | |
801 | calibra->SetMinEntries(fMinStatsVdriftLinear); // If there is less than 1000 entries in the histo: no fit | |
802 | //printf("Fill PE Array\n"); | |
803 | fAliTRDCalibraExbAltFit->FillPEArray(); | |
804 | //printf("AliTRDCalibraFit\n"); | |
805 | calibra->AnalyseExbAltFit(fAliTRDCalibraExbAltFit); | |
806 | //printf("After\n"); | |
807 | ||
808 | //Int_t nbtg = 540; | |
809 | Int_t nbfit = calibra->GetNumberFit(); | |
810 | Int_t nbE = calibra->GetNumberEnt(); | |
811 | ||
812 | ||
813 | Bool_t ok = kFALSE; | |
814 | // enough statistics | |
815 | if ((nbfit >= 0.5*nbE) && (nbE > 30)) { | |
816 | // create the cal objects | |
817 | //calibra->RemoveOutliers(1,kTRUE); | |
818 | calibra->PutMeanValueOtherVectorFit2(1,kTRUE); | |
819 | // | |
820 | TObjArray object = calibra->GetVectorFit2(); | |
821 | AliTRDCalDet *calDetLorentz = calibra->CreateDetObjectExbAlt(&object); | |
822 | TH1F *coefLorentzAngle = calDetLorentz->MakeHisto1DAsFunctionOfDet(); | |
823 | //if(!calDetLorentz) printf("No lorentz created\n"); | |
824 | // Put them in the array | |
825 | fCalibObjects->AddAt(calDetLorentz,kExbAlt); | |
826 | fPlots->AddAt(coefLorentzAngle,kExbAlt); | |
827 | // | |
828 | ok = kTRUE; | |
829 | } | |
830 | ||
831 | calibra->ResetVectorFit(); | |
832 | ||
833 | return ok; | |
834 | ||
835 | } | |
836 | //________________________________________________________________________________________________________________ | |
837 | ||
b3fcfd96 | 838 | Bool_t AliTRDPreprocessorOffline::AnalyzePRF(){ |
839 | // | |
840 | // Analyze PRF - produce the calibration objects | |
841 | // | |
842 | ||
843 | AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance(); | |
a2a4ec8e | 844 | calibra->SetMinEntries(fMinStatsPRF); // If there is less than 1000 entries in the histo: no fit |
b3fcfd96 | 845 | calibra->AnalysePRFMarianFit(fPRF2d); |
846 | ||
847 | Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2)) | |
848 | + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2)); | |
849 | Int_t nbfit = calibra->GetNumberFit(); | |
850 | Int_t nbE = calibra->GetNumberEnt(); | |
851 | ||
852 | ||
853 | Bool_t ok = kFALSE; | |
854 | // enough statistics | |
855 | if ((nbtg > 0) && | |
856 | (nbfit >= 0.95*nbE) && (nbE > 30)) { | |
857 | // create the cal objects | |
858 | TObjArray object = calibra->GetVectorFit(); | |
859 | AliTRDCalPad *calPadPRF = (AliTRDCalPad*) calibra->CreatePadObjectPRF(&object); | |
860 | TH1F *coefPRF = calPadPRF->MakeHisto1D(); | |
861 | // Put them in the array | |
862 | fCalibObjects->AddAt(calPadPRF,kPRF); | |
863 | fPlots->AddAt(coefPRF,kPRF); | |
864 | // | |
865 | ok = kTRUE; | |
866 | } | |
867 | ||
868 | calibra->ResetVectorFit(); | |
869 | ||
870 | return ok; | |
871 | ||
872 | } | |
a5dcf618 | 873 | |
874 | //_____________________________________________________________________________ | |
875 | Bool_t AliTRDPreprocessorOffline::AnalyzeChamberStatus() | |
876 | { | |
4c865c34 | 877 | // |
a5dcf618 | 878 | // Produce AliTRDCalChamberStatus out of calibration results |
4c865c34 | 879 | // |
b61b92a0 | 880 | |
83d0cc79 | 881 | // set up AliTRDCalibChamberStatus |
882 | AliTRDCalibChamberStatus *ChamberStatus = new AliTRDCalibChamberStatus(); | |
883 | ChamberStatus->SetSparseI(fSparse); | |
884 | ChamberStatus->AnalyseHisto(fMinStatsChamberStatus); | |
885 | // get AliTRDCalChamberStatus | |
886 | AliTRDCalChamberStatus *CalChamberStatus = ChamberStatus->GetCalChamberStatus(); | |
01239968 | 887 | |
a5dcf618 | 888 | // get calibration objects |
889 | AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain); | |
890 | AliTRDCalDet *calDetVDrift = (AliTRDCalDet *) fCalibObjects->At(kVdriftLinear); | |
840ec79d | 891 | AliTRDCalDet *calDetExB = (AliTRDCalDet *) fCalibObjects->At(kLorentzLinear); |
a5dcf618 | 892 | |
893 | // Check | |
83d0cc79 | 894 | if((!calDetGain) || (!calDetVDrift) || (!fCH2d) || (!calDetExB) || (!CalChamberStatus)) return kFALSE; |
a5dcf618 | 895 | |
896 | // Gain | |
83d0cc79 | 897 | Double_t gainmean = calDetGain->GetMean(); |
a5dcf618 | 898 | Double_t vdriftmean = calDetVDrift->GetMean(); |
83d0cc79 | 899 | Double_t exbmean = calDetExB->GetMean(); |
a5dcf618 | 900 | |
83d0cc79 | 901 | Double_t gainrms = calDetGain->GetRMSRobust(); |
902 | Double_t vdriftrms = calDetVDrift->GetRMSRobust(); | |
903 | Double_t exbrms = calDetExB->GetRMSRobust(); | |
a5dcf618 | 904 | |
905 | //printf("Gain mean: %f, rms: %f\n",gainmean,gainrms); | |
906 | //printf("Vdrift mean: %f, rms: %f\n",vdriftmean,vdriftrms); | |
840ec79d | 907 | //printf("ExB mean: %f, rms: %f\n",exbmean,exbrms); |
908 | ||
a5dcf618 | 909 | // Check |
840ec79d | 910 | if((TMath::Abs(gainrms) < 0.001) || (TMath::Abs(vdriftrms) < 0.001) || (TMath::Abs(exbrms) < 0.0000001)) return kFALSE; |
01239968 | 911 | |
a5dcf618 | 912 | // mask chambers with empty gain entries |
913 | //Int_t counter = 0; | |
914 | for (Int_t idet = 0; idet < 540; idet++) { | |
01239968 | 915 | |
a5dcf618 | 916 | // ch2d |
917 | TH1I *projch = (TH1I *) fCH2d->ProjectionX("projch",idet+1,idet+1,(Option_t *)"e"); | |
918 | Double_t entries = projch->GetEntries(); | |
83d0cc79 | 919 | //printf("Number of entries %f for det %d\n",entries,idet); |
4c865c34 | 920 | |
a5dcf618 | 921 | // gain |
922 | Double_t gain = calDetGain->GetValue(idet); | |
4c865c34 | 923 | |
a5dcf618 | 924 | // vdrift |
925 | Double_t vdrift = calDetVDrift->GetValue(idet); | |
4c865c34 | 926 | |
840ec79d | 927 | // exb |
928 | Double_t exb = calDetExB->GetValue(idet); | |
929 | ||
4c865c34 | 930 | |
83d0cc79 | 931 | if( (entries<50 && !CalChamberStatus->IsNoData(idet)) || |
932 | TMath::Abs(gainmean-gain) > (fRMSBadCalibratedGain*gainrms) || | |
933 | TMath::Abs(vdriftmean-vdrift) > (fRMSBadCalibratedVdrift*vdriftrms) || | |
934 | TMath::Abs(exbmean-exb) > (fRMSBadCalibratedExB*exbrms) ) { | |
a5dcf618 | 935 | |
936 | //printf(" chamber det %03d masked \n",idet); | |
937 | //printf(" gainmean %f and gain %f, gainrms %f \n",gainmean,gain,gainrms); | |
938 | //printf(" vdriftmean %f and vdrift %f, vdriftrms %f \n",vdriftmean,vdrift,vdriftrms); | |
840ec79d | 939 | //printf(" exbmean %f and exb %f, exbrms %f \n",exbmean,exb,exbrms); |
940 | ||
83d0cc79 | 941 | CalChamberStatus->SetStatus(idet,AliTRDCalChamberStatus::kBadCalibrated); |
a5dcf618 | 942 | //counter++; |
943 | } | |
4c865c34 | 944 | |
b61b92a0 | 945 | delete projch; |
946 | ||
a5dcf618 | 947 | } |
b3fcfd96 | 948 | |
83d0cc79 | 949 | // Security |
950 | for(Int_t sm=0; sm < 18; sm++) { | |
951 | Int_t smnodata = 0; | |
952 | Int_t smbadcalib = 0; | |
953 | for(Int_t det = 0; det < 30; det++){ | |
954 | Int_t detector = sm*30+det; | |
955 | if(CalChamberStatus->IsNoData(detector)) smnodata++; | |
956 | else { | |
957 | if(CalChamberStatus->IsBadCalibrated(detector)) smbadcalib++; | |
958 | } | |
959 | } | |
960 | fNoData[sm] = smnodata; | |
961 | fBadCalib[sm]= smbadcalib; | |
962 | //printf("No Data %d, bad calibrated %d for %d\n",fNoData[sm],fBadCalib[sm],sm); | |
963 | } | |
964 | ||
965 | // Security | |
966 | // for(Int_t sm=0; sm < 18; sm++) { | |
967 | // Int_t counter = 0; | |
968 | // for(Int_t det = 0; det < 30; det++){ | |
969 | // Int_t detector = sm*30+det; | |
970 | // if(CalChamberStatus->IsBadCalibrated(detector)) counter++; | |
971 | // } | |
972 | // if(counter >= 20) { | |
973 | // for(Int_t det = 0; det < 30; det++){ | |
974 | // Int_t detector = sm*30+det; | |
975 | // CalChamberStatus->SetStatus(detector,AliTRDCalChamberStatus::kGood); | |
976 | // } | |
977 | // } | |
978 | // } | |
b3fcfd96 | 979 | |
a5dcf618 | 980 | fCalibObjects->AddAt(CalChamberStatus,kChamberStatus); |
981 | return kTRUE; | |
b3fcfd96 | 982 | |
a5dcf618 | 983 | } |
b3fcfd96 | 984 | |
b3fcfd96 | 985 | |
a5dcf618 | 986 | //________________________________________________________________________________________________ |
987 | void AliTRDPreprocessorOffline::CorrectFromDetGainUsed() { | |
988 | // | |
989 | // Correct from the gas gain used afterwards | |
990 | // | |
991 | AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain); | |
992 | if(!calDetGain) return; | |
b3fcfd96 | 993 | |
a5dcf618 | 994 | // Calculate mean |
995 | Double_t mean = 0.0; | |
996 | Int_t nbdet = 0; | |
b3fcfd96 | 997 | |
a5dcf618 | 998 | for(Int_t det = 0; det < 540; det++) { |
00d203b6 | 999 | |
a5dcf618 | 1000 | Float_t gaininit = fCalDetGainUsed->GetValue(det); |
1001 | Float_t gainout = calDetGain->GetValue(det); | |
00d203b6 | 1002 | |
00d203b6 | 1003 | |
a5dcf618 | 1004 | if(TMath::Abs(gainout-1.0) > 0.000001) { |
1005 | mean += (gaininit*gainout); | |
1006 | nbdet++; | |
1007 | } | |
1008 | } | |
1009 | if(nbdet > 0) mean = mean/nbdet; | |
00d203b6 | 1010 | |
a5dcf618 | 1011 | for(Int_t det = 0; det < 540; det++) { |
00d203b6 | 1012 | |
a5dcf618 | 1013 | Float_t gaininit = fCalDetGainUsed->GetValue(det); |
1014 | Float_t gainout = calDetGain->GetValue(det); | |
1015 | ||
ba1aa7a7 | 1016 | if(TMath::Abs(gainout-1.0) > 0.000001) { |
1017 | Double_t newgain = gaininit*gainout; | |
1018 | if(newgain < 0.1) newgain = 0.1; | |
58e2b572 | 1019 | if(newgain > 1.9) newgain = 1.9; |
ba1aa7a7 | 1020 | calDetGain->SetValue(det,newgain); |
1021 | } | |
1022 | else { | |
1023 | Double_t newgain = mean; | |
1024 | if(newgain < 0.1) newgain = 0.1; | |
58e2b572 | 1025 | if(newgain > 1.9) newgain = 1.9; |
ba1aa7a7 | 1026 | calDetGain->SetValue(det,newgain); |
1027 | } | |
a5dcf618 | 1028 | } |
1029 | ||
1030 | ||
1031 | } | |
1032 | //________________________________________________________________________________________________ | |
1033 | void AliTRDPreprocessorOffline::CorrectFromDetVdriftUsed() { | |
1034 | // | |
1035 | // Correct from the drift velocity | |
1036 | // | |
1037 | ||
1038 | //printf("Correct for vdrift\n"); | |
1039 | ||
1040 | AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain); | |
1041 | if(!calDetGain) return; | |
1042 | ||
1043 | Int_t detVdrift = kVdriftPHDet; | |
1044 | if(fMethodSecond) detVdrift = kVdriftLinear; | |
1045 | AliTRDCalDet *calDetVdrift = (AliTRDCalDet *) fCalibObjects->At(detVdrift); | |
1046 | if(!calDetVdrift) return; | |
1047 | ||
1048 | // Calculate mean | |
8dac2af1 | 1049 | if(!fNotEnoughStatisticsForTheVdriftLinear) { |
54f2ff1c | 1050 | for(Int_t det = 0; det < 540; det++) { |
1051 | ||
1052 | Float_t vdriftinit = fCalDetVdriftUsed->GetValue(det); | |
1053 | Float_t vdriftout = calDetVdrift->GetValue(det); | |
1054 | ||
1055 | Float_t gain = calDetGain->GetValue(det); | |
1056 | if(vdriftout > 0.0) gain = gain*vdriftinit/vdriftout; | |
1057 | if(gain < 0.1) gain = 0.1; | |
1058 | if(gain > 1.9) gain = 1.9; | |
1059 | calDetGain->SetValue(det,gain); | |
1060 | } | |
a5dcf618 | 1061 | } |
54f2ff1c | 1062 | else { |
1063 | ||
1064 | Float_t vdriftinit = fCalDetVdriftUsed->CalcMean(kFALSE); | |
1065 | Float_t vdriftout = calDetVdrift->CalcMean(kFALSE); | |
1066 | Float_t factorcorrectif = 1.0; | |
1067 | if(vdriftout > 0.0) factorcorrectif = vdriftinit/vdriftout; | |
1068 | for(Int_t det = 0; det < 540; det++) { | |
1069 | Float_t gain = calDetGain->GetValue(det); | |
1070 | gain = gain*factorcorrectif; | |
1071 | if(gain < 0.1) gain = 0.1; | |
1072 | if(gain > 1.9) gain = 1.9; | |
1073 | calDetGain->SetValue(det,gain); | |
1074 | } | |
1075 | ||
1076 | } | |
1077 | ||
a5dcf618 | 1078 | } |
54f2ff1c | 1079 | //_________________________________________________________________________________________________________________ |
a5dcf618 | 1080 | void AliTRDPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){ |
1081 | // | |
1082 | // Update OCDB entry | |
1083 | // | |
1084 | ||
1085 | AliCDBMetaData *metaData= new AliCDBMetaData(); | |
1086 | metaData->SetObjectClassName("AliTRDCalDet"); | |
1087 | metaData->SetResponsible("Raphaelle Bailhache"); | |
1088 | metaData->SetBeamPeriod(1); | |
1089 | ||
1090 | AliCDBId id1("TRD/Calib/ChamberGainFactor", startRunNumber, endRunNumber); | |
1091 | AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath); | |
1092 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain); | |
1093 | if(calDet) gStorage->Put(calDet, id1, metaData); | |
1094 | ||
1095 | ||
eec29d56 | 1096 | } |
1097 | //___________________________________________________________________________________________________________________ | |
1098 | void AliTRDPreprocessorOffline::UpdateOCDBExB(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){ | |
1099 | // | |
1100 | // Update OCDB entry | |
1101 | // | |
1102 | ||
1103 | Int_t detExB = kLorentzLinear; | |
1104 | if(!fMethodSecond) return; | |
1105 | ||
1106 | //printf("Pass\n"); | |
1107 | ||
1108 | AliCDBMetaData *metaData= new AliCDBMetaData(); | |
1109 | metaData->SetObjectClassName("AliTRDCalDet"); | |
1110 | metaData->SetResponsible("Raphaelle Bailhache"); | |
1111 | metaData->SetBeamPeriod(1); | |
1112 | ||
1113 | AliCDBId id1("TRD/Calib/ChamberExB", startRunNumber, endRunNumber); | |
1114 | AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath); | |
1115 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detExB); | |
1116 | if(calDet) gStorage->Put(calDet, id1, metaData); | |
1117 | //if(!calDet) printf("No caldet\n"); | |
1118 | ||
a5dcf618 | 1119 | } |
a0bb5615 | 1120 | //___________________________________________________________________________________________________________________ |
1121 | void AliTRDPreprocessorOffline::UpdateOCDBExBAlt(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){ | |
1122 | // | |
1123 | // Update OCDB entry | |
1124 | // | |
1125 | ||
1126 | Int_t detExB = kExbAlt; | |
1127 | if(!fMethodSecond) return; | |
1128 | ||
1129 | //printf("Pass\n"); | |
1130 | ||
1131 | AliCDBMetaData *metaData= new AliCDBMetaData(); | |
1132 | metaData->SetObjectClassName("AliTRDCalDet"); | |
1133 | metaData->SetResponsible("Raphaelle Bailhache"); | |
1134 | metaData->SetBeamPeriod(1); | |
1135 | ||
1136 | AliCDBId id1("TRD/Calib/ChamberExBAlt", startRunNumber, endRunNumber); | |
1137 | AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath); | |
1138 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detExB); | |
1139 | if(calDet) gStorage->Put(calDet, id1, metaData); | |
1140 | //if(!calDet) printf("No caldet\n"); | |
1141 | ||
1142 | } | |
a5dcf618 | 1143 | //___________________________________________________________________________________________________________________ |
1144 | void AliTRDPreprocessorOffline::UpdateOCDBVdrift(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){ | |
1145 | // | |
1146 | // Update OCDB entry | |
1147 | // | |
1148 | ||
1149 | Int_t detVdrift = kVdriftPHDet; | |
1150 | ||
1151 | if(fMethodSecond) detVdrift = kVdriftLinear; | |
1152 | ||
1153 | AliCDBMetaData *metaData= new AliCDBMetaData(); | |
1154 | metaData->SetObjectClassName("AliTRDCalDet"); | |
1155 | metaData->SetResponsible("Raphaelle Bailhache"); | |
1156 | metaData->SetBeamPeriod(1); | |
1157 | ||
1158 | AliCDBId id1("TRD/Calib/ChamberVdrift", startRunNumber, endRunNumber); | |
1159 | AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath); | |
1160 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift); | |
1161 | if(calDet) gStorage->Put(calDet, id1, metaData); | |
1162 | ||
1163 | // | |
1164 | ||
1165 | if(!fMethodSecond) { | |
1166 | ||
1167 | AliCDBMetaData *metaDataPad= new AliCDBMetaData(); | |
1168 | metaDataPad->SetObjectClassName("AliTRDCalPad"); | |
1169 | metaDataPad->SetResponsible("Raphaelle Bailhache"); | |
1170 | metaDataPad->SetBeamPeriod(1); | |
1171 | ||
1172 | AliCDBId id1Pad("TRD/Calib/LocalVdrift", startRunNumber, endRunNumber); | |
1173 | AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad); | |
1174 | if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad); | |
1175 | ||
1176 | } | |
1177 | ||
1178 | } | |
1179 | //________________________________________________________________________________________________________________________ | |
1180 | void AliTRDPreprocessorOffline::UpdateOCDBT0(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){ | |
1181 | // | |
1182 | // Update OCDB entry | |
1183 | // | |
1184 | ||
1185 | AliCDBMetaData *metaData= new AliCDBMetaData(); | |
1186 | metaData->SetObjectClassName("AliTRDCalDet"); | |
1187 | metaData->SetResponsible("Raphaelle Bailhache"); | |
1188 | metaData->SetBeamPeriod(1); | |
1189 | ||
1190 | AliCDBId id1("TRD/Calib/ChamberT0", startRunNumber, endRunNumber); | |
1191 | AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath); | |
1192 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet); | |
1193 | if(calDet) gStorage->Put(calDet, id1, metaData); | |
1194 | ||
1195 | // | |
1196 | ||
1197 | AliCDBMetaData *metaDataPad= new AliCDBMetaData(); | |
1198 | metaDataPad->SetObjectClassName("AliTRDCalPad"); | |
1199 | metaDataPad->SetResponsible("Raphaelle Bailhache"); | |
1200 | metaDataPad->SetBeamPeriod(1); | |
1201 | ||
1202 | AliCDBId id1Pad("TRD/Calib/LocalT0", startRunNumber, endRunNumber); | |
1203 | AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad); | |
1204 | if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad); | |
1205 | ||
1206 | ||
1207 | ||
1208 | } | |
1209 | //_________________________________________________________________________________________________________________ | |
1210 | void AliTRDPreprocessorOffline::UpdateOCDBPRF(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){ | |
1211 | // | |
1212 | // Update OCDB entry | |
1213 | // | |
1214 | ||
1215 | AliCDBMetaData *metaData= new AliCDBMetaData(); | |
1216 | metaData->SetObjectClassName("AliTRDCalPad"); | |
1217 | metaData->SetResponsible("Raphaelle Bailhache"); | |
1218 | metaData->SetBeamPeriod(1); | |
1219 | ||
1220 | AliCDBId id1("TRD/Calib/PRFWidth", startRunNumber, endRunNumber); | |
1221 | AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath); | |
1222 | AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF); | |
1223 | if(calPad) gStorage->Put(calPad, id1, metaData); | |
1224 | ||
1225 | ||
1226 | } | |
1227 | //_________________________________________________________________________________________________________________ | |
1228 | void AliTRDPreprocessorOffline::UpdateOCDBChamberStatus(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){ | |
1229 | // | |
1230 | // Update OCDB entry | |
1231 | // | |
1232 | ||
1233 | AliCDBMetaData *metaData= new AliCDBMetaData(); | |
1234 | metaData->SetObjectClassName("AliTRDCalChamberStatus"); | |
1235 | metaData->SetResponsible("Raphaelle Bailhache"); | |
1236 | metaData->SetBeamPeriod(1); | |
1237 | ||
1238 | AliCDBId id1("TRD/Calib/ChamberStatus", startRunNumber, endRunNumber); | |
1239 | AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath); | |
1240 | AliTRDCalChamberStatus *calChamberStatus = (AliTRDCalChamberStatus *) fCalibObjects->At(kChamberStatus); | |
1241 | if(calChamberStatus) gStorage->Put(calChamberStatus, id1, metaData); | |
1242 | ||
1243 | ||
1244 | } | |
1245 | //__________________________________________________________________________________________________________________________ | |
54f2ff1c | 1246 | Bool_t AliTRDPreprocessorOffline::ValidateGain() { |
a5dcf618 | 1247 | // |
1248 | // Validate OCDB entry | |
1249 | // | |
1250 | ||
1251 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain); | |
1252 | if(calDet) { | |
1253 | Double_t mean = calDet->GetMean(); | |
1254 | Double_t rms = calDet->GetRMSRobust(); | |
1255 | if((mean > 0.2) && (mean < 1.4) && (rms < 0.5)) return kTRUE; | |
1256 | //if((mean > 0.2) && (mean < 1.4)) return kTRUE; | |
54f2ff1c | 1257 | else { |
82b413fd | 1258 | fStatusPos = fStatusPos | kGainErrorRange; |
54f2ff1c | 1259 | return kFALSE; |
1260 | } | |
a5dcf618 | 1261 | } |
1262 | else return kFALSE; | |
54f2ff1c | 1263 | |
a5dcf618 | 1264 | |
1265 | ||
1266 | } | |
1267 | //__________________________________________________________________________________________________________________________ | |
1268 | Bool_t AliTRDPreprocessorOffline::ValidateVdrift(){ | |
1269 | // | |
1270 | // Update OCDB entry | |
1271 | // | |
1272 | ||
1273 | Int_t detVdrift = kVdriftPHDet; | |
1274 | Bool_t ok = kTRUE; | |
1275 | ||
1276 | if(fMethodSecond) detVdrift = kVdriftLinear; | |
1277 | ||
1278 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift); | |
1279 | if(calDet) { | |
1280 | Double_t mean = calDet->GetMean(); | |
1281 | Double_t rms = calDet->GetRMSRobust(); | |
1282 | //printf("Vdrift::mean %f, rms %f\n",mean,rms); | |
54f2ff1c | 1283 | if(!((mean > 1.0) && (mean < 2.0) && (rms < 0.5))) { |
82b413fd | 1284 | fStatusPos = fStatusPos | kVdriftErrorRange; |
54f2ff1c | 1285 | ok = kFALSE; |
1286 | } | |
a5dcf618 | 1287 | } |
1288 | else return kFALSE; | |
1289 | ||
1290 | if(!fMethodSecond) { | |
1291 | AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad); | |
1292 | if(calPad) { | |
1293 | Double_t mean = calPad->GetMean(); | |
1294 | Double_t rms = calPad->GetRMS(); | |
1295 | //printf("Vdrift::meanpad %f, rmspad %f\n",mean,rms); | |
54f2ff1c | 1296 | if(!((mean > 0.9) && (mean < 1.1) && (rms < 0.6))) { |
82b413fd | 1297 | fStatusPos = fStatusPos | kVdriftErrorRange; |
54f2ff1c | 1298 | ok = kFALSE; |
1299 | } | |
a5dcf618 | 1300 | } |
1301 | else return kFALSE; | |
1302 | } | |
1303 | ||
1304 | return ok; | |
1305 | ||
840ec79d | 1306 | } |
1307 | //__________________________________________________________________________________________________________________________ | |
1308 | Bool_t AliTRDPreprocessorOffline::ValidateExB(){ | |
1309 | // | |
1310 | // Update OCDB entry | |
1311 | // | |
1312 | ||
1313 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kLorentzLinear); | |
1314 | if(calDet) { | |
1315 | Double_t mean = calDet->GetMean(); | |
1316 | Double_t rms = calDet->GetRMSRobust(); | |
1317 | //printf("Vdrift::mean %f, rms %f\n",mean,rms); | |
1318 | if(!((mean > -1.0) && (mean < 1.0) && (rms < 0.5))) { | |
1319 | fStatusNeg = fStatusNeg | kExBErrorRange; | |
1320 | return kFALSE; | |
1321 | } | |
1322 | else return kTRUE; | |
1323 | } | |
1324 | else return kFALSE; | |
1325 | ||
a5dcf618 | 1326 | } |
1327 | //__________________________________________________________________________________________________________________________ | |
1328 | Bool_t AliTRDPreprocessorOffline::ValidateT0(){ | |
1329 | // | |
1330 | // Update OCDB entry | |
1331 | // | |
1332 | ||
1333 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet); | |
1334 | AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad); | |
1335 | if(calDet && calPad) { | |
1336 | Double_t meandet = calDet->GetMean(); | |
1337 | Double_t rmsdet = calDet->GetRMSRobust(); | |
1338 | Double_t meanpad = calPad->GetMean(); | |
1339 | //Double_t rmspad = calPad->GetRMS(); | |
1340 | //printf("T0::minimum %f, rmsdet %f,meanpad %f, rmspad %f\n",meandet,rmsdet,meanpad,rmspad); | |
1341 | if((meandet > -1.5) && (meandet < 5.0) && (rmsdet < 4.0) && (meanpad < 5.0) && (meanpad > -0.5)) return kTRUE; | |
54f2ff1c | 1342 | else { |
82b413fd | 1343 | fStatusPos = fStatusPos | kTimeOffsetErrorRange; |
54f2ff1c | 1344 | return kFALSE; |
1345 | } | |
a5dcf618 | 1346 | } |
1347 | else return kFALSE; | |
1348 | ||
1349 | } | |
1350 | //__________________________________________________________________________________________________________________________ | |
1351 | Bool_t AliTRDPreprocessorOffline::ValidatePRF() const{ | |
1352 | // | |
1353 | // Update OCDB entry | |
1354 | // | |
1355 | ||
1356 | AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF); | |
1357 | if(calPad) { | |
1358 | Double_t meanpad = calPad->GetMean(); | |
1359 | Double_t rmspad = calPad->GetRMS(); | |
1360 | //printf("PRF::meanpad %f, rmspad %f\n",meanpad,rmspad); | |
1361 | if((meanpad < 1.0) && (rmspad < 0.8)) return kTRUE; | |
1362 | else return kFALSE; | |
1363 | } | |
1364 | else return kFALSE; | |
1365 | ||
1366 | ||
1367 | } | |
1368 | //__________________________________________________________________________________________________________________________ | |
82b413fd | 1369 | Bool_t AliTRDPreprocessorOffline::ValidateChamberStatus(){ |
00d203b6 | 1370 | // |
1371 | // Update OCDB entry | |
1372 | // | |
1373 | ||
a5dcf618 | 1374 | AliTRDCalChamberStatus *calChamberStatus = (AliTRDCalChamberStatus *) fCalibObjects->At(kChamberStatus); |
1375 | if(calChamberStatus) { | |
83d0cc79 | 1376 | |
1377 | Int_t detectornodata = 0; | |
1378 | Int_t detectorbadcalib = 0; | |
1379 | ||
1380 | for(Int_t sm=0; sm < 18; sm++) { | |
1381 | //printf("%d chambers w/o data in sm %d\n",fNoData[sm],sm); | |
1382 | //printf("%d bad calibrated chambers in sm %d\n",fBadCalib[sm],sm); | |
1383 | if(fNoData[sm] != 30) detectornodata += fNoData[sm]; | |
1384 | detectorbadcalib+=fBadCalib[sm]; | |
a5dcf618 | 1385 | } |
83d0cc79 | 1386 | //printf("Number of chambers w/o data %d\n",detectornodata); |
1387 | //printf("Number of chambers bad calibrated %d\n",detectorbadcalib); | |
1388 | ||
1389 | if((detectornodata > fNoDataValidate) || | |
1390 | (detectorbadcalib > fBadCalibValidate)){ | |
82b413fd | 1391 | fStatusPos = fStatusPos | kChamberStatusErrorRange; |
1392 | return kFALSE; | |
1393 | } | |
83d0cc79 | 1394 | return kTRUE; |
00d203b6 | 1395 | } |
1396 | else return kFALSE; | |
83d0cc79 | 1397 | |
b3fcfd96 | 1398 | } |
4c865c34 | 1399 | //_____________________________________________________________________________ |
1400 | Int_t AliTRDPreprocessorOffline::GetVersion(TString name) const | |
1401 | { | |
1402 | // | |
1403 | // Get version from the title | |
1404 | // | |
1405 | ||
1406 | // Some patterns | |
1407 | const Char_t *version = "Ver"; | |
1408 | if(!strstr(name.Data(),version)) return -1; | |
f558cb62 | 1409 | const Char_t *after = "Subver"; |
1410 | if(!strstr(name.Data(),after)) return -1; | |
1411 | ||
4c865c34 | 1412 | for(Int_t ver = 0; ver < 999999999; ver++) { |
1413 | ||
1414 | TString vertry(version); | |
1415 | vertry += ver; | |
f558cb62 | 1416 | vertry += after; |
4c865c34 | 1417 | |
1418 | //printf("vertry %s and name %s\n",vertry.Data(),name.Data()); | |
1419 | ||
1420 | if(strstr(name.Data(),vertry.Data())) return ver; | |
1421 | ||
1422 | } | |
1423 | ||
1424 | return -1; | |
1425 | ||
1426 | } | |
1427 | ||
1428 | //_____________________________________________________________________________ | |
1429 | Int_t AliTRDPreprocessorOffline::GetSubVersion(TString name) const | |
1430 | { | |
1431 | // | |
1432 | // Get subversion from the title | |
1433 | // | |
1434 | ||
1435 | // Some patterns | |
1436 | const Char_t *subversion = "Subver"; | |
1437 | if(!strstr(name.Data(),subversion)) return -1; | |
f558cb62 | 1438 | const Char_t *after = "FirstRun"; |
1439 | if(!strstr(name.Data(),after)) { | |
1440 | after = "Nz"; | |
1441 | } | |
1442 | if(!strstr(name.Data(),after)) return -1; | |
1443 | ||
4c865c34 | 1444 | |
1445 | for(Int_t ver = 0; ver < 999999999; ver++) { | |
1446 | ||
1447 | TString vertry(subversion); | |
1448 | vertry += ver; | |
f558cb62 | 1449 | vertry += after; |
ca7e6e64 | 1450 | |
1451 | //printf("vertry %s and name %s\n",vertry.Data(),name.Data()); | |
1452 | ||
1453 | if(strstr(name.Data(),vertry.Data())) return ver; | |
1454 | ||
1455 | } | |
1456 | ||
1457 | return -1; | |
1458 | ||
1459 | } | |
1460 | ||
1461 | //_____________________________________________________________________________ | |
1462 | Int_t AliTRDPreprocessorOffline::GetFirstRun(TString name) const | |
1463 | { | |
1464 | // | |
1465 | // Get first run from the title | |
1466 | // | |
1467 | ||
1468 | // Some patterns | |
1469 | const Char_t *firstrun = "FirstRun"; | |
1470 | if(!strstr(name.Data(),firstrun)) return -1; | |
f558cb62 | 1471 | const Char_t *after = "Nz"; |
1472 | if(!strstr(name.Data(),after)) return -1; | |
1473 | ||
ca7e6e64 | 1474 | |
1475 | for(Int_t ver = 0; ver < 999999999; ver++) { | |
1476 | ||
1477 | TString vertry(firstrun); | |
1478 | vertry += ver; | |
f558cb62 | 1479 | vertry += after; |
4c865c34 | 1480 | |
1481 | //printf("vertry %s and name %s\n",vertry.Data(),name.Data()); | |
1482 | ||
1483 | if(strstr(name.Data(),vertry.Data())) return ver; | |
1484 | ||
1485 | } | |
1486 | ||
1487 | return -1; | |
1488 | ||
1489 | } | |
82b413fd | 1490 | //_____________________________________________________________________________ |
1491 | Bool_t AliTRDPreprocessorOffline::CheckStatus(Int_t status, Int_t bitMask) const | |
1492 | { | |
1493 | // | |
1494 | // Checks the status | |
1495 | // | |
1496 | ||
1497 | return (status & bitMask) ? kTRUE : kFALSE; | |
1498 | ||
1499 | } | |
1500 | //_____________________________________________________________________________ | |
1501 | Int_t AliTRDPreprocessorOffline::GetStatus() const | |
1502 | { | |
1503 | // | |
1504 | // Checks the status | |
1505 | // fStatusPos: errors | |
1506 | // fStatusNeg: only info | |
1507 | // | |
1508 | ||
1509 | if(fStatusPos > 0) return fStatusPos; | |
1510 | else return (-TMath::Abs(fStatusNeg)); | |
1511 | ||
1512 | } | |
1513 | //_____________________________________________________________________________ | |
1514 | void AliTRDPreprocessorOffline::PrintStatus() const | |
1515 | { | |
1516 | // | |
1517 | // Do Summary | |
1518 | // | |
1519 | ||
1520 | AliInfo(Form("The error status is %d",fStatusPos)); | |
1521 | AliInfo(Form("IsGainErrorOld? %d",(Int_t)IsGainErrorOld())); | |
1522 | AliInfo(Form("IsVdriftErrorOld? %d",(Int_t)IsVdriftErrorOld())); | |
1523 | AliInfo(Form("IsGainErrorRange? %d",(Int_t)IsGainErrorRange())); | |
1524 | AliInfo(Form("IsVdriftErrorRange? %d",(Int_t)IsVdriftErrorRange())); | |
1525 | AliInfo(Form("IsTimeOffsetErrorRange? %d",(Int_t)IsTimeOffsetErrorRange())); | |
1526 | AliInfo(Form("IsChamberStatusErrorRange? %d",(Int_t)IsChamberStatusErrorRange())); | |
1527 | ||
1528 | ||
1529 | AliInfo(Form("The info status is %d",fStatusNeg)); | |
1530 | AliInfo(Form("IsGainNotEnoughStatsButFill? %d",(Int_t)IsGainNotEnoughStatsButFill())); | |
1531 | AliInfo(Form("IsVdriftNotEnoughStatsButFill? %d",(Int_t)IsVdriftNotEnoughStatsButFill())); | |
1532 | AliInfo(Form("IsGainNotEnoughStatsNotFill? %d",(Int_t)IsGainNotEnoughStatsNotFill())); | |
1533 | AliInfo(Form("IsVdriftNotEnoughStatsNotFill? %d",(Int_t)IsVdriftNotEnoughStatsNotFill())); | |
1534 | AliInfo(Form("IsTimeOffsetNotEnoughStatsNotFill? %d",(Int_t)IsTimeOffsetNotEnoughStatsNotFill())); | |
840ec79d | 1535 | |
1536 | AliInfo(Form("IsExBErrorRange? %d",(Int_t)IsExBErrorRange())); | |
1537 | AliInfo(Form("IsExBErrorOld? %d",(Int_t)IsExBErrorOld())); | |
82b413fd | 1538 | |
1539 | } | |
840ec79d | 1540 | //___________________________________________________________________________________ |
1541 | void AliTRDPreprocessorOffline::SetCalDetVdrift(AliTRDCalDet *calDetVdriftUsed) | |
1542 | { | |
1543 | ||
1544 | fCalDetVdriftUsed = calDetVdriftUsed; | |
1545 | ||
1546 | fCalDetExBUsed = new AliTRDCalDet("lorentz angle tan","lorentz angle tan (detector value)"); | |
1547 | for(Int_t k = 0; k < 540; k++){ | |
1548 | fCalDetExBUsed->SetValue(k,AliTRDCommonParam::Instance()->GetOmegaTau(fCalDetVdriftUsed->GetValue(k))); | |
1549 | //printf("Set the exb object for detector %d, vdrift %f and exb %f\n",k,fCalDetVdriftUsed->GetValue(k),fCalDetExBUsed->GetValue(k)); | |
1550 | } | |
1551 | ||
1552 | }; | |
83d0cc79 | 1553 | //___________________________________________________________________________________ |
1554 | Bool_t AliTRDPreprocessorOffline::SetCalDetGain(Int_t runNumber, Int_t version, Int_t subversion) | |
1555 | { | |
1556 | // | |
1557 | // Set the fCalDetGainUsed | |
1558 | // | |
1559 | ||
1560 | if((version == 0) && (subversion == 0)) return kFALSE; | |
82b413fd | 1561 | |
83d0cc79 | 1562 | AliCDBEntry *entry = AliCDBManager::Instance()->Get("TRD/Calib/ChamberGainFactor",runNumber, version, subversion); |
1563 | if(!entry) { | |
1564 | AliError("Found no entry\n"); | |
1565 | fStatusPos = fStatusPos | kGainErrorOld; | |
1566 | return kFALSE; | |
1567 | } | |
1568 | //const AliCDBId id = entry->GetId(); | |
1569 | //version = id.GetVersion(); | |
1570 | //subversion = id.GetSubVersion(); | |
1571 | //printf("Found version %d and subversion %d for vdrift\n",version,subversion); | |
1572 | AliTRDCalDet* calDet = (AliTRDCalDet *)entry->GetObject(); | |
1573 | if(calDet) fCalDetGainUsed = calDet; | |
1574 | else { | |
1575 | fStatusPos = fStatusPos | kGainErrorOld; | |
1576 | return kFALSE; | |
1577 | } | |
1578 | ||
1579 | return kTRUE; | |
1580 | ||
1581 | } | |
1582 | //___________________________________________________________________________________ | |
1583 | Bool_t AliTRDPreprocessorOffline::SetCalDetVdriftExB(Int_t runNumber, Int_t versionv, Int_t subversionv, Int_t versionexb, Int_t subversionexb) | |
1584 | { | |
1585 | // | |
1586 | // Set the fCalDetVdriftUsed and fCalDetExBUsed | |
1587 | // | |
1588 | ||
1589 | if((versionv == 0) && (subversionv == 0)) return kFALSE; | |
1590 | ||
1591 | AliCDBEntry *entry = AliCDBManager::Instance()->Get("TRD/Calib/ChamberVdrift",runNumber, versionv, subversionv); | |
1592 | if(!entry) { | |
1593 | AliError("Found no entry\n"); | |
1594 | fStatusPos = fStatusPos | kVdriftErrorOld; | |
1595 | return kFALSE; | |
1596 | } | |
1597 | AliTRDCalDet* calDet = (AliTRDCalDet *)entry->GetObject(); | |
1598 | if(calDet) fCalDetVdriftUsed = calDet; | |
1599 | else { | |
1600 | fStatusPos = fStatusPos | kVdriftErrorOld; | |
1601 | return kFALSE; | |
1602 | } | |
1603 | ||
1604 | // ExB object | |
1605 | ||
1606 | if((versionexb == 0) && (subversionexb == 0)) { | |
1607 | ||
1608 | fCalDetExBUsed = new AliTRDCalDet("lorentz angle tan","lorentz angle tan (detector value)"); | |
1609 | for(Int_t k = 0; k < 540; k++){ | |
1610 | fCalDetExBUsed->SetValue(k,AliTRDCommonParam::Instance()->GetOmegaTau(fCalDetVdriftUsed->GetValue(k))); | |
1611 | //printf("Nothing found: set the exb object for detector %d, vdrift %f and exb %f\n",k,fCalDetVdriftUsed->GetValue(k),fCalDetExBUsed->GetValue(k)); | |
1612 | } | |
1613 | } | |
1614 | else { | |
1615 | ||
1616 | entry = 0x0; | |
1617 | entry = AliCDBManager::Instance()->Get("TRD/Calib/ChamberExB",runNumber, versionexb, subversionexb); | |
1618 | if(!entry) { | |
1619 | //printf("Found no entry\n"); | |
1620 | fStatusPos = fStatusPos | kExBErrorOld; | |
1621 | return kFALSE; | |
1622 | } | |
1623 | AliTRDCalDet* calDetexb = (AliTRDCalDet *)entry->GetObject(); | |
1624 | if(!calDetexb) { | |
1625 | fStatusPos = fStatusPos | kExBErrorOld; | |
1626 | return kFALSE; | |
1627 | } | |
1628 | ||
1629 | Double_t meanexb = calDetexb->GetMean(); | |
1630 | //printf("Mean value %f\n",meanexb); | |
1631 | if((meanexb > 70) || (fNoExBUsedInReco)) { | |
1632 | fCalDetExBUsed = new AliTRDCalDet("lorentz angle tan","lorentz angle tan (detector value)"); | |
1633 | for(Int_t k = 0; k < 540; k++){ | |
1634 | fCalDetExBUsed->SetValue(k,AliTRDCommonParam::Instance()->GetOmegaTau(fCalDetVdriftUsed->GetValue(k))); | |
1635 | //printf("Found but: set the exb object for detector %d, vdrift %f and exb %f\n",k,fCalDetVdriftUsed->GetValue(k),fCalDetExBUsed->GetValue(k)); | |
1636 | } | |
1637 | } | |
1638 | else { | |
1639 | fCalDetExBUsed = calDetexb; | |
1640 | } | |
1641 | ||
1642 | } | |
1643 | ||
1644 | ||
1645 | return kTRUE; | |
1646 | ||
1647 | } | |
01239968 | 1648 |