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 // Gettimg 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
337 TClonesArray &aTOFTrack = *fTracks;
338 for (Int_t i=0; i<fNseeds; i++) {
340 AliESDtrack *t =(AliESDtrack*)fSeeds->At(i);
341 if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
343 AliTOFtrack *track = new AliTOFtrack(*t); // New
344 Float_t x = (Float_t)track->GetX(); //New
346 // TRD 'good' tracks, already propagated at 371 cm
347 if ( ( (t->GetStatus()&AliESDtrack::kTRDout)!=0 ) &&
348 ( x >= AliTOFGeometry::Rmin() ) ) {
349 if ( track->PropagateToInnerTOF() ) {
351 AliDebug(1,Form(" TRD propagated track till rho = %fcm."
352 " And then the track has been propagated till rho = %fcm.",
353 x, (Float_t)track->GetX()));
355 track->SetSeedIndex(i);
356 t->UpdateTrackParams(track,AliESDtrack::kTOFin);
357 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
364 // Propagate the rest of TPCbp
366 if ( track->PropagateToInnerTOF() ) {
368 AliDebug(1,Form(" TPC propagated track till rho = %fcm."
369 " And then the track has been propagated till rho = %fcm.",
370 x, (Float_t)track->GetX()));
372 track->SetSeedIndex(i);
373 t->UpdateTrackParams(track,AliESDtrack::kTOFin);
374 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
382 AliInfo(Form("Number of TOF seeds = %d (Type 1 = %d, Type 2 = %d)",fNseedsTOF,seedsTOF1,seedsTOF2));
384 // Sort according uncertainties on track position
388 //_________________________________________________________________________
389 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
392 // Parameters used/regulating the reconstruction
394 static Float_t corrLen=0.75;
395 static Float_t detDepth=15.3;
396 static Float_t padDepth=0.5;
398 const Float_t kSpeedOfLight= 2.99792458e-2; // speed of light [cm/ps]
399 const Float_t kTimeOffset = 0.; // time offset for tracking algorithm [ps]
401 Float_t dY=AliTOFGeometry::XPad();
402 Float_t dZ=AliTOFGeometry::ZPad();
404 Float_t sensRadius = fkRecoParam->GetSensRadius();
405 Float_t stepSize = fkRecoParam->GetStepSize();
406 Float_t scaleFact = fkRecoParam->GetWindowScaleFact();
407 Float_t dyMax=fkRecoParam->GetWindowSizeMaxY();
408 Float_t dzMax=fkRecoParam->GetWindowSizeMaxZ();
409 Float_t dCut=fkRecoParam->GetDistanceCut();
410 Double_t maxChi2=fkRecoParam->GetMaxChi2TRD();
411 Bool_t timeWalkCorr = fkRecoParam->GetTimeWalkCorr();
413 AliDebug(1,"++++++++++++++TOF Reconstruction Parameters:++++++++++++ \n");
414 AliDebug(1,Form("TOF sens radius: %f",sensRadius));
415 AliDebug(1,Form("TOF step size: %f",stepSize));
416 AliDebug(1,Form("TOF Window scale factor: %f",scaleFact));
417 AliDebug(1,Form("TOF Window max dy: %f",dyMax));
418 AliDebug(1,Form("TOF Window max dz: %f",dzMax));
419 AliDebug(1,Form("TOF distance Cut: %f",dCut));
420 AliDebug(1,Form("TOF Max Chi2: %f",maxChi2));
421 AliDebug(1,Form("Time Walk Correction? : %d",timeWalkCorr));
424 //Match ESD tracks to clusters in TOF
426 // Get the number of propagation steps
427 Int_t nSteps=(Int_t)(detDepth/stepSize);
429 //PH Arrays (moved outside of the loop)
430 Float_t * trackPos[4];
431 for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
432 Int_t * clind = new Int_t[fN];
435 const Int_t kNclusterMax = 1000; // related to fN value
436 TGeoHMatrix global[kNclusterMax];
439 for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
441 fTOFtrackPoints->Delete();
443 for (Int_t ii=0; ii<kNclusterMax; ii++)
445 AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
446 AliESDtrack *t =(AliESDtrack*)fSeeds->At(track->GetSeedIndex());
447 //if ( t->GetTOFsignal()>0. ) continue;
448 if ( (t->GetStatus()&AliESDtrack::kTOFout)!=0 ) continue;
449 AliTOFtrack *trackTOFin = new AliTOFtrack(*track);
451 // Determine a window around the track
453 trackTOFin->GetExternalParameters(x,par);
455 trackTOFin->GetExternalCovariance(cov);
457 if (cov[0]<0. || cov[2]<0.) {
458 AliWarning(Form("Very strange track (%d)! At least one of its covariance matrix diagonal elements is negative!",iseed));
465 ((5*TMath::Sqrt(TMath::Abs(cov[0])) + 0.5*dY + 2.5*TMath::Abs(par[2]))/sensRadius);
468 (5*TMath::Sqrt(TMath::Abs(cov[2])) + 0.5*dZ + 2.5*TMath::Abs(par[3]));
470 Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
471 if (phi<-TMath::Pi())phi+=2*TMath::Pi();
472 if (phi>=TMath::Pi())phi-=2*TMath::Pi();
475 //upper limit on window's size.
476 if (dz> dzMax) dz=dzMax;
477 if (dphi*sensRadius> dyMax) dphi=dyMax/sensRadius;
480 // find the clusters in the window of the track
482 for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
484 if (nc>=kNclusterMax) {
485 AliWarning("No more matchable clusters can be stored! Please, increase the corresponding vectors size.");
489 AliTOFcluster *c=fClusters[k];
490 if (c->GetZ() > z+dz) break;
491 if (c->IsUsed()) continue;
492 if (!c->GetStatus()) {
493 AliDebug(1,"Cluster in channel declared bad!");
494 continue; // skip bad channels as declared in OCDB
497 Double_t dph=TMath::Abs(c->GetPhi()-phi);
498 if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
499 if (TMath::Abs(dph)>dphi) continue;
501 Double_t yc=(c->GetPhi() - trackTOFin->GetAlpha())*c->GetR();
502 Double_t p[2]={yc, c->GetZ()};
503 Double_t cov2[3]= {dY*dY/12., 0., dZ*dZ/12.};
504 if (trackTOFin->AliExternalTrackParam::GetPredictedChi2(p,cov2) > maxChi2)continue;
509 ind[0]=c->GetDetInd(0);
510 ind[1]=c->GetDetInd(1);
511 ind[2]=c->GetDetInd(2);
512 ind[3]=c->GetDetInd(3);
513 ind[4]=c->GetDetInd(4);
514 fGeom->GetVolumePath(ind,path);
515 gGeoManager->cd(path);
516 global[nc] = *gGeoManager->GetCurrentMatrix();
520 AliDebug(1,Form(" Number of matchable TOF clusters for the track number %d: %d",iseed,nc));
522 //start fine propagation
524 Int_t nStepsDone = 0;
525 for( Int_t istep=0; istep<nSteps; istep++){
527 Float_t xs=AliTOFGeometry::RinTOF()+istep*0.1;
528 Double_t ymax=xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
531 Double_t ysect=trackTOFin->GetYat(xs,skip);
534 if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
537 } else if (ysect <-ymax) {
538 if (!trackTOFin->Rotate(-AliTOFGeometry::GetAlpha())) {
543 if(!trackTOFin->PropagateTo(xs)) {
549 // store the running point (Globalrf) - fine propagation
552 trackTOFin->GetXYZ(r);
553 trackPos[0][istep]= (Float_t) r[0];
554 trackPos[1][istep]= (Float_t) r[1];
555 trackPos[2][istep]= (Float_t) r[2];
556 trackPos[3][istep]= trackTOFin->GetIntegratedLength();
561 Bool_t accept = kFALSE;
562 Bool_t isInside = kFALSE;
563 for (Int_t istep=0; istep<nStepsDone; istep++) {
565 Float_t ctrackPos[3];
566 ctrackPos[0] = trackPos[0][istep];
567 ctrackPos[1] = trackPos[1][istep];
568 ctrackPos[2] = trackPos[2][istep];
570 //now see whether the track matches any of the TOF clusters
574 for (Int_t i=0; i<nc; i++) {
575 isInside = fGeom->IsInsideThePad((TGeoHMatrix*)(&global[i]),ctrackPos,dist3d);
578 Float_t yLoc = dist3d[1];
579 Float_t rLoc = TMath::Sqrt(dist3d[0]*dist3d[0]+dist3d[2]*dist3d[2]);
580 accept = (TMath::Abs(yLoc)<padDepth*0.5 && rLoc<dCut);
581 AliDebug(2," I am in the case mLastStep==kTRUE ");
588 fTOFtrackPoints->AddLast(new AliTOFtrackPoint(clind[i],
589 TMath::Sqrt(dist3d[0]*dist3d[0] + dist3d[1]*dist3d[1] + dist3d[2]*dist3d[2]),
590 dist3d[2], dist3d[0],
591 AliTOFGeometry::RinTOF()+istep*0.1,trackPos[3][istep]));
593 AliDebug(2,Form(" dist3dLoc[0] = %f, dist3dLoc[1] = %f, dist3dLoc[2] = %f ",dist3d[0],dist3d[1],dist3d[2]));
595 if(accept &&!mLastStep)break;
598 } //end for on the clusters
599 if(accept &&!mLastStep)break;
600 } //end for on the steps
602 AliDebug(1,Form(" Number of track points for the track number %d: %d",iseed,nfound));
613 // now choose the cluster to be matched with the track.
618 Float_t mindist=1000.;
621 for (Int_t iclus= 0; iclus<nfound;iclus++){
622 AliTOFtrackPoint *matchableTOFcluster = (AliTOFtrackPoint*)fTOFtrackPoints->At(iclus);
623 if (matchableTOFcluster->Distance()<mindist) {
624 mindist = matchableTOFcluster->Distance();
625 mindistZ = matchableTOFcluster->DistanceZ(); // Z distance in the
630 mindistY = matchableTOFcluster->DistanceY(); // Y distance in the
635 xpos = matchableTOFcluster->PropRadius();
636 idclus = matchableTOFcluster->Index();
637 recL = matchableTOFcluster->Length() + corrLen*0.5;
639 } // loop on found TOF track points
642 AliDebug(1,Form("Reconstructed track %d doesn't match any TOF cluster", iseed));
647 AliTOFcluster *c=fClusters[idclus];
649 AliDebug(2, Form("%7d %7d %10d %10d %10d %10d %7d",
652 TMath::Abs(trackTOFin->GetLabel()),
653 c->GetLabel(0), c->GetLabel(1), c->GetLabel(2),
658 // Track length correction for matching Step 2
661 Float_t rc = TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
662 Float_t rt = TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
663 +trackPos[1][70]*trackPos[1][70]
664 +trackPos[2][70]*trackPos[2][70]);
666 recL=trackPos[3][70]+dlt;
670 (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
672 (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
674 (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
678 AliDebug(2,Form(" track label good %5d",trackTOFin->GetLabel()));
684 AliDebug(2,Form(" track label bad %5d",trackTOFin->GetLabel()));
690 // Store quantities to be used in the TOF Calibration
691 Float_t tToT=AliTOFGeometry::ToTBinWidth()*c->GetToT()*1E-3; // in ns
692 t->SetTOFsignalToT(tToT);
693 Float_t rawTime=AliTOFGeometry::TdcBinWidth()*c->GetTDCRAW()+kTimeOffset; // RAW time,in ps
694 t->SetTOFsignalRaw(rawTime);
695 t->SetTOFsignalDz(mindistZ);
696 t->SetTOFsignalDx(mindistY);
697 t->SetTOFDeltaBC(c->GetDeltaBC());
698 t->SetTOFL0L1(c->GetL0L1Latency());
700 Float_t info[10] = {mindist,mindistY,mindistZ,
701 0.,0.,0.,0.,0.,0.,0.};
703 AliDebug(2,Form(" distance=%f; residual in the pad reference frame: dX=%f, dZ=%f", info[0],info[1],info[2]));
707 ind[0]=c->GetDetInd(0);
708 ind[1]=c->GetDetInd(1);
709 ind[2]=c->GetDetInd(2);
710 ind[3]=c->GetDetInd(3);
711 ind[4]=c->GetDetInd(4);
712 Int_t calindex = AliTOFGeometry::GetIndex(ind);
713 t->SetTOFCalChannel(calindex);
715 // keep track of the track labels in the matched cluster
717 tlab[0]=c->GetLabel(0);
718 tlab[1]=c->GetLabel(1);
719 tlab[2]=c->GetLabel(2);
720 AliDebug(2,Form(" tdc time of the matched track %6d = ",c->GetTDC()));
721 Double_t tof=AliTOFGeometry::TdcBinWidth()*c->GetTDC()+kTimeOffset; // in ps
722 AliDebug(2,Form(" tof time of the matched track: %f = ",tof));
723 Double_t tofcorr=tof;
724 if(timeWalkCorr)tofcorr=CorrectTimeWalk(mindistZ,tof);
725 AliDebug(2,Form(" tof time of the matched track, after TW corr: %f = ",tofcorr));
726 //Set TOF time signal and pointer to the matched cluster
727 t->SetTOFsignal(tofcorr);
728 t->SetTOFcluster(idclus); // pointing to the recPoints tree
730 AliDebug(2,Form(" Setting TOF raw time: %f, z distance: %f corrected time: %f ",rawTime,mindistZ,tofcorr));
733 Double_t time[AliPID::kSPECIES]; t->GetIntegratedTimes(time); // in ps
734 Double_t mom=t->GetP();
735 AliDebug(2,Form(" Momentum for track %d -> %f", iseed,mom));
736 for (Int_t j=0;j<AliPID::kSPECIES;j++) {
737 Double_t mass=AliPID::ParticleMass(j);
738 time[j]+=(recL-trackPos[3][0])/kSpeedOfLight*TMath::Sqrt(mom*mom+mass*mass)/mom;
741 AliTOFtrack *trackTOFout = new AliTOFtrack(*t);
742 trackTOFout->PropagateTo(xpos);
744 // Fill the track residual histograms.
745 FillResiduals(trackTOFout,c,kFALSE);
747 t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
748 t->SetIntegratedLength(recL);
749 t->SetIntegratedTimes(time);
750 t->SetTOFLabel(tlab);
753 // Fill Reco-QA histos for Reconstruction
754 fHRecNClus->Fill(nc);
755 fHRecDist->Fill(mindist);
757 fHRecSigYVsP->Fill(mom,TMath::Sqrt(cov[0]));
759 fHRecSigYVsP->Fill(mom,-TMath::Sqrt(-cov[0]));
761 fHRecSigZVsP->Fill(mom,TMath::Sqrt(cov[2]));
763 fHRecSigZVsP->Fill(mom,-TMath::Sqrt(-cov[2]));
764 fHRecSigYVsPWin->Fill(mom,dphi*sensRadius);
765 fHRecSigZVsPWin->Fill(mom,dz);
767 // Fill Tree for on-the-fly offline Calibration
769 if ( !((t->GetStatus() & AliESDtrack::kTIME)==0 ) ) {
781 for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
785 //_________________________________________________________________________
786 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
787 //--------------------------------------------------------------------
788 //This function loads the TOF clusters
789 //--------------------------------------------------------------------
791 Int_t npadX = AliTOFGeometry::NpadX();
792 Int_t npadZ = AliTOFGeometry::NpadZ();
793 Int_t nStripA = AliTOFGeometry::NStripA();
794 Int_t nStripB = AliTOFGeometry::NStripB();
795 Int_t nStripC = AliTOFGeometry::NStripC();
797 TBranch *branch=cTree->GetBranch("TOF");
799 AliError("can't get the branch with the TOF clusters !");
803 static TClonesArray dummy("AliTOFcluster",10000);
805 TClonesArray *clusters=&dummy;
806 branch->SetAddress(&clusters);
809 Int_t nc=clusters->GetEntriesFast();
810 fHDigNClus->Fill(nc);
812 AliInfo(Form("Number of clusters: %d",nc));
814 for (Int_t i=0; i<nc; i++) {
815 AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
816 //PH fClusters[i]=new AliTOFcluster(*c); fN++;
817 fClusters[i]=c; fN++;
819 // Fill Digits QA histos
821 Int_t isector = c->GetDetInd(0);
822 Int_t iplate = c->GetDetInd(1);
823 Int_t istrip = c->GetDetInd(2);
824 Int_t ipadX = c->GetDetInd(4);
825 Int_t ipadZ = c->GetDetInd(3);
827 Float_t time =(AliTOFGeometry::TdcBinWidth()*c->GetTDC())*1E-3; // in ns
828 Float_t tot = (AliTOFGeometry::TdcBinWidth()*c->GetToT())*1E-3;//in ns
830 Int_t stripOffset = 0;
836 stripOffset = nStripC;
839 stripOffset = nStripC+nStripB;
842 stripOffset = nStripC+nStripB+nStripA;
845 stripOffset = nStripC+nStripB+nStripA+nStripB;
848 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
851 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
852 Int_t phiindex=npadX*isector+ipadX+1;
853 fHDigClusMap->Fill(zindex,phiindex);
854 fHDigClusTime->Fill(time);
855 fHDigClusToT->Fill(tot);
862 //_________________________________________________________________________
863 void AliTOFtracker::UnloadClusters() {
864 //--------------------------------------------------------------------
865 //This function unloads TOF clusters
866 //--------------------------------------------------------------------
867 for (Int_t i=0; i<fN; i++) {
868 //PH delete fClusters[i];
874 //_________________________________________________________________________
875 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
876 //--------------------------------------------------------------------
877 // This function returns the index of the nearest cluster
878 //--------------------------------------------------------------------
880 if (z <= fClusters[0]->GetZ()) return 0;
881 if (z > fClusters[fN-1]->GetZ()) return fN;
882 Int_t b=0, e=fN-1, m=(b+e)/2;
883 for (; b<e; m=(b+e)/2) {
884 if (z > fClusters[m]->GetZ()) b=m+1;
890 //_________________________________________________________________________
891 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
893 // Get track space point with index i
894 // Coordinates are in the global system
895 AliTOFcluster *cl = fClusters[index];
897 xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
898 xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
900 Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
901 Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
902 Float_t tiltangle = AliTOFGeometry::GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
903 Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
904 Float_t sigmay2 = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
905 Float_t sigmaz2 = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
907 cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
908 cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
909 cov[2] = -cosphi*sinth*costh*sigmaz2;
910 cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
911 cov[4] = -sinphi*sinth*costh*sigmaz2;
912 cov[5] = costh*costh*sigmaz2;
913 p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
915 // Detector numbering scheme
916 Int_t nSector = AliTOFGeometry::NSectors();
917 Int_t nPlate = AliTOFGeometry::NPlates();
918 Int_t nStripA = AliTOFGeometry::NStripA();
919 Int_t nStripB = AliTOFGeometry::NStripB();
920 Int_t nStripC = AliTOFGeometry::NStripC();
922 Int_t isector = cl->GetDetInd(0);
923 if (isector >= nSector)
924 AliError(Form("Wrong sector number in TOF (%d) !",isector));
925 Int_t iplate = cl->GetDetInd(1);
926 if (iplate >= nPlate)
927 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
928 Int_t istrip = cl->GetDetInd(2);
930 Int_t stripOffset = 0;
936 stripOffset = nStripC;
939 stripOffset = nStripC+nStripB;
942 stripOffset = nStripC+nStripB+nStripA;
945 stripOffset = nStripC+nStripB+nStripA+nStripB;
948 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
952 Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
955 UShort_t volid = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,idet);
956 p.SetVolumeID((UShort_t)volid);
959 //_________________________________________________________________________
960 void AliTOFtracker::InitCheckHists() {
962 //Init histos for Digits/Reco QA and Calibration
965 TDirectory *dir = gDirectory;
966 TFile *logFileTOF = 0;
968 TSeqCollection *list = gROOT->GetListOfFiles();
969 int n = list->GetEntries();
970 Bool_t isThere=kFALSE;
971 for(int i=0; i<n; i++) {
972 logFileTOF = (TFile*)list->At(i);
973 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
979 if(!isThere)logFileTOF = new TFile( "TOFQA.root","RECREATE");
982 fCalTree = new TTree("CalTree", "Tree for TOF calibration");
983 fCalTree->Branch("TOFchannelindex",&fIch,"iTOFch/I");
984 fCalTree->Branch("ToT",&fToT,"TOFToT/F");
985 fCalTree->Branch("TOFtime",&fTime,"TOFtime/F");
986 fCalTree->Branch("PionExpTime",&fExpTimePi,"PiExpTime/F");
987 fCalTree->Branch("KaonExpTime",&fExpTimeKa,"KaExpTime/F");
988 fCalTree->Branch("ProtonExpTime",&fExpTimePr,"PrExpTime/F");
991 fHDigClusMap = new TH2F("TOFDig_ClusMap", "",182,0.5,182.5,864, 0.5,864.5);
992 fHDigNClus = new TH1F("TOFDig_NClus", "",200,0.5,200.5);
993 fHDigClusTime = new TH1F("TOFDig_ClusTime", "",2000,0.,200.);
994 fHDigClusToT = new TH1F("TOFDig_ClusToT", "",500,0.,100);
997 fHRecNClus =new TH1F("TOFRec_NClusW", "",50,0.5,50.5);
998 fHRecDist=new TH1F("TOFRec_Dist", "",50,0.5,10.5);
999 fHRecSigYVsP=new TH2F("TOFDig_SigYVsP", "",40,0.,4.,100, 0.,5.);
1000 fHRecSigZVsP=new TH2F("TOFDig_SigZVsP", "",40,0.,4.,100, 0.,5.);
1001 fHRecSigYVsPWin=new TH2F("TOFDig_SigYVsPWin", "",40,0.,4.,100, 0.,50.);
1002 fHRecSigZVsPWin=new TH2F("TOFDig_SigZVsPWin", "",40,0.,4.,100, 0.,50.);
1008 //_________________________________________________________________________
1009 void AliTOFtracker::SaveCheckHists() {
1011 //write histos for Digits/Reco QA and Calibration
1013 TDirectory *dir = gDirectory;
1014 TFile *logFileTOF = 0;
1016 TSeqCollection *list = gROOT->GetListOfFiles();
1017 int n = list->GetEntries();
1018 Bool_t isThere=kFALSE;
1019 for(int i=0; i<n; i++) {
1020 logFileTOF = (TFile*)list->At(i);
1021 if (strstr(logFileTOF->GetName(), "TOFQA.root")){
1028 AliError(Form("File TOFQA.root not found!! not wring histograms...."));
1032 fHDigClusMap->Write(fHDigClusMap->GetName(), TObject::kOverwrite);
1033 fHDigNClus->Write(fHDigNClus->GetName(), TObject::kOverwrite);
1034 fHDigClusTime->Write(fHDigClusTime->GetName(), TObject::kOverwrite);
1035 fHDigClusToT->Write(fHDigClusToT->GetName(), TObject::kOverwrite);
1036 fHRecNClus->Write(fHRecNClus->GetName(), TObject::kOverwrite);
1037 fHRecDist->Write(fHRecDist->GetName(), TObject::kOverwrite);
1038 fHRecSigYVsP->Write(fHRecSigYVsP->GetName(), TObject::kOverwrite);
1039 fHRecSigZVsP->Write(fHRecSigZVsP->GetName(), TObject::kOverwrite);
1040 fHRecSigYVsPWin->Write(fHRecSigYVsPWin->GetName(), TObject::kOverwrite);
1041 fHRecSigZVsPWin->Write(fHRecSigZVsPWin->GetName(), TObject::kOverwrite);
1042 fCalTree->Write(fCalTree->GetName(),TObject::kOverwrite);
1043 logFileTOF->Flush();
1047 //_________________________________________________________________________
1048 Float_t AliTOFtracker::CorrectTimeWalk( Float_t dist, Float_t tof) const {
1050 //dummy, for the moment
1052 if(dist<AliTOFGeometry::ZPad()*0.5){
1054 //place here the actual correction
1060 //_________________________________________________________________________
1062 void AliTOFtracker::FillClusterArray(TObjArray* arr) const
1065 // Returns the TOF cluster array
1071 for (Int_t i=0; i<fN; ++i) arr->Add(fClusters[i]);
1074 //_________________________________________________________________________