BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
|
Integrator for combined dipole and quadrupolar field. More...
#include <BDSIntegratorDipoleQuadrupole.hh>
Public Member Functions | |
BDSIntegratorDipoleQuadrupole (BDSMagnetStrength const *strength, G4double brhoIn, G4Mag_EqRhs *eqOfMIn, G4double minimumRadiusOfCurvatureIn, const BDSParticleDefinition *designParticle, const G4double &tiltIn) | |
virtual void | Stepper (const G4double y[6], const G4double dydx[6], const G4double h, G4double yOut[6], G4double yErr[6]) |
BDSIntegratorDipoleQuadrupole & | operator= (const BDSIntegratorDipoleQuadrupole &)=delete |
Assignment and copy constructor not implemented nor used. | |
BDSIntegratorDipoleQuadrupole (BDSIntegratorDipoleQuadrupole &)=delete | |
Assignment and copy constructor not implemented nor used. | |
![]() | |
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 | OneStep (const G4ThreeVector &posIn, const G4ThreeVector &momIn, const G4ThreeVector &momUIn, const G4double &h, const G4double &fcof, G4ThreeVector &posOut, G4ThreeVector &momOut, const G4double rho, const G4double beta, const G4double deltaEnergy) const |
![]() | |
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. | |
Private Member Functions | |
BDSIntegratorDipoleQuadrupole ()=delete | |
Private default constructor to enforce use of supplied constructor. | |
Private Attributes | |
const G4double | nominalBRho |
Cached magnet property, nominal magnetic rigidity. | |
BDSMagUsualEqRhs * | eq |
BDSIM's eqRhs class to give access to particle properties. | |
const G4double | bPrime |
Cached magnet property, B field gradient for calculating K1. | |
const G4double | nominalBeta |
Cached nominal relativistic beta of the nominal beam particle. | |
G4double | nominalRho |
Cached magnet property, nominal bending radius. | |
const G4double | nominalField |
Cached magnet property, nominal field strength. | |
G4double | fieldRatio |
Ratio of supplied field to nominal field. Needed for over/underpowered magnets. | |
const G4double | nominalEnergy |
Nominal beam energy. | |
const G4double | nominalMass |
Primary particle mass. Needed for recalculating nominal energy with scaling. | |
const G4double | nominalCharge |
Nominal beam charge. | |
G4ThreeVector | unitField |
Cache of the unit field direction. | |
const G4double | fieldArcLength |
Cache of the field arc length. | |
const G4double | nominalAngle |
Cache of the field angle. | |
G4double | angleForCL |
Angle used for curvilinear transforms. | |
G4double | tilt |
Tilt offset transform for field. | |
const G4double | scaling |
Cache field scaling factor. | |
G4bool | isScaled |
Cache if field is scaled. | |
BDSIntegratorDipoleRodrigues2 * | dipole |
Backup integrator. | |
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 |
![]() | |
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 |
![]() | |
static G4Navigator * | auxNavigator = new G4Navigator() |
static G4Navigator * | auxNavigatorCL = new G4Navigator() |
static G4Navigator * | auxNavigatorCLB = new G4Navigator() |
Integrator for combined dipole and quadrupolar field.
Definition at line 39 of file BDSIntegratorDipoleQuadrupole.hh.
BDSIntegratorDipoleQuadrupole::BDSIntegratorDipoleQuadrupole | ( | BDSMagnetStrength const * | strength, |
G4double | brhoIn, | ||
G4Mag_EqRhs * | eqOfMIn, | ||
G4double | minimumRadiusOfCurvatureIn, | ||
const BDSParticleDefinition * | designParticle, | ||
const G4double & | tiltIn | ||
) |
Definition at line 42 of file BDSIntegratorDipoleQuadrupole.cc.
|
virtual |
Definition at line 80 of file BDSIntegratorDipoleQuadrupole.cc.
|
protected |
Calculate a single step in curvilinear coordinates using dipole quadrupole matrix. Unit momentum is provided as an argument becuase it is already calculated in the Stepper method.
Definition at line 261 of file BDSIntegratorDipoleQuadrupole.cc.
References bPrime, and nominalBRho.
Referenced by Stepper().
|
virtual |
Check if the quadrupole has finite strength and use drift if not. If finite strength, convert to local curvilinear coordiantes and check for paraxial approximation. If paraxial, use thick quadrupole matrix for transport, else use the G4ClassicalRK4 backup stepper.
Definition at line 85 of file BDSIntegratorDipoleQuadrupole.cc.
References BDSIntegratorDrift::AdvanceDriftMag(), angleForCL, BDSIntegratorMag::backupStepperMomLimit, dipole, eq, BDSIntegratorMag::eqOfM, fieldArcLength, fieldRatio, BDSAuxiliaryNavigator::GlobalToCurvilinear(), BDS::IsFinite(), BDS::IsFiniteStrength(), isScaled, nominalAngle, nominalBeta, nominalBRho, nominalCharge, nominalEnergy, nominalMass, BDSIntegratorMag::nominalMatrixRelativeMomCut, nominalRho, OneStep(), BDSStep::PostStepPoint(), BDSStep::PreStepPoint(), BDSIntegratorDipoleRodrigues2::RadiusOfHelix(), scaling, BDSIntegratorMag::SetDistChord(), BDSIntegratorDipoleRodrigues2::SingleStep(), BDSIntegratorDipoleRodrigues2::Stepper(), tilt, BDSMagUsualEqRhs::TotalEnergy(), unitField, and BDSIntegratorMag::zeroStrength.
|
private |
Angle used for curvilinear transforms.
Definition at line 100 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Cached magnet property, B field gradient for calculating K1.
Definition at line 89 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by OneStep().
|
private |
Backup integrator.
Definition at line 105 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
BDSIM's eqRhs class to give access to particle properties.
Definition at line 88 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Cache of the field arc length.
Definition at line 98 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Ratio of supplied field to nominal field. Needed for over/underpowered magnets.
Definition at line 93 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Cache if field is scaled.
Definition at line 103 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Cache of the field angle.
Definition at line 99 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Cached nominal relativistic beta of the nominal beam particle.
Definition at line 90 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Cached magnet property, nominal magnetic rigidity.
Definition at line 87 of file BDSIntegratorDipoleQuadrupole.hh.
|
private |
Nominal beam charge.
Definition at line 96 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Nominal beam energy.
Definition at line 94 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Cached magnet property, nominal field strength.
Definition at line 92 of file BDSIntegratorDipoleQuadrupole.hh.
|
private |
Primary particle mass. Needed for recalculating nominal energy with scaling.
Definition at line 95 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Cached magnet property, nominal bending radius.
Definition at line 91 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Cache field scaling factor.
Definition at line 102 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Tilt offset transform for field.
Definition at line 101 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().
|
private |
Cache of the unit field direction.
Definition at line 97 of file BDSIntegratorDipoleQuadrupole.hh.
Referenced by Stepper().