]>
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" | |
49 | #include "TProfile2D.h" | |
50 | #include "AliTRDCalDet.h" | |
51 | #include "AliTRDCalPad.h" | |
52 | #include "AliCDBMetaData.h" | |
53 | #include "AliCDBId.h" | |
54 | #include "AliCDBManager.h" | |
55 | #include "AliCDBStorage.h" | |
56 | #include "AliTRDCalibraMode.h" | |
57 | #include "AliTRDCalibraFit.h" | |
58 | #include "AliTRDCalibraVdriftLinearFit.h" | |
59 | #include "AliTRDPreprocessorOffline.h" | |
81a5aeca | 60 | #include "AliTRDCalChamberStatus.h" |
b3fcfd96 | 61 | |
62 | ||
63 | ClassImp(AliTRDPreprocessorOffline) | |
64 | ||
65 | AliTRDPreprocessorOffline::AliTRDPreprocessorOffline(): | |
01239968 | 66 | TNamed("TPCPreprocessorOffline","TPCPreprocessorOffline"), |
b3fcfd96 | 67 | fMethodSecond(kTRUE), |
01239968 | 68 | fNameList("TRDCalib"), |
69 | fCalDetGainUsed(0x0), | |
4c865c34 | 70 | fCalDetVdriftUsed(0x0), |
b3fcfd96 | 71 | fCH2d(0x0), |
72 | fPH2d(0x0), | |
73 | fPRF2d(0x0), | |
74 | fAliTRDCalibraVdriftLinearFit(0x0), | |
75 | fNEvents(0x0), | |
76 | fAbsoluteGain(0x0), | |
77 | fPlots(new TObjArray(8)), | |
4c865c34 | 78 | fCalibObjects(new TObjArray(8)), |
79 | fVersionGainUsed(0), | |
80 | fSubVersionGainUsed(0), | |
81 | fVersionVdriftUsed(0), | |
00d203b6 | 82 | fSubVersionVdriftUsed(0), |
83 | fSwitchOnValidation(kTRUE), | |
84 | fVdriftValidated(kFALSE), | |
a2a4ec8e | 85 | fT0Validated(kFALSE), |
86 | fMinStatsVdriftT0PH(800*20), | |
87 | fMinStatsVdriftLinear(800), | |
88 | fMinStatsGain(800), | |
89 | fMinStatsPRF(600) | |
b3fcfd96 | 90 | { |
91 | // | |
92 | // default constructor | |
93 | // | |
94 | } | |
00d203b6 | 95 | //_________________________________________________________________________________________________________________ |
b3fcfd96 | 96 | AliTRDPreprocessorOffline::~AliTRDPreprocessorOffline() { |
97 | // | |
98 | // Destructor | |
99 | // | |
100 | ||
01239968 | 101 | if(fCalDetGainUsed) delete fCalDetGainUsed; |
4c865c34 | 102 | if(fCalDetVdriftUsed) delete fCalDetVdriftUsed; |
b3fcfd96 | 103 | if(fCH2d) delete fCH2d; |
104 | if(fPH2d) delete fPH2d; | |
105 | if(fPRF2d) delete fPRF2d; | |
106 | if(fAliTRDCalibraVdriftLinearFit) delete fAliTRDCalibraVdriftLinearFit; | |
107 | if(fNEvents) delete fNEvents; | |
108 | if(fAbsoluteGain) delete fAbsoluteGain; | |
109 | if(fPlots) delete fPlots; | |
110 | if(fCalibObjects) delete fCalibObjects; | |
111 | ||
112 | } | |
00d203b6 | 113 | //___________________________________________________________________________________________________________________ |
b3fcfd96 | 114 | |
115 | void AliTRDPreprocessorOffline::CalibVdriftT0(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){ | |
116 | // | |
117 | // make calibration of the drift velocity | |
118 | // Input parameters: | |
119 | // file - the location of input file | |
120 | // startRunNumber, endRunNumber - run validity period | |
121 | // ocdbStorage - path to the OCDB storage | |
122 | // - if empty - local storage 'pwd' uesed | |
123 | if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB"; | |
124 | // | |
125 | // 1. Initialization | |
126 | // | |
008817a3 | 127 | fVdriftValidated = kTRUE; |
128 | fT0Validated = kTRUE; | |
b3fcfd96 | 129 | // |
130 | // 2. extraction of the information | |
131 | // | |
132 | if(ReadVdriftT0Global(file)) AnalyzeVdriftT0(); | |
a5dcf618 | 133 | //printf("Finish PH\n"); |
b3fcfd96 | 134 | if(ReadVdriftLinearFitGlobal(file)) AnalyzeVdriftLinearFit(); |
135 | // | |
136 | // 3. Append QA plots | |
137 | // | |
138 | //MakeDefaultPlots(fVdriftArray,fVdriftArray); | |
139 | // | |
140 | // | |
00d203b6 | 141 | // 4. validate OCDB entries |
b3fcfd96 | 142 | // |
00d203b6 | 143 | if(fSwitchOnValidation==kTRUE && ValidateVdrift()==kFALSE) { |
144 | AliError("TRD vdrift OCDB parameters out of range!"); | |
145 | fVdriftValidated = kFALSE; | |
146 | } | |
147 | if(fSwitchOnValidation==kTRUE && ValidateT0()==kFALSE) { | |
148 | AliError("TRD t0 OCDB parameters out of range!"); | |
149 | fT0Validated = kFALSE; | |
150 | } | |
b3fcfd96 | 151 | // |
00d203b6 | 152 | // 5. update of OCDB |
153 | // | |
154 | // | |
155 | if(fVdriftValidated) UpdateOCDBVdrift(startRunNumber,endRunNumber,ocdbStorage); | |
156 | if(fT0Validated) UpdateOCDBT0(startRunNumber,endRunNumber,ocdbStorage); | |
157 | ||
b3fcfd96 | 158 | } |
00d203b6 | 159 | //_________________________________________________________________________________________________________________ |
b3fcfd96 | 160 | |
161 | void AliTRDPreprocessorOffline::CalibGain(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){ | |
162 | // | |
163 | // make calibration of the drift velocity | |
164 | // Input parameters: | |
165 | // file - the location of input file | |
166 | // startRunNumber, endRunNumber - run validity period | |
167 | // ocdbStorage - path to the OCDB storage | |
168 | // - if empty - local storage 'pwd' uesed | |
169 | if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB"; | |
170 | // | |
171 | // 1. Initialization | |
172 | if(!ReadGainGlobal(file)) return; | |
173 | // | |
174 | // | |
175 | // 2. extraction of the information | |
176 | // | |
177 | AnalyzeGain(); | |
01239968 | 178 | if(fCalDetGainUsed) CorrectFromDetGainUsed(); |
4c865c34 | 179 | if(fCalDetVdriftUsed) CorrectFromDetVdriftUsed(); |
b3fcfd96 | 180 | // |
181 | // 3. Append QA plots | |
182 | // | |
183 | //MakeDefaultPlots(fVdriftArray,fVdriftArray); | |
184 | // | |
185 | // | |
00d203b6 | 186 | // 4. validate OCDB entries |
187 | // | |
188 | if(fSwitchOnValidation==kTRUE && ValidateGain()==kFALSE) { | |
189 | AliError("TRD gain OCDB parameters out of range!"); | |
190 | return; | |
191 | } | |
b3fcfd96 | 192 | // |
00d203b6 | 193 | // 5. update of OCDB |
b3fcfd96 | 194 | // |
00d203b6 | 195 | // |
196 | if((!fCalDetVdriftUsed) || (fCalDetVdriftUsed && fVdriftValidated)) UpdateOCDBGain(startRunNumber,endRunNumber,ocdbStorage); | |
197 | ||
b3fcfd96 | 198 | |
199 | } | |
00d203b6 | 200 | //________________________________________________________________________________________________________________ |
b3fcfd96 | 201 | |
202 | void AliTRDPreprocessorOffline::CalibPRF(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){ | |
203 | // | |
204 | // make calibration of the drift velocity | |
205 | // Input parameters: | |
206 | // file - the location of input file | |
207 | // startRunNumber, endRunNumber - run validity period | |
208 | // ocdbStorage - path to the OCDB storage | |
209 | // - if empty - local storage 'pwd' uesed | |
210 | if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB"; | |
211 | // | |
212 | // 1. Initialization | |
213 | if(!ReadPRFGlobal(file)) return; | |
214 | // | |
215 | // | |
216 | // 2. extraction of the information | |
217 | // | |
218 | AnalyzePRF(); | |
219 | // | |
220 | // 3. Append QA plots | |
221 | // | |
222 | //MakeDefaultPlots(fVdriftArray,fVdriftArray); | |
223 | // | |
224 | // | |
00d203b6 | 225 | // |
226 | // 4. validate OCDB entries | |
227 | // | |
228 | if(fSwitchOnValidation==kTRUE && ValidatePRF()==kFALSE) { | |
229 | AliError("TRD prf OCDB parameters out of range!"); | |
230 | return; | |
231 | } | |
232 | // | |
233 | // 5. update of OCDB | |
b3fcfd96 | 234 | // |
235 | // | |
236 | UpdateOCDBPRF(startRunNumber,endRunNumber,ocdbStorage); | |
237 | ||
a5dcf618 | 238 | } |
239 | //________________________________________________________________________________________________________________ | |
240 | ||
241 | void AliTRDPreprocessorOffline::CalibChamberStatus(Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){ | |
242 | // | |
243 | // make calibration of the chamber status | |
244 | // Input parameters: | |
245 | // startRunNumber, endRunNumber - run validity period | |
246 | // ocdbStorage - path to the OCDB storage | |
247 | // - if empty - local storage 'pwd' uesed | |
248 | if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB"; | |
249 | // | |
250 | // | |
251 | // 2. extraction of the information | |
252 | // | |
253 | if(!AnalyzeChamberStatus()) return; | |
254 | // | |
255 | // 3. Append QA plots | |
256 | // | |
257 | //MakeDefaultPlots(fVdriftArray,fVdriftArray); | |
258 | // | |
259 | // | |
260 | // | |
261 | // 4. validate OCDB entries | |
262 | // | |
263 | if(fSwitchOnValidation==kTRUE && ValidateChamberStatus()==kFALSE) { | |
264 | AliError("TRD Chamber status OCDB parameters not ok!"); | |
265 | return; | |
266 | } | |
267 | // | |
268 | // 5. update of OCDB | |
269 | // | |
270 | // | |
271 | UpdateOCDBChamberStatus(startRunNumber,endRunNumber,ocdbStorage); | |
272 | ||
b3fcfd96 | 273 | } |
00d203b6 | 274 | //______________________________________________________________________________________________________ |
4c865c34 | 275 | Bool_t AliTRDPreprocessorOffline::Init(const Char_t* fileName){ |
276 | // | |
277 | // read the calibration used during the reconstruction | |
278 | // | |
279 | ||
280 | if(ReadVdriftT0Global(fileName)) { | |
281 | ||
282 | TString nameph = fPH2d->GetTitle(); | |
283 | fVersionVdriftUsed = GetVersion(nameph); | |
284 | fSubVersionVdriftUsed = GetSubVersion(nameph); | |
285 | ||
286 | //printf("Found Version %d, Subversion %d for vdrift\n",fVersionVdriftUsed,fSubVersionVdriftUsed); | |
287 | ||
288 | } | |
289 | ||
290 | if(ReadGainGlobal(fileName)) { | |
291 | ||
292 | TString namech = fCH2d->GetTitle(); | |
293 | fVersionGainUsed = GetVersion(namech); | |
294 | fSubVersionGainUsed = GetSubVersion(namech); | |
295 | ||
296 | //printf("Found Version %d, Subversion %d for gain\n",fVersionGainUsed,fSubVersionGainUsed); | |
297 | ||
298 | } | |
299 | ||
300 | return kTRUE; | |
301 | ||
302 | } | |
00d203b6 | 303 | //___________________________________________________________________________________________________________________ |
b3fcfd96 | 304 | |
305 | Bool_t AliTRDPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){ | |
306 | // | |
307 | // read calibration entries from file | |
308 | // | |
4c865c34 | 309 | if(fCH2d) return kTRUE; |
b3fcfd96 | 310 | TFile fcalib(fileName); |
01239968 | 311 | TObjArray * array = (TObjArray*)fcalib.Get(fNameList); |
b3fcfd96 | 312 | if (array){ |
313 | TH2I *ch2d = (TH2I *) array->FindObject("CH2d"); | |
314 | if(!ch2d) return kFALSE; | |
315 | fCH2d = (TH2I*)ch2d->Clone(); | |
316 | //fNEvents = (TH1I *) array->FindObject("NEvents"); | |
317 | //fAbsoluteGain = (TH2F *) array->FindObject("AbsoluteGain"); | |
318 | }else{ | |
319 | TH2I *ch2d = (TH2I *) fcalib.Get("CH2d"); | |
320 | if(!ch2d) return kFALSE; | |
321 | fCH2d = (TH2I*)ch2d->Clone(); | |
322 | //fNEvents = (TH1I *) fcalib.Get("NEvents"); | |
323 | //fAbsoluteGain = (TH2F *) fcalib.Get("AbsoluteGain"); | |
324 | } | |
325 | fCH2d->SetDirectory(0); | |
326 | //printf("title of CH2d %s\n",fCH2d->GetTitle()); | |
327 | ||
328 | return kTRUE; | |
329 | ||
330 | } | |
00d203b6 | 331 | //_________________________________________________________________________________________________________________ |
b3fcfd96 | 332 | |
333 | Bool_t AliTRDPreprocessorOffline::ReadVdriftT0Global(const Char_t* fileName){ | |
334 | // | |
335 | // read calibration entries from file | |
336 | // | |
4c865c34 | 337 | if(fPH2d) return kTRUE; |
b3fcfd96 | 338 | TFile fcalib(fileName); |
01239968 | 339 | TObjArray * array = (TObjArray*)fcalib.Get(fNameList); |
b3fcfd96 | 340 | if (array){ |
341 | TProfile2D *ph2d = (TProfile2D *) array->FindObject("PH2d"); | |
342 | if(!ph2d) return kFALSE; | |
343 | fPH2d = (TProfile2D*)ph2d->Clone(); | |
344 | //fNEvents = (TH1I *) array->FindObject("NEvents"); | |
345 | }else{ | |
346 | TProfile2D *ph2d = (TProfile2D *) fcalib.Get("PH2d"); | |
347 | if(!ph2d) return kFALSE; | |
348 | fPH2d = (TProfile2D*)ph2d->Clone(); | |
349 | //fNEvents = (TH1I *) fcalib.Get("NEvents"); | |
350 | } | |
351 | fPH2d->SetDirectory(0); | |
352 | //printf("title of PH2d %s\n",fPH2d->GetTitle()); | |
353 | ||
354 | return kTRUE; | |
355 | ||
356 | } | |
00d203b6 | 357 | //___________________________________________________________________________________________________________________ |
b3fcfd96 | 358 | |
359 | Bool_t AliTRDPreprocessorOffline::ReadVdriftLinearFitGlobal(const Char_t* fileName){ | |
360 | // | |
361 | // read calibration entries from file | |
362 | // | |
4c865c34 | 363 | if(fAliTRDCalibraVdriftLinearFit) return kTRUE; |
b3fcfd96 | 364 | TFile fcalib(fileName); |
01239968 | 365 | TObjArray * array = (TObjArray*)fcalib.Get(fNameList); |
b3fcfd96 | 366 | if (array){ |
367 | fAliTRDCalibraVdriftLinearFit = (AliTRDCalibraVdriftLinearFit *) array->FindObject("AliTRDCalibraVdriftLinearFit"); | |
368 | //fNEvents = (TH1I *) array->FindObject("NEvents"); | |
369 | }else{ | |
370 | fAliTRDCalibraVdriftLinearFit = (AliTRDCalibraVdriftLinearFit *) fcalib.Get("AliTRDCalibraVdriftLinearFit"); | |
371 | //fNEvents = (TH1I *) fcalib.Get("NEvents"); | |
372 | } | |
373 | if(!fAliTRDCalibraVdriftLinearFit) { | |
374 | //printf("No AliTRDCalibraVdriftLinearFit\n"); | |
375 | return kFALSE; | |
376 | } | |
377 | return kTRUE; | |
378 | ||
379 | } | |
00d203b6 | 380 | //_____________________________________________________________________________________________________________ |
b3fcfd96 | 381 | |
382 | Bool_t AliTRDPreprocessorOffline::ReadPRFGlobal(const Char_t* fileName){ | |
383 | // | |
384 | // read calibration entries from file | |
385 | // | |
4c865c34 | 386 | if(fPRF2d) return kTRUE; |
b3fcfd96 | 387 | TFile fcalib(fileName); |
01239968 | 388 | TObjArray * array = (TObjArray*)fcalib.Get(fNameList); |
b3fcfd96 | 389 | if (array){ |
390 | TProfile2D *prf2d = (TProfile2D *) array->FindObject("PRF2d"); | |
391 | if(!prf2d) return kFALSE; | |
392 | fPRF2d = (TProfile2D*)prf2d->Clone(); | |
393 | //fNEvents = (TH1I *) array->FindObject("NEvents"); | |
394 | }else{ | |
395 | TProfile2D *prf2d = (TProfile2D *) fcalib.Get("PRF2d"); | |
396 | if(!prf2d) return kFALSE; | |
397 | fPRF2d = (TProfile2D*)prf2d->Clone(); | |
398 | //fNEvents = (TH1I *) fcalib.Get("NEvents"); | |
399 | } | |
400 | fPRF2d->SetDirectory(0); | |
401 | //printf("title of PRF2d %s\n",fPRF2d->GetTitle()); | |
402 | ||
403 | return kTRUE; | |
404 | ||
405 | } | |
00d203b6 | 406 | //__________________________________________________________________________________________________________ |
b3fcfd96 | 407 | |
408 | Bool_t AliTRDPreprocessorOffline::AnalyzeGain(){ | |
409 | // | |
410 | // Analyze gain - produce the calibration objects | |
411 | // | |
412 | ||
413 | AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance(); | |
a2a4ec8e | 414 | calibra->SetMinEntries(fMinStatsGain); // If there is less than 1000 entries in the histo: no fit |
b3fcfd96 | 415 | calibra->AnalyseCH(fCH2d); |
416 | ||
417 | Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(0)) | |
418 | + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(0)); | |
419 | Int_t nbfit = calibra->GetNumberFit(); | |
420 | Int_t nbE = calibra->GetNumberEnt(); | |
421 | ||
422 | ||
423 | Bool_t ok = kFALSE; | |
4c865c34 | 424 | Bool_t meanother = kFALSE; |
b3fcfd96 | 425 | // enough statistics |
426 | if ((nbtg > 0) && | |
427 | (nbfit >= 0.5*nbE) && (nbE > 30)) { | |
428 | // create the cal objects | |
4c865c34 | 429 | if(!fCalDetGainUsed) { |
430 | calibra->PutMeanValueOtherVectorFit(1,kTRUE); | |
431 | meanother = kTRUE; | |
432 | } | |
b3fcfd96 | 433 | TObjArray object = calibra->GetVectorFit(); |
4c865c34 | 434 | AliTRDCalDet *calDetGain = calibra->CreateDetObjectGain(&object,meanother); |
b3fcfd96 | 435 | TH1F *coefGain = calDetGain->MakeHisto1DAsFunctionOfDet(); |
436 | // Put them in the array | |
437 | fCalibObjects->AddAt(calDetGain,kGain); | |
438 | fPlots->AddAt(coefGain,kGain); | |
439 | // | |
440 | ok = kTRUE; | |
441 | } | |
442 | ||
443 | calibra->ResetVectorFit(); | |
444 | ||
445 | return ok; | |
446 | ||
447 | } | |
00d203b6 | 448 | //_____________________________________________________________________________________________________ |
b3fcfd96 | 449 | Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftT0(){ |
450 | // | |
451 | // Analyze VdriftT0 - produce the calibration objects | |
452 | // | |
453 | ||
454 | AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance(); | |
a2a4ec8e | 455 | calibra->SetMinEntries(fMinStatsVdriftT0PH); // If there is less than 1000 entries in the histo: no fit |
b3fcfd96 | 456 | calibra->AnalysePH(fPH2d); |
457 | ||
458 | Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1)) | |
459 | + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1)); | |
460 | Int_t nbfit = calibra->GetNumberFit(); | |
461 | Int_t nbfitSuccess = calibra->GetNumberFitSuccess(); | |
462 | Int_t nbE = calibra->GetNumberEnt(); | |
463 | ||
464 | //printf("nbtg %d, nbfit %d, nbE %d, nbfitSuccess %d\n",nbtg,nbfit,nbE,nbfitSuccess); | |
465 | ||
466 | Bool_t ok = kFALSE; | |
467 | if ((nbtg > 0) && | |
468 | (nbfit >= 0.5*nbE) && (nbE > 30) && (nbfitSuccess > 30)) { | |
469 | //printf("Pass the cut for VdriftT0\n"); | |
470 | // create the cal objects | |
01239968 | 471 | calibra->RemoveOutliers(1,kFALSE); |
472 | calibra->PutMeanValueOtherVectorFit(1,kFALSE); | |
473 | calibra->RemoveOutliers2(kFALSE); | |
474 | calibra->PutMeanValueOtherVectorFit2(1,kFALSE); | |
475 | // | |
b3fcfd96 | 476 | TObjArray object = calibra->GetVectorFit(); |
477 | AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object); | |
478 | TH1F *coefVdriftPH = calDetVdrift->MakeHisto1DAsFunctionOfDet(); | |
479 | AliTRDCalPad *calPadVdrift = (AliTRDCalPad *)calibra->CreatePadObjectVdrift(&object,calDetVdrift); | |
480 | TH1F *coefPadVdrift = calPadVdrift->MakeHisto1D(); | |
481 | object = calibra->GetVectorFit2(); | |
482 | AliTRDCalDet *calDetT0 = calibra->CreateDetObjectT0(&object); | |
483 | TH1F *coefT0 = calDetT0->MakeHisto1DAsFunctionOfDet(); | |
484 | AliTRDCalPad *calPadT0 = (AliTRDCalPad *)calibra->CreatePadObjectT0(&object,calDetT0); | |
485 | TH1F *coefPadT0 = calPadT0->MakeHisto1D(); | |
486 | // Put them in the array | |
487 | fCalibObjects->AddAt(calDetT0,kT0PHDet); | |
488 | fCalibObjects->AddAt(calDetVdrift,kVdriftPHDet); | |
489 | fCalibObjects->AddAt(calPadT0,kT0PHPad); | |
490 | fCalibObjects->AddAt(calPadVdrift,kVdriftPHPad); | |
491 | fPlots->AddAt(coefVdriftPH,kVdriftPHDet); | |
492 | fPlots->AddAt(coefT0,kT0PHDet); | |
493 | fPlots->AddAt(coefPadVdrift,kVdriftPHPad); | |
494 | fPlots->AddAt(coefPadT0,kT0PHPad); | |
495 | // | |
496 | ok = kTRUE; | |
497 | } | |
498 | calibra->ResetVectorFit(); | |
499 | ||
500 | return ok; | |
501 | ||
502 | } | |
00d203b6 | 503 | //____________________________________________________________________________________________________________________ |
b3fcfd96 | 504 | Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftLinearFit(){ |
505 | // | |
506 | // Analyze vdrift linear fit - produce the calibration objects | |
507 | // | |
508 | ||
a5dcf618 | 509 | //printf("Analyse linear fit\n"); |
510 | ||
b3fcfd96 | 511 | AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance(); |
a2a4ec8e | 512 | calibra->SetMinEntries(fMinStatsVdriftLinear); // If there is less than 1000 entries in the histo: no fit |
a5dcf618 | 513 | //printf("Fill PE Array\n"); |
01239968 | 514 | fAliTRDCalibraVdriftLinearFit->FillPEArray(); |
a5dcf618 | 515 | //printf("AliTRDCalibraFit\n"); |
b3fcfd96 | 516 | calibra->AnalyseLinearFitters(fAliTRDCalibraVdriftLinearFit); |
a5dcf618 | 517 | //printf("After\n"); |
b3fcfd96 | 518 | |
519 | //Int_t nbtg = 540; | |
520 | Int_t nbfit = calibra->GetNumberFit(); | |
521 | Int_t nbE = calibra->GetNumberEnt(); | |
522 | ||
523 | ||
524 | Bool_t ok = kFALSE; | |
525 | // enough statistics | |
526 | if ((nbfit >= 0.5*nbE) && (nbE > 30)) { | |
527 | // create the cal objects | |
01239968 | 528 | //calibra->RemoveOutliers(1,kTRUE); |
529 | calibra->PutMeanValueOtherVectorFit(1,kTRUE); | |
530 | //calibra->RemoveOutliers2(kTRUE); | |
531 | calibra->PutMeanValueOtherVectorFit2(1,kTRUE); | |
532 | // | |
b3fcfd96 | 533 | TObjArray object = calibra->GetVectorFit(); |
534 | AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object,kTRUE); | |
535 | TH1F *coefDriftLinear = calDetVdrift->MakeHisto1DAsFunctionOfDet(); | |
536 | object = calibra->GetVectorFit2(); | |
537 | AliTRDCalDet *calDetLorentz = calibra->CreateDetObjectLorentzAngle(&object); | |
538 | TH1F *coefLorentzAngle = calDetLorentz->MakeHisto1DAsFunctionOfDet(); | |
539 | // Put them in the array | |
540 | fCalibObjects->AddAt(calDetVdrift,kVdriftLinear); | |
541 | fCalibObjects->AddAt(calDetLorentz,kLorentzLinear); | |
542 | fPlots->AddAt(coefDriftLinear,kVdriftLinear); | |
543 | fPlots->AddAt(coefLorentzAngle,kLorentzLinear); | |
544 | // | |
545 | ok = kTRUE; | |
546 | } | |
547 | ||
548 | calibra->ResetVectorFit(); | |
549 | ||
550 | return ok; | |
551 | ||
552 | } | |
00d203b6 | 553 | //________________________________________________________________________________________________________________ |
b3fcfd96 | 554 | |
555 | Bool_t AliTRDPreprocessorOffline::AnalyzePRF(){ | |
556 | // | |
557 | // Analyze PRF - produce the calibration objects | |
558 | // | |
559 | ||
560 | AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance(); | |
a2a4ec8e | 561 | calibra->SetMinEntries(fMinStatsPRF); // If there is less than 1000 entries in the histo: no fit |
b3fcfd96 | 562 | calibra->AnalysePRFMarianFit(fPRF2d); |
563 | ||
564 | Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2)) | |
565 | + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2)); | |
566 | Int_t nbfit = calibra->GetNumberFit(); | |
567 | Int_t nbE = calibra->GetNumberEnt(); | |
568 | ||
569 | ||
570 | Bool_t ok = kFALSE; | |
571 | // enough statistics | |
572 | if ((nbtg > 0) && | |
573 | (nbfit >= 0.95*nbE) && (nbE > 30)) { | |
574 | // create the cal objects | |
575 | TObjArray object = calibra->GetVectorFit(); | |
576 | AliTRDCalPad *calPadPRF = (AliTRDCalPad*) calibra->CreatePadObjectPRF(&object); | |
577 | TH1F *coefPRF = calPadPRF->MakeHisto1D(); | |
578 | // Put them in the array | |
579 | fCalibObjects->AddAt(calPadPRF,kPRF); | |
580 | fPlots->AddAt(coefPRF,kPRF); | |
581 | // | |
582 | ok = kTRUE; | |
583 | } | |
584 | ||
585 | calibra->ResetVectorFit(); | |
586 | ||
587 | return ok; | |
588 | ||
589 | } | |
a5dcf618 | 590 | |
591 | //_____________________________________________________________________________ | |
592 | Bool_t AliTRDPreprocessorOffline::AnalyzeChamberStatus() | |
593 | { | |
4c865c34 | 594 | // |
a5dcf618 | 595 | // Produce AliTRDCalChamberStatus out of calibration results |
4c865c34 | 596 | // |
b61b92a0 | 597 | |
a5dcf618 | 598 | // set up AliTRDCalChamberStatus |
599 | AliTRDCalChamberStatus *CalChamberStatus = new AliTRDCalChamberStatus(); | |
600 | for(Int_t det = 0; det < 540; det++) CalChamberStatus->SetStatus(det,1); | |
01239968 | 601 | |
a5dcf618 | 602 | // get calibration objects |
603 | AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain); | |
604 | AliTRDCalDet *calDetVDrift = (AliTRDCalDet *) fCalibObjects->At(kVdriftLinear); | |
605 | ||
606 | // Check | |
b61b92a0 | 607 | if((!calDetGain) || (!calDetVDrift) || (!fCH2d)) return kFALSE; |
a5dcf618 | 608 | |
609 | // Gain | |
610 | Double_t gainmean = calDetGain->GetMean(); | |
611 | Double_t vdriftmean = calDetVDrift->GetMean(); | |
612 | ||
613 | Double_t gainrms = calDetGain->GetRMSRobust(); | |
614 | Double_t vdriftrms = calDetVDrift->GetRMSRobust(); | |
615 | ||
616 | //printf("Gain mean: %f, rms: %f\n",gainmean,gainrms); | |
617 | //printf("Vdrift mean: %f, rms: %f\n",vdriftmean,vdriftrms); | |
4c865c34 | 618 | |
a5dcf618 | 619 | // Check |
620 | if((TMath::Abs(gainrms) < 0.001) || (TMath::Abs(vdriftrms) < 0.001)) return kFALSE; | |
01239968 | 621 | |
a5dcf618 | 622 | // mask chambers with empty gain entries |
623 | //Int_t counter = 0; | |
624 | for (Int_t idet = 0; idet < 540; idet++) { | |
01239968 | 625 | |
a5dcf618 | 626 | // ch2d |
627 | TH1I *projch = (TH1I *) fCH2d->ProjectionX("projch",idet+1,idet+1,(Option_t *)"e"); | |
628 | Double_t entries = projch->GetEntries(); | |
4c865c34 | 629 | |
a5dcf618 | 630 | // gain |
631 | Double_t gain = calDetGain->GetValue(idet); | |
4c865c34 | 632 | |
a5dcf618 | 633 | // vdrift |
634 | Double_t vdrift = calDetVDrift->GetValue(idet); | |
4c865c34 | 635 | |
4c865c34 | 636 | |
a5dcf618 | 637 | if(entries<=0.5 || |
638 | TMath::Abs(gainmean-gain) > (15.0*gainrms) || | |
639 | TMath::Abs(vdriftmean-vdrift) > (15.0*vdriftrms)) { | |
640 | ||
641 | //printf(" chamber det %03d masked \n",idet); | |
642 | //printf(" gainmean %f and gain %f, gainrms %f \n",gainmean,gain,gainrms); | |
643 | //printf(" vdriftmean %f and vdrift %f, vdriftrms %f \n",vdriftmean,vdrift,vdriftrms); | |
644 | CalChamberStatus->SetStatus(idet,AliTRDCalChamberStatus::kMasked); | |
645 | //counter++; | |
646 | } | |
4c865c34 | 647 | |
a5dcf618 | 648 | /* |
649 | // installed supermodules+1 -> abort | |
650 | if(counter > (7+1)*30) { | |
651 | printf("ERROR: more than one SM to be masked!! \n Abort...\n"); | |
652 | if(projch) delete projch; | |
653 | return 0x0; | |
654 | } | |
655 | */ | |
4c865c34 | 656 | |
b61b92a0 | 657 | delete projch; |
658 | ||
a5dcf618 | 659 | } |
b3fcfd96 | 660 | |
a5dcf618 | 661 | // Security |
662 | for(Int_t sm=0; sm < 18; sm++) { | |
663 | Int_t counter = 0; | |
664 | for(Int_t det = 0; det < 30; det++){ | |
665 | Int_t detector = sm*30+det; | |
666 | if(CalChamberStatus->IsMasked(detector)) counter++; | |
667 | } | |
668 | if(counter >= 10) { | |
669 | for(Int_t det = 0; det < 30; det++){ | |
670 | Int_t detector = sm*30+det; | |
671 | CalChamberStatus->SetStatus(detector,AliTRDCalChamberStatus::kInstalled); | |
672 | } | |
673 | } | |
674 | } | |
b3fcfd96 | 675 | |
a5dcf618 | 676 | fCalibObjects->AddAt(CalChamberStatus,kChamberStatus); |
677 | return kTRUE; | |
b3fcfd96 | 678 | |
a5dcf618 | 679 | } |
b3fcfd96 | 680 | |
b3fcfd96 | 681 | |
a5dcf618 | 682 | //________________________________________________________________________________________________ |
683 | void AliTRDPreprocessorOffline::CorrectFromDetGainUsed() { | |
684 | // | |
685 | // Correct from the gas gain used afterwards | |
686 | // | |
687 | AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain); | |
688 | if(!calDetGain) return; | |
b3fcfd96 | 689 | |
a5dcf618 | 690 | // Calculate mean |
691 | Double_t mean = 0.0; | |
692 | Int_t nbdet = 0; | |
b3fcfd96 | 693 | |
a5dcf618 | 694 | for(Int_t det = 0; det < 540; det++) { |
00d203b6 | 695 | |
a5dcf618 | 696 | Float_t gaininit = fCalDetGainUsed->GetValue(det); |
697 | Float_t gainout = calDetGain->GetValue(det); | |
00d203b6 | 698 | |
00d203b6 | 699 | |
a5dcf618 | 700 | if(TMath::Abs(gainout-1.0) > 0.000001) { |
701 | mean += (gaininit*gainout); | |
702 | nbdet++; | |
703 | } | |
704 | } | |
705 | if(nbdet > 0) mean = mean/nbdet; | |
00d203b6 | 706 | |
a5dcf618 | 707 | for(Int_t det = 0; det < 540; det++) { |
00d203b6 | 708 | |
a5dcf618 | 709 | Float_t gaininit = fCalDetGainUsed->GetValue(det); |
710 | Float_t gainout = calDetGain->GetValue(det); | |
711 | ||
712 | if(TMath::Abs(gainout-1.0) > 0.000001) calDetGain->SetValue(det,gaininit*gainout); | |
713 | else calDetGain->SetValue(det,mean); | |
714 | } | |
715 | ||
716 | ||
717 | } | |
718 | //________________________________________________________________________________________________ | |
719 | void AliTRDPreprocessorOffline::CorrectFromDetVdriftUsed() { | |
720 | // | |
721 | // Correct from the drift velocity | |
722 | // | |
723 | ||
724 | //printf("Correct for vdrift\n"); | |
725 | ||
726 | AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain); | |
727 | if(!calDetGain) return; | |
728 | ||
729 | Int_t detVdrift = kVdriftPHDet; | |
730 | if(fMethodSecond) detVdrift = kVdriftLinear; | |
731 | AliTRDCalDet *calDetVdrift = (AliTRDCalDet *) fCalibObjects->At(detVdrift); | |
732 | if(!calDetVdrift) return; | |
733 | ||
734 | // Calculate mean | |
735 | ||
736 | for(Int_t det = 0; det < 540; det++) { | |
737 | ||
738 | Float_t vdriftinit = fCalDetVdriftUsed->GetValue(det); | |
739 | Float_t vdriftout = calDetVdrift->GetValue(det); | |
740 | ||
741 | Float_t gain = calDetGain->GetValue(det); | |
742 | if(vdriftout > 0.0) gain = gain*vdriftinit/vdriftout; | |
743 | calDetGain->SetValue(det,gain); | |
744 | ||
745 | ||
746 | } | |
747 | ||
748 | } | |
749 | //_________________________________________________________________________________________________________________ | |
750 | void AliTRDPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){ | |
751 | // | |
752 | // Update OCDB entry | |
753 | // | |
754 | ||
755 | AliCDBMetaData *metaData= new AliCDBMetaData(); | |
756 | metaData->SetObjectClassName("AliTRDCalDet"); | |
757 | metaData->SetResponsible("Raphaelle Bailhache"); | |
758 | metaData->SetBeamPeriod(1); | |
759 | ||
760 | AliCDBId id1("TRD/Calib/ChamberGainFactor", startRunNumber, endRunNumber); | |
761 | AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath); | |
762 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain); | |
763 | if(calDet) gStorage->Put(calDet, id1, metaData); | |
764 | ||
765 | ||
766 | } | |
767 | //___________________________________________________________________________________________________________________ | |
768 | void AliTRDPreprocessorOffline::UpdateOCDBVdrift(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){ | |
769 | // | |
770 | // Update OCDB entry | |
771 | // | |
772 | ||
773 | Int_t detVdrift = kVdriftPHDet; | |
774 | ||
775 | if(fMethodSecond) detVdrift = kVdriftLinear; | |
776 | ||
777 | AliCDBMetaData *metaData= new AliCDBMetaData(); | |
778 | metaData->SetObjectClassName("AliTRDCalDet"); | |
779 | metaData->SetResponsible("Raphaelle Bailhache"); | |
780 | metaData->SetBeamPeriod(1); | |
781 | ||
782 | AliCDBId id1("TRD/Calib/ChamberVdrift", startRunNumber, endRunNumber); | |
783 | AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath); | |
784 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift); | |
785 | if(calDet) gStorage->Put(calDet, id1, metaData); | |
786 | ||
787 | // | |
788 | ||
789 | if(!fMethodSecond) { | |
790 | ||
791 | AliCDBMetaData *metaDataPad= new AliCDBMetaData(); | |
792 | metaDataPad->SetObjectClassName("AliTRDCalPad"); | |
793 | metaDataPad->SetResponsible("Raphaelle Bailhache"); | |
794 | metaDataPad->SetBeamPeriod(1); | |
795 | ||
796 | AliCDBId id1Pad("TRD/Calib/LocalVdrift", startRunNumber, endRunNumber); | |
797 | AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad); | |
798 | if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad); | |
799 | ||
800 | } | |
801 | ||
802 | } | |
803 | //________________________________________________________________________________________________________________________ | |
804 | void AliTRDPreprocessorOffline::UpdateOCDBT0(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){ | |
805 | // | |
806 | // Update OCDB entry | |
807 | // | |
808 | ||
809 | AliCDBMetaData *metaData= new AliCDBMetaData(); | |
810 | metaData->SetObjectClassName("AliTRDCalDet"); | |
811 | metaData->SetResponsible("Raphaelle Bailhache"); | |
812 | metaData->SetBeamPeriod(1); | |
813 | ||
814 | AliCDBId id1("TRD/Calib/ChamberT0", startRunNumber, endRunNumber); | |
815 | AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath); | |
816 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet); | |
817 | if(calDet) gStorage->Put(calDet, id1, metaData); | |
818 | ||
819 | // | |
820 | ||
821 | AliCDBMetaData *metaDataPad= new AliCDBMetaData(); | |
822 | metaDataPad->SetObjectClassName("AliTRDCalPad"); | |
823 | metaDataPad->SetResponsible("Raphaelle Bailhache"); | |
824 | metaDataPad->SetBeamPeriod(1); | |
825 | ||
826 | AliCDBId id1Pad("TRD/Calib/LocalT0", startRunNumber, endRunNumber); | |
827 | AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad); | |
828 | if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad); | |
829 | ||
830 | ||
831 | ||
832 | } | |
833 | //_________________________________________________________________________________________________________________ | |
834 | void AliTRDPreprocessorOffline::UpdateOCDBPRF(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){ | |
835 | // | |
836 | // Update OCDB entry | |
837 | // | |
838 | ||
839 | AliCDBMetaData *metaData= new AliCDBMetaData(); | |
840 | metaData->SetObjectClassName("AliTRDCalPad"); | |
841 | metaData->SetResponsible("Raphaelle Bailhache"); | |
842 | metaData->SetBeamPeriod(1); | |
843 | ||
844 | AliCDBId id1("TRD/Calib/PRFWidth", startRunNumber, endRunNumber); | |
845 | AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath); | |
846 | AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF); | |
847 | if(calPad) gStorage->Put(calPad, id1, metaData); | |
848 | ||
849 | ||
850 | } | |
851 | //_________________________________________________________________________________________________________________ | |
852 | void AliTRDPreprocessorOffline::UpdateOCDBChamberStatus(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){ | |
853 | // | |
854 | // Update OCDB entry | |
855 | // | |
856 | ||
857 | AliCDBMetaData *metaData= new AliCDBMetaData(); | |
858 | metaData->SetObjectClassName("AliTRDCalChamberStatus"); | |
859 | metaData->SetResponsible("Raphaelle Bailhache"); | |
860 | metaData->SetBeamPeriod(1); | |
861 | ||
862 | AliCDBId id1("TRD/Calib/ChamberStatus", startRunNumber, endRunNumber); | |
863 | AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath); | |
864 | AliTRDCalChamberStatus *calChamberStatus = (AliTRDCalChamberStatus *) fCalibObjects->At(kChamberStatus); | |
865 | if(calChamberStatus) gStorage->Put(calChamberStatus, id1, metaData); | |
866 | ||
867 | ||
868 | } | |
869 | //__________________________________________________________________________________________________________________________ | |
870 | Bool_t AliTRDPreprocessorOffline::ValidateGain() const { | |
871 | // | |
872 | // Validate OCDB entry | |
873 | // | |
874 | ||
875 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain); | |
876 | if(calDet) { | |
877 | Double_t mean = calDet->GetMean(); | |
878 | Double_t rms = calDet->GetRMSRobust(); | |
879 | if((mean > 0.2) && (mean < 1.4) && (rms < 0.5)) return kTRUE; | |
880 | //if((mean > 0.2) && (mean < 1.4)) return kTRUE; | |
881 | else return kFALSE; | |
882 | } | |
883 | else return kFALSE; | |
884 | ||
885 | ||
886 | } | |
887 | //__________________________________________________________________________________________________________________________ | |
888 | Bool_t AliTRDPreprocessorOffline::ValidateVdrift(){ | |
889 | // | |
890 | // Update OCDB entry | |
891 | // | |
892 | ||
893 | Int_t detVdrift = kVdriftPHDet; | |
894 | Bool_t ok = kTRUE; | |
895 | ||
896 | if(fMethodSecond) detVdrift = kVdriftLinear; | |
897 | ||
898 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift); | |
899 | if(calDet) { | |
900 | Double_t mean = calDet->GetMean(); | |
901 | Double_t rms = calDet->GetRMSRobust(); | |
902 | //printf("Vdrift::mean %f, rms %f\n",mean,rms); | |
903 | if(!((mean > 1.0) && (mean < 2.0) && (rms < 0.5))) ok = kFALSE; | |
904 | } | |
905 | else return kFALSE; | |
906 | ||
907 | if(!fMethodSecond) { | |
908 | AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad); | |
909 | if(calPad) { | |
910 | Double_t mean = calPad->GetMean(); | |
911 | Double_t rms = calPad->GetRMS(); | |
912 | //printf("Vdrift::meanpad %f, rmspad %f\n",mean,rms); | |
913 | if(!((mean > 0.9) && (mean < 1.1) && (rms < 0.6))) ok = kFALSE; | |
914 | } | |
915 | else return kFALSE; | |
916 | } | |
917 | ||
918 | return ok; | |
919 | ||
920 | } | |
921 | //__________________________________________________________________________________________________________________________ | |
922 | Bool_t AliTRDPreprocessorOffline::ValidateT0(){ | |
923 | // | |
924 | // Update OCDB entry | |
925 | // | |
926 | ||
927 | AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet); | |
928 | AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad); | |
929 | if(calDet && calPad) { | |
930 | Double_t meandet = calDet->GetMean(); | |
931 | Double_t rmsdet = calDet->GetRMSRobust(); | |
932 | Double_t meanpad = calPad->GetMean(); | |
933 | //Double_t rmspad = calPad->GetRMS(); | |
934 | //printf("T0::minimum %f, rmsdet %f,meanpad %f, rmspad %f\n",meandet,rmsdet,meanpad,rmspad); | |
935 | if((meandet > -1.5) && (meandet < 5.0) && (rmsdet < 4.0) && (meanpad < 5.0) && (meanpad > -0.5)) return kTRUE; | |
936 | else return kFALSE; | |
937 | } | |
938 | else return kFALSE; | |
939 | ||
940 | } | |
941 | //__________________________________________________________________________________________________________________________ | |
942 | Bool_t AliTRDPreprocessorOffline::ValidatePRF() const{ | |
943 | // | |
944 | // Update OCDB entry | |
945 | // | |
946 | ||
947 | AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF); | |
948 | if(calPad) { | |
949 | Double_t meanpad = calPad->GetMean(); | |
950 | Double_t rmspad = calPad->GetRMS(); | |
951 | //printf("PRF::meanpad %f, rmspad %f\n",meanpad,rmspad); | |
952 | if((meanpad < 1.0) && (rmspad < 0.8)) return kTRUE; | |
953 | else return kFALSE; | |
954 | } | |
955 | else return kFALSE; | |
956 | ||
957 | ||
958 | } | |
959 | //__________________________________________________________________________________________________________________________ | |
960 | Bool_t AliTRDPreprocessorOffline::ValidateChamberStatus() const{ | |
00d203b6 | 961 | // |
962 | // Update OCDB entry | |
963 | // | |
964 | ||
a5dcf618 | 965 | AliTRDCalChamberStatus *calChamberStatus = (AliTRDCalChamberStatus *) fCalibObjects->At(kChamberStatus); |
966 | if(calChamberStatus) { | |
967 | Int_t detectormasked = 0; | |
968 | for(Int_t det = 0; det < 540; det++) { | |
969 | if(calChamberStatus->IsMasked(det)) detectormasked++; | |
970 | } | |
971 | //printf("Number of chambers masked %d\n",detectormasked); | |
972 | if(detectormasked > 29) return kFALSE; | |
973 | else return kTRUE; | |
00d203b6 | 974 | } |
975 | else return kFALSE; | |
976 | ||
b3fcfd96 | 977 | } |
4c865c34 | 978 | //_____________________________________________________________________________ |
979 | Int_t AliTRDPreprocessorOffline::GetVersion(TString name) const | |
980 | { | |
981 | // | |
982 | // Get version from the title | |
983 | // | |
984 | ||
985 | // Some patterns | |
986 | const Char_t *version = "Ver"; | |
987 | if(!strstr(name.Data(),version)) return -1; | |
988 | ||
989 | for(Int_t ver = 0; ver < 999999999; ver++) { | |
990 | ||
991 | TString vertry(version); | |
992 | vertry += ver; | |
993 | vertry += "Subver"; | |
994 | ||
995 | //printf("vertry %s and name %s\n",vertry.Data(),name.Data()); | |
996 | ||
997 | if(strstr(name.Data(),vertry.Data())) return ver; | |
998 | ||
999 | } | |
1000 | ||
1001 | return -1; | |
1002 | ||
1003 | } | |
1004 | ||
1005 | //_____________________________________________________________________________ | |
1006 | Int_t AliTRDPreprocessorOffline::GetSubVersion(TString name) const | |
1007 | { | |
1008 | // | |
1009 | // Get subversion from the title | |
1010 | // | |
1011 | ||
1012 | // Some patterns | |
1013 | const Char_t *subversion = "Subver"; | |
1014 | if(!strstr(name.Data(),subversion)) return -1; | |
1015 | ||
1016 | for(Int_t ver = 0; ver < 999999999; ver++) { | |
1017 | ||
1018 | TString vertry(subversion); | |
1019 | vertry += ver; | |
1020 | vertry += "Nz"; | |
1021 | ||
1022 | //printf("vertry %s and name %s\n",vertry.Data(),name.Data()); | |
1023 | ||
1024 | if(strstr(name.Data(),vertry.Data())) return ver; | |
1025 | ||
1026 | } | |
1027 | ||
1028 | return -1; | |
1029 | ||
1030 | } | |
01239968 | 1031 |