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

Common functionality to BDSIM integrators. More...

#include <BDSIntegratorMag.hh>

Inheritance diagram for BDSIntegratorMag:
Inheritance graph
Collaboration diagram for BDSIntegratorMag:
Collaboration graph

Public Member Functions

 BDSIntegratorMag (G4Mag_EqRhs *eqOfMIn, G4int nVariablesIn)
 
virtual G4double DistChord () const
 Estimate maximum distance of curved solution and chord. More...
 
virtual G4int IntegratorOrder () const
 Geant4 requires that the integrator order must be supplied by the derived class. More...
 
BDSIntegratorMagoperator= (const BDSIntegratorMag &)=delete
 Assignment and copy constructor not implemented nor used.
 
 BDSIntegratorMag (BDSIntegratorMag &)=delete
 Assignment and copy constructor not implemented nor used.
 
- 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
 
- Public Member Functions inherited from BDSAuxiliaryNavigator
G4VPhysicalVolume * LocateGlobalPointAndSetup (const G4ThreeVector &point, const G4ThreeVector *direction=nullptr, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true, G4bool useCurvilinear=true) const
 A wrapper for the underlying static navigator instance located within this class. More...
 
G4VPhysicalVolume * LocateGlobalPointAndSetup (G4Step const *const step, G4bool useCurvilinear=true) const
 
BDSStep ConvertToLocal (G4Step const *const step, G4bool useCurvilinear=true) const
 
BDSStep ConvertToLocal (const G4ThreeVector &globalPosition, const G4ThreeVector &globalDirection, const G4double stepLength=0, const G4bool useCurvilinear=true, const G4double marginLength=1) const
 
BDSStep ConvertToGlobalStep (const G4ThreeVector &localPosition, const G4ThreeVector &localDirection, const G4bool useCurvilinear=true) const
 
G4ThreeVector ConvertToLocal (const G4double globalPoint[3], const G4bool useCurvilinear=true) const
 
G4ThreeVector ConvertToLocal (const G4ThreeVector &globalPosition, const G4bool useCurvilinear=true) const
 Vector version - see notes above. More...
 
G4ThreeVector ConvertToLocalNoSetup (const G4ThreeVector &globalPosition, const G4bool useCurvilinear=true) const
 Similar to above function but does NOT initialise the transforms. More...
 
G4ThreeVector ConvertAxisToLocal (const G4ThreeVector &globalAxis, const G4bool useCurvilinear=true) const
 
G4ThreeVector ConvertAxisToLocal (const G4double globalPoint[3], const G4double globalAxis[3], const G4bool useCurvilinear=true) const
 
G4ThreeVector ConvertAxisToLocal (const G4ThreeVector &globalPoint, const G4ThreeVector &globalAxis, const G4bool useCurvilinear=true) const
 Vector version. More...
 
G4ThreeVector ConvertAxisToGlobal (const G4ThreeVector &localAxis, const G4bool useCurvilinear=true) const
 
std::pair< G4ThreeVector, G4ThreeVector > ConvertAxisToGlobal (const std::pair< G4ThreeVector, G4ThreeVector > &localAxis, const G4bool useCurvilinear=true) const
 
G4ThreeVector ConvertToGlobal (const G4ThreeVector &localPosition, const G4bool useCurvilinear=true) const
 
G4ThreeVector ConvertAxisToGlobal (const G4ThreeVector &globalPosition, const G4ThreeVector &localAxis, const G4bool useCurvilinear=true) const
 
G4ThreeVector ConvertToGlobal (const G4ThreeVector &globalPosition, const G4ThreeVector &localPosition, const G4bool useCurvilinear=true) const
 
BDSStep GlobalToCurvilinear (const G4double fieldArcLength, const G4ThreeVector &unitField, const G4double angle, const G4ThreeVector &position, const G4ThreeVector &unitMomentum, const G4double h, const G4bool useCurvilinearWorld, const G4double FCof, const G4double tilt=0)
 
BDSStep GlobalToCurvilinear (const G4ThreeVector &position, const G4ThreeVector &unitMomentum, const G4double h, const G4bool useCurvilinearWorld)
 
BDSStep CurvilinearToGlobal (const G4ThreeVector &localPosition, const G4ThreeVector &localMomentum, const G4bool useCurvilinearWorld)
 
BDSStep CurvilinearToGlobal (const G4double fieldArcLength, const G4ThreeVector &unitField, const G4double angle, const G4ThreeVector &CLPosition, const G4ThreeVector &CLMomentum, const G4bool useCurvilinearWorld, const G4double FCof, const G4double tilt=0)
 

Static Public Attributes

static G4double thinElementLength = -1
 
static G4double nominalMatrixRelativeMomCut = -1
 
static G4bool currentTrackIsPrimary = false
 

Protected Member Functions

void ConvertToGlobal (const G4ThreeVector &localPos, const G4ThreeVector &localMom, G4double yOut[], G4double yErr[], const G4double momScaling=1.0)
 scaling of momentum in case localMom is a unit vector More...
 
void SetDistChord (G4double distChordIn)
 Setter for distChord to private member. More...
 

Protected Attributes

G4Mag_EqRhs * eqOfM
 
const G4int nVariables
 Cache of the number of variables. More...
 
G4MagIntegratorStepper * backupStepper
 
G4bool zeroStrength
 
G4double backupStepperMomLimit
 
- Protected Attributes inherited from BDSAuxiliaryNavigator
G4AffineTransform globalToLocal
 
G4AffineTransform localToGlobal
 
G4AffineTransform globalToLocalCL
 
G4AffineTransform localToGlobalCL
 
G4bool bridgeVolumeWasUsed
 

Private Member Functions

 BDSIntegratorMag ()=delete
 Private default constructor to force use of specific constructor.
 

Private Attributes

G4double distChordPrivate
 Variable used to record the distance from the chord calculated during the step. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from BDSAuxiliaryNavigator
static void AttachWorldVolumeToNavigator (G4VPhysicalVolume *worldPVIn)
 Setup the navigator w.r.t. to a world volume - typically real world. More...
 
static void AttachWorldVolumeToNavigatorCL (G4VPhysicalVolume *curvilinearWorldPVIn)
 
static void RegisterCurvilinearBridgeWorld (G4VPhysicalVolume *curvilinearBridgeWorldPVIn)
 
static void ResetNavigatorStates ()
 
- Static Protected Attributes inherited from BDSAuxiliaryNavigator
static G4Navigator * auxNavigator = new G4Navigator()
 
static G4Navigator * auxNavigatorCL = new G4Navigator()
 
static G4Navigator * auxNavigatorCLB = new G4Navigator()
 

Detailed Description

Common functionality to BDSIM integrators.

This provides a general 4th order Runge Kutta integrator that can be used by the derived class if the coordinates are non-paraxial for example.

The derived class must also satisfy G4MagIntegratorStepper.hh's virtual method Stepper.

Author
Laurie Nevay

Definition at line 43 of file BDSIntegratorMag.hh.

Constructor & Destructor Documentation

◆ BDSIntegratorMag()

BDSIntegratorMag::BDSIntegratorMag ( G4Mag_EqRhs *  eqOfMIn,
G4int  nVariablesIn 
)

Definition at line 32 of file BDSIntegratorMag.cc.

◆ ~BDSIntegratorMag()

BDSIntegratorMag::~BDSIntegratorMag ( )
virtual

Definition at line 50 of file BDSIntegratorMag.cc.

Member Function Documentation

◆ ConvertToGlobal()

void BDSIntegratorMag::ConvertToGlobal ( const G4ThreeVector &  localPos,
const G4ThreeVector &  localMom,
G4double  yOut[],
G4double  yErr[],
const G4double  momScaling = 1.0 
)
protected

scaling of momentum in case localMom is a unit vector

Convert final local position and direction to global frame. Allow

Definition at line 55 of file BDSIntegratorMag.cc.

References BDSAuxiliaryNavigator::ConvertToGlobalStep(), BDSStep::PostStepPoint(), and BDSStep::PreStepPoint().

Referenced by BDSIntegratorDipoleRodrigues::AdvanceHelix(), BDSIntegratorDecapole::AdvanceHelix(), BDSIntegratorOctupole::AdvanceHelix(), BDSIntegratorSextupole::AdvanceHelix(), BDSIntegratorKickerThin::Stepper(), and BDSIntegratorMultipoleThin::Stepper().

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

◆ DistChord()

virtual G4double BDSIntegratorMag::DistChord ( ) const
inlinevirtual

Estimate maximum distance of curved solution and chord.

Definition at line 56 of file BDSIntegratorMag.hh.

References distChordPrivate.

◆ IntegratorOrder()

virtual G4int BDSIntegratorMag::IntegratorOrder ( ) const
inlinevirtual

Geant4 requires that the integrator order must be supplied by the derived class.

Definition at line 59 of file BDSIntegratorMag.hh.

◆ SetDistChord()

void BDSIntegratorMag::SetDistChord ( G4double  distChordIn)
inlineprotected

Field Documentation

◆ backupStepper

G4MagIntegratorStepper* BDSIntegratorMag::backupStepper
protected

General integrator that can be used as a backup if the particle momentum is outside the (transverse) momentum range applicable for the integration scheme used by the derived integrator.

Definition at line 97 of file BDSIntegratorMag.hh.

Referenced by BDSIntegratorDipoleRodrigues::AdvanceHelix(), BDSIntegratorQuadrupole::Stepper(), BDSIntegratorSolenoid::Stepper(), and BDSIntegratorEulerOld::Stepper().

◆ backupStepperMomLimit

G4double BDSIntegratorMag::backupStepperMomLimit
protected

Cache of fractional unit momentum limit above which the matrix routines in the derived integrators revert to a backup integrator to ensure tracking validity.

Definition at line 105 of file BDSIntegratorMag.hh.

Referenced by BDSIntegratorDipoleRodrigues::AdvanceHelix(), BDSIntegratorDipoleQuadrupole::Stepper(), BDSIntegratorQuadrupole::Stepper(), BDSIntegratorSolenoid::Stepper(), BDSIntegratorKickerThin::Stepper(), and BDSIntegratorMultipoleThin::Stepper().

◆ currentTrackIsPrimary

G4bool BDSIntegratorMag::currentTrackIsPrimary = false
static

This static variable is updated by BDSFieldManager that marks each track as primary or not here. This variable is used throughout our integrators for magnetic fields which inherit this class.

Definition at line 72 of file BDSIntegratorMag.hh.

Referenced by BDSFieldManager::ConfigureForTrack(), and BDSTrackingAction::PreUserTrackingAction().

◆ distChordPrivate

G4double BDSIntegratorMag::distChordPrivate
private

Variable used to record the distance from the chord calculated during the step.

Definition at line 112 of file BDSIntegratorMag.hh.

Referenced by DistChord(), and SetDistChord().

◆ eqOfM

G4Mag_EqRhs* BDSIntegratorMag::eqOfM
protected

◆ nominalMatrixRelativeMomCut

G4double BDSIntegratorMag::nominalMatrixRelativeMomCut = -1
static

Cache of the fraction of the momentum outside which don't use a matrix as it's just not feasible.

Definition at line 67 of file BDSIntegratorMag.hh.

Referenced by BDSIntegratorDipoleQuadrupole::Stepper().

◆ nVariables

const G4int BDSIntegratorMag::nVariables
protected

Cache of the number of variables.

Definition at line 92 of file BDSIntegratorMag.hh.

Referenced by BDSIntegratorDipoleRodrigues::Stepper(), and BDSIntegratorEulerOld::Stepper().

◆ thinElementLength

G4double BDSIntegratorMag::thinElementLength = -1
static

Cache of thin element length to know maximum possible length scale step for coordinate lookup.

Definition at line 63 of file BDSIntegratorMag.hh.

Referenced by BDSIntegratorKickerThin::Stepper(), and BDSIntegratorMultipoleThin::Stepper().

◆ zeroStrength

G4bool BDSIntegratorMag::zeroStrength
protected

Whether a magnet has a finite strength or not. Can be set in the constructor for zero strength elements and then a drift routine is used before anything else.

Definition at line 101 of file BDSIntegratorMag.hh.

Referenced by BDSIntegratorDipoleRodrigues::AdvanceHelix(), BDSIntegratorDipoleQuadrupole::Stepper(), BDSIntegratorQuadrupole::Stepper(), BDSIntegratorSolenoid::Stepper(), BDSIntegratorEulerOld::Stepper(), BDSIntegratorKickerThin::Stepper(), BDSIntegratorMultipoleThin::Stepper(), and BDSIntegratorEuler::Stepper().


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