]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFtracker.cxx
Add required AliVParticle functionality (Markus)
[u/mrichter/AliRoot.git] / TOF / AliTOFtracker.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 // AliTOFtracker 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 <TGeoManager.h>
34 #include <TTree.h>
35 #include <TFile.h>
36 #include <TH2F.h>
37
38 #include "AliAlignObj.h"
39 #include "AliESDtrack.h"
40 #include "AliESDEvent.h"
41 #include "AliLog.h"
42 #include "AliTrackPointArray.h"
43
44 #include "AliTOFcalib.h"
45 #include "AliTOFRecoParam.h"
46 #include "AliTOFcluster.h"
47 #include "AliTOFGeometryV5.h"
48 #include "AliTOFtracker.h"
49 #include "AliTOFtrack.h"
50
51 extern TGeoManager *gGeoManager;
52 extern TROOT *gROOT;
53
54 ClassImp(AliTOFtracker)
55
56 //_____________________________________________________________________________
57 AliTOFtracker::AliTOFtracker():
58   fRecoParam(0x0),
59   fGeom(0x0),
60   fPid(0x0),
61   fN(0),
62   fNseeds(0),
63   fNseedsTOF(0),
64   fngoodmatch(0),
65   fnbadmatch(0),
66   fnunmatch(0),
67   fnmatch(0),
68   fTracks(0x0),
69   fSeeds(0x0),
70   fHDigClusMap(0x0),
71   fHDigNClus(0x0),
72   fHDigClusTime(0x0),
73   fHDigClusToT(0x0),
74   fHRecNClus(0x0),
75   fHRecDist(0x0),
76   fHRecSigYVsP(0x0),
77   fHRecSigZVsP(0x0),
78   fHRecSigYVsPWin(0x0),
79   fHRecSigZVsPWin(0x0),
80   fCalTree(0x0),
81   fIch(-1),
82   fToT(-1.),
83   fTime(-1.),
84   fExpTimePi(-1.),
85   fExpTimeKa(-1.),
86   fExpTimePr(-1.)
87  { 
88   //AliTOFtracker main Ctor
89    
90    // Gettimg the geometry 
91    fGeom=new AliTOFGeometryV5();
92    // Read the reconstruction parameters from the OCDB
93    AliTOFcalib *calib = new AliTOFcalib(fGeom);
94    fRecoParam = (AliTOFRecoParam*)calib->ReadRecParFromCDB("TOF/Calib",-1);
95    if(fRecoParam->GetApplyPbPbCuts())fRecoParam=fRecoParam->GetPbPbparam();
96    Double_t parPID[2];   
97    parPID[0]=fRecoParam->GetTimeResolution();
98    parPID[1]=fRecoParam->GetTimeNSigma();
99    fPid=new AliTOFpidESD(parPID);
100    InitCheckHists();
101
102 }
103 //_____________________________________________________________________________
104 AliTOFtracker::AliTOFtracker(const AliTOFtracker &t):
105   AliTracker(),
106   fRecoParam(0x0),
107   fGeom(0x0),
108   fPid(0x0),
109   fN(0),
110   fNseeds(0),
111   fNseedsTOF(0),
112   fngoodmatch(0),
113   fnbadmatch(0),
114   fnunmatch(0),
115   fnmatch(0),
116   fTracks(0x0),
117   fSeeds(0x0),
118   fHDigClusMap(0x0),
119   fHDigNClus(0x0),
120   fHDigClusTime(0x0),
121   fHDigClusToT(0x0),
122   fHRecNClus(0x0),
123   fHRecDist(0x0),
124   fHRecSigYVsP(0x0),
125   fHRecSigZVsP(0x0),
126   fHRecSigYVsPWin(0x0),
127   fHRecSigZVsPWin(0x0),
128   fCalTree(0x0),
129   fIch(-1),
130   fToT(-1.),
131   fTime(-1.),
132   fExpTimePi(-1.),
133   fExpTimeKa(-1.),
134   fExpTimePr(-1.)
135  { 
136   //AliTOFtracker copy Ctor
137
138   fNseeds=t.fNseeds;
139   fNseeds=t.fNseeds;
140   fNseedsTOF=t.fNseedsTOF;
141   fngoodmatch=t.fngoodmatch;
142   fnbadmatch=t.fnbadmatch;
143   fnunmatch=t.fnunmatch;
144   fnmatch=t.fnmatch;
145   fRecoParam=t.fRecoParam;
146   fGeom=t.fGeom;
147   fPid=t.fPid;
148   fSeeds=t.fSeeds;
149   fTracks=t.fTracks;
150   fN=t.fN;
151 }
152
153 //_____________________________________________________________________________
154 AliTOFtracker& AliTOFtracker::operator=(const AliTOFtracker &t)
155
156   //AliTOFtracker assignment operator
157
158   this->fNseeds=t.fNseeds;
159   this->fNseedsTOF=t.fNseedsTOF;
160   this->fngoodmatch=t.fngoodmatch;
161   this->fnbadmatch=t.fnbadmatch;
162   this->fnunmatch=t.fnunmatch;
163   this->fnmatch=t.fnmatch;
164   this->fRecoParam = t.fRecoParam;
165   this->fGeom = t.fGeom;
166   this->fPid = t.fPid;
167   this->fSeeds=t.fSeeds;
168   this->fTracks=t.fTracks;
169   this->fN=t.fN;
170   return *this;
171
172 }
173 //_____________________________________________________________________________
174 AliTOFtracker::~AliTOFtracker() {
175   //
176   // Dtor
177   //
178
179   SaveCheckHists();
180
181   delete fRecoParam; 
182   delete fGeom; 
183   delete fPid; 
184   delete fHDigClusMap;
185   delete fHDigNClus;
186   delete fHDigClusTime;
187   delete fHDigClusToT;
188   delete fHRecNClus;
189   delete fHRecDist;
190   delete fHRecSigYVsP;
191   delete fHRecSigZVsP;
192   delete fHRecSigYVsPWin;
193   delete fHRecSigZVsPWin;
194   delete fCalTree;
195 }
196 //_____________________________________________________________________________
197 Int_t AliTOFtracker::PropagateBack(AliESDEvent* event) {
198   //
199   // Gets seeds from ESD event and Match with TOF Clusters
200   //
201
202
203   //Initialise some counters
204
205   fNseeds=0;
206   fNseedsTOF=0;
207   fngoodmatch=0;
208   fnbadmatch=0;
209   fnunmatch=0;
210   fnmatch=0;
211
212   Int_t ntrk=event->GetNumberOfTracks();
213   fNseeds = ntrk;
214   fSeeds= new TClonesArray("AliESDtrack",ntrk);
215   TClonesArray &aESDTrack = *fSeeds;
216
217
218   //Load ESD tracks into a local Array of ESD Seeds
219
220   for (Int_t i=0; i<fNseeds; i++) {
221     AliESDtrack *t=event->GetTrack(i);
222     new(aESDTrack[i]) AliESDtrack(*t);
223   }
224
225   //Prepare ESD tracks candidates for TOF Matching
226   CollectESD();
227
228   //First Step with Strict Matching Criterion
229   MatchTracks(kFALSE);
230
231   //Second Step with Looser Matching Criterion
232   MatchTracks(kTRUE);
233
234   AliInfo(Form("Number of matched tracks: %d",fnmatch));
235   AliInfo(Form("Number of good matched tracks: %d",fngoodmatch));
236   AliInfo(Form("Number of bad  matched tracks: %d",fnbadmatch));
237
238   //Update the matched ESD tracks
239
240   for (Int_t i=0; i<ntrk; i++) {
241     AliESDtrack *t=event->GetTrack(i);
242     AliESDtrack *seed =(AliESDtrack*)fSeeds->UncheckedAt(i);
243     if(seed->GetTOFsignal()>0){
244       Float_t info[10];
245       t->SetTOFsignal(seed->GetTOFsignal());
246       t->SetTOFcluster(seed->GetTOFcluster());
247       t->SetTOFsignalToT(seed->GetTOFsignalToT());
248       t->SetTOFsignalRaw(seed->GetTOFsignalRaw());
249       t->SetTOFsignalDz(seed->GetTOFsignalDz());
250       t->SetTOFCalChannel(seed->GetTOFCalChannel());
251       Int_t tlab[3]; seed->GetTOFLabel(tlab);    
252       t->SetTOFLabel(tlab);
253       AliTOFtrack *track = new AliTOFtrack(*seed); 
254       t->UpdateTrackParams(track,AliESDtrack::kTOFout);   
255       delete track;
256     }
257   }
258
259   //Handle Time Zero information
260
261   Double_t timeZero=0.;
262   Double_t timeZeroMax=99999.;
263   Bool_t usetimeZero     = fRecoParam->UseTimeZero();
264   Bool_t timeZeroFromT0  = fRecoParam->GetTimeZerofromT0();
265   Bool_t timeZeroFromTOF = fRecoParam->GetTimeZerofromTOF();
266
267   AliDebug(1,Form("Use Time Zero?: %d",usetimeZero));
268   AliDebug(1,Form("Time Zero from T0? : %d",timeZeroFromT0));
269   AliDebug(1,Form("Time Zero From TOF? : %d",timeZeroFromTOF));
270
271   if(usetimeZero){
272     if(timeZeroFromT0){
273       timeZero=GetTimeZerofromT0(event); 
274     }
275     if(timeZeroFromTOF && (timeZero>timeZeroMax || !timeZeroFromT0)){
276       timeZero=GetTimeZerofromTOF(event); 
277     }
278   }
279   AliDebug(2,Form("time Zero used in PID: %f",timeZero));
280   //Make TOF PID
281   fPid->MakePID(event,timeZero);
282
283   if (fSeeds) {
284     fSeeds->Delete();
285     delete fSeeds;
286     fSeeds = 0x0;
287   }
288   if (fTracks) {
289     fTracks->Delete();
290     delete fTracks;
291     fTracks = 0x0;
292   }
293   return 0;
294   
295 }
296 //_________________________________________________________________________
297 void AliTOFtracker::CollectESD() {
298    //prepare the set of ESD tracks to be matched to clusters in TOF
299
300   Int_t seedsTOF1=0;
301   Int_t seedsTOF2=0;
302  
303   fTracks= new TClonesArray("AliTOFtrack");
304   TClonesArray &aTOFTrack = *fTracks;
305   for (Int_t i=0; i<fNseeds; i++) {
306
307     AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(i);
308     if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
309
310     // TRD 'good' tracks, already propagated at 371 cm
311
312     AliTOFtrack *track = new AliTOFtrack(*t); // New
313     Double_t x = track->GetX(); //New
314
315     if (((t->GetStatus()&AliESDtrack::kTRDout)!=0 ) && 
316          ( x >= fGeom->RinTOF()) ){
317       track->SetSeedIndex(i);
318       t->UpdateTrackParams(track,AliESDtrack::kTOFout);    
319       new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
320       fNseedsTOF++;
321       seedsTOF1++;
322       delete track;
323     }
324
325     // Propagate the rest of TPCbp  
326
327     else {
328       if(track->PropagateToInnerTOF()){ 
329         track->SetSeedIndex(i);
330         t->UpdateTrackParams(track,AliESDtrack::kTOFout);    
331         new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
332         fNseedsTOF++;
333         seedsTOF2++;
334       }
335       delete track;
336     }
337   }
338
339   AliInfo(Form("Number of TOF seeds %i",fNseedsTOF));
340   AliInfo(Form("Number of TOF seeds Type 1 %i",seedsTOF1));
341   AliInfo(Form("Number of TOF seeds Type 2 %i",seedsTOF2));
342
343   // Sort according uncertainties on track position 
344   fTracks->Sort();
345
346 }
347 //_________________________________________________________________________
348 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
349
350
351   // Parameters used/regulating the reconstruction
352
353   static Float_t corrLen=0.75;
354   static Float_t detDepth=15.3;
355   static Float_t padDepth=0.5;
356
357   Float_t dY=fGeom->XPad(); 
358   Float_t dZ=fGeom->ZPad(); 
359
360   Float_t sensRadius = fRecoParam->GetSensRadius();
361   Float_t stepSize   = fRecoParam->GetStepSize();
362   Float_t scaleFact   = fRecoParam->GetWindowScaleFact();
363   Float_t dyMax=fRecoParam->GetWindowSizeMaxY(); 
364   Float_t dzMax=fRecoParam->GetWindowSizeMaxZ();
365   Float_t dCut=fRecoParam->GetDistanceCut();
366   Double_t maxChi2=fRecoParam->GetMaxChi2();
367   Bool_t timeWalkCorr    = fRecoParam->GetTimeWalkCorr();
368   if(!mLastStep){
369     AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++ \n");
370     AliDebug(1,Form("TOF sens radius: %f",sensRadius));
371     AliDebug(1,Form("TOF step size: %f",stepSize));
372     AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
373     AliDebug(1,Form("TOF Window max dy: %f",dyMax));
374     AliDebug(1,Form("TOF Window max dz: %f",dzMax));
375     AliDebug(1,Form("TOF distance Cut: %f",dCut));
376     AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
377     AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));   
378   }
379   //Match ESD tracks to clusters in TOF
380
381
382   // Get the number of propagation steps
383
384   Int_t nSteps=(Int_t)(detDepth/stepSize);
385
386   AliTOFcalib *calib = new AliTOFcalib(fGeom);
387
388   //PH Arrays (moved outside of the loop)
389
390   Float_t * trackPos[4];
391   for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
392   Int_t * clind[6];
393   for (Int_t ii=0;ii<6;ii++) clind[ii] = new Int_t[fN];
394   
395   for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
396
397     AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
398     AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(track->GetSeedIndex());
399     if(t->GetTOFsignal()>0. ) continue;
400     AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
401
402     // Some init
403
404     Int_t         index[10000];
405     Float_t        dist[10000];
406     Float_t       distZ[10000];
407     Float_t       cxpos[10000];
408     Float_t       crecL[10000];
409     TGeoHMatrix   global[1000];
410      
411     // Determine a window around the track
412
413     Double_t x,par[5]; 
414     trackTOFin->GetExternalParameters(x,par);
415     Double_t cov[15]; 
416     trackTOFin->GetExternalCovariance(cov);
417
418     Double_t dphi=
419       scaleFact*
420       ((5*TMath::Sqrt(cov[0]) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius); 
421     Double_t dz=
422        scaleFact*
423        (5*TMath::Sqrt(cov[2]) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
424
425     Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
426     if (phi<-TMath::Pi())phi+=2*TMath::Pi();
427     if (phi>=TMath::Pi())phi-=2*TMath::Pi();
428     Double_t z=par[1];   
429
430     //upper limit on window's size.
431
432     if(dz> dzMax) dz=dzMax;
433     if(dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
434
435
436     Int_t nc=0;
437
438     // find the clusters in the window of the track
439
440     for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
441
442       AliTOFcluster *c=fClusters[k];
443       if (c->GetZ() > z+dz) break;
444       if (c->IsUsed()) continue;
445
446       if (!c->GetStatus()) continue; // skip bad channels as declared in OCDB
447
448       Double_t dph=TMath::Abs(c->GetPhi()-phi);
449       if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
450       if (TMath::Abs(dph)>dphi) continue;
451
452       Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
453       Double_t p[2]={yc, c->GetZ()};
454       Double_t cov[3]= {dY*dY/12., 0., dZ*dZ/12.};
455       if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov) > maxChi2)continue;
456
457       clind[0][nc] = c->GetDetInd(0);
458       clind[1][nc] = c->GetDetInd(1);
459       clind[2][nc] = c->GetDetInd(2);
460       clind[3][nc] = c->GetDetInd(3);
461       clind[4][nc] = c->GetDetInd(4);
462       clind[5][nc] = k;      
463       Char_t path[100];
464       Int_t ind[5];
465       ind[0]=clind[0][nc];
466       ind[1]=clind[1][nc];
467       ind[2]=clind[2][nc];
468       ind[3]=clind[3][nc];
469       ind[4]=clind[4][nc];
470       fGeom->GetVolumePath(ind,path);
471       gGeoManager->cd(path);
472       global[nc] = *gGeoManager->GetCurrentMatrix();
473       nc++;
474     }
475
476     //start fine propagation 
477
478     Int_t nStepsDone = 0;
479     for( Int_t istep=0; istep<nSteps; istep++){ 
480
481       Float_t xs=fGeom->RinTOF()+istep*0.1;
482       Double_t ymax=xs*TMath::Tan(0.5*fGeom->GetAlpha());
483
484       Bool_t skip=kFALSE;
485       Double_t ysect=trackTOFin->GetYat(xs,skip);
486       if (skip) break;
487       if (ysect > ymax) {
488         if (!trackTOFin->Rotate(fGeom->GetAlpha())) {
489           break;
490         }
491       } else if (ysect <-ymax) {
492         if (!trackTOFin->Rotate(fGeom->GetAlpha())) {
493           break;
494         }
495       }
496
497       if(!trackTOFin->PropagateTo(xs)) {
498         break;
499       }
500
501       nStepsDone++;
502
503       // store the running point (Globalrf) - fine propagation     
504
505       Double_t r[3];
506       trackTOFin->GetXYZ(r);
507       trackPos[0][istep]= (Float_t) r[0];
508       trackPos[1][istep]= (Float_t) r[1];
509       trackPos[2][istep]= (Float_t) r[2];   
510       trackPos[3][istep]= trackTOFin->GetIntegratedLength();
511     }
512
513
514     Int_t nfound = 0;
515     Bool_t accept = kFALSE;
516     Bool_t isInside =kFALSE;
517     for (Int_t istep=0; istep<nStepsDone; istep++) {
518
519       Float_t ctrackPos[3];     
520
521       ctrackPos[0]= trackPos[0][istep];
522       ctrackPos[1]= trackPos[1][istep];
523       ctrackPos[2]= trackPos[2][istep];
524
525       //now see whether the track matches any of the TOF clusters            
526
527       Float_t dist3d[3];
528       accept=kFALSE;
529       for (Int_t i=0; i<nc; i++){
530         Int_t cind[5];
531         cind[0]= clind[0][i];
532         cind[1]= clind[1][i];
533         cind[2]= clind[2][i];
534         cind[3]= clind[3][i];
535         cind[4]= clind[4][i];
536         isInside=fGeom->IsInsideThePad(global[i],ctrackPos,dist3d);
537
538         if( mLastStep){
539           Float_t xLoc=dist3d[0];
540           Float_t rLoc=TMath::Sqrt(dist3d[1]*dist3d[1]+dist3d[2]*dist3d[2]);
541           accept = (TMath::Abs(xLoc)<padDepth*0.5 && rLoc<dCut);
542         }
543         else{
544           accept = isInside;
545         }
546         if(accept){
547           dist[nfound]=TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[1]*dist3d[1]+dist3d[2]*dist3d[2]);
548           distZ[nfound]=dist3d[2];
549           crecL[nfound]=trackPos[3][istep];
550           index[nfound]=clind[5][i]; // store cluster id            
551           cxpos[nfound]=fGeom->RinTOF()+istep*0.1; //store prop.radius
552           nfound++;
553           if(accept &&!mLastStep)break;
554         }//end if accept
555       } //end for on the clusters
556       if(accept &&!mLastStep)break;
557     } //end for on the steps     
558
559
560
561     if (nfound == 0 ) {
562       fnunmatch++;
563       delete trackTOFin;
564       continue;
565     }
566     
567     fnmatch++;
568
569     // now choose the cluster to be matched with the track.
570
571     Int_t idclus=0;
572     Float_t  recL = 0.;
573     Float_t  xpos=0.;
574     Float_t  mindist=1000.;
575     Float_t  mindistZ=0.;
576     for (Int_t iclus= 0; iclus<nfound;iclus++){
577       if (dist[iclus]< mindist){
578         mindist = dist[iclus];
579         mindistZ = distZ[iclus];
580         xpos = cxpos[iclus];
581         idclus =index[iclus]; 
582         recL=crecL[iclus]+corrLen*0.5;
583       }
584     }
585
586
587     AliTOFcluster *c=fClusters[idclus];
588     c->Use(); 
589
590     // Track length correction for matching Step 2 
591
592     if(mLastStep){
593       Float_t rc=TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
594       Float_t rt=TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
595                              +trackPos[1][70]*trackPos[1][70]
596                              +trackPos[2][70]*trackPos[2][70]);
597       Float_t dlt=rc-rt;      
598       recL=trackPos[3][70]+dlt;
599     }    
600
601     if (
602         (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
603         ||
604         (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
605         ||
606         (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
607         ) {
608       fngoodmatch++;
609
610        AliDebug(2,Form(" track label good %5i",trackTOFin->GetLabel()));
611
612     }
613     else{
614       fnbadmatch++;
615
616       AliDebug(2,Form(" track label  bad %5i",trackTOFin->GetLabel()));
617
618     }
619
620     delete trackTOFin;
621
622     //  Store quantities to be used in the TOF Calibration
623     Float_t tToT=fGeom->ToTBinWidth()*c->GetToT()*1E-3; // in ns
624     t->SetTOFsignalToT(tToT);
625     Float_t rawTime=fGeom->TdcBinWidth()*c->GetTDCRAW()+32; // RAW time,in ps
626     t->SetTOFsignalRaw(rawTime);
627     t->SetTOFsignalDz(mindistZ);
628     AliDebug(2,Form(" Setting TOF raw time: %f  z distance: %f time: %f = ",rawTime,mindistZ));    
629     Int_t ind[5];
630     ind[0]=c->GetDetInd(0);
631     ind[1]=c->GetDetInd(1);
632     ind[2]=c->GetDetInd(2);
633     ind[3]=c->GetDetInd(3);
634     ind[4]=c->GetDetInd(4);
635     Int_t calindex = calib->GetIndex(ind);
636     t->SetTOFCalChannel(calindex);
637
638     // keep track of the track labels in the matched cluster
639     Int_t tlab[3];
640     tlab[0]=c->GetLabel(0);
641     tlab[1]=c->GetLabel(1);
642     tlab[2]=c->GetLabel(2);
643     AliDebug(2,Form(" tdc time of the matched track %i = ",c->GetTDC()));    
644     Double_t tof=fGeom->TdcBinWidth()*c->GetTDC()+32; // in ps
645     AliDebug(2,Form(" tof time of the matched track: %f = ",tof));    
646     Double_t tofcorr=tof;
647     if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
648     AliDebug(2,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));    
649     //Set TOF time signal and pointer to the matched cluster
650     t->SetTOFsignal(tofcorr);
651     t->SetTOFcluster(idclus); // pointing to the recPoints tree
652
653     //Tracking info
654     Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time);
655     Double_t mom=t->GetP();
656     for(Int_t j=0;j<AliPID::kSPECIES;j++){
657       Double_t mass=AliPID::ParticleMass(j);
658       time[j]+=(recL-trackPos[3][0])/3e-2*TMath::Sqrt(mom*mom+mass*mass)/mom;
659     }
660
661     AliTOFtrack *trackTOFout = new AliTOFtrack(*t); 
662     trackTOFout->PropagateTo(xpos);
663     t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);    
664     t->SetIntegratedLength(recL);
665     t->SetIntegratedTimes(time);
666     t->SetTOFLabel(tlab);
667
668  
669     // Fill Reco-QA histos for Reconstruction
670     fHRecNClus->Fill(nc);
671     fHRecDist->Fill(mindist);
672     fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
673     fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
674     fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
675     fHRecSigZVsPWin->Fill(mom,dz);
676
677     // Fill Tree for on-the-fly offline Calibration
678
679     if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 )){    
680       fIch=calindex;
681       fToT=tToT;
682       fTime=rawTime;
683       fExpTimePi=time[2];
684       fExpTimeKa=time[3];
685       fExpTimePr=time[4];
686       fCalTree->Fill();
687     }
688     delete trackTOFout;
689   }
690   for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
691   for (Int_t ii=0;ii<6;ii++) delete [] clind[ii];
692   delete calib;
693 }
694 //_________________________________________________________________________
695 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
696   //--------------------------------------------------------------------
697   //This function loads the TOF clusters
698   //--------------------------------------------------------------------
699
700   Int_t npadX = fGeom->NpadX();
701   Int_t npadZ = fGeom->NpadZ();
702   Int_t nStripA = fGeom->NStripA();
703   Int_t nStripB = fGeom->NStripB();
704   Int_t nStripC = fGeom->NStripC();
705
706   TBranch *branch=cTree->GetBranch("TOF");
707   if (!branch) { 
708     AliError("can't get the branch with the TOF clusters !");
709     return 1;
710   }
711
712   TClonesArray dummy("AliTOFcluster",10000), *clusters=&dummy;
713   branch->SetAddress(&clusters);
714
715   cTree->GetEvent(0);
716   Int_t nc=clusters->GetEntriesFast();
717   fHDigNClus->Fill(nc);
718
719   AliInfo(Form("Number of clusters: %d",nc));
720
721   for (Int_t i=0; i<nc; i++) {
722     AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
723     fClusters[i]=new AliTOFcluster(*c); fN++;
724
725   // Fill Digits QA histos
726  
727     Int_t isector = c->GetDetInd(0);
728     Int_t iplate = c->GetDetInd(1);
729     Int_t istrip = c->GetDetInd(2);
730     Int_t ipadX = c->GetDetInd(4);
731     Int_t ipadZ = c->GetDetInd(3);
732
733     Float_t time =(fGeom->TdcBinWidth()*c->GetTDC())*1E-3; // in ns
734     Float_t tot = (fGeom->TdcBinWidth()*c->GetToT())*1E-3;//in ns
735  
736     Int_t stripOffset = 0;
737     switch (iplate) {
738     case 0:
739       stripOffset = 0;
740       break;
741     case 1:
742       stripOffset = nStripC;
743       break;
744     case 2:
745       stripOffset = nStripC+nStripB;
746       break;
747     case 3:
748       stripOffset = nStripC+nStripB+nStripA;
749       break;
750     case 4:
751       stripOffset = nStripC+nStripB+nStripA+nStripB;
752       break;
753     default:
754       AliError(Form("Wrong plate number in TOF (%d) !",iplate));
755       break;
756     };
757     Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
758     Int_t phiindex=npadX*isector+ipadX+1;
759     fHDigClusMap->Fill(zindex,phiindex);
760     fHDigClusTime->Fill(time);
761     fHDigClusToT->Fill(tot);
762
763   }
764
765
766   return 0;
767 }
768 //_________________________________________________________________________
769 void AliTOFtracker::UnloadClusters() {
770   //--------------------------------------------------------------------
771   //This function unloads TOF clusters
772   //--------------------------------------------------------------------
773   for (Int_t i=0; i<fN; i++) {
774     delete fClusters[i];
775     fClusters[i] = 0x0;
776   }
777   fN=0;
778 }
779
780 //_________________________________________________________________________
781 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
782   //--------------------------------------------------------------------
783   // This function returns the index of the nearest cluster 
784   //--------------------------------------------------------------------
785   if (fN==0) return 0;
786   if (z <= fClusters[0]->GetZ()) return 0;
787   if (z > fClusters[fN-1]->GetZ()) return fN;
788   Int_t b=0, e=fN-1, m=(b+e)/2;
789   for (; b<e; m=(b+e)/2) {
790     if (z > fClusters[m]->GetZ()) b=m+1;
791     else e=m; 
792   }
793   return m;
794 }
795
796 //_________________________________________________________________________
797 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
798 {
799   // Get track space point with index i
800   // Coordinates are in the global system
801   AliTOFcluster *cl = fClusters[index];
802   Float_t xyz[3];
803   xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
804   xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
805   xyz[2] = cl->GetZ();
806   Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
807   Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
808   Float_t tiltangle = fGeom->GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
809   Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
810   Float_t sigmay2 = fGeom->XPad()*fGeom->XPad()/12.;
811   Float_t sigmaz2 = fGeom->ZPad()*fGeom->ZPad()/12.;
812   Float_t cov[6];
813   cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
814   cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
815   cov[2] = -cosphi*sinth*costh*sigmaz2;
816   cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
817   cov[4] = -sinphi*sinth*costh*sigmaz2;
818   cov[5] = costh*costh*sigmaz2;
819   p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
820
821   // Detector numbering scheme
822   Int_t nSector = fGeom->NSectors();
823   Int_t nPlate  = fGeom->NPlates();
824   Int_t nStripA = fGeom->NStripA();
825   Int_t nStripB = fGeom->NStripB();
826   Int_t nStripC = fGeom->NStripC();
827
828   Int_t isector = cl->GetDetInd(0);
829   if (isector >= nSector)
830     AliError(Form("Wrong sector number in TOF (%d) !",isector));
831   Int_t iplate = cl->GetDetInd(1);
832   if (iplate >= nPlate)
833     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
834   Int_t istrip = cl->GetDetInd(2);
835
836   Int_t stripOffset = 0;
837   switch (iplate) {
838   case 0:
839     stripOffset = 0;
840     break;
841   case 1:
842     stripOffset = nStripC;
843     break;
844   case 2:
845     stripOffset = nStripC+nStripB;
846     break;
847   case 3:
848     stripOffset = nStripC+nStripB+nStripA;
849     break;
850   case 4:
851     stripOffset = nStripC+nStripB+nStripA+nStripB;
852     break;
853   default:
854     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
855     break;
856   };
857
858   Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
859                stripOffset +
860                istrip;
861   UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
862   p.SetVolumeID((UShort_t)volid);
863   return kTRUE;
864 }
865 //_________________________________________________________________________
866 void AliTOFtracker::InitCheckHists() {
867
868   //Init histos for Digits/Reco QA and Calibration
869
870
871   fCalTree = new TTree("CalTree", "Tree for TOF calibration");
872   fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
873   fCalTree->Branch("ToT",&fToT,"TOFToT/F");
874   fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
875   fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
876   fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
877   fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
878
879   //Digits "QA" 
880   fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);  
881   fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);  
882   fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);  
883   fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);  
884
885   //Reco "QA"
886   fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
887   fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
888   fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
889   fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
890   fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
891   fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
892 }
893
894 //_________________________________________________________________________
895 void AliTOFtracker::SaveCheckHists() {
896
897   //write histos for Digits/Reco QA and Calibration
898
899   TDirectory *dir = gDirectory;
900   TFile *logFile = 0;
901   TFile *logFileTOF = 0;
902
903   TSeqCollection *list = gROOT->GetListOfFiles();
904   int n = list->GetEntries();
905   for(int i=0; i<n; i++) {
906     logFile = (TFile*)list->At(i);
907     if (strstr(logFile->GetName(), "AliESDs.root")) break;
908   }
909
910   Bool_t isThere=kFALSE;
911   for(int i=0; i<n; i++) {
912     logFileTOF = (TFile*)list->At(i);
913     if (strstr(logFileTOF->GetName(), "TOFQA.root")){
914       isThere=kTRUE;
915       break;
916     } 
917   }
918    
919   logFile->cd();
920   fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
921   fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
922   fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
923   fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
924   fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
925   fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
926   fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
927   fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
928   fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
929   fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
930   //fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
931   logFile->Flush();  
932
933   if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
934   logFileTOF->cd(); 
935   fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
936   fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
937   fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
938   fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
939   fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
940   fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
941   fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
942   fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
943   fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
944   fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
945   fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
946   logFileTOF->Flush();  
947
948   dir->cd();
949   }
950 //_________________________________________________________________________
951 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) {
952
953   //dummy, for the moment
954   Float_t tofcorr=0.;
955   if(dist<fGeom->ZPad()*0.5){
956     tofcorr=tof;
957     //place here the actual correction
958   }else{
959     tofcorr=tof; 
960   } 
961   return tofcorr;
962 }
963 //_________________________________________________________________________
964 Float_t AliTOFtracker::GetTimeZerofromT0(AliESDEvent *event) const {
965
966   //Returns TimeZero as measured by T0 detector
967
968   return event->GetT0();
969 }
970 //_________________________________________________________________________
971 Float_t AliTOFtracker::GetTimeZerofromTOF(AliESDEvent * /*event*/) const {
972
973   //dummy, for the moment. T0 algorithm using tracks on TOF
974   {
975     //place T0 algo here...
976   }
977   return 0.;
978 }