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