BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSTunnelFactoryRectAboveGround.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 "BDSTunnelFactoryBase.hh"
20#include "BDSTunnelFactoryRectAboveGround.hh"
21#include "BDSTunnelInfo.hh"
22#include "BDSUtilities.hh"
23
24#include "globals.hh" // geant4 globals / types
25#include "G4Box.hh"
26#include "G4CutTubs.hh"
27#include "G4IntersectionSolid.hh"
28#include "G4LogicalVolume.hh"
29#include "G4SubtractionSolid.hh"
30#include "G4ThreeVector.hh"
31#include "G4UnionSolid.hh"
32#include "G4VSolid.hh"
33
34#include <cmath> // sin, cos, fabs
35#include <set>
36#include <utility> // for std::pair
37
39
40BDSTunnelFactoryRectAboveGround::BDSTunnelFactoryRectAboveGround()
41{
42 slabSolid = nullptr;
43 slabXHalfWidth = 0;
44 slabYHalfWidth = 1*CLHEP::m;
45 containerXRadius = 0;
46 containerYRadius = 0;
47}
48
50{
52 containerXRadius = 0;
53 containerYRadius = 0;
54 slabXHalfWidth = 0;
55}
56
58 G4double length,
59 G4double tunnelThickness,
60 G4double tunnelSoilThickness,
61 G4Material* tunnelMaterial,
62 G4Material* tunnelSoilMaterial,
63 G4bool /*tunnelFloor*/,
64 G4double /*tunnelFloorOffset*/,
65 G4double tunnel1,
66 G4double tunnel2,
67 G4bool visible)
68{
69 // tidy up things from previous usage if any - base class method
70 CleanUp();
71
72 // test input parameters - set global options as default if not specified
73 TestInputParameters(length, tunnelThickness, tunnelSoilThickness, tunnelMaterial,
74 tunnelSoilMaterial, tunnel1, tunnel2);
75
76 // build the solids
77 // work out how wide the slab needs to be
78 slabXHalfWidth = std::max(tunnel1, tunnel2) * 2;
79
80 // lengthSafety ensures clean union with floor plane
81 G4double yDisp = tunnel2 + slabYHalfWidth + lengthSafety;
82
83 G4ThreeVector slabDisplacement = G4ThreeVector(0,-yDisp,0);
84
85 slabSolid = new G4Box(name + "_slab_solid", // name
86 slabXHalfWidth, // x half width
87 slabYHalfWidth, // y half width
88 0.5*length - lengthSafety); // z half width
89
90 G4VSolid* tunnelOuterSolid = new G4Box(name + "_tunnel_outer_solid", // name
91 tunnel1 + tunnelThickness, // x radius
92 tunnel2 + tunnelThickness, // y radius
93 0.5*length - lengthSafety); // z half length (to fit in container)
94
95 G4VSolid* tunnelRectSolid = new G4UnionSolid(name + "_tunnel_outer_solid", // name,
96 tunnelOuterSolid,
97 slabSolid,
98 0,
99 slabDisplacement);
100
101 G4VSolid* tunnelInnerSolid = new G4Box(name + "_tunnel_outer_solid", // name
102 tunnel1 + lengthSafety, // x radius
103 tunnel2 + lengthSafety, // y radius
104 length); // z half length - long for unambiguous subtraction
105
106 tunnelSolid = new G4SubtractionSolid(name + "_tunnel_solid",// name
107 tunnelRectSolid, // this
108 tunnelInnerSolid); // minus this
109
110 // register solids
111 allSolids.insert(tunnelOuterSolid);
112 allSolids.insert(tunnelRectSolid);
113 allSolids.insert(tunnelInnerSolid);
114 allSolids.insert(tunnelRectSolid);
115
116 containerSolid = BuildContainerStraight(name, length, tunnel1, tunnel2,
117 tunnelThickness, slabDisplacement);
118
119 CommonConstruction(name, tunnelMaterial, tunnelSoilMaterial, length, containerXRadius, containerYRadius, visible);
120
121 return tunnelSection;
122}
123
125 G4double length,
126 G4ThreeVector inputFace,
127 G4ThreeVector outputFace,
128 G4double tunnelThickness,
129 G4double tunnelSoilThickness,
130 G4Material* tunnelMaterial,
131 G4Material* tunnelSoilMaterial,
132 G4bool /*tunnelFloor*/,
133 G4double /*tunnelFloorOffset*/,
134 G4double tunnel1,
135 G4double tunnel2,
136 G4bool visible)
137{
138 // tidy up things from previous usage if any - base class method
139 CleanUp();
140
141 // test input parameters - set global options as default if not specified
142 TestInputParameters(length, tunnelThickness, tunnelSoilThickness, tunnelMaterial,
143 tunnelSoilMaterial, tunnel1, tunnel2);
144
145 // build the solids
146 // create an intersection cut tubs to get the faces - make it bigger than everything else
147 // then make elliptical solids longer than they need to be
148 G4double intersectionRadius = ( std::max(tunnel1,tunnel2) + tunnelThickness + tunnelSoilThickness ) * 3;
149 G4VSolid* faceSolid = new G4CutTubs(name + "_face_intersection_solid", // name
150 0, // inner radius
151 intersectionRadius, // outer radius
152 length*0.5 - lengthSafety, // z half length
153 0, // start angle
154 CLHEP::twopi, // sweep angle
155 inputFace, // input face normal vector
156 outputFace); // output face normal vector
157
158 // build the solids
159 // work out how wide the slab needs to be
160 slabXHalfWidth = std::max(tunnel1, tunnel2) * 2;
161
162 // lengthSafety ensures clean union with floor plane
163 G4double yDisp = tunnel2 + slabYHalfWidth + lengthSafety;
164
165 G4ThreeVector slabDisplacement = G4ThreeVector(0,-yDisp,0);
166
167 slabSolid = new G4Box(name + "_slab_solid", // name
168 slabXHalfWidth, // x half width
169 slabYHalfWidth, // y half width
170 1.2*length - lengthSafety); // z half width (long for uabig sub.)
171
172
173 // tunnel
174 G4VSolid* tunnelOuterSolid = new G4Box(name + "_tunnel_outer_solid", // name
175 tunnel1 + tunnelThickness, // x radius
176 tunnel2 + tunnelThickness, // y radius
177 length); // z half length
178
179 G4VSolid* tunnelRectSolid = new G4UnionSolid(name + "_tunnel_outer_solid", // name,
180 tunnelOuterSolid,
181 slabSolid,
182 0,
183 slabDisplacement);
184
185 G4VSolid* tunnelInnerSolid = new G4Box(name + "_tunnel_outer_solid", // name
186 tunnel1 + lengthSafety, // x radius
187 tunnel2 + lengthSafety, // y radius
188 length*1.5); // z half length - long for unambiguous subtraction
189
190 G4VSolid* tunnelSolidUnAngled = new G4SubtractionSolid(name + "_tunnel_square_solid", // name
191 tunnelRectSolid, // this
192 tunnelInnerSolid); // minus this
193
194 // register solids
195 allSolids.insert(tunnelOuterSolid);
196 allSolids.insert(tunnelRectSolid);
197 allSolids.insert(tunnelInnerSolid);
198 allSolids.insert(tunnelSolidUnAngled);
199
200 // cut off the faces with the angled face solid
201 tunnelSolid = new G4IntersectionSolid(name + "_tunnel_solid", // name
202 tunnelSolidUnAngled,
203 faceSolid);
204
205 G4VSolid* straightContainer = BuildContainerStraight(name, length*2.0, tunnel1, tunnel2,
206 tunnelThickness, slabDisplacement);
207
208 containerSolid = new G4IntersectionSolid(name + "_tunnel_container_solid", // name
209 straightContainer,
210 faceSolid);
211
212 CommonConstruction(name, tunnelMaterial, tunnelSoilMaterial, length, containerXRadius,
213 containerYRadius, visible);
214
215 return tunnelSection;
216}
217
220 G4double& tunnelThickness,
221 G4double& tunnelSoilThickness,
222 G4Material*& tunnelMaterial,
223 G4Material*& tunnelSoilMaterial,
224 G4double& tunnel1,
225 G4double& tunnel2)
226{
227 CommontTestInputParameters(length, tunnelThickness, tunnelSoilThickness, tunnelMaterial, tunnelSoilMaterial);
228
229 if (!BDS::IsFinite(tunnel1))
230 {tunnel1 = defaultModel->aper1;}
231
232 if (!BDS::IsFinite(tunnel2))
233 {tunnel2 = defaultModel->aper2;}
234}
235
236G4VSolid* BDSTunnelFactoryRectAboveGround::BuildContainerStraight(G4String name,
237 G4double lengthIn,
238 G4double tunnel1,
239 G4double tunnel2,
240 G4double tunnelThickness,
241 G4ThreeVector slabDisplacement)
242{
243 G4VSolid* containerSolidL; // local
244
245 containerXRadius = slabXHalfWidth + lengthSafety;
246 containerYRadius = tunnel2 + tunnelThickness + lengthSafety;
247
248 // have to do a subtraction
249 G4VSolid* tunnelContainerOuterSlab = new G4Box(name + "_cont_slab_solid", // name
250 slabXHalfWidth + lengthSafety, // x half width
251 slabYHalfWidth + lengthSafety, // y half width
252 0.5*lengthIn); // z half width
253 G4VSolid* tunnelContainerOuterTunnel = new G4Box(name + "_cont_tunnel_solid", // name
254 tunnel1 + tunnelThickness + lengthSafety, // x radius
255 tunnel2 + tunnelThickness + lengthSafety, // y radius
256 0.5*lengthIn); // z half length
257 G4VSolid* tunnelContainerOuter = new G4UnionSolid(name + "_cont_outer_solid",
258 tunnelContainerOuterTunnel,
259 tunnelContainerOuterSlab,
260 0,
261 slabDisplacement);
262
263 G4VSolid* tunnelContainerInner = new G4Box(name + "_tunnel_cont_solid_inner", // name
264 tunnel1, // x radius
265 tunnel2, // y radius
266 lengthIn);
267
268 // regsiter solids
269 allSolids.insert(tunnelContainerOuterSlab);
270 allSolids.insert(tunnelContainerOuterTunnel);
271 allSolids.insert(tunnelContainerOuter);
272 allSolids.insert(tunnelContainerInner);
273
274 containerSolidL = new G4SubtractionSolid(name + "_tunnel_cont_solid", // name
275 tunnelContainerOuter, // this
276 tunnelContainerInner); // minus this
277
278 return containerSolidL;
279}
const G4String name
Const protected member variable that may not be changed by derived classes.
G4double lengthSafety
Cache of global constants variable.
void CommontTestInputParameters(G4double &length, G4double &tunnelThickness, G4double &tunnelSoilThickness, G4Material *&tunnelMaterial, G4Material *&tunnelSoilMaterial)
virtual BDSTunnelSection * CreateTunnelSectionAngled(G4String name, G4double length, G4ThreeVector inputFace, G4ThreeVector outputFace, G4double tunnelThickness, G4double tunnelSoilThickness, G4Material *tunnelMaterial, G4Material *tunnelSoilMaterial, G4bool tunnelFloor, G4double tunnelFloorOffset, G4double tunnel1, G4double tunnel2, G4bool visible)
Create a tunnel section with an angled input and output face.
void TestInputParameters(G4double &length, G4double &tunnelThickness, G4double &tunnelSoilThickness, G4Material *&tunnelMaterial, G4Material *&tunnelSoildMaterial, G4double &tunnel1, G4double &tunnel2)
functions below here are private to this particular factory
virtual BDSTunnelSection * CreateTunnelSection(G4String name, G4double length, G4double tunnelThickness, G4double tunnelSoilThickness, G4Material *tunnelMaterial, G4Material *tunnelSoilMaterial, G4bool tunnelFloor, G4double tunnelFloorOffset, G4double tunnel1, G4double tunnel2, G4bool visible)
Create a tunnel section with flat input and output faces.
G4double aper2
Tunnel aperture / shape parameter 2.
G4double aper1
Tunnel aperture / shape parameter 1.
Class that represents a section of tunnel.
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())