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