]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDPreprocessorOffline.cxx
Fixes for memory leaks (Chiara)
[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   TObjArray * array = (TObjArray*)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   TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
536   if (array){
537     TH2I *ch2d = (TH2I *) array->FindObject("CH2d");
538     if(!ch2d) return kFALSE;
539     fCH2d = (TH2I*)ch2d->Clone();
540     //fNEvents = (TH1I *) array->FindObject("NEvents");
541     //fAbsoluteGain = (TH2F *) array->FindObject("AbsoluteGain");
542     delete array;
543   }else{
544     TH2I *ch2d = (TH2I *) fcalib.Get("CH2d");
545     if(!ch2d) return kFALSE;
546     fCH2d = (TH2I*)ch2d->Clone();
547     //fNEvents = (TH1I *) fcalib.Get("NEvents");
548     //fAbsoluteGain = (TH2F *) fcalib.Get("AbsoluteGain");
549   }
550   fCH2d->SetDirectory(0);
551   //printf("title of CH2d %s\n",fCH2d->GetTitle());
552
553   return kTRUE;
554   
555 }
556 //_________________________________________________________________________________________________________________
557
558 Bool_t AliTRDPreprocessorOffline::ReadVdriftT0Global(const Char_t* fileName){
559   //
560   // read calibration entries from file
561   // 
562   if(fPH2d) return kTRUE;
563   TFile fcalib(fileName);
564   TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
565   if (array){
566     TProfile2D *ph2d = (TProfile2D *) array->FindObject("PH2d");
567     if(!ph2d) return kFALSE;
568     fPH2d = (TProfile2D*)ph2d->Clone();
569     //fNEvents = (TH1I *) array->FindObject("NEvents");
570     delete array;
571   }else{
572     TProfile2D *ph2d = (TProfile2D *) fcalib.Get("PH2d");
573     if(!ph2d) return kFALSE;
574     fPH2d = (TProfile2D*)ph2d->Clone();
575     //fNEvents = (TH1I *) fcalib.Get("NEvents");
576   }
577   fPH2d->SetDirectory(0);
578   //printf("title of PH2d %s\n",fPH2d->GetTitle());
579   
580   return kTRUE;
581   
582 }
583 //___________________________________________________________________________________________________________________
584
585 Bool_t AliTRDPreprocessorOffline::ReadVdriftLinearFitGlobal(const Char_t* fileName){
586   //
587   // read calibration entries from file
588   // 
589   if(fAliTRDCalibraVdriftLinearFit) return kTRUE;
590   TFile fcalib(fileName);
591   TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
592   if (array){
593     fAliTRDCalibraVdriftLinearFit = (AliTRDCalibraVdriftLinearFit *) array->FindObject("AliTRDCalibraVdriftLinearFit");
594     //fNEvents = (TH1I *) array->FindObject("NEvents");
595     delete array;
596   }else{
597     fAliTRDCalibraVdriftLinearFit = (AliTRDCalibraVdriftLinearFit *) fcalib.Get("AliTRDCalibraVdriftLinearFit");
598     //fNEvents = (TH1I *) fcalib.Get("NEvents");
599   }
600   if(!fAliTRDCalibraVdriftLinearFit) {
601     //printf("No AliTRDCalibraVdriftLinearFit\n");
602     return kFALSE;
603   }
604   return kTRUE;
605   
606 }
607 //_____________________________________________________________________________________________________________
608 Bool_t AliTRDPreprocessorOffline::ReadExbAltFitGlobal(const Char_t* fileName){
609   //
610   // read calibration entries from file
611   // 
612   if(fAliTRDCalibraExbAltFit) return kTRUE;
613   TFile fcalib(fileName);
614   TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
615   if (array){
616     fAliTRDCalibraExbAltFit = (AliTRDCalibraExbAltFit *) array->FindObject("AliTRDCalibraExbAltFit");
617     //fNEvents = (TH1I *) array->FindObject("NEvents");
618     delete array;
619   }else{
620     fAliTRDCalibraExbAltFit = (AliTRDCalibraExbAltFit *) fcalib.Get("AliTRDCalibraExbAltFit");
621     //fNEvents = (TH1I *) fcalib.Get("NEvents");
622   }
623   if(!fAliTRDCalibraExbAltFit) {
624     //printf("No AliTRDCalibraExbAltFit\n");
625     return kFALSE;
626   }
627   return kTRUE;
628   
629 }
630 //_____________________________________________________________________________________________________________
631
632 Bool_t AliTRDPreprocessorOffline::ReadPRFGlobal(const Char_t* fileName){
633   //
634   // read calibration entries from file
635   // 
636   if(fPRF2d) return kTRUE;
637   TFile fcalib(fileName);
638   TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
639   if (array){
640     TProfile2D *prf2d = (TProfile2D *) array->FindObject("PRF2d");
641     if(!prf2d) return kFALSE;
642     fPRF2d = (TProfile2D*)prf2d->Clone();
643     //fNEvents = (TH1I *) array->FindObject("NEvents");
644   }else{
645     TProfile2D *prf2d = (TProfile2D *) fcalib.Get("PRF2d");
646     if(!prf2d) return kFALSE;
647     fPRF2d = (TProfile2D*)prf2d->Clone();
648     //fNEvents = (TH1I *) fcalib.Get("NEvents");
649   }
650   fPRF2d->SetDirectory(0);
651   //printf("title of PRF2d %s\n",fPRF2d->GetTitle());
652   
653   return kTRUE;
654
655 }
656 //__________________________________________________________________________________________________________
657
658 Bool_t AliTRDPreprocessorOffline::AnalyzeGain(){
659   //
660   // Analyze gain - produce the calibration objects
661   //
662
663   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
664   calibra->ChooseMethod(fMethodeGain);
665   calibra->SetBeginFitCharge(fBeginFitCharge);
666   calibra->SetFitOutliersChargeLow(fOutliersFitChargeLow);
667   calibra->SetFitOutliersChargeHigh(fOutliersFitChargeHigh);
668   calibra->SetMinEntries(fMinStatsGain); // If there is less than 1000 entries in the histo: no fit
669   calibra->AnalyseCH(fCH2d);
670
671   Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(0))
672     + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(0));
673   Int_t nbfit       = calibra->GetNumberFit();
674   Int_t nbE         = calibra->GetNumberEnt();
675
676
677   Bool_t ok = kFALSE;
678   Bool_t meanother = kFALSE;
679   // enough statistics
680   if ((nbtg >                  0) && 
681       (nbfit        >= 0.5*nbE) && (nbE > 30)) {
682     // create the cal objects
683     if(!fBackCorrectGain) {
684       calibra->PutMeanValueOtherVectorFit(1,kTRUE);
685       meanother = kTRUE;
686     }
687     TObjArray object           = calibra->GetVectorFit();
688     AliTRDCalDet *calDetGain   = calibra->CreateDetObjectGain(&object,meanother);
689     TH1F *coefGain  = calDetGain->MakeHisto1DAsFunctionOfDet();
690     // Put them in the array
691     fCalibObjects->AddAt(calDetGain,kGain);
692     fPlots->AddAt(coefGain,kGain);
693     //
694     ok = kTRUE;
695   }
696   else {
697     fNotEnoughStatisticsForTheGain = kTRUE;
698     Int_t minStatsGain = fMinStatsGain*30;
699     calibra->SetMinEntries(minStatsGain); // Because we do it for all, we increase this
700     Double_t gainoverallnotnormalized =  calibra->AnalyseCHAllTogether(fCH2d);
701     if(fCalDetGainUsed && (gainoverallnotnormalized > 0.0)) {
702       AliTRDCalDet *calDetGain = new AliTRDCalDet(*fCalDetGainUsed);
703       Double_t oldmean = fCalDetGainUsed->CalcMean(kFALSE);
704       //printf("oldmean %f\n",oldmean);
705       if(oldmean > 0.0)  {
706         Double_t scalefactor = calibra->GetScaleFactorGain();
707         //printf("Correction factor %f\n",gainoverallnotnormalized*scalefactor);
708         calDetGain->Multiply(gainoverallnotnormalized*scalefactor/oldmean);
709         //printf("newmean %f\n",calDetGain->CalcMean(kFALSE));
710         TH1F *coefGain  = calDetGain->MakeHisto1DAsFunctionOfDet();
711         fCalibObjects->AddAt(calDetGain,kGain);
712         fPlots->AddAt(coefGain,kGain);
713         // 
714         ok = kTRUE;
715         fStatusNeg = fStatusNeg | kGainNotEnoughStatsButFill;
716       }
717       else {
718         fStatusPos = fStatusPos | kGainErrorOld;
719       }      
720     }
721     else {
722       if(gainoverallnotnormalized <= 0.0) fStatusNeg = fStatusNeg | kGainNotEnoughStatsNotFill;
723       if(!fCalDetGainUsed) fStatusPos = fStatusPos | kGainErrorOld;
724     }
725   }
726   
727   calibra->ResetVectorFit();
728   
729   return ok;
730   
731 }
732 //_____________________________________________________________________________________________________
733 Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftT0(){
734   //
735   // Analyze VdriftT0 - produce the calibration objects
736   //
737
738   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
739   calibra->SetT0Shift0(fT0Shift0);
740   calibra->SetT0Shift1(fT0Shift1);
741   calibra->SetMinEntries(fMinStatsVdriftT0PH); // If there is less than 1000 entries in the histo: no fit
742   calibra->AnalysePH(fPH2d); 
743   //calibra->SetDebugLevel(2);
744
745   Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1))
746     + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1));
747   Int_t nbfit       = calibra->GetNumberFit();
748   Int_t nbfitSuccess = calibra->GetNumberFitSuccess();
749   Int_t nbE         = calibra->GetNumberEnt();
750
751   //printf("nbtg %d, nbfit %d, nbE %d, nbfitSuccess %d\n",nbtg,nbfit,nbE,nbfitSuccess);
752
753   Bool_t ok = kFALSE;
754   if ((nbtg >                  0) && 
755       (nbfit        >= 0.5*nbE) && (nbE > 30) && (nbfitSuccess > 30)) {
756     //printf("nbtg %d, nbfit %d, nbE %d, nbfitSucess %d\n",nbtg,nbfit,nbE,nbfitSuccess);
757     //printf("Pass the cut for VdriftT0\n");
758     // create the cal objects
759     calibra->RemoveOutliers(1,kFALSE);
760     calibra->PutMeanValueOtherVectorFit(1,kFALSE);
761     calibra->RemoveOutliers2(kFALSE);
762     calibra->PutMeanValueOtherVectorFit2(1,kFALSE);
763     //
764     TObjArray object = calibra->GetVectorFit();
765     AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object);
766     TH1F *coefVdriftPH  = calDetVdrift->MakeHisto1DAsFunctionOfDet();
767     AliTRDCalPad *calPadVdrift = (AliTRDCalPad *)calibra->CreatePadObjectVdrift(&object,calDetVdrift);
768     TH1F *coefPadVdrift   = calPadVdrift->MakeHisto1D();
769     object       = calibra->GetVectorFit2();
770     AliTRDCalDet *calDetT0  = calibra->CreateDetObjectT0(&object);
771     TH1F *coefT0  = calDetT0->MakeHisto1DAsFunctionOfDet();
772     AliTRDCalPad *calPadT0 = (AliTRDCalPad *)calibra->CreatePadObjectT0(&object,calDetT0);
773     TH1F *coefPadT0  = calPadT0->MakeHisto1D();
774     // Put them in the array
775     fCalibObjects->AddAt(calDetT0,kT0PHDet);
776     fCalibObjects->AddAt(calDetVdrift,kVdriftPHDet);
777     fCalibObjects->AddAt(calPadT0,kT0PHPad);
778     fCalibObjects->AddAt(calPadVdrift,kVdriftPHPad);
779     fPlots->AddAt(coefVdriftPH,kVdriftPHDet);
780     fPlots->AddAt(coefT0,kT0PHDet);
781     fPlots->AddAt(coefPadVdrift,kVdriftPHPad);
782     fPlots->AddAt(coefPadT0,kT0PHPad);
783     //
784     ok = kTRUE;
785   }
786   else {
787     //printf("Not enough stats timeoffset\n");
788     fStatusNeg = fStatusNeg | kTimeOffsetNotEnoughStatsNotFill;
789   }
790   calibra->ResetVectorFit();
791  
792   return ok;
793   
794 }
795 //____________________________________________________________________________________________________________________
796 Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftLinearFit(){
797   //
798   // Analyze vdrift linear fit - produce the calibration objects
799   //
800
801   //printf("Analyse linear fit\n");
802
803   
804   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
805   calibra->SetCalDetVdriftExB(fCalDetVdriftUsed,fCalDetExBUsed);
806   calibra->SetMinEntries(fMinStatsVdriftLinear); // If there is less than 1000 entries in the histo: no fit
807   //printf("The mean stat is by %d for VdriftLinear\n",fMinStatsVdriftLinear);
808   //fAliTRDCalibraVdriftLinearFit->SetSeeDetector(0);
809   //fAliTRDCalibraVdriftLinearFit->SetDebugLevel(1);
810   //printf("Fill PE Array\n");
811   fAliTRDCalibraVdriftLinearFit->SetRobustFit(fRobustFitDriftVelocity);
812   fAliTRDCalibraVdriftLinearFit->SetMinNumberOfPointsForFit(fMinNbOfPointVdriftFit);
813   if(!fAlternativeVdrfitFit)
814     fAliTRDCalibraVdriftLinearFit->FillPEArray();
815   else
816     fAliTRDCalibraVdriftLinearFit->FillPEArray2();
817   //printf("AliTRDCalibraFit\n");
818   calibra->AnalyseLinearFitters(fAliTRDCalibraVdriftLinearFit);
819   //printf("After\n");
820
821   //Int_t nbtg        = 540;
822   Int_t nbfit       = calibra->GetNumberFit();
823   Int_t nbE         = calibra->GetNumberEnt();
824
825   
826   Bool_t ok = kFALSE;
827   // enough statistics
828   if ((nbfit        >= 0.5*nbE) && (nbE > 30)) {
829     // create the cal objects
830     //calibra->RemoveOutliers(1,kTRUE);
831     calibra->PutMeanValueOtherVectorFit(1,kTRUE);
832     //calibra->RemoveOutliers2(kTRUE);
833     calibra->PutMeanValueOtherVectorFit2(1,kTRUE);
834     //
835     TObjArray object  = calibra->GetVectorFit();
836     AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object,kTRUE);
837     TH1F *coefDriftLinear      = calDetVdrift->MakeHisto1DAsFunctionOfDet();
838     object                     = calibra->GetVectorFit2();
839     AliTRDCalDet *calDetLorentz = calibra->CreateDetObjectLorentzAngle(&object);
840     TH1F *coefLorentzAngle = calDetLorentz->MakeHisto1DAsFunctionOfDet();
841     //if(!calDetLorentz) printf("No lorentz created\n");
842     // Put them in the array
843     fCalibObjects->AddAt(calDetVdrift,kVdriftLinear);
844     fCalibObjects->AddAt(calDetLorentz,kLorentzLinear);
845     fPlots->AddAt(coefDriftLinear,kVdriftLinear);
846     fPlots->AddAt(coefLorentzAngle,kLorentzLinear);
847     //
848     ok = kTRUE;
849   }
850   else {
851     fNotEnoughStatisticsForTheVdriftLinear = kTRUE;
852     Int_t minNumberOfEntriesForAll = fMinStatsVdriftLinear*30;
853     calibra->SetMinEntries(minNumberOfEntriesForAll); // Because we do it for all, we increase this
854     Double_t vdriftoverall = -100.0;
855     Double_t exboverall = 100.0;
856     calibra->AnalyseLinearFittersAllTogether(fAliTRDCalibraVdriftLinearFit,vdriftoverall,exboverall);
857     //printf("Found mean vdrift %f and exb %f\n",vdriftoverall,exboverall);
858     if(fCalDetVdriftUsed && (vdriftoverall > 0.0) && (exboverall < 70.0)) {
859       AliTRDCalDet *calDetVdrift = new AliTRDCalDet(*fCalDetVdriftUsed);
860       AliTRDCalDet *calDetLorentz = new AliTRDCalDet(*fCalDetExBUsed);
861       Double_t oldmeanvdrift = fCalDetVdriftUsed->CalcMean(kFALSE);
862       Double_t oldmeanexb = fCalDetExBUsed->CalcMean(kFALSE);
863       //printf("oldmean %f\n",oldmean);
864       if((oldmeanvdrift > 0.0) && (oldmeanexb < 70.0))  {
865         //printf("Correction factor %f\n",vdriftoverall);
866         calDetVdrift->Multiply(vdriftoverall/oldmeanvdrift);
867         if(TMath::Abs(oldmeanexb) > 0.0001) calDetLorentz->Multiply(exboverall/oldmeanexb);
868         //printf("newmean %f\n",calDetVdrift->CalcMean(kFALSE));
869         TH1F *coefDriftLinear  = calDetVdrift->MakeHisto1DAsFunctionOfDet();
870         TH1F *coefLorentzAngle = calDetLorentz->MakeHisto1DAsFunctionOfDet();
871         // Put them in the array
872         fCalibObjects->AddAt(calDetVdrift,kVdriftLinear);
873         fCalibObjects->AddAt(calDetLorentz,kLorentzLinear);
874         fPlots->AddAt(coefDriftLinear,kVdriftLinear);
875         fPlots->AddAt(coefLorentzAngle,kLorentzLinear);
876         // 
877         ok = kTRUE;
878         fStatusNeg = fStatusNeg | kVdriftNotEnoughStatsButFill;
879       }
880       else {
881         if(oldmeanvdrift) fStatusPos = fStatusPos | kVdriftErrorOld;
882         if(oldmeanexb) fStatusPos = fStatusPos | kExBErrorOld;
883       }      
884     }
885     else {
886       if((vdriftoverall <= 0.0) && (exboverall > 70.0)) fStatusNeg = fStatusNeg | kVdriftNotEnoughStatsNotFill;
887       if(!fCalDetVdriftUsed) fStatusPos = fStatusPos | kVdriftErrorOld;
888       if(!fCalDetExBUsed) fStatusPos = fStatusPos | kExBErrorOld;
889     }
890   }
891   
892   calibra->ResetVectorFit();
893   
894   return ok;
895   
896 }
897 //________________________________________________________________________________________________________________
898
899 Bool_t AliTRDPreprocessorOffline::AnalyzeExbAltFit(){
900   //
901   // Analyze vdrift linear fit - produce the calibration objects
902   //
903
904   //printf("Analyse linear fit\n");
905
906   
907   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
908   calibra->SetMinEntries(fMinStatsVdriftLinear); // If there is less than 1000 entries in the histo: no fit
909   //printf("Fill PE Array\n");
910   fAliTRDCalibraExbAltFit->SetRobustFit(fRobustFitExbAlt);
911   if(!fAlternativeExbAltFit)
912     fAliTRDCalibraExbAltFit->FillPEArray();
913   else
914     fAliTRDCalibraExbAltFit->FillPEArray2();
915   //printf("AliTRDCalibraFit\n");
916   calibra->AnalyseExbAltFit(fAliTRDCalibraExbAltFit);
917   //printf("After\n");
918
919   //Int_t nbtg        = 540;
920   Int_t nbfit       = calibra->GetNumberFit();
921   Int_t nbE         = calibra->GetNumberEnt();
922
923   
924   Bool_t ok = kFALSE;
925   // enough statistics
926   if ((nbfit        >= 0.5*nbE) && (nbE > 30)) {
927     // create the cal objects
928     //calibra->RemoveOutliers(1,kTRUE);
929     calibra->PutMeanValueOtherVectorFit2(1,kTRUE);
930     //
931     TObjArray object  = calibra->GetVectorFit2();
932     AliTRDCalDet *calDetLorentz = calibra->CreateDetObjectExbAlt(&object);
933     TH1F *coefLorentzAngle = calDetLorentz->MakeHisto1DAsFunctionOfDet();
934     //if(!calDetLorentz) printf("No lorentz created\n");
935     // Put them in the array
936     fCalibObjects->AddAt(calDetLorentz,kExbAlt);
937     fPlots->AddAt(coefLorentzAngle,kExbAlt);
938     //
939     ok = kTRUE;
940   }
941   
942   calibra->ResetVectorFit();
943   
944   return ok;
945   
946 }
947 //________________________________________________________________________________________________________________
948
949 Bool_t AliTRDPreprocessorOffline::AnalyzePRF(){
950   //
951   // Analyze PRF - produce the calibration objects
952   //
953
954   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
955   calibra->SetMinEntries(fMinStatsPRF); // If there is less than 1000 entries in the histo: no fit
956   calibra->AnalysePRFMarianFit(fPRF2d);
957
958   Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2))
959     + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2));
960   Int_t nbfit       = calibra->GetNumberFit();
961   Int_t nbE         = calibra->GetNumberEnt();
962
963   
964   Bool_t ok = kFALSE;
965   // enough statistics
966   if ((nbtg >                  0) && 
967       (nbfit        >= 0.95*nbE) && (nbE > 30)) {
968     // create the cal objects
969     TObjArray object  = calibra->GetVectorFit();
970     AliTRDCalPad *calPadPRF = (AliTRDCalPad*) calibra->CreatePadObjectPRF(&object);
971     TH1F *coefPRF           = calPadPRF->MakeHisto1D();
972     // Put them in the array
973     fCalibObjects->AddAt(calPadPRF,kPRF);
974     fPlots->AddAt(coefPRF,kPRF);
975     //
976     ok = kTRUE;
977   }
978   
979   calibra->ResetVectorFit();
980   
981   return ok;
982   
983 }
984
985 //_____________________________________________________________________________
986 Bool_t AliTRDPreprocessorOffline::AnalyzePHQ(Int_t startRunNumber)
987 {
988   //
989   //Produce PHQ calibration results
990   //
991   TList *lout = 0x0;
992   TTreeSRedirector *calibStream = 0x0;
993   if(IsDebugPHQon()){
994     lout = new TList;
995     lout->SetOwner();
996
997     calibStream = new TTreeSRedirector(Form("TRDCalibStream_%010d.root", startRunNumber));
998   }
999
1000   for(Int_t iter=0; iter<AliTRDdEdxCalibUtils::GetHistArray()->GetSize(); iter++){
1001     THnBase *hi = (THnBase*) AliTRDdEdxCalibUtils::GetHistAt(iter);
1002     TObjArray *obji = AliTRDdEdxCalibUtils::HistToObj(hi, startRunNumber, lout, calibStream);
1003     //printf("test analyze %s\n", obji->GetName());
1004     AliTRDdEdxCalibUtils::GetObjArray()->AddAt(obji, iter);
1005   }
1006
1007   fCalibObjects->AddAt(AliTRDdEdxCalibUtils::GetObjArray(), kPHQ);
1008
1009   if(lout){
1010     TFile *fout=new TFile(Form("TRDCalibList_%010d.root", startRunNumber),"recreate");
1011     fout->cd();
1012     lout->Write();
1013     fout->Save();
1014     fout->Close();
1015     delete fout;
1016   }
1017   delete calibStream;
1018   delete lout;
1019
1020   return kTRUE;
1021 }
1022
1023 //_____________________________________________________________________________
1024 Bool_t AliTRDPreprocessorOffline::AnalyzeChamberStatus()
1025 {
1026   //
1027   // Produce AliTRDCalChamberStatus out of calibration results
1028   //
1029   
1030   // set up AliTRDCalibChamberStatus
1031   AliTRDCalibChamberStatus *chamberStatus = new AliTRDCalibChamberStatus();
1032   chamberStatus->SetSparseI(fSparse);
1033   chamberStatus->AnalyseHisto(fMinStatsChamberStatus, fMinSingleStatsChamberStatus);
1034   // get AliTRDCalChamberStatus
1035   AliTRDCalChamberStatus *calChamberStatus = chamberStatus->GetCalChamberStatus();
1036
1037   // get calibration objects
1038   AliTRDCalDet *calDetGain   = (AliTRDCalDet *) fCalibObjects->At(kGain);
1039   AliTRDCalDet *calDetVDrift = (AliTRDCalDet *) fCalibObjects->At(kVdriftLinear);
1040   AliTRDCalDet *calDetExB    = (AliTRDCalDet *) fCalibObjects->At(kLorentzLinear);
1041
1042   // Check
1043   if((!calDetGain) || (!calDetVDrift) || (!fCH2d) || (!calDetExB) || (!calChamberStatus)) return kFALSE;
1044
1045   // Gain
1046   Double_t gainmean   = calDetGain->GetMean();
1047   Double_t vdriftmean = calDetVDrift->GetMean();
1048   Double_t exbmean    = calDetExB->GetMean();
1049
1050   Double_t gainrms    = calDetGain->GetRMSRobust();
1051   Double_t vdriftrms  = calDetVDrift->GetRMSRobust();
1052   Double_t exbrms     = calDetExB->GetRMSRobust();
1053
1054   //printf("Gain mean: %f, rms: %f\n",gainmean,gainrms);
1055   //printf("Vdrift mean: %f, rms: %f\n",vdriftmean,vdriftrms);
1056   //printf("ExB mean: %f, rms: %f\n",exbmean,exbrms);
1057
1058   // Check
1059   if((TMath::Abs(gainrms) < 0.001) || (TMath::Abs(vdriftrms) < 0.001) || (TMath::Abs(exbrms) < 0.0000001)) return kFALSE;
1060
1061   // Take mean each SM
1062   Double_t *gainmeanSM = new Double_t[18];
1063   Double_t *vdriftmeanSM = new Double_t[18];
1064   Double_t *exbmeanSM = new Double_t[18];
1065   //Double_t *t0meanSM = new Double_t[18];
1066   for(Int_t sm=0; sm< 18; sm++) {
1067     gainmeanSM[sm] = calDetGain->GetMeanSM(kFALSE,sm);
1068     vdriftmeanSM[sm] = calDetVDrift->GetMeanSM(kFALSE,sm);
1069     exbmeanSM[sm] = calDetExB->GetMeanSM(kFALSE,sm);
1070     //t0meanSM[sm] = calDetGain->GetMeanSM(kFALSE);
1071   }
1072
1073
1074   // mask chambers with empty gain entries
1075   //Int_t counter = 0;
1076   for (Int_t idet = 0; idet < 540; idet++) {
1077
1078     // ch2d
1079     TH1I *projch =  (TH1I *) fCH2d->ProjectionX("projch",idet+1,idet+1,(Option_t *)"e");
1080     Double_t entries = projch->GetEntries();
1081     //printf("Number of entries %f for det %d\n",entries,idet);
1082
1083     // sm number
1084     Int_t smnumber = (Int_t) idet/30;
1085
1086     // gain
1087     Double_t gain = calDetGain->GetValue(idet);
1088
1089     // vdrift
1090     Double_t vdrift = calDetVDrift->GetValue(idet);
1091
1092     // exb
1093     Double_t exb = calDetExB->GetValue(idet);
1094
1095
1096     if( (entries<50 && !calChamberStatus->IsNoData(idet))  ||
1097         TMath::Abs(gainmean-gain) > (fRMSBadCalibratedGain*gainrms)          ||
1098         TMath::Abs(vdriftmean-vdrift) > (fRMSBadCalibratedVdrift*vdriftrms)    ||
1099         TMath::Abs(exbmean-exb) > (fRMSBadCalibratedExB*exbrms) ) {
1100      
1101       //printf(" chamber det %03d masked \n",idet);
1102       //printf(" gainmean %f and gain %f, gainrms %f \n",gainmean,gain,gainrms);
1103       //printf(" vdriftmean %f and vdrift %f, vdriftrms %f \n",vdriftmean,vdrift,vdriftrms);
1104       //printf(" exbmean %f and exb %f, exbrms %f \n",exbmean,exb,exbrms);
1105       
1106       calChamberStatus->SetStatus(idet,AliTRDCalChamberStatus::kBadCalibrated);
1107       //counter++;
1108     }
1109
1110     if(TMath::Abs(gainmeanSM[smnumber]-gain) < 0.000001  ||
1111        TMath::Abs(vdriftmeanSM[smnumber]-vdrift) < 0.000001 ||
1112        TMath::Abs(exbmeanSM[smnumber]-exb) < 0.000001) {
1113       
1114       //printf(" chamber det %03d notcalibrated sm %d \n",idet,smnumber);
1115       //printf(" gainmeanSM %f and gain %f\n",gainmeanSM[smnumber],gain);
1116       //printf(" vdriftmeanSM %f and vdrift %f \n",vdriftmeanSM[smnumber],vdrift);
1117       //printf(" exbmeanSM %f and exb %f \n",exbmeanSM[smnumber],exb);
1118
1119       calChamberStatus->SetStatus(idet,AliTRDCalChamberStatus::kNotCalibrated);
1120     }
1121
1122
1123     delete projch;
1124     
1125    }
1126
1127   // Security
1128   for(Int_t sm=0; sm < 18; sm++) {
1129     Int_t smnodata   = 0;
1130     Int_t smbadcalib = 0;
1131     for(Int_t det = 0; det < 30; det++){
1132       Int_t detector = sm*30+det;
1133       if(calChamberStatus->IsNoData(detector)) smnodata++;
1134       else {
1135         if(calChamberStatus->IsBadCalibrated(detector)) smbadcalib++;
1136       }
1137     }
1138     fNoData[sm]  = smnodata;
1139     fBadCalib[sm]= smbadcalib;
1140     //printf("No Data %d, bad calibrated %d for %d\n",fNoData[sm],fBadCalib[sm],sm);
1141   }
1142
1143   // delete
1144   delete []gainmeanSM;
1145   delete []vdriftmeanSM;
1146   delete []exbmeanSM;
1147   
1148   // Security
1149   //   for(Int_t sm=0; sm < 18; sm++) {
1150   //     Int_t counter = 0;
1151   //     for(Int_t det = 0; det < 30; det++){
1152   //       Int_t detector = sm*30+det;
1153   //       if(calChamberStatus->IsBadCalibrated(detector)) counter++;
1154   //     }
1155   //     if(counter >= 20) {
1156   //       for(Int_t det = 0; det < 30; det++){
1157   //    Int_t detector = sm*30+det;
1158   //    calChamberStatus->SetStatus(detector,AliTRDCalChamberStatus::kGood);
1159   //       }
1160   //     }
1161   //  }
1162
1163    fCalibObjects->AddAt(calChamberStatus,kChamberStatus);
1164    return kTRUE;
1165
1166  }
1167
1168
1169  //________________________________________________________________________________________________
1170  void AliTRDPreprocessorOffline::CorrectFromDetGainUsed() {
1171    //
1172    // Correct from the gas gain used afterwards
1173    //
1174    AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
1175    if(!calDetGain) return;
1176
1177    // Calculate mean
1178    Double_t mean = 0.0;
1179    Int_t nbdet = 0;
1180
1181    for(Int_t det = 0; det < 540; det++) {
1182
1183      Float_t gaininit = fCalDetGainUsed->GetValue(det);
1184      Float_t gainout = calDetGain->GetValue(det);
1185
1186
1187      if(TMath::Abs(gainout-1.0) > 0.000001) {
1188        mean += (gaininit*gainout);
1189        nbdet++;
1190      }  
1191    }
1192    if(nbdet > 0) mean = mean/nbdet;
1193
1194    for(Int_t det = 0; det < 540; det++) {
1195
1196      Float_t gaininit = fCalDetGainUsed->GetValue(det);
1197      Float_t gainout = calDetGain->GetValue(det);
1198
1199      if(TMath::Abs(gainout-1.0) > 0.000001) {
1200        Double_t newgain = gaininit*gainout;
1201        if(newgain < 0.1) newgain = 0.1;
1202        if(newgain > 1.9) newgain = 1.9;
1203        calDetGain->SetValue(det,newgain);
1204      }
1205      else {
1206        Double_t newgain = mean;
1207        if(newgain < 0.1) newgain = 0.1;
1208        if(newgain > 1.9) newgain = 1.9;
1209        calDetGain->SetValue(det,newgain);
1210      }
1211    }
1212
1213
1214  }
1215  //________________________________________________________________________________________________
1216  void AliTRDPreprocessorOffline::CorrectFromDetVdriftUsed() {
1217    //
1218    // Correct from the drift velocity
1219    //
1220
1221    //printf("Correct for vdrift\n");
1222
1223    AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
1224    if(!calDetGain) return;
1225
1226    Int_t detVdrift = kVdriftPHDet;
1227    if(fMethodSecond) detVdrift = kVdriftLinear;
1228    AliTRDCalDet *calDetVdrift = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
1229    if(!calDetVdrift) return;
1230
1231    // Calculate mean
1232    if(!fNotEnoughStatisticsForTheVdriftLinear) {
1233      for(Int_t det = 0; det < 540; det++) {
1234        
1235        Float_t vdriftinit = fCalDetVdriftUsed->GetValue(det);
1236        Float_t vdriftout = calDetVdrift->GetValue(det);
1237        
1238        Float_t gain = calDetGain->GetValue(det);
1239        if(vdriftout > 0.0) gain = gain*vdriftinit/vdriftout;
1240        if(gain < 0.1) gain = 0.1;
1241        if(gain > 1.9) gain = 1.9;
1242        calDetGain->SetValue(det,gain);
1243      }
1244    }
1245    else {
1246      
1247      Float_t vdriftinit = fCalDetVdriftUsed->CalcMean(kFALSE);
1248      Float_t vdriftout = calDetVdrift->CalcMean(kFALSE);
1249      Float_t factorcorrectif = 1.0;
1250      if(vdriftout > 0.0) factorcorrectif = vdriftinit/vdriftout;
1251      for(Int_t det = 0; det < 540; det++) {
1252        Float_t gain = calDetGain->GetValue(det);
1253        gain = gain*factorcorrectif;
1254        if(gain < 0.1) gain = 0.1;
1255        if(gain > 1.9) gain = 1.9;
1256        calDetGain->SetValue(det,gain);
1257      }
1258      
1259    }
1260    
1261  }
1262 //_________________________________________________________________________________________________________________
1263  void AliTRDPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, AliCDBStorage *storage){
1264    //
1265    // Update OCDB entry
1266    //
1267
1268    AliCDBMetaData *metaData= new AliCDBMetaData();
1269    metaData->SetObjectClassName("AliTRDCalDet");
1270    metaData->SetResponsible("Raphaelle Bailhache");
1271    metaData->AddDateToComment();
1272    metaData->SetBeamPeriod(1);
1273
1274    AliCDBId id1("TRD/Calib/ChamberGainFactor", startRunNumber, endRunNumber);
1275    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain);
1276    if(calDet) storage->Put(calDet, id1, metaData);
1277
1278
1279  }
1280  //___________________________________________________________________________________________________________________
1281  void AliTRDPreprocessorOffline::UpdateOCDBExB(Int_t startRunNumber, Int_t endRunNumber, AliCDBStorage *storage){
1282    //
1283    // Update OCDB entry
1284    //
1285
1286    Int_t detExB = kLorentzLinear;
1287    if(!fMethodSecond) return;
1288
1289    //printf("Pass\n");
1290
1291    AliCDBMetaData *metaData= new AliCDBMetaData();
1292    metaData->SetObjectClassName("AliTRDCalDet");
1293    metaData->SetResponsible("Raphaelle Bailhache");
1294    metaData->AddDateToComment();
1295    metaData->SetBeamPeriod(1);
1296
1297    AliCDBId id1("TRD/Calib/ChamberExB", startRunNumber, endRunNumber);
1298    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detExB);
1299    if(calDet) storage->Put(calDet, id1, metaData);
1300    //if(!calDet) printf("No caldet\n");
1301
1302  }
1303 //___________________________________________________________________________________________________________________
1304 void AliTRDPreprocessorOffline::UpdateOCDBExBAlt(Int_t startRunNumber, Int_t endRunNumber, AliCDBStorage *storage){
1305   //
1306   // Update OCDB entry
1307   //
1308
1309   Int_t detExB = kExbAlt;
1310   if(!fMethodSecond) return;
1311
1312   //printf("Pass\n");
1313
1314   AliCDBMetaData *metaData= new AliCDBMetaData();
1315   metaData->SetObjectClassName("AliTRDCalDet");
1316   metaData->SetResponsible("Theo Rascanu");
1317   metaData->AddDateToComment();
1318   metaData->SetBeamPeriod(1);
1319
1320   AliCDBId id1("TRD/Calib/ChamberExBAlt", startRunNumber, endRunNumber);
1321   AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detExB);
1322   if(calDet) storage->Put(calDet, id1, metaData);
1323   //if(!calDet) printf("No caldet\n");
1324
1325 }
1326  //___________________________________________________________________________________________________________________
1327  void AliTRDPreprocessorOffline::UpdateOCDBVdrift(Int_t startRunNumber, Int_t endRunNumber, AliCDBStorage* storage){
1328    //
1329    // Update OCDB entry
1330    //
1331
1332    Int_t detVdrift = kVdriftPHDet;
1333
1334    if(fMethodSecond) detVdrift = kVdriftLinear;
1335
1336    AliCDBMetaData *metaData= new AliCDBMetaData();
1337    metaData->SetObjectClassName("AliTRDCalDet");
1338    metaData->SetResponsible("Raphaelle Bailhache");
1339    metaData->AddDateToComment();
1340    metaData->SetBeamPeriod(1);
1341
1342    AliCDBId id1("TRD/Calib/ChamberVdrift", startRunNumber, endRunNumber);
1343    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
1344    if(calDet) storage->Put(calDet, id1, metaData);
1345
1346    //
1347
1348    if(!fMethodSecond) {
1349
1350      AliCDBMetaData *metaDataPad= new AliCDBMetaData();
1351      metaDataPad->SetObjectClassName("AliTRDCalPad");
1352      metaDataPad->SetResponsible("Raphaelle Bailhache");
1353      metaDataPad->AddDateToComment();
1354      metaDataPad->SetBeamPeriod(1);
1355
1356      AliCDBId id1Pad("TRD/Calib/LocalVdrift", startRunNumber, endRunNumber);
1357      AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad);
1358      if(calPad) storage->Put(calPad, id1Pad, metaDataPad);
1359
1360    }
1361
1362  }
1363  //________________________________________________________________________________________________________________________
1364  void AliTRDPreprocessorOffline::UpdateOCDBT0(Int_t startRunNumber, Int_t endRunNumber, AliCDBStorage* storage){
1365    //
1366    // Update OCDB entry
1367    //
1368
1369    AliCDBMetaData *metaData= new AliCDBMetaData();
1370    metaData->SetObjectClassName("AliTRDCalDet");
1371    metaData->SetResponsible("Raphaelle Bailhache");
1372    metaData->AddDateToComment();
1373    metaData->SetBeamPeriod(1);
1374
1375    AliCDBId id1("TRD/Calib/ChamberT0", startRunNumber, endRunNumber);
1376    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet);
1377    if(calDet) storage->Put(calDet, id1, metaData);
1378
1379    //
1380
1381    AliCDBMetaData *metaDataPad= new AliCDBMetaData();
1382    metaDataPad->SetObjectClassName("AliTRDCalPad");
1383    metaDataPad->SetResponsible("Raphaelle Bailhache");
1384    metaDataPad->AddDateToComment();
1385    metaDataPad->SetBeamPeriod(1);
1386
1387    AliCDBId id1Pad("TRD/Calib/LocalT0", startRunNumber, endRunNumber);
1388    AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad);
1389    if(calPad) storage->Put(calPad, id1Pad, metaDataPad);
1390
1391
1392
1393  }
1394  //_________________________________________________________________________________________________________________
1395  void AliTRDPreprocessorOffline::UpdateOCDBPRF(Int_t startRunNumber, Int_t endRunNumber, AliCDBStorage* storage){
1396    //
1397    // Update OCDB entry
1398    //
1399
1400    AliCDBMetaData *metaData= new AliCDBMetaData();
1401    metaData->SetObjectClassName("AliTRDCalPad");
1402    metaData->SetResponsible("Raphaelle Bailhache");
1403    metaData->AddDateToComment();
1404    metaData->SetBeamPeriod(1);
1405
1406
1407    AliCDBId id1("TRD/Calib/PRFWidth", startRunNumber, endRunNumber);
1408    AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF);
1409    if(calPad) storage->Put(calPad, id1, metaData);
1410
1411
1412  }
1413 //_________________________________________________________________________________________________________________
1414 void AliTRDPreprocessorOffline::UpdateOCDBPHQ(Int_t startRunNumber, Int_t endRunNumber, AliCDBStorage* storage)
1415 {
1416   //
1417   // Update OCDB entry
1418   //
1419   AliCDBMetaData *metaData= new AliCDBMetaData();
1420   metaData->SetObjectClassName("TObjArray");
1421   metaData->SetResponsible("Raphaelle Bailhache and Xianguo Lu");
1422   metaData->AddDateToComment();
1423   metaData->SetBeamPeriod(1);
1424
1425   AliCDBId id1("TRD/Calib/PHQ", startRunNumber, endRunNumber);
1426   TObjArray *cobj = (TObjArray *) fCalibObjects->At(kPHQ);
1427   if(cobj){
1428     //cobj->Print();
1429     storage->Put(cobj, id1, metaData);
1430   }
1431 }
1432
1433  //_________________________________________________________________________________________________________________
1434  void AliTRDPreprocessorOffline::UpdateOCDBChamberStatus(Int_t startRunNumber, Int_t endRunNumber, AliCDBStorage* storage){
1435    //
1436    // Update OCDB entry
1437    //
1438
1439    AliCDBMetaData *metaData= new AliCDBMetaData();
1440    metaData->SetObjectClassName("AliTRDCalChamberStatus");
1441    metaData->SetResponsible("Raphaelle Bailhache and Julian Book");
1442    metaData->AddDateToComment();
1443    metaData->SetBeamPeriod(1);
1444
1445    AliCDBId id1("TRD/Calib/ChamberStatus", startRunNumber, endRunNumber);
1446    AliTRDCalChamberStatus *calChamberStatus = (AliTRDCalChamberStatus *) fCalibObjects->At(kChamberStatus);
1447    if(calChamberStatus) storage->Put(calChamberStatus, id1, metaData);
1448
1449
1450  }
1451  //__________________________________________________________________________________________________________________________
1452  Bool_t AliTRDPreprocessorOffline::ValidateGain() {
1453    //
1454    // Validate OCDB entry
1455    //
1456
1457    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain);
1458    if(calDet) {
1459      Double_t mean = calDet->GetMean();
1460      Double_t rms = calDet->GetRMSRobust();
1461      if((mean > 0.2) && (mean < 1.4) && (rms < 0.5)) return kTRUE;
1462      //if((mean > 0.2) && (mean < 1.4)) return kTRUE;
1463      else {
1464        fStatusPos = fStatusPos | kGainErrorRange;
1465        return kFALSE;
1466      }
1467    }
1468    else return kFALSE;
1469    
1470
1471
1472  }
1473  //__________________________________________________________________________________________________________________________
1474  Bool_t AliTRDPreprocessorOffline::ValidateVdrift(){
1475    //
1476    // Update OCDB entry
1477    //
1478
1479    Int_t detVdrift = kVdriftPHDet;
1480    Bool_t ok = kTRUE;
1481
1482    if(fMethodSecond) detVdrift = kVdriftLinear;
1483
1484    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
1485    if(calDet) {
1486      Double_t mean = calDet->GetMean();
1487      Double_t rms = calDet->GetRMSRobust();
1488      //printf("Vdrift::mean %f, rms %f\n",mean,rms);
1489      if(!((mean > 1.0) && (mean < 2.0) && (rms < 0.5))) {
1490        fStatusPos = fStatusPos | kVdriftErrorRange;
1491        ok = kFALSE;
1492      }
1493    }
1494    else return kFALSE; 
1495
1496    if(!fMethodSecond) {
1497      AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad);
1498      if(calPad) {
1499        Double_t mean = calPad->GetMean();
1500        Double_t rms = calPad->GetRMS();
1501        //printf("Vdrift::meanpad %f, rmspad %f\n",mean,rms);
1502        if(!((mean > 0.9) && (mean < 1.1) && (rms < 0.6))) {
1503          fStatusPos = fStatusPos | kVdriftErrorRange;
1504          ok = kFALSE;
1505        }
1506      }
1507      else return kFALSE;
1508    }
1509
1510    return ok;
1511
1512  }
1513  //__________________________________________________________________________________________________________________________
1514  Bool_t AliTRDPreprocessorOffline::ValidateExB(){
1515    //
1516    // Update OCDB entry
1517    //
1518
1519    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kLorentzLinear);
1520    if(calDet) {
1521      Double_t mean = calDet->GetMean();
1522      Double_t rms = calDet->GetRMSRobust();
1523      //printf("Vdrift::mean %f, rms %f\n",mean,rms);
1524      if(!((mean > -1.0) && (mean < 1.0) && (rms < 0.5))) {
1525        fStatusNeg = fStatusNeg | kExBErrorRange;
1526        return kFALSE;
1527      }
1528      else return kTRUE;
1529    }
1530    else return kFALSE; 
1531    
1532  }
1533  //__________________________________________________________________________________________________________________________
1534  Bool_t AliTRDPreprocessorOffline::ValidateT0(){
1535    //
1536    // Update OCDB entry
1537    //
1538
1539    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet);
1540    AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad);
1541    if(calDet && calPad) {
1542      Double_t meandet = calDet->GetMean();
1543      Double_t rmsdet = calDet->GetRMSRobust();
1544      Double_t meanpad = calPad->GetMean();
1545      //Double_t rmspad = calPad->GetRMS();
1546      printf("T0::meandet %f, rmsdet %f,meanpad %f\n",meandet,rmsdet,meanpad);
1547      if((meandet >   fMinTimeOffsetValidate) && (meandet < 5.0) && (rmsdet < 4.0) && (meanpad < 5.0) && (meanpad > -0.5)) return kTRUE;
1548      else {
1549        fStatusPos = fStatusPos | kTimeOffsetErrorRange;
1550        return kFALSE;
1551      }
1552    }
1553    else return kFALSE;
1554
1555  }
1556  //__________________________________________________________________________________________________________________________
1557  Bool_t AliTRDPreprocessorOffline::ValidatePRF() const{
1558    //
1559    // Update OCDB entry
1560    //
1561
1562    AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF);
1563    if(calPad) {
1564      Double_t meanpad = calPad->GetMean();
1565      Double_t rmspad = calPad->GetRMS();
1566      //printf("PRF::meanpad %f, rmspad %f\n",meanpad,rmspad);
1567      if((meanpad < 1.0) && (rmspad < 0.8)) return kTRUE;
1568      else return kFALSE;
1569    }
1570    else return kFALSE;
1571
1572
1573  }
1574  //__________________________________________________________________________________________________________________________
1575 Bool_t AliTRDPreprocessorOffline::ValidateChamberStatus(){
1576   //
1577   // Update OCDB entry
1578   //
1579   
1580   AliTRDCalChamberStatus *calChamberStatus = (AliTRDCalChamberStatus *) fCalibObjects->At(kChamberStatus);
1581   if(calChamberStatus) {
1582
1583     Int_t detectornodata   = 0;
1584     Int_t detectorbadcalib = 0;
1585     
1586     for(Int_t sm=0; sm < 18; sm++) {
1587       //printf("%d chambers w/o data in sm %d\n",fNoData[sm],sm);
1588       //printf("%d bad calibrated chambers in sm %d\n",fBadCalib[sm],sm);
1589       if(fNoData[sm] != 30) detectornodata += fNoData[sm];
1590       detectorbadcalib+=fBadCalib[sm];
1591     }
1592     //printf("Number of chambers w/o data %d\n",detectornodata);
1593     //printf("Number of chambers bad calibrated %d\n",detectorbadcalib);
1594
1595     if((detectornodata > fNoDataValidate) ||
1596        (detectorbadcalib > fBadCalibValidate)){
1597       fStatusPos = fStatusPos | kChamberStatusErrorRange;
1598       return kFALSE;
1599     }
1600     return kTRUE;
1601   }
1602   else return kFALSE;
1603   
1604 }
1605 //_____________________________________________________________________________
1606 Int_t AliTRDPreprocessorOffline::GetVersion(TString name) const
1607 {
1608   //
1609   // Get version from the title
1610   //
1611   
1612   // Some patterns
1613   const Char_t *version = "Ver";
1614   if(!strstr(name.Data(),version)) return -1;
1615   const Char_t *after = "Subver";  
1616   if(!strstr(name.Data(),after)) return -1;
1617
1618   for(Int_t ver = 0; ver < 999999999; ver++) {
1619
1620     TString vertry(version);
1621     vertry += ver;
1622     vertry += after;
1623
1624     //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
1625
1626     if(strstr(name.Data(),vertry.Data())) return ver;
1627     
1628   }
1629   
1630   return -1;
1631
1632 }
1633
1634 //_____________________________________________________________________________
1635 Int_t AliTRDPreprocessorOffline::GetSubVersion(TString name) const
1636 {
1637   //
1638   // Get subversion from the title
1639   //
1640   
1641   // Some patterns
1642   const Char_t *subversion = "Subver";
1643   if(!strstr(name.Data(),subversion)) return -1;
1644   const Char_t *after = "FirstRun";
1645   if(!strstr(name.Data(),after)) {
1646     after = "Nz";
1647   }
1648   if(!strstr(name.Data(),after)) return -1;
1649
1650   
1651   for(Int_t ver = 0; ver < 999999999; ver++) {
1652     
1653     TString vertry(subversion);
1654     vertry += ver;
1655     vertry += after;
1656
1657     //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
1658
1659     if(strstr(name.Data(),vertry.Data())) return ver;
1660     
1661   }
1662   
1663   return -1;
1664
1665 }
1666
1667 //_____________________________________________________________________________
1668 Int_t AliTRDPreprocessorOffline::GetFirstRun(TString name) const
1669 {
1670   //
1671   // Get first run from the title
1672   //
1673   
1674   // Some patterns
1675   const Char_t *firstrun = "FirstRun";
1676   if(!strstr(name.Data(),firstrun)) return -1;
1677   const Char_t *after = "Nz";  
1678   if(!strstr(name.Data(),after)) return -1;
1679   
1680   
1681   for(Int_t ver = 0; ver < 999999999; ver++) {
1682
1683     TString vertry(firstrun);
1684     vertry += ver;
1685     vertry += after;
1686
1687     //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
1688
1689     if(strstr(name.Data(),vertry.Data())) return ver;
1690     
1691   }
1692   
1693   return -1;
1694
1695 }
1696 //_____________________________________________________________________________
1697 Bool_t AliTRDPreprocessorOffline::CheckStatus(Int_t status, Int_t bitMask) const
1698 {
1699   //
1700   // Checks the status
1701   //
1702
1703   return (status & bitMask) ? kTRUE : kFALSE;
1704   
1705 }
1706 //_____________________________________________________________________________
1707 Int_t AliTRDPreprocessorOffline::GetStatus() const
1708 {
1709   //
1710   // Checks the status
1711   // fStatusPos: errors
1712   // fStatusNeg: only info
1713   //
1714
1715   if(fStatusPos > 0) return fStatusPos;
1716   else return (-TMath::Abs(fStatusNeg));
1717   
1718 }
1719 //_____________________________________________________________________________
1720 void AliTRDPreprocessorOffline::PrintStatus() const
1721 {
1722   //
1723   // Do Summary
1724   //
1725
1726   AliInfo(Form("The error status is %d",fStatusPos));
1727   AliInfo(Form("IsGainErrorOld? %d",(Int_t)IsGainErrorOld()));
1728   AliInfo(Form("IsVdriftErrorOld? %d",(Int_t)IsVdriftErrorOld()));
1729   AliInfo(Form("IsGainErrorRange? %d",(Int_t)IsGainErrorRange()));
1730   AliInfo(Form("IsVdriftErrorRange? %d",(Int_t)IsVdriftErrorRange()));
1731   AliInfo(Form("IsTimeOffsetErrorRange? %d",(Int_t)IsTimeOffsetErrorRange()));
1732   AliInfo(Form("IsChamberStatusErrorRange? %d",(Int_t)IsChamberStatusErrorRange()));
1733
1734  
1735   AliInfo(Form("The info status is %d",fStatusNeg));
1736   AliInfo(Form("IsGainNotEnoughStatsButFill? %d",(Int_t)IsGainNotEnoughStatsButFill()));
1737   AliInfo(Form("IsVdriftNotEnoughStatsButFill? %d",(Int_t)IsVdriftNotEnoughStatsButFill()));
1738   AliInfo(Form("IsGainNotEnoughStatsNotFill? %d",(Int_t)IsGainNotEnoughStatsNotFill()));
1739   AliInfo(Form("IsVdriftNotEnoughStatsNotFill? %d",(Int_t)IsVdriftNotEnoughStatsNotFill()));
1740   AliInfo(Form("IsTimeOffsetNotEnoughStatsNotFill? %d",(Int_t)IsTimeOffsetNotEnoughStatsNotFill()));
1741
1742   AliInfo(Form("IsExBErrorRange? %d",(Int_t)IsExBErrorRange()));
1743   AliInfo(Form("IsExBErrorOld? %d",(Int_t)IsExBErrorOld()));
1744   
1745 }
1746 //___________________________________________________________________________________
1747 void AliTRDPreprocessorOffline::SetCalDetVdrift(AliTRDCalDet *calDetVdriftUsed) 
1748 {
1749
1750   fCalDetVdriftUsed = calDetVdriftUsed;
1751
1752   fCalDetExBUsed = new AliTRDCalDet("lorentz angle tan","lorentz angle tan (detector value)");
1753   for(Int_t k = 0; k < 540; k++){
1754     fCalDetExBUsed->SetValue(k,AliTRDCommonParam::Instance()->GetOmegaTau(fCalDetVdriftUsed->GetValue(k)));
1755     //printf("Set the exb object for detector %d, vdrift %f and exb %f\n",k,fCalDetVdriftUsed->GetValue(k),fCalDetExBUsed->GetValue(k));
1756   }
1757   
1758 };
1759 //___________________________________________________________________________________
1760 Bool_t AliTRDPreprocessorOffline::SetCalDetGain(Int_t runNumber, Int_t version, Int_t subversion) 
1761 {
1762   //
1763   // Set the fCalDetGainUsed
1764   //
1765
1766   if((version == 0) && (subversion == 0)) return kFALSE;
1767
1768   AliCDBEntry *entry = AliCDBManager::Instance()->Get("TRD/Calib/ChamberGainFactor",runNumber, version, subversion);
1769   if(!entry) {
1770     AliError("Found no entry\n");
1771     fStatusPos = fStatusPos | kGainErrorOld;
1772     return kFALSE;
1773   }
1774   //const AliCDBId id = entry->GetId();
1775   //version = id.GetVersion();
1776   //subversion = id.GetSubVersion();
1777   //printf("Found version %d and subversion %d for vdrift\n",version,subversion);
1778   AliTRDCalDet* calDet = (AliTRDCalDet *)entry->GetObject();
1779   if(calDet) fCalDetGainUsed = calDet;
1780   else {
1781     fStatusPos = fStatusPos | kGainErrorOld;
1782     return kFALSE;
1783   }
1784   
1785   return kTRUE;
1786
1787 }
1788 //___________________________________________________________________________________
1789 Bool_t AliTRDPreprocessorOffline::SetCalDetVdriftExB(Int_t runNumber, Int_t versionv, Int_t subversionv, Int_t versionexb, Int_t subversionexb) 
1790 {
1791   //
1792   // Set the fCalDetVdriftUsed and fCalDetExBUsed
1793   //
1794
1795   if((versionv == 0) && (subversionv == 0)) return kFALSE;
1796
1797   AliCDBEntry *entry = AliCDBManager::Instance()->Get("TRD/Calib/ChamberVdrift",runNumber, versionv, subversionv);
1798   if(!entry) {
1799     AliError("Found no entry\n");
1800     fStatusPos = fStatusPos | kVdriftErrorOld;
1801     return kFALSE;
1802   }
1803   AliTRDCalDet* calDet = (AliTRDCalDet *)entry->GetObject();
1804   if(calDet) fCalDetVdriftUsed = calDet;
1805   else {
1806     fStatusPos = fStatusPos | kVdriftErrorOld;
1807     return kFALSE;
1808   }
1809
1810   // ExB object
1811
1812   if((versionexb == 0) && (subversionexb == 0)) {
1813     
1814     fCalDetExBUsed = new AliTRDCalDet("lorentz angle tan","lorentz angle tan (detector value)");
1815     for(Int_t k = 0; k < 540; k++){
1816       fCalDetExBUsed->SetValue(k,AliTRDCommonParam::Instance()->GetOmegaTau(fCalDetVdriftUsed->GetValue(k)));
1817       //printf("Nothing found: set the exb object for detector %d, vdrift %f and exb %f\n",k,fCalDetVdriftUsed->GetValue(k),fCalDetExBUsed->GetValue(k));
1818     }
1819   }
1820   else {
1821
1822     entry = 0x0;
1823     entry = AliCDBManager::Instance()->Get("TRD/Calib/ChamberExB",runNumber, versionexb, subversionexb);
1824     if(!entry) {
1825       //printf("Found no entry\n");
1826       fStatusPos = fStatusPos | kExBErrorOld;   
1827       return kFALSE;
1828     }
1829     AliTRDCalDet* calDetexb = (AliTRDCalDet *)entry->GetObject();
1830     if(!calDetexb) {
1831       fStatusPos = fStatusPos | kExBErrorOld;   
1832       return kFALSE;
1833     }
1834     
1835     Double_t meanexb = calDetexb->GetMean();
1836     //printf("Mean value %f\n",meanexb);
1837     if((meanexb > 70) || (fNoExBUsedInReco)) {
1838       fCalDetExBUsed = new AliTRDCalDet("lorentz angle tan","lorentz angle tan (detector value)");
1839       for(Int_t k = 0; k < 540; k++){
1840         fCalDetExBUsed->SetValue(k,AliTRDCommonParam::Instance()->GetOmegaTau(fCalDetVdriftUsed->GetValue(k)));
1841         //printf("Found but: set the exb object for detector %d, vdrift %f and exb %f\n",k,fCalDetVdriftUsed->GetValue(k),fCalDetExBUsed->GetValue(k));
1842       }
1843     }
1844     else {
1845       fCalDetExBUsed = calDetexb;
1846     }
1847     
1848   }
1849
1850   
1851   return kTRUE;
1852
1853 }
1854