BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSGeometryInspector.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 "BDSExtent.hh"
20#include "BDSGeometryInspector.hh"
21#include "BDSGlobalConstants.hh"
22
23#include "globals.hh"
24#include "G4AffineTransform.hh"
25#include "G4Box.hh"
26#include "G4CutTubs.hh"
27#include "G4DisplacedSolid.hh"
28#include "G4EllipticalTube.hh"
29#include "G4IntersectionSolid.hh"
30#include "G4SubtractionSolid.hh"
31#include "G4Tubs.hh"
32#include "G4UnionSolid.hh"
33#include "G4Version.hh"
34#include "G4VSolid.hh"
35
36#include <utility>
37
38// for Geant4.10.2 and below we have a different algorithm as some accessors are not available
39#if G4VERSION_NUMBER < 1030
40#include "G4VisExtent.hh"
41#include <algorithm>
42#include <cmath>
43#endif
44
45std::pair<BDSExtent, BDSExtent> BDS::DetermineExtents(const G4VSolid* solid)
46{
47 G4String className = solid->GetEntityType();
48 if (className == "G4DisplacedSolid")
49 {return BDS::InspectDisplacedSolid(solid);}
50 else if (className == "G4SubtractionSolid")
51 {return BDS::InspectSubtractionSolid(solid);}
52 else if (className == "G4IntersectionSolid")
53 {return BDS::InspectIntersectionSolid(solid);}
54 else if (className == "G4UnionSolid")
55 {return BDS::InspectUnionSolid(solid);}
56 else if (className == "G4Box")
57 {return BDS::InspectBox(solid);}
58 else if (className == "G4Tubs")
59 {return BDS::InspectTubs(solid);}
60 else if (className == "G4CutTubs")
61 {return BDS::InspectCutTubs(solid);}
62 else if (className == "G4EllipticalTube")
63 {return BDS::InspectEllipticalTube(solid);}
64 else
65 {
66 G4cout << "BDS::DetermineExtents> Solid of type \"" << className << "\" is not supported" << G4endl;
67 return std::make_pair(BDSExtent(), BDSExtent());
68 }
69}
70
71std::pair<BDSExtent, BDSExtent> BDS::InspectDisplacedSolid(const G4VSolid* solidIn)
72{
73 const G4DisplacedSolid* solid = dynamic_cast<const G4DisplacedSolid*>(solidIn);
74 if (!solid)
75 {return std::make_pair(BDSExtent(), BDSExtent());}
76
77 const G4VSolid* constituentSolid = solid->GetConstituentMovedSolid();
78
79 auto ext = BDS::DetermineExtents(constituentSolid);
80 G4AffineTransform transform = solid->GetTransform();
81
82 // only going to shift x,y without any rotation.
83 G4ThreeVector translation = transform.NetTranslation();
84 BDSExtent outer = ext.first.Translate(translation.x(), translation.y(), 0);
85 BDSExtent inner = ext.second.Translate(translation.x(), translation.y(), 0);
86
87 return std::make_pair(outer, inner);
88}
89
90std::pair<BDSExtent, BDSExtent> BDS::InspectSubtractionSolid(const G4VSolid* solidIn)
91{
92 const G4SubtractionSolid* solid = dynamic_cast<const G4SubtractionSolid*>(solidIn);
93 if (!solid)
94 {return std::make_pair(BDSExtent(), BDSExtent());}
95
96 const G4VSolid* a = solid->GetConstituentSolid(0);
97 const G4VSolid* b = solid->GetConstituentSolid(1);
98
99 auto extA = BDS::DetermineExtents(a);
100 auto extB = BDS::DetermineExtents(b);
101
102 BDSExtent outer = extA.first.TransverselyLessThan(extB.first) ? extB.first : extA.first;
103
104 // Here we make an ASSUMPTION that the b solid is the interior one - convention.
105 // Also it's the outer extent of that solid.
106 BDSExtent inner = extB.first;
107
108 return std::make_pair(outer, inner);
109}
110
111std::pair<BDSExtent, BDSExtent> BDS::InspectIntersectionSolid(const G4VSolid* solidIn)
112{
113 const G4IntersectionSolid* solid = dynamic_cast<const G4IntersectionSolid*>(solidIn);
114 if (!solid)
115 {return std::make_pair(BDSExtent(), BDSExtent());}
116
117 const G4VSolid* a = solid->GetConstituentSolid(0);
118 const G4VSolid* b = solid->GetConstituentSolid(1);
119
120 auto extA = BDS::DetermineExtents(a);
121 auto extB = BDS::DetermineExtents(b);
122
123 BDSExtent outer = extA.first.TransverselyLessThan(extB.first) ? extB.first : extA.first;
124 BDSExtent inner = extA.first.TransverselyLessThan(extB.first) ? extA.second : extB.second;
125
126 return std::make_pair(outer, inner);
127}
128
129std::pair<BDSExtent, BDSExtent> BDS::InspectUnionSolid(const G4VSolid* solidIn)
130{
131 const G4UnionSolid* solid = dynamic_cast<const G4UnionSolid*>(solidIn);
132 if (!solid)
133 {return std::make_pair(BDSExtent(), BDSExtent());}
134
135 const G4VSolid* a = solid->GetConstituentSolid(0);
136 const G4VSolid* b = solid->GetConstituentSolid(1);
137
138 auto extA = BDS::DetermineExtents(a);
139 auto extB = BDS::DetermineExtents(b);
140
141 BDSExtent outer = extA.first.TransverselyLessThan(extB.first) ? extB.first : extA.first;
142 BDSExtent inner = extA.first.TransverselyLessThan(extB.first) ? extA.second : extB.second;
143
144 return std::make_pair(outer, inner);
145}
146
147std::pair<BDSExtent, BDSExtent> BDS::InspectBox(const G4VSolid* solidIn)
148{
149 const G4Box* solid = dynamic_cast<const G4Box*>(solidIn);
150 if (!solid)
151 {return std::make_pair(BDSExtent(), BDSExtent());}
152
153 G4double dx = solid->GetXHalfLength();
154 G4double dy = solid->GetYHalfLength();
155 G4double dz = solid->GetZHalfLength();
156
157 // accurate along z, but margin in x,y
158 G4double lsl = BDSGlobalConstants::Instance()->LengthSafetyLarge();
159 BDSExtent outer(dx + lsl, dy + lsl, dz); // symmetric +/-
160 BDSExtent inner = BDSExtent();
161 return std::make_pair(outer, inner);
162}
163
164std::pair<BDSExtent, BDSExtent> BDS::InspectTubs(const G4VSolid* solidIn)
165{
166 const G4Tubs* solid = dynamic_cast<const G4Tubs*>(solidIn);
167 if (!solid)
168 {return std::make_pair(BDSExtent(), BDSExtent());}
169
170 G4double innerR = solid->GetInnerRadius();
171 G4double outerR = solid->GetOuterRadius();
172 G4double zHalfL = solid->GetZHalfLength();
173
174 G4double lsl = BDSGlobalConstants::Instance()->LengthSafetyLarge();
175 BDSExtent outer(outerR + lsl, outerR + lsl, zHalfL);
176 BDSExtent inner(innerR, innerR, zHalfL);
177 return std::make_pair(outer, inner);
178}
179
180std::pair<BDSExtent, BDSExtent> BDS::InspectCutTubs(const G4VSolid* solidIn)
181{
182 G4ThreeVector zminV = G4ThreeVector();
183 G4ThreeVector zmaxV = G4ThreeVector();
184#if G4VERSION_NUMBER > 1039
185 solidIn->BoundingLimits(zminV, zmaxV);
186#elif G4VERSION_NUMBER > 1029
187 solidIn->Extent(zminV, zmaxV);
188#else
189 G4VisExtent v = solidIn->GetExtent();
190 zmaxV = G4ThreeVector(std::max(std::abs(v.GetXmin()),std::abs(v.GetXmax())),
191 std::max(std::abs(v.GetYmin()),std::abs(v.GetYmax())),
192 std::max(std::abs(v.GetZmin()),std::abs(v.GetZmax())));
193#endif
194
195 const G4CutTubs* solid = dynamic_cast<const G4CutTubs*>(solidIn);
196 if (!solid)
197 {return std::make_pair(BDSExtent(zmaxV), BDSExtent());}
198
199 G4double innerR = solid->GetInnerRadius();
200 G4double outerR = solid->GetOuterRadius();
201 G4double lsl = BDSGlobalConstants::Instance()->LengthSafetyLarge();
202 outerR += lsl;
203
204 BDSExtent outer(-outerR, outerR,
205 -outerR, outerR,
206 -zmaxV.z(), zmaxV.z());
207 BDSExtent inner(-innerR, innerR,
208 -innerR, innerR,
209 -zmaxV.z(), zmaxV.z());
210 return std::make_pair(outer, inner);
211}
212
213std::pair<BDSExtent, BDSExtent> BDS::InspectEllipticalTube(const G4VSolid* solidIn)
214{
215 const G4EllipticalTube* solid = dynamic_cast<const G4EllipticalTube*>(solidIn);
216 if (!solid)
217 {return std::make_pair(BDSExtent(), BDSExtent());}
218
219 G4double dX = solid->GetDx();
220 G4double dY = solid->GetDy();
221 G4double dZ = solid->GetDz();
222
223 G4double lsl = BDSGlobalConstants::Instance()->LengthSafetyLarge();
224 BDSExtent outer(dX + lsl, dY + lsl, dZ);
225 BDSExtent inner(0, 0, dZ);
226 return std::make_pair(outer, inner);
227}
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
BDSExtent Translate(const G4ThreeVector &offset) const
Provide a new copy of this extent with an offset applied.
Definition: BDSExtent.hh:125
G4bool TransverselyLessThan(const BDSExtent &r) const
Comparison operator for x,y only. Ignores z (length).
Definition: BDSExtent.hh:191
static BDSGlobalConstants * Instance()
Access method.
std::pair< BDSExtent, BDSExtent > InspectIntersectionSolid(const G4VSolid *solidIn)
Inspect a G4IntersectionSolid.
std::pair< BDSExtent, BDSExtent > InspectEllipticalTube(const G4VSolid *solidIn)
Inspect a G4EllipticalTube.
std::pair< BDSExtent, BDSExtent > InspectCutTubs(const G4VSolid *solidIn)
Inspect a G4CutTubs.
std::pair< BDSExtent, BDSExtent > InspectSubtractionSolid(const G4VSolid *solidIn)
std::pair< BDSExtent, BDSExtent > DetermineExtents(const G4VSolid *solid)
std::pair< BDSExtent, BDSExtent > InspectBox(const G4VSolid *solidIn)
Inspect a G4Box.
std::pair< BDSExtent, BDSExtent > InspectTubs(const G4VSolid *solidIn)
Inspect a G4Tubs.
std::pair< BDSExtent, BDSExtent > InspectDisplacedSolid(const G4VSolid *solidIn)
Inspect a G4DisplacedSolid.
std::pair< BDSExtent, BDSExtent > InspectUnionSolid(const G4VSolid *solidIn)
Inspect a G4UnionSolid.