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