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