BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSComponentFactory.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 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;
55class BDSTiltOffset;
56
77{
78public:
79 explicit BDSComponentFactory(const BDSParticleDefinition* designParticleIn,
80 BDSComponentFactoryUser* userComponentFactoryIn = nullptr,
81 G4bool usualPrintOut = true);
83
90 GMAD::Element const* prevElementIn,
91 GMAD::Element const* nextElementIn,
92 G4double currentArcLengthIn = 0);
93
97
100 BDSAcceleratorComponent* CreateTeleporter(const G4double teleporterLength,
101 const G4double teleporterWidth,
102 const G4Transform3D& transformIn);
103
106
108 static G4Transform3D CreateFieldTransform(const BDSTiltOffset* tiltOffset);
109
112 static G4Transform3D CreateFieldTransform(GMAD::Element const* el);
113
117 const G4ThreeVector& inputFaceNormal = G4ThreeVector(0,0,-1),
118 const G4ThreeVector& outputFaceNormal = G4ThreeVector(0,0,1));
119
123 G4double angleIn,
124 G4double angleOut);
125
128 static G4bool YokeOnLeft(const GMAD::Element* el,
129 const BDSMagnetStrength* st);
130
132 static G4double PrepareHorizontalWidth(GMAD::Element const* el,
133 G4double defaultHorizontalWidth = -1);
134
136 static G4double ScalingFieldOuter(const GMAD::Element* ele);
137
140 const BDSFieldType& fieldType,
141 const BDSBeamPipeInfo* bpInfo,
142 const BDSMagnetOuterInfo* outerInfo,
143 const G4Transform3D& fieldTransform,
144 const BDSIntegratorSet* integratorSetIn,
145 G4double brhoIn,
146 G4double outerFieldScaling = 1.0,
147 BDSModulatorInfo* modulatorInfo = nullptr);
148
152 static BDSMagnetOuterInfo* PrepareMagnetOuterInfo(const G4String& elementNameIn,
153 const GMAD::Element* el,
154 const BDSMagnetStrength* st,
155 const BDSBeamPipeInfo* beamPipe,
156 G4double defaultHorizontalWidth = -1,
157 G4double defaultVHRatio = 1.0,
158 G4double defaultCoilWidthFraction = -1,
159 G4double defaultCoilHeightFraction = -1);
160
164
168 static BDSMagnetOuterInfo* PrepareMagnetOuterInfo(const G4String& elementNameIn,
169 const GMAD::Element* el,
170 const G4double angleIn,
171 const G4double angleOut,
172 const BDSBeamPipeInfo* beamPipe,
173 const G4bool yokeOnLeft = false,
174 G4double defaultHorizontalWidth = -1,
175 G4double defaultVHRatio = -1,
176 G4double defaultCoilWidthFraction = -1,
177 G4double defaultCoilHeightFraction = -1);
178
180 static G4Colour* PrepareColour(GMAD::Element const* element);
181
183 static G4Material* PrepareMaterial(GMAD::Element const* element,
184 const G4String& defaultMaterialName);
185
187 static G4Material* PrepareMaterial(GMAD::Element const* element);
188
191 static void CheckBendLengthAngleWidthCombo(G4double arcLength,
192 G4double angle,
193 G4double horizontalWidth,
194 const G4String& name = "not given");
195
197 static void PoleFaceRotationsNotTooLarge(const GMAD::Element* el,
198 G4double maxAngle = 0.5*CLHEP::halfpi);
199
202 static G4double EFieldFromElement(GMAD::Element const* el,
203 G4double cavityLength);
204
207 BDSCrystalInfo* PrepareCrystalInfo(const G4String& crystalName) const;
208
209private:
212
214 G4double brho;
215 G4double beta0;
217 G4double lengthSafety;
220 G4bool yokeFields;
222
225
227 inline void SetBeta0(BDSMagnetStrength* stIn) const {(*stIn)["beta0"] = beta0;}
228
230 GMAD::Element const* element = nullptr;
232 GMAD::Element const* prevElement = nullptr;
234 GMAD::Element const* nextElement = nullptr;
235
237 enum class KickerType {horizontal, vertical, general};
238
240 enum class RFFieldDirection {x, y, z};
241
242 BDSAcceleratorComponent* CreateDrift(G4double angleIn, G4double angleOut);
243 BDSAcceleratorComponent* CreateRF(RFFieldDirection direction);
244 BDSAcceleratorComponent* CreateSBend();
245 BDSAcceleratorComponent* CreateRBend();
246 BDSAcceleratorComponent* CreateKicker(KickerType type);
247 BDSAcceleratorComponent* CreateQuad();
248 BDSAcceleratorComponent* CreateSextupole();
249 BDSAcceleratorComponent* CreateOctupole();
250 BDSAcceleratorComponent* CreateDecapole();
251 BDSAcceleratorComponent* CreateMultipole();
252 BDSAcceleratorComponent* CreateThinMultipole(G4double angleIn);
253 BDSAcceleratorComponent* CreateElement();
254 BDSAcceleratorComponent* CreateSolenoid();
255 BDSAcceleratorComponent* CreateParallelTransporter();
256 BDSAcceleratorComponent* CreateRectangularCollimator();
257 BDSAcceleratorComponent* CreateTarget();
258 BDSAcceleratorComponent* CreateEllipticalCollimator();
259 BDSAcceleratorComponent* CreateJawCollimator();
260 BDSAcceleratorComponent* CreateMuonSpoiler();
261 BDSAcceleratorComponent* CreateShield();
262 BDSAcceleratorComponent* CreateDegrader();
263 BDSAcceleratorComponent* CreateGap();
264 BDSAcceleratorComponent* CreateCrystalCollimator();
265 BDSAcceleratorComponent* CreateLaser();
266 BDSAcceleratorComponent* CreateScreen();
267 BDSAcceleratorComponent* CreateTransform3D();
268 BDSAcceleratorComponent* CreateWireScanner();
269 BDSAcceleratorComponent* CreateRMatrix();
270 BDSAcceleratorComponent* CreateThinRMatrix(G4double angleIn,
271 const G4String& name);
272 BDSAcceleratorComponent* CreateThinRMatrix(G4double angleIn,
273 BDSMagnetStrength* stIn,
274 const G4String& name,
275 BDSIntegratorType intType = BDSIntegratorType::rmatrixthin,
276 BDSFieldType fieldType = BDSFieldType::rmatrix,
277 G4double beamPipeRadius = 0,
278 BDSModulatorInfo* fieldModulator = nullptr);
279 BDSAcceleratorComponent* CreateUndulator();
280 BDSAcceleratorComponent* CreateDump();
281#ifdef USE_DICOM
282 BDSAcceleratorComponent* CreateCT();
283#endif
284 BDSAcceleratorComponent* CreateCavityFringe(G4double angleIn,
285 BDSMagnetStrength* stIn,
286 const G4String& name,
287 G4double irisRadius,
288 BDSModulatorInfo* fieldModulator = nullptr);
289
290#ifdef USE_AWAKE
291 BDSAcceleratorComponent* CreateAwakeScreen();
292 BDSAcceleratorComponent* CreateAwakeSpectrometer();
293#endif
294
298 BDSFieldType fieldType,
299 BDSMagnetType magnetType,
300 G4double angle = 0.0,
301 const G4String& nameSuffix = "") const;
302
305 G4bool printWarning = true);
306
308 G4Material* PrepareVacuumMaterial(GMAD::Element const* el) const;
309
312 void PrepareCavityModels();
313
315 void PrepareColours();
316
318 void PrepareCrystals();
319
326 G4double frequency) const;
327
331 G4double frequency) const;
332
336 BDSFieldType fieldType,
337 G4double cavityLength,
338 BDSMagnetStrength*& fringeIn,
339 BDSMagnetStrength*& fringeOut) const;
340
344 void SetFieldDefinitions(GMAD::Element const* el,
345 BDSAcceleratorComponent* component) const;
346
350 BDSFieldInfo* info) const;
351
354
357
359 std::map<G4String, BDSCavityInfo*> cavityInfos;
360
362 std::map<G4String, BDSCrystalInfo*> crystalInfos;
363
366
369
372 G4double FieldFromAngle(const G4double angle,
373 const G4double arcLength) const;
374
377 G4double AngleFromField(const G4double field,
378 const G4double arcLength) const;
379
383 G4double& angle,
384 G4double& field) const;
385
391 G4double& arcLength,
392 G4double& chordLength,
393 G4double& field,
394 G4double& angle) const;
395
397 G4double BendAngle(const GMAD::Element* el) const;
398
404 G4double OutgoingFaceAngle(const GMAD::Element* el) const;
405
411 G4double IncomingFaceAngle(const GMAD::Element* el) const;
412
415 G4double elementArcLength) const;
416
419 BDSModulatorInfo* ModulatorDefinition(const GMAD::Element* el, G4bool inDevelopment=false) const; // TBC
420
422 void INDEVELOPMENTERROR() const;
423
426 void GetKickValue(G4double& hkick,
427 G4double& vkick,
428 const KickerType type) const;
429
435 std::map<G4String, G4int> modifiedElements;
436
438 G4String elementName;
439
442 static G4bool coloursInitialised;
443};
444#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
void INDEVELOPMENTERROR() const
TBC - remove when modulators are implemented fully.
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)
BDSModulatorInfo * defaultModulator
Default modulator for all components.
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
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, BDSModulatorInfo *modulatorInfo=nullptr)
Prepare the field definition for the yoke of a magnet.
G4double beta0
Cache of relativistic 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 * PrepareMagnetStrengthForRMatrix(GMAD::Element const *el) const
Prepare magnet strength for rmatrix.
G4double currentArcLength
Updated each time CreateComponent is called - supplied from outside. Only here to pass around all fun...
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
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.
BDSMagnetStrength * PrepareCavityStrength(GMAD::Element const *el, BDSFieldType fieldType, G4double cavityLength, BDSMagnetStrength *&fringeIn, BDSMagnetStrength *&fringeOut) const
G4String elementName
Variable used to pass around the possibly modified name of an element.
void AddSynchronousTimeInformation(BDSMagnetStrength *st, G4double elementArcLength) const
Update the BDSMagnetStrength key synchronousT0 with the time at the centre of the element.
BDSAcceleratorComponent * CreateComponent(GMAD::Element const *elementIn, GMAD::Element const *prevElementIn, GMAD::Element const *nextElementIn, G4double currentArcLengthIn=0)
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.
void SetModulatorDefinition(GMAD::Element const *el, BDSFieldInfo *info) const
const BDSIntegratorSet * integratorSet
Local copy of reference to integrator set to use.
BDSModulatorInfo * ModulatorDefinition(const GMAD::Element *el, G4bool inDevelopment=false) const
RFFieldDirection
Private enum for RF cavity principle accelerating direction.
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)
static G4double EFieldFromElement(GMAD::Element const *el, G4double cavityLength)
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.
static G4double ScalingFieldOuter(const GMAD::Element *ele)
Get the scaling factor for a particular outer field depending on the global and individual setting.
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:66
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
Holder class for all information required to describe a modulator.
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