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