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