BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSFieldFactory.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef BDSFIELDFACTORY_H
20#define BDSFIELDFACTORY_H
21
22#include "BDSFieldType.hh"
23#include "BDSInterpolatorType.hh"
24#include "BDSMagnetType.hh"
25
26#include "globals.hh"
27
28#include <map>
29#include <vector>
30
31namespace GMAD
32{
33 class Field;
34 class Modulator;
35}
36
37class BDSFieldE;
38class BDSFieldInfo;
39class BDSFieldMag;
40class BDSFieldObjects;
42class BDSModulator;
46
47class G4EquationOfMotion;
48class G4MagIntegratorStepper;
49class G4Mag_EqRhs;
50
77{
78public:
80 static BDSFieldFactory* Instance();
81
83 static void SetDesignParticle(const BDSParticleDefinition* designParticleIn)
84 {designParticle = designParticleIn;}
87
89
92 const BDSMagnetStrength* scalingStrength = nullptr,
93 const G4String& scalingKey = "none");
94
98 BDSFieldInfo* GetDefinition(const G4String& name) const;
99
103 BDSModulatorInfo* GetModulatorDefinition(const G4String& modulatorName) const;
104
106 static BDSInterpolatorType DefaultInterpolatorType(G4int numberOfDimensions);
107
108 static G4double CalculateGlobalPhase(G4double oscillatorFrequency,
109 G4double tOffsetIn);
110
111 static G4double CalculateGlobalPhase(const BDSModulatorInfo& modulatorInfo,
112 const BDSFieldInfo& fieldInfo);
113
114private:
117 const BDSMagnetStrength* scalingStrength = nullptr,
118 const G4String& scalingKey = "none");
119
122
125
128
131 const BDSMagnetStrength* scalingStrength = nullptr,
132 const G4String& scalingKey = "none");
133
136
139 G4MagIntegratorStepper* CreateIntegratorMag(const BDSFieldInfo& info,
140 G4Mag_EqRhs* eqOfM,
141 const BDSMagnetStrength* strength);
142
145 G4MagIntegratorStepper* CreateIntegratorEM(const BDSFieldInfo& info,
146 G4EquationOfMotion* eqOfM);
147
151 G4MagIntegratorStepper* CreateIntegratorE(const BDSFieldInfo& info,
152 G4EquationOfMotion* eqOfM);
153
157
160
163
167
169 G4double GetOuterScaling(const BDSMagnetStrength* st) const;
170
172 BDSModulator* CreateModulator(const BDSModulatorInfo* modulatorRecipe,
173 const BDSFieldInfo& info) const;
174
177
180
182 void PrepareFieldDefinitions(const std::vector<GMAD::Field>& definitions,
183 G4double defaultBRho);
184
186 void PrepareModulatorDefinitions(const std::vector<GMAD::Modulator>& definitions);
187
189 G4double ConvertToDoubleWithException(const G4String& value,
190 const G4String& parameterNameForError) const;
191
195 const G4String& fieldParameters,
196 G4double& poleTipRadius) const;
197
199 std::map<G4String, BDSFieldInfo*> parserDefinitions;
200 std::map<G4String, BDSModulatorInfo*> parserModulatorDefinitions;
201
204
207
208 G4bool useOldMultipoleOuterFields;
209};
210#endif
Interface for BDSIM electric fields that may or may not be local.
Definition: BDSFieldE.hh:39
Factory that produces fields and their associated objects.
G4double GetOuterScaling(const BDSMagnetStrength *st) const
Return the parameter "outerScaling" from strength st, but default to 1.
static void SetDesignParticle(const BDSParticleDefinition *designParticleIn)
Update the internal cache of the rigidity.
void PrepareFieldDefinitions(const std::vector< GMAD::Field > &definitions, G4double defaultBRho)
Prepare all required definitions that can be used dynamically.
void PrepareModulatorDefinitions(const std::vector< GMAD::Modulator > &definitions)
Prepare all required modulator definitions that can be used dynamically.
static void SetPrimaryGeneratorAction(BDSPrimaryGeneratorAction *pgaIn)
Update the internal cache of the primary generator action.
BDSFieldObjects * CreateFieldE(const BDSFieldInfo &info)
Create an electric field.
BDSFieldObjects * CreateFieldIrregular(const BDSFieldInfo &info)
Create an irregular (special) field.
BDSFieldObjects * CreateRMatrix(const BDSFieldInfo &info)
Create special rmatrix 'field' that applies an rmatrix.
void PrepareFieldStrengthFromParameters(BDSMagnetStrength *st, const G4String &fieldParameters, G4double &poleTipRadius) const
static BDSPrimaryGeneratorAction * primaryGeneratorAction
Cache of primary generator action.
BDSModulator * CreateModulator(const BDSModulatorInfo *modulatorRecipe, const BDSFieldInfo &info) const
Create the necessary modulator.
BDSFieldE * CreateFieldERaw(const BDSFieldInfo &info)
Creat just the electric field object.
G4double ConvertToDoubleWithException(const G4String &value, const G4String &parameterNameForError) const
Convert the string 'value' to a double. Throw an exception including the parameterNameForError if it ...
static BDSFieldFactory * Instance()
Public accessor method for singleton pattern.
BDSFieldObjects * CreateField(const BDSFieldInfo &info, const BDSMagnetStrength *scalingStrength=nullptr, const G4String &scalingKey="none")
Main interface to field factory.
BDSFieldObjects * CreateParallelTransport(const BDSFieldInfo &info)
static BDSFieldFactory * instance
Instance - singleton pattern.
BDSFieldInfo * GetDefinition(const G4String &name) const
BDSFieldObjects * CreateCavityFringe(const BDSFieldInfo &info)
Create special rf cavity fringe 'field' that applies an rmatrix.
std::map< G4String, BDSFieldInfo * > parserDefinitions
BDSFieldInfo definitions prepare from parser vector of definitions.
G4MagIntegratorStepper * CreateIntegratorMag(const BDSFieldInfo &info, G4Mag_EqRhs *eqOfM, const BDSMagnetStrength *strength)
G4MagIntegratorStepper * CreateIntegratorE(const BDSFieldInfo &info, G4EquationOfMotion *eqOfM)
BDSFieldObjects * CreateFieldMag(const BDSFieldInfo &info, const BDSMagnetStrength *scalingStrength=nullptr, const G4String &scalingKey="none")
Create a purely magnetic field.
static BDSInterpolatorType DefaultInterpolatorType(G4int numberOfDimensions)
Suggest a default interpolator.
BDSFieldObjects * CreateFieldEM(const BDSFieldInfo &info)
Create a general EM field.
BDSFieldMag * CreateFieldMagRaw(const BDSFieldInfo &info, const BDSMagnetStrength *scalingStrength=nullptr, const G4String &scalingKey="none")
Creat just the magnetic field object.
BDSFieldFactory()
Private default constructor as singleton class.
BDSModulatorInfo * GetModulatorDefinition(const G4String &modulatorName) const
BDSFieldObjects * CreateTeleporter(const BDSFieldInfo &info)
G4MagIntegratorStepper * CreateIntegratorEM(const BDSFieldInfo &info, G4EquationOfMotion *eqOfM)
static const BDSParticleDefinition * designParticle
Cache of design particle for fields.
All info required to build complete field of any type.
Definition: BDSFieldInfo.hh:66
Interface for static magnetic fields that may or may not be local.
Definition: BDSFieldMag.hh:39
A holder for all the Geant4 field related objects.
Efficient storage of magnet strengths.
Holder class for all information required to describe a modulator.
Base class for a modulator.
Definition: BDSModulator.hh:37
Wrapper for particle definition.
Generates primary particle vertices using BDSBunch.
Parser namespace for GMAD language. Combination of Geant4 and MAD.