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