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
93 // Gettimg the geometry
94 fGeom=new AliTOFGeometryV5();
95 // Read the reconstruction parameters from the OCDB
96 AliTOFcalib *calib = new AliTOFcalib(fGeom);
97 fRecoParam = (AliTOFRecoParam*)calib->ReadRecParFromCDB("TOF/Calib",-1);
98 if(!fRecoParam) {AliFatal("Exiting, no Reconstruction Parameters object found!!!");exit(0);}
99 if(fRecoParam->GetApplyPbPbCuts())fRecoParam=fRecoParam->GetPbPbparam();
101 parPID[0]=fRecoParam->GetTimeResolution();
102 parPID[1]=fRecoParam->GetTimeNSigma();
103 AliDebug(2,Form("TOF PID pars: Sigma= %f Range= %f",parPID[0],parPID[1]));
104 fPid=new AliTOFpidESD(parPID);
108 //_____________________________________________________________________________
109 AliTOFtracker::AliTOFtracker(const AliTOFtracker &t):
131 fHRecSigYVsPWin(0x0),
132 fHRecSigZVsPWin(0x0),
141 //AliTOFtracker copy Ctor
145 fNseedsTOF=t.fNseedsTOF;
146 fngoodmatch=t.fngoodmatch;
147 fnbadmatch=t.fnbadmatch;
148 fnunmatch=t.fnunmatch;
150 fRecoParam=t.fRecoParam;
158 //_____________________________________________________________________________
159 AliTOFtracker& AliTOFtracker::operator=(const AliTOFtracker &t)
161 //AliTOFtracker assignment operator
163 this->fNseeds=t.fNseeds;
164 this->fNseedsTOF=t.fNseedsTOF;
165 this->fngoodmatch=t.fngoodmatch;
166 this->fnbadmatch=t.fnbadmatch;
167 this->fnunmatch=t.fnunmatch;
168 this->fnmatch=t.fnmatch;
169 this->fRecoParam = t.fRecoParam;
170 this->fGeom = t.fGeom;
172 this->fSeeds=t.fSeeds;
173 this->fTracks=t.fTracks;
178 //_____________________________________________________________________________
179 AliTOFtracker::~AliTOFtracker() {
191 delete fHDigClusTime;
197 delete fHRecSigYVsPWin;
198 delete fHRecSigZVsPWin;
201 //_____________________________________________________________________________
202 Int_t AliTOFtracker::PropagateBack(AliESD* event) {
204 // Gets seeds from ESD event and Match with TOF Clusters
208 //Initialise some counters
217 Int_t ntrk=event->GetNumberOfTracks();
219 fSeeds= new TClonesArray("AliESDtrack",ntrk);
220 TClonesArray &aESDTrack = *fSeeds;
223 //Load ESD tracks into a local Array of ESD Seeds
225 for (Int_t i=0; i<fNseeds; i++) {
226 AliESDtrack *t=event->GetTrack(i);
227 new(aESDTrack[i]) AliESDtrack(*t);
230 //Prepare ESD tracks candidates for TOF Matching
233 //First Step with Strict Matching Criterion
236 //Second Step with Looser Matching Criterion
239 AliInfo(Form("Number of matched tracks: %d",fnmatch));
240 AliInfo(Form("Number of good matched tracks: %d",fngoodmatch));
241 AliInfo(Form("Number of bad matched tracks: %d",fnbadmatch));
243 //Update the matched ESD tracks
245 for (Int_t i=0; i<ntrk; i++) {
246 AliESDtrack *t=event->GetTrack(i);
247 AliESDtrack *seed =(AliESDtrack*)fSeeds->UncheckedAt(i);
248 if(seed->GetTOFsignal()>0){
249 t->SetTOFsignal(seed->GetTOFsignal());
250 t->SetTOFcluster(seed->GetTOFcluster());
251 t->SetTOFsignalToT(seed->GetTOFsignalToT());
252 t->SetTOFCalChannel(seed->GetTOFCalChannel());
253 Int_t tlab[3]; seed->GetTOFLabel(tlab);
254 t->SetTOFLabel(tlab);
255 AliTOFtrack *track = new AliTOFtrack(*seed);
256 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
261 //Handle Time Zero information
263 Double_t timeZero=0.;
264 Double_t timeZeroMax=99999.;
265 Bool_t usetimeZero = fRecoParam->UseTimeZero();
266 Bool_t timeZeroFromT0 = fRecoParam->GetTimeZerofromT0();
267 Bool_t timeZeroFromTOF = fRecoParam->GetTimeZerofromTOF();
269 AliDebug(2,Form("Use Time Zero?: %d",usetimeZero));
270 AliDebug(2,Form("Time Zero from T0? : %d",timeZeroFromT0));
271 AliDebug(2,Form("Time Zero From TOF? : %d",timeZeroFromTOF));
275 timeZero=GetTimeZerofromT0(event);
277 if(timeZeroFromTOF && (timeZero>timeZeroMax || !timeZeroFromT0)){
278 timeZero=GetTimeZerofromTOF(event);
281 AliDebug(2,Form("time Zero used in PID: %f",timeZero));
283 fPid->MakePID(event,timeZero);
298 //_________________________________________________________________________
299 void AliTOFtracker::CollectESD() {
300 //prepare the set of ESD tracks to be matched to clusters in TOF
302 fTracks= new TClonesArray("AliTOFtrack");
303 TClonesArray &aTOFTrack = *fTracks;
304 for (Int_t i=0; i<fNseeds; i++) {
306 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(i);
307 if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
309 // TRD good tracks, already propagated at 371 cm
311 AliTOFtrack *track = new AliTOFtrack(*t); // New
312 Double_t x = track->GetX(); //New
314 if (((t->GetStatus()&AliESDtrack::kTRDout)!=0 ) &&
315 ( x >= fGeom->RinTOF()) ){
316 track->SetSeedIndex(i);
317 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
318 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
323 // Propagate the rest of TPCbp
326 if(track->PropagateToInnerTOF()){
327 track->SetSeedIndex(i);
328 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
329 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
336 AliInfo(Form("Number of TOF seedds %i",fNseedsTOF));
338 // Sort according uncertainties on track position
342 //_________________________________________________________________________
343 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
346 // Parameters used/regulating the reconstruction
348 static Float_t corrLen=0.75;
349 static Float_t detDepth=15.3;
351 Float_t dY=fGeom->XPad();
352 Float_t dZ=fGeom->ZPad();
354 Float_t sensRadius = fRecoParam->GetSensRadius();
355 Float_t stepSize = fRecoParam->GetStepSize();
356 Float_t scaleFact = fRecoParam->GetWindowScaleFact();
357 Float_t dyMax=fRecoParam->GetWindowSizeMaxY();
358 Float_t dzMax=fRecoParam->GetWindowSizeMaxZ();
359 Float_t dCut=fRecoParam->GetDistanceCut();
360 Double_t maxChi2=fRecoParam->GetMaxChi2();
362 Bool_t timeWalkCorr = fRecoParam->GetTimeWalkCorr();
364 AliDebug(2,Form("Time Walk Correction? : %d",timeWalkCorr));
366 AliDebug(2,"TOF RecPars pars: \n");
367 AliDebug(2,Form("TOF sens radius: %f",sensRadius));
368 AliDebug(2,Form("TOF step size: %f",stepSize));
369 AliDebug(2,Form("TOF Window scale factor: %f",scaleFact));
370 AliDebug(2,Form("TOF Window max dy: %f",dyMax));
371 AliDebug(2,Form("TOF Window max dz: %f",dzMax));
372 AliDebug(2,Form("TOF distance Cut: %f",dCut));
373 AliDebug(2,Form("TOF Max Chi2: %f",maxChi2));
375 //Match ESD tracks to clusters in TOF
378 // Get the number of propagation steps
380 Int_t nSteps=(Int_t)(detDepth/stepSize);
382 AliTOFcalib *calib = new AliTOFcalib(fGeom);
384 //PH Arrays (moved outside of the loop)
386 Float_t * trackPos[4];
387 for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
389 for (Int_t ii=0;ii<6;ii++) clind[ii] = new Int_t[fN];
391 for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
393 AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
394 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(track->GetSeedIndex());
395 if(t->GetTOFsignal()>0. ) continue;
396 AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
402 Float_t distZ[10000];
403 Float_t cxpos[10000];
404 Float_t crecL[10000];
405 TGeoHMatrix global[1000];
407 // Determine a window around the track
410 trackTOFin->GetExternalParameters(x,par);
412 trackTOFin->GetExternalCovariance(cov);
416 ((5*TMath::Sqrt(cov[0]) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius);
419 (5*TMath::Sqrt(cov[2]) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
421 Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
422 if (phi<-TMath::Pi())phi+=2*TMath::Pi();
423 if (phi>=TMath::Pi())phi-=2*TMath::Pi();
426 //upper limit on window's size.
428 if(dz> dzMax) dz=dzMax;
429 if(dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
434 // find the clusters in the window of the track
436 for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
438 AliTOFcluster *c=fClusters[k];
439 if (c->GetZ() > z+dz) break;
440 if (c->IsUsed()) continue;
442 if (!c->GetStatus()) continue; // skip bad channels as declared in OCDB
444 Double_t dph=TMath::Abs(c->GetPhi()-phi);
445 if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
446 if (TMath::Abs(dph)>dphi) continue;
448 Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
449 Double_t p[2]={yc, c->GetZ()};
450 Double_t cov[3]= {dY*dY/12., 0., dZ*dZ/12.};
451 if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov) > maxChi2)continue;
453 clind[0][nc] = c->GetDetInd(0);
454 clind[1][nc] = c->GetDetInd(1);
455 clind[2][nc] = c->GetDetInd(2);
456 clind[3][nc] = c->GetDetInd(3);
457 clind[4][nc] = c->GetDetInd(4);
466 fGeom->GetVolumePath(ind,path);
467 gGeoManager->cd(path);
468 global[nc] = *gGeoManager->GetCurrentMatrix();
472 //start fine propagation
474 Int_t nStepsDone = 0;
475 for( Int_t istep=0; istep<nSteps; istep++){
477 Float_t xs=fGeom->RinTOF()+istep*0.1;
478 Double_t ymax=xs*TMath::Tan(0.5*fGeom->GetAlpha());
481 Double_t ysect=trackTOFin->GetYat(xs,skip);
484 if (!trackTOFin->Rotate(fGeom->GetAlpha())) {
487 } else if (ysect <-ymax) {
488 if (!trackTOFin->Rotate(fGeom->GetAlpha())) {
493 if(!trackTOFin->PropagateTo(xs)) {
499 // store the running point (Globalrf) - fine propagation
502 trackTOFin->GetXYZ(r);
503 trackPos[0][istep]= (Float_t) r[0];
504 trackPos[1][istep]= (Float_t) r[1];
505 trackPos[2][istep]= (Float_t) r[2];
506 trackPos[3][istep]= trackTOFin->GetIntegratedLength();
511 for (Int_t istep=0; istep<nStepsDone; istep++) {
513 Bool_t isInside =kFALSE;
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
522 for (Int_t i=0; i<nc; i++){
524 cind[0]= clind[0][i];
525 cind[1]= clind[1][i];
526 cind[2]= clind[2][i];
527 cind[3]= clind[3][i];
528 cind[4]= clind[4][i];
529 Bool_t accept = kFALSE;
530 if( mLastStep)accept = (fGeom->DistanceToPad(cind,global[i],ctrackPos)<dCut);
531 if(!mLastStep)accept = (fGeom->IsInsideThePad(cind,global[i],ctrackPos));
533 if(!mLastStep)isInside=kTRUE;
535 dist[nfound]=fGeom->DistanceToPad(cind,global[i],ctrackPos,dist3d);
536 distZ[nfound]=dist3d[2];
537 crecL[nfound]=trackPos[3][istep];
538 index[nfound]=clind[5][i]; // store cluster id
539 cxpos[nfound]=fGeom->RinTOF()+istep*0.1; //store prop.radius
543 } //end for on the clusters
547 } //end for on the steps
559 // now choose the cluster to be matched with the track.
564 Float_t mindist=1000.;
566 for (Int_t iclus= 0; iclus<nfound;iclus++){
567 if (dist[iclus]< mindist){
568 mindist = dist[iclus];
569 mindistZ = distZ[iclus];
571 idclus =index[iclus];
572 recL=crecL[iclus]+corrLen*0.5;
576 AliTOFcluster *c=fClusters[idclus];
579 // Track length correction for matching Step 2
582 Float_t rc=TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
583 Float_t rt=TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
584 +trackPos[1][70]*trackPos[1][70]
585 +trackPos[2][70]*trackPos[2][70]);
587 recL=trackPos[3][70]+dlt;
591 (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
593 (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
595 (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
599 AliDebug(2,Form(" track label good %5i",trackTOFin->GetLabel()));
605 AliDebug(2,Form(" track label bad %5i",trackTOFin->GetLabel()));
611 // Store quantities to be used in the TOF Calibration
612 Float_t tToT=fGeom->TdcBinWidth()*c->GetToT()*1E-3; // in ns
613 t->SetTOFsignalToT(tToT);
615 ind[0]=c->GetDetInd(0);
616 ind[1]=c->GetDetInd(1);
617 ind[2]=c->GetDetInd(2);
618 ind[3]=c->GetDetInd(3);
619 ind[4]=c->GetDetInd(4);
620 Int_t calindex = calib->GetIndex(ind);
621 t->SetTOFCalChannel(calindex);
623 // keep track of the track labels in the matched cluster
625 tlab[0]=c->GetLabel(0);
626 tlab[1]=c->GetLabel(1);
627 tlab[2]=c->GetLabel(2);
629 Double_t tof=fGeom->TdcBinWidth()*c->GetTDC()+32; // in ps
630 if(timeWalkCorr)tof=CorrectTimeWalk(mindistZ,tof);
631 t->SetTOFsignal(tof);
632 t->SetTOFcluster(idclus); // pointing to the recPoints tree
633 Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time);
634 Double_t mom=t->GetP();
635 for(Int_t j=0;j<AliPID::kSPECIES;j++){
636 Double_t mass=AliPID::ParticleMass(j);
637 time[j]+=(recL-trackPos[3][0])/3e-2*TMath::Sqrt(mom*mom+mass*mass)/mom;
640 AliTOFtrack *trackTOFout = new AliTOFtrack(*t);
641 trackTOFout->PropagateTo(xpos);
642 t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
643 t->SetIntegratedLength(recL);
644 t->SetIntegratedTimes(time);
645 t->SetTOFLabel(tlab);
646 // Fill Reco-QA histos for Reconstruction
647 fHRecNClus->Fill(nc);
648 fHRecDist->Fill(mindist);
649 fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
650 fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
651 fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
652 fHRecSigZVsPWin->Fill(mom,dz);
654 // Fill Tree for on-the-fly offline Calibration
656 if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 )){
657 Float_t rawtime=fGeom->TdcBinWidth()*c->GetTDCRAW()+32; // RAW time,in ps
668 for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
669 for (Int_t ii=0;ii<6;ii++) delete [] clind[ii];
672 //_________________________________________________________________________
673 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
674 //--------------------------------------------------------------------
675 //This function loads the TOF clusters
676 //--------------------------------------------------------------------
678 Int_t npadX = fGeom->NpadX();
679 Int_t npadZ = fGeom->NpadZ();
680 Int_t nStripA = fGeom->NStripA();
681 Int_t nStripB = fGeom->NStripB();
682 Int_t nStripC = fGeom->NStripC();
684 TBranch *branch=cTree->GetBranch("TOF");
686 AliError("can't get the branch with the TOF clusters !");
690 TClonesArray dummy("AliTOFcluster",10000), *clusters=&dummy;
691 branch->SetAddress(&clusters);
694 Int_t nc=clusters->GetEntriesFast();
695 fHDigNClus->Fill(nc);
697 AliInfo(Form("Number of clusters: %d",nc));
699 for (Int_t i=0; i<nc; i++) {
700 AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
701 fClusters[i]=new AliTOFcluster(*c); fN++;
703 // Fill Digits QA histos
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);
711 Float_t time =(fGeom->TdcBinWidth()*c->GetTDC())*1E-3; // in ns
712 Float_t tot = (fGeom->TdcBinWidth()*c->GetToT())*1E-3;//in ns
714 Int_t stripOffset = 0;
720 stripOffset = nStripC;
723 stripOffset = nStripC+nStripB;
726 stripOffset = nStripC+nStripB+nStripA;
729 stripOffset = nStripC+nStripB+nStripA+nStripB;
732 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
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);
746 //_________________________________________________________________________
747 void AliTOFtracker::UnloadClusters() {
748 //--------------------------------------------------------------------
749 //This function unloads TOF clusters
750 //--------------------------------------------------------------------
751 for (Int_t i=0; i<fN; i++) {
758 //_________________________________________________________________________
759 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
760 //--------------------------------------------------------------------
761 // This function returns the index of the nearest cluster
762 //--------------------------------------------------------------------
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;
774 //_________________________________________________________________________
775 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
777 // Get track space point with index i
778 // Coordinates are in the global system
779 AliTOFcluster *cl = fClusters[index];
781 xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
782 xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
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 = fGeom->GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
787 Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
788 Float_t sigmay2 = fGeom->XPad()*fGeom->XPad()/12.;
789 Float_t sigmaz2 = fGeom->ZPad()*fGeom->ZPad()/12.;
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);
799 // Detector numbering scheme
800 Int_t nSector = fGeom->NSectors();
801 Int_t nPlate = fGeom->NPlates();
802 Int_t nStripA = fGeom->NStripA();
803 Int_t nStripB = fGeom->NStripB();
804 Int_t nStripC = fGeom->NStripC();
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);
814 Int_t stripOffset = 0;
820 stripOffset = nStripC;
823 stripOffset = nStripC+nStripB;
826 stripOffset = nStripC+nStripB+nStripA;
829 stripOffset = nStripC+nStripB+nStripA+nStripB;
832 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
836 Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
839 UShort_t volid = AliAlignObj::LayerToVolUID(AliAlignObj::kTOF,idet);
840 p.SetVolumeID((UShort_t)volid);
843 //_________________________________________________________________________
844 void AliTOFtracker::InitCheckHists() {
846 //Init histos for Digits/Reco QA and Calibration
849 fCalTree = new TTree("CalTree", "Tree for TOF calibration");
850 fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
851 fCalTree->Branch("ToT",&fToT,"TOFToT/F");
852 fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
853 fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
854 fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
855 fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
858 fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);
859 fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);
860 fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);
861 fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);
864 fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
865 fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
866 fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
867 fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
868 fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
869 fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
872 //_________________________________________________________________________
873 void AliTOFtracker::SaveCheckHists() {
875 //write histos for Digits/Reco QA and Calibration
877 TDirectory *dir = gDirectory;
879 TFile *logFileTOF = 0;
881 TSeqCollection *list = gROOT->GetListOfFiles();
882 int n = list->GetEntries();
883 for(int i=0; i<n; i++) {
884 logFile = (TFile*)list->At(i);
885 if (strstr(logFile->GetName(), "AliESDs.root")) break;
888 Bool_t isThere=kFALSE;
889 for(int i=0; i<n; i++) {
890 logFileTOF = (TFile*)list->At(i);
891 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
898 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
899 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
900 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
901 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
902 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
903 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
904 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
905 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
906 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
907 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
908 fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
911 if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
913 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
914 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
915 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
916 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
917 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
918 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
919 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
920 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
921 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
922 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
923 fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
928 //_________________________________________________________________________
929 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) {
931 //dummy, for the moment
933 if(dist<fGeom->ZPad()*0.5){
935 //place here the actual correction
941 //_________________________________________________________________________
942 Float_t AliTOFtracker::GetTimeZerofromT0(AliESD *event) {
944 //Returns TimeZero as measured by T0 detector
946 return event->GetT0();
948 //_________________________________________________________________________
949 Float_t AliTOFtracker::GetTimeZerofromTOF(AliESD *event) {
951 //dummy, for the moment. T0 algorithm using tracks on TOF
953 //place T0 algo here...