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

#include <G4Tet.hh>

+ Inheritance diagram for G4Tet:

Public Member Functions

 G4Tet (const G4String &pName, const G4ThreeVector &anchor, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3, G4bool *degeneracyFlag=nullptr)
 
 ~G4Tet () override
 
void SetVertices (const G4ThreeVector &anchor, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3, G4bool *degeneracyFlag=nullptr)
 
void GetVertices (G4ThreeVector &anchor, G4ThreeVector &p1, G4ThreeVector &p2, G4ThreeVector &p3) const
 
std::vector< G4ThreeVectorGetVertices () const
 
void PrintWarnings (G4bool)
 
G4bool CheckDegeneracy (const G4ThreeVector &p0, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3) const
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
 
void SetBoundingLimits (const G4ThreeVector &pMin, const G4ThreeVector &pMax)
 
void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const override
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const override
 
EInside Inside (const G4ThreeVector &p) const override
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const override
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const override
 
G4double DistanceToIn (const G4ThreeVector &p) const override
 
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
 
G4double DistanceToOut (const G4ThreeVector &p) const override
 
G4GeometryType GetEntityType () const override
 
G4bool IsFaceted () const override
 
G4VSolidClone () const override
 
std::ostream & StreamInfo (std::ostream &os) const override
 
G4double GetCubicVolume () override
 
G4double GetSurfaceArea () override
 
G4ThreeVector GetPointOnSurface () const override
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const override
 
G4VisExtent GetExtent () const override
 
G4PolyhedronCreatePolyhedron () const override
 
G4PolyhedronGetPolyhedron () const override
 
 G4Tet (__void__ &)
 
 G4Tet (const G4Tet &rhs)
 
G4Tetoperator= (const G4Tet &rhs)
 
- Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
 
virtual ~G4VSolid ()
 
G4bool operator== (const G4VSolid &s) const
 
G4String GetName () const
 
void SetName (const G4String &name)
 
G4double GetTolerance () const
 
virtual G4int GetNumOfConstituents () const
 
void DumpInfo () const
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
 G4VSolid (__void__ &)
 
 G4VSolid (const G4VSolid &rhs)
 
G4VSolidoperator= (const G4VSolid &rhs)
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 

Additional Inherited Members

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

Detailed Description

Definition at line 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 != nullptr)
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
93 halfTolerance = 0.5 * kCarTolerance;
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)
bool G4bool
Definition G4Types.hh:86
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:306

Referenced by Clone(), G4Tet(), and operator=().

◆ ~G4Tet()

G4Tet::~G4Tet ( )
override

Definition at line 113 of file G4Tet.cc.

114{
115 delete fpPolyhedron; fpPolyhedron = nullptr;
116}

◆ 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{
125 halfTolerance = rhs.halfTolerance;
126 fCubicVolume = rhs.fCubicVolume;
127 fSurfaceArea = rhs.fSurfaceArea;
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

Member Function Documentation

◆ BoundingLimits()

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

Reimplemented from G4VSolid.

Definition at line 360 of file G4Tet.cc.

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

Referenced by CalculateExtent().

◆ CalculateExtent()

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

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
CLHEP::Hep3Vector G4ThreeVector
std::vector< G4ThreeVector > GetVertices() const
Definition G4Tet.cc:301
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
Definition G4Tet.cc:360

◆ 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}
double G4double
Definition G4Types.hh:83

Referenced by G4Tet(), and SetVertices().

◆ Clone()

G4VSolid * G4Tet::Clone ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 603 of file G4Tet.cc.

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

◆ ComputeDimensions()

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

Reimplemented from G4VSolid.

Definition at line 313 of file G4Tet.cc.

316{
317}

◆ CreatePolyhedron()

G4Polyhedron * G4Tet::CreatePolyhedron ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 698 of file G4Tet.cc.

699{
700 // Check orientation of vertices
701 G4ThreeVector v1 = fVertex[1] - fVertex[0];
702 G4ThreeVector v2 = fVertex[2] - fVertex[0];
703 G4ThreeVector v3 = fVertex[3] - fVertex[0];
704 G4bool invert = v1.cross(v2).dot(v3) < 0.;
705 G4int k2 = (invert) ? 3 : 2;
706 G4int k3 = (invert) ? 2 : 3;
707
708 // Set coordinates of vertices
709 G4double xyz[4][3];
710 for (G4int i = 0; i < 3; ++i)
711 {
712 xyz[0][i] = fVertex[0][i];
713 xyz[1][i] = fVertex[1][i];
714 xyz[2][i] = fVertex[k2][i];
715 xyz[3][i] = fVertex[k3][i];
716 }
717
718 // Create polyhedron
719 G4int faces[4][4] = { {1,3,2,0}, {1,4,3,0}, {1,2,4,0}, {2,3,4,0} };
720 auto ph = new G4Polyhedron;
721 ph->createPolyhedron(4,4,xyz,faces);
722
723 return ph;
724}
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const

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

◆ DescribeYourselfTo()

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

Implements G4VSolid.

Definition at line 678 of file G4Tet.cc.

679{
680 scene.AddSolid (*this);
681}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

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

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}

◆ DistanceToIn() [2/2]

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

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}
#define DBL_MAX
Definition templates.hh:62

◆ DistanceToOut() [1/2]

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

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}

◆ DistanceToOut() [2/2]

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

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] = (G4int)(tmp > 0) * i;
542 nside += (G4int)(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}

◆ GetCubicVolume()

G4double G4Tet::GetCubicVolume ( )
overridevirtual

Reimplemented from G4VSolid.

Definition at line 660 of file G4Tet.cc.

661{
662 return fCubicVolume;
663}

◆ GetEntityType()

G4GeometryType G4Tet::GetEntityType ( ) const
overridevirtual

Implements G4VSolid.

Definition at line 585 of file G4Tet.cc.

586{
587 return {"G4Tet"};
588}

Referenced by StreamInfo().

◆ GetExtent()

G4VisExtent G4Tet::GetExtent ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 687 of file G4Tet.cc.

688{
689 return { fBmin.x(), fBmax.x(),
690 fBmin.y(), fBmax.y(),
691 fBmin.z(), fBmax.z() };
692}

◆ GetPointOnSurface()

G4ThreeVector G4Tet::GetPointOnSurface ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 633 of file G4Tet.cc.

634{
635 constexpr G4int iface[4][3] = { {0,1,2}, {0,2,3}, {0,3,1}, {1,2,3} };
636
637 // Select face
638 G4double select = fSurfaceArea*G4QuickRand();
639 G4int i = 0;
640 i += (G4int)(select > fArea[0]);
641 i += (G4int)(select > fArea[0] + fArea[1]);
642 i += (G4int)(select > fArea[0] + fArea[1] + fArea[2]);
643
644 // Set selected triangle
645 G4ThreeVector p0 = fVertex[iface[i][0]];
646 G4ThreeVector e1 = fVertex[iface[i][1]] - p0;
647 G4ThreeVector e2 = fVertex[iface[i][2]] - p0;
648
649 // Return random point
650 G4double r1 = G4QuickRand();
651 G4double r2 = G4QuickRand();
652 return (r1 + r2 > 1.) ?
653 p0 + e1*(1. - r1) + e2*(1. - r2) : p0 + e1*r1 + e2*r2;
654}
G4double G4QuickRand()

◆ GetPolyhedron()

G4Polyhedron * G4Tet::GetPolyhedron ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 730 of file G4Tet.cc.

731{
732 if (fpPolyhedron == nullptr ||
733 fRebuildPolyhedron ||
734 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
735 fpPolyhedron->GetNumberOfRotationSteps())
736 {
737 G4AutoLock l(&polyhedronMutex);
738 delete fpPolyhedron;
739 fpPolyhedron = CreatePolyhedron();
740 fRebuildPolyhedron = false;
741 l.unlock();
742 }
743 return fpPolyhedron;
744}
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4Polyhedron * CreatePolyhedron() const override
Definition G4Tet.cc:698

◆ GetSurfaceArea()

G4double G4Tet::GetSurfaceArea ( )
overridevirtual

Reimplemented from G4VSolid.

Definition at line 669 of file G4Tet.cc.

670{
671 return fSurfaceArea;
672}

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

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}

Referenced by G4GDMLWriteSolids::TetWrite().

◆ Inside()

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

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

◆ IsFaceted()

G4bool G4Tet::IsFaceted ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 594 of file G4Tet.cc.

595{
596 return true;
597}

◆ 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 //
152 halfTolerance = rhs.halfTolerance;
153 fCubicVolume = rhs.fCubicVolume;
154 fSurfaceArea = rhs.fSurfaceArea;
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

◆ PrintWarnings()

void G4Tet::PrintWarnings ( G4bool )
inline

Definition at line 82 of file G4Tet.hh.

82{};

◆ 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] = (G4int)(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]) != 0 ? " is outside\n" : "\n")
346 << " p1 " << fVertex[1] << ((iout[1]) != 0 ? " is outside\n" : "\n")
347 << " p2 " << fVertex[2] << ((iout[2]) != 0 ? " is outside\n" : "\n")
348 << " p3 " << fVertex[3] << ((iout[3]) != 0 ? " is outside" : "");
349 G4Exception("G4Tet::SetBoundingLimits()", "GeomSolids0002",
350 FatalException, message);
351 }
352 fBmin = pMin;
353 fBmax = pMax;
354}
double z() const
double x() const
double y() const

◆ 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 != nullptr)
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}

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 612 of file G4Tet.cc.

613{
614 G4long oldprc = os.precision(16);
615 os << "-----------------------------------------------------------\n"
616 << " *** Dump for solid - " << GetName() << " ***\n"
617 << " ===================================================\n"
618 << " Solid type: " << GetEntityType() << "\n"
619 << " Parameters: \n"
620 << " anchor: " << fVertex[0]/mm << " mm\n"
621 << " p1 : " << fVertex[1]/mm << " mm\n"
622 << " p2 : " << fVertex[2]/mm << " mm\n"
623 << " p3 : " << fVertex[3]/mm << " mm\n"
624 << "-----------------------------------------------------------\n";
625 os.precision(oldprc);
626 return os;
627}
long G4long
Definition G4Types.hh:87
G4GeometryType GetEntityType() const override
Definition G4Tet.cc:585

◆ SurfaceNormal()

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

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] = (G4double)(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 G4long 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
Hep3Vector unit() const
void DumpInfo() const

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