]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONTrackParam.cxx
Add pdfset_alice and structm_alice handle calls to the nuclear structure
[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
16/*
17$Log$
956019b6 18Revision 1.4 2000/07/03 07:53:31 morsch
19Double declaration problem on HP solved.
20
88962f0c 21Revision 1.3 2000/06/30 10:15:48 gosset
22Changes to EventReconstructor...:
23precision fit with multiple Coulomb scattering;
24extrapolation to vertex with Branson correction in absorber (JPC)
25
04b5ea16 26Revision 1.2 2000/06/15 07:58:49 morsch
27Code from MUON-dev joined
28
a9e2aefa 29Revision 1.1.2.3 2000/06/09 21:03:09 morsch
30Make includes consistent with new file structure.
31
32Revision 1.1.2.2 2000/06/09 12:58:05 gosset
33Removed comment beginnings in Log sections of .cxx files
34Suppressed most violations of coding rules
35
36Revision 1.1.2.1 2000/06/07 14:44:53 gosset
37Addition of files for track reconstruction in C++
38*/
39
40//__________________________________________________________________________
41//
42// Track parameters in ALICE dimuon spectrometer
43//__________________________________________________________________________
44
45#include <iostream.h>
46
47#include "AliCallf77.h"
48#include "AliMUON.h"
49#include "AliMUONHitForRec.h"
50#include "AliMUONSegment.h"
51#include "AliMUONTrackParam.h"
52#include "AliMUONChamber.h"
53#include "AliRun.h"
54
55ClassImp(AliMUONTrackParam) // Class implementation in ROOT context
56
57#ifndef WIN32
58# define reco_ghelix reco_ghelix_
59#else
60# define reco_ghelix RECO_GHELIX
61#endif
62
63extern "C"
64{
65void type_of_call reco_ghelix(Double_t &Charge, Double_t &StepLength, Double_t *VGeant3, Double_t *VGeant3New);
66}
67
68// Inline functions for Get and Set: inline removed because it does not work !!!!
69Double_t AliMUONTrackParam::GetInverseBendingMomentum(void) {
70 // Get fInverseBendingMomentum
71 return fInverseBendingMomentum;}
72void AliMUONTrackParam::SetInverseBendingMomentum(Double_t InverseBendingMomentum) {
73 // Set fInverseBendingMomentum
74 fInverseBendingMomentum = InverseBendingMomentum;}
75Double_t AliMUONTrackParam::GetBendingSlope(void) {
76 // Get fBendingSlope
77 return fBendingSlope;}
78void AliMUONTrackParam::SetBendingSlope(Double_t BendingSlope) {
79 // Set fBendingSlope
80 fBendingSlope = BendingSlope;}
81Double_t AliMUONTrackParam::GetNonBendingSlope(void) {
82 // Get fNonBendingSlope
83 return fNonBendingSlope;}
84void AliMUONTrackParam::SetNonBendingSlope(Double_t NonBendingSlope) {
85 // Set fNonBendingSlope
86 fNonBendingSlope = NonBendingSlope;}
87Double_t AliMUONTrackParam::GetZ(void) {
88 // Get fZ
89 return fZ;}
90void AliMUONTrackParam::SetZ(Double_t Z) {
91 // Set fZ
92 fZ = Z;}
93Double_t AliMUONTrackParam::GetBendingCoor(void) {
94 // Get fBendingCoor
95 return fBendingCoor;}
96void AliMUONTrackParam::SetBendingCoor(Double_t BendingCoor) {
97 // Set fBendingCoor
98 fBendingCoor = BendingCoor;}
99Double_t AliMUONTrackParam::GetNonBendingCoor(void) {
100 // Get fNonBendingCoor
101 return fNonBendingCoor;}
102void AliMUONTrackParam::SetNonBendingCoor(Double_t NonBendingCoor) {
103 // Set fNonBendingCoor
104 fNonBendingCoor = NonBendingCoor;}
105
106 //__________________________________________________________________________
107void AliMUONTrackParam::ExtrapToZ(Double_t Z)
108{
109 // Track parameter extrapolation to the plane at "Z".
110 // On return, the track parameters resulting from the extrapolation
111 // replace the current track parameters.
112 // Use "reco_ghelix" which should be replaced by something else !!!!
113 if (this->fZ == Z) return; // nothing to be done if same Z
114 Double_t forwardBackward; // +1 if forward, -1 if backward
115 if (Z > this->fZ) forwardBackward = 1.0;
116 else forwardBackward = -1.0;
117 Double_t temp, vGeant3[7], vGeant3New[7]; // 7 in parameter ????
118 Int_t iGeant3, stepNumber;
119 Int_t maxStepNumber = 5000; // in parameter ????
120 // For safety: return kTRUE or kFALSE ????
121 // Parameter vector for calling GHELIX in Geant3
122 SetGeant3Parameters(vGeant3, forwardBackward);
123 // For use of reco_ghelix...: invert X and Y, PX/PTOT and PY/PTOT !!!!
124 temp = vGeant3[0]; vGeant3[0] = vGeant3[1]; vGeant3[1] = temp;
125 temp = vGeant3[3]; vGeant3[3] = vGeant3[4]; vGeant3[4] = temp;
956019b6 126 // sign of charge (sign of fInverseBendingMomentum if forward motion)
127 // must be also changed if backward extrapolation
128 Double_t chargeExtrap = forwardBackward *
129 TMath::Sign(Double_t(1.0), this->fInverseBendingMomentum);
a9e2aefa 130 Double_t stepLength = 6.0; // in parameter ????
131 // Extrapolation loop
132 stepNumber = 0;
133 while (((forwardBackward * (vGeant3[2] - Z)) <= 0.0) &&
134 (stepNumber < maxStepNumber)) {
135 stepNumber++;
136 // call Geant3 "ghelix" subroutine through a copy in "reco_muon.F":
137 // the true function should be called, but how ???? and remove prototyping ...
956019b6 138 reco_ghelix(chargeExtrap, stepLength, vGeant3, vGeant3New);
a9e2aefa 139 if ((forwardBackward * (vGeant3New[2] - Z)) > 0.0) break; // one is beyond Z
140 // better use TArray ????
141 for (iGeant3 = 0; iGeant3 < 7; iGeant3++)
142 {vGeant3[iGeant3] = vGeant3New[iGeant3];}
143 }
144 // check maxStepNumber ????
145 // For use of reco_ghelix...:
146 // invert back X and Y, PX/PTOT and PY/PTOT, both for vGeant3 and vGeant3New !!!!
147 temp = vGeant3[0]; vGeant3[0] = vGeant3[1]; vGeant3[1] = temp;
148 temp = vGeant3New[0]; vGeant3New[0] = vGeant3New[1]; vGeant3New[1] = temp;
149 temp = vGeant3[3]; vGeant3[3] = vGeant3[4]; vGeant3[4] = temp;
150 temp = vGeant3New[3]; vGeant3New[3] = vGeant3New[4]; vGeant3New[4] = temp;
151 // Interpolation back to exact Z (2nd order)
152 // should be in function ???? using TArray ????
153 Double_t dZ12 = vGeant3New[2] - vGeant3[2]; // 1->2
154 Double_t dZ1i = Z - vGeant3[2]; // 1-i
155 Double_t dZi2 = vGeant3New[2] - Z; // i->2
156 Double_t xPrime = (vGeant3New[0] - vGeant3[0]) / dZ12;
157 Double_t xSecond =
158 ((vGeant3New[3] / vGeant3New[5]) - (vGeant3[3] / vGeant3[5])) / dZ12;
159 Double_t yPrime = (vGeant3New[1] - vGeant3[1]) / dZ12;
160 Double_t ySecond =
161 ((vGeant3New[4] / vGeant3New[5]) - (vGeant3[4] / vGeant3[5])) / dZ12;
162 vGeant3[0] = vGeant3[0] + xPrime * dZ1i - 0.5 * xSecond * dZ1i * dZi2; // X
163 vGeant3[1] = vGeant3[1] + yPrime * dZ1i - 0.5 * ySecond * dZ1i * dZi2; // Y
164 vGeant3[2] = Z; // Z
165 Double_t xPrimeI = xPrime - 0.5 * xSecond * (dZi2 - dZ1i);
166 Double_t yPrimeI = yPrime - 0.5 * ySecond * (dZi2 - dZ1i);
956019b6 167 // (PX, PY, PZ)/PTOT assuming forward motion
a9e2aefa 168 vGeant3[5] =
169 1.0 / TMath::Sqrt(1.0 + xPrimeI * xPrimeI + yPrimeI * yPrimeI); // PZ/PTOT
170 vGeant3[3] = xPrimeI * vGeant3[5]; // PX/PTOT
171 vGeant3[4] = yPrimeI * vGeant3[5]; // PY/PTOT
956019b6 172 // Track parameters from Geant3 parameters,
173 // with charge back for forward motion
174 GetFromGeant3Parameters(vGeant3, chargeExtrap * forwardBackward);
a9e2aefa 175}
176
177 //__________________________________________________________________________
178void AliMUONTrackParam::SetGeant3Parameters(Double_t *VGeant3, Double_t ForwardBackward)
179{
180 // Set vector of Geant3 parameters pointed to by "VGeant3"
181 // from track parameters in current AliMUONTrackParam.
182 // Since AliMUONTrackParam is only geometry, one uses "ForwardBackward"
183 // to know whether the particle is going forward (+1) or backward (-1).
184 VGeant3[0] = this->fNonBendingCoor; // X
185 VGeant3[1] = this->fBendingCoor; // Y
186 VGeant3[2] = this->fZ; // Z
187 Double_t pYZ = TMath::Abs(1.0 / this->fInverseBendingMomentum);
188 Double_t pZ =
189 pYZ / TMath::Sqrt(1.0 + this->fBendingSlope * this->fBendingSlope);
190 VGeant3[6] =
191 TMath::Sqrt(pYZ * pYZ +
192 pZ * pZ * this->fNonBendingSlope * this->fNonBendingSlope); // PTOT
193 VGeant3[5] = ForwardBackward * pZ / VGeant3[6]; // PZ/PTOT
194 VGeant3[3] = this->fNonBendingSlope * VGeant3[5]; // PX/PTOT
195 VGeant3[4] = this->fBendingSlope * VGeant3[5]; // PY/PTOT
196}
197
198 //__________________________________________________________________________
199void AliMUONTrackParam::GetFromGeant3Parameters(Double_t *VGeant3, Double_t Charge)
200{
201 // Get track parameters in current AliMUONTrackParam
956019b6 202 // from Geant3 parameters pointed to by "VGeant3",
203 // assumed to be calculated for forward motion in Z.
a9e2aefa 204 // "InverseBendingMomentum" is signed with "Charge".
205 this->fNonBendingCoor = VGeant3[0]; // X
206 this->fBendingCoor = VGeant3[1]; // Y
207 this->fZ = VGeant3[2]; // Z
208 Double_t pYZ = VGeant3[6] * TMath::Sqrt(1.0 - VGeant3[3] * VGeant3[3]);
209 this->fInverseBendingMomentum = Charge / pYZ;
210 this->fBendingSlope = VGeant3[4] / VGeant3[5];
211 this->fNonBendingSlope = VGeant3[3] / VGeant3[5];
212}
213
214 //__________________________________________________________________________
215void AliMUONTrackParam::ExtrapToStation(Int_t Station, AliMUONTrackParam *TrackParam)
216{
217 // Track parameters extrapolated from current track parameters ("this")
218 // to both chambers of the station(0..) "Station"
219 // are returned in the array (dimension 2) of track parameters
220 // pointed to by "TrackParam" (index 0 and 1 for first and second chambers).
221 Double_t extZ[2], z1, z2;
222 Int_t i1, i2;
223 AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
224 // range of Station to be checked ????
225 z1 = (&(pMUON->Chamber(2 * Station)))->Z(); // Z of first chamber
226 z2 = (&(pMUON->Chamber(2 * Station + 1)))->Z(); // Z of second chamber
227 // First and second Z to extrapolate at
228 if ((z1 > this->fZ) && (z2 > this->fZ)) {i1 = 0; i2 = 1;}
229 else if ((z1 < this->fZ) && (z2 < this->fZ)) {i1 = 1; i2 = 0;}
230 else {
231 cout << "ERROR in AliMUONTrackParam::CreateExtrapSegmentInStation" << endl;
232 cout << "Starting Z (" << this->fZ << ") in between z1 (" << z1 <<
233 ") and z2 (" << z2 << ") of station(0..) " << Station << endl;
234 }
235 extZ[i1] = z1;
236 extZ[i2] = z2;
237 // copy of track parameters
238 TrackParam[i1] = *this;
239 // first extrapolation
240 (&(TrackParam[i1]))->ExtrapToZ(extZ[0]);
241 TrackParam[i2] = TrackParam[i1];
242 // second extrapolation
243 (&(TrackParam[i2]))->ExtrapToZ(extZ[1]);
244 return;
245}
246
04b5ea16 247 //__________________________________________________________________________
248void AliMUONTrackParam::ExtrapToVertex()
249{
250 // Extrapolation to the vertex.
251 // Returns the track parameters resulting from the extrapolation,
252 // in the current TrackParam.
956019b6 253 // Changes parameters according to Branson correction through the absorber
04b5ea16 254
255 Double_t zAbsorber = 503.0; // to be coherent with the Geant absorber geometry !!!!
256 // Extrapolates track parameters upstream to the "Z" end of the front absorber
257 ExtrapToZ(zAbsorber);
258 // Makes Branson correction (multiple scattering + energy loss)
259 BransonCorrection();
260}
261
262 //__________________________________________________________________________
263void AliMUONTrackParam::BransonCorrection()
264{
265 // Branson correction of track parameters
266 // the entry parameters have to be calculated at the end of the absorber
267 Double_t zEndAbsorber, zBP, xBP, yBP;
268 Double_t pYZ, pX, pY, pZ, pTotal, xEndAbsorber, yEndAbsorber, radiusEndAbsorber2, pT, theta;
269 Int_t sign;
270 // Would it be possible to calculate all that from Geant configuration ????
956019b6 271 // and to get the Branson parameters from a function in ABSO module ????
272 // with an eventual contribution from other detectors like START ????
04b5ea16 273 // Radiation lengths outer part theta > 3 degres
274 static Double_t x01[9] = { 18.8, // C (cm)
275 10.397, // Concrete (cm)
276 0.56, // Plomb (cm)
277 47.26, // Polyethylene (cm)
278 0.56, // Plomb (cm)
279 47.26, // Polyethylene (cm)
280 0.56, // Plomb (cm)
281 47.26, // Polyethylene (cm)
282 0.56 }; // Plomb (cm)
283 // inner part theta < 3 degres
284 static Double_t x02[3] = { 18.8, // C (cm)
285 10.397, // Concrete (cm)
286 0.35 }; // W (cm)
287 // z positions of the materials inside the absober outer part theta > 3 degres
288 static Double_t z1[10] = { 90, 315, 467, 472, 477, 482, 487, 492, 497, 502 };
289 // inner part theta < 3 degres
290 static Double_t z2[4] = { 90, 315, 467, 503 };
291 static Bool_t first = kTRUE;
292 static Double_t zBP1, zBP2, rLimit;
293 // Calculates z positions of the Branson's planes: zBP1 for outer part and zBP2 for inner part (only at the first call)
294 if (first) {
295 first = kFALSE;
296 Double_t aNBP = 0.0;
297 Double_t aDBP = 0.0;
88962f0c 298 Int_t iBound;
299
300 for (iBound = 0; iBound < 9; iBound++) {
04b5ea16 301 aNBP = aNBP +
302 (z1[iBound+1] * z1[iBound+1] * z1[iBound+1] -
303 z1[iBound] * z1[iBound] * z1[iBound] ) / x01[iBound];
304 aDBP = aDBP +
305 (z1[iBound+1] * z1[iBound+1] - z1[iBound] * z1[iBound] ) / x01[iBound];
306 }
307 zBP1 = (2.0 * aNBP) / (3.0 * aDBP);
308 aNBP = 0.0;
309 aDBP = 0.0;
88962f0c 310 for (iBound = 0; iBound < 3; iBound++) {
04b5ea16 311 aNBP = aNBP +
312 (z2[iBound+1] * z2[iBound+1] * z2[iBound+1] -
313 z2[iBound] * z2[iBound ] * z2[iBound] ) / x02[iBound];
314 aDBP = aDBP +
315 (z2[iBound+1] * z2[iBound+1] - z2[iBound] * z2[iBound]) / x02[iBound];
316 }
317 zBP2 = (2.0 * aNBP) / (3.0 * aDBP);
318 rLimit = z2[3] * TMath::Tan(3.0 * (TMath::Pi()) / 180.);
319 }
320
321 pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
322 sign = 1;
323 if (fInverseBendingMomentum < 0) sign = -1;
324 pZ = pYZ / (TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope));
325 pX = pZ * fNonBendingSlope;
326 pY = pZ * fBendingSlope;
327 pTotal = TMath::Sqrt(pYZ *pYZ + pX * pX);
328 xEndAbsorber = fNonBendingCoor;
329 yEndAbsorber = fBendingCoor;
330 radiusEndAbsorber2 = xEndAbsorber * xEndAbsorber + yEndAbsorber * yEndAbsorber;
331
332 if (radiusEndAbsorber2 > rLimit*rLimit) {
333 zEndAbsorber = z1[9];
334 zBP = zBP1;
335 } else {
956019b6 336 zEndAbsorber = z2[3];
04b5ea16 337 zBP = zBP2;
338 }
339
340 xBP = xEndAbsorber - (pX / pZ) * (zEndAbsorber - zBP);
341 yBP = yEndAbsorber - (pY / pZ) * (zEndAbsorber - zBP);
342
343 // new parameters after Branson and energy loss corrections
344 pZ = pTotal * zBP / TMath::Sqrt(xBP * xBP + yBP * yBP + zBP * zBP);
345 pX = pZ * xBP / zBP;
346 pY = pZ * yBP / zBP;
347 fBendingSlope = pY / pZ;
348 fNonBendingSlope = pX / pZ;
349
350 pT = TMath::Sqrt(pX * pX + pY * pY);
351 theta = TMath::ATan2(pT, pZ);
352 pTotal =
353 TotalMomentumEnergyLoss(rLimit, pTotal, theta, xEndAbsorber, yEndAbsorber);
354
355 fInverseBendingMomentum = (sign / pTotal) *
356 TMath::Sqrt(1.0 +
357 fBendingSlope * fBendingSlope +
358 fNonBendingSlope * fNonBendingSlope) /
359 TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope);
360
361 // vertex position at (0,0,0)
362 // should be taken from vertex measurement ???
363 fBendingCoor = 0.0;
364 fNonBendingCoor = 0;
365 fZ= 0;
366}
367
368 //__________________________________________________________________________
369Double_t AliMUONTrackParam::TotalMomentumEnergyLoss(Double_t rLimit, Double_t pTotal, Double_t theta, Double_t xEndAbsorber, Double_t yEndAbsorber)
370{
371 // Returns the total momentum corrected from energy loss in the front absorber
372 Double_t deltaP, pTotalCorrected;
373
374 Double_t radiusEndAbsorber2 =
375 xEndAbsorber *xEndAbsorber + yEndAbsorber * yEndAbsorber;
376 // Parametrization to be redone according to change of absorber material ????
956019b6 377 // See remark in function BransonCorrection !!!!
04b5ea16 378 // The name is not so good, and there are many arguments !!!!
379 if (radiusEndAbsorber2 < rLimit * rLimit) {
380 if (pTotal < 15) {
381 deltaP = 2.737 + 0.0494 * pTotal - 0.001123 * pTotal * pTotal;
382 } else {
383 deltaP = 3.0643 + 0.01346 *pTotal;
384 }
385 } else {
386 if (pTotal < 15) {
387 deltaP = 2.1380 + 0.0351 * pTotal - 0.000853 * pTotal * pTotal;
388 } else {
389 deltaP = 2.407 + 0.00702 * pTotal;
390 }
391 }
392 pTotalCorrected = pTotal + deltaP / TMath::Cos(theta);
393 return pTotalCorrected;
394}
395