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