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

#include <G4Tet.hh>

Inheritance diagram for G4Tet:
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
 
G4bool CheckDegeneracy (const G4ThreeVector &p0, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3) 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
 
 G4Tet (__void__ &)
 
 G4Tet (const G4String &pName, const G4ThreeVector &anchor, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3, G4bool *degeneracyFlag=nullptr)
 
 G4Tet (const G4Tet &rhs)
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
G4double GetCubicVolume ()
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
G4GeometryType GetEntityType () const
 
G4VisExtent GetExtent () const
 
G4String GetName () const
 
G4ThreeVector GetPointOnSurface () const
 
G4PolyhedronGetPolyhedron () const
 
G4double GetSurfaceArea ()
 
G4double GetTolerance () const
 
std::vector< G4ThreeVectorGetVertices () const
 
void GetVertices (G4ThreeVector &anchor, G4ThreeVector &p1, G4ThreeVector &p2, G4ThreeVector &p3) const
 
EInside Inside (const G4ThreeVector &p) const
 
G4Tetoperator= (const G4Tet &rhs)
 
G4bool operator== (const G4VSolid &s) const
 
void PrintWarnings (G4bool)
 
void SetBoundingLimits (const G4ThreeVector &pMin, const G4ThreeVector &pMax)
 
void SetName (const G4String &name)
 
void SetVertices (const G4ThreeVector &anchor, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3, G4bool *degeneracyFlag=nullptr)
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
virtual ~G4Tet ()
 

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 ClipPolygonToSimpleLimits (G4ThreeVectorList &pPolygon, G4ThreeVectorList &outputPolygon, const G4VoxelLimits &pVoxelLimit) const
 
void Initialize (const G4ThreeVector &p0, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3)
 

Private Attributes

G4double fArea [4] = {0}
 
G4ThreeVector fBmax
 
G4ThreeVector fBmin
 
G4double fCubicVolume = 0
 
G4double fDist [4] = {0}
 
G4ThreeVector fNormal [4]
 
G4PolyhedronfpPolyhedron = nullptr
 
G4bool fRebuildPolyhedron = false
 
G4String fshapeName
 
G4double fSurfaceArea = 0
 
G4ThreeVector fVertex [4]
 
G4double halfTolerance = 0
 

Detailed Description

Definition at line 55 of file G4Tet.hh.

Constructor & Destructor Documentation

◆ G4Tet() [1/3]

G4Tet::G4Tet ( const G4String pName,
const G4ThreeVector anchor,
const G4ThreeVector p1,
const G4ThreeVector p2,
const G4ThreeVector p3,
G4bool degeneracyFlag = nullptr 
)

Definition at line 66 of file G4Tet.cc.

71 : G4VSolid(pName)
72{
73 // Check for degeneracy
74 G4bool degenerate = CheckDegeneracy(p0, p1, p2, p3);
75 if (degeneracyFlag)
76 {
77 *degeneracyFlag = degenerate;
78 }
79 else if (degenerate)
80 {
81 std::ostringstream message;
82 message << "Degenerate tetrahedron: " << GetName() << " !\n"
83 << " anchor: " << p0 << "\n"
84 << " p1 : " << p1 << "\n"
85 << " p2 : " << p2 << "\n"
86 << " p3 : " << p3 << "\n"
87 << " volume: "
88 << std::abs((p1 - p0).cross(p2 - p0).dot(p3 - p0))/6.;
89 G4Exception("G4Tet::G4Tet()", "GeomSolids0002", FatalException, message);
90 }
91
92 // Define surface thickness
94
95 // Set data members
96 Initialize(p0, p1, p2, p3);
97}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
bool G4bool
Definition: G4Types.hh:86
G4double halfTolerance
Definition: G4Tet.hh:161
void Initialize(const G4ThreeVector &p0, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3)
Definition: G4Tet.cc:202
G4bool CheckDegeneracy(const G4ThreeVector &p0, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3) const
Definition: G4Tet.cc:173
G4String GetName() const
G4VSolid(const G4String &name)
Definition: G4VSolid.cc:57
G4double kCarTolerance
Definition: G4VSolid.hh:299

References CheckDegeneracy(), FatalException, G4Exception(), G4VSolid::GetName(), halfTolerance, Initialize(), and G4VSolid::kCarTolerance.

Referenced by Clone().

◆ ~G4Tet()

G4Tet::~G4Tet ( )
virtual

Definition at line 113 of file G4Tet.cc.

114{
115 delete fpPolyhedron; fpPolyhedron = nullptr;
116}
G4Polyhedron * fpPolyhedron
Definition: G4Tet.hh:165

References fpPolyhedron.

◆ G4Tet() [2/3]

G4Tet::G4Tet ( __void__ &  a)

Definition at line 104 of file G4Tet.cc.

105 : G4VSolid(a)
106{
107}

◆ G4Tet() [3/3]

G4Tet::G4Tet ( const G4Tet rhs)

Definition at line 122 of file G4Tet.cc.

123 : G4VSolid(rhs)
124{
128 for (G4int i = 0; i < 4; ++i) { fVertex[i] = rhs.fVertex[i]; }
129 for (G4int i = 0; i < 4; ++i) { fNormal[i] = rhs.fNormal[i]; }
130 for (G4int i = 0; i < 4; ++i) { fDist[i] = rhs.fDist[i]; }
131 for (G4int i = 0; i < 4; ++i) { fArea[i] = rhs.fArea[i]; }
132 fBmin = rhs.fBmin;
133 fBmax = rhs.fBmax;
134}
int G4int
Definition: G4Types.hh:85
G4ThreeVector fVertex[4]
Definition: G4Tet.hh:167
G4double fCubicVolume
Definition: G4Tet.hh:162
G4ThreeVector fBmax
Definition: G4Tet.hh:171
G4ThreeVector fNormal[4]
Definition: G4Tet.hh:168
G4ThreeVector fBmin
Definition: G4Tet.hh:171
G4double fArea[4]
Definition: G4Tet.hh:170
G4double fDist[4]
Definition: G4Tet.hh:169
G4double fSurfaceArea
Definition: G4Tet.hh:163

References fArea, fBmax, fBmin, fCubicVolume, fDist, fNormal, fSurfaceArea, fVertex, and halfTolerance.

Member Function Documentation

◆ ApproxSurfaceNormal()

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

Definition at line 471 of file G4Tet.cc.

472{
473 G4double dist = -DBL_MAX;
474 G4int iside = 0;
475 for (G4int i = 0; i < 4; ++i)
476 {
477 G4double d = fNormal[i].dot(p) - fDist[i];
478 if (d > dist) { dist = d; iside = i; }
479 }
480 return fNormal[iside];
481}
double G4double
Definition: G4Types.hh:83
double dot(const Hep3Vector &) const
#define DBL_MAX
Definition: templates.hh:62

References DBL_MAX, CLHEP::Hep3Vector::dot(), fDist, and fNormal.

Referenced by SurfaceNormal().

◆ BoundingLimits()

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

Reimplemented from G4VSolid.

Definition at line 360 of file G4Tet.cc.

361{
362 pMin = fBmin;
363 pMax = fBmax;
364}
static const G4double pMax
static const G4double pMin

References fBmax, fBmin, 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}
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 G4Tet::CalculateExtent ( const EAxis  pAxis,
const G4VoxelLimits pVoxelLimit,
const G4AffineTransform pTransform,
G4double pmin,
G4double pmax 
) const
virtual

Implements G4VSolid.

Definition at line 370 of file G4Tet.cc.

374{
375 G4ThreeVector bmin, bmax;
376
377 // Check bounding box (bbox)
378 //
379 BoundingLimits(bmin,bmax);
380 G4BoundingEnvelope bbox(bmin,bmax);
381
382 // Use simple bounding-box to help in the case of complex 3D meshes
383 //
384 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
385
386#if 0
387 // Precise extent computation (disabled by default for this shape)
388 //
389 G4bool exist;
390 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
391 {
392 return exist = (pMin < pMax) ? true : false;
393 }
394
395 // Set bounding envelope (benv) and calculate extent
396 //
397 std::vector<G4ThreeVector> vec = GetVertices();
398
399 G4ThreeVectorList anchor(1);
400 anchor[0].set(vec[0].x(),vec[0].y(),vec[0].z());
401
402 G4ThreeVectorList base(3);
403 base[0].set(vec[1].x(),vec[1].y(),vec[1].z());
404 base[1].set(vec[2].x(),vec[2].y(),vec[2].z());
405 base[2].set(vec[3].x(),vec[3].y(),vec[3].z());
406
407 std::vector<const G4ThreeVectorList *> polygons(2);
408 polygons[0] = &anchor;
409 polygons[1] = &base;
410
411 G4BoundingEnvelope benv(bmin,bmax,polygons);
412 return exists = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
413#endif
414}
std::vector< G4ThreeVector > G4ThreeVectorList
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Tet.cc:360
std::vector< G4ThreeVector > GetVertices() const
Definition: G4Tet.cc:301

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

◆ CheckDegeneracy()

G4bool G4Tet::CheckDegeneracy ( const G4ThreeVector p0,
const G4ThreeVector p1,
const G4ThreeVector p2,
const G4ThreeVector p3 
) const

Definition at line 173 of file G4Tet.cc.

177{
178 G4double hmin = 4. * kCarTolerance; // degeneracy tolerance
179
180 // Calculate volume
181 G4double vol = std::abs((p1 - p0).cross(p2 - p0).dot(p3 - p0));
182
183 // Calculate face areas squared
184 G4double ss[4];
185 ss[0] = ((p1 - p0).cross(p2 - p0)).mag2();
186 ss[1] = ((p2 - p0).cross(p3 - p0)).mag2();
187 ss[2] = ((p3 - p0).cross(p1 - p0)).mag2();
188 ss[3] = ((p2 - p1).cross(p3 - p1)).mag2();
189
190 // Find face with max area
191 G4int k = 0;
192 for (G4int i = 1; i < 4; ++i) { if (ss[i] > ss[k]) k = i; }
193
194 // Check: vol^2 / s^2 <= hmin^2
195 return (vol*vol <= ss[k]*hmin*hmin);
196}

References G4VSolid::kCarTolerance.

Referenced by G4Tet(), and SetVertices().

◆ 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}
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 * G4Tet::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 594 of file G4Tet.cc.

595{
596 return new G4Tet(*this);
597}
G4Tet(const G4String &pName, const G4ThreeVector &anchor, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3, G4bool *degeneracyFlag=nullptr)
Definition: G4Tet.cc:66

References G4Tet().

◆ ComputeDimensions()

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

Reimplemented from G4VSolid.

Definition at line 313 of file G4Tet.cc.

316{
317}

◆ CreatePolyhedron()

G4Polyhedron * G4Tet::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 687 of file G4Tet.cc.

688{
689 // Check orientation of vertices
690 G4ThreeVector v1 = fVertex[1] - fVertex[0];
691 G4ThreeVector v2 = fVertex[2] - fVertex[0];
692 G4ThreeVector v3 = fVertex[3] - fVertex[0];
693 G4bool invert = v1.cross(v2).dot(v3) < 0.;
694 G4int k2 = (invert) ? 3 : 2;
695 G4int k3 = (invert) ? 2 : 3;
696
697 // Set coordinates of vertices
698 G4double xyz[4][3];
699 for (G4int i = 0; i < 3; ++i)
700 {
701 xyz[0][i] = fVertex[0][i];
702 xyz[1][i] = fVertex[1][i];
703 xyz[2][i] = fVertex[k2][i];
704 xyz[3][i] = fVertex[k3][i];
705 }
706
707 // Create polyhedron
708 G4int faces[4][4] = { {1,3,2,0}, {1,4,3,0}, {1,2,4,0}, {2,3,4,0} };
709 G4Polyhedron* ph = new G4Polyhedron;
710 ph->createPolyhedron(4,4,xyz,faces);
711
712 return ph;
713}
Hep3Vector cross(const Hep3Vector &) const
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])

References HepPolyhedron::createPolyhedron(), CLHEP::Hep3Vector::cross(), CLHEP::Hep3Vector::dot(), and fVertex.

Referenced by G4ArrowModel::G4ArrowModel(), and GetPolyhedron().

◆ DescribeYourselfTo()

void G4Tet::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 667 of file G4Tet.cc.

668{
669 scene.AddSolid (*this);
670}
virtual void AddSolid(const G4Box &)=0

References G4VGraphicsScene::AddSolid().

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 515 of file G4Tet.cc.

516{
517 G4double dd[4];
518 for (G4int i = 0; i < 4; ++i) { dd[i] = fNormal[i].dot(p) - fDist[i]; }
519
520 G4double dist = std::max(std::max(std::max(dd[0], dd[1]), dd[2]), dd[3]);
521 return (dist > 0.) ? dist : 0.;
522}
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References CLHEP::Hep3Vector::dot(), fDist, fNormal, and G4INCL::Math::max().

◆ DistanceToIn() [2/2]

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

Implements G4VSolid.

Definition at line 488 of file G4Tet.cc.

490{
491 G4double tin = -DBL_MAX, tout = DBL_MAX;
492 for (G4int i = 0; i < 4; ++i)
493 {
494 G4double cosa = fNormal[i].dot(v);
495 G4double dist = fNormal[i].dot(p) - fDist[i];
496 if (dist >= -halfTolerance)
497 {
498 if (cosa >= 0.) { return kInfinity; }
499 tin = std::max(tin, -dist/cosa);
500 }
501 else if (cosa > 0.)
502 {
503 tout = std::min(tout, -dist/cosa);
504 }
505 }
506
507 return (tout - tin <= halfTolerance) ?
508 kInfinity : ((tin < halfTolerance) ? 0. : tin);
509}
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References DBL_MAX, CLHEP::Hep3Vector::dot(), fDist, fNormal, halfTolerance, kInfinity, G4INCL::Math::max(), and G4INCL::Math::min().

◆ DistanceToOut() [1/2]

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

Implements G4VSolid.

Definition at line 572 of file G4Tet.cc.

573{
574 G4double dd[4];
575 for (G4int i = 0; i < 4; ++i) { dd[i] = fDist[i] - fNormal[i].dot(p); }
576
577 G4double dist = std::min(std::min(std::min(dd[0], dd[1]), dd[2]), dd[3]);
578 return (dist > 0.) ? dist : 0.;
579}

References CLHEP::Hep3Vector::dot(), fDist, fNormal, and G4INCL::Math::min().

◆ DistanceToOut() [2/2]

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

Implements G4VSolid.

Definition at line 528 of file G4Tet.cc.

533{
534 // Calculate distances and cosines
535 G4double cosa[4], dist[4];
536 G4int ind[4] = {0}, nside = 0;
537 for (G4int i = 0; i < 4; ++i)
538 {
539 G4double tmp = fNormal[i].dot(v);
540 cosa[i] = tmp;
541 ind[nside] = (tmp > 0) * i;
542 nside += (tmp > 0);
543 dist[i] = fNormal[i].dot(p) - fDist[i];
544 }
545
546 // Find intersection (in most of cases nside == 1)
547 G4double tout = DBL_MAX;
548 G4int iside = 0;
549 for (G4int i = 0; i < nside; ++i)
550 {
551 G4int k = ind[i];
552 // Check: leaving the surface
553 if (dist[k] >= -halfTolerance) { tout = 0.; iside = k; break; }
554 // Compute distance to intersection
555 G4double tmp = -dist[k]/cosa[k];
556 if (tmp < tout) { tout = tmp; iside = k; }
557 }
558
559 // Set normal, if required, and return distance to out
560 if (calcNorm)
561 {
562 *validNorm = true;
563 *n = fNormal[iside];
564 }
565 return tout;
566}

References DBL_MAX, CLHEP::Hep3Vector::dot(), fDist, fNormal, halfTolerance, and CLHEP::detail::n.

◆ 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(), G4Ellipsoid::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(), G4Ellipsoid::SurfaceNormal(), G4EllipticalCone::SurfaceNormal(), G4EllipticalTube::SurfaceNormal(), G4ExtrudedSolid::SurfaceNormal(), and 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
CLHEP::Hep3Vector G4ThreeVector
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 z() const
double x() const
double y() 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 G4Tet::GetCubicVolume ( )
virtual

Reimplemented from G4VSolid.

Definition at line 649 of file G4Tet.cc.

650{
651 return fCubicVolume;
652}

References fCubicVolume.

◆ 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; }

◆ GetEntityType()

G4GeometryType G4Tet::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 585 of file G4Tet.cc.

586{
587 return G4String("G4Tet");
588}

Referenced by StreamInfo().

◆ GetExtent()

G4VisExtent G4Tet::GetExtent ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 676 of file G4Tet.cc.

677{
678 return G4VisExtent(fBmin.x(), fBmax.x(),
679 fBmin.y(), fBmax.y(),
680 fBmin.z(), fBmax.z());
681}

References fBmax, fBmin, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ 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(), G4Ellipsoid::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(), G4Ellipsoid::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(), 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(), SetBoundingLimits(), G4Polycone::SetOriginalParameters(), G4Polyhedra::SetOriginalParameters(), G4TessellatedSolid::SetSolidClosed(), 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(), G4Ellipsoid::StreamInfo(), G4EllipticalCone::StreamInfo(), G4EllipticalTube::StreamInfo(), G4ExtrudedSolid::StreamInfo(), G4GenericPolycone::StreamInfo(), G4GenericTrap::StreamInfo(), G4Hype::StreamInfo(), G4Paraboloid::StreamInfo(), G4Polycone::StreamInfo(), G4Polyhedra::StreamInfo(), G4TessellatedSolid::StreamInfo(), 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(), G4Ellipsoid::SurfaceNormal(), G4EllipticalCone::SurfaceNormal(), G4EllipticalTube::SurfaceNormal(), G4ExtrudedSolid::SurfaceNormal(), 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 G4Tet::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 624 of file G4Tet.cc.

625{
626 constexpr G4int iface[4][3] = { {0,1,2}, {0,2,3}, {0,3,1}, {1,2,3} };
627
628 // Select face
630 G4int i = 0;
631 for ( ; i < 4; ++i) { if ((select -= fArea[i]) <= 0.) break; }
632
633 // Set selected triangle
634 G4ThreeVector p0 = fVertex[iface[i][0]];
635 G4ThreeVector e1 = fVertex[iface[i][1]] - p0;
636 G4ThreeVector e2 = fVertex[iface[i][2]] - p0;
637
638 // Return random point
639 G4double r1 = G4QuickRand();
640 G4double r2 = G4QuickRand();
641 return (r1 + r2 > 1.) ?
642 p0 + e1*(1. - r1) + e2*(1. - r2) : p0 + e1*r1 + e2*r2;
643}
static const G4double e1[44]
static const G4double e2[44]

References e1, e2, fArea, fSurfaceArea, fVertex, and G4QuickRand().

◆ GetPolyhedron()

G4Polyhedron * G4Tet::GetPolyhedron ( ) const
virtual

◆ GetSurfaceArea()

G4double G4Tet::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 658 of file G4Tet.cc.

659{
660 return fSurfaceArea;
661}

References fSurfaceArea.

◆ GetTolerance()

G4double G4VSolid::GetTolerance ( ) const
inlineinherited

◆ GetVertices() [1/2]

std::vector< G4ThreeVector > G4Tet::GetVertices ( ) const

Definition at line 301 of file G4Tet.cc.

302{
303 std::vector<G4ThreeVector> vertices(4);
304 for (G4int i = 0; i < 4; ++i) { vertices[i] = fVertex[i]; }
305 return vertices;
306}

References fVertex.

Referenced by CalculateExtent().

◆ GetVertices() [2/2]

void G4Tet::GetVertices ( G4ThreeVector anchor,
G4ThreeVector p1,
G4ThreeVector p2,
G4ThreeVector p3 
) const

Definition at line 286 of file G4Tet.cc.

290{
291 p0 = fVertex[0];
292 p1 = fVertex[1];
293 p2 = fVertex[2];
294 p3 = fVertex[3];
295}

References fVertex.

Referenced by G4GDMLWriteSolids::TetWrite().

◆ Initialize()

void G4Tet::Initialize ( const G4ThreeVector p0,
const G4ThreeVector p1,
const G4ThreeVector p2,
const G4ThreeVector p3 
)
private

Definition at line 202 of file G4Tet.cc.

206{
207 // Set vertices
208 fVertex[0] = p0;
209 fVertex[1] = p1;
210 fVertex[2] = p2;
211 fVertex[3] = p3;
212
213 G4ThreeVector norm[4];
214 norm[0] = (p2 - p0).cross(p1 - p0);
215 norm[1] = (p3 - p0).cross(p2 - p0);
216 norm[2] = (p1 - p0).cross(p3 - p0);
217 norm[3] = (p2 - p1).cross(p3 - p1);
218 G4double volume = norm[0].dot(p3 - p0);
219 if (volume > 0.)
220 {
221 for (G4int i = 0; i < 4; ++i) { norm[i] = -norm[i]; }
222 }
223
224 // Set normals to face planes
225 for (G4int i = 0; i < 4; ++i) { fNormal[i] = norm[i].unit(); }
226
227 // Set distances to planes
228 for (G4int i = 0; i < 3; ++i) { fDist[i] = fNormal[i].dot(p0); }
229 fDist[3] = fNormal[3].dot(p1);
230
231 // Set face areas
232 for (G4int i = 0; i < 4; ++i) { fArea[i] = 0.5*norm[i].mag(); }
233
234 // Set bounding box
235 for (G4int i = 0; i < 3; ++i)
236 {
237 fBmin[i] = std::min(std::min(std::min(p0[i], p1[i]), p2[i]), p3[i]);
238 fBmax[i] = std::max(std::max(std::max(p0[i], p1[i]), p2[i]), p3[i]);
239 }
240
241 // Set volume and surface area
242 fCubicVolume = std::abs(volume)/6.;
243 fSurfaceArea = fArea[0] + fArea[1] + fArea[2] + fArea[3];
244}
Hep3Vector unit() const
double mag() const

References CLHEP::Hep3Vector::dot(), fArea, fBmax, fBmin, fCubicVolume, fDist, fNormal, fSurfaceArea, fVertex, CLHEP::Hep3Vector::mag(), G4INCL::Math::max(), G4INCL::Math::min(), and CLHEP::Hep3Vector::unit().

Referenced by G4Tet(), and SetVertices().

◆ Inside()

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

Implements G4VSolid.

Definition at line 420 of file G4Tet.cc.

421{
422 G4double dd[4];
423 for (G4int i = 0; i < 4; ++i) { dd[i] = fNormal[i].dot(p) - fDist[i]; }
424
425 G4double dist = std::max(std::max(std::max(dd[0], dd[1]), dd[2]), dd[3]);
426 return (dist <= -halfTolerance) ?
427 kInside : ((dist <= halfTolerance) ? kSurface : kOutside);
428}
@ kSurface
Definition: geomdefs.hh:69

References CLHEP::Hep3Vector::dot(), fDist, fNormal, halfTolerance, kInside, kOutside, kSurface, and G4INCL::Math::max().

◆ operator=()

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

Definition at line 140 of file G4Tet.cc.

141{
142 // Check assignment to self
143 //
144 if (this == &rhs) { return *this; }
145
146 // Copy base class data
147 //
149
150 // Copy data
151 //
155 for (G4int i = 0; i < 4; ++i) { fVertex[i] = rhs.fVertex[i]; }
156 for (G4int i = 0; i < 4; ++i) { fNormal[i] = rhs.fNormal[i]; }
157 for (G4int i = 0; i < 4; ++i) { fDist[i] = rhs.fDist[i]; }
158 for (G4int i = 0; i < 4; ++i) { fArea[i] = rhs.fArea[i]; }
159 fBmin = rhs.fBmin;
160 fBmax = rhs.fBmax;
161 fRebuildPolyhedron = false;
162 delete fpPolyhedron; fpPolyhedron = nullptr;
163
164 return *this;
165}
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:107

References fArea, fBmax, fBmin, fCubicVolume, fDist, fNormal, fpPolyhedron, fRebuildPolyhedron, fSurfaceArea, fVertex, halfTolerance, and G4VSolid::operator=().

◆ operator==()

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

◆ PrintWarnings()

void G4Tet::PrintWarnings ( G4bool  )
inline

Definition at line 86 of file G4Tet.hh.

86{};

◆ SetBoundingLimits()

void G4Tet::SetBoundingLimits ( const G4ThreeVector pMin,
const G4ThreeVector pMax 
)

Definition at line 323 of file G4Tet.cc.

325{
326 G4int iout[4] = { 0, 0, 0, 0 };
327 for (G4int i = 0; i < 4; ++i)
328 {
329 iout[i] = (fVertex[i].x() < pMin.x() ||
330 fVertex[i].y() < pMin.y() ||
331 fVertex[i].z() < pMin.z() ||
332 fVertex[i].x() > pMax.x() ||
333 fVertex[i].y() > pMax.y() ||
334 fVertex[i].z() > pMax.z());
335 }
336 if (iout[0] + iout[1] + iout[2] + iout[3] != 0)
337 {
338 std::ostringstream message;
339 message << "Attempt to set bounding box that does not encapsulate solid: "
340 << GetName() << " !\n"
341 << " Specified bounding box limits:\n"
342 << " pmin: " << pMin << "\n"
343 << " pmax: " << pMax << "\n"
344 << " Tetrahedron vertices:\n"
345 << " anchor " << fVertex[0] << ((iout[0]) ? " is outside\n" : "\n")
346 << " p1 " << fVertex[1] << ((iout[1]) ? " is outside\n" : "\n")
347 << " p2 " << fVertex[2] << ((iout[2]) ? " is outside\n" : "\n")
348 << " p3 " << fVertex[3] << ((iout[3]) ? " is outside" : "");
349 G4Exception("G4Tet::SetBoundingLimits()", "GeomSolids0002",
350 FatalException, message);
351 }
352 fBmin = pMin;
353 fBmax = pMax;
354}

References FatalException, fBmax, fBmin, fVertex, G4Exception(), G4VSolid::GetName(), pMax, pMin, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ SetName()

void G4VSolid::SetName ( const G4String name)
inherited

Definition at line 127 of file G4VSolid.cc.

128{
131}
void SetMapValid(G4bool val)
Definition: G4SolidStore.hh:76
static G4SolidStore * GetInstance()
G4String fshapeName
Definition: G4VSolid.hh:312
const char * name(G4int ptype)

References G4VSolid::fshapeName, G4SolidStore::GetInstance(), G4InuclParticleNames::name(), and G4SolidStore::SetMapValid().

Referenced by export_G4VSolid(), G4MultiUnion::G4MultiUnion(), and G4GDMLRead::StripNames().

◆ SetVertices()

void G4Tet::SetVertices ( const G4ThreeVector anchor,
const G4ThreeVector p1,
const G4ThreeVector p2,
const G4ThreeVector p3,
G4bool degeneracyFlag = nullptr 
)

Definition at line 250 of file G4Tet.cc.

254{
255 // Check for degeneracy
256 G4bool degenerate = CheckDegeneracy(p0, p1, p2, p3);
257 if (degeneracyFlag)
258 {
259 *degeneracyFlag = degenerate;
260 }
261 else if (degenerate)
262 {
263 std::ostringstream message;
264 message << "Degenerate tetrahedron is not permitted: " << GetName() << " !\n"
265 << " anchor: " << p0 << "\n"
266 << " p1 : " << p1 << "\n"
267 << " p2 : " << p2 << "\n"
268 << " p3 : " << p3 << "\n"
269 << " volume: "
270 << std::abs((p1 - p0).cross(p2 - p0).dot(p3 - p0))/6.;
271 G4Exception("G4Tet::SetVertices()", "GeomSolids0002",
272 FatalException, message);
273 }
274
275 // Set data members
276 Initialize(p0, p1, p2, p3);
277
278 // Set flag to rebuild polyhedron
279 fRebuildPolyhedron = true;
280}

References CheckDegeneracy(), FatalException, fRebuildPolyhedron, G4Exception(), G4VSolid::GetName(), and Initialize().

◆ StreamInfo()

std::ostream & G4Tet::StreamInfo ( std::ostream &  os) const
virtual

Implements G4VSolid.

Definition at line 603 of file G4Tet.cc.

604{
605 G4int oldprc = os.precision(16);
606 os << "-----------------------------------------------------------\n"
607 << " *** Dump for solid - " << GetName() << " ***\n"
608 << " ===================================================\n"
609 << " Solid type: " << GetEntityType() << "\n"
610 << " Parameters: \n"
611 << " anchor: " << fVertex[0]/mm << " mm\n"
612 << " p1 : " << fVertex[1]/mm << " mm\n"
613 << " p2 : " << fVertex[2]/mm << " mm\n"
614 << " p3 : " << fVertex[3]/mm << " mm\n"
615 << "-----------------------------------------------------------\n";
616 os.precision(oldprc);
617 return os;
618}
static constexpr double mm
Definition: G4SIunits.hh:95
G4GeometryType GetEntityType() const
Definition: G4Tet.cc:585

References fVertex, GetEntityType(), G4VSolid::GetName(), and mm.

◆ SurfaceNormal()

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

Implements G4VSolid.

Definition at line 434 of file G4Tet.cc.

435{
436 G4double k[4];
437 for (G4int i = 0; i < 4; ++i)
438 {
439 k[i] = std::abs(fNormal[i].dot(p) - fDist[i]) <= halfTolerance;
440 }
441 G4double nsurf = k[0] + k[1] + k[2] + k[3];
442 G4ThreeVector norm =
443 k[0]*fNormal[0] + k[1]*fNormal[1] + k[2]*fNormal[2] + k[3]*fNormal[3];
444
445 if (nsurf == 1.) return norm;
446 else if (nsurf > 1.) return norm.unit(); // edge or vertex
447 {
448#ifdef G4SPECSDEBUG
449 std::ostringstream message;
450 G4int oldprc = message.precision(16);
451 message << "Point p is not on surface (!?) of solid: "
452 << GetName() << "\n";
453 message << "Position:\n";
454 message << " p.x() = " << p.x()/mm << " mm\n";
455 message << " p.y() = " << p.y()/mm << " mm\n";
456 message << " p.z() = " << p.z()/mm << " mm";
457 G4cout.precision(oldprc);
458 G4Exception("G4Tet::SurfaceNormal(p)", "GeomSolids1002",
459 JustWarning, message );
460 DumpInfo();
461#endif
462 return ApproxSurfaceNormal(p);
463 }
464}
@ JustWarning
G4GLOB_DLL std::ostream G4cout
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tet.cc:471
void DumpInfo() const

References ApproxSurfaceNormal(), G4VSolid::DumpInfo(), fDist, fNormal, G4cout, G4Exception(), G4VSolid::GetName(), halfTolerance, JustWarning, mm, CLHEP::Hep3Vector::unit(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Field Documentation

◆ fArea

G4double G4Tet::fArea[4] = {0}
private

Definition at line 170 of file G4Tet.hh.

Referenced by G4Tet(), GetPointOnSurface(), Initialize(), and operator=().

◆ fBmax

G4ThreeVector G4Tet::fBmax
private

Definition at line 171 of file G4Tet.hh.

Referenced by BoundingLimits(), G4Tet(), GetExtent(), Initialize(), operator=(), and SetBoundingLimits().

◆ fBmin

G4ThreeVector G4Tet::fBmin
private

Definition at line 171 of file G4Tet.hh.

Referenced by BoundingLimits(), G4Tet(), GetExtent(), Initialize(), operator=(), and SetBoundingLimits().

◆ fCubicVolume

G4double G4Tet::fCubicVolume = 0
private

Definition at line 162 of file G4Tet.hh.

Referenced by G4Tet(), GetCubicVolume(), Initialize(), and operator=().

◆ fDist

G4double G4Tet::fDist[4] = {0}
private

◆ fNormal

G4ThreeVector G4Tet::fNormal[4]
private

◆ fpPolyhedron

G4Polyhedron* G4Tet::fpPolyhedron = nullptr
mutableprivate

Definition at line 165 of file G4Tet.hh.

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

◆ fRebuildPolyhedron

G4bool G4Tet::fRebuildPolyhedron = false
mutableprivate

Definition at line 164 of file G4Tet.hh.

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

◆ fshapeName

G4String G4VSolid::fshapeName
privateinherited

Definition at line 312 of file G4VSolid.hh.

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

◆ fSurfaceArea

G4double G4Tet::fSurfaceArea = 0
private

Definition at line 163 of file G4Tet.hh.

Referenced by G4Tet(), GetPointOnSurface(), GetSurfaceArea(), Initialize(), and operator=().

◆ fVertex

G4ThreeVector G4Tet::fVertex[4]
private

◆ halfTolerance

G4double G4Tet::halfTolerance = 0
private

Definition at line 161 of file G4Tet.hh.

Referenced by DistanceToIn(), DistanceToOut(), G4Tet(), Inside(), operator=(), and SurfaceNormal().

◆ kCarTolerance

G4double G4VSolid::kCarTolerance
protectedinherited

Definition at line 299 of file G4VSolid.hh.

Referenced by G4TessellatedSolid::AddFacet(), G4Polycone::CalculateExtent(), G4Polyhedra::CalculateExtent(), CheckDegeneracy(), G4Para::CheckParameters(), G4Trd::CheckParameters(), G4Ellipsoid::CheckParameters(), G4EllipticalTube::CheckParameters(), G4GenericTrap::ComputeIsTwisted(), G4Polyhedra::Create(), G4GenericPolycone::Create(), G4Polycone::Create(), G4CutTubs::CreatePolyhedron(), G4TessellatedSolid::CreateVertexList(), G4VCSGfaceted::DistanceTo(), G4Sphere::DistanceToIn(), G4Ellipsoid::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(), 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: