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