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