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