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