AliRoot Core  3dc7879 (3dc7879)
AliMUONTrackExtrap Class Reference

Track parameters in ALICE dimuon spectrometer. More...

#include <AliMUONTrackExtrap.h>

Inheritance diagram for AliMUONTrackExtrap:

Public Member Functions

 AliMUONTrackExtrap ()
 Constructor. More...
 
virtual ~AliMUONTrackExtrap ()
 Destructor. More...
 

Static Public Member Functions

static void SetField ()
 
static Bool_t IsFieldON ()
 return kTRUE if the field is switched ON More...
 
static Double_t GetImpactParamFromBendingMomentum (Double_t bendingMomentum)
 
static Double_t GetBendingMomentumFromImpactParam (Double_t impactParam)
 
static void LinearExtrapToZ (AliMUONTrackParam *trackParam, Double_t zEnd)
 
static void LinearExtrapToZCov (AliMUONTrackParam *trackParam, Double_t zEnd, Bool_t updatePropagator=kFALSE)
 
static Bool_t ExtrapToZ (AliMUONTrackParam *trackParam, Double_t zEnd)
 
static Bool_t ExtrapToZCov (AliMUONTrackParam *trackParam, Double_t zEnd, Bool_t updatePropagator=kFALSE)
 
static void ExtrapToVertex (AliMUONTrackParam *trackParam, Double_t xVtx, Double_t yVtx, Double_t zVtx, Double_t errXVtx, Double_t errYVtx)
 
static void ExtrapToVertexWithoutELoss (AliMUONTrackParam *trackParam, Double_t xVtx, Double_t yVtx, Double_t zVtx, Double_t errXVtx, Double_t errYVtx)
 
static void ExtrapToVertexWithoutBranson (AliMUONTrackParam *trackParam, Double_t zVtx)
 
static void ExtrapToVertexUncorrected (AliMUONTrackParam *trackParam, Double_t zVtx)
 
static Double_t TotalMomentumEnergyLoss (AliMUONTrackParam *trackParam, Double_t xVtx, Double_t yVtx, Double_t zVtx)
 
static Double_t GetMCSAngle2 (const AliMUONTrackParam &param, Double_t dZ, Double_t x0)
 
static void AddMCSEffect (AliMUONTrackParam *param, Double_t dZ, Double_t x0)
 
static Bool_t ExtrapOneStepRungekutta (Double_t charge, Double_t step, const Double_t *vect, Double_t *vout)
 

Private Member Functions

 AliMUONTrackExtrap (const AliMUONTrackExtrap &trackExtrap)
 Not implemented. More...
 
AliMUONTrackExtrapoperator= (const AliMUONTrackExtrap &trackExtrap)
 Not implemented. More...
 

Static Private Member Functions

static Bool_t ExtrapToZHelix (AliMUONTrackParam *trackParam, Double_t Z)
 
static Bool_t ExtrapToZRungekutta (AliMUONTrackParam *trackParam, Double_t Z)
 
static void ConvertTrackParamForExtrap (AliMUONTrackParam *trackParam, Double_t forwardBackward, Double_t *v3)
 
static void RecoverTrackParam (Double_t *v3, Double_t Charge, AliMUONTrackParam *trackParam)
 
static void ExtrapToVertex (AliMUONTrackParam *trackParam, Double_t xVtx, Double_t yVtx, Double_t zVtx, Double_t errXVtx, Double_t errYVtx, Bool_t correctForMCS, Bool_t correctForEnergyLoss)
 
static void AddMCSEffectInAbsorber (AliMUONTrackParam *trackParam, Double_t signedPathLength, Double_t f0, Double_t f1, Double_t f2)
 
static void CorrectMCSEffectInAbsorber (AliMUONTrackParam *param, Double_t xVtx, Double_t yVtx, Double_t zVtx, Double_t errXVtx, Double_t errYVtx, Double_t absZBeg, Double_t pathLength, Double_t f0, Double_t f1, Double_t f2)
 
static void CorrectELossEffectInAbsorber (AliMUONTrackParam *param, Double_t eLoss, Double_t sigmaELoss2)
 
static Bool_t GetAbsorberCorrectionParam (Double_t trackXYZIn[3], Double_t trackXYZOut[3], Double_t pTotal, Double_t &pathLength, Double_t &f0, Double_t &f1, Double_t &f2, Double_t &meanRho, Double_t &totalELoss, Double_t &sigmaELoss2)
 
static Double_t BetheBloch (Double_t pTotal, Double_t pathLength, Double_t rho, Double_t atomicZ, Double_t atomicZoverA)
 
static Double_t EnergyLossFluctuation (Double_t pTotal, Double_t pathLength, Double_t rho, Double_t atomicZoverA)
 
static void Cov2CovP (const TMatrixD &param, TMatrixD &cov)
 
static void CovP2Cov (const TMatrixD &param, TMatrixD &cov)
 
static void ExtrapOneStepHelix (Double_t charge, Double_t step, const Double_t *vect, Double_t *vout)
 
static void ExtrapOneStepHelix3 (Double_t field, Double_t step, const Double_t *vect, Double_t *vout)
 

Static Private Attributes

static const Double_t fgkSimpleBPosition = 0.5 * (AliMUONConstants::CoilZ() + AliMUONConstants::YokeZ())
 ! position of the dipole More...
 
static const Double_t fgkSimpleBLength = 0.5 * (AliMUONConstants::CoilL() + AliMUONConstants::YokeL())
 ! length of the dipole More...
 
static Double_t fgSimpleBValue = 0.
 ! magnetic field value at the centre More...
 
static Bool_t fgFieldON = kFALSE
 ! kTRUE if the field is switched ON More...
 
static const Bool_t fgkUseHelix = kFALSE
 ! Tell whether to use Helix or not (default is Runge-Kutta) More...
 
static const Int_t fgkMaxStepNumber = 5000
 ! Maximum number of steps for track extrapolation More...
 
static const Double_t fgkHelixStepLength = 6.
 ! Step lenght for track extrapolation (used in Helix) More...
 
static const Double_t fgkRungeKuttaMaxResidue = 0.002
 ! Maximal distance (in Z) to destination to stop the track extrapolation (used in Runge-Kutta) More...
 

Detailed Description

Track parameters in ALICE dimuon spectrometer.

Tools for track extrapolation in ALICE dimuon spectrometer

Definition at line 23 of file AliMUONTrackExtrap.h.

Constructor & Destructor Documentation

AliMUONTrackExtrap::AliMUONTrackExtrap ( )
inline

Constructor.

Definition at line 27 of file AliMUONTrackExtrap.h.

virtual AliMUONTrackExtrap::~AliMUONTrackExtrap ( )
inlinevirtual

Destructor.

Definition at line 29 of file AliMUONTrackExtrap.h.

AliMUONTrackExtrap::AliMUONTrackExtrap ( const AliMUONTrackExtrap trackExtrap)
private

Not implemented.

Member Function Documentation

void AliMUONTrackExtrap::AddMCSEffect ( AliMUONTrackParam param,
Double_t  dZ,
Double_t  x0 
)
static

Add to the track parameter covariances the effects of multiple Coulomb scattering through a material of thickness "Abs(dZ)" and of radiation length "x0" assuming linear propagation and using the small angle approximation. dZ = zOut - zIn (sign is important) and "param" is assumed to be given zOut. If x0 <= 0., assume dZ = pathLength/x0 and consider the material thickness as negligible.

Definition at line 711 of file AliMUONTrackExtrap.cxx.

Referenced by AliMUONTrackHitPattern::ApplyMCSCorrections(), AliMFTTracker::Clusters2Tracks(), AliMUONVTrackReconstructor::FollowLinearTrackInChamber(), AliMUONVTrackReconstructor::FollowLinearTrackInStation(), AliMUONTrackReconstructor::FollowTrackInChamber(), AliMUONTrackReconstructorK::FollowTrackInChamber(), AliMUONTrackReconstructor::FollowTrackInStation(), AliMUONTrackReconstructorK::FollowTrackInStation(), AliMUONTriggerTrackToTrackerClusters::GenerateClusters(), IsFieldON(), AliMUONTrackReconstructorK::RetracePartialTrack(), and AliMUONTrack::UpdateCovTrackParamAtCluster().

void AliMUONTrackExtrap::AddMCSEffectInAbsorber ( AliMUONTrackParam trackParam,
Double_t  signedPathLength,
Double_t  f0,
Double_t  f1,
Double_t  f2 
)
staticprivate

Add to the track parameter covariances the effects of multiple Coulomb scattering signedPathLength must have the sign of (zOut - zIn) where all other parameters are assumed to be given at zOut.

Definition at line 422 of file AliMUONTrackExtrap.cxx.

Referenced by CorrectMCSEffectInAbsorber(), and ExtrapToVertex().

Double_t AliMUONTrackExtrap::BetheBloch ( Double_t  pTotal,
Double_t  pathLength,
Double_t  rho,
Double_t  atomicZ,
Double_t  atomicZoverA 
)
staticprivate

Returns the mean total momentum energy loss of muon with total momentum='pTotal' in the absorber layer of lenght='pathLength', density='rho', A='atomicA' and Z='atomicZ'

Definition at line 953 of file AliMUONTrackExtrap.cxx.

Referenced by GetAbsorberCorrectionParam().

void AliMUONTrackExtrap::ConvertTrackParamForExtrap ( AliMUONTrackParam trackParam,
Double_t  forwardBackward,
Double_t *  v3 
)
staticprivate

Set vector of Geant3 parameters pointed to by "v3" from track parameters in trackParam. Since AliMUONTrackParam is only geometry, one uses "forwardBackward" to know whether the particle is going forward (+1) or backward (-1).

Definition at line 309 of file AliMUONTrackExtrap.cxx.

Referenced by ExtrapToZHelix(), and ExtrapToZRungekutta().

void AliMUONTrackExtrap::CorrectELossEffectInAbsorber ( AliMUONTrackParam param,
Double_t  eLoss,
Double_t  sigmaELoss2 
)
staticprivate

Correct parameters for energy loss and add energy loss fluctuation effect to covariances

Definition at line 537 of file AliMUONTrackExtrap.cxx.

Referenced by ExtrapToVertex().

void AliMUONTrackExtrap::CorrectMCSEffectInAbsorber ( AliMUONTrackParam param,
Double_t  xVtx,
Double_t  yVtx,
Double_t  zVtx,
Double_t  errXVtx,
Double_t  errYVtx,
Double_t  absZBeg,
Double_t  pathLength,
Double_t  f0,
Double_t  f1,
Double_t  f2 
)
staticprivate

Correct parameters and corresponding covariances using Branson correction

  • input param are parameters and covariances at the end of absorber
  • output param are parameters and covariances at vertex Absorber correction parameters are supposed to be calculated at the current track z-position

Definition at line 466 of file AliMUONTrackExtrap.cxx.

Referenced by ExtrapToVertex().

void AliMUONTrackExtrap::Cov2CovP ( const TMatrixD &  param,
TMatrixD &  cov 
)
staticprivate

change coordinate system: (X, SlopeX, Y, SlopeY, q/Pyz) -> (X, SlopeX, Y, SlopeY, q*PTot) parameters (param) are given in the (X, SlopeX, Y, SlopeY, q/Pyz) coordinate system

Definition at line 987 of file AliMUONTrackExtrap.cxx.

Referenced by CorrectELossEffectInAbsorber(), and CorrectMCSEffectInAbsorber().

void AliMUONTrackExtrap::CovP2Cov ( const TMatrixD &  param,
TMatrixD &  cov 
)
staticprivate

change coordinate system: (X, SlopeX, Y, SlopeY, q*PTot) -> (X, SlopeX, Y, SlopeY, q/Pyz) parameters (param) are given in the (X, SlopeX, Y, SlopeY, q/Pyz) coordinate system

Definition at line 1010 of file AliMUONTrackExtrap.cxx.

Referenced by CorrectELossEffectInAbsorber(), and CorrectMCSEffectInAbsorber().

Double_t AliMUONTrackExtrap::EnergyLossFluctuation ( Double_t  pTotal,
Double_t  pathLength,
Double_t  rho,
Double_t  atomicZoverA 
)
staticprivate

Returns the total momentum energy loss fluctuation of muon with total momentum='pTotal' in the absorber layer of lenght='pathLength', density='rho', A='atomicA' and Z='atomicZ'

Definition at line 968 of file AliMUONTrackExtrap.cxx.

Referenced by GetAbsorberCorrectionParam().

void AliMUONTrackExtrap::ExtrapOneStepHelix ( Double_t  charge,
Double_t  step,
const Double_t *  vect,
Double_t *  vout 
)
staticprivate
   ******************************************************************
   *                                                                *
   *  Performs the tracking of one step in a magnetic field         *
   *  The trajectory is assumed to be a helix in a constant field   *
   *  taken at the mid point of the step.                           *
   *  Parameters:                                                   *
   *   input                                                        *
   *     STEP =arc length of the step asked                         *
   *     VECT =input vector (position,direction cos and momentum)   *
   *     CHARGE=  electric charge of the particle                   *
   *   output                                                       *
   *     VOUT = same as VECT after completion of the step           *
   *                                                                *
   *    ==>Called by : USER, GUSWIM                               *
   *       Author    m.hansroul  *********                          *
   *       modified  s.egli, s.v.levonian                           *
   *       modified  v.perevoztchikov
   *                                                                *
   ******************************************************************

Definition at line 1033 of file AliMUONTrackExtrap.cxx.

Referenced by ExtrapToZHelix().

void AliMUONTrackExtrap::ExtrapOneStepHelix3 ( Double_t  field,
Double_t  step,
const Double_t *  vect,
Double_t *  vout 
)
staticprivate
******************************************************************
*                *
* Tracking routine in a constant field oriented    *
* along axis 3             *
* Tracking is performed with a conventional    *
* helix step method          *
*                *
*    ==>Called by : USER, GUSWIM         *
* Authors    R.Brun, M.Hansroul  *********     *
* Rewritten  V.Perevoztchikov
*                *
******************************************************************

Definition at line 1152 of file AliMUONTrackExtrap.cxx.

Referenced by ExtrapOneStepHelix().

Bool_t AliMUONTrackExtrap::ExtrapOneStepRungekutta ( Double_t  charge,
Double_t  step,
const Double_t *  vect,
Double_t *  vout 
)
static
******************************************************************
*                *
*  Runge-Kutta method for tracking a particle through a magnetic *
*  field. Uses Nystroem algorithm (See Handbook Nat. Bur. of   *
*  Standards, procedure 25.5.20)         *
*                *
*  Input parameters            *
* CHARGE    Particle charge        *
* STEP    Step size          *
* VECT    Initial co-ords,direction cosines,momentum   *
*  Output parameters             *
* VOUT    Output co-ords,direction cosines,momentum  *
*  User routine called             *
* CALL GUFLD(X,F)            *
*                *
*    ==>Called by : USER, GUSWIM         *
* Authors    R.Brun, M.Hansroul  *********     *
*      V.Perevoztchikov (CUT STEP implementation)  *
*                *
*                *
******************************************************************

Definition at line 1230 of file AliMUONTrackExtrap.cxx.

Referenced by ExtrapToZRungekutta(), and IsFieldON().

void AliMUONTrackExtrap::ExtrapToVertex ( AliMUONTrackParam trackParam,
Double_t  xVtx,
Double_t  yVtx,
Double_t  zVtx,
Double_t  errXVtx,
Double_t  errYVtx 
)
static

Extrapolate track parameters to vertex, corrected for multiple scattering and energy loss effects Add branson correction resolution and energy loss fluctuation to parameter covariances

Definition at line 882 of file AliMUONTrackExtrap.cxx.

Referenced by AliMFTTracker::Clusters2Tracks(), ExtrapToVertexUncorrected(), ExtrapToVertexWithoutBranson(), ExtrapToVertexWithoutELoss(), AliMuonForwardTrackPair::GetMassWithoutMFT(), IsFieldON(), MUONefficiency(), MUONmassPlot(), and AliMUONESDInterface::MUONToESD().

void AliMUONTrackExtrap::ExtrapToVertex ( AliMUONTrackParam trackParam,
Double_t  xVtx,
Double_t  yVtx,
Double_t  zVtx,
Double_t  errXVtx,
Double_t  errYVtx,
Bool_t  correctForMCS,
Bool_t  correctForEnergyLoss 
)
staticprivate

Main method for extrapolation to the vertex: Returns the track parameters and covariances resulting from the extrapolation of the current trackParam Changes parameters and covariances according to multiple scattering and energy loss corrections: if correctForMCS=kTRUE: compute parameters using Branson correction and add correction resolution to covariances if correctForMCS=kFALSE: add parameter dispersion due to MCS in parameter covariances if correctForEnergyLoss=kTRUE: correct parameters for energy loss and add energy loss fluctuation to covariances if correctForEnergyLoss=kFALSE: do nothing about energy loss

Definition at line 766 of file AliMUONTrackExtrap.cxx.

void AliMUONTrackExtrap::ExtrapToVertexUncorrected ( AliMUONTrackParam trackParam,
Double_t  zVtx 
)
static

Extrapolate track parameters to vertex without multiple scattering and energy loss corrections Add dispersion due to multiple scattering to parameter covariances

Definition at line 910 of file AliMUONTrackExtrap.cxx.

Referenced by IsFieldON(), and AliMUONTrackReconstructor::SetVertexErrXY2ForFit().

void AliMUONTrackExtrap::ExtrapToVertexWithoutBranson ( AliMUONTrackParam trackParam,
Double_t  zVtx 
)
static

Extrapolate track parameters to vertex, corrected for energy loss effects only Add dispersion due to multiple scattering and energy loss fluctuation to parameter covariances

Definition at line 902 of file AliMUONTrackExtrap.cxx.

Referenced by AliMFTTracker::Clusters2Tracks(), IsFieldON(), MUONRecoCheck(), and AliMUONESDInterface::MUONToESD().

void AliMUONTrackExtrap::ExtrapToVertexWithoutELoss ( AliMUONTrackParam trackParam,
Double_t  xVtx,
Double_t  yVtx,
Double_t  zVtx,
Double_t  errXVtx,
Double_t  errYVtx 
)
static

Extrapolate track parameters to vertex, corrected for multiple scattering effects only Add branson correction resolution to parameter covariances

Definition at line 892 of file AliMUONTrackExtrap.cxx.

Referenced by IsFieldON().

Bool_t AliMUONTrackExtrap::ExtrapToZHelix ( AliMUONTrackParam trackParam,
Double_t  Z 
)
staticprivate

Track parameter extrapolation to the plane at "Z" using Helix algorithm. On return, the track parameters resulting from the extrapolation are updated in trackParam.

Definition at line 172 of file AliMUONTrackExtrap.cxx.

Referenced by ExtrapToZ().

Bool_t AliMUONTrackExtrap::ExtrapToZRungekutta ( AliMUONTrackParam trackParam,
Double_t  Z 
)
staticprivate

Track parameter extrapolation to the plane at "Z" using Rungekutta algorithm. On return, the track parameters resulting from the extrapolation are updated in trackParam.

Definition at line 226 of file AliMUONTrackExtrap.cxx.

Referenced by ExtrapToZ().

Bool_t AliMUONTrackExtrap::GetAbsorberCorrectionParam ( Double_t  trackXYZIn[3],
Double_t  trackXYZOut[3],
Double_t  pTotal,
Double_t &  pathLength,
Double_t &  f0,
Double_t &  f1,
Double_t &  f2,
Double_t &  meanRho,
Double_t &  totalELoss,
Double_t &  sigmaELoss2 
)
staticprivate

Parameters used to correct for Multiple Coulomb Scattering and energy loss in absorber Calculated assuming a linear propagation from trackXYZIn to trackXYZOut (order is important)

Definition at line 568 of file AliMUONTrackExtrap.cxx.

Referenced by ExtrapToVertex(), and TotalMomentumEnergyLoss().

Double_t AliMUONTrackExtrap::GetBendingMomentumFromImpactParam ( Double_t  impactParam)
static

Returns signed bending momentum in bending plane (GeV/c), the sign being the sign of the charge for particles moving forward in Z, from the impact parameter "ImpactParam" at vertex in bending plane (cm), using simple values for dipole magnetic field.

Definition at line 85 of file AliMUONTrackExtrap.cxx.

Referenced by AliMUONTrack::AliMUONTrack(), AliMUONTriggerTrackToTrackerClusters::GenerateClusters(), IsFieldON(), AliMUONVTrackReconstructor::MakeSegmentsBetweenChambers(), and AliMUONTrackReconstructorK::RetraceTrack().

Double_t AliMUONTrackExtrap::GetImpactParamFromBendingMomentum ( Double_t  bendingMomentum)
static

Returns impact parameter at vertex in bending plane (cm), from the signed bending momentum "BendingMomentum" in bending plane (GeV/c), using simple values for dipole magnetic field. The sign of "BendingMomentum" is the sign of the charge.

Definition at line 69 of file AliMUONTrackExtrap.cxx.

Referenced by IsFieldON().

Double_t AliMUONTrackExtrap::GetMCSAngle2 ( const AliMUONTrackParam param,
Double_t  dZ,
Double_t  x0 
)
static

Return the angular dispersion square due to multiple Coulomb scattering through a material of thickness "dZ" and of radiation length "x0" assuming linear propagation and using the small angle approximation.

Definition at line 689 of file AliMUONTrackExtrap.cxx.

Referenced by AliMUONVTrackReconstructor::AliMUONVTrackReconstructor(), AliMUONTrack::ComputeMCSCovariances(), AliMUONVTrackReconstructor::IsAcceptable(), and IsFieldON().

static Bool_t AliMUONTrackExtrap::IsFieldON ( )
inlinestatic
void AliMUONTrackExtrap::LinearExtrapToZ ( AliMUONTrackParam trackParam,
Double_t  zEnd 
)
static

Track parameters linearly extrapolated to the plane at "zEnd". On return, results from the extrapolation are updated in trackParam.

Definition at line 107 of file AliMUONTrackExtrap.cxx.

Referenced by ExtrapToZ(), AliMUONTrackHitPattern::FindMatchingPads(), AliMUONTrackHitPattern::GetTrackParamAtChamber(), IsFieldON(), and LinearExtrapToZCov().

void AliMUONTrackExtrap::LinearExtrapToZCov ( AliMUONTrackParam trackParam,
Double_t  zEnd,
Bool_t  updatePropagator = kFALSE 
)
static
AliMUONTrackExtrap& AliMUONTrackExtrap::operator= ( const AliMUONTrackExtrap trackExtrap)
private

Not implemented.

void AliMUONTrackExtrap::RecoverTrackParam ( Double_t *  v3,
Double_t  Charge,
AliMUONTrackParam trackParam 
)
staticprivate

Set track parameters in trackParam from Geant3 parameters pointed to by "v3", assumed to be calculated for forward motion in Z. "InverseBendingMomentum" is signed with "charge".

Definition at line 326 of file AliMUONTrackExtrap.cxx.

Referenced by ExtrapToZHelix(), and ExtrapToZRungekutta().

Double_t AliMUONTrackExtrap::TotalMomentumEnergyLoss ( AliMUONTrackParam trackParam,
Double_t  xVtx,
Double_t  yVtx,
Double_t  zVtx 
)
static

Calculate the total momentum energy loss in-between the track position and the vertex assuming a linear propagation

Definition at line 918 of file AliMUONTrackExtrap.cxx.

Referenced by IsFieldON().

Member Data Documentation

Bool_t AliMUONTrackExtrap::fgFieldON = kFALSE
staticprivate

! kTRUE if the field is switched ON

Definition at line 83 of file AliMUONTrackExtrap.h.

Referenced by AddMCSEffect(), AddMCSEffectInAbsorber(), ExtrapToZ(), ExtrapToZCov(), GetBendingMomentumFromImpactParam(), and IsFieldON().

const Double_t AliMUONTrackExtrap::fgkHelixStepLength = 6.
staticprivate

! Step lenght for track extrapolation (used in Helix)

Definition at line 86 of file AliMUONTrackExtrap.h.

Referenced by ExtrapToZHelix().

const Int_t AliMUONTrackExtrap::fgkMaxStepNumber = 5000
staticprivate

! Maximum number of steps for track extrapolation

Definition at line 85 of file AliMUONTrackExtrap.h.

Referenced by ExtrapToZHelix(), and ExtrapToZRungekutta().

const Double_t AliMUONTrackExtrap::fgkRungeKuttaMaxResidue = 0.002
staticprivate

! Maximal distance (in Z) to destination to stop the track extrapolation (used in Runge-Kutta)

Definition at line 87 of file AliMUONTrackExtrap.h.

Referenced by ExtrapToZRungekutta().

const Double_t AliMUONTrackExtrap::fgkSimpleBLength = 0.5 * (AliMUONConstants::CoilL() + AliMUONConstants::YokeL())
staticprivate

! length of the dipole

Definition at line 81 of file AliMUONTrackExtrap.h.

Referenced by GetBendingMomentumFromImpactParam(), and GetImpactParamFromBendingMomentum().

const Double_t AliMUONTrackExtrap::fgkSimpleBPosition = 0.5 * (AliMUONConstants::CoilZ() + AliMUONConstants::YokeZ())
staticprivate

! position of the dipole

Definition at line 80 of file AliMUONTrackExtrap.h.

Referenced by GetBendingMomentumFromImpactParam(), and GetImpactParamFromBendingMomentum().

const Bool_t AliMUONTrackExtrap::fgkUseHelix = kFALSE
staticprivate

! Tell whether to use Helix or not (default is Runge-Kutta)

Definition at line 84 of file AliMUONTrackExtrap.h.

Referenced by ExtrapToZ().

Double_t AliMUONTrackExtrap::fgSimpleBValue = 0.
staticprivate

! magnetic field value at the centre

Definition at line 82 of file AliMUONTrackExtrap.h.

Referenced by ExtrapToZRungekutta(), GetBendingMomentumFromImpactParam(), and GetImpactParamFromBendingMomentum().


The documentation for this class was generated from the following files: