Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4Ellipsoid Class Reference

#include <G4Ellipsoid.hh>

Inheritance diagram for G4Ellipsoid:
G4VSolid

Public Member Functions

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

Protected Member Functions

void CalculateClippedPolygonExtent (G4ThreeVectorList &pPolygon, 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 ClipCrossSection (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

G4double kCarTolerance
 

Private Member Functions

G4ThreeVector ApproxSurfaceNormal (const G4ThreeVector &p) const
 
void CheckParameters ()
 
void ClipPolygonToSimpleLimits (G4ThreeVectorList &pPolygon, G4ThreeVectorList &outputPolygon, const G4VoxelLimits &pVoxelLimit) const
 
G4double LateralSurfaceArea () const
 

Private Attributes

G4double fCubicVolume = 0.0
 
G4double fDx
 
G4double fDy
 
G4double fDz
 
G4double fLateralArea = 0.0
 
G4PolyhedronfpPolyhedron = nullptr
 
G4double fQ1
 
G4double fQ2
 
G4double fR
 
G4bool fRebuildPolyhedron = false
 
G4double fRsph
 
G4String fshapeName
 
G4double fSurfaceArea = 0.0
 
G4double fSx
 
G4double fSy
 
G4double fSz
 
G4double fXmax
 
G4double fYmax
 
G4double fZBottomCut
 
G4double fZDimCut
 
G4double fZMidCut
 
G4double fZTopCut
 
G4double halfTolerance
 

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{
77}
G4double fZTopCut
Definition: G4Ellipsoid.hh:161
G4double fDx
Definition: G4Ellipsoid.hh:157
G4double fDy
Definition: G4Ellipsoid.hh:158
G4double fZBottomCut
Definition: G4Ellipsoid.hh:160
G4double fDz
Definition: G4Ellipsoid.hh:159
void CheckParameters()
Definition: G4Ellipsoid.cc:167
G4VSolid(const G4String &name)
Definition: G4VSolid.cc:57
const char * name(G4int ptype)

References CheckParameters().

Referenced by Clone().

◆ ~G4Ellipsoid()

G4Ellipsoid::~G4Ellipsoid ( )
virtual

Definition at line 93 of file G4Ellipsoid.cc.

94{
95 delete fpPolyhedron; fpPolyhedron = nullptr;
96}
G4Polyhedron * fpPolyhedron
Definition: G4Ellipsoid.hh:181

References fpPolyhedron.

◆ 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),
107 fXmax(rhs.fXmax), fYmax(rhs.fYmax),
108 fRsph(rhs.fRsph), fR(rhs.fR),
109 fSx(rhs.fSx), fSy(rhs.fSy), fSz(rhs.fSz),
111 fQ1(rhs.fQ1), fQ2(rhs.fQ2),
115{
116}
G4double fSurfaceArea
Definition: G4Ellipsoid.hh:178
G4double fRsph
Definition: G4Ellipsoid.hh:167
G4double fCubicVolume
Definition: G4Ellipsoid.hh:177
G4double fZDimCut
Definition: G4Ellipsoid.hh:173
G4double fSx
Definition: G4Ellipsoid.hh:169
G4double fYmax
Definition: G4Ellipsoid.hh:166
G4double halfTolerance
Definition: G4Ellipsoid.hh:164
G4double fR
Definition: G4Ellipsoid.hh:168
G4double fSz
Definition: G4Ellipsoid.hh:171
G4double fSy
Definition: G4Ellipsoid.hh:170
G4double fQ2
Definition: G4Ellipsoid.hh:175
G4double fQ1
Definition: G4Ellipsoid.hh:174
G4double fLateralArea
Definition: G4Ellipsoid.hh:179
G4double fZMidCut
Definition: G4Ellipsoid.hh:172
G4double fXmax
Definition: G4Ellipsoid.hh:165

Member Function Documentation

◆ ApproxSurfaceNormal()

G4ThreeVector G4Ellipsoid::ApproxSurfaceNormal ( const G4ThreeVector p) const
private

Definition at line 366 of file G4Ellipsoid.cc.

367{
368 G4double x = p.x() * fSx;
369 G4double y = p.y() * fSy;
370 G4double z = p.z() * fSz;
371 G4double rr = x * x + y * y + z * z;
372 G4double distZ = std::abs(z - fZMidCut) - fZDimCut;
373 G4double distR = std::sqrt(rr) - fR;
374 if (distR > distZ && rr > 0.) // distR > distZ is correct!
375 return G4ThreeVector(x*fSx, y*fSy, z*fSz).unit();
376 else
377 return G4ThreeVector(0., 0., std::copysign(1., z - fZMidCut));
378}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
double z() const
Hep3Vector unit() const
double x() const
double y() const

References fR, fSx, fSy, fSz, fZDimCut, fZMidCut, CLHEP::Hep3Vector::unit(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by DistanceToOut(), and SurfaceNormal().

◆ BoundingLimits()

void G4Ellipsoid::BoundingLimits ( G4ThreeVector pMin,
G4ThreeVector pMax 
) const
virtual

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}
static const G4double pMax
static const G4double pMin

References fXmax, fYmax, fZBottomCut, fZTopCut, pMax, and pMin.

Referenced by CalculateExtent().

◆ CalculateClippedPolygonExtent()

void G4VSolid::CalculateClippedPolygonExtent ( G4ThreeVectorList pPolygon,
const G4VoxelLimits pVoxelLimit,
const EAxis  pAxis,
G4double pMin,
G4double pMax 
) const
protectedinherited

Definition at line 489 of file G4VSolid.cc.

494{
495 G4int noLeft,i;
496 G4double component;
497
498 ClipPolygon(pPolygon,pVoxelLimit,pAxis);
499 noLeft = pPolygon.size();
500
501 if ( noLeft )
502 {
503 for (i=0; i<noLeft; ++i)
504 {
505 component = pPolygon[i].operator()(pAxis);
506
507 if (component < pMin)
508 {
509 pMin = component;
510 }
511 if (component > pMax)
512 {
513 pMax = component;
514 }
515 }
516 }
517}
int G4int
Definition: G4Types.hh:85
void ClipPolygon(G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
Definition: G4VSolid.cc:539

References G4VSolid::ClipPolygon(), pMax, and pMin.

Referenced by G4VSolid::ClipBetweenSections(), and G4VSolid::ClipCrossSection().

◆ CalculateExtent()

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

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
Definition: G4Ellipsoid.cc:264

References BoundingLimits(), G4BoundingEnvelope::CalculateExtent(), pMax, and pMin.

◆ CheckParameters()

void G4Ellipsoid::CheckParameters ( )
private

Definition at line 167 of file G4Ellipsoid.cc.

168{
169 halfTolerance = 0.5 * kCarTolerance; // half tolerance
170 G4double dmin = 2. * kCarTolerance;
171
172 // Check dimensions
173 //
174 if (fDx < dmin || fDy < dmin || fDz < dmin)
175 {
176 std::ostringstream message;
177 message << "Invalid (too small or negative) dimensions for Solid: "
178 << GetName() << "\n"
179 << " semi-axis x: " << fDx << "\n"
180 << " semi-axis y: " << fDy << "\n"
181 << " semi-axis z: " << fDz;
182 G4Exception("G4Ellipsoid::CheckParameters()", "GeomSolids0002",
183 FatalException, message);
184 }
185 G4double A = fDx;
186 G4double B = fDy;
187 G4double C = fDz;
188
189 // Check cuts
190 //
191 if (fZBottomCut == 0. && fZTopCut == 0.)
192 {
193 fZBottomCut = -C;
194 fZTopCut = C;
195 }
196 if (fZBottomCut >= C || fZTopCut <= -C || fZBottomCut >= fZTopCut)
197 {
198 std::ostringstream message;
199 message << "Invalid Z cuts for Solid: "
200 << GetName() << "\n"
201 << " bottom cut: " << fZBottomCut << "\n"
202 << " top cut: " << fZTopCut;
203 G4Exception("G4Ellipsoid::CheckParameters()", "GeomSolids0002",
204 FatalException, message);
205
206 }
209
210 // Set extent in x and y
211 fXmax = A;
212 fYmax = B;
213 if (fZBottomCut > 0.)
214 {
215 G4double ratio = fZBottomCut / C;
216 G4double scale = std::sqrt((1. - ratio) * (1 + ratio));
217 fXmax *= scale;
218 fYmax *= scale;
219 }
220 if (fZTopCut < 0.)
221 {
222 G4double ratio = fZTopCut / C;
223 G4double scale = std::sqrt((1. - ratio) * (1 + ratio));
224 fXmax *= scale;
225 fYmax *= scale;
226 }
227
228 // Set scale factors
229 fRsph = std::max(std::max(A, B), C); // bounding sphere
230 fR = std::min(std::min(A, B), C); // radius of sphere after scaling
231 fSx = fR / A; // X scale factor
232 fSy = fR / B; // Y scale factor
233 fSz = fR / C; // Z scale factor
234
235 // Scaled cuts
236 fZMidCut = 0.5 * (fZTopCut + fZBottomCut) * fSz; // middle position
237 fZDimCut = 0.5 * (fZTopCut - fZBottomCut) * fSz; // half distance
238
239 // Coefficients for approximation of distance: Q1 * (x^2 + y^2 + z^2) - Q2
240 fQ1 = 0.5 / fR;
241 fQ2 = 0.5 * fR + halfTolerance * halfTolerance * fQ1;
242
243 fCubicVolume = 0.; // volume
244 fSurfaceArea = 0.; // surface area
245 fLateralArea = 0.; // lateral surface area
246}
G4double C(G4double temp)
G4double B(G4double temperature)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
const G4double A[17]
G4String GetName() const
G4double kCarTolerance
Definition: G4VSolid.hh:299
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

References A, B(), C(), FatalException, fCubicVolume, fDx, fDy, fDz, fLateralArea, fQ1, fQ2, fR, fRsph, fSurfaceArea, fSx, fSy, fSz, fXmax, fYmax, fZBottomCut, fZDimCut, fZMidCut, fZTopCut, G4Exception(), G4VSolid::GetName(), halfTolerance, G4VSolid::kCarTolerance, G4INCL::Math::max(), and G4INCL::Math::min().

Referenced by G4Ellipsoid().

◆ ClipBetweenSections()

void G4VSolid::ClipBetweenSections ( G4ThreeVectorList pVertices,
const G4int  pSectionIndex,
const G4VoxelLimits pVoxelLimit,
const EAxis  pAxis,
G4double pMin,
G4double pMax 
) const
protectedinherited

Definition at line 444 of file G4VSolid.cc.

449{
450 G4ThreeVectorList polygon;
451 polygon.reserve(4);
452 polygon.push_back((*pVertices)[pSectionIndex]);
453 polygon.push_back((*pVertices)[pSectionIndex+4]);
454 polygon.push_back((*pVertices)[pSectionIndex+5]);
455 polygon.push_back((*pVertices)[pSectionIndex+1]);
456 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
457 polygon.clear();
458
459 polygon.push_back((*pVertices)[pSectionIndex+1]);
460 polygon.push_back((*pVertices)[pSectionIndex+5]);
461 polygon.push_back((*pVertices)[pSectionIndex+6]);
462 polygon.push_back((*pVertices)[pSectionIndex+2]);
463 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
464 polygon.clear();
465
466 polygon.push_back((*pVertices)[pSectionIndex+2]);
467 polygon.push_back((*pVertices)[pSectionIndex+6]);
468 polygon.push_back((*pVertices)[pSectionIndex+7]);
469 polygon.push_back((*pVertices)[pSectionIndex+3]);
470 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
471 polygon.clear();
472
473 polygon.push_back((*pVertices)[pSectionIndex+3]);
474 polygon.push_back((*pVertices)[pSectionIndex+7]);
475 polygon.push_back((*pVertices)[pSectionIndex+4]);
476 polygon.push_back((*pVertices)[pSectionIndex]);
477 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
478 return;
479}
std::vector< G4ThreeVector > G4ThreeVectorList
void CalculateClippedPolygonExtent(G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:489

References G4VSolid::CalculateClippedPolygonExtent(), pMax, and pMin.

◆ ClipCrossSection()

void G4VSolid::ClipCrossSection ( G4ThreeVectorList pVertices,
const G4int  pSectionIndex,
const G4VoxelLimits pVoxelLimit,
const EAxis  pAxis,
G4double pMin,
G4double pMax 
) const
protectedinherited

Definition at line 414 of file G4VSolid.cc.

419{
420
421 G4ThreeVectorList polygon;
422 polygon.reserve(4);
423 polygon.push_back((*pVertices)[pSectionIndex]);
424 polygon.push_back((*pVertices)[pSectionIndex+1]);
425 polygon.push_back((*pVertices)[pSectionIndex+2]);
426 polygon.push_back((*pVertices)[pSectionIndex+3]);
427 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
428 return;
429}

References G4VSolid::CalculateClippedPolygonExtent(), pMax, and pMin.

◆ ClipPolygon()

void G4VSolid::ClipPolygon ( G4ThreeVectorList pPolygon,
const G4VoxelLimits pVoxelLimit,
const EAxis  pAxis 
) const
protectedinherited

Definition at line 539 of file G4VSolid.cc.

542{
543 G4ThreeVectorList outputPolygon;
544
545 if ( pVoxelLimit.IsLimited() )
546 {
547 if (pVoxelLimit.IsXLimited() ) // && pAxis != kXAxis)
548 {
549 G4VoxelLimits simpleLimit1;
550 simpleLimit1.AddLimit(kXAxis,pVoxelLimit.GetMinXExtent(),kInfinity);
551 ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
552
553 pPolygon.clear();
554
555 if ( !outputPolygon.size() ) return;
556
557 G4VoxelLimits simpleLimit2;
558 simpleLimit2.AddLimit(kXAxis,-kInfinity,pVoxelLimit.GetMaxXExtent());
559 ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
560
561 if ( !pPolygon.size() ) return;
562 else outputPolygon.clear();
563 }
564 if ( pVoxelLimit.IsYLimited() ) // && pAxis != kYAxis)
565 {
566 G4VoxelLimits simpleLimit1;
567 simpleLimit1.AddLimit(kYAxis,pVoxelLimit.GetMinYExtent(),kInfinity);
568 ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
569
570 // Must always clear pPolygon - for clip to simpleLimit2 and in case of
571 // early exit
572
573 pPolygon.clear();
574
575 if ( !outputPolygon.size() ) return;
576
577 G4VoxelLimits simpleLimit2;
578 simpleLimit2.AddLimit(kYAxis,-kInfinity,pVoxelLimit.GetMaxYExtent());
579 ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
580
581 if ( !pPolygon.size() ) return;
582 else outputPolygon.clear();
583 }
584 if ( pVoxelLimit.IsZLimited() ) // && pAxis != kZAxis)
585 {
586 G4VoxelLimits simpleLimit1;
587 simpleLimit1.AddLimit(kZAxis,pVoxelLimit.GetMinZExtent(),kInfinity);
588 ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
589
590 // Must always clear pPolygon - for clip to simpleLimit2 and in case of
591 // early exit
592
593 pPolygon.clear();
594
595 if ( !outputPolygon.size() ) return;
596
597 G4VoxelLimits simpleLimit2;
598 simpleLimit2.AddLimit(kZAxis,-kInfinity,pVoxelLimit.GetMaxZExtent());
599 ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
600
601 // Return after final clip - no cleanup
602 }
603 }
604}
void ClipPolygonToSimpleLimits(G4ThreeVectorList &pPolygon, G4ThreeVectorList &outputPolygon, const G4VoxelLimits &pVoxelLimit) const
Definition: G4VSolid.cc:612
G4bool IsYLimited() const
G4double GetMinZExtent() const
void AddLimit(const EAxis pAxis, const G4double pMin, const G4double pMax)
G4bool IsXLimited() const
G4double GetMaxYExtent() const
G4double GetMaxZExtent() const
G4double GetMinYExtent() const
G4double GetMinXExtent() const
G4bool IsZLimited() const
G4bool IsLimited() const
G4double GetMaxXExtent() const
@ kYAxis
Definition: geomdefs.hh:56
@ kXAxis
Definition: geomdefs.hh:55
@ kZAxis
Definition: geomdefs.hh:57
static const G4double kInfinity
Definition: geomdefs.hh:41

References G4VoxelLimits::AddLimit(), G4VSolid::ClipPolygonToSimpleLimits(), G4VoxelLimits::GetMaxXExtent(), G4VoxelLimits::GetMaxYExtent(), G4VoxelLimits::GetMaxZExtent(), G4VoxelLimits::GetMinXExtent(), G4VoxelLimits::GetMinYExtent(), G4VoxelLimits::GetMinZExtent(), G4VoxelLimits::IsLimited(), G4VoxelLimits::IsXLimited(), G4VoxelLimits::IsYLimited(), G4VoxelLimits::IsZLimited(), kInfinity, kXAxis, kYAxis, and kZAxis.

Referenced by G4VSolid::CalculateClippedPolygonExtent().

◆ ClipPolygonToSimpleLimits()

void G4VSolid::ClipPolygonToSimpleLimits ( G4ThreeVectorList pPolygon,
G4ThreeVectorList outputPolygon,
const G4VoxelLimits pVoxelLimit 
) const
privateinherited

Definition at line 612 of file G4VSolid.cc.

615{
616 G4int i;
617 G4int noVertices=pPolygon.size();
618 G4ThreeVector vEnd,vStart;
619
620 for (i = 0 ; i < noVertices ; ++i )
621 {
622 vStart = pPolygon[i];
623 if ( i == noVertices-1 ) vEnd = pPolygon[0];
624 else vEnd = pPolygon[i+1];
625
626 if ( pVoxelLimit.Inside(vStart) )
627 {
628 if (pVoxelLimit.Inside(vEnd))
629 {
630 // vStart and vEnd inside -> output end point
631 //
632 outputPolygon.push_back(vEnd);
633 }
634 else
635 {
636 // vStart inside, vEnd outside -> output crossing point
637 //
638 pVoxelLimit.ClipToLimits(vStart,vEnd);
639 outputPolygon.push_back(vEnd);
640 }
641 }
642 else
643 {
644 if (pVoxelLimit.Inside(vEnd))
645 {
646 // vStart outside, vEnd inside -> output inside section
647 //
648 pVoxelLimit.ClipToLimits(vStart,vEnd);
649 outputPolygon.push_back(vStart);
650 outputPolygon.push_back(vEnd);
651 }
652 else // Both point outside -> no output
653 {
654 // outputPolygon.push_back(vStart);
655 // outputPolygon.push_back(vEnd);
656 }
657 }
658 }
659}
G4bool ClipToLimits(G4ThreeVector &pStart, G4ThreeVector &pEnd) const
G4bool Inside(const G4ThreeVector &pVec) const

References G4VoxelLimits::ClipToLimits(), and G4VoxelLimits::Inside().

Referenced by G4VSolid::ClipPolygon().

◆ Clone()

G4VSolid * G4Ellipsoid::Clone ( ) const
virtual

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.)
Definition: G4Ellipsoid.cc:67

References G4Ellipsoid().

◆ ComputeDimensions()

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

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

References G4VPVParameterisation::ComputeDimensions(), and CLHEP::detail::n.

◆ CreatePolyhedron()

G4Polyhedron * G4Ellipsoid::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 877 of file G4Ellipsoid.cc.

References fDx, fDy, fDz, fZBottomCut, and fZTopCut.

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

void G4Ellipsoid::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 859 of file G4Ellipsoid.cc.

860{
861 scene.AddSolid(*this);
862}
virtual void AddSolid(const G4Box &)=0

References G4VGraphicsScene::AddSolid().

◆ DistanceToIn() [1/2]

G4double G4Ellipsoid::DistanceToIn ( const G4ThreeVector p) const
virtual

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}

References fR, fSx, fSy, fSz, fXmax, fYmax, fZBottomCut, fZTopCut, G4INCL::Math::max(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ DistanceToIn() [2/2]

G4double G4Ellipsoid::DistanceToIn ( const G4ThreeVector p,
const G4ThreeVector v 
) const
virtual

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 D(G4double temp)
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Ellipsoid.cc:384
#define EPS
#define DBL_MAX
Definition: templates.hh:62

References A, B(), C(), D(), DBL_MAX, DistanceToIn(), EPS, fQ1, fQ2, fR, fRsph, fSx, fSy, fSz, fXmax, fYmax, fZBottomCut, fZDimCut, fZMidCut, fZTopCut, halfTolerance, G4VSolid::kCarTolerance, kInfinity, G4INCL::Math::max(), G4INCL::Math::min(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by DistanceToIn().

◆ DistanceToOut() [1/2]

G4double G4Ellipsoid::DistanceToOut ( const G4ThreeVector p) const
virtual

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}

References fR, fSx, fSy, fSz, fZBottomCut, fZTopCut, G4INCL::Math::min(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ DistanceToOut() [2/2]

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

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 G4int 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;
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
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Ellipsoid.cc:366
void DumpInfo() const
#define DBL_EPSILON
Definition: templates.hh:66

References A, ApproxSurfaceNormal(), B(), C(), D(), DBL_EPSILON, DBL_MAX, G4VSolid::DumpInfo(), EPS, fQ1, fQ2, fR, fSx, fSy, fSz, fZDimCut, fZMidCut, G4cout, G4endl, G4Exception(), G4VSolid::GetName(), halfTolerance, JustWarning, G4INCL::Math::max(), G4INCL::Math::min(), CLHEP::detail::n, CLHEP::Hep3Vector::unit(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ DumpInfo()

void G4VSolid::DumpInfo ( ) const
inlineinherited

Referenced by G4Cons::ApproxSurfaceNormal(), G4CutTubs::ApproxSurfaceNormal(), G4Sphere::ApproxSurfaceNormal(), G4Torus::ApproxSurfaceNormal(), G4Tubs::ApproxSurfaceNormal(), G4ReflectedSolid::BoundingLimits(), G4DisplacedSolid::BoundingLimits(), G4IntersectionSolid::BoundingLimits(), G4ScaledSolid::BoundingLimits(), G4SubtractionSolid::BoundingLimits(), G4UnionSolid::BoundingLimits(), G4Box::BoundingLimits(), G4Cons::BoundingLimits(), G4CutTubs::BoundingLimits(), G4Orb::BoundingLimits(), G4Para::BoundingLimits(), G4Sphere::BoundingLimits(), G4Torus::BoundingLimits(), G4Trap::BoundingLimits(), G4Trd::BoundingLimits(), G4Tubs::BoundingLimits(), G4EllipticalCone::BoundingLimits(), G4ExtrudedSolid::BoundingLimits(), G4GenericPolycone::BoundingLimits(), G4GenericTrap::BoundingLimits(), G4Hype::BoundingLimits(), G4Paraboloid::BoundingLimits(), G4Polycone::BoundingLimits(), G4Polyhedra::BoundingLimits(), G4TessellatedSolid::BoundingLimits(), G4TwistedTubs::BoundingLimits(), G4ParameterisationBoxX::ComputeDimensions(), G4ParameterisationBoxY::ComputeDimensions(), G4ParameterisationBoxZ::ComputeDimensions(), G4ParameterisationConsRho::ComputeDimensions(), G4ParameterisationConsPhi::ComputeDimensions(), G4ParameterisationConsZ::ComputeDimensions(), G4ParameterisationParaX::ComputeDimensions(), G4ParameterisationParaY::ComputeDimensions(), G4ParameterisationParaZ::ComputeDimensions(), G4ParameterisationPolyconeRho::ComputeDimensions(), G4ParameterisationPolyconePhi::ComputeDimensions(), G4ParameterisationPolyconeZ::ComputeDimensions(), G4ParameterisationPolyhedraRho::ComputeDimensions(), G4ParameterisationPolyhedraPhi::ComputeDimensions(), G4ParameterisationPolyhedraZ::ComputeDimensions(), G4ParameterisationTrdX::ComputeDimensions(), G4ParameterisationTrdY::ComputeDimensions(), G4ParameterisationTrdZ::ComputeDimensions(), G4ParameterisationTubsRho::ComputeDimensions(), G4ParameterisationTubsPhi::ComputeDimensions(), G4ParameterisationTubsZ::ComputeDimensions(), G4ReflectedSolid::ComputeDimensions(), G4DisplacedSolid::ComputeDimensions(), G4ScaledSolid::ComputeDimensions(), G4ParameterisedNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4DisplacedSolid::CreatePolyhedron(), G4ScaledSolid::CreatePolyhedron(), G4SubtractionSolid::DistanceToIn(), G4Box::DistanceToOut(), G4Orb::DistanceToOut(), G4Para::DistanceToOut(), G4Trap::DistanceToOut(), G4Trd::DistanceToOut(), G4Paraboloid::DistanceToOut(), G4VTwistedFaceted::DistanceToOut(), G4Cons::DistanceToOut(), G4CutTubs::DistanceToOut(), G4Sphere::DistanceToOut(), G4Torus::DistanceToOut(), G4Tubs::DistanceToOut(), DistanceToOut(), G4EllipticalCone::DistanceToOut(), G4EllipticalTube::DistanceToOut(), G4GenericTrap::DistanceToOut(), export_G4VSolid(), G4Polycone::G4Polycone(), G4Polyhedra::G4Polyhedra(), G4BooleanSolid::GetConstituentSolid(), G4NavigationLogger::PostComputeStepLog(), G4Box::SurfaceNormal(), G4Para::SurfaceNormal(), G4Trap::SurfaceNormal(), G4Trd::SurfaceNormal(), SurfaceNormal(), G4EllipticalCone::SurfaceNormal(), G4EllipticalTube::SurfaceNormal(), G4ExtrudedSolid::SurfaceNormal(), and G4Tet::SurfaceNormal().

◆ EstimateCubicVolume()

G4double G4VSolid::EstimateCubicVolume ( G4int  nStat,
G4double  epsilon 
) const
inherited

Definition at line 203 of file G4VSolid.cc.

204{
205 G4int iInside=0;
206 G4double px,py,pz,minX,maxX,minY,maxY,minZ,maxZ,volume,halfepsilon;
208 EInside in;
209
210 // values needed for CalculateExtent signature
211
212 G4VoxelLimits limit; // Unlimited
213 G4AffineTransform origin;
214
215 // min max extents of pSolid along X,Y,Z
216
217 CalculateExtent(kXAxis,limit,origin,minX,maxX);
218 CalculateExtent(kYAxis,limit,origin,minY,maxY);
219 CalculateExtent(kZAxis,limit,origin,minZ,maxZ);
220
221 // limits
222
223 if(nStat < 100) nStat = 100;
224 if(epsilon > 0.01) epsilon = 0.01;
225 halfepsilon = 0.5*epsilon;
226
227 for(auto i = 0; i < nStat; ++i )
228 {
229 px = minX-halfepsilon+(maxX-minX+epsilon)*G4QuickRand();
230 py = minY-halfepsilon+(maxY-minY+epsilon)*G4QuickRand();
231 pz = minZ-halfepsilon+(maxZ-minZ+epsilon)*G4QuickRand();
232 p = G4ThreeVector(px,py,pz);
233 in = Inside(p);
234 if(in != kOutside) ++iInside;
235 }
236 volume = (maxX-minX+epsilon)*(maxY-minY+epsilon)
237 * (maxZ-minZ+epsilon)*iInside/nStat;
238 return volume;
239}
G4double epsilon(G4double density, G4double temperature)
static const G4int maxZ
G4double G4QuickRand()
Definition: G4QuickRand.hh:34
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
virtual EInside Inside(const G4ThreeVector &p) const =0
EInside
Definition: geomdefs.hh:67
@ kOutside
Definition: geomdefs.hh:68

References G4VSolid::CalculateExtent(), epsilon(), G4QuickRand(), G4VSolid::Inside(), kOutside, kXAxis, kYAxis, kZAxis, and maxZ.

Referenced by G4VSolid::GetCubicVolume(), G4BooleanSolid::GetCubicVolume(), and G4VCSGfaceted::GetCubicVolume().

◆ EstimateSurfaceArea()

G4double G4VSolid::EstimateSurfaceArea ( G4int  nStat,
G4double  ell 
) const
inherited

Definition at line 265 of file G4VSolid.cc.

266{
267 static const G4double s2 = 1./std::sqrt(2.);
268 static const G4double s3 = 1./std::sqrt(3.);
269 static const G4ThreeVector directions[64] =
270 {
271 G4ThreeVector( 0, 0, 0), G4ThreeVector( -1, 0, 0), // ( , , ) ( -, , )
272 G4ThreeVector( 1, 0, 0), G4ThreeVector( -1, 0, 0), // ( +, , ) (-+, , )
273 G4ThreeVector( 0, -1, 0), G4ThreeVector(-s2,-s2, 0), // ( , -, ) ( -, -, )
274 G4ThreeVector( s2, -s2, 0), G4ThreeVector( 0, -1, 0), // ( +, -, ) (-+, -, )
275
276 G4ThreeVector( 0, 1, 0), G4ThreeVector( -s2, s2, 0), // ( , +, ) ( -, +, )
277 G4ThreeVector( s2, s2, 0), G4ThreeVector( 0, 1, 0), // ( +, +, ) (-+, +, )
278 G4ThreeVector( 0, -1, 0), G4ThreeVector( -1, 0, 0), // ( ,-+, ) ( -,-+, )
279 G4ThreeVector( 1, 0, 0), G4ThreeVector( -1, 0, 0), // ( +,-+, ) (-+,-+, )
280
281 G4ThreeVector( 0, 0, -1), G4ThreeVector(-s2, 0,-s2), // ( , , -) ( -, , -)
282 G4ThreeVector( s2, 0,-s2), G4ThreeVector( 0, 0, -1), // ( +, , -) (-+, , -)
283 G4ThreeVector( 0,-s2,-s2), G4ThreeVector(-s3,-s3,-s3), // ( , -, -) ( -, -, -)
284 G4ThreeVector( s3,-s3,-s3), G4ThreeVector( 0,-s2,-s2), // ( +, -, -) (-+, -, -)
285
286 G4ThreeVector( 0, s2,-s2), G4ThreeVector(-s3, s3,-s3), // ( , +, -) ( -, +, -)
287 G4ThreeVector( s3, s3,-s3), G4ThreeVector( 0, s2,-s2), // ( +, +, -) (-+, +, -)
288 G4ThreeVector( 0, 0, -1), G4ThreeVector(-s2, 0,-s2), // ( ,-+, -) ( -,-+, -)
289 G4ThreeVector( s2, 0,-s2), G4ThreeVector( 0, 0, -1), // ( +,-+, -) (-+,-+, -)
290
291 G4ThreeVector( 0, 0, 1), G4ThreeVector(-s2, 0, s2), // ( , , +) ( -, , +)
292 G4ThreeVector( s2, 0, s2), G4ThreeVector( 0, 0, 1), // ( +, , +) (-+, , +)
293 G4ThreeVector( 0,-s2, s2), G4ThreeVector(-s3,-s3, s3), // ( , -, +) ( -, -, +)
294 G4ThreeVector( s3,-s3, s3), G4ThreeVector( 0,-s2, s2), // ( +, -, +) (-+, -, +)
295
296 G4ThreeVector( 0, s2, s2), G4ThreeVector(-s3, s3, s3), // ( , +, +) ( -, +, +)
297 G4ThreeVector( s3, s3, s3), G4ThreeVector( 0, s2, s2), // ( +, +, +) (-+, +, +)
298 G4ThreeVector( 0, 0, 1), G4ThreeVector(-s2, 0, s2), // ( ,-+, +) ( -,-+, +)
299 G4ThreeVector( s2, 0, s2), G4ThreeVector( 0, 0, 1), // ( +,-+, +) (-+,-+, +)
300
301 G4ThreeVector( 0, 0, -1), G4ThreeVector( -1, 0, 0), // ( , ,-+) ( -, ,-+)
302 G4ThreeVector( 1, 0, 0), G4ThreeVector( -1, 0, 0), // ( +, ,-+) (-+, ,-+)
303 G4ThreeVector( 0, -1, 0), G4ThreeVector(-s2,-s2, 0), // ( , -,-+) ( -, -,-+)
304 G4ThreeVector( s2, -s2, 0), G4ThreeVector( 0, -1, 0), // ( +, -,-+) (-+, -,-+)
305
306 G4ThreeVector( 0, 1, 0), G4ThreeVector( -s2, s2, 0), // ( , +,-+) ( -, +,-+)
307 G4ThreeVector( s2, s2, 0), G4ThreeVector( 0, 1, 0), // ( +, +,-+) (-+, +,-+)
308 G4ThreeVector( 0, -1, 0), G4ThreeVector( -1, 0, 0), // ( ,-+,-+) ( -,-+,-+)
309 G4ThreeVector( 1, 0, 0), G4ThreeVector( -1, 0, 0), // ( +,-+,-+) (-+,-+,-+)
310 };
311
312 G4ThreeVector bmin, bmax;
313 BoundingLimits(bmin, bmax);
314
315 G4double dX = bmax.x() - bmin.x();
316 G4double dY = bmax.y() - bmin.y();
317 G4double dZ = bmax.z() - bmin.z();
318
319 // Define statistics and shell thickness
320 //
321 G4int npoints = (nstat < 1000) ? 1000 : nstat;
322 G4double coeff = 0.5 / std::cbrt(G4double(npoints));
323 G4double eps = (ell > 0) ? ell : coeff * std::min(std::min(dX, dY), dZ);
324 G4double del = 1.8 * eps; // shold be more than sqrt(3.)
325
326 G4double minX = bmin.x() - eps;
327 G4double minY = bmin.y() - eps;
328 G4double minZ = bmin.z() - eps;
329
330 G4double dd = 2. * eps;
331 dX += dd;
332 dY += dd;
333 dZ += dd;
334
335 // Calculate surface area
336 //
337 G4int icount = 0;
338 for(auto i = 0; i < npoints; ++i)
339 {
340 G4double px = minX + dX*G4QuickRand();
341 G4double py = minY + dY*G4QuickRand();
342 G4double pz = minZ + dZ*G4QuickRand();
343 G4ThreeVector p = G4ThreeVector(px, py, pz);
344 EInside in = Inside(p);
345 G4double dist = 0;
346 if (in == kInside)
347 {
348 if (DistanceToOut(p) >= eps) continue;
349 G4int icase = 0;
350 if (Inside(G4ThreeVector(px-del, py, pz)) != kInside) icase += 1;
351 if (Inside(G4ThreeVector(px+del, py, pz)) != kInside) icase += 2;
352 if (Inside(G4ThreeVector(px, py-del, pz)) != kInside) icase += 4;
353 if (Inside(G4ThreeVector(px, py+del, pz)) != kInside) icase += 8;
354 if (Inside(G4ThreeVector(px, py, pz-del)) != kInside) icase += 16;
355 if (Inside(G4ThreeVector(px, py, pz+del)) != kInside) icase += 32;
356 if (icase == 0) continue;
357 G4ThreeVector v = directions[icase];
358 dist = DistanceToOut(p, v);
359 G4ThreeVector n = SurfaceNormal(p + v*dist);
360 dist *= v.dot(n);
361 }
362 else if (in == kOutside)
363 {
364 if (DistanceToIn(p) >= eps) continue;
365 G4int icase = 0;
366 if (Inside(G4ThreeVector(px-del, py, pz)) != kOutside) icase += 1;
367 if (Inside(G4ThreeVector(px+del, py, pz)) != kOutside) icase += 2;
368 if (Inside(G4ThreeVector(px, py-del, pz)) != kOutside) icase += 4;
369 if (Inside(G4ThreeVector(px, py+del, pz)) != kOutside) icase += 8;
370 if (Inside(G4ThreeVector(px, py, pz-del)) != kOutside) icase += 16;
371 if (Inside(G4ThreeVector(px, py, pz+del)) != kOutside) icase += 32;
372 if (icase == 0) continue;
373 G4ThreeVector v = directions[icase];
374 dist = DistanceToIn(p, v);
375 if (dist == kInfinity) continue;
376 G4ThreeVector n = SurfaceNormal(p + v*dist);
377 dist *= -(v.dot(n));
378 }
379 if (dist < eps) ++icount;
380 }
381 return dX*dY*dZ*icount/npoints/dd;
382}
static const G4double eps
double dot(const Hep3Vector &) const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4VSolid.cc:665
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
@ kInside
Definition: geomdefs.hh:70

References G4VSolid::BoundingLimits(), G4VSolid::DistanceToIn(), G4VSolid::DistanceToOut(), CLHEP::Hep3Vector::dot(), eps, G4QuickRand(), G4VSolid::Inside(), kInfinity, kInside, kOutside, G4INCL::Math::min(), CLHEP::detail::n, G4VSolid::SurfaceNormal(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by G4VSolid::GetSurfaceArea(), G4MultiUnion::GetSurfaceArea(), and G4VCSGfaceted::GetSurfaceArea().

◆ GetConstituentSolid() [1/2]

G4VSolid * G4VSolid::GetConstituentSolid ( G4int  no)
virtualinherited

Reimplemented in G4BooleanSolid.

Definition at line 170 of file G4VSolid.cc.

171{ return nullptr; }

◆ GetConstituentSolid() [2/2]

const G4VSolid * G4VSolid::GetConstituentSolid ( G4int  no) const
virtualinherited

Reimplemented in G4BooleanSolid.

Definition at line 167 of file G4VSolid.cc.

168{ return nullptr; }

Referenced by G4BooleanSolid::StackPolyhedron().

◆ GetCubicVolume()

G4double G4Ellipsoid::GetCubicVolume ( )
virtual

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}
static constexpr double pi
Definition: SystemOfUnits.h:55

References fCubicVolume, fDx, fDy, fDz, fZBottomCut, fZTopCut, and CLHEP::pi.

◆ GetDisplacedSolidPtr() [1/2]

G4DisplacedSolid * G4VSolid::GetDisplacedSolidPtr ( )
virtualinherited

Reimplemented in G4DisplacedSolid.

Definition at line 176 of file G4VSolid.cc.

177{ return nullptr; }

◆ GetDisplacedSolidPtr() [2/2]

const G4DisplacedSolid * G4VSolid::GetDisplacedSolidPtr ( ) const
virtualinherited

Reimplemented in G4DisplacedSolid.

Definition at line 173 of file G4VSolid.cc.

174{ return nullptr; }

◆ 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
virtual

Implements G4VSolid.

Definition at line 642 of file G4Ellipsoid.cc.

643{
644 return G4String("G4Ellipsoid");
645}

Referenced by StreamInfo().

◆ GetExtent()

G4VisExtent G4Ellipsoid::GetExtent ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 868 of file G4Ellipsoid.cc.

References fXmax, fYmax, fZBottomCut, and fZTopCut.

◆ GetName()

G4String G4VSolid::GetName ( ) const
inlineinherited

Referenced by G4GMocrenFileSceneHandler::AddDetector(), G4HepRepFileSceneHandler::AddHepRepInstance(), G4GMocrenFileSceneHandler::AddPrimitive(), G4HepRepFileSceneHandler::AddSolid(), G4GMocrenFileSceneHandler::AddSolid(), G4VtkSceneHandler::AddSolid(), G4GDMLWriteSolids::AddSolid(), G4NavigationLogger::AlongComputeStepLog(), G4GDMLWriteSolids::BooleanWrite(), G4ReflectedSolid::BoundingLimits(), G4DisplacedSolid::BoundingLimits(), G4IntersectionSolid::BoundingLimits(), G4ScaledSolid::BoundingLimits(), G4SubtractionSolid::BoundingLimits(), G4UnionSolid::BoundingLimits(), G4Box::BoundingLimits(), G4Cons::BoundingLimits(), G4CutTubs::BoundingLimits(), G4Orb::BoundingLimits(), G4Para::BoundingLimits(), G4Sphere::BoundingLimits(), G4Torus::BoundingLimits(), G4Trap::BoundingLimits(), G4Trd::BoundingLimits(), G4Tubs::BoundingLimits(), G4EllipticalCone::BoundingLimits(), G4ExtrudedSolid::BoundingLimits(), G4GenericPolycone::BoundingLimits(), G4GenericTrap::BoundingLimits(), G4Hype::BoundingLimits(), G4Paraboloid::BoundingLimits(), G4Polycone::BoundingLimits(), G4Polyhedra::BoundingLimits(), G4TessellatedSolid::BoundingLimits(), G4TwistedTubs::BoundingLimits(), G4GDMLWriteSolids::BoxWrite(), G4ExtrudedSolid::CalculateExtent(), G4GenericPolycone::CalculateExtent(), G4Polycone::CalculateExtent(), G4Polyhedra::CalculateExtent(), G4NavigationLogger::CheckDaughterEntryPoint(), G4VDivisionParameterisation::CheckNDivAndWidth(), G4VDivisionParameterisation::CheckOffset(), G4GenericTrap::CheckOrder(), G4Para::CheckParameters(), G4Trap::CheckParameters(), G4Trd::CheckParameters(), CheckParameters(), G4EllipticalTube::CheckParameters(), G4ParameterisationPolyconeRho::CheckParametersValidity(), G4ParameterisationPolyconeZ::CheckParametersValidity(), G4ParameterisationPolyhedraRho::CheckParametersValidity(), G4ParameterisationPolyhedraPhi::CheckParametersValidity(), G4ParameterisationPolyhedraZ::CheckParametersValidity(), G4PhantomParameterisation::CheckVoxelsFillContainer(), G4GenericTrap::ComputeIsTwisted(), G4VoxelNavigation::ComputeSafety(), G4VoxelSafety::ComputeSafety(), G4NavigationLogger::ComputeSafetyLog(), G4ParameterisedNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4GDMLWriteSolids::ConeWrite(), G4Polyhedra::Create(), G4GenericPolycone::Create(), G4Polycone::Create(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G4ReflectedSolid::CreatePolyhedron(), G4ReflectionFactory::CreateReflectedLV(), G4GenericTrap::CreateTessellatedSolid(), G4GDMLWriteSolids::CutTubeWrite(), G4SolidStore::DeRegister(), G4PhysicalVolumeModel::DescribeSolid(), G4SubtractionSolid::DistanceToIn(), G4Paraboloid::DistanceToIn(), G4TessellatedSolid::DistanceToIn(), G4Box::DistanceToOut(), G4Orb::DistanceToOut(), G4Para::DistanceToOut(), G4Trap::DistanceToOut(), G4Trd::DistanceToOut(), G4EllipticalCone::DistanceToOut(), G4TessellatedSolid::DistanceToOut(), DistanceToOut(), G4EllipticalTube::DistanceToOut(), G4tgbGeometryDumper::DumpMultiUnionVolume(), G4tgbGeometryDumper::DumpScaledVolume(), G4tgbGeometryDumper::DumpSolid(), G4GDMLWriteSolids::ElconeWrite(), G4GDMLWriteSolids::EllipsoidWrite(), G4GDMLWriteSolids::EltubeWrite(), G4PVDivision::ErrorInAxis(), G4ReplicatedSlice::ErrorInAxis(), export_G4VSolid(), G4Box::G4Box(), G4Cons::G4Cons(), G4CutTubs::G4CutTubs(), G4EllipticalCone::G4EllipticalCone(), G4Hype::G4Hype(), G4Para::G4Para(), G4Paraboloid::G4Paraboloid(), G4Polycone::G4Polycone(), G4Polyhedra::G4Polyhedra(), G4Sphere::G4Sphere(), G4Tet::G4Tet(), G4Trap::G4Trap(), G4Tubs::G4Tubs(), G4VParameterisationCons::G4VParameterisationCons(), G4VParameterisationPara::G4VParameterisationPara(), G4VParameterisationPolycone::G4VParameterisationPolycone(), G4VParameterisationPolyhedra::G4VParameterisationPolyhedra(), G4VParameterisationTrd::G4VParameterisationTrd(), G4VTwistedFaceted::G4VTwistedFaceted(), G4GDMLWriteSolids::GenericPolyconeWrite(), G4GDMLWriteSolids::GenTrapWrite(), G4Navigator::GetGlobalExitNormal(), G4Navigator::GetLocalExitNormal(), G4ITNavigator1::GetLocalExitNormal(), G4ITNavigator2::GetLocalExitNormal(), G4BooleanSolid::GetPointOnSurface(), G4PhantomParameterisation::GetReplicaNo(), G4GDMLWriteSolids::HypeWrite(), G4TessellatedSolid::InsideNoVoxels(), G4TessellatedSolid::InsideVoxels(), G4ITNavigator1::LocateGlobalPointAndSetup(), G4ITNavigator2::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointAndSetup(), G4GenericTrap::MakeDownFacet(), G4Trap::MakePlanes(), G4GenericTrap::MakeUpFacet(), G4GDMLWriteSolids::MultiUnionWrite(), G4GDMLWriteSolids::OrbWrite(), G4GDMLWriteSolids::ParaboloidWrite(), G4GDMLWriteParamvol::ParametersWrite(), G4GDMLWriteSolids::ParaWrite(), G4GDMLWriteSolids::PolyconeWrite(), G4GDMLWriteSolids::PolyhedraWrite(), G4NavigationLogger::PostComputeStepLog(), G4NavigationLogger::PreComputeStepLog(), G4NavigationLogger::PrintDaughterLog(), G4PseudoScene::ProcessVolume(), G4SolidStore::Register(), G4tgbVolumeMgr::RegisterMe(), G4NavigationLogger::ReportOutsideMother(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4VSceneHandler::RequestPrimitives(), G4GenericPolycone::Reset(), G4Polyhedra::Reset(), G4VoxelSafety::SafetyForVoxelNode(), G4GDMLWriteSolids::ScaledWrite(), G4Torus::SetAllParameters(), G4Tet::SetBoundingLimits(), G4Polycone::SetOriginalParameters(), G4Polyhedra::SetOriginalParameters(), G4TessellatedSolid::SetSolidClosed(), G4Tet::SetVertices(), G4Box::SetXHalfLength(), G4Box::SetYHalfLength(), G4Box::SetZHalfLength(), G4GDMLWriteSolids::SphereWrite(), G4BooleanSolid::StackPolyhedron(), G4ReflectedSolid::StreamInfo(), G4BooleanSolid::StreamInfo(), G4DisplacedSolid::StreamInfo(), G4MultiUnion::StreamInfo(), G4ScaledSolid::StreamInfo(), G4Box::StreamInfo(), G4Cons::StreamInfo(), G4CSGSolid::StreamInfo(), G4CutTubs::StreamInfo(), G4Orb::StreamInfo(), G4Para::StreamInfo(), G4Sphere::StreamInfo(), G4Torus::StreamInfo(), G4Trap::StreamInfo(), G4Trd::StreamInfo(), G4Tubs::StreamInfo(), StreamInfo(), G4EllipticalCone::StreamInfo(), G4EllipticalTube::StreamInfo(), G4ExtrudedSolid::StreamInfo(), G4GenericPolycone::StreamInfo(), G4GenericTrap::StreamInfo(), G4Hype::StreamInfo(), G4Paraboloid::StreamInfo(), G4Polycone::StreamInfo(), G4Polyhedra::StreamInfo(), G4TessellatedSolid::StreamInfo(), G4Tet::StreamInfo(), G4TwistedBox::StreamInfo(), G4TwistedTrap::StreamInfo(), G4TwistedTrd::StreamInfo(), G4TwistedTubs::StreamInfo(), G4VCSGfaceted::StreamInfo(), G4VTwistedFaceted::StreamInfo(), G4GDMLRead::StripNames(), SubstractSolids(), G4UnionSolid::SurfaceNormal(), G4Box::SurfaceNormal(), G4Para::SurfaceNormal(), G4Trap::SurfaceNormal(), G4Trd::SurfaceNormal(), SurfaceNormal(), G4EllipticalCone::SurfaceNormal(), G4EllipticalTube::SurfaceNormal(), G4ExtrudedSolid::SurfaceNormal(), G4Tet::SurfaceNormal(), G4GDMLWriteSolids::TessellatedWrite(), G4GDMLWriteSolids::TetWrite(), G4GDMLWriteSolids::TorusWrite(), G4GDMLWriteSolids::TrapWrite(), G4GDMLWriteStructure::TraverseVolumeTree(), G4GDMLWriteSolids::TrdWrite(), G4GDMLWriteSolids::TubeWrite(), G4GDMLWriteSolids::TwistedboxWrite(), G4GDMLWriteSolids::TwistedtrapWrite(), G4GDMLWriteSolids::TwistedtrdWrite(), G4GDMLWriteSolids::TwistedtubsWrite(), G4PhysicalVolumeModel::VisitGeometryAndGetVisReps(), and G4GDMLWriteSolids::XtruWrite().

◆ GetPointOnSurface()

G4ThreeVector G4Ellipsoid::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 781 of file G4Ellipsoid.cc.

782{
783 G4double A = GetDx();
784 G4double B = GetDy();
785 G4double C = GetDz();
786 G4double Zbot = GetZBottomCut();
787 G4double Ztop = GetZTopCut();
788
789 // Calculate cut areas
790 G4double Hbot = 1. + Zbot / C;
791 G4double Htop = 1. - Ztop / C;
792 G4double piAB = CLHEP::pi * A * B;
793 G4double Sbot = piAB * Hbot * (2. - Hbot);
794 G4double Stop = piAB * Htop * (2. - Htop);
795
796 // Get area of lateral surface
797 if (fLateralArea == 0.)
798 {
801 l.unlock();
802 }
803 G4double Slat = fLateralArea;
804
805 // Select surface (0 - bottom cut, 1 - lateral surface, 2 - top cut)
806 G4double select = (Sbot + Slat + Stop) * G4QuickRand();
807 G4int k = 0;
808 if (select > Sbot) k = 1;
809 if (select > Sbot + Slat) k = 2;
810
811 // Pick random point on selected surface (rejection sampling)
813 switch (k)
814 {
815 case 0: // bootom z-cut
816 {
817 G4double scale = std::sqrt(Hbot * (2. - Hbot));
818 G4TwoVector rho = G4RandomPointInEllipse(A * scale, B * scale);
819 p.set(rho.x(), rho.y(), Zbot);
820 break;
821 }
822 case 1: // lateral surface
823 {
824 G4double x, y, z;
825 G4double mu_max = std::max(std::max(A * B, A * C), B * C);
826 for (G4int i = 0; i < 1000; ++i)
827 {
828 // generate random point on unit sphere
829 z = (Zbot + (Ztop - Zbot) * G4QuickRand()) / C;
830 G4double rho = std::sqrt((1. + z) * (1. - z));
832 x = rho * std::cos(phi);
833 y = rho * std::sin(phi);
834 // check acceptance
835 G4double xbc = x * B * C;
836 G4double yac = y * A * C;
837 G4double zab = z * A * B;
838 G4double mu = std::sqrt(xbc * xbc + yac * yac + zab * zab);
839 if (mu_max * G4QuickRand() <= mu) break;
840 }
841 p.set(A * x, B * y, C * z);
842 break;
843 }
844 case 2: // top z-cut
845 {
846 G4double scale = std::sqrt(Htop * (2. - Htop));
847 G4TwoVector rho = G4RandomPointInEllipse(A * scale, B * scale);
848 p.set(rho.x(), rho.y(), Ztop);
849 break;
850 }
851 }
852 return p;
853}
G4TwoVector G4RandomPointInEllipse(G4double a, G4double b)
double x() const
double y() const
void set(double x, double y, double z)
G4double GetDx() const
G4double LateralSurfaceArea() const
Definition: G4Ellipsoid.cc:706
G4double GetDy() const
G4double GetZTopCut() const
G4double GetZBottomCut() const
G4double GetDz() const
static constexpr double twopi
Definition: SystemOfUnits.h:56

References A, B(), C(), fLateralArea, G4QuickRand(), G4RandomPointInEllipse(), GetDx(), GetDy(), GetDz(), GetZBottomCut(), GetZTopCut(), anonymous_namespace{G4Ellipsoid.cc}::lateralareaMutex, LateralSurfaceArea(), G4INCL::Math::max(), CLHEP::pi, CLHEP::Hep3Vector::set(), CLHEP::twopi, G4TemplateAutoLock< _Mutex_t >::unlock(), CLHEP::Hep2Vector::x(), and CLHEP::Hep2Vector::y().

◆ GetPolyhedron()

G4Polyhedron * G4Ellipsoid::GetPolyhedron ( ) const
virtual

◆ GetSemiAxisMax()

G4double G4Ellipsoid::GetSemiAxisMax ( G4int  i) const
inline

◆ GetSurfaceArea()

G4double G4Ellipsoid::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 757 of file G4Ellipsoid.cc.

758{
759 if (fSurfaceArea == 0.)
760 {
761 G4double piAB = CLHEP::pi * fDx * fDy;
763 if (fZBottomCut > -fDz)
764 {
765 G4double hbot = 1. + fZBottomCut / fDz;
766 fSurfaceArea += piAB * hbot * (2. - hbot);
767 }
768 if (fZTopCut < fDz)
769 {
770 G4double htop = 1. - fZTopCut / fDz;
771 fSurfaceArea += piAB * htop * (2. - htop);
772 }
773 }
774 return fSurfaceArea;
775}

References fDx, fDy, fDz, fSurfaceArea, fZBottomCut, fZTopCut, LateralSurfaceArea(), and CLHEP::pi.

◆ GetTolerance()

G4double G4VSolid::GetTolerance ( ) const
inlineinherited

◆ GetZBottomCut()

G4double G4Ellipsoid::GetZBottomCut ( ) const
inline

◆ GetZTopCut()

G4double G4Ellipsoid::GetZTopCut ( ) const
inline

◆ Inside()

EInside G4Ellipsoid::Inside ( const G4ThreeVector p) const
virtual

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}
@ kSurface
Definition: geomdefs.hh:69

References fQ1, fQ2, fSx, fSy, fSz, fZDimCut, fZMidCut, halfTolerance, kInside, kOutside, kSurface, G4INCL::Math::max(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ LateralSurfaceArea()

G4double G4Ellipsoid::LateralSurfaceArea ( ) const
private

Definition at line 706 of file G4Ellipsoid.cc.

707{
708 const G4int Nphi = 100;
709 const G4int Nz = 200;
710 G4double rho[Nz + 1];
711
712 // Set array of rho
713 G4double zbot = fZBottomCut / fDz;
714 G4double ztop = fZTopCut / fDz;
715 G4double dz = (ztop - zbot) / Nz;
716 for (G4int iz = 0; iz < Nz; ++iz)
717 {
718 G4double z = zbot + iz * dz;
719 rho[iz] = std::sqrt((1. + z) * (1. - z));
720 }
721 rho[Nz] = std::sqrt((1. + ztop) * (1. - ztop));
722
723 // Compute area
724 zbot = fZBottomCut;
725 ztop = fZTopCut;
726 dz = (ztop - zbot) / Nz;
727 G4double area = 0.;
728 G4double dphi = CLHEP::halfpi / Nphi;
729 for (G4int iphi = 0; iphi < Nphi; ++iphi)
730 {
731 G4double phi1 = iphi * dphi;
732 G4double phi2 = (iphi == Nphi - 1) ? CLHEP::halfpi : phi1 + dphi;
733 G4double cos1 = std::cos(phi1) * fDx;
734 G4double cos2 = std::cos(phi2) * fDx;
735 G4double sin1 = std::sin(phi1) * fDy;
736 G4double sin2 = std::sin(phi2) * fDy;
737 for (G4int iz = 0; iz < Nz; ++iz)
738 {
739 G4double z1 = zbot + iz * dz;
740 G4double z2 = (iz == Nz - 1) ? ztop : z1 + dz;
741 G4double rho1 = rho[iz];
742 G4double rho2 = rho[iz + 1];
743 G4ThreeVector p1(rho1 * cos1, rho1 * sin1, z1);
744 G4ThreeVector p2(rho1 * cos2, rho1 * sin2, z1);
745 G4ThreeVector p3(rho2 * cos1, rho2 * sin1, z2);
746 G4ThreeVector p4(rho2 * cos2, rho2 * sin2, z2);
747 area += ((p4 - p1).cross(p3 - p2)).mag();
748 }
749 }
750 return 2. * area;
751}
static constexpr double halfpi
Definition: SystemOfUnits.h:57

References fDx, fDy, fDz, fZBottomCut, fZTopCut, and CLHEP::halfpi.

Referenced by GetPointOnSurface(), and GetSurfaceArea().

◆ 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;
138 fZTopCut = rhs.fZTopCut;
139
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
156
157 fRebuildPolyhedron = false;
158 delete fpPolyhedron; fpPolyhedron = nullptr;
159
160 return *this;
161}
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:107

References fCubicVolume, fDx, fDy, fDz, fLateralArea, fpPolyhedron, fQ1, fQ2, fR, fRebuildPolyhedron, fRsph, fSurfaceArea, fSx, fSy, fSz, fXmax, fYmax, fZBottomCut, fZDimCut, fZMidCut, fZTopCut, halfTolerance, and G4VSolid::operator=().

◆ operator==()

G4bool G4VSolid::operator== ( const G4VSolid s) const
inlineinherited

◆ SetName()

void G4VSolid::SetName ( const G4String name)
inherited

◆ 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
virtual

Implements G4VSolid.

Definition at line 660 of file G4Ellipsoid.cc.

661{
662 G4int 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}
static constexpr double mm
Definition: G4SIunits.hh:95
G4GeometryType GetEntityType() const
Definition: G4Ellipsoid.cc:642

References GetDx(), GetDy(), GetDz(), GetEntityType(), G4VSolid::GetName(), GetZBottomCut(), GetZTopCut(), and mm.

◆ SurfaceNormal()

G4ThreeVector G4Ellipsoid::SurfaceNormal ( const G4ThreeVector p) const
virtual

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 G4int 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}

References ApproxSurfaceNormal(), G4VSolid::DumpInfo(), fQ1, fQ2, fSx, fSy, fSz, fZDimCut, fZMidCut, G4cout, G4Exception(), G4VSolid::GetName(), halfTolerance, JustWarning, mm, CLHEP::Hep3Vector::setZ(), CLHEP::Hep3Vector::unit(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Field Documentation

◆ fCubicVolume

G4double G4Ellipsoid::fCubicVolume = 0.0
private

Definition at line 177 of file G4Ellipsoid.hh.

Referenced by CheckParameters(), GetCubicVolume(), and operator=().

◆ fDx

G4double G4Ellipsoid::fDx
private

◆ fDy

G4double G4Ellipsoid::fDy
private

◆ fDz

G4double G4Ellipsoid::fDz
private

◆ fLateralArea

G4double G4Ellipsoid::fLateralArea = 0.0
mutableprivate

Definition at line 179 of file G4Ellipsoid.hh.

Referenced by CheckParameters(), GetPointOnSurface(), and operator=().

◆ fpPolyhedron

G4Polyhedron* G4Ellipsoid::fpPolyhedron = nullptr
mutableprivate

Definition at line 181 of file G4Ellipsoid.hh.

Referenced by GetPolyhedron(), operator=(), and ~G4Ellipsoid().

◆ fQ1

G4double G4Ellipsoid::fQ1
private

◆ fQ2

G4double G4Ellipsoid::fQ2
private

◆ fR

G4double G4Ellipsoid::fR
private

◆ fRebuildPolyhedron

G4bool G4Ellipsoid::fRebuildPolyhedron = false
mutableprivate

Definition at line 180 of file G4Ellipsoid.hh.

Referenced by GetPolyhedron(), and operator=().

◆ fRsph

G4double G4Ellipsoid::fRsph
private

Definition at line 167 of file G4Ellipsoid.hh.

Referenced by CheckParameters(), DistanceToIn(), and operator=().

◆ fshapeName

G4String G4VSolid::fshapeName
privateinherited

Definition at line 312 of file G4VSolid.hh.

Referenced by G4VSolid::operator=(), and G4VSolid::SetName().

◆ fSurfaceArea

G4double G4Ellipsoid::fSurfaceArea = 0.0
private

Definition at line 178 of file G4Ellipsoid.hh.

Referenced by CheckParameters(), GetSurfaceArea(), and operator=().

◆ fSx

G4double G4Ellipsoid::fSx
private

◆ fSy

G4double G4Ellipsoid::fSy
private

◆ fSz

G4double G4Ellipsoid::fSz
private

◆ fXmax

G4double G4Ellipsoid::fXmax
private

Definition at line 165 of file G4Ellipsoid.hh.

Referenced by BoundingLimits(), CheckParameters(), DistanceToIn(), GetExtent(), and operator=().

◆ fYmax

G4double G4Ellipsoid::fYmax
private

Definition at line 166 of file G4Ellipsoid.hh.

Referenced by BoundingLimits(), CheckParameters(), DistanceToIn(), GetExtent(), and operator=().

◆ fZBottomCut

G4double G4Ellipsoid::fZBottomCut
private

◆ fZDimCut

G4double G4Ellipsoid::fZDimCut
private

◆ fZMidCut

G4double G4Ellipsoid::fZMidCut
private

◆ fZTopCut

G4double G4Ellipsoid::fZTopCut
private

◆ halfTolerance

G4double G4Ellipsoid::halfTolerance
private

◆ kCarTolerance

G4double G4VSolid::kCarTolerance
protectedinherited

Definition at line 299 of file G4VSolid.hh.

Referenced by G4TessellatedSolid::AddFacet(), G4Polycone::CalculateExtent(), G4Polyhedra::CalculateExtent(), G4Tet::CheckDegeneracy(), G4Para::CheckParameters(), G4Trd::CheckParameters(), CheckParameters(), G4EllipticalTube::CheckParameters(), G4GenericTrap::ComputeIsTwisted(), G4Polyhedra::Create(), G4GenericPolycone::Create(), G4Polycone::Create(), G4CutTubs::CreatePolyhedron(), G4TessellatedSolid::CreateVertexList(), G4VCSGfaceted::DistanceTo(), G4Sphere::DistanceToIn(), DistanceToIn(), G4Hype::DistanceToIn(), G4Paraboloid::DistanceToIn(), G4VCSGfaceted::DistanceToIn(), G4TessellatedSolid::DistanceToInCore(), G4Cons::DistanceToOut(), G4CutTubs::DistanceToOut(), G4Sphere::DistanceToOut(), G4Torus::DistanceToOut(), G4Tubs::DistanceToOut(), G4GenericTrap::DistanceToOut(), G4Hype::DistanceToOut(), G4Paraboloid::DistanceToOut(), G4VCSGfaceted::DistanceToOut(), G4TessellatedSolid::DistanceToOutCandidates(), G4TessellatedSolid::DistanceToOutCore(), G4TessellatedSolid::DistanceToOutNoVoxels(), G4GenericTrap::DistToPlane(), G4GenericTrap::DistToTriangle(), G4Box::G4Box(), G4Cons::G4Cons(), G4CutTubs::G4CutTubs(), G4EllipticalCone::G4EllipticalCone(), G4ExtrudedSolid::G4ExtrudedSolid(), G4GenericTrap::G4GenericTrap(), G4Hype::G4Hype(), G4Para::G4Para(), G4Sphere::G4Sphere(), G4Tet::G4Tet(), G4Trap::G4Trap(), G4Tubs::G4Tubs(), G4UnionSolid::G4UnionSolid(), G4VSolid::G4VSolid(), G4VTwistedFaceted::G4VTwistedFaceted(), G4GenericPolycone::GetPointOnSurface(), G4Polycone::GetPointOnSurface(), G4UnionSolid::Init(), G4Orb::Initialize(), G4TessellatedSolid::Initialize(), G4SubtractionSolid::Inside(), G4Hype::Inside(), G4Paraboloid::Inside(), G4VCSGfaceted::Inside(), G4VTwistedFaceted::Inside(), G4TessellatedSolid::InsideNoVoxels(), G4GenericTrap::InsidePolygone(), G4TessellatedSolid::InsideVoxels(), G4CutTubs::IsCrossingCutPlanes(), G4GenericTrap::IsSegCrossingZ(), G4Trap::MakePlane(), G4GenericTrap::NormalToPlane(), G4VSolid::operator=(), G4TessellatedSolid::SafetyFromInside(), G4TessellatedSolid::SafetyFromOutside(), G4Torus::SetAllParameters(), G4Polycone::SetOriginalParameters(), G4Polyhedra::SetOriginalParameters(), G4Box::SetXHalfLength(), G4Box::SetYHalfLength(), G4Box::SetZHalfLength(), G4Torus::SurfaceNormal(), G4GenericTrap::SurfaceNormal(), and G4Paraboloid::SurfaceNormal().


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