Geant4 11.4.0
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
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// Copy constructor
59
60G4USphere::G4USphere(const G4USphere& rhs)
61 : Base_t(rhs)
62{
63}
64
65//////////////////////////////////////////////////////////////////////////
66//
67// Assignment operator
68
69G4USphere& G4USphere::operator = (const G4USphere& rhs)
70{
71 // Check assignment to self
72 //
73 if (this == &rhs) { return *this; }
74
75 // Copy base class data
76 //
77 Base_t::operator=(rhs);
78
79 return *this;
80}
81
82//////////////////////////////////////////////////////////////////////////
83//
84// Accessors & modifiers
85
86G4double G4USphere::GetInnerRadius() const
87{
88 return Base_t::GetInnerRadius();
89}
90G4double G4USphere::GetOuterRadius() const
91{
92 return Base_t::GetOuterRadius();
93}
94G4double G4USphere::GetStartPhiAngle() const
95{
96 return Base_t::GetStartPhiAngle();
97}
98G4double G4USphere::GetDeltaPhiAngle() const
99{
100 return Base_t::GetDeltaPhiAngle();
101}
102G4double G4USphere::GetStartThetaAngle() const
103{
104 return Base_t::GetStartThetaAngle();
105}
106G4double G4USphere::GetDeltaThetaAngle() const
107{
108 return Base_t::GetDeltaThetaAngle();
109}
110G4double G4USphere::GetSinStartPhi() const
111{
112 return Base_t::GetSinSPhi();
113}
114G4double G4USphere::GetCosStartPhi() const
115{
116 return Base_t::GetCosSPhi();
117}
118G4double G4USphere::GetSinEndPhi() const
119{
120 return Base_t::GetSinEPhi();
121}
122G4double G4USphere::GetCosEndPhi() const
123{
124 return Base_t::GetCosEPhi();
125}
126G4double G4USphere::GetSinStartTheta() const
127{
128 return Base_t::GetSinSTheta();
129}
130G4double G4USphere::GetCosStartTheta() const
131{
132 return Base_t::GetCosSTheta();
133}
134G4double G4USphere::GetSinEndTheta() const
135{
136 return Base_t::GetSinETheta();
137}
138G4double G4USphere::GetCosEndTheta() const
139{
140 return Base_t::GetCosETheta();
141}
142
143void G4USphere::SetInnerRadius(G4double newRMin)
144{
145 Base_t::SetInnerRadius(newRMin);
146 fRebuildPolyhedron = true;
147}
148void G4USphere::SetOuterRadius(G4double newRmax)
149{
150 Base_t::SetOuterRadius(newRmax);
151 fRebuildPolyhedron = true;
152}
153void G4USphere::SetStartPhiAngle(G4double newSphi, G4bool trig)
154{
155 Base_t::SetStartPhiAngle(newSphi, trig);
156 fRebuildPolyhedron = true;
157}
158void G4USphere::SetDeltaPhiAngle(G4double newDphi)
159{
160 Base_t::SetDeltaPhiAngle(newDphi);
161 fRebuildPolyhedron = true;
162}
163void G4USphere::SetStartThetaAngle(G4double newSTheta)
164{
165 Base_t::SetStartThetaAngle(newSTheta);
166 fRebuildPolyhedron = true;
167}
168void G4USphere::SetDeltaThetaAngle(G4double newDTheta)
169{
170 Base_t::SetDeltaThetaAngle(newDTheta);
171 fRebuildPolyhedron = true;
172}
173
174//////////////////////////////////////////////////////////////////////////
175//
176// Dispatch to parameterisation for replication mechanism dimension
177// computation & modification.
178
179void G4USphere::ComputeDimensions( G4VPVParameterisation* p,
180 const G4int n,
181 const G4VPhysicalVolume* pRep)
182{
183 p->ComputeDimensions(*(G4Sphere*)this,n,pRep);
184}
185
186//////////////////////////////////////////////////////////////////////////
187//
188// Make a clone of the object
189
190G4VSolid* G4USphere::Clone() const
191{
192 return new G4USphere(*this);
193}
194
195//////////////////////////////////////////////////////////////////////////
196//
197// Get bounding box
198
199void G4USphere::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
200{
201 static G4bool checkBBox = true;
202
203 G4double rmin = GetInnerRadius();
204 G4double rmax = GetOuterRadius();
205
206 // Find bounding box
207 //
208 if (GetDeltaThetaAngle() >= pi && GetDeltaPhiAngle() >= twopi)
209 {
210 pMin.set(-rmax,-rmax,-rmax);
211 pMax.set( rmax, rmax, rmax);
212 }
213 else
214 {
215 G4double sinStart = GetSinStartTheta();
216 G4double cosStart = GetCosStartTheta();
217 G4double sinEnd = GetSinEndTheta();
218 G4double cosEnd = GetCosEndTheta();
219
220 G4double stheta = GetStartThetaAngle();
221 G4double etheta = stheta + GetDeltaThetaAngle();
222 G4double rhomin = rmin*std::min(sinStart,sinEnd);
223 G4double rhomax = rmax;
224 if (stheta > halfpi) rhomax = rmax*sinStart;
225 if (etheta < halfpi) rhomax = rmax*sinEnd;
226
227 G4TwoVector xymin,xymax;
228 G4GeomTools::DiskExtent(rhomin,rhomax,
229 GetSinStartPhi(),GetCosStartPhi(),
230 GetSinEndPhi(),GetCosEndPhi(),
231 xymin,xymax);
232
233 G4double zmin = std::min(rmin*cosEnd,rmax*cosEnd);
234 G4double zmax = std::max(rmin*cosStart,rmax*cosStart);
235 pMin.set(xymin.x(),xymin.y(),zmin);
236 pMax.set(xymax.x(),xymax.y(),zmax);
237 }
238
239 // Check correctness of the bounding box
240 //
241 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
242 {
243 std::ostringstream message;
244 message << "Bad bounding box (min >= max) for solid: "
245 << GetName() << " !"
246 << "\npMin = " << pMin
247 << "\npMax = " << pMax;
248 G4Exception("G4USphere::BoundingLimits()", "GeomMgt0001",
249 JustWarning, message);
250 StreamInfo(G4cout);
251 }
252
253 // Check consistency of bounding boxes
254 //
255 if (checkBBox)
256 {
257 U3Vector vmin, vmax;
258 Extent(vmin,vmax);
259 if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
260 std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
261 std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
262 std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
263 std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
264 std::abs(pMax.z()-vmax.z()) > kCarTolerance)
265 {
266 std::ostringstream message;
267 message << "Inconsistency in bounding boxes for solid: "
268 << GetName() << " !"
269 << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
270 << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
271 G4Exception("G4USphere::BoundingLimits()", "GeomMgt0001",
272 JustWarning, message);
273 checkBBox = false;
274 }
275 }
276}
277
278//////////////////////////////////////////////////////////////////////////
279//
280// Calculate extent under transform and specified limit
281
282G4bool G4USphere::CalculateExtent(const EAxis pAxis,
283 const G4VoxelLimits& pVoxelLimit,
284 const G4AffineTransform& pTransform,
285 G4double& pMin, G4double& pMax) const
286{
287 G4ThreeVector bmin, bmax;
288
289 // Get bounding box
290 BoundingLimits(bmin,bmax);
291
292 // Find extent
293 G4BoundingEnvelope bbox(bmin,bmax);
294 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
295}
296
297//////////////////////////////////////////////////////////////////////////
298//
299// Create polyhedron for visualization
300
301G4Polyhedron* G4USphere::CreatePolyhedron() const
302{
303 return new G4PolyhedronSphere(GetInnerRadius(),
304 GetOuterRadius(),
305 GetStartPhiAngle(),
306 GetDeltaPhiAngle(),
307 GetStartThetaAngle(),
308 GetDeltaThetaAngle());
309}
310
311#endif // G4GEOM_USE_USOLIDS
const G4double kCarTolerance
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
CLHEP::Hep2Vector G4TwoVector
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)
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4BoundingEnvelope is a helper class to facilitate calculation of the extent of a solid within the li...
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
G4Sphere is, in the general case, a section of a spherical shell, between specified phi and theta ang...
Definition G4Sphere.hh:89
G4VPVParameterisation ia an abstract base class for Parameterisation, able to compute the transformat...
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4VSolid is an abstract base class for solids, physical shapes that can be tracked through....
Definition G4VSolid.hh:80
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
EAxis
Definition geomdefs.hh:54