]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFtracker.cxx
Coding conventions, simplifications
[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
199     if ( (seed->GetStatus()&AliESDtrack::kTOFin)!=0 ) {
200       t->SetStatus(AliESDtrack::kTOFin);
201       //if(seed->GetTOFsignal()>0){
202       if ( (seed->GetStatus()&AliESDtrack::kTOFout)!=0 ) {
203         t->SetStatus(AliESDtrack::kTOFout);
204         t->SetTOFsignal(seed->GetTOFsignal());
205         t->SetTOFcluster(seed->GetTOFcluster());
206         t->SetTOFsignalToT(seed->GetTOFsignalToT());
207         t->SetTOFsignalRaw(seed->GetTOFsignalRaw());
208         t->SetTOFsignalDz(seed->GetTOFsignalDz());
209         t->SetTOFCalChannel(seed->GetTOFCalChannel());
210         Int_t tlab[3]; seed->GetTOFLabel(tlab);    
211         t->SetTOFLabel(tlab);
212
213         AliTOFtrack *track = new AliTOFtrack(*seed);
214         t->UpdateTrackParams(track,AliESDtrack::kTOFout); // to be checked - AdC
215         delete track;
216         Double_t time[10]; t->GetIntegratedTimes(time);
217         AliDebug(2,Form(" %6d  %f %f %f %f %6d %3d %f  %f %f %f %f %f",
218                         i,
219                         t->GetTOFsignalRaw(),
220                         t->GetTOFsignal(),
221                         t->GetTOFsignalToT(),
222                         t->GetTOFsignalDz(),
223                         t->GetTOFCalChannel(),
224                         t->GetTOFcluster(),
225                         t->GetIntegratedLength(),
226                         time[0], time[1], time[2], time[3], time[4]
227                         )
228                  );
229       }
230     }
231   }
232
233   //Handle Time Zero information
234
235   Double_t timeZero=0.;
236   Double_t timeZeroMax=99999.;
237   Bool_t usetimeZero     = fRecoParam->UseTimeZero();
238   Bool_t timeZeroFromT0  = fRecoParam->GetTimeZerofromT0();
239   Bool_t timeZeroFromTOF = fRecoParam->GetTimeZerofromTOF();
240
241   AliDebug(2,Form("Use Time Zero?: %d",usetimeZero));
242   AliDebug(2,Form("Time Zero from T0? : %d",timeZeroFromT0));
243   AliDebug(2,Form("Time Zero From TOF? : %d",timeZeroFromTOF));
244
245   if(usetimeZero){
246     if(timeZeroFromT0){
247       timeZero=GetTimeZerofromT0(event); 
248     }
249     if(timeZeroFromTOF && (timeZero>timeZeroMax || !timeZeroFromT0)){
250       timeZero=GetTimeZerofromTOF(event); 
251     }
252   }
253   AliDebug(2,Form("time Zero used in PID: %f",timeZero));
254   //Make TOF PID
255   fPid->MakePID(event,timeZero);
256
257   fSeeds->Clear();
258   fTracks->Clear();
259   return 0;
260   
261 }
262 //_________________________________________________________________________
263 void AliTOFtracker::CollectESD() {
264    //prepare the set of ESD tracks to be matched to clusters in TOF
265
266   Int_t seedsTOF1=0;
267   Int_t seedsTOF2=0;
268  
269   TClonesArray &aTOFTrack = *fTracks;
270   for (Int_t i=0; i<fNseeds; i++) {
271
272     AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(i);
273     if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
274
275     AliTOFtrack *track = new AliTOFtrack(*t); // New
276     Float_t x = (Float_t)track->GetX(); //New
277
278     // TRD 'good' tracks, already propagated at 371 cm
279     if ( ( (t->GetStatus()&AliESDtrack::kTRDout)!=0 ) &&
280          ( x >= AliTOFGeometry::Rmin() ) ) {
281       if  ( track->PropagateToInnerTOF() ) {
282
283         AliDebug(1,Form(" TRD propagated track till rho = %fcm."
284                         " And then the track has been propagated till rho = %fcm.",
285                         x, (Float_t)track->GetX()));
286
287         track->SetSeedIndex(i);
288         t->UpdateTrackParams(track,AliESDtrack::kTOFin);
289         new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
290         fNseedsTOF++;
291         seedsTOF1++;
292       }
293       delete track;
294     }
295
296     // Propagate the rest of TPCbp  
297     else {
298       if ( track->PropagateToInnerTOF() ) { 
299
300         AliDebug(1,Form(" TPC propagated track till rho = %fcm."
301                         " And then the track has been propagated till rho = %fcm.",
302                         x, (Float_t)track->GetX()));
303
304         track->SetSeedIndex(i);
305         t->UpdateTrackParams(track,AliESDtrack::kTOFin);
306         new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
307         fNseedsTOF++;
308         seedsTOF2++;
309       }
310       delete track;
311     }
312   }
313
314   AliInfo(Form("Number of TOF seeds %d",fNseedsTOF));
315   AliInfo(Form("Number of TOF seeds Type 1 %d",seedsTOF1));
316   AliInfo(Form("Number of TOF seeds Type 2 %d",seedsTOF2));
317
318   // Sort according uncertainties on track position 
319   fTracks->Sort();
320
321 }
322 //_________________________________________________________________________
323 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
324
325
326   // Parameters used/regulating the reconstruction
327
328   static Float_t corrLen=0.75;
329   static Float_t detDepth=15.3;
330   static Float_t padDepth=0.5;
331
332   const Float_t kSpeedOfLight= 2.99792458e-2; // speed of light [cm/ps]
333   const Float_t kTimeOffset = 32.; // time offset for tracking algorithm [ps]
334
335   Float_t dY=AliTOFGeometry::XPad(); 
336   Float_t dZ=AliTOFGeometry::ZPad(); 
337
338   Float_t sensRadius = fRecoParam->GetSensRadius();
339   Float_t stepSize   = fRecoParam->GetStepSize();
340   Float_t scaleFact   = fRecoParam->GetWindowScaleFact();
341   Float_t dyMax=fRecoParam->GetWindowSizeMaxY(); 
342   Float_t dzMax=fRecoParam->GetWindowSizeMaxZ();
343   Float_t dCut=fRecoParam->GetDistanceCut();
344   Double_t maxChi2=fRecoParam->GetMaxChi2TRD();
345   Bool_t timeWalkCorr    = fRecoParam->GetTimeWalkCorr();
346   if(!mLastStep){
347     AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++ \n");
348     AliDebug(1,Form("TOF sens radius: %f",sensRadius));
349     AliDebug(1,Form("TOF step size: %f",stepSize));
350     AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
351     AliDebug(1,Form("TOF Window max dy: %f",dyMax));
352     AliDebug(1,Form("TOF Window max dz: %f",dzMax));
353     AliDebug(1,Form("TOF distance Cut: %f",dCut));
354     AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
355     AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));   
356   }
357
358   //Match ESD tracks to clusters in TOF
359
360   // Get the number of propagation steps
361   Int_t nSteps=(Int_t)(detDepth/stepSize);
362
363   //PH Arrays (moved outside of the loop)
364   Float_t * trackPos[4];
365   for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
366   Int_t * clind = new Int_t[fN];
367   
368   // Some init
369   const Int_t kNfoundMax = 10000; // related to nSteps value
370   Int_t         index[kNfoundMax];
371   Float_t        dist[kNfoundMax];
372   Float_t       distZ[kNfoundMax];
373   Float_t       cxpos[kNfoundMax];
374   Float_t       crecL[kNfoundMax];
375   const Int_t kNclusterMax = 1000; // related to fN value
376   TGeoHMatrix   global[kNclusterMax];
377      
378   //The matching loop
379   for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
380
381     for (Int_t ii=0; ii<kNfoundMax; ii++) {
382       index[ii] = -1;
383       dist[ii] = 9999.;
384       distZ[ii] = 9999.;
385       cxpos[ii] = 9999.;
386       crecL[ii] = 0.;
387     }
388     for (Int_t ii=0; ii<kNclusterMax; ii++)
389       global[ii] = 0x0;
390
391     AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
392     AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(track->GetSeedIndex());
393     //if ( t->GetTOFsignal()>0. ) continue;
394     if ( (t->GetStatus()&AliESDtrack::kTOFout)!=0 ) continue;
395     AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
396
397     // Determine a window around the track
398     Double_t x,par[5]; 
399     trackTOFin->GetExternalParameters(x,par);
400     Double_t cov[15]; 
401     trackTOFin->GetExternalCovariance(cov);
402
403     if (cov[0]<0. || cov[2]<0.) {
404       AliWarning(Form("Very strange track (%d)! At least one of its covariance matrix diagonal elements is negative!",iseed));
405       delete trackTOFin;
406       continue;
407     }
408
409     Double_t dphi=
410       scaleFact*
411       ((5*TMath::Sqrt(cov[0]) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius); 
412     Double_t dz=
413        scaleFact*
414        (5*TMath::Sqrt(cov[2]) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
415
416     Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
417     if (phi<-TMath::Pi())phi+=2*TMath::Pi();
418     if (phi>=TMath::Pi())phi-=2*TMath::Pi();
419     Double_t z=par[1];   
420
421     //upper limit on window's size.
422
423     if(dz> dzMax) dz=dzMax;
424     if(dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
425
426
427     Int_t nc=0;
428
429     // find the clusters in the window of the track
430
431     for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
432
433       if (nc>=kNclusterMax) {
434         AliWarning("No more matchable clusters can be stored! Please, increase the corresponding vectors size.");
435         break;
436       }
437
438       AliTOFcluster *c=fClusters[k];
439       if (c->GetZ() > z+dz) break;
440       if (c->IsUsed()) continue;
441
442       if (!c->GetStatus()) {
443         AliDebug(1,"Cluster in channel declared bad!");
444         continue; // skip bad channels as declared in OCDB
445       }
446
447       Double_t dph=TMath::Abs(c->GetPhi()-phi);
448       if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
449       if (TMath::Abs(dph)>dphi) continue;
450
451       Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
452       Double_t p[2]={yc, c->GetZ()};
453       Double_t cov2[3]= {dY*dY/12., 0., dZ*dZ/12.};
454       if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov2) > maxChi2)continue;
455
456       clind[nc] = k;      
457       Char_t path[100];
458       Int_t ind[5];
459       ind[0]=c->GetDetInd(0);
460       ind[1]=c->GetDetInd(1);
461       ind[2]=c->GetDetInd(2);
462       ind[3]=c->GetDetInd(3);
463       ind[4]=c->GetDetInd(4);
464       fGeom->GetVolumePath(ind,path);
465       gGeoManager->cd(path);
466       global[nc] = *gGeoManager->GetCurrentMatrix();
467       nc++;
468     }
469
470     AliDebug(1,Form(" Number of matchable TOF clusters for the track number %d: %d",iseed,nc));
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=AliTOFGeometry::RinTOF()+istep*0.1;
478       Double_t ymax=xs*TMath::Tan(0.5*AliTOFGeometry::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(AliTOFGeometry::GetAlpha())) {
485           break;
486         }
487       } else if (ysect <-ymax) {
488         if (!trackTOFin->Rotate(-AliTOFGeometry::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     Bool_t accept = kFALSE;
512     Bool_t isInside = kFALSE;
513     for (Int_t istep=0; istep<nStepsDone; istep++) {
514
515       if (nfound>=kNfoundMax) {
516         AliWarning("No more track positions can be stored! Please, increase the corresponding vectors size.");
517         break;
518       }
519
520       Float_t ctrackPos[3];     
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         isInside = fGeom->IsInsideThePad(global[i],ctrackPos,dist3d);
531
532         if ( mLastStep ) {
533           Float_t yLoc = dist3d[1];
534           Float_t rLoc = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[2]*dist3d[2]);
535           accept = (TMath::Abs(yLoc)<padDepth*0.5 && rLoc<dCut);
536           AliDebug(2," I am in the case mLastStep==kTRUE ");
537         }
538         else {
539           accept = isInside;
540         }
541         if (accept) {
542
543           dist[nfound] = TMath::Sqrt(dist3d[0]*dist3d[0] +
544                                      dist3d[1]*dist3d[1] +
545                                      dist3d[2]*dist3d[2]);
546           AliDebug(2,Form(" dist3dLoc[0] = %f, dist3dLoc[1] = %f, dist3dLoc[2] = %f ",
547                           dist3d[0],dist3d[1],dist3d[2]));
548           distZ[nfound] = dist3d[2]; // Z distance in the RF of the
549                                      // hit pad closest to the
550                                      // reconstructed track
551
552           AliDebug(2,Form("   dist3dLoc[2] = %f --- distZ[%d] = %f",
553                           dist3d[2],nfound,distZ[nfound]));
554
555           crecL[nfound] = trackPos[3][istep];
556           index[nfound] = clind[i]; // store cluster id             
557           cxpos[nfound] = AliTOFGeometry::RinTOF()+istep*0.1; //store prop.radius
558           nfound++;
559           if(accept &&!mLastStep)break;
560         }//end if accept
561
562       } //end for on the clusters
563       if(accept &&!mLastStep)break;
564     } //end for on the steps     
565
566     AliDebug(1,Form(" Number of track points for the track number %d: %d",iseed,nfound));
567
568
569     if (nfound == 0 ) {
570       fnunmatch++;
571       delete trackTOFin;
572       continue;
573     }
574     
575     fnmatch++;
576
577     // now choose the cluster to be matched with the track.
578
579     Int_t idclus=-1;
580     Float_t  recL = 0.;
581     Float_t  xpos=0.;
582     Float_t  mindist=1000.;
583     Float_t  mindistZ=0.;
584     for (Int_t iclus= 0; iclus<nfound;iclus++){
585       if (dist[iclus]< mindist){
586         mindist = dist[iclus];
587         mindistZ = distZ[iclus]; // Z distance in the RF of the hit
588                                  // pad closest to the reconstructed
589                                  // track
590         xpos = cxpos[iclus];
591         idclus =index[iclus]; 
592         recL=crecL[iclus]+corrLen*0.5;
593       }
594     }
595
596
597     AliTOFcluster *c=fClusters[idclus];
598     /*
599     Float_t tiltangle = AliTOFGeometry::GetAngles(c->GetDetInd(1),c->GetDetInd(2))*TMath::DegToRad();
600     Float_t localCheck=-mindistZ;
601     localCheck/=TMath::Cos(tiltangle); // Z (tracking/ALICE RF) component of
602                                        // the distance between the
603                                        // reconstructed track and the
604                                        // TOF closest cluster
605     */
606     AliDebug(2, Form("%7d     %7d     %10d     %10d  %10d  %10d      %7d",
607                      iseed,
608                      fnmatch-1,
609                      TMath::Abs(trackTOFin->GetLabel()),
610                      c->GetLabel(0), c->GetLabel(1), c->GetLabel(2),
611                      idclus)); // AdC
612
613     c->Use(); 
614
615     // Track length correction for matching Step 2 
616
617     if(mLastStep){
618       Float_t rc=TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
619       Float_t rt=TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
620                              +trackPos[1][70]*trackPos[1][70]
621                              +trackPos[2][70]*trackPos[2][70]);
622       Float_t dlt=rc-rt;      
623       recL=trackPos[3][70]+dlt;
624     }    
625
626     if (
627         (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
628         ||
629         (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
630         ||
631         (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
632         ) {
633       fngoodmatch++;
634
635        AliDebug(2,Form(" track label good %5d",trackTOFin->GetLabel()));
636
637     }
638     else {
639       fnbadmatch++;
640
641       AliDebug(2,Form(" track label  bad %5d",trackTOFin->GetLabel()));
642
643     }
644
645     delete trackTOFin;
646
647     //  Store quantities to be used in the TOF Calibration
648     Float_t tToT=AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3; // in ns
649     t->SetTOFsignalToT(tToT);
650     Float_t rawTime=AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()+kTimeOffset; // RAW time,in ps
651     t->SetTOFsignalRaw(rawTime);
652     t->SetTOFsignalDz(mindistZ);
653
654     Int_t ind[5];
655     ind[0]=c->GetDetInd(0);
656     ind[1]=c->GetDetInd(1);
657     ind[2]=c->GetDetInd(2);
658     ind[3]=c->GetDetInd(3);
659     ind[4]=c->GetDetInd(4);
660     Int_t calindex = AliTOFGeometry::GetIndex(ind);
661     t->SetTOFCalChannel(calindex);
662
663     // keep track of the track labels in the matched cluster
664     Int_t tlab[3];
665     tlab[0]=c->GetLabel(0);
666     tlab[1]=c->GetLabel(1);
667     tlab[2]=c->GetLabel(2);
668     AliDebug(2,Form(" tdc time of the matched track %6d = ",c->GetTDC()));    
669     Double_t tof=AliTOFGeometry::TdcBinWidth()*c->GetTDC()+kTimeOffset; // in ps
670     AliDebug(2,Form(" tof time of the matched track: %f = ",tof));    
671     Double_t tofcorr=tof;
672     if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
673     AliDebug(2,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));    
674     //Set TOF time signal and pointer to the matched cluster
675     t->SetTOFsignal(tofcorr);
676     t->SetTOFcluster(idclus); // pointing to the recPoints tree
677
678     AliDebug(2,Form(" Setting TOF raw time: %f, z distance: %f  corrected time: %f ",rawTime,mindistZ,tofcorr));
679
680     //Tracking info
681     Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time); // in ps
682     Double_t mom=t->GetP();
683     AliDebug(2,Form(" Momentum for track %d -> %f", iseed,mom));
684     for (Int_t j=0;j<AliPID::kSPECIES;j++) {
685       Double_t mass=AliPID::ParticleMass(j);
686       time[j]+=(recL-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(mom*mom+mass*mass)/mom;
687     }
688
689     AliTOFtrack *trackTOFout = new AliTOFtrack(*t); 
690     trackTOFout->PropagateTo(xpos);
691
692     // Fill the track residual histograms.
693     FillResiduals(trackTOFout,c,kFALSE);
694
695     t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
696     t->SetIntegratedLength(recL);
697     t->SetIntegratedTimes(time);
698     t->SetTOFLabel(tlab);
699
700  
701     // Fill Reco-QA histos for Reconstruction
702     fHRecNClus->Fill(nc);
703     fHRecDist->Fill(mindist);
704     fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
705     fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
706     fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
707     fHRecSigZVsPWin->Fill(mom,dz);
708
709     // Fill Tree for on-the-fly offline Calibration
710
711     if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 ) ) {
712       fIch=calindex;
713       fToT=tToT;
714       fTime=rawTime;
715       fExpTimePi=time[2];
716       fExpTimeKa=time[3];
717       fExpTimePr=time[4];
718       fCalTree->Fill();
719     }
720     delete trackTOFout;
721   }
722   for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
723   delete [] clind;
724  
725 }
726 //_________________________________________________________________________
727 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
728   //--------------------------------------------------------------------
729   //This function loads the TOF clusters
730   //--------------------------------------------------------------------
731
732   Int_t npadX = AliTOFGeometry::NpadX();
733   Int_t npadZ = AliTOFGeometry::NpadZ();
734   Int_t nStripA = AliTOFGeometry::NStripA();
735   Int_t nStripB = AliTOFGeometry::NStripB();
736   Int_t nStripC = AliTOFGeometry::NStripC();
737
738   TBranch *branch=cTree->GetBranch("TOF");
739   if (!branch) { 
740     AliError("can't get the branch with the TOF clusters !");
741     return 1;
742   }
743
744   static TClonesArray dummy("AliTOFcluster",10000);
745   dummy.Clear();
746   TClonesArray *clusters=&dummy;
747   branch->SetAddress(&clusters);
748
749   cTree->GetEvent(0);
750   Int_t nc=clusters->GetEntriesFast();
751   fHDigNClus->Fill(nc);
752
753   AliInfo(Form("Number of clusters: %d",nc));
754
755   for (Int_t i=0; i<nc; i++) {
756     AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
757     fClusters[i]=new AliTOFcluster(*c); fN++;
758
759   // Fill Digits QA histos
760  
761     Int_t isector = c->GetDetInd(0);
762     Int_t iplate = c->GetDetInd(1);
763     Int_t istrip = c->GetDetInd(2);
764     Int_t ipadX = c->GetDetInd(4);
765     Int_t ipadZ = c->GetDetInd(3);
766
767     Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
768     Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
769  
770     Int_t stripOffset = 0;
771     switch (iplate) {
772     case 0:
773       stripOffset = 0;
774       break;
775     case 1:
776       stripOffset = nStripC;
777       break;
778     case 2:
779       stripOffset = nStripC+nStripB;
780       break;
781     case 3:
782       stripOffset = nStripC+nStripB+nStripA;
783       break;
784     case 4:
785       stripOffset = nStripC+nStripB+nStripA+nStripB;
786       break;
787     default:
788       AliError(Form("Wrong plate number in TOF (%d) !",iplate));
789       break;
790     };
791     Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
792     Int_t phiindex=npadX*isector+ipadX+1;
793     fHDigClusMap->Fill(zindex,phiindex);
794     fHDigClusTime->Fill(time);
795     fHDigClusToT->Fill(tot);
796
797   }
798
799
800   return 0;
801 }
802 //_________________________________________________________________________
803 void AliTOFtracker::UnloadClusters() {
804   //--------------------------------------------------------------------
805   //This function unloads TOF clusters
806   //--------------------------------------------------------------------
807   for (Int_t i=0; i<fN; i++) {
808     delete fClusters[i];
809     fClusters[i] = 0x0;
810   }
811   fN=0;
812 }
813
814 //_________________________________________________________________________
815 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
816   //--------------------------------------------------------------------
817   // This function returns the index of the nearest cluster 
818   //--------------------------------------------------------------------
819   if (fN==0) return 0;
820   if (z <= fClusters[0]->GetZ()) return 0;
821   if (z > fClusters[fN-1]->GetZ()) return fN;
822   Int_t b=0, e=fN-1, m=(b+e)/2;
823   for (; b<e; m=(b+e)/2) {
824     if (z > fClusters[m]->GetZ()) b=m+1;
825     else e=m; 
826   }
827   return m;
828 }
829
830 //_________________________________________________________________________
831 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
832 {
833   // Get track space point with index i
834   // Coordinates are in the global system
835   AliTOFcluster *cl = fClusters[index];
836   Float_t xyz[3];
837   xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
838   xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
839   xyz[2] = cl->GetZ();
840   Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
841   Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
842   Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
843   Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
844   Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
845   Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
846   Float_t cov[6];
847   cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
848   cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
849   cov[2] = -cosphi*sinth*costh*sigmaz2;
850   cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
851   cov[4] = -sinphi*sinth*costh*sigmaz2;
852   cov[5] = costh*costh*sigmaz2;
853   p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
854
855   // Detector numbering scheme
856   Int_t nSector = AliTOFGeometry::NSectors();
857   Int_t nPlate  = AliTOFGeometry::NPlates();
858   Int_t nStripA = AliTOFGeometry::NStripA();
859   Int_t nStripB = AliTOFGeometry::NStripB();
860   Int_t nStripC = AliTOFGeometry::NStripC();
861
862   Int_t isector = cl->GetDetInd(0);
863   if (isector >= nSector)
864     AliError(Form("Wrong sector number in TOF (%d) !",isector));
865   Int_t iplate = cl->GetDetInd(1);
866   if (iplate >= nPlate)
867     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
868   Int_t istrip = cl->GetDetInd(2);
869
870   Int_t stripOffset = 0;
871   switch (iplate) {
872   case 0:
873     stripOffset = 0;
874     break;
875   case 1:
876     stripOffset = nStripC;
877     break;
878   case 2:
879     stripOffset = nStripC+nStripB;
880     break;
881   case 3:
882     stripOffset = nStripC+nStripB+nStripA;
883     break;
884   case 4:
885     stripOffset = nStripC+nStripB+nStripA+nStripB;
886     break;
887   default:
888     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
889     break;
890   };
891
892   Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
893                stripOffset +
894                istrip;
895   UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
896   p.SetVolumeID((UShort_t)volid);
897   return kTRUE;
898 }
899 //_________________________________________________________________________
900 void AliTOFtracker::InitCheckHists() {
901
902   //Init histos for Digits/Reco QA and Calibration
903
904
905   TDirectory *dir = gDirectory;
906   TFile *logFileTOF = 0;
907
908   TSeqCollection *list = gROOT->GetListOfFiles();
909   int n = list->GetEntries();
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   if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
920   logFileTOF->cd(); 
921
922   fCalTree = new TTree("CalTree", "Tree for TOF calibration");
923   fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
924   fCalTree->Branch("ToT",&fToT,"TOFToT/F");
925   fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
926   fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
927   fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
928   fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
929
930   //Digits "QA" 
931   fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);  
932   fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);  
933   fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);  
934   fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);  
935
936   //Reco "QA"
937   fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
938   fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
939   fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
940   fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
941   fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
942   fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
943
944   dir->cd();
945
946 }
947
948 //_________________________________________________________________________
949 void AliTOFtracker::SaveCheckHists() {
950
951   //write histos for Digits/Reco QA and Calibration
952
953   TDirectory *dir = gDirectory;
954   TFile *logFileTOF = 0;
955
956   TSeqCollection *list = gROOT->GetListOfFiles();
957   int n = list->GetEntries();
958   Bool_t isThere=kFALSE;
959   for(int i=0; i<n; i++) {
960     logFileTOF = (TFile*)list->At(i);
961     if (strstr(logFileTOF->GetName(), "TOFQA.root")){
962       isThere=kTRUE;
963       break;
964     } 
965   }
966    
967   if(!isThere) {
968           AliError(Form("File TOFQA.root not found!! not wring histograms...."));
969           return;
970   }
971   logFileTOF->cd(); 
972   fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
973   fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
974   fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
975   fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
976   fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
977   fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
978   fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
979   fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
980   fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
981   fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
982   fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
983   logFileTOF->Flush();  
984
985   dir->cd();
986   }
987 //_________________________________________________________________________
988 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) {
989
990   //dummy, for the moment
991   Float_t tofcorr=0.;
992   if(dist<AliTOFGeometry::ZPad()*0.5){
993     tofcorr=tof;
994     //place here the actual correction
995   }else{
996     tofcorr=tof; 
997   } 
998   return tofcorr;
999 }
1000 //_________________________________________________________________________
1001 Float_t AliTOFtracker::GetTimeZerofromT0(AliESDEvent *event) const {
1002
1003   //Returns TimeZero as measured by T0 detector
1004
1005   return event->GetT0();
1006 }
1007 //_________________________________________________________________________
1008 Float_t AliTOFtracker::GetTimeZerofromTOF(AliESDEvent * /*event*/) const {
1009
1010   //dummy, for the moment. T0 algorithm using tracks on TOF
1011   {
1012     //place T0 algo here...
1013   }
1014   return 0.;
1015 }
1016 //_________________________________________________________________________
1017
1018 void AliTOFtracker::FillClusterArray(TObjArray* arr) const
1019 {
1020   //
1021   // Returns the TOF cluster array
1022   //
1023
1024   if (fN==0)
1025     arr = 0x0;
1026   else
1027     for (Int_t i=0; i<fN; ++i) arr->Add(fClusters[i]);
1028
1029 }
1030 //_________________________________________________________________________
1031