BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSBendBuilder.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 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;
35
36namespace GMAD
37{
38 struct Element;
39}
40
41namespace BDS
42{
50 BDSAcceleratorComponent* BuildSBendLine(const G4String& elementName,
51 const GMAD::Element* element,
53 G4double brho,
54 const BDSIntegratorSet* integratorSet,
55 G4double incomingFaceAngle,
56 G4double outgoingFaceAngle,
57 G4bool buildFringeFields,
58 const GMAD::Element* prevElement,
59 const GMAD::Element* nextElement,
60 BDSModulatorInfo* fieldModulator = nullptr);
61
65 BDSLine* BuildRBendLine(const G4String& elementName,
66 const GMAD::Element* element,
67 const GMAD::Element* prevElement,
68 const GMAD::Element* nextElement,
69 G4double brho,
71 const BDSIntegratorSet* integratorSet,
72 G4double incomingFaceAngle,
73 G4double outgoingFaceAngle,
74 G4bool buildFringeFields,
75 BDSModulatorInfo* fieldModulator = nullptr);
76
79 G4int CalculateNSBendSegments(G4double length,
80 G4double angle,
81 G4double incomingFaceAngle = 0,
82 G4double outgoingFaceAngle = 0,
83 G4double aperturePrecision = 1.0);
84
88 G4double angleIn,
89 G4double angleOut,
90 const G4String& name,
92 G4double brho,
93 const BDSIntegratorSet* integratorSet,
94 BDSFieldType dipoleFieldType,
95 BDSModulatorInfo* fieldModulator = nullptr);
96
99 const G4String& name,
100 G4double arcLength,
101 G4double angle,
102 G4double angleIn,
103 G4double angleOut,
104 const BDSMagnetStrength* strength,
105 G4double brho,
106 const BDSIntegratorSet* integratorSet,
107 G4bool yokeOnLeft,
108 const BDSFieldInfo* outerFieldIn,
109 BDSModulatorInfo* fieldModulator = nullptr);
110
111 void UpdateSegmentAngles(G4int index,
112 G4int nSBends,
113 G4double semiAngle,
114 G4double incomingFaceAngle,
115 G4double outgoingFaceAngle,
116 G4double& segmentAngleIn,
117 G4double& segmentAngleOut);
118
119 BDSMagnetStrength* GetFringeMagnetStrength(const GMAD::Element* element,
120 const BDSMagnetStrength* st,
121 G4double fringeAngle,
122 G4double e1,
123 G4double e2,
124 G4double fintx,
125 G4bool entranceOrExit);
126
130 const GMAD::Element* element);
131
133 G4bool ZeroStrengthDipole(const BDSMagnetStrength* st);
134}
135
136#endif
Abstract class that represents a component of an accelerator.
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.
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
Holder class for all information required to describe a modulator.
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 * 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, BDSModulatorInfo *fieldModulator=nullptr)
Function to return a single sector bend section.
BDSMagnet * BuildDipoleFringe(const GMAD::Element *element, G4double angleIn, G4double angleOut, const G4String &name, BDSMagnetStrength *st, G4double brho, const BDSIntegratorSet *integratorSet, BDSFieldType dipoleFieldType, BDSModulatorInfo *fieldModulator=nullptr)
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, BDSModulatorInfo *fieldModulator=nullptr)
BDSIntegratorType GetDipoleIntegratorType(const BDSIntegratorSet *integratorSet, const GMAD::Element *element)
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, BDSModulatorInfo *fieldModulator=nullptr)
Parser namespace for GMAD language. Combination of Geant4 and MAD.
Element class.
Definition: element.h:43