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 **************************************************************************/
15 /* $Id: AliTOFT0maker.cxx,v 1.8 2010/01/19 16:32:20 noferini Exp $ */
17 /////////////////////////////////////////////////////////////////////////////
19 // This class contains the basic functions for the time zero //
20 // evaluation with TOF detector informations. //
21 // Use case in an analysis task: //
23 // Create the object in the task constructor (fTOFmaker is a private var) //
24 // fTOFmaker = new AliTOFT0maker(); //
25 // fTOFmaker->SetTimeResolution(130.0); // if you want set the TOF res //
26 // 115 ps is the TOF default resolution value //
28 // Use the RemakePID method in the task::Exec //
29 // Double_t* calcolot0; //
30 // calcolot0=fTOFmaker->RemakePID(fESD); //
31 // //calcolot0[0] = calculated event time //
32 // //calcolot0[1] = event time time resolution //
33 // //calcolot0[2] = average event time for the current fill //
34 // //calcolot0[3] = tracks at TOF //
35 // //calcolot0[4] = calculated event time (only TOF) //
36 // //calcolot0[5] = event time time resolution (only TOF) //
37 // //calcolot0[6] = sigma t0 fill //
38 // //calcolot0[7] = tracks at TOF really used in tht algorithm //
40 // Let consider that: //
41 // - the PIF is automatically recalculated with the event time subtrction //
43 /////////////////////////////////////////////////////////////////////////////
45 #include "AliTOFT0v1.h"
46 #include "AliTOFT0maker.h"
47 #include "AliTOFcalibHisto.h"
49 #include "AliESDpid.h"
50 #include "AliESDEvent.h"
54 ClassImp(AliTOFT0maker)
56 //____________________________________________________________________________
57 AliTOFT0maker::AliTOFT0maker():
58 fCalib(new AliTOFcalibHisto()),
74 // fCalib->SetCalibParFileName("./AliTOFcalibPar.LHC10b.7000GeV.20100405.root");
75 fCalib->LoadCalibPar();
77 if(AliPID::ParticleMass(0) == 0) new AliPID();
81 //____________________________________________________________________________
82 AliTOFT0maker::AliTOFT0maker(const AliTOFT0maker & t) :
88 fESDswitch(t.fESDswitch),
89 fTimeResolution(t.fTimeResolution),
91 fHmapChannel(t.fHmapChannel),
97 //____________________________________________________________________________
98 AliTOFT0maker& AliTOFT0maker::operator=(const AliTOFT0maker &t)
107 fESDswitch = t.fESDswitch;
108 fTimeResolution = t.fTimeResolution;
109 fT0sigma = t.fT0sigma;
113 //____________________________________________________________________________
114 AliTOFT0maker::~AliTOFT0maker()
117 if(fCalib) delete fCalib;
119 //____________________________________________________________________________
120 Double_t* AliTOFT0maker::RemakePID(AliESDEvent *esd,Double_t t0time,Double_t t0sigma){
122 // Remake TOF PID probabilities
127 if(fKmask) ApplyMask(esd);
129 AliTOFT0v1* t0maker=new AliTOFT0v1(esd);
130 // t0maker->SetCalib(fCalib);
131 t0maker->SetTimeResolution(fTimeResolution*1e-12);
134 TakeTimeRawCorrection(esd);
137 t0tof=t0maker->DefineT0("all");
139 Float_t lT0Current=0.;
142 // Int_t nrun = esd->GetRunNumber();
143 Double_t t0fill = GetT0Fill();
146 Float_t sigmaFill = (t0fill - Int_t(t0fill))*1000;
147 if(sigmaFill < 0) sigmaFill += 1000;
149 fCalculated[0]=-1000*t0tof[0]; // best t0
150 fCalculated[1]=1000*t0tof[1]; // sigma best t0
151 fCalculated[2] = t0fill; //t0 fill
152 fCalculated[3] = t0tof[2]; // n TOF tracks
153 fCalculated[4]=-1000*t0tof[0]; // TOF t0
154 fCalculated[5]=1000*t0tof[1]; // TOF t0 sigma
155 fCalculated[6]=sigmaFill; // sigma t0 fill
156 fCalculated[7] = t0tof[3]; // n TOF tracks used for T0
158 if(fCalculated[1] < sigmaFill){
160 fT0fill[fiT0] = fCalculated[0];
161 fT0sigmaTOF[fiT0] = fCalculated[1];
165 else if(TMath::Abs(fCalculated[0] - t0fill) < 500){
166 fT0fill[fiT0] = fCalculated[0];
167 fT0sigmaTOF[fiT0] = fCalculated[1];
172 // printf("%i - %i) %f\n",fiT0,fnT0,t0fill);
176 if(fiT0 > fgkNmaxT0step-1) fiT0=0;
182 fCalculated[2] = t0fill; //t0 fill
185 if(fCalculated[1] < sigmaFill && TMath::Abs(fCalculated[0] - t0fill) < 500){
186 fT0sigma=fCalculated[1];
187 lT0Current=fCalculated[0];
190 fCalculated[4] = t0fill;
191 fCalculated[5] = sigmaFill;
194 if(fCalculated[1] < 1 || fT0sigma > sigmaFill){
196 fCalculated[4] = t0fill;
197 fCalculated[5] = sigmaFill;
202 Double_t w1 = 1./t0sigma/t0sigma;
203 Double_t w2 = 1./fCalculated[1]/fCalculated[1];
205 Double_t wtot = w1+w2;
207 lT0Current = (w1*t0time + w2*fCalculated[0]) / wtot;
208 fT0sigma = TMath::Sqrt(1./wtot);
216 if(fT0sigma < sigmaFill && TMath::Abs(lT0Current - t0fill) < 500){
217 fCalculated[1]=fT0sigma;
218 fCalculated[0]=lT0Current;
221 if(fT0sigma >= 1000 || fNoTOFT0){
223 fT0sigma = sigmaFill;
225 fCalculated[0] = t0fill;
226 fCalculated[1] = sigmaFill;
231 RemakeTOFpid(esd,lT0Current);
235 //____________________________________________________________________________
236 void AliTOFT0maker::TakeTimeRawCorrection(AliESDEvent * const esd){
238 // Take raw corrections for time measurements
241 Int_t ntracks = esd->GetNumberOfTracks();
244 AliESDtrack *t=esd->GetTrack(ntracks);
246 if ((t->GetStatus()&AliESDtrack::kTOFout)==0) continue;
248 Double_t time=t->GetTOFsignalRaw();
250 Double_t tot = t->GetTOFsignalToT();
251 Int_t chan = t->GetTOFCalChannel();
252 Double_t corr = fCalib->GetFullCorrection(chan,tot) - fCalib->GetCorrection(AliTOFcalibHisto::kTimeSlewingCorr,chan,0);
255 t->SetTOFsignal(time);
259 //____________________________________________________________________________
260 void AliTOFT0maker::RemakeTOFpid(AliESDEvent *esd,Float_t timezero){
262 // Recalculate TOF PID probabilities
266 pidESD.GetTOFResponse().SetTimeResolution(TMath::Sqrt(fT0sigma*fT0sigma + fTimeResolution*fTimeResolution));
267 pidESD.MakePID(esd,kFALSE,timezero);
270 //____________________________________________________________________________
271 Double_t AliTOFT0maker::GetT0Fill() const {
273 // Return T0 of filling
280 if(n >10 && n <= 20) n = 10;
285 if(n > fgkNmaxT0step) n = fgkNmaxT0step;
289 Double_t lT0sigmaav=0;
291 for(Int_t i=0;i<n;i++){
293 lT0sigmaav += fT0sigmaTOF[fiT0];
294 lT0avErr+=fT0fill[i]*fT0fill[i];
296 lT0avErr -= lT0av*lT0av/n;
299 lT0avErr = TMath::Sqrt(TMath::Max(lT0avErr/(n-1) - lT0sigmaav*lT0sigmaav,0.00001));
302 if(lT0avErr > 300) lT0avErr = 300;
304 lT0av = Int_t(lT0av) + lT0avErr/1000.;
312 //____________________________________________________________________________
313 void AliTOFT0maker::LoadChannelMap(char *filename){
314 // Load the histo with the channel off map
315 TFile *f= new TFile(filename);
317 printf("Cannot open the channel map file (%s)\n",filename);
321 fHmapChannel = (TH1F *) f->Get("hChEnabled");
324 printf("Cannot laod the channel map histo (from %s)\n",filename);
329 //____________________________________________________________________________
330 void AliTOFT0maker::ApplyMask(AliESDEvent * const esd){
331 // Switch off the disable channel
333 printf("Channel Map is not available\n");
337 Int_t ntracks = esd->GetNumberOfTracks();
340 AliESDtrack *t=esd->GetTrack(ntracks);
342 if ((t->GetStatus()&AliESDtrack::kTOFout)==0) continue;
344 Int_t chan = t->GetTOFCalChannel();
346 if(fHmapChannel->GetBinContent(chan) < 0.01){
347 t->ResetStatus(AliESDtrack::kTOFout);
352 //____________________________________________________________________________