]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFtrackerV2.cxx
Adding set of AliESDEvent in the track for TOF matching
[u/mrichter/AliRoot.git] / TOF / AliTOFtrackerV2.cxx
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43 /**************************************************************************
44  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
45  *                                                                        *
46  * Author: The ALICE Off-line Project.                                    *
47  * Contributors are mentioned in the code where appropriate.              *
48  *                                                                        *
49  * Permission to use, copy, modify and distribute this software and its   *
50  * documentation strictly for non-commercial purposes is hereby granted   *
51  * without fee, provided that the above copyright notice appears in all   *
52  * copies and that both the copyright notice and this permission notice   *
53  * appear in the supporting documentation. The authors make no claims     *
54  * about the suitability of this software for any purpose. It is          *
55  * provided "as is" without express or implied warranty.                  *
56  **************************************************************************/
57
58 //--------------------------------------------------------------------//
59 //                                                                    //
60 // AliTOFtrackerV2 Class                                              //
61 // Task: Perform association of the ESD tracks to TOF Clusters        //
62 // and Update ESD track with associated TOF Cluster parameters        //
63 //                                                                    //
64 // -- Authors : S. Arcelli, C. Zampolli (Bologna University and INFN) //
65 // -- Contacts: Annalisa.De.Caro@cern.ch                              //
66 // --         : Chiara.Zampolli@bo.infn.it                            //
67 // --         : Silvia.Arcelli@bo.infn.it                             //
68 //                                                                    //
69 //--------------------------------------------------------------------//
70
71 #include <Rtypes.h>
72 #include <TROOT.h>
73
74 #include <TClonesArray.h>
75 #include <TObjArray.h>
76 #include <TGeoManager.h>
77 #include <TTree.h>
78
79 #include "AliGeomManager.h"
80 #include "AliESDtrack.h"
81 #include "AliESDEvent.h"
82 #include "AliESDpid.h"
83 #include "AliESDTOFcluster.h"
84 #include "AliLog.h"
85 #include "AliTrackPointArray.h"
86 #include "AliCDBManager.h"
87
88 #include "AliTOFRecoParam.h"
89 #include "AliTOFReconstructor.h"
90 #include "AliTOFcluster.h"
91 #include "AliTOFGeometry.h"
92 #include "AliTOFtrackerV2.h"
93 #include "AliTOFtrack.h"
94
95 extern TGeoManager *gGeoManager;
96
97 ClassImp(AliTOFtrackerV2)
98
99 //_____________________________________________________________________________
100 AliTOFtrackerV2::AliTOFtrackerV2():
101   fkRecoParam(0x0),
102   fGeom(0x0),
103   fN(0),
104   fNseeds(0),
105   fNseedsTOF(0),
106   fnunmatch(0),
107   fnmatch(0),
108   fTracks(new TClonesArray("AliTOFtrack")),
109   fSeeds(new TObjArray(100)),
110   fClusters(0x0)
111 {
112   //AliTOFtrackerV2 main Ctor
113
114   // Getting the geometry
115   fGeom = new AliTOFGeometry();
116
117 }
118 //_____________________________________________________________________________
119 AliTOFtrackerV2::~AliTOFtrackerV2() {
120   //
121   // Dtor
122   //
123
124   if(!(AliCDBManager::Instance()->GetCacheFlag())){
125     delete fkRecoParam;
126   }
127   delete fGeom; 
128   if (fTracks){
129     fTracks->Delete();
130     delete fTracks;
131     fTracks=0x0;
132   }
133   if (fSeeds){
134     fSeeds->Delete();
135     delete fSeeds;
136     fSeeds=0x0;
137   }
138
139   if(fClusters){
140     delete[] fClusters;
141     fClusters = NULL;
142   }
143
144 }
145 //_____________________________________________________________________________
146 void AliTOFtrackerV2::GetPidSettings(AliESDpid *esdPID) {
147   // 
148   // Sets TOF resolution from RecoParams
149   //
150   if (fkRecoParam)
151     esdPID->GetTOFResponse().SetTimeResolution(fkRecoParam->GetTimeResolution());
152   else
153     AliWarning("fkRecoParam not yet set; cannot set PID settings");
154
155 //_____________________________________________________________________________
156 Int_t AliTOFtrackerV2::PropagateBack(AliESDEvent * const event) {
157   //
158   // Gets seeds from ESD event and Match with TOF Clusters
159   //
160
161   //Update the matched ESD tracks
162   // needed in case of call of TOF info before of the selection of matching and in case of no clusters available at all
163   if(fN==0)
164     event->SetTOFcluster(1,fClusters); 
165   else
166     event->SetTOFcluster(fN,fClusters);
167
168   if (fN==0) {
169     AliInfo("No TOF recPoints to be matched with reconstructed tracks");
170     return 0;
171   }
172
173   // initialize RecoParam for current event
174   AliDebug(1,"Initializing params for TOF");
175
176   fkRecoParam = AliTOFReconstructor::GetRecoParam();  // instantiate reco param from STEER...
177
178   if (fkRecoParam == 0x0) { 
179     AliFatal("No Reco Param found for TOF!!!");
180   }
181
182   //Initialise some counters
183
184   fNseeds=0;
185   fNseedsTOF=0;
186   fnunmatch=0;
187   fnmatch=0;
188
189   Int_t ntrk=event->GetNumberOfTracks();
190   fNseeds = ntrk;
191
192   //Load ESD tracks into a local Array of ESD Seeds
193   for (Int_t i=0; i<fNseeds; i++){
194     fSeeds->AddLast(event->GetTrack(i));
195     event->GetTrack(i)->SetESDEvent(event);
196   }
197   //Prepare ESD tracks candidates for TOF Matching
198   CollectESD();
199
200   if (fNseeds==0 || fNseedsTOF==0) {
201     AliInfo("No seeds to try TOF match");
202     return 0;
203   }
204
205   // clusterize before of matching
206   Clusterize();
207
208   //Second Step with Looser Matching Criterion
209   MatchTracks();
210
211   AliInfo(Form("Number of matched tracks = %d",fnmatch));
212
213   for (Int_t i=0; i<ntrk; i++) {
214     AliESDtrack *t=event->GetTrack(i);
215     //t->SetESDEvent(event);
216     AliESDtrack *seed =(AliESDtrack*)fSeeds->At(i);
217
218     if ( (seed->GetStatus()&AliESDtrack::kTOFin)!=0 ) {
219       t->SetStatus(AliESDtrack::kTOFin);
220       if ( (seed->GetStatus()&AliESDtrack::kTOFout)!=0 ) {
221         t->SetStatus(AliESDtrack::kTOFout);
222         //t->SetTOFclusterArray(seed->GetNTOFclusters(),seed->GetTOFclusterArray());
223         t->SortTOFcluster();
224
225         // Make attention, please:
226         //      AliESDtrack::fTOFInfo array does not be stored in the AliESDs.root file
227         //      it is there only for a check during the reconstruction step.
228         Float_t info[10]; seed->GetTOFInfo(info);
229         t->SetTOFInfo(info);
230         AliDebug(3,Form(" distance=%f; residual in the pad reference frame: dX=%f, dZ=%f", info[0],info[1],info[2]));
231
232         /*
233         Double_t alphaA = (Double_t)t->GetAlpha();
234         Double_t xA = (Double_t)t->GetX();
235         Double_t yA = (Double_t)t->GetY();
236         Double_t zA = (Double_t)t->GetZ();
237         Double_t p1A = (Double_t)t->GetSnp();
238         Double_t p2A = (Double_t)t->GetTgl();
239         Double_t p3A = (Double_t)t->GetSigned1Pt();
240         const Double_t *covA = (Double_t*)t->GetCovariance();
241
242         // Check done:
243         //       by calling the AliESDtrack::UpdateTrackParams,
244         //       the current track parameters are changed
245         //       and it could cause refit problems.
246         //       We need to update only the following track parameters:
247         //            the track length and expected times.
248         //       Removed AliESDtrack::UpdateTrackParams call
249         //       Called AliESDtrack::SetIntegratedTimes(...) and
250         //       AliESDtrack::SetIntegratedLength() routines.
251
252         AliTOFtrack *track = new AliTOFtrack(*seed);
253         t->UpdateTrackParams(track,AliESDtrack::kTOFout); // to be checked - AdC
254         delete track;
255         Double_t time[AliPID::kSPECIESC]; t->GetIntegratedTimes(time);
256         Double_t alphaB = (Double_t)t->GetAlpha();
257         Double_t xB = (Double_t)t->GetX();
258         Double_t yB = (Double_t)t->GetY();
259         Double_t zB = (Double_t)t->GetZ();
260         Double_t p1B = (Double_t)t->GetSnp();
261         Double_t p2B = (Double_t)t->GetTgl();
262         Double_t p3B = (Double_t)t->GetSigned1Pt();
263         const Double_t *covB = (Double_t*)t->GetCovariance();
264         AliDebug(2,"Track params -now(before)-:");
265         AliDebug(2,Form("    X: %f(%f), Y: %f(%f), Z: %f(%f) --- alpha: %f(%f)",
266                         xB,xA,
267                         yB,yA,
268                         zB,zA,
269                         alphaB,alphaA));
270         AliDebug(2,Form("    p1: %f(%f), p2: %f(%f), p3: %f(%f)",
271                         p1B,p1A,
272                         p2B,p2A,
273                         p3B,p3A));
274         AliDebug(2,Form("    cov1: %f(%f), cov2: %f(%f), cov3: %f(%f)"
275                         " cov4: %f(%f), cov5: %f(%f), cov6: %f(%f)"
276                         " cov7: %f(%f), cov8: %f(%f), cov9: %f(%f)"
277                         " cov10: %f(%f), cov11: %f(%f), cov12: %f(%f)"
278                         " cov13: %f(%f), cov14: %f(%f), cov15: %f(%f)",
279                         covB[0],covA[0],
280                         covB[1],covA[1],
281                         covB[2],covA[2],
282                         covB[3],covA[3],
283                         covB[4],covA[4],
284                         covB[5],covA[5],
285                         covB[6],covA[6],
286                         covB[7],covA[7],
287                         covB[8],covA[8],
288                         covB[9],covA[9],
289                         covB[10],covA[10],
290                         covB[11],covA[11],
291                         covB[12],covA[12],
292                         covB[13],covA[13],
293                         covB[14],covA[14]
294                         ));
295         */
296         Double_t time[AliPID::kSPECIESC]; t->GetIntegratedTimes(time);
297         AliDebug(2,Form(" TOF params: %6d  %f %f %f %f %f %6d %3d %f",
298                         i,
299                         t->GetTOFsignalRaw(),t->GetTOFsignal(),t->GetTOFsignalToT(),
300                         t->GetTOFsignalDz(),t->GetTOFsignalDx(),t->GetTOFCalChannel(),
301                         t->GetTOFcluster(),t->GetIntegratedLength()));
302         AliDebug(2,Form("  %f %f %f %f %f %f %f %f %f",
303                         time[0], time[1], time[2], time[3], time[4], time[5], time[6], time[7], time[8]));
304       }
305     }
306   }
307   
308   fSeeds->Clear();
309   fTracks->Clear();
310   
311   AliInfo(Form("Number of cluster to be checked = %d",fN));
312   if(fN){
313     Int_t *matchmap = new Int_t[fN];
314     event->SetTOFcluster(fN,fClusters,matchmap);
315     for (Int_t i=0; i<ntrk; i++) { // remapping after TOF matching selection
316       AliESDtrack *t=event->GetTrack(i);
317       t->ReMapTOFcluster(fN,matchmap);
318     }
319
320     delete[] matchmap;
321   }
322   
323   
324
325   return 0;
326   
327 }
328 //_________________________________________________________________________
329 void AliTOFtrackerV2::CollectESD() {
330    //prepare the set of ESD tracks to be matched to clusters in TOF
331
332   Int_t seedsTOF1=0;
333   Int_t seedsTOF3=0;
334   Int_t seedsTOF2=0;
335  
336   TClonesArray &aTOFTrack = *fTracks;
337   for (Int_t i=0; i<fNseeds; i++) {
338
339     AliESDtrack *t =(AliESDtrack*)fSeeds->At(i);
340     if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
341
342     AliTOFtrack *track = new AliTOFtrack(*t); // New
343     Float_t x = (Float_t)track->GetX(); //New
344
345     // TRD 'good' tracks
346     if ( ( (t->GetStatus()&AliESDtrack::kTRDout)!=0 ) ) {
347
348       AliDebug(1,Form(" Before propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
349
350       // TRD 'good' tracks, already propagated at 371 cm
351       if ( x >= AliTOFGeometry::Rmin() ) {
352
353         if ( track->PropagateToInnerTOF() ) {
354
355           AliDebug(1,Form(" TRD propagated track till rho = %fcm."
356                           " And then the track has been propagated till rho = %fcm.",
357                           x, (Float_t)track->GetX()));
358
359           track->SetSeedIndex(i);
360           t->UpdateTrackParams(track,AliESDtrack::kTOFin);
361           new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
362           fNseedsTOF++;
363           seedsTOF1++;
364
365           AliDebug(1,Form(" After propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
366         }
367         delete track;
368
369       }
370       else { // TRD 'good' tracks, propagated rho<371cm
371
372         if  ( track->PropagateToInnerTOF() ) {
373
374           AliDebug(1,Form(" TRD propagated track till rho = %fcm."
375                           " And then the track has been propagated till rho = %fcm.",
376                           x, (Float_t)track->GetX()));
377
378           track->SetSeedIndex(i);
379           t->UpdateTrackParams(track,AliESDtrack::kTOFin);
380           new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
381           fNseedsTOF++;
382           seedsTOF3++;
383
384           AliDebug(1,Form(" After propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
385         }
386         delete track;
387
388       }
389       //delete track;
390     }
391
392     else { // Propagate the rest of TPCbp
393
394       AliDebug(1,Form(" Before propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
395
396       if ( track->PropagateToInnerTOF() ) { 
397
398         AliDebug(1,Form(" TPC propagated track till rho = %fcm."
399                         " And then the track has been propagated till rho = %fcm.",
400                         x, (Float_t)track->GetX()));
401
402         track->SetSeedIndex(i);
403         t->UpdateTrackParams(track,AliESDtrack::kTOFin);
404         new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
405         fNseedsTOF++;
406         seedsTOF2++;
407
408         AliDebug(1,Form(" After propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
409       }
410       delete track;
411     }
412   }
413
414   AliInfo(Form("Number of TOF seeds = %d (kTRDout371 = %d, kTRDoutLess371 = %d, !kTRDout = %d)",fNseedsTOF,seedsTOF1,seedsTOF3,seedsTOF2));
415
416   // Sort according uncertainties on track position 
417   fTracks->Sort();
418
419 }
420
421 //_________________________________________________________________________
422 void AliTOFtrackerV2::MatchTracks() {
423   //
424   //Match ESD tracks to clusters in TOF
425   //
426
427   // Parameters used/regulating the reconstruction
428   static Float_t detDepth=18.;
429   static Float_t padDepth=0.5;
430
431   const Float_t kSpeedOfLight= 2.99792458e-2; // speed of light [cm/ps]
432
433   Float_t dY=AliTOFGeometry::XPad(); 
434   Float_t dZ=AliTOFGeometry::ZPad(); 
435
436   Float_t sensRadius = fkRecoParam->GetSensRadius();
437   Float_t stepSize   = fkRecoParam->GetStepSize();
438   Float_t scaleFact  = fkRecoParam->GetWindowScaleFact();
439   Float_t dyMax=fkRecoParam->GetWindowSizeMaxY(); 
440   Float_t dzMax=fkRecoParam->GetWindowSizeMaxZ();
441   Float_t dCut=10.;//fkRecoParam->GetDistanceCut(); // This is to be loaded by OCDB. It should be 10cm always.
442   Double_t maxChi2=fkRecoParam->GetMaxChi2TRD();
443   Bool_t timeWalkCorr = fkRecoParam->GetTimeWalkCorr();
444   AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++++");
445   AliDebug(1,Form("TOF sens radius: %f",sensRadius));
446   AliDebug(1,Form("TOF step size: %f",stepSize));
447   AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
448   AliDebug(1,Form("TOF Window max dy: %f",dyMax));
449   AliDebug(1,Form("TOF Window max dz: %f",dzMax));
450   AliDebug(1,Form("TOF distance Cut: %f",dCut));
451   AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
452   AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));   
453
454   //Match ESD tracks to clusters in TOF
455
456   // Get the number of propagation steps
457   Int_t nSteps=(Int_t)(detDepth/stepSize);
458   AliDebug(1,Form(" Number of steps to be done %d",nSteps));
459
460   AliDebug(1,"++++++++++++++++++++++++++++++++++++++++++++++++++++++++++");
461
462   //PH Arrays (moved outside of the loop)
463   Float_t * trackPos[4];
464   for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
465   Int_t * clind = new Int_t[fN];
466   
467   // Some init
468   const Int_t kNclusterMax = 1000; // related to fN value
469   TGeoHMatrix global[kNclusterMax];
470
471   //The matching loop
472   for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
473
474     for (Int_t ii=0; ii<fN; ii++) clind[ii]=-1;
475     for (Int_t ii=0; ii<kNclusterMax; ii++) global[ii] = 0x0;
476     for (Int_t ii=0; ii<4; ii++)
477       for (Int_t jj=0; jj<nSteps; jj++) trackPos[ii][jj]=0.;
478
479     AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
480     AliESDtrack *t =(AliESDtrack*)fSeeds->At(track->GetSeedIndex());
481     AliTOFtrack *trackTOFin = new AliTOFtrack(*track);
482
483     Double_t timesOr[AliPID::kSPECIESC]; t->GetIntegratedTimes(timesOr); // in ps
484
485     // Determine a window around the track
486     Double_t x,par[5]; 
487     trackTOFin->GetExternalParameters(x,par);
488     Double_t cov[15]; 
489     trackTOFin->GetExternalCovariance(cov);
490
491     if (cov[0]<0. || cov[2]<0.) {
492       AliWarning(Form("Very strange track (%d)! At least one of its covariance matrix diagonal elements is negative!",iseed));
493       delete trackTOFin;
494       continue;
495     }
496
497     Double_t dphi=
498       scaleFact*
499       ((5*TMath::Sqrt(TMath::Abs(cov[0])) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius); 
500     Double_t dz=
501        scaleFact*
502        (5*TMath::Sqrt(TMath::Abs(cov[2])) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
503
504     Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
505     if (phi<-TMath::Pi())phi+=2*TMath::Pi();
506     if (phi>=TMath::Pi())phi-=2*TMath::Pi();
507     Double_t z=par[1];   
508
509     //upper limit on window's size.
510     if (dz> dzMax) dz=dzMax;
511     if (dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
512
513
514     // find the clusters in the window of the track
515     Int_t nc=0;
516     for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
517
518       if (nc>=kNclusterMax) {
519         AliWarning("No more matchable clusters can be stored! Please, increase the corresponding vectors size.");
520         break;
521       }
522
523       AliESDTOFcluster *c=&(fClusters[k]);
524       if (c->GetZ() > z+dz) break;
525       if (!c->GetStatus()) {
526         AliDebug(1,"Cluster in channel declared bad!");
527         continue; // skip bad channels as declared in OCDB
528       }
529
530       Double_t dph=TMath::Abs(c->GetPhi()-phi);
531       if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
532       if (TMath::Abs(dph)>dphi) continue;
533
534       Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
535       Double_t p[2]={yc, c->GetZ()};
536       Double_t cov2[3]= {dY*dY/12., 0., dZ*dZ/12.};
537       if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov2) > maxChi2)continue;
538
539       clind[nc] = k;      
540       Char_t path[200];
541       Int_t ind[5]; fGeom->GetVolumeIndices(c->GetTOFchannel(),ind);
542       fGeom->GetVolumePath(ind,path);
543       gGeoManager->cd(path);
544       global[nc] = *gGeoManager->GetCurrentMatrix();
545       nc++;
546     }
547
548     if (nc == 0 ) {
549       AliDebug(1,Form("No available clusters for the track number %d",iseed));
550       fnunmatch++;
551       delete trackTOFin;
552       continue;
553     }
554
555     AliDebug(1,Form(" Number of available TOF clusters for the track number %d: %d",iseed,nc));
556
557     //start fine propagation 
558
559     Double_t *times[AliPID::kSPECIESC];
560     for(Int_t isp=0;isp < AliPID::kSPECIESC;isp++){
561       times[isp] = new Double_t[nSteps];
562     }
563
564     Int_t nStepsDone = 0;
565     for( Int_t istep=0; istep<nSteps; istep++){ 
566       
567       // First of all, propagate the track...
568       Float_t xs = AliTOFGeometry::RinTOF()+istep*stepSize;
569       if (!(trackTOFin->PropagateTo(xs))) break;
570
571       //  ...and then, if necessary, rotate the track
572       Double_t ymax = xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
573       Double_t ysect = trackTOFin->GetY();
574       if (ysect > ymax) {
575         if (!(trackTOFin->Rotate(AliTOFGeometry::GetAlpha()))) break;
576       } else if (ysect <-ymax) {
577         if (!(trackTOFin->Rotate(-AliTOFGeometry::GetAlpha()))) break;
578       }
579
580       Double_t mom = trackTOFin->P();
581
582       if(istep == 0){
583         for(Int_t isp=0;isp<AliPID::kSPECIESC;isp++){
584           Double_t mass=AliPID::ParticleMass(isp);
585           Double_t momz = mom*AliPID::ParticleCharge(isp);
586           times[isp][nStepsDone] = stepSize/kSpeedOfLight*TMath::Sqrt(momz*momz+mass*mass)/momz;
587         }
588       }
589       else{
590         for(Int_t isp=0;isp<AliPID::kSPECIESC;isp++){
591           Double_t mass=AliPID::ParticleMass(isp);
592           Double_t momz = mom*AliPID::ParticleCharge(isp);
593           times[isp][nStepsDone] = times[isp][nStepsDone-1] + (trackTOFin->GetIntegratedLength()-trackPos[3][nStepsDone-1])/kSpeedOfLight*TMath::Sqrt(momz*momz+mass*mass)/momz;
594         }
595       }
596
597       // store the running point (Globalrf) - fine propagation     
598
599       Double_t r[3]; trackTOFin->GetXYZ(r);
600       trackPos[0][nStepsDone]= (Float_t) r[0];
601       trackPos[1][nStepsDone]= (Float_t) r[1];
602       trackPos[2][nStepsDone]= (Float_t) r[2];   
603       trackPos[3][nStepsDone]= trackTOFin->GetIntegratedLength();
604
605       nStepsDone++;
606       AliDebug(3,Form(" current step %d (%d) - nStepsDone=%d",istep,nSteps,nStepsDone));
607     }
608
609     if ( nStepsDone == 0 ) {
610       AliDebug(1,Form(" No track points for track number %d",iseed));
611       fnunmatch++;
612       delete trackTOFin;
613       continue;
614     }
615
616     AliDebug(3,Form(" Number of steps done for the track number %d: %d",iseed,nStepsDone));
617
618     Int_t *isClusterMatchable = NULL;
619     if(nc){
620       isClusterMatchable = new Int_t[nc];
621       for (Int_t i=0; i<nc; i++) isClusterMatchable[i] = kFALSE;                
622     }
623
624     Int_t nfound = 0;
625     Bool_t accept = kFALSE;
626     Bool_t isInside = kFALSE;
627     for (Int_t istep=0; istep<nStepsDone; istep++) {
628
629       Bool_t gotInsideCluster = kFALSE;
630       Int_t trackInsideCluster = -1;
631
632       Float_t ctrackPos[3];     
633       ctrackPos[0] = trackPos[0][istep];
634       ctrackPos[1] = trackPos[1][istep];
635       ctrackPos[2] = trackPos[2][istep];
636
637       //now see whether the track matches any of the TOF clusters            
638
639       Float_t dist3d[3]={0.,0.,0.};
640       accept = kFALSE;
641      
642       for (Int_t i=0; i<nc; i++) {
643
644         // ***** NEW *****
645         /* check whether track was inside another cluster
646          * and in case inhibit this cluster.
647          * this will allow to only go on and add track points for
648          * that cluster where the track got inside first */
649         if (gotInsideCluster && trackInsideCluster != i) {
650           AliDebug(3,Form(" A - istep=%d ~ %d %d ~ nfound=%d",istep,trackInsideCluster,i,nfound));
651           continue;
652         }
653         AliDebug(3,Form(" B - istep=%d ~ %d %d ~ nfound=%d",istep,trackInsideCluster,i,nfound));
654
655         /* check whether track is inside this cluster */
656         for (Int_t hh=0; hh<3; hh++) dist3d[hh]=0.;
657         isInside = fGeom->IsInsideThePad((TGeoHMatrix*)(&global[i]),ctrackPos,dist3d);
658
659         // ***** NEW *****
660         /* if track is inside this cluster set flags which will then
661          * inhibit to add track points for the other clusters */
662         if (isInside) {
663           gotInsideCluster = kTRUE;
664           trackInsideCluster = i;
665         }
666
667         Float_t yLoc = dist3d[1];
668         Float_t rLoc = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[2]*dist3d[2]);
669         accept = (TMath::Abs(yLoc)<padDepth*0.5 && rLoc<dCut);
670
671         //***** NEW *****
672         /* add point everytime that:
673          * - the track is inside the cluster
674          * - the track got inside the cluster, even when it eventually exited the cluster
675          * - the tracks is within dCut from the cluster
676          */
677         if (accept || isInside || gotInsideCluster) {
678
679           Double_t timesCurrent[AliPID::kSPECIESC];
680           AliDebug(3,Form(" Momentum for track %d -> %f", iseed,t->P()));
681           for (Int_t j=0;j<AliPID::kSPECIESC;j++) {
682             timesCurrent[j] = timesOr[j] + times[j][istep];
683           }
684
685           if (TMath::Abs(dist3d[1])<stepSize && !isClusterMatchable[i]) {
686             isClusterMatchable[i] = kTRUE;
687             fClusters[clind[i]].Update(t->GetID(),dist3d[1],dist3d[0],dist3d[2],trackPos[3][istep],timesCurrent);//x,y,z -> tracking RF
688             t->AddTOFcluster(clind[i]);
689             t->SetStatus(AliESDtrack::kTOFout);
690           }
691           AliDebug(2,Form(" dist3dLoc[0] = %f, dist3dLoc[1] = %f, dist3dLoc[2] = %f ",dist3d[0],dist3d[1],dist3d[2]));
692
693           nfound++;
694
695           AliDebug(3,Form(" C - istep=%d ~ %d %d ~ nfound=%d",istep,trackInsideCluster,i,nfound));
696         
697           // ***** NEW *****
698         }//end if accept
699         
700       } //end for on the clusters
701     } //end for on the steps     
702     if(nc) delete[] isClusterMatchable;
703
704     for(Int_t isp=0;isp < AliPID::kSPECIESC;isp++){
705       delete[] times[isp];
706     }
707
708
709     if (nfound == 0 ) {
710       AliDebug(1,Form(" No matchable track points for the track number %d",iseed));
711       fnunmatch++;
712       delete trackTOFin;
713       continue;
714     }
715
716     AliDebug(1,Form(" Number of track points for the track number %d: %d",iseed,nfound));
717
718     Int_t nMatchedClusters = t->GetNTOFclusters();
719  
720     if (nMatchedClusters==0) {
721       AliDebug(1,Form("Reconstructed track %d doesn't match any TOF cluster", iseed));
722       fnunmatch++;
723       delete trackTOFin;
724       continue;
725     }
726
727     AliDebug(1,Form(" %d - matched (%d)",track->GetSeedIndex()/*iseed*/,nMatchedClusters));
728
729     fnmatch++;
730
731     /*
732     AliTOFcluster cTOF = AliTOFcluster(volIdClus,
733     (Float_t)posClus[0],(Float_t)posClus[1],(Float_t)posClus[2],
734     (Float_t)covClus[0],(Float_t)covClus[1],(Float_t)covClus[2],
735     (Float_t)covClus[3],(Float_t)covClus[4],(Float_t)covClus[5],
736     tofLabels,volIndices,parClu,kTRUE,index[i]);
737
738     // Fill the track residual histograms.
739     FillResiduals(trackTOFin,c,kFALSE);
740     */
741
742     delete trackTOFin;
743
744   } // loop on fSeeds
745
746   for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
747   delete [] clind;
748  
749 }
750 //_________________________________________________________________________
751 Int_t AliTOFtrackerV2::LoadClusters(TTree *cTree) {
752   //--------------------------------------------------------------------
753   //This function loads the TOF clusters
754   //--------------------------------------------------------------------
755
756   TBranch *branch=cTree->GetBranch("TOF");
757   if (!branch) { 
758     AliError("can't get the branch with the TOF clusters !");
759     return 1;
760   }
761
762   static TClonesArray dummy("AliTOFcluster",10000);
763   dummy.Clear();
764   TClonesArray *clusters=&dummy;
765   branch->SetAddress(&clusters);
766
767   cTree->GetEvent(0);
768   fN=clusters->GetEntriesFast();
769   AliInfo(Form("Number of clusters: %d",fN));
770
771   if(fClusters){
772     delete[] fClusters;
773     fClusters = NULL;
774   }
775
776   if(fN)
777     fClusters = new AliESDTOFcluster[fN];
778   else{
779     fClusters = new AliESDTOFcluster[1];
780     fN = 1;
781     return 0;
782   }
783
784   for (Int_t i=0; i<fN; i++) {
785     AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
786     Int_t ind[5];
787     ind[0]=c->GetDetInd(0);
788     ind[1]=c->GetDetInd(1);
789     ind[2]=c->GetDetInd(2);
790     ind[3]=c->GetDetInd(3);
791     ind[4]=c->GetDetInd(4);
792     Int_t calindex = AliTOFGeometry::GetIndex(ind);
793     Int_t tofLabels[3]={c->GetLabel(0),c->GetLabel(1),c->GetLabel(2)};
794     AliESDTOFcluster esdTOFclus(i,calindex,
795                                 AliTOFGeometry::TdcBinWidth()*c->GetTDC()/*ps*/,
796                                 AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()/*ps*/,
797                                 AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3/*ns*/,
798                                 tofLabels,
799                                 c->GetDeltaBC(),c->GetL0L1Latency(),
800                                 c->GetStatus(),c->GetZ(),c->GetPhi(),c->GetR());
801
802     fClusters[i] = esdTOFclus;
803
804   }
805
806   return 0;
807 }
808 //_________________________________________________________________________
809 void AliTOFtrackerV2::UnloadClusters() {
810   //--------------------------------------------------------------------
811   //This function unloads TOF clusters
812   //--------------------------------------------------------------------
813
814   // don't delete TOF clusters here because they should be written
815 }
816
817 //_________________________________________________________________________
818 Int_t AliTOFtrackerV2::FindClusterIndex(Double_t z) const {
819   //--------------------------------------------------------------------
820   // This function returns the index of the nearest cluster 
821   //--------------------------------------------------------------------
822   if (fN==0) return 0;
823   if (z <= fClusters[0].GetZ()) return 0;
824   if (z > fClusters[fN-1].GetZ()) return fN;
825   Int_t b=0, e=fN-1, m=(b+e)/2;
826   for (; b<e; m=(b+e)/2) {
827     if (z > fClusters[m].GetZ()) b=m+1;
828
829     else e=m; 
830   }
831   return m;
832 }
833
834 //_________________________________________________________________________
835 Bool_t AliTOFtrackerV2::GetTrackPoint(Int_t index, AliTrackPoint& p) const
836 {
837   // Get track space point with index i
838   // Coordinates are in the global system
839   AliESDTOFcluster *cl = &(fClusters[index]);
840   Float_t xyz[3];
841   xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
842   xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
843   xyz[2] = cl->GetZ();
844   Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
845   Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
846   Int_t tofChannel=cl->GetTOFchannel();
847   Int_t ind[5]; fGeom->GetVolumeIndices(tofChannel,ind);
848   Float_t tiltangle = AliTOFGeometry::GetAngles(ind[1],ind[2])*TMath::DegToRad();
849   Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
850   Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
851   Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
852   Float_t cov[6];
853   cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
854   cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
855   cov[2] = -cosphi*sinth*costh*sigmaz2;
856   cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
857   cov[4] = -sinphi*sinth*costh*sigmaz2;
858   cov[5] = costh*costh*sigmaz2;
859   p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
860
861   // Detector numbering scheme
862   Int_t nSector = AliTOFGeometry::NSectors();
863   Int_t nPlate  = AliTOFGeometry::NPlates();
864   Int_t nStripA = AliTOFGeometry::NStripA();
865   Int_t nStripB = AliTOFGeometry::NStripB();
866   Int_t nStripC = AliTOFGeometry::NStripC();
867
868   Int_t isector = ind[0];//cl->GetDetInd(0);
869   if (isector >= nSector)
870     AliError(Form("Wrong sector number in TOF (%d) !",isector));
871   Int_t iplate = ind[1];//cl->GetDetInd(1);
872   if (iplate >= nPlate)
873     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
874   Int_t istrip = ind[2];//cl->GetDetInd(2);
875
876   Int_t stripOffset = 0;
877   switch (iplate) {
878   case 0:
879     stripOffset = 0;
880     break;
881   case 1:
882     stripOffset = nStripC;
883     break;
884   case 2:
885     stripOffset = nStripC+nStripB;
886     break;
887   case 3:
888     stripOffset = nStripC+nStripB+nStripA;
889     break;
890   case 4:
891     stripOffset = nStripC+nStripB+nStripA+nStripB;
892     break;
893   default:
894     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
895     break;
896   };
897
898   Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
899                stripOffset +
900                istrip;
901   UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
902   p.SetVolumeID((UShort_t)volid);
903   return kTRUE;
904 }
905 //_________________________________________________________________________
906
907 void AliTOFtrackerV2::FillClusterArray(TObjArray* arr) const
908 {
909   //
910   // Returns the TOF cluster array
911   //
912
913   if (fN==0)
914     arr = 0x0;
915   else
916     for (Int_t i=0; i<fN; ++i) arr->Add(&(fClusters[i]));
917
918 }
919 //_________________________________________________________________________
920 Float_t AliTOFtrackerV2::CorrectTimeWalk( Float_t dist, Float_t tof) const {
921
922   //dummy, for the moment
923   Float_t tofcorr=0.;
924   if(dist<AliTOFGeometry::ZPad()*0.5){
925     tofcorr=tof;
926     //place here the actual correction
927   }else{
928     tofcorr=tof; 
929   } 
930   return tofcorr;
931 }
932 //_________________________________________________________________________
933 void AliTOFtrackerV2::Clusterize(){
934   Int_t detId[5];
935   for(Int_t i=0; i < fN-1;i++){
936     AliESDTOFcluster *c1=&(fClusters[i]);
937     if(!c1->GetStatus()) continue;
938
939     Int_t chan1 = c1->GetTOFchannel();
940     AliTOFGeometry::GetVolumeIndices(chan1, detId); // Get volume index from channel index
941
942     Int_t ieta = detId[2]/*strip*/*2 + detId[3]/*pad Z*/;
943     if(detId[1]/*module*/ == 0) ieta += 0;
944     else if(detId[1] == 1) ieta += 38;
945     else if(detId[1] == 2) ieta += 76;
946     else if(detId[1] == 3) ieta += 106;
947     else if(detId[1] == 4) ieta += 144;
948     Int_t iphi = detId[0]/*phi sector*/*48 + detId[4]/*pad x*/;
949     
950
951     for(Int_t j=i+1; j < fN;j++){
952       AliESDTOFcluster *c2=&(fClusters[j]);
953       if(!c2->GetStatus()) continue;
954
955       Int_t chan2 = c2->GetTOFchannel();
956
957       // check if the two TOF hits are in the same strip
958       if(chan1/96 != chan2/96) continue;
959
960       AliTOFGeometry::GetVolumeIndices(chan2, detId); // Get volume index from channel index
961       Int_t ieta2 = detId[2]/*strip*/*2 + detId[3]/*pad Z*/;
962       if(detId[1]/*module*/ == 0) ieta2 += 0;
963       else if(detId[1] == 1) ieta2 += 38;
964       else if(detId[1] == 2) ieta2 += 76;
965       else if(detId[1] == 3) ieta2 += 106;
966       else if(detId[1] == 4) ieta2 += 144;
967       Int_t iphi2 = detId[0]/*phi sector*/*48 + detId[4]/*pad x*/;
968       
969       // check if the fired pad are close in space
970       if(TMath::Abs(iphi-iphi2)>1 || TMath::Abs(ieta-ieta2)>1) continue;
971
972       // check if the TOF time are close enough to be merged
973       if(TMath::Abs(c1->GetTime() - c2->GetTime()) > 500/*in ps*/) continue;
974
975       // merge them
976       Int_t label[3] = {c2->GetLabel(0),c2->GetLabel(1),c2->GetLabel(2)};
977       fClusters[i].AddTOFhit(c2->GetClusterIndex(),chan2,c2->GetTime(),c2->GetTimeRaw(),c2->GetTOT(),label,
978                                        c2->GetDeltaBC(),c2->GetL0L1Latency(),1,c2->GetZ(),c2->GetPhi(),c2->GetR());
979       
980       c2->SetStatus(0); // only the merged one should be used
981       j = fN; // cluster "i" merged go to the next one ("i+1")
982     }
983   }
984
985   // second step of clusterization
986   for(Int_t i=0; i < fN-1;i++){
987     AliESDTOFcluster *c1=&(fClusters[i]);
988     if(!c1->GetStatus()) continue;
989
990     Int_t chan1 = c1->GetTOFchannel(0);
991     AliTOFGeometry::GetVolumeIndices(chan1, detId); // Get volume index from channel index
992
993     Int_t ieta = detId[2]/*strip*/*2 + detId[3]/*pad Z*/;
994     if(detId[1]/*module*/ == 0) ieta += 0;
995     else if(detId[1] == 1) ieta += 38;
996     else if(detId[1] == 2) ieta += 76;
997     else if(detId[1] == 3) ieta += 106;
998     else if(detId[1] == 4) ieta += 144;
999     Int_t iphi = detId[0]/*phi sector*/*48 + detId[4]/*pad x*/;
1000
1001     Int_t ieta2,iphi2,chan2=chan1;
1002     if(c1->GetNTOFhits() > 1){
1003       chan2 = c1->GetTOFchannel(1);
1004       AliTOFGeometry::GetVolumeIndices(chan2, detId); // Get volume index from channel index
1005
1006       ieta2 = detId[2]/*strip*/*2 + detId[3]/*pad Z*/;
1007       if(detId[1]/*module*/ == 0) ieta2 += 0;
1008       else if(detId[1] == 1) ieta2 += 38;
1009       else if(detId[1] == 2) ieta2 += 76;
1010       else if(detId[1] == 3) ieta2 += 106;
1011       else if(detId[1] == 4) ieta2 += 144;
1012       iphi2 = detId[0]/*phi sector*/*48 + detId[4]/*pad x*/;
1013     }
1014     else{
1015       iphi2=iphi;
1016       ieta2=ieta;
1017     }
1018
1019     for(Int_t j=i+1; j < i;j++){
1020       AliESDTOFcluster *c2=&(fClusters[j]);
1021       if(!c2->GetStatus()) continue;
1022
1023       Int_t chan3 = c2->GetTOFchannel();
1024
1025       // check if the two TOF hits are in the same strip
1026       if(chan1/96 != chan3/96) continue;
1027
1028       AliTOFGeometry::GetVolumeIndices(chan3, detId); // Get volume index from channel index
1029       Int_t ieta3 = detId[2]/*strip*/*2 + detId[3]/*pad Z*/;
1030       if(detId[1]/*module*/ == 0) ieta3 += 0;
1031       else if(detId[1] == 1) ieta3 += 38;
1032       else if(detId[1] == 2) ieta3 += 76;
1033       else if(detId[1] == 3) ieta3 += 106;
1034       else if(detId[1] == 4) ieta3 += 144;
1035       Int_t iphi3 = detId[0]/*phi sector*/*48 + detId[4]/*pad x*/;
1036       
1037       // check if the fired pad are close in space
1038       if((TMath::Abs(iphi-iphi3)>1 && TMath::Abs(iphi2-iphi3)>1) || (TMath::Abs(ieta-ieta3)>1 && TMath::Abs(ieta2-ieta3)>1)) 
1039 continue;
1040       
1041       // check if the TOF time are close enough to be merged
1042       if(TMath::Abs(c1->GetTime() - c2->GetTime()) > 500/*in ps*/) continue;
1043       
1044       // merge them
1045       Int_t label[3] = {c2->GetLabel(0),c2->GetLabel(1),c2->GetLabel(2)};
1046       fClusters[i].AddTOFhit(c2->GetClusterIndex(),chan2,c2->GetTime(),c2->GetTimeRaw(),c2->GetTOT(),label,
1047                                        c2->GetDeltaBC(),c2->GetL0L1Latency(),1,c2->GetZ(),c2->GetPhi(),c2->GetR());
1048
1049       if(c2->GetNTOFhits() > 1){ // in case also the second cluster has two hits
1050         Int_t label2[3] = {c2->GetLabel(0,1),c2->GetLabel(1,1),c2->GetLabel(2,1)};
1051         fClusters[i].AddTOFhit(c2->GetClusterIndex(1),c2->GetTOFchannel(1),c2->GetTime(1),c2->GetTimeRaw(1),
1052                                          c2->GetTOT(1),label2,c2->GetDeltaBC(2),c2->GetL0L1Latency(2),1,c2->GetZ(),
1053                                          c2->GetPhi(),c2->GetR());
1054       }
1055
1056       c1->SetStatus(0); // only the merged one should be used
1057       j = fN; // cluster "i" merged go to the next one ("i+1")
1058     }
1059   }  
1060 }
1061