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
142 AliDebug(1,"Initializing params for TOF");
144 fRecoParam = AliTOFReconstructor::GetRecoParam(); // instantiate reco param from STEER...
146 if (fRecoParam == 0x0) {
147 AliFatal("No Reco Param found for TOF!!!");
149 //fRecoParam->Dump();
150 //if(fRecoParam->GetApplyPbPbCuts())fRecoParam=fRecoParam->GetPbPbparam();
151 //fRecoParam->PrintParameters();
154 parPID[0]=fRecoParam->GetTimeResolution();
155 parPID[1]=fRecoParam->GetTimeNSigma();
156 fPid=new AliTOFpidESD(parPID);
158 //Initialise some counters
167 Int_t ntrk=event->GetNumberOfTracks();
169 TClonesArray &aESDTrack = *fSeeds;
172 //Load ESD tracks into a local Array of ESD Seeds
174 for (Int_t i=0; i<fNseeds; i++) {
175 AliESDtrack *t=event->GetTrack(i);
176 new(aESDTrack[i]) AliESDtrack(*t);
179 //Prepare ESD tracks candidates for TOF Matching
182 //First Step with Strict Matching Criterion
185 //Second Step with Looser Matching Criterion
188 AliInfo(Form("Number of matched tracks = %d (good = %d, bad = %d)",fnmatch,fngoodmatch,fnbadmatch));
190 //Update the matched ESD tracks
192 for (Int_t i=0; i<ntrk; i++) {
193 AliESDtrack *t=event->GetTrack(i);
194 AliESDtrack *seed =(AliESDtrack*)fSeeds->UncheckedAt(i);
196 if ( (seed->GetStatus()&AliESDtrack::kTOFin)!=0 ) {
197 t->SetStatus(AliESDtrack::kTOFin);
198 //if(seed->GetTOFsignal()>0){
199 if ( (seed->GetStatus()&AliESDtrack::kTOFout)!=0 ) {
200 t->SetStatus(AliESDtrack::kTOFout);
201 t->SetTOFsignal(seed->GetTOFsignal());
202 t->SetTOFcluster(seed->GetTOFcluster());
203 t->SetTOFsignalToT(seed->GetTOFsignalToT());
204 t->SetTOFsignalRaw(seed->GetTOFsignalRaw());
205 t->SetTOFsignalDz(seed->GetTOFsignalDz());
206 t->SetTOFsignalDx(seed->GetTOFsignalDx());
207 t->SetTOFCalChannel(seed->GetTOFCalChannel());
208 Int_t tlab[3]; seed->GetTOFLabel(tlab);
209 t->SetTOFLabel(tlab);
211 Double_t alphaA = dynamic_cast<AliExternalTrackParam*>(t)->GetAlpha();
212 Double_t xA = dynamic_cast<AliExternalTrackParam*>(t)->GetX();
213 Double_t yA = dynamic_cast<AliExternalTrackParam*>(t)->GetY();
214 Double_t zA = dynamic_cast<AliExternalTrackParam*>(t)->GetZ();
215 Double_t p1A = dynamic_cast<AliExternalTrackParam*>(t)->GetSnp();
216 Double_t p2A = dynamic_cast<AliExternalTrackParam*>(t)->GetTgl();
217 Double_t p3A = dynamic_cast<AliExternalTrackParam*>(t)->GetSigned1Pt();
218 const Double_t *covA = dynamic_cast<AliExternalTrackParam*>(t)->GetCovariance();
220 // Make attention, please:
221 // AliESDtrack::fTOFInfo array does not be stored in the AliESDs.root file
222 // it is there only for a check during the reconstruction step.
223 Float_t info[10]; seed->GetTOFInfo(info);
225 AliDebug(3,Form(" distance=%f; residual in the pad reference frame: dX=%f, dZ=%f", info[0],info[1],info[2]));
228 // by calling the AliESDtrack::UpdateTrackParams,
229 // the current track parameters are changed
230 // and it could cause refit problems.
231 // We need to update only the following track parameters:
232 // the track length and expected times.
233 // Removed AliESDtrack::UpdateTrackParams call
234 // Called AliESDtrack::SetIntegratedTimes(...) and
235 // AliESDtrack::SetIntegratedLength() routines.
237 AliTOFtrack *track = new AliTOFtrack(*seed);
238 t->UpdateTrackParams(track,AliESDtrack::kTOFout); // to be checked - AdC
240 Double_t time[10]; t->GetIntegratedTimes(time);
243 Double_t time[10]; seed->GetIntegratedTimes(time);
244 t->SetIntegratedTimes(time);
246 Double_t length = seed->GetIntegratedLength();
247 t->SetIntegratedLength(length);
249 Double_t alphaB = dynamic_cast<AliExternalTrackParam*>(t)->GetAlpha();
250 Double_t xB = dynamic_cast<AliExternalTrackParam*>(t)->GetX();
251 Double_t yB = dynamic_cast<AliExternalTrackParam*>(t)->GetY();
252 Double_t zB = dynamic_cast<AliExternalTrackParam*>(t)->GetZ();
253 Double_t p1B = dynamic_cast<AliExternalTrackParam*>(t)->GetSnp();
254 Double_t p2B = dynamic_cast<AliExternalTrackParam*>(t)->GetTgl();
255 Double_t p3B = dynamic_cast<AliExternalTrackParam*>(t)->GetSigned1Pt();
256 const Double_t *covB = dynamic_cast<AliExternalTrackParam*>(t)->GetCovariance();
257 AliDebug(2,"Track params -now(before)-:");
258 AliDebug(2,Form(" X: %f(%f), Y: %f(%f), Z: %f(%f) --- alpha: %f(%f)",
263 AliDebug(2,Form(" p1: %f(%f), p2: %f(%f), p3: %f(%f)",
267 AliDebug(2,Form(" cov1: %f(%f), cov2: %f(%f), cov3: %f(%f)"
268 " cov4: %f(%f), cov5: %f(%f), cov6: %f(%f)"
269 " cov7: %f(%f), cov8: %f(%f), cov9: %f(%f)"
270 " cov10: %f(%f), cov11: %f(%f), cov12: %f(%f)"
271 " cov13: %f(%f), cov14: %f(%f), cov15: %f(%f)",
288 AliDebug(3,Form(" TOF params: %6d %f %f %f %f %f %6d %3d %f %f %f %f %f %f",
290 t->GetTOFsignalRaw(),
292 t->GetTOFsignalToT(),
295 t->GetTOFCalChannel(),
297 t->GetIntegratedLength(),
298 time[0], time[1], time[2], time[3], time[4]
305 //Handle Time Zero information
307 Double_t timeZero=0.;
308 Double_t timeZeroMax=99999.;
309 Bool_t usetimeZero = fRecoParam->UseTimeZero();
310 Bool_t timeZeroFromT0 = fRecoParam->GetTimeZerofromT0();
311 Bool_t timeZeroFromTOF = fRecoParam->GetTimeZerofromTOF();
313 AliDebug(2,Form("Use Time Zero?: %d",usetimeZero));
314 AliDebug(2,Form("Time Zero from T0? : %d",timeZeroFromT0));
315 AliDebug(2,Form("Time Zero From TOF? : %d",timeZeroFromTOF));
319 timeZero=GetTimeZerofromT0(event);
321 if(timeZeroFromTOF && (timeZero>timeZeroMax || !timeZeroFromT0)){
322 timeZero=GetTimeZerofromTOF(event);
325 AliDebug(2,Form("time Zero used in PID: %f",timeZero));
327 fPid->MakePID(event,timeZero);
334 //_________________________________________________________________________
335 void AliTOFtracker::CollectESD() {
336 //prepare the set of ESD tracks to be matched to clusters in TOF
341 TClonesArray &aTOFTrack = *fTracks;
342 for (Int_t i=0; i<fNseeds; i++) {
344 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(i);
345 if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
347 AliTOFtrack *track = new AliTOFtrack(*t); // New
348 Float_t x = (Float_t)track->GetX(); //New
350 // TRD 'good' tracks, already propagated at 371 cm
351 if ( ( (t->GetStatus()&AliESDtrack::kTRDout)!=0 ) &&
352 ( x >= AliTOFGeometry::Rmin() ) ) {
353 if ( track->PropagateToInnerTOF() ) {
355 AliDebug(1,Form(" TRD propagated track till rho = %fcm."
356 " And then the track has been propagated till rho = %fcm.",
357 x, (Float_t)track->GetX()));
359 track->SetSeedIndex(i);
360 t->UpdateTrackParams(track,AliESDtrack::kTOFin);
361 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
368 // Propagate the rest of TPCbp
370 if ( track->PropagateToInnerTOF() ) {
372 AliDebug(1,Form(" TPC propagated track till rho = %fcm."
373 " And then the track has been propagated till rho = %fcm.",
374 x, (Float_t)track->GetX()));
376 track->SetSeedIndex(i);
377 t->UpdateTrackParams(track,AliESDtrack::kTOFin);
378 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
386 AliInfo(Form("Number of TOF seeds = %d (Type 1 = %d, Type 2 = %d)",fNseedsTOF,seedsTOF1,seedsTOF2));
388 // Sort according uncertainties on track position
392 //_________________________________________________________________________
393 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
396 // Parameters used/regulating the reconstruction
398 static Float_t corrLen=0.75;
399 static Float_t detDepth=15.3;
400 static Float_t padDepth=0.5;
402 const Float_t kSpeedOfLight= 2.99792458e-2; // speed of light [cm/ps]
403 const Float_t kTimeOffset = 32.; // time offset for tracking algorithm [ps]
405 Float_t dY=AliTOFGeometry::XPad();
406 Float_t dZ=AliTOFGeometry::ZPad();
408 Float_t sensRadius = fRecoParam->GetSensRadius();
409 Float_t stepSize = fRecoParam->GetStepSize();
410 Float_t scaleFact = fRecoParam->GetWindowScaleFact();
411 Float_t dyMax=fRecoParam->GetWindowSizeMaxY();
412 Float_t dzMax=fRecoParam->GetWindowSizeMaxZ();
413 Float_t dCut=fRecoParam->GetDistanceCut();
414 Double_t maxChi2=fRecoParam->GetMaxChi2TRD();
415 Bool_t timeWalkCorr = fRecoParam->GetTimeWalkCorr();
417 AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++ \n");
418 AliDebug(1,Form("TOF sens radius: %f",sensRadius));
419 AliDebug(1,Form("TOF step size: %f",stepSize));
420 AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
421 AliDebug(1,Form("TOF Window max dy: %f",dyMax));
422 AliDebug(1,Form("TOF Window max dz: %f",dzMax));
423 AliDebug(1,Form("TOF distance Cut: %f",dCut));
424 AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
425 AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));
428 //Match ESD tracks to clusters in TOF
430 // Get the number of propagation steps
431 Int_t nSteps=(Int_t)(detDepth/stepSize);
433 //PH Arrays (moved outside of the loop)
434 Float_t * trackPos[4];
435 for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
436 Int_t * clind = new Int_t[fN];
439 const Int_t kNfoundMax = 10000; // related to nSteps value
440 Int_t index[kNfoundMax];
441 Float_t dist[kNfoundMax];
442 Float_t distZ[kNfoundMax];
443 Float_t distY[kNfoundMax]; // delta(rhoXphi) // AdC
444 Float_t cxpos[kNfoundMax];
445 Float_t crecL[kNfoundMax];
446 const Int_t kNclusterMax = 1000; // related to fN value
447 TGeoHMatrix global[kNclusterMax];
450 for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
452 for (Int_t ii=0; ii<kNfoundMax; ii++) {
460 for (Int_t ii=0; ii<kNclusterMax; ii++)
463 AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
464 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(track->GetSeedIndex());
465 //if ( t->GetTOFsignal()>0. ) continue;
466 if ( (t->GetStatus()&AliESDtrack::kTOFout)!=0 ) continue;
467 AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
469 // Determine a window around the track
471 trackTOFin->GetExternalParameters(x,par);
473 trackTOFin->GetExternalCovariance(cov);
475 if (cov[0]<0. || cov[2]<0.) {
476 AliWarning(Form("Very strange track (%d)! At least one of its covariance matrix diagonal elements is negative!",iseed));
483 ((5*TMath::Sqrt(TMath::Abs(cov[0])) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius);
486 (5*TMath::Sqrt(TMath::Abs(cov[2])) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
488 Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
489 if (phi<-TMath::Pi())phi+=2*TMath::Pi();
490 if (phi>=TMath::Pi())phi-=2*TMath::Pi();
493 //upper limit on window's size.
495 if(dz> dzMax) dz=dzMax;
496 if(dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
501 // find the clusters in the window of the track
503 for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
505 if (nc>=kNclusterMax) {
506 AliWarning("No more matchable clusters can be stored! Please, increase the corresponding vectors size.");
510 AliTOFcluster *c=fClusters[k];
511 if (c->GetZ() > z+dz) break;
512 if (c->IsUsed()) continue;
514 if (!c->GetStatus()) {
515 AliDebug(1,"Cluster in channel declared bad!");
516 continue; // skip bad channels as declared in OCDB
519 Double_t dph=TMath::Abs(c->GetPhi()-phi);
520 if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
521 if (TMath::Abs(dph)>dphi) continue;
523 Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
524 Double_t p[2]={yc, c->GetZ()};
525 Double_t cov2[3]= {dY*dY/12., 0., dZ*dZ/12.};
526 if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov2) > maxChi2)continue;
531 ind[0]=c->GetDetInd(0);
532 ind[1]=c->GetDetInd(1);
533 ind[2]=c->GetDetInd(2);
534 ind[3]=c->GetDetInd(3);
535 ind[4]=c->GetDetInd(4);
536 fGeom->GetVolumePath(ind,path);
537 gGeoManager->cd(path);
538 global[nc] = *gGeoManager->GetCurrentMatrix();
542 AliDebug(1,Form(" Number of matchable TOF clusters for the track number %d: %d",iseed,nc));
544 //start fine propagation
546 Int_t nStepsDone = 0;
547 for( Int_t istep=0; istep<nSteps; istep++){
549 Float_t xs=AliTOFGeometry::RinTOF()+istep*0.1;
550 Double_t ymax=xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
553 Double_t ysect=trackTOFin->GetYat(xs,skip);
556 if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
559 } else if (ysect <-ymax) {
560 if (!trackTOFin->Rotate(-AliTOFGeometry::GetAlpha())) {
565 if(!trackTOFin->PropagateTo(xs)) {
571 // store the running point (Globalrf) - fine propagation
574 trackTOFin->GetXYZ(r);
575 trackPos[0][istep]= (Float_t) r[0];
576 trackPos[1][istep]= (Float_t) r[1];
577 trackPos[2][istep]= (Float_t) r[2];
578 trackPos[3][istep]= trackTOFin->GetIntegratedLength();
583 Bool_t accept = kFALSE;
584 Bool_t isInside = kFALSE;
585 for (Int_t istep=0; istep<nStepsDone; istep++) {
587 if (nfound>=kNfoundMax) {
588 AliWarning("No more track positions can be stored! Please, increase the corresponding vectors size.");
592 Float_t ctrackPos[3];
593 ctrackPos[0] = trackPos[0][istep];
594 ctrackPos[1] = trackPos[1][istep];
595 ctrackPos[2] = trackPos[2][istep];
597 //now see whether the track matches any of the TOF clusters
601 for (Int_t i=0; i<nc; i++) {
602 isInside = fGeom->IsInsideThePad(global[i],ctrackPos,dist3d);
605 Float_t yLoc = dist3d[1];
606 Float_t rLoc = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[2]*dist3d[2]);
607 accept = (TMath::Abs(yLoc)<padDepth*0.5 && rLoc<dCut);
608 AliDebug(2," I am in the case mLastStep==kTRUE ");
615 dist[nfound] = TMath::Sqrt(dist3d[0]*dist3d[0] +
616 dist3d[1]*dist3d[1] +
617 dist3d[2]*dist3d[2]);
618 AliDebug(2,Form(" dist3dLoc[0] = %f, dist3dLoc[1] = %f, dist3dLoc[2] = %f ",
619 dist3d[0],dist3d[1],dist3d[2]));
620 distZ[nfound] = dist3d[2]; // Z distance in the RF of the
621 // hit pad closest to the
622 // reconstructed track
623 distY[nfound] = dist3d[0]; // X distance in the RF of the
624 // hit pad closest to the
625 // reconstructed track
626 // It corresponds to Y coordinate
630 AliDebug(2,Form(" dist3dLoc[0] = %f --- distY[%d] = %f",
631 dist3d[0],nfound,distY[nfound]));
632 AliDebug(2,Form(" dist3dLoc[2] = %f --- distZ[%d] = %f",
633 dist3d[2],nfound,distZ[nfound]));
636 crecL[nfound] = trackPos[3][istep];
637 index[nfound] = clind[i]; // store cluster id
638 cxpos[nfound] = AliTOFGeometry::RinTOF()+istep*0.1; //store prop.radius
640 if(accept &&!mLastStep)break;
643 } //end for on the clusters
644 if(accept &&!mLastStep)break;
645 } //end for on the steps
647 AliDebug(1,Form(" Number of track points for the track number %d: %d",iseed,nfound));
658 // now choose the cluster to be matched with the track.
663 Float_t mindist=1000.;
666 for (Int_t iclus= 0; iclus<nfound;iclus++){
667 if (dist[iclus]< mindist){
668 mindist = dist[iclus];
669 mindistZ = distZ[iclus]; // Z distance in the RF of the hit
670 // pad closest to the reconstructed
672 mindistY = distY[iclus]; // Y distance in the RF of the hit
673 // pad closest to the reconstructed
676 idclus =index[iclus];
677 recL=crecL[iclus]+corrLen*0.5;
682 AliTOFcluster *c=fClusters[idclus];
684 Float_t tiltangle = AliTOFGeometry::GetAngles(c->GetDetInd(1),c->GetDetInd(2))*TMath::DegToRad();
685 Float_t localCheck=-mindistZ;
686 localCheck/=TMath::Cos(tiltangle); // Z (tracking/ALICE RF) component of
687 // the distance between the
688 // reconstructed track and the
689 // TOF closest cluster
691 AliDebug(2, Form("%7d %7d %10d %10d %10d %10d %7d",
694 TMath::Abs(trackTOFin->GetLabel()),
695 c->GetLabel(0), c->GetLabel(1), c->GetLabel(2),
700 // Track length correction for matching Step 2
703 Float_t rc=TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
704 Float_t rt=TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
705 +trackPos[1][70]*trackPos[1][70]
706 +trackPos[2][70]*trackPos[2][70]);
708 recL=trackPos[3][70]+dlt;
712 (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
714 (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
716 (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
720 AliDebug(2,Form(" track label good %5d",trackTOFin->GetLabel()));
726 AliDebug(2,Form(" track label bad %5d",trackTOFin->GetLabel()));
732 // Store quantities to be used in the TOF Calibration
733 Float_t tToT=AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3; // in ns
734 t->SetTOFsignalToT(tToT);
735 Float_t rawTime=AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()+kTimeOffset; // RAW time,in ps
736 t->SetTOFsignalRaw(rawTime);
737 t->SetTOFsignalDz(mindistZ);
738 t->SetTOFsignalDx(mindistY);
740 Float_t info[10] = {mindist,mindistY,mindistZ,
741 0.,0.,0.,0.,0.,0.,0.};
743 AliDebug(2,Form(" distance=%f; residual in the pad reference frame: dX=%f, dZ=%f", info[0],info[1],info[2]));
747 ind[0]=c->GetDetInd(0);
748 ind[1]=c->GetDetInd(1);
749 ind[2]=c->GetDetInd(2);
750 ind[3]=c->GetDetInd(3);
751 ind[4]=c->GetDetInd(4);
752 Int_t calindex = AliTOFGeometry::GetIndex(ind);
753 t->SetTOFCalChannel(calindex);
755 // keep track of the track labels in the matched cluster
757 tlab[0]=c->GetLabel(0);
758 tlab[1]=c->GetLabel(1);
759 tlab[2]=c->GetLabel(2);
760 AliDebug(2,Form(" tdc time of the matched track %6d = ",c->GetTDC()));
761 Double_t tof=AliTOFGeometry::TdcBinWidth()*c->GetTDC()+kTimeOffset; // in ps
762 AliDebug(2,Form(" tof time of the matched track: %f = ",tof));
763 Double_t tofcorr=tof;
764 if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
765 AliDebug(2,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));
766 //Set TOF time signal and pointer to the matched cluster
767 t->SetTOFsignal(tofcorr);
768 t->SetTOFcluster(idclus); // pointing to the recPoints tree
770 AliDebug(2,Form(" Setting TOF raw time: %f, z distance: %f corrected time: %f ",rawTime,mindistZ,tofcorr));
773 Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time); // in ps
774 Double_t mom=t->GetP();
775 AliDebug(2,Form(" Momentum for track %d -> %f", iseed,mom));
776 for (Int_t j=0;j<AliPID::kSPECIES;j++) {
777 Double_t mass=AliPID::ParticleMass(j);
778 time[j]+=(recL-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(mom*mom+mass*mass)/mom;
781 AliTOFtrack *trackTOFout = new AliTOFtrack(*t);
782 trackTOFout->PropagateTo(xpos);
784 // Fill the track residual histograms.
785 FillResiduals(trackTOFout,c,kFALSE);
787 t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
788 t->SetIntegratedLength(recL);
789 t->SetIntegratedTimes(time);
790 t->SetTOFLabel(tlab);
793 // Fill Reco-QA histos for Reconstruction
794 fHRecNClus->Fill(nc);
795 fHRecDist->Fill(mindist);
796 fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
797 fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
798 fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
799 fHRecSigZVsPWin->Fill(mom,dz);
801 // Fill Tree for on-the-fly offline Calibration
803 if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 ) ) {
814 for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
818 //_________________________________________________________________________
819 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
820 //--------------------------------------------------------------------
821 //This function loads the TOF clusters
822 //--------------------------------------------------------------------
824 Int_t npadX = AliTOFGeometry::NpadX();
825 Int_t npadZ = AliTOFGeometry::NpadZ();
826 Int_t nStripA = AliTOFGeometry::NStripA();
827 Int_t nStripB = AliTOFGeometry::NStripB();
828 Int_t nStripC = AliTOFGeometry::NStripC();
830 TBranch *branch=cTree->GetBranch("TOF");
832 AliError("can't get the branch with the TOF clusters !");
836 static TClonesArray dummy("AliTOFcluster",10000);
838 TClonesArray *clusters=&dummy;
839 branch->SetAddress(&clusters);
842 Int_t nc=clusters->GetEntriesFast();
843 fHDigNClus->Fill(nc);
845 AliInfo(Form("Number of clusters: %d",nc));
847 for (Int_t i=0; i<nc; i++) {
848 AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
849 //PH fClusters[i]=new AliTOFcluster(*c); fN++;
850 fClusters[i]=c; fN++;
852 // Fill Digits QA histos
854 Int_t isector = c->GetDetInd(0);
855 Int_t iplate = c->GetDetInd(1);
856 Int_t istrip = c->GetDetInd(2);
857 Int_t ipadX = c->GetDetInd(4);
858 Int_t ipadZ = c->GetDetInd(3);
860 Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
861 Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
863 Int_t stripOffset = 0;
869 stripOffset = nStripC;
872 stripOffset = nStripC+nStripB;
875 stripOffset = nStripC+nStripB+nStripA;
878 stripOffset = nStripC+nStripB+nStripA+nStripB;
881 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
884 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
885 Int_t phiindex=npadX*isector+ipadX+1;
886 fHDigClusMap->Fill(zindex,phiindex);
887 fHDigClusTime->Fill(time);
888 fHDigClusToT->Fill(tot);
895 //_________________________________________________________________________
896 void AliTOFtracker::UnloadClusters() {
897 //--------------------------------------------------------------------
898 //This function unloads TOF clusters
899 //--------------------------------------------------------------------
900 for (Int_t i=0; i<fN; i++) {
901 //PH delete fClusters[i];
907 //_________________________________________________________________________
908 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
909 //--------------------------------------------------------------------
910 // This function returns the index of the nearest cluster
911 //--------------------------------------------------------------------
913 if (z <= fClusters[0]->GetZ()) return 0;
914 if (z > fClusters[fN-1]->GetZ()) return fN;
915 Int_t b=0, e=fN-1, m=(b+e)/2;
916 for (; b<e; m=(b+e)/2) {
917 if (z > fClusters[m]->GetZ()) b=m+1;
923 //_________________________________________________________________________
924 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
926 // Get track space point with index i
927 // Coordinates are in the global system
928 AliTOFcluster *cl = fClusters[index];
930 xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
931 xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
933 Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
934 Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
935 Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
936 Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
937 Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
938 Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
940 cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
941 cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
942 cov[2] = -cosphi*sinth*costh*sigmaz2;
943 cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
944 cov[4] = -sinphi*sinth*costh*sigmaz2;
945 cov[5] = costh*costh*sigmaz2;
946 p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
948 // Detector numbering scheme
949 Int_t nSector = AliTOFGeometry::NSectors();
950 Int_t nPlate = AliTOFGeometry::NPlates();
951 Int_t nStripA = AliTOFGeometry::NStripA();
952 Int_t nStripB = AliTOFGeometry::NStripB();
953 Int_t nStripC = AliTOFGeometry::NStripC();
955 Int_t isector = cl->GetDetInd(0);
956 if (isector >= nSector)
957 AliError(Form("Wrong sector number in TOF (%d) !",isector));
958 Int_t iplate = cl->GetDetInd(1);
959 if (iplate >= nPlate)
960 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
961 Int_t istrip = cl->GetDetInd(2);
963 Int_t stripOffset = 0;
969 stripOffset = nStripC;
972 stripOffset = nStripC+nStripB;
975 stripOffset = nStripC+nStripB+nStripA;
978 stripOffset = nStripC+nStripB+nStripA+nStripB;
981 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
985 Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
988 UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
989 p.SetVolumeID((UShort_t)volid);
992 //_________________________________________________________________________
993 void AliTOFtracker::InitCheckHists() {
995 //Init histos for Digits/Reco QA and Calibration
998 TDirectory *dir = gDirectory;
999 TFile *logFileTOF = 0;
1001 TSeqCollection *list = gROOT->GetListOfFiles();
1002 int n = list->GetEntries();
1003 Bool_t isThere=kFALSE;
1004 for(int i=0; i<n; i++) {
1005 logFileTOF = (TFile*)list->At(i);
1006 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
1012 if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
1015 fCalTree = new TTree("CalTree", "Tree for TOF calibration");
1016 fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
1017 fCalTree->Branch("ToT",&fToT,"TOFToT/F");
1018 fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
1019 fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
1020 fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
1021 fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
1024 fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);
1025 fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);
1026 fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);
1027 fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);
1030 fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
1031 fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
1032 fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
1033 fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
1034 fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
1035 fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
1041 //_________________________________________________________________________
1042 void AliTOFtracker::SaveCheckHists() {
1044 //write histos for Digits/Reco QA and Calibration
1046 TDirectory *dir = gDirectory;
1047 TFile *logFileTOF = 0;
1049 TSeqCollection *list = gROOT->GetListOfFiles();
1050 int n = list->GetEntries();
1051 Bool_t isThere=kFALSE;
1052 for(int i=0; i<n; i++) {
1053 logFileTOF = (TFile*)list->At(i);
1054 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
1061 AliError(Form("File TOFQA.root not found!! not wring histograms...."));
1065 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
1066 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
1067 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
1068 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
1069 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
1070 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
1071 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
1072 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
1073 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
1074 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
1075 fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
1076 logFileTOF->Flush();
1080 //_________________________________________________________________________
1081 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) {
1083 //dummy, for the moment
1085 if(dist<AliTOFGeometry::ZPad()*0.5){
1087 //place here the actual correction
1093 //_________________________________________________________________________
1094 Float_t AliTOFtracker::GetTimeZerofromT0(AliESDEvent *event) const {
1096 //Returns TimeZero as measured by T0 detector
1098 return event->GetT0();
1100 //_________________________________________________________________________
1101 Float_t AliTOFtracker::GetTimeZerofromTOF(AliESDEvent * /*event*/) const {
1103 //dummy, for the moment. T0 algorithm using tracks on TOF
1105 //place T0 algo here...
1109 //_________________________________________________________________________
1111 void AliTOFtracker::FillClusterArray(TObjArray* arr) const
1114 // Returns the TOF cluster array
1120 for (Int_t i=0; i<fN; ++i) arr->Add(fClusters[i]);
1123 //_________________________________________________________________________