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;
57 //extern TROOT *gROOT;
60 ClassImp(AliTOFtracker)
62 //_____________________________________________________________________________
63 AliTOFtracker::AliTOFtracker():
73 fTracks(new TClonesArray("AliTOFtrack")),
74 fSeeds(new TObjArray(100)),
75 fTOFtrackPoints(new TObjArray(10)),
95 //AliTOFtracker main Ctor
97 for (Int_t ii=0; ii<kMaxCluster; ii++) fClusters[ii]=0x0;
99 // Getting the geometry
100 fGeom = new AliTOFGeometry();
105 //_____________________________________________________________________________
106 AliTOFtracker::~AliTOFtracker() {
113 if(!(AliCDBManager::Instance()->GetCacheFlag())){
119 delete fHDigClusTime;
125 delete fHRecSigYVsPWin;
126 delete fHRecSigZVsPWin;
138 if (fTOFtrackPoints){
139 fTOFtrackPoints->Delete();
140 delete fTOFtrackPoints;
144 for (Int_t ii=0; ii<kMaxCluster; ii++)
145 if (fClusters[ii]) fClusters[ii]->Delete();
147 for(Int_t i=0; i< 20000;i++){
149 delete fClusterESD[i];
150 fClusterESD[i] = NULL;
156 //_____________________________________________________________________________
157 void AliTOFtracker::GetPidSettings(AliESDpid *esdPID) {
159 // Sets TOF resolution from RecoParams
162 esdPID->GetTOFResponse().SetTimeResolution(fkRecoParam->GetTimeResolution());
164 AliWarning("fkRecoParam not yet set; cannot set PID settings");
166 //_____________________________________________________________________________
167 Int_t AliTOFtracker::PropagateBack(AliESDEvent * const event) {
169 // Gets seeds from ESD event and Match with TOF Clusters
173 event->SetTOFcluster(1,fClusterESD);
175 if (fNTOFmatched==0) {
176 AliInfo("No TOF recPoints to be matched with reconstructed tracks");
180 // initialize RecoParam for current event
181 AliDebug(1,"Initializing params for TOF");
183 fkRecoParam = AliTOFReconstructor::GetRecoParam(); // instantiate reco param from STEER...
185 if (fkRecoParam == 0x0) {
186 AliFatal("No Reco Param found for TOF!!!");
188 //fkRecoParam->Dump();
189 //if(fkRecoParam->GetApplyPbPbCuts())fkRecoParam=fkRecoParam->GetPbPbparam();
190 //fkRecoParam->PrintParameters();
192 //Initialise some counters
201 Int_t ntrk=event->GetNumberOfTracks();
205 //Load ESD tracks into a local Array of ESD Seeds
206 for (Int_t i=0; i<fNseeds; i++)
207 fSeeds->AddLast(event->GetTrack(i));
209 //Prepare ESD tracks candidates for TOF Matching
212 if (fNseeds==0 || fNseedsTOF==0) {
213 event->SetTOFcluster(1,fClusterESD);
214 AliInfo("No seeds to try TOF match");
218 //First Step with Strict Matching Criterion
221 //Second Step with Looser Matching Criterion
224 //Third Step without kTOFout flag (just to update clusters)
227 event->SetTOFcluster(fNTOFmatched,fClusterESD);
229 if (fNTOFmatched==0) {
230 AliInfo("No TOF recPoints to be matched with reconstructed tracks");
234 AliInfo(Form("Number of matched tracks = %d (good = %d, bad = %d)",fnmatch,fngoodmatch,fnbadmatch));
236 //Update the matched ESD tracks
238 for (Int_t i=0; i<ntrk; i++) {
239 AliESDtrack *t=event->GetTrack(i);
240 AliESDtrack *seed =(AliESDtrack*)fSeeds->At(i);
242 if ( (seed->GetStatus()&AliESDtrack::kTOFin)!=0 ) {
243 t->SetStatus(AliESDtrack::kTOFin);
244 //if(seed->GetTOFsignal()>0){
245 if ( (seed->GetStatus()&AliESDtrack::kTOFout)!=0 ) {
246 t->SetStatus(AliESDtrack::kTOFout);
247 t->SetTOFsignal(seed->GetTOFsignal());
248 t->SetTOFcluster(seed->GetTOFcluster());
249 t->SetTOFsignalToT(seed->GetTOFsignalToT());
250 t->SetTOFsignalRaw(seed->GetTOFsignalRaw());
251 t->SetTOFsignalDz(seed->GetTOFsignalDz());
252 t->SetTOFsignalDx(seed->GetTOFsignalDx());
253 t->SetTOFDeltaBC(seed->GetTOFDeltaBC());
254 t->SetTOFL0L1(seed->GetTOFL0L1());
255 t->SetTOFCalChannel(seed->GetTOFCalChannel());
256 Int_t tlab[3]; seed->GetTOFLabel(tlab);
257 t->SetTOFLabel(tlab);
259 Double_t alphaA = (Double_t)t->GetAlpha();
260 Double_t xA = (Double_t)t->GetX();
261 Double_t yA = (Double_t)t->GetY();
262 Double_t zA = (Double_t)t->GetZ();
263 Double_t p1A = (Double_t)t->GetSnp();
264 Double_t p2A = (Double_t)t->GetTgl();
265 Double_t p3A = (Double_t)t->GetSigned1Pt();
266 const Double_t *covA = (Double_t*)t->GetCovariance();
268 // Make attention, please:
269 // AliESDtrack::fTOFInfo array does not be stored in the AliESDs.root file
270 // it is there only for a check during the reconstruction step.
271 Float_t info[10]; seed->GetTOFInfo(info);
273 AliDebug(3,Form(" distance=%f; residual in the pad reference frame: dX=%f, dZ=%f", info[0],info[1],info[2]));
276 // by calling the AliESDtrack::UpdateTrackParams,
277 // the current track parameters are changed
278 // and it could cause refit problems.
279 // We need to update only the following track parameters:
280 // the track length and expected times.
281 // Removed AliESDtrack::UpdateTrackParams call
282 // Called AliESDtrack::SetIntegratedTimes(...) and
283 // AliESDtrack::SetIntegratedLength() routines.
285 AliTOFtrack *track = new AliTOFtrack(*seed);
286 t->UpdateTrackParams(track,AliESDtrack::kTOFout); // to be checked - AdC
288 Double_t time[AliPID::kSPECIESC]; t->GetIntegratedTimes(time);
291 Double_t time[AliPID::kSPECIESC]; seed->GetIntegratedTimes(time);
292 t->SetIntegratedTimes(time);
294 Double_t length = seed->GetIntegratedLength();
295 t->SetIntegratedLength(length);
297 Double_t alphaB = (Double_t)t->GetAlpha();
298 Double_t xB = (Double_t)t->GetX();
299 Double_t yB = (Double_t)t->GetY();
300 Double_t zB = (Double_t)t->GetZ();
301 Double_t p1B = (Double_t)t->GetSnp();
302 Double_t p2B = (Double_t)t->GetTgl();
303 Double_t p3B = (Double_t)t->GetSigned1Pt();
304 const Double_t *covB = (Double_t*)t->GetCovariance();
305 AliDebug(3,"Track params -now(before)-:");
306 AliDebug(3,Form(" X: %f(%f), Y: %f(%f), Z: %f(%f) --- alpha: %f(%f)",
311 AliDebug(3,Form(" p1: %f(%f), p2: %f(%f), p3: %f(%f)",
315 AliDebug(3,Form(" cov1: %f(%f), cov2: %f(%f), cov3: %f(%f)"
316 " cov4: %f(%f), cov5: %f(%f), cov6: %f(%f)"
317 " cov7: %f(%f), cov8: %f(%f), cov9: %f(%f)"
318 " cov10: %f(%f), cov11: %f(%f), cov12: %f(%f)"
319 " cov13: %f(%f), cov14: %f(%f), cov15: %f(%f)",
336 AliDebug(2,Form(" TOF params: %6d %f %f %f %f %f %6d %3d %f",
338 t->GetTOFsignalRaw(),
340 t->GetTOFsignalToT(),
343 t->GetTOFCalChannel(),
345 t->GetIntegratedLength()));
346 AliDebug(2,Form(" %f %f %f %f %f %f %f %f %f",
347 time[0], time[1], time[2], time[3], time[4], time[5], time[6], time[7], time[8]));
353 Int_t *matchmap = new Int_t[fNTOFmatched];
354 event->SetTOFcluster(fNTOFmatched,fClusterESD,matchmap);
355 for (Int_t i=0; i<ntrk; i++) { // remapping after TOF matching selection
356 AliESDtrack *t=event->GetTrack(i);
357 t->ReMapTOFcluster(fNTOFmatched,matchmap);
364 // Now done in AliESDpid
365 // fPid->MakePID(event,timeZero);
373 //_________________________________________________________________________
374 void AliTOFtracker::CollectESD() {
375 //prepare the set of ESD tracks to be matched to clusters in TOF
381 TClonesArray &aTOFTrack = *fTracks;
382 for (Int_t i=0; i<fNseeds; i++) {
384 AliESDtrack *t =(AliESDtrack*)fSeeds->At(i);
385 if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
387 AliTOFtrack *track = new AliTOFtrack(*t); // New
388 Float_t x = (Float_t)track->GetX(); //New
391 if ( ( (t->GetStatus()&AliESDtrack::kTRDout)!=0 ) ) {
393 AliDebug(1,Form(" Before propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
395 // TRD 'good' tracks, already propagated at 371 cm
396 if( x >= AliTOFGeometry::Rmin() ) {
398 if ( track->PropagateToInnerTOF() ) {
400 AliDebug(1,Form(" TRD 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()));
415 else { // TRD 'good' tracks, propagated rho<371cm
417 if ( track->PropagateToInnerTOF() ) {
419 AliDebug(1,Form(" TRD propagated track till rho = %fcm."
420 " And then the track has been propagated till rho = %fcm.",
421 x, (Float_t)track->GetX()));
423 track->SetSeedIndex(i);
424 t->UpdateTrackParams(track,AliESDtrack::kTOFin);
425 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
429 AliDebug(1,Form(" After propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
437 else { // Propagate the rest of TPCbp
439 AliDebug(1,Form(" Before propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
441 if ( track->PropagateToInnerTOF() ) {
443 AliDebug(1,Form(" TPC propagated track till rho = %fcm."
444 " And then the track has been propagated till rho = %fcm.",
445 x, (Float_t)track->GetX()));
447 track->SetSeedIndex(i);
448 t->UpdateTrackParams(track,AliESDtrack::kTOFin);
449 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
453 AliDebug(1,Form(" After propagation till inner TOF radius, ESDtrackLength=%f, TOFtrackLength=%f",t->GetIntegratedLength(),track->GetIntegratedLength()));
459 AliInfo(Form("Number of TOF seeds = %d (kTRDout371 = %d, kTRDoutLess371 = %d, !kTRDout = %d)",fNseedsTOF,seedsTOF1,seedsTOF3,seedsTOF2));
461 // Sort according uncertainties on track position
466 //_________________________________________________________________________
467 void AliTOFtracker::MatchTracks( Int_t mLastStep){
469 // Parameters used/regulating the reconstruction
471 //static Float_t corrLen=0.;//0.75;
472 static Float_t detDepth=18.;
473 static Float_t padDepth=0.5;
475 const Float_t kSpeedOfLight= 2.99792458e-2; // speed of light [cm/ps]
476 const Float_t kTimeOffset = 0.; // time offset for tracking algorithm [ps]
478 Float_t dY=AliTOFGeometry::XPad();
479 Float_t dZ=AliTOFGeometry::ZPad();
481 Float_t sensRadius = fkRecoParam->GetSensRadius();
482 Float_t stepSize = fkRecoParam->GetStepSize();
483 Float_t scaleFact = fkRecoParam->GetWindowScaleFact();
484 Float_t dyMax=fkRecoParam->GetWindowSizeMaxY();
485 Float_t dzMax=fkRecoParam->GetWindowSizeMaxZ();
486 Float_t dCut=fkRecoParam->GetDistanceCut();
487 if (dCut==3. && fNseedsTOF<=10) {
489 AliInfo(Form("Matching window=%f, since low multiplicity event (fNseedsTOF=%d)",
496 Double_t maxChi2=fkRecoParam->GetMaxChi2TRD();
497 Bool_t timeWalkCorr = fkRecoParam->GetTimeWalkCorr();
499 AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++");
500 AliDebug(1,Form("TOF sens radius: %f",sensRadius));
501 AliDebug(1,Form("TOF step size: %f",stepSize));
502 AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
503 AliDebug(1,Form("TOF Window max dy: %f",dyMax));
504 AliDebug(1,Form("TOF Window max dz: %f",dzMax));
505 AliDebug(1,Form("TOF distance Cut: %f",dCut));
506 AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
507 AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));
510 //Match ESD tracks to clusters in TOF
512 // Get the number of propagation steps
513 Int_t nSteps=(Int_t)(detDepth/stepSize);
514 AliDebug(1,Form(" Number of steps to be done %d",nSteps));
516 AliDebug(1,"++++++++++++++++++++++++++++++++++++++++++++++++++++++++");
518 //PH Arrays (moved outside of the loop)
519 Float_t * trackPos[4];
520 for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
521 Int_t * clind = new Int_t[fN];
524 const Int_t kNclusterMax = 1000; // related to fN value
525 TGeoHMatrix global[kNclusterMax];
528 for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
530 fTOFtrackPoints->Delete();
532 for (Int_t ii=0; ii<kNclusterMax; ii++)
534 AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
535 AliESDtrack *t =(AliESDtrack*)fSeeds->At(track->GetSeedIndex());
536 //if ( t->GetTOFsignal()>0. ) continue;
537 if ( ((t->GetStatus()&AliESDtrack::kTOFout)!=0 ) && mLastStep < 2) continue;
538 AliTOFtrack *trackTOFin = new AliTOFtrack(*track);
540 // Determine a window around the track
542 trackTOFin->GetExternalParameters(x,par);
544 trackTOFin->GetExternalCovariance(cov);
546 if (cov[0]<0. || cov[2]<0.) {
547 AliWarning(Form("Very strange track (%d)! At least one of its covariance matrix diagonal elements is negative!",iseed));
554 ((5*TMath::Sqrt(TMath::Abs(cov[0])) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius);
557 (5*TMath::Sqrt(TMath::Abs(cov[2])) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
559 Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
560 if (phi<-TMath::Pi())phi+=2*TMath::Pi();
561 if (phi>=TMath::Pi())phi-=2*TMath::Pi();
564 //upper limit on window's size.
565 if (dz> dzMax) dz=dzMax;
566 if (dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
569 // find the clusters in the window of the track
571 for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
573 if (nc>=kNclusterMax) {
574 AliWarning("No more matchable clusters can be stored! Please, increase the corresponding vectors size.");
578 AliTOFcluster *c=fClusters[k];
579 if (c->GetZ() > z+dz) break;
580 if (c->IsUsed() && mLastStep < 2) continue;
581 if (!c->GetStatus()) {
582 AliDebug(1,"Cluster in channel declared bad!");
583 continue; // skip bad channels as declared in OCDB
586 Double_t dph=TMath::Abs(c->GetPhi()-phi);
587 if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
588 if (TMath::Abs(dph)>dphi) continue;
590 Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
591 Double_t p[2]={yc, c->GetZ()};
592 Double_t cov2[3]= {dY*dY/12., 0., dZ*dZ/12.};
593 if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov2) > maxChi2)continue;
598 ind[0]=c->GetDetInd(0);
599 ind[1]=c->GetDetInd(1);
600 ind[2]=c->GetDetInd(2);
601 ind[3]=c->GetDetInd(3);
602 ind[4]=c->GetDetInd(4);
603 fGeom->GetVolumePath(ind,path);
604 gGeoManager->cd(path);
605 global[nc] = *gGeoManager->GetCurrentMatrix();
610 AliDebug(1,Form("No available clusters for the track number %d",iseed));
616 AliDebug(1,Form(" Number of available TOF clusters for the track number %d: %d",iseed,nc));
618 //start fine propagation
620 Int_t nStepsDone = 0;
621 for( Int_t istep=0; istep<nSteps; istep++){
623 // First of all, propagate the track...
624 Float_t xs = AliTOFGeometry::RinTOF()+istep*stepSize;
625 if (!(trackTOFin->PropagateTo(xs))) break;
627 // ...and then, if necessary, rotate the track
628 Double_t ymax = xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
629 Double_t ysect = trackTOFin->GetY();
631 if (!(trackTOFin->Rotate(AliTOFGeometry::GetAlpha()))) break;
632 } else if (ysect <-ymax) {
633 if (!(trackTOFin->Rotate(-AliTOFGeometry::GetAlpha()))) break;
637 AliDebug(3,Form(" current step %d (%d) - nStepsDone=%d",istep,nSteps,nStepsDone));
639 // store the running point (Globalrf) - fine propagation
642 trackTOFin->GetXYZ(r);
643 trackPos[0][nStepsDone-1]= (Float_t) r[0];
644 trackPos[1][nStepsDone-1]= (Float_t) r[1];
645 trackPos[2][nStepsDone-1]= (Float_t) r[2];
646 trackPos[3][nStepsDone-1]= trackTOFin->GetIntegratedLength();
656 Bool_t accept = kFALSE;
657 Bool_t isInside = kFALSE;
658 for (Int_t istep=0; istep<nStepsDone; istep++) {
660 Float_t ctrackPos[3];
661 ctrackPos[0] = trackPos[0][istep];
662 ctrackPos[1] = trackPos[1][istep];
663 ctrackPos[2] = trackPos[2][istep];
665 //now see whether the track matches any of the TOF clusters
669 for (Int_t i=0; i<nc; i++) {
670 isInside = fGeom->IsInsideThePad((TGeoHMatrix*)(&global[i]),ctrackPos,dist3d);
673 Float_t yLoc = dist3d[1];
674 Float_t rLoc = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[2]*dist3d[2]);
675 accept = (TMath::Abs(yLoc)<padDepth*0.5 && rLoc<dCut);
676 AliDebug(3," I am in the case mLastStep==kTRUE ");
683 fTOFtrackPoints->AddLast(new AliTOFtrackPoint(clind[i],
684 TMath::Sqrt(dist3d[0]*dist3d[0] + dist3d[1]*dist3d[1] + dist3d[2]*dist3d[2]),
685 dist3d[2],dist3d[0],dist3d[1],
686 AliTOFGeometry::RinTOF()+istep*stepSize,trackPos[3][istep]));
688 AliDebug(3,Form(" dist3dLoc[0] = %f, dist3dLoc[1] = %f, dist3dLoc[2] = %f ",dist3d[0],dist3d[1],dist3d[2]));
690 if(accept &&!mLastStep)break;
693 } //end for on the clusters
694 if(accept &&!mLastStep)break;
695 } //end for on the steps
702 if ( nStepsDone == 0 ) {
703 AliDebug(1,Form(" No track points for the track number %d",iseed));
709 AliDebug(2,Form(" Number of steps done for the track number %d: %d",iseed,nStepsDone));
715 Int_t *isClusterMatchable = NULL;
717 isClusterMatchable = new Int_t[nc];
718 for (Int_t i=0; i<nc; i++) isClusterMatchable[i] = kFALSE;
722 Bool_t accept = kFALSE;
723 Bool_t isInside = kFALSE;
724 for (Int_t istep=0; istep<nStepsDone; istep++) {
726 Bool_t gotInsideCluster = kFALSE;
727 Int_t trackInsideCluster = -1;
729 Float_t ctrackPos[3];
730 ctrackPos[0] = trackPos[0][istep];
731 ctrackPos[1] = trackPos[1][istep];
732 ctrackPos[2] = trackPos[2][istep];
734 //now see whether the track matches any of the TOF clusters
736 Float_t dist3d[3]={0.,0.,0.};
738 for (Int_t i=0; i<nc; i++) {
741 /* check whether track was inside another cluster
742 * and in case inhibit this cluster.
743 * this will allow to only go on and add track points for
744 * that cluster where the track got inside first */
745 if (gotInsideCluster && trackInsideCluster != i) {
746 AliDebug(3,Form(" A - istep=%d ~ %d %d ~ nfound=%d",istep,trackInsideCluster,i,nfound));
749 AliDebug(3,Form(" B - istep=%d ~ %d %d ~ nfound=%d",istep,trackInsideCluster,i,nfound));
751 /* check whether track is inside this cluster */
752 for (Int_t hh=0; hh<3; hh++) dist3d[hh]=0.;
753 isInside = fGeom->IsInsideThePad((TGeoHMatrix*)(&global[i]),ctrackPos,dist3d);
756 /* if track is inside this cluster set flags which will then
757 * inhibit to add track points for the other clusters */
759 gotInsideCluster = kTRUE;
760 trackInsideCluster = i;
764 Float_t yLoc = dist3d[1];
765 Float_t rLoc = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[2]*dist3d[2]);
766 accept = (TMath::Abs(yLoc)<padDepth*0.5 && rLoc<dCut);
767 AliDebug(3," I am in the case mLastStep==kTRUE ");
771 /* add point everytime that:
772 * - the track is inside the cluster
773 * - the track got inside the cluster, even when it eventually exited the cluster
774 * - the tracks is within dCut from the cluster
776 if (accept || isInside || gotInsideCluster) {
778 fTOFtrackPoints->AddLast(new AliTOFtrackPoint(clind[i],
779 TMath::Sqrt(dist3d[0]*dist3d[0] + dist3d[1]*dist3d[1] + dist3d[2]*dist3d[2]),
780 dist3d[2],dist3d[0],dist3d[1],
781 AliTOFGeometry::RinTOF()+istep*stepSize,trackPos[3][istep]));
783 AliDebug(2,Form(" dist3dLoc[0] = %f, dist3dLoc[1] = %f, dist3dLoc[2] = %f ",dist3d[0],dist3d[1],dist3d[2]));
786 AliDebug(3,Form(" C - istep=%d ~ %d %d ~ nfound=%d",istep,trackInsideCluster,i,nfound));
789 if(clind[i] < 20000 && mLastStep==2 && !isClusterMatchable[i]){ // add TOF clusters to the track
790 isClusterMatchable[i] = kTRUE;
793 Double_t mom=t->GetP();
794 AliDebug(3,Form(" Momentum for track %d -> %f", iseed,mom));
795 Double_t time[AliPID::kSPECIESC];
796 // read from old structure (the one used by TPC in reco)
797 for(Int_t isp=0;isp<AliPID::kSPECIESC;isp++){
798 time[isp] = t->GetIntegratedTimesOld(isp); // in ps
799 Double_t mass=AliPID::ParticleMass(isp);
800 Double_t momz = mom*AliPID::ParticleCharge(isp);
801 time[isp]+=(trackPos[3][istep]-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(momz*momz+mass*mass)/momz;
802 //time[isp]+=(trackPos[3][istep]-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(mom*mom+mass*mass)/mom;
805 if(!fClusterESD[clind[i]]->Update(t->GetID(),dist3d[1],dist3d[0],dist3d[2],trackPos[3][istep],time))//x,y,z -> tracking RF
806 t->AddTOFcluster(clind[i]);
810 /* do not break loop in any case
811 * if the track got inside a cluster all other clusters
813 // if(accept &&!mLastStep)break;
817 } //end for on the clusters
820 /* do not break loop in any case
821 * if the track got inside a cluster all other clusters
822 * are inhibited but we want to go on adding track points */
823 // if(accept &&!mLastStep)break;
825 } //end for on the steps
826 if(nc) delete[] isClusterMatchable;
829 AliDebug(1,Form("No track points for the track number %d",iseed));
835 AliDebug(1,Form(" Number of track points for the track number %d: %d",iseed,nfound));
837 // now choose the cluster to be matched with the track.
842 Float_t mindist=1000.;
845 Float_t mindistX=stepSize;
846 for (Int_t iclus= 0; iclus<nfound;iclus++) {
847 AliTOFtrackPoint *matchableTOFcluster = (AliTOFtrackPoint*)fTOFtrackPoints->At(iclus);
848 //if ( matchableTOFcluster->Distance()<mindist ) {
849 if ( TMath::Abs(matchableTOFcluster->DistanceX())<TMath::Abs(mindistX) &&
850 TMath::Abs(matchableTOFcluster->DistanceX())<=stepSize ) {
851 mindist = matchableTOFcluster->Distance();
852 mindistZ = matchableTOFcluster->DistanceZ(); // Z distance in the
857 mindistY = matchableTOFcluster->DistanceY(); // Y distance in the
862 mindistX = matchableTOFcluster->DistanceX(); // X distance in the
867 xpos = matchableTOFcluster->PropRadius();
868 idclus = matchableTOFcluster->Index();
869 recL = matchableTOFcluster->Length();// + corrLen*0.5;
871 AliDebug(2,Form(" %d(%d) --- %f (%f, %f, %f), step=%f -- idclus=%d --- seed=%d, trackId=%d, trackLab=%d", iclus,nfound,
872 mindist,mindistX,mindistY,mindistZ,stepSize,idclus,iseed,track->GetSeedIndex(),track->GetLabel()));
875 } // loop on found TOF track points
877 if (TMath::Abs(mindistX)>stepSize && idclus!=-1) {
878 AliInfo(Form(" %d - not matched --- but idclus=%d, trackId=%d, trackLab=%d",iseed,
879 idclus,track->GetSeedIndex(),track->GetLabel()));
884 AliDebug(1,Form("Reconstructed track %d doesn't match any TOF cluster", iseed));
890 AliDebug(1,Form(" %d - matched",iseed));
894 AliTOFcluster *c=fClusters[idclus];
896 AliDebug(3, Form("%7d %7d %10d %10d %10d %10d %7d",
899 TMath::Abs(trackTOFin->GetLabel()),
900 c->GetLabel(0), c->GetLabel(1), c->GetLabel(2),
905 // Track length correction for matching Step 2
908 Float_t rc = TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
909 Float_t rt = TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
910 +trackPos[1][70]*trackPos[1][70]
911 +trackPos[2][70]*trackPos[2][70]);
913 recL=trackPos[3][70]+dlt;
917 (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
919 (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
921 (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
925 AliDebug(2,Form(" track label good %5d",trackTOFin->GetLabel()));
931 AliDebug(2,Form(" track label bad %5d",trackTOFin->GetLabel()));
937 // Store quantities to be used in the TOF Calibration
938 Float_t tToT=AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3; // in ns
939 t->SetTOFsignalToT(tToT);
940 Float_t rawTime=AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()+kTimeOffset; // RAW time,in ps
941 t->SetTOFsignalRaw(rawTime);
942 t->SetTOFsignalDz(mindistZ);
943 t->SetTOFsignalDx(mindistY);
944 t->SetTOFDeltaBC(c->GetDeltaBC());
945 t->SetTOFL0L1(c->GetL0L1Latency());
947 Float_t info[10] = {mindist,mindistY,mindistZ,
948 0.,0.,0.,0.,0.,0.,0.};
950 AliDebug(3,Form(" distance=%f; residual in the pad reference frame: dX=%f, dZ=%f", info[0],info[1],info[2]));
954 ind[0]=c->GetDetInd(0);
955 ind[1]=c->GetDetInd(1);
956 ind[2]=c->GetDetInd(2);
957 ind[3]=c->GetDetInd(3);
958 ind[4]=c->GetDetInd(4);
959 Int_t calindex = AliTOFGeometry::GetIndex(ind);
960 t->SetTOFCalChannel(calindex);
962 // keep track of the track labels in the matched cluster
964 tlab[0]=c->GetLabel(0);
965 tlab[1]=c->GetLabel(1);
966 tlab[2]=c->GetLabel(2);
967 AliDebug(3,Form(" tdc time of the matched track %6d = ",c->GetTDC()));
968 Double_t tof=AliTOFGeometry::TdcBinWidth()*c->GetTDC()+kTimeOffset; // in ps
969 AliDebug(3,Form(" tof time of the matched track: %f = ",tof));
970 Double_t tofcorr=tof;
971 if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
972 AliDebug(3,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));
973 //Set TOF time signal and pointer to the matched cluster
974 t->SetTOFsignal(tofcorr);
975 t->SetTOFcluster(idclus); // pointing to the recPoints tree
977 AliDebug(3,Form(" Setting TOF raw time: %f, z distance: %f corrected time: %f ",rawTime,mindistZ,tofcorr));
980 Double_t time[AliPID::kSPECIESC];
981 // read from old structure (the one used by TPC in reco)
982 for(Int_t isp=0;isp<AliPID::kSPECIESC;isp++){
983 time[isp] = t->GetIntegratedTimesOld(isp); // in ps
985 Double_t mom=t->GetP();
986 AliDebug(3,Form(" Momentum for track %d -> %f", iseed,mom));
987 for (Int_t j=0;j<AliPID::kSPECIESC;j++) {
988 Double_t mass=AliPID::ParticleMass(j);
989 Double_t momz = mom*AliPID::ParticleCharge(j);
990 time[j]+=(recL-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(momz*momz+mass*mass)/momz;
991 //time[j]+=(recL-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(mom*mom+mass*mass)/mom;
994 AliTOFtrack *trackTOFout = new AliTOFtrack(*t);
995 if (!(trackTOFout->PropagateTo(xpos))) {
1000 // If necessary, rotate the track
1001 Double_t yATxposMax=xpos*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
1002 Double_t yATxpos=trackTOFout->GetY();
1003 if (yATxpos > yATxposMax) {
1004 if (!(trackTOFout->Rotate(AliTOFGeometry::GetAlpha()))) {
1008 } else if (yATxpos < -yATxposMax) {
1009 if (!(trackTOFout->Rotate(-AliTOFGeometry::GetAlpha()))) {
1015 // Fill the track residual histograms and update track only if in the first two step (0 and 1)
1017 FillResiduals(trackTOFout,c,kFALSE);
1019 t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
1021 // don't update old structure with TOF info
1022 // t->SetIntegratedLength(recL);
1023 // t->SetIntegratedTimes(time);
1024 // t->SetTOFLabel(tlab);
1026 // add tof cluster to the track also for step 2
1028 fClusterESD[idclus]->Update(t->GetID(),mindistY,mindist,mindistZ,recL,time);//x,y,z -> tracking RF
1030 t->AddTOFcluster(idclus);
1033 AliInfo("Too many TOF clusters matched with tracks (> 20000)");
1037 // Fill Reco-QA histos for Reconstruction
1038 fHRecNClus->Fill(nc);
1039 fHRecDist->Fill(mindist);
1041 fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
1043 fHRecSigYVsP->Fill(mom,-TMath::Sqrt(-cov[0]));
1045 fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
1047 fHRecSigZVsP->Fill(mom,-TMath::Sqrt(-cov[2]));
1048 fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
1049 fHRecSigZVsPWin->Fill(mom,dz);
1051 // Fill Tree for on-the-fly offline Calibration
1053 if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 ) ) {
1066 for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
1070 //_________________________________________________________________________
1071 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
1072 //--------------------------------------------------------------------
1073 //This function loads the TOF clusters
1074 //--------------------------------------------------------------------
1076 Int_t npadX = AliTOFGeometry::NpadX();
1077 Int_t npadZ = AliTOFGeometry::NpadZ();
1078 Int_t nStripA = AliTOFGeometry::NStripA();
1079 Int_t nStripB = AliTOFGeometry::NStripB();
1080 Int_t nStripC = AliTOFGeometry::NStripC();
1082 TBranch *branch=cTree->GetBranch("TOF");
1084 AliError("can't get the branch with the TOF clusters !");
1088 static TClonesArray dummy("AliTOFcluster",10000);
1090 TClonesArray *clusters=&dummy;
1091 branch->SetAddress(&clusters);
1094 Int_t nc=clusters->GetEntriesFast();
1095 fHDigNClus->Fill(nc);
1097 AliInfo(Form("Number of clusters: %d",nc));
1101 for(Int_t i=0; i< 20000;i++){
1103 delete fClusterESD[i];
1104 fClusterESD[i] = NULL;
1108 for (Int_t i=0; i<nc; i++) {
1109 AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
1110 //PH fClusters[i]=new AliTOFcluster(*c); fN++;
1111 fClusters[i]=c; fN++;
1113 // Fill Digits QA histos
1115 Int_t isector = c->GetDetInd(0);
1116 Int_t iplate = c->GetDetInd(1);
1117 Int_t istrip = c->GetDetInd(2);
1118 Int_t ipadX = c->GetDetInd(4);
1119 Int_t ipadZ = c->GetDetInd(3);
1121 Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
1122 Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
1131 Int_t calindex = AliTOFGeometry::GetIndex(ind);
1132 Int_t tofLabels[3]={c->GetLabel(0),c->GetLabel(1),c->GetLabel(2)};
1134 Int_t stripOffset = 0;
1140 stripOffset = nStripC;
1143 stripOffset = nStripC+nStripB;
1146 stripOffset = nStripC+nStripB+nStripA;
1149 stripOffset = nStripC+nStripB+nStripA+nStripB;
1152 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1155 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
1156 Int_t phiindex=npadX*isector+ipadX+1;
1157 fHDigClusMap->Fill(zindex,phiindex);
1158 fHDigClusTime->Fill(time);
1159 fHDigClusToT->Fill(tot);
1161 if(fNTOFmatched < 20000){
1162 fClusterESD[fNTOFmatched] = new AliESDTOFcluster(i,calindex,
1163 AliTOFGeometry::TdcBinWidth()*c->GetTDC()/*ps*/,
1164 AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()/*ps*/,
1165 AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3/*ns*/,
1167 c->GetDeltaBC(),c->GetL0L1Latency(),
1168 c->GetStatus(),c->GetZ(),c->GetPhi(),c->GetR());
1175 if(fNTOFmatched == 0)
1176 fClusterESD[0] = new AliESDTOFcluster();
1180 //_________________________________________________________________________
1181 void AliTOFtracker::UnloadClusters() {
1182 //--------------------------------------------------------------------
1183 //This function unloads TOF clusters
1184 //--------------------------------------------------------------------
1185 for (Int_t i=0; i<fN; i++) {
1186 //PH delete fClusters[i];
1192 //_________________________________________________________________________
1193 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
1194 //--------------------------------------------------------------------
1195 // This function returns the index of the nearest cluster
1196 //--------------------------------------------------------------------
1197 if (fN==0) return 0;
1198 if (z <= fClusters[0]->GetZ()) return 0;
1199 if (z > fClusters[fN-1]->GetZ()) return fN;
1200 Int_t b=0, e=fN-1, m=(b+e)/2;
1201 for (; b<e; m=(b+e)/2) {
1202 if (z > fClusters[m]->GetZ()) b=m+1;
1208 //_________________________________________________________________________
1209 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
1211 // Get track space point with index i
1212 // Coordinates are in the global system
1213 AliTOFcluster *cl = fClusters[index];
1215 xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
1216 xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
1217 xyz[2] = cl->GetZ();
1218 Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
1219 Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
1220 Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
1221 Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
1222 Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
1223 Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
1225 cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
1226 cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
1227 cov[2] = -cosphi*sinth*costh*sigmaz2;
1228 cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
1229 cov[4] = -sinphi*sinth*costh*sigmaz2;
1230 cov[5] = costh*costh*sigmaz2;
1231 p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
1233 // Detector numbering scheme
1234 Int_t nSector = AliTOFGeometry::NSectors();
1235 Int_t nPlate = AliTOFGeometry::NPlates();
1236 Int_t nStripA = AliTOFGeometry::NStripA();
1237 Int_t nStripB = AliTOFGeometry::NStripB();
1238 Int_t nStripC = AliTOFGeometry::NStripC();
1240 Int_t isector = cl->GetDetInd(0);
1241 if (isector >= nSector)
1242 AliError(Form("Wrong sector number in TOF (%d) !",isector));
1243 Int_t iplate = cl->GetDetInd(1);
1244 if (iplate >= nPlate)
1245 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1246 Int_t istrip = cl->GetDetInd(2);
1248 Int_t stripOffset = 0;
1254 stripOffset = nStripC;
1257 stripOffset = nStripC+nStripB;
1260 stripOffset = nStripC+nStripB+nStripA;
1263 stripOffset = nStripC+nStripB+nStripA+nStripB;
1266 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1270 Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
1273 UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
1274 p.SetVolumeID((UShort_t)volid);
1277 //_________________________________________________________________________
1278 void AliTOFtracker::InitCheckHists() {
1280 //Init histos for Digits/Reco QA and Calibration
1283 TDirectory *dir = gDirectory;
1284 TFile *logFileTOF = 0;
1286 TSeqCollection *list = gROOT->GetListOfFiles();
1287 int n = list->GetEntries();
1288 Bool_t isThere=kFALSE;
1289 for(int i=0; i<n; i++) {
1290 logFileTOF = (TFile*)list->At(i);
1291 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
1297 if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
1300 fCalTree = new TTree("CalTree", "Tree for TOF calibration");
1301 fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
1302 fCalTree->Branch("ToT",&fToT,"TOFToT/F");
1303 fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
1304 fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
1305 fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
1306 fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
1309 fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);
1310 fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);
1311 fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);
1312 fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);
1315 fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
1316 fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
1317 fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
1318 fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
1319 fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
1320 fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
1326 //_________________________________________________________________________
1327 void AliTOFtracker::SaveCheckHists() {
1329 //write histos for Digits/Reco QA and Calibration
1331 TDirectory *dir = gDirectory;
1332 TFile *logFileTOF = 0;
1334 TSeqCollection *list = gROOT->GetListOfFiles();
1335 int n = list->GetEntries();
1336 Bool_t isThere=kFALSE;
1337 for(int i=0; i<n; i++) {
1338 logFileTOF = (TFile*)list->At(i);
1339 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
1346 AliError(Form("File TOFQA.root not found!! not wring histograms...."));
1350 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
1351 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
1352 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
1353 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
1354 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
1355 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
1356 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
1357 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
1358 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
1359 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
1360 fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
1361 logFileTOF->Flush();
1365 //_________________________________________________________________________
1366 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) const {
1368 //dummy, for the moment
1370 if(dist<AliTOFGeometry::ZPad()*0.5){
1372 //place here the actual correction
1378 //_________________________________________________________________________
1380 void AliTOFtracker::FillClusterArray(TObjArray* arr) const
1383 // Returns the TOF cluster array
1389 for (Int_t i=0; i<fN; ++i) arr->Add(fClusters[i]);
1392 //_________________________________________________________________________