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 <TSeqCollection.h>
33 #include <TClonesArray.h>
34 #include <TGeoManager.h>
39 #include "AliGeomManager.h"
40 #include "AliESDtrack.h"
41 #include "AliESDEvent.h"
43 #include "AliTrackPointArray.h"
44 #include "AliCDBManager.h"
46 #include "AliTOFcalib.h"
47 #include "AliTOFpidESD.h"
48 #include "AliTOFRecoParam.h"
49 #include "AliTOFcluster.h"
50 #include "AliTOFGeometry.h"
51 #include "AliTOFtracker.h"
52 #include "AliTOFtrack.h"
54 extern TGeoManager *gGeoManager;
57 ClassImp(AliTOFtracker)
59 //_____________________________________________________________________________
60 AliTOFtracker::AliTOFtracker():
71 fTracks(new TClonesArray("AliTOFtrack")),
72 fSeeds(new TClonesArray("AliESDtrack")),
91 //AliTOFtracker main Ctor
93 // Gettimg the geometry
94 fGeom= new AliTOFGeometry();
95 // Read the reconstruction parameters from the OCDB
96 AliTOFcalib* calib=new AliTOFcalib();
97 fRecoParam = (AliTOFRecoParam*)calib->ReadRecParFromCDB("TOF/Calib",-1);
98 if(fRecoParam->GetApplyPbPbCuts())fRecoParam=fRecoParam->GetPbPbparam();
100 parPID[0]=fRecoParam->GetTimeResolution();
101 parPID[1]=fRecoParam->GetTimeNSigma();
102 fPid=new AliTOFpidESD(parPID);
106 //_____________________________________________________________________________
107 AliTOFtracker::AliTOFtracker(const AliTOFtracker &t):
119 fTracks(new TClonesArray("AliTOFtrack")),
120 fSeeds(new TClonesArray("AliESDtrack")),
129 fHRecSigYVsPWin(0x0),
130 fHRecSigZVsPWin(0x0),
139 //AliTOFtracker copy Ctor
143 fNseedsTOF=t.fNseedsTOF;
144 fngoodmatch=t.fngoodmatch;
145 fnbadmatch=t.fnbadmatch;
146 fnunmatch=t.fnunmatch;
148 fRecoParam=t.fRecoParam;
155 //_____________________________________________________________________________
156 AliTOFtracker& AliTOFtracker::operator=(const AliTOFtracker &t)
158 //AliTOFtracker assignment operator
160 this->fNseeds=t.fNseeds;
161 this->fNseedsTOF=t.fNseedsTOF;
162 this->fngoodmatch=t.fngoodmatch;
163 this->fnbadmatch=t.fnbadmatch;
164 this->fnunmatch=t.fnunmatch;
165 this->fnmatch=t.fnmatch;
166 this->fRecoParam = t.fRecoParam;
169 this->fSeeds=t.fSeeds;
170 this->fTracks=t.fTracks;
175 //_____________________________________________________________________________
176 AliTOFtracker::~AliTOFtracker() {
183 if(!(AliCDBManager::Instance()->GetCacheFlag())){
190 delete fHDigClusTime;
196 delete fHRecSigYVsPWin;
197 delete fHRecSigZVsPWin;
211 //_____________________________________________________________________________
212 Int_t AliTOFtracker::PropagateBack(AliESDEvent* event) {
214 // Gets seeds from ESD event and Match with TOF Clusters
218 //Initialise some counters
227 Int_t ntrk=event->GetNumberOfTracks();
229 TClonesArray &aESDTrack = *fSeeds;
232 //Load ESD tracks into a local Array of ESD Seeds
234 for (Int_t i=0; i<fNseeds; i++) {
235 AliESDtrack *t=event->GetTrack(i);
236 new(aESDTrack[i]) AliESDtrack(*t);
239 //Prepare ESD tracks candidates for TOF Matching
242 //First Step with Strict Matching Criterion
245 //Second Step with Looser Matching Criterion
248 AliInfo(Form("Number of matched tracks: %d",fnmatch));
249 AliInfo(Form("Number of good matched tracks: %d",fngoodmatch));
250 AliInfo(Form("Number of bad matched tracks: %d",fnbadmatch));
252 //Update the matched ESD tracks
254 for (Int_t i=0; i<ntrk; i++) {
255 AliESDtrack *t=event->GetTrack(i);
256 AliESDtrack *seed =(AliESDtrack*)fSeeds->UncheckedAt(i);
257 if(seed->GetTOFsignal()>0){
258 t->SetTOFsignal(seed->GetTOFsignal());
259 t->SetTOFcluster(seed->GetTOFcluster());
260 t->SetTOFsignalToT(seed->GetTOFsignalToT());
261 t->SetTOFsignalRaw(seed->GetTOFsignalRaw());
262 t->SetTOFsignalDz(seed->GetTOFsignalDz());
263 t->SetTOFCalChannel(seed->GetTOFCalChannel());
264 Int_t tlab[3]; seed->GetTOFLabel(tlab);
265 t->SetTOFLabel(tlab);
266 AliTOFtrack *track = new AliTOFtrack(*seed);
267 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
272 //Handle Time Zero information
274 Double_t timeZero=0.;
275 Double_t timeZeroMax=99999.;
276 Bool_t usetimeZero = fRecoParam->UseTimeZero();
277 Bool_t timeZeroFromT0 = fRecoParam->GetTimeZerofromT0();
278 Bool_t timeZeroFromTOF = fRecoParam->GetTimeZerofromTOF();
280 AliDebug(1,Form("Use Time Zero?: %d",usetimeZero));
281 AliDebug(1,Form("Time Zero from T0? : %d",timeZeroFromT0));
282 AliDebug(1,Form("Time Zero From TOF? : %d",timeZeroFromTOF));
286 timeZero=GetTimeZerofromT0(event);
288 if(timeZeroFromTOF && (timeZero>timeZeroMax || !timeZeroFromT0)){
289 timeZero=GetTimeZerofromTOF(event);
292 AliDebug(2,Form("time Zero used in PID: %f",timeZero));
294 fPid->MakePID(event,timeZero);
301 //_________________________________________________________________________
302 void AliTOFtracker::CollectESD() {
303 //prepare the set of ESD tracks to be matched to clusters in TOF
308 TClonesArray &aTOFTrack = *fTracks;
309 for (Int_t i=0; i<fNseeds; i++) {
311 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(i);
312 if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
314 // TRD 'good' tracks, already propagated at 372 cm
316 AliTOFtrack *track = new AliTOFtrack(*t); // New
317 Double_t x = track->GetX(); //New
318 if (((t->GetStatus()&AliESDtrack::kTRDout)!=0 ) &&
319 ( x >= AliTOFGeometry::RinTOF()) ){
320 track->SetSeedIndex(i);
321 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
322 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
328 // Propagate the rest of TPCbp
331 if(track->PropagateToInnerTOF()){
332 track->SetSeedIndex(i);
333 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
334 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
342 AliInfo(Form("Number of TOF seeds %i",fNseedsTOF));
343 AliInfo(Form("Number of TOF seeds Type 1 %i",seedsTOF1));
344 AliInfo(Form("Number of TOF seeds Type 2 %i",seedsTOF2));
346 // Sort according uncertainties on track position
350 //_________________________________________________________________________
351 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
354 // Parameters used/regulating the reconstruction
356 static Float_t corrLen=0.75;
357 static Float_t detDepth=15.3;
358 static Float_t padDepth=0.5;
360 Float_t dY=AliTOFGeometry::XPad();
361 Float_t dZ=AliTOFGeometry::ZPad();
363 Float_t sensRadius = fRecoParam->GetSensRadius();
364 Float_t stepSize = fRecoParam->GetStepSize();
365 Float_t scaleFact = fRecoParam->GetWindowScaleFact();
366 Float_t dyMax=fRecoParam->GetWindowSizeMaxY();
367 Float_t dzMax=fRecoParam->GetWindowSizeMaxZ();
368 Float_t dCut=fRecoParam->GetDistanceCut();
369 Double_t maxChi2=fRecoParam->GetMaxChi2TRD();
370 Bool_t timeWalkCorr = fRecoParam->GetTimeWalkCorr();
372 AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++ \n");
373 AliDebug(1,Form("TOF sens radius: %f",sensRadius));
374 AliDebug(1,Form("TOF step size: %f",stepSize));
375 AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
376 AliDebug(1,Form("TOF Window max dy: %f",dyMax));
377 AliDebug(1,Form("TOF Window max dz: %f",dzMax));
378 AliDebug(1,Form("TOF distance Cut: %f",dCut));
379 AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
380 AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));
382 //Match ESD tracks to clusters in TOF
385 // Get the number of propagation steps
387 Int_t nSteps=(Int_t)(detDepth/stepSize);
389 //PH Arrays (moved outside of the loop)
391 Float_t * trackPos[4];
392 for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
393 Int_t * clind = 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()) {
447 AliDebug(1,"Cluster in channel declared bad!");
448 continue; // skip bad channels as declared in OCDB
451 Double_t dph=TMath::Abs(c->GetPhi()-phi);
452 if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
453 if (TMath::Abs(dph)>dphi) continue;
455 Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
456 Double_t p[2]={yc, c->GetZ()};
457 Double_t cov2[3]= {dY*dY/12., 0., dZ*dZ/12.};
458 if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov2) > maxChi2)continue;
463 ind[0]=c->GetDetInd(0);
464 ind[1]=c->GetDetInd(1);
465 ind[2]=c->GetDetInd(2);
466 ind[3]=c->GetDetInd(3);
467 ind[4]=c->GetDetInd(4);
468 fGeom->GetVolumePath(ind,path);
469 gGeoManager->cd(path);
470 global[nc] = *gGeoManager->GetCurrentMatrix();
474 //start fine propagation
476 Int_t nStepsDone = 0;
477 for( Int_t istep=0; istep<nSteps; istep++){
479 Float_t xs=AliTOFGeometry::RinTOF()+istep*0.1;
480 Double_t ymax=xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
483 Double_t ysect=trackTOFin->GetYat(xs,skip);
486 if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
489 } else if (ysect <-ymax) {
490 if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
495 if(!trackTOFin->PropagateTo(xs)) {
501 // store the running point (Globalrf) - fine propagation
504 trackTOFin->GetXYZ(r);
505 trackPos[0][istep]= (Float_t) r[0];
506 trackPos[1][istep]= (Float_t) r[1];
507 trackPos[2][istep]= (Float_t) r[2];
508 trackPos[3][istep]= trackTOFin->GetIntegratedLength();
513 Bool_t accept = kFALSE;
514 Bool_t isInside =kFALSE;
515 for (Int_t istep=0; istep<nStepsDone; istep++) {
517 Float_t ctrackPos[3];
519 ctrackPos[0]= trackPos[0][istep];
520 ctrackPos[1]= trackPos[1][istep];
521 ctrackPos[2]= trackPos[2][istep];
523 //now see whether the track matches any of the TOF clusters
527 for (Int_t i=0; i<nc; i++){
528 isInside=fGeom->IsInsideThePad(global[i],ctrackPos,dist3d);
531 Float_t xLoc=dist3d[0];
532 Float_t rLoc=TMath::Sqrt(dist3d[1]*dist3d[1]+dist3d[2]*dist3d[2]);
533 accept = (TMath::Abs(xLoc)<padDepth*0.5 && rLoc<dCut);
539 dist[nfound]=TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[1]*dist3d[1]+dist3d[2]*dist3d[2]);
540 distZ[nfound]=dist3d[2];
541 crecL[nfound]=trackPos[3][istep];
542 index[nfound]=clind[i]; // store cluster id
543 cxpos[nfound]=AliTOFGeometry::RinTOF()+istep*0.1; //store prop.radius
545 if(accept &&!mLastStep)break;
547 } //end for on the clusters
548 if(accept &&!mLastStep)break;
549 } //end for on the steps
561 // now choose the cluster to be matched with the track.
566 Float_t mindist=1000.;
568 for (Int_t iclus= 0; iclus<nfound;iclus++){
569 if (dist[iclus]< mindist){
570 mindist = dist[iclus];
571 mindistZ = distZ[iclus];
573 idclus =index[iclus];
574 recL=crecL[iclus]+corrLen*0.5;
579 AliTOFcluster *c=fClusters[idclus];
581 AliDebug(2, Form("%7i %7i %10i %10i %10i %10i %7i",
584 TMath::Abs(trackTOFin->GetLabel()),
585 c->GetLabel(0), c->GetLabel(1), c->GetLabel(2),
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=AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3; // in ns
624 t->SetTOFsignalToT(tToT);
625 Float_t rawTime=AliTOFGeometry::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 = AliTOFGeometry::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=AliTOFGeometry::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];
694 //_________________________________________________________________________
695 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
696 //--------------------------------------------------------------------
697 //This function loads the TOF clusters
698 //--------------------------------------------------------------------
700 Int_t npadX = AliTOFGeometry::NpadX();
701 Int_t npadZ = AliTOFGeometry::NpadZ();
702 Int_t nStripA = AliTOFGeometry::NStripA();
703 Int_t nStripB = AliTOFGeometry::NStripB();
704 Int_t nStripC = AliTOFGeometry::NStripC();
706 TBranch *branch=cTree->GetBranch("TOF");
708 AliError("can't get the branch with the TOF clusters !");
712 static TClonesArray dummy("AliTOFcluster",10000);
714 TClonesArray *clusters=&dummy;
715 branch->SetAddress(&clusters);
718 Int_t nc=clusters->GetEntriesFast();
719 fHDigNClus->Fill(nc);
721 AliInfo(Form("Number of clusters: %d",nc));
723 for (Int_t i=0; i<nc; i++) {
724 AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
725 fClusters[i]=new AliTOFcluster(*c); fN++;
727 // Fill Digits QA histos
729 Int_t isector = c->GetDetInd(0);
730 Int_t iplate = c->GetDetInd(1);
731 Int_t istrip = c->GetDetInd(2);
732 Int_t ipadX = c->GetDetInd(4);
733 Int_t ipadZ = c->GetDetInd(3);
735 Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
736 Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
738 Int_t stripOffset = 0;
744 stripOffset = nStripC;
747 stripOffset = nStripC+nStripB;
750 stripOffset = nStripC+nStripB+nStripA;
753 stripOffset = nStripC+nStripB+nStripA+nStripB;
756 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
759 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
760 Int_t phiindex=npadX*isector+ipadX+1;
761 fHDigClusMap->Fill(zindex,phiindex);
762 fHDigClusTime->Fill(time);
763 fHDigClusToT->Fill(tot);
770 //_________________________________________________________________________
771 void AliTOFtracker::UnloadClusters() {
772 //--------------------------------------------------------------------
773 //This function unloads TOF clusters
774 //--------------------------------------------------------------------
775 for (Int_t i=0; i<fN; i++) {
782 //_________________________________________________________________________
783 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
784 //--------------------------------------------------------------------
785 // This function returns the index of the nearest cluster
786 //--------------------------------------------------------------------
788 if (z <= fClusters[0]->GetZ()) return 0;
789 if (z > fClusters[fN-1]->GetZ()) return fN;
790 Int_t b=0, e=fN-1, m=(b+e)/2;
791 for (; b<e; m=(b+e)/2) {
792 if (z > fClusters[m]->GetZ()) b=m+1;
798 //_________________________________________________________________________
799 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
801 // Get track space point with index i
802 // Coordinates are in the global system
803 AliTOFcluster *cl = fClusters[index];
805 xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
806 xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
808 Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
809 Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
810 Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
811 Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
812 Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
813 Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
815 cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
816 cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
817 cov[2] = -cosphi*sinth*costh*sigmaz2;
818 cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
819 cov[4] = -sinphi*sinth*costh*sigmaz2;
820 cov[5] = costh*costh*sigmaz2;
821 p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
823 // Detector numbering scheme
824 Int_t nSector = AliTOFGeometry::NSectors();
825 Int_t nPlate = AliTOFGeometry::NPlates();
826 Int_t nStripA = AliTOFGeometry::NStripA();
827 Int_t nStripB = AliTOFGeometry::NStripB();
828 Int_t nStripC = AliTOFGeometry::NStripC();
830 Int_t isector = cl->GetDetInd(0);
831 if (isector >= nSector)
832 AliError(Form("Wrong sector number in TOF (%d) !",isector));
833 Int_t iplate = cl->GetDetInd(1);
834 if (iplate >= nPlate)
835 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
836 Int_t istrip = cl->GetDetInd(2);
838 Int_t stripOffset = 0;
844 stripOffset = nStripC;
847 stripOffset = nStripC+nStripB;
850 stripOffset = nStripC+nStripB+nStripA;
853 stripOffset = nStripC+nStripB+nStripA+nStripB;
856 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
860 Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
863 UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
864 p.SetVolumeID((UShort_t)volid);
867 //_________________________________________________________________________
868 void AliTOFtracker::InitCheckHists() {
870 //Init histos for Digits/Reco QA and Calibration
873 TDirectory *dir = gDirectory;
874 TFile *logFileTOF = 0;
876 TSeqCollection *list = gROOT->GetListOfFiles();
877 int n = list->GetEntries();
878 Bool_t isThere=kFALSE;
879 for(int i=0; i<n; i++) {
880 logFileTOF = (TFile*)list->At(i);
881 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
887 if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
890 fCalTree = new TTree("CalTree", "Tree for TOF calibration");
891 fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
892 fCalTree->Branch("ToT",&fToT,"TOFToT/F");
893 fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
894 fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
895 fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
896 fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
899 fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);
900 fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);
901 fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);
902 fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);
905 fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
906 fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
907 fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
908 fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
909 fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
910 fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
916 //_________________________________________________________________________
917 void AliTOFtracker::SaveCheckHists() {
919 //write histos for Digits/Reco QA and Calibration
921 TDirectory *dir = gDirectory;
922 TFile *logFileTOF = 0;
924 TSeqCollection *list = gROOT->GetListOfFiles();
925 int n = list->GetEntries();
926 Bool_t isThere=kFALSE;
927 for(int i=0; i<n; i++) {
928 logFileTOF = (TFile*)list->At(i);
929 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
936 AliError(Form("File TOFQA.root not found!! not wring histograms...."));
940 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
941 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
942 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
943 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
944 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
945 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
946 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
947 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
948 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
949 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
950 fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
955 //_________________________________________________________________________
956 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) {
958 //dummy, for the moment
960 if(dist<AliTOFGeometry::ZPad()*0.5){
962 //place here the actual correction
968 //_________________________________________________________________________
969 Float_t AliTOFtracker::GetTimeZerofromT0(AliESDEvent *event) const {
971 //Returns TimeZero as measured by T0 detector
973 return event->GetT0();
975 //_________________________________________________________________________
976 Float_t AliTOFtracker::GetTimeZerofromTOF(AliESDEvent * /*event*/) const {
978 //dummy, for the moment. T0 algorithm using tracks on TOF
980 //place T0 algo here...
984 //_________________________________________________________________________
986 void AliTOFtracker::FillClusterArray(TObjArray* arr) const
989 // Returns the TOF cluster array
995 for (Int_t i=0; i<fN; ++i) arr->Add(fClusters[i]);