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 <TObjArray.h>
35 #include <TGeoManager.h>
40 #include "AliGeomManager.h"
41 #include "AliESDtrack.h"
42 #include "AliESDEvent.h"
43 #include "AliESDpid.h"
45 #include "AliTrackPointArray.h"
46 #include "AliCDBManager.h"
48 //#include "AliTOFpidESD.h"
49 #include "AliTOFRecoParam.h"
50 #include "AliTOFReconstructor.h"
51 #include "AliTOFcluster.h"
52 #include "AliTOFGeometry.h"
53 #include "AliTOFtracker.h"
54 #include "AliTOFtrack.h"
56 extern TGeoManager *gGeoManager;
59 ClassImp(AliTOFtracker)
61 //_____________________________________________________________________________
62 AliTOFtracker::AliTOFtracker():
72 fTracks(new TClonesArray("AliTOFtrack")),
73 fSeeds(new TObjArray(15000)),
92 //AliTOFtracker main Ctor
94 // Gettimg the geometry
95 fGeom= new AliTOFGeometry();
100 //_____________________________________________________________________________
101 AliTOFtracker::~AliTOFtracker() {
108 if(!(AliCDBManager::Instance()->GetCacheFlag())){
114 delete fHDigClusTime;
120 delete fHRecSigYVsPWin;
121 delete fHRecSigZVsPWin;
135 //_____________________________________________________________________________
136 void AliTOFtracker::GetPidSettings(AliESDpid *esdPID) {
138 // Sets TOF resolution from RecoParams
141 esdPID->GetTOFResponse().SetTimeResolution(fkRecoParam->GetTimeResolution());
143 AliWarning("fkRecoParam not yet set; cannot set PID settings");
145 //_____________________________________________________________________________
146 Int_t AliTOFtracker::PropagateBack(AliESDEvent * const event) {
148 // Gets seeds from ESD event and Match with TOF Clusters
151 // initialize RecoParam for current event
152 AliDebug(1,"Initializing params for TOF");
154 fkRecoParam = AliTOFReconstructor::GetRecoParam(); // instantiate reco param from STEER...
156 if (fkRecoParam == 0x0) {
157 AliFatal("No Reco Param found for TOF!!!");
159 //fkRecoParam->Dump();
160 //if(fkRecoParam->GetApplyPbPbCuts())fkRecoParam=fkRecoParam->GetPbPbparam();
161 //fkRecoParam->PrintParameters();
163 //Initialise some counters
172 Int_t ntrk=event->GetNumberOfTracks();
176 //Load ESD tracks into a local Array of ESD Seeds
178 fSeeds = new TObjArray(fNseeds);
179 for (Int_t i=0; i<fNseeds; i++)
180 fSeeds->AddLast(event->GetTrack(i));
182 //Prepare ESD tracks candidates for TOF Matching
185 //First Step with Strict Matching Criterion
188 //Second Step with Looser Matching Criterion
191 AliInfo(Form("Number of matched tracks = %d (good = %d, bad = %d)",fnmatch,fngoodmatch,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->At(i);
199 if ( (seed->GetStatus()&AliESDtrack::kTOFin)!=0 ) {
200 t->SetStatus(AliESDtrack::kTOFin);
201 //if(seed->GetTOFsignal()>0){
202 if ( (seed->GetStatus()&AliESDtrack::kTOFout)!=0 ) {
203 t->SetStatus(AliESDtrack::kTOFout);
204 t->SetTOFsignal(seed->GetTOFsignal());
205 t->SetTOFcluster(seed->GetTOFcluster());
206 t->SetTOFsignalToT(seed->GetTOFsignalToT());
207 t->SetTOFsignalRaw(seed->GetTOFsignalRaw());
208 t->SetTOFsignalDz(seed->GetTOFsignalDz());
209 t->SetTOFsignalDx(seed->GetTOFsignalDx());
210 t->SetTOFDeltaBC(seed->GetTOFDeltaBC());
211 t->SetTOFL0L1(seed->GetTOFL0L1());
212 t->SetTOFCalChannel(seed->GetTOFCalChannel());
213 Int_t tlab[3]; seed->GetTOFLabel(tlab);
214 t->SetTOFLabel(tlab);
216 Double_t alphaA = (Double_t)t->GetAlpha();
217 Double_t xA = (Double_t)t->GetX();
218 Double_t yA = (Double_t)t->GetY();
219 Double_t zA = (Double_t)t->GetZ();
220 Double_t p1A = (Double_t)t->GetSnp();
221 Double_t p2A = (Double_t)t->GetTgl();
222 Double_t p3A = (Double_t)t->GetSigned1Pt();
223 const Double_t *covA = (Double_t*)t->GetCovariance();
225 // Make attention, please:
226 // AliESDtrack::fTOFInfo array does not be stored in the AliESDs.root file
227 // it is there only for a check during the reconstruction step.
228 Float_t info[10]; seed->GetTOFInfo(info);
230 AliDebug(3,Form(" distance=%f; residual in the pad reference frame: dX=%f, dZ=%f", info[0],info[1],info[2]));
233 // by calling the AliESDtrack::UpdateTrackParams,
234 // the current track parameters are changed
235 // and it could cause refit problems.
236 // We need to update only the following track parameters:
237 // the track length and expected times.
238 // Removed AliESDtrack::UpdateTrackParams call
239 // Called AliESDtrack::SetIntegratedTimes(...) and
240 // AliESDtrack::SetIntegratedLength() routines.
242 AliTOFtrack *track = new AliTOFtrack(*seed);
243 t->UpdateTrackParams(track,AliESDtrack::kTOFout); // to be checked - AdC
245 Double_t time[10]; t->GetIntegratedTimes(time);
248 Double_t time[10]; seed->GetIntegratedTimes(time);
249 t->SetIntegratedTimes(time);
251 Double_t length = seed->GetIntegratedLength();
252 t->SetIntegratedLength(length);
254 Double_t alphaB = (Double_t)t->GetAlpha();
255 Double_t xB = (Double_t)t->GetX();
256 Double_t yB = (Double_t)t->GetY();
257 Double_t zB = (Double_t)t->GetZ();
258 Double_t p1B = (Double_t)t->GetSnp();
259 Double_t p2B = (Double_t)t->GetTgl();
260 Double_t p3B = (Double_t)t->GetSigned1Pt();
261 const Double_t *covB = (Double_t*)t->GetCovariance();
262 AliDebug(2,"Track params -now(before)-:");
263 AliDebug(2,Form(" X: %f(%f), Y: %f(%f), Z: %f(%f) --- alpha: %f(%f)",
268 AliDebug(2,Form(" p1: %f(%f), p2: %f(%f), p3: %f(%f)",
272 AliDebug(2,Form(" cov1: %f(%f), cov2: %f(%f), cov3: %f(%f)"
273 " cov4: %f(%f), cov5: %f(%f), cov6: %f(%f)"
274 " cov7: %f(%f), cov8: %f(%f), cov9: %f(%f)"
275 " cov10: %f(%f), cov11: %f(%f), cov12: %f(%f)"
276 " cov13: %f(%f), cov14: %f(%f), cov15: %f(%f)",
293 AliDebug(3,Form(" TOF params: %6d %f %f %f %f %f %6d %3d %f %f %f %f %f %f",
295 t->GetTOFsignalRaw(),
297 t->GetTOFsignalToT(),
300 t->GetTOFCalChannel(),
302 t->GetIntegratedLength(),
303 time[0], time[1], time[2], time[3], time[4]
311 // Now done in AliESDpid
312 // fPid->MakePID(event,timeZero);
314 fSeeds->Clear(); delete fSeeds; fSeeds=0;
319 //_________________________________________________________________________
320 void AliTOFtracker::CollectESD() {
321 //prepare the set of ESD tracks to be matched to clusters in TOF
326 TClonesArray &aTOFTrack = *fTracks;
327 for (Int_t i=0; i<fNseeds; i++) {
329 AliESDtrack *t =(AliESDtrack*)fSeeds->At(i);
330 if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
332 AliTOFtrack *track = new AliTOFtrack(*t); // New
333 Float_t x = (Float_t)track->GetX(); //New
335 // TRD 'good' tracks, already propagated at 371 cm
336 if ( ( (t->GetStatus()&AliESDtrack::kTRDout)!=0 ) &&
337 ( x >= AliTOFGeometry::Rmin() ) ) {
338 if ( track->PropagateToInnerTOF() ) {
340 AliDebug(1,Form(" TRD propagated track till rho = %fcm."
341 " And then the track has been propagated till rho = %fcm.",
342 x, (Float_t)track->GetX()));
344 track->SetSeedIndex(i);
345 t->UpdateTrackParams(track,AliESDtrack::kTOFin);
346 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
353 // Propagate the rest of TPCbp
355 if ( track->PropagateToInnerTOF() ) {
357 AliDebug(1,Form(" TPC propagated track till rho = %fcm."
358 " And then the track has been propagated till rho = %fcm.",
359 x, (Float_t)track->GetX()));
361 track->SetSeedIndex(i);
362 t->UpdateTrackParams(track,AliESDtrack::kTOFin);
363 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
371 AliInfo(Form("Number of TOF seeds = %d (Type 1 = %d, Type 2 = %d)",fNseedsTOF,seedsTOF1,seedsTOF2));
373 // Sort according uncertainties on track position
377 //_________________________________________________________________________
378 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
381 // Parameters used/regulating the reconstruction
383 static Float_t corrLen=0.75;
384 static Float_t detDepth=15.3;
385 static Float_t padDepth=0.5;
387 const Float_t kSpeedOfLight= 2.99792458e-2; // speed of light [cm/ps]
388 const Float_t kTimeOffset = 32.; // time offset for tracking algorithm [ps]
390 Float_t dY=AliTOFGeometry::XPad();
391 Float_t dZ=AliTOFGeometry::ZPad();
393 Float_t sensRadius = fkRecoParam->GetSensRadius();
394 Float_t stepSize = fkRecoParam->GetStepSize();
395 Float_t scaleFact = fkRecoParam->GetWindowScaleFact();
396 Float_t dyMax=fkRecoParam->GetWindowSizeMaxY();
397 Float_t dzMax=fkRecoParam->GetWindowSizeMaxZ();
398 Float_t dCut=fkRecoParam->GetDistanceCut();
399 Double_t maxChi2=fkRecoParam->GetMaxChi2TRD();
400 Bool_t timeWalkCorr = fkRecoParam->GetTimeWalkCorr();
402 AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++ \n");
403 AliDebug(1,Form("TOF sens radius: %f",sensRadius));
404 AliDebug(1,Form("TOF step size: %f",stepSize));
405 AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
406 AliDebug(1,Form("TOF Window max dy: %f",dyMax));
407 AliDebug(1,Form("TOF Window max dz: %f",dzMax));
408 AliDebug(1,Form("TOF distance Cut: %f",dCut));
409 AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
410 AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));
413 //Match ESD tracks to clusters in TOF
415 // Get the number of propagation steps
416 Int_t nSteps=(Int_t)(detDepth/stepSize);
418 //PH Arrays (moved outside of the loop)
419 Float_t * trackPos[4];
420 for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
421 Int_t * clind = new Int_t[fN];
424 const Int_t kNfoundMax = 10000; // related to nSteps value
425 Int_t index[kNfoundMax];
426 Float_t dist[kNfoundMax];
427 Float_t distZ[kNfoundMax];
428 Float_t distY[kNfoundMax]; // delta(rhoXphi) // AdC
429 Float_t cxpos[kNfoundMax];
430 Float_t crecL[kNfoundMax];
431 const Int_t kNclusterMax = 1000; // related to fN value
432 TGeoHMatrix global[kNclusterMax];
435 for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
437 for (Int_t ii=0; ii<kNfoundMax; ii++) {
445 for (Int_t ii=0; ii<kNclusterMax; ii++)
448 AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
449 AliESDtrack *t =(AliESDtrack*)fSeeds->At(track->GetSeedIndex());
450 //if ( t->GetTOFsignal()>0. ) continue;
451 if ( (t->GetStatus()&AliESDtrack::kTOFout)!=0 ) continue;
452 AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
454 // Determine a window around the track
456 trackTOFin->GetExternalParameters(x,par);
458 trackTOFin->GetExternalCovariance(cov);
460 if (cov[0]<0. || cov[2]<0.) {
461 AliWarning(Form("Very strange track (%d)! At least one of its covariance matrix diagonal elements is negative!",iseed));
468 ((5*TMath::Sqrt(TMath::Abs(cov[0])) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius);
471 (5*TMath::Sqrt(TMath::Abs(cov[2])) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
473 Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
474 if (phi<-TMath::Pi())phi+=2*TMath::Pi();
475 if (phi>=TMath::Pi())phi-=2*TMath::Pi();
478 //upper limit on window's size.
480 if(dz> dzMax) dz=dzMax;
481 if(dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
486 // find the clusters in the window of the track
488 for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
490 if (nc>=kNclusterMax) {
491 AliWarning("No more matchable clusters can be stored! Please, increase the corresponding vectors size.");
495 AliTOFcluster *c=fClusters[k];
496 if (c->GetZ() > z+dz) break;
497 if (c->IsUsed()) continue;
499 if (!c->GetStatus()) {
500 AliDebug(1,"Cluster in channel declared bad!");
501 continue; // skip bad channels as declared in OCDB
504 Double_t dph=TMath::Abs(c->GetPhi()-phi);
505 if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
506 if (TMath::Abs(dph)>dphi) continue;
508 Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
509 Double_t p[2]={yc, c->GetZ()};
510 Double_t cov2[3]= {dY*dY/12., 0., dZ*dZ/12.};
511 if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov2) > maxChi2)continue;
516 ind[0]=c->GetDetInd(0);
517 ind[1]=c->GetDetInd(1);
518 ind[2]=c->GetDetInd(2);
519 ind[3]=c->GetDetInd(3);
520 ind[4]=c->GetDetInd(4);
521 fGeom->GetVolumePath(ind,path);
522 gGeoManager->cd(path);
523 global[nc] = *gGeoManager->GetCurrentMatrix();
527 AliDebug(1,Form(" Number of matchable TOF clusters for the track number %d: %d",iseed,nc));
529 //start fine propagation
531 Int_t nStepsDone = 0;
532 for( Int_t istep=0; istep<nSteps; istep++){
534 Float_t xs=AliTOFGeometry::RinTOF()+istep*0.1;
535 Double_t ymax=xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
538 Double_t ysect=trackTOFin->GetYat(xs,skip);
541 if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
544 } else if (ysect <-ymax) {
545 if (!trackTOFin->Rotate(-AliTOFGeometry::GetAlpha())) {
550 if(!trackTOFin->PropagateTo(xs)) {
556 // store the running point (Globalrf) - fine propagation
559 trackTOFin->GetXYZ(r);
560 trackPos[0][istep]= (Float_t) r[0];
561 trackPos[1][istep]= (Float_t) r[1];
562 trackPos[2][istep]= (Float_t) r[2];
563 trackPos[3][istep]= trackTOFin->GetIntegratedLength();
568 Bool_t accept = kFALSE;
569 Bool_t isInside = kFALSE;
570 for (Int_t istep=0; istep<nStepsDone; istep++) {
572 if (nfound>=kNfoundMax) {
573 AliWarning("No more track positions can be stored! Please, increase the corresponding vectors size.");
577 Float_t ctrackPos[3];
578 ctrackPos[0] = trackPos[0][istep];
579 ctrackPos[1] = trackPos[1][istep];
580 ctrackPos[2] = trackPos[2][istep];
582 //now see whether the track matches any of the TOF clusters
586 for (Int_t i=0; i<nc; i++) {
587 isInside = fGeom->IsInsideThePad(global[i],ctrackPos,dist3d);
590 Float_t yLoc = dist3d[1];
591 Float_t rLoc = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[2]*dist3d[2]);
592 accept = (TMath::Abs(yLoc)<padDepth*0.5 && rLoc<dCut);
593 AliDebug(2," I am in the case mLastStep==kTRUE ");
600 dist[nfound] = TMath::Sqrt(dist3d[0]*dist3d[0] +
601 dist3d[1]*dist3d[1] +
602 dist3d[2]*dist3d[2]);
603 AliDebug(2,Form(" dist3dLoc[0] = %f, dist3dLoc[1] = %f, dist3dLoc[2] = %f ",
604 dist3d[0],dist3d[1],dist3d[2]));
605 distZ[nfound] = dist3d[2]; // Z distance in the RF of the
606 // hit pad closest to the
607 // reconstructed track
608 distY[nfound] = dist3d[0]; // X distance in the RF of the
609 // hit pad closest to the
610 // reconstructed track
611 // It corresponds to Y coordinate
615 AliDebug(2,Form(" dist3dLoc[0] = %f --- distY[%d] = %f",
616 dist3d[0],nfound,distY[nfound]));
617 AliDebug(2,Form(" dist3dLoc[2] = %f --- distZ[%d] = %f",
618 dist3d[2],nfound,distZ[nfound]));
621 crecL[nfound] = trackPos[3][istep];
622 index[nfound] = clind[i]; // store cluster id
623 cxpos[nfound] = AliTOFGeometry::RinTOF()+istep*0.1; //store prop.radius
625 if(accept &&!mLastStep)break;
628 } //end for on the clusters
629 if(accept &&!mLastStep)break;
630 } //end for on the steps
632 AliDebug(1,Form(" Number of track points for the track number %d: %d",iseed,nfound));
643 // now choose the cluster to be matched with the track.
648 Float_t mindist=1000.;
651 for (Int_t iclus= 0; iclus<nfound;iclus++){
652 if (dist[iclus]< mindist){
653 mindist = dist[iclus];
654 mindistZ = distZ[iclus]; // Z distance in the RF of the hit
655 // pad closest to the reconstructed
657 mindistY = distY[iclus]; // Y distance in the RF of the hit
658 // pad closest to the reconstructed
661 idclus =index[iclus];
662 recL=crecL[iclus]+corrLen*0.5;
667 AliTOFcluster *c=fClusters[idclus];
669 Float_t tiltangle = AliTOFGeometry::GetAngles(c->GetDetInd(1),c->GetDetInd(2))*TMath::DegToRad();
670 Float_t localCheck=-mindistZ;
671 localCheck/=TMath::Cos(tiltangle); // Z (tracking/ALICE RF) component of
672 // the distance between the
673 // reconstructed track and the
674 // TOF closest cluster
676 AliDebug(2, Form("%7d %7d %10d %10d %10d %10d %7d",
679 TMath::Abs(trackTOFin->GetLabel()),
680 c->GetLabel(0), c->GetLabel(1), c->GetLabel(2),
685 // Track length correction for matching Step 2
688 Float_t rc=TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
689 Float_t rt=TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
690 +trackPos[1][70]*trackPos[1][70]
691 +trackPos[2][70]*trackPos[2][70]);
693 recL=trackPos[3][70]+dlt;
697 (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
699 (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
701 (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
705 AliDebug(2,Form(" track label good %5d",trackTOFin->GetLabel()));
711 AliDebug(2,Form(" track label bad %5d",trackTOFin->GetLabel()));
717 // Store quantities to be used in the TOF Calibration
718 Float_t tToT=AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3; // in ns
719 t->SetTOFsignalToT(tToT);
720 Float_t rawTime=AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()+kTimeOffset; // RAW time,in ps
721 t->SetTOFsignalRaw(rawTime);
722 t->SetTOFsignalDz(mindistZ);
723 t->SetTOFsignalDx(mindistY);
724 t->SetTOFDeltaBC(c->GetDeltaBC());
725 t->SetTOFL0L1(c->GetL0L1Latency());
727 Float_t info[10] = {mindist,mindistY,mindistZ,
728 0.,0.,0.,0.,0.,0.,0.};
730 AliDebug(2,Form(" distance=%f; residual in the pad reference frame: dX=%f, dZ=%f", info[0],info[1],info[2]));
734 ind[0]=c->GetDetInd(0);
735 ind[1]=c->GetDetInd(1);
736 ind[2]=c->GetDetInd(2);
737 ind[3]=c->GetDetInd(3);
738 ind[4]=c->GetDetInd(4);
739 Int_t calindex = AliTOFGeometry::GetIndex(ind);
740 t->SetTOFCalChannel(calindex);
742 // keep track of the track labels in the matched cluster
744 tlab[0]=c->GetLabel(0);
745 tlab[1]=c->GetLabel(1);
746 tlab[2]=c->GetLabel(2);
747 AliDebug(2,Form(" tdc time of the matched track %6d = ",c->GetTDC()));
748 Double_t tof=AliTOFGeometry::TdcBinWidth()*c->GetTDC()+kTimeOffset; // in ps
749 AliDebug(2,Form(" tof time of the matched track: %f = ",tof));
750 Double_t tofcorr=tof;
751 if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
752 AliDebug(2,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));
753 //Set TOF time signal and pointer to the matched cluster
754 t->SetTOFsignal(tofcorr);
755 t->SetTOFcluster(idclus); // pointing to the recPoints tree
757 AliDebug(2,Form(" Setting TOF raw time: %f, z distance: %f corrected time: %f ",rawTime,mindistZ,tofcorr));
760 Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time); // in ps
761 Double_t mom=t->GetP();
762 AliDebug(2,Form(" Momentum for track %d -> %f", iseed,mom));
763 for (Int_t j=0;j<AliPID::kSPECIES;j++) {
764 Double_t mass=AliPID::ParticleMass(j);
765 time[j]+=(recL-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(mom*mom+mass*mass)/mom;
768 AliTOFtrack *trackTOFout = new AliTOFtrack(*t);
769 trackTOFout->PropagateTo(xpos);
771 // Fill the track residual histograms.
772 FillResiduals(trackTOFout,c,kFALSE);
774 t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
775 t->SetIntegratedLength(recL);
776 t->SetIntegratedTimes(time);
777 t->SetTOFLabel(tlab);
780 // Fill Reco-QA histos for Reconstruction
781 fHRecNClus->Fill(nc);
782 fHRecDist->Fill(mindist);
783 fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
784 fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
785 fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
786 fHRecSigZVsPWin->Fill(mom,dz);
788 // Fill Tree for on-the-fly offline Calibration
790 if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 ) ) {
801 for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
805 //_________________________________________________________________________
806 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
807 //--------------------------------------------------------------------
808 //This function loads the TOF clusters
809 //--------------------------------------------------------------------
811 Int_t npadX = AliTOFGeometry::NpadX();
812 Int_t npadZ = AliTOFGeometry::NpadZ();
813 Int_t nStripA = AliTOFGeometry::NStripA();
814 Int_t nStripB = AliTOFGeometry::NStripB();
815 Int_t nStripC = AliTOFGeometry::NStripC();
817 TBranch *branch=cTree->GetBranch("TOF");
819 AliError("can't get the branch with the TOF clusters !");
823 static TClonesArray dummy("AliTOFcluster",10000);
825 TClonesArray *clusters=&dummy;
826 branch->SetAddress(&clusters);
829 Int_t nc=clusters->GetEntriesFast();
830 fHDigNClus->Fill(nc);
832 AliInfo(Form("Number of clusters: %d",nc));
834 for (Int_t i=0; i<nc; i++) {
835 AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
836 //PH fClusters[i]=new AliTOFcluster(*c); fN++;
837 fClusters[i]=c; fN++;
839 // Fill Digits QA histos
841 Int_t isector = c->GetDetInd(0);
842 Int_t iplate = c->GetDetInd(1);
843 Int_t istrip = c->GetDetInd(2);
844 Int_t ipadX = c->GetDetInd(4);
845 Int_t ipadZ = c->GetDetInd(3);
847 Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
848 Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
850 Int_t stripOffset = 0;
856 stripOffset = nStripC;
859 stripOffset = nStripC+nStripB;
862 stripOffset = nStripC+nStripB+nStripA;
865 stripOffset = nStripC+nStripB+nStripA+nStripB;
868 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
871 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
872 Int_t phiindex=npadX*isector+ipadX+1;
873 fHDigClusMap->Fill(zindex,phiindex);
874 fHDigClusTime->Fill(time);
875 fHDigClusToT->Fill(tot);
882 //_________________________________________________________________________
883 void AliTOFtracker::UnloadClusters() {
884 //--------------------------------------------------------------------
885 //This function unloads TOF clusters
886 //--------------------------------------------------------------------
887 for (Int_t i=0; i<fN; i++) {
888 //PH delete fClusters[i];
894 //_________________________________________________________________________
895 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
896 //--------------------------------------------------------------------
897 // This function returns the index of the nearest cluster
898 //--------------------------------------------------------------------
900 if (z <= fClusters[0]->GetZ()) return 0;
901 if (z > fClusters[fN-1]->GetZ()) return fN;
902 Int_t b=0, e=fN-1, m=(b+e)/2;
903 for (; b<e; m=(b+e)/2) {
904 if (z > fClusters[m]->GetZ()) b=m+1;
910 //_________________________________________________________________________
911 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
913 // Get track space point with index i
914 // Coordinates are in the global system
915 AliTOFcluster *cl = fClusters[index];
917 xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
918 xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
920 Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
921 Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
922 Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
923 Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
924 Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
925 Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
927 cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
928 cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
929 cov[2] = -cosphi*sinth*costh*sigmaz2;
930 cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
931 cov[4] = -sinphi*sinth*costh*sigmaz2;
932 cov[5] = costh*costh*sigmaz2;
933 p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
935 // Detector numbering scheme
936 Int_t nSector = AliTOFGeometry::NSectors();
937 Int_t nPlate = AliTOFGeometry::NPlates();
938 Int_t nStripA = AliTOFGeometry::NStripA();
939 Int_t nStripB = AliTOFGeometry::NStripB();
940 Int_t nStripC = AliTOFGeometry::NStripC();
942 Int_t isector = cl->GetDetInd(0);
943 if (isector >= nSector)
944 AliError(Form("Wrong sector number in TOF (%d) !",isector));
945 Int_t iplate = cl->GetDetInd(1);
946 if (iplate >= nPlate)
947 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
948 Int_t istrip = cl->GetDetInd(2);
950 Int_t stripOffset = 0;
956 stripOffset = nStripC;
959 stripOffset = nStripC+nStripB;
962 stripOffset = nStripC+nStripB+nStripA;
965 stripOffset = nStripC+nStripB+nStripA+nStripB;
968 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
972 Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
975 UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
976 p.SetVolumeID((UShort_t)volid);
979 //_________________________________________________________________________
980 void AliTOFtracker::InitCheckHists() {
982 //Init histos for Digits/Reco QA and Calibration
985 TDirectory *dir = gDirectory;
986 TFile *logFileTOF = 0;
988 TSeqCollection *list = gROOT->GetListOfFiles();
989 int n = list->GetEntries();
990 Bool_t isThere=kFALSE;
991 for(int i=0; i<n; i++) {
992 logFileTOF = (TFile*)list->At(i);
993 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
999 if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
1002 fCalTree = new TTree("CalTree", "Tree for TOF calibration");
1003 fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
1004 fCalTree->Branch("ToT",&fToT,"TOFToT/F");
1005 fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
1006 fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
1007 fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
1008 fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
1011 fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);
1012 fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);
1013 fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);
1014 fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);
1017 fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
1018 fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
1019 fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
1020 fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
1021 fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
1022 fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
1028 //_________________________________________________________________________
1029 void AliTOFtracker::SaveCheckHists() {
1031 //write histos for Digits/Reco QA and Calibration
1033 TDirectory *dir = gDirectory;
1034 TFile *logFileTOF = 0;
1036 TSeqCollection *list = gROOT->GetListOfFiles();
1037 int n = list->GetEntries();
1038 Bool_t isThere=kFALSE;
1039 for(int i=0; i<n; i++) {
1040 logFileTOF = (TFile*)list->At(i);
1041 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
1048 AliError(Form("File TOFQA.root not found!! not wring histograms...."));
1052 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
1053 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
1054 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
1055 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
1056 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
1057 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
1058 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
1059 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
1060 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
1061 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
1062 fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
1063 logFileTOF->Flush();
1067 //_________________________________________________________________________
1068 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) const {
1070 //dummy, for the moment
1072 if(dist<AliTOFGeometry::ZPad()*0.5){
1074 //place here the actual correction
1080 //_________________________________________________________________________
1082 void AliTOFtracker::FillClusterArray(TObjArray* arr) const
1085 // Returns the TOF cluster array
1091 for (Int_t i=0; i<fN; ++i) arr->Add(fClusters[i]);
1094 //_________________________________________________________________________