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