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