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 | Static Private Attributes
BDSFieldFactory Class Reference

Factory that produces fields and their associated objects. More...

#include <BDSFieldFactory.hh>

Collaboration diagram for BDSFieldFactory:
Collaboration graph

Public Member Functions

BDSFieldObjectsCreateField (const BDSFieldInfo &info, const BDSMagnetStrength *scalingStrength=nullptr, const G4String &scalingKey="none")
 Main interface to field factory. More...
 
BDSFieldInfoGetDefinition (const G4String &name) const
 

Static Public Member Functions

static BDSFieldFactoryInstance ()
 Public accessor method for singleton pattern. More...
 
static void SetDesignParticle (const BDSParticleDefinition *designParticleIn)
 Update the internal cache of the rigidity. More...
 
static void SetPrimaryGeneratorAction (BDSPrimaryGeneratorAction *pgaIn)
 Update the internal cache of the primary generator action. More...
 
static BDSInterpolatorType DefaultInterpolatorType (G4int numberOfDimensions)
 Suggest a default interpolator. More...
 

Private Member Functions

BDSFieldObjectsCreateFieldMag (const BDSFieldInfo &info, const BDSMagnetStrength *scalingStrength=nullptr, const G4String &scalingKey="none")
 Create a purely magnetic field. More...
 
BDSFieldObjectsCreateFieldEM (const BDSFieldInfo &info)
 Create a general EM field. More...
 
BDSFieldObjectsCreateFieldE (const BDSFieldInfo &info)
 Create an electric field. More...
 
BDSFieldObjectsCreateFieldIrregular (const BDSFieldInfo &info)
 Create an irregular (special) field. More...
 
BDSFieldMagCreateFieldMagRaw (const BDSFieldInfo &info, const BDSMagnetStrength *scalingStrength=nullptr, const G4String &scalingKey="none")
 Creat just the magnetic field object. More...
 
BDSFieldECreateFieldERaw (const BDSFieldInfo &info)
 Creat just the electric field object. More...
 
G4MagIntegratorStepper * CreateIntegratorMag (const BDSFieldInfo &info, G4Mag_EqRhs *eqOfM, const BDSMagnetStrength *strength)
 
G4MagIntegratorStepper * CreateIntegratorEM (const BDSFieldInfo &info, G4EquationOfMotion *eqOfM)
 
G4MagIntegratorStepper * CreateIntegratorE (const BDSFieldInfo &info, G4EquationOfMotion *eqOfM)
 
BDSFieldObjectsCreateTeleporter (const BDSFieldInfo &info)
 
BDSFieldObjectsCreateRMatrix (const BDSFieldInfo &info)
 Create special rmatrix 'field' that applies an rmatrix. More...
 
BDSFieldObjectsCreateCavityFringe (const BDSFieldInfo &info)
 Create special rf cavity fringe 'field' that applies an rmatrix. More...
 
BDSFieldObjectsCreateParallelTransport (const BDSFieldInfo &info)
 
G4double GetOuterScaling (const BDSMagnetStrength *st) const
 Return the parameter "outerScaling" from strength st, but default to 1. More...
 
 BDSFieldFactory ()
 Private default constructor as singleton class. More...
 
void PrepareFieldDefinitions (const std::vector< GMAD::Field > &definitions, G4double defaultBRho)
 Prepare all required definitions that can be used dynamically. More...
 
G4double ConvertToDoubleWithException (const G4String &value, const G4String &parameterNameForError) const
 Convert the string 'value' to a double. Throw an exception including the parameterNameForError if it doesn't work. More...
 
void PrepareFieldStrengthFromParameters (BDSMagnetStrength *st, const G4String &fieldParameters, G4double &poleTipRadius) const
 

Private Attributes

std::map< G4String, BDSFieldInfo * > parserDefinitions
 BDSFieldInfo definitions prepare from parser vector of definitions. More...
 
G4bool useOldMultipoleOuterFields
 

Static Private Attributes

static BDSFieldFactoryinstance = nullptr
 Instance - singleton pattern. More...
 
static const BDSParticleDefinitiondesignParticle = nullptr
 Cache of design particle for fields. More...
 
static BDSPrimaryGeneratorActionprimaryGeneratorAction = nullptr
 Cache of primary generator action. More...
 

Detailed Description

Factory that produces fields and their associated objects.

Field objects are created according to a BDSFieldType and the associated and required Geant4 objects to properly implement a field. These are packaged together in one object. This factory does not own any of its products. Construction follows in this order:

1 field 2 equation of motion (based on field object) 3 integrator 4 chord finder 5 field manager 6 package it up

This also makes use of BDSParser singleton class to create a series of BDSFieldInfo field specifications as defined by the parser.

Author
Laurie Nevay

Definition at line 69 of file BDSFieldFactory.hh.

Constructor & Destructor Documentation

◆ ~BDSFieldFactory()

BDSFieldFactory::~BDSFieldFactory ( )

Definition at line 168 of file BDSFieldFactory.cc.

◆ BDSFieldFactory()

BDSFieldFactory::BDSFieldFactory ( )
private

Private default constructor as singleton class.

Definition at line 158 of file BDSFieldFactory.cc.

References BDSParticleDefinition::BRho(), designParticle, BDSGlobalConstants::Instance(), BDSParser::Instance(), and PrepareFieldDefinitions().

Referenced by Instance().

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

Member Function Documentation

◆ ConvertToDoubleWithException()

G4double BDSFieldFactory::ConvertToDoubleWithException ( const G4String &  value,
const G4String &  parameterNameForError 
) const
private

Convert the string 'value' to a double. Throw an exception including the parameterNameForError if it doesn't work.

Definition at line 358 of file BDSFieldFactory.cc.

Referenced by PrepareFieldStrengthFromParameters().

Here is the caller graph for this function:

◆ CreateCavityFringe()

BDSFieldObjects * BDSFieldFactory::CreateCavityFringe ( const BDSFieldInfo info)
private

Create special rf cavity fringe 'field' that applies an rmatrix.

Definition at line 1133 of file BDSFieldFactory.cc.

References BDSFieldInfo::BeamPipeRadius(), and BDSFieldInfo::MagnetStrength().

Referenced by CreateFieldIrregular().

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

◆ CreateField()

BDSFieldObjects * BDSFieldFactory::CreateField ( const BDSFieldInfo info,
const BDSMagnetStrength scalingStrength = nullptr,
const G4String &  scalingKey = "none" 
)

Main interface to field factory.

Definition at line 426 of file BDSFieldFactory.cc.

References CreateFieldE(), CreateFieldEM(), CreateFieldIrregular(), CreateFieldMag(), BDS::DetermineFieldClassType(), BDSFieldInfo::FieldType(), and BDSTypeSafeEnum< def, inner >::underlying().

Here is the call graph for this function:

◆ CreateFieldE()

BDSFieldObjects * BDSFieldFactory::CreateFieldE ( const BDSFieldInfo info)
private

Create an electric field.

Definition at line 814 of file BDSFieldFactory.cc.

References CreateFieldERaw(), CreateIntegratorE(), BDSFieldInfo::ProvideGlobal(), and BDSFieldInfo::UsePlacementWorldTransform().

Referenced by CreateField().

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

◆ CreateFieldEM()

BDSFieldObjects * BDSFieldFactory::CreateFieldEM ( const BDSFieldInfo info)
private

◆ CreateFieldERaw()

BDSFieldE * BDSFieldFactory::CreateFieldERaw ( const BDSFieldInfo info)
private

◆ CreateFieldIrregular()

BDSFieldObjects * BDSFieldFactory::CreateFieldIrregular ( const BDSFieldInfo info)
private

Create an irregular (special) field.

Definition at line 889 of file BDSFieldFactory.cc.

References CreateCavityFringe(), CreateParallelTransport(), CreateRMatrix(), CreateTeleporter(), BDSFieldInfo::FieldType(), and BDSTypeSafeEnum< def, inner >::underlying().

Referenced by CreateField().

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

◆ CreateFieldMag()

BDSFieldObjects * BDSFieldFactory::CreateFieldMag ( const BDSFieldInfo info,
const BDSMagnetStrength scalingStrength = nullptr,
const G4String &  scalingKey = "none" 
)
private

Create a purely magnetic field.

Definition at line 476 of file BDSFieldFactory.cc.

References CreateFieldMagRaw(), CreateIntegratorMag(), BDSFieldInfo::MagnetStrength(), BDSFieldInfo::ProvideGlobal(), and BDSFieldInfo::UsePlacementWorldTransform().

Referenced by CreateField().

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

◆ CreateFieldMagRaw()

BDSFieldMag * BDSFieldFactory::CreateFieldMagRaw ( const BDSFieldInfo info,
const BDSMagnetStrength scalingStrength = nullptr,
const G4String &  scalingKey = "none" 
)
private

◆ CreateIntegratorE()

G4MagIntegratorStepper * BDSFieldFactory::CreateIntegratorE ( const BDSFieldInfo info,
G4EquationOfMotion *  eqOfM 
)
private

Create an integrator for a general E field. Same ones as EM but keep this method for clarity as Geant4 unclear - only based on their examples. examples/extended/field/field02/src/F02ElectricFieldSetup.cc

Definition at line 1093 of file BDSFieldFactory.cc.

References CreateIntegratorEM().

Referenced by CreateFieldE().

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

◆ CreateIntegratorEM()

G4MagIntegratorStepper * BDSFieldFactory::CreateIntegratorEM ( const BDSFieldInfo info,
G4EquationOfMotion *  eqOfM 
)
private

Create an integrator for a general EM field. As it's a general field, this takes a G4EquationOfMotion* equation of motion instance.

Definition at line 997 of file BDSFieldFactory.cc.

References BDSFieldInfo::IntegratorType(), and BDSTypeSafeEnum< def, inner >::underlying().

Referenced by CreateFieldEM(), CreateIntegratorE(), and CreateIntegratorMag().

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

◆ CreateIntegratorMag()

G4MagIntegratorStepper * BDSFieldFactory::CreateIntegratorMag ( const BDSFieldInfo info,
G4Mag_EqRhs *  eqOfM,
const BDSMagnetStrength strength 
)
private

Create a purely magnetic integrator. As it's purely magnetic, this requires a G4Mag_EqRhs* equation of motion instance.

Definition at line 909 of file BDSFieldFactory.cc.

References BDSFieldInfo::BeamPipeRadius(), BDSFieldInfo::BRho(), CreateIntegratorEM(), designParticle, BDSGlobalConstants::Instance(), BDSFieldInfo::IntegratorType(), BDSFieldInfo::Tilt(), and BDSTypeSafeEnum< def, inner >::underlying().

Referenced by CreateFieldMag().

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

◆ CreateParallelTransport()

BDSFieldObjects * BDSFieldFactory::CreateParallelTransport ( const BDSFieldInfo info)
private

Create special parallel transport 'field' that applies a parallel transport along beam line.

Definition at line 1143 of file BDSFieldFactory.cc.

Referenced by CreateFieldIrregular().

Here is the caller graph for this function:

◆ CreateRMatrix()

BDSFieldObjects * BDSFieldFactory::CreateRMatrix ( const BDSFieldInfo info)
private

Create special rmatrix 'field' that applies an rmatrix.

Definition at line 1123 of file BDSFieldFactory.cc.

References BDSFieldInfo::BeamPipeRadius(), and BDSFieldInfo::MagnetStrength().

Referenced by CreateFieldIrregular().

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

◆ CreateTeleporter()

BDSFieldObjects * BDSFieldFactory::CreateTeleporter ( const BDSFieldInfo info)
private

Create a special teleporter 'field' that shifts particles at the end of rings to match up correctly.

Definition at line 1099 of file BDSFieldFactory.cc.

References designParticle, BDSGlobalConstants::Instance(), BDSFieldInfo::MagnetStrength(), primaryGeneratorAction, BDSPrimaryGeneratorAction::RegisterPTCOneTurnMap(), and BDSFieldInfo::TransformComplete().

Referenced by CreateFieldIrregular().

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

◆ DefaultInterpolatorType()

BDSInterpolatorType BDSFieldFactory::DefaultInterpolatorType ( G4int  numberOfDimensions)
static

Suggest a default interpolator.

Definition at line 457 of file BDSFieldFactory.cc.

Referenced by PrepareFieldDefinitions().

Here is the caller graph for this function:

◆ GetDefinition()

BDSFieldInfo * BDSFieldFactory::GetDefinition ( const G4String &  name) const

Return a BDSFieldInfo instance from the parser definitions. Will exit if no matching field definition found.

Definition at line 407 of file BDSFieldFactory.cc.

References parserDefinitions.

Referenced by CreateFieldERaw(), and CreateFieldMagRaw().

Here is the caller graph for this function:

◆ GetOuterScaling()

G4double BDSFieldFactory::GetOuterScaling ( const BDSMagnetStrength st) const
private

Return the parameter "outerScaling" from strength st, but default to 1.

Definition at line 1153 of file BDSFieldFactory.cc.

References BDSMagnetStrength::KeyHasBeenSet().

Referenced by CreateFieldMagRaw().

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

◆ Instance()

BDSFieldFactory * BDSFieldFactory::Instance ( )
static

Public accessor method for singleton pattern.

Definition at line 151 of file BDSFieldFactory.cc.

References BDSFieldFactory(), and instance.

Referenced by BDS::BuildPlacementGeometry(), BDSComponentFactory::SetFieldDefinitions(), BDSIM::~BDSIM(), and BDSIMLink::~BDSIMLink().

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

◆ PrepareFieldDefinitions()

void BDSFieldFactory::PrepareFieldDefinitions ( const std::vector< GMAD::Field > &  definitions,
G4double  defaultBRho 
)
private

◆ PrepareFieldStrengthFromParameters()

void BDSFieldFactory::PrepareFieldStrengthFromParameters ( BDSMagnetStrength st,
const G4String &  fieldParameters,
G4double &  poleTipRadius 
) const
private

Fill an instance of BDSMagnetStrength with parameters as defined in a string "fieldParameters" that is assumed to be a space-delimited set of parameter=value strings.

Definition at line 374 of file BDSFieldFactory.cc.

References BDSMagnetStrength::AllKeys(), ConvertToDoubleWithException(), BDS::GetUserParametersMap(), BDSMagnetStrength::Unit(), and BDSMagnetStrength::ValidKey().

Referenced by PrepareFieldDefinitions().

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

◆ SetDesignParticle()

static void BDSFieldFactory::SetDesignParticle ( const BDSParticleDefinition designParticleIn)
inlinestatic

Update the internal cache of the rigidity.

Definition at line 76 of file BDSFieldFactory.hh.

References designParticle.

Referenced by BDSIM::Initialise().

Here is the caller graph for this function:

◆ SetPrimaryGeneratorAction()

static void BDSFieldFactory::SetPrimaryGeneratorAction ( BDSPrimaryGeneratorAction pgaIn)
inlinestatic

Update the internal cache of the primary generator action.

Definition at line 79 of file BDSFieldFactory.hh.

References primaryGeneratorAction.

Referenced by BDSIM::Initialise().

Here is the caller graph for this function:

Field Documentation

◆ designParticle

const BDSParticleDefinition * BDSFieldFactory::designParticle = nullptr
staticprivate

Cache of design particle for fields.

Definition at line 176 of file BDSFieldFactory.hh.

Referenced by BDSFieldFactory(), CreateIntegratorMag(), CreateTeleporter(), and SetDesignParticle().

◆ instance

BDSFieldFactory * BDSFieldFactory::instance = nullptr
staticprivate

Instance - singleton pattern.

Definition at line 156 of file BDSFieldFactory.hh.

Referenced by Instance().

◆ parserDefinitions

std::map<G4String, BDSFieldInfo*> BDSFieldFactory::parserDefinitions
private

BDSFieldInfo definitions prepare from parser vector of definitions.

Definition at line 173 of file BDSFieldFactory.hh.

Referenced by GetDefinition(), and PrepareFieldDefinitions().

◆ primaryGeneratorAction

BDSPrimaryGeneratorAction * BDSFieldFactory::primaryGeneratorAction = nullptr
staticprivate

Cache of primary generator action.

Definition at line 179 of file BDSFieldFactory.hh.

Referenced by CreateTeleporter(), and SetPrimaryGeneratorAction().

◆ useOldMultipoleOuterFields

G4bool BDSFieldFactory::useOldMultipoleOuterFields
private

Definition at line 181 of file BDSFieldFactory.hh.


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