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