File indexing completed on 2025-11-05 09:42:42
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039 #ifndef G4TET_HH
0040 #define G4TET_HH
0041
0042 #include "G4GeomTypes.hh"
0043
0044 #if defined(G4GEOM_USE_USOLIDS)
0045 #define G4GEOM_USE_UTET 1
0046 #endif
0047
0048 #if defined(G4GEOM_USE_UTET)
0049 #define G4UTet G4Tet
0050 #include "G4UTet.hh"
0051 #else
0052
0053 #include "G4VSolid.hh"
0054
0055 class G4Tet : public G4VSolid
0056 {
0057 public:
0058
0059 G4Tet(const G4String& pName,
0060 const G4ThreeVector& anchor,
0061 const G4ThreeVector& p1,
0062 const G4ThreeVector& p2,
0063 const G4ThreeVector& p3,
0064 G4bool* degeneracyFlag = nullptr);
0065
0066 ~G4Tet() override;
0067
0068 void SetVertices(const G4ThreeVector& anchor,
0069 const G4ThreeVector& p1,
0070 const G4ThreeVector& p2,
0071 const G4ThreeVector& p3,
0072 G4bool* degeneracyFlag = nullptr);
0073
0074
0075 void GetVertices(G4ThreeVector& anchor,
0076 G4ThreeVector& p1,
0077 G4ThreeVector& p2,
0078 G4ThreeVector& p3) const;
0079 std::vector<G4ThreeVector> GetVertices() const;
0080
0081
0082 inline void PrintWarnings(G4bool) {};
0083
0084
0085 G4bool CheckDegeneracy(const G4ThreeVector& p0,
0086 const G4ThreeVector& p1,
0087 const G4ThreeVector& p2,
0088 const G4ThreeVector& p3) const;
0089
0090
0091 void ComputeDimensions(G4VPVParameterisation* p,
0092 const G4int n,
0093 const G4VPhysicalVolume* pRep) override;
0094
0095 void SetBoundingLimits(const G4ThreeVector& pMin, const G4ThreeVector& pMax);
0096 void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const override;
0097 G4bool CalculateExtent(const EAxis pAxis,
0098 const G4VoxelLimits& pVoxelLimit,
0099 const G4AffineTransform& pTransform,
0100 G4double& pmin, G4double& pmax) const override;
0101
0102 EInside Inside(const G4ThreeVector& p) const override;
0103 G4ThreeVector SurfaceNormal( const G4ThreeVector& p) const override;
0104 G4double DistanceToIn(const G4ThreeVector& p,
0105 const G4ThreeVector& v) const override;
0106 G4double DistanceToIn(const G4ThreeVector& p) const override;
0107 G4double DistanceToOut(const G4ThreeVector& p,
0108 const G4ThreeVector& v,
0109 const G4bool calcNorm = false,
0110 G4bool* validNorm = nullptr,
0111 G4ThreeVector* n = nullptr) const override;
0112 G4double DistanceToOut(const G4ThreeVector& p) const override;
0113
0114 G4GeometryType GetEntityType() const override;
0115
0116 G4bool IsFaceted () const override;
0117
0118 G4VSolid* Clone() const override;
0119
0120 std::ostream& StreamInfo(std::ostream& os) const override;
0121
0122 G4double GetCubicVolume() override;
0123 G4double GetSurfaceArea() override;
0124
0125 G4ThreeVector GetPointOnSurface() const override;
0126
0127
0128 void DescribeYourselfTo (G4VGraphicsScene& scene) const override;
0129 G4VisExtent GetExtent () const override;
0130 G4Polyhedron* CreatePolyhedron () const override;
0131 G4Polyhedron* GetPolyhedron () const override;
0132
0133
0134
0135
0136 G4Tet(__void__&);
0137
0138
0139 G4Tet(const G4Tet& rhs);
0140
0141
0142 G4Tet& operator=(const G4Tet& rhs);
0143
0144 private:
0145
0146
0147 void Initialize(const G4ThreeVector& p0,
0148 const G4ThreeVector& p1,
0149 const G4ThreeVector& p2,
0150 const G4ThreeVector& p3);
0151
0152
0153 G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector& p) const;
0154
0155 private:
0156
0157 G4double halfTolerance = 0;
0158 G4double fCubicVolume = 0;
0159 G4double fSurfaceArea = 0;
0160 mutable G4bool fRebuildPolyhedron = false;
0161 mutable G4Polyhedron* fpPolyhedron = nullptr;
0162
0163 G4ThreeVector fVertex[4];
0164 G4ThreeVector fNormal[4];
0165 G4double fDist[4] = {0};
0166 G4double fArea[4] = {0};
0167 G4ThreeVector fBmin, fBmax;
0168 };
0169
0170 #endif
0171
0172 #endif