Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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
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
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
67//
68// Destructor
69
70G4USphere::~G4USphere()
71{
72}
73
75//
76// Copy constructor
77
78G4USphere::G4USphere(const G4USphere& rhs)
79 : Base_t(rhs)
80{
81}
82
84//
85// Assignment operator
86
87G4USphere& G4USphere::operator = (const G4USphere& rhs)
88{
89 // Check assignment to self
90 //
91 if (this == &rhs) { return *this; }
92
93 // Copy base class data
94 //
95 Base_t::operator=(rhs);
96
97 return *this;
98}
99
101//
102// Accessors & modifiers
103
104G4double G4USphere::GetInnerRadius() const
105{
106 return Base_t::GetInnerRadius();
107}
108G4double G4USphere::GetOuterRadius() const
109{
110 return Base_t::GetOuterRadius();
111}
112G4double G4USphere::GetStartPhiAngle() const
113{
114 return Base_t::GetStartPhiAngle();
115}
116G4double G4USphere::GetDeltaPhiAngle() const
117{
118 return Base_t::GetDeltaPhiAngle();
119}
120G4double G4USphere::GetStartThetaAngle() const
121{
122 return Base_t::GetStartThetaAngle();
123}
124G4double G4USphere::GetDeltaThetaAngle() const
125{
126 return Base_t::GetDeltaThetaAngle();
127}
128G4double G4USphere::GetSinStartPhi() const
129{
130 return Base_t::GetSinSPhi();
131}
132G4double G4USphere::GetCosStartPhi() const
133{
134 return Base_t::GetCosSPhi();
135}
136G4double G4USphere::GetSinEndPhi() const
137{
138 return Base_t::GetSinEPhi();
139}
140G4double G4USphere::GetCosEndPhi() const
141{
142 return Base_t::GetCosEPhi();
143}
144G4double G4USphere::GetSinStartTheta() const
145{
146 return Base_t::GetSinSTheta();
147}
148G4double G4USphere::GetCosStartTheta() const
149{
150 return Base_t::GetCosSTheta();
151}
152G4double G4USphere::GetSinEndTheta() const
153{
154 return Base_t::GetSinETheta();
155}
156G4double G4USphere::GetCosEndTheta() const
157{
158 return Base_t::GetCosETheta();
159}
160
161void G4USphere::SetInnerRadius(G4double newRMin)
162{
163 Base_t::SetInnerRadius(newRMin);
164 fRebuildPolyhedron = true;
165}
166void G4USphere::SetOuterRadius(G4double newRmax)
167{
168 Base_t::SetOuterRadius(newRmax);
169 fRebuildPolyhedron = true;
170}
171void G4USphere::SetStartPhiAngle(G4double newSphi, G4bool trig)
172{
173 Base_t::SetStartPhiAngle(newSphi, trig);
174 fRebuildPolyhedron = true;
175}
176void G4USphere::SetDeltaPhiAngle(G4double newDphi)
177{
178 Base_t::SetDeltaPhiAngle(newDphi);
179 fRebuildPolyhedron = true;
180}
181void G4USphere::SetStartThetaAngle(G4double newSTheta)
182{
183 Base_t::SetStartThetaAngle(newSTheta);
184 fRebuildPolyhedron = true;
185}
186void G4USphere::SetDeltaThetaAngle(G4double newDTheta)
187{
188 Base_t::SetDeltaThetaAngle(newDTheta);
189 fRebuildPolyhedron = true;
190}
191
193//
194// Dispatch to parameterisation for replication mechanism dimension
195// computation & modification.
196
197void G4USphere::ComputeDimensions( G4VPVParameterisation* p,
198 const G4int n,
199 const G4VPhysicalVolume* pRep)
200{
201 p->ComputeDimensions(*(G4Sphere*)this,n,pRep);
202}
203
205//
206// Make a clone of the object
207
208G4VSolid* G4USphere::Clone() const
209{
210 return new G4USphere(*this);
211}
212
214//
215// Get bounding box
216
217void G4USphere::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
218{
219 static G4bool checkBBox = true;
220
221 G4double rmin = GetInnerRadius();
222 G4double rmax = GetOuterRadius();
223
224 // Find bounding box
225 //
226 if (GetDeltaThetaAngle() >= pi && GetDeltaPhiAngle() >= twopi)
227 {
228 pMin.set(-rmax,-rmax,-rmax);
229 pMax.set( rmax, rmax, rmax);
230 }
231 else
232 {
233 G4double sinStart = GetSinStartTheta();
234 G4double cosStart = GetCosStartTheta();
235 G4double sinEnd = GetSinEndTheta();
236 G4double cosEnd = GetCosEndTheta();
237
238 G4double stheta = GetStartThetaAngle();
239 G4double etheta = stheta + GetDeltaThetaAngle();
240 G4double rhomin = rmin*std::min(sinStart,sinEnd);
241 G4double rhomax = rmax;
242 if (stheta > halfpi) rhomax = rmax*sinStart;
243 if (etheta < halfpi) rhomax = rmax*sinEnd;
244
245 G4TwoVector xymin,xymax;
246 G4GeomTools::DiskExtent(rhomin,rhomax,
247 GetSinStartPhi(),GetCosStartPhi(),
248 GetSinEndPhi(),GetCosEndPhi(),
249 xymin,xymax);
250
251 G4double zmin = std::min(rmin*cosEnd,rmax*cosEnd);
252 G4double zmax = std::max(rmin*cosStart,rmax*cosStart);
253 pMin.set(xymin.x(),xymin.y(),zmin);
254 pMax.set(xymax.x(),xymax.y(),zmax);
255 }
256
257 // Check correctness of the bounding box
258 //
259 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
260 {
261 std::ostringstream message;
262 message << "Bad bounding box (min >= max) for solid: "
263 << GetName() << " !"
264 << "\npMin = " << pMin
265 << "\npMax = " << pMax;
266 G4Exception("G4USphere::BoundingLimits()", "GeomMgt0001",
267 JustWarning, message);
268 StreamInfo(G4cout);
269 }
270
271 // Check consistency of bounding boxes
272 //
273 if (checkBBox)
274 {
275 U3Vector vmin, vmax;
276 Extent(vmin,vmax);
277 if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
278 std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
279 std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
280 std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
281 std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
282 std::abs(pMax.z()-vmax.z()) > kCarTolerance)
283 {
284 std::ostringstream message;
285 message << "Inconsistency in bounding boxes for solid: "
286 << GetName() << " !"
287 << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
288 << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
289 G4Exception("G4USphere::BoundingLimits()", "GeomMgt0001",
290 JustWarning, message);
291 checkBBox = false;
292 }
293 }
294}
295
297//
298// Calculate extent under transform and specified limit
299
300G4bool G4USphere::CalculateExtent(const EAxis pAxis,
301 const G4VoxelLimits& pVoxelLimit,
302 const G4AffineTransform& pTransform,
303 G4double& pMin, G4double& pMax) const
304{
305 G4ThreeVector bmin, bmax;
306
307 // Get bounding box
308 BoundingLimits(bmin,bmax);
309
310 // Find extent
311 G4BoundingEnvelope bbox(bmin,bmax);
312 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
313}
314
316//
317// Create polyhedron for visualization
318
319G4Polyhedron* G4USphere::CreatePolyhedron() const
320{
321 return new G4PolyhedronSphere(GetInnerRadius(),
322 GetOuterRadius(),
323 GetStartPhiAngle(),
324 GetDeltaPhiAngle(),
325 GetStartThetaAngle(),
326 GetDeltaThetaAngle());
327}
328
329#endif // G4GEOM_USE_USOLIDS
const G4double kCarTolerance
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static const G4double pMax
static const G4double pMin
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double pi
Definition: G4SIunits.hh:55
static constexpr double halfpi
Definition: G4SIunits.hh:57
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
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
Definition: G4GeomTools.cc:390
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
EAxis
Definition: geomdefs.hh:54
Definition: DoubConv.h:17
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments