BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSDegrader.cc
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#include "BDSAcceleratorComponent.hh"
20#include "BDSDebug.hh"
21#include "BDSDegrader.hh"
22#include "BDSException.hh"
23#include "BDSSDType.hh"
24
25#include "G4Box.hh"
26#include "G4ExtrudedSolid.hh"
27#include "G4LogicalVolume.hh"
28#include "G4PVPlacement.hh"
29#include "G4RotationMatrix.hh"
30#include "G4TwoVector.hh"
31#include "G4VisAttributes.hh"
32#include "G4VSolid.hh"
33
34#include "globals.hh" // geant4 globals / types
35#include <string>
36#include <vector>
37
38#include "CLHEP/Units/SystemOfUnits.h"
39
40BDSDegrader::BDSDegrader(const G4String& nameIn,
41 G4double lengthIn,
42 G4double horizontalWidthIn,
43 G4int numberWedgesIn,
44 G4double wedgeLengthIn,
45 G4double degraderHeightIn,
46 G4double degraderOffsetIn,
47 G4double baseWidthIn,
48 G4Material* degraderMaterialIn,
49 G4Material* vacuumMaterialIn,
50 G4Colour* colourIn):
51 BDSAcceleratorComponent(nameIn, lengthIn, 0, "degrader"),
52 horizontalWidth(horizontalWidthIn),
53 numberWedges(numberWedgesIn),
54 wedgeLength(wedgeLengthIn),
55 degraderHeight(degraderHeightIn),
56 degraderOffset(degraderOffsetIn),
57 baseWidth(baseWidthIn),
58 degraderMaterial(degraderMaterialIn),
59 vacuumMaterial(vacuumMaterialIn),
60 colour(colourIn)
61{;}
62
63BDSDegrader::~BDSDegrader()
64{;}
65
67{
68 // Input Checks
69 if (horizontalWidth <= 0)
70 {throw BDSException(__METHOD_NAME__,"option \"horizontalWidth\" is not defined or must be greater than 0 for element \"" + name + "\"");}
71
72 if (numberWedges < 1)
73 {throw BDSException(__METHOD_NAME__, "option \"numberWedges\" is not defined or must be greater than 0 for element \"" + name + "\"");}
74
75 if (wedgeLength <= 0)
76 {throw BDSException(__METHOD_NAME__, "option \"wedgeLength\" is not defined or must be greater than 0 for element \"" + name + "\"");}
77
78 if (degraderHeight <= 0)
79 {throw BDSException(__METHOD_NAME__, "option \"degraderHeight\" is not defined or must be greater than 0 for element \"" + name + "\"");}
80
81 if (degraderHeight > (0.5*horizontalWidth))
82 {throw BDSException(__METHOD_NAME__, "option \"degraderHeight\" must be less than 0.5 times \"horizontalWidth\" for element \"" + name + "\"");}
83
84 if ((degraderOffset + baseWidth) > 0.5*horizontalWidth)
85 {throw BDSException(__METHOD_NAME__, "option \"degraderOffset\" must be less than (0.5 times \"horizontalWidth\", minus aper1) for element \"" + name + "\"");}
86
87 // adjust wedge offset to account for added base width
88 degraderOffset += baseWidth;
89
90 // ensure the container is wide enough
91 G4double minimumWidth = wedgeLength;
92 G4double containerWidth = degraderOffset;
93 if (degraderOffset < minimumWidth)
94 {containerWidth = minimumWidth;}
95
96 containerSolid = new G4Box(name + "_container_solid",
97 containerWidth + lengthSafety,
98 degraderHeight*0.5,
99 chordLength*0.5);
100
101 containerLogicalVolume = new G4LogicalVolume(containerSolid,
102 vacuumMaterial,
103 name + "_container_lv");
104}
105
107{
109
110 // width of the wedge at its base
111 G4double wedgeWidth = chordLength/numberWedges - lengthSafety;
112
113 // additional width to make the additional base section triangular
114 // reduce by 10nm to remove annoying angstrom sized overlaps with fewer lengthSafety's
115 G4double addedBasewidth = baseWidth* wedgeWidth/(2.0*wedgeLength) - 10*lengthSafety;
116
117 std::vector<G4TwoVector> fullWedgeVertices; // vertex co-ordinates
118 std::vector<G4TwoVector> halfWedgeVertices; // vertex co-ordinates
119
120 // vertices of full wedge plus base
121 fullWedgeVertices.push_back(G4TwoVector(0.5*wedgeWidth + addedBasewidth, 0));
122 fullWedgeVertices.push_back(G4TwoVector(0, wedgeLength+baseWidth));
123 fullWedgeVertices.push_back(G4TwoVector(-0.5*wedgeWidth - addedBasewidth, 0));
124
125 // vertices of half wedge plus base
126 halfWedgeVertices.push_back(G4TwoVector(0, 0));
127 halfWedgeVertices.push_back(G4TwoVector(0.5*wedgeWidth + addedBasewidth, 0));
128 halfWedgeVertices.push_back(G4TwoVector(0, wedgeLength+baseWidth));
129
130 // full wedge Solid and logical Volume
131 G4ExtrudedSolid* fullWedge = new G4ExtrudedSolid(name + "_fullwedge_solid",
132 fullWedgeVertices,
133 degraderHeight*0.5 - lengthSafety,
134 G4TwoVector(),1, G4TwoVector(), 1);
135 RegisterSolid(fullWedge);
136 G4LogicalVolume* fullWedgeLV = new G4LogicalVolume(fullWedge, // solid
137 degraderMaterial, // material
138 name + "_fullwedge_lv"); // name
139 RegisterLogicalVolume(fullWedgeLV);
140
141 // half wedge Solid and logical Volume
142 G4ExtrudedSolid* halfWedge = new G4ExtrudedSolid(name + "_halfwedge_solid",
143 halfWedgeVertices,
144 degraderHeight*0.5 - lengthSafety,
145 G4TwoVector(),1, G4TwoVector(), 1);
146 RegisterSolid(halfWedge);
147 G4LogicalVolume* halfWedgeLV = new G4LogicalVolume(halfWedge, // solid
148 degraderMaterial, // material
149 name + "_halfwedge_lv"); // name
150 RegisterLogicalVolume(halfWedgeLV);
151
152 // Wedge color
153 G4VisAttributes* degraderVisAttr = new G4VisAttributes(colour);
154 fullWedgeLV->SetVisAttributes(degraderVisAttr);
155 halfWedgeLV->SetVisAttributes(degraderVisAttr);
156 RegisterVisAttributes(degraderVisAttr);
157
158 // Rotation of wedges. Left taken to be +VE x direction, right is -VE x direction.
159 G4RotationMatrix* rightRot = new G4RotationMatrix();
160 rightRot->rotateY(-CLHEP::halfpi);
161 rightRot->rotateX(-CLHEP::halfpi);
162 RegisterRotationMatrix(rightRot);
163
164 G4RotationMatrix* leftRot = new G4RotationMatrix();
165 leftRot->rotateY(CLHEP::halfpi);
166 leftRot->rotateX(-CLHEP::halfpi);
167 RegisterRotationMatrix(leftRot);
168
169 // seperate rotation for first half wedge in even number of wedges.
170 G4RotationMatrix* altRightRot = new G4RotationMatrix();
171 altRightRot->rotateY(CLHEP::halfpi);
172 altRightRot->rotateX(CLHEP::halfpi);
173 RegisterRotationMatrix(altRightRot);
174
175 // loop over number of wedges and place
176 for (G4int i = 0; i < (numberWedges + 1); i++)
177 {
178 G4String wedgeName = name + "_" + std::to_string(i)+"_pv";
179 if (IsEven(numberWedges))
180 {
181 if (i == 0)
182 {PlaceWedge(true, lengthSafety, wedgeName, halfWedgeLV, altRightRot);}
183 else if (i==numberWedges)
184 {PlaceWedge(true, arcLength-lengthSafety, wedgeName, halfWedgeLV, rightRot);}
185 else if(IsEven(i))
186 {PlaceWedge(true, i*wedgeWidth, wedgeName, fullWedgeLV, rightRot);}
187 else if (IsOdd(i))
188 {PlaceWedge(false, i*wedgeWidth, wedgeName, fullWedgeLV, leftRot);}
189 }
190 else if (IsOdd(numberWedges))
191 {
192 if (i == 0)
193 {PlaceWedge(false, lengthSafety, wedgeName, halfWedgeLV, leftRot);}
194 else if (i == numberWedges)
195 {PlaceWedge(true, arcLength-lengthSafety, wedgeName, halfWedgeLV, rightRot);}
196 else if (IsEven(i))
197 {PlaceWedge(false, i*wedgeWidth, wedgeName, fullWedgeLV, leftRot);}
198 else if (IsOdd(i))
199 {PlaceWedge(true, i*wedgeWidth, wedgeName, fullWedgeLV, rightRot);}
200 }
201 }
202
203 if (sensitiveOuter)
204 {
205 RegisterSensitiveVolume(fullWedgeLV, BDSSDType::energydep);
206 RegisterSensitiveVolume(halfWedgeLV, BDSSDType::energydep);
207 }
208}
209
210void BDSDegrader::PlaceWedge(G4bool placeRight,
211 G4double zOffset,
212 const G4String& placementName,
213 G4LogicalVolume* lv,
214 G4RotationMatrix* rotation)
215{
216 G4double xOffset = -degraderOffset; // default is right placement
217 if (!placeRight)
218 {xOffset = degraderOffset;}
219
220 G4VPhysicalVolume* wedgePV = new G4PVPlacement(rotation,
221 G4ThreeVector(xOffset, 0, -0.5*arcLength + zOffset),
222 lv,
223 placementName,
224 containerLogicalVolume,
225 false,
226 0,
228 RegisterPhysicalVolume(wedgePV);
229}
Abstract class that represents a component of an accelerator.
const G4double arcLength
Const protected member variable that may not be changed by derived classes.
const G4String name
Const protected member variable that may not be changed by derived classes.
static G4bool sensitiveOuter
Useful variable often used in construction.
static G4double lengthSafety
Useful variable often used in construction.
static G4bool checkOverlaps
Useful variable often used in construction.
G4double chordLength
Protected member variable that can be modified by derived classes.
void PlaceWedge(G4bool placeRight, G4double zOffset, const G4String &placementName, G4LogicalVolume *lv, G4RotationMatrix *rotation)
Definition: BDSDegrader.cc:210
virtual void Build()
Definition: BDSDegrader.cc:106
BDSDegrader()=delete
Private default constructor to force the use of the supplied one.
virtual void BuildContainerLogicalVolume()
Definition: BDSDegrader.cc:66
General exception with possible name of object and message.
Definition: BDSException.hh:35
void RegisterRotationMatrix(G4RotationMatrix *rotationMatrix)
void RegisterLogicalVolume(G4LogicalVolume *logicalVolume)
void RegisterPhysicalVolume(G4VPhysicalVolume *physicalVolume)
void RegisterVisAttributes(G4VisAttributes *visAttribute)
void RegisterSolid(G4VSolid *solid)
void RegisterSensitiveVolume(G4LogicalVolume *sensitiveVolume, BDSSDType sensitivityType)