BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSUndulator.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 "BDSAcceleratorComponent.hh"
20#include "BDSBeamPipe.hh"
21#include "BDSBeamPipeFactory.hh"
22#include "BDSBeamPipeInfo.hh"
23#include "BDSColours.hh"
24#include "BDSDebug.hh"
25#include "BDSException.hh"
26#include "BDSFieldBuilder.hh"
27#include "BDSFieldInfo.hh"
28#include "BDSGlobalConstants.hh"
29#include "BDSMaterials.hh"
30#include "BDSSDType.hh"
31#include "BDSUndulator.hh"
32#include "BDSUtilities.hh"
33#include "BDSWarning.hh"
34
35#include "globals.hh" // geant4 globals / types
36#include "G4Box.hh"
37#include "G4PVPlacement.hh"
38#include "G4VisAttributes.hh"
39
40#include <cmath>
41
42BDSUndulator::BDSUndulator(const G4String& nameIn,
43 G4double lengthIn,
44 G4double periodIn,
45 G4double undulatorMagnetHeightIn,
46 G4double horizontalWidthIn,
47 G4double undulatorGapIn,
48 BDSBeamPipeInfo* beamPipeInfoIn,
49 BDSFieldInfo* vacuumFieldInfoIn,
50 BDSFieldInfo* outerFieldInfoIn,
51 const G4String& materialIn):
52 BDSAcceleratorComponent(nameIn, lengthIn, 0, "undulator", beamPipeInfoIn),
53 vacuumFieldInfo(vacuumFieldInfoIn),
54 outerFieldInfo(outerFieldInfoIn),
55 undulatorPeriod(periodIn),
56 horizontalWidth(horizontalWidthIn),
57 undulatorMagnetHeight(undulatorMagnetHeightIn),
58 undulatorGap(undulatorGapIn),
59 numMagnets(0)
60{
61 if (materialIn.empty())
62 {
63 BDS::Warning(__METHOD_NAME__, "element \"" + name + "\" no material set for undulator magnet - using iron");
64 material = "iron";
65 }
66 else
67 {material = materialIn;}
68
69 if (vacuumFieldInfo)
70 {vacuumFieldInfo->SetBeamPipeRadius(beamPipeInfoIn->IndicativeRadius());}
71}
72
73BDSUndulator::~BDSUndulator()
74{;}
75
77{
78 // input Checks
80 if (!BDS::IsFinite(undulatorPeriod))
81 {throw BDSException(__METHOD_NAME__, "undulator period is 0, period must be finite.");}
82
83 // check if the undulator period is an integer factor of the element length
84 if (BDS::IsFinite(std::fmod(chordLength, undulatorPeriod)))
85 {throw BDSException(__METHOD_NAME__, "undulator length \"arcLength\" does not divide into an integer number of\n undulator periods (length \"undulatorPeriod\".");}
86
87 // can now cast num magnets to integer as above check should catch if it isnt an integer.
88 numMagnets = (G4int) 2*chordLength/undulatorPeriod;
89
90 G4double beampipeThickness = BDSGlobalConstants::Instance()->DefaultBeamPipeModel()->beamPipeThickness;
92 {
93 G4cout << __METHOD_NAME__ << "\"undulatorGap\" = 0 -> using 2x beam pipe height." << G4endl;
94 undulatorGap = 2*(bp.DY() + 2*beampipeThickness);
95 }
96 if (undulatorGap < (bp.DY() + 2*beampipeThickness + lengthSafetyLarge))
97 {throw BDSException(__METHOD_NAME__, "\"undulatorGap\" for element: \"" + name + "\" smaller than beam pipe aperture.");}
98
100 {throw BDSException(__METHOD_NAME__, "\"undulatorGap\" for element: \"" + name + "\" larger than horizontalWidth.");}
101
103 {
104 // update single magnet box height in case of undulator gap change.
106 }
108 {throw BDSException(__METHOD_NAME__, "\"undulatorMagnetHeight\" larger than 0.5*horizontalWidth in component " + name);}
110 {
111 throw BDSException(__METHOD_NAME__, "Total undulator height (2*undulatorMagnetHeight + undulatorGap) is larger than horizontalWidth in component " + name);
112 }
113
114 G4double halfWidth = 0.5 * (horizontalWidth + lengthSafetyLarge);
115 containerSolid = new G4Box(name + "_container_solid",
116 halfWidth,
117 halfWidth,
118 chordLength*0.5);
119
120 containerLogicalVolume = new G4LogicalVolume(containerSolid,
122 name + "_container_lv");
123
124 BDSExtent ext = BDSExtent(2 * halfWidth, 2 * halfWidth, chordLength);
125 SetExtent(ext);
126}
127
129{
131
133 BDSBeamPipe* pipe = factory->CreateBeamPipe(name, chordLength, beamPipeInfo);
134 RegisterDaughter(pipe);
135
136 G4double singleMagnetLength = (undulatorPeriod * 0.5) - lengthSafetyLarge;
137
138 // magnet geometry
139 G4Box* magnet = new G4Box(name + "_single_magnet_solid",
142 0.5*singleMagnetLength);
143 RegisterSolid(magnet);
144
145 G4Material* materialBox = BDSMaterials::Instance()->GetMaterial(material);
146
147 G4LogicalVolume* lowerBoxLV = new G4LogicalVolume(magnet,
148 materialBox,
149 name + "_lower_box_lv");
150
151 G4LogicalVolume* upperBoxLV = new G4LogicalVolume(magnet,
152 materialBox,
153 name + "_upper_box_lv");
154 RegisterLogicalVolume(lowerBoxLV);
155 RegisterLogicalVolume(upperBoxLV);
156 if (sensitiveOuter)
157 {
158 RegisterSensitiveVolume(lowerBoxLV, BDSSDType::energydep);
159 RegisterSensitiveVolume(upperBoxLV, BDSSDType::energydep);
160 }
161
162 // colour
163 G4VisAttributes* lowerBoxcolour = new G4VisAttributes(*BDSColours::Instance()->GetColour("quadrupole"));
164 lowerBoxLV->SetVisAttributes(lowerBoxcolour);
165 RegisterVisAttributes(lowerBoxcolour);
166
167 G4VisAttributes* upperBoxcolour = new G4VisAttributes(*BDSColours::Instance()->GetColour("sectorbend"));
168 upperBoxLV->SetVisAttributes(upperBoxcolour);
169 RegisterVisAttributes(upperBoxcolour);
170
171 G4double verticalOffset = 0.5 * (undulatorGap + undulatorMagnetHeight);
172 // place upper and lower magnets in a loop
173 for (G4int i = 1; i <= numMagnets; i++)
174 {
175 G4bool sign = BDS::IsFinite(std::fmod(i, 2));
176 G4LogicalVolume* uVol = sign ? upperBoxLV : lowerBoxLV;
177 G4LogicalVolume* lVol = !sign ? upperBoxLV : lowerBoxLV;
178
179 G4ThreeVector upperBoxPos(0, verticalOffset, (0.5*chordLength - undulatorPeriod/4.0) - ((i-1) *undulatorPeriod/2.0));
180 G4ThreeVector lowerBoxPos(0, -verticalOffset, (0.5*chordLength - undulatorPeriod/4.0) - ((i-1) *undulatorPeriod/2.0));
181
182 G4PVPlacement* upperBoxPV = new G4PVPlacement(nullptr, // rotation
183 upperBoxPos, // position
184 uVol, // logical volume
185 name + "_upper_pos_" + std::to_string(i) + "_pv", // name
186 containerLogicalVolume, // mother volume
187 false, // no boolean operation
188 i, // copy number
190
191 G4PVPlacement* lowerBoxPV= new G4PVPlacement(nullptr,
192 lowerBoxPos,
193 lVol,
194 name + "_lower_pos_" + std::to_string(i) + "_pv",
195 containerLogicalVolume,
196 false,
197 i,
199 RegisterPhysicalVolume(upperBoxPV);
200 RegisterPhysicalVolume(lowerBoxPV);
201 }
202
203 // place beam pipe volume
207 G4PVPlacement* bpPV = new G4PVPlacement(nullptr,
208 G4ThreeVector(),
210 name+"_beampipe_pv",
211 containerLogicalVolume,
212 false,
213 0,
216
219 true);
220
221 if (outerFieldInfo)
222 {
223 // Attach to the container but don't propagate to daughter volumes. This ensures
224 // any gap between the beam pipe and the outer also has a field.
226 containerLogicalVolume,
227 false,
228 vacuumFieldInfo->MagnetStrength(),
229 "field");
230 }
231}
Abstract class that represents a component of an accelerator.
void SetAcceleratorVacuumLogicalVolume(G4LogicalVolume *accVacLVIn)
const G4String name
Const protected member variable that may not be changed by derived classes.
static G4Material * emptyMaterial
Useful variable often used in construction.
static G4bool sensitiveOuter
Useful variable often used in construction.
virtual void SetOutputFaceNormal(const G4ThreeVector &output)
Allow updating of face normals. Virtual so derived class may apply it to daughters.
static G4bool checkOverlaps
Useful variable often used in construction.
virtual void SetInputFaceNormal(const G4ThreeVector &input)
Allow updating of face normals. Virtual so derived class may apply it to daughters.
G4double chordLength
Protected member variable that can be modified by derived classes.
BDSBeamPipeInfo * beamPipeInfo
Optional beam pipe recipe that is written out to the survey if it exists.
The main interface for using the beam pipe factories.
static BDSBeamPipeFactory * Instance()
Singleton accessor.
Holder class for all information required to describe a beam pipe model.
BDSExtent Extent() const
G4double beamPipeThickness
Public member for direct access.
G4double IndicativeRadius() const
Return an indicative extent of the beam pipe - typically the maximum of x or y extent.
A holder class for a piece of beam pipe geometry.
Definition: BDSBeamPipe.hh:45
G4LogicalVolume * GetVacuumLogicalVolume() const
Access the vacuum volume to set fields and limits.
Definition: BDSBeamPipe.hh:65
G4ThreeVector InputFaceNormal() const
Accessor.
Definition: BDSBeamPipe.hh:73
G4ThreeVector OutputFaceNormal() const
Accessor.
Definition: BDSBeamPipe.hh:74
static BDSColours * Instance()
singleton pattern
Definition: BDSColours.cc:33
General exception with possible name of object and message.
Definition: BDSException.hh:35
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
G4double DY() const
The difference in a dimension.
Definition: BDSExtent.hh:84
void RegisterFieldForConstruction(const BDSFieldInfo *info, G4LogicalVolume *logicalVolume, const G4bool propagateToDaughters=false, const BDSMagnetStrength *magnetStrengthForScaling=nullptr, const G4String &scalingKey="none")
static BDSFieldBuilder * Instance()
Singleton pattern accessor.
All info required to build complete field of any type.
Definition: BDSFieldInfo.hh:65
BDSMagnetStrength * MagnetStrength() const
Accessor.
G4LogicalVolume * GetContainerLogicalVolume() const
Accessor - see member for more info.
void RegisterDaughter(BDSGeometryComponent *anotherComponent)
void RegisterLogicalVolume(G4LogicalVolume *logicalVolume)
void RegisterPhysicalVolume(G4VPhysicalVolume *physicalVolume)
void SetExtent(const BDSExtent &extIn)
Set extent.
void RegisterVisAttributes(G4VisAttributes *visAttribute)
void RegisterSolid(G4VSolid *solid)
void RegisterSensitiveVolume(G4LogicalVolume *sensitiveVolume, BDSSDType sensitivityType)
static BDSGlobalConstants * Instance()
Access method.
static BDSMaterials * Instance()
Singleton pattern access.
Definition: BDSMaterials.cc:38
G4Material * GetMaterial(G4String material) const
Get material by name.
G4double undulatorGap
Full undulator gap.
Definition: BDSUndulator.hh:64
const G4double horizontalWidth
Element width (and height)
Definition: BDSUndulator.hh:62
G4double undulatorMagnetHeight
Full magnet box height.
Definition: BDSUndulator.hh:63
G4String material
Undulator magnet material.
Definition: BDSUndulator.hh:66
G4int numMagnets
Total number of magnets (1 undulator period is 2 magnets)
Definition: BDSUndulator.hh:65
virtual void BuildContainerLogicalVolume()
Definition: BDSUndulator.cc:76
virtual void Build()
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())