BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
|
Find multipole components of fieldmaps by numerical differentiation. More...
#include <BDSFieldMagGradient.hh>
Public Member Functions | |
BDSMagnetStrength * | CalculateMultipoles (BDSFieldMag *BField, G4int order, G4double Brho=4) |
Find the Magnet strengths of a multipole field to nth order. | |
G4double | GetBy (BDSFieldMag *BField, G4double x, G4double y=0) const |
Query a point on x axis (by default) for By. | |
std::vector< G4double > | PrepareValues (BDSFieldMag *field, G4int order, G4double centreX, G4double h, G4int ¢reIndex) const |
std::vector< std::vector< G4double > > | PrepareSkewValues (BDSFieldMag *field, G4int order, G4double centreX, G4double h, G4int ¢reIndex) const |
G4double | Derivative (const std::vector< G4double > &data, const G4int order, const G4int startIndex, const G4double h) const |
Find multipole components of fieldmaps by numerical differentiation.
Not optimised right now. Uses the Central difference method with each level taking the previous as the function (i.e. 1st order takes b, 2nd order takes 1st). Step size works right now but dropping it too low introduces large, inconsistent errors.
Definition at line 42 of file BDSFieldMagGradient.hh.
BDSFieldMagGradient::BDSFieldMagGradient | ( | ) |
Definition at line 34 of file BDSFieldMagGradient.cc.
BDSFieldMagGradient::~BDSFieldMagGradient | ( | ) |
Definition at line 37 of file BDSFieldMagGradient.cc.
BDSMagnetStrength * BDSFieldMagGradient::CalculateMultipoles | ( | BDSFieldMag * | BField, |
G4int | order, | ||
G4double | Brho = 4 |
||
) |
Find the Magnet strengths of a multipole field to nth order.
Definition at line 48 of file BDSFieldMagGradient.cc.
References Derivative(), PrepareSkewValues(), and PrepareValues().
Referenced by BDSFieldLoader::LoadMagField().
G4double BDSFieldMagGradient::Derivative | ( | const std::vector< G4double > & | data, |
const G4int | order, | ||
const G4int | startIndex, | ||
const G4double | h | ||
) | const |
Returns 'order' derivative of data about startIndex in data with assumed spacing of samples of h (in Geant4 units). Derivative returned is in SI T/(m^order).
Definition at line 107 of file BDSFieldMagGradient.cc.
References Derivative().
Referenced by CalculateMultipoles(), and Derivative().
G4double BDSFieldMagGradient::GetBy | ( | BDSFieldMag * | BField, |
G4double | x, | ||
G4double | y = 0 |
||
) | const |
Query a point on x axis (by default) for By.
Definition at line 40 of file BDSFieldMagGradient.cc.
References BDSFieldMag::GetField().
Referenced by PrepareValues().
std::vector< std::vector< G4double > > BDSFieldMagGradient::PrepareSkewValues | ( | BDSFieldMag * | field, |
G4int | order, | ||
G4double | centreX, | ||
G4double | h, | ||
G4int & | centreIndex | ||
) | const |
Similar to PrepareValues but for skew fields. For each order, rotate field by CLHEP::pi / (2*order + 2) - ie quadrupole = order 1 -> pi/4. For each rotation provide the vector of values to calculate gradients on.
Definition at line 90 of file BDSFieldMagGradient.cc.
References PrepareValues().
Referenced by CalculateMultipoles().
std::vector< G4double > BDSFieldMagGradient::PrepareValues | ( | BDSFieldMag * | field, |
G4int | order, | ||
G4double | centreX, | ||
G4double | h, | ||
G4int & | centreIndex | ||
) | const |
Prepare a vector of values corresponding to +- (2*order) in units of h about (spatial coordinate) centreX. centreIndex is passed by reference to write out which index represents the centre of the vector with x = centreX.
Definition at line 74 of file BDSFieldMagGradient.cc.
References GetBy().
Referenced by CalculateMultipoles(), and PrepareSkewValues().