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