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

Class that provides the magnetic field due to a cylinder of current. More...

#include <BDSFieldMagSolenoidSheet.hh>

Inheritance diagram for BDSFieldMagSolenoidSheet:
Inheritance graph
Collaboration diagram for BDSFieldMagSolenoidSheet:
Collaboration graph

Public Member Functions

 BDSFieldMagSolenoidSheet (BDSMagnetStrength const *strength, G4double radiusIn)
 
 BDSFieldMagSolenoidSheet (G4double fullLength, G4double sheetRadius, G4double B0)
 More reasonable constructor for the internal parameterisation. More...
 
virtual G4ThreeVector GetField (const G4ThreeVector &position, const G4double t=0) const
 Calculate the field value. More...
 
- Public Member Functions inherited from BDSFieldMag
 BDSFieldMag ()
 
 BDSFieldMag (G4Transform3D transformIn)
 
virtual G4ThreeVector GetField (const G4ThreeVector &position, const G4double t=0) const =0
 
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. More...
 
virtual void SetTransform (const G4Transform3D &transformIn)
 
G4bool FiniteStrength () const
 Accessor. More...
 

Static Public Member Functions

static G4double CEL (G4double kc, G4double p, G4double c, G4double s, G4int nIterationLimit=1000)
 Generalised Complete Elliptical Integral. More...
 

Private Member Functions

G4double OnAxisBz (G4double zp, G4double zm) const
 

Private Attributes

G4double a
 
G4double halfLength
 
G4double B0
 
G4double spatialLimit
 
G4double normalisation
 

Additional Inherited Members

- Protected Attributes inherited from BDSFieldMag
G4bool finiteStrength
 Flag to cache whether finite nor not. More...
 
G4Transform3D transform
 Transform to apply for the field relative to the local coordinates of the geometry. More...
 

Detailed Description

Class that provides the magnetic field due to a cylinder of current.

This follows the parameterisation and uses the algorithm for the generalised complete elliptical integral as described in:

Cylindrical Magnets and ideal Solenoids, N. Derby and S. Olbert, American Journal of Physics 78, 229 (2010); https://doi.org/10.1119/1.3256157 and also at https://arxiv.org/abs/0909.3880.

The field is calculated in cylindrical coordinates. A complete description is in the manual.

Author
Laurie Nevay

Definition at line 44 of file BDSFieldMagSolenoidSheet.hh.

Constructor & Destructor Documentation

◆ BDSFieldMagSolenoidSheet() [1/2]

BDSFieldMagSolenoidSheet::BDSFieldMagSolenoidSheet ( BDSMagnetStrength const *  strength,
G4double  radiusIn 
)

This constructor uses the "field" and "length" parameters from the BDSMagnetStrength instance and forwards to the next constructor.

Definition at line 33 of file BDSFieldMagSolenoidSheet.cc.

◆ BDSFieldMagSolenoidSheet() [2/2]

BDSFieldMagSolenoidSheet::BDSFieldMagSolenoidSheet ( G4double  fullLength,
G4double  sheetRadius,
G4double  B0 
)

More reasonable constructor for the internal parameterisation.

Definition at line 38 of file BDSFieldMagSolenoidSheet.cc.

References BDSFieldMag::finiteStrength, BDS::IsFinite(), and OnAxisBz().

Here is the call graph for this function:

◆ ~BDSFieldMagSolenoidSheet()

virtual BDSFieldMagSolenoidSheet::~BDSFieldMagSolenoidSheet ( )
inlinevirtual

Definition at line 57 of file BDSFieldMagSolenoidSheet.hh.

Member Function Documentation

◆ CEL()

G4double BDSFieldMagSolenoidSheet::CEL ( G4double  kc,
G4double  p,
G4double  c,
G4double  s,
G4int  nIterationLimit = 1000 
)
static

Generalised Complete Elliptical Integral.

Definition at line 119 of file BDSFieldMagSolenoidSheet.cc.

References BDS::IsFinite().

Referenced by GetField().

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

◆ GetField()

G4ThreeVector BDSFieldMagSolenoidSheet::GetField ( const G4ThreeVector &  position,
const G4double  t = 0 
) const
virtual

Calculate the field value.

Implements BDSFieldMag.

Definition at line 56 of file BDSFieldMagSolenoidSheet.cc.

References CEL(), and OnAxisBz().

Here is the call graph for this function:

◆ OnAxisBz()

G4double BDSFieldMagSolenoidSheet::OnAxisBz ( G4double  zp,
G4double  zm 
) const
private

Approximation for rho=0 Bz field. Brho=0 by definition. zp and zm are z+halfLength and z-halfLength. Returns Bz.

Definition at line 177 of file BDSFieldMagSolenoidSheet.cc.

Referenced by BDSFieldMagSolenoidSheet(), and GetField().

Here is the caller graph for this function:

Field Documentation

◆ a

G4double BDSFieldMagSolenoidSheet::a
private

Definition at line 75 of file BDSFieldMagSolenoidSheet.hh.

◆ B0

G4double BDSFieldMagSolenoidSheet::B0
private

Definition at line 77 of file BDSFieldMagSolenoidSheet.hh.

◆ halfLength

G4double BDSFieldMagSolenoidSheet::halfLength
private

Definition at line 76 of file BDSFieldMagSolenoidSheet.hh.

◆ normalisation

G4double BDSFieldMagSolenoidSheet::normalisation
private

Definition at line 79 of file BDSFieldMagSolenoidSheet.hh.

◆ spatialLimit

G4double BDSFieldMagSolenoidSheet::spatialLimit
private

Definition at line 78 of file BDSFieldMagSolenoidSheet.hh.


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