Geant4 9.6.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, G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4, G4bool *degeneracyFlag=0)
 
virtual ~G4Tet ()
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
 
EInside Inside (const G4ThreeVector &p) const
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
G4double DistanceToIn (const G4ThreeVector &p) const
 
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
 
G4double DistanceToOut (const G4ThreeVector &p) const
 
G4double GetCubicVolume ()
 
G4double GetSurfaceArea ()
 
G4GeometryType GetEntityType () const
 
G4VSolidClone () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4ThreeVector GetPointOnSurface () const
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
G4VisExtent GetExtent () const
 
G4PolyhedronCreatePolyhedron () const
 
G4NURBSCreateNURBS () const
 
G4PolyhedronGetPolyhedron () const
 
 G4Tet (__void__ &)
 
 G4Tet (const G4Tet &rhs)
 
G4Tetoperator= (const G4Tet &rhs)
 
const char * CVSHeaderVers ()
 
const char * CVSFileVers ()
 
void PrintWarnings (G4bool flag)
 
std::vector< G4ThreeVectorGetVertices () const
 
- 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 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
 
virtual G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p) const =0
 
virtual void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
virtual G4double GetCubicVolume ()
 
virtual G4double GetSurfaceArea ()
 
virtual G4GeometryType GetEntityType () const =0
 
virtual G4ThreeVector GetPointOnSurface () const
 
virtual G4VSolidClone () const
 
virtual std::ostream & StreamInfo (std::ostream &os) const =0
 
void DumpInfo () const
 
virtual void DescribeYourselfTo (G4VGraphicsScene &scene) const =0
 
virtual G4VisExtent GetExtent () const
 
virtual G4PolyhedronCreatePolyhedron () const
 
virtual G4NURBSCreateNURBS () const
 
virtual G4PolyhedronGetPolyhedron () 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)
 

Static Public Member Functions

static G4bool CheckDegeneracy (G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4)
 

Protected Member Functions

G4ThreeVectorListCreateRotatedVertices (const G4AffineTransform &pTransform) const
 
- 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
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 

Additional Inherited Members

- Protected Attributes inherited from G4VSolid
G4double kCarTolerance
 

Detailed Description

Definition at line 56 of file G4Tet.hh.

Constructor & Destructor Documentation

◆ G4Tet() [1/3]

G4Tet::G4Tet ( const G4String pName,
G4ThreeVector  anchor,
G4ThreeVector  p2,
G4ThreeVector  p3,
G4ThreeVector  p4,
G4bool degeneracyFlag = 0 
)

Definition at line 89 of file G4Tet.cc.

94 : G4VSolid(pName), fpPolyhedron(0), warningFlag(0)
95{
96 // fV<x><y> is vector from vertex <y> to vertex <x>
97 //
98 G4ThreeVector fV21=p2-anchor;
99 G4ThreeVector fV31=p3-anchor;
100 G4ThreeVector fV41=p4-anchor;
101
102 // make sure this is a correctly oriented set of points for the tetrahedron
103 //
104 G4double signed_vol=fV21.cross(fV31).dot(fV41);
105 if(signed_vol<0.0)
106 {
107 G4ThreeVector temp(p4);
108 p4=p3;
109 p3=temp;
110 temp=fV41;
111 fV41=fV31;
112 fV31=temp;
113 }
114 fCubicVolume = std::fabs(signed_vol) / 6.;
115
116 G4ThreeVector fV24=p2-p4;
117 G4ThreeVector fV43=p4-p3;
118 G4ThreeVector fV32=p3-p2;
119
120 fXMin=std::min(std::min(std::min(anchor.x(), p2.x()),p3.x()),p4.x());
121 fXMax=std::max(std::max(std::max(anchor.x(), p2.x()),p3.x()),p4.x());
122 fYMin=std::min(std::min(std::min(anchor.y(), p2.y()),p3.y()),p4.y());
123 fYMax=std::max(std::max(std::max(anchor.y(), p2.y()),p3.y()),p4.y());
124 fZMin=std::min(std::min(std::min(anchor.z(), p2.z()),p3.z()),p4.z());
125 fZMax=std::max(std::max(std::max(anchor.z(), p2.z()),p3.z()),p4.z());
126
127 fDx=(fXMax-fXMin)*0.5; fDy=(fYMax-fYMin)*0.5; fDz=(fZMax-fZMin)*0.5;
128
129 fMiddle=G4ThreeVector(fXMax+fXMin, fYMax+fYMin, fZMax+fZMin)*0.5;
130 fMaxSize=std::max(std::max(std::max((anchor-fMiddle).mag(),
131 (p2-fMiddle).mag()),
132 (p3-fMiddle).mag()),
133 (p4-fMiddle).mag());
134
135 G4bool degenerate=std::fabs(signed_vol) < 1e-9*fMaxSize*fMaxSize*fMaxSize;
136
137 if(degeneracyFlag) *degeneracyFlag=degenerate;
138 else if (degenerate)
139 {
140 G4Exception("G4Tet::G4Tet()", "GeomSolids0002", FatalException,
141 "Degenerate tetrahedron not allowed.");
142 }
143
144 fTol=1e-9*(std::fabs(fXMin)+std::fabs(fXMax)+std::fabs(fYMin)
145 +std::fabs(fYMax)+std::fabs(fZMin)+std::fabs(fZMax));
146 //fTol=kCarTolerance;
147
148 fAnchor=anchor;
149 fP2=p2;
150 fP3=p3;
151 fP4=p4;
152
153 G4ThreeVector fCenter123=(anchor+p2+p3)*(1.0/3.0); // face center
154 G4ThreeVector fCenter134=(anchor+p4+p3)*(1.0/3.0);
155 G4ThreeVector fCenter142=(anchor+p4+p2)*(1.0/3.0);
156 G4ThreeVector fCenter234=(p2+p3+p4)*(1.0/3.0);
157
158 // compute area of each triangular face by cross product
159 // and sum for total surface area
160
161 G4ThreeVector normal123=fV31.cross(fV21);
162 G4ThreeVector normal134=fV41.cross(fV31);
163 G4ThreeVector normal142=fV21.cross(fV41);
164 G4ThreeVector normal234=fV32.cross(fV43);
165
166 fSurfaceArea=(
167 normal123.mag()+
168 normal134.mag()+
169 normal142.mag()+
170 normal234.mag()
171 )/2.0;
172
173 fNormal123=normal123.unit();
174 fNormal134=normal134.unit();
175 fNormal142=normal142.unit();
176 fNormal234=normal234.unit();
177
178 fCdotN123=fCenter123.dot(fNormal123);
179 fCdotN134=fCenter134.dot(fNormal134);
180 fCdotN142=fCenter142.dot(fNormal142);
181 fCdotN234=fCenter234.dot(fNormal234);
182}
@ FatalException
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
Hep3Vector unit() const
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const
double mag() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ ~G4Tet()

G4Tet::~G4Tet ( )
virtual

Definition at line 204 of file G4Tet.cc.

205{
206 delete fpPolyhedron;
207}

◆ G4Tet() [2/3]

G4Tet::G4Tet ( __void__ &  a)

Definition at line 189 of file G4Tet.cc.

190 : G4VSolid(a), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0),
191 fAnchor(0,0,0), fP2(0,0,0), fP3(0,0,0), fP4(0,0,0), fMiddle(0,0,0),
192 fNormal123(0,0,0), fNormal142(0,0,0), fNormal134(0,0,0),
193 fNormal234(0,0,0), warningFlag(0),
194 fCdotN123(0.), fCdotN142(0.), fCdotN134(0.), fCdotN234(0.),
195 fXMin(0.), fXMax(0.), fYMin(0.), fYMax(0.), fZMin(0.), fZMax(0.),
196 fDx(0.), fDy(0.), fDz(0.), fTol(0.), fMaxSize(0.)
197{
198}

◆ G4Tet() [3/3]

G4Tet::G4Tet ( const G4Tet rhs)

Definition at line 213 of file G4Tet.cc.

214 : G4VSolid(rhs),
215 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
216 fpPolyhedron(0), fAnchor(rhs.fAnchor),
217 fP2(rhs.fP2), fP3(rhs.fP3), fP4(rhs.fP4), fMiddle(rhs.fMiddle),
218 fNormal123(rhs.fNormal123), fNormal142(rhs.fNormal142),
219 fNormal134(rhs.fNormal134), fNormal234(rhs.fNormal234),
220 warningFlag(rhs.warningFlag), fCdotN123(rhs.fCdotN123),
221 fCdotN142(rhs.fCdotN142), fCdotN134(rhs.fCdotN134),
222 fCdotN234(rhs.fCdotN234), fXMin(rhs.fXMin), fXMax(rhs.fXMax),
223 fYMin(rhs.fYMin), fYMax(rhs.fYMax), fZMin(rhs.fZMin), fZMax(rhs.fZMax),
224 fDx(rhs.fDx), fDy(rhs.fDy), fDz(rhs.fDz), fTol(rhs.fTol),
225 fMaxSize(rhs.fMaxSize)
226{
227}

Member Function Documentation

◆ CalculateExtent()

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

Implements G4VSolid.

Definition at line 291 of file G4Tet.cc.

295{
296 G4double xMin,xMax;
297 G4double yMin,yMax;
298 G4double zMin,zMax;
299
300 if (pTransform.IsRotated())
301 {
302 G4ThreeVector pp0=pTransform.TransformPoint(fAnchor);
303 G4ThreeVector pp1=pTransform.TransformPoint(fP2);
304 G4ThreeVector pp2=pTransform.TransformPoint(fP3);
305 G4ThreeVector pp3=pTransform.TransformPoint(fP4);
306
307 xMin = std::min(std::min(std::min(pp0.x(), pp1.x()),pp2.x()),pp3.x());
308 xMax = std::max(std::max(std::max(pp0.x(), pp1.x()),pp2.x()),pp3.x());
309 yMin = std::min(std::min(std::min(pp0.y(), pp1.y()),pp2.y()),pp3.y());
310 yMax = std::max(std::max(std::max(pp0.y(), pp1.y()),pp2.y()),pp3.y());
311 zMin = std::min(std::min(std::min(pp0.z(), pp1.z()),pp2.z()),pp3.z());
312 zMax = std::max(std::max(std::max(pp0.z(), pp1.z()),pp2.z()),pp3.z());
313
314 }
315 else
316 {
317 G4double xoffset = pTransform.NetTranslation().x() ;
318 xMin = xoffset + fXMin;
319 xMax = xoffset + fXMax;
320 G4double yoffset = pTransform.NetTranslation().y() ;
321 yMin = yoffset + fYMin;
322 yMax = yoffset + fYMax;
323 G4double zoffset = pTransform.NetTranslation().z() ;
324 zMin = zoffset + fZMin;
325 zMax = zoffset + fZMax;
326 }
327
328 if (pVoxelLimit.IsXLimited())
329 {
330 if ( (xMin > pVoxelLimit.GetMaxXExtent()+fTol) ||
331 (xMax < pVoxelLimit.GetMinXExtent()-fTol) ) { return false; }
332 else
333 {
334 xMin = std::max(xMin, pVoxelLimit.GetMinXExtent());
335 xMax = std::min(xMax, pVoxelLimit.GetMaxXExtent());
336 }
337 }
338
339 if (pVoxelLimit.IsYLimited())
340 {
341 if ( (yMin > pVoxelLimit.GetMaxYExtent()+fTol) ||
342 (yMax < pVoxelLimit.GetMinYExtent()-fTol) ) { return false; }
343 else
344 {
345 yMin = std::max(yMin, pVoxelLimit.GetMinYExtent());
346 yMax = std::min(yMax, pVoxelLimit.GetMaxYExtent());
347 }
348 }
349
350 if (pVoxelLimit.IsZLimited())
351 {
352 if ( (zMin > pVoxelLimit.GetMaxZExtent()+fTol) ||
353 (zMax < pVoxelLimit.GetMinZExtent()-fTol) ) { return false; }
354 else
355 {
356 zMin = std::max(zMin, pVoxelLimit.GetMinZExtent());
357 zMax = std::min(zMax, pVoxelLimit.GetMaxZExtent());
358 }
359 }
360
361 switch (pAxis)
362 {
363 case kXAxis:
364 pMin=xMin;
365 pMax=xMax;
366 break;
367 case kYAxis:
368 pMin=yMin;
369 pMax=yMax;
370 break;
371 case kZAxis:
372 pMin=zMin;
373 pMax=zMax;
374 break;
375 default:
376 break;
377 }
378
379 return true;
380}
double z() const
double x() const
double y() const
G4bool IsRotated() const
G4ThreeVector NetTranslation() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4bool IsYLimited() const
G4double GetMinZExtent() const
G4bool IsXLimited() const
G4double GetMaxYExtent() const
G4double GetMaxZExtent() const
G4double GetMinYExtent() const
G4double GetMinXExtent() const
G4bool IsZLimited() const
G4double GetMaxXExtent() const
@ kYAxis
Definition: geomdefs.hh:54
@ kXAxis
Definition: geomdefs.hh:54
@ kZAxis
Definition: geomdefs.hh:54

◆ CheckDegeneracy()

G4bool G4Tet::CheckDegeneracy ( G4ThreeVector  anchor,
G4ThreeVector  p2,
G4ThreeVector  p3,
G4ThreeVector  p4 
)
static

Definition at line 265 of file G4Tet.cc.

269{
270 G4bool result;
271 G4Tet *object=new G4Tet("temp",anchor,p2,p3,p4,&result);
272 delete object;
273 return result;
274}
Definition: G4Tet.hh:57

◆ Clone()

G4VSolid * G4Tet::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 659 of file G4Tet.cc.

660{
661 return new G4Tet(*this);
662}

◆ ComputeDimensions()

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

Reimplemented from G4VSolid.

Definition at line 281 of file G4Tet.cc.

284{
285}

◆ CreateNURBS()

G4NURBS * G4Tet::CreateNURBS ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 811 of file G4Tet.cc.

812{
813 return new G4NURBSbox (fDx, fDy, fDz);
814}

◆ CreatePolyhedron()

G4Polyhedron * G4Tet::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 792 of file G4Tet.cc.

793{
795 G4double xyz[4][3];
796 static G4int faces[4][4]={{1,3,2,0},{1,4,3,0},{1,2,4,0},{2,3,4,0}};
797 xyz[0][0]=fAnchor.x(); xyz[0][1]=fAnchor.y(); xyz[0][2]=fAnchor.z();
798 xyz[1][0]=fP2.x(); xyz[1][1]=fP2.y(); xyz[1][2]=fP2.z();
799 xyz[2][0]=fP3.x(); xyz[2][1]=fP3.y(); xyz[2][2]=fP3.z();
800 xyz[3][0]=fP4.x(); xyz[3][1]=fP4.y(); xyz[3][2]=fP4.z();
801
802 ph->createPolyhedron(4,4,xyz,faces);
803
804 return ph;
805}
int G4int
Definition: G4Types.hh:66
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])

Referenced by GetPolyhedron().

◆ CreateRotatedVertices()

G4ThreeVectorList * G4Tet::CreateRotatedVertices ( const G4AffineTransform pTransform) const
protected

Definition at line 619 of file G4Tet.cc.

620{
621 G4ThreeVectorList* vertices = new G4ThreeVectorList();
622
623 if (vertices)
624 {
625 vertices->reserve(4);
626 G4ThreeVector vertex0(fAnchor);
627 G4ThreeVector vertex1(fP2);
628 G4ThreeVector vertex2(fP3);
629 G4ThreeVector vertex3(fP4);
630
631 vertices->push_back(pTransform.TransformPoint(vertex0));
632 vertices->push_back(pTransform.TransformPoint(vertex1));
633 vertices->push_back(pTransform.TransformPoint(vertex2));
634 vertices->push_back(pTransform.TransformPoint(vertex3));
635 }
636 else
637 {
638 DumpInfo();
639 G4Exception("G4Tet::CreateRotatedVertices()",
640 "GeomSolids0003", FatalException,
641 "Error in allocation of vertices. Out of memory !");
642 }
643 return vertices;
644}
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:85
void DumpInfo() const

◆ CVSFileVers()

const char * G4Tet::CVSFileVers ( )
inline

Definition at line 126 of file G4Tet.hh.

127 { return CVSVers; }

◆ CVSHeaderVers()

const char * G4Tet::CVSHeaderVers ( )
inline

Definition at line 124 of file G4Tet.hh.

125 { return "$Id: G4Tet.hh,v 1.11 2010-10-20 08:54:18 gcosmo Exp $"; }

◆ DescribeYourselfTo()

void G4Tet::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 774 of file G4Tet.cc.

775{
776 scene.AddSolid (*this);
777}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 517 of file G4Tet.cc.

518{
519 G4double dd=(p-fMiddle).mag() - fMaxSize - fTol;
520 return std::max(0.0, dd);
521}

◆ DistanceToIn() [2/2]

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

Implements G4VSolid.

Definition at line 436 of file G4Tet.cc.

438{
439 G4ThreeVector vu(v.unit()), hp;
440 G4double vdotn, t, tmin=kInfinity;
441
442 G4double extraDistance=10.0*fTol; // a little ways into the solid
443
444 vdotn=-vu.dot(fNormal123);
445 if(vdotn > 1e-12)
446 { // this is a candidate face, since it is pointing at us
447 t=(p.dot(fNormal123)-fCdotN123)/vdotn; // # distance to intersection
448 if( (t>=-fTol) && (t<tmin) )
449 { // if not true, we're going away from this face or it's not close
450 hp=p+vu*(t+extraDistance); // a little beyond point of intersection
451 if ( ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
452 ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
453 ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
454 {
455 tmin=t;
456 }
457 }
458 }
459
460 vdotn=-vu.dot(fNormal134);
461 if(vdotn > 1e-12)
462 { // # this is a candidate face, since it is pointing at us
463 t=(p.dot(fNormal134)-fCdotN134)/vdotn; // # distance to intersection
464 if( (t>=-fTol) && (t<tmin) )
465 { // if not true, we're going away from this face
466 hp=p+vu*(t+extraDistance); // a little beyond point of intersection
467 if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
468 ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
469 ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
470 {
471 tmin=t;
472 }
473 }
474 }
475
476 vdotn=-vu.dot(fNormal142);
477 if(vdotn > 1e-12)
478 { // # this is a candidate face, since it is pointing at us
479 t=(p.dot(fNormal142)-fCdotN142)/vdotn; // # distance to intersection
480 if( (t>=-fTol) && (t<tmin) )
481 { // if not true, we're going away from this face
482 hp=p+vu*(t+extraDistance); // a little beyond point of intersection
483 if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
484 ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
485 ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
486 {
487 tmin=t;
488 }
489 }
490 }
491
492 vdotn=-vu.dot(fNormal234);
493 if(vdotn > 1e-12)
494 { // # this is a candidate face, since it is pointing at us
495 t=(p.dot(fNormal234)-fCdotN234)/vdotn; // # distance to intersection
496 if( (t>=-fTol) && (t<tmin) )
497 { // if not true, we're going away from this face
498 hp=p+vu*(t+extraDistance); // a little beyond point of intersection
499 if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
500 ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
501 ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) )
502 {
503 tmin=t;
504 }
505 }
506 }
507
508 return std::max(0.0,tmin);
509}

◆ DistanceToOut() [1/2]

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

Implements G4VSolid.

Definition at line 598 of file G4Tet.cc.

599{
600 G4double t1,t2,t3,t4;
601 t1=fCdotN123-p.dot(fNormal123); // distance to plane, positive if inside
602 t2=fCdotN134-p.dot(fNormal134); // distance to plane
603 t3=fCdotN142-p.dot(fNormal142); // distance to plane
604 t4=fCdotN234-p.dot(fNormal234); // distance to plane
605
606 // if any one of these is negative, we are outside,
607 // so return zero in that case
608
609 G4double tmin=std::min(std::min(std::min(t1,t2),t3),t4);
610 return (tmin < fTol)? 0:tmin;
611}

◆ DistanceToOut() [2/2]

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

Implements G4VSolid.

Definition at line 529 of file G4Tet.cc.

532{
533 G4ThreeVector vu(v.unit());
534 G4double t1=kInfinity,t2=kInfinity,t3=kInfinity,t4=kInfinity, vdotn, tt;
535
536 vdotn=vu.dot(fNormal123);
537 if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
538 {
539 t1=(fCdotN123-p.dot(fNormal123))/vdotn; // # distance to intersection
540 }
541
542 vdotn=vu.dot(fNormal134);
543 if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
544 {
545 t2=(fCdotN134-p.dot(fNormal134))/vdotn; // # distance to intersection
546 }
547
548 vdotn=vu.dot(fNormal142);
549 if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
550 {
551 t3=(fCdotN142-p.dot(fNormal142))/vdotn; // # distance to intersection
552 }
553
554 vdotn=vu.dot(fNormal234);
555 if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
556 {
557 t4=(fCdotN234-p.dot(fNormal234))/vdotn; // # distance to intersection
558 }
559
560 tt=std::min(std::min(std::min(t1,t2),t3),t4);
561
562 if (warningFlag && (tt == kInfinity || tt < -fTol))
563 {
564 DumpInfo();
565 std::ostringstream message;
566 message << "No good intersection found or already outside!?" << G4endl
567 << "p = " << p / mm << "mm" << G4endl
568 << "v = " << v << G4endl
569 << "t1, t2, t3, t4 (mm) "
570 << t1/mm << ", " << t2/mm << ", " << t3/mm << ", " << t4/mm;
571 G4Exception("G4Tet::DistanceToOut(p,v,...)", "GeomSolids1002",
572 JustWarning, message);
573 if(validNorm)
574 {
575 *validNorm=false; // flag normal as meaningless
576 }
577 }
578 else if(calcNorm && n)
579 {
580 static G4ThreeVector normal;
581 if(tt==t1) { normal=fNormal123; }
582 else if (tt==t2) { normal=fNormal134; }
583 else if (tt==t3) { normal=fNormal142; }
584 else if (tt==t4) { normal=fNormal234; }
585 n=&normal;
586 if(validNorm) { *validNorm=true; }
587 }
588
589 return std::max(tt,0.0); // avoid tt<0.0 by a tiny bit
590 // if we are right on a face
591}
@ JustWarning
#define G4endl
Definition: G4ios.hh:52

◆ GetCubicVolume()

G4double G4Tet::GetCubicVolume ( )
virtual

Reimplemented from G4VSolid.

Definition at line 754 of file G4Tet.cc.

755{
756 return fCubicVolume;
757}

◆ GetEntityType()

G4GeometryType G4Tet::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 650 of file G4Tet.cc.

651{
652 return G4String("G4Tet");
653}

◆ GetExtent()

G4VisExtent G4Tet::GetExtent ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 783 of file G4Tet.cc.

784{
785 return G4VisExtent (fXMin, fXMax, fYMin, fYMax, fZMin, fZMax);
786}

◆ GetPointOnSurface()

G4ThreeVector G4Tet::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 718 of file G4Tet.cc.

719{
720 G4double chose,aOne,aTwo,aThree,aFour;
721 G4ThreeVector p1, p2, p3, p4;
722
723 p1 = GetPointOnFace(fAnchor,fP2,fP3,aOne);
724 p2 = GetPointOnFace(fAnchor,fP4,fP3,aTwo);
725 p3 = GetPointOnFace(fAnchor,fP4,fP2,aThree);
726 p4 = GetPointOnFace(fP4,fP3,fP2,aFour);
727
728 chose = RandFlat::shoot(0.,aOne+aTwo+aThree+aFour);
729 if( (chose>=0.) && (chose <aOne) ) {return p1;}
730 else if( (chose>=aOne) && (chose < aOne+aTwo) ) {return p2;}
731 else if( (chose>=aOne+aTwo) && (chose<aOne+aTwo+aThree) ) {return p3;}
732 return p4;
733}
static double shoot()
Definition: RandFlat.cc:59

◆ GetPolyhedron()

G4Polyhedron * G4Tet::GetPolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 820 of file G4Tet.cc.

821{
822 if (!fpPolyhedron ||
824 fpPolyhedron->GetNumberOfRotationSteps())
825 {
826 delete fpPolyhedron;
827 fpPolyhedron = CreatePolyhedron();
828 }
829 return fpPolyhedron;
830}
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tet.cc:792
static G4int GetNumberOfRotationSteps()

◆ GetSurfaceArea()

G4double G4Tet::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 763 of file G4Tet.cc.

764{
765 return fSurfaceArea;
766}

◆ GetVertices()

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

Definition at line 739 of file G4Tet.cc.

740{
741 std::vector<G4ThreeVector> vertices(4);
742 vertices[0] = fAnchor;
743 vertices[1] = fP2;
744 vertices[2] = fP3;
745 vertices[3] = fP4;
746
747 return vertices;
748}

Referenced by G4GDMLWriteSolids::TetWrite().

◆ Inside()

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

Implements G4VSolid.

Definition at line 386 of file G4Tet.cc.

387{
388 G4double r123, r134, r142, r234;
389
390 // this is written to allow if-statement truncation so the outside test
391 // (where most of the world is) can fail very quickly and efficiently
392
393 if ( (r123=p.dot(fNormal123)-fCdotN123) > fTol ||
394 (r134=p.dot(fNormal134)-fCdotN134) > fTol ||
395 (r142=p.dot(fNormal142)-fCdotN142) > fTol ||
396 (r234=p.dot(fNormal234)-fCdotN234) > fTol )
397 {
398 return kOutside; // at least one is out!
399 }
400 else if( (r123 < -fTol)&&(r134 < -fTol)&&(r142 < -fTol)&&(r234 < -fTol) )
401 {
402 return kInside; // all are definitively inside
403 }
404 else
405 {
406 return kSurface; // too close to tell
407 }
408}
@ kInside
Definition: geomdefs.hh:58
@ kOutside
Definition: geomdefs.hh:58
@ kSurface
Definition: geomdefs.hh:58

◆ operator=()

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

Definition at line 234 of file G4Tet.cc.

235{
236 // Check assignment to self
237 //
238 if (this == &rhs) { return *this; }
239
240 // Copy base class data
241 //
243
244 // Copy data
245 //
246 fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
247 fpPolyhedron = 0; fAnchor = rhs.fAnchor;
248 fP2 = rhs.fP2; fP3 = rhs.fP3; fP4 = rhs.fP4; fMiddle = rhs.fMiddle;
249 fNormal123 = rhs.fNormal123; fNormal142 = rhs.fNormal142;
250 fNormal134 = rhs.fNormal134; fNormal234 = rhs.fNormal234;
251 warningFlag = rhs.warningFlag; fCdotN123 = rhs.fCdotN123;
252 fCdotN142 = rhs.fCdotN142; fCdotN134 = rhs.fCdotN134;
253 fCdotN234 = rhs.fCdotN234; fXMin = rhs.fXMin; fXMax = rhs.fXMax;
254 fYMin = rhs.fYMin; fYMax = rhs.fYMax; fZMin = rhs.fZMin; fZMax = rhs.fZMax;
255 fDx = rhs.fDx; fDy = rhs.fDy; fDz = rhs.fDz; fTol = rhs.fTol;
256 fMaxSize = rhs.fMaxSize;
257
258 return *this;
259}
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110

◆ PrintWarnings()

void G4Tet::PrintWarnings ( G4bool  flag)
inline

Definition at line 128 of file G4Tet.hh.

129 { warningFlag=flag; }

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 668 of file G4Tet.cc.

669{
670 G4int oldprc = os.precision(16);
671 os << "-----------------------------------------------------------\n"
672 << " *** Dump for solid - " << GetName() << " ***\n"
673 << " ===================================================\n"
674 << " Solid type: G4Tet\n"
675 << " Parameters: \n"
676 << " anchor: " << fAnchor/mm << " mm \n"
677 << " p2: " << fP2/mm << " mm \n"
678 << " p3: " << fP3/mm << " mm \n"
679 << " p4: " << fP4/mm << " mm \n"
680 << " normal123: " << fNormal123 << " \n"
681 << " normal134: " << fNormal134 << " \n"
682 << " normal142: " << fNormal142 << " \n"
683 << " normal234: " << fNormal234 << " \n"
684 << "-----------------------------------------------------------\n";
685 os.precision(oldprc);
686
687 return os;
688}
G4String GetName() const

◆ SurfaceNormal()

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

Implements G4VSolid.

Definition at line 417 of file G4Tet.cc.

418{
419 G4double r123=std::fabs(p.dot(fNormal123)-fCdotN123);
420 G4double r134=std::fabs(p.dot(fNormal134)-fCdotN134);
421 G4double r142=std::fabs(p.dot(fNormal142)-fCdotN142);
422 G4double r234=std::fabs(p.dot(fNormal234)-fCdotN234);
423
424 if( (r123<=r134) && (r123<=r142) && (r123<=r234) ) { return fNormal123; }
425 else if ( (r134<=r142) && (r134<=r234) ) { return fNormal134; }
426 else if (r142 <= r234) { return fNormal142; }
427 return fNormal234;
428}

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