BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSBeamPipeFactoryCircular.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 "BDSBeamPipeFactoryCircular.hh"
20#include "BDSBeamPipe.hh"
21#include "BDSExtent.hh"
22#include "BDSUtilities.hh"
23
24#include "globals.hh" // geant4 globals / types
25#include "G4CutTubs.hh"
26#include "G4LogicalVolume.hh"
27#include "G4SubtractionSolid.hh"
28#include "G4ThreeVector.hh"
29#include "G4Tubs.hh"
30#include "G4VSolid.hh"
31
32#include <cmath>
33#include <set>
34#include <utility> // for std::pair
35
36BDSBeamPipeFactoryCircular::BDSBeamPipeFactoryCircular()
37{;}
38
40 G4double lengthIn,
41 G4double aper1In,
42 G4double /*aper2In*/,
43 G4double /*aper3In*/,
44 G4double /*aper4In*/,
45 G4Material* vacuumMaterialIn,
46 G4double beamPipeThicknessIn,
47 G4Material* beamPipeMaterialIn,
48 const G4String& /*pointsFileIn*/,
49 const G4String& /*pointsUnitIn*/)
50{
51 // clean up after last usage
52 CleanUp();
53
54 // build the solids
55 vacuumSolid = new G4Tubs(nameIn + "_vacuum_solid", // name
56 0, // inner radius
57 aper1In, // outer radius
58 lengthIn*0.5 - lengthSafety, // half length
59 0, // rotation start angle
60 CLHEP::twopi); // rotation finish angle
61
62 beamPipeSolid = new G4Tubs(nameIn + "_pipe_solid", // name
63 aper1In + lengthSafetyLarge, // inner radius + length safety to avoid overlaps
64 aper1In + lengthSafetyLarge + beamPipeThicknessIn, // outer radius
65 lengthIn*0.5 - lengthSafety, // half length
66 0, // rotation start angle
67 CLHEP::twopi); // rotation finish angle
68
69 G4double containerRadius = aper1In + beamPipeThicknessIn + 2*lengthSafetyLarge;
70 containerSolid = new G4Tubs(nameIn + "_container_solid", // name
71 0, // inner radius
72 containerRadius, // outer radius
73 lengthIn*0.5, // half length
74 0, // rotation start angle
75 CLHEP::twopi); // rotation finish angle
76
77 return CommonFinalConstruction(nameIn, vacuumMaterialIn, beamPipeMaterialIn,
78 lengthIn, containerRadius);
79}
80
82 G4double lengthIn,
83 const G4ThreeVector& inputFaceNormalIn,
84 const G4ThreeVector& outputFaceNormalIn,
85 G4double aper1In,
86 G4double /*aper2In*/,
87 G4double /*aper3In*/,
88 G4double /*aper4In */,
89 G4Material* vacuumMaterialIn,
90 G4double beamPipeThicknessIn,
91 G4Material* beamPipeMaterialIn,
92 const G4String& /*pointsFileIn*/,
93 const G4String& /*pointsUnitIn*/)
94{
95 // clean up after last usage
96 CleanUp();
97 inputFaceNormal = inputFaceNormalIn;
98 outputFaceNormal = outputFaceNormalIn;
99
100 G4double containerRadius = 0;
101 CreateGeneralAngledSolids(nameIn, lengthIn, aper1In, beamPipeThicknessIn,
102 inputFaceNormal, outputFaceNormal, containerRadius);
103
104 return CommonFinalConstruction(nameIn, vacuumMaterialIn, beamPipeMaterialIn,
105 lengthIn, containerRadius);
106}
107
111 G4Material* vacuumMaterialIn,
112 G4Material* beamPipeMaterialIn,
113 G4double lengthIn,
114 G4double containerRadiusIn)
115{
116 // prepare a longer container subtraction solid
117 // doesn't have to be angled as it's only used for transverse subtraction
118 containerSubtractionSolid = new G4Tubs(nameIn + "_container_sub_solid",// name
119 0, // inner radius
120 containerRadiusIn, // outer radius
121 lengthIn*4, // full length for unambiguous subtraction
122 0, // rotation start angle
123 CLHEP::twopi); // rotation finish angle
124
125 BDSBeamPipeFactoryBase::CommonConstruction(nameIn, vacuumMaterialIn,
126 beamPipeMaterialIn, lengthIn);
127
128 // record extents
129 BDSExtent ext = BDSExtent(containerRadiusIn, containerRadiusIn, lengthIn*0.5);
130
131 // true for containerIsCircular - true for this factory
132 return BDSBeamPipeFactoryBase::BuildBeamPipeAndRegisterVolumes(ext, containerRadiusIn, true);
133}
134
138 G4double lengthIn,
139 G4double aper1In,
140 G4double beamPipeThicknessIn,
141 const G4ThreeVector& inputfaceIn,
142 const G4ThreeVector& outputfaceIn,
143 G4double& containerRadius)
144{
145 // long length for unambiguous boolean - ensure no gaps in beam pipe geometry
146 G4double angledVolumeLength = BDS::CalculateSafeAngledVolumeLength(inputfaceIn, outputfaceIn, lengthIn, aper1In);
147
148 G4double extraWidth = lengthSafetyLarge + beamPipeThicknessIn;
149 containerRadius = aper1In + extraWidth + lengthSafetyLarge;
150
151 // check faces of angled container volume don't intersect - if it can be built, remaining angled volumes can be built
152 CheckAngledVolumeCanBeBuilt(lengthIn, inputfaceIn, outputfaceIn, containerRadius, nameIn);
153
154 containerSolid = new G4CutTubs(nameIn + "_container_solid", // name
155 0, // inner radius
156 containerRadius, // outer radius
157 lengthIn*0.5, // half length - no -length safety!
158 0, // rotation start angle
159 CLHEP::twopi, // rotation finish angle
160 inputfaceIn, // input face normal
161 outputfaceIn); // rotation finish angle
162
163 // build the solids
164 vacuumSolid = new G4CutTubs(nameIn + "_vacuum_solid", // name
165 0, // inner radius
166 aper1In, // outer radius
167 lengthIn*0.5 - lengthSafety, // half length
168 0, // rotation start angle
169 CLHEP::twopi, // rotation finish angle
170 inputfaceIn, // input face normal
171 outputfaceIn ); // output face normal
172
173 // beampipesolid created as subtraction since direct G4CutTubs creation
174 // created scattering in sector bends. not really understood but likely
175 // fault in G4CutTubs
176 G4VSolid* inner = new G4CutTubs(nameIn + "_pipe_inner_solid", // name
177 0, // inner radius
178 aper1In + lengthSafetyLarge, // outer radius
179 angledVolumeLength, // long length!
180 0, // rotation start angle
181 CLHEP::twopi, // rotation finish angle
182 inputfaceIn, // input face normal
183 outputfaceIn); // output face normal
184
185 G4VSolid* outer = new G4CutTubs(nameIn + "_pipe_outer_solid", // name
186 0, // inner radius + length safety to avoid overlaps
187 aper1In + extraWidth, // outer radius
188 lengthIn*0.5 - lengthSafety, // half length
189 0, // rotation start angle
190 CLHEP::twopi, // rotation finish angle
191 inputfaceIn, // input face normal
192 outputfaceIn); // output face normal
193 allSolids.insert(inner);
194 allSolids.insert(outer);
195
196 beamPipeSolid = new G4SubtractionSolid(nameIn + "_pipe_solid",
197 outer,
198 inner);
199}
G4ThreeVector outputFaceNormal
For recording the face normals in the finished pipe component.
G4ThreeVector inputFaceNormal
For recording the face normals in the finished pipe component.
void CommonConstruction(const G4String &nameIn, G4Material *vacuumMaterialIn, G4Material *beamPipeMaterialIn, G4double length)
finalise beampipe construction
BDSBeamPipe * BuildBeamPipeAndRegisterVolumes(BDSExtent extent, G4double containerRadius, G4bool containerIsCircular=false)
build beampipe and register logical volumes
static void CheckAngledVolumeCanBeBuilt(G4double length, const G4ThreeVector &inputfaceAngle, const G4ThreeVector &outputfaceAngle, G4double horizontalWidth, const G4String &name)
check if a beam pipe volume with angled faces can be constructed
G4VSolid * containerSubtractionSolid
Longer (in length) version of container solid for unambiguous subtraction.
virtual BDSBeamPipe * CreateBeamPipe(const G4String &nameIn, G4double lengthIn, G4double aper1=0, G4double aper2=0, G4double aper3=0, G4double aper4=0, G4Material *vacuumMaterialIn=nullptr, G4double beamPipeThicknessIn=0, G4Material *beamPipeMaterialIn=nullptr, const G4String &pointsFileIn="", const G4String &pointsUnitIn="")
create a flat ended beampipe
void CreateGeneralAngledSolids(const G4String &nameIn, G4double lengthIn, G4double aper1In, G4double beamPipeThicknessIn, const G4ThreeVector &inputfaceIn, const G4ThreeVector &outputfaceIn, G4double &containerRadius)
Pass containerRadius by reference to update.
BDSBeamPipe * CommonFinalConstruction(const G4String &nameIn, G4Material *vacuumMaterialIn, G4Material *beamPipeMaterialIn, G4double lengthIn, G4double containerRadiusIn)
A holder class for a piece of beam pipe geometry.
Definition: BDSBeamPipe.hh:45
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
G4double lengthSafety
Cache of global constants variable.
G4double lengthSafetyLarge
Cache of global constants variable.
G4double CalculateSafeAngledVolumeLength(G4double angleIn, G4double angleOut, G4double length, G4double containerWidth, G4double containerHeight=0)
Calculate safe length of an angled volume so it fills the length of its container.