Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4USphere.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// Implementation for G4USphere wrapper class
27//
28// 13.09.13 G.Cosmo, CERN/PH
29// --------------------------------------------------------------------
30
31#include "G4Sphere.hh"
32#include "G4USphere.hh"
33
34#if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
35
36#include "G4GeomTools.hh"
37#include "G4AffineTransform.hh"
39#include "G4BoundingEnvelope.hh"
40
41using namespace CLHEP;
42
43////////////////////////////////////////////////////////////////////////
44//
45// constructor - check parameters, convert angles so 0<sphi+dpshi<=2_PI
46// - note if pDPhi>2PI then reset to 2PI
47
48G4USphere::G4USphere( const G4String& pName,
49 G4double pRmin, G4double pRmax,
50 G4double pSPhi, G4double pDPhi,
51 G4double pSTheta, G4double pDTheta )
52 : Base_t(pName, pRmin, pRmax, pSPhi, pDPhi, pSTheta, pDTheta)
53{
54}
55
56///////////////////////////////////////////////////////////////////////
57//
58// Fake default constructor - sets only member data and allocates memory
59// for usage restricted to object persistency.
60//
61G4USphere::G4USphere( __void__& a )
62 : Base_t(a)
63{
64}
65
66/////////////////////////////////////////////////////////////////////
67//
68// Destructor
69
70G4USphere::~G4USphere() = default;
71
72//////////////////////////////////////////////////////////////////////////
73//
74// Copy constructor
75
76G4USphere::G4USphere(const G4USphere& rhs)
77 : Base_t(rhs)
78{
79}
80
81//////////////////////////////////////////////////////////////////////////
82//
83// Assignment operator
84
85G4USphere& G4USphere::operator = (const G4USphere& rhs)
86{
87 // Check assignment to self
88 //
89 if (this == &rhs) { return *this; }
90
91 // Copy base class data
92 //
93 Base_t::operator=(rhs);
94
95 return *this;
96}
97
98//////////////////////////////////////////////////////////////////////////
99//
100// Accessors & modifiers
101
102G4double G4USphere::GetInnerRadius() const
103{
104 return Base_t::GetInnerRadius();
105}
106G4double G4USphere::GetOuterRadius() const
107{
108 return Base_t::GetOuterRadius();
109}
110G4double G4USphere::GetStartPhiAngle() const
111{
112 return Base_t::GetStartPhiAngle();
113}
114G4double G4USphere::GetDeltaPhiAngle() const
115{
116 return Base_t::GetDeltaPhiAngle();
117}
118G4double G4USphere::GetStartThetaAngle() const
119{
120 return Base_t::GetStartThetaAngle();
121}
122G4double G4USphere::GetDeltaThetaAngle() const
123{
124 return Base_t::GetDeltaThetaAngle();
125}
126G4double G4USphere::GetSinStartPhi() const
127{
128 return Base_t::GetSinSPhi();
129}
130G4double G4USphere::GetCosStartPhi() const
131{
132 return Base_t::GetCosSPhi();
133}
134G4double G4USphere::GetSinEndPhi() const
135{
136 return Base_t::GetSinEPhi();
137}
138G4double G4USphere::GetCosEndPhi() const
139{
140 return Base_t::GetCosEPhi();
141}
142G4double G4USphere::GetSinStartTheta() const
143{
144 return Base_t::GetSinSTheta();
145}
146G4double G4USphere::GetCosStartTheta() const
147{
148 return Base_t::GetCosSTheta();
149}
150G4double G4USphere::GetSinEndTheta() const
151{
152 return Base_t::GetSinETheta();
153}
154G4double G4USphere::GetCosEndTheta() const
155{
156 return Base_t::GetCosETheta();
157}
158
159void G4USphere::SetInnerRadius(G4double newRMin)
160{
161 Base_t::SetInnerRadius(newRMin);
162 fRebuildPolyhedron = true;
163}
164void G4USphere::SetOuterRadius(G4double newRmax)
165{
166 Base_t::SetOuterRadius(newRmax);
167 fRebuildPolyhedron = true;
168}
169void G4USphere::SetStartPhiAngle(G4double newSphi, G4bool trig)
170{
171 Base_t::SetStartPhiAngle(newSphi, trig);
172 fRebuildPolyhedron = true;
173}
174void G4USphere::SetDeltaPhiAngle(G4double newDphi)
175{
176 Base_t::SetDeltaPhiAngle(newDphi);
177 fRebuildPolyhedron = true;
178}
179void G4USphere::SetStartThetaAngle(G4double newSTheta)
180{
181 Base_t::SetStartThetaAngle(newSTheta);
182 fRebuildPolyhedron = true;
183}
184void G4USphere::SetDeltaThetaAngle(G4double newDTheta)
185{
186 Base_t::SetDeltaThetaAngle(newDTheta);
187 fRebuildPolyhedron = true;
188}
189
190//////////////////////////////////////////////////////////////////////////
191//
192// Dispatch to parameterisation for replication mechanism dimension
193// computation & modification.
194
195void G4USphere::ComputeDimensions( G4VPVParameterisation* p,
196 const G4int n,
197 const G4VPhysicalVolume* pRep)
198{
199 p->ComputeDimensions(*(G4Sphere*)this,n,pRep);
200}
201
202//////////////////////////////////////////////////////////////////////////
203//
204// Make a clone of the object
205
206G4VSolid* G4USphere::Clone() const
207{
208 return new G4USphere(*this);
209}
210
211//////////////////////////////////////////////////////////////////////////
212//
213// Get bounding box
214
215void G4USphere::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
216{
217 static G4bool checkBBox = true;
218
219 G4double rmin = GetInnerRadius();
220 G4double rmax = GetOuterRadius();
221
222 // Find bounding box
223 //
224 if (GetDeltaThetaAngle() >= pi && GetDeltaPhiAngle() >= twopi)
225 {
226 pMin.set(-rmax,-rmax,-rmax);
227 pMax.set( rmax, rmax, rmax);
228 }
229 else
230 {
231 G4double sinStart = GetSinStartTheta();
232 G4double cosStart = GetCosStartTheta();
233 G4double sinEnd = GetSinEndTheta();
234 G4double cosEnd = GetCosEndTheta();
235
236 G4double stheta = GetStartThetaAngle();
237 G4double etheta = stheta + GetDeltaThetaAngle();
238 G4double rhomin = rmin*std::min(sinStart,sinEnd);
239 G4double rhomax = rmax;
240 if (stheta > halfpi) rhomax = rmax*sinStart;
241 if (etheta < halfpi) rhomax = rmax*sinEnd;
242
243 G4TwoVector xymin,xymax;
244 G4GeomTools::DiskExtent(rhomin,rhomax,
245 GetSinStartPhi(),GetCosStartPhi(),
246 GetSinEndPhi(),GetCosEndPhi(),
247 xymin,xymax);
248
249 G4double zmin = std::min(rmin*cosEnd,rmax*cosEnd);
250 G4double zmax = std::max(rmin*cosStart,rmax*cosStart);
251 pMin.set(xymin.x(),xymin.y(),zmin);
252 pMax.set(xymax.x(),xymax.y(),zmax);
253 }
254
255 // Check correctness of the bounding box
256 //
257 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
258 {
259 std::ostringstream message;
260 message << "Bad bounding box (min >= max) for solid: "
261 << GetName() << " !"
262 << "\npMin = " << pMin
263 << "\npMax = " << pMax;
264 G4Exception("G4USphere::BoundingLimits()", "GeomMgt0001",
265 JustWarning, message);
266 StreamInfo(G4cout);
267 }
268
269 // Check consistency of bounding boxes
270 //
271 if (checkBBox)
272 {
273 U3Vector vmin, vmax;
274 Extent(vmin,vmax);
275 if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
276 std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
277 std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
278 std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
279 std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
280 std::abs(pMax.z()-vmax.z()) > kCarTolerance)
281 {
282 std::ostringstream message;
283 message << "Inconsistency in bounding boxes for solid: "
284 << GetName() << " !"
285 << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
286 << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
287 G4Exception("G4USphere::BoundingLimits()", "GeomMgt0001",
288 JustWarning, message);
289 checkBBox = false;
290 }
291 }
292}
293
294//////////////////////////////////////////////////////////////////////////
295//
296// Calculate extent under transform and specified limit
297
298G4bool G4USphere::CalculateExtent(const EAxis pAxis,
299 const G4VoxelLimits& pVoxelLimit,
300 const G4AffineTransform& pTransform,
301 G4double& pMin, G4double& pMax) const
302{
303 G4ThreeVector bmin, bmax;
304
305 // Get bounding box
306 BoundingLimits(bmin,bmax);
307
308 // Find extent
309 G4BoundingEnvelope bbox(bmin,bmax);
310 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
311}
312
313//////////////////////////////////////////////////////////////////////////
314//
315// Create polyhedron for visualization
316
317G4Polyhedron* G4USphere::CreatePolyhedron() const
318{
319 return new G4PolyhedronSphere(GetInnerRadius(),
320 GetOuterRadius(),
321 GetStartPhiAngle(),
322 GetDeltaPhiAngle(),
323 GetStartThetaAngle(),
324 GetDeltaThetaAngle());
325}
326
327#endif // G4GEOM_USE_USOLIDS
const G4double kCarTolerance
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
double x() const
double y() const
double z() const
double x() const
double y() const
void set(double x, double y, double z)
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
EAxis
Definition geomdefs.hh:54