]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFtrackerV1.cxx
AliESDtrack::fTOFsignalDz corresponds to Z_track in the matched TOF pad RF. As a...
[u/mrichter/AliRoot.git] / TOF / AliTOFtrackerV1.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 // AliTOFtrackerV1 Class                                              //
19 // Task: Perform association of the ESD tracks to TOF Clusters        //
20 // and Update ESD track with associated TOF Cluster parameters        //
21 //                                                                    //
22 // -- Authors : S. Arcelli, C. Zampolli (Bologna University and INFN) //
23 // -- Contacts: Annalisa.De.Caro@cern.ch                              //
24 // --         : Chiara.Zampolli@bo.infn.it                            //
25 // --         : Silvia.Arcelli@bo.infn.it                             //
26 //                                                                    //
27 //--------------------------------------------------------------------//
28
29 #include <Rtypes.h>
30 #include <TROOT.h>
31
32 #include <TClonesArray.h>
33 #include <TTree.h>
34 #include <TFile.h>
35 #include <TH1F.h>
36 #include <TH2F.h>
37 #include <TSeqCollection.h>
38
39 #include "AliESDtrack.h"
40 #include "AliESDEvent.h"
41 #include "AliLog.h"
42 #include "AliTrackPointArray.h"
43 #include "AliGeomManager.h"
44 #include "AliCDBManager.h"
45
46 #include "AliTOFRecoParam.h"
47 #include "AliTOFReconstructor.h"
48 #include "AliTOFcluster.h"
49 #include "AliTOFGeometry.h"
50 #include "AliTOFtrackerV1.h"
51 #include "AliTOFtrack.h"
52 #include "AliTOFpidESD.h"
53
54 extern TROOT *gROOT;
55
56 ClassImp(AliTOFtrackerV1)
57
58 //_____________________________________________________________________________
59 AliTOFtrackerV1::AliTOFtrackerV1():
60   fRecoParam(0x0),
61   fPid(0x0),
62   fN(0),
63   fNseeds(0),
64   fNseedsTOF(0),
65   fngoodmatch(0),
66   fnbadmatch(0),
67   fnunmatch(0),
68   fnmatch(0),
69   fTracks(new TClonesArray("AliTOFtrack")),
70   fSeeds(new TClonesArray("AliESDtrack")),
71   fHDigClusMap(0x0),
72   fHDigNClus(0x0),
73   fHDigClusTime(0x0),
74   fHDigClusToT(0x0),
75   fHRecNClus(0x0),
76   fHRecChi2(0x0),
77   fHRecDistZ(0x0),
78   fHRecSigYVsP(0x0),
79   fHRecSigZVsP(0x0),
80   fHRecSigYVsPWin(0x0),
81   fHRecSigZVsPWin(0x0)
82  { 
83   //AliTOFtrackerV1 main Ctor
84
85    InitCheckHists();
86
87 }
88 //_____________________________________________________________________________
89 AliTOFtrackerV1::~AliTOFtrackerV1() {
90   //
91   // Dtor
92   //
93
94   SaveCheckHists();
95
96   if(!(AliCDBManager::Instance()->GetCacheFlag())){
97     delete fRecoParam;
98   }
99   delete fPid; 
100   delete fHDigClusMap;
101   delete fHDigNClus;
102   delete fHDigClusTime;
103   delete fHDigClusToT;
104   delete fHRecNClus;
105   delete fHRecChi2;
106   delete fHRecDistZ;
107   delete fHRecSigYVsP;
108   delete fHRecSigZVsP;
109   delete fHRecSigYVsPWin;
110   delete fHRecSigZVsPWin;
111   if (fTracks){
112     fTracks->Delete();
113     delete fTracks;
114     fTracks=0x0;
115   }
116   if (fSeeds){
117     fSeeds->Delete();
118     delete fSeeds;
119     fSeeds=0x0;
120   }
121 }
122 //_____________________________________________________________________________
123 Int_t AliTOFtrackerV1::PropagateBack(AliESDEvent* event) {
124   //
125   // Gets seeds from ESD event and Match with TOF Clusters
126   //
127
128   // initialize RecoParam for current event
129
130   AliInfo("Initializing params for TOF... ");
131
132   fRecoParam = AliTOFReconstructor::GetRecoParam();  // instantiate reco param from STEER...
133
134   if (fRecoParam == 0x0) { 
135     AliFatal("No Reco Param found for TOF!!!");
136   }
137   //fRecoParam->Dump();
138   //if(fRecoParam->GetApplyPbPbCuts())fRecoParam=fRecoParam->GetPbPbparam();
139   //fRecoParam->PrintParameters();
140
141   Double_t parPID[2];   
142   parPID[0]=fRecoParam->GetTimeResolution();
143   parPID[1]=fRecoParam->GetTimeNSigma();
144   fPid=new AliTOFpidESD(parPID);
145
146   //Initialise some counters
147
148   fNseeds=0;
149   fNseedsTOF=0;
150   fngoodmatch=0;
151   fnbadmatch=0;
152   fnunmatch=0;
153   fnmatch=0;
154
155   Int_t ntrk=event->GetNumberOfTracks();
156   fNseeds = ntrk;
157   TClonesArray &aESDTrack = *fSeeds;
158
159
160   //Load ESD tracks into a local Array of ESD Seeds
161
162   for (Int_t i=0; i<fNseeds; i++) {
163     AliESDtrack *t=event->GetTrack(i);
164     new(aESDTrack[i]) AliESDtrack(*t);
165   }
166
167   //Prepare ESD tracks candidates for TOF Matching
168   CollectESD();
169
170   //Matching Step
171   MatchTracks();
172
173   AliInfo(Form("Number of matched tracks: %d",fnmatch));
174   AliInfo(Form("Number of good matched tracks: %d",fngoodmatch));
175   AliInfo(Form("Number of bad  matched tracks: %d",fnbadmatch));
176
177   //Update the matched ESD tracks
178
179   for (Int_t i=0; i<ntrk; i++) {
180     AliESDtrack *t=event->GetTrack(i);
181     AliESDtrack *seed =(AliESDtrack*)fSeeds->UncheckedAt(i);
182     if(seed->GetTOFsignal()>0){
183       t->SetTOFsignal(seed->GetTOFsignal());
184       t->SetTOFcluster(seed->GetTOFcluster());
185       t->SetTOFsignalToT(seed->GetTOFsignalToT());
186       t->SetTOFsignalRaw(seed->GetTOFsignalRaw());
187       t->SetTOFsignalDz(seed->GetTOFsignalDz());
188       t->SetTOFCalChannel(seed->GetTOFCalChannel());
189       Int_t tlab[3]; seed->GetTOFLabel(tlab);    
190       t->SetTOFLabel(tlab);
191       AliTOFtrack *track = new AliTOFtrack(*seed); 
192       t->UpdateTrackParams(track,AliESDtrack::kTOFout);   
193
194       Double_t time[10]; t->GetIntegratedTimes(time);
195       AliDebug(1, Form(" %d %f %f %f %f %f", i,
196                        time[0], time[1], time[2], time[3], time[4]));
197
198       delete track;
199     }
200   }
201
202   //Handle Time Zero information
203
204   Double_t timeZero=0.;
205   Double_t timeZeroMax=99999.;
206   Bool_t usetimeZero     = fRecoParam->UseTimeZero();
207   Bool_t timeZeroFromT0  = fRecoParam->GetTimeZerofromT0();
208   Bool_t timeZeroFromTOF = fRecoParam->GetTimeZerofromTOF();
209
210   AliDebug(1,Form("Use Time Zero?: %d",usetimeZero));
211   AliDebug(1,Form("Time Zero from T0? : %d",timeZeroFromT0));
212   AliDebug(1,Form("Time Zero From TOF? : %d",timeZeroFromTOF));
213
214   if(usetimeZero){
215     if(timeZeroFromT0){
216       timeZero=GetTimeZerofromT0(event); 
217     }
218     if(timeZeroFromTOF && (timeZero>timeZeroMax || !timeZeroFromT0)){
219       timeZero=GetTimeZerofromTOF(event); 
220     }
221   }
222   AliDebug(2,Form("time Zero used in PID: %f",timeZero));
223   //Make TOF PID
224   fPid->MakePID(event,timeZero);
225
226   fSeeds->Clear();
227   fTracks->Clear();
228   return 0;
229   
230 }
231 //_________________________________________________________________________
232 void AliTOFtrackerV1::CollectESD() {
233    //prepare the set of ESD tracks to be matched to clusters in TOF
234
235   Int_t seedsTOF1=0;
236   Int_t seedsTOF2=0;
237  
238   TClonesArray &aTOFTrack = *fTracks;
239   for (Int_t i=0; i<fNseeds; i++) {
240
241     AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(i);
242     if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
243
244     // TRD 'good' tracks, already propagated at 372 cm
245
246     AliTOFtrack *track = new AliTOFtrack(*t); // New
247     Double_t x = track->GetX(); //New
248
249     if (((t->GetStatus()&AliESDtrack::kTRDout)!=0 ) && 
250          ( x >= AliTOFGeometry::RinTOF()) ){
251       track->SetSeedIndex(i);
252       t->UpdateTrackParams(track,AliESDtrack::kTOFout);    
253       new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
254       fNseedsTOF++;
255       seedsTOF1++;
256       delete track;
257     }
258
259     // Propagate the rest of TPCbp  
260
261     else {
262       if(track->PropagateToInnerTOF()){ 
263         track->SetSeedIndex(i);
264         t->UpdateTrackParams(track,AliESDtrack::kTOFout);    
265         new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
266         fNseedsTOF++;
267         seedsTOF2++;
268       }
269       delete track;
270     }
271   }
272
273   AliInfo(Form("Number of TOF seeds %d",fNseedsTOF));
274   AliInfo(Form("Number of TOF seeds Type 1 %d",seedsTOF1));
275   AliInfo(Form("Number of TOF seeds Type 2 %d",seedsTOF2));
276
277   // Sort according uncertainties on track position 
278   fTracks->Sort();
279
280 }
281 //_________________________________________________________________________
282 void AliTOFtrackerV1::MatchTracks( ){
283   //
284   //Match ESD tracks to clusters in TOF
285   //
286
287
288   // Parameters regulating the reconstruction
289   Float_t dY=AliTOFGeometry::XPad(); 
290   Float_t dZ=AliTOFGeometry::ZPad(); 
291
292   const Float_t kTimeOffset = 32.; // time offset for tracking algorithm [ps]
293
294   const Int_t kncmax = 100;
295   Float_t sensRadius = fRecoParam->GetSensRadius();
296   Float_t scaleFact   = fRecoParam->GetWindowScaleFact();
297   Float_t dyMax=fRecoParam->GetWindowSizeMaxY(); 
298   Float_t dzMax=fRecoParam->GetWindowSizeMaxZ();
299   Double_t maxChi2=fRecoParam->GetMaxChi2();
300   Bool_t timeWalkCorr    = fRecoParam->GetTimeWalkCorr();
301   AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++ \n");
302   AliDebug(1,Form("TOF sens radius: %f",sensRadius));
303   AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
304   AliDebug(1,Form("TOF Window max dy: %f",dyMax));
305   AliDebug(1,Form("TOF Window max dz: %f",dzMax));
306   AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
307   AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));   
308
309
310   //The matching loop
311
312   for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
313
314     AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
315     AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(track->GetSeedIndex());
316     if(t->GetTOFsignal()>0. ) continue;
317     AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
318      
319     // Determine a window around the track
320     Double_t x,par[5]; trackTOFin->GetExternalParameters(x,par);
321     Double_t cov[15]; trackTOFin->GetExternalCovariance(cov);
322
323     Double_t z    = par[1];   
324     Double_t dz   =  scaleFact*3.*TMath::Sqrt(cov[2]+dZ*dZ/12.);
325     Double_t dphi =  scaleFact*3.*TMath::Sqrt(cov[0]+dY*dY/12.)/sensRadius; 
326
327     Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
328     if (phi<-TMath::Pi())phi+=2*TMath::Pi();
329     if (phi>=TMath::Pi())phi-=2*TMath::Pi();
330
331     //upper limit on window's size.
332     if(dz> dzMax) dz=dzMax;
333     if(dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
334
335
336     // find the clusters inside the selected window 
337     Int_t nc=0;
338     AliTOFcluster *clusters[kncmax]; // pointers to the clusters in the window
339     Int_t index[kncmax];//to keep track of the cluster index
340     for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {  
341       AliTOFcluster *c=fClusters[k];
342       //      if(nc>kncmax)break; /* R+ fix (buffer overflow) */
343       if(nc>=kncmax)break; /* R+ fix (buffer overflow protection) */
344       if(c->GetZ() > z+dz) break;
345       if(c->IsUsed()) continue;      
346       if(!c->GetStatus()) {
347               AliDebug(1,"Cluster in channel declared bad!");
348               continue; // skip bad channels as declared in OCDB  
349       }
350       Float_t xyz[3]; c->GetGlobalXYZ(xyz);
351       Double_t clPhi=TMath::ATan2(xyz[1],xyz[0]);
352       Double_t dph=TMath::Abs(clPhi-phi);
353       if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
354       if (TMath::Abs(dph)>dphi) continue;
355       clusters[nc]=c;
356       index[nc] = k;      
357       nc++;  
358     }
359
360     //start propagation: go to the average TOF pad middle plane at ~379.5 cm
361
362     Float_t  xTOF = sensRadius;
363     Double_t ymax = xTOF*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
364     Bool_t skip = kFALSE;
365     Double_t ysect = trackTOFin->GetYat(xTOF,skip);
366     if (skip) break;
367     if (ysect > ymax) {
368       if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
369         break;
370       }
371     } else if (ysect <-ymax) {
372       if (!trackTOFin->Rotate(-AliTOFGeometry::GetAlpha())) {
373         break;
374       }
375     }
376     if(!trackTOFin->PropagateTo(xTOF)) {
377       break;
378     }
379
380
381     AliTOFcluster *bestCluster=0;
382     Double_t bestChi2=maxChi2; 
383     Int_t idclus=-1;
384     //    for (Int_t i=0; i<nc; i++){ /* R+ fix (unsafe) */
385     for (Int_t i=0; i<nc && i<kncmax; i++){ /* R+ fix (buffer overflow protection) */
386       AliTOFcluster *c=clusters[i];  // one of the preselected clusters     
387       Double_t chi2=trackTOFin->GetPredictedChi2((AliCluster3D*)c); 
388       if (chi2 >= bestChi2) continue;
389       bestChi2=chi2;
390       bestCluster=c;
391       idclus=index[i];
392     }
393     
394     if (!bestCluster) {  // no matching , go to the next track 
395       fnunmatch++;
396       delete trackTOFin;
397       continue;
398     }
399
400     fnmatch++;
401
402     AliDebug(2, Form("%7i     %7i     %10i     %10i  %10i  %10i      %7i",
403                      iseed,
404                      fnmatch-1,
405                      TMath::Abs(trackTOFin->GetLabel()),
406                      bestCluster->GetLabel(0), 
407                      bestCluster->GetLabel(1), 
408                      bestCluster->GetLabel(2),
409                      idclus)); // AdC
410
411     bestCluster->Use(); 
412     if (
413         (bestCluster->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
414         ||
415         (bestCluster->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
416         ||
417         (bestCluster->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
418         ) {
419       fngoodmatch++;
420        AliDebug(2,Form(" track label good %5d",trackTOFin->GetLabel()));
421
422     }
423     else{
424       fnbadmatch++;
425       AliDebug(2,Form(" track label bad %5d",trackTOFin->GetLabel()));
426     }
427
428     //Propagate the track to the best matched cluster
429     trackTOFin->PropagateTo(bestCluster);
430
431     // Fill the track residual histograms.
432     FillResiduals(trackTOFin,bestCluster,kFALSE);
433
434     //now take the local distance in Z from the pad center for time walk correction
435     Float_t tiltangle = AliTOFGeometry::GetAngles(bestCluster->GetDetInd(1),bestCluster->GetDetInd(2))*TMath::DegToRad();
436     Double_t dzTW=trackTOFin->GetZ()-bestCluster->GetZ(); // in cm - in the ALICE RF -
437     dzTW/=TMath::Cos(tiltangle); // from ALICE/tracking RF to pad RF (1)
438     dzTW=-dzTW; // from ALICE/tracking RF to pad RF (2)
439     if (tiltangle!=0.) AliDebug(2,Form(" rho_track = %f --- rho_cluster = %f ",trackTOFin->GetX(),bestCluster->GetX()));
440
441     //update the ESD track and delete the TOFtrack
442     t->UpdateTrackParams(trackTOFin,AliESDtrack::kTOFout);    
443     delete trackTOFin;
444
445     //  Store quantities to be used in the TOF Calibration
446     Float_t tToT=AliTOFGeometry::ToTBinWidth()*bestCluster->GetToT()*1E-3; // in ns
447     t->SetTOFsignalToT(tToT);
448     Float_t rawTime=AliTOFGeometry::TdcBinWidth()*bestCluster->GetTDCRAW()+kTimeOffset; // RAW time,in ps
449     t->SetTOFsignalRaw(rawTime);
450     t->SetTOFsignalDz(dzTW);
451
452     Int_t ind[5];
453     ind[0]=bestCluster->GetDetInd(0);
454     ind[1]=bestCluster->GetDetInd(1);
455     ind[2]=bestCluster->GetDetInd(2);
456     ind[3]=bestCluster->GetDetInd(3);
457     ind[4]=bestCluster->GetDetInd(4);
458     Int_t calindex = AliTOFGeometry::GetIndex(ind);
459     t->SetTOFCalChannel(calindex);
460
461     // keep track of the track labels in the matched cluster
462     Int_t tlab[3];
463     tlab[0]=bestCluster->GetLabel(0);
464     tlab[1]=bestCluster->GetLabel(1);
465     tlab[2]=bestCluster->GetLabel(2);
466     AliDebug(2,Form(" tdc time of the matched track %6d = ",bestCluster->GetTDC()));    
467     Double_t tof=AliTOFGeometry::TdcBinWidth()*bestCluster->GetTDC()+kTimeOffset; // in ps
468     AliDebug(2,Form(" tof time of the matched track: %f = ",tof));
469     Double_t tofcorr=tof;
470     if(timeWalkCorr)tofcorr=CorrectTimeWalk(dzTW,tof);
471     AliDebug(2,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));    
472     //Set TOF time signal and pointer to the matched cluster
473     t->SetTOFsignal(tofcorr);
474     t->SetTOFcluster(idclus); // pointing to the recPoints tree
475     t->SetTOFLabel(tlab);
476
477     AliDebug(2,Form(" Setting TOF raw time: %f  z distance: %f  corrected time: %f",rawTime,dzTW,tofcorr));
478
479     Double_t mom=t->GetP();
480     AliDebug(1,Form(" Momentum for track %d -> %f", iseed,mom));
481     // Fill Reco-QA histos for Reconstruction
482     fHRecNClus->Fill(nc);
483     fHRecChi2->Fill(bestChi2);
484     fHRecDistZ->Fill(dzTW);
485     fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
486     fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
487     fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
488     fHRecSigZVsPWin->Fill(mom,dz);
489
490     // Fill Tree for on-the-fly offline Calibration
491     // no longer there - all info is in the ESDs now
492
493   }
494
495 }
496 //_________________________________________________________________________
497 Int_t AliTOFtrackerV1::LoadClusters(TTree *cTree) {
498   //--------------------------------------------------------------------
499   //This function loads the TOF clusters
500   //--------------------------------------------------------------------
501
502   Int_t npadX = AliTOFGeometry::NpadX();
503   Int_t npadZ = AliTOFGeometry::NpadZ();
504   Int_t nStripA = AliTOFGeometry::NStripA();
505   Int_t nStripB = AliTOFGeometry::NStripB();
506   Int_t nStripC = AliTOFGeometry::NStripC();
507
508   TBranch *branch=cTree->GetBranch("TOF");
509   if (!branch) { 
510     AliError("can't get the branch with the TOF clusters !");
511     return 1;
512   }
513
514   static TClonesArray dummy("AliTOFcluster",10000);
515   dummy.Clear();
516   TClonesArray *clusters=&dummy;
517   branch->SetAddress(&clusters);
518
519   cTree->GetEvent(0);
520   Int_t nc=clusters->GetEntriesFast();
521   fHDigNClus->Fill(nc);
522
523   AliInfo(Form("Number of clusters: %d",nc));
524
525   for (Int_t i=0; i<nc; i++) {
526     AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
527     fClusters[i]=new AliTOFcluster(*c); fN++;
528
529   // Fill Digits QA histos
530  
531     Int_t isector = c->GetDetInd(0);
532     Int_t iplate = c->GetDetInd(1);
533     Int_t istrip = c->GetDetInd(2);
534     Int_t ipadX = c->GetDetInd(4);
535     Int_t ipadZ = c->GetDetInd(3);
536
537     Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
538     Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
539  
540     Int_t stripOffset = 0;
541     switch (iplate) {
542     case 0:
543       stripOffset = 0;
544       break;
545     case 1:
546       stripOffset = nStripC;
547       break;
548     case 2:
549       stripOffset = nStripC+nStripB;
550       break;
551     case 3:
552       stripOffset = nStripC+nStripB+nStripA;
553       break;
554     case 4:
555       stripOffset = nStripC+nStripB+nStripA+nStripB;
556       break;
557     default:
558       AliError(Form("Wrong plate number in TOF (%d) !",iplate));
559       break;
560     };
561     Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
562     Int_t phiindex=npadX*isector+ipadX+1;
563     fHDigClusMap->Fill(zindex,phiindex);
564     fHDigClusTime->Fill(time);
565     fHDigClusToT->Fill(tot);
566   }
567
568
569   return 0;
570 }
571 //_________________________________________________________________________
572 void AliTOFtrackerV1::UnloadClusters() {
573   //--------------------------------------------------------------------
574   //This function unloads TOF clusters
575   //--------------------------------------------------------------------
576   for (Int_t i=0; i<fN; i++) {
577     delete fClusters[i];
578     fClusters[i] = 0x0;
579   }
580   fN=0;
581 }
582
583 //_________________________________________________________________________
584 Int_t AliTOFtrackerV1::FindClusterIndex(Double_t z) const {
585   //--------------------------------------------------------------------
586   // This function returns the index of the nearest cluster 
587   //--------------------------------------------------------------------
588   //MOD
589   //Here we need to get the Z in the tracking system
590
591   if (fN==0) return 0;
592   if (z <= fClusters[0]->GetZ()) return 0;
593   if (z > fClusters[fN-1]->GetZ()) return fN;
594   Int_t b=0, e=fN-1, m=(b+e)/2;
595   for (; b<e; m=(b+e)/2) {
596     if (z > fClusters[m]->GetZ()) b=m+1;
597     else e=m; 
598   }
599   return m;
600 }
601
602 //_________________________________________________________________________
603 Bool_t AliTOFtrackerV1::GetTrackPoint(Int_t index, AliTrackPoint& p) const
604 {
605   // Get track space point with index i
606   // Coordinates are in the global system
607   AliTOFcluster *cl = fClusters[index];
608   Float_t xyz[3];
609   cl->GetGlobalXYZ(xyz);
610   Float_t phi=TMath::ATan2(xyz[1],xyz[0]);
611   Float_t phiangle = (Int_t(phi*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
612   Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
613   Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
614   Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
615   Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
616   Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
617   Float_t cov[6];
618   cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
619   cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
620   cov[2] = -cosphi*sinth*costh*sigmaz2;
621   cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
622   cov[4] = -sinphi*sinth*costh*sigmaz2;
623   cov[5] = costh*costh*sigmaz2;
624   p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
625
626   // Detector numbering scheme
627   Int_t nSector = AliTOFGeometry::NSectors();
628   Int_t nPlate  = AliTOFGeometry::NPlates();
629   Int_t nStripA = AliTOFGeometry::NStripA();
630   Int_t nStripB = AliTOFGeometry::NStripB();
631   Int_t nStripC = AliTOFGeometry::NStripC();
632
633   Int_t isector = cl->GetDetInd(0);
634   if (isector >= nSector)
635     AliError(Form("Wrong sector number in TOF (%d) !",isector));
636   Int_t iplate = cl->GetDetInd(1);
637   if (iplate >= nPlate)
638     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
639   Int_t istrip = cl->GetDetInd(2);
640
641   Int_t stripOffset = 0;
642   switch (iplate) {
643   case 0:
644     stripOffset = 0;
645     break;
646   case 1:
647     stripOffset = nStripC;
648     break;
649   case 2:
650     stripOffset = nStripC+nStripB;
651     break;
652   case 3:
653     stripOffset = nStripC+nStripB+nStripA;
654     break;
655   case 4:
656     stripOffset = nStripC+nStripB+nStripA+nStripB;
657     break;
658   default:
659     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
660     break;
661   };
662
663   Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
664                stripOffset +
665                istrip;
666   UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
667   p.SetVolumeID((UShort_t)volid);
668   return kTRUE;
669 }
670 //_________________________________________________________________________
671 void AliTOFtrackerV1::InitCheckHists() {
672
673   //Init histos for Digits/Reco QA and Calibration
674
675   TDirectory *dir = gDirectory;
676   TFile *logFileTOF = 0;
677
678   TSeqCollection *list = gROOT->GetListOfFiles();
679   int n = list->GetEntries();
680   Bool_t isThere=kFALSE;
681   for(int i=0; i<n; i++) {
682     logFileTOF = (TFile*)list->At(i);
683     if (strstr(logFileTOF->GetName(), "TOFQA.root")){
684       isThere=kTRUE;
685       break;
686     } 
687   }
688
689   if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
690   logFileTOF->cd(); 
691
692   //Digits "QA" 
693   fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);  
694   fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);  
695   fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);  
696   fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);  
697
698   //Reco "QA"
699   fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
700   fHRecDistZ=new TH1F("TOFRec_DistZ", "",50,0.5,10.5);
701   fHRecChi2=new TH1F("TOFRec_Chi2", "",100,0.,10.);
702   fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
703   fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
704   fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
705   fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
706
707   dir->cd();
708
709 }
710
711 //_________________________________________________________________________
712 void AliTOFtrackerV1::SaveCheckHists() {
713
714   //write histos for Digits/Reco QA and Calibration
715
716   TDirectory *dir = gDirectory;
717   TFile *logFile = 0;
718   TFile *logFileTOF = 0;
719
720   TSeqCollection *list = gROOT->GetListOfFiles();
721   int n = list->GetEntries();
722   for(int i=0; i<n; i++) {
723     logFile = (TFile*)list->At(i);
724     if (strstr(logFile->GetName(), "AliESDs.root")) break;
725   }
726
727   Bool_t isThere=kFALSE;
728   for(int i=0; i<n; i++) {
729     logFileTOF = (TFile*)list->At(i);
730     if (strstr(logFileTOF->GetName(), "TOFQA.root")){
731       isThere=kTRUE;
732       break;
733     } 
734   }
735    
736   if(!isThere) {
737           AliError(Form("File TOFQA.root not found!! not wring histograms...."));
738           return;
739   }
740   logFile->cd();
741   fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
742   fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
743   fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
744   fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
745   fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
746   fHRecChi2->Write(fHRecChi2->GetName(), TObject::kOverwrite);
747   fHRecDistZ->Write(fHRecDistZ->GetName(), TObject::kOverwrite);
748   fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
749   fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
750   fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
751   fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
752   logFile->Flush();  
753
754   dir->cd();
755
756   }
757 //_________________________________________________________________________
758 Float_t AliTOFtrackerV1::CorrectTimeWalk( Float_t dist, Float_t tof) {
759
760   //dummy, for the moment
761   Float_t tofcorr=0.;
762   if(dist<AliTOFGeometry::ZPad()*0.5){
763     tofcorr=tof;
764     //place here the actual correction
765   }else{
766     tofcorr=tof; 
767   } 
768   return tofcorr;
769 }
770 //_________________________________________________________________________
771 Float_t AliTOFtrackerV1::GetTimeZerofromT0(AliESDEvent *event) const {
772
773   //Returns TimeZero as measured by T0 detector
774
775   return event->GetT0();
776 }
777 //_________________________________________________________________________
778 Float_t AliTOFtrackerV1::GetTimeZerofromTOF(AliESDEvent * /*event*/) const {
779
780   //dummy, for the moment. T0 algorithm using tracks on TOF
781   {
782     //place T0 algo here...
783   }
784   return 0.;
785 }
786 //_________________________________________________________________________
787
788 void AliTOFtrackerV1::FillClusterArray(TObjArray* arr) const
789 {
790   //
791   // Returns the TOF cluster array
792   //
793
794   if (fN==0)
795     arr = 0x0;
796   else
797     for (Int_t i=0; i<fN; ++i) arr->Add(fClusters[i]);
798
799 }