BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldInfo.cc
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#include "BDSArrayReflectionType.hh"
20#include "BDSFieldInfo.hh"
21#include "BDSFieldType.hh"
22#include "BDSIntegratorType.hh"
23#include "BDSInterpolatorType.hh"
24#include "BDSMagnetStrength.hh"
25#include "BDSUtilities.hh"
26
27#include "globals.hh" // geant4 types / globals
28#include "G4RotationMatrix.hh"
29#include "G4ThreeVector.hh"
30#include "G4Track.hh"
31#include "G4Transform3D.hh"
32#include "G4UserLimits.hh"
33
34#include <algorithm>
35#include <ostream>
36
37G4UserLimits* BDSFieldInfo::defaultUL = nullptr;
38
40 fieldType(BDSFieldType::none),
41 brho(0),
42 integratorType(BDSIntegratorType::none),
43 magnetStrength(nullptr),
44 provideGlobalTransform(false),
45 transform(nullptr),
46 magneticFieldFilePath(""),
47 magneticFieldFormat(BDSFieldFormat::none),
48 magneticInterpolatorType(BDSInterpolatorType::nearest3d),
49 magneticArrayReflectionTypeSet(BDSArrayReflectionTypeSet()),
50 electricFieldFilePath(""),
51 electricFieldFormat(BDSFieldFormat::none),
52 electricInterpolatorType(BDSInterpolatorType::nearest3d),
53 electricArrayReflectionTypeSet(BDSArrayReflectionTypeSet()),
54 cacheTransforms(true),
55 eScaling(1.0),
56 bScaling(1.0),
57 timeOffset(0),
58 autoScale(false),
59 stepLimit(nullptr),
60 poleTipRadius(1),
61 beamPipeRadius(0),
62 chordStepMinimum(-1),
63 tilt(0),
64 secondFieldOnLeft(false),
65 magneticSubFieldName(""),
66 electricSubFieldName(""),
67 usePlacementWorldTransform(false),
68 transformBeamline(nullptr),
69 nameOfParserDefinition("")
70{;}
71
73 G4double brhoIn,
74 BDSIntegratorType integratorTypeIn,
75 BDSMagnetStrength* magnetStrengthIn,
76 G4bool provideGlobalTransformIn,
77 const G4Transform3D& transformIn,
78 const G4String& magneticFieldFilePathIn,
79 BDSFieldFormat magneticFieldFormatIn,
80 BDSInterpolatorType magneticInterpolatorTypeIn,
81 const G4String& electricFieldFilePathIn,
82 BDSFieldFormat electricFieldFormatIn,
83 BDSInterpolatorType electricInterpolatorTypeIn,
84 G4bool cacheTransformsIn,
85 G4double eScalingIn,
86 G4double bScalingIn,
87 G4double timeOffsetIn,
88 G4bool autoScaleIn,
89 G4UserLimits* stepLimitIn,
90 G4double poleTipRadiusIn,
91 G4double beamPipeRadiusIn,
92 G4bool leftIn,
93 const G4String& magneticSubFieldNameIn,
94 const G4String& electricSubFieldNameIn):
95 fieldType(fieldTypeIn),
96 brho(brhoIn),
97 integratorType(integratorTypeIn),
98 magnetStrength(magnetStrengthIn),
99 provideGlobalTransform(provideGlobalTransformIn),
100 transform(nullptr),
101 magneticFieldFilePath(magneticFieldFilePathIn),
102 magneticFieldFormat(magneticFieldFormatIn),
103 magneticInterpolatorType(magneticInterpolatorTypeIn),
104 magneticArrayReflectionTypeSet(BDSArrayReflectionTypeSet()),
105 electricFieldFilePath(electricFieldFilePathIn),
106 electricFieldFormat(electricFieldFormatIn),
107 electricInterpolatorType(electricInterpolatorTypeIn),
108 electricArrayReflectionTypeSet(BDSArrayReflectionTypeSet()),
109 cacheTransforms(cacheTransformsIn),
110 eScaling(eScalingIn),
111 bScaling(bScalingIn),
112 timeOffset(timeOffsetIn),
113 autoScale(autoScaleIn),
114 stepLimit(stepLimitIn),
115 poleTipRadius(poleTipRadiusIn),
116 beamPipeRadius(beamPipeRadiusIn),
117 chordStepMinimum(-1),
118 secondFieldOnLeft(leftIn),
119 magneticSubFieldName(magneticSubFieldNameIn),
120 electricSubFieldName(electricSubFieldNameIn),
121 usePlacementWorldTransform(false),
122 transformBeamline(nullptr),
123 nameOfParserDefinition("")
124{
125 if (transformIn != G4Transform3D::Identity)
126 {transform = new G4Transform3D(transformIn);}
127
128 // back calculate tilt angle from field transform
129 G4ThreeVector unitY(0,1,0);
130 G4ThreeVector unitYR = unitY.transform(transformIn.getRotation());
131 tilt = -(CLHEP::halfpi - unitYR.getPhi()); // halfpi offset for which axis we choose by convention
132}
133
134BDSFieldInfo::~BDSFieldInfo()
135{
136 delete magnetStrength;
137 delete transform;
138 delete stepLimit;
139 delete transformBeamline;
140}
141
143 fieldType(other.fieldType),
144 brho(other.brho),
145 integratorType(other.integratorType),
146 provideGlobalTransform(other.provideGlobalTransform),
147 transform(nullptr),
148 magneticFieldFilePath(other.magneticFieldFilePath),
149 magneticFieldFormat(other.magneticFieldFormat),
150 magneticInterpolatorType(other.magneticInterpolatorType),
151 magneticArrayReflectionTypeSet(other.magneticArrayReflectionTypeSet),
152 electricFieldFilePath(other.electricFieldFilePath),
153 electricFieldFormat(other.electricFieldFormat),
154 electricInterpolatorType(other.electricInterpolatorType),
155 electricArrayReflectionTypeSet(other.electricArrayReflectionTypeSet),
156 cacheTransforms(other.cacheTransforms),
157 eScaling(other.eScaling),
158 bScaling(other.bScaling),
159 timeOffset(other.timeOffset),
160 autoScale(other.autoScale),
161 poleTipRadius(other.poleTipRadius),
162 beamPipeRadius(other.beamPipeRadius),
163 chordStepMinimum(other.chordStepMinimum),
164 tilt(other.tilt),
165 secondFieldOnLeft(other.secondFieldOnLeft),
166 magneticSubFieldName(other.magneticSubFieldName),
167 electricSubFieldName(other.electricSubFieldName),
168 usePlacementWorldTransform(other.usePlacementWorldTransform),
169 transformBeamline(nullptr),
170 nameOfParserDefinition(other.nameOfParserDefinition)
171{
172 if (other.transform)
173 {transform = new G4Transform3D(*other.transform);}
174
175 if (other.magnetStrength)
176 {magnetStrength = new BDSMagnetStrength(*other.magnetStrength);}
177 else
178 {magnetStrength = nullptr;} // also nullptr
179
180 if (other.stepLimit)
181 {stepLimit = new G4UserLimits(*other.stepLimit);}
182 else
183 {stepLimit = nullptr;}
184
185 if (other.transformBeamline)
186 {transformBeamline = new G4Transform3D(*other.transformBeamline);}
187}
188
189void BDSFieldInfo::SetUserLimits(G4UserLimits* userLimitsIn)
190{
191 if (stepLimit != defaultUL)
192 {delete stepLimit;} // shouldn't delete global default step limit!
193 stepLimit = userLimitsIn;
194}
195
197 G4bool warn) const
198{
199 if (stepLimit && (stepLimit != defaultUL))
200 {
201 G4UserLimits* old = stepLimit;
202 stepLimit = BDS::CreateUserLimits(stepLimit, maximumStepSize, 1.0);
203 if (stepLimit == old)
204 {warn = false;} // no change and warning would print out wrong number
205 if ((stepLimit != old) && (old != defaultUL))
206 {delete old;}
207 }
208 else
209 {stepLimit = new G4UserLimits(maximumStepSize);}
210
211 if (warn)
212 {
213 G4cout << "Reducing maximum step size of field definition \"" << nameOfParserDefinition
214 << "\" to " << maximumStepSize << " mm " << G4endl;
215 }
216}
217
218std::ostream& operator<< (std::ostream& out, BDSFieldInfo const& info)
219{
220 out << "Parser definition name: \"" << info.nameOfParserDefinition << "\"" << G4endl;
221 out << "Field type: " << info.fieldType << G4endl;
222 out << "Rigidity: " << info.brho << G4endl;
223 out << "Integrator: " << info.integratorType << G4endl;
224 out << "Global transform? " << info.provideGlobalTransform << G4endl;
225 out << "B map file: " << info.magneticFieldFilePath << G4endl;
226 out << "B map file format: " << info.magneticFieldFormat << G4endl;
227 out << "B interpolator " << info.magneticInterpolatorType << G4endl;
228 out << "B array reflection: " << info.magneticArrayReflectionTypeSet << G4endl;
229 out << "E map file: " << info.electricFieldFilePath << G4endl;
230 out << "E map file format: " << info.electricFieldFormat << G4endl;
231 out << "E interpolator " << info.electricInterpolatorType << G4endl;
232 out << "E array reflection: " << info.electricArrayReflectionTypeSet << G4endl;
233 out << "Transform caching: " << info.cacheTransforms << G4endl;
234 out << "E Scaling: " << info.eScaling << G4endl;
235 out << "B Scaling: " << info.bScaling << G4endl;
236 out << "t offset " << info.timeOffset << G4endl;
237 out << "Auto scale " << info.autoScale << G4endl;
238 out << "Pole tip radius: " << info.poleTipRadius << G4endl;
239 out << "Beam pipe radius: " << info.beamPipeRadius << G4endl;
240 out << "Chord Step Min: " << info.chordStepMinimum << G4endl;
241 out << "Tilt: " << info.tilt << G4endl;
242 out << "Second field on left " << info.secondFieldOnLeft << G4endl;
243 out << "Magnetic Sub Field " << info.magneticSubFieldName << G4endl;
244 out << "Electric Sub Field " << info.electricSubFieldName << G4endl;
245 out << "Use Placement World Transform " << info.usePlacementWorldTransform << G4endl;
246 if (info.magnetStrength)
247 {out << "Magnet strength: " << *(info.magnetStrength) << G4endl;}
248 if (info.stepLimit)
249 {
250 G4Track t = G4Track(); // dummy track
251 G4double maxStep = info.stepLimit->GetMaxAllowedStep(t);
252 out << "Step limit: " << maxStep << G4endl;
253 }
254 return out;
255}
256
257void BDSFieldInfo::Translate(const G4ThreeVector& translationIn)
258{
259 if (!transform)
260 {transform = new G4Transform3D();}
261 G4RotationMatrix rm = transform->getRotation();
262 G4ThreeVector translation = transform->getTranslation();
263 translation += translationIn;
264 G4Transform3D* newTransform = new G4Transform3D(rm, translation);
265 delete transform;
266 transform = newTransform;
267}
268
269G4Transform3D BDSFieldInfo::Transform() const
270{
271 return transform ? *transform : G4Transform3D::Identity;
272}
273
275{
276 return transformBeamline ? *transformBeamline : G4Transform3D::Identity;
277}
278
280{
281 return Transform() * TransformBeamline();
282}
283
284void BDSFieldInfo::SetTransform(const G4Transform3D& transformIn)
285{
286 delete transform;
287 transform = new G4Transform3D(transformIn);
288}
289
290void BDSFieldInfo::SetTransformBeamline(const G4Transform3D& transformIn)
291{
292 delete transformBeamline;
293 transformBeamline = new G4Transform3D(transformIn);
294}
All info required to build complete field of any type.
Definition: BDSFieldInfo.hh:65
G4double beamPipeRadius
Optional radius of beam pipe.
void UpdateUserLimitsLengthMaximumStepSize(G4double maximumStepSize, G4bool warn=false) const
static G4UserLimits * defaultUL
Cache of default user limits.
G4Transform3D * transformBeamline
Transform from curvilinear frame to this field - ie beam line bit only.
G4Transform3D TransformComplete() const
Compound transform of field + beam line transform.
void SetTransformBeamline(const G4Transform3D &transformIn)
Set the beam line transform.
G4double poleTipRadius
Radius at which point the field will be scaled to.
G4double tilt
Cache of tilt of field.
void SetTransform(const G4Transform3D &transformIn)
Set the field definition transform.
G4Transform3D Transform() const
Transform for the field definition only.
G4bool secondFieldOnLeft
Flag for case of two-beam field - if not left, it's right.
void Translate(const G4ThreeVector &translationIn)
G4Transform3D TransformBeamline() const
Transform from the curvilinear coordinates to the beam line component.
void SetUserLimits(G4UserLimits *userLimitsIn)
Delete and replace the user limits which this class owns (only if not default ul).
G4Transform3D * transform
Transform w.r.t. solid field will be attached to.
BDSFieldInfo()
Default constructor for zero field effectively.
Definition: BDSFieldInfo.cc:39
Efficient storage of magnet strengths.
G4UserLimits * CreateUserLimits(G4UserLimits *defaultUL, G4double length, G4double fraction=1.6)