BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
|
Stepper that calculates trajectory through uniform magnetic field. More...
#include <BDSIntegratorDipoleRodrigues.hh>
Public Member Functions | |
BDSIntegratorDipoleRodrigues (BDSMagnetStrength const *strength, G4double brho, G4Mag_EqRhs *eqOfMIn) | |
virtual void | Stepper (const G4double yIn[], const G4double dydx[], const G4double h, G4double yOut[], G4double yErr[]) |
![]() | |
BDSIntegratorMag (G4Mag_EqRhs *eqOfMIn, G4int nVariablesIn) | |
virtual G4double | DistChord () const |
Estimate maximum distance of curved solution and chord. | |
virtual G4int | IntegratorOrder () const |
Geant4 requires that the integrator order must be supplied by the derived class. | |
BDSIntegratorMag & | operator= (const BDSIntegratorMag &)=delete |
Assignment and copy constructor not implemented nor used. | |
BDSIntegratorMag (BDSIntegratorMag &)=delete | |
Assignment and copy constructor not implemented nor used. | |
![]() | |
void | AdvanceDriftMag (const G4double yIn[], const G4double h, G4double yOut[], G4double yErr[]) const |
Error array [6] all 0. | |
void | AdvanceDriftMag (const G4double yIn[], const G4ThreeVector &unitMomentum, const G4double h, G4double yOut[], G4double yErr[]) const |
![]() | |
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. | |
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. | |
G4ThreeVector | ConvertToLocalNoSetup (const G4ThreeVector &globalPosition, const G4bool useCurvilinear=true) const |
Similar to above function but does NOT initialise the transforms. | |
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. | |
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 Member Functions | |
void | AdvanceHelix (const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[], G4double yErr[]) |
Calculate the new particle coordinates. | |
std::pair< G4ThreeVector, G4ThreeVector > | UpdatePandR (G4double rho, G4double h, G4ThreeVector localPos, G4ThreeVector localMomUnit) |
![]() | |
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 | |
void | SetDistChord (G4double distChordIn) |
Setter for distChord to private member. | |
Protected Attributes | |
G4double | cOverGeV |
Scaling factor in brho calculation. | |
const G4double | angle |
Angle that the dipole induces in the reference trajectory. | |
const G4double | length |
Arc length of the magnetic field. | |
G4double | bField |
Uniform magnetic field in global Y direction. | |
BDSMagnetStrength const * | strength |
Magnet strength object needed by curvilinear transforms. | |
![]() | |
G4Mag_EqRhs * | eqOfM |
const G4int | nVariables |
Cache of the number of variables. | |
G4MagIntegratorStepper * | backupStepper |
G4bool | zeroStrength |
G4double | backupStepperMomLimit |
![]() | |
G4AffineTransform | globalToLocal |
G4AffineTransform | localToGlobal |
G4AffineTransform | globalToLocalCL |
G4AffineTransform | localToGlobalCL |
G4bool | bridgeVolumeWasUsed |
Private Attributes | |
G4ThreeVector | yInitial |
Data stored in order to find the chord. | |
G4ThreeVector | yMidPoint |
G4ThreeVector | yFinal |
const G4double | minimumRadiusOfCurvature |
Minimum tolerable radius of curvature - used to prevent spiraling particles. | |
Additional Inherited Members | |
![]() | |
static void | AttachWorldVolumeToNavigator (G4VPhysicalVolume *worldPVIn) |
Setup the navigator w.r.t. to a world volume - typically real world. | |
static void | AttachWorldVolumeToNavigatorCL (G4VPhysicalVolume *curvilinearWorldPVIn) |
static void | RegisterCurvilinearBridgeWorld (G4VPhysicalVolume *curvilinearBridgeWorldPVIn) |
static void | ResetNavigatorStates () |
![]() | |
static G4double | thinElementLength = -1 |
static G4double | nominalMatrixRelativeMomCut = -1 |
static G4bool | currentTrackIsPrimary = false |
![]() | |
static G4Navigator * | auxNavigator = new G4Navigator() |
static G4Navigator * | auxNavigatorCL = new G4Navigator() |
static G4Navigator * | auxNavigatorCLB = new G4Navigator() |
Stepper that calculates trajectory through uniform magnetic field.
This calculates the coordinates of a particle moving through a dipole field. This ignores the field value supplied from Geant4 and uses its own field and quadrupole gradient from the BDSMagnetStrength instance supplied w.r.t. the rigidity also supplied.
Originally part of BDSIM by many authors.
Definition at line 43 of file BDSIntegratorDipoleRodrigues.hh.
BDSIntegratorDipoleRodrigues::BDSIntegratorDipoleRodrigues | ( | BDSMagnetStrength const * | strength, |
G4double | brho, | ||
G4Mag_EqRhs * | eqOfMIn | ||
) |
Definition at line 35 of file BDSIntegratorDipoleRodrigues.cc.
|
inlinevirtual |
Definition at line 50 of file BDSIntegratorDipoleRodrigues.hh.
|
protected |
Calculate the new particle coordinates.
Definition at line 52 of file BDSIntegratorDipoleRodrigues.cc.
References BDSIntegratorDrift::AdvanceDriftMag(), BDSIntegratorMag::backupStepper, BDSIntegratorMag::backupStepperMomLimit, bField, BDSIntegratorMag::ConvertToGlobal(), BDSAuxiliaryNavigator::ConvertToLocal(), cOverGeV, BDSIntegratorMag::eqOfM, BDS::IsFiniteStrength(), minimumRadiusOfCurvature, BDSStep::PostStepPoint(), BDSStep::PreStepPoint(), BDSIntegratorMag::SetDistChord(), and BDSIntegratorMag::zeroStrength.
Referenced by Stepper().
|
virtual |
Stepper for this integrator. Calculates the new coordinates of a particle through a uniform magnetic field.
Definition at line 135 of file BDSIntegratorDipoleRodrigues.cc.
References AdvanceHelix(), and BDSIntegratorMag::nVariables.
|
protected |
Definition at line 148 of file BDSIntegratorDipoleRodrigues.cc.
|
protected |
Angle that the dipole induces in the reference trajectory.
Definition at line 72 of file BDSIntegratorDipoleRodrigues.hh.
|
protected |
Uniform magnetic field in global Y direction.
Definition at line 78 of file BDSIntegratorDipoleRodrigues.hh.
Referenced by AdvanceHelix().
|
protected |
Scaling factor in brho calculation.
Definition at line 69 of file BDSIntegratorDipoleRodrigues.hh.
Referenced by AdvanceHelix().
|
protected |
Arc length of the magnetic field.
Definition at line 75 of file BDSIntegratorDipoleRodrigues.hh.
|
private |
Minimum tolerable radius of curvature - used to prevent spiraling particles.
Definition at line 93 of file BDSIntegratorDipoleRodrigues.hh.
Referenced by AdvanceHelix().
|
protected |
Magnet strength object needed by curvilinear transforms.
Definition at line 81 of file BDSIntegratorDipoleRodrigues.hh.
|
private |
Definition at line 90 of file BDSIntegratorDipoleRodrigues.hh.
|
private |
Data stored in order to find the chord.
Definition at line 90 of file BDSIntegratorDipoleRodrigues.hh.
|
private |
Definition at line 90 of file BDSIntegratorDipoleRodrigues.hh.