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"
40 #include "AliESDEvent.h"
42 #include "AliTrackPointArray.h"
44 #include "AliTOFcalib.h"
45 #include "AliTOFRecoParam.h"
46 #include "AliTOFcluster.h"
47 #include "AliTOFGeometry.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 AliTOFGeometry();
92 // Read the reconstruction parameters from the OCDB
93 AliTOFcalib* calib=new AliTOFcalib();
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;
152 //_____________________________________________________________________________
153 AliTOFtracker& AliTOFtracker::operator=(const AliTOFtracker &t)
155 //AliTOFtracker assignment operator
157 this->fNseeds=t.fNseeds;
158 this->fNseedsTOF=t.fNseedsTOF;
159 this->fngoodmatch=t.fngoodmatch;
160 this->fnbadmatch=t.fnbadmatch;
161 this->fnunmatch=t.fnunmatch;
162 this->fnmatch=t.fnmatch;
163 this->fRecoParam = t.fRecoParam;
166 this->fSeeds=t.fSeeds;
167 this->fTracks=t.fTracks;
172 //_____________________________________________________________________________
173 AliTOFtracker::~AliTOFtracker() {
185 delete fHDigClusTime;
191 delete fHRecSigYVsPWin;
192 delete fHRecSigZVsPWin;
195 //_____________________________________________________________________________
196 Int_t AliTOFtracker::PropagateBack(AliESDEvent* event) {
198 // Gets seeds from ESD event and Match with TOF Clusters
202 //Initialise some counters
211 Int_t ntrk=event->GetNumberOfTracks();
213 fSeeds= new TClonesArray("AliESDtrack",ntrk);
214 TClonesArray &aESDTrack = *fSeeds;
217 //Load ESD tracks into a local Array of ESD Seeds
219 for (Int_t i=0; i<fNseeds; i++) {
220 AliESDtrack *t=event->GetTrack(i);
221 new(aESDTrack[i]) AliESDtrack(*t);
224 //Prepare ESD tracks candidates for TOF Matching
227 //First Step with Strict Matching Criterion
230 //Second Step with Looser Matching Criterion
233 AliInfo(Form("Number of matched tracks: %d",fnmatch));
234 AliInfo(Form("Number of good matched tracks: %d",fngoodmatch));
235 AliInfo(Form("Number of bad matched tracks: %d",fnbadmatch));
237 //Update the matched ESD tracks
239 for (Int_t i=0; i<ntrk; i++) {
240 AliESDtrack *t=event->GetTrack(i);
241 AliESDtrack *seed =(AliESDtrack*)fSeeds->UncheckedAt(i);
242 if(seed->GetTOFsignal()>0){
243 t->SetTOFsignal(seed->GetTOFsignal());
244 t->SetTOFcluster(seed->GetTOFcluster());
245 t->SetTOFsignalToT(seed->GetTOFsignalToT());
246 t->SetTOFsignalRaw(seed->GetTOFsignalRaw());
247 t->SetTOFsignalDz(seed->GetTOFsignalDz());
248 t->SetTOFCalChannel(seed->GetTOFCalChannel());
249 Int_t tlab[3]; seed->GetTOFLabel(tlab);
250 t->SetTOFLabel(tlab);
251 AliTOFtrack *track = new AliTOFtrack(*seed);
252 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
257 //Handle Time Zero information
259 Double_t timeZero=0.;
260 Double_t timeZeroMax=99999.;
261 Bool_t usetimeZero = fRecoParam->UseTimeZero();
262 Bool_t timeZeroFromT0 = fRecoParam->GetTimeZerofromT0();
263 Bool_t timeZeroFromTOF = fRecoParam->GetTimeZerofromTOF();
265 AliDebug(1,Form("Use Time Zero?: %d",usetimeZero));
266 AliDebug(1,Form("Time Zero from T0? : %d",timeZeroFromT0));
267 AliDebug(1,Form("Time Zero From TOF? : %d",timeZeroFromTOF));
271 timeZero=GetTimeZerofromT0(event);
273 if(timeZeroFromTOF && (timeZero>timeZeroMax || !timeZeroFromT0)){
274 timeZero=GetTimeZerofromTOF(event);
277 AliDebug(2,Form("time Zero used in PID: %f",timeZero));
279 fPid->MakePID(event,timeZero);
294 //_________________________________________________________________________
295 void AliTOFtracker::CollectESD() {
296 //prepare the set of ESD tracks to be matched to clusters in TOF
301 fTracks= new TClonesArray("AliTOFtrack");
302 TClonesArray &aTOFTrack = *fTracks;
303 for (Int_t i=0; i<fNseeds; i++) {
305 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(i);
306 if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
308 // TRD 'good' tracks, already propagated at 371 cm
310 AliTOFtrack *track = new AliTOFtrack(*t); // New
311 Double_t x = track->GetX(); //New
312 if (((t->GetStatus()&AliESDtrack::kTRDout)!=0 ) &&
313 ( x >= AliTOFGeometry::RinTOF()) ){
314 track->SetSeedIndex(i);
315 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
316 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
322 // Propagate the rest of TPCbp
325 if(track->PropagateToInnerTOF()){
326 track->SetSeedIndex(i);
327 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
328 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
336 AliInfo(Form("Number of TOF seeds %i",fNseedsTOF));
337 AliInfo(Form("Number of TOF seeds Type 1 %i",seedsTOF1));
338 AliInfo(Form("Number of TOF seeds Type 2 %i",seedsTOF2));
340 // Sort according uncertainties on track position
344 //_________________________________________________________________________
345 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
348 // Parameters used/regulating the reconstruction
350 static Float_t corrLen=0.75;
351 static Float_t detDepth=15.3;
352 static Float_t padDepth=0.5;
354 Float_t dY=AliTOFGeometry::XPad();
355 Float_t dZ=AliTOFGeometry::ZPad();
357 Float_t sensRadius = fRecoParam->GetSensRadius();
358 Float_t stepSize = fRecoParam->GetStepSize();
359 Float_t scaleFact = fRecoParam->GetWindowScaleFact();
360 Float_t dyMax=fRecoParam->GetWindowSizeMaxY();
361 Float_t dzMax=fRecoParam->GetWindowSizeMaxZ();
362 Float_t dCut=fRecoParam->GetDistanceCut();
363 Double_t maxChi2=fRecoParam->GetMaxChi2TRD();
364 Bool_t timeWalkCorr = fRecoParam->GetTimeWalkCorr();
366 AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++ \n");
367 AliDebug(1,Form("TOF sens radius: %f",sensRadius));
368 AliDebug(1,Form("TOF step size: %f",stepSize));
369 AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
370 AliDebug(1,Form("TOF Window max dy: %f",dyMax));
371 AliDebug(1,Form("TOF Window max dz: %f",dzMax));
372 AliDebug(1,Form("TOF distance Cut: %f",dCut));
373 AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
374 AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));
376 //Match ESD tracks to clusters in TOF
379 // Get the number of propagation steps
381 Int_t nSteps=(Int_t)(detDepth/stepSize);
383 //PH Arrays (moved outside of the loop)
385 Float_t * trackPos[4];
386 for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
388 for (Int_t ii=0;ii<6;ii++) clind[ii] = new Int_t[fN];
390 for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
392 AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
393 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(track->GetSeedIndex());
394 if(t->GetTOFsignal()>0. ) continue;
395 AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
401 Float_t distZ[10000];
402 Float_t cxpos[10000];
403 Float_t crecL[10000];
404 TGeoHMatrix global[1000];
406 // Determine a window around the track
409 trackTOFin->GetExternalParameters(x,par);
411 trackTOFin->GetExternalCovariance(cov);
415 ((5*TMath::Sqrt(cov[0]) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius);
418 (5*TMath::Sqrt(cov[2]) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
420 Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
421 if (phi<-TMath::Pi())phi+=2*TMath::Pi();
422 if (phi>=TMath::Pi())phi-=2*TMath::Pi();
425 //upper limit on window's size.
427 if(dz> dzMax) dz=dzMax;
428 if(dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
433 // find the clusters in the window of the track
435 for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
437 AliTOFcluster *c=fClusters[k];
438 if (c->GetZ() > z+dz) break;
439 if (c->IsUsed()) continue;
441 if (!c->GetStatus()) continue; // skip bad channels as declared in OCDB
443 Double_t dph=TMath::Abs(c->GetPhi()-phi);
444 if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
445 if (TMath::Abs(dph)>dphi) continue;
447 Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
448 Double_t p[2]={yc, c->GetZ()};
449 Double_t cov[3]= {dY*dY/12., 0., dZ*dZ/12.};
450 if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov) > maxChi2)continue;
452 clind[0][nc] = c->GetDetInd(0);
453 clind[1][nc] = c->GetDetInd(1);
454 clind[2][nc] = c->GetDetInd(2);
455 clind[3][nc] = c->GetDetInd(3);
456 clind[4][nc] = c->GetDetInd(4);
465 fGeom->GetVolumePath(ind,path);
466 gGeoManager->cd(path);
467 global[nc] = *gGeoManager->GetCurrentMatrix();
471 //start fine propagation
473 Int_t nStepsDone = 0;
474 for( Int_t istep=0; istep<nSteps; istep++){
476 Float_t xs=AliTOFGeometry::RinTOF()+istep*0.1;
477 Double_t ymax=xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
480 Double_t ysect=trackTOFin->GetYat(xs,skip);
483 if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
486 } else if (ysect <-ymax) {
487 if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
492 if(!trackTOFin->PropagateTo(xs)) {
498 // store the running point (Globalrf) - fine propagation
501 trackTOFin->GetXYZ(r);
502 trackPos[0][istep]= (Float_t) r[0];
503 trackPos[1][istep]= (Float_t) r[1];
504 trackPos[2][istep]= (Float_t) r[2];
505 trackPos[3][istep]= trackTOFin->GetIntegratedLength();
510 Bool_t accept = kFALSE;
511 Bool_t isInside =kFALSE;
512 for (Int_t istep=0; istep<nStepsDone; istep++) {
514 Float_t ctrackPos[3];
516 ctrackPos[0]= trackPos[0][istep];
517 ctrackPos[1]= trackPos[1][istep];
518 ctrackPos[2]= trackPos[2][istep];
520 //now see whether the track matches any of the TOF clusters
524 for (Int_t i=0; i<nc; i++){
526 cind[0]= clind[0][i];
527 cind[1]= clind[1][i];
528 cind[2]= clind[2][i];
529 cind[3]= clind[3][i];
530 cind[4]= clind[4][i];
531 isInside=fGeom->IsInsideThePad(global[i],ctrackPos,dist3d);
534 Float_t xLoc=dist3d[0];
535 Float_t rLoc=TMath::Sqrt(dist3d[1]*dist3d[1]+dist3d[2]*dist3d[2]);
536 accept = (TMath::Abs(xLoc)<padDepth*0.5 && rLoc<dCut);
542 dist[nfound]=TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[1]*dist3d[1]+dist3d[2]*dist3d[2]);
543 distZ[nfound]=dist3d[2];
544 crecL[nfound]=trackPos[3][istep];
545 index[nfound]=clind[5][i]; // store cluster id
546 cxpos[nfound]=AliTOFGeometry::RinTOF()+istep*0.1; //store prop.radius
548 if(accept &&!mLastStep)break;
550 } //end for on the clusters
551 if(accept &&!mLastStep)break;
552 } //end for on the steps
564 // now choose the cluster to be matched with the track.
569 Float_t mindist=1000.;
571 for (Int_t iclus= 0; iclus<nfound;iclus++){
572 if (dist[iclus]< mindist){
573 mindist = dist[iclus];
574 mindistZ = distZ[iclus];
576 idclus =index[iclus];
577 recL=crecL[iclus]+corrLen*0.5;
582 AliTOFcluster *c=fClusters[idclus];
585 // Track length correction for matching Step 2
588 Float_t rc=TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
589 Float_t rt=TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
590 +trackPos[1][70]*trackPos[1][70]
591 +trackPos[2][70]*trackPos[2][70]);
593 recL=trackPos[3][70]+dlt;
597 (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
599 (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
601 (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
605 AliDebug(2,Form(" track label good %5i",trackTOFin->GetLabel()));
611 AliDebug(2,Form(" track label bad %5i",trackTOFin->GetLabel()));
617 // Store quantities to be used in the TOF Calibration
618 Float_t tToT=AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3; // in ns
619 t->SetTOFsignalToT(tToT);
620 Float_t rawTime=AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()+32; // RAW time,in ps
621 t->SetTOFsignalRaw(rawTime);
622 t->SetTOFsignalDz(mindistZ);
623 AliDebug(2,Form(" Setting TOF raw time: %f, z distance: %f time: %f = ",rawTime,mindistZ));
625 ind[0]=c->GetDetInd(0);
626 ind[1]=c->GetDetInd(1);
627 ind[2]=c->GetDetInd(2);
628 ind[3]=c->GetDetInd(3);
629 ind[4]=c->GetDetInd(4);
630 Int_t calindex = AliTOFGeometry::GetIndex(ind);
631 t->SetTOFCalChannel(calindex);
633 // keep track of the track labels in the matched cluster
635 tlab[0]=c->GetLabel(0);
636 tlab[1]=c->GetLabel(1);
637 tlab[2]=c->GetLabel(2);
638 AliDebug(2,Form(" tdc time of the matched track %i = ",c->GetTDC()));
639 Double_t tof=AliTOFGeometry::TdcBinWidth()*c->GetTDC()+32; // in ps
640 AliDebug(2,Form(" tof time of the matched track: %f = ",tof));
641 Double_t tofcorr=tof;
642 if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
643 AliDebug(2,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));
644 //Set TOF time signal and pointer to the matched cluster
645 t->SetTOFsignal(tofcorr);
646 t->SetTOFcluster(idclus); // pointing to the recPoints tree
649 Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time);
650 Double_t mom=t->GetP();
651 for(Int_t j=0;j<AliPID::kSPECIES;j++){
652 Double_t mass=AliPID::ParticleMass(j);
653 time[j]+=(recL-trackPos[3][0])/3e-2*TMath::Sqrt(mom*mom+mass*mass)/mom;
656 AliTOFtrack *trackTOFout = new AliTOFtrack(*t);
657 trackTOFout->PropagateTo(xpos);
658 t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
659 t->SetIntegratedLength(recL);
660 t->SetIntegratedTimes(time);
661 t->SetTOFLabel(tlab);
664 // Fill Reco-QA histos for Reconstruction
665 fHRecNClus->Fill(nc);
666 fHRecDist->Fill(mindist);
667 fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
668 fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
669 fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
670 fHRecSigZVsPWin->Fill(mom,dz);
672 // Fill Tree for on-the-fly offline Calibration
674 if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 )){
685 for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
686 for (Int_t ii=0;ii<6;ii++) delete [] clind[ii];
689 //_________________________________________________________________________
690 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
691 //--------------------------------------------------------------------
692 //This function loads the TOF clusters
693 //--------------------------------------------------------------------
695 Int_t npadX = AliTOFGeometry::NpadX();
696 Int_t npadZ = AliTOFGeometry::NpadZ();
697 Int_t nStripA = AliTOFGeometry::NStripA();
698 Int_t nStripB = AliTOFGeometry::NStripB();
699 Int_t nStripC = AliTOFGeometry::NStripC();
701 TBranch *branch=cTree->GetBranch("TOF");
703 AliError("can't get the branch with the TOF clusters !");
707 TClonesArray dummy("AliTOFcluster",10000), *clusters=&dummy;
708 branch->SetAddress(&clusters);
711 Int_t nc=clusters->GetEntriesFast();
712 fHDigNClus->Fill(nc);
714 AliInfo(Form("Number of clusters: %d",nc));
716 for (Int_t i=0; i<nc; i++) {
717 AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
718 fClusters[i]=new AliTOFcluster(*c); fN++;
720 // Fill Digits QA histos
722 Int_t isector = c->GetDetInd(0);
723 Int_t iplate = c->GetDetInd(1);
724 Int_t istrip = c->GetDetInd(2);
725 Int_t ipadX = c->GetDetInd(4);
726 Int_t ipadZ = c->GetDetInd(3);
728 Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
729 Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
731 Int_t stripOffset = 0;
737 stripOffset = nStripC;
740 stripOffset = nStripC+nStripB;
743 stripOffset = nStripC+nStripB+nStripA;
746 stripOffset = nStripC+nStripB+nStripA+nStripB;
749 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
752 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
753 Int_t phiindex=npadX*isector+ipadX+1;
754 fHDigClusMap->Fill(zindex,phiindex);
755 fHDigClusTime->Fill(time);
756 fHDigClusToT->Fill(tot);
763 //_________________________________________________________________________
764 void AliTOFtracker::UnloadClusters() {
765 //--------------------------------------------------------------------
766 //This function unloads TOF clusters
767 //--------------------------------------------------------------------
768 for (Int_t i=0; i<fN; i++) {
775 //_________________________________________________________________________
776 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
777 //--------------------------------------------------------------------
778 // This function returns the index of the nearest cluster
779 //--------------------------------------------------------------------
781 if (z <= fClusters[0]->GetZ()) return 0;
782 if (z > fClusters[fN-1]->GetZ()) return fN;
783 Int_t b=0, e=fN-1, m=(b+e)/2;
784 for (; b<e; m=(b+e)/2) {
785 if (z > fClusters[m]->GetZ()) b=m+1;
791 //_________________________________________________________________________
792 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
794 // Get track space point with index i
795 // Coordinates are in the global system
796 AliTOFcluster *cl = fClusters[index];
798 xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
799 xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
801 Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
802 Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
803 Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
804 Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
805 Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
806 Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
808 cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
809 cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
810 cov[2] = -cosphi*sinth*costh*sigmaz2;
811 cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
812 cov[4] = -sinphi*sinth*costh*sigmaz2;
813 cov[5] = costh*costh*sigmaz2;
814 p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
816 // Detector numbering scheme
817 Int_t nSector = AliTOFGeometry::NSectors();
818 Int_t nPlate = AliTOFGeometry::NPlates();
819 Int_t nStripA = AliTOFGeometry::NStripA();
820 Int_t nStripB = AliTOFGeometry::NStripB();
821 Int_t nStripC = AliTOFGeometry::NStripC();
823 Int_t isector = cl->GetDetInd(0);
824 if (isector >= nSector)
825 AliError(Form("Wrong sector number in TOF (%d) !",isector));
826 Int_t iplate = cl->GetDetInd(1);
827 if (iplate >= nPlate)
828 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
829 Int_t istrip = cl->GetDetInd(2);
831 Int_t stripOffset = 0;
837 stripOffset = nStripC;
840 stripOffset = nStripC+nStripB;
843 stripOffset = nStripC+nStripB+nStripA;
846 stripOffset = nStripC+nStripB+nStripA+nStripB;
849 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
853 Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
856 UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
857 p.SetVolumeID((UShort_t)volid);
860 //_________________________________________________________________________
861 void AliTOFtracker::InitCheckHists() {
863 //Init histos for Digits/Reco QA and Calibration
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");
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);
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.);
889 //_________________________________________________________________________
890 void AliTOFtracker::SaveCheckHists() {
892 //write histos for Digits/Reco QA and Calibration
894 TDirectory *dir = gDirectory;
896 TFile *logFileTOF = 0;
898 TSeqCollection *list = gROOT->GetListOfFiles();
899 int n = list->GetEntries();
900 for(int i=0; i<n; i++) {
901 logFile = (TFile*)list->At(i);
902 if (strstr(logFile->GetName(), "AliESDs.root")) break;
905 Bool_t isThere=kFALSE;
906 for(int i=0; i<n; i++) {
907 logFileTOF = (TFile*)list->At(i);
908 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
915 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
916 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
917 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
918 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
919 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
920 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
921 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
922 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
923 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
924 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
925 //fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
928 if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
930 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
931 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
932 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
933 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
934 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
935 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
936 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
937 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
938 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
939 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
940 fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
945 //_________________________________________________________________________
946 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) {
948 //dummy, for the moment
950 if(dist<AliTOFGeometry::ZPad()*0.5){
952 //place here the actual correction
958 //_________________________________________________________________________
959 Float_t AliTOFtracker::GetTimeZerofromT0(AliESDEvent *event) const {
961 //Returns TimeZero as measured by T0 detector
963 return event->GetT0();
965 //_________________________________________________________________________
966 Float_t AliTOFtracker::GetTimeZerofromTOF(AliESDEvent * /*event*/) const {
968 //dummy, for the moment. T0 algorithm using tracks on TOF
970 //place T0 algo here...