Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Ellipsoid Class Reference

G4Ellipsoid is an ellipsoidal solid, optionally cut at a given Z. More...

#include <G4Ellipsoid.hh>

Inheritance diagram for G4Ellipsoid:

Public Member Functions

 G4Ellipsoid (const G4String &name, G4double xSemiAxis, G4double ySemiAxis, G4double zSemiAxis, G4double zBottomCut=0., G4double zTopCut=0.)
 ~G4Ellipsoid () override
G4double GetDx () const
G4double GetDy () const
G4double GetDz () const
G4double GetSemiAxisMax (G4int i) const
G4double GetZBottomCut () const
G4double GetZTopCut () const
void SetSemiAxis (G4double x, G4double y, G4double z)
void SetZCuts (G4double newzBottomCut, G4double newzTopCut)
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const override
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const override
EInside Inside (const G4ThreeVector &p) const override
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const override
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const override
G4double DistanceToIn (const G4ThreeVector &p) const override
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
G4double DistanceToOut (const G4ThreeVector &p) const override
G4GeometryType GetEntityType () const override
G4VSolidClone () const override
std::ostream & StreamInfo (std::ostream &os) const override
G4double GetCubicVolume () override
G4double GetSurfaceArea () override
G4ThreeVector GetPointOnSurface () const override
void DescribeYourselfTo (G4VGraphicsScene &scene) const override
G4VisExtent GetExtent () const override
G4PolyhedronCreatePolyhedron () const override
G4PolyhedronGetPolyhedron () const override
 G4Ellipsoid (__void__ &)
 G4Ellipsoid (const G4Ellipsoid &rhs)
G4Ellipsoidoperator= (const G4Ellipsoid &rhs)
Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
virtual ~G4VSolid ()
 G4VSolid (const G4VSolid &rhs)
G4VSolidoperator= (const G4VSolid &rhs)
G4bool operator== (const G4VSolid &s) const
G4String GetName () const
void SetName (const G4String &name)
G4double GetTolerance () const
virtual G4int GetNumOfConstituents () const
virtual G4bool IsFaceted () const
void DumpInfo () const
virtual const G4VSolidGetConstituentSolid (G4int no) const
virtual G4VSolidGetConstituentSolid (G4int no)
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 G4VSolid (__void__ &)
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
G4double EstimateSurfaceArea (G4int nStat, G4double epsilon) const

Additional Inherited Members

Protected Member Functions inherited from G4VSolid
void CalculateClippedPolygonExtent (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
void ClipCrossSection (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
void ClipBetweenSections (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
void ClipPolygon (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
Protected Attributes inherited from G4VSolid
G4double kCarTolerance

Detailed Description

G4Ellipsoid is an ellipsoidal solid, optionally cut at a given Z.

Definition at line 66 of file G4Ellipsoid.hh.

Constructor & Destructor Documentation

◆ G4Ellipsoid() [1/3]

G4Ellipsoid::G4Ellipsoid ( const G4String & name,
G4double xSemiAxis,
G4double ySemiAxis,
G4double zSemiAxis,
G4double zBottomCut = 0.,
G4double zTopCut = 0. )

Constructs an ellipsoid, given its input parameters.

Parameters
[in]nameThe solid name.
[in]xSemiAxisSemiaxis in X.
[in]ySemiAxisSemiaxis in Y.
[in]zSemiAxisSemiaxis in Z.
[in]zBottomCutOptional lower cut plane level in Z.
[in]zTopCutOptional upper cut plane level in Z.

Definition at line 67 of file G4Ellipsoid.cc.

73 : G4VSolid(name), fDx(xSemiAxis), fDy(ySemiAxis), fDz(zSemiAxis),
74 fZBottomCut(zBottomCut), fZTopCut(zTopCut)
75{
76 CheckParameters();
77}
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59

Referenced by Clone(), G4Ellipsoid(), operator=(), and SetZCuts().

◆ ~G4Ellipsoid()

G4Ellipsoid::~G4Ellipsoid ( )
override

Destructor.

Definition at line 93 of file G4Ellipsoid.cc.

94{
95 delete fpPolyhedron; fpPolyhedron = nullptr;
96}

◆ G4Ellipsoid() [2/3]

G4Ellipsoid::G4Ellipsoid ( __void__ & a)

Fake default constructor for usage restricted to direct object persistency for clients requiring preallocation of memory for persistifiable objects.

Definition at line 84 of file G4Ellipsoid.cc.

85 : G4VSolid(a), fDx(0.), fDy(0.), fDz(0.), fZBottomCut(0.), fZTopCut(0.)
86{
87}

◆ G4Ellipsoid() [3/3]

G4Ellipsoid::G4Ellipsoid ( const G4Ellipsoid & rhs)

Copy constructor and assignment operator.

Definition at line 102 of file G4Ellipsoid.cc.

103 : G4VSolid(rhs),
104 fDx(rhs.fDx), fDy(rhs.fDy), fDz(rhs.fDz),
105 fZBottomCut(rhs.fZBottomCut), fZTopCut(rhs.fZTopCut),
106 halfTolerance(rhs.halfTolerance),
107 fXmax(rhs.fXmax), fYmax(rhs.fYmax),
108 fRsph(rhs.fRsph), fR(rhs.fR),
109 fSx(rhs.fSx), fSy(rhs.fSy), fSz(rhs.fSz),
110 fZMidCut(rhs.fZMidCut), fZDimCut(rhs.fZDimCut),
111 fQ1(rhs.fQ1), fQ2(rhs.fQ2),
112 fCubicVolume(rhs.fCubicVolume),
113 fSurfaceArea(rhs.fSurfaceArea),
114 fLateralArea(rhs.fLateralArea)
115{
116}

Member Function Documentation

◆ BoundingLimits()

void G4Ellipsoid::BoundingLimits ( G4ThreeVector & pMin,
G4ThreeVector & pMax ) const
overridevirtual

Computes the bounding limits of the solid.

Parameters
[out]pMinThe minimum bounding limit point.
[out]pMaxThe maximum bounding limit point.

Reimplemented from G4VSolid.

Definition at line 264 of file G4Ellipsoid.cc.

266{
267 pMin.set(-fXmax,-fYmax, fZBottomCut);
268 pMax.set( fXmax, fYmax, fZTopCut);
269}
void set(double x, double y, double z)

Referenced by CalculateExtent().

◆ CalculateExtent()

G4bool G4Ellipsoid::CalculateExtent ( const EAxis pAxis,
const G4VoxelLimits & pVoxelLimit,
const G4AffineTransform & pTransform,
G4double & pmin,
G4double & pmax ) const
overridevirtual

Calculates the minimum and maximum extent of the solid, when under the specified transform, and within the specified limits.

Parameters
[in]pAxisThe axis along which compute the extent.
[in]pVoxelLimitThe limiting space dictated by voxels.
[in]pTransformThe internal transformation applied to the solid.
[out]pMinThe minimum extent value.
[out]pMaxThe maximum extent value.
Returns
True if the solid is intersected by the extent region.

Implements G4VSolid.

Definition at line 276 of file G4Ellipsoid.cc.

280{
281 G4ThreeVector bmin, bmax;
282
283 // Get bounding box
284 BoundingLimits(bmin,bmax);
285
286 // Find extent
287 G4BoundingEnvelope bbox(bmin,bmax);
288 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
289}
CLHEP::Hep3Vector G4ThreeVector
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override

◆ Clone()

G4VSolid * G4Ellipsoid::Clone ( ) const
overridevirtual

Makes a clone of the object for use in multi-treading.

Returns
A pointer to the new cloned allocated solid.

Reimplemented from G4VSolid.

Definition at line 659 of file G4Ellipsoid.cc.

660{
661 return new G4Ellipsoid(*this);
662}
G4Ellipsoid(const G4String &name, G4double xSemiAxis, G4double ySemiAxis, G4double zSemiAxis, G4double zBottomCut=0., G4double zTopCut=0.)

◆ ComputeDimensions()

void G4Ellipsoid::ComputeDimensions ( G4VPVParameterisation * p,
const G4int n,
const G4VPhysicalVolume * pRep )
overridevirtual

Dispatch method for parameterisation replication mechanism and dimension computation.

Reimplemented from G4VSolid.

Definition at line 253 of file G4Ellipsoid.cc.

256{
257 p->ComputeDimensions(*this,n,pRep);
258}
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const

◆ CreatePolyhedron()

G4Polyhedron * G4Ellipsoid::CreatePolyhedron ( ) const
overridevirtual

Creates a Polyhedron used for Visualisation. It is the caller's responsibility to delete it. A null pointer means "not created".

Reimplemented from G4VSolid.

Definition at line 924 of file G4Ellipsoid.cc.

925{
926 return new G4PolyhedronEllipsoid(fDx, fDy, fDz, fZBottomCut, fZTopCut);
927}

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

void G4Ellipsoid::DescribeYourselfTo ( G4VGraphicsScene & scene) const
overridevirtual

Methods for creating graphical representations (i.e. for visualisation).

Implements G4VSolid.

Definition at line 906 of file G4Ellipsoid.cc.

907{
908 scene.AddSolid(*this);
909}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

G4double G4Ellipsoid::DistanceToIn ( const G4ThreeVector & p) const
overridevirtual

Calculates the distance to the nearest surface of a shape from an outside point. The distance can be an underestimate.

Parameters
[in]pThe point at offset p.
Returns
The safety distance to enter the shape.

Implements G4VSolid.

Definition at line 476 of file G4Ellipsoid.cc.

477{
478 G4double px = p.x();
479 G4double py = p.y();
480 G4double pz = p.z();
481
482 // Safety distance to bounding box
483 G4double distX = std::abs(px) - fXmax;
484 G4double distY = std::abs(py) - fYmax;
485 G4double distZ = std::max(pz - fZTopCut, fZBottomCut - pz);
486 G4double distB = std::max(std::max(distX, distY), distZ);
487
488 // Safety distance to lateral surface
489 G4double x = px * fSx;
490 G4double y = py * fSy;
491 G4double z = pz * fSz;
492 G4double distR = std::sqrt(x*x + y*y + z*z) - fR;
493
494 // Return safety to in
495 G4double dist = std::max(distB, distR);
496 return (dist < 0.) ? 0. : dist;
497}
double G4double
Definition G4Types.hh:83
double z() const
double x() const
double y() const

◆ DistanceToIn() [2/2]

G4double G4Ellipsoid::DistanceToIn ( const G4ThreeVector & p,
const G4ThreeVector & v ) const
overridevirtual

Returns the distance along the normalised vector 'v' to the shape, from the point at offset 'p'. If there is no intersection, returns kInfinity. The first intersection resulting from 'leaving' a surface/volume is discarded. Hence, it is tolerant of points on the surface of the shape.

Parameters
[in]pThe point at offset p.
[in]vThe normalised direction vector.
Returns
The distance to enter the shape.

Implements G4VSolid.

Definition at line 389 of file G4Ellipsoid.cc.

391{
392 G4double offset = 0.;
393 G4ThreeVector pcur = p;
394
395 // Check if point is flying away, relative to bounding box
396 //
397 G4double safex = std::abs(p.x()) - fXmax;
398 G4double safey = std::abs(p.y()) - fYmax;
399 G4double safet = p.z() - fZTopCut;
400 G4double safeb = fZBottomCut - p.z();
401
402 if (safex >= -halfTolerance && p.x() * v.x() >= 0.) { return kInfinity; }
403 if (safey >= -halfTolerance && p.y() * v.y() >= 0.) { return kInfinity; }
404 if (safet >= -halfTolerance && v.z() >= 0.) { return kInfinity; }
405 if (safeb >= -halfTolerance && v.z() <= 0.) { return kInfinity; }
406
407 // Relocate point, if required
408 //
409 G4double safe = std::max(std::max(std::max(safex, safey), safet), safeb);
410 if (safe > 32. * fRsph)
411 {
412 offset = (1. - 1.e-08) * safe - 2. * fRsph;
413 pcur += offset * v;
414 G4double dist = DistanceToIn(pcur, v);
415 return (dist == kInfinity) ? kInfinity : dist + offset;
416 }
417
418 // Scale ellipsoid to sphere
419 //
420 G4double px = pcur.x() * fSx;
421 G4double py = pcur.y() * fSy;
422 G4double pz = pcur.z() * fSz;
423 G4double vx = v.x() * fSx;
424 G4double vy = v.y() * fSy;
425 G4double vz = v.z() * fSz;
426
427 // Check if point is leaving the solid
428 //
429 G4double dzcut = fZDimCut;
430 G4double pzcut = pz - fZMidCut;
431 G4double distZ = std::abs(pzcut) - dzcut;
432 if (distZ >= -halfTolerance && pzcut * vz >= 0.) { return kInfinity; }
433
434 G4double rr = px * px + py * py + pz * pz;
435 G4double pv = px * vx + py * vy + pz * vz;
436 G4double distR = fQ1 * rr - fQ2;
437 if (distR >= -halfTolerance && pv >= 0.) { return kInfinity; }
438
439 G4double A = vx * vx + vy * vy + vz * vz;
440 G4double B = pv;
441 G4double C = rr - fR * fR;
442 G4double D = B * B - A * C;
443 // scratch^2 = R^2 - (R - halfTolerance)^2 = 2 * R * halfTolerance
444 G4double EPS = A * A * fR * kCarTolerance; // discriminant at scratching
445 if (D <= EPS) { return kInfinity; } // no intersection or scratching
446
447 // Find intersection with Z planes
448 //
449 G4double invz = (vz == 0) ? DBL_MAX : -1./vz;
450 G4double dz = std::copysign(dzcut, invz);
451 G4double tzmin = (pzcut - dz) * invz;
452 G4double tzmax = (pzcut + dz) * invz;
453
454 // Find intersection with lateral surface
455 //
456 G4double tmp = -B - std::copysign(std::sqrt(D), B);
457 G4double t1 = tmp / A;
458 G4double t2 = C / tmp;
459 G4double trmin = std::min(t1, t2);
460 G4double trmax = std::max(t1, t2);
461
462 // Return distance
463 //
464 G4double tmin = std::max(tzmin, trmin);
465 G4double tmax = std::min(tzmax, trmax);
466
467 if (tmax - tmin <= halfTolerance) { return kInfinity; } // touch or no hit
468
469 return (tmin < halfTolerance) ? offset : tmin + offset;
470}
G4double C(G4double temp)
G4double B(G4double temperature)
G4double D(G4double temp)
G4ThreadLocal T * G4GeomSplitter< T >::offset
const G4double A[17]
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
G4double kCarTolerance
Definition G4VSolid.hh:418
#define DBL_MAX
Definition templates.hh:62

Referenced by DistanceToIn().

◆ DistanceToOut() [1/2]

G4double G4Ellipsoid::DistanceToOut ( const G4ThreeVector & p) const
overridevirtual

Calculates the distance to the nearest surface of a shape from an inside point 'p'. The distance can be an underestimate.

Parameters
[in]pThe point at offset p.
Returns
The safety distance to exit the shape.

Implements G4VSolid.

Definition at line 630 of file G4Ellipsoid.cc.

631{
632 // Safety distance in z direction
633 G4double distZ = std::min(fZTopCut - p.z(), p.z() - fZBottomCut);
634
635 // Safety distance to lateral surface
636 G4double x = p.x() * fSx;
637 G4double y = p.y() * fSy;
638 G4double z = p.z() * fSz;
639 G4double distR = fR - std::sqrt(x*x + y*y + z*z);
640
641 // Return safety to out
642 G4double dist = std::min(distZ, distR);
643 return (dist < 0.) ? 0. : dist;
644}

◆ DistanceToOut() [2/2]

G4double G4Ellipsoid::DistanceToOut ( const G4ThreeVector & p,
const G4ThreeVector & v,
const G4bool calcNorm = false,
G4bool * validNorm = nullptr,
G4ThreeVector * n = nullptr ) const
overridevirtual

Returns the distance along the normalised vector 'v' to the shape, from a point at an offset 'p' inside or on the surface of the shape. Intersections with surfaces, when the point is less than Tolerance/2 from a surface must be ignored.

Parameters
[in]pThe point at offset p.
[in]vThe normalised direction vector.
[in]calcNormFlag to indicate if to calculate the normal or not.
[out]validNormFlag set to true if the solid lies entirely behind or on the exiting surface. It is set false if the solid does not lie entirely behind or on the exiting surface. 'calcNorm' must be true, otherwise it is unused.
[out]nThe exiting outwards normal vector (undefined Magnitude). 'calcNorm' must be true, otherwise it is unused.
Returns
The distance to exit the shape.

Implements G4VSolid.

Definition at line 503 of file G4Ellipsoid.cc.

508{
509 // Check if point flying away relative to Z planes
510 //
511 G4double pz = p.z() * fSz;
512 G4double vz = v.z() * fSz;
513 G4double dzcut = fZDimCut;
514 G4double pzcut = pz - fZMidCut;
515 G4double distZ = std::abs(pzcut) - dzcut;
516 if (distZ >= -halfTolerance && pzcut * vz > 0.)
517 {
518 if (calcNorm)
519 {
520 *validNorm = true;
521 n->set(0., 0., std::copysign(1., pzcut));
522 }
523 return 0.;
524 }
525
526 // Check if point is flying away relative to lateral surface
527 //
528 G4double px = p.x() * fSx;
529 G4double py = p.y() * fSy;
530 G4double vx = v.x() * fSx;
531 G4double vy = v.y() * fSy;
532 G4double rr = px * px + py * py + pz * pz;
533 G4double pv = px * vx + py * vy + pz * vz;
534 G4double distR = fQ1 * rr - fQ2;
535 if (distR >= -halfTolerance && pv > 0.)
536 {
537 if (calcNorm)
538 {
539 *validNorm = true;
540 *n = G4ThreeVector(px*fSx, py*fSy, pz*fSz).unit();
541 }
542 return 0.;
543 }
544
545 // Just in case check if point is outside (normally it should never be)
546 //
547 if (std::max(distZ, distR) > halfTolerance)
548 {
549#ifdef G4SPECSDEBUG
550 std::ostringstream message;
551 G4long oldprc = message.precision(16);
552 message << "Point p is outside (!?) of solid: "
553 << GetName() << G4endl;
554 message << "Position: " << p << G4endl;;
555 message << "Direction: " << v;
556 G4cout.precision(oldprc);
557 G4Exception("G4Ellipsoid::DistanceToOut(p,v)", "GeomSolids1002",
558 JustWarning, message );
559 DumpInfo();
560#endif
561 if (calcNorm)
562 {
563 *validNorm = true;
564 *n = ApproxSurfaceNormal(p);
565 }
566 return 0.;
567 }
568
569 // Set coefficients of quadratic equation: A t^2 + 2B t + C = 0
570 //
571 G4double A = vx * vx + vy * vy + vz * vz;
572 G4double B = pv;
573 G4double C = rr - fR * fR;
574 G4double D = B * B - A * C;
575 // It is expected that the point is located inside the sphere, so
576 // max term in the expression for discriminant is A * R^2 and
577 // max calculation error can be derived as follows:
578 // A * (1 + 2e) * R^2 * (1 + 2e) = A * R^2 + (4 * A * R^2 * e)
579 G4double EPS = 4. * A * fR * fR * DBL_EPSILON; // calculation error
580
581 if (D <= EPS) // no intersection
582 {
583 if (calcNorm)
584 {
585 *validNorm = true;
586 *n = G4ThreeVector(px*fSx, py*fSy, pz*fSz).unit();
587 }
588 return 0.;
589 }
590
591 // Find intersection with Z cuts
592 //
593 G4double tzmax = (vz == 0.)
594 ? DBL_MAX
595 : (std::copysign(dzcut, vz) - pzcut) / vz;
596
597 // Find intersection with lateral surface
598 //
599 G4double tmp = -B - std::copysign(std::sqrt(D), B);
600 G4double trmax = (tmp < 0.) ? C/tmp : tmp/A;
601
602 // Find distance and set normal, if required
603 //
604 G4double tmax = std::min(tzmax, trmax);
605 //if (tmax < halfTolerance) tmax = 0.;
606
607 if (calcNorm)
608 {
609 *validNorm = true;
610 if (tmax == tzmax)
611 {
612 G4double pznew = pz + tmax * vz;
613 n->set(0., 0., (pznew > fZMidCut) ? 1. : -1.);
614 }
615 else
616 {
617 G4double nx = (px + tmax * vx) * fSx;
618 G4double ny = (py + tmax * vy) * fSy;
619 G4double nz = (pz + tmax * vz) * fSz;
620 *n = G4ThreeVector(nx, ny, nz).unit();
621 }
622 }
623 return tmax;
624}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
long G4long
Definition G4Types.hh:87
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
G4String GetName() const
void DumpInfo() const
#define DBL_EPSILON
Definition templates.hh:66

◆ GetCubicVolume()

G4double G4Ellipsoid::GetCubicVolume ( )
overridevirtual

Returning an estimation of the solid volume (capacity) and surface area, in internal units.

Reimplemented from G4VSolid.

Definition at line 776 of file G4Ellipsoid.cc.

777{
778 if (fCubicVolume == 0)
779 {
780 G4AutoLock l(&ellipseMutex);
781 G4double piAB_3 = CLHEP::pi * fDx * fDy / 3.;
782 fCubicVolume = 4. * piAB_3 * fDz;
783 if (fZBottomCut > -fDz)
784 {
785 G4double hbot = 1. + fZBottomCut / fDz;
786 fCubicVolume -= piAB_3 * hbot * hbot * (2. * fDz - fZBottomCut);
787 }
788 if (fZTopCut < fDz)
789 {
790 G4double htop = 1. - fZTopCut / fDz;
791 fCubicVolume -= piAB_3 * htop * htop * (2. * fDz + fZTopCut);
792 }
793 l.unlock();
794 }
795 return fCubicVolume;
796}
G4TemplateAutoLock< G4Mutex > G4AutoLock

◆ GetDx()

G4double G4Ellipsoid::GetDx ( ) const
inline

Accessors.

Referenced by GetPointOnSurface(), and StreamInfo().

◆ GetDy()

G4double G4Ellipsoid::GetDy ( ) const
inline

Referenced by GetPointOnSurface(), and StreamInfo().

◆ GetDz()

G4double G4Ellipsoid::GetDz ( ) const
inline

Referenced by GetPointOnSurface(), and StreamInfo().

◆ GetEntityType()

G4GeometryType G4Ellipsoid::GetEntityType ( ) const
overridevirtual

Returns the type ID, "G4Ellipsoid" of the solid.

Implements G4VSolid.

Definition at line 650 of file G4Ellipsoid.cc.

651{
652 return {"G4Ellipsoid"};
653}

Referenced by StreamInfo().

◆ GetExtent()

G4VisExtent G4Ellipsoid::GetExtent ( ) const
overridevirtual

Provides extent (bounding box) as possible hint to the graphics view.

Reimplemented from G4VSolid.

Definition at line 915 of file G4Ellipsoid.cc.

916{
917 return { -fXmax, fXmax, -fYmax, fYmax, fZBottomCut, fZTopCut };
918}

◆ GetPointOnSurface()

G4ThreeVector G4Ellipsoid::GetPointOnSurface ( ) const
overridevirtual

Returns a random point located and uniformly distributed on the surface of the solid.

Reimplemented from G4VSolid.

Definition at line 828 of file G4Ellipsoid.cc.

829{
830 G4double A = GetDx();
831 G4double B = GetDy();
832 G4double C = GetDz();
833 G4double Zbot = GetZBottomCut();
834 G4double Ztop = GetZTopCut();
835
836 // Calculate cut areas
837 G4double invC = 1. / C;
838 G4double Hbot = 1. + Zbot * invC;
839 G4double Htop = 1. - Ztop * invC;
840 G4double piAB = CLHEP::pi * A * B;
841 G4double Sbot = piAB * Hbot * (2. - Hbot);
842 G4double Stop = piAB * Htop * (2. - Htop);
843
844 // Get area of lateral surface
845 if (fLateralArea == 0.)
846 {
847 G4AutoLock l(&lateralareaMutex);
848 fLateralArea = LateralSurfaceArea();
849 l.unlock();
850 }
851 G4double Slat = fLateralArea;
852
853 // Select surface (0 - bottom cut, 1 - lateral surface, 2 - top cut)
854 G4double select = (Sbot + Slat + Stop) * G4QuickRand();
855 G4int k = 0;
856 k += (G4int)(select > Sbot);
857 k += (G4int)(select > Sbot + Slat);
858
859 // Pick random point on selected surface
860 G4double phi = CLHEP::twopi * G4QuickRand();
861 G4double cosphi = std::cos(phi);
862 G4double sinphi = std::sin(phi);
864 switch (k)
865 {
866 case 0: // bootom z-cut, ellipse
867 {
868 G4double scale = std::sqrt(Hbot * (2. - Hbot));
869 G4double rho = scale*std::sqrt(G4QuickRand());
870 p.set(A * rho * cosphi, B * rho * sinphi, Zbot);
871 break;
872 }
873 case 1: // lateral surface (rejection sampling)
874 {
875 G4double x, y, z;
876 G4double s_max = std::max(std::max(A * B, A * C), B * C);
877 for (G4int i = 0; i < 10000; ++i)
878 {
879 // generate random point on unit sphere
880 z = (Zbot + (Ztop - Zbot) * G4QuickRand()) * invC;
881 G4double rho = std::sqrt((1. + z) * (1. - z));
882 x = rho * cosphi;
883 y = rho * sinphi;
884 // check acceptance
885 G4double ss = sqr(B * C * x) + sqr(A * C * y) + sqr(A * B * z);
886 if (sqr(s_max * G4QuickRand()) <= ss) { break; }
887 }
888 p.set(A * x, B * y, C * z);
889 break;
890 }
891 case 2: // top z-cut, ellipse
892 {
893 G4double scale = std::sqrt(Htop * (2. - Htop));
894 G4double rho = scale*std::sqrt(G4QuickRand());
895 p.set(A * rho * cosphi, B * rho * sinphi, Ztop);
896 break;
897 }
898 }
899 return p;
900}
G4double G4QuickRand(uint32_t seed=0)
int G4int
Definition G4Types.hh:85
G4double GetDx() const
G4double GetDy() const
G4double GetZTopCut() const
G4double GetZBottomCut() const
G4double GetDz() const
T sqr(const T &x)
Definition templates.hh:128

◆ GetPolyhedron()

G4Polyhedron * G4Ellipsoid::GetPolyhedron ( ) const
overridevirtual

Smart access function - creates on request and stores for future access. A null pointer means "not available".

Reimplemented from G4VSolid.

Definition at line 933 of file G4Ellipsoid.cc.

934{
935 if (fpPolyhedron == nullptr ||
936 fRebuildPolyhedron ||
937 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
938 fpPolyhedron->GetNumberOfRotationSteps())
939 {
940 G4AutoLock l(&polyhedronMutex);
941 delete fpPolyhedron;
942 fpPolyhedron = CreatePolyhedron();
943 fRebuildPolyhedron = false;
944 l.unlock();
945 }
946 return fpPolyhedron;
947}
G4Polyhedron * CreatePolyhedron() const override

◆ GetSemiAxisMax()

◆ GetSurfaceArea()

G4double G4Ellipsoid::GetSurfaceArea ( )
overridevirtual

Returns an estimation of the solid surface area in internal units. This method may be overloaded by derived classes to compute the exact geometrical quantity for solids where this is possible, or anyway to cache the computed value. Note: the computed value is NOT cached.

Reimplemented from G4VSolid.

Definition at line 802 of file G4Ellipsoid.cc.

803{
804 if (fSurfaceArea == 0.)
805 {
806 G4AutoLock l(&ellipseMutex);
807 G4double piAB = CLHEP::pi * fDx * fDy;
808 fSurfaceArea = LateralSurfaceArea();
809 if (fZBottomCut > -fDz)
810 {
811 G4double hbot = 1. + fZBottomCut / fDz;
812 fSurfaceArea += piAB * hbot * (2. - hbot);
813 }
814 if (fZTopCut < fDz)
815 {
816 G4double htop = 1. - fZTopCut / fDz;
817 fSurfaceArea += piAB * htop * (2. - htop);
818 }
819 l.unlock();
820 }
821 return fSurfaceArea;
822}

◆ GetZBottomCut()

◆ GetZTopCut()

◆ Inside()

EInside G4Ellipsoid::Inside ( const G4ThreeVector & p) const
overridevirtual

Concrete implementations of the expected query interfaces for solids, as defined in the base class G4VSolid.

Implements G4VSolid.

Definition at line 295 of file G4Ellipsoid.cc.

296{
297 G4double x = p.x() * fSx;
298 G4double y = p.y() * fSy;
299 G4double z = p.z() * fSz;
300 G4double rr = x * x + y * y + z * z;
301 G4double distZ = std::abs(z - fZMidCut) - fZDimCut;
302 G4double distR = fQ1 * rr - fQ2;
303 G4double dist = std::max(distZ, distR);
304
305 if (dist > halfTolerance) { return kOutside; }
306 return (dist > -halfTolerance) ? kSurface : kInside;
307}
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68
@ kSurface
Definition geomdefs.hh:69

◆ operator=()

G4Ellipsoid & G4Ellipsoid::operator= ( const G4Ellipsoid & rhs)

Definition at line 122 of file G4Ellipsoid.cc.

123{
124 // Check assignment to self
125 //
126 if (this == &rhs) { return *this; }
127
128 // Copy base class data
129 //
131
132 // Copy data
133 //
134 fDx = rhs.fDx;
135 fDy = rhs.fDy;
136 fDz = rhs.fDz;
137 fZBottomCut = rhs.fZBottomCut;
138 fZTopCut = rhs.fZTopCut;
139
140 halfTolerance = rhs.halfTolerance;
141 fXmax = rhs.fXmax;
142 fYmax = rhs.fYmax;
143 fRsph = rhs.fRsph;
144 fR = rhs.fR;
145 fSx = rhs.fSx;
146 fSy = rhs.fSy;
147 fSz = rhs.fSz;
148 fZMidCut = rhs.fZMidCut;
149 fZDimCut = rhs.fZDimCut;
150 fQ1 = rhs.fQ1;
151 fQ2 = rhs.fQ2;
152
153 fCubicVolume = rhs.fCubicVolume;
154 fSurfaceArea = rhs.fSurfaceArea;
155 fLateralArea = rhs.fLateralArea;
156
157 fRebuildPolyhedron = false;
158 delete fpPolyhedron; fpPolyhedron = nullptr;
159
160 return *this;
161}
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:108

◆ SetSemiAxis()

void G4Ellipsoid::SetSemiAxis ( G4double x,
G4double y,
G4double z )
inline

Modifiers.

◆ SetZCuts()

void G4Ellipsoid::SetZCuts ( G4double newzBottomCut,
G4double newzTopCut )
inline

◆ StreamInfo()

std::ostream & G4Ellipsoid::StreamInfo ( std::ostream & os) const
overridevirtual

Streams the object contents to an output stream.

Implements G4VSolid.

Definition at line 668 of file G4Ellipsoid.cc.

669{
670 G4long oldprc = os.precision(16);
671 os << "-----------------------------------------------------------\n"
672 << " *** Dump for solid - " << GetName() << " ***\n"
673 << " ===================================================\n"
674 << " Solid type: " << GetEntityType() << "\n"
675 << " Parameters: \n"
676 << " semi-axis x: " << GetDx()/mm << " mm \n"
677 << " semi-axis y: " << GetDy()/mm << " mm \n"
678 << " semi-axis z: " << GetDz()/mm << " mm \n"
679 << " lower cut in z: " << GetZBottomCut()/mm << " mm \n"
680 << " upper cut in z: " << GetZTopCut()/mm << " mm \n"
681 << "-----------------------------------------------------------\n";
682 os.precision(oldprc);
683 return os;
684}
G4GeometryType GetEntityType() const override

◆ SurfaceNormal()

G4ThreeVector G4Ellipsoid::SurfaceNormal ( const G4ThreeVector & p) const
overridevirtual

Returns the outwards pointing unit normal of the shape for the surface closest to the point at offset 'p'.

Parameters
[in]pThe point at offset p.
Returns
The outwards pointing unit normal.

Implements G4VSolid.

Definition at line 313 of file G4Ellipsoid.cc.

314{
315 G4ThreeVector norm(0., 0., 0.);
316 G4int nsurf = 0;
317
318 // Check cuts
319 G4double x = p.x() * fSx;
320 G4double y = p.y() * fSy;
321 G4double z = p.z() * fSz;
322 G4double distZ = std::abs(z - fZMidCut) - fZDimCut;
323 if (std::abs(distZ) <= halfTolerance)
324 {
325 norm.setZ(std::copysign(1., z - fZMidCut));
326 ++nsurf;
327 }
328
329 // Check lateral surface
330 G4double distR = fQ1*(x*x + y*y + z*z) - fQ2;
331 if (std::abs(distR) <= halfTolerance)
332 {
333 // normal = (p.x/a^2, p.y/b^2, p.z/c^2)
334 norm += G4ThreeVector(x*fSx, y*fSy, z*fSz).unit();
335 ++nsurf;
336 }
337
338 // Return normal
339 if (nsurf == 1)
340 {
341 return norm;
342 }
343 if (nsurf > 1)
344 {
345 return norm.unit(); // edge
346 }
347
348#ifdef G4SPECSDEBUG
349 std::ostringstream message;
350 G4long oldprc = message.precision(16);
351 message << "Point p is not on surface (!?) of solid: "
352 << GetName() << "\n";
353 message << "Position:\n";
354 message << " p.x() = " << p.x()/mm << " mm\n";
355 message << " p.y() = " << p.y()/mm << " mm\n";
356 message << " p.z() = " << p.z()/mm << " mm";
357 G4cout.precision(oldprc);
358 G4Exception("G4Ellipsoid::SurfaceNormal(p)", "GeomSolids1002",
359 JustWarning, message );
360 DumpInfo();
361#endif
362 return ApproxSurfaceNormal(p);
363}

The documentation for this class was generated from the following files: