BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
|
A base class for magnetic fields in local to be used in global coordinates. More...
#include <BDSFieldMagGlobal.hh>
Public Member Functions | |
BDSFieldMagGlobal (BDSFieldMag *fieldIn) | |
virtual G4ThreeVector | GetFieldTransformed (const G4ThreeVector &position, const G4double t) const |
virtual G4ThreeVector | GetField (const G4ThreeVector &position, const G4double t=0) const |
![]() | |
BDSFieldMag () | |
BDSFieldMag (G4Transform3D transformIn) | |
virtual G4ThreeVector | GetField (const G4ThreeVector &position, const G4double t=0) const =0 |
virtual G4bool | TimeVarying () const |
virtual void | GetFieldValue (const G4double point[4], G4double *field) const |
virtual G4ThreeVector | GetFieldTransformed (const G4ThreeVector &position, const G4double t) const |
Get the field value after applying transform for local offset. | |
virtual void | SetTransform (const G4Transform3D &transformIn) |
void | SetModulator (BDSModulator *modulatorIn) |
Set the optional modulator. | |
G4bool | FiniteStrength () const |
Accessor. | |
![]() | |
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) |
Private Attributes | |
BDSFieldMag * | field |
The field on which this is based. | |
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 () |
![]() | |
G4bool | finiteStrength |
Flag to cache whether finite nor not. | |
G4Transform3D | transform |
Transform to apply for the field relative to the local coordinates of the geometry. | |
G4bool | transformIsNotIdentity |
Cache of whether to use transform at all. | |
BDSModulator * | modulator |
Optional modulator;. | |
![]() | |
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() |
A base class for magnetic fields in local to be used in global coordinates.
This base class provides the aggregate inheritance and utility functions for magnetic fields in local coordinates to be used in global coordinates.
Constness is particularly important here as member functions are called from inside G4MagneticField::GetField function which is const.
This owns the field it wraps.
Definition at line 42 of file BDSFieldMagGlobal.hh.
|
explicit |
Definition at line 26 of file BDSFieldMagGlobal.cc.
|
virtual |
Definition at line 32 of file BDSFieldMagGlobal.cc.
|
virtual |
Apply the global to local transform, query the wrapped field object and transform this field to global coordinates before returning.
Implements BDSFieldMag.
Definition at line 46 of file BDSFieldMagGlobal.cc.
References BDSAuxiliaryNavigator::ConvertAxisToGlobal(), BDSAuxiliaryNavigator::ConvertToLocal(), field, and BDSFieldMag::GetFieldTransformed().
Referenced by GetFieldTransformed().
|
virtual |
As we use a discrete member field object, we do not need to apply the transform. Override default method and just directly call GetField().
Reimplemented from BDSFieldMag.
Definition at line 37 of file BDSFieldMagGlobal.cc.
References BDSFieldMag::finiteStrength, and GetField().
|
private |
The field on which this is based.
Definition at line 61 of file BDSFieldMagGlobal.hh.
Referenced by GetField().