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

#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 ()
 
G4bool operator== (const G4VSolid &s) const
 
G4String GetName () const
 
void SetName (const G4String &name)
 
G4double GetTolerance () const
 
void DumpInfo () const
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
 G4VSolid (__void__ &)
 
 G4VSolid (const G4VSolid &rhs)
 
G4VSolidoperator= (const G4VSolid &rhs)
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) 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

Definition at line 63 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. )

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:57

Referenced by Clone().

◆ ~G4Ellipsoid()

G4Ellipsoid::~G4Ellipsoid ( )
override

Definition at line 93 of file G4Ellipsoid.cc.

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

◆ G4Ellipsoid() [2/3]

G4Ellipsoid::G4Ellipsoid ( __void__ & a)

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)

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

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

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}
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override

◆ Clone()

G4VSolid * G4Ellipsoid::Clone ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 651 of file G4Ellipsoid.cc.

652{
653 return new G4Ellipsoid(*this);
654}
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

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

Reimplemented from G4VSolid.

Definition at line 912 of file G4Ellipsoid.cc.

913{
914 return new G4PolyhedronEllipsoid(fDx, fDy, fDz, fZBottomCut, fZTopCut);
915}

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

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

Implements G4VSolid.

Definition at line 894 of file G4Ellipsoid.cc.

895{
896 scene.AddSolid(*this);
897}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 470 of file G4Ellipsoid.cc.

471{
472 G4double px = p.x();
473 G4double py = p.y();
474 G4double pz = p.z();
475
476 // Safety distance to bounding box
477 G4double distX = std::abs(px) - fXmax;
478 G4double distY = std::abs(py) - fYmax;
479 G4double distZ = std::max(pz - fZTopCut, fZBottomCut - pz);
480 G4double distB = std::max(std::max(distX, distY), distZ);
481
482 // Safety distance to lateral surface
483 G4double x = px * fSx;
484 G4double y = py * fSy;
485 G4double z = pz * fSz;
486 G4double distR = std::sqrt(x*x + y*y + z*z) - fR;
487
488 // Return safety to in
489 G4double dist = std::max(distB, distR);
490 return (dist < 0.) ? 0. : dist;
491}
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

Implements G4VSolid.

Definition at line 384 of file G4Ellipsoid.cc.

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

Referenced by DistanceToIn().

◆ DistanceToOut() [1/2]

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

Implements G4VSolid.

Definition at line 622 of file G4Ellipsoid.cc.

623{
624 // Safety distance in z direction
625 G4double distZ = std::min(fZTopCut - p.z(), p.z() - fZBottomCut);
626
627 // Safety distance to lateral surface
628 G4double x = p.x() * fSx;
629 G4double y = p.y() * fSy;
630 G4double z = p.z() * fSz;
631 G4double distR = fR - std::sqrt(x*x + y*y + z*z);
632
633 // Return safety to out
634 G4double dist = std::min(distZ, distR);
635 return (dist < 0.) ? 0. : dist;
636}

◆ DistanceToOut() [2/2]

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

Implements G4VSolid.

Definition at line 497 of file G4Ellipsoid.cc.

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

Reimplemented from G4VSolid.

Definition at line 682 of file G4Ellipsoid.cc.

683{
684 if (fCubicVolume == 0.)
685 {
686 G4double piAB_3 = CLHEP::pi * fDx * fDy / 3.;
687 fCubicVolume = 4. * piAB_3 * fDz;
688 if (fZBottomCut > -fDz)
689 {
690 G4double hbot = 1. + fZBottomCut / fDz;
691 fCubicVolume -= piAB_3 * hbot * hbot * (2. * fDz - fZBottomCut);
692 }
693 if (fZTopCut < fDz)
694 {
695 G4double htop = 1. - fZTopCut / fDz;
696 fCubicVolume -= piAB_3 * htop * htop * (2. * fDz + fZTopCut);
697 }
698 }
699 return fCubicVolume;
700}

◆ GetDx()

G4double G4Ellipsoid::GetDx ( ) const
inline

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

Implements G4VSolid.

Definition at line 642 of file G4Ellipsoid.cc.

643{
644 return {"G4Ellipsoid"};
645}

Referenced by StreamInfo().

◆ GetExtent()

G4VisExtent G4Ellipsoid::GetExtent ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 903 of file G4Ellipsoid.cc.

904{
905 return { -fXmax, fXmax, -fYmax, fYmax, fZBottomCut, fZTopCut };
906}

◆ GetPointOnSurface()

G4ThreeVector G4Ellipsoid::GetPointOnSurface ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 816 of file G4Ellipsoid.cc.

817{
818 G4double A = GetDx();
819 G4double B = GetDy();
820 G4double C = GetDz();
821 G4double Zbot = GetZBottomCut();
822 G4double Ztop = GetZTopCut();
823
824 // Calculate cut areas
825 G4double Hbot = 1. + Zbot / C;
826 G4double Htop = 1. - Ztop / C;
827 G4double piAB = CLHEP::pi * A * B;
828 G4double Sbot = piAB * Hbot * (2. - Hbot);
829 G4double Stop = piAB * Htop * (2. - Htop);
830
831 // Get area of lateral surface
832 if (fLateralArea == 0.)
833 {
834 G4AutoLock l(&lateralareaMutex);
835 fLateralArea = LateralSurfaceArea();
836 l.unlock();
837 }
838 G4double Slat = fLateralArea;
839
840 // Select surface (0 - bottom cut, 1 - lateral surface, 2 - top cut)
841 G4double select = (Sbot + Slat + Stop) * G4QuickRand();
842 G4int k = 0;
843 if (select > Sbot) k = 1;
844 if (select > Sbot + Slat) k = 2;
845
846 // Pick random point on selected surface (rejection sampling)
848 switch (k)
849 {
850 case 0: // bootom z-cut
851 {
852 G4double scale = std::sqrt(Hbot * (2. - Hbot));
853 G4TwoVector rho = G4RandomPointInEllipse(A * scale, B * scale);
854 p.set(rho.x(), rho.y(), Zbot);
855 break;
856 }
857 case 1: // lateral surface
858 {
859 G4double x, y, z;
860 G4double mu_max = std::max(std::max(A * B, A * C), B * C);
861 for (G4int i = 0; i < 1000; ++i)
862 {
863 // generate random point on unit sphere
864 z = (Zbot + (Ztop - Zbot) * G4QuickRand()) / C;
865 G4double rho = std::sqrt((1. + z) * (1. - z));
866 G4double phi = CLHEP::twopi * G4QuickRand();
867 x = rho * std::cos(phi);
868 y = rho * std::sin(phi);
869 // check acceptance
870 G4double xbc = x * B * C;
871 G4double yac = y * A * C;
872 G4double zab = z * A * B;
873 G4double mu = std::sqrt(xbc * xbc + yac * yac + zab * zab);
874 if (mu_max * G4QuickRand() <= mu) break;
875 }
876 p.set(A * x, B * y, C * z);
877 break;
878 }
879 case 2: // top z-cut
880 {
881 G4double scale = std::sqrt(Htop * (2. - Htop));
882 G4TwoVector rho = G4RandomPointInEllipse(A * scale, B * scale);
883 p.set(rho.x(), rho.y(), Ztop);
884 break;
885 }
886 }
887 return p;
888}
G4double G4QuickRand()
G4TwoVector G4RandomPointInEllipse(G4double a, G4double b)
int G4int
Definition G4Types.hh:85
double x() const
double y() const
G4double GetDx() const
G4double GetDy() const
G4double GetZTopCut() const
G4double GetZBottomCut() const
G4double GetDz() const

◆ GetPolyhedron()

G4Polyhedron * G4Ellipsoid::GetPolyhedron ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 921 of file G4Ellipsoid.cc.

922{
923 if (fpPolyhedron == nullptr ||
924 fRebuildPolyhedron ||
926 fpPolyhedron->GetNumberOfRotationSteps())
927 {
928 G4AutoLock l(&polyhedronMutex);
929 delete fpPolyhedron;
930 fpPolyhedron = CreatePolyhedron();
931 fRebuildPolyhedron = false;
932 l.unlock();
933 }
934 return fpPolyhedron;
935}
G4Polyhedron * CreatePolyhedron() const override
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
static G4int GetNumberOfRotationSteps()

◆ GetSemiAxisMax()

◆ GetSurfaceArea()

G4double G4Ellipsoid::GetSurfaceArea ( )
overridevirtual

Reimplemented from G4VSolid.

Definition at line 792 of file G4Ellipsoid.cc.

793{
794 if (fSurfaceArea == 0.)
795 {
796 G4double piAB = CLHEP::pi * fDx * fDy;
797 fSurfaceArea = LateralSurfaceArea();
798 if (fZBottomCut > -fDz)
799 {
800 G4double hbot = 1. + fZBottomCut / fDz;
801 fSurfaceArea += piAB * hbot * (2. - hbot);
802 }
803 if (fZTopCut < fDz)
804 {
805 G4double htop = 1. - fZTopCut / fDz;
806 fSurfaceArea += piAB * htop * (2. - htop);
807 }
808 }
809 return fSurfaceArea;
810}

◆ GetZBottomCut()

◆ GetZTopCut()

◆ Inside()

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

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:107

◆ SetSemiAxis()

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

◆ SetZCuts()

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

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 660 of file G4Ellipsoid.cc.

661{
662 G4long oldprc = os.precision(16);
663 os << "-----------------------------------------------------------\n"
664 << " *** Dump for solid - " << GetName() << " ***\n"
665 << " ===================================================\n"
666 << " Solid type: " << GetEntityType() << "\n"
667 << " Parameters: \n"
668 << " semi-axis x: " << GetDx()/mm << " mm \n"
669 << " semi-axis y: " << GetDy()/mm << " mm \n"
670 << " semi-axis z: " << GetDz()/mm << " mm \n"
671 << " lower cut in z: " << GetZBottomCut()/mm << " mm \n"
672 << " upper cut in z: " << GetZTopCut()/mm << " mm \n"
673 << "-----------------------------------------------------------\n";
674 os.precision(oldprc);
675 return os;
676}
G4GeometryType GetEntityType() const override

◆ SurfaceNormal()

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

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) return norm;
340 else if (nsurf > 1) return norm.unit(); // edge
341 else
342 {
343#ifdef G4SPECSDEBUG
344 std::ostringstream message;
345 G4long oldprc = message.precision(16);
346 message << "Point p is not on surface (!?) of solid: "
347 << GetName() << "\n";
348 message << "Position:\n";
349 message << " p.x() = " << p.x()/mm << " mm\n";
350 message << " p.y() = " << p.y()/mm << " mm\n";
351 message << " p.z() = " << p.z()/mm << " mm";
352 G4cout.precision(oldprc);
353 G4Exception("G4Ellipsoid::SurfaceNormal(p)", "GeomSolids1002",
354 JustWarning, message );
355 DumpInfo();
356#endif
357 return ApproxSurfaceNormal(p);
358 }
359}

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