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 "AliAlignObj.h"
40 #include "AliGeomManager.h"
41 #include "AliESDtrack.h"
42 #include "AliESDEvent.h"
44 #include "AliTrackPointArray.h"
45 #include "AliCDBManager.h"
47 #include "AliTOFcalib.h"
48 #include "AliTOFpidESD.h"
49 #include "AliTOFRecoParam.h"
50 #include "AliTOFcluster.h"
51 #include "AliTOFGeometry.h"
52 #include "AliTOFtracker.h"
53 #include "AliTOFtrack.h"
55 extern TGeoManager *gGeoManager;
58 ClassImp(AliTOFtracker)
60 //_____________________________________________________________________________
61 AliTOFtracker::AliTOFtracker():
92 //AliTOFtracker main Ctor
94 // Gettimg the geometry
95 fGeom= new AliTOFGeometry();
96 // Read the reconstruction parameters from the OCDB
97 AliTOFcalib* calib=new AliTOFcalib();
98 fRecoParam = (AliTOFRecoParam*)calib->ReadRecParFromCDB("TOF/Calib",-1);
99 if(fRecoParam->GetApplyPbPbCuts())fRecoParam=fRecoParam->GetPbPbparam();
101 parPID[0]=fRecoParam->GetTimeResolution();
102 parPID[1]=fRecoParam->GetTimeNSigma();
103 fPid=new AliTOFpidESD(parPID);
107 //_____________________________________________________________________________
108 AliTOFtracker::AliTOFtracker(const AliTOFtracker &t):
130 fHRecSigYVsPWin(0x0),
131 fHRecSigZVsPWin(0x0),
140 //AliTOFtracker copy Ctor
144 fNseedsTOF=t.fNseedsTOF;
145 fngoodmatch=t.fngoodmatch;
146 fnbadmatch=t.fnbadmatch;
147 fnunmatch=t.fnunmatch;
149 fRecoParam=t.fRecoParam;
156 //_____________________________________________________________________________
157 AliTOFtracker& AliTOFtracker::operator=(const AliTOFtracker &t)
159 //AliTOFtracker assignment operator
161 this->fNseeds=t.fNseeds;
162 this->fNseedsTOF=t.fNseedsTOF;
163 this->fngoodmatch=t.fngoodmatch;
164 this->fnbadmatch=t.fnbadmatch;
165 this->fnunmatch=t.fnunmatch;
166 this->fnmatch=t.fnmatch;
167 this->fRecoParam = t.fRecoParam;
170 this->fSeeds=t.fSeeds;
171 this->fTracks=t.fTracks;
176 //_____________________________________________________________________________
177 AliTOFtracker::~AliTOFtracker() {
184 if(!(AliCDBManager::Instance()->GetCacheFlag())){
191 delete fHDigClusTime;
197 delete fHRecSigYVsPWin;
198 delete fHRecSigZVsPWin;
201 //_____________________________________________________________________________
202 Int_t AliTOFtracker::PropagateBack(AliESDEvent* 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->SetTOFsignalRaw(seed->GetTOFsignalRaw());
253 t->SetTOFsignalDz(seed->GetTOFsignalDz());
254 t->SetTOFCalChannel(seed->GetTOFCalChannel());
255 Int_t tlab[3]; seed->GetTOFLabel(tlab);
256 t->SetTOFLabel(tlab);
257 AliTOFtrack *track = new AliTOFtrack(*seed);
258 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
263 //Handle Time Zero information
265 Double_t timeZero=0.;
266 Double_t timeZeroMax=99999.;
267 Bool_t usetimeZero = fRecoParam->UseTimeZero();
268 Bool_t timeZeroFromT0 = fRecoParam->GetTimeZerofromT0();
269 Bool_t timeZeroFromTOF = fRecoParam->GetTimeZerofromTOF();
271 AliDebug(1,Form("Use Time Zero?: %d",usetimeZero));
272 AliDebug(1,Form("Time Zero from T0? : %d",timeZeroFromT0));
273 AliDebug(1,Form("Time Zero From TOF? : %d",timeZeroFromTOF));
277 timeZero=GetTimeZerofromT0(event);
279 if(timeZeroFromTOF && (timeZero>timeZeroMax || !timeZeroFromT0)){
280 timeZero=GetTimeZerofromTOF(event);
283 AliDebug(2,Form("time Zero used in PID: %f",timeZero));
285 fPid->MakePID(event,timeZero);
300 //_________________________________________________________________________
301 void AliTOFtracker::CollectESD() {
302 //prepare the set of ESD tracks to be matched to clusters in TOF
307 fTracks= new TClonesArray("AliTOFtrack");
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 371 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];
394 for (Int_t ii=0;ii<6;ii++) clind[ii] = new Int_t[fN];
396 for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
398 AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
399 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(track->GetSeedIndex());
400 if(t->GetTOFsignal()>0. ) continue;
401 AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
407 Float_t distZ[10000];
408 Float_t cxpos[10000];
409 Float_t crecL[10000];
410 TGeoHMatrix global[1000];
412 // Determine a window around the track
415 trackTOFin->GetExternalParameters(x,par);
417 trackTOFin->GetExternalCovariance(cov);
421 ((5*TMath::Sqrt(cov[0]) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius);
424 (5*TMath::Sqrt(cov[2]) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
426 Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
427 if (phi<-TMath::Pi())phi+=2*TMath::Pi();
428 if (phi>=TMath::Pi())phi-=2*TMath::Pi();
431 //upper limit on window's size.
433 if(dz> dzMax) dz=dzMax;
434 if(dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
439 // find the clusters in the window of the track
441 for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
443 AliTOFcluster *c=fClusters[k];
444 if (c->GetZ() > z+dz) break;
445 if (c->IsUsed()) continue;
447 if (!c->GetStatus()) continue; // skip bad channels as declared in OCDB
449 Double_t dph=TMath::Abs(c->GetPhi()-phi);
450 if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
451 if (TMath::Abs(dph)>dphi) continue;
453 Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
454 Double_t p[2]={yc, c->GetZ()};
455 Double_t cov[3]= {dY*dY/12., 0., dZ*dZ/12.};
456 if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov) > maxChi2)continue;
458 clind[0][nc] = c->GetDetInd(0);
459 clind[1][nc] = c->GetDetInd(1);
460 clind[2][nc] = c->GetDetInd(2);
461 clind[3][nc] = c->GetDetInd(3);
462 clind[4][nc] = c->GetDetInd(4);
471 fGeom->GetVolumePath(ind,path);
472 gGeoManager->cd(path);
473 global[nc] = *gGeoManager->GetCurrentMatrix();
477 //start fine propagation
479 Int_t nStepsDone = 0;
480 for( Int_t istep=0; istep<nSteps; istep++){
482 Float_t xs=AliTOFGeometry::RinTOF()+istep*0.1;
483 Double_t ymax=xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
486 Double_t ysect=trackTOFin->GetYat(xs,skip);
489 if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
492 } else if (ysect <-ymax) {
493 if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
498 if(!trackTOFin->PropagateTo(xs)) {
504 // store the running point (Globalrf) - fine propagation
507 trackTOFin->GetXYZ(r);
508 trackPos[0][istep]= (Float_t) r[0];
509 trackPos[1][istep]= (Float_t) r[1];
510 trackPos[2][istep]= (Float_t) r[2];
511 trackPos[3][istep]= trackTOFin->GetIntegratedLength();
516 Bool_t accept = kFALSE;
517 Bool_t isInside =kFALSE;
518 for (Int_t istep=0; istep<nStepsDone; istep++) {
520 Float_t ctrackPos[3];
522 ctrackPos[0]= trackPos[0][istep];
523 ctrackPos[1]= trackPos[1][istep];
524 ctrackPos[2]= trackPos[2][istep];
526 //now see whether the track matches any of the TOF clusters
530 for (Int_t i=0; i<nc; i++){
532 cind[0]= clind[0][i];
533 cind[1]= clind[1][i];
534 cind[2]= clind[2][i];
535 cind[3]= clind[3][i];
536 cind[4]= clind[4][i];
537 isInside=fGeom->IsInsideThePad(global[i],ctrackPos,dist3d);
540 Float_t xLoc=dist3d[0];
541 Float_t rLoc=TMath::Sqrt(dist3d[1]*dist3d[1]+dist3d[2]*dist3d[2]);
542 accept = (TMath::Abs(xLoc)<padDepth*0.5 && rLoc<dCut);
548 dist[nfound]=TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[1]*dist3d[1]+dist3d[2]*dist3d[2]);
549 distZ[nfound]=dist3d[2];
550 crecL[nfound]=trackPos[3][istep];
551 index[nfound]=clind[5][i]; // store cluster id
552 cxpos[nfound]=AliTOFGeometry::RinTOF()+istep*0.1; //store prop.radius
554 if(accept &&!mLastStep)break;
556 } //end for on the clusters
557 if(accept &&!mLastStep)break;
558 } //end for on the steps
570 // now choose the cluster to be matched with the track.
575 Float_t mindist=1000.;
577 for (Int_t iclus= 0; iclus<nfound;iclus++){
578 if (dist[iclus]< mindist){
579 mindist = dist[iclus];
580 mindistZ = distZ[iclus];
582 idclus =index[iclus];
583 recL=crecL[iclus]+corrLen*0.5;
588 AliTOFcluster *c=fClusters[idclus];
591 // Track length correction for matching Step 2
594 Float_t rc=TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
595 Float_t rt=TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
596 +trackPos[1][70]*trackPos[1][70]
597 +trackPos[2][70]*trackPos[2][70]);
599 recL=trackPos[3][70]+dlt;
603 (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
605 (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
607 (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
611 AliDebug(2,Form(" track label good %5i",trackTOFin->GetLabel()));
617 AliDebug(2,Form(" track label bad %5i",trackTOFin->GetLabel()));
623 // Store quantities to be used in the TOF Calibration
624 Float_t tToT=AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3; // in ns
625 t->SetTOFsignalToT(tToT);
626 Float_t rawTime=AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()+32; // RAW time,in ps
627 t->SetTOFsignalRaw(rawTime);
628 t->SetTOFsignalDz(mindistZ);
629 AliDebug(2,Form(" Setting TOF raw time: %f, z distance: %f time: %f = ",rawTime,mindistZ));
631 ind[0]=c->GetDetInd(0);
632 ind[1]=c->GetDetInd(1);
633 ind[2]=c->GetDetInd(2);
634 ind[3]=c->GetDetInd(3);
635 ind[4]=c->GetDetInd(4);
636 Int_t calindex = AliTOFGeometry::GetIndex(ind);
637 t->SetTOFCalChannel(calindex);
639 // keep track of the track labels in the matched cluster
641 tlab[0]=c->GetLabel(0);
642 tlab[1]=c->GetLabel(1);
643 tlab[2]=c->GetLabel(2);
644 AliDebug(2,Form(" tdc time of the matched track %i = ",c->GetTDC()));
645 Double_t tof=AliTOFGeometry::TdcBinWidth()*c->GetTDC()+32; // in ps
646 AliDebug(2,Form(" tof time of the matched track: %f = ",tof));
647 Double_t tofcorr=tof;
648 if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
649 AliDebug(2,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));
650 //Set TOF time signal and pointer to the matched cluster
651 t->SetTOFsignal(tofcorr);
652 t->SetTOFcluster(idclus); // pointing to the recPoints tree
655 Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time);
656 Double_t mom=t->GetP();
657 for(Int_t j=0;j<AliPID::kSPECIES;j++){
658 Double_t mass=AliPID::ParticleMass(j);
659 time[j]+=(recL-trackPos[3][0])/3e-2*TMath::Sqrt(mom*mom+mass*mass)/mom;
662 AliTOFtrack *trackTOFout = new AliTOFtrack(*t);
663 trackTOFout->PropagateTo(xpos);
664 t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
665 t->SetIntegratedLength(recL);
666 t->SetIntegratedTimes(time);
667 t->SetTOFLabel(tlab);
670 // Fill Reco-QA histos for Reconstruction
671 fHRecNClus->Fill(nc);
672 fHRecDist->Fill(mindist);
673 fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
674 fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
675 fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
676 fHRecSigZVsPWin->Fill(mom,dz);
678 // Fill Tree for on-the-fly offline Calibration
680 if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 )){
691 for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
692 for (Int_t ii=0;ii<6;ii++) delete [] clind[ii];
695 //_________________________________________________________________________
696 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
697 //--------------------------------------------------------------------
698 //This function loads the TOF clusters
699 //--------------------------------------------------------------------
701 Int_t npadX = AliTOFGeometry::NpadX();
702 Int_t npadZ = AliTOFGeometry::NpadZ();
703 Int_t nStripA = AliTOFGeometry::NStripA();
704 Int_t nStripB = AliTOFGeometry::NStripB();
705 Int_t nStripC = AliTOFGeometry::NStripC();
707 TBranch *branch=cTree->GetBranch("TOF");
709 AliError("can't get the branch with the TOF clusters !");
713 TClonesArray dummy("AliTOFcluster",10000), *clusters=&dummy;
714 branch->SetAddress(&clusters);
717 Int_t nc=clusters->GetEntriesFast();
718 fHDigNClus->Fill(nc);
720 AliInfo(Form("Number of clusters: %d",nc));
722 for (Int_t i=0; i<nc; i++) {
723 AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
724 fClusters[i]=new AliTOFcluster(*c); fN++;
726 // Fill Digits QA histos
728 Int_t isector = c->GetDetInd(0);
729 Int_t iplate = c->GetDetInd(1);
730 Int_t istrip = c->GetDetInd(2);
731 Int_t ipadX = c->GetDetInd(4);
732 Int_t ipadZ = c->GetDetInd(3);
734 Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
735 Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
737 Int_t stripOffset = 0;
743 stripOffset = nStripC;
746 stripOffset = nStripC+nStripB;
749 stripOffset = nStripC+nStripB+nStripA;
752 stripOffset = nStripC+nStripB+nStripA+nStripB;
755 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
758 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
759 Int_t phiindex=npadX*isector+ipadX+1;
760 fHDigClusMap->Fill(zindex,phiindex);
761 fHDigClusTime->Fill(time);
762 fHDigClusToT->Fill(tot);
769 //_________________________________________________________________________
770 void AliTOFtracker::UnloadClusters() {
771 //--------------------------------------------------------------------
772 //This function unloads TOF clusters
773 //--------------------------------------------------------------------
774 for (Int_t i=0; i<fN; i++) {
781 //_________________________________________________________________________
782 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
783 //--------------------------------------------------------------------
784 // This function returns the index of the nearest cluster
785 //--------------------------------------------------------------------
787 if (z <= fClusters[0]->GetZ()) return 0;
788 if (z > fClusters[fN-1]->GetZ()) return fN;
789 Int_t b=0, e=fN-1, m=(b+e)/2;
790 for (; b<e; m=(b+e)/2) {
791 if (z > fClusters[m]->GetZ()) b=m+1;
797 //_________________________________________________________________________
798 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
800 // Get track space point with index i
801 // Coordinates are in the global system
802 AliTOFcluster *cl = fClusters[index];
804 xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
805 xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
807 Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
808 Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
809 Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
810 Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
811 Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
812 Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
814 cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
815 cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
816 cov[2] = -cosphi*sinth*costh*sigmaz2;
817 cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
818 cov[4] = -sinphi*sinth*costh*sigmaz2;
819 cov[5] = costh*costh*sigmaz2;
820 p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
822 // Detector numbering scheme
823 Int_t nSector = AliTOFGeometry::NSectors();
824 Int_t nPlate = AliTOFGeometry::NPlates();
825 Int_t nStripA = AliTOFGeometry::NStripA();
826 Int_t nStripB = AliTOFGeometry::NStripB();
827 Int_t nStripC = AliTOFGeometry::NStripC();
829 Int_t isector = cl->GetDetInd(0);
830 if (isector >= nSector)
831 AliError(Form("Wrong sector number in TOF (%d) !",isector));
832 Int_t iplate = cl->GetDetInd(1);
833 if (iplate >= nPlate)
834 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
835 Int_t istrip = cl->GetDetInd(2);
837 Int_t stripOffset = 0;
843 stripOffset = nStripC;
846 stripOffset = nStripC+nStripB;
849 stripOffset = nStripC+nStripB+nStripA;
852 stripOffset = nStripC+nStripB+nStripA+nStripB;
855 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
859 Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
862 UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
863 p.SetVolumeID((UShort_t)volid);
866 //_________________________________________________________________________
867 void AliTOFtracker::InitCheckHists() {
869 //Init histos for Digits/Reco QA and Calibration
872 fCalTree = new TTree("CalTree", "Tree for TOF calibration");
873 fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
874 fCalTree->Branch("ToT",&fToT,"TOFToT/F");
875 fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
876 fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
877 fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
878 fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
881 fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);
882 fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);
883 fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);
884 fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);
887 fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
888 fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
889 fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
890 fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
891 fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
892 fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
895 //_________________________________________________________________________
896 void AliTOFtracker::SaveCheckHists() {
898 //write histos for Digits/Reco QA and Calibration
900 TDirectory *dir = gDirectory;
902 TFile *logFileTOF = 0;
904 TSeqCollection *list = gROOT->GetListOfFiles();
905 int n = list->GetEntries();
906 for(int i=0; i<n; i++) {
907 logFile = (TFile*)list->At(i);
908 if (strstr(logFile->GetName(), "AliESDs.root")) break;
911 Bool_t isThere=kFALSE;
912 for(int i=0; i<n; i++) {
913 logFileTOF = (TFile*)list->At(i);
914 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
921 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
922 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
923 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
924 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
925 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
926 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
927 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
928 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
929 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
930 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
931 //fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
934 if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
936 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
937 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
938 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
939 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
940 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
941 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
942 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
943 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
944 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
945 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
946 fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
951 //_________________________________________________________________________
952 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) {
954 //dummy, for the moment
956 if(dist<AliTOFGeometry::ZPad()*0.5){
958 //place here the actual correction
964 //_________________________________________________________________________
965 Float_t AliTOFtracker::GetTimeZerofromT0(AliESDEvent *event) const {
967 //Returns TimeZero as measured by T0 detector
969 return event->GetT0();
971 //_________________________________________________________________________
972 Float_t AliTOFtracker::GetTimeZerofromTOF(AliESDEvent * /*event*/) const {
974 //dummy, for the moment. T0 algorithm using tracks on TOF
976 //place T0 algo here...