BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
BDSIntegratorDipoleRodrigues2 Class Reference

Exact helix through pure dipole field. More...

#include <BDSIntegratorDipoleRodrigues2.hh>

Inheritance diagram for BDSIntegratorDipoleRodrigues2:
Inheritance graph
Collaboration diagram for BDSIntegratorDipoleRodrigues2:
Collaboration graph

Public Member Functions

 BDSIntegratorDipoleRodrigues2 (G4Mag_EqRhs *eqOfMIn, G4double minimumRadiusOfCurvature)
 
virtual void DumbStepper (const G4double yIn[6], G4ThreeVector field, G4double stepLength, G4double yOut[6])
 
void SingleStep (const G4double yIn[6], const G4double &h, G4double yOut[6])
 
virtual void Stepper (const G4double yIn[6], const G4double dydx[], G4double h, G4double yOut[6], G4double yErr[6])
 Output error array. More...
 
virtual G4int IntegratorOrder () const
 
void AdvanceHelixForSpiralling (const G4double yIn[6], const G4ThreeVector &field, const G4double &h, G4double yOut[6], G4double yErr[6])
 
G4double RadiusOfHelix () const
 Public accessor for protected variable in base class. More...
 
- Public Member Functions inherited from BDSIntegratorDrift
void AdvanceDriftMag (const G4double yIn[], const G4double h, G4double yOut[], G4double yErr[]) const
 Error array [6] all 0. More...
 
void AdvanceDriftMag (const G4double yIn[], const G4ThreeVector &unitMomentum, const G4double h, G4double yOut[], G4double yErr[]) const
 

Protected Member Functions

void FudgeDistChordToZero ()
 

Protected Attributes

G4Mag_EqRhs * eqOfM
 Cache of equation of motion. This class does not own it. More...
 

Private Member Functions

 BDSIntegratorDipoleRodrigues2 ()=delete
 Private default constructor to force use of provided one.
 

Private Attributes

G4double minimumRadiusOfCurvature
 

Detailed Description

Exact helix through pure dipole field.

Integrator based on G4MagHelicalStepper that provides required tracking througha pure magnetic field, but is a pure virtual class.

This class implements the required Stepper method that uses the provided AdvanceHelix method from the base class.

Unlike other BDSIM integrators this works in global coordinates and does not require the transforms to curvilinear coordinates.

For low momentum particles (|p| < 40 MeV) we use the G4ClassicalRK4 integrator instead of the one for the pure magnetic field as it is more robust.

Note, some things that could be passed by reference are not purely for correct overload of Geant4 function.

Author
Laurie Nevay

Definition at line 53 of file BDSIntegratorDipoleRodrigues2.hh.

Constructor & Destructor Documentation

◆ BDSIntegratorDipoleRodrigues2()

BDSIntegratorDipoleRodrigues2::BDSIntegratorDipoleRodrigues2 ( G4Mag_EqRhs *  eqOfMIn,
G4double  minimumRadiusOfCurvature 
)

Definition at line 32 of file BDSIntegratorDipoleRodrigues2.cc.

◆ ~BDSIntegratorDipoleRodrigues2()

virtual BDSIntegratorDipoleRodrigues2::~BDSIntegratorDipoleRodrigues2 ( )
inlinevirtual

Definition at line 59 of file BDSIntegratorDipoleRodrigues2.hh.

Member Function Documentation

◆ AdvanceHelixForSpiralling()

void BDSIntegratorDipoleRodrigues2::AdvanceHelixForSpiralling ( const G4double  yIn[6],
const G4ThreeVector &  field,
const G4double &  h,
G4double  yOut[6],
G4double  yErr[6] 
)

Variation of AdvanceHelix specifically to deal with particles that are likely to be spiralling in the magnetic field.

Definition at line 156 of file BDSIntegratorDipoleRodrigues2.cc.

References BDSIntegratorDrift::AdvanceDriftMag(), FudgeDistChordToZero(), and BDS::IsFinite().

Referenced by Stepper().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ DumbStepper()

void BDSIntegratorDipoleRodrigues2::DumbStepper ( const G4double  yIn[6],
G4ThreeVector  field,
G4double  stepLength,
G4double  yOut[6] 
)
virtual

Required to be provided by base class, but apparently should never be called by the driver. Simply calls AdvanceHelix.

Definition at line 39 of file BDSIntegratorDipoleRodrigues2.cc.

◆ FudgeDistChordToZero()

void BDSIntegratorDipoleRodrigues2::FudgeDistChordToZero ( )
protected

DistChord() is non-virtual function in base class so set Ang and RadHelix appropriately such that DistChord() will return 0.

Definition at line 196 of file BDSIntegratorDipoleRodrigues2.cc.

Referenced by AdvanceHelixForSpiralling(), BDSIntegratorDipoleFringe::BaseStepper(), and Stepper().

Here is the caller graph for this function:

◆ IntegratorOrder()

virtual G4int BDSIntegratorDipoleRodrigues2::IntegratorOrder ( ) const
inlinevirtual

Definition at line 84 of file BDSIntegratorDipoleRodrigues2.hh.

◆ RadiusOfHelix()

G4double BDSIntegratorDipoleRodrigues2::RadiusOfHelix ( ) const
inline

Public accessor for protected variable in base class.

Definition at line 95 of file BDSIntegratorDipoleRodrigues2.hh.

Referenced by BDSIntegratorDipoleQuadrupole::Stepper().

Here is the caller graph for this function:

◆ SingleStep()

void BDSIntegratorDipoleRodrigues2::SingleStep ( const G4double  yIn[6],
const G4double &  h,
G4double  yOut[6] 
)

Perform a single step. Useful for other integrators that want to use this one to estimate the radius of curvature in the field, in which case, no need for usual two step calculation, nor error calculation.

Definition at line 47 of file BDSIntegratorDipoleRodrigues2.cc.

References eqOfM.

Referenced by BDSIntegratorDipoleQuadrupole::Stepper().

Here is the caller graph for this function:

◆ Stepper()

void BDSIntegratorDipoleRodrigues2::Stepper ( const G4double  yIn[6],
const G4double  dydx[],
G4double  h,
G4double  yOut[6],
G4double  yErr[6] 
)
virtual

Output error array.

Calculate output coordinates. Decide if particle is spiralling or not. Nominally calculate two half steps and compare to one full step for error estimation.

Parameters
yInInput coordinates x,y,z,px,py,pz,t
dydxPartial differentials at yInput
hLength of trajectory to calculate
yOutOutput array

Definition at line 59 of file BDSIntegratorDipoleRodrigues2.cc.

References BDSIntegratorDrift::AdvanceDriftMag(), AdvanceHelixForSpiralling(), eqOfM, FudgeDistChordToZero(), BDS::IsFiniteStrength(), and minimumRadiusOfCurvature.

Referenced by BDSIntegratorDipoleFringe::BaseStepper(), and BDSIntegratorDipoleQuadrupole::Stepper().

Here is the call graph for this function:
Here is the caller graph for this function:

Field Documentation

◆ eqOfM

G4Mag_EqRhs* BDSIntegratorDipoleRodrigues2::eqOfM
protected

Cache of equation of motion. This class does not own it.

Definition at line 103 of file BDSIntegratorDipoleRodrigues2.hh.

Referenced by SingleStep(), BDSIntegratorDipoleFringe::Stepper(), Stepper(), and BDSIntegratorDipoleFringeScaling::Stepper().

◆ minimumRadiusOfCurvature

G4double BDSIntegratorDipoleRodrigues2::minimumRadiusOfCurvature
private

The minimum tolerable radius of curvature before we decide the particle is spiralling and should be treated differently.

Definition at line 111 of file BDSIntegratorDipoleRodrigues2.hh.

Referenced by Stepper().


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