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;
60 ClassImp(AliTOFtracker)
62 //_____________________________________________________________________________
63 AliTOFtracker::AliTOFtracker():
73 fTracks(new TClonesArray("AliTOFtrack")),
74 fSeeds(new TObjArray(100)),
75 fTOFtrackPoints(new TObjArray(10)),
94 //AliTOFtracker main Ctor
96 for (Int_t ii=0; ii<kMaxCluster; ii++) fClusters[ii]=0x0;
98 // Getting the geometry
99 fGeom = new AliTOFGeometry();
104 //_____________________________________________________________________________
105 AliTOFtracker::~AliTOFtracker() {
112 if(!(AliCDBManager::Instance()->GetCacheFlag())){
118 delete fHDigClusTime;
124 delete fHRecSigYVsPWin;
125 delete fHRecSigZVsPWin;
137 if (fTOFtrackPoints){
138 fTOFtrackPoints->Delete();
139 delete fTOFtrackPoints;
143 for (Int_t ii=0; ii<kMaxCluster; ii++)
144 if (fClusters[ii]) fClusters[ii]->Delete();
147 //_____________________________________________________________________________
148 void AliTOFtracker::GetPidSettings(AliESDpid *esdPID) {
150 // Sets TOF resolution from RecoParams
153 esdPID->GetTOFResponse().SetTimeResolution(fkRecoParam->GetTimeResolution());
155 AliWarning("fkRecoParam not yet set; cannot set PID settings");
157 //_____________________________________________________________________________
158 Int_t AliTOFtracker::PropagateBack(AliESDEvent * const event) {
160 // Gets seeds from ESD event and Match with TOF Clusters
163 // initialize RecoParam for current event
164 AliDebug(1,"Initializing params for TOF");
166 fkRecoParam = AliTOFReconstructor::GetRecoParam(); // instantiate reco param from STEER...
168 if (fkRecoParam == 0x0) {
169 AliFatal("No Reco Param found for TOF!!!");
171 //fkRecoParam->Dump();
172 //if(fkRecoParam->GetApplyPbPbCuts())fkRecoParam=fkRecoParam->GetPbPbparam();
173 //fkRecoParam->PrintParameters();
175 //Initialise some counters
184 Int_t ntrk=event->GetNumberOfTracks();
188 //Load ESD tracks into a local Array of ESD Seeds
189 for (Int_t i=0; i<fNseeds; i++)
190 fSeeds->AddLast(event->GetTrack(i));
192 //Prepare ESD tracks candidates for TOF Matching
195 //First Step with Strict Matching Criterion
198 //Second Step with Looser Matching Criterion
201 AliInfo(Form("Number of matched tracks = %d (good = %d, bad = %d)",fnmatch,fngoodmatch,fnbadmatch));
203 //Update the matched ESD tracks
205 for (Int_t i=0; i<ntrk; i++) {
206 AliESDtrack *t=event->GetTrack(i);
207 AliESDtrack *seed =(AliESDtrack*)fSeeds->At(i);
209 if ( (seed->GetStatus()&AliESDtrack::kTOFin)!=0 ) {
210 t->SetStatus(AliESDtrack::kTOFin);
211 //if(seed->GetTOFsignal()>0){
212 if ( (seed->GetStatus()&AliESDtrack::kTOFout)!=0 ) {
213 t->SetStatus(AliESDtrack::kTOFout);
214 t->SetTOFsignal(seed->GetTOFsignal());
215 t->SetTOFcluster(seed->GetTOFcluster());
216 t->SetTOFsignalToT(seed->GetTOFsignalToT());
217 t->SetTOFsignalRaw(seed->GetTOFsignalRaw());
218 t->SetTOFsignalDz(seed->GetTOFsignalDz());
219 t->SetTOFsignalDx(seed->GetTOFsignalDx());
220 t->SetTOFDeltaBC(seed->GetTOFDeltaBC());
221 t->SetTOFL0L1(seed->GetTOFL0L1());
222 t->SetTOFCalChannel(seed->GetTOFCalChannel());
223 Int_t tlab[3]; seed->GetTOFLabel(tlab);
224 t->SetTOFLabel(tlab);
226 Double_t alphaA = (Double_t)t->GetAlpha();
227 Double_t xA = (Double_t)t->GetX();
228 Double_t yA = (Double_t)t->GetY();
229 Double_t zA = (Double_t)t->GetZ();
230 Double_t p1A = (Double_t)t->GetSnp();
231 Double_t p2A = (Double_t)t->GetTgl();
232 Double_t p3A = (Double_t)t->GetSigned1Pt();
233 const Double_t *covA = (Double_t*)t->GetCovariance();
235 // Make attention, please:
236 // AliESDtrack::fTOFInfo array does not be stored in the AliESDs.root file
237 // it is there only for a check during the reconstruction step.
238 Float_t info[10]; seed->GetTOFInfo(info);
240 AliDebug(3,Form(" distance=%f; residual in the pad reference frame: dX=%f, dZ=%f", info[0],info[1],info[2]));
243 // by calling the AliESDtrack::UpdateTrackParams,
244 // the current track parameters are changed
245 // and it could cause refit problems.
246 // We need to update only the following track parameters:
247 // the track length and expected times.
248 // Removed AliESDtrack::UpdateTrackParams call
249 // Called AliESDtrack::SetIntegratedTimes(...) and
250 // AliESDtrack::SetIntegratedLength() routines.
252 AliTOFtrack *track = new AliTOFtrack(*seed);
253 t->UpdateTrackParams(track,AliESDtrack::kTOFout); // to be checked - AdC
255 Double_t time[10]; t->GetIntegratedTimes(time);
258 Double_t time[10]; seed->GetIntegratedTimes(time);
259 t->SetIntegratedTimes(time);
261 Double_t length = seed->GetIntegratedLength();
262 t->SetIntegratedLength(length);
264 Double_t alphaB = (Double_t)t->GetAlpha();
265 Double_t xB = (Double_t)t->GetX();
266 Double_t yB = (Double_t)t->GetY();
267 Double_t zB = (Double_t)t->GetZ();
268 Double_t p1B = (Double_t)t->GetSnp();
269 Double_t p2B = (Double_t)t->GetTgl();
270 Double_t p3B = (Double_t)t->GetSigned1Pt();
271 const Double_t *covB = (Double_t*)t->GetCovariance();
272 AliDebug(2,"Track params -now(before)-:");
273 AliDebug(2,Form(" X: %f(%f), Y: %f(%f), Z: %f(%f) --- alpha: %f(%f)",
278 AliDebug(2,Form(" p1: %f(%f), p2: %f(%f), p3: %f(%f)",
282 AliDebug(2,Form(" cov1: %f(%f), cov2: %f(%f), cov3: %f(%f)"
283 " cov4: %f(%f), cov5: %f(%f), cov6: %f(%f)"
284 " cov7: %f(%f), cov8: %f(%f), cov9: %f(%f)"
285 " cov10: %f(%f), cov11: %f(%f), cov12: %f(%f)"
286 " cov13: %f(%f), cov14: %f(%f), cov15: %f(%f)",
303 AliDebug(3,Form(" TOF params: %6d %f %f %f %f %f %6d %3d %f %f %f %f %f %f",
305 t->GetTOFsignalRaw(),
307 t->GetTOFsignalToT(),
310 t->GetTOFCalChannel(),
312 t->GetIntegratedLength(),
313 time[0], time[1], time[2], time[3], time[4]
321 // Now done in AliESDpid
322 // fPid->MakePID(event,timeZero);
330 //_________________________________________________________________________
331 void AliTOFtracker::CollectESD() {
332 //prepare the set of ESD tracks to be matched to clusters in TOF
338 TClonesArray &aTOFTrack = *fTracks;
339 for (Int_t i=0; i<fNseeds; i++) {
341 AliESDtrack *t =(AliESDtrack*)fSeeds->At(i);
342 if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
344 AliTOFtrack *track = new AliTOFtrack(*t); // New
345 Float_t x = (Float_t)track->GetX(); //New
348 if ( ( (t->GetStatus()&AliESDtrack::kTRDout)!=0 ) ) {
350 AliDebug(1,Form(" Before propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
352 // TRD 'good' tracks, already propagated at 371 cm
353 if( x >= AliTOFGeometry::Rmin() ) {
355 if ( track->PropagateToInnerTOF() ) {
357 AliDebug(1,Form(" TRD 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);
367 AliDebug(1,Form(" After propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
372 else { // TRD 'good' tracks, propagated rho<371cm
374 if ( track->PropagateToInnerTOF() ) {
376 AliDebug(1,Form(" TRD propagated track till rho = %fcm."
377 " And then the track has been propagated till rho = %fcm.",
378 x, (Float_t)track->GetX()));
380 track->SetSeedIndex(i);
381 t->UpdateTrackParams(track,AliESDtrack::kTOFin);
382 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
386 AliDebug(1,Form(" After propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
394 else { // Propagate the rest of TPCbp
396 AliDebug(1,Form(" Before propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
398 if ( track->PropagateToInnerTOF() ) {
400 AliDebug(1,Form(" TPC propagated track till rho = %fcm."
401 " And then the track has been propagated till rho = %fcm.",
402 x, (Float_t)track->GetX()));
404 track->SetSeedIndex(i);
405 t->UpdateTrackParams(track,AliESDtrack::kTOFin);
406 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
410 AliDebug(1,Form(" After propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
416 AliInfo(Form("Number of TOF seeds = %d (kTRDout371 = %d, kTRDoutLess371 = %d, !kTRDout = %d)",fNseedsTOF,seedsTOF1,seedsTOF3,seedsTOF2));
418 // Sort according uncertainties on track position
423 //_________________________________________________________________________
424 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
426 // Parameters used/regulating the reconstruction
428 //static Float_t corrLen=0.;//0.75;
429 static Float_t detDepth=18.;
430 static Float_t padDepth=0.5;
432 const Float_t kSpeedOfLight= 2.99792458e-2; // speed of light [cm/ps]
433 const Float_t kTimeOffset = 0.; // time offset for tracking algorithm [ps]
435 Float_t dY=AliTOFGeometry::XPad();
436 Float_t dZ=AliTOFGeometry::ZPad();
438 Float_t sensRadius = fkRecoParam->GetSensRadius();
439 Float_t stepSize = fkRecoParam->GetStepSize();
440 Float_t scaleFact = fkRecoParam->GetWindowScaleFact();
441 Float_t dyMax=fkRecoParam->GetWindowSizeMaxY();
442 Float_t dzMax=fkRecoParam->GetWindowSizeMaxZ();
443 Float_t dCut=fkRecoParam->GetDistanceCut();
444 if (dCut==3. && fNseedsTOF<=10) {
446 AliInfo(Form("Matching window=%f, since low multiplicity event (fNseedsTOF=%d)",
449 Double_t maxChi2=fkRecoParam->GetMaxChi2TRD();
450 Bool_t timeWalkCorr = fkRecoParam->GetTimeWalkCorr();
452 AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++ \n");
453 AliDebug(1,Form("TOF sens radius: %f",sensRadius));
454 AliDebug(1,Form("TOF step size: %f",stepSize));
455 AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
456 AliDebug(1,Form("TOF Window max dy: %f",dyMax));
457 AliDebug(1,Form("TOF Window max dz: %f",dzMax));
458 AliDebug(1,Form("TOF distance Cut: %f",dCut));
459 AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
460 AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));
463 //Match ESD tracks to clusters in TOF
465 // Get the number of propagation steps
466 Int_t nSteps=(Int_t)(detDepth/stepSize);
467 AliDebug(1,Form(" Number of steps to be done %d",nSteps));
469 //PH Arrays (moved outside of the loop)
470 Float_t * trackPos[4];
471 for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
472 Int_t * clind = new Int_t[fN];
475 const Int_t kNclusterMax = 1000; // related to fN value
476 TGeoHMatrix global[kNclusterMax];
479 for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
481 fTOFtrackPoints->Delete();
483 for (Int_t ii=0; ii<kNclusterMax; ii++)
485 AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
486 AliESDtrack *t =(AliESDtrack*)fSeeds->At(track->GetSeedIndex());
487 //if ( t->GetTOFsignal()>0. ) continue;
488 if ( (t->GetStatus()&AliESDtrack::kTOFout)!=0 ) continue;
489 AliTOFtrack *trackTOFin = new AliTOFtrack(*track);
491 // Determine a window around the track
493 trackTOFin->GetExternalParameters(x,par);
495 trackTOFin->GetExternalCovariance(cov);
497 if (cov[0]<0. || cov[2]<0.) {
498 AliWarning(Form("Very strange track (%d)! At least one of its covariance matrix diagonal elements is negative!",iseed));
505 ((5*TMath::Sqrt(TMath::Abs(cov[0])) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius);
508 (5*TMath::Sqrt(TMath::Abs(cov[2])) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
510 Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
511 if (phi<-TMath::Pi())phi+=2*TMath::Pi();
512 if (phi>=TMath::Pi())phi-=2*TMath::Pi();
515 //upper limit on window's size.
516 if (dz> dzMax) dz=dzMax;
517 if (dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
520 // find the clusters in the window of the track
522 for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
524 if (nc>=kNclusterMax) {
525 AliWarning("No more matchable clusters can be stored! Please, increase the corresponding vectors size.");
529 AliTOFcluster *c=fClusters[k];
530 if (c->GetZ() > z+dz) break;
531 if (c->IsUsed()) continue;
532 if (!c->GetStatus()) {
533 AliDebug(1,"Cluster in channel declared bad!");
534 continue; // skip bad channels as declared in OCDB
537 Double_t dph=TMath::Abs(c->GetPhi()-phi);
538 if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
539 if (TMath::Abs(dph)>dphi) continue;
541 Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
542 Double_t p[2]={yc, c->GetZ()};
543 Double_t cov2[3]= {dY*dY/12., 0., dZ*dZ/12.};
544 if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov2) > maxChi2)continue;
549 ind[0]=c->GetDetInd(0);
550 ind[1]=c->GetDetInd(1);
551 ind[2]=c->GetDetInd(2);
552 ind[3]=c->GetDetInd(3);
553 ind[4]=c->GetDetInd(4);
554 fGeom->GetVolumePath(ind,path);
555 gGeoManager->cd(path);
556 global[nc] = *gGeoManager->GetCurrentMatrix();
560 AliDebug(1,Form(" Number of matchable TOF clusters for the track number %d: %d",iseed,nc));
568 //start fine propagation
570 Int_t nStepsDone = 0;
571 for( Int_t istep=0; istep<nSteps; istep++){
573 // First of all, propagate the track...
574 Float_t xs = AliTOFGeometry::RinTOF()+istep*stepSize;
575 if (!(trackTOFin->PropagateTo(xs))) break;
577 // ...and then, if necessary, rotate the track
578 Double_t ymax = xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
579 Double_t ysect = trackTOFin->GetY();
581 if (!(trackTOFin->Rotate(AliTOFGeometry::GetAlpha()))) break;
582 } else if (ysect <-ymax) {
583 if (!(trackTOFin->Rotate(-AliTOFGeometry::GetAlpha()))) break;
587 AliDebug(2,Form(" current step %d (%d) - nStepsDone=%d",istep,nSteps,nStepsDone));
589 // store the running point (Globalrf) - fine propagation
592 trackTOFin->GetXYZ(r);
593 trackPos[0][nStepsDone-1]= (Float_t) r[0];
594 trackPos[1][nStepsDone-1]= (Float_t) r[1];
595 trackPos[2][nStepsDone-1]= (Float_t) r[2];
596 trackPos[3][nStepsDone-1]= trackTOFin->GetIntegratedLength();
606 Bool_t accept = kFALSE;
607 Bool_t isInside = kFALSE;
608 for (Int_t istep=0; istep<nStepsDone; istep++) {
610 Float_t ctrackPos[3];
611 ctrackPos[0] = trackPos[0][istep];
612 ctrackPos[1] = trackPos[1][istep];
613 ctrackPos[2] = trackPos[2][istep];
615 //now see whether the track matches any of the TOF clusters
619 for (Int_t i=0; i<nc; i++) {
620 isInside = fGeom->IsInsideThePad((TGeoHMatrix*)(&global[i]),ctrackPos,dist3d);
623 Float_t yLoc = dist3d[1];
624 Float_t rLoc = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[2]*dist3d[2]);
625 accept = (TMath::Abs(yLoc)<padDepth*0.5 && rLoc<dCut);
626 AliDebug(2," I am in the case mLastStep==kTRUE ");
633 fTOFtrackPoints->AddLast(new AliTOFtrackPoint(clind[i],
634 TMath::Sqrt(dist3d[0]*dist3d[0] + dist3d[1]*dist3d[1] + dist3d[2]*dist3d[2]),
635 dist3d[2],dist3d[0],dist3d[1],
636 AliTOFGeometry::RinTOF()+istep*stepSize,trackPos[3][istep]));
638 AliDebug(2,Form(" dist3dLoc[0] = %f, dist3dLoc[1] = %f, dist3dLoc[2] = %f ",dist3d[0],dist3d[1],dist3d[2]));
640 if(accept &&!mLastStep)break;
643 } //end for on the clusters
644 if(accept &&!mLastStep)break;
645 } //end for on the steps
652 AliDebug(1,Form(" Number of steps done for the track number %d: %d",iseed,nStepsDone));
654 if ( nStepsDone == 0 ) {
665 Bool_t accept = kFALSE;
666 Bool_t isInside = kFALSE;
667 for (Int_t istep=0; istep<nStepsDone; istep++) {
669 Bool_t gotInsideCluster = kFALSE;
670 Int_t trackInsideCluster = -1;
672 Float_t ctrackPos[3];
673 ctrackPos[0] = trackPos[0][istep];
674 ctrackPos[1] = trackPos[1][istep];
675 ctrackPos[2] = trackPos[2][istep];
677 //now see whether the track matches any of the TOF clusters
679 Float_t dist3d[3]={0.,0.,0.};
681 for (Int_t i=0; i<nc; i++) {
684 /* check whether track was inside another cluster
685 * and in case inhibit this cluster.
686 * this will allow to only go on and add track points for
687 * that cluster where the track got inside first */
688 if (gotInsideCluster && trackInsideCluster != i) {
689 AliDebug(2,Form(" A - istep=%d ~ %d %d ~ nfound=%d",istep,trackInsideCluster,i,nfound));
692 AliDebug(2,Form(" B - istep=%d ~ %d %d ~ nfound=%d",istep,trackInsideCluster,i,nfound));
694 /* check whether track is inside this cluster */
695 for (Int_t hh=0; hh<3; hh++) dist3d[hh]=0.;
696 isInside = fGeom->IsInsideThePad((TGeoHMatrix*)(&global[i]),ctrackPos,dist3d);
699 /* if track is inside this cluster set flags which will then
700 * inhibit to add track points for the other clusters */
702 gotInsideCluster = kTRUE;
703 trackInsideCluster = i;
707 Float_t yLoc = dist3d[1];
708 Float_t rLoc = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[2]*dist3d[2]);
709 accept = (TMath::Abs(yLoc)<padDepth*0.5 && rLoc<dCut);
710 AliDebug(2," I am in the case mLastStep==kTRUE ");
714 /* add point everytime that:
715 * - the track is inside the cluster
716 * - the track got inside the cluster, even when it eventually exited the cluster
717 * - the tracks is within dCut from the cluster
719 if (accept || isInside || gotInsideCluster) {
721 fTOFtrackPoints->AddLast(new AliTOFtrackPoint(clind[i],
722 TMath::Sqrt(dist3d[0]*dist3d[0] + dist3d[1]*dist3d[1] + dist3d[2]*dist3d[2]),
723 dist3d[2],dist3d[0],dist3d[1],
724 AliTOFGeometry::RinTOF()+istep*stepSize,trackPos[3][istep]));
726 AliDebug(2,Form(" dist3dLoc[0] = %f, dist3dLoc[1] = %f, dist3dLoc[2] = %f ",dist3d[0],dist3d[1],dist3d[2]));
729 AliDebug(2,Form(" C - istep=%d ~ %d %d ~ nfound=%d",istep,trackInsideCluster,i,nfound));
732 /* do not break loop in any case
733 * if the track got inside a cluster all other clusters
735 // if(accept &&!mLastStep)break;
739 } //end for on the clusters
742 /* do not break loop in any case
743 * if the track got inside a cluster all other clusters
744 * are inhibited but we want to go on adding track points */
745 // if(accept &&!mLastStep)break;
747 } //end for on the steps
750 AliDebug(1,Form(" Number of track points for the track number %d: %d",iseed,nfound));
758 // now choose the cluster to be matched with the track.
763 Float_t mindist=1000.;
766 Float_t mindistX=stepSize;
767 for (Int_t iclus= 0; iclus<nfound;iclus++) {
768 AliTOFtrackPoint *matchableTOFcluster = (AliTOFtrackPoint*)fTOFtrackPoints->At(iclus);
769 //if ( matchableTOFcluster->Distance()<mindist ) {
770 if ( TMath::Abs(matchableTOFcluster->DistanceX())<TMath::Abs(mindistX) &&
771 TMath::Abs(matchableTOFcluster->DistanceX())<=stepSize ) {
772 mindist = matchableTOFcluster->Distance();
773 mindistZ = matchableTOFcluster->DistanceZ(); // Z distance in the
778 mindistY = matchableTOFcluster->DistanceY(); // Y distance in the
783 mindistX = matchableTOFcluster->DistanceX(); // X distance in the
788 xpos = matchableTOFcluster->PropRadius();
789 idclus = matchableTOFcluster->Index();
790 recL = matchableTOFcluster->Length();// + corrLen*0.5;
792 AliDebug(1,Form(" %d(%d) --- %f (%f, %f, %f), step=%f -- idclus=%d --- seed=%d, trackId=%d, trackLab=%d", iclus,nfound,
793 mindist,mindistX,mindistY,mindistZ,stepSize,idclus,iseed,track->GetSeedIndex(),track->GetLabel()));
796 } // loop on found TOF track points
798 if (TMath::Abs(mindistX)>stepSize && idclus!=-1) {
799 AliInfo(Form("--------Not matched --- but idclus=%d, trackId=%d, trackLab=%d",
800 idclus,track->GetSeedIndex(),track->GetLabel()));
805 AliDebug(1,Form("Reconstructed track %d doesn't match any TOF cluster", iseed));
811 AliDebug(1,"--------Matched");
815 AliTOFcluster *c=fClusters[idclus];
817 AliDebug(2, Form("%7d %7d %10d %10d %10d %10d %7d",
820 TMath::Abs(trackTOFin->GetLabel()),
821 c->GetLabel(0), c->GetLabel(1), c->GetLabel(2),
826 // Track length correction for matching Step 2
829 Float_t rc = TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
830 Float_t rt = TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
831 +trackPos[1][70]*trackPos[1][70]
832 +trackPos[2][70]*trackPos[2][70]);
834 recL=trackPos[3][70]+dlt;
838 (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
840 (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
842 (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
846 AliDebug(2,Form(" track label good %5d",trackTOFin->GetLabel()));
852 AliDebug(2,Form(" track label bad %5d",trackTOFin->GetLabel()));
858 // Store quantities to be used in the TOF Calibration
859 Float_t tToT=AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3; // in ns
860 t->SetTOFsignalToT(tToT);
861 Float_t rawTime=AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()+kTimeOffset; // RAW time,in ps
862 t->SetTOFsignalRaw(rawTime);
863 t->SetTOFsignalDz(mindistZ);
864 t->SetTOFsignalDx(mindistY);
865 t->SetTOFDeltaBC(c->GetDeltaBC());
866 t->SetTOFL0L1(c->GetL0L1Latency());
868 Float_t info[10] = {mindist,mindistY,mindistZ,
869 0.,0.,0.,0.,0.,0.,0.};
871 AliDebug(2,Form(" distance=%f; residual in the pad reference frame: dX=%f, dZ=%f", info[0],info[1],info[2]));
875 ind[0]=c->GetDetInd(0);
876 ind[1]=c->GetDetInd(1);
877 ind[2]=c->GetDetInd(2);
878 ind[3]=c->GetDetInd(3);
879 ind[4]=c->GetDetInd(4);
880 Int_t calindex = AliTOFGeometry::GetIndex(ind);
881 t->SetTOFCalChannel(calindex);
883 // keep track of the track labels in the matched cluster
885 tlab[0]=c->GetLabel(0);
886 tlab[1]=c->GetLabel(1);
887 tlab[2]=c->GetLabel(2);
888 AliDebug(2,Form(" tdc time of the matched track %6d = ",c->GetTDC()));
889 Double_t tof=AliTOFGeometry::TdcBinWidth()*c->GetTDC()+kTimeOffset; // in ps
890 AliDebug(2,Form(" tof time of the matched track: %f = ",tof));
891 Double_t tofcorr=tof;
892 if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
893 AliDebug(2,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));
894 //Set TOF time signal and pointer to the matched cluster
895 t->SetTOFsignal(tofcorr);
896 t->SetTOFcluster(idclus); // pointing to the recPoints tree
898 AliDebug(2,Form(" Setting TOF raw time: %f, z distance: %f corrected time: %f ",rawTime,mindistZ,tofcorr));
901 Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time); // in ps
902 Double_t mom=t->GetP();
903 AliDebug(2,Form(" Momentum for track %d -> %f", iseed,mom));
904 for (Int_t j=0;j<AliPID::kSPECIES;j++) {
905 Double_t mass=AliPID::ParticleMass(j);
906 time[j]+=(recL-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(mom*mom+mass*mass)/mom;
909 AliTOFtrack *trackTOFout = new AliTOFtrack(*t);
910 if (!(trackTOFout->PropagateTo(xpos))) {
915 // If necessary, rotate the track
916 Double_t yATxposMax=xpos*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
917 Double_t yATxpos=trackTOFout->GetY();
918 if (yATxpos > yATxposMax) {
919 if (!(trackTOFout->Rotate(AliTOFGeometry::GetAlpha()))) {
923 } else if (yATxpos < -yATxposMax) {
924 if (!(trackTOFout->Rotate(-AliTOFGeometry::GetAlpha()))) {
930 // Fill the track residual histograms.
931 FillResiduals(trackTOFout,c,kFALSE);
933 t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
934 t->SetIntegratedLength(recL);
935 t->SetIntegratedTimes(time);
936 t->SetTOFLabel(tlab);
939 // Fill Reco-QA histos for Reconstruction
940 fHRecNClus->Fill(nc);
941 fHRecDist->Fill(mindist);
943 fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
945 fHRecSigYVsP->Fill(mom,-TMath::Sqrt(-cov[0]));
947 fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
949 fHRecSigZVsP->Fill(mom,-TMath::Sqrt(-cov[2]));
950 fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
951 fHRecSigZVsPWin->Fill(mom,dz);
953 // Fill Tree for on-the-fly offline Calibration
955 if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 ) ) {
967 for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
971 //_________________________________________________________________________
972 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
973 //--------------------------------------------------------------------
974 //This function loads the TOF clusters
975 //--------------------------------------------------------------------
977 Int_t npadX = AliTOFGeometry::NpadX();
978 Int_t npadZ = AliTOFGeometry::NpadZ();
979 Int_t nStripA = AliTOFGeometry::NStripA();
980 Int_t nStripB = AliTOFGeometry::NStripB();
981 Int_t nStripC = AliTOFGeometry::NStripC();
983 TBranch *branch=cTree->GetBranch("TOF");
985 AliError("can't get the branch with the TOF clusters !");
989 static TClonesArray dummy("AliTOFcluster",10000);
991 TClonesArray *clusters=&dummy;
992 branch->SetAddress(&clusters);
995 Int_t nc=clusters->GetEntriesFast();
996 fHDigNClus->Fill(nc);
998 AliInfo(Form("Number of clusters: %d",nc));
1000 for (Int_t i=0; i<nc; i++) {
1001 AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
1002 //PH fClusters[i]=new AliTOFcluster(*c); fN++;
1003 fClusters[i]=c; fN++;
1005 // Fill Digits QA histos
1007 Int_t isector = c->GetDetInd(0);
1008 Int_t iplate = c->GetDetInd(1);
1009 Int_t istrip = c->GetDetInd(2);
1010 Int_t ipadX = c->GetDetInd(4);
1011 Int_t ipadZ = c->GetDetInd(3);
1013 Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
1014 Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
1016 Int_t stripOffset = 0;
1022 stripOffset = nStripC;
1025 stripOffset = nStripC+nStripB;
1028 stripOffset = nStripC+nStripB+nStripA;
1031 stripOffset = nStripC+nStripB+nStripA+nStripB;
1034 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1037 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
1038 Int_t phiindex=npadX*isector+ipadX+1;
1039 fHDigClusMap->Fill(zindex,phiindex);
1040 fHDigClusTime->Fill(time);
1041 fHDigClusToT->Fill(tot);
1048 //_________________________________________________________________________
1049 void AliTOFtracker::UnloadClusters() {
1050 //--------------------------------------------------------------------
1051 //This function unloads TOF clusters
1052 //--------------------------------------------------------------------
1053 for (Int_t i=0; i<fN; i++) {
1054 //PH delete fClusters[i];
1060 //_________________________________________________________________________
1061 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
1062 //--------------------------------------------------------------------
1063 // This function returns the index of the nearest cluster
1064 //--------------------------------------------------------------------
1065 if (fN==0) return 0;
1066 if (z <= fClusters[0]->GetZ()) return 0;
1067 if (z > fClusters[fN-1]->GetZ()) return fN;
1068 Int_t b=0, e=fN-1, m=(b+e)/2;
1069 for (; b<e; m=(b+e)/2) {
1070 if (z > fClusters[m]->GetZ()) b=m+1;
1076 //_________________________________________________________________________
1077 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
1079 // Get track space point with index i
1080 // Coordinates are in the global system
1081 AliTOFcluster *cl = fClusters[index];
1083 xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
1084 xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
1085 xyz[2] = cl->GetZ();
1086 Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
1087 Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
1088 Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
1089 Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
1090 Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
1091 Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
1093 cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
1094 cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
1095 cov[2] = -cosphi*sinth*costh*sigmaz2;
1096 cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
1097 cov[4] = -sinphi*sinth*costh*sigmaz2;
1098 cov[5] = costh*costh*sigmaz2;
1099 p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
1101 // Detector numbering scheme
1102 Int_t nSector = AliTOFGeometry::NSectors();
1103 Int_t nPlate = AliTOFGeometry::NPlates();
1104 Int_t nStripA = AliTOFGeometry::NStripA();
1105 Int_t nStripB = AliTOFGeometry::NStripB();
1106 Int_t nStripC = AliTOFGeometry::NStripC();
1108 Int_t isector = cl->GetDetInd(0);
1109 if (isector >= nSector)
1110 AliError(Form("Wrong sector number in TOF (%d) !",isector));
1111 Int_t iplate = cl->GetDetInd(1);
1112 if (iplate >= nPlate)
1113 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1114 Int_t istrip = cl->GetDetInd(2);
1116 Int_t stripOffset = 0;
1122 stripOffset = nStripC;
1125 stripOffset = nStripC+nStripB;
1128 stripOffset = nStripC+nStripB+nStripA;
1131 stripOffset = nStripC+nStripB+nStripA+nStripB;
1134 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1138 Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
1141 UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
1142 p.SetVolumeID((UShort_t)volid);
1145 //_________________________________________________________________________
1146 void AliTOFtracker::InitCheckHists() {
1148 //Init histos for Digits/Reco QA and Calibration
1151 TDirectory *dir = gDirectory;
1152 TFile *logFileTOF = 0;
1154 TSeqCollection *list = gROOT->GetListOfFiles();
1155 int n = list->GetEntries();
1156 Bool_t isThere=kFALSE;
1157 for(int i=0; i<n; i++) {
1158 logFileTOF = (TFile*)list->At(i);
1159 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
1165 if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
1168 fCalTree = new TTree("CalTree", "Tree for TOF calibration");
1169 fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
1170 fCalTree->Branch("ToT",&fToT,"TOFToT/F");
1171 fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
1172 fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
1173 fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
1174 fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
1177 fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);
1178 fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);
1179 fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);
1180 fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);
1183 fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
1184 fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
1185 fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
1186 fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
1187 fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
1188 fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
1194 //_________________________________________________________________________
1195 void AliTOFtracker::SaveCheckHists() {
1197 //write histos for Digits/Reco QA and Calibration
1199 TDirectory *dir = gDirectory;
1200 TFile *logFileTOF = 0;
1202 TSeqCollection *list = gROOT->GetListOfFiles();
1203 int n = list->GetEntries();
1204 Bool_t isThere=kFALSE;
1205 for(int i=0; i<n; i++) {
1206 logFileTOF = (TFile*)list->At(i);
1207 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
1214 AliError(Form("File TOFQA.root not found!! not wring histograms...."));
1218 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
1219 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
1220 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
1221 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
1222 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
1223 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
1224 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
1225 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
1226 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
1227 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
1228 fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
1229 logFileTOF->Flush();
1233 //_________________________________________________________________________
1234 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) const {
1236 //dummy, for the moment
1238 if(dist<AliTOFGeometry::ZPad()*0.5){
1240 //place here the actual correction
1246 //_________________________________________________________________________
1248 void AliTOFtracker::FillClusterArray(TObjArray* arr) const
1251 // Returns the TOF cluster array
1257 for (Int_t i=0; i<fN; ++i) arr->Add(fClusters[i]);
1260 //_________________________________________________________________________