Geant4  10.02.p03
G4Tet Class Reference

#include <G4Tet.hh>

Inheritance diagram for G4Tet:
Collaboration 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
 
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
 
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
 

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
 

Private Member Functions

G4ThreeVector GetPointOnFace (G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4double &area) const
 

Private Attributes

G4double fCubicVolume
 
G4double fSurfaceArea
 
G4bool fRebuildPolyhedron
 
G4PolyhedronfpPolyhedron
 
G4ThreeVector fAnchor
 
G4ThreeVector fP2
 
G4ThreeVector fP3
 
G4ThreeVector fP4
 
G4ThreeVector fMiddle
 
G4ThreeVector fNormal123
 
G4ThreeVector fNormal142
 
G4ThreeVector fNormal134
 
G4ThreeVector fNormal234
 
G4bool warningFlag
 
G4double fCdotN123
 
G4double fCdotN142
 
G4double fCdotN134
 
G4double fCdotN234
 
G4double fXMin
 
G4double fXMax
 
G4double fYMin
 
G4double fYMax
 
G4double fZMin
 
G4double fZMax
 
G4double fDx
 
G4double fDy
 
G4double fDz
 
G4double fTol
 
G4double fMaxSize
 

Static Private Attributes

static const char CVSVers [] ="$Id: G4Tet.cc 102297 2017-01-20 13:33:54Z gcosmo $"
 

Additional Inherited Members

- Protected Attributes inherited from G4VSolid
G4double kCarTolerance
 

Detailed Description

Definition at line 65 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 96 of file G4Tet.cc.

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

◆ ~G4Tet()

G4Tet::~G4Tet ( )
virtual

Definition at line 212 of file G4Tet.cc.

213 {
214  delete fpPolyhedron; fpPolyhedron = 0;
215 }
G4Polyhedron * fpPolyhedron
Definition: G4Tet.hh:157

◆ G4Tet() [2/3]

G4Tet::G4Tet ( __void__ &  a)

Definition at line 196 of file G4Tet.cc.

197  : G4VSolid(a), fCubicVolume(0.), fSurfaceArea(0.),
198  fRebuildPolyhedron(false), fpPolyhedron(0),
199  fAnchor(0,0,0), fP2(0,0,0), fP3(0,0,0), fP4(0,0,0), fMiddle(0,0,0),
200  fNormal123(0,0,0), fNormal142(0,0,0), fNormal134(0,0,0),
201  fNormal234(0,0,0), warningFlag(0),
202  fCdotN123(0.), fCdotN142(0.), fCdotN134(0.), fCdotN234(0.),
203  fXMin(0.), fXMax(0.), fYMin(0.), fYMax(0.), fZMin(0.), fZMax(0.),
204  fDx(0.), fDy(0.), fDz(0.), fTol(0.), fMaxSize(0.)
205 {
206 }
G4ThreeVector fP4
Definition: G4Tet.hh:165
G4double fTol
Definition: G4Tet.hh:172
G4double fYMin
Definition: G4Tet.hh:171
G4double fXMax
Definition: G4Tet.hh:171
G4ThreeVector fNormal142
Definition: G4Tet.hh:166
G4bool fRebuildPolyhedron
Definition: G4Tet.hh:156
G4double fDz
Definition: G4Tet.hh:172
G4ThreeVector fP2
Definition: G4Tet.hh:165
G4double fYMax
Definition: G4Tet.hh:171
G4ThreeVector fNormal123
Definition: G4Tet.hh:166
G4double fCdotN142
Definition: G4Tet.hh:170
G4ThreeVector fP3
Definition: G4Tet.hh:165
G4double fCdotN134
Definition: G4Tet.hh:170
G4double fMaxSize
Definition: G4Tet.hh:172
G4ThreeVector fAnchor
Definition: G4Tet.hh:165
G4double fZMin
Definition: G4Tet.hh:171
G4double fXMin
Definition: G4Tet.hh:171
G4ThreeVector fMiddle
Definition: G4Tet.hh:165
G4double fSurfaceArea
Definition: G4Tet.hh:154
G4double fCdotN234
Definition: G4Tet.hh:170
G4double fDy
Definition: G4Tet.hh:172
G4ThreeVector fNormal134
Definition: G4Tet.hh:166
G4double fCdotN123
Definition: G4Tet.hh:170
G4double fCubicVolume
Definition: G4Tet.hh:154
G4double fDx
Definition: G4Tet.hh:172
G4VSolid(const G4String &name)
Definition: G4VSolid.cc:60
G4ThreeVector fNormal234
Definition: G4Tet.hh:166
G4Polyhedron * fpPolyhedron
Definition: G4Tet.hh:157
G4bool warningFlag
Definition: G4Tet.hh:168
G4double fZMax
Definition: G4Tet.hh:171

◆ G4Tet() [3/3]

G4Tet::G4Tet ( const G4Tet rhs)

Definition at line 221 of file G4Tet.cc.

222  : G4VSolid(rhs),
225  fP2(rhs.fP2), fP3(rhs.fP3), fP4(rhs.fP4), fMiddle(rhs.fMiddle),
230  fCdotN234(rhs.fCdotN234), fXMin(rhs.fXMin), fXMax(rhs.fXMax),
231  fYMin(rhs.fYMin), fYMax(rhs.fYMax), fZMin(rhs.fZMin), fZMax(rhs.fZMax),
232  fDx(rhs.fDx), fDy(rhs.fDy), fDz(rhs.fDz), fTol(rhs.fTol),
233  fMaxSize(rhs.fMaxSize)
234 {
235 }
G4ThreeVector fP4
Definition: G4Tet.hh:165
G4double fTol
Definition: G4Tet.hh:172
G4double fYMin
Definition: G4Tet.hh:171
G4double fXMax
Definition: G4Tet.hh:171
G4ThreeVector fNormal142
Definition: G4Tet.hh:166
G4bool fRebuildPolyhedron
Definition: G4Tet.hh:156
G4double fDz
Definition: G4Tet.hh:172
G4ThreeVector fP2
Definition: G4Tet.hh:165
G4double fYMax
Definition: G4Tet.hh:171
G4ThreeVector fNormal123
Definition: G4Tet.hh:166
G4double fCdotN142
Definition: G4Tet.hh:170
G4ThreeVector fP3
Definition: G4Tet.hh:165
G4double fCdotN134
Definition: G4Tet.hh:170
G4double fMaxSize
Definition: G4Tet.hh:172
G4ThreeVector fAnchor
Definition: G4Tet.hh:165
G4double fZMin
Definition: G4Tet.hh:171
G4double fXMin
Definition: G4Tet.hh:171
G4ThreeVector fMiddle
Definition: G4Tet.hh:165
G4double fSurfaceArea
Definition: G4Tet.hh:154
G4double fCdotN234
Definition: G4Tet.hh:170
G4double fDy
Definition: G4Tet.hh:172
G4ThreeVector fNormal134
Definition: G4Tet.hh:166
G4double fCdotN123
Definition: G4Tet.hh:170
G4double fCubicVolume
Definition: G4Tet.hh:154
G4double fDx
Definition: G4Tet.hh:172
G4VSolid(const G4String &name)
Definition: G4VSolid.cc:60
G4ThreeVector fNormal234
Definition: G4Tet.hh:166
G4Polyhedron * fpPolyhedron
Definition: G4Tet.hh:157
G4bool warningFlag
Definition: G4Tet.hh:168
G4double fZMax
Definition: G4Tet.hh:171

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 301 of file G4Tet.cc.

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

◆ CheckDegeneracy()

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

Definition at line 275 of file G4Tet.cc.

279 {
280  G4bool result;
281  G4Tet *object=new G4Tet("temp",anchor,p2,p3,p4,&result);
282  delete object;
283  return result;
284 }
Definition: G4Tet.hh:65
bool G4bool
Definition: G4Types.hh:79
G4Tet(const G4String &pName, G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4, G4bool *degeneracyFlag=0)
Definition: G4Tet.cc:96
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Clone()

G4VSolid * G4Tet::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 710 of file G4Tet.cc.

711 {
712  return new G4Tet(*this);
713 }
G4Tet(const G4String &pName, G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4, G4bool *degeneracyFlag=0)
Definition: G4Tet.cc:96
Here is the call graph for this function:

◆ ComputeDimensions()

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

Reimplemented from G4VSolid.

Definition at line 291 of file G4Tet.cc.

294 {
295 }

◆ CreatePolyhedron()

G4Polyhedron * G4Tet::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 843 of file G4Tet.cc.

844 {
845  G4Polyhedron *ph=new G4Polyhedron;
846  G4double xyz[4][3];
847  const G4int faces[4][4]={{1,3,2,0},{1,4,3,0},{1,2,4,0},{2,3,4,0}};
848  xyz[0][0]=fAnchor.x(); xyz[0][1]=fAnchor.y(); xyz[0][2]=fAnchor.z();
849  xyz[1][0]=fP2.x(); xyz[1][1]=fP2.y(); xyz[1][2]=fP2.z();
850  xyz[2][0]=fP3.x(); xyz[2][1]=fP3.y(); xyz[2][2]=fP3.z();
851  xyz[3][0]=fP4.x(); xyz[3][1]=fP4.y(); xyz[3][2]=fP4.z();
852 
853  ph->createPolyhedron(4,4,xyz,faces);
854 
855  return ph;
856 }
G4ThreeVector fP4
Definition: G4Tet.hh:165
G4ThreeVector fP2
Definition: G4Tet.hh:165
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
G4ThreeVector fP3
Definition: G4Tet.hh:165
int G4int
Definition: G4Types.hh:78
G4ThreeVector fAnchor
Definition: G4Tet.hh:165
double x() const
double y() const
double z() const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CreateRotatedVertices()

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

Definition at line 670 of file G4Tet.cc.

671 {
672  G4ThreeVectorList* vertices = new G4ThreeVectorList();
673 
674  if (vertices)
675  {
676  vertices->reserve(4);
677  G4ThreeVector vertex0(fAnchor);
678  G4ThreeVector vertex1(fP2);
679  G4ThreeVector vertex2(fP3);
680  G4ThreeVector vertex3(fP4);
681 
682  vertices->push_back(pTransform.TransformPoint(vertex0));
683  vertices->push_back(pTransform.TransformPoint(vertex1));
684  vertices->push_back(pTransform.TransformPoint(vertex2));
685  vertices->push_back(pTransform.TransformPoint(vertex3));
686  }
687  else
688  {
689  DumpInfo();
690  G4Exception("G4Tet::CreateRotatedVertices()",
691  "GeomSolids0003", FatalException,
692  "Error in allocation of vertices. Out of memory !");
693  }
694  return vertices;
695 }
G4ThreeVector fP4
Definition: G4Tet.hh:165
G4ThreeVector fP2
Definition: G4Tet.hh:165
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector fP3
Definition: G4Tet.hh:165
void DumpInfo() const
G4ThreeVector fAnchor
Definition: G4Tet.hh:165
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CVSFileVers()

const char* G4Tet::CVSFileVers ( )
inline

Definition at line 134 of file G4Tet.hh.

135  { return CVSVers; }
static const char CVSVers[]
Definition: G4Tet.hh:161

◆ CVSHeaderVers()

const char* G4Tet::CVSHeaderVers ( )
inline

Definition at line 132 of file G4Tet.hh.

133  { return "$Id: G4Tet.hh 83572 2014-09-01 15:23:27Z gcosmo $"; }

◆ DescribeYourselfTo()

void G4Tet::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 825 of file G4Tet.cc.

826 {
827  scene.AddSolid (*this);
828 }
virtual void AddSolid(const G4Box &)=0
Here is the call graph for this function:

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 487 of file G4Tet.cc.

489 {
490  G4ThreeVector vu(v.unit()), hp;
491  G4double vdotn, t, tmin=kInfinity;
492 
493  G4double extraDistance=10.0*fTol; // a little ways into the solid
494 
495  vdotn=-vu.dot(fNormal123);
496  if(vdotn > 1e-12)
497  { // this is a candidate face, since it is pointing at us
498  t=(p.dot(fNormal123)-fCdotN123)/vdotn; // # distance to intersection
499  if( (t>=-fTol) && (t<tmin) )
500  { // if not true, we're going away from this face or it's not close
501  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
502  if ( ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
503  ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
504  ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
505  {
506  tmin=t;
507  }
508  }
509  }
510 
511  vdotn=-vu.dot(fNormal134);
512  if(vdotn > 1e-12)
513  { // # this is a candidate face, since it is pointing at us
514  t=(p.dot(fNormal134)-fCdotN134)/vdotn; // # distance to intersection
515  if( (t>=-fTol) && (t<tmin) )
516  { // if not true, we're going away from this face
517  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
518  if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
519  ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
520  ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
521  {
522  tmin=t;
523  }
524  }
525  }
526 
527  vdotn=-vu.dot(fNormal142);
528  if(vdotn > 1e-12)
529  { // # this is a candidate face, since it is pointing at us
530  t=(p.dot(fNormal142)-fCdotN142)/vdotn; // # distance to intersection
531  if( (t>=-fTol) && (t<tmin) )
532  { // if not true, we're going away from this face
533  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
534  if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
535  ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
536  ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
537  {
538  tmin=t;
539  }
540  }
541  }
542 
543  vdotn=-vu.dot(fNormal234);
544  if(vdotn > 1e-12)
545  { // # this is a candidate face, since it is pointing at us
546  t=(p.dot(fNormal234)-fCdotN234)/vdotn; // # distance to intersection
547  if( (t>=-fTol) && (t<tmin) )
548  { // if not true, we're going away from this face
549  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
550  if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
551  ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
552  ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) )
553  {
554  tmin=t;
555  }
556  }
557  }
558 
559  return std::max(0.0,tmin);
560 }
G4double fTol
Definition: G4Tet.hh:172
G4ThreeVector fNormal142
Definition: G4Tet.hh:166
static const G4double kInfinity
Definition: geomdefs.hh:42
G4ThreeVector fNormal123
Definition: G4Tet.hh:166
G4double fCdotN142
Definition: G4Tet.hh:170
G4double fCdotN134
Definition: G4Tet.hh:170
Hep3Vector unit() const
double dot(const Hep3Vector &) const
G4double fCdotN234
Definition: G4Tet.hh:170
G4ThreeVector fNormal134
Definition: G4Tet.hh:166
G4double fCdotN123
Definition: G4Tet.hh:170
G4ThreeVector fNormal234
Definition: G4Tet.hh:166
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ DistanceToIn() [2/2]

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

Implements G4VSolid.

Definition at line 568 of file G4Tet.cc.

569 {
570  G4double dd=(p-fMiddle).mag() - fMaxSize - fTol;
571  return std::max(0.0, dd);
572 }
G4double fTol
Definition: G4Tet.hh:172
G4double fMaxSize
Definition: G4Tet.hh:172
G4ThreeVector fMiddle
Definition: G4Tet.hh:165
double G4double
Definition: G4Types.hh:76

◆ DistanceToOut() [1/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 580 of file G4Tet.cc.

583 {
584  G4ThreeVector vu(v.unit());
586 
587  vdotn=vu.dot(fNormal123);
588  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
589  {
590  t1=(fCdotN123-p.dot(fNormal123))/vdotn; // # distance to intersection
591  }
592 
593  vdotn=vu.dot(fNormal134);
594  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
595  {
596  t2=(fCdotN134-p.dot(fNormal134))/vdotn; // # distance to intersection
597  }
598 
599  vdotn=vu.dot(fNormal142);
600  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
601  {
602  t3=(fCdotN142-p.dot(fNormal142))/vdotn; // # distance to intersection
603  }
604 
605  vdotn=vu.dot(fNormal234);
606  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
607  {
608  t4=(fCdotN234-p.dot(fNormal234))/vdotn; // # distance to intersection
609  }
610 
611  tt=std::min(std::min(std::min(t1,t2),t3),t4);
612 
613  if (warningFlag && (tt == kInfinity || tt < -fTol))
614  {
615  DumpInfo();
616  std::ostringstream message;
617  message << "No good intersection found or already outside!?" << G4endl
618  << "p = " << p / mm << "mm" << G4endl
619  << "v = " << v << G4endl
620  << "t1, t2, t3, t4 (mm) "
621  << t1/mm << ", " << t2/mm << ", " << t3/mm << ", " << t4/mm;
622  G4Exception("G4Tet::DistanceToOut(p,v,...)", "GeomSolids1002",
623  JustWarning, message);
624  if(validNorm)
625  {
626  *validNorm=false; // flag normal as meaningless
627  }
628  }
629  else if(calcNorm && n)
630  {
632  if(tt==t1) { normal=fNormal123; }
633  else if (tt==t2) { normal=fNormal134; }
634  else if (tt==t3) { normal=fNormal142; }
635  else if (tt==t4) { normal=fNormal234; }
636  *n=normal;
637  if(validNorm) { *validNorm=true; }
638  }
639 
640  return std::max(tt,0.0); // avoid tt<0.0 by a tiny bit
641  // if we are right on a face
642 }
G4double fTol
Definition: G4Tet.hh:172
TTree * t1
Definition: plottest35.C:26
G4ThreeVector fNormal142
Definition: G4Tet.hh:166
static const G4double kInfinity
Definition: geomdefs.hh:42
G4ThreeVector fNormal123
Definition: G4Tet.hh:166
G4double fCdotN142
Definition: G4Tet.hh:170
G4double fCdotN134
Definition: G4Tet.hh:170
void DumpInfo() const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
Hep3Vector unit() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double dot(const Hep3Vector &) const
G4double fCdotN234
Definition: G4Tet.hh:170
G4ThreeVector fNormal134
Definition: G4Tet.hh:166
G4double fCdotN123
Definition: G4Tet.hh:170
TTree * t2
Definition: plottest35.C:36
#define G4endl
Definition: G4ios.hh:61
G4ThreeVector fNormal234
Definition: G4Tet.hh:166
double G4double
Definition: G4Types.hh:76
static const double mm
Definition: G4SIunits.hh:114
G4bool warningFlag
Definition: G4Tet.hh:168
Here is the call graph for this function:

◆ DistanceToOut() [2/2]

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

Implements G4VSolid.

Definition at line 649 of file G4Tet.cc.

650 {
651  G4double t1,t2,t3,t4;
652  t1=fCdotN123-p.dot(fNormal123); // distance to plane, positive if inside
653  t2=fCdotN134-p.dot(fNormal134); // distance to plane
654  t3=fCdotN142-p.dot(fNormal142); // distance to plane
655  t4=fCdotN234-p.dot(fNormal234); // distance to plane
656 
657  // if any one of these is negative, we are outside,
658  // so return zero in that case
659 
660  G4double tmin=std::min(std::min(std::min(t1,t2),t3),t4);
661  return (tmin < fTol)? 0:tmin;
662 }
G4double fTol
Definition: G4Tet.hh:172
TTree * t1
Definition: plottest35.C:26
G4ThreeVector fNormal142
Definition: G4Tet.hh:166
G4ThreeVector fNormal123
Definition: G4Tet.hh:166
G4double fCdotN142
Definition: G4Tet.hh:170
G4double fCdotN134
Definition: G4Tet.hh:170
double dot(const Hep3Vector &) const
G4double fCdotN234
Definition: G4Tet.hh:170
G4ThreeVector fNormal134
Definition: G4Tet.hh:166
G4double fCdotN123
Definition: G4Tet.hh:170
TTree * t2
Definition: plottest35.C:36
G4ThreeVector fNormal234
Definition: G4Tet.hh:166
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetCubicVolume()

G4double G4Tet::GetCubicVolume ( )
virtual

Reimplemented from G4VSolid.

Definition at line 805 of file G4Tet.cc.

806 {
807  return fCubicVolume;
808 }
G4double fCubicVolume
Definition: G4Tet.hh:154

◆ GetEntityType()

G4GeometryType G4Tet::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 701 of file G4Tet.cc.

◆ GetExtent()

G4VisExtent G4Tet::GetExtent ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 834 of file G4Tet.cc.

835 {
836  return G4VisExtent (fXMin, fXMax, fYMin, fYMax, fZMin, fZMax);
837 }
G4double fYMin
Definition: G4Tet.hh:171
G4double fXMax
Definition: G4Tet.hh:171
G4double fYMax
Definition: G4Tet.hh:171
G4double fZMin
Definition: G4Tet.hh:171
G4double fXMin
Definition: G4Tet.hh:171
G4double fZMax
Definition: G4Tet.hh:171

◆ GetPointOnFace()

G4ThreeVector G4Tet::GetPointOnFace ( G4ThreeVector  p1,
G4ThreeVector  p2,
G4ThreeVector  p3,
G4double area 
) const
private

Definition at line 748 of file G4Tet.cc.

750 {
751  G4double lambda1,lambda2;
752  G4ThreeVector v, w;
753 
754  v = p3 - p1;
755  w = p1 - p2;
756 
757  lambda1 = G4RandFlat::shoot(0.,1.);
758  lambda2 = G4RandFlat::shoot(0.,lambda1);
759 
760  area = 0.5*(v.cross(w)).mag();
761 
762  return (p2 + lambda1*w + lambda2*v);
763 }
ThreeVector shoot(const G4int Ap, const G4int Af)
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetPointOnSurface()

G4ThreeVector G4Tet::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 769 of file G4Tet.cc.

770 {
771  G4double chose,aOne,aTwo,aThree,aFour;
772  G4ThreeVector p1, p2, p3, p4;
773 
774  p1 = GetPointOnFace(fAnchor,fP2,fP3,aOne);
775  p2 = GetPointOnFace(fAnchor,fP4,fP3,aTwo);
776  p3 = GetPointOnFace(fAnchor,fP4,fP2,aThree);
777  p4 = GetPointOnFace(fP4,fP3,fP2,aFour);
778 
779  chose = G4RandFlat::shoot(0.,aOne+aTwo+aThree+aFour);
780  if( (chose>=0.) && (chose <aOne) ) {return p1;}
781  else if( (chose>=aOne) && (chose < aOne+aTwo) ) {return p2;}
782  else if( (chose>=aOne+aTwo) && (chose<aOne+aTwo+aThree) ) {return p3;}
783  return p4;
784 }
G4ThreeVector fP4
Definition: G4Tet.hh:165
G4ThreeVector GetPointOnFace(G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4double &area) const
Definition: G4Tet.cc:748
ThreeVector shoot(const G4int Ap, const G4int Af)
G4ThreeVector fP2
Definition: G4Tet.hh:165
G4ThreeVector fP3
Definition: G4Tet.hh:165
G4ThreeVector fAnchor
Definition: G4Tet.hh:165
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetPolyhedron()

G4Polyhedron * G4Tet::GetPolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 862 of file G4Tet.cc.

863 {
864  if (!fpPolyhedron ||
868  {
869  G4AutoLock l(&polyhedronMutex);
870  delete fpPolyhedron;
872  fRebuildPolyhedron = false;
873  l.unlock();
874  }
875  return fpPolyhedron;
876 }
G4bool fRebuildPolyhedron
Definition: G4Tet.hh:156
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tet.cc:843
static G4int GetNumberOfRotationSteps()
G4Polyhedron * fpPolyhedron
Definition: G4Tet.hh:157
Here is the call graph for this function:

◆ GetSurfaceArea()

G4double G4Tet::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 814 of file G4Tet.cc.

815 {
816  return fSurfaceArea;
817 }
G4double fSurfaceArea
Definition: G4Tet.hh:154

◆ GetVertices()

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

Definition at line 790 of file G4Tet.cc.

791 {
792  std::vector<G4ThreeVector> vertices(4);
793  vertices[0] = fAnchor;
794  vertices[1] = fP2;
795  vertices[2] = fP3;
796  vertices[3] = fP4;
797 
798  return vertices;
799 }
G4ThreeVector fP4
Definition: G4Tet.hh:165
G4ThreeVector fP2
Definition: G4Tet.hh:165
G4ThreeVector fP3
Definition: G4Tet.hh:165
G4ThreeVector fAnchor
Definition: G4Tet.hh:165
Here is the caller graph for this function:

◆ Inside()

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

Implements G4VSolid.

Definition at line 396 of file G4Tet.cc.

397 {
398  G4double r123, r134, r142, r234;
399 
400  // this is written to allow if-statement truncation so the outside test
401  // (where most of the world is) can fail very quickly and efficiently
402 
403  if ( (r123=p.dot(fNormal123)-fCdotN123) > fTol ||
404  (r134=p.dot(fNormal134)-fCdotN134) > fTol ||
405  (r142=p.dot(fNormal142)-fCdotN142) > fTol ||
406  (r234=p.dot(fNormal234)-fCdotN234) > fTol )
407  {
408  return kOutside; // at least one is out!
409  }
410  else if( (r123 < -fTol)&&(r134 < -fTol)&&(r142 < -fTol)&&(r234 < -fTol) )
411  {
412  return kInside; // all are definitively inside
413  }
414  else
415  {
416  return kSurface; // too close to tell
417  }
418 }
G4double fTol
Definition: G4Tet.hh:172
G4ThreeVector fNormal142
Definition: G4Tet.hh:166
G4ThreeVector fNormal123
Definition: G4Tet.hh:166
G4double fCdotN142
Definition: G4Tet.hh:170
G4double fCdotN134
Definition: G4Tet.hh:170
double dot(const Hep3Vector &) const
G4double fCdotN234
Definition: G4Tet.hh:170
G4ThreeVector fNormal134
Definition: G4Tet.hh:166
G4double fCdotN123
Definition: G4Tet.hh:170
G4ThreeVector fNormal234
Definition: G4Tet.hh:166
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ operator=()

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

Definition at line 242 of file G4Tet.cc.

243 {
244  // Check assignment to self
245  //
246  if (this == &rhs) { return *this; }
247 
248  // Copy base class data
249  //
250  G4VSolid::operator=(rhs);
251 
252  // Copy data
253  //
255  fAnchor = rhs.fAnchor;
256  fP2 = rhs.fP2; fP3 = rhs.fP3; fP4 = rhs.fP4; fMiddle = rhs.fMiddle;
260  fCdotN142 = rhs.fCdotN142; fCdotN134 = rhs.fCdotN134;
261  fCdotN234 = rhs.fCdotN234; fXMin = rhs.fXMin; fXMax = rhs.fXMax;
262  fYMin = rhs.fYMin; fYMax = rhs.fYMax; fZMin = rhs.fZMin; fZMax = rhs.fZMax;
263  fDx = rhs.fDx; fDy = rhs.fDy; fDz = rhs.fDz; fTol = rhs.fTol;
264  fMaxSize = rhs.fMaxSize;
265  fRebuildPolyhedron = false;
266  delete fpPolyhedron; fpPolyhedron = 0;
267 
268  return *this;
269 }
G4ThreeVector fP4
Definition: G4Tet.hh:165
G4double fTol
Definition: G4Tet.hh:172
G4double fYMin
Definition: G4Tet.hh:171
G4double fXMax
Definition: G4Tet.hh:171
G4ThreeVector fNormal142
Definition: G4Tet.hh:166
G4bool fRebuildPolyhedron
Definition: G4Tet.hh:156
G4double fDz
Definition: G4Tet.hh:172
G4ThreeVector fP2
Definition: G4Tet.hh:165
G4double fYMax
Definition: G4Tet.hh:171
G4ThreeVector fNormal123
Definition: G4Tet.hh:166
G4double fCdotN142
Definition: G4Tet.hh:170
G4ThreeVector fP3
Definition: G4Tet.hh:165
G4double fCdotN134
Definition: G4Tet.hh:170
G4double fMaxSize
Definition: G4Tet.hh:172
G4ThreeVector fAnchor
Definition: G4Tet.hh:165
G4double fZMin
Definition: G4Tet.hh:171
G4double fXMin
Definition: G4Tet.hh:171
G4ThreeVector fMiddle
Definition: G4Tet.hh:165
G4double fSurfaceArea
Definition: G4Tet.hh:154
G4double fCdotN234
Definition: G4Tet.hh:170
G4double fDy
Definition: G4Tet.hh:172
G4ThreeVector fNormal134
Definition: G4Tet.hh:166
G4double fCdotN123
Definition: G4Tet.hh:170
G4double fCubicVolume
Definition: G4Tet.hh:154
G4double fDx
Definition: G4Tet.hh:172
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
G4ThreeVector fNormal234
Definition: G4Tet.hh:166
G4Polyhedron * fpPolyhedron
Definition: G4Tet.hh:157
G4bool warningFlag
Definition: G4Tet.hh:168
G4double fZMax
Definition: G4Tet.hh:171
Here is the call graph for this function:

◆ PrintWarnings()

void G4Tet::PrintWarnings ( G4bool  flag)
inline

Definition at line 136 of file G4Tet.hh.

137  { warningFlag=flag; }
G4bool warningFlag
Definition: G4Tet.hh:168
Here is the call graph for this function:

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 719 of file G4Tet.cc.

720 {
721  G4int oldprc = os.precision(16);
722  os << "-----------------------------------------------------------\n"
723  << " *** Dump for solid - " << GetName() << " ***\n"
724  << " ===================================================\n"
725  << " Solid type: G4Tet\n"
726  << " Parameters: \n"
727  << " anchor: " << fAnchor/mm << " mm \n"
728  << " p2: " << fP2/mm << " mm \n"
729  << " p3: " << fP3/mm << " mm \n"
730  << " p4: " << fP4/mm << " mm \n"
731  << " normal123: " << fNormal123 << " \n"
732  << " normal134: " << fNormal134 << " \n"
733  << " normal142: " << fNormal142 << " \n"
734  << " normal234: " << fNormal234 << " \n"
735  << "-----------------------------------------------------------\n";
736  os.precision(oldprc);
737 
738  return os;
739 }
G4ThreeVector fP4
Definition: G4Tet.hh:165
G4ThreeVector fNormal142
Definition: G4Tet.hh:166
G4ThreeVector fP2
Definition: G4Tet.hh:165
G4ThreeVector fNormal123
Definition: G4Tet.hh:166
G4ThreeVector fP3
Definition: G4Tet.hh:165
int G4int
Definition: G4Types.hh:78
G4String GetName() const
G4ThreeVector fAnchor
Definition: G4Tet.hh:165
G4ThreeVector fNormal134
Definition: G4Tet.hh:166
G4ThreeVector fNormal234
Definition: G4Tet.hh:166
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:

◆ SurfaceNormal()

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

Implements G4VSolid.

Definition at line 427 of file G4Tet.cc.

428 {
429  G4double r123=std::fabs(p.dot(fNormal123)-fCdotN123);
430  G4double r134=std::fabs(p.dot(fNormal134)-fCdotN134);
431  G4double r142=std::fabs(p.dot(fNormal142)-fCdotN142);
432  G4double r234=std::fabs(p.dot(fNormal234)-fCdotN234);
433 
434  const G4double delta = 0.5*kCarTolerance;
435  G4ThreeVector sumnorm(0., 0., 0.);
436  G4int noSurfaces=0;
437 
438  if (r123 <= delta)
439  {
440  noSurfaces ++;
441  sumnorm= fNormal123;
442  }
443 
444  if (r134 <= delta)
445  {
446  noSurfaces ++;
447  sumnorm += fNormal134;
448  }
449 
450  if (r142 <= delta)
451  {
452  noSurfaces ++;
453  sumnorm += fNormal142;
454  }
455  if (r234 <= delta)
456  {
457  noSurfaces ++;
458  sumnorm += fNormal234;
459  }
460 
461  if( noSurfaces > 0 )
462  {
463  if( noSurfaces == 1 )
464  {
465  return sumnorm;
466  }
467  else
468  {
469  return sumnorm.unit();
470  }
471  }
472  else // Approximative Surface Normal
473  {
474 
475  if( (r123<=r134) && (r123<=r142) && (r123<=r234) ) { return fNormal123; }
476  else if ( (r134<=r142) && (r134<=r234) ) { return fNormal134; }
477  else if (r142 <= r234) { return fNormal142; }
478  return fNormal234;
479  }
480 }
G4ThreeVector fNormal142
Definition: G4Tet.hh:166
G4ThreeVector fNormal123
Definition: G4Tet.hh:166
G4double fCdotN142
Definition: G4Tet.hh:170
G4double fCdotN134
Definition: G4Tet.hh:170
int G4int
Definition: G4Types.hh:78
Hep3Vector unit() const
double dot(const Hep3Vector &) const
G4double fCdotN234
Definition: G4Tet.hh:170
G4ThreeVector fNormal134
Definition: G4Tet.hh:166
G4double fCdotN123
Definition: G4Tet.hh:170
G4double kCarTolerance
Definition: G4VSolid.hh:304
G4ThreeVector fNormal234
Definition: G4Tet.hh:166
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

Member Data Documentation

◆ CVSVers

const char G4Tet::CVSVers ="$Id: G4Tet.cc 102297 2017-01-20 13:33:54Z gcosmo $"
staticprivate

Definition at line 161 of file G4Tet.hh.

◆ fAnchor

G4ThreeVector G4Tet::fAnchor
private

Definition at line 165 of file G4Tet.hh.

◆ fCdotN123

G4double G4Tet::fCdotN123
private

Definition at line 170 of file G4Tet.hh.

◆ fCdotN134

G4double G4Tet::fCdotN134
private

Definition at line 170 of file G4Tet.hh.

◆ fCdotN142

G4double G4Tet::fCdotN142
private

Definition at line 170 of file G4Tet.hh.

◆ fCdotN234

G4double G4Tet::fCdotN234
private

Definition at line 170 of file G4Tet.hh.

◆ fCubicVolume

G4double G4Tet::fCubicVolume
private

Definition at line 154 of file G4Tet.hh.

◆ fDx

G4double G4Tet::fDx
private

Definition at line 172 of file G4Tet.hh.

◆ fDy

G4double G4Tet::fDy
private

Definition at line 172 of file G4Tet.hh.

◆ fDz

G4double G4Tet::fDz
private

Definition at line 172 of file G4Tet.hh.

◆ fMaxSize

G4double G4Tet::fMaxSize
private

Definition at line 172 of file G4Tet.hh.

◆ fMiddle

G4ThreeVector G4Tet::fMiddle
private

Definition at line 165 of file G4Tet.hh.

◆ fNormal123

G4ThreeVector G4Tet::fNormal123
private

Definition at line 166 of file G4Tet.hh.

◆ fNormal134

G4ThreeVector G4Tet::fNormal134
private

Definition at line 166 of file G4Tet.hh.

◆ fNormal142

G4ThreeVector G4Tet::fNormal142
private

Definition at line 166 of file G4Tet.hh.

◆ fNormal234

G4ThreeVector G4Tet::fNormal234
private

Definition at line 166 of file G4Tet.hh.

◆ fP2

G4ThreeVector G4Tet::fP2
private

Definition at line 165 of file G4Tet.hh.

◆ fP3

G4ThreeVector G4Tet::fP3
private

Definition at line 165 of file G4Tet.hh.

◆ fP4

G4ThreeVector G4Tet::fP4
private

Definition at line 165 of file G4Tet.hh.

◆ fpPolyhedron

G4Polyhedron* G4Tet::fpPolyhedron
mutableprivate

Definition at line 157 of file G4Tet.hh.

◆ fRebuildPolyhedron

G4bool G4Tet::fRebuildPolyhedron
mutableprivate

Definition at line 156 of file G4Tet.hh.

◆ fSurfaceArea

G4double G4Tet::fSurfaceArea
private

Definition at line 154 of file G4Tet.hh.

◆ fTol

G4double G4Tet::fTol
private

Definition at line 172 of file G4Tet.hh.

◆ fXMax

G4double G4Tet::fXMax
private

Definition at line 171 of file G4Tet.hh.

◆ fXMin

G4double G4Tet::fXMin
private

Definition at line 171 of file G4Tet.hh.

◆ fYMax

G4double G4Tet::fYMax
private

Definition at line 171 of file G4Tet.hh.

◆ fYMin

G4double G4Tet::fYMin
private

Definition at line 171 of file G4Tet.hh.

◆ fZMax

G4double G4Tet::fZMax
private

Definition at line 171 of file G4Tet.hh.

◆ fZMin

G4double G4Tet::fZMin
private

Definition at line 171 of file G4Tet.hh.

◆ warningFlag

G4bool G4Tet::warningFlag
private

Definition at line 168 of file G4Tet.hh.


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