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