1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 ///////////////////////////////////////////////////////////////////////////////
19 // TOF tender: - load updated calibrations (for TOF and T0)
20 // - set tofHeader if missing in ESDs (2010 data)
23 // Contacts: Pietro.Antonioli@bo.infn.it //
24 // Francesco.Noferini@bo.infn.it //
25 ///////////////////////////////////////////////////////////////////////////////
32 #include <AliESDEvent.h>
33 #include <AliESDtrack.h>
34 #include <AliESDInputHandler.h>
35 #include <AliAnalysisManager.h>
36 #include <AliESDpid.h>
37 #include <AliTender.h>
39 #include <AliTOFcalib.h>
40 #include <AliTOFT0maker.h>
42 #include <AliGeomManager.h>
43 #include <AliCDBManager.h>
44 #include <AliCDBEntry.h>
46 #include <AliOADBContainer.h>
47 #include <AliTOFPIDParams.h>
49 #include <AliT0CalibSeasonTimeShift.h>
51 #include "AliTOFTenderSupply.h"
53 ClassImp(AliTOFTenderSupply)
55 Float_t AliTOFTenderSupply::fgT0Aresolution = 75.;
56 Float_t AliTOFTenderSupply::fgT0Cresolution = 65.;
58 AliTOFTenderSupply::AliTOFTenderSupply() :
61 fTenderNoAction(kTRUE),
63 fCorrectExpTimes(kTRUE),
64 fCorrectTRDBug(kFALSE),
66 fT0DetectorAdjust(kFALSE),
68 fAutomaticSettings(kTRUE),
71 fForceCorrectTRDBug(kFALSE),
75 fT0IntercalibrationShift(0),
77 fIsEnteringInTRD(kFALSE),
79 fIsComingOutTRD(kFALSE),
81 fRhoTRDin(288.38), // cm
82 fRhoTRDout(366.38), // cm
99 //_____________________________________________________
100 AliTOFTenderSupply::AliTOFTenderSupply(const char *name, const AliTender *tender) :
101 AliTenderSupply(name,tender),
103 fTenderNoAction(kTRUE),
105 fCorrectExpTimes(kTRUE),
106 fCorrectTRDBug(kFALSE),
107 fLHC10dPatch(kFALSE),
108 fT0DetectorAdjust(kFALSE),
110 fAutomaticSettings(kTRUE),
113 fForceCorrectTRDBug(kFALSE),
117 fT0IntercalibrationShift(0),
119 fIsEnteringInTRD(kFALSE),
121 fIsComingOutTRD(kFALSE),
123 fRhoTRDin(288.38), // cm
124 fRhoTRDout(366.38), // cm
140 //_____________________________________________________
141 void AliTOFTenderSupply::Init()
144 fTenderNoAction = kFALSE;
145 // Initialise TOF tender (this is called at each detected run change)
146 AliLog::SetClassDebugLevel("AliTOFTenderSupply",10);
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
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");
160 if (fAutomaticSettings) {
162 if (fUserRecoPass == 0) DetectRecoPass();
163 else fRecoPass = fUserRecoPass;
166 fTenderNoAction = kTRUE;
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;}
172 fT0IntercalibrationShift = 0;
173 fT0DetectorAdjust=kTRUE;
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;}
179 fT0IntercalibrationShift = 0;
180 fT0DetectorAdjust=kFALSE;
182 else if (run>=122195&&run<=126437) { //period="LHC10D";
183 fCorrectExpTimes=kFALSE;
185 fT0IntercalibrationShift = 0;
186 fT0DetectorAdjust=kTRUE;
188 else if (run>=127719&&run<=130850) { //period="LHC10E";
189 fCorrectExpTimes=kFALSE;
191 fT0IntercalibrationShift = 30.;
192 fT0DetectorAdjust=kTRUE;
194 else if (run>=133004&&run<=135029) { //period="LHC10F";
195 fTenderNoAction=kTRUE;
197 else if (run>=135654&&run<=136377) { //period="LHC10G";
198 fTenderNoAction=kTRUE;
200 else if (run>=136851&&run<=139517) { //period="LHC10H" - pass2;
201 fCorrectExpTimes=kFALSE;
203 fT0IntercalibrationShift = 0.;
204 fT0DetectorAdjust=kTRUE;
206 else if (run>=139699) { //period="LHC11A";
207 fTenderNoAction=kTRUE;
211 if (fTenderNoAction) {
212 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(" |---------------------------------------------------------------------------|");
224 // Load from OADB TOF resolution
225 LoadTOFPIDParams(run);
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();
231 fESDpid=new AliESDpid;
232 fTender->GetESDhandler()->SetESDpid(fESDpid);
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)
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);
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("|******************************************************|");
264 //_____________________________________________________
265 void AliTOFTenderSupply::ProcessEvent()
268 // Use updated calibrations for TOF and T0, reapply PID information
269 // For MC: timeZero sampling and additional smearing for T0
271 if (fDebugLevel > 1) AliInfo("process event");
273 AliESDEvent *event=fTender->GetEvent();
275 if (fDebugLevel > 1) AliInfo("event read");
279 if (fTender->RunChanged()){
282 if (fTenderNoAction) return;
283 Int_t versionNumber = GetOCDBVersion(fTender->GetRun());
284 fTOFCalib->SetRunParamsSpecificVersion(versionNumber);
285 fTOFCalib->Init(fTender->GetRun());
287 if(event->GetT0TOF()){ // read T0 detector correction from OCDB
289 if (fT0DetectorAdjust) {
290 AliCDBManager* ocdbMan = AliCDBManager::Instance();
291 ocdbMan->SetRun(fTender->GetRun());
292 AliCDBEntry *entry = ocdbMan->Get("T0/Calib/TimeAdjust/");
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;
302 for (Int_t i=0;i<4;i++) fT0shift[i]=0;
303 AliWarning("TofTender no T0 entry found T0shift set to 0");
306 for (Int_t i=0;i<4;i++) fT0shift[i]=0;
311 if (fTenderNoAction) return;
313 fTOFCalib->CalibrateESD(event); //recalculate TOF signal (no harm for MC, see settings inside init)
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
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
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));
326 // event by event TO detector treatment
327 if(event->GetT0TOF()){ // protection: we adjust T0 only if it is there....
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.);
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]);
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);
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);
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);
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));
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.);
367 if (fDebugLevel > 1) Printf(" TofTender: T0 time (FINAL) %f %f %f",event->GetT0TOF(0),event->GetT0TOF(1),event->GetT0TOF(2));
369 //compute timeZero of the event via TOF-TO
370 fTOFT0maker->ComputeT0TOF(event);
371 fTOFT0maker->WriteInESD(event);
373 // set preferred startTime: this is now done via AliPIDResponseTask
374 fESDpid->SetTOFResponse(event, (AliESDpid::EStartTimeType_t)fTOFPIDParams->GetStartTimeMethod());
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);
392 //_____________________________________________________
393 void AliTOFTenderSupply::RecomputeTExp(AliESDEvent *event) const
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);
410 //_____________________________________________________
411 void AliTOFTenderSupply::RecomputeTExp(AliESDtrack *track) const
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
420 Double_t texp[AliPID::kSPECIES];
421 if (!track || !(track->GetStatus() & AliESDtrack::kTOFout)) return;
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);
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);
443 //______________________________________________________________________________
444 void AliTOFTenderSupply::DetectRecoPass()
447 // Detect reconstruction information
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;
458 TTree *tree= (TTree*)inputHandler->GetTree();
459 TFile *file= (TFile*)tree->GetCurrentFile();
462 AliFatal("Current file not found");
466 //find pass from file name (UGLY, but not stored in ESD... )
467 TString fileName(file->GetName());
468 if (fileName.Contains("pass1") ) {
470 } else if (fileName.Contains("pass2") ) {
472 } else if (fileName.Contains("pass3") ) {
474 } else if (fileName.Contains("pass4") ) {
476 } else if (fileName.Contains("pass5") ) {
478 } else if (fileName.Contains("pass6") ) {
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 ---------- ");
489 //______________________________________________________________________________
490 void AliTOFTenderSupply::InitGeom()
493 if (fGeomSet == kTRUE) return;
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);
502 AliCDBEntry *entry = (AliCDBEntry*)man->Get("GRP/Geometry/Data");
504 AliGeomManager::LoadGeometry();
505 AliGeomManager::ApplyAlignObjsFromCDB("ITS TPC TRD TOF");
506 // fCDBkey = man->SetLock(kTRUE, fCDBkey);
507 // Printf("\n \n ----- Geometry loaded ------ \n \n");
514 //______________________________________________________________________________
515 void AliTOFTenderSupply::FixTRDBug(AliESDEvent* event)
517 // recompute texp fixing wrong dE/dx from TRD (LHC10b,c pass3)
521 if (fGeomSet == kFALSE) InitGeom();
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);
533 //_____________________________________________________
534 void AliTOFTenderSupply::FixTRDBug(AliESDtrack *track)
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;
548 fIsEnteringInTRD=kFALSE;
550 fIsComingOutTRD=kFALSE;
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]);
562 for (Int_t jj=0; jj<5; jj++) expectedTimes[jj]+=correctionTimes[jj];
563 track->SetIntegratedTimes(expectedTimes);
568 //________________________________________________________________________
569 void AliTOFTenderSupply::FindTRDFix(AliESDtrack *track,Double_t *corrections)
572 Double_t pT = track->Pt();
573 ULong_t status=track->GetStatus();
574 Bool_t isTRDout = (status & AliVTrack::kTRDout)==AliVTrack::kTRDout;
576 Double_t length = 0.;
578 Double_t xyzIN[3]={0.,0.,0.};
579 fIsEnteringInTRD = track->GetXYZAt(fRhoTRDin,fMagField,xyzIN);
581 Double_t xyzOUT[3]={0.,0.,0.};
582 fIsComingOutTRD = track->GetXYZAt(fRhoTRDout,fMagField,xyzOUT);
584 if (fIsEnteringInTRD && fIsComingOutTRD) {
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
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
601 if (fInTRD || fOutTRD) {
603 if ( ( fInTRD && fOutTRD ) || ( fInTRD && !fOutTRD ) ) {
604 length = EstimateLengthInTRD1(track);
605 } else if ( !fInTRD && fOutTRD ) {
606 length = EstimateLengthInTRD2(track);
609 } else { // ( !fInTRD && !fOutTRD )
611 length = EstimateLengthOutTRD(track);
616 // Printf("estimated length in TRD %f [isTRDout %d]",length,isTRDout);
617 CorrectDeltaTimes(pT,length,isTRDout,corrections);
621 //________________________________________________________________________
622 void AliTOFTenderSupply::CorrectDeltaTimes(Double_t pT,
625 Double_t *corrections)
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);
636 //________________________________________________________________________
637 Double_t AliTOFTenderSupply::CorrectExpectedPionTime(Double_t pT,
641 // correction for expected time for pions
645 // Printf("Flags Ent In ComingOut Out %d %d %d %d",fIsEnteringInTRD,fInTRD,fIsComingOutTRD,fOutTRD);
646 if (!fIsEnteringInTRD || !fIsComingOutTRD) { // zone 5
648 Float_t p[2]={0.,0.};
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.;
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.;
682 delta = p[0]+p[1]*pT;
683 // Printf("Pion time: %f %f %f %f",p[0],p[1],length,delta);
687 Float_t p[2] = {0.,0.};
689 if ( fInTRD && fOutTRD) { // zone 1
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.;
701 } else { // !isTRDout
703 p[0] = 20.; p[1] = 0.;
707 } else if (!fInTRD && !fOutTRD) { // zone 2
709 p[0] = 0.; p[1] = 0.;
711 } else if ( fInTRD && !fOutTRD) { // zone 3
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.;
725 } else { // !isTRDout
727 p[0] = 0.; p[1] = 0.;
731 } else if (!fInTRD && fOutTRD) { // zone 4
736 p[0] = 0.; p[1] = 0.;
737 } else {//if (length>=80.)
738 p[0] = 10.; p[1] = 0.;
741 } else { // !isTRDout
744 p[0] = 0.; p[1] = 0.;
745 } else {//if (length>=30.)
746 p[0] = 6.; p[1] = 0.;
753 delta = p[0]+p[1]*length;
754 // Printf("Pion time: %f %f %f %f",p[0],p[1],length,delta);
762 //________________________________________________________________________
763 Double_t AliTOFTenderSupply::CorrectExpectedKaonTime(Double_t pT,
767 // correction for expected time for kaons
770 // Printf("Flags Ent In ComingOut Out %d %d %d %d",fIsEnteringInTRD,fInTRD,fIsComingOutTRD,fOutTRD);
772 if (!fIsEnteringInTRD || !fIsComingOutTRD) { // zone 5
774 Float_t p[2]={0.,0.};
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.;
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.;
812 delta = p[0]+p[1]*pT;
813 // Printf("Kaon time: %f %f %f %f",p[0],p[1],length,delta);
817 Float_t p[2] = {0.,0.};
819 if ( fInTRD && fOutTRD) { // zone 1
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.;
831 } else { // !isTRDout
833 p[0] = 40.; p[1] = 0.;
837 } else if (!fInTRD && !fOutTRD) { // zone 2
839 p[0] = 0.; p[1] = 0.;
841 } else if ( fInTRD && !fOutTRD) { // zone 3
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.;
853 } else { // !isTRDout
855 p[0] = 0.; p[1] = 0.;
859 } else if (!fInTRD && fOutTRD) { // zone 4
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.;
871 } else { // !isTRDout
873 p[0] = 0.; p[1] = 0.;
879 delta = p[0]+p[1]*length;
880 // Printf("Kaon time: %f %f %f %f",p[0],p[1],length,delta);
888 //________________________________________________________________________
889 Double_t AliTOFTenderSupply::CorrectExpectedProtonTime(Double_t pT,
893 // correction for expected time for protons
896 // Printf("Flags Ent In ComingOut Out %d %d %d %d",fIsEnteringInTRD,fInTRD,fIsComingOutTRD,fOutTRD);
898 if (!fIsEnteringInTRD || !fIsComingOutTRD) { // zone 5
899 Float_t p[2]={0.,0.};
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.;
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.;
936 delta = p[0]+p[1]*pT;
937 // Printf("Proton time: %f %f %f %f",p[0],p[1],length,delta);
941 Float_t p[2] = {0.,0.};
943 if ( fInTRD && fOutTRD) { // zone 1
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.;
955 } else { // !isTRDout
957 p[0] = 80.; p[1] = 0.;
961 } else if (!fInTRD && !fOutTRD) { // zone 2
964 p[0] = 0.; p[1] = 0.;
967 p[0] = 0.; p[1] = 0.;
968 } else if (length>=125. && length<180.) {
969 p[0] = -132.; p[1] = 1.3;
971 p[0] = 100.; p[1] = 0.;
976 } else if ( fInTRD && !fOutTRD) { // zone 3
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.;
988 } else { // !isTRDout
991 p[0] = 140.; p[1] = -4.5;
993 p[0] = 0.; p[1] = 0.;
998 } else if (!fInTRD && fOutTRD) { // zone 4
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.;
1010 } else { // !isTRDout
1013 p[0] = 0.; p[1] = 0.;
1014 } else if (length>= 75.5 && length<90.) {
1015 p[0] = -830.; p[1] = 11.;
1017 p[0] = 160.; p[1] = 0.;
1024 delta = p[0]+p[1]*length;
1025 // Printf("Proton time: %f %f %f %f",p[0],p[1],length,delta);
1033 //________________________________________________________________________
1034 Double_t AliTOFTenderSupply::EstimateLengthInTRD1(AliESDtrack *track)
1037 Double_t xyz0[3]={0.,0.,0.};
1038 Bool_t stayInTRD = track->GetXYZAt(fRhoTRDin,fMagField,xyz0);
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.) );
1046 Double_t trackLengthInTRD = 0.;
1049 Double_t b[3];track->GetBxByBz(b);
1051 Double_t xyz1[3]={0.,0.,0.};
1052 Double_t rho = fRhoTRDin;
1053 while (stayInTRD && rho<=fRhoTRDout) {
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.) );
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;
1070 for (Int_t ii=0; ii<3; ii++) xyz0[ii]=xyz1[ii];
1073 return trackLengthInTRD;
1077 //________________________________________________________________________
1078 Double_t AliTOFTenderSupply::EstimateLengthInTRD2(AliESDtrack *track)
1081 Double_t xyz0[3]={0.,0.,0.};
1082 Bool_t stayInTRD = track->GetXYZAt(fRhoTRDout,fMagField,xyz0);
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.) );
1090 Double_t trackLengthInTRD = 0.;
1093 Double_t b[3];track->GetBxByBz(b);
1095 Double_t xyz1[3]={0.,0.,0.};
1096 Double_t rho = fRhoTRDout;
1097 while (stayInTRD && rho>=fRhoTRDin) {
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.) );
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;
1114 for (Int_t ii=0; ii<3; ii++) xyz0[ii]=xyz1[ii];
1117 return trackLengthInTRD;
1121 //________________________________________________________________________
1122 Double_t AliTOFTenderSupply::EstimateLengthOutTRD(AliESDtrack *track)
1125 Double_t xyz0[3]={0.,0.,0.};
1126 Bool_t stayInTRD = track->GetXYZAt(fRhoTRDin,fMagField,xyz0);
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.) );
1134 Double_t trackLengthInTRD = 0.;
1137 Double_t b[3];track->GetBxByBz(b);
1139 Double_t xyz1[3]={0.,0.,0.};
1140 Double_t rho = fRhoTRDin;
1141 while (stayInTRD && rho<=fRhoTRDout) {
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.) );
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;
1158 for (Int_t ii=0; ii<3; ii++) xyz0[ii]=xyz1[ii];
1161 return trackLengthInTRD;
1165 //________________________________________________________________________
1166 Int_t AliTOFTenderSupply::GetOCDBVersion(Int_t runNo)
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;
1183 //__________________________________________________________________________
1184 void AliTOFTenderSupply::LoadTOFPIDParams(Int_t runNumber)
1187 // Load the TOF pid params from the OADB
1190 if (fTOFPIDParams) delete fTOFPIDParams;
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"));
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);