]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDPreprocessorOffline.cxx
Protection against non existing online gain table (Theo)
[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     Double_t gainoverallnotnormalized =  calibra->AnalyseCHAllTogether(fCH2d);
457     if(fCalDetGainUsed && (gainoverallnotnormalized > 0.0)) {
458       AliTRDCalDet *calDetGain = new AliTRDCalDet(*fCalDetGainUsed);
459       Double_t oldmean = fCalDetGainUsed->CalcMean(kFALSE);
460       //printf("oldmean %f\n",oldmean);
461       if(oldmean > 0.0)  {
462         Double_t scalefactor = calibra->GetScaleFactorGain();
463         //printf("Correction factor %f\n",gainoverallnotnormalized*scalefactor);
464         calDetGain->Multiply(gainoverallnotnormalized*scalefactor/oldmean);
465         //printf("newmean %f\n",calDetGain->CalcMean(kFALSE));
466         TH1F *coefGain  = calDetGain->MakeHisto1DAsFunctionOfDet();
467         fCalibObjects->AddAt(calDetGain,kGain);
468         fPlots->AddAt(coefGain,kGain);
469         // 
470         ok = kTRUE;
471         fStatus += 1000000;
472       }
473       else {
474         fStatus += 1000;
475       }      
476     }
477     else {
478       fStatus += 1000;
479     }
480   }
481   
482   calibra->ResetVectorFit();
483   
484   return ok;
485   
486 }
487 //_____________________________________________________________________________________________________
488 Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftT0(){
489   //
490   // Analyze VdriftT0 - produce the calibration objects
491   //
492
493   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
494   calibra->SetMinEntries(fMinStatsVdriftT0PH); // If there is less than 1000 entries in the histo: no fit
495   calibra->AnalysePH(fPH2d);
496
497   Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1))
498     + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1));
499   Int_t nbfit       = calibra->GetNumberFit();
500   Int_t nbfitSuccess = calibra->GetNumberFitSuccess();
501   Int_t nbE         = calibra->GetNumberEnt();
502
503   //printf("nbtg %d, nbfit %d, nbE %d, nbfitSuccess %d\n",nbtg,nbfit,nbE,nbfitSuccess);
504
505   Bool_t ok = kFALSE;
506   if ((nbtg >                  0) && 
507       (nbfit        >= 0.5*nbE) && (nbE > 30) && (nbfitSuccess > 30)) {
508     //printf("Pass the cut for VdriftT0\n");
509     // create the cal objects
510     calibra->RemoveOutliers(1,kFALSE);
511     calibra->PutMeanValueOtherVectorFit(1,kFALSE);
512     calibra->RemoveOutliers2(kFALSE);
513     calibra->PutMeanValueOtherVectorFit2(1,kFALSE);
514     //
515     TObjArray object = calibra->GetVectorFit();
516     AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object);
517     TH1F *coefVdriftPH  = calDetVdrift->MakeHisto1DAsFunctionOfDet();
518     AliTRDCalPad *calPadVdrift = (AliTRDCalPad *)calibra->CreatePadObjectVdrift(&object,calDetVdrift);
519     TH1F *coefPadVdrift   = calPadVdrift->MakeHisto1D();
520     object       = calibra->GetVectorFit2();
521     AliTRDCalDet *calDetT0  = calibra->CreateDetObjectT0(&object);
522     TH1F *coefT0  = calDetT0->MakeHisto1DAsFunctionOfDet();
523     AliTRDCalPad *calPadT0 = (AliTRDCalPad *)calibra->CreatePadObjectT0(&object,calDetT0);
524     TH1F *coefPadT0  = calPadT0->MakeHisto1D();
525     // Put them in the array
526     fCalibObjects->AddAt(calDetT0,kT0PHDet);
527     fCalibObjects->AddAt(calDetVdrift,kVdriftPHDet);
528     fCalibObjects->AddAt(calPadT0,kT0PHPad);
529     fCalibObjects->AddAt(calPadVdrift,kVdriftPHPad);
530     fPlots->AddAt(coefVdriftPH,kVdriftPHDet);
531     fPlots->AddAt(coefT0,kT0PHDet);
532     fPlots->AddAt(coefPadVdrift,kVdriftPHPad);
533     fPlots->AddAt(coefPadT0,kT0PHPad);
534     //
535     ok = kTRUE;
536   }
537   else {
538     fStatus += 100;
539   }
540   calibra->ResetVectorFit();
541  
542   return ok;
543   
544 }
545 //____________________________________________________________________________________________________________________
546 Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftLinearFit(){
547   //
548   // Analyze vdrift linear fit - produce the calibration objects
549   //
550
551   //printf("Analyse linear fit\n");
552
553   
554   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
555   calibra->SetMinEntries(fMinStatsVdriftLinear); // If there is less than 1000 entries in the histo: no fit
556   //printf("Fill PE Array\n");
557   fAliTRDCalibraVdriftLinearFit->FillPEArray();
558   //printf("AliTRDCalibraFit\n");
559   calibra->AnalyseLinearFitters(fAliTRDCalibraVdriftLinearFit);
560   //printf("After\n");
561
562   //Int_t nbtg        = 540;
563   Int_t nbfit       = calibra->GetNumberFit();
564   Int_t nbE         = calibra->GetNumberEnt();
565
566   
567   Bool_t ok = kFALSE;
568   // enough statistics
569   if ((nbfit        >= 0.5*nbE) && (nbE > 30)) {
570     // create the cal objects
571     //calibra->RemoveOutliers(1,kTRUE);
572     calibra->PutMeanValueOtherVectorFit(1,kTRUE);
573     //calibra->RemoveOutliers2(kTRUE);
574     calibra->PutMeanValueOtherVectorFit2(1,kTRUE);
575     //
576     TObjArray object  = calibra->GetVectorFit();
577     AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object,kTRUE);
578     TH1F *coefDriftLinear      = calDetVdrift->MakeHisto1DAsFunctionOfDet();
579     object                     = calibra->GetVectorFit2();
580     AliTRDCalDet *calDetLorentz = calibra->CreateDetObjectLorentzAngle(&object);
581     TH1F *coefLorentzAngle = calDetLorentz->MakeHisto1DAsFunctionOfDet();
582     //if(!calDetLorentz) printf("No lorentz created\n");
583     // Put them in the array
584     fCalibObjects->AddAt(calDetVdrift,kVdriftLinear);
585     fCalibObjects->AddAt(calDetLorentz,kLorentzLinear);
586     fPlots->AddAt(coefDriftLinear,kVdriftLinear);
587     fPlots->AddAt(coefLorentzAngle,kLorentzLinear);
588     //
589     ok = kTRUE;
590   }
591   else {
592     fNotEnoughStatisticsForTheVdriftLinear = kTRUE;
593     Double_t vdriftoverall =  calibra->AnalyseLinearFittersAllTogether(fAliTRDCalibraVdriftLinearFit);
594     if(fCalDetVdriftUsed && (vdriftoverall > 0.0)) {
595       AliTRDCalDet *calDetVdrift = new AliTRDCalDet(*fCalDetVdriftUsed);
596       Double_t oldmean = fCalDetVdriftUsed->CalcMean(kFALSE);
597       //printf("oldmean %f\n",oldmean);
598       if(oldmean > 0.0)  {
599         //printf("Correction factor %f\n",vdriftoverall);
600         calDetVdrift->Multiply(vdriftoverall/oldmean);
601         //printf("newmean %f\n",calDetVdrift->CalcMean(kFALSE));
602         TH1F *coefDriftLinear  = calDetVdrift->MakeHisto1DAsFunctionOfDet();
603         // Put them in the array
604         fCalibObjects->AddAt(calDetVdrift,kVdriftLinear);
605         fPlots->AddAt(coefDriftLinear,kVdriftLinear);
606         // 
607         ok = kTRUE;
608         fStatus += 10000;
609       }
610       else fStatus += 1;      
611     }
612     else fStatus += 1;
613   }
614   
615   calibra->ResetVectorFit();
616   
617   return ok;
618   
619 }
620 //________________________________________________________________________________________________________________
621
622 Bool_t AliTRDPreprocessorOffline::AnalyzePRF(){
623   //
624   // Analyze PRF - produce the calibration objects
625   //
626
627   AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
628   calibra->SetMinEntries(fMinStatsPRF); // If there is less than 1000 entries in the histo: no fit
629   calibra->AnalysePRFMarianFit(fPRF2d);
630
631   Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2))
632     + 6*  18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2));
633   Int_t nbfit       = calibra->GetNumberFit();
634   Int_t nbE         = calibra->GetNumberEnt();
635
636   
637   Bool_t ok = kFALSE;
638   // enough statistics
639   if ((nbtg >                  0) && 
640       (nbfit        >= 0.95*nbE) && (nbE > 30)) {
641     // create the cal objects
642     TObjArray object  = calibra->GetVectorFit();
643     AliTRDCalPad *calPadPRF = (AliTRDCalPad*) calibra->CreatePadObjectPRF(&object);
644     TH1F *coefPRF           = calPadPRF->MakeHisto1D();
645     // Put them in the array
646     fCalibObjects->AddAt(calPadPRF,kPRF);
647     fPlots->AddAt(coefPRF,kPRF);
648     //
649     ok = kTRUE;
650   }
651   
652   calibra->ResetVectorFit();
653   
654   return ok;
655   
656 }
657
658 //_____________________________________________________________________________
659 Bool_t AliTRDPreprocessorOffline::AnalyzeChamberStatus()
660 {
661   //
662   // Produce AliTRDCalChamberStatus out of calibration results
663   //
664   
665   // set up AliTRDCalChamberStatus
666   AliTRDCalChamberStatus *CalChamberStatus = new AliTRDCalChamberStatus();
667   for(Int_t det = 0; det < 540; det++) CalChamberStatus->SetStatus(det,1);
668
669   // get calibration objects
670   AliTRDCalDet *calDetGain   = (AliTRDCalDet *) fCalibObjects->At(kGain);
671   AliTRDCalDet *calDetVDrift = (AliTRDCalDet *) fCalibObjects->At(kVdriftLinear);
672
673   // Check
674   if((!calDetGain) || (!calDetVDrift) || (!fCH2d)) return kFALSE;
675
676   // Gain
677   Double_t gainmean = calDetGain->GetMean();
678   Double_t vdriftmean = calDetVDrift->GetMean();
679
680   Double_t gainrms = calDetGain->GetRMSRobust();
681   Double_t vdriftrms = calDetVDrift->GetRMSRobust();
682
683   //printf("Gain mean: %f, rms: %f\n",gainmean,gainrms);
684   //printf("Vdrift mean: %f, rms: %f\n",vdriftmean,vdriftrms);
685   
686   // Check
687   if((TMath::Abs(gainrms) < 0.001) || (TMath::Abs(vdriftrms) < 0.001)) return kFALSE;
688
689   // mask chambers with empty gain entries
690   //Int_t counter = 0;
691   for (Int_t idet = 0; idet < 540; idet++) {
692
693     // ch2d
694     TH1I *projch =  (TH1I *) fCH2d->ProjectionX("projch",idet+1,idet+1,(Option_t *)"e");
695     Double_t entries = projch->GetEntries();
696
697     // gain
698     Double_t gain = calDetGain->GetValue(idet);
699
700     // vdrift
701     Double_t vdrift = calDetVDrift->GetValue(idet);
702
703
704     if(entries<=0.5 ||
705        TMath::Abs(gainmean-gain) > (15.0*gainrms) ||
706        TMath::Abs(vdriftmean-vdrift) > (15.0*vdriftrms)) {
707      
708       //printf(" chamber det %03d masked \n",idet);
709       //printf(" gainmean %f and gain %f, gainrms %f \n",gainmean,gain,gainrms);
710       //printf(" vdriftmean %f and vdrift %f, vdriftrms %f \n",vdriftmean,vdrift,vdriftrms);
711       CalChamberStatus->SetStatus(idet,AliTRDCalChamberStatus::kMasked);
712       //counter++;
713     }
714
715      /*
716      // installed supermodules+1 -> abort
717      if(counter > (7+1)*30) {
718        printf("ERROR: more than one SM to be masked!! \n Abort...\n");
719        if(projch) delete projch;
720        return 0x0;
721      }
722      */
723
724     delete projch;
725     
726    }
727
728    // Security
729    for(Int_t sm=0; sm < 18; sm++) {
730      Int_t counter = 0;
731      for(Int_t det = 0; det < 30; det++){
732        Int_t detector = sm*30+det;
733        if(CalChamberStatus->IsMasked(detector)) counter++;
734      }
735      if(counter >= 10) {
736        for(Int_t det = 0; det < 30; det++){
737          Int_t detector = sm*30+det;
738          CalChamberStatus->SetStatus(detector,AliTRDCalChamberStatus::kInstalled);
739        }
740      }
741    }
742
743    fCalibObjects->AddAt(CalChamberStatus,kChamberStatus);
744    return kTRUE;
745
746  }
747
748
749  //________________________________________________________________________________________________
750  void AliTRDPreprocessorOffline::CorrectFromDetGainUsed() {
751    //
752    // Correct from the gas gain used afterwards
753    //
754    AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
755    if(!calDetGain) return;
756
757    // Calculate mean
758    Double_t mean = 0.0;
759    Int_t nbdet = 0;
760
761    for(Int_t det = 0; det < 540; det++) {
762
763      Float_t gaininit = fCalDetGainUsed->GetValue(det);
764      Float_t gainout = calDetGain->GetValue(det);
765
766
767      if(TMath::Abs(gainout-1.0) > 0.000001) {
768        mean += (gaininit*gainout);
769        nbdet++;
770      }  
771    }
772    if(nbdet > 0) mean = mean/nbdet;
773
774    for(Int_t det = 0; det < 540; det++) {
775
776      Float_t gaininit = fCalDetGainUsed->GetValue(det);
777      Float_t gainout = calDetGain->GetValue(det);
778
779      if(TMath::Abs(gainout-1.0) > 0.000001) {
780        Double_t newgain = gaininit*gainout;
781        if(newgain < 0.1) newgain = 0.1;
782        if(newgain > 1.9) newgain = 1.9;
783        calDetGain->SetValue(det,newgain);
784      }
785      else {
786        Double_t newgain = mean;
787        if(newgain < 0.1) newgain = 0.1;
788        if(newgain > 1.9) newgain = 1.9;
789        calDetGain->SetValue(det,newgain);
790      }
791    }
792
793
794  }
795  //________________________________________________________________________________________________
796  void AliTRDPreprocessorOffline::CorrectFromDetVdriftUsed() {
797    //
798    // Correct from the drift velocity
799    //
800
801    //printf("Correct for vdrift\n");
802
803    AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
804    if(!calDetGain) return;
805
806    Int_t detVdrift = kVdriftPHDet;
807    if(fMethodSecond) detVdrift = kVdriftLinear;
808    AliTRDCalDet *calDetVdrift = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
809    if(!calDetVdrift) return;
810
811    // Calculate mean
812    if(!fNotEnoughStatisticsForTheGain) {
813      for(Int_t det = 0; det < 540; det++) {
814        
815        Float_t vdriftinit = fCalDetVdriftUsed->GetValue(det);
816        Float_t vdriftout = calDetVdrift->GetValue(det);
817        
818        Float_t gain = calDetGain->GetValue(det);
819        if(vdriftout > 0.0) gain = gain*vdriftinit/vdriftout;
820        if(gain < 0.1) gain = 0.1;
821        if(gain > 1.9) gain = 1.9;
822        calDetGain->SetValue(det,gain);
823      }
824    }
825    else {
826      
827      Float_t vdriftinit = fCalDetVdriftUsed->CalcMean(kFALSE);
828      Float_t vdriftout = calDetVdrift->CalcMean(kFALSE);
829      Float_t factorcorrectif = 1.0;
830      if(vdriftout > 0.0) factorcorrectif = vdriftinit/vdriftout;
831      for(Int_t det = 0; det < 540; det++) {
832        Float_t gain = calDetGain->GetValue(det);
833        gain = gain*factorcorrectif;
834        if(gain < 0.1) gain = 0.1;
835        if(gain > 1.9) gain = 1.9;
836        calDetGain->SetValue(det,gain);
837      }
838      
839    }
840    
841  }
842 //_________________________________________________________________________________________________________________
843  void AliTRDPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
844    //
845    // Update OCDB entry
846    //
847
848    AliCDBMetaData *metaData= new AliCDBMetaData();
849    metaData->SetObjectClassName("AliTRDCalDet");
850    metaData->SetResponsible("Raphaelle Bailhache");
851    metaData->SetBeamPeriod(1);
852
853    AliCDBId id1("TRD/Calib/ChamberGainFactor", startRunNumber, endRunNumber);
854    AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
855    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain);
856    if(calDet) gStorage->Put(calDet, id1, metaData);
857
858
859  }
860  //___________________________________________________________________________________________________________________
861  void AliTRDPreprocessorOffline::UpdateOCDBExB(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
862    //
863    // Update OCDB entry
864    //
865
866    Int_t detExB = kLorentzLinear;
867    if(!fMethodSecond) return;
868
869    //printf("Pass\n");
870
871    AliCDBMetaData *metaData= new AliCDBMetaData();
872    metaData->SetObjectClassName("AliTRDCalDet");
873    metaData->SetResponsible("Raphaelle Bailhache");
874    metaData->SetBeamPeriod(1);
875
876    AliCDBId id1("TRD/Calib/ChamberExB", startRunNumber, endRunNumber);
877    AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
878    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detExB);
879    if(calDet) gStorage->Put(calDet, id1, metaData);
880    //if(!calDet) printf("No caldet\n");
881
882  }
883  //___________________________________________________________________________________________________________________
884  void AliTRDPreprocessorOffline::UpdateOCDBVdrift(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
885    //
886    // Update OCDB entry
887    //
888
889    Int_t detVdrift = kVdriftPHDet;
890
891    if(fMethodSecond) detVdrift = kVdriftLinear;
892
893    AliCDBMetaData *metaData= new AliCDBMetaData();
894    metaData->SetObjectClassName("AliTRDCalDet");
895    metaData->SetResponsible("Raphaelle Bailhache");
896    metaData->SetBeamPeriod(1);
897
898    AliCDBId id1("TRD/Calib/ChamberVdrift", startRunNumber, endRunNumber);
899    AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
900    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
901    if(calDet) gStorage->Put(calDet, id1, metaData);
902
903    //
904
905    if(!fMethodSecond) {
906
907      AliCDBMetaData *metaDataPad= new AliCDBMetaData();
908      metaDataPad->SetObjectClassName("AliTRDCalPad");
909      metaDataPad->SetResponsible("Raphaelle Bailhache");
910      metaDataPad->SetBeamPeriod(1);
911
912      AliCDBId id1Pad("TRD/Calib/LocalVdrift", startRunNumber, endRunNumber);
913      AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad);
914      if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad);
915
916    }
917
918  }
919  //________________________________________________________________________________________________________________________
920  void AliTRDPreprocessorOffline::UpdateOCDBT0(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
921    //
922    // Update OCDB entry
923    //
924
925    AliCDBMetaData *metaData= new AliCDBMetaData();
926    metaData->SetObjectClassName("AliTRDCalDet");
927    metaData->SetResponsible("Raphaelle Bailhache");
928    metaData->SetBeamPeriod(1);
929
930    AliCDBId id1("TRD/Calib/ChamberT0", startRunNumber, endRunNumber);
931    AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
932    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet);
933    if(calDet) gStorage->Put(calDet, id1, metaData);
934
935    //
936
937    AliCDBMetaData *metaDataPad= new AliCDBMetaData();
938    metaDataPad->SetObjectClassName("AliTRDCalPad");
939    metaDataPad->SetResponsible("Raphaelle Bailhache");
940    metaDataPad->SetBeamPeriod(1);
941
942    AliCDBId id1Pad("TRD/Calib/LocalT0", startRunNumber, endRunNumber);
943    AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad);
944    if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad);
945
946
947
948  }
949  //_________________________________________________________________________________________________________________
950  void AliTRDPreprocessorOffline::UpdateOCDBPRF(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
951    //
952    // Update OCDB entry
953    //
954
955    AliCDBMetaData *metaData= new AliCDBMetaData();
956    metaData->SetObjectClassName("AliTRDCalPad");
957    metaData->SetResponsible("Raphaelle Bailhache");
958    metaData->SetBeamPeriod(1);
959
960    AliCDBId id1("TRD/Calib/PRFWidth", startRunNumber, endRunNumber);
961    AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
962    AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF);
963    if(calPad) gStorage->Put(calPad, id1, metaData);
964
965
966  }
967  //_________________________________________________________________________________________________________________
968  void AliTRDPreprocessorOffline::UpdateOCDBChamberStatus(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
969    //
970    // Update OCDB entry
971    //
972
973    AliCDBMetaData *metaData= new AliCDBMetaData();
974    metaData->SetObjectClassName("AliTRDCalChamberStatus");
975    metaData->SetResponsible("Raphaelle Bailhache");
976    metaData->SetBeamPeriod(1);
977
978    AliCDBId id1("TRD/Calib/ChamberStatus", startRunNumber, endRunNumber);
979    AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
980    AliTRDCalChamberStatus *calChamberStatus = (AliTRDCalChamberStatus *) fCalibObjects->At(kChamberStatus);
981    if(calChamberStatus) gStorage->Put(calChamberStatus, id1, metaData);
982
983
984  }
985  //__________________________________________________________________________________________________________________________
986  Bool_t AliTRDPreprocessorOffline::ValidateGain() {
987    //
988    // Validate OCDB entry
989    //
990
991    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain);
992    if(calDet) {
993      Double_t mean = calDet->GetMean();
994      Double_t rms = calDet->GetRMSRobust();
995      if((mean > 0.2) && (mean < 1.4) && (rms < 0.5)) return kTRUE;
996      //if((mean > 0.2) && (mean < 1.4)) return kTRUE;
997      else {
998        fStatus += 1000000000;
999        return kFALSE;
1000      }
1001    }
1002    else return kFALSE;
1003    
1004
1005
1006  }
1007  //__________________________________________________________________________________________________________________________
1008  Bool_t AliTRDPreprocessorOffline::ValidateVdrift(){
1009    //
1010    // Update OCDB entry
1011    //
1012
1013    Int_t detVdrift = kVdriftPHDet;
1014    Bool_t ok = kTRUE;
1015
1016    if(fMethodSecond) detVdrift = kVdriftLinear;
1017
1018    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
1019    if(calDet) {
1020      Double_t mean = calDet->GetMean();
1021      Double_t rms = calDet->GetRMSRobust();
1022      //printf("Vdrift::mean %f, rms %f\n",mean,rms);
1023      if(!((mean > 1.0) && (mean < 2.0) && (rms < 0.5))) {
1024        fStatus += 10000000;
1025        ok = kFALSE;
1026      }
1027    }
1028    else return kFALSE; 
1029
1030    if(!fMethodSecond) {
1031      AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad);
1032      if(calPad) {
1033        Double_t mean = calPad->GetMean();
1034        Double_t rms = calPad->GetRMS();
1035        //printf("Vdrift::meanpad %f, rmspad %f\n",mean,rms);
1036        if(!((mean > 0.9) && (mean < 1.1) && (rms < 0.6))) {
1037          fStatus += 10000000;
1038          ok = kFALSE;
1039        }
1040      }
1041      else return kFALSE;
1042    }
1043
1044    return ok;
1045
1046  }
1047  //__________________________________________________________________________________________________________________________
1048  Bool_t AliTRDPreprocessorOffline::ValidateT0(){
1049    //
1050    // Update OCDB entry
1051    //
1052
1053    AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet);
1054    AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad);
1055    if(calDet && calPad) {
1056      Double_t meandet = calDet->GetMean();
1057      Double_t rmsdet = calDet->GetRMSRobust();
1058      Double_t meanpad = calPad->GetMean();
1059      //Double_t rmspad = calPad->GetRMS();
1060      //printf("T0::minimum %f, rmsdet %f,meanpad %f, rmspad %f\n",meandet,rmsdet,meanpad,rmspad);
1061      if((meandet > -1.5) && (meandet < 5.0) && (rmsdet < 4.0) && (meanpad < 5.0) && (meanpad > -0.5)) return kTRUE;
1062      else {
1063        fStatus += 100000000;
1064        return kFALSE;
1065      }
1066    }
1067    else return kFALSE;
1068
1069  }
1070  //__________________________________________________________________________________________________________________________
1071  Bool_t AliTRDPreprocessorOffline::ValidatePRF() const{
1072    //
1073    // Update OCDB entry
1074    //
1075
1076    AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF);
1077    if(calPad) {
1078      Double_t meanpad = calPad->GetMean();
1079      Double_t rmspad = calPad->GetRMS();
1080      //printf("PRF::meanpad %f, rmspad %f\n",meanpad,rmspad);
1081      if((meanpad < 1.0) && (rmspad < 0.8)) return kTRUE;
1082      else return kFALSE;
1083    }
1084    else return kFALSE;
1085
1086
1087  }
1088  //__________________________________________________________________________________________________________________________
1089 Bool_t AliTRDPreprocessorOffline::ValidateChamberStatus() const{
1090   //
1091   // Update OCDB entry
1092   //
1093   
1094   AliTRDCalChamberStatus *calChamberStatus = (AliTRDCalChamberStatus *) fCalibObjects->At(kChamberStatus);
1095   if(calChamberStatus) {
1096     Int_t detectormasked = 0;
1097     for(Int_t det = 0; det < 540; det++) {
1098       if(calChamberStatus->IsMasked(det)) detectormasked++;
1099     }
1100     //printf("Number of chambers masked %d\n",detectormasked);
1101     if(detectormasked > 40) return kFALSE;
1102     else return kTRUE;
1103   }
1104   else return kFALSE;
1105  
1106 }
1107 //_____________________________________________________________________________
1108 Int_t AliTRDPreprocessorOffline::GetVersion(TString name) const
1109 {
1110   //
1111   // Get version from the title
1112   //
1113   
1114   // Some patterns
1115   const Char_t *version = "Ver";
1116   if(!strstr(name.Data(),version)) return -1;
1117   const Char_t *after = "Subver";  
1118   if(!strstr(name.Data(),after)) return -1;
1119
1120   for(Int_t ver = 0; ver < 999999999; ver++) {
1121
1122     TString vertry(version);
1123     vertry += ver;
1124     vertry += after;
1125
1126     //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
1127
1128     if(strstr(name.Data(),vertry.Data())) return ver;
1129     
1130   }
1131   
1132   return -1;
1133
1134 }
1135
1136 //_____________________________________________________________________________
1137 Int_t AliTRDPreprocessorOffline::GetSubVersion(TString name) const
1138 {
1139   //
1140   // Get subversion from the title
1141   //
1142   
1143   // Some patterns
1144   const Char_t *subversion = "Subver";
1145   if(!strstr(name.Data(),subversion)) return -1;
1146   const Char_t *after = "FirstRun";
1147   if(!strstr(name.Data(),after)) {
1148     after = "Nz";
1149   }
1150   if(!strstr(name.Data(),after)) return -1;
1151
1152   
1153   for(Int_t ver = 0; ver < 999999999; ver++) {
1154     
1155     TString vertry(subversion);
1156     vertry += ver;
1157     vertry += after;
1158
1159     //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
1160
1161     if(strstr(name.Data(),vertry.Data())) return ver;
1162     
1163   }
1164   
1165   return -1;
1166
1167 }
1168
1169 //_____________________________________________________________________________
1170 Int_t AliTRDPreprocessorOffline::GetFirstRun(TString name) const
1171 {
1172   //
1173   // Get first run from the title
1174   //
1175   
1176   // Some patterns
1177   const Char_t *firstrun = "FirstRun";
1178   if(!strstr(name.Data(),firstrun)) return -1;
1179   const Char_t *after = "Nz";  
1180   if(!strstr(name.Data(),after)) return -1;
1181   
1182   
1183   for(Int_t ver = 0; ver < 999999999; ver++) {
1184
1185     TString vertry(firstrun);
1186     vertry += ver;
1187     vertry += after;
1188
1189     //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
1190
1191     if(strstr(name.Data(),vertry.Data())) return ver;
1192     
1193   }
1194   
1195   return -1;
1196
1197 }
1198