]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONTrackParam.cxx
Bug in AliMUONHit constructor for fPz and fPx calculation (Arthur, Bruce and Ivana)
[u/mrichter/AliRoot.git] / MUON / AliMUONTrackParam.cxx
CommitLineData
a9e2aefa 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
88cb7938 16/* $Id$ */
a9e2aefa 17
3831f268 18///////////////////////////////////////////////////
19//
20// Track parameters
21// in
22// ALICE
23// dimuon
24// spectrometer
a9e2aefa 25//
3831f268 26///////////////////////////////////////////////////
a9e2aefa 27
70479d0e 28#include <Riostream.h>
a9e2aefa 29
30#include "AliCallf77.h"
3831f268 31#include "AliMUON.h"
a9e2aefa 32#include "AliMUONTrackParam.h"
3831f268 33#include "AliMUONChamber.h"
a9e2aefa 34#include "AliRun.h"
94de3818 35#include "AliMagF.h"
a9e2aefa 36
37ClassImp(AliMUONTrackParam) // Class implementation in ROOT context
38
a6f03ddb 39 // A few calls in Fortran or from Fortran (extrap.F).
40 // Needed, instead of calls to Geant subroutines,
41 // because double precision is necessary for track fit converging with Minuit.
42 // The "extrap" functions should be translated into C++ ????
a9e2aefa 43#ifndef WIN32
a6f03ddb 44# define extrap_onestep_helix extrap_onestep_helix_
45# define extrap_onestep_helix3 extrap_onestep_helix3_
46# define extrap_onestep_rungekutta extrap_onestep_rungekutta_
47# define gufld_double gufld_double_
a9e2aefa 48#else
a6f03ddb 49# define extrap_onestep_helix EXTRAP_ONESTEP_HELIX
50# define extrap_onestep_helix3 EXTRAP_ONESTEP_HELIX3
51# define extrap_onestep_rungekutta EXTRAP_ONESTEP_RUNGEKUTTA
52# define gufld_double GUFLD_DOUBLE
a9e2aefa 53#endif
54
a6f03ddb 55extern "C" {
56 void type_of_call extrap_onestep_helix
57 (Double_t &Charge, Double_t &StepLength, Double_t *VGeant3, Double_t *VGeant3New);
58
59 void type_of_call extrap_onestep_helix3
60 (Double_t &Field, Double_t &StepLength, Double_t *VGeant3, Double_t *VGeant3New);
61
62 void type_of_call extrap_onestep_rungekutta
63 (Double_t &Charge, Double_t &StepLength, Double_t *VGeant3, Double_t *VGeant3New);
64
65 void type_of_call gufld_double(Double_t *Position, Double_t *Field) {
66 // interface to "gAlice->Field()->Field" for arguments in double precision
67 Float_t x[3], b[3];
68 x[0] = Position[0]; x[1] = Position[1]; x[2] = Position[2];
69 gAlice->Field()->Field(x, b);
70 Field[0] = b[0]; Field[1] = b[1]; Field[2] = b[2];
71 }
a9e2aefa 72}
30178c30 73
61adb9bd 74 //_________________________________________________________________________
30178c30 75AliMUONTrackParam::AliMUONTrackParam()
76 : TObject()
77{
78// Constructor
79
80 fInverseBendingMomentum = 0;
81 fBendingSlope = 0;
82 fNonBendingSlope = 0;
83 fZ = 0;
84 fBendingCoor = 0;
85 fNonBendingCoor = 0;
86}
61adb9bd 87
30178c30 88 //_________________________________________________________________________
89AliMUONTrackParam&
90AliMUONTrackParam::operator=(const AliMUONTrackParam& theMUONTrackParam)
61adb9bd 91{
30178c30 92 if (this == &theMUONTrackParam)
61adb9bd 93 return *this;
94
30178c30 95 // base class assignement
96 TObject::operator=(theMUONTrackParam);
97
98 fInverseBendingMomentum = theMUONTrackParam.fInverseBendingMomentum;
99 fBendingSlope = theMUONTrackParam.fBendingSlope;
100 fNonBendingSlope = theMUONTrackParam.fNonBendingSlope;
101 fZ = theMUONTrackParam.fZ;
102 fBendingCoor = theMUONTrackParam.fBendingCoor;
103 fNonBendingCoor = theMUONTrackParam.fNonBendingCoor;
61adb9bd 104
105 return *this;
106}
107 //_________________________________________________________________________
30178c30 108AliMUONTrackParam::AliMUONTrackParam(const AliMUONTrackParam& theMUONTrackParam)
109 : TObject(theMUONTrackParam)
61adb9bd 110{
30178c30 111 fInverseBendingMomentum = theMUONTrackParam.fInverseBendingMomentum;
112 fBendingSlope = theMUONTrackParam.fBendingSlope;
113 fNonBendingSlope = theMUONTrackParam.fNonBendingSlope;
114 fZ = theMUONTrackParam.fZ;
115 fBendingCoor = theMUONTrackParam.fBendingCoor;
116 fNonBendingCoor = theMUONTrackParam.fNonBendingCoor;
61adb9bd 117}
a9e2aefa 118
a9e2aefa 119 //__________________________________________________________________________
120void AliMUONTrackParam::ExtrapToZ(Double_t Z)
121{
122 // Track parameter extrapolation to the plane at "Z".
123 // On return, the track parameters resulting from the extrapolation
124 // replace the current track parameters.
a9e2aefa 125 if (this->fZ == Z) return; // nothing to be done if same Z
126 Double_t forwardBackward; // +1 if forward, -1 if backward
5b64e914 127 if (Z < this->fZ) forwardBackward = 1.0; // spectro. z<0
a9e2aefa 128 else forwardBackward = -1.0;
a6f03ddb 129 Double_t vGeant3[7], vGeant3New[7]; // 7 in parameter ????
a9e2aefa 130 Int_t iGeant3, stepNumber;
131 Int_t maxStepNumber = 5000; // in parameter ????
132 // For safety: return kTRUE or kFALSE ????
a6f03ddb 133 // Parameter vector for calling EXTRAP_ONESTEP
a9e2aefa 134 SetGeant3Parameters(vGeant3, forwardBackward);
956019b6 135 // sign of charge (sign of fInverseBendingMomentum if forward motion)
a6f03ddb 136 // must be changed if backward extrapolation
956019b6 137 Double_t chargeExtrap = forwardBackward *
138 TMath::Sign(Double_t(1.0), this->fInverseBendingMomentum);
a9e2aefa 139 Double_t stepLength = 6.0; // in parameter ????
140 // Extrapolation loop
141 stepNumber = 0;
5b64e914 142 while (((-forwardBackward * (vGeant3[2] - Z)) <= 0.0) && // spectro. z<0
a9e2aefa 143 (stepNumber < maxStepNumber)) {
144 stepNumber++;
a6f03ddb 145 // Option for switching between helix and Runge-Kutta ????
146 // extrap_onestep_rungekutta(chargeExtrap, stepLength, vGeant3, vGeant3New);
147 extrap_onestep_helix(chargeExtrap, stepLength, vGeant3, vGeant3New);
5b64e914 148 if ((-forwardBackward * (vGeant3New[2] - Z)) > 0.0) break; // one is beyond Z spectro. z<0
a9e2aefa 149 // better use TArray ????
150 for (iGeant3 = 0; iGeant3 < 7; iGeant3++)
151 {vGeant3[iGeant3] = vGeant3New[iGeant3];}
152 }
153 // check maxStepNumber ????
a9e2aefa 154 // Interpolation back to exact Z (2nd order)
155 // should be in function ???? using TArray ????
156 Double_t dZ12 = vGeant3New[2] - vGeant3[2]; // 1->2
157 Double_t dZ1i = Z - vGeant3[2]; // 1-i
158 Double_t dZi2 = vGeant3New[2] - Z; // i->2
159 Double_t xPrime = (vGeant3New[0] - vGeant3[0]) / dZ12;
160 Double_t xSecond =
161 ((vGeant3New[3] / vGeant3New[5]) - (vGeant3[3] / vGeant3[5])) / dZ12;
162 Double_t yPrime = (vGeant3New[1] - vGeant3[1]) / dZ12;
163 Double_t ySecond =
164 ((vGeant3New[4] / vGeant3New[5]) - (vGeant3[4] / vGeant3[5])) / dZ12;
165 vGeant3[0] = vGeant3[0] + xPrime * dZ1i - 0.5 * xSecond * dZ1i * dZi2; // X
166 vGeant3[1] = vGeant3[1] + yPrime * dZ1i - 0.5 * ySecond * dZ1i * dZi2; // Y
167 vGeant3[2] = Z; // Z
168 Double_t xPrimeI = xPrime - 0.5 * xSecond * (dZi2 - dZ1i);
169 Double_t yPrimeI = yPrime - 0.5 * ySecond * (dZi2 - dZ1i);
956019b6 170 // (PX, PY, PZ)/PTOT assuming forward motion
a9e2aefa 171 vGeant3[5] =
172 1.0 / TMath::Sqrt(1.0 + xPrimeI * xPrimeI + yPrimeI * yPrimeI); // PZ/PTOT
173 vGeant3[3] = xPrimeI * vGeant3[5]; // PX/PTOT
174 vGeant3[4] = yPrimeI * vGeant3[5]; // PY/PTOT
956019b6 175 // Track parameters from Geant3 parameters,
176 // with charge back for forward motion
177 GetFromGeant3Parameters(vGeant3, chargeExtrap * forwardBackward);
a9e2aefa 178}
179
180 //__________________________________________________________________________
181void AliMUONTrackParam::SetGeant3Parameters(Double_t *VGeant3, Double_t ForwardBackward)
182{
183 // Set vector of Geant3 parameters pointed to by "VGeant3"
184 // from track parameters in current AliMUONTrackParam.
185 // Since AliMUONTrackParam is only geometry, one uses "ForwardBackward"
186 // to know whether the particle is going forward (+1) or backward (-1).
187 VGeant3[0] = this->fNonBendingCoor; // X
188 VGeant3[1] = this->fBendingCoor; // Y
189 VGeant3[2] = this->fZ; // Z
190 Double_t pYZ = TMath::Abs(1.0 / this->fInverseBendingMomentum);
191 Double_t pZ =
192 pYZ / TMath::Sqrt(1.0 + this->fBendingSlope * this->fBendingSlope);
193 VGeant3[6] =
194 TMath::Sqrt(pYZ * pYZ +
195 pZ * pZ * this->fNonBendingSlope * this->fNonBendingSlope); // PTOT
5b64e914 196 VGeant3[5] = -ForwardBackward * pZ / VGeant3[6]; // PZ/PTOT spectro. z<0
a9e2aefa 197 VGeant3[3] = this->fNonBendingSlope * VGeant3[5]; // PX/PTOT
198 VGeant3[4] = this->fBendingSlope * VGeant3[5]; // PY/PTOT
199}
200
201 //__________________________________________________________________________
202void AliMUONTrackParam::GetFromGeant3Parameters(Double_t *VGeant3, Double_t Charge)
203{
204 // Get track parameters in current AliMUONTrackParam
956019b6 205 // from Geant3 parameters pointed to by "VGeant3",
206 // assumed to be calculated for forward motion in Z.
a9e2aefa 207 // "InverseBendingMomentum" is signed with "Charge".
208 this->fNonBendingCoor = VGeant3[0]; // X
209 this->fBendingCoor = VGeant3[1]; // Y
210 this->fZ = VGeant3[2]; // Z
211 Double_t pYZ = VGeant3[6] * TMath::Sqrt(1.0 - VGeant3[3] * VGeant3[3]);
212 this->fInverseBendingMomentum = Charge / pYZ;
213 this->fBendingSlope = VGeant3[4] / VGeant3[5];
214 this->fNonBendingSlope = VGeant3[3] / VGeant3[5];
215}
216
217 //__________________________________________________________________________
218void AliMUONTrackParam::ExtrapToStation(Int_t Station, AliMUONTrackParam *TrackParam)
219{
220 // Track parameters extrapolated from current track parameters ("this")
221 // to both chambers of the station(0..) "Station"
222 // are returned in the array (dimension 2) of track parameters
223 // pointed to by "TrackParam" (index 0 and 1 for first and second chambers).
224 Double_t extZ[2], z1, z2;
ecfa008b 225 Int_t i1 = -1, i2 = -1; // = -1 to avoid compilation warnings
a9e2aefa 226 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
227 // range of Station to be checked ????
228 z1 = (&(pMUON->Chamber(2 * Station)))->Z(); // Z of first chamber
229 z2 = (&(pMUON->Chamber(2 * Station + 1)))->Z(); // Z of second chamber
230 // First and second Z to extrapolate at
231 if ((z1 > this->fZ) && (z2 > this->fZ)) {i1 = 0; i2 = 1;}
232 else if ((z1 < this->fZ) && (z2 < this->fZ)) {i1 = 1; i2 = 0;}
233 else {
234 cout << "ERROR in AliMUONTrackParam::CreateExtrapSegmentInStation" << endl;
235 cout << "Starting Z (" << this->fZ << ") in between z1 (" << z1 <<
236 ") and z2 (" << z2 << ") of station(0..) " << Station << endl;
237 }
238 extZ[i1] = z1;
239 extZ[i2] = z2;
240 // copy of track parameters
241 TrackParam[i1] = *this;
242 // first extrapolation
243 (&(TrackParam[i1]))->ExtrapToZ(extZ[0]);
244 TrackParam[i2] = TrackParam[i1];
245 // second extrapolation
246 (&(TrackParam[i2]))->ExtrapToZ(extZ[1]);
247 return;
248}
249
04b5ea16 250 //__________________________________________________________________________
251void AliMUONTrackParam::ExtrapToVertex()
252{
253 // Extrapolation to the vertex.
254 // Returns the track parameters resulting from the extrapolation,
255 // in the current TrackParam.
956019b6 256 // Changes parameters according to Branson correction through the absorber
04b5ea16 257
5b64e914 258 Double_t zAbsorber = -503.0; // to be coherent with the Geant absorber geometry !!!!
259 // spectro. (z<0)
04b5ea16 260 // Extrapolates track parameters upstream to the "Z" end of the front absorber
b45fd22b 261 ExtrapToZ(zAbsorber); // !!!
5b64e914 262 // Makes Branson correction (multiple scattering + energy loss)
04b5ea16 263 BransonCorrection();
5b64e914 264 // Makes a simple magnetic field correction through the absorber
b45fd22b 265 FieldCorrection(zAbsorber);
04b5ea16 266}
267
43af2cb6 268
269// Keep this version for future developments
04b5ea16 270 //__________________________________________________________________________
43af2cb6 271// void AliMUONTrackParam::BransonCorrection()
272// {
273// // Branson correction of track parameters
274// // the entry parameters have to be calculated at the end of the absorber
275// Double_t zEndAbsorber, zBP, xBP, yBP;
276// Double_t pYZ, pX, pY, pZ, pTotal, xEndAbsorber, yEndAbsorber, radiusEndAbsorber2, pT, theta;
277// Int_t sign;
278// // Would it be possible to calculate all that from Geant configuration ????
279// // and to get the Branson parameters from a function in ABSO module ????
280// // with an eventual contribution from other detectors like START ????
281// // Radiation lengths outer part theta > 3 degres
282// static Double_t x01[9] = { 18.8, // C (cm)
283// 10.397, // Concrete (cm)
284// 0.56, // Plomb (cm)
285// 47.26, // Polyethylene (cm)
286// 0.56, // Plomb (cm)
287// 47.26, // Polyethylene (cm)
288// 0.56, // Plomb (cm)
289// 47.26, // Polyethylene (cm)
290// 0.56 }; // Plomb (cm)
291// // inner part theta < 3 degres
292// static Double_t x02[3] = { 18.8, // C (cm)
293// 10.397, // Concrete (cm)
294// 0.35 }; // W (cm)
295// // z positions of the materials inside the absober outer part theta > 3 degres
296// static Double_t z1[10] = { 90, 315, 467, 472, 477, 482, 487, 492, 497, 502 };
297// // inner part theta < 3 degres
298// static Double_t z2[4] = { 90, 315, 467, 503 };
299// static Bool_t first = kTRUE;
300// static Double_t zBP1, zBP2, rLimit;
301// // Calculates z positions of the Branson's planes: zBP1 for outer part and zBP2 for inner part (only at the first call)
302// if (first) {
303// first = kFALSE;
304// Double_t aNBP = 0.0;
305// Double_t aDBP = 0.0;
306// Int_t iBound;
307
308// for (iBound = 0; iBound < 9; iBound++) {
309// aNBP = aNBP +
310// (z1[iBound+1] * z1[iBound+1] * z1[iBound+1] -
311// z1[iBound] * z1[iBound] * z1[iBound] ) / x01[iBound];
312// aDBP = aDBP +
313// (z1[iBound+1] * z1[iBound+1] - z1[iBound] * z1[iBound] ) / x01[iBound];
314// }
315// zBP1 = (2.0 * aNBP) / (3.0 * aDBP);
316// aNBP = 0.0;
317// aDBP = 0.0;
318// for (iBound = 0; iBound < 3; iBound++) {
319// aNBP = aNBP +
320// (z2[iBound+1] * z2[iBound+1] * z2[iBound+1] -
321// z2[iBound] * z2[iBound ] * z2[iBound] ) / x02[iBound];
322// aDBP = aDBP +
323// (z2[iBound+1] * z2[iBound+1] - z2[iBound] * z2[iBound]) / x02[iBound];
324// }
325// zBP2 = (2.0 * aNBP) / (3.0 * aDBP);
326// rLimit = z2[3] * TMath::Tan(3.0 * (TMath::Pi()) / 180.);
327// }
328
329// pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
330// sign = 1;
331// if (fInverseBendingMomentum < 0) sign = -1;
332// pZ = pYZ / (TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope));
333// pX = pZ * fNonBendingSlope;
334// pY = pZ * fBendingSlope;
335// pTotal = TMath::Sqrt(pYZ *pYZ + pX * pX);
336// xEndAbsorber = fNonBendingCoor;
337// yEndAbsorber = fBendingCoor;
338// radiusEndAbsorber2 = xEndAbsorber * xEndAbsorber + yEndAbsorber * yEndAbsorber;
339
340// if (radiusEndAbsorber2 > rLimit*rLimit) {
341// zEndAbsorber = z1[9];
342// zBP = zBP1;
343// } else {
344// zEndAbsorber = z2[3];
345// zBP = zBP2;
346// }
347
348// xBP = xEndAbsorber - (pX / pZ) * (zEndAbsorber - zBP);
349// yBP = yEndAbsorber - (pY / pZ) * (zEndAbsorber - zBP);
350
351// // new parameters after Branson and energy loss corrections
352// pZ = pTotal * zBP / TMath::Sqrt(xBP * xBP + yBP * yBP + zBP * zBP);
353// pX = pZ * xBP / zBP;
354// pY = pZ * yBP / zBP;
355// fBendingSlope = pY / pZ;
356// fNonBendingSlope = pX / pZ;
357
358// pT = TMath::Sqrt(pX * pX + pY * pY);
359// theta = TMath::ATan2(pT, pZ);
360// pTotal =
361// TotalMomentumEnergyLoss(rLimit, pTotal, theta, xEndAbsorber, yEndAbsorber);
362
363// fInverseBendingMomentum = (sign / pTotal) *
364// TMath::Sqrt(1.0 +
365// fBendingSlope * fBendingSlope +
366// fNonBendingSlope * fNonBendingSlope) /
367// TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope);
368
369// // vertex position at (0,0,0)
370// // should be taken from vertex measurement ???
371// fBendingCoor = 0.0;
372// fNonBendingCoor = 0;
373// fZ= 0;
374// }
375
04b5ea16 376void AliMUONTrackParam::BransonCorrection()
377{
378 // Branson correction of track parameters
379 // the entry parameters have to be calculated at the end of the absorber
43af2cb6 380 // simplified version: the z positions of Branson's planes are no longer calculated
381 // but are given as inputs. One can use the macros MUONTestAbso.C and DrawTestAbso.C
382 // to test this correction.
04b5ea16 383 // Would it be possible to calculate all that from Geant configuration ????
956019b6 384 // and to get the Branson parameters from a function in ABSO module ????
385 // with an eventual contribution from other detectors like START ????
43af2cb6 386 Double_t zBP, xBP, yBP;
387 Double_t pYZ, pX, pY, pZ, pTotal, xEndAbsorber, yEndAbsorber, radiusEndAbsorber2, pT, theta;
388 Int_t sign;
04b5ea16 389 static Bool_t first = kTRUE;
b45fd22b 390 static Double_t zBP1, zBP2, rLimit, thetaLimit, zEndAbsorber;
43af2cb6 391 // zBP1 for outer part and zBP2 for inner part (only at the first call)
04b5ea16 392 if (first) {
393 first = kFALSE;
43af2cb6 394
5b64e914 395 zEndAbsorber = -503; // spectro (z<0)
b45fd22b 396 thetaLimit = 3.0 * (TMath::Pi()) / 180.;
5b64e914 397 rLimit = TMath::Abs(zEndAbsorber) * TMath::Tan(thetaLimit);
398 zBP1 = -450; // values close to those calculated with EvalAbso.C
399 zBP2 = -480;
04b5ea16 400 }
401
402 pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
403 sign = 1;
404 if (fInverseBendingMomentum < 0) sign = -1;
5b64e914 405 pZ = -pYZ / (TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope)); // spectro (z<0)
04b5ea16 406 pX = pZ * fNonBendingSlope;
407 pY = pZ * fBendingSlope;
408 pTotal = TMath::Sqrt(pYZ *pYZ + pX * pX);
409 xEndAbsorber = fNonBendingCoor;
410 yEndAbsorber = fBendingCoor;
411 radiusEndAbsorber2 = xEndAbsorber * xEndAbsorber + yEndAbsorber * yEndAbsorber;
412
413 if (radiusEndAbsorber2 > rLimit*rLimit) {
04b5ea16 414 zBP = zBP1;
415 } else {
04b5ea16 416 zBP = zBP2;
417 }
418
419 xBP = xEndAbsorber - (pX / pZ) * (zEndAbsorber - zBP);
420 yBP = yEndAbsorber - (pY / pZ) * (zEndAbsorber - zBP);
421
422 // new parameters after Branson and energy loss corrections
b45fd22b 423// Float_t zSmear = zBP - gRandom->Gaus(0.,2.); // !!! possible smearing of Z vertex position
424 Float_t zSmear = zBP;
425
426 pZ = pTotal * zSmear / TMath::Sqrt(xBP * xBP + yBP * yBP + zSmear * zSmear);
427 pX = pZ * xBP / zSmear;
428 pY = pZ * yBP / zSmear;
04b5ea16 429 fBendingSlope = pY / pZ;
430 fNonBendingSlope = pX / pZ;
5b64e914 431
04b5ea16 432
433 pT = TMath::Sqrt(pX * pX + pY * pY);
5b64e914 434 theta = TMath::ATan2(pT, TMath::Abs(pZ));
b45fd22b 435 pTotal = TotalMomentumEnergyLoss(thetaLimit, pTotal, theta);
04b5ea16 436
437 fInverseBendingMomentum = (sign / pTotal) *
438 TMath::Sqrt(1.0 +
439 fBendingSlope * fBendingSlope +
440 fNonBendingSlope * fNonBendingSlope) /
441 TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope);
442
443 // vertex position at (0,0,0)
444 // should be taken from vertex measurement ???
445 fBendingCoor = 0.0;
446 fNonBendingCoor = 0;
447 fZ= 0;
448}
b45fd22b 449
04b5ea16 450 //__________________________________________________________________________
b45fd22b 451Double_t AliMUONTrackParam::TotalMomentumEnergyLoss(Double_t thetaLimit, Double_t pTotal, Double_t theta)
04b5ea16 452{
453 // Returns the total momentum corrected from energy loss in the front absorber
43af2cb6 454 // One can use the macros MUONTestAbso.C and DrawTestAbso.C
455 // to test this correction.
b45fd22b 456 // Momentum energy loss behaviour evaluated with the simulation of single muons (april 2002)
04b5ea16 457 Double_t deltaP, pTotalCorrected;
458
b45fd22b 459 // Parametrization to be redone according to change of absorber material ????
956019b6 460 // See remark in function BransonCorrection !!!!
04b5ea16 461 // The name is not so good, and there are many arguments !!!!
b45fd22b 462 if (theta < thetaLimit ) {
463 if (pTotal < 20) {
464 deltaP = 2.5938 + 0.0570 * pTotal - 0.001151 * pTotal * pTotal;
04b5ea16 465 } else {
b45fd22b 466 deltaP = 3.0714 + 0.011767 *pTotal;
04b5ea16 467 }
468 } else {
b45fd22b 469 if (pTotal < 20) {
470 deltaP = 2.1207 + 0.05478 * pTotal - 0.00145079 * pTotal * pTotal;
04b5ea16 471 } else {
b45fd22b 472 deltaP = 2.6069 + 0.0051705 * pTotal;
04b5ea16 473 }
474 }
475 pTotalCorrected = pTotal + deltaP / TMath::Cos(theta);
476 return pTotalCorrected;
477}
478
b45fd22b 479 //__________________________________________________________________________
480void AliMUONTrackParam::FieldCorrection(Double_t Z)
481{
482 //
483 // Correction of the effect of the magnetic field in the absorber
484 // Assume a constant field along Z axis.
485
486 Float_t b[3],x[3];
487 Double_t bZ;
488 Double_t pYZ,pX,pY,pZ,pT;
489 Double_t pXNew,pYNew;
490 Double_t c;
491
492 pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
493 c = TMath::Sign(1.0,fInverseBendingMomentum); // particle charge
494
5b64e914 495 pZ = -pYZ / (TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope)); // spectro. (z<0)
b45fd22b 496 pX = pZ * fNonBendingSlope;
497 pY = pZ * fBendingSlope;
498 pT = TMath::Sqrt(pX*pX+pY*pY);
499
5b64e914 500 if (TMath::Abs(pZ) <= 0) return;
b45fd22b 501 x[2] = Z/2;
502 x[0] = x[2]*fNonBendingSlope;
503 x[1] = x[2]*fBendingSlope;
504
505 // Take magn. field value at position x.
506 gAlice->Field()->Field(x, b);
507 bZ = b[2];
508
509 // Transverse momentum rotation
510 // Parameterized with the study of DeltaPhi = phiReco - phiGen as a function of pZ.
5b64e914 511 Double_t phiShift = c*0.436*0.0003*bZ*Z/pZ;
b45fd22b 512 // Rotate momentum around Z axis.
513 pXNew = pX*TMath::Cos(phiShift) - pY*TMath::Sin(phiShift);
514 pYNew = pX*TMath::Sin(phiShift) + pY*TMath::Cos(phiShift);
515
516 fBendingSlope = pYNew / pZ;
517 fNonBendingSlope = pXNew / pZ;
518
519 fInverseBendingMomentum = c / TMath::Sqrt(pYNew*pYNew+pZ*pZ);
520
521}