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