]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFtracker.cxx
Reduced QA output (Yves)
[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 dist3d[3];
477       accept = kFALSE;
478       for (Int_t i=0; i<nc; i++) {
479         isInside = fGeom->IsInsideThePad(global[i],ctrackPos,dist3d);
480
481         if( mLastStep){
482           Float_t xLoc = dist3d[0];
483           Float_t rLoc = TMath::Sqrt(dist3d[1]*dist3d[1]+dist3d[2]*dist3d[2]);
484           accept = (TMath::Abs(xLoc)<padDepth*0.5 && rLoc<dCut);
485         }
486         else {
487           accept = isInside;
488         }
489         if (accept) {
490           dist[nfound]  = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[1]*dist3d[1]+dist3d[2]*dist3d[2]);
491           distZ[nfound] = dist3d[2];
492           crecL[nfound] = trackPos[3][istep];
493           index[nfound] = clind[i]; // store cluster id             
494           cxpos[nfound] = AliTOFGeometry::RinTOF()+istep*0.1; //store prop.radius
495           nfound++;
496           if(accept &&!mLastStep)break;
497         }//end if accept
498
499       } //end for on the clusters
500       if(accept &&!mLastStep)break;
501     } //end for on the steps     
502
503
504
505     if (nfound == 0 ) {
506       fnunmatch++;
507       delete trackTOFin;
508       continue;
509     }
510     
511     fnmatch++;
512
513     // now choose the cluster to be matched with the track.
514
515     Int_t idclus=-1;
516     Float_t  recL = 0.;
517     Float_t  xpos=0.;
518     Float_t  mindist=1000.;
519     Float_t  mindistZ=0.;
520     for (Int_t iclus= 0; iclus<nfound;iclus++){
521       if (dist[iclus]< mindist){
522         mindist = dist[iclus];
523         mindistZ = distZ[iclus];
524         xpos = cxpos[iclus];
525         idclus =index[iclus]; 
526         recL=crecL[iclus]+corrLen*0.5;
527       }
528     }
529
530
531     AliTOFcluster *c=fClusters[idclus];
532
533     AliDebug(2, Form("%7i     %7i     %10i     %10i  %10i  %10i      %7i",
534                      iseed,
535                      fnmatch-1,
536                      TMath::Abs(trackTOFin->GetLabel()),
537                      c->GetLabel(0), c->GetLabel(1), c->GetLabel(2),
538                      idclus)); // AdC
539
540     c->Use(); 
541
542     // Track length correction for matching Step 2 
543
544     if(mLastStep){
545       Float_t rc=TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
546       Float_t rt=TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
547                              +trackPos[1][70]*trackPos[1][70]
548                              +trackPos[2][70]*trackPos[2][70]);
549       Float_t dlt=rc-rt;      
550       recL=trackPos[3][70]+dlt;
551     }    
552
553     if (
554         (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
555         ||
556         (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
557         ||
558         (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
559         ) {
560       fngoodmatch++;
561
562        AliDebug(2,Form(" track label good %5i",trackTOFin->GetLabel()));
563
564     }
565     else{
566       fnbadmatch++;
567
568       AliDebug(2,Form(" track label  bad %5i",trackTOFin->GetLabel()));
569
570     }
571
572     delete trackTOFin;
573
574     //  Store quantities to be used in the TOF Calibration
575     Float_t tToT=AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3; // in ns
576     t->SetTOFsignalToT(tToT);
577     Float_t rawTime=AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()+32; // RAW time,in ps
578     t->SetTOFsignalRaw(rawTime);
579     t->SetTOFsignalDz(mindistZ);
580     AliDebug(2,Form(" Setting TOF raw time: %f, z distance: %f time: %f = ",rawTime,mindistZ));    
581     Int_t ind[5];
582     ind[0]=c->GetDetInd(0);
583     ind[1]=c->GetDetInd(1);
584     ind[2]=c->GetDetInd(2);
585     ind[3]=c->GetDetInd(3);
586     ind[4]=c->GetDetInd(4);
587     Int_t calindex = AliTOFGeometry::GetIndex(ind);
588     t->SetTOFCalChannel(calindex);
589
590     // keep track of the track labels in the matched cluster
591     Int_t tlab[3];
592     tlab[0]=c->GetLabel(0);
593     tlab[1]=c->GetLabel(1);
594     tlab[2]=c->GetLabel(2);
595     AliDebug(2,Form(" tdc time of the matched track %i = ",c->GetTDC()));    
596     Double_t tof=AliTOFGeometry::TdcBinWidth()*c->GetTDC()+32; // in ps
597     AliDebug(2,Form(" tof time of the matched track: %f = ",tof));    
598     Double_t tofcorr=tof;
599     if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
600     AliDebug(2,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));    
601     //Set TOF time signal and pointer to the matched cluster
602     t->SetTOFsignal(tofcorr);
603     t->SetTOFcluster(idclus); // pointing to the recPoints tree
604
605     //Tracking info
606     Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time);
607     Double_t mom=t->GetP();
608     AliDebug(1,Form("Momentum for track %d -> %f", iseed,mom));
609     for(Int_t j=0;j<AliPID::kSPECIES;j++){
610       Double_t mass=AliPID::ParticleMass(j);
611       time[j]+=(recL-trackPos[3][0])/2.99792458e-2*TMath::Sqrt(mom*mom+mass*mass)/mom;
612     }
613
614     AliTOFtrack *trackTOFout = new AliTOFtrack(*t); 
615     trackTOFout->PropagateTo(xpos);
616
617     // Fill the track residual histograms.
618     FillResiduals(trackTOFout,c,kFALSE);
619
620     t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);    
621     t->SetIntegratedLength(recL);
622     t->SetIntegratedTimes(time);
623     t->SetTOFLabel(tlab);
624
625  
626     // Fill Reco-QA histos for Reconstruction
627     fHRecNClus->Fill(nc);
628     fHRecDist->Fill(mindist);
629     fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
630     fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
631     fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
632     fHRecSigZVsPWin->Fill(mom,dz);
633
634     // Fill Tree for on-the-fly offline Calibration
635
636     if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 )){    
637       fIch=calindex;
638       fToT=tToT;
639       fTime=rawTime;
640       fExpTimePi=time[2];
641       fExpTimeKa=time[3];
642       fExpTimePr=time[4];
643       fCalTree->Fill();
644     }
645     delete trackTOFout;
646   }
647   for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
648   delete [] clind;
649  
650 }
651 //_________________________________________________________________________
652 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
653   //--------------------------------------------------------------------
654   //This function loads the TOF clusters
655   //--------------------------------------------------------------------
656
657   Int_t npadX = AliTOFGeometry::NpadX();
658   Int_t npadZ = AliTOFGeometry::NpadZ();
659   Int_t nStripA = AliTOFGeometry::NStripA();
660   Int_t nStripB = AliTOFGeometry::NStripB();
661   Int_t nStripC = AliTOFGeometry::NStripC();
662
663   TBranch *branch=cTree->GetBranch("TOF");
664   if (!branch) { 
665     AliError("can't get the branch with the TOF clusters !");
666     return 1;
667   }
668
669   static TClonesArray dummy("AliTOFcluster",10000);
670   dummy.Clear();
671   TClonesArray *clusters=&dummy;
672   branch->SetAddress(&clusters);
673
674   cTree->GetEvent(0);
675   Int_t nc=clusters->GetEntriesFast();
676   fHDigNClus->Fill(nc);
677
678   AliInfo(Form("Number of clusters: %d",nc));
679
680   for (Int_t i=0; i<nc; i++) {
681     AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
682     fClusters[i]=new AliTOFcluster(*c); fN++;
683
684   // Fill Digits QA histos
685  
686     Int_t isector = c->GetDetInd(0);
687     Int_t iplate = c->GetDetInd(1);
688     Int_t istrip = c->GetDetInd(2);
689     Int_t ipadX = c->GetDetInd(4);
690     Int_t ipadZ = c->GetDetInd(3);
691
692     Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
693     Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
694  
695     Int_t stripOffset = 0;
696     switch (iplate) {
697     case 0:
698       stripOffset = 0;
699       break;
700     case 1:
701       stripOffset = nStripC;
702       break;
703     case 2:
704       stripOffset = nStripC+nStripB;
705       break;
706     case 3:
707       stripOffset = nStripC+nStripB+nStripA;
708       break;
709     case 4:
710       stripOffset = nStripC+nStripB+nStripA+nStripB;
711       break;
712     default:
713       AliError(Form("Wrong plate number in TOF (%d) !",iplate));
714       break;
715     };
716     Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
717     Int_t phiindex=npadX*isector+ipadX+1;
718     fHDigClusMap->Fill(zindex,phiindex);
719     fHDigClusTime->Fill(time);
720     fHDigClusToT->Fill(tot);
721
722   }
723
724
725   return 0;
726 }
727 //_________________________________________________________________________
728 void AliTOFtracker::UnloadClusters() {
729   //--------------------------------------------------------------------
730   //This function unloads TOF clusters
731   //--------------------------------------------------------------------
732   for (Int_t i=0; i<fN; i++) {
733     delete fClusters[i];
734     fClusters[i] = 0x0;
735   }
736   fN=0;
737 }
738
739 //_________________________________________________________________________
740 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
741   //--------------------------------------------------------------------
742   // This function returns the index of the nearest cluster 
743   //--------------------------------------------------------------------
744   if (fN==0) return 0;
745   if (z <= fClusters[0]->GetZ()) return 0;
746   if (z > fClusters[fN-1]->GetZ()) return fN;
747   Int_t b=0, e=fN-1, m=(b+e)/2;
748   for (; b<e; m=(b+e)/2) {
749     if (z > fClusters[m]->GetZ()) b=m+1;
750     else e=m; 
751   }
752   return m;
753 }
754
755 //_________________________________________________________________________
756 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
757 {
758   // Get track space point with index i
759   // Coordinates are in the global system
760   AliTOFcluster *cl = fClusters[index];
761   Float_t xyz[3];
762   xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
763   xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
764   xyz[2] = cl->GetZ();
765   Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
766   Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
767   Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
768   Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
769   Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
770   Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
771   Float_t cov[6];
772   cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
773   cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
774   cov[2] = -cosphi*sinth*costh*sigmaz2;
775   cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
776   cov[4] = -sinphi*sinth*costh*sigmaz2;
777   cov[5] = costh*costh*sigmaz2;
778   p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
779
780   // Detector numbering scheme
781   Int_t nSector = AliTOFGeometry::NSectors();
782   Int_t nPlate  = AliTOFGeometry::NPlates();
783   Int_t nStripA = AliTOFGeometry::NStripA();
784   Int_t nStripB = AliTOFGeometry::NStripB();
785   Int_t nStripC = AliTOFGeometry::NStripC();
786
787   Int_t isector = cl->GetDetInd(0);
788   if (isector >= nSector)
789     AliError(Form("Wrong sector number in TOF (%d) !",isector));
790   Int_t iplate = cl->GetDetInd(1);
791   if (iplate >= nPlate)
792     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
793   Int_t istrip = cl->GetDetInd(2);
794
795   Int_t stripOffset = 0;
796   switch (iplate) {
797   case 0:
798     stripOffset = 0;
799     break;
800   case 1:
801     stripOffset = nStripC;
802     break;
803   case 2:
804     stripOffset = nStripC+nStripB;
805     break;
806   case 3:
807     stripOffset = nStripC+nStripB+nStripA;
808     break;
809   case 4:
810     stripOffset = nStripC+nStripB+nStripA+nStripB;
811     break;
812   default:
813     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
814     break;
815   };
816
817   Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
818                stripOffset +
819                istrip;
820   UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
821   p.SetVolumeID((UShort_t)volid);
822   return kTRUE;
823 }
824 //_________________________________________________________________________
825 void AliTOFtracker::InitCheckHists() {
826
827   //Init histos for Digits/Reco QA and Calibration
828
829
830   TDirectory *dir = gDirectory;
831   TFile *logFileTOF = 0;
832
833   TSeqCollection *list = gROOT->GetListOfFiles();
834   int n = list->GetEntries();
835   Bool_t isThere=kFALSE;
836   for(int i=0; i<n; i++) {
837     logFileTOF = (TFile*)list->At(i);
838     if (strstr(logFileTOF->GetName(), "TOFQA.root")){
839       isThere=kTRUE;
840       break;
841     } 
842   }
843
844   if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
845   logFileTOF->cd(); 
846
847   fCalTree = new TTree("CalTree", "Tree for TOF calibration");
848   fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
849   fCalTree->Branch("ToT",&fToT,"TOFToT/F");
850   fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
851   fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
852   fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
853   fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
854
855   //Digits "QA" 
856   fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);  
857   fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);  
858   fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);  
859   fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);  
860
861   //Reco "QA"
862   fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
863   fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
864   fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
865   fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
866   fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
867   fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
868
869   dir->cd();
870
871 }
872
873 //_________________________________________________________________________
874 void AliTOFtracker::SaveCheckHists() {
875
876   //write histos for Digits/Reco QA and Calibration
877
878   TDirectory *dir = gDirectory;
879   TFile *logFileTOF = 0;
880
881   TSeqCollection *list = gROOT->GetListOfFiles();
882   int n = list->GetEntries();
883   Bool_t isThere=kFALSE;
884   for(int i=0; i<n; i++) {
885     logFileTOF = (TFile*)list->At(i);
886     if (strstr(logFileTOF->GetName(), "TOFQA.root")){
887       isThere=kTRUE;
888       break;
889     } 
890   }
891    
892   if(!isThere) {
893           AliError(Form("File TOFQA.root not found!! not wring histograms...."));
894           return;
895   }
896   logFileTOF->cd(); 
897   fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
898   fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
899   fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
900   fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
901   fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
902   fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
903   fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
904   fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
905   fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
906   fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
907   fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
908   logFileTOF->Flush();  
909
910   dir->cd();
911   }
912 //_________________________________________________________________________
913 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) {
914
915   //dummy, for the moment
916   Float_t tofcorr=0.;
917   if(dist<AliTOFGeometry::ZPad()*0.5){
918     tofcorr=tof;
919     //place here the actual correction
920   }else{
921     tofcorr=tof; 
922   } 
923   return tofcorr;
924 }
925 //_________________________________________________________________________
926 Float_t AliTOFtracker::GetTimeZerofromT0(AliESDEvent *event) const {
927
928   //Returns TimeZero as measured by T0 detector
929
930   return event->GetT0();
931 }
932 //_________________________________________________________________________
933 Float_t AliTOFtracker::GetTimeZerofromTOF(AliESDEvent * /*event*/) const {
934
935   //dummy, for the moment. T0 algorithm using tracks on TOF
936   {
937     //place T0 algo here...
938   }
939   return 0.;
940 }
941 //_________________________________________________________________________
942
943 void AliTOFtracker::FillClusterArray(TObjArray* arr) const
944 {
945   //
946   // Returns the TOF cluster array
947   //
948
949   if (fN==0)
950     arr = 0x0;
951   else
952     for (Int_t i=0; i<fN; ++i) arr->Add(fClusters[i]);
953
954 }
955 //_________________________________________________________________________
956
957 void AliTOFtracker::PadRS2TrackingRS(Float_t *ctrackPos, Float_t *differenceT)
958 {
959   //
960   // To convert the 3D distance ctrackPos, referred to the hit pad,
961   // into the 3D distance differenceT, referred to the tracking
962   // reference system.
963   //
964
965   for (Int_t ii=0; ii<3; ii++) differenceT[ii] = 0.;
966
967   AliDebug(1,Form(" track position in ALICE global Ref. frame -> %f, %f, %f",
968                   ctrackPos[0],ctrackPos[1],ctrackPos[2]));
969
970   Int_t detId[5] = {
971     fGeom->GetSector(ctrackPos),
972     fGeom->GetPlate(ctrackPos),
973     fGeom->GetStrip(ctrackPos),
974     fGeom->GetPadZ(ctrackPos),
975     fGeom->GetPadX(ctrackPos)};
976   UShort_t alignableStripIndex = fGeom->GetAliSensVolIndex(detId[0],detId[1],detId[2]);
977   AliDebug(1,Form(" sector = %2d, plate = %1d, strip = %2d ---> stripIndex = %4d",
978                   detId[0], detId[1], detId[2], alignableStripIndex));
979
980   // pad centre coordinates in the strip ref. frame
981   Double_t padCentreL[3] = {(detId[4]-AliTOFGeometry::NpadX()/2)*AliTOFGeometry::XPad()+AliTOFGeometry::XPad()/2.,
982                             0.,
983                             (detId[3]-AliTOFGeometry::NpadZ()/2)*AliTOFGeometry::XPad()+AliTOFGeometry::XPad()/2.};
984   // pad centre coordinates in the strip tracking frame
985   Double_t padCentreT[3] = {0., 0., 0.};
986   TGeoHMatrix l2t = *AliGeomManager::GetTracking2LocalMatrix(alignableStripIndex);
987   l2t.MasterToLocal(padCentreL,padCentreT);
988
989
990   Char_t path[100];
991   // pad centre coordinates in its ref. frame
992   Double_t padCentreL2[3] = {0., 0., 0.};
993   // pad centre coordinates in the ALICE global ref. frame
994   Double_t padCentreG[3] = {0., 0., 0.};
995   fGeom->GetVolumePath(detId,path);
996   gGeoManager->cd(path);
997   TGeoHMatrix g2l = *gGeoManager->GetCurrentMatrix();
998   TGeoHMatrix l2g = g2l.Inverse();
999   l2g.MasterToLocal(padCentreL2,padCentreG);
1000
1001
1002   Char_t path2[100];
1003   // strip centre coordinates in its ref. frame
1004   Double_t stripCentreL[3] = {0., 0., 0.};
1005   // strip centre coordinates in the ALICE global ref. frame
1006   Double_t stripCentreG[3] = {0., 0., 0.};
1007   fGeom->GetVolumePath(detId[0],detId[1],detId[2],path2);
1008   gGeoManager->cd(path2);
1009   TGeoHMatrix g2lb = *gGeoManager->GetCurrentMatrix();
1010   TGeoHMatrix l2gb = g2lb.Inverse();
1011   l2gb.MasterToLocal(stripCentreL,stripCentreG);
1012
1013   TGeoHMatrix g2t = 0;
1014   AliGeomManager::GetTrackingMatrix(alignableStripIndex, g2t);
1015
1016   // track position in the ALICE global ref. frame
1017   Double_t posG[3];
1018   for (Int_t ii=0; ii<3; ii++) posG[ii] = (Double_t)ctrackPos[ii];
1019
1020   // strip centre coordinates in the tracking ref. frame
1021   Double_t stripCentreT[3] = {0., 0., 0.};
1022   // track position in the tracking ref. frame
1023   Double_t posT[3] = {0., 0., 0.};
1024   g2t.MasterToLocal(posG,posT);
1025   g2t.MasterToLocal(stripCentreG,stripCentreT);
1026
1027   for (Int_t ii=0; ii<3; ii++)
1028     AliDebug(1,Form(" track position in ALICE global and tracking Ref. frames -> posG[%d] = %f --- posT[%d] = %f",
1029                     ii, posG[ii], ii, posT[ii]));
1030   for (Int_t ii=0; ii<3; ii++)
1031     AliDebug(1,Form(" pad centre coordinates in its, the ALICE global and tracking Ref. frames -> padCentreL[%d] = %f --- padCentreG[%d] = %f --- padCentreT[%d] = %f",
1032                     ii, padCentreL[ii],
1033                     ii, padCentreG[ii],
1034                     ii, padCentreT[ii]));
1035   for (Int_t ii=0; ii<3; ii++)
1036     AliDebug(1,Form(" strip centre coordinates in its, the ALICE global and tracking Ref. frames -> stripCentreL[%d] = %f --- stripCentreG[%d] = %f --- stripCentreT[%d] = %f",
1037                     ii, stripCentreL[ii],
1038                     ii, stripCentreG[ii],
1039                     ii, stripCentreT[ii]));
1040   for (Int_t ii=0; ii<3; ii++)
1041     AliDebug(1,Form(" difference between the track position and the pad centre in the tracking Ref. frame -> posT[%d]-padCentreT[%d] = %f",
1042                     ii,ii,
1043                     posT[ii]-padCentreT[ii]));
1044
1045   for (Int_t ii=0; ii<3; ii++) differenceT[ii] = (Float_t)(posT[ii]-padCentreT[ii]);
1046
1047 }