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