]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFtrackerMI.cxx
Shadowed declarations: suppression
[u/mrichter/AliRoot.git] / TOF / AliTOFtrackerMI.cxx
CommitLineData
d88fbf15 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
0e46b9ae 15
16//-----------------------------------------------------------------//
17// //
18// AliTOFtracker Class //
19// Task: Perform association of the ESD tracks to TOF Clusters //
20// and Update ESD track with associated TOF Cluster parameters //
21// //
22//-----------------------------------------------------------------//
d88fbf15 23
24#include <Rtypes.h>
0e46b9ae 25
d88fbf15 26#include "TClonesArray.h"
5c7c93fa 27//#include "TGeoManager.h"
0e46b9ae 28#include "TTree.h"
d88fbf15 29#include "TTreeStream.h"
0e46b9ae 30
5c7c93fa 31//#include "AliRun.h"
af885e0f 32#include "AliESDEvent.h"
0e46b9ae 33#include "AliESDtrack.h"
34
e0ddb533 35#include "AliTOFRecoParam.h"
a533f541 36#include "AliTOFcalib.h"
0e46b9ae 37#include "AliTOFcluster.h"
96f01799 38#include "AliTOFGeometry.h"
0e46b9ae 39#include "AliTOFtrackerMI.h"
40#include "AliTOFtrack.h"
5c7c93fa 41#include "AliTOFpidESD.h"
42
43class TGeoManager;
0e46b9ae 44
45extern TGeoManager *gGeoManager;
d88fbf15 46
47ClassImp(AliTOFtrackerMI)
48
49//_____________________________________________________________________________
e0ddb533 50AliTOFtrackerMI::AliTOFtrackerMI():
51 fRecoParam(0x0),
52 fGeom(0x0),
53 fPid(0x0),
58d8d9a3 54 fN(0),
55 fNseeds(0),
56 fNseedsTOF(0),
57 fngoodmatch(0),
58 fnbadmatch(0),
59 fnunmatch(0),
60 fnmatch(0),
8dacd1bb 61 fR(379.),
58d8d9a3 62 fTOFHeigth(15.3),
63 fdCut(3.),
64 fDx(1.5),
65 fDy(0),
66 fDz(0),
b42a5ea1 67 fTracks(new TClonesArray("AliTOFtrack")),
68 fSeeds(new TClonesArray("AliESDtrack")),
58d8d9a3 69 fDebugStreamer(0x0)
70 {
d88fbf15 71 //AliTOFtrackerMI main Ctor
72
e0ddb533 73 fRecoParam=new AliTOFRecoParam();
96f01799 74 fGeom=new AliTOFGeometry();
e0ddb533 75 Double_t parPID[2];
76 parPID[0]=fRecoParam->GetTimeResolution();
77 parPID[1]=fRecoParam->GetTimeNSigma();
78 fPid=new AliTOFpidESD(parPID);
96f01799 79 fDy=AliTOFGeometry::XPad();
80 fDz=AliTOFGeometry::ZPad();
e0ddb533 81 fDebugStreamer = new TTreeSRedirector("TOFdebug.root");
d88fbf15 82}
83//_____________________________________________________________________________
58d8d9a3 84AliTOFtrackerMI::AliTOFtrackerMI(const AliTOFtrackerMI &t):
85 AliTracker(),
e0ddb533 86 fRecoParam(0x0),
58d8d9a3 87 fGeom(0x0),
e0ddb533 88 fPid(0x0),
58d8d9a3 89 fN(0),
90 fNseeds(0),
91 fNseedsTOF(0),
92 fngoodmatch(0),
93 fnbadmatch(0),
94 fnunmatch(0),
95 fnmatch(0),
8dacd1bb 96 fR(379.),
58d8d9a3 97 fTOFHeigth(15.3),
98 fdCut(3.),
99 fDx(1.5),
100 fDy(0),
101 fDz(0),
b42a5ea1 102 fTracks(new TClonesArray("AliTOFtrack")),
103 fSeeds(new TClonesArray("AliESDtrack")),
58d8d9a3 104 fDebugStreamer(0x0)
105 {
d88fbf15 106 //AliTOFtrackerMI copy Ctor
107
d88fbf15 108 fNseeds=t.fNseeds;
109 fNseedsTOF=t.fNseedsTOF;
110 fngoodmatch=t.fngoodmatch;
111 fnbadmatch=t.fnbadmatch;
112 fnunmatch=t.fnunmatch;
113 fnmatch=t.fnmatch;
e0ddb533 114 fRecoParam = t.fRecoParam;
96f01799 115 fGeom=t.fGeom;
e0ddb533 116 fPid = t.fPid;
d88fbf15 117 fR=t.fR;
118 fTOFHeigth=t.fTOFHeigth;
119 fdCut=t.fdCut;
120 fDy=t.fDy;
121 fDz=t.fDz;
122 fDx=1.5;
123 fSeeds=t.fSeeds;
124 fTracks=t.fTracks;
125 fN=t.fN;
126}
127
128//_________________________________________________________________________________
7aeeaf38 129AliTOFtrackerMI& AliTOFtrackerMI::operator=(const AliTOFtrackerMI &t)
130{
131 //AliTOFtrackerMI assignment operator
132
7aeeaf38 133 this->fNseeds=t.fNseeds;
134 this->fNseedsTOF=t.fNseedsTOF;
135 this->fngoodmatch=t.fngoodmatch;
136 this->fnbadmatch=t.fnbadmatch;
137 this->fnunmatch=t.fnunmatch;
138 this->fnmatch=t.fnmatch;
e0ddb533 139 this->fRecoParam = t.fRecoParam;
7aeeaf38 140 this->fGeom = t.fGeom;
e0ddb533 141 this->fPid = t.fPid;
7aeeaf38 142 this->fR=t.fR;
143 this->fTOFHeigth=t.fTOFHeigth;
144 this->fdCut=t.fdCut;
145 this->fDy=t.fDy;
146 this->fDz=t.fDz;
147 this->fDx=t.fDx;
148 this->fSeeds=t.fSeeds;
149 this->fTracks=t.fTracks;
150 this->fN=t.fN;
151 return *this;
152
153}
154
155//_____________________________________________________________________________
d88fbf15 156AliTOFtrackerMI::~AliTOFtrackerMI(){
157 //
158 //
159 //
160 if (fDebugStreamer) {
161 //fDebugStreamer->Close();
162 delete fDebugStreamer;
163 }
e0ddb533 164 delete fRecoParam;
165 delete fGeom;
166 delete fPid;
b42a5ea1 167 if (fTracks){
168 fTracks->Delete();
169 delete fTracks;
170 fTracks=0x0;
171 }
172 if (fSeeds){
173 fSeeds->Delete();
174 delete fSeeds;
175 fSeeds=0x0;
176 }
d88fbf15 177}
178
d88fbf15 179//_____________________________________________________________________________
af885e0f 180Int_t AliTOFtrackerMI::PropagateBack(AliESDEvent* event) {
d88fbf15 181 //
182 // Gets seeds from ESD event and Match with TOF Clusters
183 //
184
185
186 //Initialise some counters
187
188 fNseeds=0;
189 fNseedsTOF=0;
190 fngoodmatch=0;
191 fnbadmatch=0;
192 fnunmatch=0;
193 fnmatch=0;
194
195 Int_t ntrk=event->GetNumberOfTracks();
196 fNseeds = ntrk;
d88fbf15 197 TClonesArray &aESDTrack = *fSeeds;
198
199
200 //Load ESD tracks into a local Array of ESD Seeds
201
202 for (Int_t i=0; i<fNseeds; i++) {
203 AliESDtrack *t=event->GetTrack(i);
204 new(aESDTrack[i]) AliESDtrack(*t);
205 }
206
207 //Prepare ESD tracks candidates for TOF Matching
208 CollectESD();
209
210 //First Step with Strict Matching Criterion
211 //MatchTracks(kFALSE);
212
213 //Second Step with Looser Matching Criterion
214 //MatchTracks(kTRUE);
215 MatchTracksMI(kFALSE); // assign track to clusters
216 MatchTracksMI(kTRUE); // assign clusters to esd
217
218 Info("PropagateBack","Number of matched tracks: %d",fnmatch);
219 Info("PropagateBack","Number of good matched tracks: %d",fngoodmatch);
220 Info("PropagateBack","Number of bad matched tracks: %d",fnbadmatch);
221
222 //Update the matched ESD tracks
223
224 for (Int_t i=0; i<ntrk; i++) {
225 AliESDtrack *t=event->GetTrack(i);
226 AliESDtrack *seed =(AliESDtrack*)fSeeds->UncheckedAt(i);
227 if(seed->GetTOFsignal()>0){
228 t->SetTOFsignal(seed->GetTOFsignal());
229 t->SetTOFcluster(seed->GetTOFcluster());
230 Int_t tlab[3];
231 seed->GetTOFLabel(tlab);
232 t->SetTOFLabel(tlab);
233 AliTOFtrack *track = new AliTOFtrack(*seed);
234 Float_t info[10];
235 Double_t times[10];
236 seed->GetTOFInfo(info);
237 seed->GetIntegratedTimes(times);
238 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
239 t->SetIntegratedLength(seed->GetIntegratedLength());
240 t->SetIntegratedTimes(times);
a533f541 241 t->SetTOFsignalToT(seed->GetTOFsignalToT());
242 t->SetTOFCalChannel(seed->GetTOFCalChannel());
d88fbf15 243 //
244 t->SetTOFInfo(info);
245 delete track;
246 }
247 }
248
249
250 //Make TOF PID
e0ddb533 251 fPid->MakePID(event);
d88fbf15 252
b42a5ea1 253 fSeeds->Clear();
254 fTracks->Clear();
d88fbf15 255 return 0;
256
257}
258//_________________________________________________________________________
259void AliTOFtrackerMI::CollectESD() {
260 //prepare the set of ESD tracks to be matched to clusters in TOF
261
d88fbf15 262 TClonesArray &aTOFTrack = *fTracks;
263 Int_t c0=0;
264 Int_t c1=0;
265 for (Int_t i=0; i<fNseeds; i++) {
266
267 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(i);
268 if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
269
8dacd1bb 270 // TRD good tracks, already propagated at 372 cm
d88fbf15 271
272 AliTOFtrack *track = new AliTOFtrack(*t); // New
273 Double_t x = track->GetX(); //New
274
275 if (((t->GetStatus()&AliESDtrack::kTRDout)!=0 ) &&
96f01799 276 ( x >= AliTOFGeometry::RinTOF()) ){
d88fbf15 277 track->SetSeedIndex(i);
278 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
279 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
280 fNseedsTOF++;
281 c0++;
282 delete track;
283 }
284
285 // Propagate the rest of TPCbp
286
287 else {
e0ddb533 288 if(track->PropagateToInnerTOF()){ // temporary solution
d88fbf15 289 track->SetSeedIndex(i);
290 t->UpdateTrackParams(track,AliESDtrack::kTOFout);
291 new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
292 fNseedsTOF++;
293 c1++;
294 }
295 delete track;
296 }
297 }
298 //
299 //
300 printf("TRD\tOn\t%d\tOff\t%d\n",c0,c1);
301
302 // Sort according uncertainties on track position
303 fTracks->Sort();
304
305}
306
307
308
309
310
311
312
313//
314//
315//_________________________________________________________________________
316void AliTOFtrackerMI::MatchTracks( Bool_t /*mLastStep*/){
317 return;
318}
319//
320//
321//_________________________________________________________________________
322void AliTOFtrackerMI::MatchTracksMI(Bool_t mLastStep){
323
324 //Match ESD tracks to clusters in TOF
325 const Float_t kTofOffset = 26; // time offset
326 const Float_t kMinQuality = -6.; // minimal quality
327 const Float_t kMaxQualityD = 1.; // max delta quality if cluster used
328 const Float_t kForbiddenR = 0.1; // minimal PID according TPC
329
330 static const Double_t kMasses[]={
331 0.000511, 0.105658, 0.139570, 0.493677, 0.938272, 1.875613
332 };
333
334 Int_t nSteps=(Int_t)(fTOFHeigth/0.1);
335
96f01799 336 AliTOFcalib *calib = new AliTOFcalib();
a533f541 337
d88fbf15 338 //PH Arrays (moved outside of the loop)
339 Float_t * trackPos[4];
340 for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
341 Int_t * clind[6];
342 for (Int_t ii=0;ii<6;ii++) clind[ii] = new Int_t[fN];
343
344 // Some init
345
346 Int_t index[1000];
3c609b5c 347 Float_t quality[1000];
d88fbf15 348 Float_t dist3D[1000][6];
349 Double_t times[1000][6];
350 Float_t mintimedist[1000];
351 Float_t likelihood[1000];
352 Float_t length[1000];
353 AliTOFcluster *clusters[1000];
3c609b5c 354 Double_t tpcpid[5];
d88fbf15 355 dist3D[0][0]=1;
356
357 for (Int_t i=0; i<fNseedsTOF; i++) {
358
359 AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(i);
360 AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(track->GetSeedIndex());
361 Bool_t hasTime = ( (t->GetStatus()& AliESDtrack::kTIME)>0) ? kTRUE:kFALSE; // did we integrate time
362 Float_t trdquality = t->GetTRDQuality();
363 //
364 // Normalize tpc pid
365 //
366 t->GetTPCpid(tpcpid);
367 Double_t sumpid=0;
368 for (Int_t ipid=0;ipid<5;ipid++){
369 sumpid+=tpcpid[ipid];
370 }
371 for (Int_t ipid=0;ipid<5;ipid++){
372 if (sumpid>0) tpcpid[ipid]/=sumpid;
373 else{
374 tpcpid[ipid]=0.2;
375 }
376 }
377
378 if (trdquality<0) continue; // no chance
379 //
380 AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
381 //
382 //propagat track to the middle of TOF
383 //
8dacd1bb 384 Float_t xs = 379.2; // should be defined in the TOF geometry
96f01799 385 Double_t ymax=xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
d3c7bfac 386 Bool_t skip=kFALSE;
d88fbf15 387 Double_t ysect=trackTOFin->GetYat(xs,skip);
388 if (skip){
8dacd1bb 389 xs = 373.;
96f01799 390 ymax=xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
d88fbf15 391 ysect=trackTOFin->GetYat(xs,skip);
392 }
393 if (ysect > ymax) {
96f01799 394 if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
d88fbf15 395 continue;
396 }
397 } else if (ysect <-ymax) {
96f01799 398 if (!trackTOFin->Rotate(-AliTOFGeometry::GetAlpha())) {
d88fbf15 399 continue;
400 }
401 }
402 if(!trackTOFin->PropagateTo(xs)) {
403 continue;
404 }
405 //
406 // Determine a window around the track
407 //
408 Double_t x,par[5];
409 trackTOFin->GetExternalParameters(x,par);
410 Double_t cov[15];
411 trackTOFin->GetExternalCovariance(cov);
412 Float_t scalefact=3.;
413 Double_t dphi=
414 scalefact*
415 ((5*TMath::Sqrt(cov[0]) + 3.*fDy +10.*TMath::Abs(par[2]))/fR);
416 Double_t dz=
417 scalefact*
418 (5*TMath::Sqrt(cov[2]) + 3.*fDz +10.*TMath::Abs(par[3]));
419
420 Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
421 if (phi<-TMath::Pi())phi+=2*TMath::Pi();
422 if (phi>=TMath::Pi())phi-=2*TMath::Pi();
423 Double_t z=par[1];
424
425 Int_t nc =0;
426 Int_t nfound =0;
427 // find the clusters in the window of the track
428
429 for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
430 AliTOFcluster *c=fClusters[k];
431 if (c->GetZ() > z+dz) break;
432 // if (c->IsUsed()) continue;
433
434 Double_t dph=TMath::Abs(c->GetPhi()-phi);
435 if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
436 if (TMath::Abs(dph)>dphi) continue;
3c609b5c 437
d88fbf15 438 clind[0][nc] = c->GetDetInd(0);
439 clind[1][nc] = c->GetDetInd(1);
440 clind[2][nc] = c->GetDetInd(2);
441 clind[3][nc] = c->GetDetInd(3);
442 clind[4][nc] = c->GetDetInd(4);
3c609b5c 443 clind[5][nc] = k;
d88fbf15 444 nc++;
445 }
446
447 //
448 // select close clusters
449 //
450 Double_t mom=t->GetP();
451 // Bool_t dump = kTRUE;
452 for (Int_t icl=0; icl<nc; icl++){
453 Float_t distances[5];
454 if (nfound>=1000) break;
455 index[nfound]=clind[5][icl];
456 AliTOFcluster *cluster = fClusters[clind[5][icl]];
457 GetLinearDistances(trackTOFin,cluster, distances);
458 dist3D[nfound][0] = distances[4];
459 dist3D[nfound][1] = distances[1];
460 dist3D[nfound][2] = distances[2];
461 // cut on distance
462 if (TMath::Abs(dist3D[nfound][1])>20 || TMath::Abs(dist3D[nfound][2])>20) continue;
463 //
464 GetLikelihood(distances[1],distances[2],cov,trackTOFin, dist3D[nfound][3], dist3D[nfound][4]);
465 //
466 // cut on likelihood
467 if (dist3D[nfound][3]*dist3D[nfound][4]<0.00000000000001) continue; // log likelihood
468 if (TMath::Log(dist3D[nfound][3]*dist3D[nfound][4])<-9) continue; // log likelihood
469 //
470 clusters[nfound] = cluster;
471 //
472 //length and TOF updates
473 trackTOFin->GetIntegratedTimes(times[nfound]);
474 length[nfound] = trackTOFin->GetIntegratedLength();
475 length[nfound]+=distances[4];
476 mintimedist[nfound]=1000;
96f01799 477 Double_t tof2=AliTOFGeometry::TdcBinWidth()*cluster->GetTDC()+kTofOffset; // in ps
d88fbf15 478 // Float_t tgamma = TMath::Sqrt(cluster->GetR()*cluster->GetR()+cluster->GetZ()*cluster->GetZ())/0.03; //time for "primary" gamma
479 //if (trackTOFin->GetPt()<0.7 && TMath::Abs(tgamma-tof2)<350) continue; // gamma conversion candidate - TEMPORARY
480 for(Int_t j=0;j<=5;j++){
481 Double_t mass=kMasses[j];
482 times[nfound][j]+=distances[4]/3e-2*TMath::Sqrt(mom*mom+mass*mass)/mom; // add time distance
483 if ( TMath::Abs(times[nfound][j]-tof2)<mintimedist[nfound] && tpcpid[j]>kForbiddenR){
484 mintimedist[nfound]=TMath::Abs(times[nfound][j]-tof2);
485 }
486 }
487 //
488 Float_t liketime = TMath::Exp(-mintimedist[nfound]/90.)+0.25*TMath::Exp(-mintimedist[nfound]/180.);
489 if (!hasTime) liketime=0.2;
490 likelihood[nfound] = TMath::Log(dist3D[nfound][3]*dist3D[nfound][4]*liketime);
491
492 if (TMath::Log(dist3D[nfound][3]*dist3D[nfound][4])<-1){
493 if (likelihood[nfound]<-9.) continue;
494 }
495 //
496 nfound++;
497 }
498 if (nfound == 0 ) {
499 fnunmatch++;
500 delete trackTOFin;
501 continue;
502 }
503 //
504 //choose the best cluster
505 //
3c609b5c 506 //Float_t quality[1000];
507 //Int_t index[1000];
508 for (Int_t kk=0; kk<1000; kk++) quality[kk]=0;
d88fbf15 509 //
510 AliTOFcluster * cgold=0;
511 Int_t igold =-1;
512 for (Int_t icl=0;icl<nfound;icl++){
513 quality[icl] = dist3D[icl][3]*dist3D[icl][4];
514 }
515 TMath::Sort(nfound,quality,index,kTRUE);
516 igold = index[0];
517 cgold = clusters[igold];
518 if (nfound>1 &&likelihood[index[1]]>likelihood[index[0]]){
519 if ( quality[index[0]]<quality[index[1]]+0.5){
520 igold = index[1];
521 cgold = clusters[igold];
522 }
523 }
524 //
525 //
526 Float_t qualityGold = TMath::Log(0.0000001+(quality[igold]*(0.1+TMath::Exp(-mintimedist[igold]/80.))*(0.2+trdquality)));
527 if (!mLastStep){
528 if (cgold->GetQuality()<qualityGold) cgold->SetQuality(qualityGold);
529 continue;
530 }
531 //
532 if (mLastStep){
533 //signed better cluster
534 if (cgold->GetQuality()>qualityGold+kMaxQualityD) continue;
535 if (2.*qualityGold-cgold->GetQuality()<kMinQuality) continue;
536 }
537
538 Int_t inonfake=-1;
539
540 for (Int_t icl=0;icl<nfound;icl++){
541 if (
542 (clusters[index[icl]]->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
543 ||
544 (clusters[index[icl]]->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
545 ||
546 (clusters[index[icl]]->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
547 ) {
548 inonfake = icl;
549 break;
550 }
551 }
552 fnmatch++;;
553 if (inonfake==0) fngoodmatch++;
554 else{
555 fnbadmatch++;
556 }
557
558 Int_t tlab[3];
559 tlab[0]=cgold->GetLabel(0);
560 tlab[1]=cgold->GetLabel(1);
561 tlab[2]=cgold->GetLabel(2);
562 // Double_t tof2=25.*cgold->GetTDC()-350; // in ps
96f01799 563 Double_t tof2=AliTOFGeometry::TdcBinWidth()*cgold->GetTDC()+kTofOffset; // in ps
d88fbf15 564 Float_t tgamma = TMath::Sqrt(cgold->GetR()*cgold->GetR()+cgold->GetZ()*cgold->GetZ())/0.03;
565 Float_t info[11]={dist3D[igold][0],dist3D[igold][1],dist3D[igold][2],dist3D[igold][3],dist3D[igold][4],mintimedist[igold],
566 -1,tgamma, qualityGold,cgold->GetQuality(),0};
567 // GetLinearDistances(trackTOFin,cgold,&info[6]);
568 if (inonfake>=0){
569 info[6] = inonfake;
570 // info[7] = mintimedist[index[inonfake]];
571 }
572 //
a533f541 573 // Store quantities to be used for TOF Calibration
7aeeaf38 574 Float_t tToT=cgold->GetToT(); // in ps
575 t->SetTOFsignalToT(tToT);
a533f541 576 Int_t ind[5];
577 ind[0]=cgold->GetDetInd(0);
578 ind[1]=cgold->GetDetInd(1);
579 ind[2]=cgold->GetDetInd(2);
580 ind[3]=cgold->GetDetInd(3);
581 ind[4]=cgold->GetDetInd(4);
96f01799 582 Int_t calindex = AliTOFGeometry::GetIndex(ind);
a533f541 583 t->SetTOFCalChannel(calindex);
584
d88fbf15 585 t->SetTOFInfo(info);
586 t->SetTOFsignal(tof2);
587 t->SetTOFcluster(cgold->GetIndex());
588 AliTOFtrack *trackTOFout = new AliTOFtrack(*t);
8dacd1bb 589 trackTOFout->PropagateTo(379.);
d88fbf15 590 t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);
591 t->SetIntegratedLength(length[igold]);
592 t->SetIntegratedTimes(times[igold]);
593 t->SetTOFLabel(tlab);
594 //
595 delete trackTOFin;
596 delete trackTOFout;
597 //
598 }
599 //
600 //
601 //
602 for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
603 for (Int_t ii=0;ii<6;ii++) delete [] clind[ii];
a533f541 604 delete calib;
d88fbf15 605}
606
607
608
609
610// //_________________________________________________________________________
611// Int_t AliTOFtrackerMI::LoadClusters(TTree *dTree) {
612// //--------------------------------------------------------------------
613// //This function loads the TOF clusters
614// //--------------------------------------------------------------------
615
616// TBranch *branch=dTree->GetBranch("TOF");
617// if (!branch) {
0e46b9ae 618// AliError(" can't get the branch with the TOF digits !");
d88fbf15 619// return 1;
620// }
621
622// TClonesArray dummy("AliTOFdigit",10000), *digits=&dummy;
623// branch->SetAddress(&digits);
624
625// dTree->GetEvent(0);
626// Int_t nd=digits->GetEntriesFast();
627// Info("LoadClusters","number of digits: %d",nd);
628
629// for (Int_t i=0; i<nd; i++) {
630// AliTOFdigit *d=(AliTOFdigit*)digits->UncheckedAt(i);
631// Int_t dig[5]; Float_t g[3];
632// dig[0]=d->GetSector();
633// dig[1]=d->GetPlate();
634// dig[2]=d->GetStrip();
635// dig[3]=d->GetPadz();
636// dig[4]=d->GetPadx();
637
638// fGeom->GetPos(dig,g);
639
640// Double_t h[5];
641// h[0]=TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
642// h[1]=TMath::ATan2(g[1],g[0]); h[2]=g[2];
643// h[3]=d->GetTdc(); h[4]=d->GetAdc();
644
645// AliTOFcluster *cl=new AliTOFcluster(h,d->GetTracks(),dig,i);
646// InsertCluster(cl);
647// }
648
649// return 0;
650// }
651// //_________________________________________________________________________
652// void AliTOFtrackerMI::UnloadClusters() {
653// //--------------------------------------------------------------------
654// //This function unloads TOF clusters
655// //--------------------------------------------------------------------
656// for (Int_t i=0; i<fN; i++) delete fClusters[i];
657// fN=0;
658// }
659
660
661
662//_________________________________________________________________________
663Int_t AliTOFtrackerMI::LoadClusters(TTree *cTree) {
664 //--------------------------------------------------------------------
665 //This function loads the TOF clusters
666 //--------------------------------------------------------------------
667
668 TBranch *branch=cTree->GetBranch("TOF");
669 if (!branch) {
670 AliError("can't get the branch with the TOF clusters !");
671 return 1;
672 }
673
b42a5ea1 674 static TClonesArray dummy("AliTOFcluster",10000);
675 dummy.Clear();
676 TClonesArray *clusters=&dummy;
d88fbf15 677 branch->SetAddress(&clusters);
678
679 cTree->GetEvent(0);
680 Int_t nc=clusters->GetEntriesFast();
681 AliInfo(Form("Number of clusters: %d",nc));
682
683 for (Int_t i=0; i<nc; i++) {
684 AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
685 /*
686 for (Int_t jj=0; jj<5; jj++) dig[jj]=c->GetDetInd(jj);
687
688 Double_t h[5];
689 h[0]=c->GetR();
690 h[1]=c->GetPhi();
691 h[2]=c->GetZ();
692 h[3]=c->GetTDC();
693 h[4]=c->GetADC();
694
695 Int_t indexDig[3];
696 for (Int_t jj=0; jj<3; jj++) indexDig[jj] = c->GetLabel(jj);
697
698 AliTOFcluster *cl=new AliTOFcluster(h,c->GetTracks(),dig,i);
699 */
700
701 // fClusters[i]=c; fN++;
702 fClusters[i]=new AliTOFcluster(*c); fN++;
703
704 //AliInfo(Form("%4i %4i %f %f %f %f %f %2i %1i %2i %1i %2i",i, fClusters[i]->GetIndex(),fClusters[i]->GetZ(),fClusters[i]->GetR(),fClusters[i]->GetPhi(), fClusters[i]->GetTDC(),fClusters[i]->GetADC(),fClusters[i]->GetDetInd(0),fClusters[i]->GetDetInd(1),fClusters[i]->GetDetInd(2),fClusters[i]->GetDetInd(3),fClusters[i]->GetDetInd(4)));
705 //AliInfo(Form("%i %f",i, fClusters[i]->GetZ()));
706 }
707
708 //AliInfo(Form("Number of clusters: %d",fN));
709
710 return 0;
711}
712//_________________________________________________________________________
713void AliTOFtrackerMI::UnloadClusters() {
714 //--------------------------------------------------------------------
715 //This function unloads TOF clusters
716 //--------------------------------------------------------------------
717 for (Int_t i=0; i<fN; i++) {
718 delete fClusters[i];
719 fClusters[i] = 0x0;
720 }
721 fN=0;
722}
723
724
725
726
727//_________________________________________________________________________
728Int_t AliTOFtrackerMI::InsertCluster(AliTOFcluster *c) {
729 //--------------------------------------------------------------------
730 //This function adds a cluster to the array of clusters sorted in Z
731 //--------------------------------------------------------------------
732 if (fN==kMaxCluster) {
0e46b9ae 733 AliError("Too many clusters !");
d88fbf15 734 return 1;
735 }
736
737 if (fN==0) {fClusters[fN++]=c; return 0;}
738 Int_t i=FindClusterIndex(c->GetZ());
739 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTOFcluster*));
740 fClusters[i]=c; fN++;
741
742 return 0;
743}
744
745//_________________________________________________________________________
746Int_t AliTOFtrackerMI::FindClusterIndex(Double_t z) const {
747 //--------------------------------------------------------------------
748 // This function returns the index of the nearest cluster
749 //--------------------------------------------------------------------
750 if (fN==0) return 0;
751 if (z <= fClusters[0]->GetZ()) return 0;
752 if (z > fClusters[fN-1]->GetZ()) return fN;
753 Int_t b=0, e=fN-1, m=(b+e)/2;
754 for (; b<e; m=(b+e)/2) {
755 if (z > fClusters[m]->GetZ()) b=m+1;
756 else e=m;
757 }
758 return m;
759}
760
761
762
763Float_t AliTOFtrackerMI::GetLinearDistances(AliTOFtrack * track, AliTOFcluster *cluster, Float_t distances[5])
764{
765 //
766 // calclates distance between cluster and track
767 // use linear aproximation
768 //
769 const Float_t kRaddeg = 180/3.14159265358979312;
770 //
771 // Float_t tiltangle = fGeom->GetAngles(cluster->fdetIndex[1],cluster->fdetIndex[2])/kRaddeg; //tiltangle
772 Int_t cind[5];
773 cind[0]= cluster->GetDetInd(0);
774 cind[1]= cluster->GetDetInd(1);
775 cind[2]= cluster->GetDetInd(2);
776 cind[3]= cluster->GetDetInd(3);
777 cind[4]= cluster->GetDetInd(4);
96f01799 778 Float_t tiltangle = AliTOFGeometry::GetAngles(cluster->GetDetInd(1),cluster->GetDetInd(2))/kRaddeg; //tiltangle
d88fbf15 779
780 Float_t cpos[3]; //cluster position
781 Float_t cpos0[3]; //cluster position
782 // fGeom->GetPos(cluster->fdetIndex,cpos);
783 //fGeom->GetPos(cluster->fdetIndex,cpos0);
784 //
785 fGeom->GetPos(cind,cpos);
786 fGeom->GetPos(cind,cpos0);
787
788 Float_t phi = TMath::ATan2(cpos[1],cpos[0]);
789 if(phi<0) phi=2.*TMath::Pi()+phi;
790 // Get the local angle in the sector philoc
791 Float_t phiangle = (Int_t (phi*kRaddeg/20.) + 0.5)*20./kRaddeg;
792 //
793 Double_t v0[3];
794 Double_t dir[3];
6c94f330 795 track->GetXYZ(v0);
796 track->GetPxPyPz(dir);
d88fbf15 797 dir[0]/=track->GetP();
798 dir[1]/=track->GetP();
799 dir[2]/=track->GetP();
800 //
801 //
802 //rotate 0
803 Float_t sinphi = TMath::Sin(phiangle);
804 Float_t cosphi = TMath::Cos(phiangle);
805 Float_t sinth = TMath::Sin(tiltangle);
806 Float_t costh = TMath::Cos(tiltangle);
807 //
808 Float_t temp;
809 temp = cpos[0]*cosphi+cpos[1]*sinphi;
810 cpos[1] = -cpos[0]*sinphi+cpos[1]*cosphi;
811 cpos[0] = temp;
812 temp = v0[0]*cosphi+v0[1]*sinphi;
813 v0[1] = -v0[0]*sinphi+v0[1]*cosphi;
814 v0[0] = temp;
815 //
816 temp = cpos[0]*costh+cpos[2]*sinth;
817 cpos[2] = -cpos[0]*sinth+cpos[2]*costh;
818 cpos[0] = temp;
819 temp = v0[0]*costh+v0[2]*sinth;
820 v0[2] = -v0[0]*sinth+v0[2]*costh;
821 v0[0] = temp;
822 //
823 //
824 //rotate direction vector
825 //
826 temp = dir[0]*cosphi+dir[1]*sinphi;
827 dir[1] = -dir[0]*sinphi+dir[1]*cosphi;
828 dir[0] = temp;
829 //
830 temp = dir[0]*costh+dir[2]*sinth;
831 dir[2] = -dir[0]*sinth+dir[2]*costh;
832 dir[0] = temp;
833 //
834 Float_t v3[3];
835 Float_t k = (cpos[0]-v0[0])/dir[0];
836 v3[0] = v0[0]+k*dir[0];
837 v3[1] = v0[1]+k*dir[1];
838 v3[2] = v0[2]+k*dir[2];
839 //
840 distances[0] = v3[0]-cpos[0];
841 distances[1] = v3[1]-cpos[1];
842 distances[2] = v3[2]-cpos[2];
843 distances[3] = TMath::Sqrt( distances[0]*distances[0]+distances[1]*distances[1]+distances[2]*distances[2]); //distance
844 distances[4] = k; //length
845
846 //
847 // Debuging part of the matching
848 //
849 if (track->GetLabel()==cluster->GetLabel(0) ||track->GetLabel()==cluster->GetLabel(1) ){
850 TTreeSRedirector& cstream = *fDebugStreamer;
851 Float_t tdc = cluster->GetTDC();
852 cstream<<"Tracks"<<
853 "TOF.="<<track<<
854 "Cx="<<cpos0[0]<<
855 "Cy="<<cpos0[1]<<
856 "Cz="<<cpos0[2]<<
857 "Dist="<<k<<
858 "Dist0="<<distances[0]<<
859 "Dist1="<<distances[1]<<
860 "Dist2="<<distances[2]<<
861 "TDC="<<tdc<<
862 "\n";
863 }
864 return distances[3];
865}
866
867
868
869void AliTOFtrackerMI::GetLikelihood(Float_t dy, Float_t dz, const Double_t *cov, AliTOFtrack * /*track*/, Float_t & py, Float_t &pz)
870{
871 //
872 // get likelihood - track covariance taken
873 // 75 % of gauss with expected sigma
874 // 25 % of gauss with extended sigma
875
876 Double_t kMaxSigmaY = 0.6; // ~ 90% of TRD tracks
877 Double_t kMaxSigmaZ = 1.2; // ~ 90% of TRD tracks
878 Double_t kMeanSigmaY = 0.25; // mean TRD sigma
879 Double_t kMeanSigmaZ = 0.5; // mean TRD sigma
880
881
882 Float_t normwidth, normd, p0,p1;
883 Float_t sigmay = TMath::Max(TMath::Sqrt(cov[0]+kMeanSigmaY*kMeanSigmaY),kMaxSigmaY);
884 Float_t sigmaz = TMath::Max(TMath::Sqrt(cov[2]+kMeanSigmaZ*kMeanSigmaZ),kMaxSigmaZ);
885
886 py=0;
887 pz=0;
888 //
889 // py calculation - 75% admixture of original sigma - 25% tails
890 //
891 normwidth = fDy/sigmay;
892 normd = dy/sigmay;
893 p0 = 0.5*(1+TMath::Erf(normd-normwidth*0.5));
894 p1 = 0.5*(1+TMath::Erf(normd+normwidth*0.5));
895 py+= 0.75*(p1-p0);
896 //
897 normwidth = fDy/(3.*sigmay);
898 normd = dy/(3.*sigmay);
899 p0 = 0.5*(1+TMath::Erf(normd-normwidth*0.5));
900 p1 = 0.5*(1+TMath::Erf(normd+normwidth*0.5));
901 py+= 0.25*(p1-p0);
902 //
903 // pz calculation - 75% admixture of original sigma - 25% tails
904 //
905 normwidth = fDz/sigmaz;
906 normd = dz/sigmaz;
907 p0 = 0.5*(1+TMath::Erf(normd-normwidth*0.5));
908 p1 = 0.5*(1+TMath::Erf(normd+normwidth*0.5));
909 pz+= 0.75*(p1-p0);
910 //
911 normwidth = fDz/(3.*sigmaz);
912 normd = dz/(3.*sigmaz);
913 p0 = 0.5*(1+TMath::Erf(normd-normwidth*0.5));
914 p1 = 0.5*(1+TMath::Erf(normd+normwidth*0.5));
915 pz+= 0.25*(p1-p0);
916}