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

Integrator that ignores the field and uses the analytical solution for a dipole kick. More...

#include <BDSIntegratorDipoleFringe.hh>

Inheritance diagram for BDSIntegratorDipoleFringe:
Inheritance graph
Collaboration diagram for BDSIntegratorDipoleFringe:
Collaboration graph

Public Member Functions

 BDSIntegratorDipoleFringe (BDSMagnetStrength const *strength, G4double brhoIn, G4Mag_EqRhs *eqOfMIn, G4double minimumRadiusOfCurvature, const G4double &tiltIn=0)
 
virtual void Stepper (const G4double yIn[6], const G4double dydx[6], const G4double h, G4double yOut[6], G4double yErr[6])
 The stepper for integration. Calls base class stepper. More...
 
void BaseStepper (const G4double yIn[6], const G4double dydx[6], const G4double &h, G4double yOut[6], G4double yErr[6], const G4double &fcof, const G4double &momScaling)
 
void OneStep (const G4ThreeVector &posIn, const G4ThreeVector &momUIn, G4ThreeVector &posOut, G4ThreeVector &momOut, const G4double &bendingRadius) const
 
void MultipoleStep (const G4double yIn[6], G4double yMultipoleOut[6], const G4double &h)
 
G4double GetPolefaceAngle ()
 
G4double GetFringeCorr ()
 
G4double GetSecondFringeCorr ()
 
G4double GetPolefaceCurv ()
 
BDSIntegratorDipoleFringeoperator= (const BDSIntegratorDipoleFringe &)=delete
 Assignment and copy constructor not implemented nor used.
 
 BDSIntegratorDipoleFringe (BDSIntegratorDipoleFringe &)=delete
 Assignment and copy constructor not implemented nor used.
 
- Public Member Functions inherited from BDSIntegratorDipoleRodrigues2
 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
 
- 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)
 

Protected Attributes

G4double polefaceAngle
 Poleface rotation angle. More...
 
G4double fringeCorr
 Fringe field correction term. More...
 
G4double secondFringeCorr
 Second fringe field correction term. More...
 
G4double polefaceCurvature
 Poleface curvature. More...
 
G4double rho
 Nominal magnet bending radius. More...
 
G4ThreeVector unitField
 Cache of the unit field direction. More...
 
const G4double fieldArcLength
 Cache of the field arc length. More...
 
const G4double fieldAngle
 Cache of the field angle. More...
 
const G4double tilt
 
const G4bool finiteTilt
 
G4double bx
 
G4double by
 
G4bool zeroStrength
 
BDSIntegratorMultipoleThinmultipoleIntegrator
 
G4bool isEntrance
 store if fringe is entrance or exit More...
 
G4double backupStepperMomLimit
 
- Protected Attributes inherited from BDSIntegratorDipoleRodrigues2
G4Mag_EqRhs * eqOfM
 Cache of equation of motion. This class does not own it. More...
 
- Protected Attributes inherited from BDSAuxiliaryNavigator
G4AffineTransform globalToLocal
 
G4AffineTransform localToGlobal
 
G4AffineTransform globalToLocalCL
 
G4AffineTransform localToGlobalCL
 
G4bool bridgeVolumeWasUsed
 

Static Protected Attributes

static G4double thinElementLength = -1
 
- Static Protected Attributes inherited from BDSAuxiliaryNavigator
static G4Navigator * auxNavigator = new G4Navigator()
 
static G4Navigator * auxNavigatorCL = new G4Navigator()
 
static G4Navigator * auxNavigatorCLB = new G4Navigator()
 

Private Member Functions

 BDSIntegratorDipoleFringe ()=delete
 Private default constructor to enforce use of supplied constructor.
 

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 ()
 
- Protected Member Functions inherited from BDSIntegratorDipoleRodrigues2
void FudgeDistChordToZero ()
 

Detailed Description

Integrator that ignores the field and uses the analytical solution for a dipole kick.

Author
Will Shields

Definition at line 47 of file BDSIntegratorDipoleFringe.hh.

Constructor & Destructor Documentation

◆ BDSIntegratorDipoleFringe()

BDSIntegratorDipoleFringe::BDSIntegratorDipoleFringe ( BDSMagnetStrength const *  strength,
G4double  brhoIn,
G4Mag_EqRhs *  eqOfMIn,
G4double  minimumRadiusOfCurvature,
const G4double &  tiltIn = 0 
)

Definition at line 34 of file BDSIntegratorDipoleFringe.cc.

◆ ~BDSIntegratorDipoleFringe()

BDSIntegratorDipoleFringe::~BDSIntegratorDipoleFringe ( )
virtual

Definition at line 102 of file BDSIntegratorDipoleFringe.cc.

Member Function Documentation

◆ BaseStepper()

void BDSIntegratorDipoleFringe::BaseStepper ( const G4double  yIn[6],
const G4double  dydx[6],
const G4double &  h,
G4double  yOut[6],
G4double  yErr[6],
const G4double &  fcof,
const G4double &  momScaling 
)

The stepper for integration. Particle charge and momentum scaling are provided as arguments as they are calculated in the derived classes. Uses BDSIntegratorDipole2::Stepper and then adds a kick in yp in curvilinear coordinates.

Definition at line 119 of file BDSIntegratorDipoleFringe.cc.

References BDSIntegratorDrift::AdvanceDriftMag(), backupStepperMomLimit, BDSAuxiliaryNavigator::ConvertAxisToGlobal(), BDSIntegratorDipoleRodrigues2::FudgeDistChordToZero(), BDSAuxiliaryNavigator::GlobalToCurvilinear(), isEntrance, BDS::IsFinite(), BDS::IsFiniteStrength(), MultipoleStep(), OneStep(), BDSStep::PostStepPoint(), BDSStep::PreStepPoint(), rho, BDSStep::rotateZ(), BDSIntegratorDipoleRodrigues2::Stepper(), thinElementLength, and zeroStrength.

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

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

◆ GetFringeCorr()

G4double BDSIntegratorDipoleFringe::GetFringeCorr ( )
inline

Definition at line 100 of file BDSIntegratorDipoleFringe.hh.

◆ GetPolefaceAngle()

G4double BDSIntegratorDipoleFringe::GetPolefaceAngle ( )
inline

Getter functions for poleface and fringe variables. Values need to be read in at least the thin kicker integrator, but should not be overwritten.

Definition at line 99 of file BDSIntegratorDipoleFringe.hh.

References polefaceAngle.

◆ GetPolefaceCurv()

G4double BDSIntegratorDipoleFringe::GetPolefaceCurv ( )
inline

Definition at line 102 of file BDSIntegratorDipoleFringe.hh.

◆ GetSecondFringeCorr()

G4double BDSIntegratorDipoleFringe::GetSecondFringeCorr ( )
inline

Definition at line 101 of file BDSIntegratorDipoleFringe.hh.

◆ MultipoleStep()

void BDSIntegratorDipoleFringe::MultipoleStep ( const G4double  yIn[6],
G4double  yMultipoleOut[6],
const G4double &  h 
)

Calculate the thin multipole kick to represent the dipole poleface curvature effect. Step length is passed in as it is needed by the transforms.

Definition at line 334 of file BDSIntegratorDipoleFringe.cc.

References BDSAuxiliaryNavigator::ConvertAxisToGlobal(), BDSAuxiliaryNavigator::GlobalToCurvilinear(), BDSIntegratorMultipoleThin::OneStep(), BDSStep::PostStepPoint(), and BDSStep::PreStepPoint().

Referenced by BaseStepper().

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

◆ OneStep()

void BDSIntegratorDipoleFringe::OneStep ( const G4ThreeVector &  posIn,
const G4ThreeVector &  momUIn,
G4ThreeVector &  posOut,
G4ThreeVector &  momOut,
const G4double &  bendingRadius 
) const

Calculate a single step using dipole fringe field matrix. Unit momentum, momentum magnitude, and normalised bending radius are provided as arguments because they already calculated in the BaseStepper method.

Definition at line 285 of file BDSIntegratorDipoleFringe.cc.

References fringeCorr, polefaceAngle, and secondFringeCorr.

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

Here is the caller graph for this function:

◆ Stepper()

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

The stepper for integration. Calls base class stepper.

Definition at line 107 of file BDSIntegratorDipoleFringe.cc.

References BaseStepper(), and BDSIntegratorDipoleRodrigues2::eqOfM.

Here is the call graph for this function:

Field Documentation

◆ backupStepperMomLimit

G4double BDSIntegratorDipoleFringe::backupStepperMomLimit
protected

Fractional momentum limit above which the thin matrix fringe kick accuracy becomes questionable, so advance with the dipolerodrigues2 stepper to ensure correct dipole tracking length is applied

Definition at line 138 of file BDSIntegratorDipoleFringe.hh.

Referenced by BaseStepper().

◆ bx

G4double BDSIntegratorDipoleFringe::bx
protected

Definition at line 122 of file BDSIntegratorDipoleFringe.hh.

◆ by

G4double BDSIntegratorDipoleFringe::by
protected

Definition at line 123 of file BDSIntegratorDipoleFringe.hh.

◆ fieldAngle

const G4double BDSIntegratorDipoleFringe::fieldAngle
protected

Cache of the field angle.

Definition at line 118 of file BDSIntegratorDipoleFringe.hh.

◆ fieldArcLength

const G4double BDSIntegratorDipoleFringe::fieldArcLength
protected

Cache of the field arc length.

Definition at line 117 of file BDSIntegratorDipoleFringe.hh.

◆ finiteTilt

const G4bool BDSIntegratorDipoleFringe::finiteTilt
protected

Definition at line 121 of file BDSIntegratorDipoleFringe.hh.

◆ fringeCorr

G4double BDSIntegratorDipoleFringe::fringeCorr
protected

Fringe field correction term.

Definition at line 108 of file BDSIntegratorDipoleFringe.hh.

Referenced by OneStep().

◆ isEntrance

G4bool BDSIntegratorDipoleFringe::isEntrance
protected

store if fringe is entrance or exit

Definition at line 134 of file BDSIntegratorDipoleFringe.hh.

Referenced by BaseStepper().

◆ multipoleIntegrator

BDSIntegratorMultipoleThin* BDSIntegratorDipoleFringe::multipoleIntegrator
protected

Definition at line 132 of file BDSIntegratorDipoleFringe.hh.

◆ polefaceAngle

G4double BDSIntegratorDipoleFringe::polefaceAngle
protected

Poleface rotation angle.

Definition at line 106 of file BDSIntegratorDipoleFringe.hh.

Referenced by GetPolefaceAngle(), and OneStep().

◆ polefaceCurvature

G4double BDSIntegratorDipoleFringe::polefaceCurvature
protected

Poleface curvature.

Definition at line 112 of file BDSIntegratorDipoleFringe.hh.

◆ rho

G4double BDSIntegratorDipoleFringe::rho
protected

Nominal magnet bending radius.

Definition at line 114 of file BDSIntegratorDipoleFringe.hh.

Referenced by BaseStepper().

◆ secondFringeCorr

G4double BDSIntegratorDipoleFringe::secondFringeCorr
protected

Second fringe field correction term.

Definition at line 110 of file BDSIntegratorDipoleFringe.hh.

Referenced by OneStep().

◆ thinElementLength

G4double BDSIntegratorDipoleFringe::thinElementLength = -1
staticprotected

Cache of thin element length from global constants. Initialised via check on unphysical -1 value as global constants doesn't exist at compile time.

Definition at line 130 of file BDSIntegratorDipoleFringe.hh.

Referenced by BaseStepper().

◆ tilt

const G4double BDSIntegratorDipoleFringe::tilt
protected

Definition at line 120 of file BDSIntegratorDipoleFringe.hh.

◆ unitField

G4ThreeVector BDSIntegratorDipoleFringe::unitField
protected

Cache of the unit field direction.

Definition at line 116 of file BDSIntegratorDipoleFringe.hh.

◆ zeroStrength

G4bool BDSIntegratorDipoleFringe::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 126 of file BDSIntegratorDipoleFringe.hh.

Referenced by BaseStepper().


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