BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldLoader.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
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 BDSFIELDLOADER_H
20#define BDSFIELDLOADER_H
21
22#include "BDSArrayReflectionType.hh"
23#include "BDSInterpolatorType.hh"
24#include "G4String.hh"
25#include "G4Transform3D.hh"
26
27#include <array>
28#include <set>
29
34class BDSArrayInfo;
37class BDSFieldInfo;
38class BDSFieldMag;
47
63{
64public:
66 static BDSFieldLoader* Instance();
67
69
70 void DeleteArrays();
71
74 const BDSMagnetStrength* scalingStrength = nullptr,
75 const G4String& scalingKey = "none");
76
79
82
83private:
86
89
91 static void BFilePathOK(const BDSFieldInfo& info);
92 static void EFilePathOK(const BDSFieldInfo& info);
94
96 BDSArray1DCoords* Get1DCached(const G4String& filePath);
97 BDSArray2DCoords* Get2DCached(const G4String& filePath);
98 BDSArray3DCoords* Get3DCached(const G4String& filePath);
99 BDSArray4DCoords* Get4DCached(const G4String& filePath);
101
103 BDSArray2DCoords* LoadPoissonMag2D(const G4String& filePath);
104 BDSArray1DCoords* LoadBDSIM1D(const G4String& filePath);
105 BDSArray2DCoords* LoadBDSIM2D(const G4String& filePath);
106 BDSArray3DCoords* LoadBDSIM3D(const G4String& filePath);
107 BDSArray4DCoords* LoadBDSIM4D(const G4String& filePath);
109
112 void CreateOperators(const BDSArrayReflectionTypeSet* reflectionTypes,
113 const BDSArray4DCoords* existingArray,
114 BDSArrayOperatorIndex*& indexOperator,
115 BDSArrayOperatorValue*& valueOperator) const;
116
121 const std::array<G4bool, 4>& operatesOnXYZT,
122 G4double tolerance=0.05) const;
123
126 G4bool NeedToProvideTransform(const BDSArrayReflectionTypeSet* reflectionTypes) const;
127
128 BDSArray1DCoords* CreateArrayReflected(BDSArray1DCoords* existingArray,
129 const BDSArrayReflectionTypeSet* reflectionType) const;
130 BDSArray2DCoords* CreateArrayReflected(BDSArray2DCoords* existingArray,
131 const BDSArrayReflectionTypeSet* reflectionType) const;
132 BDSArray3DCoords* CreateArrayReflected(BDSArray3DCoords* existingArray,
133 const BDSArrayReflectionTypeSet* reflectionType) const;
134 BDSArray4DCoords* CreateArrayReflected(BDSArray4DCoords* existingArray,
135 const BDSArrayReflectionTypeSet* reflectionType) const;
136
139 BDSInterpolatorType interpolatorType) const;
140
143 BDSInterpolatorType interpolatorType) const;
144
147 BDSInterpolatorType interpolatorType) const;
148
151 BDSInterpolatorType interpolatorType) const;
152
154 BDSFieldMagInterpolated* LoadBDSIM1DB(const G4String& filePath,
155 BDSInterpolatorType interpolatorType,
156 const G4Transform3D& transform,
157 G4double bScaling,
158 const BDSArrayReflectionTypeSet* reflection = nullptr);
159
161 BDSFieldMagInterpolated* LoadBDSIM2DB(const G4String& filePath,
162 BDSInterpolatorType interpolatorType,
163 const G4Transform3D& transform,
164 G4double bScaling,
165 const BDSArrayReflectionTypeSet* reflection = nullptr);
166
168 BDSFieldMagInterpolated* LoadBDSIM3DB(const G4String& filePath,
169 BDSInterpolatorType interpolatorType,
170 const G4Transform3D& transform,
171 G4double bScaling,
172 const BDSArrayReflectionTypeSet* reflection = nullptr);
173
175 BDSFieldMagInterpolated* LoadBDSIM4DB(const G4String& filePath,
176 BDSInterpolatorType interpolatorType,
177 const G4Transform3D& transform,
178 G4double bScaling,
179 const BDSArrayReflectionTypeSet* reflection = nullptr);
180
182 BDSFieldMagInterpolated* LoadPoissonSuperFishB(const G4String& filePath,
183 BDSInterpolatorType interpolatorType,
184 const G4Transform3D& transform,
185 G4double bScaling,
186 const BDSArrayReflectionTypeSet* reflection = nullptr);
187
190 BDSFieldMagInterpolated* LoadPoissonSuperFishBQuad(const G4String& filePath,
191 BDSInterpolatorType interpolatorType,
192 const G4Transform3D& transform,
193 G4double bScaling,
194 const BDSArrayReflectionTypeSet* reflection = nullptr);
195
199 BDSInterpolatorType interpolatorType,
200 const G4Transform3D& transform,
201 G4double bScaling,
202 const BDSArrayReflectionTypeSet* reflection = nullptr);
203
205 BDSFieldEInterpolated* LoadBDSIM1DE(const G4String& filePath,
206 BDSInterpolatorType interpolatorType,
207 const G4Transform3D& transform,
208 G4double eScaling,
209 const BDSArrayReflectionTypeSet* reflection = nullptr);
210
212 BDSFieldEInterpolated* LoadBDSIM2DE(const G4String& filePath,
213 BDSInterpolatorType interpolatorType,
214 const G4Transform3D& transform,
215 G4double eScaling,
216 const BDSArrayReflectionTypeSet* reflection = nullptr);
217
219 BDSFieldEInterpolated* LoadBDSIM3DE(const G4String& filePath,
220 BDSInterpolatorType interpolatorType,
221 const G4Transform3D& transform,
222 G4double eScaling,
223 const BDSArrayReflectionTypeSet* reflection = nullptr);
224
226 BDSFieldEInterpolated* LoadBDSIM4DE(const G4String& filePath,
227 BDSInterpolatorType interpolatorType,
228 const G4Transform3D& transform,
229 G4double eScaling,
230 const BDSArrayReflectionTypeSet* reflection = nullptr);
231
233 BDSFieldEMInterpolated* LoadBDSIM1DEM(const G4String& eFilePath,
234 const G4String& bFilePath,
235 BDSInterpolatorType eInterpolatorType,
236 BDSInterpolatorType bInterpolatorType,
237 const G4Transform3D& transform,
238 G4double eScaling,
239 G4double bScaling,
240 const BDSArrayReflectionTypeSet* eReflection = nullptr,
241 const BDSArrayReflectionTypeSet* bReflection = nullptr);
242
244 BDSFieldEMInterpolated* LoadBDSIM2DEM(const G4String& eFilePath,
245 const G4String& bFilePath,
246 BDSInterpolatorType eInterpolatorType,
247 BDSInterpolatorType bInterpolatorType,
248 const G4Transform3D& transform,
249 G4double eScaling,
250 G4double bScaling,
251 const BDSArrayReflectionTypeSet* eReflection = nullptr,
252 const BDSArrayReflectionTypeSet* bReflection = nullptr);
253
255 BDSFieldEMInterpolated* LoadBDSIM3DEM(const G4String& eFilePath,
256 const G4String& bFilePath,
257 BDSInterpolatorType eInterpolatorType,
258 BDSInterpolatorType bInterpolatorType,
259 const G4Transform3D& transform,
260 G4double eScaling,
261 G4double bScaling,
262 const BDSArrayReflectionTypeSet* eReflection = nullptr,
263 const BDSArrayReflectionTypeSet* bReflection = nullptr);
264
266 BDSFieldEMInterpolated* LoadBDSIM4DEM(const G4String& eFilePath,
267 const G4String& bFilePath,
268 BDSInterpolatorType eInterpolatorType,
269 BDSInterpolatorType bInterpolatorType,
270 const G4Transform3D& transform,
271 G4double eScaling,
272 G4double bScaling,
273 const BDSArrayReflectionTypeSet* eReflection = nullptr,
274 const BDSArrayReflectionTypeSet* bReflection = nullptr);
275
277 std::map<G4String, BDSArray1DCoords*> arrays1d;
278 std::map<G4String, BDSArray2DCoords*> arrays2d;
279 std::map<G4String, BDSArray3DCoords*> arrays3d;
280 std::map<G4String, BDSArray4DCoords*> arrays4d;
282};
283
284#endif
1D array with spatial mapping derived from BDSArray4DCoords.
2D array with spatial mapping derived from BDSArray4DCoords.
3D array with spatial mapping derived from BDSArray4DCoords.
Overlay of 4D array that provides uniform only spatial coordinate mapping.
Simple holder of information about an array.
Definition: BDSArrayInfo.hh:35
Interface for modifying by reference array indices.
Interface for modifying field values.
Class to provide scaling and a base class pointer for interpolator fields.
Class to provide scaling and a base class pointer for interpolator fields.
All info required to build complete field of any type.
Definition: BDSFieldInfo.hh:65
A loader for various field map formats.
BDSFieldEInterpolated * LoadBDSIM3DE(const G4String &filePath, BDSInterpolatorType interpolatorType, const G4Transform3D &transform, G4double eScaling, const BDSArrayReflectionTypeSet *reflection=nullptr)
Load a 3D BDSIM format electric field.
static void EFilePathOK(const BDSFieldInfo &info)
Check file path isn't empty and throw exception if it is.
BDSFieldEMInterpolated * LoadBDSIM3DEM(const G4String &eFilePath, const G4String &bFilePath, BDSInterpolatorType eInterpolatorType, BDSInterpolatorType bInterpolatorType, const G4Transform3D &transform, G4double eScaling, G4double bScaling, const BDSArrayReflectionTypeSet *eReflection=nullptr, const BDSArrayReflectionTypeSet *bReflection=nullptr)
Load a 3D BDSIM format electro-magnetic field.
BDSInterpolator3D * CreateInterpolator3D(BDSArray3DCoords *array, BDSInterpolatorType interpolatorType) const
Create the appropriate 3D interpolator for an array.
std::map< G4String, BDSArray1DCoords * > arrays1d
Map of cached field map array.
BDSFieldEMInterpolated * LoadBDSIM4DEM(const G4String &eFilePath, const G4String &bFilePath, BDSInterpolatorType eInterpolatorType, BDSInterpolatorType bInterpolatorType, const G4Transform3D &transform, G4double eScaling, G4double bScaling, const BDSArrayReflectionTypeSet *eReflection=nullptr, const BDSArrayReflectionTypeSet *bReflection=nullptr)
Load a 4D BDSIM format electro-magnetic field.
std::map< G4String, BDSArray2DCoords * > arrays2d
Map of cached field map array.
BDSArray2DCoords * LoadPoissonMag2D(const G4String &filePath)
Utility function to use the right templated loader class (gz or normal).
BDSArray3DCoords * Get3DCached(const G4String &filePath)
Return the cached array if there is one - may return nullptr.
BDSArray2DCoords * LoadBDSIM2D(const G4String &filePath)
Utility function to use the right templated loader class (gz or normal).
BDSFieldEMInterpolated * LoadEMField(const BDSFieldInfo &info)
Main interface to load an electro-magnetic field.
BDSInterpolator1D * CreateInterpolator1D(BDSArray1DCoords *array, BDSInterpolatorType interpolatorType) const
Create the appropriate 1D interpolator for an array.
BDSFieldEInterpolated * LoadBDSIM2DE(const G4String &filePath, BDSInterpolatorType interpolatorType, const G4Transform3D &transform, G4double eScaling, const BDSArrayReflectionTypeSet *reflection=nullptr)
Load a 2D BDSIM format electric field.
std::map< G4String, BDSArray3DCoords * > arrays3d
Map of cached field map array.
void ReportIfProblemWithReflection(const BDSArrayInfo &info, const std::array< G4bool, 4 > &operatesOnXYZT, G4double tolerance=0.05) const
BDSInterpolator2D * CreateInterpolator2D(BDSArray2DCoords *array, BDSInterpolatorType interpolatorType) const
Create the appropriate 2D interpolator for an array.
BDSFieldEInterpolated * LoadBDSIM4DE(const G4String &filePath, BDSInterpolatorType interpolatorType, const G4Transform3D &transform, G4double eScaling, const BDSArrayReflectionTypeSet *reflection=nullptr)
Load a 4D BDSIM format electric field.
BDSFieldMagInterpolated * LoadPoissonSuperFishB(const G4String &filePath, BDSInterpolatorType interpolatorType, const G4Transform3D &transform, G4double bScaling, const BDSArrayReflectionTypeSet *reflection=nullptr)
Load a 2D poisson superfish B field map.
BDSFieldLoader()
Private default constructor as singleton.
BDSFieldMagInterpolated * LoadPoissonSuperFishBQuad(const G4String &filePath, BDSInterpolatorType interpolatorType, const G4Transform3D &transform, G4double bScaling, const BDSArrayReflectionTypeSet *reflection=nullptr)
BDSArray1DCoords * Get1DCached(const G4String &filePath)
Return the cached array if there is one - may return nullptr.
BDSFieldMagInterpolated * LoadBDSIM1DB(const G4String &filePath, BDSInterpolatorType interpolatorType, const G4Transform3D &transform, G4double bScaling, const BDSArrayReflectionTypeSet *reflection=nullptr)
Load a 1D BDSIM format magnetic field.
BDSFieldEMInterpolated * LoadBDSIM1DEM(const G4String &eFilePath, const G4String &bFilePath, BDSInterpolatorType eInterpolatorType, BDSInterpolatorType bInterpolatorType, const G4Transform3D &transform, G4double eScaling, G4double bScaling, const BDSArrayReflectionTypeSet *eReflection=nullptr, const BDSArrayReflectionTypeSet *bReflection=nullptr)
Load a 1D BDSIM format electro-magnetic field.
BDSArray3DCoords * LoadBDSIM3D(const G4String &filePath)
Utility function to use the right templated loader class (gz or normal).
BDSFieldMagInterpolated * LoadBDSIM4DB(const G4String &filePath, BDSInterpolatorType interpolatorType, const G4Transform3D &transform, G4double bScaling, const BDSArrayReflectionTypeSet *reflection=nullptr)
Load a 4D BDSIM format magnetic field.
BDSFieldMagInterpolated * LoadBDSIM2DB(const G4String &filePath, BDSInterpolatorType interpolatorType, const G4Transform3D &transform, G4double bScaling, const BDSArrayReflectionTypeSet *reflection=nullptr)
Load a 2D BDSIM format magnetic field.
BDSFieldMagInterpolated * LoadMagField(const BDSFieldInfo &info, const BDSMagnetStrength *scalingStrength=nullptr, const G4String &scalingKey="none")
Main interface to load a magnetic field.
static BDSFieldLoader * instance
Singleton instance.
static BDSFieldLoader * Instance()
Singleton accessor.
BDSFieldEMInterpolated * LoadBDSIM2DEM(const G4String &eFilePath, const G4String &bFilePath, BDSInterpolatorType eInterpolatorType, BDSInterpolatorType bInterpolatorType, const G4Transform3D &transform, G4double eScaling, G4double bScaling, const BDSArrayReflectionTypeSet *eReflection=nullptr, const BDSArrayReflectionTypeSet *bReflection=nullptr)
Load a 2D BDSIM format electro-magnetic field.
BDSArray4DCoords * Get4DCached(const G4String &filePath)
Return the cached array if there is one - may return nullptr.
BDSFieldMagInterpolated * LoadPoissonSuperFishBDipole(const G4String &filePath, BDSInterpolatorType interpolatorType, const G4Transform3D &transform, G4double bScaling, const BDSArrayReflectionTypeSet *reflection=nullptr)
BDSInterpolator4D * CreateInterpolator4D(BDSArray4DCoords *array, BDSInterpolatorType interpolatorType) const
Create the appropriate 4D interpolator for an array.
BDSFieldMagInterpolated * LoadBDSIM3DB(const G4String &filePath, BDSInterpolatorType interpolatorType, const G4Transform3D &transform, G4double bScaling, const BDSArrayReflectionTypeSet *reflection=nullptr)
Load a 3D BDSIM format magnetic field.
BDSFieldEInterpolated * LoadEField(const BDSFieldInfo &info)
Main interface to load an electric field.
G4bool NeedToProvideTransform(const BDSArrayReflectionTypeSet *reflectionTypes) const
BDSFieldEInterpolated * LoadBDSIM1DE(const G4String &filePath, BDSInterpolatorType interpolatorType, const G4Transform3D &transform, G4double eScaling, const BDSArrayReflectionTypeSet *reflection=nullptr)
Load a 1D BDSIM format electric field.
void CreateOperators(const BDSArrayReflectionTypeSet *reflectionTypes, const BDSArray4DCoords *existingArray, BDSArrayOperatorIndex *&indexOperator, BDSArrayOperatorValue *&valueOperator) const
static void BFilePathOK(const BDSFieldInfo &info)
Check file path isn't empty and throw exception if it is.
BDSArray1DCoords * LoadBDSIM1D(const G4String &filePath)
Utility function to use the right templated loader class (gz or normal).
BDSArray4DCoords * LoadBDSIM4D(const G4String &filePath)
Utility function to use the right templated loader class (gz or normal).
std::map< G4String, BDSArray4DCoords * > arrays4d
Map of cached field map array.
BDSArray2DCoords * Get2DCached(const G4String &filePath)
Return the cached array if there is one - may return nullptr.
Class to provide scaling and a base class pointer for interpolator fields.
Interface for static magnetic fields that may or may not be local.
Definition: BDSFieldMag.hh:37
Interface for all 1D interpolators.
Interface for all 2D interpolators.
Interface for all 3D interpolators.
Interface for all 4D interpolators.
Efficient storage of magnet strengths.