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 "AliTOFpidESD.h"
47 #include "AliTOFRecoParam.h"
48 #include "AliTOFReconstructor.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();
99 //_____________________________________________________________________________
100 AliTOFtracker::~AliTOFtracker() {
107 if(!(AliCDBManager::Instance()->GetCacheFlag())){
114 delete fHDigClusTime;
120 delete fHRecSigYVsPWin;
121 delete fHRecSigZVsPWin;
135 //_____________________________________________________________________________
136 Int_t AliTOFtracker::PropagateBack(AliESDEvent* event) {
138 // Gets seeds from ESD event and Match with TOF Clusters
141 // initialize RecoParam for current event
143 AliInfo("Initializing params for TOF... ");
145 fRecoParam = AliTOFReconstructor::GetRecoParam(); // instantiate reco param from STEER...
147 if (fRecoParam == 0x0) {
148 AliFatal("No Reco Param found for TOF!!!");
150 //fRecoParam->Dump();
151 //if(fRecoParam->GetApplyPbPbCuts())fRecoParam=fRecoParam->GetPbPbparam();
152 //fRecoParam->PrintParameters();
155 parPID[0]=fRecoParam->GetTimeResolution();
156 parPID[1]=fRecoParam->GetTimeNSigma();
157 fPid=new AliTOFpidESD(parPID);
159 //Initialise some counters
168 Int_t ntrk=event->GetNumberOfTracks();
170 TClonesArray &aESDTrack = *fSeeds;
173 //Load ESD tracks into a local Array of ESD Seeds
175 for (Int_t i=0; i<fNseeds; i++) {
176 AliESDtrack *t=event->GetTrack(i);
177 new(aESDTrack[i]) AliESDtrack(*t);
180 //Prepare ESD tracks candidates for TOF Matching
183 //First Step with Strict Matching Criterion
186 //Second Step with Looser Matching Criterion
189 AliInfo(Form("Number of matched tracks: %d",fnmatch));
190 AliInfo(Form("Number of good matched tracks: %d",fngoodmatch));
191 AliInfo(Form("Number of bad matched tracks: %d",fnbadmatch));
193 //Update the matched ESD tracks
195 for (Int_t i=0; i<ntrk; i++) {
196 AliESDtrack *t=event->GetTrack(i);
197 AliESDtrack *seed =(AliESDtrack*)fSeeds->UncheckedAt(i);
198 if(seed->GetTOFsignal()>0){
199 t->SetTOFsignal(seed->GetTOFsignal());
200 t->SetTOFcluster(seed->GetTOFcluster());
201 t->SetTOFsignalToT(seed->GetTOFsignalToT());
202 t->SetTOFsignalRaw(seed->GetTOFsignalRaw());
203 t->SetTOFsignalDz(seed->GetTOFsignalDz());
204 t->SetTOFCalChannel(seed->GetTOFCalChannel());
205 Int_t tlab[3]; seed->GetTOFLabel(tlab);
206 t->SetTOFLabel(tlab);
207 AliTOFtrack *track = new AliTOFtrack(*seed);
208 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
210 Double_t time[10]; t->GetIntegratedTimes(time);
211 AliDebug(1,Form("%d %f %f %f %f %f",i,
223 //Handle Time Zero information
225 Double_t timeZero=0.;
226 Double_t timeZeroMax=99999.;
227 Bool_t usetimeZero = fRecoParam->UseTimeZero();
228 Bool_t timeZeroFromT0 = fRecoParam->GetTimeZerofromT0();
229 Bool_t timeZeroFromTOF = fRecoParam->GetTimeZerofromTOF();
231 AliDebug(1,Form("Use Time Zero?: %d",usetimeZero));
232 AliDebug(1,Form("Time Zero from T0? : %d",timeZeroFromT0));
233 AliDebug(1,Form("Time Zero From TOF? : %d",timeZeroFromTOF));
237 timeZero=GetTimeZerofromT0(event);
239 if(timeZeroFromTOF && (timeZero>timeZeroMax || !timeZeroFromT0)){
240 timeZero=GetTimeZerofromTOF(event);
243 AliDebug(2,Form("time Zero used in PID: %f",timeZero));
245 fPid->MakePID(event,timeZero);
252 //_________________________________________________________________________
253 void AliTOFtracker::CollectESD() {
254 //prepare the set of ESD tracks to be matched to clusters in TOF
259 TClonesArray &aTOFTrack = *fTracks;
260 for (Int_t i=0; i<fNseeds; i++) {
262 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(i);
263 if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
265 // TRD 'good' tracks, already propagated at 372 cm
267 AliTOFtrack *track = new AliTOFtrack(*t); // New
268 Double_t x = track->GetX(); //New
269 if ( ( (t->GetStatus()&AliESDtrack::kTRDout)!=0 ) &&
270 ( x >= AliTOFGeometry::RinTOF() ) ) {
271 track->SetSeedIndex(i);
272 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
273 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
279 // Propagate the rest of TPCbp
282 if(track->PropagateToInnerTOF()){
283 track->SetSeedIndex(i);
284 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
285 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
293 AliInfo(Form("Number of TOF seeds %d",fNseedsTOF));
294 AliInfo(Form("Number of TOF seeds Type 1 %d",seedsTOF1));
295 AliInfo(Form("Number of TOF seeds Type 2 %d",seedsTOF2));
297 // Sort according uncertainties on track position
301 //_________________________________________________________________________
302 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
305 // Parameters used/regulating the reconstruction
307 static Float_t corrLen=0.75;
308 static Float_t detDepth=15.3;
309 static Float_t padDepth=0.5;
311 const Float_t kSpeedOfLight= 2.99792458e-2; // speed of light [cm/ps]
312 const Float_t kTimeOffset = 32.; // time offset for tracking algorithm [ps]
314 Float_t dY=AliTOFGeometry::XPad();
315 Float_t dZ=AliTOFGeometry::ZPad();
317 Float_t sensRadius = fRecoParam->GetSensRadius();
318 Float_t stepSize = fRecoParam->GetStepSize();
319 Float_t scaleFact = fRecoParam->GetWindowScaleFact();
320 Float_t dyMax=fRecoParam->GetWindowSizeMaxY();
321 Float_t dzMax=fRecoParam->GetWindowSizeMaxZ();
322 Float_t dCut=fRecoParam->GetDistanceCut();
323 Double_t maxChi2=fRecoParam->GetMaxChi2TRD();
324 Bool_t timeWalkCorr = fRecoParam->GetTimeWalkCorr();
326 AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++ \n");
327 AliDebug(1,Form("TOF sens radius: %f",sensRadius));
328 AliDebug(1,Form("TOF step size: %f",stepSize));
329 AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
330 AliDebug(1,Form("TOF Window max dy: %f",dyMax));
331 AliDebug(1,Form("TOF Window max dz: %f",dzMax));
332 AliDebug(1,Form("TOF distance Cut: %f",dCut));
333 AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
334 AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));
336 //Match ESD tracks to clusters in TOF
339 // Get the number of propagation steps
341 Int_t nSteps=(Int_t)(detDepth/stepSize);
343 //PH Arrays (moved outside of the loop)
345 Float_t * trackPos[4];
346 for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
347 Int_t * clind = new Int_t[fN];
349 for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
351 AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
352 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(track->GetSeedIndex());
353 if(t->GetTOFsignal()>0. ) continue;
354 AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
360 Float_t distZ[10000];
361 Float_t cxpos[10000];
362 Float_t crecL[10000];
363 TGeoHMatrix global[1000];
365 // Determine a window around the track
368 trackTOFin->GetExternalParameters(x,par);
370 trackTOFin->GetExternalCovariance(cov);
374 ((5*TMath::Sqrt(cov[0]) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius);
377 (5*TMath::Sqrt(cov[2]) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
379 Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
380 if (phi<-TMath::Pi())phi+=2*TMath::Pi();
381 if (phi>=TMath::Pi())phi-=2*TMath::Pi();
384 //upper limit on window's size.
386 if(dz> dzMax) dz=dzMax;
387 if(dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
392 // find the clusters in the window of the track
394 for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
396 AliTOFcluster *c=fClusters[k];
397 if (c->GetZ() > z+dz) break;
398 if (c->IsUsed()) continue;
400 if (!c->GetStatus()) {
401 AliDebug(1,"Cluster in channel declared bad!");
402 continue; // skip bad channels as declared in OCDB
405 Double_t dph=TMath::Abs(c->GetPhi()-phi);
406 if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
407 if (TMath::Abs(dph)>dphi) continue;
409 Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
410 Double_t p[2]={yc, c->GetZ()};
411 Double_t cov2[3]= {dY*dY/12., 0., dZ*dZ/12.};
412 if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov2) > maxChi2)continue;
417 ind[0]=c->GetDetInd(0);
418 ind[1]=c->GetDetInd(1);
419 ind[2]=c->GetDetInd(2);
420 ind[3]=c->GetDetInd(3);
421 ind[4]=c->GetDetInd(4);
422 fGeom->GetVolumePath(ind,path);
423 gGeoManager->cd(path);
424 global[nc] = *gGeoManager->GetCurrentMatrix();
428 //start fine propagation
430 Int_t nStepsDone = 0;
431 for( Int_t istep=0; istep<nSteps; istep++){
433 Float_t xs=AliTOFGeometry::RinTOF()+istep*0.1;
434 Double_t ymax=xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
437 Double_t ysect=trackTOFin->GetYat(xs,skip);
440 if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
443 } else if (ysect <-ymax) {
444 if (!trackTOFin->Rotate(-AliTOFGeometry::GetAlpha())) {
449 if(!trackTOFin->PropagateTo(xs)) {
455 // store the running point (Globalrf) - fine propagation
458 trackTOFin->GetXYZ(r);
459 trackPos[0][istep]= (Float_t) r[0];
460 trackPos[1][istep]= (Float_t) r[1];
461 trackPos[2][istep]= (Float_t) r[2];
462 trackPos[3][istep]= trackTOFin->GetIntegratedLength();
467 Bool_t accept = kFALSE;
468 Bool_t isInside = kFALSE;
469 for (Int_t istep=0; istep<nStepsDone; istep++) {
471 Float_t ctrackPos[3];
473 ctrackPos[0] = trackPos[0][istep];
474 ctrackPos[1] = trackPos[1][istep];
475 ctrackPos[2] = trackPos[2][istep];
477 //now see whether the track matches any of the TOF clusters
481 for (Int_t i=0; i<nc; i++) {
482 isInside = fGeom->IsInsideThePad(global[i],ctrackPos,dist3d);
485 Float_t yLoc = dist3d[1];
486 Float_t rLoc = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[2]*dist3d[2]);
487 accept = (TMath::Abs(yLoc)<padDepth*0.5 && rLoc<dCut);
488 AliDebug(2," I am in the case mLastStep==kTRUE ");
495 dist[nfound] = TMath::Sqrt(dist3d[0]*dist3d[0] +
496 dist3d[1]*dist3d[1] +
497 dist3d[2]*dist3d[2]);
498 AliDebug(1,Form(" dist3dLoc[0] = %f, dist3dLoc[1] = %f, dist3dLoc[2] = %f ",
499 dist3d[0],dist3d[1],dist3d[2]));
500 distZ[nfound] = dist3d[2]; // Z distance in the RF of the
501 // hit pad closest to the
502 // reconstructed track
504 AliDebug(1,Form(" dist3dLoc[2] = %f --- distZ[%d] = %f",
505 dist3d[2],nfound,distZ[nfound]));
507 crecL[nfound] = trackPos[3][istep];
508 index[nfound] = clind[i]; // store cluster id
509 cxpos[nfound] = AliTOFGeometry::RinTOF()+istep*0.1; //store prop.radius
511 if(accept &&!mLastStep)break;
514 } //end for on the clusters
515 if(accept &&!mLastStep)break;
516 } //end for on the steps
528 // now choose the cluster to be matched with the track.
533 Float_t mindist=1000.;
535 for (Int_t iclus= 0; iclus<nfound;iclus++){
536 if (dist[iclus]< mindist){
537 mindist = dist[iclus];
538 mindistZ = distZ[iclus]; // Z distance in the RF of the hit
539 // pad closest to the reconstructed
542 idclus =index[iclus];
543 recL=crecL[iclus]+corrLen*0.5;
548 AliTOFcluster *c=fClusters[idclus];
550 Float_t tiltangle = AliTOFGeometry::GetAngles(c->GetDetInd(1),c->GetDetInd(2))*TMath::DegToRad();
551 Float_t localCheck=-mindistZ;
552 localCheck/=TMath::Cos(tiltangle); // Z (tracking/ALICE RF) component of
553 // the distance between the
554 // reconstructed track and the
555 // TOF closest cluster
557 AliDebug(2, Form("%7d %7d %10d %10d %10d %10d %7d",
560 TMath::Abs(trackTOFin->GetLabel()),
561 c->GetLabel(0), c->GetLabel(1), c->GetLabel(2),
566 // Track length correction for matching Step 2
569 Float_t rc=TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
570 Float_t rt=TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
571 +trackPos[1][70]*trackPos[1][70]
572 +trackPos[2][70]*trackPos[2][70]);
574 recL=trackPos[3][70]+dlt;
578 (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
580 (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
582 (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
586 AliDebug(2,Form(" track label good %5d",trackTOFin->GetLabel()));
592 AliDebug(2,Form(" track label bad %5d",trackTOFin->GetLabel()));
598 // Store quantities to be used in the TOF Calibration
599 Float_t tToT=AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3; // in ns
600 t->SetTOFsignalToT(tToT);
601 Float_t rawTime=AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()+kTimeOffset; // RAW time,in ps
602 t->SetTOFsignalRaw(rawTime);
603 t->SetTOFsignalDz(mindistZ);
606 ind[0]=c->GetDetInd(0);
607 ind[1]=c->GetDetInd(1);
608 ind[2]=c->GetDetInd(2);
609 ind[3]=c->GetDetInd(3);
610 ind[4]=c->GetDetInd(4);
611 Int_t calindex = AliTOFGeometry::GetIndex(ind);
612 t->SetTOFCalChannel(calindex);
614 // keep track of the track labels in the matched cluster
616 tlab[0]=c->GetLabel(0);
617 tlab[1]=c->GetLabel(1);
618 tlab[2]=c->GetLabel(2);
619 AliDebug(2,Form(" tdc time of the matched track %6d = ",c->GetTDC()));
620 Double_t tof=AliTOFGeometry::TdcBinWidth()*c->GetTDC()+kTimeOffset; // in ps
621 AliDebug(2,Form(" tof time of the matched track: %f = ",tof));
622 Double_t tofcorr=tof;
623 if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
624 AliDebug(2,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));
625 //Set TOF time signal and pointer to the matched cluster
626 t->SetTOFsignal(tofcorr);
627 t->SetTOFcluster(idclus); // pointing to the recPoints tree
629 AliDebug(2,Form(" Setting TOF raw time: %f, z distance: %f corrected time: %f ",rawTime,mindistZ,tofcorr));
632 Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time); // in ps
633 Double_t mom=t->GetP();
634 AliDebug(1,Form(" Momentum for track %d -> %f", iseed,mom));
635 for(Int_t j=0;j<AliPID::kSPECIES;j++){
636 Double_t mass=AliPID::ParticleMass(j);
637 time[j]+=(recL-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(mom*mom+mass*mass)/mom;
640 AliTOFtrack *trackTOFout = new AliTOFtrack(*t);
641 trackTOFout->PropagateTo(xpos);
643 // Fill the track residual histograms.
644 FillResiduals(trackTOFout,c,kFALSE);
646 t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
647 t->SetIntegratedLength(recL);
648 t->SetIntegratedTimes(time);
649 t->SetTOFLabel(tlab);
652 // Fill Reco-QA histos for Reconstruction
653 fHRecNClus->Fill(nc);
654 fHRecDist->Fill(mindist);
655 fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
656 fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
657 fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
658 fHRecSigZVsPWin->Fill(mom,dz);
660 // Fill Tree for on-the-fly offline Calibration
662 if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 )){
673 for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
677 //_________________________________________________________________________
678 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
679 //--------------------------------------------------------------------
680 //This function loads the TOF clusters
681 //--------------------------------------------------------------------
683 Int_t npadX = AliTOFGeometry::NpadX();
684 Int_t npadZ = AliTOFGeometry::NpadZ();
685 Int_t nStripA = AliTOFGeometry::NStripA();
686 Int_t nStripB = AliTOFGeometry::NStripB();
687 Int_t nStripC = AliTOFGeometry::NStripC();
689 TBranch *branch=cTree->GetBranch("TOF");
691 AliError("can't get the branch with the TOF clusters !");
695 static TClonesArray dummy("AliTOFcluster",10000);
697 TClonesArray *clusters=&dummy;
698 branch->SetAddress(&clusters);
701 Int_t nc=clusters->GetEntriesFast();
702 fHDigNClus->Fill(nc);
704 AliInfo(Form("Number of clusters: %d",nc));
706 for (Int_t i=0; i<nc; i++) {
707 AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
708 fClusters[i]=new AliTOFcluster(*c); fN++;
710 // Fill Digits QA histos
712 Int_t isector = c->GetDetInd(0);
713 Int_t iplate = c->GetDetInd(1);
714 Int_t istrip = c->GetDetInd(2);
715 Int_t ipadX = c->GetDetInd(4);
716 Int_t ipadZ = c->GetDetInd(3);
718 Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
719 Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
721 Int_t stripOffset = 0;
727 stripOffset = nStripC;
730 stripOffset = nStripC+nStripB;
733 stripOffset = nStripC+nStripB+nStripA;
736 stripOffset = nStripC+nStripB+nStripA+nStripB;
739 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
742 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
743 Int_t phiindex=npadX*isector+ipadX+1;
744 fHDigClusMap->Fill(zindex,phiindex);
745 fHDigClusTime->Fill(time);
746 fHDigClusToT->Fill(tot);
753 //_________________________________________________________________________
754 void AliTOFtracker::UnloadClusters() {
755 //--------------------------------------------------------------------
756 //This function unloads TOF clusters
757 //--------------------------------------------------------------------
758 for (Int_t i=0; i<fN; i++) {
765 //_________________________________________________________________________
766 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
767 //--------------------------------------------------------------------
768 // This function returns the index of the nearest cluster
769 //--------------------------------------------------------------------
771 if (z <= fClusters[0]->GetZ()) return 0;
772 if (z > fClusters[fN-1]->GetZ()) return fN;
773 Int_t b=0, e=fN-1, m=(b+e)/2;
774 for (; b<e; m=(b+e)/2) {
775 if (z > fClusters[m]->GetZ()) b=m+1;
781 //_________________________________________________________________________
782 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
784 // Get track space point with index i
785 // Coordinates are in the global system
786 AliTOFcluster *cl = fClusters[index];
788 xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
789 xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
791 Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
792 Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
793 Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
794 Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
795 Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
796 Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
798 cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
799 cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
800 cov[2] = -cosphi*sinth*costh*sigmaz2;
801 cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
802 cov[4] = -sinphi*sinth*costh*sigmaz2;
803 cov[5] = costh*costh*sigmaz2;
804 p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
806 // Detector numbering scheme
807 Int_t nSector = AliTOFGeometry::NSectors();
808 Int_t nPlate = AliTOFGeometry::NPlates();
809 Int_t nStripA = AliTOFGeometry::NStripA();
810 Int_t nStripB = AliTOFGeometry::NStripB();
811 Int_t nStripC = AliTOFGeometry::NStripC();
813 Int_t isector = cl->GetDetInd(0);
814 if (isector >= nSector)
815 AliError(Form("Wrong sector number in TOF (%d) !",isector));
816 Int_t iplate = cl->GetDetInd(1);
817 if (iplate >= nPlate)
818 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
819 Int_t istrip = cl->GetDetInd(2);
821 Int_t stripOffset = 0;
827 stripOffset = nStripC;
830 stripOffset = nStripC+nStripB;
833 stripOffset = nStripC+nStripB+nStripA;
836 stripOffset = nStripC+nStripB+nStripA+nStripB;
839 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
843 Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
846 UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
847 p.SetVolumeID((UShort_t)volid);
850 //_________________________________________________________________________
851 void AliTOFtracker::InitCheckHists() {
853 //Init histos for Digits/Reco QA and Calibration
856 TDirectory *dir = gDirectory;
857 TFile *logFileTOF = 0;
859 TSeqCollection *list = gROOT->GetListOfFiles();
860 int n = list->GetEntries();
861 Bool_t isThere=kFALSE;
862 for(int i=0; i<n; i++) {
863 logFileTOF = (TFile*)list->At(i);
864 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
870 if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
873 fCalTree = new TTree("CalTree", "Tree for TOF calibration");
874 fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
875 fCalTree->Branch("ToT",&fToT,"TOFToT/F");
876 fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
877 fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
878 fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
879 fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
882 fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);
883 fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);
884 fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);
885 fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);
888 fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
889 fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
890 fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
891 fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
892 fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
893 fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
899 //_________________________________________________________________________
900 void AliTOFtracker::SaveCheckHists() {
902 //write histos for Digits/Reco QA and Calibration
904 TDirectory *dir = gDirectory;
905 TFile *logFileTOF = 0;
907 TSeqCollection *list = gROOT->GetListOfFiles();
908 int n = list->GetEntries();
909 Bool_t isThere=kFALSE;
910 for(int i=0; i<n; i++) {
911 logFileTOF = (TFile*)list->At(i);
912 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
919 AliError(Form("File TOFQA.root not found!! not wring histograms...."));
923 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
924 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
925 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
926 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
927 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
928 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
929 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
930 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
931 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
932 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
933 fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
938 //_________________________________________________________________________
939 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) {
941 //dummy, for the moment
943 if(dist<AliTOFGeometry::ZPad()*0.5){
945 //place here the actual correction
951 //_________________________________________________________________________
952 Float_t AliTOFtracker::GetTimeZerofromT0(AliESDEvent *event) const {
954 //Returns TimeZero as measured by T0 detector
956 return event->GetT0();
958 //_________________________________________________________________________
959 Float_t AliTOFtracker::GetTimeZerofromTOF(AliESDEvent * /*event*/) const {
961 //dummy, for the moment. T0 algorithm using tracks on TOF
963 //place T0 algo here...
967 //_________________________________________________________________________
969 void AliTOFtracker::FillClusterArray(TObjArray* arr) const
972 // Returns the TOF cluster array
978 for (Int_t i=0; i<fN; ++i) arr->Add(fClusters[i]);
981 //_________________________________________________________________________