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