]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFtracker.cxx
Coding conventions.
[u/mrichter/AliRoot.git] / TOF / AliTOFtracker.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 // AliTOFtracker 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 <TGeoManager.h>
34 #include <TTree.h>
35 #include <TFile.h>
36 #include <TH2F.h>
37
38 #include "AliAlignObj.h"
39 #include "AliESDtrack.h"
40 #include "AliESD.h"
41 #include "AliLog.h"
42 #include "AliTrackPointArray.h"
43
44 #include "AliTOFcalib.h"
45 #include "AliTOFRecoParam.h"
46 #include "AliTOFcluster.h"
47 #include "AliTOFGeometryV5.h"
48 #include "AliTOFtracker.h"
49 #include "AliTOFtrack.h"
50
51 extern TGeoManager *gGeoManager;
52 extern TROOT *gROOT;
53
54 ClassImp(AliTOFtracker)
55
56 //_____________________________________________________________________________
57 AliTOFtracker::AliTOFtracker():
58   fRecoParam(0x0),
59   fGeom(0x0),
60   fPid(0x0),
61   fN(0),
62   fNseeds(0),
63   fNseedsTOF(0),
64   fngoodmatch(0),
65   fnbadmatch(0),
66   fnunmatch(0),
67   fnmatch(0),
68   fTracks(0x0),
69   fSeeds(0x0),
70   fHDigClusMap(0x0),
71   fHDigNClus(0x0),
72   fHDigClusTime(0x0),
73   fHDigClusToT(0x0),
74   fHRecNClus(0x0),
75   fHRecDist(0x0),
76   fHRecSigYVsP(0x0),
77   fHRecSigZVsP(0x0),
78   fHRecSigYVsPWin(0x0),
79   fHRecSigZVsPWin(0x0),
80   fCalTree(0x0),
81   fIch(-1),
82   fToT(-1.),
83   fTime(-1.),
84   fExpTimePi(-1.),
85   fExpTimeKa(-1.),
86   fExpTimePr(-1.)
87  { 
88   //AliTOFtracker main Ctor
89    
90    // Gettimg the geometry 
91    fGeom=new AliTOFGeometryV5();
92    // Read the reconstruction parameters from the OCDB
93    AliTOFcalib *calib = new AliTOFcalib(fGeom);
94    fRecoParam = (AliTOFRecoParam*)calib->ReadRecParFromCDB("TOF/Calib",-1);
95    if(fRecoParam->GetApplyPbPbCuts())fRecoParam=fRecoParam->GetPbPbparam();
96    Double_t parPID[2];   
97    parPID[0]=fRecoParam->GetTimeResolution();
98    parPID[1]=fRecoParam->GetTimeNSigma();
99    fPid=new AliTOFpidESD(parPID);
100    InitCheckHists();
101
102 }
103 //_____________________________________________________________________________
104 AliTOFtracker::AliTOFtracker(const AliTOFtracker &t):
105   AliTracker(),
106   fRecoParam(0x0),
107   fGeom(0x0),
108   fPid(0x0),
109   fN(0),
110   fNseeds(0),
111   fNseedsTOF(0),
112   fngoodmatch(0),
113   fnbadmatch(0),
114   fnunmatch(0),
115   fnmatch(0),
116   fTracks(0x0),
117   fSeeds(0x0),
118   fHDigClusMap(0x0),
119   fHDigNClus(0x0),
120   fHDigClusTime(0x0),
121   fHDigClusToT(0x0),
122   fHRecNClus(0x0),
123   fHRecDist(0x0),
124   fHRecSigYVsP(0x0),
125   fHRecSigZVsP(0x0),
126   fHRecSigYVsPWin(0x0),
127   fHRecSigZVsPWin(0x0),
128   fCalTree(0x0),
129   fIch(-1),
130   fToT(-1.),
131   fTime(-1.),
132   fExpTimePi(-1.),
133   fExpTimeKa(-1.),
134   fExpTimePr(-1.)
135  { 
136   //AliTOFtracker copy Ctor
137
138   fNseeds=t.fNseeds;
139   fNseeds=t.fNseeds;
140   fNseedsTOF=t.fNseedsTOF;
141   fngoodmatch=t.fngoodmatch;
142   fnbadmatch=t.fnbadmatch;
143   fnunmatch=t.fnunmatch;
144   fnmatch=t.fnmatch;
145   fRecoParam=t.fRecoParam;
146   fGeom=t.fGeom;
147   fPid=t.fPid;
148   fSeeds=t.fSeeds;
149   fTracks=t.fTracks;
150   fN=t.fN;
151 }
152
153 //_____________________________________________________________________________
154 AliTOFtracker& AliTOFtracker::operator=(const AliTOFtracker &t)
155
156   //AliTOFtracker assignment operator
157
158   this->fNseeds=t.fNseeds;
159   this->fNseedsTOF=t.fNseedsTOF;
160   this->fngoodmatch=t.fngoodmatch;
161   this->fnbadmatch=t.fnbadmatch;
162   this->fnunmatch=t.fnunmatch;
163   this->fnmatch=t.fnmatch;
164   this->fRecoParam = t.fRecoParam;
165   this->fGeom = t.fGeom;
166   this->fPid = t.fPid;
167   this->fSeeds=t.fSeeds;
168   this->fTracks=t.fTracks;
169   this->fN=t.fN;
170   return *this;
171
172 }
173 //_____________________________________________________________________________
174 AliTOFtracker::~AliTOFtracker() {
175   //
176   // Dtor
177   //
178
179   SaveCheckHists();
180
181   delete fRecoParam; 
182   delete fGeom; 
183   delete fPid; 
184   delete fHDigClusMap;
185   delete fHDigNClus;
186   delete fHDigClusTime;
187   delete fHDigClusToT;
188   delete fHRecNClus;
189   delete fHRecDist;
190   delete fHRecSigYVsP;
191   delete fHRecSigZVsP;
192   delete fHRecSigYVsPWin;
193   delete fHRecSigZVsPWin;
194   delete fCalTree;
195 }
196 //_____________________________________________________________________________
197 Int_t AliTOFtracker::PropagateBack(AliESD* event) {
198   //
199   // Gets seeds from ESD event and Match with TOF Clusters
200   //
201
202
203   //Initialise some counters
204
205   fNseeds=0;
206   fNseedsTOF=0;
207   fngoodmatch=0;
208   fnbadmatch=0;
209   fnunmatch=0;
210   fnmatch=0;
211
212   Int_t ntrk=event->GetNumberOfTracks();
213   fNseeds = ntrk;
214   fSeeds= new TClonesArray("AliESDtrack",ntrk);
215   TClonesArray &aESDTrack = *fSeeds;
216
217
218   //Load ESD tracks into a local Array of ESD Seeds
219
220   for (Int_t i=0; i<fNseeds; i++) {
221     AliESDtrack *t=event->GetTrack(i);
222     new(aESDTrack[i]) AliESDtrack(*t);
223   }
224
225   //Prepare ESD tracks candidates for TOF Matching
226   CollectESD();
227
228   //First Step with Strict Matching Criterion
229   MatchTracks(kFALSE);
230
231   //Second Step with Looser Matching Criterion
232   MatchTracks(kTRUE);
233
234   AliInfo(Form("Number of matched tracks: %d",fnmatch));
235   AliInfo(Form("Number of good matched tracks: %d",fngoodmatch));
236   AliInfo(Form("Number of bad  matched tracks: %d",fnbadmatch));
237
238   //Update the matched ESD tracks
239
240   for (Int_t i=0; i<ntrk; i++) {
241     AliESDtrack *t=event->GetTrack(i);
242     AliESDtrack *seed =(AliESDtrack*)fSeeds->UncheckedAt(i);
243     if(seed->GetTOFsignal()>0){
244       t->SetTOFsignal(seed->GetTOFsignal());
245       t->SetTOFcluster(seed->GetTOFcluster());
246       t->SetTOFsignalToT(seed->GetTOFsignalToT());
247       t->SetTOFCalChannel(seed->GetTOFCalChannel());
248       Int_t tlab[3]; seed->GetTOFLabel(tlab);    
249       t->SetTOFLabel(tlab);
250       AliTOFtrack *track = new AliTOFtrack(*seed); 
251       t->UpdateTrackParams(track,AliESDtrack::kTOFout);   
252       delete track;
253     }
254   }
255
256   //Handle Time Zero information
257
258   Double_t timeZero=0.;
259   Double_t timeZeroMax=99999.;
260   Bool_t usetimeZero     = fRecoParam->UseTimeZero();
261   Bool_t timeZeroFromT0  = fRecoParam->GetTimeZerofromT0();
262   Bool_t timeZeroFromTOF = fRecoParam->GetTimeZerofromTOF();
263
264   AliDebug(1,Form("Use Time Zero?: %d",usetimeZero));
265   AliDebug(1,Form("Time Zero from T0? : %d",timeZeroFromT0));
266   AliDebug(1,Form("Time Zero From TOF? : %d",timeZeroFromTOF));
267
268   if(usetimeZero){
269     if(timeZeroFromT0){
270       timeZero=GetTimeZerofromT0(event); 
271     }
272     if(timeZeroFromTOF && (timeZero>timeZeroMax || !timeZeroFromT0)){
273       timeZero=GetTimeZerofromTOF(event); 
274     }
275   }
276   AliDebug(2,Form("time Zero used in PID: %f",timeZero));
277   //Make TOF PID
278   fPid->MakePID(event,timeZero);
279
280   if (fSeeds) {
281     fSeeds->Delete();
282     delete fSeeds;
283     fSeeds = 0x0;
284   }
285   if (fTracks) {
286     fTracks->Delete();
287     delete fTracks;
288     fTracks = 0x0;
289   }
290   return 0;
291   
292 }
293 //_________________________________________________________________________
294 void AliTOFtracker::CollectESD() {
295    //prepare the set of ESD tracks to be matched to clusters in TOF
296  
297   fTracks= new TClonesArray("AliTOFtrack");
298   TClonesArray &aTOFTrack = *fTracks;
299   for (Int_t i=0; i<fNseeds; i++) {
300
301     AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(i);
302     if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
303
304     // TRD good tracks, already propagated at 371 cm
305
306     AliTOFtrack *track = new AliTOFtrack(*t); // New
307     Double_t x = track->GetX(); //New
308
309     if (((t->GetStatus()&AliESDtrack::kTRDout)!=0 ) && 
310          ( x >= fGeom->RinTOF()) ){
311       track->SetSeedIndex(i);
312       t->UpdateTrackParams(track,AliESDtrack::kTOFout);    
313       new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
314       fNseedsTOF++;
315       delete track;
316     }
317
318     // Propagate the rest of TPCbp  
319
320     else {
321       if(track->PropagateToInnerTOF()){ 
322         track->SetSeedIndex(i);
323         t->UpdateTrackParams(track,AliESDtrack::kTOFout);    
324         new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
325         fNseedsTOF++;
326       }
327       delete track;
328     }
329   }
330
331   AliInfo(Form("Number of TOF seedds %i",fNseedsTOF));
332
333   // Sort according uncertainties on track position 
334   fTracks->Sort();
335
336 }
337 //_________________________________________________________________________
338 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
339
340
341   // Parameters used/regulating the reconstruction
342
343   static Float_t corrLen=0.75;
344   static Float_t detDepth=15.3;
345
346   Float_t dY=fGeom->XPad(); 
347   Float_t dZ=fGeom->ZPad(); 
348
349   Float_t sensRadius = fRecoParam->GetSensRadius();
350   Float_t stepSize   = fRecoParam->GetStepSize();
351   Float_t scaleFact   = fRecoParam->GetWindowScaleFact();
352   Float_t dyMax=fRecoParam->GetWindowSizeMaxY(); 
353   Float_t dzMax=fRecoParam->GetWindowSizeMaxZ();
354   Float_t dCut=fRecoParam->GetDistanceCut();
355   Double_t maxChi2=fRecoParam->GetMaxChi2();
356   Bool_t timeWalkCorr    = fRecoParam->GetTimeWalkCorr();
357   if(!mLastStep){
358     AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++ \n");
359     AliDebug(1,Form("TOF sens radius: %f",sensRadius));
360     AliDebug(1,Form("TOF step size: %f",stepSize));
361     AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
362     AliDebug(1,Form("TOF Window max dy: %f",dyMax));
363     AliDebug(1,Form("TOF Window max dz: %f",dzMax));
364     AliDebug(1,Form("TOF distance Cut: %f",dCut));
365     AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
366     AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));   
367   }
368   //Match ESD tracks to clusters in TOF
369
370
371   // Get the number of propagation steps
372
373   Int_t nSteps=(Int_t)(detDepth/stepSize);
374
375   AliTOFcalib *calib = new AliTOFcalib(fGeom);
376
377   //PH Arrays (moved outside of the loop)
378
379   Float_t * trackPos[4];
380   for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
381   Int_t * clind[6];
382   for (Int_t ii=0;ii<6;ii++) clind[ii] = new Int_t[fN];
383   
384   for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
385
386     AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
387     AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(track->GetSeedIndex());
388     if(t->GetTOFsignal()>0. ) continue;
389     AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
390
391     // Some init
392
393     Int_t         index[10000];
394     Float_t        dist[10000];
395     Float_t       distZ[10000];
396     Float_t       cxpos[10000];
397     Float_t       crecL[10000];
398     TGeoHMatrix   global[1000];
399      
400     // Determine a window around the track
401
402     Double_t x,par[5]; 
403     trackTOFin->GetExternalParameters(x,par);
404     Double_t cov[15]; 
405     trackTOFin->GetExternalCovariance(cov);
406
407     Double_t dphi=
408       scaleFact*
409       ((5*TMath::Sqrt(cov[0]) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius); 
410     Double_t dz=
411        scaleFact*
412        (5*TMath::Sqrt(cov[2]) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
413
414     Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
415     if (phi<-TMath::Pi())phi+=2*TMath::Pi();
416     if (phi>=TMath::Pi())phi-=2*TMath::Pi();
417     Double_t z=par[1];   
418
419     //upper limit on window's size.
420
421     if(dz> dzMax) dz=dzMax;
422     if(dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
423
424
425     Int_t nc=0;
426
427     // find the clusters in the window of the track
428
429     for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
430
431       AliTOFcluster *c=fClusters[k];
432       if (c->GetZ() > z+dz) break;
433       if (c->IsUsed()) continue;
434
435       if (!c->GetStatus()) continue; // skip bad channels as declared in OCDB
436
437       Double_t dph=TMath::Abs(c->GetPhi()-phi);
438       if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
439       if (TMath::Abs(dph)>dphi) continue;
440
441       Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
442       Double_t p[2]={yc, c->GetZ()};
443       Double_t cov[3]= {dY*dY/12., 0., dZ*dZ/12.};
444       if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov) > maxChi2)continue;
445
446       clind[0][nc] = c->GetDetInd(0);
447       clind[1][nc] = c->GetDetInd(1);
448       clind[2][nc] = c->GetDetInd(2);
449       clind[3][nc] = c->GetDetInd(3);
450       clind[4][nc] = c->GetDetInd(4);
451       clind[5][nc] = k;      
452       Char_t path[100];
453       Int_t ind[5];
454       ind[0]=clind[0][nc];
455       ind[1]=clind[1][nc];
456       ind[2]=clind[2][nc];
457       ind[3]=clind[3][nc];
458       ind[4]=clind[4][nc];
459       fGeom->GetVolumePath(ind,path);
460       gGeoManager->cd(path);
461       global[nc] = *gGeoManager->GetCurrentMatrix();
462       nc++;
463     }
464
465     //start fine propagation 
466
467     Int_t nStepsDone = 0;
468     for( Int_t istep=0; istep<nSteps; istep++){ 
469
470       Float_t xs=fGeom->RinTOF()+istep*0.1;
471       Double_t ymax=xs*TMath::Tan(0.5*fGeom->GetAlpha());
472
473       Bool_t skip=kFALSE;
474       Double_t ysect=trackTOFin->GetYat(xs,skip);
475       if (skip) break;
476       if (ysect > ymax) {
477         if (!trackTOFin->Rotate(fGeom->GetAlpha())) {
478           break;
479         }
480       } else if (ysect <-ymax) {
481         if (!trackTOFin->Rotate(fGeom->GetAlpha())) {
482           break;
483         }
484       }
485
486       if(!trackTOFin->PropagateTo(xs)) {
487         break;
488       }
489
490       nStepsDone++;
491
492       // store the running point (Globalrf) - fine propagation     
493
494       Double_t r[3];
495       trackTOFin->GetXYZ(r);
496       trackPos[0][istep]= (Float_t) r[0];
497       trackPos[1][istep]= (Float_t) r[1];
498       trackPos[2][istep]= (Float_t) r[2];   
499       trackPos[3][istep]= trackTOFin->GetIntegratedLength();
500     }
501
502
503     Int_t nfound = 0;
504     for (Int_t istep=0; istep<nStepsDone; istep++) {
505
506       Bool_t isInside =kFALSE;
507       Float_t ctrackPos[3];     
508
509       ctrackPos[0]= trackPos[0][istep];
510       ctrackPos[1]= trackPos[1][istep];
511       ctrackPos[2]= trackPos[2][istep];
512
513       //now see whether the track matches any of the TOF clusters            
514
515       for (Int_t i=0; i<nc; i++){
516         Int_t cind[5];
517         cind[0]= clind[0][i];
518         cind[1]= clind[1][i];
519         cind[2]= clind[2][i];
520         cind[3]= clind[3][i];
521         cind[4]= clind[4][i];
522         Bool_t accept = kFALSE;
523         if( mLastStep)accept = (fGeom->DistanceToPad(cind,global[i],ctrackPos)<dCut);
524         if(!mLastStep)accept = (fGeom->IsInsideThePad(cind,global[i],ctrackPos));
525         if(accept){
526           if(!mLastStep)isInside=kTRUE;
527           Float_t dist3d[3];
528           dist[nfound]=fGeom->DistanceToPad(cind,global[i],ctrackPos,dist3d);
529           distZ[nfound]=dist3d[2];
530           crecL[nfound]=trackPos[3][istep];
531           index[nfound]=clind[5][i]; // store cluster id            
532           cxpos[nfound]=fGeom->RinTOF()+istep*0.1; //store prop.radius
533           nfound++;
534           if(isInside)break;
535         }//end if accept
536       } //end for on the clusters
537
538
539       if(isInside)break;
540     } //end for on the steps     
541
542
543
544     if (nfound == 0 ) {
545       fnunmatch++;
546       delete trackTOFin;
547       continue;
548     }
549     
550     fnmatch++;
551
552     // now choose the cluster to be matched with the track.
553
554     Int_t idclus=0;
555     Float_t  recL = 0.;
556     Float_t  xpos=0.;
557     Float_t  mindist=1000.;
558     Float_t  mindistZ=0.;
559     for (Int_t iclus= 0; iclus<nfound;iclus++){
560       if (dist[iclus]< mindist){
561         mindist = dist[iclus];
562         mindistZ = distZ[iclus];
563         xpos = cxpos[iclus];
564         idclus =index[iclus]; 
565         recL=crecL[iclus]+corrLen*0.5;
566       }
567     }
568
569     AliTOFcluster *c=fClusters[idclus];
570     c->Use(); 
571
572     // Track length correction for matching Step 2 
573
574     if(mLastStep){
575       Float_t rc=TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
576       Float_t rt=TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
577                              +trackPos[1][70]*trackPos[1][70]
578                              +trackPos[2][70]*trackPos[2][70]);
579       Float_t dlt=rc-rt;      
580       recL=trackPos[3][70]+dlt;
581     }    
582
583     if (
584         (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
585         ||
586         (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
587         ||
588         (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
589         ) {
590       fngoodmatch++;
591
592        AliDebug(2,Form(" track label good %5i",trackTOFin->GetLabel()));
593
594     }
595     else{
596       fnbadmatch++;
597
598       AliDebug(2,Form(" track label  bad %5i",trackTOFin->GetLabel()));
599
600     }
601
602     delete trackTOFin;
603
604     //  Store quantities to be used in the TOF Calibration
605     Float_t tToT=fGeom->ToTBinWidth()*c->GetToT()*1E-3; // in ns
606     t->SetTOFsignalToT(tToT);
607     Int_t ind[5];
608     ind[0]=c->GetDetInd(0);
609     ind[1]=c->GetDetInd(1);
610     ind[2]=c->GetDetInd(2);
611     ind[3]=c->GetDetInd(3);
612     ind[4]=c->GetDetInd(4);
613     Int_t calindex = calib->GetIndex(ind);
614     t->SetTOFCalChannel(calindex);
615
616     // keep track of the track labels in the matched cluster
617     Int_t tlab[3];
618     tlab[0]=c->GetLabel(0);
619     tlab[1]=c->GetLabel(1);
620     tlab[2]=c->GetLabel(2);
621     
622     Double_t tof=fGeom->TdcBinWidth()*c->GetTDC()+32; // in ps
623     if(timeWalkCorr)tof=CorrectTimeWalk(mindistZ,tof);
624     t->SetTOFsignal(tof);
625     t->SetTOFcluster(idclus); // pointing to the recPoints tree
626     Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time);
627     Double_t mom=t->GetP();
628     for(Int_t j=0;j<AliPID::kSPECIES;j++){
629       Double_t mass=AliPID::ParticleMass(j);
630       time[j]+=(recL-trackPos[3][0])/3e-2*TMath::Sqrt(mom*mom+mass*mass)/mom;
631     }
632
633     AliTOFtrack *trackTOFout = new AliTOFtrack(*t); 
634     trackTOFout->PropagateTo(xpos);
635     t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);    
636     t->SetIntegratedLength(recL);
637     t->SetIntegratedTimes(time);
638     t->SetTOFLabel(tlab);
639     // Fill Reco-QA histos for Reconstruction
640     fHRecNClus->Fill(nc);
641     fHRecDist->Fill(mindist);
642     fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
643     fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
644     fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
645     fHRecSigZVsPWin->Fill(mom,dz);
646
647     // Fill Tree for on-the-fly offline Calibration
648
649     if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 )){    
650       Float_t rawtime=fGeom->TdcBinWidth()*c->GetTDCRAW()+32; // RAW time,in ps
651       fIch=calindex;
652       fToT=tToT;
653       fTime=rawtime;
654       fExpTimePi=time[2];
655       fExpTimeKa=time[3];
656       fExpTimePr=time[4];
657       fCalTree->Fill();
658     }
659     delete trackTOFout;
660   }
661   for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
662   for (Int_t ii=0;ii<6;ii++) delete [] clind[ii];
663   delete calib;
664 }
665 //_________________________________________________________________________
666 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
667   //--------------------------------------------------------------------
668   //This function loads the TOF clusters
669   //--------------------------------------------------------------------
670
671   Int_t npadX = fGeom->NpadX();
672   Int_t npadZ = fGeom->NpadZ();
673   Int_t nStripA = fGeom->NStripA();
674   Int_t nStripB = fGeom->NStripB();
675   Int_t nStripC = fGeom->NStripC();
676
677   TBranch *branch=cTree->GetBranch("TOF");
678   if (!branch) { 
679     AliError("can't get the branch with the TOF clusters !");
680     return 1;
681   }
682
683   TClonesArray dummy("AliTOFcluster",10000), *clusters=&dummy;
684   branch->SetAddress(&clusters);
685
686   cTree->GetEvent(0);
687   Int_t nc=clusters->GetEntriesFast();
688   fHDigNClus->Fill(nc);
689
690   AliInfo(Form("Number of clusters: %d",nc));
691
692   for (Int_t i=0; i<nc; i++) {
693     AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
694     fClusters[i]=new AliTOFcluster(*c); fN++;
695
696   // Fill Digits QA histos
697  
698     Int_t isector = c->GetDetInd(0);
699     Int_t iplate = c->GetDetInd(1);
700     Int_t istrip = c->GetDetInd(2);
701     Int_t ipadX = c->GetDetInd(4);
702     Int_t ipadZ = c->GetDetInd(3);
703
704     Float_t time =(fGeom->TdcBinWidth()*c->GetTDC())*1E-3; // in ns
705     Float_t tot = (fGeom->TdcBinWidth()*c->GetToT())*1E-3;//in ns
706  
707     Int_t stripOffset = 0;
708     switch (iplate) {
709     case 0:
710       stripOffset = 0;
711       break;
712     case 1:
713       stripOffset = nStripC;
714       break;
715     case 2:
716       stripOffset = nStripC+nStripB;
717       break;
718     case 3:
719       stripOffset = nStripC+nStripB+nStripA;
720       break;
721     case 4:
722       stripOffset = nStripC+nStripB+nStripA+nStripB;
723       break;
724     default:
725       AliError(Form("Wrong plate number in TOF (%d) !",iplate));
726       break;
727     };
728     Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
729     Int_t phiindex=npadX*isector+ipadX+1;
730     fHDigClusMap->Fill(zindex,phiindex);
731     fHDigClusTime->Fill(time);
732     fHDigClusToT->Fill(tot);
733
734   }
735
736
737   return 0;
738 }
739 //_________________________________________________________________________
740 void AliTOFtracker::UnloadClusters() {
741   //--------------------------------------------------------------------
742   //This function unloads TOF clusters
743   //--------------------------------------------------------------------
744   for (Int_t i=0; i<fN; i++) {
745     delete fClusters[i];
746     fClusters[i] = 0x0;
747   }
748   fN=0;
749 }
750
751 //_________________________________________________________________________
752 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
753   //--------------------------------------------------------------------
754   // This function returns the index of the nearest cluster 
755   //--------------------------------------------------------------------
756   if (fN==0) return 0;
757   if (z <= fClusters[0]->GetZ()) return 0;
758   if (z > fClusters[fN-1]->GetZ()) return fN;
759   Int_t b=0, e=fN-1, m=(b+e)/2;
760   for (; b<e; m=(b+e)/2) {
761     if (z > fClusters[m]->GetZ()) b=m+1;
762     else e=m; 
763   }
764   return m;
765 }
766
767 //_________________________________________________________________________
768 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
769 {
770   // Get track space point with index i
771   // Coordinates are in the global system
772   AliTOFcluster *cl = fClusters[index];
773   Float_t xyz[3];
774   xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
775   xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
776   xyz[2] = cl->GetZ();
777   Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
778   Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
779   Float_t tiltangle = fGeom->GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
780   Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
781   Float_t sigmay2 = fGeom->XPad()*fGeom->XPad()/12.;
782   Float_t sigmaz2 = fGeom->ZPad()*fGeom->ZPad()/12.;
783   Float_t cov[6];
784   cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
785   cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
786   cov[2] = -cosphi*sinth*costh*sigmaz2;
787   cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
788   cov[4] = -sinphi*sinth*costh*sigmaz2;
789   cov[5] = costh*costh*sigmaz2;
790   p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
791
792   // Detector numbering scheme
793   Int_t nSector = fGeom->NSectors();
794   Int_t nPlate  = fGeom->NPlates();
795   Int_t nStripA = fGeom->NStripA();
796   Int_t nStripB = fGeom->NStripB();
797   Int_t nStripC = fGeom->NStripC();
798
799   Int_t isector = cl->GetDetInd(0);
800   if (isector >= nSector)
801     AliError(Form("Wrong sector number in TOF (%d) !",isector));
802   Int_t iplate = cl->GetDetInd(1);
803   if (iplate >= nPlate)
804     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
805   Int_t istrip = cl->GetDetInd(2);
806
807   Int_t stripOffset = 0;
808   switch (iplate) {
809   case 0:
810     stripOffset = 0;
811     break;
812   case 1:
813     stripOffset = nStripC;
814     break;
815   case 2:
816     stripOffset = nStripC+nStripB;
817     break;
818   case 3:
819     stripOffset = nStripC+nStripB+nStripA;
820     break;
821   case 4:
822     stripOffset = nStripC+nStripB+nStripA+nStripB;
823     break;
824   default:
825     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
826     break;
827   };
828
829   Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
830                stripOffset +
831                istrip;
832   UShort_t volid = AliAlignObj::LayerToVolUID(AliAlignObj::kTOF,idet);
833   p.SetVolumeID((UShort_t)volid);
834   return kTRUE;
835 }
836 //_________________________________________________________________________
837 void AliTOFtracker::InitCheckHists() {
838
839   //Init histos for Digits/Reco QA and Calibration
840
841
842   fCalTree = new TTree("CalTree", "Tree for TOF calibration");
843   fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
844   fCalTree->Branch("ToT",&fToT,"TOFToT/F");
845   fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
846   fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
847   fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
848   fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
849
850   //Digits "QA" 
851   fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);  
852   fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);  
853   fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);  
854   fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);  
855
856   //Reco "QA"
857   fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
858   fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
859   fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
860   fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
861   fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
862   fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
863 }
864
865 //_________________________________________________________________________
866 void AliTOFtracker::SaveCheckHists() {
867
868   //write histos for Digits/Reco QA and Calibration
869
870   TDirectory *dir = gDirectory;
871   TFile *logFile = 0;
872   TFile *logFileTOF = 0;
873
874   TSeqCollection *list = gROOT->GetListOfFiles();
875   int n = list->GetEntries();
876   for(int i=0; i<n; i++) {
877     logFile = (TFile*)list->At(i);
878     if (strstr(logFile->GetName(), "AliESDs.root")) break;
879   }
880
881   Bool_t isThere=kFALSE;
882   for(int i=0; i<n; i++) {
883     logFileTOF = (TFile*)list->At(i);
884     if (strstr(logFileTOF->GetName(), "TOFQA.root")){
885       isThere=kTRUE;
886       break;
887     } 
888   }
889    
890   logFile->cd();
891   fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
892   fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
893   fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
894   fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
895   fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
896   fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
897   fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
898   fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
899   fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
900   fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
901   fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
902   logFile->Flush();  
903
904   if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
905   logFileTOF->cd(); 
906   fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
907   fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
908   fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
909   fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
910   fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
911   fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
912   fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
913   fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
914   fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
915   fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
916   fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
917   logFileTOF->Flush();  
918
919   dir->cd();
920   }
921 //_________________________________________________________________________
922 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) {
923
924   //dummy, for the moment
925   Float_t tofcorr=0.;
926   if(dist<fGeom->ZPad()*0.5){
927     tofcorr=tof;
928     //place here the actual correction
929   }else{
930     tofcorr=tof; 
931   } 
932   return tofcorr;
933 }
934 //_________________________________________________________________________
935 Float_t AliTOFtracker::GetTimeZerofromT0(AliESD *event) {
936
937   //Returns TimeZero as measured by T0 detector
938
939   return event->GetT0();
940 }
941 //_________________________________________________________________________
942 Float_t AliTOFtracker::GetTimeZerofromTOF(AliESD *event) {
943
944   //dummy, for the moment. T0 algorithm using tracks on TOF
945   {
946     //place T0 algo here...
947   }
948   return 0.;
949 }