]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONTrackExtrap.cxx
Adding CreateIterator(void) and GetNeighbours() pure virtual methods,
[u/mrichter/AliRoot.git] / MUON / AliMUONTrackExtrap.cxx
CommitLineData
c04e3238 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//
20// Tools
21// for
22// track
23// extrapolation
24// in
25// ALICE
26// dimuon
27// spectrometer
28//
29///////////////////////////////////////////////////
30
31#include <Riostream.h>
208f139e 32#include <TMatrixD.h>
c04e3238 33
34#include "AliMUONTrackExtrap.h"
35#include "AliMUONTrackParam.h"
36#include "AliMUONConstants.h"
37#include "AliMagF.h"
38#include "AliLog.h"
39#include "AliTracker.h"
40
41ClassImp(AliMUONTrackExtrap) // Class implementation in ROOT context
42
43const AliMagF* AliMUONTrackExtrap::fgkField = 0x0;
208f139e 44const Int_t AliMUONTrackExtrap::fgkMaxStepNumber = 5000;
45const Double_t AliMUONTrackExtrap::fgkStepLength = 6.;
46
47 //__________________________________________________________________________
48Double_t AliMUONTrackExtrap::GetImpactParamFromBendingMomentum(Double_t bendingMomentum)
49{
50 /// Returns impact parameter at vertex in bending plane (cm),
51 /// from the signed bending momentum "BendingMomentum" in bending plane (GeV/c),
52 /// using simple values for dipole magnetic field.
53 /// The sign of "BendingMomentum" is the sign of the charge.
54
55 if (bendingMomentum == 0.) return 1.e10;
56
57 Double_t simpleBPosition = 0.5 * (AliMUONConstants::CoilZ() + AliMUONConstants::YokeZ());
58 Double_t simpleBLength = 0.5 * (AliMUONConstants::CoilL() + AliMUONConstants::YokeL());
59 Float_t b[3], x[3] = {0.,0.,(Float_t) simpleBPosition};
60 if (fgkField) fgkField->Field(x,b);
61 else {
62 cout<<"F-AliMUONTrackExtrap::GetField: fgkField = 0x0"<<endl;
63 exit(-1);
64 }
65 Double_t simpleBValue = (Double_t) b[0];
66
67 return (-0.0003 * simpleBValue * simpleBLength * simpleBPosition / bendingMomentum);
68}
69
70 //__________________________________________________________________________
71Double_t AliMUONTrackExtrap::GetBendingMomentumFromImpactParam(Double_t impactParam)
72{
73 /// Returns signed bending momentum in bending plane (GeV/c),
74 /// the sign being the sign of the charge for particles moving forward in Z,
75 /// from the impact parameter "ImpactParam" at vertex in bending plane (cm),
76 /// using simple values for dipole magnetic field.
77
78 if (impactParam == 0.) return 1.e10;
79
80 Double_t simpleBPosition = 0.5 * (AliMUONConstants::CoilZ() + AliMUONConstants::YokeZ());
81 Double_t simpleBLength = 0.5 * (AliMUONConstants::CoilL() + AliMUONConstants::YokeL());
82 Float_t b[3], x[3] = {0.,0.,(Float_t) simpleBPosition};
83 if (fgkField) fgkField->Field(x,b);
84 else {
85 cout<<"F-AliMUONTrackExtrap::GetField: fgkField = 0x0"<<endl;
86 exit(-1);
87 }
88 Double_t simpleBValue = (Double_t) b[0];
89
90 return (-0.0003 * simpleBValue * simpleBLength * simpleBPosition / impactParam);
91}
c04e3238 92
93 //__________________________________________________________________________
94void AliMUONTrackExtrap::ExtrapToZ(AliMUONTrackParam* trackParam, Double_t zEnd)
95{
96 /// Track parameter extrapolation to the plane at "Z".
97 /// On return, the track parameters resulting from the extrapolation are updated in trackParam.
98 if (trackParam->GetZ() == zEnd) return; // nothing to be done if same Z
99 Double_t forwardBackward; // +1 if forward, -1 if backward
100 if (zEnd < trackParam->GetZ()) forwardBackward = 1.0; // spectro. z<0
101 else forwardBackward = -1.0;
dade8580 102 Double_t v3[7], v3New[7]; // 7 in parameter ????
103 Int_t i3, stepNumber;
c04e3238 104 // For safety: return kTRUE or kFALSE ????
105 // Parameter vector for calling EXTRAP_ONESTEP
dade8580 106 ConvertTrackParamForExtrap(trackParam, v3, forwardBackward);
c04e3238 107 // sign of charge (sign of fInverseBendingMomentum if forward motion)
108 // must be changed if backward extrapolation
208f139e 109 Double_t chargeExtrap = forwardBackward * TMath::Sign(Double_t(1.0), trackParam->GetInverseBendingMomentum());
c04e3238 110 // Extrapolation loop
111 stepNumber = 0;
208f139e 112 while (((-forwardBackward * (v3[2] - zEnd)) <= 0.0) && (stepNumber < fgkMaxStepNumber)) { // spectro. z<0
c04e3238 113 stepNumber++;
114 // Option for switching between helix and Runge-Kutta ????
208f139e 115 //ExtrapOneStepRungekutta(chargeExtrap, fgkStepLength, v3, v3New);
116 ExtrapOneStepHelix(chargeExtrap, fgkStepLength, v3, v3New);
dade8580 117 if ((-forwardBackward * (v3New[2] - zEnd)) > 0.0) break; // one is beyond Z spectro. z<0
c04e3238 118 // better use TArray ????
208f139e 119 for (i3 = 0; i3 < 7; i3++) {v3[i3] = v3New[i3];}
c04e3238 120 }
208f139e 121 // check fgkMaxStepNumber ????
c04e3238 122 // Interpolation back to exact Z (2nd order)
123 // should be in function ???? using TArray ????
dade8580 124 Double_t dZ12 = v3New[2] - v3[2]; // 1->2
c04e3238 125 if (TMath::Abs(dZ12) > 0) {
dade8580 126 Double_t dZ1i = zEnd - v3[2]; // 1-i
127 Double_t dZi2 = v3New[2] - zEnd; // i->2
128 Double_t xPrime = (v3New[0] - v3[0]) / dZ12;
129 Double_t xSecond = ((v3New[3] / v3New[5]) - (v3[3] / v3[5])) / dZ12;
130 Double_t yPrime = (v3New[1] - v3[1]) / dZ12;
131 Double_t ySecond = ((v3New[4] / v3New[5]) - (v3[4] / v3[5])) / dZ12;
132 v3[0] = v3[0] + xPrime * dZ1i - 0.5 * xSecond * dZ1i * dZi2; // X
133 v3[1] = v3[1] + yPrime * dZ1i - 0.5 * ySecond * dZ1i * dZi2; // Y
134 v3[2] = zEnd; // Z
c04e3238 135 Double_t xPrimeI = xPrime - 0.5 * xSecond * (dZi2 - dZ1i);
136 Double_t yPrimeI = yPrime - 0.5 * ySecond * (dZi2 - dZ1i);
137 // (PX, PY, PZ)/PTOT assuming forward motion
208f139e 138 v3[5] = 1.0 / TMath::Sqrt(1.0 + xPrimeI * xPrimeI + yPrimeI * yPrimeI); // PZ/PTOT
dade8580 139 v3[3] = xPrimeI * v3[5]; // PX/PTOT
140 v3[4] = yPrimeI * v3[5]; // PY/PTOT
c04e3238 141 } else {
142 cout<<"W-AliMUONTrackExtrap::ExtrapToZ: Extrap. to Z not reached, Z = "<<zEnd<<endl;
143 }
dade8580 144 // Track parameters from 3 parameters,
c04e3238 145 // with charge back for forward motion
dade8580 146 RecoverTrackParam(v3, chargeExtrap * forwardBackward, trackParam);
c04e3238 147}
148
149 //__________________________________________________________________________
dade8580 150void AliMUONTrackExtrap::ConvertTrackParamForExtrap(AliMUONTrackParam* trackParam, Double_t *v3, Double_t forwardBackward)
c04e3238 151{
dade8580 152 /// Set vector of Geant3 parameters pointed to by "v3" from track parameters in trackParam.
c04e3238 153 /// Since AliMUONTrackParam is only geometry, one uses "forwardBackward"
154 /// to know whether the particle is going forward (+1) or backward (-1).
dade8580 155 v3[0] = trackParam->GetNonBendingCoor(); // X
156 v3[1] = trackParam->GetBendingCoor(); // Y
157 v3[2] = trackParam->GetZ(); // Z
c04e3238 158 Double_t pYZ = TMath::Abs(1.0 / trackParam->GetInverseBendingMomentum());
159 Double_t pZ = pYZ / TMath::Sqrt(1.0 + trackParam->GetBendingSlope() * trackParam->GetBendingSlope());
dade8580 160 v3[6] = TMath::Sqrt(pYZ * pYZ + pZ * pZ * trackParam->GetNonBendingSlope() * trackParam->GetNonBendingSlope()); // PTOT
161 v3[5] = -forwardBackward * pZ / v3[6]; // PZ/PTOT spectro. z<0
162 v3[3] = trackParam->GetNonBendingSlope() * v3[5]; // PX/PTOT
163 v3[4] = trackParam->GetBendingSlope() * v3[5]; // PY/PTOT
c04e3238 164}
165
166 //__________________________________________________________________________
dade8580 167void AliMUONTrackExtrap::RecoverTrackParam(Double_t *v3, Double_t charge, AliMUONTrackParam* trackParam)
c04e3238 168{
dade8580 169 /// Set track parameters in trackParam from Geant3 parameters pointed to by "v3",
c04e3238 170 /// assumed to be calculated for forward motion in Z.
171 /// "InverseBendingMomentum" is signed with "charge".
dade8580 172 trackParam->SetNonBendingCoor(v3[0]); // X
173 trackParam->SetBendingCoor(v3[1]); // Y
174 trackParam->SetZ(v3[2]); // Z
175 Double_t pYZ = v3[6] * TMath::Sqrt(1.0 - v3[3] * v3[3]);
c04e3238 176 trackParam->SetInverseBendingMomentum(charge/pYZ);
dade8580 177 trackParam->SetBendingSlope(v3[4]/v3[5]);
178 trackParam->SetNonBendingSlope(v3[3]/v3[5]);
208f139e 179}
180
181 //__________________________________________________________________________
182void AliMUONTrackExtrap::ExtrapToZCov(AliMUONTrackParam* trackParam, Double_t zEnd)
183{
184 /// Track parameters and their covariances extrapolated to the plane at "zEnd".
185 /// On return, results from the extrapolation are updated in trackParam.
186
187 if (trackParam->GetZ() == zEnd) return; // nothing to be done if same z
188
189 // Save the actual track parameters
190 AliMUONTrackParam trackParamSave(*trackParam);
191 Double_t nonBendingCoor = trackParamSave.GetNonBendingCoor();
192 Double_t nonBendingSlope = trackParamSave.GetNonBendingSlope();
193 Double_t bendingCoor = trackParamSave.GetBendingCoor();
194 Double_t bendingSlope = trackParamSave.GetBendingSlope();
195 Double_t inverseBendingMomentum = trackParamSave.GetInverseBendingMomentum();
196 Double_t zBegin = trackParamSave.GetZ();
197
198 // Extrapolate track parameters to "zEnd"
199 ExtrapToZ(trackParam,zEnd);
200 Double_t extrapNonBendingCoor = trackParam->GetNonBendingCoor();
201 Double_t extrapNonBendingSlope = trackParam->GetNonBendingSlope();
202 Double_t extrapBendingCoor = trackParam->GetBendingCoor();
203 Double_t extrapBendingSlope = trackParam->GetBendingSlope();
204 Double_t extrapInverseBendingMomentum = trackParam->GetInverseBendingMomentum();
205
206 // Get the pointer to the parameter covariance matrix
207 if (!trackParam->CovariancesExist()) {
208 //cout<<"W-AliMUONTrackExtrap::ExtrapToZCov: track parameter covariance matrix does not exist"<<endl;
209 //cout<<" -> nothing to extrapolate !!"<<endl;
210 return;
211 }
212 TMatrixD* paramCov = trackParam->GetCovariances();
213
214 // Calculate the jacobian related to the track parameters extrapolation to "zEnd"
215 TMatrixD jacob(5,5);
216 jacob = 0.;
217 Double_t dParam[5];
218 for (Int_t i=0; i<5; i++) {
219 // Skip jacobian calculation for parameters with no associated error
220 if ((*paramCov)(i,i) == 0.) continue;
221 // Small variation of parameter i only
222 for (Int_t j=0; j<5; j++) {
223 if (j==i) {
224 dParam[j] = TMath::Sqrt((*paramCov)(i,i));
225 if (j == 4) dParam[j] *= TMath::Sign(1.,-inverseBendingMomentum); // variation always in the same direction
226 } else dParam[j] = 0.;
227 }
228 // Set new parameters
229 trackParamSave.SetNonBendingCoor (nonBendingCoor + dParam[0]);
230 trackParamSave.SetNonBendingSlope (nonBendingSlope + dParam[1]);
231 trackParamSave.SetBendingCoor (bendingCoor + dParam[2]);
232 trackParamSave.SetBendingSlope (bendingSlope + dParam[3]);
233 trackParamSave.SetInverseBendingMomentum(inverseBendingMomentum + dParam[4]);
234 trackParamSave.SetZ (zBegin);
235 // Extrapolate new track parameters to "zEnd"
236 ExtrapToZ(&trackParamSave,zEnd);
237 // Calculate the jacobian
238 jacob(0,i) = (trackParamSave.GetNonBendingCoor() - extrapNonBendingCoor ) / dParam[i];
239 jacob(1,i) = (trackParamSave.GetNonBendingSlope() - extrapNonBendingSlope ) / dParam[i];
240 jacob(2,i) = (trackParamSave.GetBendingCoor() - extrapBendingCoor ) / dParam[i];
241 jacob(3,i) = (trackParamSave.GetBendingSlope() - extrapBendingSlope ) / dParam[i];
242 jacob(4,i) = (trackParamSave.GetInverseBendingMomentum() - extrapInverseBendingMomentum) / dParam[i];
243 }
244
245 // Extrapolate track parameter covariances to "zEnd"
246 TMatrixD tmp((*paramCov),TMatrixD::kMultTranspose,jacob);
247 (*paramCov) = TMatrixD(jacob,TMatrixD::kMult,tmp);
248
c04e3238 249}
250
251 //__________________________________________________________________________
252void AliMUONTrackExtrap::ExtrapToStation(AliMUONTrackParam* trackParamIn, Int_t station, AliMUONTrackParam *trackParamOut)
253{
254 /// Track parameters extrapolated from "trackParamIn" to both chambers of the station(0..) "station"
255 /// are returned in the array (dimension 2) of track parameters pointed to by "TrackParamOut"
256 /// (index 0 and 1 for first and second chambers).
257 Double_t extZ[2], z1, z2;
258 Int_t i1 = -1, i2 = -1; // = -1 to avoid compilation warnings
259 // range of station to be checked ????
260 z1 = AliMUONConstants::DefaultChamberZ(2 * station);
261 z2 = AliMUONConstants::DefaultChamberZ(2 * station + 1);
262 // First and second Z to extrapolate at
263 if ((z1 > trackParamIn->GetZ()) && (z2 > trackParamIn->GetZ())) {i1 = 0; i2 = 1;}
264 else if ((z1 < trackParamIn->GetZ()) && (z2 < trackParamIn->GetZ())) {i1 = 1; i2 = 0;}
265 else {
208f139e 266 cout<<"E-AliMUONTrackExtrap::ExtrapToStation: Starting Z ("<<trackParamIn->GetZ()
c04e3238 267 <<") in between z1 ("<<z1<<") and z2 ("<<z2<<") of station(0..)"<<station<<endl;
268 exit(-1);
269 }
270 extZ[i1] = z1;
271 extZ[i2] = z2;
272 // copy of track parameters
273 trackParamOut[i1] = *trackParamIn;
274 // first extrapolation
275 ExtrapToZ(&(trackParamOut[i1]),extZ[0]);
276 trackParamOut[i2] = trackParamOut[i1];
277 // second extrapolation
278 ExtrapToZ(&(trackParamOut[i2]),extZ[1]);
279 return;
208f139e 280}
281
282 //__________________________________________________________________________
283void AliMUONTrackExtrap::ExtrapToVertexUncorrected(AliMUONTrackParam* trackParam, Double_t zVtx)
284{
285 /// Extrapolation to the vertex (at the z position "zVtx") without Branson and Field correction.
286 /// Returns the track parameters resulting from the extrapolation in the current TrackParam.
287 /// Include multiple Coulomb scattering effects in trackParam covariances.
288
289 if (trackParam->GetZ() > zVtx) { // spectro. (z<0)
290 cout<<"W-AliMUONTrackExtrap::ExtrapToVertexUncorrected: Starting Z ("<<trackParam->GetZ()
291 <<") upstream the vertex (zVtx = "<<zVtx<<")"<<endl;
292 exit(-1);
293 }
294
295 if (zVtx < AliMUONConstants::ZAbsorberEnd()) { // spectro. (z<0)
296 cout<<"W-AliMUONTrackExtrap::ExtrapToVertexUncorrected: Ending Z ("<<zVtx
297 <<") downstream the front absorber (zAbsorberEnd = "<<AliMUONConstants::ZAbsorberEnd()<<")"<<endl;
298
299 ExtrapToZCov(trackParam,zVtx);
300 return;
301 }
302
303 // First Extrapolates track parameters upstream to the "Z" end of the front absorber
304 if (trackParam->GetZ() < AliMUONConstants::ZAbsorberEnd()) { // spectro. (z<0)
305 ExtrapToZCov(trackParam,AliMUONConstants::ZAbsorberEnd());
306 } else {
307 cout<<"W-AliMUONTrackExtrap::ExtrapToVertexUncorrected: Starting Z ("<<trackParam->GetZ()
308 <<") upstream or inside the front absorber (zAbsorberEnd = "<<AliMUONConstants::ZAbsorberEnd()<<")"<<endl;
309 }
310
311 // Then go through all the absorber layers
312 Double_t tan3 = TMath::Tan(3./180.*TMath::Pi());
313 Double_t r0Norm, x0, z, zElement, dZ, nonBendingCoor, bendingCoor;
314 for (Int_t iElement=AliMUONConstants::NAbsorberElements()-1; iElement>=0; iElement--) {
315 zElement = AliMUONConstants::ZAbsorberElement(iElement);
316 z = trackParam->GetZ();
317 if (z > zElement) continue; // spectro. (z<0)
318 nonBendingCoor = trackParam->GetNonBendingCoor();
319 bendingCoor = trackParam->GetBendingCoor();
320 r0Norm = nonBendingCoor * nonBendingCoor + bendingCoor * bendingCoor;
321 r0Norm = TMath::Sqrt(r0Norm) / TMath::Abs(trackParam->GetZ()) / tan3;
322 if (r0Norm > 1.) x0 = AliMUONConstants::X0AbsorberOut(iElement); // outer part of the absorber
323 else x0 = AliMUONConstants::X0AbsorberIn(iElement); // inner part of the absorber
324
325 if (zVtx > zElement) {
326 ExtrapToZCov(trackParam,zElement); // extrapolate to absorber element "iElement"
327 dZ = zElement - z;
328 AddMCSEffectInTrackParamCov(trackParam,dZ,x0); // include MCS effect in covariances
329 } else {
330 ExtrapToZCov(trackParam,zVtx); // extrapolate to zVtx
331 dZ = zVtx - z;
332 AddMCSEffectInTrackParamCov(trackParam,dZ,x0); // include MCS effect in covariances
333 break;
334 }
335 }
336
337 // finally go to the vertex
338 ExtrapToZCov(trackParam,zVtx);
339
340}
341
342 //__________________________________________________________________________
343void AliMUONTrackExtrap::AddMCSEffectInTrackParamCov(AliMUONTrackParam *param, Double_t dZ, Double_t x0)
344{
345 /// Add to the track parameter covariances the effects of multiple Coulomb scattering
346 /// through a material of thickness "dZ" and of radiation length "x0"
347 /// assuming linear propagation and using the small angle approximation.
348
349 Double_t bendingSlope = param->GetBendingSlope();
350 Double_t nonBendingSlope = param->GetNonBendingSlope();
351 Double_t inverseTotalMomentum2 = param->GetInverseBendingMomentum() * param->GetInverseBendingMomentum() *
352 (1.0 + bendingSlope * bendingSlope) /
353 (1.0 + bendingSlope *bendingSlope + nonBendingSlope * nonBendingSlope);
354 // Path length in the material
355 Double_t pathLength = TMath::Abs(dZ) * TMath::Sqrt(1.0 + bendingSlope*bendingSlope + nonBendingSlope*nonBendingSlope);
356 Double_t pathLength2 = pathLength * pathLength;
357 // relativistic velocity
358 Double_t velo = 1.;
359 // Angular dispersion square of the track (variance) in a plane perpendicular to the trajectory
360 Double_t theta02 = 0.0136 / velo * (1 + 0.038 * TMath::Log(pathLength/x0));
361 theta02 *= theta02 * inverseTotalMomentum2 * pathLength / x0;
362
363 // Add effects of multiple Coulomb scattering in track parameter covariances
364 TMatrixD* paramCov = param->GetCovariances();
365 Double_t varCoor = pathLength2 * theta02 / 3.;
366 Double_t varSlop = theta02;
367 Double_t covCorrSlope = pathLength * theta02 / 2.;
368 // Non bending plane
369 (*paramCov)(0,0) += varCoor; (*paramCov)(0,1) += covCorrSlope;
370 (*paramCov)(1,0) += covCorrSlope; (*paramCov)(1,1) += varSlop;
371 // Bending plane
372 (*paramCov)(2,2) += varCoor; (*paramCov)(2,3) += covCorrSlope;
373 (*paramCov)(3,2) += covCorrSlope; (*paramCov)(3,3) += varSlop;
374
c04e3238 375}
376
377 //__________________________________________________________________________
378void AliMUONTrackExtrap::ExtrapToVertex(AliMUONTrackParam* trackParam, Double_t xVtx, Double_t yVtx, Double_t zVtx)
379{
380 /// Extrapolation to the vertex.
381 /// Returns the track parameters resulting from the extrapolation in the current TrackParam.
382 /// Changes parameters according to Branson correction through the absorber
383
c04e3238 384 // Extrapolates track parameters upstream to the "Z" end of the front absorber
208f139e 385 ExtrapToZ(trackParam,AliMUONConstants::ZAbsorberEnd()); // !!!
c04e3238 386 // Makes Branson correction (multiple scattering + energy loss)
387 BransonCorrection(trackParam,xVtx,yVtx,zVtx);
388 // Makes a simple magnetic field correction through the absorber
208f139e 389 FieldCorrection(trackParam,AliMUONConstants::ZAbsorberEnd());
c04e3238 390}
391
392
393// Keep this version for future developments
394 //__________________________________________________________________________
395// void AliMUONTrackExtrap::BransonCorrection(AliMUONTrackParam* trackParam)
396// {
397// /// Branson correction of track parameters
398// // the entry parameters have to be calculated at the end of the absorber
399// Double_t zEndAbsorber, zBP, xBP, yBP;
400// Double_t pYZ, pX, pY, pZ, pTotal, xEndAbsorber, yEndAbsorber, radiusEndAbsorber2, pT, theta;
401// Int_t sign;
402// // Would it be possible to calculate all that from Geant configuration ????
403// // and to get the Branson parameters from a function in ABSO module ????
404// // with an eventual contribution from other detectors like START ????
405// // Radiation lengths outer part theta > 3 degres
406// static Double_t x01[9] = { 18.8, // C (cm)
407// 10.397, // Concrete (cm)
408// 0.56, // Plomb (cm)
409// 47.26, // Polyethylene (cm)
410// 0.56, // Plomb (cm)
411// 47.26, // Polyethylene (cm)
412// 0.56, // Plomb (cm)
413// 47.26, // Polyethylene (cm)
414// 0.56 }; // Plomb (cm)
415// // inner part theta < 3 degres
416// static Double_t x02[3] = { 18.8, // C (cm)
417// 10.397, // Concrete (cm)
418// 0.35 }; // W (cm)
419// // z positions of the materials inside the absober outer part theta > 3 degres
420// static Double_t z1[10] = { 90, 315, 467, 472, 477, 482, 487, 492, 497, 502 };
421// // inner part theta < 3 degres
422// static Double_t z2[4] = { 90, 315, 467, 503 };
423// static Bool_t first = kTRUE;
424// static Double_t zBP1, zBP2, rLimit;
425// // Calculates z positions of the Branson's planes: zBP1 for outer part and zBP2 for inner part (only at the first call)
426// if (first) {
427// first = kFALSE;
428// Double_t aNBP = 0.0;
429// Double_t aDBP = 0.0;
430// Int_t iBound;
431//
432// for (iBound = 0; iBound < 9; iBound++) {
433// aNBP = aNBP +
434// (z1[iBound+1] * z1[iBound+1] * z1[iBound+1] -
435// z1[iBound] * z1[iBound] * z1[iBound] ) / x01[iBound];
436// aDBP = aDBP +
437// (z1[iBound+1] * z1[iBound+1] - z1[iBound] * z1[iBound] ) / x01[iBound];
438// }
439// zBP1 = (2.0 * aNBP) / (3.0 * aDBP);
440// aNBP = 0.0;
441// aDBP = 0.0;
442// for (iBound = 0; iBound < 3; iBound++) {
443// aNBP = aNBP +
444// (z2[iBound+1] * z2[iBound+1] * z2[iBound+1] -
445// z2[iBound] * z2[iBound ] * z2[iBound] ) / x02[iBound];
446// aDBP = aDBP +
447// (z2[iBound+1] * z2[iBound+1] - z2[iBound] * z2[iBound]) / x02[iBound];
448// }
449// zBP2 = (2.0 * aNBP) / (3.0 * aDBP);
450// rLimit = z2[3] * TMath::Tan(3.0 * (TMath::Pi()) / 180.);
451// }
452//
453// pYZ = TMath::Abs(1.0 / trackParam->GetInverseBendingMomentum());
454// sign = 1;
455// if (trackParam->GetInverseBendingMomentum() < 0) sign = -1;
456// pZ = pYZ / (TMath::Sqrt(1.0 + trackParam->GetBendingSlope() * trackParam->GetBendingSlope()));
457// pX = pZ * trackParam->GetNonBendingSlope();
458// pY = pZ * trackParam->GetBendingSlope();
459// pTotal = TMath::Sqrt(pYZ *pYZ + pX * pX);
460// xEndAbsorber = trackParam->GetNonBendingCoor();
461// yEndAbsorber = trackParam->GetBendingCoor();
462// radiusEndAbsorber2 = xEndAbsorber * xEndAbsorber + yEndAbsorber * yEndAbsorber;
463//
464// if (radiusEndAbsorber2 > rLimit*rLimit) {
465// zEndAbsorber = z1[9];
466// zBP = zBP1;
467// } else {
468// zEndAbsorber = z2[3];
469// zBP = zBP2;
470// }
471//
472// xBP = xEndAbsorber - (pX / pZ) * (zEndAbsorber - zBP);
473// yBP = yEndAbsorber - (pY / pZ) * (zEndAbsorber - zBP);
474//
475// // new parameters after Branson and energy loss corrections
476// pZ = pTotal * zBP / TMath::Sqrt(xBP * xBP + yBP * yBP + zBP * zBP);
477// pX = pZ * xBP / zBP;
478// pY = pZ * yBP / zBP;
479// trackParam->SetBendingSlope(pY/pZ);
480// trackParam->SetNonBendingSlope(pX/pZ);
481//
482// pT = TMath::Sqrt(pX * pX + pY * pY);
483// theta = TMath::ATan2(pT, pZ);
484// pTotal = TotalMomentumEnergyLoss(rLimit, pTotal, theta, xEndAbsorber, yEndAbsorber);
485//
486// trackParam->SetInverseBendingMomentum((sign / pTotal) *
487// TMath::Sqrt(1.0 +
488// trackParam->GetBendingSlope() * trackParam->GetBendingSlope() +
489// trackParam->GetNonBendingSlope() * trackParam->GetNonBendingSlope()) /
490// TMath::Sqrt(1.0 + trackParam->GetBendingSlope() * trackParam->GetBendingSlope()));
491//
492// // vertex position at (0,0,0)
493// // should be taken from vertex measurement ???
494// trackParam->SetBendingCoor(0.);
495// trackParam->SetNonBendingCoor(0.);
496// trackParam->SetZ(0.);
497// }
498
499void AliMUONTrackExtrap::BransonCorrection(AliMUONTrackParam* trackParam, Double_t xVtx, Double_t yVtx, Double_t zVtx)
500{
501 /// Branson correction of track parameters
502 // the entry parameters have to be calculated at the end of the absorber
503 // simplified version: the z positions of Branson's planes are no longer calculated
504 // but are given as inputs. One can use the macros MUONTestAbso.C and DrawTestAbso.C
505 // to test this correction.
506 // Would it be possible to calculate all that from Geant configuration ????
507 // and to get the Branson parameters from a function in ABSO module ????
508 // with an eventual contribution from other detectors like START ????
509 // change to take into account the vertex postition (real, reconstruct,....)
510
511 Double_t zBP, xBP, yBP;
512 Double_t pYZ, pX, pY, pZ, pTotal, xEndAbsorber, yEndAbsorber, radiusEndAbsorber2, pT, theta;
513 Int_t sign;
514 static Bool_t first = kTRUE;
208f139e 515 static Double_t zBP1, zBP2, rLimit, thetaLimit;
c04e3238 516 // zBP1 for outer part and zBP2 for inner part (only at the first call)
517 if (first) {
518 first = kFALSE;
519
c04e3238 520 thetaLimit = 3.0 * (TMath::Pi()) / 180.;
208f139e 521 rLimit = TMath::Abs(AliMUONConstants::ZAbsorberEnd()) * TMath::Tan(thetaLimit);
c04e3238 522 zBP1 = -450; // values close to those calculated with EvalAbso.C
523 zBP2 = -480;
524 }
525
526 pYZ = TMath::Abs(1.0 / trackParam->GetInverseBendingMomentum());
527 sign = 1;
528 if (trackParam->GetInverseBendingMomentum() < 0) sign = -1;
529 pZ = trackParam->Pz();
530 pX = trackParam->Px();
531 pY = trackParam->Py();
532 pTotal = TMath::Sqrt(pYZ *pYZ + pX * pX);
533 xEndAbsorber = trackParam->GetNonBendingCoor();
534 yEndAbsorber = trackParam->GetBendingCoor();
535 radiusEndAbsorber2 = xEndAbsorber * xEndAbsorber + yEndAbsorber * yEndAbsorber;
536
537 if (radiusEndAbsorber2 > rLimit*rLimit) {
538 zBP = zBP1;
539 } else {
540 zBP = zBP2;
541 }
542
208f139e 543 xBP = xEndAbsorber - (pX / pZ) * (AliMUONConstants::ZAbsorberEnd() - zBP);
544 yBP = yEndAbsorber - (pY / pZ) * (AliMUONConstants::ZAbsorberEnd() - zBP);
c04e3238 545
546 // new parameters after Branson and energy loss corrections
547// Float_t zSmear = zBP - gRandom->Gaus(0.,2.); // !!! possible smearing of Z vertex position
548
549 Float_t zSmear = zBP ;
550
551 pZ = pTotal * (zSmear-zVtx) / TMath::Sqrt((xBP-xVtx) * (xBP-xVtx) + (yBP-yVtx) * (yBP-yVtx) +( zSmear-zVtx) * (zSmear-zVtx) );
552 pX = pZ * (xBP - xVtx)/ (zSmear-zVtx);
553 pY = pZ * (yBP - yVtx) / (zSmear-zVtx);
554 trackParam->SetBendingSlope(pY/pZ);
555 trackParam->SetNonBendingSlope(pX/pZ);
556
557
558 pT = TMath::Sqrt(pX * pX + pY * pY);
559 theta = TMath::ATan2(pT, TMath::Abs(pZ));
560 pTotal = TotalMomentumEnergyLoss(thetaLimit, pTotal, theta);
561
562 trackParam->SetInverseBendingMomentum((sign / pTotal) *
563 TMath::Sqrt(1.0 +
564 trackParam->GetBendingSlope() * trackParam->GetBendingSlope() +
565 trackParam->GetNonBendingSlope() * trackParam->GetNonBendingSlope()) /
566 TMath::Sqrt(1.0 + trackParam->GetBendingSlope() * trackParam->GetBendingSlope()));
567
568 // vertex position at (0,0,0)
569 // should be taken from vertex measurement ???
570
571 trackParam->SetBendingCoor(xVtx);
572 trackParam->SetNonBendingCoor(yVtx);
573 trackParam->SetZ(zVtx);
574
575}
576
577 //__________________________________________________________________________
578Double_t AliMUONTrackExtrap::TotalMomentumEnergyLoss(Double_t thetaLimit, Double_t pTotal, Double_t theta)
579{
580 /// Returns the total momentum corrected from energy loss in the front absorber
581 // One can use the macros MUONTestAbso.C and DrawTestAbso.C
582 // to test this correction.
583 // Momentum energy loss behaviour evaluated with the simulation of single muons (april 2002)
584 Double_t deltaP, pTotalCorrected;
585
586 // Parametrization to be redone according to change of absorber material ????
587 // See remark in function BransonCorrection !!!!
588 // The name is not so good, and there are many arguments !!!!
589 if (theta < thetaLimit ) {
590 if (pTotal < 20) {
591 deltaP = 2.5938 + 0.0570 * pTotal - 0.001151 * pTotal * pTotal;
592 } else {
593 deltaP = 3.0714 + 0.011767 *pTotal;
594 }
595 deltaP *= 0.75; // AZ
596 } else {
597 if (pTotal < 20) {
598 deltaP = 2.1207 + 0.05478 * pTotal - 0.00145079 * pTotal * pTotal;
599 } else {
600 deltaP = 2.6069 + 0.0051705 * pTotal;
601 }
602 deltaP *= 0.9; // AZ
603 }
604 pTotalCorrected = pTotal + deltaP / TMath::Cos(theta);
605 return pTotalCorrected;
606}
607
608 //__________________________________________________________________________
609void AliMUONTrackExtrap::FieldCorrection(AliMUONTrackParam *trackParam, Double_t zEnd)
610{
611 /// Correction of the effect of the magnetic field in the absorber
612 // Assume a constant field along Z axis.
613 Float_t b[3],x[3];
614 Double_t bZ;
615 Double_t pYZ,pX,pY,pZ,pT;
616 Double_t pXNew,pYNew;
617 Double_t c;
618
619 pYZ = TMath::Abs(1.0 / trackParam->GetInverseBendingMomentum());
620 c = TMath::Sign(1.0,trackParam->GetInverseBendingMomentum()); // particle charge
621
622 pZ = trackParam->Pz();
623 pX = trackParam->Px();
624 pY = trackParam->Py();
625 pT = TMath::Sqrt(pX*pX+pY*pY);
626
627 if (TMath::Abs(pZ) <= 0) return;
628 x[2] = zEnd/2;
629 x[0] = x[2]*trackParam->GetNonBendingSlope();
630 x[1] = x[2]*trackParam->GetBendingSlope();
631
632 // Take magn. field value at position x.
633 if (fgkField) fgkField->Field(x,b);
634 else {
635 cout<<"F-AliMUONTrackExtrap::FieldCorrection: fgkField = 0x0"<<endl;
636 exit(-1);
637 }
638 bZ = b[2];
639
640 // Transverse momentum rotation
641 // Parameterized with the study of DeltaPhi = phiReco - phiGen as a function of pZ.
642 Double_t phiShift = c*0.436*0.0003*bZ*zEnd/pZ;
643 // Rotate momentum around Z axis.
644 pXNew = pX*TMath::Cos(phiShift) - pY*TMath::Sin(phiShift);
645 pYNew = pX*TMath::Sin(phiShift) + pY*TMath::Cos(phiShift);
646
647 trackParam->SetBendingSlope(pYNew/pZ);
648 trackParam->SetNonBendingSlope(pXNew/pZ);
649
650 trackParam->SetInverseBendingMomentum(c/TMath::Sqrt(pYNew*pYNew+pZ*pZ));
651
652}
653
654 //__________________________________________________________________________
655void AliMUONTrackExtrap::ExtrapOneStepHelix(Double_t charge, Double_t step, Double_t *vect, Double_t *vout)
656{
657/// ******************************************************************
658/// * *
659/// * Performs the tracking of one step in a magnetic field *
660/// * The trajectory is assumed to be a helix in a constant field *
661/// * taken at the mid point of the step. *
662/// * Parameters: *
663/// * input *
664/// * STEP =arc length of the step asked *
665/// * VECT =input vector (position,direction cos and momentum) *
666/// * CHARGE= electric charge of the particle *
667/// * output *
668/// * VOUT = same as VECT after completion of the step *
669/// * *
670/// * ==>Called by : <USER>, GUSWIM *
671/// * Author m.hansroul ********* *
672/// * modified s.egli, s.v.levonian *
673/// * modified v.perevoztchikov
674/// * *
675/// ******************************************************************
676
677// modif: everything in double precision
678
679 Double_t xyz[3], h[4], hxp[3];
680 Double_t h2xy, hp, rho, tet;
681 Double_t sint, sintt, tsint, cos1t;
682 Double_t f1, f2, f3, f4, f5, f6;
683
684 const Int_t kix = 0;
685 const Int_t kiy = 1;
686 const Int_t kiz = 2;
687 const Int_t kipx = 3;
688 const Int_t kipy = 4;
689 const Int_t kipz = 5;
690 const Int_t kipp = 6;
691
692 const Double_t kec = 2.9979251e-4;
693 //
694 // ------------------------------------------------------------------
695 //
696 // units are kgauss,centimeters,gev/c
697 //
698 vout[kipp] = vect[kipp];
699 if (TMath::Abs(charge) < 0.00001) {
700 for (Int_t i = 0; i < 3; i++) {
701 vout[i] = vect[i] + step * vect[i+3];
702 vout[i+3] = vect[i+3];
703 }
704 return;
705 }
706 xyz[0] = vect[kix] + 0.5 * step * vect[kipx];
707 xyz[1] = vect[kiy] + 0.5 * step * vect[kipy];
708 xyz[2] = vect[kiz] + 0.5 * step * vect[kipz];
709
710 //cmodif: call gufld (xyz, h) changed into:
711 GetField (xyz, h);
712
713 h2xy = h[0]*h[0] + h[1]*h[1];
714 h[3] = h[2]*h[2]+ h2xy;
715 if (h[3] < 1.e-12) {
716 for (Int_t i = 0; i < 3; i++) {
717 vout[i] = vect[i] + step * vect[i+3];
718 vout[i+3] = vect[i+3];
719 }
720 return;
721 }
722 if (h2xy < 1.e-12*h[3]) {
723 ExtrapOneStepHelix3(charge*h[2], step, vect, vout);
724 return;
725 }
726 h[3] = TMath::Sqrt(h[3]);
727 h[0] /= h[3];
728 h[1] /= h[3];
729 h[2] /= h[3];
730 h[3] *= kec;
731
732 hxp[0] = h[1]*vect[kipz] - h[2]*vect[kipy];
733 hxp[1] = h[2]*vect[kipx] - h[0]*vect[kipz];
734 hxp[2] = h[0]*vect[kipy] - h[1]*vect[kipx];
735
736 hp = h[0]*vect[kipx] + h[1]*vect[kipy] + h[2]*vect[kipz];
737
738 rho = -charge*h[3]/vect[kipp];
739 tet = rho * step;
740
741 if (TMath::Abs(tet) > 0.15) {
742 sint = TMath::Sin(tet);
743 sintt = (sint/tet);
744 tsint = (tet-sint)/tet;
745 cos1t = 2.*(TMath::Sin(0.5*tet))*(TMath::Sin(0.5*tet))/tet;
746 } else {
747 tsint = tet*tet/36.;
748 sintt = (1. - tsint);
749 sint = tet*sintt;
750 cos1t = 0.5*tet;
751 }
752
753 f1 = step * sintt;
754 f2 = step * cos1t;
755 f3 = step * tsint * hp;
756 f4 = -tet*cos1t;
757 f5 = sint;
758 f6 = tet * cos1t * hp;
759
760 vout[kix] = vect[kix] + f1*vect[kipx] + f2*hxp[0] + f3*h[0];
761 vout[kiy] = vect[kiy] + f1*vect[kipy] + f2*hxp[1] + f3*h[1];
762 vout[kiz] = vect[kiz] + f1*vect[kipz] + f2*hxp[2] + f3*h[2];
763
764 vout[kipx] = vect[kipx] + f4*vect[kipx] + f5*hxp[0] + f6*h[0];
765 vout[kipy] = vect[kipy] + f4*vect[kipy] + f5*hxp[1] + f6*h[1];
766 vout[kipz] = vect[kipz] + f4*vect[kipz] + f5*hxp[2] + f6*h[2];
767
768 return;
769}
770
771 //__________________________________________________________________________
772void AliMUONTrackExtrap::ExtrapOneStepHelix3(Double_t field, Double_t step, Double_t *vect, Double_t *vout)
773{
774/// ******************************************************************
775/// * *
776/// * Tracking routine in a constant field oriented *
777/// * along axis 3 *
778/// * Tracking is performed with a conventional *
779/// * helix step method *
780/// * *
781/// * ==>Called by : <USER>, GUSWIM *
782/// * Authors R.Brun, M.Hansroul ********* *
783/// * Rewritten V.Perevoztchikov
784/// * *
785/// ******************************************************************
786
787 Double_t hxp[3];
788 Double_t h4, hp, rho, tet;
789 Double_t sint, sintt, tsint, cos1t;
790 Double_t f1, f2, f3, f4, f5, f6;
791
792 const Int_t kix = 0;
793 const Int_t kiy = 1;
794 const Int_t kiz = 2;
795 const Int_t kipx = 3;
796 const Int_t kipy = 4;
797 const Int_t kipz = 5;
798 const Int_t kipp = 6;
799
800 const Double_t kec = 2.9979251e-4;
801
802//
803// ------------------------------------------------------------------
804//
805// units are kgauss,centimeters,gev/c
806//
807 vout[kipp] = vect[kipp];
808 h4 = field * kec;
809
810 hxp[0] = - vect[kipy];
811 hxp[1] = + vect[kipx];
812
813 hp = vect[kipz];
814
815 rho = -h4/vect[kipp];
816 tet = rho * step;
817 if (TMath::Abs(tet) > 0.15) {
818 sint = TMath::Sin(tet);
819 sintt = (sint/tet);
820 tsint = (tet-sint)/tet;
821 cos1t = 2.* TMath::Sin(0.5*tet) * TMath::Sin(0.5*tet)/tet;
822 } else {
823 tsint = tet*tet/36.;
824 sintt = (1. - tsint);
825 sint = tet*sintt;
826 cos1t = 0.5*tet;
827 }
828
829 f1 = step * sintt;
830 f2 = step * cos1t;
831 f3 = step * tsint * hp;
832 f4 = -tet*cos1t;
833 f5 = sint;
834 f6 = tet * cos1t * hp;
835
836 vout[kix] = vect[kix] + f1*vect[kipx] + f2*hxp[0];
837 vout[kiy] = vect[kiy] + f1*vect[kipy] + f2*hxp[1];
838 vout[kiz] = vect[kiz] + f1*vect[kipz] + f3;
839
840 vout[kipx] = vect[kipx] + f4*vect[kipx] + f5*hxp[0];
841 vout[kipy] = vect[kipy] + f4*vect[kipy] + f5*hxp[1];
842 vout[kipz] = vect[kipz] + f4*vect[kipz] + f6;
843
844 return;
845}
846 //__________________________________________________________________________
847void AliMUONTrackExtrap::ExtrapOneStepRungekutta(Double_t charge, Double_t step, Double_t* vect, Double_t* vout)
848{
849/// ******************************************************************
850/// * *
851/// * Runge-Kutta method for tracking a particle through a magnetic *
852/// * field. Uses Nystroem algorithm (See Handbook Nat. Bur. of *
853/// * Standards, procedure 25.5.20) *
854/// * *
855/// * Input parameters *
856/// * CHARGE Particle charge *
857/// * STEP Step size *
858/// * VECT Initial co-ords,direction cosines,momentum *
859/// * Output parameters *
860/// * VOUT Output co-ords,direction cosines,momentum *
861/// * User routine called *
862/// * CALL GUFLD(X,F) *
863/// * *
864/// * ==>Called by : <USER>, GUSWIM *
865/// * Authors R.Brun, M.Hansroul ********* *
866/// * V.Perevoztchikov (CUT STEP implementation) *
867/// * *
868/// * *
869/// ******************************************************************
870
871 Double_t h2, h4, f[4];
872 Double_t xyzt[3], a, b, c, ph,ph2;
873 Double_t secxs[4],secys[4],seczs[4],hxp[3];
874 Double_t g1, g2, g3, g4, g5, g6, ang2, dxt, dyt, dzt;
875 Double_t est, at, bt, ct, cba;
876 Double_t f1, f2, f3, f4, rho, tet, hnorm, hp, rho1, sint, cost;
877
878 Double_t x;
879 Double_t y;
880 Double_t z;
881
882 Double_t xt;
883 Double_t yt;
884 Double_t zt;
885
886 Double_t maxit = 1992;
887 Double_t maxcut = 11;
888
889 const Double_t kdlt = 1e-4;
890 const Double_t kdlt32 = kdlt/32.;
891 const Double_t kthird = 1./3.;
892 const Double_t khalf = 0.5;
893 const Double_t kec = 2.9979251e-4;
894
895 const Double_t kpisqua = 9.86960440109;
896 const Int_t kix = 0;
897 const Int_t kiy = 1;
898 const Int_t kiz = 2;
899 const Int_t kipx = 3;
900 const Int_t kipy = 4;
901 const Int_t kipz = 5;
902
903 // *.
904 // *. ------------------------------------------------------------------
905 // *.
906 // * this constant is for units cm,gev/c and kgauss
907 // *
908 Int_t iter = 0;
909 Int_t ncut = 0;
910 for(Int_t j = 0; j < 7; j++)
911 vout[j] = vect[j];
912
913 Double_t pinv = kec * charge / vect[6];
914 Double_t tl = 0.;
915 Double_t h = step;
916 Double_t rest;
917
918
919 do {
920 rest = step - tl;
921 if (TMath::Abs(h) > TMath::Abs(rest)) h = rest;
922 //cmodif: call gufld(vout,f) changed into:
923
924 GetField(vout,f);
925
926 // *
927 // * start of integration
928 // *
929 x = vout[0];
930 y = vout[1];
931 z = vout[2];
932 a = vout[3];
933 b = vout[4];
934 c = vout[5];
935
936 h2 = khalf * h;
937 h4 = khalf * h2;
938 ph = pinv * h;
939 ph2 = khalf * ph;
940 secxs[0] = (b * f[2] - c * f[1]) * ph2;
941 secys[0] = (c * f[0] - a * f[2]) * ph2;
942 seczs[0] = (a * f[1] - b * f[0]) * ph2;
943 ang2 = (secxs[0]*secxs[0] + secys[0]*secys[0] + seczs[0]*seczs[0]);
944 if (ang2 > kpisqua) break;
945
946 dxt = h2 * a + h4 * secxs[0];
947 dyt = h2 * b + h4 * secys[0];
948 dzt = h2 * c + h4 * seczs[0];
949 xt = x + dxt;
950 yt = y + dyt;
951 zt = z + dzt;
952 // *
953 // * second intermediate point
954 // *
955
956 est = TMath::Abs(dxt) + TMath::Abs(dyt) + TMath::Abs(dzt);
957 if (est > h) {
958 if (ncut++ > maxcut) break;
959 h *= khalf;
960 continue;
961 }
962
963 xyzt[0] = xt;
964 xyzt[1] = yt;
965 xyzt[2] = zt;
966
967 //cmodif: call gufld(xyzt,f) changed into:
968 GetField(xyzt,f);
969
970 at = a + secxs[0];
971 bt = b + secys[0];
972 ct = c + seczs[0];
973
974 secxs[1] = (bt * f[2] - ct * f[1]) * ph2;
975 secys[1] = (ct * f[0] - at * f[2]) * ph2;
976 seczs[1] = (at * f[1] - bt * f[0]) * ph2;
977 at = a + secxs[1];
978 bt = b + secys[1];
979 ct = c + seczs[1];
980 secxs[2] = (bt * f[2] - ct * f[1]) * ph2;
981 secys[2] = (ct * f[0] - at * f[2]) * ph2;
982 seczs[2] = (at * f[1] - bt * f[0]) * ph2;
983 dxt = h * (a + secxs[2]);
984 dyt = h * (b + secys[2]);
985 dzt = h * (c + seczs[2]);
986 xt = x + dxt;
987 yt = y + dyt;
988 zt = z + dzt;
989 at = a + 2.*secxs[2];
990 bt = b + 2.*secys[2];
991 ct = c + 2.*seczs[2];
992
993 est = TMath::Abs(dxt)+TMath::Abs(dyt)+TMath::Abs(dzt);
994 if (est > 2.*TMath::Abs(h)) {
995 if (ncut++ > maxcut) break;
996 h *= khalf;
997 continue;
998 }
999
1000 xyzt[0] = xt;
1001 xyzt[1] = yt;
1002 xyzt[2] = zt;
1003
1004 //cmodif: call gufld(xyzt,f) changed into:
1005 GetField(xyzt,f);
1006
1007 z = z + (c + (seczs[0] + seczs[1] + seczs[2]) * kthird) * h;
1008 y = y + (b + (secys[0] + secys[1] + secys[2]) * kthird) * h;
1009 x = x + (a + (secxs[0] + secxs[1] + secxs[2]) * kthird) * h;
1010
1011 secxs[3] = (bt*f[2] - ct*f[1])* ph2;
1012 secys[3] = (ct*f[0] - at*f[2])* ph2;
1013 seczs[3] = (at*f[1] - bt*f[0])* ph2;
1014 a = a+(secxs[0]+secxs[3]+2. * (secxs[1]+secxs[2])) * kthird;
1015 b = b+(secys[0]+secys[3]+2. * (secys[1]+secys[2])) * kthird;
1016 c = c+(seczs[0]+seczs[3]+2. * (seczs[1]+seczs[2])) * kthird;
1017
1018 est = TMath::Abs(secxs[0]+secxs[3] - (secxs[1]+secxs[2]))
1019 + TMath::Abs(secys[0]+secys[3] - (secys[1]+secys[2]))
1020 + TMath::Abs(seczs[0]+seczs[3] - (seczs[1]+seczs[2]));
1021
1022 if (est > kdlt && TMath::Abs(h) > 1.e-4) {
1023 if (ncut++ > maxcut) break;
1024 h *= khalf;
1025 continue;
1026 }
1027
1028 ncut = 0;
1029 // * if too many iterations, go to helix
1030 if (iter++ > maxit) break;
1031
1032 tl += h;
1033 if (est < kdlt32)
1034 h *= 2.;
1035 cba = 1./ TMath::Sqrt(a*a + b*b + c*c);
1036 vout[0] = x;
1037 vout[1] = y;
1038 vout[2] = z;
1039 vout[3] = cba*a;
1040 vout[4] = cba*b;
1041 vout[5] = cba*c;
1042 rest = step - tl;
1043 if (step < 0.) rest = -rest;
1044 if (rest < 1.e-5*TMath::Abs(step)) return;
1045
1046 } while(1);
1047
1048 // angle too big, use helix
1049
1050 f1 = f[0];
1051 f2 = f[1];
1052 f3 = f[2];
1053 f4 = TMath::Sqrt(f1*f1+f2*f2+f3*f3);
1054 rho = -f4*pinv;
1055 tet = rho * step;
1056
1057 hnorm = 1./f4;
1058 f1 = f1*hnorm;
1059 f2 = f2*hnorm;
1060 f3 = f3*hnorm;
1061
1062 hxp[0] = f2*vect[kipz] - f3*vect[kipy];
1063 hxp[1] = f3*vect[kipx] - f1*vect[kipz];
1064 hxp[2] = f1*vect[kipy] - f2*vect[kipx];
1065
1066 hp = f1*vect[kipx] + f2*vect[kipy] + f3*vect[kipz];
1067
1068 rho1 = 1./rho;
1069 sint = TMath::Sin(tet);
1070 cost = 2.*TMath::Sin(khalf*tet)*TMath::Sin(khalf*tet);
1071
1072 g1 = sint*rho1;
1073 g2 = cost*rho1;
1074 g3 = (tet-sint) * hp*rho1;
1075 g4 = -cost;
1076 g5 = sint;
1077 g6 = cost * hp;
1078
1079 vout[kix] = vect[kix] + g1*vect[kipx] + g2*hxp[0] + g3*f1;
1080 vout[kiy] = vect[kiy] + g1*vect[kipy] + g2*hxp[1] + g3*f2;
1081 vout[kiz] = vect[kiz] + g1*vect[kipz] + g2*hxp[2] + g3*f3;
1082
1083 vout[kipx] = vect[kipx] + g4*vect[kipx] + g5*hxp[0] + g6*f1;
1084 vout[kipy] = vect[kipy] + g4*vect[kipy] + g5*hxp[1] + g6*f2;
1085 vout[kipz] = vect[kipz] + g4*vect[kipz] + g5*hxp[2] + g6*f3;
1086
1087 return;
1088}
1089//___________________________________________________________
1090 void AliMUONTrackExtrap::GetField(Double_t *Position, Double_t *Field)
1091{
1092 /// interface for arguments in double precision (Why ? ChF)
1093 Float_t x[3], b[3];
1094
1095 x[0] = Position[0]; x[1] = Position[1]; x[2] = Position[2];
1096
1097 if (fgkField) fgkField->Field(x,b);
1098 else {
1099 cout<<"F-AliMUONTrackExtrap::GetField: fgkField = 0x0"<<endl;
1100 exit(-1);
1101 }
1102
1103 Field[0] = b[0]; Field[1] = b[1]; Field[2] = b[2];
1104
1105 return;
1106}
1107