Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Mag_UsualEqRhs Class Reference

#include <G4Mag_UsualEqRhs.hh>

+ Inheritance diagram for G4Mag_UsualEqRhs:

Public Member Functions

 G4Mag_UsualEqRhs (G4MagneticField *MagField)
 
virtual ~G4Mag_UsualEqRhs ()
 
void EvaluateRhsGivenB (const G4double y[], const G4double B[3], G4double dydx[]) const
 
virtual void SetChargeMomentumMass (G4ChargeState particleCharge, G4double MomentumXc, G4double mass)
 
- Public Member Functions inherited from G4Mag_EqRhs
 G4Mag_EqRhs (G4MagneticField *magField)
 
virtual ~G4Mag_EqRhs ()
 
virtual void EvaluateRhsGivenB (const G4double y[], const G4double B[3], G4double dydx[]) const =0
 
G4double FCof () const
 
virtual void SetChargeMomentumMass (G4ChargeState particleCharge, G4double MomentumXc, G4double mass)
 
- Public Member Functions inherited from G4EquationOfMotion
 G4EquationOfMotion (G4Field *Field)
 
virtual ~G4EquationOfMotion ()
 
virtual void EvaluateRhsGivenB (const G4double y[], const G4double B[3], G4double dydx[]) const =0
 
virtual void SetChargeMomentumMass (G4ChargeState particleCharge, G4double MomentumXc, G4double MassXc2)=0
 
void RightHandSide (const G4double y[], G4double dydx[]) const
 
void EvaluateRhsReturnB (const G4double y[], G4double dydx[], G4double Field[]) const
 
void GetFieldValue (const G4double Point[4], G4double Field[]) const
 
const G4FieldGetFieldObj () const
 
G4FieldGetFieldObj ()
 
void SetFieldObj (G4Field *pField)
 

Detailed Description

Definition at line 45 of file G4Mag_UsualEqRhs.hh.

Constructor & Destructor Documentation

◆ G4Mag_UsualEqRhs()

G4Mag_UsualEqRhs::G4Mag_UsualEqRhs ( G4MagneticField MagField)

Definition at line 36 of file G4Mag_UsualEqRhs.cc.

37 : G4Mag_EqRhs( MagField )
38{
39}

◆ ~G4Mag_UsualEqRhs()

G4Mag_UsualEqRhs::~G4Mag_UsualEqRhs ( )
virtual

Definition at line 41 of file G4Mag_UsualEqRhs.cc.

42{
43}

Member Function Documentation

◆ EvaluateRhsGivenB()

void G4Mag_UsualEqRhs::EvaluateRhsGivenB ( const G4double  y[],
const G4double  B[3],
G4double  dydx[] 
) const
virtual

Implements G4Mag_EqRhs.

Definition at line 46 of file G4Mag_UsualEqRhs.cc.

49{
50 G4double momentum_mag_square = y[3]*y[3] + y[4]*y[4] + y[5]*y[5];
51 G4double inv_momentum_magnitude = 1.0 / std::sqrt( momentum_mag_square );
52
53 G4double cof = FCof()*inv_momentum_magnitude;
54
55 dydx[0] = y[3]*inv_momentum_magnitude; // (d/ds)x = Vx/V
56 dydx[1] = y[4]*inv_momentum_magnitude; // (d/ds)y = Vy/V
57 dydx[2] = y[5]*inv_momentum_magnitude; // (d/ds)z = Vz/V
58
59 dydx[3] = cof*(y[4]*B[2] - y[5]*B[1]) ; // Ax = a*(Vy*Bz - Vz*By)
60 dydx[4] = cof*(y[5]*B[0] - y[3]*B[2]) ; // Ay = a*(Vz*Bx - Vx*Bz)
61 dydx[5] = cof*(y[3]*B[1] - y[4]*B[0]) ; // Az = a*(Vx*By - Vy*Bx)
62
63 return;
64}
double B(double temperature)
double G4double
Definition: G4Types.hh:83
G4double FCof() const
Definition: G4Mag_EqRhs.hh:62

Referenced by G4ErrorMag_UsualEqRhs::EvaluateRhsGivenB().

◆ SetChargeMomentumMass()

void G4Mag_UsualEqRhs::SetChargeMomentumMass ( G4ChargeState  particleCharge,
G4double  MomentumXc,
G4double  mass 
)
virtual

Reimplemented from G4Mag_EqRhs.

Definition at line 67 of file G4Mag_UsualEqRhs.cc.

71{
72 G4Mag_EqRhs::SetChargeMomentumMass( particleCharge, MomentumXc, mass);
73}
virtual void SetChargeMomentumMass(G4ChargeState particleCharge, G4double MomentumXc, G4double mass)
Definition: G4Mag_EqRhs.cc:49

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