]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFtracker.cxx
Some additional protections
[u/mrichter/AliRoot.git] / TOF / AliTOFtracker.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 // -- 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                             //
26 //                                                                    //
27 //--------------------------------------------------------------------//
28
29 #include "Rtypes.h"
30
31 #include "TClonesArray.h"
32 #include "TGeoManager.h"
33 #include "TTree.h"
34
35 #include "AliAlignObj.h"
36 #include "AliESDtrack.h"
37 #include "AliESD.h"
38 #include "AliLog.h"
39 #include "AliTrackPointArray.h"
40
41 #include "AliTOFcalib.h"
42 #include "AliTOFcluster.h"
43 #include "AliTOFGeometry.h"
44 #include "AliTOFtracker.h"
45 #include "AliTOFtrack.h"
46
47 extern TGeoManager *gGeoManager;
48
49 ClassImp(AliTOFtracker)
50
51 //_____________________________________________________________________________
52 AliTOFtracker::AliTOFtracker(AliTOFGeometry * geom, Double_t parPID[2]) { 
53   //AliTOFtracker main Ctor
54
55   //fHoles=true;
56   fNseeds=0;
57   fNseedsTOF=0;
58   fngoodmatch=0;
59   fnbadmatch=0;
60   fnunmatch=0;
61   fnmatch=0;
62   fGeom = geom;
63   fTOFpid = new AliTOFpidESD(parPID);
64   fR=378.; 
65   fTOFHeigth=15.3;  
66   fdCut=3.; 
67   fDy=AliTOFGeometry::XPad(); 
68   fDz=AliTOFGeometry::ZPad(); 
69   fDx=1.5; 
70   fDzMax=35.; 
71   fDyMax=50.; 
72   fSeeds=0x0;
73   fTracks=0x0;
74   fN=0;
75   fHoles = fGeom->GetHoles();
76 }
77 //_____________________________________________________________________________
78 AliTOFtracker::AliTOFtracker(const AliTOFtracker &t):AliTracker() { 
79   //AliTOFtracker copy Ctor
80
81   fHoles=t.fHoles;
82   fNseeds=t.fNseeds;
83   fNseedsTOF=t.fNseedsTOF;
84   fngoodmatch=t.fngoodmatch;
85   fnbadmatch=t.fnbadmatch;
86   fnunmatch=t.fnunmatch;
87   fnmatch=t.fnmatch;
88   fGeom = t.fGeom;
89   fTOFpid = t.fTOFpid;
90   fR=t.fR; 
91   fTOFHeigth=t.fTOFHeigth;  
92   fdCut=t.fdCut; 
93   fDy=t.fDy; 
94   fDz=t.fDz; 
95   fDx=t.fDx; 
96   fDzMax=t.fDzMax; 
97   fDyMax=t.fDyMax; 
98   fSeeds=t.fSeeds;
99   fTracks=t.fTracks;
100   fN=t.fN;
101 }
102
103 //_____________________________________________________________________________
104 AliTOFtracker& AliTOFtracker::operator=(const AliTOFtracker &t)
105
106   //AliTOFtracker assignment operator
107
108   this->fHoles=t.fHoles;
109   this->fNseeds=t.fNseeds;
110   this->fNseedsTOF=t.fNseedsTOF;
111   this->fngoodmatch=t.fngoodmatch;
112   this->fnbadmatch=t.fnbadmatch;
113   this->fnunmatch=t.fnunmatch;
114   this->fnmatch=t.fnmatch;
115   this->fGeom = t.fGeom;
116   this->fTOFpid = t.fTOFpid;
117   this->fR=t.fR; 
118   this->fTOFHeigth=t.fTOFHeigth;  
119   this->fdCut=t.fdCut; 
120   this->fDy=t.fDy; 
121   this->fDz=t.fDz; 
122   this->fDx=t.fDx; 
123   this->fDzMax=t.fDzMax; 
124   this->fDyMax=t.fDyMax; 
125   this->fSeeds=t.fSeeds;
126   this->fTracks=t.fTracks;
127   this->fN=t.fN;
128   return *this;
129
130 }
131
132 //_____________________________________________________________________________
133 Int_t AliTOFtracker::PropagateBack(AliESD* event) {
134   //
135   // Gets seeds from ESD event and Match with TOF Clusters
136   //
137
138
139   //Initialise some counters
140
141   fNseeds=0;
142   fNseedsTOF=0;
143   fngoodmatch=0;
144   fnbadmatch=0;
145   fnunmatch=0;
146   fnmatch=0;
147
148   Int_t ntrk=event->GetNumberOfTracks();
149   fNseeds = ntrk;
150   fSeeds= new TClonesArray("AliESDtrack",ntrk);
151   TClonesArray &aESDTrack = *fSeeds;
152
153
154   //Load ESD tracks into a local Array of ESD Seeds
155
156   for (Int_t i=0; i<fNseeds; i++) {
157     AliESDtrack *t=event->GetTrack(i);
158     new(aESDTrack[i]) AliESDtrack(*t);
159   }
160
161   //Prepare ESD tracks candidates for TOF Matching
162   CollectESD();
163
164   //First Step with Strict Matching Criterion
165   MatchTracks(kFALSE);
166   /*
167   for (Int_t ijk=0; ijk<fN; ijk++) {
168     AliInfo(Form("%4i %4i  %f %f %f  %f %f   %2i %1i %2i %1i %2i",ijk, fClusters[ijk]->GetIndex(),fClusters[ijk]->GetZ(),fClusters[ijk]->GetR(),fClusters[ijk]->GetPhi(), fClusters[ijk]->GetTDC(),fClusters[ijk]->GetADC(),fClusters[ijk]->GetDetInd(0),fClusters[ijk]->GetDetInd(1),fClusters[ijk]->GetDetInd(2),fClusters[ijk]->GetDetInd(3),fClusters[ijk]->GetDetInd(4)));
169   }
170   */
171
172   //Second Step with Looser Matching Criterion
173   MatchTracks(kTRUE);
174
175   AliInfo(Form("Number of matched tracks: %d",fnmatch));
176   AliInfo(Form("Number of good matched tracks: %d",fngoodmatch));
177   AliInfo(Form("Number of bad  matched tracks: %d",fnbadmatch));
178
179   //Update the matched ESD tracks
180
181   for (Int_t i=0; i<ntrk; i++) {
182     AliESDtrack *t=event->GetTrack(i);
183     AliESDtrack *seed =(AliESDtrack*)fSeeds->UncheckedAt(i);
184     if(seed->GetTOFsignal()>0){
185       t->SetTOFsignal(seed->GetTOFsignal());
186       t->SetTOFcluster(seed->GetTOFcluster());
187       t->SetTOFsignalToT(seed->GetTOFsignalToT());
188       t->SetTOFCalChannel(seed->GetTOFCalChannel());
189       AliTOFtrack *track = new AliTOFtrack(*seed); 
190       t->UpdateTrackParams(track,AliESDtrack::kTOFout);   
191       delete track;
192     }
193   }
194
195
196   //Make TOF PID
197   fTOFpid->MakePID(event);
198
199   if (fSeeds) {
200     fSeeds->Delete();
201     delete fSeeds;
202     fSeeds = 0x0;
203   }
204   if (fTracks) {
205     fTracks->Delete();
206     delete fTracks;
207     fTracks = 0x0;
208   }
209   return 0;
210   
211 }
212 //_________________________________________________________________________
213 void AliTOFtracker::CollectESD() {
214    //prepare the set of ESD tracks to be matched to clusters in TOF
215  
216   fTracks= new TClonesArray("AliTOFtrack");
217   TClonesArray &aTOFTrack = *fTracks;
218   for (Int_t i=0; i<fNseeds; i++) {
219
220     AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(i);
221     if ((t->GetStatus()&AliESDtrack::kTPCout)==0)continue;
222
223     // TRD good tracks, already propagated at 371 cm
224
225     AliTOFtrack *track = new AliTOFtrack(*t); // New
226     Double_t x = track->GetX(); //New
227
228     if (((t->GetStatus()&AliESDtrack::kTRDout)!=0 ) && 
229          ( x >= fGeom->RinTOF()) ){
230       track->SetSeedIndex(i);
231       t->UpdateTrackParams(track,AliESDtrack::kTOFout);    
232       new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
233       fNseedsTOF++;
234       delete track;
235     }
236
237     // Propagate the rest of TPCbp  
238
239     else {
240       if(track->PropagateToInnerTOF(fHoles)){ // temporary solution
241         //      if(track->PropagateToInnerTOF(fGeom->GetHoles())){
242         track->SetSeedIndex(i);
243         t->UpdateTrackParams(track,AliESDtrack::kTOFout);    
244         new(aTOFTrack[fNseedsTOF]) AliTOFtrack(*track);
245         fNseedsTOF++;
246       }
247       delete track;
248     }
249   }
250
251   AliInfo(Form("Number of TOF seedds %i",fNseedsTOF));
252
253   // Sort according uncertainties on track position 
254   fTracks->Sort();
255
256 }
257 //_________________________________________________________________________
258 void AliTOFtracker::MatchTracks( Bool_t mLastStep){
259
260   //Match ESD tracks to clusters in TOF
261
262
263   Int_t nSteps=(Int_t)(fTOFHeigth/0.1);
264
265   AliTOFcalib *calib = new AliTOFcalib(fGeom);
266   //PH Arrays (moved outside of the loop)
267   Float_t * trackPos[4];
268   for (Int_t ii=0; ii<4; ii++) trackPos[ii] = new Float_t[nSteps];
269   Int_t * clind[6];
270   for (Int_t ii=0;ii<6;ii++) clind[ii] = new Int_t[fN];
271   
272   for (Int_t iseed=0; iseed<fNseedsTOF; iseed++) {
273
274     AliTOFtrack *track =(AliTOFtrack*)fTracks->UncheckedAt(iseed);
275     AliESDtrack *t =(AliESDtrack*)fSeeds->UncheckedAt(track->GetSeedIndex());
276     if(t->GetTOFsignal()>0. ) continue;
277     AliTOFtrack *trackTOFin =new AliTOFtrack(*track);
278
279     // Some init
280
281     Int_t         index[10000];
282     Float_t        dist[10000];
283     Float_t       cxpos[10000];
284     Float_t       crecL[10000];
285     TGeoHMatrix   global[1000];
286      
287     // Determine a window around the track
288
289     Double_t x,par[5]; 
290     trackTOFin->GetExternalParameters(x,par);
291     Double_t cov[15]; 
292     trackTOFin->GetExternalCovariance(cov);
293     Float_t scalefact=3.;    
294     Double_t dphi=
295       scalefact*
296       ((5*TMath::Sqrt(cov[0]) + 0.5*fDy + 2.5*TMath::Abs(par[2]))/fR); 
297     Double_t dz=
298       scalefact*
299       (5*TMath::Sqrt(cov[2]) + 0.5*fDz + 2.5*TMath::Abs(par[3]));
300     
301     Double_t phi=TMath::ATan2(par[0],x) + trackTOFin->GetAlpha();
302     if (phi<-TMath::Pi())phi+=2*TMath::Pi();
303     if (phi>=TMath::Pi())phi-=2*TMath::Pi();
304     Double_t z=par[1];   
305
306     //upper limit on window's size.
307
308     if(dz> fDzMax) dz=fDzMax;
309     if(dphi*fR>fDyMax) dphi=fDyMax/fR;
310
311
312     Int_t nc=0;
313
314     // find the clusters in the window of the track
315
316     for (Int_t k=FindClusterIndex(z-dz); k<fN; k++) {
317
318       AliTOFcluster *c=fClusters[k];
319       if (c->GetZ() > z+dz) break;
320       if (c->IsUsed()) continue;
321       
322       //AliInfo(Form(" fClusters[k]->GetZ() (%f) z-dz (%f)   %4i ", fClusters[k]->GetZ(), z-dz, k));
323
324       Double_t dph=TMath::Abs(c->GetPhi()-phi);
325       if (dph>TMath::Pi()) dph-=2.*TMath::Pi();
326       if (TMath::Abs(dph)>dphi) continue;
327     
328       clind[0][nc] = c->GetDetInd(0);
329       clind[1][nc] = c->GetDetInd(1);
330       clind[2][nc] = c->GetDetInd(2);
331       clind[3][nc] = c->GetDetInd(3);
332       clind[4][nc] = c->GetDetInd(4);
333       clind[5][nc] = k;      
334       Char_t path[100];
335       Int_t ind[5];
336       ind[0]=clind[0][nc];
337       ind[1]=clind[1][nc];
338       ind[2]=clind[2][nc];
339       ind[3]=clind[3][nc];
340       ind[4]=clind[4][nc];
341       fGeom->GetVolumePath(ind,path);
342       gGeoManager->cd(path);
343       global[nc] = *gGeoManager->GetCurrentMatrix();
344       nc++;
345     }
346
347     //if (nc) AliInfo(Form("seed for TOF %4i and number of clusters in the track window %4i (cluster index %4i)     %4i",i,nc, clind[5][0], fN));
348
349     //start fine propagation 
350
351     Int_t nStepsDone = 0;
352     for( Int_t istep=0; istep<nSteps; istep++){ 
353
354       Float_t xs=fGeom->RinTOF()+istep*0.1;
355       Double_t ymax=xs*TMath::Tan(0.5*AliTOFGeometry::GetAlpha());
356
357       Bool_t skip=kFALSE;
358       Double_t ysect=trackTOFin->GetYat(xs,skip);
359       if(skip)break;
360       if (ysect > ymax) {
361         if (!trackTOFin->Rotate(AliTOFGeometry::GetAlpha())) {
362           break;
363         }
364       } else if (ysect <-ymax) {
365         if (!trackTOFin->Rotate(-AliTOFGeometry::GetAlpha())) {
366           break;
367         }
368       }
369
370       if(!trackTOFin->PropagateTo(xs)) {
371         break;
372       }
373
374       nStepsDone++;
375
376       // store the running point (Globalrf) - fine propagation     
377
378       Double_t x,y,z;
379       trackTOFin->GetGlobalXYZ(x,y,z);
380       trackPos[0][istep]= (Float_t) x;
381       trackPos[1][istep]= (Float_t) y;
382       trackPos[2][istep]= (Float_t) z;   
383       trackPos[3][istep]= trackTOFin->GetIntegratedLength();
384     }
385
386
387     Int_t nfound = 0;
388     for (Int_t istep=0; istep<nStepsDone; istep++) {
389
390       Bool_t isInside =kFALSE;
391       Float_t ctrackPos[3];     
392
393       ctrackPos[0]= trackPos[0][istep];
394       ctrackPos[1]= trackPos[1][istep];
395       ctrackPos[2]= trackPos[2][istep];
396
397       //now see whether the track matches any of the TOF clusters            
398
399       for (Int_t i=0; i<nc; i++){
400         Int_t cind[5];
401         cind[0]= clind[0][i];
402         cind[1]= clind[1][i];
403         cind[2]= clind[2][i];
404         cind[3]= clind[3][i];
405         cind[4]= clind[4][i];
406         Bool_t accept = kFALSE;
407         if( mLastStep)accept = (fGeom->DistanceToPad(cind,global[i],ctrackPos)<fdCut);
408         if(!mLastStep)accept = (fGeom->IsInsideThePad(cind,global[i],ctrackPos));
409         if(accept){
410           if(!mLastStep)isInside=kTRUE;
411           dist[nfound]=fGeom->DistanceToPad(cind,global[i],ctrackPos);
412           crecL[nfound]=trackPos[3][istep];
413           index[nfound]=clind[5][i]; // store cluster id            
414           cxpos[nfound]=fGeom->RinTOF()+istep*0.1; //store prop.radius
415           nfound++;
416           if(isInside)break;
417         }//end if accept
418       } //end for on the clusters
419
420
421       if(isInside)break;
422     } //end for on the steps     
423
424
425
426     if (nfound == 0 ) {
427       fnunmatch++;
428       delete trackTOFin;
429       continue;
430     }
431     
432     fnmatch++;
433
434     // now choose the cluster to be matched with the track.
435
436     Int_t idclus=0;
437     Float_t  recL = 0.;
438     Float_t  xpos=0.;
439     Float_t  mindist=1000.;
440     for (Int_t iclus= 0; iclus<nfound;iclus++){
441       if (dist[iclus]< mindist){
442         mindist = dist[iclus];
443         xpos = cxpos[iclus];
444         idclus =index[iclus]; 
445         recL=crecL[iclus]+fDx*0.5;
446       }
447     }
448
449     AliTOFcluster *c=fClusters[idclus];
450     c->Use(); //AliInfo(Form("I am using the cluster"));
451
452     // Track length correction for matching Step 2 
453
454     if(mLastStep){
455       Float_t rc=TMath::Sqrt(c->GetR()*c->GetR() + c->GetZ()*c->GetZ());
456       Float_t rt=TMath::Sqrt(trackPos[0][70]*trackPos[0][70]
457                              +trackPos[1][70]*trackPos[1][70]
458                              +trackPos[2][70]*trackPos[2][70]);
459       Float_t dlt=rc-rt;      
460       recL=trackPos[3][70]+dlt;
461     }    
462
463     if (
464         (c->GetLabel(0)==TMath::Abs(trackTOFin->GetLabel()))
465         ||
466         (c->GetLabel(1)==TMath::Abs(trackTOFin->GetLabel()))
467         ||
468         (c->GetLabel(2)==TMath::Abs(trackTOFin->GetLabel()))
469         ) {
470       fngoodmatch++;
471
472       //AliInfo(Form(" track label good %5i",trackTOFin->GetLabel()));
473
474     }
475     else{
476       fnbadmatch++;
477
478       //AliInfo(Form(" track label  bad %5i",trackTOFin->GetLabel()));
479
480     }
481
482     delete trackTOFin;
483
484     //  Store quantities to be used in the TOF Calibration
485     Float_t tToT=c->GetToT(); // in ps
486     t->SetTOFsignalToT(tToT);
487     Int_t ind[5];
488     ind[0]=c->GetDetInd(0);
489     ind[1]=c->GetDetInd(1);
490     ind[2]=c->GetDetInd(2);
491     ind[3]=c->GetDetInd(3);
492     ind[4]=c->GetDetInd(4);
493     Int_t calindex = calib->GetIndex(ind);
494     t->SetTOFCalChannel(calindex);
495     
496     Double_t tof=AliTOFGeometry::TdcBinWidth()*c->GetTDC()+32; // in ps
497     t->SetTOFsignal(tof);
498     //t->SetTOFcluster(c->GetIndex()); // pointing to the digits tree
499     t->SetTOFcluster(idclus); // pointing to the recPoints tree
500     Double_t time[10]; t->GetIntegratedTimes(time);
501     Double_t mom=t->GetP();
502     for(Int_t j=0;j<=AliPID::kSPECIES;j++){
503       Double_t mass=AliPID::ParticleMass(j);
504       time[j]+=(recL-trackPos[3][0])/3e-2*TMath::Sqrt(mom*mom+mass*mass)/mom;
505     }
506
507     AliTOFtrack *trackTOFout = new AliTOFtrack(*t); 
508     trackTOFout->PropagateTo(xpos);
509     t->UpdateTrackParams(trackTOFout,AliESDtrack::kTOFout);    
510     t->SetIntegratedLength(recL);
511     t->SetIntegratedTimes(time);
512
513     delete trackTOFout;
514   }
515   for (Int_t ii=0; ii<4; ii++) delete [] trackPos[ii];
516   for (Int_t ii=0;ii<6;ii++) delete [] clind[ii];
517   delete calib;
518 }
519 //_________________________________________________________________________
520 Int_t AliTOFtracker::LoadClusters(TTree *cTree) {
521   //--------------------------------------------------------------------
522   //This function loads the TOF clusters
523   //--------------------------------------------------------------------
524
525   TBranch *branch=cTree->GetBranch("TOF");
526   if (!branch) { 
527     AliError("can't get the branch with the TOF clusters !");
528     return 1;
529   }
530
531   TClonesArray dummy("AliTOFcluster",10000), *clusters=&dummy;
532   branch->SetAddress(&clusters);
533
534   cTree->GetEvent(0);
535   Int_t nc=clusters->GetEntriesFast();
536   AliInfo(Form("Number of clusters: %d",nc));
537
538   for (Int_t i=0; i<nc; i++) {
539     AliTOFcluster *c=(AliTOFcluster*)clusters->UncheckedAt(i);
540     fClusters[i]=new AliTOFcluster(*c); fN++;
541     //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)));
542     //AliInfo(Form("%i %f",i, fClusters[i]->GetZ()));
543   }
544
545   //AliInfo(Form("Number of clusters: %d",fN));
546
547   return 0;
548 }
549 //_________________________________________________________________________
550 void AliTOFtracker::UnloadClusters() {
551   //--------------------------------------------------------------------
552   //This function unloads TOF clusters
553   //--------------------------------------------------------------------
554   for (Int_t i=0; i<fN; i++) {
555     delete fClusters[i];
556     fClusters[i] = 0x0;
557   }
558   fN=0;
559 }
560
561 //_________________________________________________________________________
562 Int_t AliTOFtracker::FindClusterIndex(Double_t z) const {
563   //--------------------------------------------------------------------
564   // This function returns the index of the nearest cluster 
565   //--------------------------------------------------------------------
566   if (fN==0) return 0;
567   if (z <= fClusters[0]->GetZ()) return 0;
568   if (z > fClusters[fN-1]->GetZ()) return fN;
569   Int_t b=0, e=fN-1, m=(b+e)/2;
570   for (; b<e; m=(b+e)/2) {
571     if (z > fClusters[m]->GetZ()) b=m+1;
572     else e=m; 
573   }
574   return m;
575 }
576
577 //_________________________________________________________________________
578 Bool_t AliTOFtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
579 {
580   // Get track space point with index i
581   // Coordinates are in the global system
582   AliTOFcluster *cl = fClusters[index];
583   Float_t xyz[3];
584   xyz[0] = cl->GetR()*TMath::Cos(cl->GetPhi());
585   xyz[1] = cl->GetR()*TMath::Sin(cl->GetPhi());
586   xyz[2] = cl->GetZ();
587   Float_t phiangle = (Int_t(cl->GetPhi()*TMath::RadToDeg()/20.)+0.5)*20.*TMath::DegToRad();
588   Float_t sinphi = TMath::Sin(phiangle), cosphi = TMath::Cos(phiangle);
589   Float_t tiltangle = fGeom->GetAngles(cl->GetDetInd(1),cl->GetDetInd(2))*TMath::DegToRad();
590   Float_t sinth = TMath::Sin(tiltangle), costh = TMath::Cos(tiltangle);
591   Float_t sigmay2 = fGeom->XPad()*fGeom->XPad()/12.;
592   Float_t sigmaz2 = fGeom->ZPad()*fGeom->ZPad()/12.;
593   Float_t cov[6];
594   cov[0] = sinphi*sinphi*sigmay2 + cosphi*cosphi*sinth*sinth*sigmaz2;
595   cov[1] = -sinphi*cosphi*sigmay2 + sinphi*cosphi*sinth*sinth*sigmaz2;
596   cov[2] = -cosphi*sinth*costh*sigmaz2;
597   cov[3] = cosphi*cosphi*sigmay2 + sinphi*sinphi*sinth*sinth*sigmaz2;
598   cov[4] = -sinphi*sinth*costh*sigmaz2;
599   cov[5] = costh*costh*sigmaz2;
600   p.SetXYZ(xyz[0],xyz[1],xyz[2],cov);
601
602   // Detector numbering scheme
603   Int_t nSector = fGeom->NSectors();
604   Int_t nPlate  = fGeom->NPlates();
605   Int_t nStripA = fGeom->NStripA();
606   Int_t nStripB = fGeom->NStripB();
607   Int_t nStripC = fGeom->NStripC();
608
609   Int_t isector = cl->GetDetInd(0);
610   if (isector >= nSector)
611     AliError(Form("Wrong sector number in TOF (%d) !",isector));
612   Int_t iplate = cl->GetDetInd(1);
613   if (iplate >= nPlate)
614     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
615   Int_t istrip = cl->GetDetInd(2);
616
617   Int_t stripOffset = 0;
618   switch (iplate) {
619   case 0:
620     stripOffset = 0;
621     break;
622   case 1:
623     stripOffset = nStripC;
624     break;
625   case 2:
626     stripOffset = nStripC+nStripB;
627     break;
628   case 3:
629     stripOffset = nStripC+nStripB+nStripA;
630     break;
631   case 4:
632     stripOffset = nStripC+nStripB+nStripA+nStripB;
633     break;
634   default:
635     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
636     break;
637   };
638
639   Int_t idet = (2*(nStripC+nStripB)+nStripA)*isector +
640                stripOffset +
641                istrip;
642   UShort_t volid = AliAlignObj::LayerToVolUID(AliAlignObj::kTOF,idet);
643   p.SetVolumeID((UShort_t)volid);
644   return kTRUE;
645 }