]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/TenderSupplies/AliTOFTenderSupply.cxx
4b5101f73408b36e295a6acff787143a21de8b0d
[u/mrichter/AliRoot.git] / ANALYSIS / TenderSupplies / AliTOFTenderSupply.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 // TOF tender: - load updated calibrations (for TOF and T0)
20 //             - set tofHeader if missing in ESDs (2010 data)
21 //             - re-apply PID 
22 //             
23 // Contacts: Pietro.Antonioli@bo.infn.it                                     //
24 //           Francesco.Noferini@bo.infn.it                                   //
25 ///////////////////////////////////////////////////////////////////////////////
26 #include <TFile.h>
27 #include <TChain.h>
28
29 #include <TMath.h>
30 #include <TRandom.h>
31 #include <AliLog.h>
32 #include <AliESDEvent.h>
33 #include <AliESDtrack.h>
34 #include <AliESDInputHandler.h>
35 #include <AliAnalysisManager.h>
36 #include <AliESDpid.h>
37 #include <AliTender.h>
38
39 #include <AliTOFcalib.h>
40 #include <AliTOFT0maker.h>
41
42 #include <AliGeomManager.h>
43 #include <AliCDBManager.h>
44 #include <AliCDBEntry.h>
45
46 #include <AliOADBContainer.h>
47 #include <AliTOFPIDParams.h>
48
49 #include <AliT0CalibSeasonTimeShift.h>
50
51 #include "AliTOFTenderSupply.h"
52
53 ClassImp(AliTOFTenderSupply)
54
55 Float_t AliTOFTenderSupply::fgT0Aresolution = 75.;
56 Float_t AliTOFTenderSupply::fgT0Cresolution = 65.;
57
58 AliTOFTenderSupply::AliTOFTenderSupply() :
59   AliTenderSupply(),
60   fESDpid(0x0),
61   fTenderNoAction(kTRUE),
62   fIsMC(kFALSE),
63   fCorrectExpTimes(kTRUE),
64   fCorrectTRDBug(kFALSE),
65   fLHC10dPatch(kFALSE),
66   fT0DetectorAdjust(kFALSE),
67   fDebugLevel(0),
68   fAutomaticSettings(kTRUE),
69   fRecoPass(0),
70   fUserRecoPass(0),
71   fForceCorrectTRDBug(kFALSE),
72   fTOFPIDParams(0x0),
73   fTOFCalib(0x0),
74   fTOFT0maker(0x0),
75   fT0IntercalibrationShift(0),
76   fGeomSet(kFALSE),
77   fIsEnteringInTRD(kFALSE),
78   fInTRD(kFALSE),
79   fIsComingOutTRD(kFALSE),
80   fOutTRD(kFALSE),
81   fRhoTRDin(288.38), // cm
82   fRhoTRDout(366.38), // cm
83   fStep(0.5),
84   fMagField(0.),
85   fCDBkey(0)
86
87
88
89 {
90   //
91   // default ctor
92   //
93   fT0shift[0] = 0;
94   fT0shift[1] = 0;
95   fT0shift[2] = 0;
96   fT0shift[3] = 0;
97 }
98
99 //_____________________________________________________
100 AliTOFTenderSupply::AliTOFTenderSupply(const char *name, const AliTender *tender) :
101   AliTenderSupply(name,tender),
102   fESDpid(0x0),
103   fTenderNoAction(kTRUE),
104   fIsMC(kFALSE),
105   fCorrectExpTimes(kTRUE),
106   fCorrectTRDBug(kFALSE),
107   fLHC10dPatch(kFALSE),
108   fT0DetectorAdjust(kFALSE),
109   fDebugLevel(0),
110   fAutomaticSettings(kTRUE),
111   fRecoPass(0),
112   fUserRecoPass(0),
113   fForceCorrectTRDBug(kFALSE),
114   fTOFPIDParams(0x0),
115   fTOFCalib(0x0),
116   fTOFT0maker(0x0),
117   fT0IntercalibrationShift(0),
118   fGeomSet(kFALSE),
119   fIsEnteringInTRD(kFALSE),
120   fInTRD(kFALSE),
121   fIsComingOutTRD(kFALSE),
122   fOutTRD(kFALSE),
123   fRhoTRDin(288.38), // cm
124   fRhoTRDout(366.38), // cm
125   fStep(0.5),
126   fMagField(0.),
127   fCDBkey(0)
128  
129 {
130   //
131   // named ctor
132   //
133
134   fT0shift[0] = 0;
135   fT0shift[1] = 0;
136   fT0shift[2] = 0;
137   fT0shift[3] = 0;
138 }
139
140 //_____________________________________________________
141 void AliTOFTenderSupply::Init()
142 {
143
144   fTenderNoAction = kFALSE;
145   // Initialise TOF tender (this is called at each detected run change)
146   AliLog::SetClassDebugLevel("AliTOFTenderSupply",10); 
147
148   // Setup PID object, check for MC, set AliTOFcalib and TOFT0 maker conf
149   Int_t run = fTender->GetRun();
150   if (run == 0) return;                // to skip first init, when we don't have yet a run number
151
152   fGeomSet=kFALSE;                         
153   // Even if the user didn't set fIsMC, we force it on if we find the MC handler  
154   AliAnalysisManager *mgr=AliAnalysisManager::GetAnalysisManager();
155   if (mgr->GetMCtruthEventHandler() && !(fIsMC) ) {
156     AliWarning("This ESD is MC, fIsMC found OFF: fIsMC turned ON");
157     fIsMC=kTRUE;
158   }
159
160   if (fAutomaticSettings) {
161     if (!fIsMC) {
162       if (fUserRecoPass == 0) DetectRecoPass();
163       else fRecoPass = fUserRecoPass;
164     }
165     if (run<114737) {
166       fTenderNoAction = kTRUE;
167     }
168     else if (run>=114737&&run<=117223) {      //period="LHC10B";
169       if (fRecoPass == 2) {fCorrectExpTimes=kTRUE; fCorrectTRDBug=kFALSE;}
170       else if (fRecoPass == 3) {fCorrectExpTimes=kFALSE; fCorrectTRDBug=kTRUE;}
171       fLHC10dPatch=kFALSE;
172       fT0IntercalibrationShift = 0;
173       fT0DetectorAdjust=kTRUE;
174     }
175     else if (run>=118503&&run<=121040) { //period="LHC10C";
176       if (fRecoPass == 2) {fCorrectExpTimes=kTRUE; fCorrectTRDBug=kFALSE;}
177       else if (fRecoPass == 3) {fCorrectExpTimes=kFALSE; fCorrectTRDBug=kTRUE;}
178       fLHC10dPatch=kFALSE;
179       fT0IntercalibrationShift = 0;
180       fT0DetectorAdjust=kFALSE;
181     }
182     else if (run>=122195&&run<=126437) { //period="LHC10D";
183       fCorrectExpTimes=kFALSE;
184       fLHC10dPatch=kTRUE;
185       fT0IntercalibrationShift = 0;
186       fT0DetectorAdjust=kTRUE;
187     }
188     else if (run>=127719&&run<=130850) { //period="LHC10E";
189       fCorrectExpTimes=kFALSE;
190       fLHC10dPatch=kFALSE;
191       fT0IntercalibrationShift = 30.;
192       fT0DetectorAdjust=kTRUE;
193     }
194     else if (run>=133004&&run<=135029) { //period="LHC10F";
195       fTenderNoAction=kTRUE;
196     }
197     else if (run>=135654&&run<=136377) { //period="LHC10G";
198       fTenderNoAction=kTRUE;
199     }
200     else if (run>=136851&&run<=139517) { //period="LHC10H" - pass2;
201       fCorrectExpTimes=kFALSE;
202       fLHC10dPatch=kFALSE;                
203       fT0IntercalibrationShift = 0.;
204       fT0DetectorAdjust=kTRUE;
205     }
206     else if (run>=139699) {              //period="LHC11A";
207       fTenderNoAction=kTRUE;
208     }
209   }
210
211   if (fTenderNoAction) {
212     AliInfo(" |---------------------------------------------------------------------------|");
213     AliInfo(" |                                                                           |");
214     AliInfo(Form(" |  TOF tender is not supported for run %d                           |",run));
215     AliInfo(" |  TOF tender will do nothing.                                              |");
216     AliInfo(" |  Check TOF tender usage for run/periods at:                               |");
217     AliInfo(" |  https://twiki.cern.ch/twiki/bin/view/ALICE/TOF.                          |");
218     AliInfo(" |---------------------------------------------------------------------------|");
219     AliInfo(" ");
220     return;
221   }
222
223
224   // Load from OADB TOF resolution
225   LoadTOFPIDParams(run);
226
227   // Check if another tender wagon already created the esd pid object
228   // if not we create it and set it to the ESD input handler
229   fESDpid=fTender->GetESDhandler()->GetESDpid();
230   if (!fESDpid) {
231     fESDpid=new AliESDpid;
232     fTender->GetESDhandler()->SetESDpid(fESDpid);
233   }
234
235
236   // Configure TOF calibration class
237   if (!fTOFCalib)fTOFCalib=new AliTOFcalib();  // create if needed
238   fTOFCalib->SetRemoveMeanT0(!(fIsMC));        // must be kFALSE on MC (default is kTRUE)
239   fTOFCalib->SetCalibrateTOFsignal(!(fIsMC));  // must be kFALSE on MC (no new calibration) (default is kTRUE)
240   fTOFCalib->SetCorrectTExp(fCorrectExpTimes); // apply a fine tuning on the expected times at low momenta
241                                                // (this is done for LHC10b/LHC10c pass2)
242
243   // Configure TOFT0 maker class
244   //  if (!fTOFT0maker) fTOFT0maker = new AliTOFT0maker(fESDpid,fTOFCalib); // create if needed
245   if (!fTOFT0maker) fTOFT0maker = new AliTOFT0maker(fESDpid); // without passing AliTOFCalib it uses the diamond
246   fTOFT0maker->SetTimeResolution(fTOFPIDParams->GetTOFresolution());     // set TOF resolution for the PID
247   fTOFT0maker->SetTOFT0algorithm(2);
248
249   AliInfo("|******************************************************|");
250   AliInfo(Form("|    Alice TOF Tender Initialisation (Run %d)  |",fTender->GetRun()));
251   AliInfo("|    Settings:                                         |");
252   AliInfo(Form("|    Correct Exp Times              :  %d               |",fCorrectExpTimes));
253   AliInfo(Form("|    Correct TRD Bug                :  %d               |",fCorrectTRDBug));
254   AliInfo(Form("|    LHC10d patch                   :  %d               |",fLHC10dPatch));
255   AliInfo(Form("|    TOF resolution for TOFT0 maker :  %5.2f (ps)     |",fTOFPIDParams->GetTOFresolution()));
256   AliInfo(Form("|    MC flag                        :  %d               |",fIsMC));
257   AliInfo(Form("|    T0 detector offsets applied    :  %d               |",fT0DetectorAdjust));
258   AliInfo(Form("|    TOF/T0 intecalibration shift   :  %5.2f (ps)     |",fT0IntercalibrationShift));
259   AliInfo("|******************************************************|");
260
261
262 }
263
264 //_____________________________________________________
265 void AliTOFTenderSupply::ProcessEvent()
266 {
267   //
268   // Use updated calibrations for TOF and T0, reapply PID information
269   // For MC: timeZero sampling and additional smearing for T0
270
271   if (fDebugLevel > 1) AliInfo("process event");
272
273   AliESDEvent *event=fTender->GetEvent();
274   if (!event) return;
275   if (fDebugLevel > 1) AliInfo("event read");
276
277
278     
279   if (fTender->RunChanged()){ 
280
281     Init();
282     if (fTenderNoAction) return;            
283     Int_t versionNumber = GetOCDBVersion(fTender->GetRun());
284     fTOFCalib->SetRunParamsSpecificVersion(versionNumber);
285     fTOFCalib->Init(fTender->GetRun());
286     
287     if(event->GetT0TOF()){ // read T0 detector correction from OCDB
288       // OCDB instance
289       if (fT0DetectorAdjust) {
290         AliCDBManager* ocdbMan = AliCDBManager::Instance();
291         ocdbMan->SetRun(fTender->GetRun());    
292         AliCDBEntry *entry = ocdbMan->Get("T0/Calib/TimeAdjust/");
293         if(entry) {
294           AliT0CalibSeasonTimeShift *clb = (AliT0CalibSeasonTimeShift*) entry->GetObject();
295           Float_t *t0means= clb->GetT0Means();
296           //      Float_t *t0sigmas = clb->GetT0Sigmas();
297           fT0shift[0] = t0means[0] + fT0IntercalibrationShift;
298           fT0shift[1] = t0means[1] + fT0IntercalibrationShift;
299           fT0shift[2] = t0means[2] + fT0IntercalibrationShift;
300           fT0shift[3] = t0means[3] + fT0IntercalibrationShift;
301         } else {
302           for (Int_t i=0;i<4;i++) fT0shift[i]=0;
303           AliWarning("TofTender no T0 entry found T0shift set to 0");
304         }
305       } else {
306         for (Int_t i=0;i<4;i++) fT0shift[i]=0;
307       }
308     }
309   }
310
311   if (fTenderNoAction) return;
312
313   fTOFCalib->CalibrateESD(event);   //recalculate TOF signal (no harm for MC, see settings inside init)
314
315
316   // patches for various reconstruction bugs
317   if (fLHC10dPatch && !(fIsMC)) RecomputeTExp(event);   // LHC10d pass2: fake full TRD geometry
318   if ( (fCorrectTRDBug && !(fIsMC)) || (fForceCorrectTRDBug)) FixTRDBug(event);     // LHC10b,c pass3: wrong TRD dE/dx 
319
320   Double_t startTime = 0.;
321   if (fIsMC) startTime = fTOFCalib->TuneForMC(event,fTOFPIDParams->GetTOFresolution());   // this is for old MC when we didn't jitter startTime in MC
322
323   if (fDebugLevel > 1) Printf(" TofTender: startTime %f",startTime);
324   if (fDebugLevel > 1) Printf(" TofTender: T0 time (orig) %f %f %f",event->GetT0TOF(0),event->GetT0TOF(1),event->GetT0TOF(2));
325
326   // event by event TO detector treatment  
327   if(event->GetT0TOF()){   // protection: we adjust T0 only if it is there....
328
329     if (event->GetT0TOF(0) == 0) event->SetT0TOF(0, 9999999.); // in case no information we set to unknown
330     if (event->GetT0TOF(1) == 0) event->SetT0TOF(1, 99999.);
331     if (event->GetT0TOF(2) == 0) event->SetT0TOF(2, 99999.);
332
333     if(!fIsMC){   // data: apply shifts to align around Zero
334       event->SetT0TOF(0,event->GetT0TOF(0) - fT0shift[0]);
335       event->SetT0TOF(1,event->GetT0TOF(1) - fT0shift[1]);
336       event->SetT0TOF(2,event->GetT0TOF(2) - fT0shift[2]);
337     } else {
338       // MC: add smearing for realistic T0A and T0C resolution
339       Double_t defResolutionT0A = 33.;   // in future we will get this from ESDrun data structure or via OCDB
340       Double_t defResolutionT0C = 30.;   // for the moment we don't trust them
341       if ( (fgT0Aresolution > defResolutionT0A) && (event->GetT0TOF(1)<90000.) ) { // add smearing only if signal is there
342         Double_t addedSmearingT0A = TMath::Sqrt(fgT0Aresolution*fgT0Aresolution - defResolutionT0A*defResolutionT0A);
343         Double_t smearingT0A = gRandom->Gaus(0.,addedSmearingT0A);
344         event->SetT0TOF(1,event->GetT0TOF(1) + smearingT0A);
345       }
346       if ( (fgT0Cresolution > defResolutionT0C) && (event->GetT0TOF(2)<90000.) ) { // add smearing only if signal is there
347         Double_t addedSmearingT0C = TMath::Sqrt(fgT0Cresolution*fgT0Cresolution - defResolutionT0C*defResolutionT0C);
348         Double_t smearingT0C = gRandom->Gaus(0.,addedSmearingT0C);
349         event->SetT0TOF(2,event->GetT0TOF(2) + smearingT0C);
350       }
351       if (event->GetT0TOF(0)<90000.) { // we recompute the AND only if it is already there...
352         Double_t smearedT0AC = (event->GetT0TOF(1)+event->GetT0TOF(2))/2.;
353         event->SetT0TOF(0,smearedT0AC); 
354       }
355       if (fDebugLevel > 1) Printf(" TofTender: T0 time (postSmear) %f %f %f",event->GetT0TOF(0),event->GetT0TOF(1),event->GetT0TOF(2));
356       // add finally the timeZero offset also to the T0 detector information
357       event->SetT0TOF(0,event->GetT0TOF(0) + startTime);
358       event->SetT0TOF(1,event->GetT0TOF(1) + startTime);
359       event->SetT0TOF(2,event->GetT0TOF(2) + startTime);  
360       if (fDebugLevel > 1) Printf(" TofTender: T0 time (postStart) %f %f %f",event->GetT0TOF(0),event->GetT0TOF(1),event->GetT0TOF(2));
361     }
362     // after shifts adjust (data) or smearing+offset (MC) we 'clean' to default if signals not there 
363     if(event->GetT0TOF(0) > 900000) event->SetT0TOF(0, 999999.);
364     if(event->GetT0TOF(1) > 90000)  event->SetT0TOF(1, 99999.);
365     if(event->GetT0TOF(2) > 90000)  event->SetT0TOF(2, 99999.);
366   }
367   if (fDebugLevel > 1) Printf(" TofTender: T0 time (FINAL) %f %f %f",event->GetT0TOF(0),event->GetT0TOF(1),event->GetT0TOF(2));
368   
369   //compute timeZero of the event via TOF-TO
370   fTOFT0maker->ComputeT0TOF(event);
371   fTOFT0maker->WriteInESD(event);
372
373   //  set preferred startTime: this is now done via AliPIDResponseTask
374   fESDpid->SetTOFResponse(event, (AliESDpid::EStartTimeType_t)fTOFPIDParams->GetStartTimeMethod());
375
376   // recalculate PID probabilities
377   // this is for safety, especially if the user doesn't attach a PID tender after TOF tender  
378   Int_t ntracks=event->GetNumberOfTracks();
379   //  AliESDtrack *track = NULL;
380   //  Float_t tzeroTrack = 0;
381   for(Int_t itrack = 0; itrack < ntracks; itrack++){
382     //    track = event->GetTrack(itrack);
383     //    tzeroTrack = fESDpid->GetTOFResponse().GetStartTime(track->P());
384     //    Printf("================> Track # %d mom: %f tzeroTrack %f",itrack,track->P(),tzeroTrack);
385     fESDpid->MakeTOFPID(event->GetTrack(itrack),0);   
386   }
387   
388   
389 }
390
391
392 //_____________________________________________________
393 void AliTOFTenderSupply::RecomputeTExp(AliESDEvent *event) const
394 {
395   /*
396    * calibrate TExp
397    */
398
399   
400   /* loop over tracks */
401   AliESDtrack *track = NULL;
402   for (Int_t itrk = 0; itrk < event->GetNumberOfTracks(); itrk++) {
403     /* get track and calibrate */
404     track = event->GetTrack(itrk);
405     RecomputeTExp(track);
406   }
407   
408 }
409
410 //_____________________________________________________
411 void AliTOFTenderSupply::RecomputeTExp(AliESDtrack *track) const
412 {
413   /*** 
414        THIS METHOD IS BASED ON THEORETICAL EXPECTED TIME COMPUTED
415        USING AVERAGE MOMENTUM BETWEEN INNER/OUTER TRACK PARAMS 
416        IT IS A ROUGH APPROXIMATION APPLIED TO FIX LHC10d-pass2 DATA
417        WHERE A WRONG GEOMETRY (FULL TRD) WAS INSERTED
418   ***/
419
420   Double_t texp[AliPID::kSPECIES];
421   if (!track || !(track->GetStatus() & AliESDtrack::kTOFout)) return;
422
423
424   /* get track params */
425   Float_t l = track->GetIntegratedLength();
426   Float_t p = track->P();
427   if (track->GetInnerParam() && track->GetOuterParam()) {
428     Float_t pin = track->GetInnerParam()->P();
429     Float_t pout = track->GetOuterParam()->P();
430     p = 0.5 * (pin + pout);
431   }
432   /* loop over particle types and compute expected time */
433   for (Int_t ipart = 0; ipart < AliPID::kSPECIES; ipart++)
434     texp[ipart] = GetExpTimeTh(AliPID::ParticleMass(ipart), p, l) - 37.; 
435   // 37 is a final semiempirical offset to further adjust (calibrations were
436   // done with "standard" integratedTimes)
437   /* set integrated times */
438   track->SetIntegratedTimes(texp);
439
440 }
441
442
443 //______________________________________________________________________________
444 void AliTOFTenderSupply::DetectRecoPass()
445 {
446   //
447   // Detect reconstruction information
448   //
449   
450   //reset information
451   fRecoPass=0;
452   
453   //Get the current file to check the reconstruction pass (UGLY, but not stored in ESD... )
454   AliAnalysisManager *mgr=AliAnalysisManager::GetAnalysisManager();
455   AliVEventHandler *inputHandler=mgr->GetInputEventHandler();
456   if (!inputHandler) return;
457   
458   TTree *tree= (TTree*)inputHandler->GetTree();
459   TFile *file= (TFile*)tree->GetCurrentFile();
460   
461   if (!file) {
462     AliFatal("Current file not found");
463     return; // coverity
464   }
465   
466   //find pass from file name (UGLY, but not stored in ESD... )
467   TString fileName(file->GetName());
468   if (fileName.Contains("pass1") ) {
469     fRecoPass=1;
470   } else if (fileName.Contains("pass2") ) {
471     fRecoPass=2;
472   } else if (fileName.Contains("pass3") ) {
473     fRecoPass=3;
474   } else if (fileName.Contains("pass4") ) {
475     fRecoPass=4;
476   } else if (fileName.Contains("pass5") ) {
477     fRecoPass=5;
478   } else if (fileName.Contains("pass6") ) {
479     fRecoPass=6;
480   }
481   if (fRecoPass == 0) {
482     AliInfo(Form("From file name %s reco pass cannot be detected",fileName.Data()));
483     AliInfo("Change file name or use SetUserRecoPass method");
484     AliFatal("------------- TOF tender cannot run with reco pass unspecified, issuing FATAL error ---------- ");
485   }
486 }
487
488
489 //______________________________________________________________________________
490 void AliTOFTenderSupply::InitGeom()
491 {
492
493   if (fGeomSet == kTRUE) return;
494
495   //  Printf("\n \n ----- calling InitGeom to fix TRD Bug ----- \n \n");
496   AliCDBManager * man = AliCDBManager::Instance();
497   man->SetDefaultStorage("raw://");
498   fCDBkey = man->SetLock(kFALSE, fCDBkey);
499   Int_t run = fTender->GetRun();
500   //  Printf(" ---------> run is %d",run);
501   man->SetRun(run);
502   AliCDBEntry *entry = (AliCDBEntry*)man->Get("GRP/Geometry/Data");
503   if (entry) {
504     AliGeomManager::LoadGeometry();
505     AliGeomManager::ApplyAlignObjsFromCDB("ITS TPC TRD TOF");
506     //    fCDBkey = man->SetLock(kTRUE, fCDBkey);
507     //    Printf("\n \n ----- Geometry loaded ------ \n \n");
508   }
509   fGeomSet=kTRUE;
510
511 }
512
513
514 //______________________________________________________________________________
515 void AliTOFTenderSupply::FixTRDBug(AliESDEvent* event)
516 //
517 // recompute texp fixing wrong dE/dx from TRD (LHC10b,c pass3)
518 //
519 {
520
521   if (fGeomSet == kFALSE) InitGeom();
522
523   //  Printf("Running FixTRD bug ");
524   /* loop over tracks */
525   AliESDtrack *track = NULL;
526   for (Int_t itrk = 0; itrk < event->GetNumberOfTracks(); itrk++) {
527     track = event->GetTrack(itrk);
528     FixTRDBug(track);
529   }
530 }
531
532
533 //_____________________________________________________
534 void AliTOFTenderSupply::FixTRDBug(AliESDtrack *track)
535 {
536   // 
537   //
538   //
539
540
541     ULong_t status=track->GetStatus();
542     if (!( ( (status & AliVTrack::kITSrefit)==AliVTrack::kITSrefit ) &&
543            ( (status & AliVTrack::kTPCrefit)==AliVTrack::kTPCrefit ) &&
544            ( (status & AliVTrack::kTPCout)==AliVTrack::kTPCout ) &&
545            ( (status & AliVTrack::kTOFout)==AliVTrack::kTOFout ) &&
546            ( (status & AliVTrack::kTIME)==AliVTrack::kTIME ) ) ) return;
547
548     fIsEnteringInTRD=kFALSE;
549     fInTRD=kFALSE;
550     fIsComingOutTRD=kFALSE;
551     fOutTRD=kFALSE;
552
553     //    Printf("Track reached TOF %f",track->P());
554     Double_t correctionTimes[5] = {0.,0.,0.,0.,0.}; // to be added to the expected times
555     FindTRDFix(track, correctionTimes);
556     Double_t expectedTimes[5] = {0.,0.,0.,0.,0.}; track->GetIntegratedTimes(expectedTimes);
557     //    Printf("Exp. times: %f %f %f %f %f",
558     //     expectedTimes[0],expectedTimes[1],expectedTimes[2],expectedTimes[3],expectedTimes[4]);
559     //    Printf("Corr. times: %f %f %f %f %f",
560     //     correctionTimes[0],correctionTimes[1],correctionTimes[2],correctionTimes[3],correctionTimes[4]);
561
562     for (Int_t jj=0; jj<5; jj++) expectedTimes[jj]+=correctionTimes[jj];
563     track->SetIntegratedTimes(expectedTimes);
564   
565 }
566
567
568 //________________________________________________________________________
569 void AliTOFTenderSupply::FindTRDFix(AliESDtrack *track,Double_t *corrections)
570 {
571
572   Double_t pT = track->Pt();
573   ULong_t status=track->GetStatus();
574   Bool_t isTRDout = (status & AliVTrack::kTRDout)==AliVTrack::kTRDout;
575
576   Double_t length = 0.;
577
578   Double_t xyzIN[3]={0.,0.,0.};
579   fIsEnteringInTRD = track->GetXYZAt(fRhoTRDin,fMagField,xyzIN);
580
581   Double_t xyzOUT[3]={0.,0.,0.};
582   fIsComingOutTRD = track->GetXYZAt(fRhoTRDout,fMagField,xyzOUT);
583
584   if (fIsEnteringInTRD && fIsComingOutTRD) {
585
586
587     Double_t phiIN = TMath::Pi()+TMath::ATan2(-xyzIN[1],-xyzIN[0]);
588     phiIN *= TMath::RadToDeg();
589     fInTRD = ( (phiIN>=  0. && phiIN<= 40.) ||
590                (phiIN>=140. && phiIN<=220.) ||
591                (phiIN>=340. && phiIN<=360.) ); // TRD SMs installed @ 2010
592
593     Double_t phiOUT = TMath::Pi()+TMath::ATan2(-xyzOUT[1],-xyzOUT[0]);
594     phiOUT *= TMath::RadToDeg();
595     fOutTRD = ( (phiOUT>=  0. && phiOUT<= 40.) ||
596                 (phiOUT>=140. && phiOUT<=220.) ||
597                 (phiOUT>=340. && phiOUT<=360.) ); // TRD SMs installed @ 2010
598
599     length = 0.;
600
601     if (fInTRD || fOutTRD) {
602
603       if ( ( fInTRD && fOutTRD ) || ( fInTRD && !fOutTRD ) ) {
604         length = EstimateLengthInTRD1(track);
605       } else if ( !fInTRD && fOutTRD ) {
606         length = EstimateLengthInTRD2(track);
607       }
608
609     } else { // ( !fInTRD && !fOutTRD )
610
611       length = EstimateLengthOutTRD(track);
612
613     }
614
615   }
616   //  Printf("estimated length in TRD %f [isTRDout %d]",length,isTRDout);
617   CorrectDeltaTimes(pT,length,isTRDout,corrections);
618
619 }
620
621 //________________________________________________________________________
622 void AliTOFTenderSupply::CorrectDeltaTimes(Double_t pT,
623                                                 Double_t length,
624                                                 Bool_t flagTRDout,
625                                                 Double_t *corrections)
626 {
627
628   corrections[2] = CorrectExpectedPionTime(pT,length,flagTRDout);
629   corrections[0] = corrections[2]; // x electrons used pion corrections
630   corrections[1] = corrections[2]; // x muons used pion corrections
631   corrections[3] = CorrectExpectedKaonTime(pT,length,flagTRDout);
632   corrections[4] = CorrectExpectedProtonTime(pT,length,flagTRDout);
633
634 }
635
636 //________________________________________________________________________
637 Double_t AliTOFTenderSupply::CorrectExpectedPionTime(Double_t pT,
638                                                           Double_t length,
639                                                           Bool_t isTRDout)
640 {
641   // correction for expected time for pions
642
643   Double_t delta=0.;
644
645   //  Printf("Flags Ent In ComingOut Out %d %d %d %d",fIsEnteringInTRD,fInTRD,fIsComingOutTRD,fOutTRD);
646   if (!fIsEnteringInTRD || !fIsComingOutTRD) { // zone 5
647
648     Float_t p[2]={0.,0.};
649
650     if (isTRDout) {
651
652       if (pT<0.30) {
653         p[0] = 180.; p[1] = 0.;
654       } else if (pT>=0.30 && pT<0.35) {
655         p[0] = 740.; p[1] = -1800.;
656       } else if (pT>=0.35 && pT<0.40) {
657         p[0] = 488.; p[1] =-1080.;
658       } else if (pT>=0.40 && pT<0.45) {
659         p[0] = 179.; p[1] = -307.;
660       } else if (pT>=0.45 && pT<0.50) {
661         p[0] =  97.; p[1] = -123.;
662       } else { //if (pT>=0.50)
663         p[0] = 120.; p[1] = -172.;
664       }
665
666     } else {
667
668       if (pT<0.30) {
669         p[0] =  70.; p[1] =    0.;
670       } else if (pT>=0.30 && pT<0.35) {
671         p[0] = 339.; p[1] = -927.;
672       } else if (pT>=0.35 && pT<0.40) {
673         p[0] =  59.; p[1] = -121.;
674       } else if (pT>=0.40 && pT<0.50) {
675         p[0] =  21.; p[1] =  -24.;
676       } else { //if (pT>=0.50)
677         p[0] =  42.; p[1] =  -67.;
678       }
679
680     }
681
682     delta = p[0]+p[1]*pT;
683     //    Printf("Pion time: %f %f %f %f",p[0],p[1],length,delta);
684
685   } else {
686
687     Float_t p[2] = {0.,0.};
688
689     if ( fInTRD &&  fOutTRD) { // zone 1
690
691       if (isTRDout) {
692
693         if (length<130.) {
694           p[0] = 0.; p[1] = 0.;
695         } else if (length>=130. && length<170.) {
696           p[0] = -20.5; p[1] = 0.25;
697         } else {//if (length>=170.)
698           p[0] = 22.; p[1] = 0.;
699         }
700
701       } else { // !isTRDout
702
703         p[0] = 20.; p[1] = 0.;
704
705       }
706
707     } else if (!fInTRD && !fOutTRD) { // zone 2
708
709       p[0] = 0.; p[1] = 0.;
710
711     } else if ( fInTRD &&  !fOutTRD) { // zone 3
712
713       if (isTRDout) {
714
715         if (length< 75.) {
716           p[0] = 17.; p[1] =  0.;
717         } else if (length>= 75. && length< 95.) {
718           p[0] = 81.; p[1] = -0.85;
719         } else if (length>= 95. && length<155.) {
720           p[0] =  0.; p[1] =  0.;
721         } else {//if (length>=155.)
722           p[0] = 10.; p[1] =  0.;
723         }
724
725       } else { // !isTRDout
726
727         p[0] = 0.; p[1] = 0.;
728
729       }
730
731     } else if (!fInTRD &&  fOutTRD) { // zone 4
732
733       if (isTRDout) {
734
735         if (length<80.) {
736           p[0] =  0.; p[1] =  0.;
737         } else {//if (length>=80.)
738           p[0] = 10.; p[1] =  0.;
739         }
740
741       } else { // !isTRDout
742
743         if (length<30.) {
744           p[0] =  0.; p[1] =  0.;
745         } else {//if (length>=30.)
746           p[0] =  6.; p[1] =  0.;
747         }
748
749       }
750
751     }
752
753     delta = p[0]+p[1]*length;
754     //    Printf("Pion time: %f %f %f %f",p[0],p[1],length,delta);
755
756   }
757
758   return delta;
759
760 }
761
762 //________________________________________________________________________
763 Double_t AliTOFTenderSupply::CorrectExpectedKaonTime(Double_t pT,
764                                                           Double_t length,
765                                                           Bool_t isTRDout)
766 {
767   // correction for expected time for kaons
768
769   Double_t delta=0.;
770   //  Printf("Flags Ent In ComingOut Out %d %d %d %d",fIsEnteringInTRD,fInTRD,fIsComingOutTRD,fOutTRD);
771
772   if (!fIsEnteringInTRD || !fIsComingOutTRD) { // zone 5
773
774     Float_t p[2]={0.,0.};
775
776     if (isTRDout) {
777
778       if (pT<0.4) {
779         p[0] =  900.; p[1] = 0.;
780       } else if (pT>=0.40 && pT<0.45) {
781         p[0] = 3100.; p[1] = -6000.;
782       } else if (pT>=0.45 && pT<0.50) {
783         p[0] = 1660.; p[1] = -2800.;
784       } else if (pT>=0.50 && pT<0.55) {
785         p[0] =  860.; p[1] = -1200.;
786       } else { //if (pT>=0.55)
787         p[0] =  200.; p[1] = 0.;
788       }
789
790     } else {
791
792       if (pT<0.30) {
793         p[0] =   0.; p[1] =    0.;
794       } else if (pT>=0.30 && pT<0.32) {
795         p[0] = 570.; p[1] =    0.;
796       } else if (pT>=0.32 && pT<0.35) {
797         p[0] = 3171.; p[1] = -8133.;
798       } else if (pT>=0.35 && pT<0.40) {
799         p[0] = 1815.; p[1] = -4260.;
800       } else if (pT>=0.40 && pT<0.45) {
801         p[0] =  715.; p[1] =  -1471.;
802       } else if (pT>=0.45 && pT<0.50) {
803         p[0] =  233.; p[1] =  -407.;
804       } else if (pT>=0.50 && pT<0.55) {
805         p[0] =  408.; p[1] =  -752.;
806       } else { //if (pT>=0.55)
807         p[0] =  408.-752.*0.55; p[1] = 0.;
808       }
809
810     }
811
812     delta = p[0]+p[1]*pT;
813     //    Printf("Kaon time: %f %f %f %f",p[0],p[1],length,delta);
814
815   } else {
816
817     Float_t p[2] = {0.,0.};
818
819     if ( fInTRD &&  fOutTRD) { // zone 1
820
821       if (isTRDout) {
822
823         if (length<95.) {
824           p[0] = 20.; p[1] = 0.;
825         } else if (length>=95. && length<195.) {
826           p[0] = -24.0; p[1] = 0.10+0.0041*length;
827         } else {//if (length>=195.)
828           p[0] =  150.; p[1] = 0.;
829         }
830
831       } else { // !isTRDout
832
833         p[0] = 40.; p[1] = 0.;
834
835       }
836
837     } else if (!fInTRD && !fOutTRD) { // zone 2
838
839       p[0] = 0.; p[1] = 0.;
840
841     } else if ( fInTRD &&  !fOutTRD) { // zone 3
842
843       if (isTRDout) {
844
845         if (length< 15.) {
846           p[0] = 180.; p[1] =  0.;
847         } else if (length>= 15. && length< 55.) {
848           p[0] = 215.; p[1] = -2.5;
849         } else {//if (length>=55.)
850           p[0] = 78.; p[1] =  0.;
851         }
852
853       } else { // !isTRDout
854
855         p[0] = 0.; p[1] = 0.;
856
857       }
858
859     } else if (!fInTRD &&  fOutTRD) { // zone 4
860
861       if (isTRDout) {
862
863         if (length< 55.) {
864           p[0] =  0.; p[1] =  0.;
865         } else if (length>= 55. && length<115.) {
866           p[0] = -85.; p[1] = 1.9;
867         } else {//if (length>=115.)
868           p[0] = 100.; p[1] =  0.;
869         }
870
871       } else { // !isTRDout
872
873         p[0] =  0.; p[1] =  0.;
874
875       }
876
877     }
878
879     delta = p[0]+p[1]*length;
880     //    Printf("Kaon time: %f %f %f %f",p[0],p[1],length,delta);
881
882   }
883
884   return delta;
885
886 }
887
888 //________________________________________________________________________
889 Double_t AliTOFTenderSupply::CorrectExpectedProtonTime(Double_t pT,
890                                                             Double_t length,
891                                                             Bool_t isTRDout)
892 {
893   // correction for expected time for protons
894
895   Double_t delta=0.;
896   //  Printf("Flags Ent In ComingOut Out %d %d %d %d",fIsEnteringInTRD,fInTRD,fIsComingOutTRD,fOutTRD);
897
898   if (!fIsEnteringInTRD || !fIsComingOutTRD) { // zone 5
899     Float_t p[2]={0.,0.};
900
901
902     if (isTRDout) {
903
904       if (pT<0.375) {
905         p[0] = 1000.; p[1] = 0.;
906       } else if (pT>=0.375 && pT<0.45) {
907         p[0] = 1500.; p[1] = 0.;
908       } else if (pT>=0.45 && pT<0.50) {
909         p[0] = 4650.; p[1] = -7000.;
910       } else if (pT>=0.50 && pT<0.55) {
911         p[0] = 3150.; p[1] = -4000.;
912       } else { //if (pT>=0.55)
913         p[0] = 3150. -4000.*0.55; p[1] = 0.;
914       }
915
916     } else {
917
918       if (pT<0.32) {
919         p[0] = 2963.-5670.*0.032; p[1] = 0.;
920       } else if (pT>=0.32 && pT<0.35) {
921         p[0] = 2963.; p[1] =  -5670.;
922       } else if (pT>=0.35 && pT<0.40) {
923         p[0] = 4270.; p[1] =  -9400.;
924       } else if (pT>=0.40 && pT<0.45) {
925         p[0] = 1550.; p[1] =  -2600.;
926       } else if (pT>=0.45 && pT<0.50) {
927         p[0] = 1946.; p[1] =  -3480.;
928       } else if (pT>=0.50 && pT<0.55) {
929         p[0] = 1193.; p[1] =  -1974.;
930       } else { //if (pT>=0.55)
931         p[0] = 1193.-1974.*0.55; p[1] = 0.;
932       }
933
934     }
935
936     delta = p[0]+p[1]*pT;
937     //    Printf("Proton time: %f %f %f %f",p[0],p[1],length,delta);
938
939   } else {
940
941     Float_t p[2] = {0.,0.};
942
943     if ( fInTRD &&  fOutTRD) { // zone 1
944
945       if (isTRDout) {
946
947         if (length<90.) {
948           p[0] = 0.; p[1] = 0.;
949         } else if (length>=90. && length<200.) {
950           p[0] = 1063.; p[1] = -32.+0.30*length-0.00072*length*length;
951         } else {//if (length>=200.)
952           p[0] =  900.; p[1] = 0.;
953         }
954
955       } else { // !isTRDout
956
957         p[0] = 80.; p[1] = 0.;
958
959       }
960
961     } else if (!fInTRD && !fOutTRD) { // zone 2
962
963       if (isTRDout) {
964         p[0] = 0.; p[1] = 0.;
965       } else {
966         if (length<125.) {
967           p[0] = 0.; p[1] = 0.;
968         } else if (length>=125. && length<180.) {
969           p[0] = -132.; p[1] = 1.3;
970         } else {
971           p[0] = 100.; p[1] = 0.;
972         }
973
974       }
975
976     } else if ( fInTRD &&  !fOutTRD) { // zone 3
977
978       if (isTRDout) {
979
980         if (length< 30.) {
981           p[0] = 670.; p[1] =  0.;
982         } else if (length>= 30. && length<155.) {
983           p[0] = 944.; p[1] = -11.+0.064*length;
984         } else {//if (length>=155.)
985           p[0] = 780.; p[1] =  0.;
986         }
987
988       } else { // !isTRDout
989
990         if (length< 30.) {
991           p[0] = 140.; p[1] = -4.5;
992         } else {
993           p[0] = 0.; p[1] = 0.;
994         }
995       
996       }
997
998     } else if (!fInTRD &&  fOutTRD) { // zone 4
999
1000       if (isTRDout) {
1001
1002         if (length< 45.) {
1003           p[0] = 130.; p[1] =  0.;
1004         } else if (length>= 45. && length<120.) {
1005           p[0] = -190.; p[1] = 6.5;
1006         } else {//if (length>=120.)
1007           p[0] = 750.; p[1] =  0.;
1008         }
1009
1010       } else { // !isTRDout
1011
1012         if (length<75.5) {
1013           p[0] =    0.; p[1] =  0.;
1014         } else if (length>= 75.5 && length<90.) {
1015           p[0] = -830.; p[1] = 11.;
1016         } else {
1017           p[0] =  160.; p[1] =  0.;
1018         }
1019
1020       }
1021
1022     }
1023
1024     delta = p[0]+p[1]*length;
1025     //    Printf("Proton time: %f %f %f %f",p[0],p[1],length,delta);
1026
1027   }
1028
1029   return delta;
1030
1031 }
1032
1033 //________________________________________________________________________
1034 Double_t AliTOFTenderSupply::EstimateLengthInTRD1(AliESDtrack *track)
1035 {
1036
1037   Double_t xyz0[3]={0.,0.,0.};
1038   Bool_t stayInTRD = track->GetXYZAt(fRhoTRDin,fMagField,xyz0);
1039
1040   Double_t phi0 = TMath::Pi()+TMath::ATan2(-xyz0[1],-xyz0[0]);
1041   phi0 *= TMath::RadToDeg();
1042   stayInTRD = stayInTRD && ( (phi0>=  0. && phi0<= 40.) ||
1043                              (phi0>=140. && phi0<=220.) ||
1044                              (phi0>=340. && phi0<=360.) );
1045
1046   Double_t trackLengthInTRD = 0.;
1047   Int_t iStep=0;
1048
1049   Double_t b[3];track->GetBxByBz(b);
1050
1051   Double_t xyz1[3]={0.,0.,0.};
1052   Double_t rho = fRhoTRDin;
1053   while (stayInTRD && rho<=fRhoTRDout) {
1054     iStep++;
1055     rho += fStep;
1056
1057     for (Int_t ii=0; ii<3; ii++) xyz1[ii]=0.;
1058     stayInTRD = track->GetXYZAt(rho,fMagField,xyz1);
1059     Double_t phi1 = TMath::Pi()+TMath::ATan2(-xyz1[1],-xyz1[0]);
1060     phi1 *= TMath::RadToDeg();
1061     stayInTRD = stayInTRD && ( (phi1>=  0. && phi1<= 40.) ||
1062                                (phi1>=140. && phi1<=220.) ||
1063                                (phi1>=340. && phi1<=360.) );
1064
1065     Double_t l2  = TMath::Sqrt((xyz1[0]-xyz0[0])*(xyz1[0]-xyz0[0]) +
1066                                (xyz1[1]-xyz0[1])*(xyz1[1]-xyz0[1]) +
1067                                (xyz1[2]-xyz0[2])*(xyz1[2]-xyz0[2]));
1068     trackLengthInTRD += l2;
1069
1070     for (Int_t ii=0; ii<3; ii++) xyz0[ii]=xyz1[ii];
1071   }
1072
1073   return trackLengthInTRD;
1074
1075 }
1076
1077 //________________________________________________________________________
1078 Double_t AliTOFTenderSupply::EstimateLengthInTRD2(AliESDtrack *track)
1079 {
1080
1081   Double_t xyz0[3]={0.,0.,0.};
1082   Bool_t stayInTRD = track->GetXYZAt(fRhoTRDout,fMagField,xyz0);
1083
1084   Double_t phi0 = TMath::Pi()+TMath::ATan2(-xyz0[1],-xyz0[0]);
1085   phi0 *= TMath::RadToDeg();
1086   stayInTRD = stayInTRD && ( (phi0>=  0. && phi0<= 40.) ||
1087                              (phi0>=140. && phi0<=220.) ||
1088                              (phi0>=340. && phi0<=360.) );
1089
1090   Double_t trackLengthInTRD = 0.;
1091   Int_t iStep=0;
1092
1093   Double_t b[3];track->GetBxByBz(b);
1094
1095   Double_t xyz1[3]={0.,0.,0.};
1096   Double_t rho = fRhoTRDout;
1097   while (stayInTRD && rho>=fRhoTRDin) {
1098     iStep++;
1099     rho -= fStep;
1100
1101     for (Int_t ii=0; ii<3; ii++) xyz1[ii]=0.;
1102     stayInTRD = track->GetXYZAt(rho,fMagField,xyz1);
1103     Double_t phi1 = TMath::Pi()+TMath::ATan2(-xyz1[1],-xyz1[0]);
1104     phi1 *= TMath::RadToDeg();
1105     stayInTRD = stayInTRD && ( (phi1>=  0. && phi1<= 40.) ||
1106                                (phi1>=140. && phi1<=220.) ||
1107                                (phi1>=340. && phi1<=360.) );
1108
1109     Double_t l2  = TMath::Sqrt((xyz0[0]-xyz1[0])*(xyz0[0]-xyz1[0]) +
1110                                (xyz0[1]-xyz1[1])*(xyz0[1]-xyz1[1]) +
1111                                (xyz0[2]-xyz1[2])*(xyz0[2]-xyz1[2]));
1112     trackLengthInTRD += l2;
1113
1114     for (Int_t ii=0; ii<3; ii++) xyz0[ii]=xyz1[ii];
1115   }
1116
1117   return trackLengthInTRD;
1118
1119 }
1120
1121 //________________________________________________________________________
1122 Double_t AliTOFTenderSupply::EstimateLengthOutTRD(AliESDtrack *track)
1123 {
1124
1125   Double_t xyz0[3]={0.,0.,0.};
1126   Bool_t stayInTRD = track->GetXYZAt(fRhoTRDin,fMagField,xyz0);
1127
1128   Double_t phi0 = TMath::Pi()+TMath::ATan2(-xyz0[1],-xyz0[0]);
1129   phi0 *= TMath::RadToDeg();
1130   stayInTRD = stayInTRD && !( (phi0>=  0. && phi0<= 40.) ||
1131                               (phi0>=140. && phi0<=220.) ||
1132                               (phi0>=340. && phi0<=360.) );
1133
1134   Double_t trackLengthInTRD = 0.;
1135   Int_t iStep=0;
1136
1137   Double_t b[3];track->GetBxByBz(b);
1138
1139   Double_t xyz1[3]={0.,0.,0.};
1140   Double_t rho = fRhoTRDin;
1141   while (stayInTRD && rho<=fRhoTRDout) {
1142     iStep++;
1143     rho += fStep;
1144
1145     for (Int_t ii=0; ii<3; ii++) xyz1[ii]=0.;
1146     stayInTRD = track->GetXYZAt(rho,fMagField,xyz1);
1147     Double_t phi1 = TMath::Pi()+TMath::ATan2(-xyz1[1],-xyz1[0]);
1148     phi1 *= TMath::RadToDeg();
1149     stayInTRD = stayInTRD && !( (phi1>=  0. && phi1<= 40.) ||
1150                                 (phi1>=140. && phi1<=220.) ||
1151                                 (phi1>=340. && phi1<=360.) );
1152
1153     Double_t l2  = TMath::Sqrt((xyz1[0]-xyz0[0])*(xyz1[0]-xyz0[0]) +
1154                                (xyz1[1]-xyz0[1])*(xyz1[1]-xyz0[1]) +
1155                                (xyz1[2]-xyz0[2])*(xyz1[2]-xyz0[2]));
1156     trackLengthInTRD += l2;
1157
1158     for (Int_t ii=0; ii<3; ii++) xyz0[ii]=xyz1[ii];
1159   }
1160
1161   return trackLengthInTRD;
1162
1163 }
1164
1165 //________________________________________________________________________
1166 Int_t AliTOFTenderSupply::GetOCDBVersion(Int_t runNo)
1167 {
1168   Int_t verNo = -1;
1169   if ( (runNo>=118503 && runNo <=121040) )  { // LHC10C
1170     if (fRecoPass == 2) {                     // on pass2
1171       if (runNo >= 119159 && runNo <= 119163) verNo=3;
1172       else if (runNo >= 119837 && runNo <= 119934) verNo=4;
1173       else if (runNo >= 120067 && runNo <= 120244) verNo=4;
1174       else if (runNo >= 120503 && runNo <= 120505) verNo=4;
1175       else if (runNo >= 120616 && runNo <= 120671) verNo=4;
1176       else if (runNo >= 120741 && runNo <= 120829) verNo=4;
1177     }
1178   }
1179   return verNo;
1180 }
1181
1182
1183 //__________________________________________________________________________
1184 void AliTOFTenderSupply::LoadTOFPIDParams(Int_t runNumber)
1185 {
1186   //
1187   // Load the TOF pid params from the OADB
1188   //
1189
1190   if (fTOFPIDParams) delete fTOFPIDParams;
1191   fTOFPIDParams=0x0;
1192   
1193   TFile *oadbf = new TFile("$ALICE_ROOT/OADB/COMMON/PID/data/TOFPIDParams.root");
1194   if (oadbf && oadbf->IsOpen()) {
1195     AliInfo("Loading TOF Params from $ALICE_ROOT/OADB/COMMON/PID/data/TOFPIDParams.root");
1196     AliOADBContainer *oadbc = (AliOADBContainer *)oadbf->Get("TOFoadb");
1197     if (oadbc) fTOFPIDParams = dynamic_cast<AliTOFPIDParams *>(oadbc->GetObject(runNumber,"TOFparams"));
1198     oadbf->Close();
1199     delete oadbc;
1200   }
1201   delete oadbf;
1202
1203   if (!fTOFPIDParams) {
1204     AliError("TOFPIDParams.root not found in $ALICE_ROOT/OADB/COMMON/PID/data !!");
1205     fTOFPIDParams = new AliTOFPIDParams;
1206     fTOFPIDParams->SetTOFresolution(90.);
1207     fTOFPIDParams->SetStartTimeMethod(AliESDpid::kTOF_T0);
1208   }  
1209 }