BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSComponentFactory.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 BDSCOMPONENTFACTORY_H
20#define BDSCOMPONENTFACTORY_H
21
22#include "BDSFieldType.hh"
23#include "BDSMagnetGeometryType.hh"
24#include "BDSMagnetStrength.hh"
25#include "BDSMagnetType.hh"
26#include "BDSIntegratorType.hh"
27#include "BDSIntegratorSetType.hh"
28
29#include "globals.hh"
30#include "G4ThreeVector.hh"
31#include "G4Transform3D.hh"
32
33#include "CLHEP/Units/PhysicalConstants.h"
34
35#include <map>
36
37class G4Colour;
38class G4Material;
39
40namespace GMAD
41{
42 struct Element;
43}
45class BDSBeamPipeInfo;
46class BDSCavityInfo;
48class BDSCrystalInfo;
49class BDSFieldInfo;
51class BDSMagnet;
54class BDSTiltOffset;
55
76{
77public:
78 explicit BDSComponentFactory(const BDSParticleDefinition* designParticleIn,
79 BDSComponentFactoryUser* userComponentFactoryIn = nullptr,
80 G4bool usualPrintOut = true);
82
89 GMAD::Element const* prevElementIn,
90 GMAD::Element const* nextElementIn,
91 G4double currentArcLength = 0);
92
96
99 BDSAcceleratorComponent* CreateTeleporter(const G4double teleporterLength,
100 const G4double teleporterWidth,
101 const G4Transform3D& transformIn);
102
105
107 static G4Transform3D CreateFieldTransform(const BDSTiltOffset* tiltOffset);
108
111 static G4Transform3D CreateFieldTransform(GMAD::Element const* el);
112
116 const G4ThreeVector& inputFaceNormal = G4ThreeVector(0,0,-1),
117 const G4ThreeVector& outputFaceNormal = G4ThreeVector(0,0,1));
118
122 G4double angleIn,
123 G4double angleOut);
124
127 static G4bool YokeOnLeft(const GMAD::Element* el,
128 const BDSMagnetStrength* st);
129
131 static G4double PrepareHorizontalWidth(GMAD::Element const* el,
132 G4double defaultHorizontalWidth = -1);
133
134 static G4double ScalingFieldOuter(const GMAD::Element* ele);
135
138 const BDSFieldType& fieldType,
139 const BDSBeamPipeInfo* bpInfo,
140 const BDSMagnetOuterInfo* outerInfo,
141 const G4Transform3D& fieldTransform,
142 const BDSIntegratorSet* integratorSetIn,
143 G4double brhoIn,
144 G4double outerFieldScaling = 1.0);
145
149 static BDSMagnetOuterInfo* PrepareMagnetOuterInfo(const G4String& elementNameIn,
150 const GMAD::Element* el,
151 const BDSMagnetStrength* st,
152 const BDSBeamPipeInfo* beamPipe,
153 G4double defaultHorizontalWidth = -1,
154 G4double defaultVHRatio = 1.0,
155 G4double defaultCoilWidthFraction = -1,
156 G4double defaultCoilHeightFraction = -1);
157
161
165 static BDSMagnetOuterInfo* PrepareMagnetOuterInfo(const G4String& elementNameIn,
166 const GMAD::Element* el,
167 const G4double angleIn,
168 const G4double angleOut,
169 const BDSBeamPipeInfo* beamPipe,
170 const G4bool yokeOnLeft = false,
171 G4double defaultHorizontalWidth = -1,
172 G4double defaultVHRatio = -1,
173 G4double defaultCoilWidthFraction = -1,
174 G4double defaultCoilHeightFraction = -1);
175
177 static G4Colour* PrepareColour(GMAD::Element const* element);
178
180 static G4Material* PrepareMaterial(GMAD::Element const* element,
181 const G4String& defaultMaterialName);
182
184 static G4Material* PrepareMaterial(GMAD::Element const* element);
185
188 static void CheckBendLengthAngleWidthCombo(G4double arcLength,
189 G4double angle,
190 G4double horizontalWidth,
191 const G4String& name = "not given");
192
194 static void PoleFaceRotationsNotTooLarge(const GMAD::Element* el,
195 G4double maxAngle = 0.5*CLHEP::halfpi);
196
199 BDSCrystalInfo* PrepareCrystalInfo(const G4String& crystalName) const;
200private:
203
205 G4double brho;
206 G4double beta0;
208 G4double lengthSafety;
211 G4bool yokeFields;
212
214 inline void SetBeta0(BDSMagnetStrength* stIn) const {(*stIn)["beta0"] = beta0;}
215
217 GMAD::Element const* element = nullptr;
219 GMAD::Element const* prevElement = nullptr;
221 GMAD::Element const* nextElement = nullptr;
222
224 enum class KickerType {horizontal, vertical, general};
225
226 BDSAcceleratorComponent* CreateDrift(G4double angleIn, G4double angleOut);
227 BDSAcceleratorComponent* CreateRF(G4double currentArcLength);
228 BDSAcceleratorComponent* CreateSBend();
229 BDSAcceleratorComponent* CreateRBend();
230 BDSAcceleratorComponent* CreateKicker(KickerType type);
231 BDSAcceleratorComponent* CreateQuad();
232 BDSAcceleratorComponent* CreateSextupole();
233 BDSAcceleratorComponent* CreateOctupole();
234 BDSAcceleratorComponent* CreateDecapole();
235 BDSAcceleratorComponent* CreateMultipole();
236 BDSAcceleratorComponent* CreateThinMultipole(G4double angleIn);
237 BDSAcceleratorComponent* CreateElement();
238 BDSAcceleratorComponent* CreateSolenoid();
239 BDSAcceleratorComponent* CreateParallelTransporter();
240 BDSAcceleratorComponent* CreateRectangularCollimator();
241 BDSAcceleratorComponent* CreateTarget();
242 BDSAcceleratorComponent* CreateEllipticalCollimator();
243 BDSAcceleratorComponent* CreateJawCollimator();
244 BDSAcceleratorComponent* CreateMuonSpoiler();
245 BDSAcceleratorComponent* CreateShield();
246 BDSAcceleratorComponent* CreateDegrader();
247 BDSAcceleratorComponent* CreateGap();
248 BDSAcceleratorComponent* CreateCrystalCollimator();
249 BDSAcceleratorComponent* CreateLaser();
250 BDSAcceleratorComponent* CreateScreen();
251 BDSAcceleratorComponent* CreateTransform3D();
252 BDSAcceleratorComponent* CreateWireScanner();
253 BDSAcceleratorComponent* CreateRMatrix();
254 BDSAcceleratorComponent* CreateThinRMatrix(G4double angleIn,
255 const G4String& name);
256 BDSAcceleratorComponent* CreateThinRMatrix(G4double angleIn,
257 BDSMagnetStrength* stIn,
258 const G4String& name,
259 BDSIntegratorType intType = BDSIntegratorType::rmatrixthin,
260 BDSFieldType fieldType = BDSFieldType::rmatrix,
261 G4double beamPipeRadius = 0);
262 BDSAcceleratorComponent* CreateUndulator();
263 BDSAcceleratorComponent* CreateDump();
264#ifdef USE_DICOM
265 BDSAcceleratorComponent* CreateCT();
266#endif
267 BDSAcceleratorComponent* CreateCavityFringe(G4double angleIn,
268 BDSMagnetStrength* stIn,
269 const G4String& name,
270 G4double irisRadius);
271
272#ifdef USE_AWAKE
273 BDSAcceleratorComponent* CreateAwakeScreen();
274 BDSAcceleratorComponent* CreateAwakeSpectrometer();
275#endif
276
280 BDSFieldType fieldType,
281 BDSMagnetType magnetType,
282 G4double angle = 0.0,
283 const G4String& nameSuffix = "") const;
284
287 G4bool printWarning = true);
288
290 G4Material* PrepareVacuumMaterial(GMAD::Element const* el) const;
291
294 void PrepareCavityModels();
295
297 void PrepareColours();
298
300 void PrepareCrystals();
301
308 G4double frequency) const;
309
313 G4double frequency) const;
314
318 G4double cavityLength,
319 G4double currentArcLength,
320 BDSMagnetStrength*& fringeIn,
321 BDSMagnetStrength*& fringeOut) const;
322
326 void SetFieldDefinitions(GMAD::Element const* el,
327 BDSAcceleratorComponent* component) const;
328
331
334
336 std::map<G4String, BDSCavityInfo*> cavityInfos;
337
339 std::map<G4String, BDSCrystalInfo*> crystalInfos;
340
343
346
349 G4double FieldFromAngle(const G4double angle,
350 const G4double arcLength) const;
351
354 G4double AngleFromField(const G4double field,
355 const G4double arcLength) const;
356
360 G4double& angle,
361 G4double& field) const;
362
368 G4double& arcLength,
369 G4double& chordLength,
370 G4double& field,
371 G4double& angle) const;
372
374 G4double BendAngle(const GMAD::Element* el) const;
375
381 G4double OutgoingFaceAngle(const GMAD::Element* el) const;
382
388 G4double IncomingFaceAngle(const GMAD::Element* el) const;
389
392 void GetKickValue(G4double& hkick,
393 G4double& vkick,
394 const KickerType type) const;
395
401 std::map<G4String, G4int> modifiedElements;
402
404 G4String elementName;
405
408 static G4bool coloursInitialised;
409};
410#endif
Abstract class that represents a component of an accelerator.
Holder class for all information required to describe a beam pipe model.
Holder for all Geometrical information required to create an RF cavity.
Factory for user specified accelerator components.
Factory to produce all types of BDSAcceleratorComponents.
static G4bool YokeOnLeft(const GMAD::Element *el, const BDSMagnetStrength *st)
std::map< G4String, BDSCrystalInfo * > crystalInfos
Maps of crystal info instances by name.
G4double OutgoingFaceAngle(const GMAD::Element *el) const
G4double AngleFromField(const G4double field, const G4double arcLength) const
BDSComponentFactoryUser * userComponentFactory
User component factory if any.
KickerType
Private enum for kicker types.
G4double FieldFromAngle(const G4double angle, const G4double arcLength) const
static G4double PrepareHorizontalWidth(GMAD::Element const *el, G4double defaultHorizontalWidth=-1)
Prepare the element horizontal width in Geant4 units - if not set, use the global default.
BDSComponentFactory()=delete
No default constructor.
BDSAcceleratorComponent * CreateTeleporter(const G4double teleporterLength, const G4double teleporterWidth, const G4Transform3D &transformIn)
G4Material * PrepareVacuumMaterial(GMAD::Element const *el) const
Prepare the vacuum material from the element or resort to default in options.
void CalculateAngleAndFieldRBend(const GMAD::Element *el, G4double &arcLength, G4double &chordLength, G4double &field, G4double &angle) const
static G4Transform3D CreateFieldTransform(const BDSTiltOffset *tiltOffset)
Create a transform from a tilt offset. If nullptr, returns identity transform.
BDSMagnet * CreateMagnet(const GMAD::Element *el, BDSMagnetStrength *st, BDSFieldType fieldType, BDSMagnetType magnetType, G4double angle=0.0, const G4String &nameSuffix="") const
Helper method for common magnet construction.
BDSMagnetStrength * PrepareMagnetStrengthForMultipoles(GMAD::Element const *el) const
Prepare magnet strength for multipoles.
GMAD::Element const * element
element for storing instead of passing around
static void PoleFaceRotationsNotTooLarge(const GMAD::Element *el, G4double maxAngle=0.5 *CLHEP::halfpi)
Check whether the pole face rotation angles are too big for practical construction.
G4double thinElementLength
Length of a thin element.
BDSCavityInfo * PrepareCavityModelInfoForElement(GMAD::Element const *el, G4double frequency) const
GMAD::Element const * nextElement
element access to previous element (can be nullptr)
BDSCavityInfo * PrepareCavityModelInfo(GMAD::Element const *el, G4double frequency) const
const BDSIntegratorSetType integratorSetType
Local copy of enum of the integrator set, i.e. which one it is specifically.
G4bool yokeFields
Cache of whether to include yoke magnetic fields.
static G4bool coloursInitialised
G4double beta0
Cache of relativisitic beta for primary particle.
void PrepareColours()
Prepare all colours defined in the parser.
static BDSMagnetOuterInfo * PrepareMagnetOuterInfo(const G4String &elementNameIn, const GMAD::Element *el, const BDSMagnetStrength *st, const BDSBeamPipeInfo *beamPipe, G4double defaultHorizontalWidth=-1, G4double defaultVHRatio=1.0, G4double defaultCoilWidthFraction=-1, G4double defaultCoilHeightFraction=-1)
BDSAcceleratorComponent * CreateTerminator(G4double witdth)
static G4Colour * PrepareColour(GMAD::Element const *element)
Checks if colour is specified for element, else uses the default for that element type.
BDSMagnetStrength * PrepareCavityStrength(GMAD::Element const *el, G4double cavityLength, G4double currentArcLength, BDSMagnetStrength *&fringeIn, BDSMagnetStrength *&fringeOut) const
BDSMagnetStrength * PrepareMagnetStrengthForRMatrix(GMAD::Element const *el) const
Prepare magnet strength for rmatrix.
const BDSParticleDefinition * designParticle
Particle w.r.t. which elements are built.
G4bool includeFringeFields
Cache of whether to include fringe fields.
void SetFieldDefinitions(GMAD::Element const *el, BDSAcceleratorComponent *component) const
BDSAcceleratorComponent * CreateComponent(GMAD::Element const *elementIn, GMAD::Element const *prevElementIn, GMAD::Element const *nextElementIn, G4double currentArcLength=0)
static void CheckBendLengthAngleWidthCombo(G4double arcLength, G4double angle, G4double horizontalWidth, const G4String &name="not given")
void PrepareCrystals()
Prepare all crystals in defined the parser.
G4double lengthSafety
Length safety from global constants.
G4String elementName
Variable used to pass around the possibly modified name of an element.
std::map< G4String, G4int > modifiedElements
void GetKickValue(G4double &hkick, G4double &vkick, const KickerType type) const
G4double BendAngle(const GMAD::Element *el) const
Calculate the angle of a bend whether it's an rbend or an sbend.
G4double IncomingFaceAngle(const GMAD::Element *el) const
void CalculateAngleAndFieldSBend(GMAD::Element const *el, G4double &angle, G4double &field) const
static G4Material * PrepareMaterial(GMAD::Element const *element, const G4String &defaultMaterialName)
Checks if a material is named in Element::material, else uses the supplied default.
static BDSTiltOffset * CreateTiltOffset(GMAD::Element const *el)
Create the tilt and offset information object by inspecting the parser element.
static BDSFieldInfo * PrepareMagnetOuterFieldInfo(const BDSMagnetStrength *vacuumSt, const BDSFieldType &fieldType, const BDSBeamPipeInfo *bpInfo, const BDSMagnetOuterInfo *outerInfo, const G4Transform3D &fieldTransform, const BDSIntegratorSet *integratorSetIn, G4double brhoIn, G4double outerFieldScaling=1.0)
Prepare the field definition for the yoke of a magnet.
const BDSIntegratorSet * integratorSet
Local copy of reference to integrator set to use.
G4double brho
Rigidity in T*m (G4units) for beam particles.
void SetBeta0(BDSMagnetStrength *stIn) const
Simple setter used to add Beta0 to a strength instance.
BDSCrystalInfo * PrepareCrystalInfo(const G4String &crystalName) const
static BDSBeamPipeInfo * PrepareBeamPipeInfo(GMAD::Element const *el, const G4ThreeVector &inputFaceNormal=G4ThreeVector(0, 0,-1), const G4ThreeVector &outputFaceNormal=G4ThreeVector(0, 0, 1))
static BDSMagnetGeometryType MagnetGeometryType(const GMAD::Element *el)
std::map< G4String, BDSCavityInfo * > cavityInfos
Map of cavity model info instances by name.
G4bool HasSufficientMinimumLength(GMAD::Element const *el, G4bool printWarning=true)
Test the component length is sufficient for practical construction.
GMAD::Element const * prevElement
element access to previous element (can be nullptr)
Holder for all information required to create a crystal.
All info required to build complete field of any type.
Definition: BDSFieldInfo.hh:65
Which integrator to use for each type of magnet / field object.
Holder struct of all information required to create the outer geometry of a magnet.
Efficient storage of magnet strengths.
Abstract base class that implements features common to all magnets.
Definition: BDSMagnet.hh:45
Wrapper for particle definition.
A holder for any placement offsets and rotations for a BDSAcceleratorComponent.
Parser namespace for GMAD language. Combination of Geant4 and MAD.
Element class.
Definition: element.h:43