Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Orb.cc
Go to the documentation of this file.
1// ********************************************************************
2// * License and Disclaimer *
3// * *
4// * The Geant4 software is copyright of the Copyright Holders of *
5// * the Geant4 Collaboration. It is provided under the terms and *
6// * conditions of the Geant4 Software License, included in the file *
7// * LICENSE and available at http://cern.ch/geant4/license . These *
8// * include a list of copyright holders. *
9// * *
10// * Neither the authors of this software system, nor their employing *
11// * institutes,nor the agencies providing financial support for this *
12// * work make any representation or warranty, express or implied, *
13// * regarding this software system or assume any liability for its *
14// * use. Please see the license in the file LICENSE and URL above *
15// * for the full disclaimer and the limitation of liability. *
16// * *
17// * This code implementation is the result of the scientific and *
18// * technical work of the GEANT4 collaboration. *
19// * By using, copying, modifying or distributing the software (or *
20// * any work based on the software) you agree to acknowledge its *
21// * use in resulting scientific publications, and indicate your *
22// * acceptance of all terms of the Geant4 Software license. *
23// ********************************************************************
24//
25// Implementation for G4Orb class
26//
27// 20.08.03 V.Grichine - created
28// 08.08.17 E.Tcherniaev - complete revision, speed-up
29// --------------------------------------------------------------------
30
31#include "G4Orb.hh"
32
33#if !defined(G4GEOM_USE_UORB)
34
35#include "G4TwoVector.hh"
36#include "G4VoxelLimits.hh"
37#include "G4AffineTransform.hh"
38#include "G4BoundingEnvelope.hh"
39#include "G4QuickRand.hh"
40
42
43#include "G4VGraphicsScene.hh"
44#include "G4VisExtent.hh"
45#include "G4AutoLock.hh"
46
47namespace
48{
50}
51
52using namespace CLHEP;
53
54//////////////////////////////////////////////////////////////////////////
55//
56// Constructor
57
58G4Orb::G4Orb( const G4String& pName, G4double pRmax )
59 : G4CSGSolid(pName), fRmax(pRmax)
60{
61 Initialize();
62}
63
64//////////////////////////////////////////////////////////////////////////
65//
66// Fake default constructor - sets only member data and allocates memory
67// for usage restricted to object persistency
68
69G4Orb::G4Orb( __void__& a )
70 : G4CSGSolid(a)
71{
72}
73
74//////////////////////////////////////////////////////////////////////////
75//
76// Assignment operator
77
79{
80 // Check assignment to self
81 //
82 if (this == &rhs) { return *this; }
83
84 // Copy base class data
85 //
87
88 // Copy data
89 //
90 fRmax = rhs.fRmax;
91 halfRmaxTol = rhs.halfRmaxTol;
92 sqrRmaxPlusTol = rhs.sqrRmaxPlusTol;
93 sqrRmaxMinusTol = rhs.sqrRmaxMinusTol;
94
95 return *this;
96}
97
98//////////////////////////////////////////////////////////////////////////
99//
100// Check radius and initialize data members
101
103{
104 const G4double fEpsilon = 2.e-11; // relative tolerance of fRmax
105
106 // Check radius
107 //
108 if ( fRmax < 10*kCarTolerance )
109 {
110 G4Exception("G4Orb::Initialize()", "GeomSolids0002", FatalException,
111 "Invalid radius < 10*kCarTolerance.");
112 }
113 halfRmaxTol = 0.5 * std::max(kCarTolerance, fEpsilon*fRmax);
114 G4double rmaxPlusTol = fRmax + halfRmaxTol;
115 G4double rmaxMinusTol = fRmax - halfRmaxTol;
116 sqrRmaxPlusTol = rmaxPlusTol*rmaxPlusTol;
117 sqrRmaxMinusTol = rmaxMinusTol*rmaxMinusTol;
118}
119
120//////////////////////////////////////////////////////////////////////////
121//
122// Dispatch to parameterisation for replication mechanism dimension
123// computation & modification
124
126 const G4int n,
127 const G4VPhysicalVolume* pRep )
128{
129 p->ComputeDimensions(*this,n,pRep);
130}
131
132//////////////////////////////////////////////////////////////////////////
133//
134// Get bounding box
135
137{
138 G4double radius = GetRadius();
139 pMin.set(-radius,-radius,-radius);
140 pMax.set( radius, radius, radius);
141
142 // Check correctness of the bounding box
143 //
144 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
145 {
146 std::ostringstream message;
147 message << "Bad bounding box (min >= max) for solid: "
148 << GetName() << " !"
149 << "\npMin = " << pMin
150 << "\npMax = " << pMax;
151 G4Exception("G4Orb::BoundingLimits()", "GeomMgt0001",
152 JustWarning, message);
153 DumpInfo();
154 }
155}
156
157//////////////////////////////////////////////////////////////////////////
158//
159// Calculate extent under transform and specified limit
160
162 const G4VoxelLimits& pVoxelLimit,
163 const G4AffineTransform& pTransform,
164 G4double& pMin, G4double& pMax) const
165{
166 G4ThreeVector bmin, bmax;
167 G4bool exist;
168
169 // Get bounding box
170 BoundingLimits(bmin,bmax);
171
172 // Check bounding box
173 G4BoundingEnvelope bbox(bmin,bmax);
174#ifdef G4BBOX_EXTENT
175 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
176#endif
177 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
178 {
179 return exist = pMin < pMax;
180 }
181
182 // Find bounding envelope and calculate extent
183 //
184 static const G4int NTHETA = 8; // number of steps along Theta
185 static const G4int NPHI = 16; // number of steps along Phi
186 static const G4double sinHalfTheta = std::sin(halfpi/NTHETA);
187 static const G4double cosHalfTheta = std::cos(halfpi/NTHETA);
188 static const G4double sinHalfPhi = std::sin(pi/NPHI);
189 static const G4double cosHalfPhi = std::cos(pi/NPHI);
190 static const G4double sinStepTheta = 2.*sinHalfTheta*cosHalfTheta;
191 static const G4double cosStepTheta = 1. - 2.*sinHalfTheta*sinHalfTheta;
192 static const G4double sinStepPhi = 2.*sinHalfPhi*cosHalfPhi;
193 static const G4double cosStepPhi = 1. - 2.*sinHalfPhi*sinHalfPhi;
194
195 G4double radius = GetRadius();
196 G4double rtheta = radius/cosHalfTheta;
197 G4double rphi = rtheta/cosHalfPhi;
198
199 // set reference circle
200 G4TwoVector xy[NPHI];
201 G4double sinCurPhi = sinHalfPhi;
202 G4double cosCurPhi = cosHalfPhi;
203 for (auto & k : xy)
204 {
205 k.set(cosCurPhi,sinCurPhi);
206 G4double sinTmpPhi = sinCurPhi;
207 sinCurPhi = sinCurPhi*cosStepPhi + cosCurPhi*sinStepPhi;
208 cosCurPhi = cosCurPhi*cosStepPhi - sinTmpPhi*sinStepPhi;
209 }
210
211 // set bounding circles
212 G4ThreeVectorList circles[NTHETA];
213 for (auto & circle : circles) { circle.resize(NPHI); }
214
215 G4double sinCurTheta = sinHalfTheta;
216 G4double cosCurTheta = cosHalfTheta;
217 for (auto & circle : circles)
218 {
219 G4double z = rtheta*cosCurTheta;
220 G4double rho = rphi*sinCurTheta;
221 for (G4int k=0; k<NPHI; ++k)
222 {
223 circle[k].set(rho*xy[k].x(),rho*xy[k].y(),z);
224 }
225 G4double sinTmpTheta = sinCurTheta;
226 sinCurTheta = sinCurTheta*cosStepTheta + cosCurTheta*sinStepTheta;
227 cosCurTheta = cosCurTheta*cosStepTheta - sinTmpTheta*sinStepTheta;
228 }
229
230 // set envelope and calculate extent
231 std::vector<const G4ThreeVectorList *> polygons;
232 polygons.resize(NTHETA);
233 for (G4int i=0; i<NTHETA; ++i) { polygons[i] = &circles[i]; }
234
235 G4BoundingEnvelope benv(bmin,bmax,polygons);
236 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
237 return exist;
238}
239
240//////////////////////////////////////////////////////////////////////////
241//
242// Return whether point is inside/outside/on surface
243
245{
246 G4double rr = p.mag2();
247 if (rr > sqrRmaxPlusTol) { return kOutside; }
248 return (rr > sqrRmaxMinusTol) ? kSurface : kInside;
249}
250
251//////////////////////////////////////////////////////////////////////////
252//
253// Return unit normal of surface closest to p
254
256{
257 return (1/p.mag())*p;
258}
259
260//////////////////////////////////////////////////////////////////////////
261//
262// Calculate distance to the surface of the orb from outside
263// - return kInfinity if no intersection or
264// intersection distance <= tolerance
265
267 const G4ThreeVector& v ) const
268{
269 // Check if point is on the surface and traveling away
270 //
271 G4double rr = p.mag2();
272 G4double pv = p.dot(v);
273 if (rr >= sqrRmaxMinusTol && pv >= 0) { return kInfinity; }
274
275 // Find intersection
276 //
277 // Sphere eqn: x^2 + y^2 + z^2 = R^2
278 //
279 // => (px + t*vx)^2 + (py + t*vy)^2 + (pz + t*vz)^2 = R^2
280 // => r^2 + 2t(p.v) + t^2 = R^2
281 // => tmin = -(p.v) - Sqrt((p.v)^2 - (r^2 - R^2))
282 //
283 G4double D = pv*pv - rr + fRmax*fRmax;
284 if (D < 0) { return kInfinity; } // no intersection
285
286 G4double sqrtD = std::sqrt(D);
287 G4double dist = -pv - sqrtD;
288
289 // Avoid rounding errors due to precision issues seen on 64 bits systems.
290 // Split long distances and recompute
291 //
292 G4double Dmax = 32*fRmax;
293 if (dist > Dmax)
294 {
295 dist = dist - 1.e-8*dist - fRmax; // to stay outside after the move
296 dist += DistanceToIn(p + dist*v, v);
297 return (dist >= kInfinity) ? kInfinity : dist;
298 }
299
300 if (sqrtD*2 <= halfRmaxTol) { return kInfinity; } // touch
301
302 return (dist < halfRmaxTol) ? 0. : dist;
303}
304
305//////////////////////////////////////////////////////////////////////////
306//
307// Calculate shortest distance to the boundary from outside
308// - Return 0 if point is inside
309
311{
312 G4double dist = p.mag() - fRmax;
313 return (dist > 0) ? dist : 0.;
314}
315
316//////////////////////////////////////////////////////////////////////////
317//
318// Calculate distance to the surface of the orb from inside and
319// find normal at exit point, if required
320// - when leaving the surface, return 0
321
323 const G4ThreeVector& v,
324 const G4bool calcNorm,
325 G4bool* validNorm,
326 G4ThreeVector* n ) const
327{
328 // Check if point is on the surface and traveling away
329 //
330 G4double rr = p.mag2();
331 G4double pv = p.dot(v);
332 if (rr >= sqrRmaxMinusTol && pv > 0)
333 {
334 if (calcNorm)
335 {
336 *validNorm = true;
337 *n = p*(1./std::sqrt(rr));
338 }
339 return 0.;
340 }
341
342 // Find intersection
343 //
344 // Sphere eqn: x^2 + y^2 + z^2 = R^2
345 //
346 // => (px + t*vx)^2 + (py + t*vy)^2 + (pz + t*vz)^2 = R^2
347 // => r^2 + 2t(p.v) + t^2 = R^2
348 // => tmax = -(p.v) + Sqrt((p.v)^2 - (r^2 - R^2))
349 //
350 G4double D = pv*pv - rr + fRmax*fRmax;
351 G4double tmax = (D <= 0) ? 0. : std::sqrt(D) - pv;
352 if (tmax < halfRmaxTol) { tmax = 0.; }
353 if (calcNorm)
354 {
355 *validNorm = true;
356 G4ThreeVector ptmax = p + tmax*v;
357 *n = ptmax*(1./ptmax.mag());
358 }
359 return tmax;
360}
361
362//////////////////////////////////////////////////////////////////////////
363//
364// Calculate distance (<=actual) to closest surface of shape from inside
365
367{
368#ifdef G4CSGDEBUG
369 if( Inside(p) == kOutside )
370 {
371 std::ostringstream message;
372 G4int oldprc = message.precision(16);
373 message << "Point p is outside (!?) of solid: " << GetName() << "\n";
374 message << "Position:\n";
375 message << " p.x() = " << p.x()/mm << " mm\n";
376 message << " p.y() = " << p.y()/mm << " mm\n";
377 message << " p.z() = " << p.z()/mm << " mm";
378 G4cout.precision(oldprc);
379 G4Exception("G4Trap::DistanceToOut(p)", "GeomSolids1002",
380 JustWarning, message );
381 DumpInfo();
382 }
383#endif
384 G4double dist = fRmax - p.mag();
385 return (dist > 0) ? dist : 0.;
386}
387
388//////////////////////////////////////////////////////////////////////////
389//
390// G4EntityType
391
393{
394 return {"G4Orb"};
395}
396
397//////////////////////////////////////////////////////////////////////////
398//
399// Make a clone of the object
400
402{
403 return new G4Orb(*this);
404}
405
406//////////////////////////////////////////////////////////////////////////
407//
408// Stream object contents to an output stream
409
410std::ostream& G4Orb::StreamInfo( std::ostream& os ) const
411{
412 G4long oldprc = os.precision(16);
413 os << "-----------------------------------------------------------\n"
414 << " *** Dump for solid - " << GetName() << " ***\n"
415 << " ===================================================\n"
416 << " Solid type: G4Orb\n"
417 << " Parameters: \n"
418 << " outer radius: " << fRmax/mm << " mm \n"
419 << "-----------------------------------------------------------\n";
420 os.precision(oldprc);
421 return os;
422}
423
424//////////////////////////////////////////////////////////////////////////
425//
426// Pick random point on the surface
427
429{
430 G4double u, v, b;
431 do
432 {
433 u = 2.*G4QuickRand() - 1.;
434 v = 2.*G4QuickRand() - 1.;
435 b = sqr(u) + sqr(v);
436 } while(b > 1.);
437 G4double a = 2.*std::sqrt(1. - b);
438 return { fRmax*a*u, fRmax*a*v, fRmax*(2.*b - 1.) };
439}
440
441//////////////////////////////////////////////////////////////////////////
442//
443// Computes/returns volume capacity
444
446{
447 if (fCubicVolume == 0)
448 {
449 G4AutoLock l(&orbMutex);
450 fCubicVolume = 4*CLHEP::pi*fRmax*fRmax*fRmax/3.;
451 l.unlock();
452 }
453 return fCubicVolume;
454}
455
456//////////////////////////////////////////////////////////////////////////
457//
458// Computes/returns surface area
459
461{
462 if (fSurfaceArea == 0)
463 {
464 G4AutoLock l(&orbMutex);
465 fSurfaceArea = 4*CLHEP::pi*fRmax*fRmax;
466 l.unlock();
467 }
468 return fSurfaceArea;
469}
470//////////////////////////////////////////////////////////////////////////
471//
472// Methods for visualisation
473
475{
476 scene.AddSolid (*this);
477}
478
480{
481 return {-fRmax, fRmax, -fRmax, fRmax, -fRmax, fRmax};
482}
483
485{
486 return new G4PolyhedronSphere (0., fRmax, 0., 2*pi, 0., pi);
487}
488
489#endif
G4TemplateAutoLock< G4Mutex > G4AutoLock
std::vector< G4ThreeVector > G4ThreeVectorList
G4double D(G4double temp)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double G4QuickRand(uint32_t seed=0)
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
CLHEP::Hep3Vector G4ThreeVector
CLHEP::Hep2Vector G4TwoVector
double G4double
Definition G4Types.hh:83
long G4long
Definition G4Types.hh:87
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4String G4GeometryType
Definition G4VSolid.hh:70
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double mag2() const
double y() const
double dot(const Hep3Vector &) const
double mag() 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...
G4bool BoundingBoxVsVoxelLimits(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4double fSurfaceArea
Definition G4CSGSolid.hh:93
G4double fCubicVolume
Definition G4CSGSolid.hh:92
G4CSGSolid(const G4String &pName)
Definition G4CSGSolid.cc:49
G4CSGSolid & operator=(const G4CSGSolid &rhs)
Definition G4CSGSolid.cc:89
G4Polyhedron * CreatePolyhedron() const override
Definition G4Orb.cc:484
G4double GetSurfaceArea() override
Definition G4Orb.cc:460
void Initialize()
Definition G4Orb.cc:102
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
Definition G4Orb.cc:125
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
Definition G4Orb.cc:266
EInside Inside(const G4ThreeVector &p) const override
Definition G4Orb.cc:244
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const override
Definition G4Orb.cc:161
G4double GetCubicVolume() override
Definition G4Orb.cc:445
G4double GetRadius() const
std::ostream & StreamInfo(std::ostream &os) const override
Definition G4Orb.cc:410
G4Orb(const G4String &pName, G4double pRmax)
Definition G4Orb.cc:58
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
Definition G4Orb.cc:322
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
Definition G4Orb.cc:255
G4ThreeVector GetPointOnSurface() const override
Definition G4Orb.cc:428
G4GeometryType GetEntityType() const override
Definition G4Orb.cc:392
G4VisExtent GetExtent() const override
Definition G4Orb.cc:479
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
Definition G4Orb.cc:474
G4VSolid * Clone() const override
Definition G4Orb.cc:401
G4Orb & operator=(const G4Orb &rhs)
Definition G4Orb.cc:78
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
Definition G4Orb.cc:136
virtual void AddSolid(const G4Box &)=0
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....
G4String GetName() const
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59
void DumpInfo() const
G4double kCarTolerance
Definition G4VSolid.hh:418
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
EAxis
Definition geomdefs.hh:54
EInside
Definition geomdefs.hh:67
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68
@ kSurface
Definition geomdefs.hh:69
T sqr(const T &x)
Definition templates.hh:128