Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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)
 
void Extent (G4ThreeVector &pMin, G4ThreeVector &pMax) const
 
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)
 

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 66 of file G4Tet.hh.

Constructor & Destructor Documentation

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

Definition at line 99 of file G4Tet.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

G4Tet::~G4Tet ( )
virtual

Definition at line 215 of file G4Tet.cc.

216 {
217  delete fpPolyhedron; fpPolyhedron = 0;
218 }
G4Tet::G4Tet ( __void__ &  a)

Definition at line 199 of file G4Tet.cc.

200  : G4VSolid(a), fCubicVolume(0.), fSurfaceArea(0.),
201  fRebuildPolyhedron(false), fpPolyhedron(0),
202  fAnchor(0,0,0), fP2(0,0,0), fP3(0,0,0), fP4(0,0,0), fMiddle(0,0,0),
203  fNormal123(0,0,0), fNormal142(0,0,0), fNormal134(0,0,0),
204  fNormal234(0,0,0), warningFlag(0),
205  fCdotN123(0.), fCdotN142(0.), fCdotN134(0.), fCdotN234(0.),
206  fXMin(0.), fXMax(0.), fYMin(0.), fYMax(0.), fZMin(0.), fZMax(0.),
207  fDx(0.), fDy(0.), fDz(0.), fTol(0.), fMaxSize(0.)
208 {
209 }
G4VSolid(const G4String &name)
Definition: G4VSolid.cc:61
G4Tet::G4Tet ( const G4Tet rhs)

Definition at line 224 of file G4Tet.cc.

225  : G4VSolid(rhs),
226  fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
227  fRebuildPolyhedron(false), fpPolyhedron(0), fAnchor(rhs.fAnchor),
228  fP2(rhs.fP2), fP3(rhs.fP3), fP4(rhs.fP4), fMiddle(rhs.fMiddle),
229  fNormal123(rhs.fNormal123), fNormal142(rhs.fNormal142),
230  fNormal134(rhs.fNormal134), fNormal234(rhs.fNormal234),
231  warningFlag(rhs.warningFlag), fCdotN123(rhs.fCdotN123),
232  fCdotN142(rhs.fCdotN142), fCdotN134(rhs.fCdotN134),
233  fCdotN234(rhs.fCdotN234), fXMin(rhs.fXMin), fXMax(rhs.fXMax),
234  fYMin(rhs.fYMin), fYMax(rhs.fYMax), fZMin(rhs.fZMin), fZMax(rhs.fZMax),
235  fDx(rhs.fDx), fDy(rhs.fDy), fDz(rhs.fDz), fTol(rhs.fTol),
236  fMaxSize(rhs.fMaxSize)
237 {
238 }
G4VSolid(const G4String &name)
Definition: G4VSolid.cc:61

Member Function Documentation

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

Implements G4VSolid.

Definition at line 327 of file G4Tet.cc.

331 {
332  G4ThreeVector bmin, bmax;
333  G4bool exist;
334 
335  // Check bounding box (bbox)
336  //
337  Extent(bmin,bmax);
338  G4BoundingEnvelope bbox(bmin,bmax);
339 #ifdef G4BBOX_EXTENT
340  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
341 #endif
342  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
343  {
344  return exist = (pMin < pMax) ? true : false;
345  }
346 
347  // Set bounding envelope (benv) and calculate extent
348  //
349  std::vector<G4ThreeVector> vec = GetVertices();
350 
351  G4ThreeVectorList anchor(1);
352  anchor[0].set(vec[0].x(),vec[0].y(),vec[0].z());
353 
355  base[0].set(vec[1].x(),vec[1].y(),vec[1].z());
356  base[1].set(vec[2].x(),vec[2].y(),vec[2].z());
357  base[2].set(vec[3].x(),vec[3].y(),vec[3].z());
358 
359  std::vector<const G4ThreeVectorList *> polygons(2);
360  polygons[0] = &anchor;
361  polygons[1] = &base;
362 
363  G4BoundingEnvelope benv(bmin,bmax,polygons);
364  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
365  return exist;
366 }
const XML_Char int const XML_Char int const XML_Char * base
Definition: expat.h:331
std::vector< G4ThreeVector > GetVertices() const
Definition: G4Tet.cc:733
bool G4bool
Definition: G4Types.hh:79
void Extent(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Tet.cc:304
std::vector< G4ThreeVector > G4ThreeVectorList

Here is the call graph for this function:

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

Definition at line 278 of file G4Tet.cc.

282 {
283  G4bool result;
284  G4Tet *object=new G4Tet("temp",anchor,p2,p3,p4,&result);
285  delete object;
286  return result;
287 }
G4double G4ParticleHPJENDLHEData::G4double result
Definition: G4Tet.hh:66
bool G4bool
Definition: G4Types.hh:79
G4Tet(const G4String &pName, G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4, G4bool *degeneracyFlag=0)
Definition: G4Tet.cc:99

Here is the call graph for this function:

G4VSolid * G4Tet::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 653 of file G4Tet.cc.

654 {
655  return new G4Tet(*this);
656 }
G4Tet(const G4String &pName, G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4, G4bool *degeneracyFlag=0)
Definition: G4Tet.cc:99

Here is the call graph for this function:

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

Reimplemented from G4VSolid.

Definition at line 294 of file G4Tet.cc.

297 {
298 }
G4Polyhedron * G4Tet::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 786 of file G4Tet.cc.

787 {
788  G4Polyhedron *ph=new G4Polyhedron;
789  G4double xyz[4][3];
790  const G4int faces[4][4]={{1,3,2,0},{1,4,3,0},{1,2,4,0},{2,3,4,0}};
791  xyz[0][0]=fAnchor.x(); xyz[0][1]=fAnchor.y(); xyz[0][2]=fAnchor.z();
792  xyz[1][0]=fP2.x(); xyz[1][1]=fP2.y(); xyz[1][2]=fP2.z();
793  xyz[2][0]=fP3.x(); xyz[2][1]=fP3.y(); xyz[2][2]=fP3.z();
794  xyz[3][0]=fP4.x(); xyz[3][1]=fP4.y(); xyz[3][2]=fP4.z();
795 
796  ph->createPolyhedron(4,4,xyz,faces);
797 
798  return ph;
799 }
double x() const
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
int G4int
Definition: G4Types.hh:78
double z() const
double y() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

const char* G4Tet::CVSFileVers ( )
inline

Definition at line 137 of file G4Tet.hh.

138  { return CVSVers; }
const char* G4Tet::CVSHeaderVers ( )
inline

Definition at line 135 of file G4Tet.hh.

136  { return "$Id: G4Tet.hh 99781 2016-10-05 10:18:54Z gcosmo $"; }
void G4Tet::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 768 of file G4Tet.cc.

769 {
770  scene.AddSolid (*this);
771 }
virtual void AddSolid(const G4Box &)=0

Here is the call graph for this function:

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

Implements G4VSolid.

Definition at line 463 of file G4Tet.cc.

465 {
466  G4ThreeVector vu(v.unit()), hp;
467  G4double vdotn, t, tmin=kInfinity;
468 
469  G4double extraDistance=10.0*fTol; // a little ways into the solid
470 
471  vdotn=-vu.dot(fNormal123);
472  if(vdotn > 1e-12)
473  { // this is a candidate face, since it is pointing at us
474  t=(p.dot(fNormal123)-fCdotN123)/vdotn; // # distance to intersection
475  if( (t>=-fTol) && (t<tmin) )
476  { // if not true, we're going away from this face or it's not close
477  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
478  if ( ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
479  ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
480  ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
481  {
482  tmin=t;
483  }
484  }
485  }
486 
487  vdotn=-vu.dot(fNormal134);
488  if(vdotn > 1e-12)
489  { // # this is a candidate face, since it is pointing at us
490  t=(p.dot(fNormal134)-fCdotN134)/vdotn; // # distance to intersection
491  if( (t>=-fTol) && (t<tmin) )
492  { // if not true, we're going away from this face
493  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
494  if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
495  ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
496  ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
497  {
498  tmin=t;
499  }
500  }
501  }
502 
503  vdotn=-vu.dot(fNormal142);
504  if(vdotn > 1e-12)
505  { // # this is a candidate face, since it is pointing at us
506  t=(p.dot(fNormal142)-fCdotN142)/vdotn; // # distance to intersection
507  if( (t>=-fTol) && (t<tmin) )
508  { // if not true, we're going away from this face
509  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
510  if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
511  ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
512  ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
513  {
514  tmin=t;
515  }
516  }
517  }
518 
519  vdotn=-vu.dot(fNormal234);
520  if(vdotn > 1e-12)
521  { // # this is a candidate face, since it is pointing at us
522  t=(p.dot(fNormal234)-fCdotN234)/vdotn; // # distance to intersection
523  if( (t>=-fTol) && (t<tmin) )
524  { // if not true, we're going away from this face
525  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
526  if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
527  ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
528  ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) )
529  {
530  tmin=t;
531  }
532  }
533  }
534 
535  return std::max(0.0,tmin);
536 }
static const G4double kInfinity
Definition: geomdefs.hh:42
double dot(const Hep3Vector &) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Hep3Vector unit() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VSolid.

Definition at line 544 of file G4Tet.cc.

545 {
546  G4double dd=(p-fMiddle).mag() - fMaxSize - fTol;
547  return std::max(0.0, dd);
548 }
T max(const T t1, const T t2)
brief Return the largest of the two arguments
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

559 {
560  G4ThreeVector vu(v.unit());
561  G4double t1=kInfinity,t2=kInfinity,t3=kInfinity,t4=kInfinity, vdotn, tt;
562 
563  vdotn=vu.dot(fNormal123);
564  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
565  {
566  t1=(fCdotN123-p.dot(fNormal123))/vdotn; // # distance to intersection
567  }
568 
569  vdotn=vu.dot(fNormal134);
570  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
571  {
572  t2=(fCdotN134-p.dot(fNormal134))/vdotn; // # distance to intersection
573  }
574 
575  vdotn=vu.dot(fNormal142);
576  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
577  {
578  t3=(fCdotN142-p.dot(fNormal142))/vdotn; // # distance to intersection
579  }
580 
581  vdotn=vu.dot(fNormal234);
582  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
583  {
584  t4=(fCdotN234-p.dot(fNormal234))/vdotn; // # distance to intersection
585  }
586 
587  tt=std::min(std::min(std::min(t1,t2),t3),t4);
588 
589  if (warningFlag && (tt == kInfinity || tt < -fTol))
590  {
591  DumpInfo();
592  std::ostringstream message;
593  message << "No good intersection found or already outside!?" << G4endl
594  << "p = " << p / mm << "mm" << G4endl
595  << "v = " << v << G4endl
596  << "t1, t2, t3, t4 (mm) "
597  << t1/mm << ", " << t2/mm << ", " << t3/mm << ", " << t4/mm;
598  G4Exception("G4Tet::DistanceToOut(p,v,...)", "GeomSolids1002",
599  JustWarning, message);
600  if(validNorm)
601  {
602  *validNorm=false; // flag normal as meaningless
603  }
604  }
605  else if(calcNorm && n)
606  {
608  if(tt==t1) { normal=fNormal123; }
609  else if (tt==t2) { normal=fNormal134; }
610  else if (tt==t3) { normal=fNormal142; }
611  else if (tt==t4) { normal=fNormal234; }
612  *n=normal;
613  if(validNorm) { *validNorm=true; }
614  }
615 
616  return std::max(tt,0.0); // avoid tt<0.0 by a tiny bit
617  // if we are right on a face
618 }
static constexpr double mm
Definition: G4SIunits.hh:115
static const G4double kInfinity
Definition: geomdefs.hh:42
double dot(const Hep3Vector &) const
void DumpInfo() const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Hep3Vector unit() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VSolid.

Definition at line 625 of file G4Tet.cc.

626 {
627  G4double t1,t2,t3,t4;
628  t1=fCdotN123-p.dot(fNormal123); // distance to plane, positive if inside
629  t2=fCdotN134-p.dot(fNormal134); // distance to plane
630  t3=fCdotN142-p.dot(fNormal142); // distance to plane
631  t4=fCdotN234-p.dot(fNormal234); // distance to plane
632 
633  // if any one of these is negative, we are outside,
634  // so return zero in that case
635 
636  G4double tmin=std::min(std::min(std::min(t1,t2),t3),t4);
637  return (tmin < fTol)? 0:tmin;
638 }
double dot(const Hep3Vector &) const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Reimplemented from G4VSolid.

Definition at line 304 of file G4Tet.cc.

305 {
306  pMin.set(fXMin,fYMin,fZMin);
307  pMax.set(fXMax,fYMax,fZMax);
308 
309  // Check correctness of the bounding box
310  //
311  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
312  {
313  std::ostringstream message;
314  message << "Bad bounding box (min >= max) for solid: "
315  << GetName() << " !"
316  << "\npMin = " << pMin
317  << "\npMax = " << pMax;
318  G4Exception("G4Tet::Extent()", "GeomMgt0001", JustWarning, message);
319  DumpInfo();
320  }
321 }
void set(double x, double y, double z)
G4String GetName() const
double x() const
double z() const
void DumpInfo() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double y() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4Tet::GetCubicVolume ( )
virtual

Reimplemented from G4VSolid.

Definition at line 748 of file G4Tet.cc.

749 {
750  return fCubicVolume;
751 }
G4GeometryType G4Tet::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 644 of file G4Tet.cc.

645 {
646  return G4String("G4Tet");
647 }
G4VisExtent G4Tet::GetExtent ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 777 of file G4Tet.cc.

778 {
779  return G4VisExtent (fXMin, fXMax, fYMin, fYMax, fZMin, fZMax);
780 }
G4ThreeVector G4Tet::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 712 of file G4Tet.cc.

713 {
714  G4double chose,aOne,aTwo,aThree,aFour;
715  G4ThreeVector p1, p2, p3, p4;
716 
717  p1 = GetPointOnFace(fAnchor,fP2,fP3,aOne);
718  p2 = GetPointOnFace(fAnchor,fP4,fP3,aTwo);
719  p3 = GetPointOnFace(fAnchor,fP4,fP2,aThree);
720  p4 = GetPointOnFace(fP4,fP3,fP2,aFour);
721 
722  chose = G4RandFlat::shoot(0.,aOne+aTwo+aThree+aFour);
723  if( (chose>=0.) && (chose <aOne) ) {return p1;}
724  else if( (chose>=aOne) && (chose < aOne+aTwo) ) {return p2;}
725  else if( (chose>=aOne+aTwo) && (chose<aOne+aTwo+aThree) ) {return p3;}
726  return p4;
727 }
ThreeVector shoot(const G4int Ap, const G4int Af)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4Polyhedron * G4Tet::GetPolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 805 of file G4Tet.cc.

806 {
807  if (!fpPolyhedron ||
808  fRebuildPolyhedron ||
810  fpPolyhedron->GetNumberOfRotationSteps())
811  {
812  G4AutoLock l(&polyhedronMutex);
813  delete fpPolyhedron;
814  fpPolyhedron = CreatePolyhedron();
815  fRebuildPolyhedron = false;
816  l.unlock();
817  }
818  return fpPolyhedron;
819 }
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tet.cc:786
static G4int GetNumberOfRotationSteps()
G4int GetNumberOfRotationStepsAtTimeOfCreation() const

Here is the call graph for this function:

G4double G4Tet::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 757 of file G4Tet.cc.

758 {
759  return fSurfaceArea;
760 }
std::vector< G4ThreeVector > G4Tet::GetVertices ( ) const

Definition at line 733 of file G4Tet.cc.

734 {
735  std::vector<G4ThreeVector> vertices(4);
736  vertices[0] = fAnchor;
737  vertices[1] = fP2;
738  vertices[2] = fP3;
739  vertices[3] = fP4;
740 
741  return vertices;
742 }

Here is the caller graph for this function:

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

Implements G4VSolid.

Definition at line 372 of file G4Tet.cc.

373 {
374  G4double r123, r134, r142, r234;
375 
376  // this is written to allow if-statement truncation so the outside test
377  // (where most of the world is) can fail very quickly and efficiently
378 
379  if ( (r123=p.dot(fNormal123)-fCdotN123) > fTol ||
380  (r134=p.dot(fNormal134)-fCdotN134) > fTol ||
381  (r142=p.dot(fNormal142)-fCdotN142) > fTol ||
382  (r234=p.dot(fNormal234)-fCdotN234) > fTol )
383  {
384  return kOutside; // at least one is out!
385  }
386  else if( (r123 < -fTol)&&(r134 < -fTol)&&(r142 < -fTol)&&(r234 < -fTol) )
387  {
388  return kInside; // all are definitively inside
389  }
390  else
391  {
392  return kSurface; // too close to tell
393  }
394 }
double dot(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Definition at line 245 of file G4Tet.cc.

246 {
247  // Check assignment to self
248  //
249  if (this == &rhs) { return *this; }
250 
251  // Copy base class data
252  //
253  G4VSolid::operator=(rhs);
254 
255  // Copy data
256  //
257  fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
258  fAnchor = rhs.fAnchor;
259  fP2 = rhs.fP2; fP3 = rhs.fP3; fP4 = rhs.fP4; fMiddle = rhs.fMiddle;
260  fNormal123 = rhs.fNormal123; fNormal142 = rhs.fNormal142;
261  fNormal134 = rhs.fNormal134; fNormal234 = rhs.fNormal234;
262  warningFlag = rhs.warningFlag; fCdotN123 = rhs.fCdotN123;
263  fCdotN142 = rhs.fCdotN142; fCdotN134 = rhs.fCdotN134;
264  fCdotN234 = rhs.fCdotN234; fXMin = rhs.fXMin; fXMax = rhs.fXMax;
265  fYMin = rhs.fYMin; fYMax = rhs.fYMax; fZMin = rhs.fZMin; fZMax = rhs.fZMax;
266  fDx = rhs.fDx; fDy = rhs.fDy; fDz = rhs.fDz; fTol = rhs.fTol;
267  fMaxSize = rhs.fMaxSize;
268  fRebuildPolyhedron = false;
269  delete fpPolyhedron; fpPolyhedron = 0;
270 
271  return *this;
272 }
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:111

Here is the call graph for this function:

void G4Tet::PrintWarnings ( G4bool  flag)
inline

Definition at line 139 of file G4Tet.hh.

140  { warningFlag=flag; }
std::ostream & G4Tet::StreamInfo ( std::ostream &  os) const
virtual

Implements G4VSolid.

Definition at line 662 of file G4Tet.cc.

663 {
664  G4int oldprc = os.precision(16);
665  os << "-----------------------------------------------------------\n"
666  << " *** Dump for solid - " << GetName() << " ***\n"
667  << " ===================================================\n"
668  << " Solid type: G4Tet\n"
669  << " Parameters: \n"
670  << " anchor: " << fAnchor/mm << " mm \n"
671  << " p2: " << fP2/mm << " mm \n"
672  << " p3: " << fP3/mm << " mm \n"
673  << " p4: " << fP4/mm << " mm \n"
674  << " normal123: " << fNormal123 << " \n"
675  << " normal134: " << fNormal134 << " \n"
676  << " normal142: " << fNormal142 << " \n"
677  << " normal234: " << fNormal234 << " \n"
678  << "-----------------------------------------------------------\n";
679  os.precision(oldprc);
680 
681  return os;
682 }
G4String GetName() const
static constexpr double mm
Definition: G4SIunits.hh:115
int G4int
Definition: G4Types.hh:78

Here is the call graph for this function:

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

Implements G4VSolid.

Definition at line 403 of file G4Tet.cc.

404 {
405  G4double r123=std::fabs(p.dot(fNormal123)-fCdotN123);
406  G4double r134=std::fabs(p.dot(fNormal134)-fCdotN134);
407  G4double r142=std::fabs(p.dot(fNormal142)-fCdotN142);
408  G4double r234=std::fabs(p.dot(fNormal234)-fCdotN234);
409 
410  const G4double delta = 0.5*kCarTolerance;
411  G4ThreeVector sumnorm(0., 0., 0.);
412  G4int noSurfaces=0;
413 
414  if (r123 <= delta)
415  {
416  noSurfaces ++;
417  sumnorm= fNormal123;
418  }
419 
420  if (r134 <= delta)
421  {
422  noSurfaces ++;
423  sumnorm += fNormal134;
424  }
425 
426  if (r142 <= delta)
427  {
428  noSurfaces ++;
429  sumnorm += fNormal142;
430  }
431  if (r234 <= delta)
432  {
433  noSurfaces ++;
434  sumnorm += fNormal234;
435  }
436 
437  if( noSurfaces > 0 )
438  {
439  if( noSurfaces == 1 )
440  {
441  return sumnorm;
442  }
443  else
444  {
445  return sumnorm.unit();
446  }
447  }
448  else // Approximative Surface Normal
449  {
450 
451  if( (r123<=r134) && (r123<=r142) && (r123<=r234) ) { return fNormal123; }
452  else if ( (r134<=r142) && (r134<=r234) ) { return fNormal134; }
453  else if (r142 <= r234) { return fNormal142; }
454  return fNormal234;
455  }
456 }
double dot(const Hep3Vector &) const
int G4int
Definition: G4Types.hh:78
Hep3Vector unit() const
G4double kCarTolerance
Definition: G4VSolid.hh:307
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:


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