]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITStrackV2.cxx
Add plot for check of drift speed calibration
[u/mrichter/AliRoot.git] / ITS / AliITStrackV2.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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////
19 //                Implementation of the ITS track class
20 //
21 //          Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
22 //     dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
23 ///////////////////////////////////////////////////////////////////////////
24 #include <TMath.h>
25
26 #include "AliCluster.h"
27 #include "AliESDVertex.h"
28 #include "AliITSReconstructor.h"
29 #include "AliITStrackV2.h"
30 #include "AliTracker.h"
31
32 const Int_t AliITStrackV2::fgkWARN = 5;
33
34 ClassImp(AliITStrackV2)
35
36
37 //____________________________________________________________________________
38 AliITStrackV2::AliITStrackV2() : AliKalmanTrack(),
39   fCheckInvariant(kTRUE),
40   fdEdx(0),
41   fESDtrack(0)
42 {
43   for(Int_t i=0; i<2*AliITSgeomTGeo::kNLayers; i++) {fIndex[i]=-1; fModule[i]=-1;}
44   for(Int_t i=0; i<4; i++) fdEdxSample[i]=0;
45 }
46
47
48 //____________________________________________________________________________
49 AliITStrackV2::AliITStrackV2(AliESDtrack& t,Bool_t c) throw (const Char_t *) :
50   AliKalmanTrack(),
51   fCheckInvariant(kTRUE),
52   fdEdx(t.GetITSsignal()),
53   fESDtrack(&t)
54 {
55   //------------------------------------------------------------------
56   // Conversion ESD track -> ITS track.
57   // If c==kTRUE, create the ITS track out of the constrained params.
58   //------------------------------------------------------------------
59   const AliExternalTrackParam *par=&t;
60   if (c) {
61     par=t.GetConstrainedParam();
62     if (!par) throw "AliITStrackV2: conversion failed !\n";
63   }
64   Set(par->GetX(),par->GetAlpha(),par->GetParameter(),par->GetCovariance());
65
66   //if (!Invariant()) throw "AliITStrackV2: conversion failed !\n";
67
68   SetLabel(t.GetLabel());
69   SetMass(t.GetMass());
70   SetNumberOfClusters(t.GetITSclusters(fIndex));
71
72   if (t.GetStatus()&AliESDtrack::kTIME) {
73     StartTimeIntegral();
74     Double_t times[10]; t.GetIntegratedTimes(times); SetIntegratedTimes(times);
75     SetIntegratedLength(t.GetIntegratedLength());
76   }
77
78   for(Int_t i=0; i<4; i++) fdEdxSample[i]=0;
79 }
80
81 //____________________________________________________________________________
82 void AliITStrackV2::ResetClusters() {
83   //------------------------------------------------------------------
84   // Reset the array of attached clusters.
85   //------------------------------------------------------------------
86   for (Int_t i=0; i<2*AliITSgeomTGeo::kNLayers; i++) fIndex[i]=-1;
87   SetChi2(0.); 
88   SetNumberOfClusters(0);
89
90
91 //____________________________________________________________________________
92 void AliITStrackV2::UpdateESDtrack(ULong_t flags) const {
93   // Update track params
94   fESDtrack->UpdateTrackParams(this,flags);
95   // copy the module indices
96   for(Int_t i=0;i<12;i++) {
97     //   printf("     %d\n",GetModuleIndex(i));
98     fESDtrack->SetITSModuleIndex(i,GetModuleIndex(i));
99   }
100   // copy the 4 dedx samples
101   Double_t sdedx[4]={0.,0.,0.,0.};
102   for(Int_t i=0; i<4; i++) sdedx[i]=fdEdxSample[i];
103   fESDtrack->SetITSdEdxSamples(sdedx);
104 }
105
106 //____________________________________________________________________________
107 AliITStrackV2::AliITStrackV2(const AliITStrackV2& t) : 
108   AliKalmanTrack(t),
109   fCheckInvariant(t.fCheckInvariant),
110   fdEdx(t.fdEdx),
111   fESDtrack(t.fESDtrack) 
112 {
113   //------------------------------------------------------------------
114   //Copy constructor
115   //------------------------------------------------------------------
116   Int_t i;
117   for (i=0; i<4; i++) fdEdxSample[i]=t.fdEdxSample[i];
118   for (i=0; i<2*AliITSgeomTGeo::GetNLayers(); i++) {
119     fIndex[i]=t.fIndex[i];
120     fModule[i]=t.fModule[i];
121   }
122 }
123
124 //_____________________________________________________________________________
125 Int_t AliITStrackV2::Compare(const TObject *o) const {
126   //-----------------------------------------------------------------
127   // This function compares tracks according to the their curvature
128   //-----------------------------------------------------------------
129   AliITStrackV2 *t=(AliITStrackV2*)o;
130   //Double_t co=OneOverPt();
131   //Double_t c =OneOverPt();
132   Double_t co=t->GetSigmaY2()*t->GetSigmaZ2();
133   Double_t c =GetSigmaY2()*GetSigmaZ2();
134   if (c>co) return 1;
135   else if (c<co) return -1;
136   return 0;
137 }
138
139 //____________________________________________________________________________
140 Bool_t 
141 AliITStrackV2::PropagateToVertex(const AliESDVertex *v,Double_t d,Double_t x0) 
142 {
143   //------------------------------------------------------------------
144   //This function propagates a track to the minimal distance from the origin
145   //------------------------------------------------------------------  
146   Double_t bz=GetBz();
147   if (PropagateToDCA(v,bz,kVeryBig)) {
148     Double_t xOverX0,xTimesRho; 
149     xOverX0 = d; xTimesRho = d*x0;
150     if (CorrectForMeanMaterial(xOverX0,xTimesRho,kTRUE)) return kTRUE;
151   }
152   return kFALSE;
153 }
154
155 //____________________________________________________________________________
156 Bool_t AliITStrackV2::
157 GetGlobalXYZat(Double_t xloc, Double_t &x, Double_t &y, Double_t &z) const {
158   //------------------------------------------------------------------
159   //This function returns a track position in the global system
160   //------------------------------------------------------------------
161   Double_t r[3];
162   Bool_t rc=GetXYZAt(xloc, GetBz(), r);
163   x=r[0]; y=r[1]; z=r[2]; 
164   return rc;
165 }
166
167 //_____________________________________________________________________________
168 Double_t AliITStrackV2::GetPredictedChi2(const AliCluster *c) const {
169   //-----------------------------------------------------------------
170   // This function calculates a predicted chi2 increment.
171   //-----------------------------------------------------------------
172   Double_t p[2]={c->GetY(), c->GetZ()};
173   Double_t cov[3]={c->GetSigmaY2(), 0., c->GetSigmaZ2()};
174   return AliExternalTrackParam::GetPredictedChi2(p,cov);
175 }
176
177 //____________________________________________________________________________
178 Bool_t AliITStrackV2::PropagateTo(Double_t xk, Double_t d, Double_t x0) {
179   //------------------------------------------------------------------
180   //This function propagates a track
181   //------------------------------------------------------------------
182
183   Double_t oldX=GetX(), oldY=GetY(), oldZ=GetZ();
184   
185   //Double_t bz=GetBz();
186   //if (!AliExternalTrackParam::PropagateTo(xk,bz)) return kFALSE;
187   Double_t b[3]; GetBxByBz(b);
188   if (!AliExternalTrackParam::PropagateToBxByBz(xk,b)) return kFALSE;
189   Double_t xOverX0,xTimesRho; 
190   xOverX0 = d; xTimesRho = d*x0;
191   if (!CorrectForMeanMaterial(xOverX0,xTimesRho,kTRUE)) return kFALSE;
192
193   Double_t x=GetX(), y=GetY(), z=GetZ();
194   if (IsStartedTimeIntegral() && x>oldX) {
195     Double_t l2 = (x-oldX)*(x-oldX) + (y-oldY)*(y-oldY) + (z-oldZ)*(z-oldZ);
196     AddTimeStep(TMath::Sqrt(l2));
197   }
198
199   return kTRUE;
200 }
201
202 //____________________________________________________________________________
203 Bool_t AliITStrackV2::PropagateToTGeo(Double_t xToGo, Int_t nstep, Double_t &xOverX0, Double_t &xTimesRho, Bool_t addTime) {
204   //-------------------------------------------------------------------
205   //  Propagates the track to a reference plane x=xToGo in n steps.
206   //  These n steps are only used to take into account the curvature.
207   //  The material is calculated with TGeo. (L.Gaudichet)
208   //-------------------------------------------------------------------
209   
210   Double_t startx = GetX(), starty = GetY(), startz = GetZ();
211   Double_t sign = (startx<xToGo) ? -1.:1.;
212   Double_t step = (xToGo-startx)/TMath::Abs(nstep);
213
214   Double_t start[3], end[3], mparam[7];
215   //Double_t bz = GetBz();
216   Double_t b[3]; GetBxByBz(b);
217   Double_t bz = b[2];
218
219   Double_t x = startx;
220   
221   for (Int_t i=0; i<nstep; i++) {
222     
223     GetXYZ(start);   //starting global position
224     x += step;
225     if (!GetXYZAt(x, bz, end)) return kFALSE;
226     //if (!AliExternalTrackParam::PropagateTo(x, bz)) return kFALSE;
227     if (!AliExternalTrackParam::PropagateToBxByBz(x, b)) return kFALSE;
228     AliTracker::MeanMaterialBudget(start, end, mparam);
229     xTimesRho = sign*mparam[4]*mparam[0];
230     xOverX0   = mparam[1];
231     if (mparam[1]<900000) {
232       if (!AliExternalTrackParam::CorrectForMeanMaterial(xOverX0,
233                            xTimesRho,GetMass())) return kFALSE;
234     } else { // this happens when MeanMaterialBudget cannot cross a boundary
235       return kFALSE;
236     }
237   }
238
239   if (addTime && IsStartedTimeIntegral() && GetX()>startx) {
240     Double_t l2 = ( (GetX()-startx)*(GetX()-startx) +
241                     (GetY()-starty)*(GetY()-starty) +
242                     (GetZ()-startz)*(GetZ()-startz) );
243     AddTimeStep(TMath::Sqrt(l2));
244   }
245
246   return kTRUE;
247 }
248
249 //____________________________________________________________________________
250 Bool_t AliITStrackV2::Update(const AliCluster* c, Double_t chi2, Int_t index) 
251 {
252   //------------------------------------------------------------------
253   //This function updates track parameters
254   //------------------------------------------------------------------
255   Double_t p[2]={c->GetY(), c->GetZ()};
256   Double_t cov[3]={c->GetSigmaY2(), c->GetSigmaYZ(), c->GetSigmaZ2()};
257
258   if (!AliExternalTrackParam::Update(p,cov)) return kFALSE;
259
260   Int_t n=GetNumberOfClusters();
261   if (!Invariant()) {
262      if (n>fgkWARN) AliWarning("Wrong invariant !");
263      return kFALSE;
264   }
265
266   if (chi2<0) return kTRUE;
267
268   // fill residuals for ITS+TPC tracks 
269   if (fESDtrack) {
270     if (fESDtrack->GetStatus()&AliESDtrack::kTPCin) {
271       AliTracker::FillResiduals(this,p,cov,c->GetVolumeId());
272     }
273   }
274
275   fIndex[n]=index;
276   SetNumberOfClusters(n+1);
277   SetChi2(GetChi2()+chi2);
278
279   return kTRUE;
280 }
281
282 Bool_t AliITStrackV2::Invariant() const {
283   //------------------------------------------------------------------
284   // This function is for debugging purpose only
285   //------------------------------------------------------------------
286   if(!fCheckInvariant) return kTRUE;
287
288   Int_t n=GetNumberOfClusters();
289
290   // take into account the misalignment error
291   Float_t maxMisalErrY2=0,maxMisalErrZ2=0;
292   for (Int_t lay=0; lay<AliITSgeomTGeo::kNLayers; lay++) {
293     maxMisalErrY2 = TMath::Max(maxMisalErrY2,AliITSReconstructor::GetRecoParam()->GetClusterMisalErrorY(lay,GetBz()));
294     maxMisalErrZ2 = TMath::Max(maxMisalErrZ2,AliITSReconstructor::GetRecoParam()->GetClusterMisalErrorZ(lay,GetBz()));
295   }
296   maxMisalErrY2 *= maxMisalErrY2;
297   maxMisalErrZ2 *= maxMisalErrZ2;
298   // this is because when we reset before refitting, we multiply the
299   // matrix by 10
300   maxMisalErrY2 *= 10.; 
301   maxMisalErrZ2 *= 10.;
302
303   Double_t sP2=GetParameter()[2];
304   if (TMath::Abs(sP2) >= kAlmost1){
305      if (n>fgkWARN) Warning("Invariant","fP2=%f\n",sP2);
306      return kFALSE;
307   }
308   Double_t sC00=GetCovariance()[0];
309   if (sC00<=0 || sC00>(9.+maxMisalErrY2)) {
310      if (n>fgkWARN) Warning("Invariant","fC00=%f\n",sC00); 
311      return kFALSE;
312   }
313   Double_t sC11=GetCovariance()[2];
314   if (sC11<=0 || sC11>(9.+maxMisalErrZ2)) {
315      if (n>fgkWARN) Warning("Invariant","fC11=%f\n",sC11); 
316      return kFALSE;
317   }
318   Double_t sC22=GetCovariance()[5];
319   if (sC22<=0 || sC22>1.) {
320      if (n>fgkWARN) Warning("Invariant","fC22=%f\n",sC22); 
321      return kFALSE;
322   }
323   Double_t sC33=GetCovariance()[9];
324   if (sC33<=0 || sC33>1.) {
325      if (n>fgkWARN) Warning("Invariant","fC33=%f\n",sC33); 
326      return kFALSE;
327   }
328   Double_t sC44=GetCovariance()[14];
329   if (sC44<=0 /*|| sC44>6e-5*/) {
330      if (n>fgkWARN) Warning("Invariant","fC44=%f\n",sC44);
331      return kFALSE;
332   }
333
334   return kTRUE;
335 }
336
337 //____________________________________________________________________________
338 Bool_t AliITStrackV2::Propagate(Double_t alp,Double_t xk) {
339   //------------------------------------------------------------------
340   //This function propagates a track
341   //------------------------------------------------------------------
342   //Double_t bz=GetBz();
343   //if (!AliExternalTrackParam::Propagate(alp,xk,bz)) return kFALSE;
344   Double_t b[3]; GetBxByBz(b);
345   if (!AliExternalTrackParam::PropagateBxByBz(alp,xk,b)) return kFALSE;
346
347   if (!Invariant()) {
348     Int_t n=GetNumberOfClusters();
349     if (n>fgkWARN) AliWarning("Wrong invariant !");
350     return kFALSE;
351   }
352
353   return kTRUE;
354 }
355
356 Bool_t AliITStrackV2::MeanBudgetToPrimVertex(Double_t xyz[3], Double_t step, Double_t &d) const {
357
358   //-------------------------------------------------------------------
359   //  Get the mean material budget between the actual point and the
360   //  primary vertex. (L.Gaudichet)
361   //-------------------------------------------------------------------
362
363   Double_t cs=TMath::Cos(GetAlpha()), sn=TMath::Sin(GetAlpha());
364   Double_t vertexX = xyz[0]*cs + xyz[1]*sn;
365
366   Int_t nstep = Int_t((GetX()-vertexX)/step);
367   if (nstep<1) nstep = 1;
368   step = (GetX()-vertexX)/nstep;
369
370   //  Double_t mparam[7], densMean=0, radLength=0, length=0;
371   Double_t mparam[7];
372   Double_t p1[3], p2[3], x = GetX(), bz = GetBz();
373   GetXYZ(p1);
374
375   d=0.;
376
377   for (Int_t i=0; i<nstep; i++) {
378     x  += step;
379     if (!GetXYZAt(x, bz, p2)) return kFALSE;
380     AliTracker::MeanMaterialBudget(p1, p2, mparam);
381     if (mparam[1]>900000) return kFALSE;
382     d  += mparam[1];
383
384     p1[0] = p2[0];
385     p1[1] = p2[1];
386     p1[2] = p2[2];
387   }
388
389   return kTRUE;
390 }
391
392 Bool_t AliITStrackV2::Improve(Double_t x0,Double_t xyz[3],Double_t ers[3]) {
393   //------------------------------------------------------------------
394   //This function improves angular track parameters
395   //------------------------------------------------------------------
396   //Store the initail track parameters
397
398   Double_t x = GetX();
399   Double_t alpha = GetAlpha();
400   Double_t par[] = {GetY(),GetZ(),GetSnp(),GetTgl(),GetSigned1Pt()};
401   Double_t cov[] = {
402     GetSigmaY2(),
403     GetSigmaZY(),
404     GetSigmaZ2(),
405     GetSigmaSnpY(),
406     GetSigmaSnpZ(),
407     GetSigmaSnp2(),
408     GetSigmaTglY(),
409     GetSigmaTglZ(),
410     GetSigmaTglSnp(),
411     GetSigmaTgl2(),
412     GetSigma1PtY(),
413     GetSigma1PtZ(),
414     GetSigma1PtSnp(),
415     GetSigma1PtTgl(),
416     GetSigma1Pt2()
417   }; 
418
419
420   Double_t cs=TMath::Cos(GetAlpha()), sn=TMath::Sin(GetAlpha());
421 //Double_t xv = xyz[0]*cs + xyz[1]*sn; // vertex
422   Double_t yv =-xyz[0]*sn + xyz[1]*cs; // in the
423   Double_t zv = xyz[2];                // local frame
424
425   Double_t dy = par[0] - yv, dz = par[1] - zv;
426   Double_t r2=GetX()*GetX() + dy*dy;
427   Double_t p2=(1.+ GetTgl()*GetTgl())/(GetSigned1Pt()*GetSigned1Pt());
428   Double_t beta2=p2/(p2 + GetMass()*GetMass());
429   x0*=TMath::Sqrt((1.+ GetTgl()*GetTgl())/(1.- GetSnp()*GetSnp()));
430   Double_t theta2=14.1*14.1/(beta2*p2*1e6)*x0;
431   //Double_t theta2=1.0259e-6*14*14/28/(beta2*p2)*x0*9.36*2.33;
432
433   Double_t cnv=GetBz()*kB2C;
434   {
435     Double_t dummy = 4/r2 - GetC()*GetC();
436     if (dummy < 0) return kFALSE;
437     Double_t parp = 0.5*(GetC()*GetX() + dy*TMath::Sqrt(dummy));
438     Double_t sigma2p = theta2*(1.-GetSnp())*(1.+GetSnp())*(1. + GetTgl()*GetTgl());
439     Double_t ovSqr2 = 1./TMath::Sqrt(r2);
440     Double_t tfact = ovSqr2*(1.-dy*ovSqr2)*(1.+dy*ovSqr2);
441     sigma2p += cov[0]*tfact*tfact;
442     sigma2p += ers[1]*ers[1]/r2;
443     sigma2p += 0.25*cov[14]*cnv*cnv*GetX()*GetX();
444     Double_t eps2p=sigma2p/(cov[5] + sigma2p);
445     par[0] += cov[3]/(cov[5] + sigma2p)*(parp - GetSnp());
446     par[2]  = eps2p*GetSnp() + (1 - eps2p)*parp;
447     cov[5] *= eps2p;
448     cov[3] *= eps2p;
449   }
450   {
451     Double_t parl=0.5*GetC()*dz/TMath::ASin(0.5*GetC()*TMath::Sqrt(r2));
452     Double_t sigma2l=theta2;
453     sigma2l += cov[2]/r2 + cov[0]*dy*dy*dz*dz/(r2*r2*r2);
454     sigma2l += ers[2]*ers[2]/r2;
455     Double_t eps2l = sigma2l/(cov[9] + sigma2l);
456     par[1] += cov[7 ]/(cov[9] + sigma2l)*(parl - par[3]);
457     par[4] += cov[13]/(cov[9] + sigma2l)*(parl - par[3]);
458     par[3]  = eps2l*par[3] + (1-eps2l)*parl;
459     cov[9] *= eps2l; 
460     cov[13]*= eps2l; 
461     cov[7] *= eps2l; 
462   }
463
464   Set(x,alpha,par,cov);
465
466   if (!Invariant()) return kFALSE;
467
468   return kTRUE;
469 }
470
471 void AliITStrackV2::CookdEdx(Double_t low, Double_t up) {
472   //-----------------------------------------------------------------
473   // This function calculates dE/dX within the "low" and "up" cuts.
474   // Origin: Boris Batyunya, JINR, Boris.Batiounia@cern.ch 
475   // Updated: F. Prino 8-June-2009
476   //-----------------------------------------------------------------
477   // The cluster order is: SDD-1, SDD-2, SSD-1, SSD-2
478
479   Int_t nc=0;
480   Float_t dedx[4];
481   for (Int_t il=0; il<4; il++) { // count good (>0) dE/dx values
482     if(fdEdxSample[il]>0.){
483       dedx[nc]= fdEdxSample[il];
484       nc++;
485     }
486   }
487   if(nc<1){
488     SetdEdx(0.);
489     return;
490   }
491
492   Int_t swap; // sort in ascending order
493   do {
494     swap=0;
495     for (Int_t i=0; i<nc-1; i++) {
496       if (dedx[i]<=dedx[i+1]) continue;
497       Float_t tmp=dedx[i];
498       dedx[i]=dedx[i+1]; 
499       dedx[i+1]=tmp;
500       swap++;
501     }
502   } while (swap);
503
504
505   Double_t nl=low*nc, nu =up*nc;
506   Float_t sumamp = 0;
507   Float_t sumweight =0;
508   for (Int_t i=0; i<nc; i++) {
509     Float_t weight =1;
510     if (i<nl+0.1)     weight = TMath::Max(1.-(nl-i),0.);
511     if (i>nu-1)     weight = TMath::Max(nu-i,0.);
512     sumamp+= dedx[i]*weight;
513     sumweight+=weight;
514   }
515   SetdEdx(sumamp/sumweight);
516 }
517
518 //____________________________________________________________________________
519 Bool_t AliITStrackV2::
520 GetPhiZat(Double_t r, Double_t &phi, Double_t &z) const {
521   //------------------------------------------------------------------
522   // This function returns the global cylindrical (phi,z) of the track 
523   // position estimated at the radius r. 
524   // The track curvature is neglected.
525   //------------------------------------------------------------------
526   Double_t d=GetD(0.,0.);
527   if (TMath::Abs(d) > r) {
528     if (r>1e-1) return kFALSE;
529     r = TMath::Abs(d);
530   }
531
532   Double_t rcurr=TMath::Sqrt(GetX()*GetX() + GetY()*GetY());
533   if (TMath::Abs(d) > rcurr) return kFALSE;
534   Double_t globXYZcurr[3]; GetXYZ(globXYZcurr); 
535   Double_t phicurr=TMath::ATan2(globXYZcurr[1],globXYZcurr[0]);
536
537   if (GetX()>=0.) {
538     phi=phicurr+TMath::ASin(d/r)-TMath::ASin(d/rcurr);
539   } else {
540     phi=phicurr+TMath::ASin(d/r)+TMath::ASin(d/rcurr)-TMath::Pi();
541   }
542
543   // return a phi in [0,2pi[ 
544   if (phi<0.) phi+=2.*TMath::Pi();
545   else if (phi>=2.*TMath::Pi()) phi-=2.*TMath::Pi();
546   z=GetZ()+GetTgl()*(TMath::Sqrt((r-d)*(r+d))-TMath::Sqrt((rcurr-d)*(rcurr+d)));
547   return kTRUE;
548 }
549 //____________________________________________________________________________
550 Bool_t AliITStrackV2::
551 GetLocalXat(Double_t r,Double_t &xloc) const {
552   //------------------------------------------------------------------
553   // This function returns the local x of the track 
554   // position estimated at the radius r. 
555   // The track curvature is neglected.
556   //------------------------------------------------------------------
557   Double_t d=GetD(0.,0.);
558   if (TMath::Abs(d) > r) { 
559     if (r>1e-1) return kFALSE; 
560     r = TMath::Abs(d); 
561   } 
562
563   Double_t rcurr=TMath::Sqrt(GetX()*GetX() + GetY()*GetY());
564   Double_t globXYZcurr[3]; GetXYZ(globXYZcurr); 
565   Double_t phicurr=TMath::ATan2(globXYZcurr[1],globXYZcurr[0]);
566   Double_t phi;
567   if (GetX()>=0.) {
568     phi=phicurr+TMath::ASin(d/r)-TMath::ASin(d/rcurr);
569   } else {
570     phi=phicurr+TMath::ASin(d/r)+TMath::ASin(d/rcurr)-TMath::Pi();
571   }
572
573   xloc=r*(TMath::Cos(phi)*TMath::Cos(GetAlpha())
574          +TMath::Sin(phi)*TMath::Sin(GetAlpha())); 
575
576   return kTRUE;
577 }