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