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