BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBendBuilder.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 BDSBENDBUILDER_H
20#define BDSBENDBUILDER_H
21
22#include "globals.hh" // geant4 globals / types
23#include "G4String.hh"
24#include "BDSFieldType.hh"
25#include "BDSIntegratorSetType.hh"
26#include "BDSIntegratorType.hh"
27
29class BDSFieldInfo;
31class BDSLine;
32class BDSMagnet;
34
35namespace GMAD
36{
37 struct Element;
38}
39
40namespace BDS
41{
49 BDSAcceleratorComponent* BuildSBendLine(const G4String& elementName,
50 const GMAD::Element* element,
52 G4double brho,
53 const BDSIntegratorSet* integratorSet,
54 G4double incomingFaceAngle,
55 G4double outgoingFaceAngle,
56 G4bool buildFringeFields,
57 const GMAD::Element* prevElement,
58 const GMAD::Element* nextElement);
59
63 BDSLine* BuildRBendLine(const G4String& elementName,
64 const GMAD::Element* element,
65 const GMAD::Element* prevElement,
66 const GMAD::Element* nextElement,
67 G4double brho,
69 const BDSIntegratorSet* integratorSet,
70 G4double incomingFaceAngle,
71 G4double outgoingFaceAngle,
72 G4bool buildFringeFields);
73
76 G4int CalculateNSBendSegments(G4double length,
77 G4double angle,
78 G4double incomingFaceAngle = 0,
79 G4double outgoingFaceAngle = 0,
80 G4double aperturePrecision = 1.0);
81
85 G4double angleIn,
86 G4double angleOut,
87 const G4String& name,
89 G4double brho,
90 const BDSIntegratorSet* integratorSet,
91 BDSFieldType dipoleFieldType);
92
95 const G4String& name,
96 G4double arcLength,
97 G4double angle,
98 G4double angleIn,
99 G4double angleOut,
100 const BDSMagnetStrength* strength,
101 G4double brho,
102 const BDSIntegratorSet* integratorSet,
103 G4bool yokeOnLeft,
104 const BDSFieldInfo* outerFieldIn);
105
106 void UpdateSegmentAngles(G4int index,
107 G4int nSBends,
108 G4double semiAngle,
109 G4double incomingFaceAngle,
110 G4double outgoingFaceAngle,
111 G4double& segmentAngleIn,
112 G4double& segmentAngleOut);
113
114 BDSMagnetStrength* GetFringeMagnetStrength(const GMAD::Element* element,
115 const BDSMagnetStrength* st,
116 G4double fringeAngle,
117 G4double e1,
118 G4double e2,
119 G4double fintx,
120 G4bool entranceOrExit);
121
125 const GMAD::Element* element);
126
128 G4bool ZeroStrengthDipole(const BDSMagnetStrength* st);
129}
130
131#endif
Abstract class that represents a component of an accelerator.
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.
A class that hold multiple accelerator components.
Definition: BDSLine.hh:38
Efficient storage of magnet strengths.
Abstract base class that implements features common to all magnets.
Definition: BDSMagnet.hh:45
Return either G4Tubs or G4CutTubs depending on flat face.
G4bool ZeroStrengthDipole(const BDSMagnetStrength *st)
Return whether finite angle or field for a dipole.
G4int CalculateNSBendSegments(G4double length, G4double angle, G4double incomingFaceAngle=0, G4double outgoingFaceAngle=0, G4double aperturePrecision=1.0)
BDSMagnet * BuildDipoleFringe(const GMAD::Element *element, G4double angleIn, G4double angleOut, const G4String &name, BDSMagnetStrength *st, G4double brho, const BDSIntegratorSet *integratorSet, BDSFieldType dipoleFieldType)
BDSLine * BuildRBendLine(const G4String &elementName, const GMAD::Element *element, const GMAD::Element *prevElement, const GMAD::Element *nextElement, G4double brho, BDSMagnetStrength *st, const BDSIntegratorSet *integratorSet, G4double incomingFaceAngle, G4double outgoingFaceAngle, G4bool buildFringeFields)
BDSAcceleratorComponent * BuildSBendLine(const G4String &elementName, const GMAD::Element *element, BDSMagnetStrength *st, G4double brho, const BDSIntegratorSet *integratorSet, G4double incomingFaceAngle, G4double outgoingFaceAngle, G4bool buildFringeFields, const GMAD::Element *prevElement, const GMAD::Element *nextElement)
BDSIntegratorType GetDipoleIntegratorType(const BDSIntegratorSet *integratorSet, const GMAD::Element *element)
BDSMagnet * BuildSingleSBend(const GMAD::Element *element, const G4String &name, G4double arcLength, G4double angle, G4double angleIn, G4double angleOut, const BDSMagnetStrength *strength, G4double brho, const BDSIntegratorSet *integratorSet, G4bool yokeOnLeft, const BDSFieldInfo *outerFieldIn)
Function to return a single sector bend section.
Parser namespace for GMAD language. Combination of Geant4 and MAD.
Element class.
Definition: element.h:43