1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //--------------------------------------------------------------------//
18 // AliTOFtracker Class //
19 // Task: Perform association of the ESD tracks to TOF Clusters //
20 // and Update ESD track with associated TOF Cluster parameters //
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 //
27 //--------------------------------------------------------------------//
32 #include <TClonesArray.h>
33 #include <TGeoManager.h>
38 #include "AliAlignObj.h"
39 #include "AliESDtrack.h"
42 #include "AliTrackPointArray.h"
44 #include "AliTOFcalib.h"
45 #include "AliTOFRecoParam.h"
46 #include "AliTOFcluster.h"
47 #include "AliTOFGeometryV5.h"
48 #include "AliTOFtracker.h"
49 #include "AliTOFtrack.h"
51 extern TGeoManager *gGeoManager;
54 ClassImp(AliTOFtracker)
56 //_____________________________________________________________________________
57 AliTOFtracker::AliTOFtracker():
88 //AliTOFtracker main Ctor
90 // Gettimg the geometry
91 fGeom=new AliTOFGeometryV5();
92 // Read the reconstruction parameters from the OCDB
93 AliTOFcalib *calib = new AliTOFcalib(fGeom);
94 fRecoParam = (AliTOFRecoParam*)calib->ReadRecParFromCDB("TOF/Calib",-1);
95 if(fRecoParam->GetApplyPbPbCuts())fRecoParam=fRecoParam->GetPbPbparam();
97 parPID[0]=fRecoParam->GetTimeResolution();
98 parPID[1]=fRecoParam->GetTimeNSigma();
99 fPid=new AliTOFpidESD(parPID);
103 //_____________________________________________________________________________
104 AliTOFtracker::AliTOFtracker(const AliTOFtracker &t):
126 fHRecSigYVsPWin(0x0),
127 fHRecSigZVsPWin(0x0),
136 //AliTOFtracker copy Ctor
140 fNseedsTOF=t.fNseedsTOF;
141 fngoodmatch=t.fngoodmatch;
142 fnbadmatch=t.fnbadmatch;
143 fnunmatch=t.fnunmatch;
145 fRecoParam=t.fRecoParam;
153 //_____________________________________________________________________________
154 AliTOFtracker& AliTOFtracker::operator=(const AliTOFtracker &t)
156 //AliTOFtracker assignment operator
158 this->fNseeds=t.fNseeds;
159 this->fNseedsTOF=t.fNseedsTOF;
160 this->fngoodmatch=t.fngoodmatch;
161 this->fnbadmatch=t.fnbadmatch;
162 this->fnunmatch=t.fnunmatch;
163 this->fnmatch=t.fnmatch;
164 this->fRecoParam = t.fRecoParam;
165 this->fGeom = t.fGeom;
167 this->fSeeds=t.fSeeds;
168 this->fTracks=t.fTracks;
173 //_____________________________________________________________________________
174 AliTOFtracker::~AliTOFtracker() {
186 delete fHDigClusTime;
192 delete fHRecSigYVsPWin;
193 delete fHRecSigZVsPWin;
196 //_____________________________________________________________________________
197 Int_t AliTOFtracker::PropagateBack(AliESD* event) {
199 // Gets seeds from ESD event and Match with TOF Clusters
203 //Initialise some counters
212 Int_t ntrk=event->GetNumberOfTracks();
214 fSeeds= new TClonesArray("AliESDtrack",ntrk);
215 TClonesArray &aESDTrack = *fSeeds;
218 //Load ESD tracks into a local Array of ESD Seeds
220 for (Int_t i=0; i<fNseeds; i++) {
221 AliESDtrack *t=event->GetTrack(i);
222 new(aESDTrack[i]) AliESDtrack(*t);
225 //Prepare ESD tracks candidates for TOF Matching
228 //First Step with Strict Matching Criterion
231 //Second Step with Looser Matching Criterion
234 AliInfo(Form("Number of matched tracks: %d",fnmatch));
235 AliInfo(Form("Number of good matched tracks: %d",fngoodmatch));
236 AliInfo(Form("Number of bad matched tracks: %d",fnbadmatch));
238 //Update the matched ESD tracks
240 for (Int_t i=0; i<ntrk; i++) {
241 AliESDtrack *t=event->GetTrack(i);
242 AliESDtrack *seed =(AliESDtrack*)fSeeds->UncheckedAt(i);
243 if(seed->GetTOFsignal()>0){
245 t->SetTOFsignal(seed->GetTOFsignal());
246 t->SetTOFcluster(seed->GetTOFcluster());
247 t->SetTOFsignalToT(seed->GetTOFsignalToT());
248 t->SetTOFsignalRaw(seed->GetTOFsignalRaw());
249 t->SetTOFsignalDz(seed->GetTOFsignalDz());
250 t->SetTOFCalChannel(seed->GetTOFCalChannel());
251 Int_t tlab[3]; seed->GetTOFLabel(tlab);
252 t->SetTOFLabel(tlab);
253 AliTOFtrack *track = new AliTOFtrack(*seed);
254 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
259 //Handle Time Zero information
261 Double_t timeZero=0.;
262 Double_t timeZeroMax=99999.;
263 Bool_t usetimeZero = fRecoParam->UseTimeZero();
264 Bool_t timeZeroFromT0 = fRecoParam->GetTimeZerofromT0();
265 Bool_t timeZeroFromTOF = fRecoParam->GetTimeZerofromTOF();
267 AliDebug(1,Form("Use Time Zero?: %d",usetimeZero));
268 AliDebug(1,Form("Time Zero from T0? : %d",timeZeroFromT0));
269 AliDebug(1,Form("Time Zero From TOF? : %d",timeZeroFromTOF));
273 timeZero=GetTimeZerofromT0(event);
275 if(timeZeroFromTOF && (timeZero>timeZeroMax || !timeZeroFromT0)){
276 timeZero=GetTimeZerofromTOF(event);
279 AliDebug(2,Form("time Zero used in PID: %f",timeZero));
281 fPid->MakePID(event,timeZero);
296 //_________________________________________________________________________
297 void AliTOFtracker::CollectESD() {
298 //prepare the set of ESD tracks to be matched to clusters in TOF
303 fTracks= new TClonesArray("AliTOFtrack");
304 TClonesArray &aTOFTrack = *fTracks;
305 for (Int_t i=0; i<fNseeds; i++) {
307 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(i);
308 if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
310 // TRD 'good' tracks, already propagated at 371 cm
312 AliTOFtrack *track = new AliTOFtrack(*t); // New
313 Double_t x = track->GetX(); //New
315 if (((t->GetStatus()&AliESDtrack::kTRDout)!=0 ) &&
316 ( x >= fGeom->RinTOF()) ){
317 track->SetSeedIndex(i);
318 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
319 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
325 // Propagate the rest of TPCbp
328 if(track->PropagateToInnerTOF()){
329 track->SetSeedIndex(i);
330 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
331 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
339 AliInfo(Form("Number of TOF seeds %i",fNseedsTOF));
340 AliInfo(Form("Number of TOF seeds Type 1 %i",seedsTOF1));
341 AliInfo(Form("Number of TOF seeds Type 2 %i",seedsTOF2));
343 // Sort according uncertainties on track position
347 //_________________________________________________________________________
348 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
351 // Parameters used/regulating the reconstruction
353 static Float_t corrLen=0.75;
354 static Float_t detDepth=15.3;
355 static Float_t padDepth=0.5;
357 Float_t dY=fGeom->XPad();
358 Float_t dZ=fGeom->ZPad();
360 Float_t sensRadius = fRecoParam->GetSensRadius();
361 Float_t stepSize = fRecoParam->GetStepSize();
362 Float_t scaleFact = fRecoParam->GetWindowScaleFact();
363 Float_t dyMax=fRecoParam->GetWindowSizeMaxY();
364 Float_t dzMax=fRecoParam->GetWindowSizeMaxZ();
365 Float_t dCut=fRecoParam->GetDistanceCut();
366 Double_t maxChi2=fRecoParam->GetMaxChi2();
367 Bool_t timeWalkCorr = fRecoParam->GetTimeWalkCorr();
369 AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++ \n");
370 AliDebug(1,Form("TOF sens radius: %f",sensRadius));
371 AliDebug(1,Form("TOF step size: %f",stepSize));
372 AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
373 AliDebug(1,Form("TOF Window max dy: %f",dyMax));
374 AliDebug(1,Form("TOF Window max dz: %f",dzMax));
375 AliDebug(1,Form("TOF distance Cut: %f",dCut));
376 AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
377 AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));
379 //Match ESD tracks to clusters in TOF
382 // Get the number of propagation steps
384 Int_t nSteps=(Int_t)(detDepth/stepSize);
386 AliTOFcalib *calib = new AliTOFcalib(fGeom);
388 //PH Arrays (moved outside of the loop)
390 Float_t * trackPos[4];
391 for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
393 for (Int_t ii=0;ii<6;ii++) clind[ii] = new Int_t[fN];
395 for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
397 AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
398 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(track->GetSeedIndex());
399 if(t->GetTOFsignal()>0. ) continue;
400 AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
406 Float_t distZ[10000];
407 Float_t cxpos[10000];
408 Float_t crecL[10000];
409 TGeoHMatrix global[1000];
411 // Determine a window around the track
414 trackTOFin->GetExternalParameters(x,par);
416 trackTOFin->GetExternalCovariance(cov);
420 ((5*TMath::Sqrt(cov[0]) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius);
423 (5*TMath::Sqrt(cov[2]) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
425 Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
426 if (phi<-TMath::Pi())phi+=2*TMath::Pi();
427 if (phi>=TMath::Pi())phi-=2*TMath::Pi();
430 //upper limit on window's size.
432 if(dz> dzMax) dz=dzMax;
433 if(dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
438 // find the clusters in the window of the track
440 for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
442 AliTOFcluster *c=fClusters[k];
443 if (c->GetZ() > z+dz) break;
444 if (c->IsUsed()) continue;
446 if (!c->GetStatus()) continue; // skip bad channels as declared in OCDB
448 Double_t dph=TMath::Abs(c->GetPhi()-phi);
449 if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
450 if (TMath::Abs(dph)>dphi) continue;
452 Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
453 Double_t p[2]={yc, c->GetZ()};
454 Double_t cov[3]= {dY*dY/12., 0., dZ*dZ/12.};
455 if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov) > maxChi2)continue;
457 clind[0][nc] = c->GetDetInd(0);
458 clind[1][nc] = c->GetDetInd(1);
459 clind[2][nc] = c->GetDetInd(2);
460 clind[3][nc] = c->GetDetInd(3);
461 clind[4][nc] = c->GetDetInd(4);
470 fGeom->GetVolumePath(ind,path);
471 gGeoManager->cd(path);
472 global[nc] = *gGeoManager->GetCurrentMatrix();
476 //start fine propagation
478 Int_t nStepsDone = 0;
479 for( Int_t istep=0; istep<nSteps; istep++){
481 Float_t xs=fGeom->RinTOF()+istep*0.1;
482 Double_t ymax=xs*TMath::Tan(0.5*fGeom->GetAlpha());
485 Double_t ysect=trackTOFin->GetYat(xs,skip);
488 if (!trackTOFin->Rotate(fGeom->GetAlpha())) {
491 } else if (ysect <-ymax) {
492 if (!trackTOFin->Rotate(fGeom->GetAlpha())) {
497 if(!trackTOFin->PropagateTo(xs)) {
503 // store the running point (Globalrf) - fine propagation
506 trackTOFin->GetXYZ(r);
507 trackPos[0][istep]= (Float_t) r[0];
508 trackPos[1][istep]= (Float_t) r[1];
509 trackPos[2][istep]= (Float_t) r[2];
510 trackPos[3][istep]= trackTOFin->GetIntegratedLength();
515 Bool_t accept = kFALSE;
516 Bool_t isInside =kFALSE;
517 for (Int_t istep=0; istep<nStepsDone; istep++) {
519 Float_t ctrackPos[3];
521 ctrackPos[0]= trackPos[0][istep];
522 ctrackPos[1]= trackPos[1][istep];
523 ctrackPos[2]= trackPos[2][istep];
525 //now see whether the track matches any of the TOF clusters
529 for (Int_t i=0; i<nc; i++){
531 cind[0]= clind[0][i];
532 cind[1]= clind[1][i];
533 cind[2]= clind[2][i];
534 cind[3]= clind[3][i];
535 cind[4]= clind[4][i];
536 isInside=fGeom->IsInsideThePad(global[i],ctrackPos,dist3d);
539 Float_t xLoc=dist3d[0];
540 Float_t rLoc=TMath::Sqrt(dist3d[1]*dist3d[1]+dist3d[2]*dist3d[2]);
541 accept = (TMath::Abs(xLoc)<padDepth*0.5 && rLoc<dCut);
547 dist[nfound]=TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[1]*dist3d[1]+dist3d[2]*dist3d[2]);
548 distZ[nfound]=dist3d[2];
549 crecL[nfound]=trackPos[3][istep];
550 index[nfound]=clind[5][i]; // store cluster id
551 cxpos[nfound]=fGeom->RinTOF()+istep*0.1; //store prop.radius
553 if(accept &&!mLastStep)break;
555 } //end for on the clusters
556 if(accept &&!mLastStep)break;
557 } //end for on the steps
569 // now choose the cluster to be matched with the track.
574 Float_t mindist=1000.;
576 for (Int_t iclus= 0; iclus<nfound;iclus++){
577 if (dist[iclus]< mindist){
578 mindist = dist[iclus];
579 mindistZ = distZ[iclus];
581 idclus =index[iclus];
582 recL=crecL[iclus]+corrLen*0.5;
587 AliTOFcluster *c=fClusters[idclus];
590 // Track length correction for matching Step 2
593 Float_t rc=TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
594 Float_t rt=TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
595 +trackPos[1][70]*trackPos[1][70]
596 +trackPos[2][70]*trackPos[2][70]);
598 recL=trackPos[3][70]+dlt;
602 (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
604 (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
606 (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
610 AliDebug(2,Form(" track label good %5i",trackTOFin->GetLabel()));
616 AliDebug(2,Form(" track label bad %5i",trackTOFin->GetLabel()));
622 // Store quantities to be used in the TOF Calibration
623 Float_t tToT=fGeom->ToTBinWidth()*c->GetToT()*1E-3; // in ns
624 t->SetTOFsignalToT(tToT);
625 Float_t rawTime=fGeom->TdcBinWidth()*c->GetTDCRAW()+32; // RAW time,in ps
626 t->SetTOFsignalRaw(rawTime);
627 t->SetTOFsignalDz(mindistZ);
628 AliDebug(2,Form(" Setting TOF raw time: %f z distance: %f time: %f = ",rawTime,mindistZ));
630 ind[0]=c->GetDetInd(0);
631 ind[1]=c->GetDetInd(1);
632 ind[2]=c->GetDetInd(2);
633 ind[3]=c->GetDetInd(3);
634 ind[4]=c->GetDetInd(4);
635 Int_t calindex = calib->GetIndex(ind);
636 t->SetTOFCalChannel(calindex);
638 // keep track of the track labels in the matched cluster
640 tlab[0]=c->GetLabel(0);
641 tlab[1]=c->GetLabel(1);
642 tlab[2]=c->GetLabel(2);
643 AliDebug(2,Form(" tdc time of the matched track %i = ",c->GetTDC()));
644 Double_t tof=fGeom->TdcBinWidth()*c->GetTDC()+32; // in ps
645 AliDebug(2,Form(" tof time of the matched track: %f = ",tof));
646 Double_t tofcorr=tof;
647 if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
648 AliDebug(2,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));
649 //Set TOF time signal and pointer to the matched cluster
650 t->SetTOFsignal(tofcorr);
651 t->SetTOFcluster(idclus); // pointing to the recPoints tree
654 Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time);
655 Double_t mom=t->GetP();
656 for(Int_t j=0;j<AliPID::kSPECIES;j++){
657 Double_t mass=AliPID::ParticleMass(j);
658 time[j]+=(recL-trackPos[3][0])/3e-2*TMath::Sqrt(mom*mom+mass*mass)/mom;
661 AliTOFtrack *trackTOFout = new AliTOFtrack(*t);
662 trackTOFout->PropagateTo(xpos);
663 t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
664 t->SetIntegratedLength(recL);
665 t->SetIntegratedTimes(time);
666 t->SetTOFLabel(tlab);
669 // Fill Reco-QA histos for Reconstruction
670 fHRecNClus->Fill(nc);
671 fHRecDist->Fill(mindist);
672 fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
673 fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
674 fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
675 fHRecSigZVsPWin->Fill(mom,dz);
677 // Fill Tree for on-the-fly offline Calibration
679 if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 )){
690 for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
691 for (Int_t ii=0;ii<6;ii++) delete [] clind[ii];
694 //_________________________________________________________________________
695 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
696 //--------------------------------------------------------------------
697 //This function loads the TOF clusters
698 //--------------------------------------------------------------------
700 Int_t npadX = fGeom->NpadX();
701 Int_t npadZ = fGeom->NpadZ();
702 Int_t nStripA = fGeom->NStripA();
703 Int_t nStripB = fGeom->NStripB();
704 Int_t nStripC = fGeom->NStripC();
706 TBranch *branch=cTree->GetBranch("TOF");
708 AliError("can't get the branch with the TOF clusters !");
712 TClonesArray dummy("AliTOFcluster",10000), *clusters=&dummy;
713 branch->SetAddress(&clusters);
716 Int_t nc=clusters->GetEntriesFast();
717 fHDigNClus->Fill(nc);
719 AliInfo(Form("Number of clusters: %d",nc));
721 for (Int_t i=0; i<nc; i++) {
722 AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
723 fClusters[i]=new AliTOFcluster(*c); fN++;
725 // Fill Digits QA histos
727 Int_t isector = c->GetDetInd(0);
728 Int_t iplate = c->GetDetInd(1);
729 Int_t istrip = c->GetDetInd(2);
730 Int_t ipadX = c->GetDetInd(4);
731 Int_t ipadZ = c->GetDetInd(3);
733 Float_t time =(fGeom->TdcBinWidth()*c->GetTDC())*1E-3; // in ns
734 Float_t tot = (fGeom->TdcBinWidth()*c->GetToT())*1E-3;//in ns
736 Int_t stripOffset = 0;
742 stripOffset = nStripC;
745 stripOffset = nStripC+nStripB;
748 stripOffset = nStripC+nStripB+nStripA;
751 stripOffset = nStripC+nStripB+nStripA+nStripB;
754 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
757 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
758 Int_t phiindex=npadX*isector+ipadX+1;
759 fHDigClusMap->Fill(zindex,phiindex);
760 fHDigClusTime->Fill(time);
761 fHDigClusToT->Fill(tot);
768 //_________________________________________________________________________
769 void AliTOFtracker::UnloadClusters() {
770 //--------------------------------------------------------------------
771 //This function unloads TOF clusters
772 //--------------------------------------------------------------------
773 for (Int_t i=0; i<fN; i++) {
780 //_________________________________________________________________________
781 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
782 //--------------------------------------------------------------------
783 // This function returns the index of the nearest cluster
784 //--------------------------------------------------------------------
786 if (z <= fClusters[0]->GetZ()) return 0;
787 if (z > fClusters[fN-1]->GetZ()) return fN;
788 Int_t b=0, e=fN-1, m=(b+e)/2;
789 for (; b<e; m=(b+e)/2) {
790 if (z > fClusters[m]->GetZ()) b=m+1;
796 //_________________________________________________________________________
797 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
799 // Get track space point with index i
800 // Coordinates are in the global system
801 AliTOFcluster *cl = fClusters[index];
803 xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
804 xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
806 Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
807 Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
808 Float_t tiltangle = fGeom->GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
809 Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
810 Float_t sigmay2 = fGeom->XPad()*fGeom->XPad()/12.;
811 Float_t sigmaz2 = fGeom->ZPad()*fGeom->ZPad()/12.;
813 cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
814 cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
815 cov[2] = -cosphi*sinth*costh*sigmaz2;
816 cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
817 cov[4] = -sinphi*sinth*costh*sigmaz2;
818 cov[5] = costh*costh*sigmaz2;
819 p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
821 // Detector numbering scheme
822 Int_t nSector = fGeom->NSectors();
823 Int_t nPlate = fGeom->NPlates();
824 Int_t nStripA = fGeom->NStripA();
825 Int_t nStripB = fGeom->NStripB();
826 Int_t nStripC = fGeom->NStripC();
828 Int_t isector = cl->GetDetInd(0);
829 if (isector >= nSector)
830 AliError(Form("Wrong sector number in TOF (%d) !",isector));
831 Int_t iplate = cl->GetDetInd(1);
832 if (iplate >= nPlate)
833 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
834 Int_t istrip = cl->GetDetInd(2);
836 Int_t stripOffset = 0;
842 stripOffset = nStripC;
845 stripOffset = nStripC+nStripB;
848 stripOffset = nStripC+nStripB+nStripA;
851 stripOffset = nStripC+nStripB+nStripA+nStripB;
854 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
858 Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
861 UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
862 p.SetVolumeID((UShort_t)volid);
865 //_________________________________________________________________________
866 void AliTOFtracker::InitCheckHists() {
868 //Init histos for Digits/Reco QA and Calibration
871 fCalTree = new TTree("CalTree", "Tree for TOF calibration");
872 fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
873 fCalTree->Branch("ToT",&fToT,"TOFToT/F");
874 fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
875 fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
876 fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
877 fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
880 fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);
881 fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);
882 fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);
883 fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);
886 fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
887 fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
888 fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
889 fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
890 fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
891 fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
894 //_________________________________________________________________________
895 void AliTOFtracker::SaveCheckHists() {
897 //write histos for Digits/Reco QA and Calibration
899 TDirectory *dir = gDirectory;
901 TFile *logFileTOF = 0;
903 TSeqCollection *list = gROOT->GetListOfFiles();
904 int n = list->GetEntries();
905 for(int i=0; i<n; i++) {
906 logFile = (TFile*)list->At(i);
907 if (strstr(logFile->GetName(), "AliESDs.root")) break;
910 Bool_t isThere=kFALSE;
911 for(int i=0; i<n; i++) {
912 logFileTOF = (TFile*)list->At(i);
913 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
920 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
921 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
922 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
923 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
924 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
925 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
926 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
927 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
928 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
929 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
930 //fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
933 if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
935 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
936 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
937 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
938 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
939 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
940 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
941 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
942 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
943 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
944 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
945 fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
950 //_________________________________________________________________________
951 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) {
953 //dummy, for the moment
955 if(dist<fGeom->ZPad()*0.5){
957 //place here the actual correction
963 //_________________________________________________________________________
964 Float_t AliTOFtracker::GetTimeZerofromT0(AliESD *event) const {
966 //Returns TimeZero as measured by T0 detector
968 return event->GetT0();
970 //_________________________________________________________________________
971 Float_t AliTOFtracker::GetTimeZerofromTOF(AliESD * /*event*/) const {
973 //dummy, for the moment. T0 algorithm using tracks on TOF
975 //place T0 algo here...