]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/TenderSupplies/AliTOFTenderSupply.cxx
disabled actions on T0 data
[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=kFALSE;   // it was 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=kFALSE;     // it was kTRUE
187     }
188     else if (run>=127719&&run<=130850) { //period="LHC10E";
189       fCorrectExpTimes=kFALSE;
190       fLHC10dPatch=kFALSE;
191       fT0IntercalibrationShift = 30.;
192       fT0DetectorAdjust=kFALSE;      // it was 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=kFALSE;      // it was 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 (fT0DetectorAdjust) {
334       if(!fIsMC){   // data: apply shifts to align around Zero
335         event->SetT0TOF(0,event->GetT0TOF(0) - fT0shift[0]);
336         event->SetT0TOF(1,event->GetT0TOF(1) - fT0shift[1]);
337         event->SetT0TOF(2,event->GetT0TOF(2) - fT0shift[2]);
338       } else {
339       // MC: add smearing for realistic T0A and T0C resolution
340         Double_t defResolutionT0A = 33.;   // in future we will get this from ESDrun data structure or via OCDB
341         Double_t defResolutionT0C = 30.;   // for the moment we don't trust them
342         if ( (fgT0Aresolution > defResolutionT0A) && (event->GetT0TOF(1)<90000.) ) { // add smearing only if signal is there
343           Double_t addedSmearingT0A = TMath::Sqrt(fgT0Aresolution*fgT0Aresolution - defResolutionT0A*defResolutionT0A);
344           Double_t smearingT0A = gRandom->Gaus(0.,addedSmearingT0A);
345           event->SetT0TOF(1,event->GetT0TOF(1) + smearingT0A);
346         }
347         if ( (fgT0Cresolution > defResolutionT0C) && (event->GetT0TOF(2)<90000.) ) { // add smearing only if signal is there
348         Double_t addedSmearingT0C = TMath::Sqrt(fgT0Cresolution*fgT0Cresolution - defResolutionT0C*defResolutionT0C);
349         Double_t smearingT0C = gRandom->Gaus(0.,addedSmearingT0C);
350         event->SetT0TOF(2,event->GetT0TOF(2) + smearingT0C);
351         }
352         if (event->GetT0TOF(0)<90000.) { // we recompute the AND only if it is already there...
353           Double_t smearedT0AC = (event->GetT0TOF(1)+event->GetT0TOF(2))/2.;
354           event->SetT0TOF(0,smearedT0AC); 
355         }
356         if (fDebugLevel > 1) Printf(" TofTender: T0 time (postSmear) %f %f %f",event->GetT0TOF(0),event->GetT0TOF(1),event->GetT0TOF(2));
357         // add finally the timeZero offset also to the T0 detector information
358         event->SetT0TOF(0,event->GetT0TOF(0) + startTime);
359         event->SetT0TOF(1,event->GetT0TOF(1) + startTime);
360         event->SetT0TOF(2,event->GetT0TOF(2) + startTime);  
361         if (fDebugLevel > 1) Printf(" TofTender: T0 time (postStart) %f %f %f",event->GetT0TOF(0),event->GetT0TOF(1),event->GetT0TOF(2));
362       }
363     }
364     // after shifts adjust (data) or smearing+offset (MC) we 'clean' to default if signals not there 
365     if(event->GetT0TOF(0) > 900000) event->SetT0TOF(0, 999999.);
366     if(event->GetT0TOF(1) > 90000)  event->SetT0TOF(1, 99999.);
367     if(event->GetT0TOF(2) > 90000)  event->SetT0TOF(2, 99999.);
368   }
369   if (fDebugLevel > 1) Printf(" TofTender: T0 time (FINAL) %f %f %f",event->GetT0TOF(0),event->GetT0TOF(1),event->GetT0TOF(2));
370   
371   //compute timeZero of the event via TOF-TO
372   fTOFT0maker->ComputeT0TOF(event);
373   fTOFT0maker->WriteInESD(event);
374
375   //  set preferred startTime: this is now done via AliPIDResponseTask
376   fESDpid->SetTOFResponse(event, (AliESDpid::EStartTimeType_t)fTOFPIDParams->GetStartTimeMethod());
377
378   // recalculate PID probabilities
379   // this is for safety, especially if the user doesn't attach a PID tender after TOF tender  
380   Int_t ntracks=event->GetNumberOfTracks();
381   //  AliESDtrack *track = NULL;
382   //  Float_t tzeroTrack = 0;
383   for(Int_t itrack = 0; itrack < ntracks; itrack++){
384     //    track = event->GetTrack(itrack);
385     //    tzeroTrack = fESDpid->GetTOFResponse().GetStartTime(track->P());
386     //    Printf("================> Track # %d mom: %f tzeroTrack %f",itrack,track->P(),tzeroTrack);
387     fESDpid->MakeTOFPID(event->GetTrack(itrack),0);   
388   }
389   
390   
391 }
392
393
394 //_____________________________________________________
395 void AliTOFTenderSupply::RecomputeTExp(AliESDEvent *event) const
396 {
397   /*
398    * calibrate TExp
399    */
400
401   
402   /* loop over tracks */
403   AliESDtrack *track = NULL;
404   for (Int_t itrk = 0; itrk < event->GetNumberOfTracks(); itrk++) {
405     /* get track and calibrate */
406     track = event->GetTrack(itrk);
407     RecomputeTExp(track);
408   }
409   
410 }
411
412 //_____________________________________________________
413 void AliTOFTenderSupply::RecomputeTExp(AliESDtrack *track) const
414 {
415   /*** 
416        THIS METHOD IS BASED ON THEORETICAL EXPECTED TIME COMPUTED
417        USING AVERAGE MOMENTUM BETWEEN INNER/OUTER TRACK PARAMS 
418        IT IS A ROUGH APPROXIMATION APPLIED TO FIX LHC10d-pass2 DATA
419        WHERE A WRONG GEOMETRY (FULL TRD) WAS INSERTED
420   ***/
421
422   Double_t texp[AliPID::kSPECIES];
423   if (!track || !(track->GetStatus() & AliESDtrack::kTOFout)) return;
424
425
426   /* get track params */
427   Float_t l = track->GetIntegratedLength();
428   Float_t p = track->P();
429   if (track->GetInnerParam() && track->GetOuterParam()) {
430     Float_t pin = track->GetInnerParam()->P();
431     Float_t pout = track->GetOuterParam()->P();
432     p = 0.5 * (pin + pout);
433   }
434   /* loop over particle types and compute expected time */
435   for (Int_t ipart = 0; ipart < AliPID::kSPECIES; ipart++)
436     texp[ipart] = GetExpTimeTh(AliPID::ParticleMass(ipart), p, l) - 37.; 
437   // 37 is a final semiempirical offset to further adjust (calibrations were
438   // done with "standard" integratedTimes)
439   /* set integrated times */
440   track->SetIntegratedTimes(texp);
441
442 }
443
444
445 //______________________________________________________________________________
446 void AliTOFTenderSupply::DetectRecoPass()
447 {
448   //
449   // Detect reconstruction information
450   //
451   
452   //reset information
453   fRecoPass=0;
454   
455   //Get the current file to check the reconstruction pass (UGLY, but not stored in ESD... )
456   AliAnalysisManager *mgr=AliAnalysisManager::GetAnalysisManager();
457   AliVEventHandler *inputHandler=mgr->GetInputEventHandler();
458   if (!inputHandler) return;
459   
460   TTree *tree= (TTree*)inputHandler->GetTree();
461   TFile *file= (TFile*)tree->GetCurrentFile();
462   
463   if (!file) {
464     AliFatal("Current file not found");
465     return; // coverity
466   }
467   
468   //find pass from file name (UGLY, but not stored in ESD... )
469   TString fileName(file->GetName());
470   if (fileName.Contains("pass1") ) {
471     fRecoPass=1;
472   } else if (fileName.Contains("pass2") ) {
473     fRecoPass=2;
474   } else if (fileName.Contains("pass3") ) {
475     fRecoPass=3;
476   } else if (fileName.Contains("pass4") ) {
477     fRecoPass=4;
478   } else if (fileName.Contains("pass5") ) {
479     fRecoPass=5;
480   } else if (fileName.Contains("pass6") ) {
481     fRecoPass=6;
482   }
483   if (fRecoPass == 0) {
484     AliInfo(Form("From file name %s reco pass cannot be detected",fileName.Data()));
485     AliInfo("Change file name or use SetUserRecoPass method");
486     AliFatal("------------- TOF tender cannot run with reco pass unspecified, issuing FATAL error ---------- ");
487   }
488 }
489
490
491 //______________________________________________________________________________
492 void AliTOFTenderSupply::InitGeom()
493 {
494
495   if (fGeomSet == kTRUE) return;
496
497   //  Printf("\n \n ----- calling InitGeom to fix TRD Bug ----- \n \n");
498   AliCDBManager * man = AliCDBManager::Instance();
499   man->SetDefaultStorage("raw://");
500   fCDBkey = man->SetLock(kFALSE, fCDBkey);
501   Int_t run = fTender->GetRun();
502   //  Printf(" ---------> run is %d",run);
503   man->SetRun(run);
504   AliCDBEntry *entry = (AliCDBEntry*)man->Get("GRP/Geometry/Data");
505   if (entry) {
506     AliGeomManager::LoadGeometry();
507     AliGeomManager::ApplyAlignObjsFromCDB("ITS TPC TRD TOF");
508     //    fCDBkey = man->SetLock(kTRUE, fCDBkey);
509     //    Printf("\n \n ----- Geometry loaded ------ \n \n");
510   }
511   fGeomSet=kTRUE;
512
513 }
514
515
516 //______________________________________________________________________________
517 void AliTOFTenderSupply::FixTRDBug(AliESDEvent* event)
518 //
519 // recompute texp fixing wrong dE/dx from TRD (LHC10b,c pass3)
520 //
521 {
522
523   if (fGeomSet == kFALSE) InitGeom();
524
525   //  Printf("Running FixTRD bug ");
526   /* loop over tracks */
527   AliESDtrack *track = NULL;
528   for (Int_t itrk = 0; itrk < event->GetNumberOfTracks(); itrk++) {
529     track = event->GetTrack(itrk);
530     FixTRDBug(track);
531   }
532 }
533
534
535 //_____________________________________________________
536 void AliTOFTenderSupply::FixTRDBug(AliESDtrack *track)
537 {
538   // 
539   //
540   //
541
542
543     ULong_t status=track->GetStatus();
544     if (!( ( (status & AliVTrack::kITSrefit)==AliVTrack::kITSrefit ) &&
545            ( (status & AliVTrack::kTPCrefit)==AliVTrack::kTPCrefit ) &&
546            ( (status & AliVTrack::kTPCout)==AliVTrack::kTPCout ) &&
547            ( (status & AliVTrack::kTOFout)==AliVTrack::kTOFout ) &&
548            ( (status & AliVTrack::kTIME)==AliVTrack::kTIME ) ) ) return;
549
550     fIsEnteringInTRD=kFALSE;
551     fInTRD=kFALSE;
552     fIsComingOutTRD=kFALSE;
553     fOutTRD=kFALSE;
554
555     //    Printf("Track reached TOF %f",track->P());
556     Double_t correctionTimes[5] = {0.,0.,0.,0.,0.}; // to be added to the expected times
557     FindTRDFix(track, correctionTimes);
558     Double_t expectedTimes[5] = {0.,0.,0.,0.,0.}; track->GetIntegratedTimes(expectedTimes);
559     //    Printf("Exp. times: %f %f %f %f %f",
560     //     expectedTimes[0],expectedTimes[1],expectedTimes[2],expectedTimes[3],expectedTimes[4]);
561     //    Printf("Corr. times: %f %f %f %f %f",
562     //     correctionTimes[0],correctionTimes[1],correctionTimes[2],correctionTimes[3],correctionTimes[4]);
563
564     for (Int_t jj=0; jj<5; jj++) expectedTimes[jj]+=correctionTimes[jj];
565     track->SetIntegratedTimes(expectedTimes);
566   
567 }
568
569
570 //________________________________________________________________________
571 void AliTOFTenderSupply::FindTRDFix(AliESDtrack *track,Double_t *corrections)
572 {
573
574   Double_t pT = track->Pt();
575   ULong_t status=track->GetStatus();
576   Bool_t isTRDout = (status & AliVTrack::kTRDout)==AliVTrack::kTRDout;
577
578   Double_t length = 0.;
579
580   Double_t xyzIN[3]={0.,0.,0.};
581   fIsEnteringInTRD = track->GetXYZAt(fRhoTRDin,fMagField,xyzIN);
582
583   Double_t xyzOUT[3]={0.,0.,0.};
584   fIsComingOutTRD = track->GetXYZAt(fRhoTRDout,fMagField,xyzOUT);
585
586   if (fIsEnteringInTRD && fIsComingOutTRD) {
587
588
589     Double_t phiIN = TMath::Pi()+TMath::ATan2(-xyzIN[1],-xyzIN[0]);
590     phiIN *= TMath::RadToDeg();
591     fInTRD = ( (phiIN>=  0. && phiIN<= 40.) ||
592                (phiIN>=140. && phiIN<=220.) ||
593                (phiIN>=340. && phiIN<=360.) ); // TRD SMs installed @ 2010
594
595     Double_t phiOUT = TMath::Pi()+TMath::ATan2(-xyzOUT[1],-xyzOUT[0]);
596     phiOUT *= TMath::RadToDeg();
597     fOutTRD = ( (phiOUT>=  0. && phiOUT<= 40.) ||
598                 (phiOUT>=140. && phiOUT<=220.) ||
599                 (phiOUT>=340. && phiOUT<=360.) ); // TRD SMs installed @ 2010
600
601     length = 0.;
602
603     if (fInTRD || fOutTRD) {
604
605       if ( ( fInTRD && fOutTRD ) || ( fInTRD && !fOutTRD ) ) {
606         length = EstimateLengthInTRD1(track);
607       } else if ( !fInTRD && fOutTRD ) {
608         length = EstimateLengthInTRD2(track);
609       }
610
611     } else { // ( !fInTRD && !fOutTRD )
612
613       length = EstimateLengthOutTRD(track);
614
615     }
616
617   }
618   //  Printf("estimated length in TRD %f [isTRDout %d]",length,isTRDout);
619   CorrectDeltaTimes(pT,length,isTRDout,corrections);
620
621 }
622
623 //________________________________________________________________________
624 void AliTOFTenderSupply::CorrectDeltaTimes(Double_t pT,
625                                                 Double_t length,
626                                                 Bool_t flagTRDout,
627                                                 Double_t *corrections)
628 {
629
630   corrections[2] = CorrectExpectedPionTime(pT,length,flagTRDout);
631   corrections[0] = corrections[2]; // x electrons used pion corrections
632   corrections[1] = corrections[2]; // x muons used pion corrections
633   corrections[3] = CorrectExpectedKaonTime(pT,length,flagTRDout);
634   corrections[4] = CorrectExpectedProtonTime(pT,length,flagTRDout);
635
636 }
637
638 //________________________________________________________________________
639 Double_t AliTOFTenderSupply::CorrectExpectedPionTime(Double_t pT,
640                                                           Double_t length,
641                                                           Bool_t isTRDout)
642 {
643   // correction for expected time for pions
644
645   Double_t delta=0.;
646
647   //  Printf("Flags Ent In ComingOut Out %d %d %d %d",fIsEnteringInTRD,fInTRD,fIsComingOutTRD,fOutTRD);
648   if (!fIsEnteringInTRD || !fIsComingOutTRD) { // zone 5
649
650     Float_t p[2]={0.,0.};
651
652     if (isTRDout) {
653
654       if (pT<0.30) {
655         p[0] = 180.; p[1] = 0.;
656       } else if (pT>=0.30 && pT<0.35) {
657         p[0] = 740.; p[1] = -1800.;
658       } else if (pT>=0.35 && pT<0.40) {
659         p[0] = 488.; p[1] =-1080.;
660       } else if (pT>=0.40 && pT<0.45) {
661         p[0] = 179.; p[1] = -307.;
662       } else if (pT>=0.45 && pT<0.50) {
663         p[0] =  97.; p[1] = -123.;
664       } else { //if (pT>=0.50)
665         p[0] = 120.; p[1] = -172.;
666       }
667
668     } else {
669
670       if (pT<0.30) {
671         p[0] =  70.; p[1] =    0.;
672       } else if (pT>=0.30 && pT<0.35) {
673         p[0] = 339.; p[1] = -927.;
674       } else if (pT>=0.35 && pT<0.40) {
675         p[0] =  59.; p[1] = -121.;
676       } else if (pT>=0.40 && pT<0.50) {
677         p[0] =  21.; p[1] =  -24.;
678       } else { //if (pT>=0.50)
679         p[0] =  42.; p[1] =  -67.;
680       }
681
682     }
683
684     delta = p[0]+p[1]*pT;
685     //    Printf("Pion time: %f %f %f %f",p[0],p[1],length,delta);
686
687   } else {
688
689     Float_t p[2] = {0.,0.};
690
691     if ( fInTRD &&  fOutTRD) { // zone 1
692
693       if (isTRDout) {
694
695         if (length<130.) {
696           p[0] = 0.; p[1] = 0.;
697         } else if (length>=130. && length<170.) {
698           p[0] = -20.5; p[1] = 0.25;
699         } else {//if (length>=170.)
700           p[0] = 22.; p[1] = 0.;
701         }
702
703       } else { // !isTRDout
704
705         p[0] = 20.; p[1] = 0.;
706
707       }
708
709     } else if (!fInTRD && !fOutTRD) { // zone 2
710
711       p[0] = 0.; p[1] = 0.;
712
713     } else if ( fInTRD &&  !fOutTRD) { // zone 3
714
715       if (isTRDout) {
716
717         if (length< 75.) {
718           p[0] = 17.; p[1] =  0.;
719         } else if (length>= 75. && length< 95.) {
720           p[0] = 81.; p[1] = -0.85;
721         } else if (length>= 95. && length<155.) {
722           p[0] =  0.; p[1] =  0.;
723         } else {//if (length>=155.)
724           p[0] = 10.; p[1] =  0.;
725         }
726
727       } else { // !isTRDout
728
729         p[0] = 0.; p[1] = 0.;
730
731       }
732
733     } else if (!fInTRD &&  fOutTRD) { // zone 4
734
735       if (isTRDout) {
736
737         if (length<80.) {
738           p[0] =  0.; p[1] =  0.;
739         } else {//if (length>=80.)
740           p[0] = 10.; p[1] =  0.;
741         }
742
743       } else { // !isTRDout
744
745         if (length<30.) {
746           p[0] =  0.; p[1] =  0.;
747         } else {//if (length>=30.)
748           p[0] =  6.; p[1] =  0.;
749         }
750
751       }
752
753     }
754
755     delta = p[0]+p[1]*length;
756     //    Printf("Pion time: %f %f %f %f",p[0],p[1],length,delta);
757
758   }
759
760   return delta;
761
762 }
763
764 //________________________________________________________________________
765 Double_t AliTOFTenderSupply::CorrectExpectedKaonTime(Double_t pT,
766                                                           Double_t length,
767                                                           Bool_t isTRDout)
768 {
769   // correction for expected time for kaons
770
771   Double_t delta=0.;
772   //  Printf("Flags Ent In ComingOut Out %d %d %d %d",fIsEnteringInTRD,fInTRD,fIsComingOutTRD,fOutTRD);
773
774   if (!fIsEnteringInTRD || !fIsComingOutTRD) { // zone 5
775
776     Float_t p[2]={0.,0.};
777
778     if (isTRDout) {
779
780       if (pT<0.4) {
781         p[0] =  900.; p[1] = 0.;
782       } else if (pT>=0.40 && pT<0.45) {
783         p[0] = 3100.; p[1] = -6000.;
784       } else if (pT>=0.45 && pT<0.50) {
785         p[0] = 1660.; p[1] = -2800.;
786       } else if (pT>=0.50 && pT<0.55) {
787         p[0] =  860.; p[1] = -1200.;
788       } else { //if (pT>=0.55)
789         p[0] =  200.; p[1] = 0.;
790       }
791
792     } else {
793
794       if (pT<0.30) {
795         p[0] =   0.; p[1] =    0.;
796       } else if (pT>=0.30 && pT<0.32) {
797         p[0] = 570.; p[1] =    0.;
798       } else if (pT>=0.32 && pT<0.35) {
799         p[0] = 3171.; p[1] = -8133.;
800       } else if (pT>=0.35 && pT<0.40) {
801         p[0] = 1815.; p[1] = -4260.;
802       } else if (pT>=0.40 && pT<0.45) {
803         p[0] =  715.; p[1] =  -1471.;
804       } else if (pT>=0.45 && pT<0.50) {
805         p[0] =  233.; p[1] =  -407.;
806       } else if (pT>=0.50 && pT<0.55) {
807         p[0] =  408.; p[1] =  -752.;
808       } else { //if (pT>=0.55)
809         p[0] =  408.-752.*0.55; p[1] = 0.;
810       }
811
812     }
813
814     delta = p[0]+p[1]*pT;
815     //    Printf("Kaon time: %f %f %f %f",p[0],p[1],length,delta);
816
817   } else {
818
819     Float_t p[2] = {0.,0.};
820
821     if ( fInTRD &&  fOutTRD) { // zone 1
822
823       if (isTRDout) {
824
825         if (length<95.) {
826           p[0] = 20.; p[1] = 0.;
827         } else if (length>=95. && length<195.) {
828           p[0] = -24.0; p[1] = 0.10+0.0041*length;
829         } else {//if (length>=195.)
830           p[0] =  150.; p[1] = 0.;
831         }
832
833       } else { // !isTRDout
834
835         p[0] = 40.; p[1] = 0.;
836
837       }
838
839     } else if (!fInTRD && !fOutTRD) { // zone 2
840
841       p[0] = 0.; p[1] = 0.;
842
843     } else if ( fInTRD &&  !fOutTRD) { // zone 3
844
845       if (isTRDout) {
846
847         if (length< 15.) {
848           p[0] = 180.; p[1] =  0.;
849         } else if (length>= 15. && length< 55.) {
850           p[0] = 215.; p[1] = -2.5;
851         } else {//if (length>=55.)
852           p[0] = 78.; p[1] =  0.;
853         }
854
855       } else { // !isTRDout
856
857         p[0] = 0.; p[1] = 0.;
858
859       }
860
861     } else if (!fInTRD &&  fOutTRD) { // zone 4
862
863       if (isTRDout) {
864
865         if (length< 55.) {
866           p[0] =  0.; p[1] =  0.;
867         } else if (length>= 55. && length<115.) {
868           p[0] = -85.; p[1] = 1.9;
869         } else {//if (length>=115.)
870           p[0] = 100.; p[1] =  0.;
871         }
872
873       } else { // !isTRDout
874
875         p[0] =  0.; p[1] =  0.;
876
877       }
878
879     }
880
881     delta = p[0]+p[1]*length;
882     //    Printf("Kaon time: %f %f %f %f",p[0],p[1],length,delta);
883
884   }
885
886   return delta;
887
888 }
889
890 //________________________________________________________________________
891 Double_t AliTOFTenderSupply::CorrectExpectedProtonTime(Double_t pT,
892                                                             Double_t length,
893                                                             Bool_t isTRDout)
894 {
895   // correction for expected time for protons
896
897   Double_t delta=0.;
898   //  Printf("Flags Ent In ComingOut Out %d %d %d %d",fIsEnteringInTRD,fInTRD,fIsComingOutTRD,fOutTRD);
899
900   if (!fIsEnteringInTRD || !fIsComingOutTRD) { // zone 5
901     Float_t p[2]={0.,0.};
902
903
904     if (isTRDout) {
905
906       if (pT<0.375) {
907         p[0] = 1000.; p[1] = 0.;
908       } else if (pT>=0.375 && pT<0.45) {
909         p[0] = 1500.; p[1] = 0.;
910       } else if (pT>=0.45 && pT<0.50) {
911         p[0] = 4650.; p[1] = -7000.;
912       } else if (pT>=0.50 && pT<0.55) {
913         p[0] = 3150.; p[1] = -4000.;
914       } else { //if (pT>=0.55)
915         p[0] = 3150. -4000.*0.55; p[1] = 0.;
916       }
917
918     } else {
919
920       if (pT<0.32) {
921         p[0] = 2963.-5670.*0.032; p[1] = 0.;
922       } else if (pT>=0.32 && pT<0.35) {
923         p[0] = 2963.; p[1] =  -5670.;
924       } else if (pT>=0.35 && pT<0.40) {
925         p[0] = 4270.; p[1] =  -9400.;
926       } else if (pT>=0.40 && pT<0.45) {
927         p[0] = 1550.; p[1] =  -2600.;
928       } else if (pT>=0.45 && pT<0.50) {
929         p[0] = 1946.; p[1] =  -3480.;
930       } else if (pT>=0.50 && pT<0.55) {
931         p[0] = 1193.; p[1] =  -1974.;
932       } else { //if (pT>=0.55)
933         p[0] = 1193.-1974.*0.55; p[1] = 0.;
934       }
935
936     }
937
938     delta = p[0]+p[1]*pT;
939     //    Printf("Proton time: %f %f %f %f",p[0],p[1],length,delta);
940
941   } else {
942
943     Float_t p[2] = {0.,0.};
944
945     if ( fInTRD &&  fOutTRD) { // zone 1
946
947       if (isTRDout) {
948
949         if (length<90.) {
950           p[0] = 0.; p[1] = 0.;
951         } else if (length>=90. && length<200.) {
952           p[0] = 1063.; p[1] = -32.+0.30*length-0.00072*length*length;
953         } else {//if (length>=200.)
954           p[0] =  900.; p[1] = 0.;
955         }
956
957       } else { // !isTRDout
958
959         p[0] = 80.; p[1] = 0.;
960
961       }
962
963     } else if (!fInTRD && !fOutTRD) { // zone 2
964
965       if (isTRDout) {
966         p[0] = 0.; p[1] = 0.;
967       } else {
968         if (length<125.) {
969           p[0] = 0.; p[1] = 0.;
970         } else if (length>=125. && length<180.) {
971           p[0] = -132.; p[1] = 1.3;
972         } else {
973           p[0] = 100.; p[1] = 0.;
974         }
975
976       }
977
978     } else if ( fInTRD &&  !fOutTRD) { // zone 3
979
980       if (isTRDout) {
981
982         if (length< 30.) {
983           p[0] = 670.; p[1] =  0.;
984         } else if (length>= 30. && length<155.) {
985           p[0] = 944.; p[1] = -11.+0.064*length;
986         } else {//if (length>=155.)
987           p[0] = 780.; p[1] =  0.;
988         }
989
990       } else { // !isTRDout
991
992         if (length< 30.) {
993           p[0] = 140.; p[1] = -4.5;
994         } else {
995           p[0] = 0.; p[1] = 0.;
996         }
997       
998       }
999
1000     } else if (!fInTRD &&  fOutTRD) { // zone 4
1001
1002       if (isTRDout) {
1003
1004         if (length< 45.) {
1005           p[0] = 130.; p[1] =  0.;
1006         } else if (length>= 45. && length<120.) {
1007           p[0] = -190.; p[1] = 6.5;
1008         } else {//if (length>=120.)
1009           p[0] = 750.; p[1] =  0.;
1010         }
1011
1012       } else { // !isTRDout
1013
1014         if (length<75.5) {
1015           p[0] =    0.; p[1] =  0.;
1016         } else if (length>= 75.5 && length<90.) {
1017           p[0] = -830.; p[1] = 11.;
1018         } else {
1019           p[0] =  160.; p[1] =  0.;
1020         }
1021
1022       }
1023
1024     }
1025
1026     delta = p[0]+p[1]*length;
1027     //    Printf("Proton time: %f %f %f %f",p[0],p[1],length,delta);
1028
1029   }
1030
1031   return delta;
1032
1033 }
1034
1035 //________________________________________________________________________
1036 Double_t AliTOFTenderSupply::EstimateLengthInTRD1(AliESDtrack *track)
1037 {
1038
1039   Double_t xyz0[3]={0.,0.,0.};
1040   Bool_t stayInTRD = track->GetXYZAt(fRhoTRDin,fMagField,xyz0);
1041
1042   Double_t phi0 = TMath::Pi()+TMath::ATan2(-xyz0[1],-xyz0[0]);
1043   phi0 *= TMath::RadToDeg();
1044   stayInTRD = stayInTRD && ( (phi0>=  0. && phi0<= 40.) ||
1045                              (phi0>=140. && phi0<=220.) ||
1046                              (phi0>=340. && phi0<=360.) );
1047
1048   Double_t trackLengthInTRD = 0.;
1049   Int_t iStep=0;
1050
1051   Double_t b[3];track->GetBxByBz(b);
1052
1053   Double_t xyz1[3]={0.,0.,0.};
1054   Double_t rho = fRhoTRDin;
1055   while (stayInTRD && rho<=fRhoTRDout) {
1056     iStep++;
1057     rho += fStep;
1058
1059     for (Int_t ii=0; ii<3; ii++) xyz1[ii]=0.;
1060     stayInTRD = track->GetXYZAt(rho,fMagField,xyz1);
1061     Double_t phi1 = TMath::Pi()+TMath::ATan2(-xyz1[1],-xyz1[0]);
1062     phi1 *= TMath::RadToDeg();
1063     stayInTRD = stayInTRD && ( (phi1>=  0. && phi1<= 40.) ||
1064                                (phi1>=140. && phi1<=220.) ||
1065                                (phi1>=340. && phi1<=360.) );
1066
1067     Double_t l2  = TMath::Sqrt((xyz1[0]-xyz0[0])*(xyz1[0]-xyz0[0]) +
1068                                (xyz1[1]-xyz0[1])*(xyz1[1]-xyz0[1]) +
1069                                (xyz1[2]-xyz0[2])*(xyz1[2]-xyz0[2]));
1070     trackLengthInTRD += l2;
1071
1072     for (Int_t ii=0; ii<3; ii++) xyz0[ii]=xyz1[ii];
1073   }
1074
1075   return trackLengthInTRD;
1076
1077 }
1078
1079 //________________________________________________________________________
1080 Double_t AliTOFTenderSupply::EstimateLengthInTRD2(AliESDtrack *track)
1081 {
1082
1083   Double_t xyz0[3]={0.,0.,0.};
1084   Bool_t stayInTRD = track->GetXYZAt(fRhoTRDout,fMagField,xyz0);
1085
1086   Double_t phi0 = TMath::Pi()+TMath::ATan2(-xyz0[1],-xyz0[0]);
1087   phi0 *= TMath::RadToDeg();
1088   stayInTRD = stayInTRD && ( (phi0>=  0. && phi0<= 40.) ||
1089                              (phi0>=140. && phi0<=220.) ||
1090                              (phi0>=340. && phi0<=360.) );
1091
1092   Double_t trackLengthInTRD = 0.;
1093   Int_t iStep=0;
1094
1095   Double_t b[3];track->GetBxByBz(b);
1096
1097   Double_t xyz1[3]={0.,0.,0.};
1098   Double_t rho = fRhoTRDout;
1099   while (stayInTRD && rho>=fRhoTRDin) {
1100     iStep++;
1101     rho -= fStep;
1102
1103     for (Int_t ii=0; ii<3; ii++) xyz1[ii]=0.;
1104     stayInTRD = track->GetXYZAt(rho,fMagField,xyz1);
1105     Double_t phi1 = TMath::Pi()+TMath::ATan2(-xyz1[1],-xyz1[0]);
1106     phi1 *= TMath::RadToDeg();
1107     stayInTRD = stayInTRD && ( (phi1>=  0. && phi1<= 40.) ||
1108                                (phi1>=140. && phi1<=220.) ||
1109                                (phi1>=340. && phi1<=360.) );
1110
1111     Double_t l2  = TMath::Sqrt((xyz0[0]-xyz1[0])*(xyz0[0]-xyz1[0]) +
1112                                (xyz0[1]-xyz1[1])*(xyz0[1]-xyz1[1]) +
1113                                (xyz0[2]-xyz1[2])*(xyz0[2]-xyz1[2]));
1114     trackLengthInTRD += l2;
1115
1116     for (Int_t ii=0; ii<3; ii++) xyz0[ii]=xyz1[ii];
1117   }
1118
1119   return trackLengthInTRD;
1120
1121 }
1122
1123 //________________________________________________________________________
1124 Double_t AliTOFTenderSupply::EstimateLengthOutTRD(AliESDtrack *track)
1125 {
1126
1127   Double_t xyz0[3]={0.,0.,0.};
1128   Bool_t stayInTRD = track->GetXYZAt(fRhoTRDin,fMagField,xyz0);
1129
1130   Double_t phi0 = TMath::Pi()+TMath::ATan2(-xyz0[1],-xyz0[0]);
1131   phi0 *= TMath::RadToDeg();
1132   stayInTRD = stayInTRD && !( (phi0>=  0. && phi0<= 40.) ||
1133                               (phi0>=140. && phi0<=220.) ||
1134                               (phi0>=340. && phi0<=360.) );
1135
1136   Double_t trackLengthInTRD = 0.;
1137   Int_t iStep=0;
1138
1139   Double_t b[3];track->GetBxByBz(b);
1140
1141   Double_t xyz1[3]={0.,0.,0.};
1142   Double_t rho = fRhoTRDin;
1143   while (stayInTRD && rho<=fRhoTRDout) {
1144     iStep++;
1145     rho += fStep;
1146
1147     for (Int_t ii=0; ii<3; ii++) xyz1[ii]=0.;
1148     stayInTRD = track->GetXYZAt(rho,fMagField,xyz1);
1149     Double_t phi1 = TMath::Pi()+TMath::ATan2(-xyz1[1],-xyz1[0]);
1150     phi1 *= TMath::RadToDeg();
1151     stayInTRD = stayInTRD && !( (phi1>=  0. && phi1<= 40.) ||
1152                                 (phi1>=140. && phi1<=220.) ||
1153                                 (phi1>=340. && phi1<=360.) );
1154
1155     Double_t l2  = TMath::Sqrt((xyz1[0]-xyz0[0])*(xyz1[0]-xyz0[0]) +
1156                                (xyz1[1]-xyz0[1])*(xyz1[1]-xyz0[1]) +
1157                                (xyz1[2]-xyz0[2])*(xyz1[2]-xyz0[2]));
1158     trackLengthInTRD += l2;
1159
1160     for (Int_t ii=0; ii<3; ii++) xyz0[ii]=xyz1[ii];
1161   }
1162
1163   return trackLengthInTRD;
1164
1165 }
1166
1167 //________________________________________________________________________
1168 Int_t AliTOFTenderSupply::GetOCDBVersion(Int_t runNo)
1169 {
1170   Int_t verNo = -1;
1171   if ( (runNo>=118503 && runNo <=121040) )  { // LHC10C
1172     if (fRecoPass == 2) {                     // on pass2
1173       if (runNo >= 119159 && runNo <= 119163) verNo=3;
1174       else if (runNo >= 119837 && runNo <= 119934) verNo=4;
1175       else if (runNo >= 120067 && runNo <= 120244) verNo=4;
1176       else if (runNo >= 120503 && runNo <= 120505) verNo=4;
1177       else if (runNo >= 120616 && runNo <= 120671) verNo=4;
1178       else if (runNo >= 120741 && runNo <= 120829) verNo=4;
1179     }
1180   }
1181   return verNo;
1182 }
1183
1184
1185 //__________________________________________________________________________
1186 void AliTOFTenderSupply::LoadTOFPIDParams(Int_t runNumber)
1187 {
1188   //
1189   // Load the TOF pid params from the OADB
1190   //
1191
1192   if (fTOFPIDParams) delete fTOFPIDParams;
1193   fTOFPIDParams=0x0;
1194   
1195   TFile *oadbf = new TFile("$ALICE_ROOT/OADB/COMMON/PID/data/TOFPIDParams.root");
1196   if (oadbf && oadbf->IsOpen()) {
1197     AliInfo("Loading TOF Params from $ALICE_ROOT/OADB/COMMON/PID/data/TOFPIDParams.root");
1198     AliOADBContainer *oadbc = (AliOADBContainer *)oadbf->Get("TOFoadb");
1199     if (oadbc) fTOFPIDParams = dynamic_cast<AliTOFPIDParams *>(oadbc->GetObject(runNumber,"TOFparams"));
1200     oadbf->Close();
1201     delete oadbc;
1202   }
1203   delete oadbf;
1204
1205   if (!fTOFPIDParams) {
1206     AliError("TOFPIDParams.root not found in $ALICE_ROOT/OADB/COMMON/PID/data !!");
1207     fTOFPIDParams = new AliTOFPIDParams;
1208     fTOFPIDParams->SetTOFresolution(90.);
1209     fTOFPIDParams->SetStartTimeMethod(AliESDpid::kTOF_T0);
1210   }  
1211 }