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

#include <G4GenericTrap.hh>

+ Inheritance diagram for G4GenericTrap:

Public Member Functions

 G4GenericTrap (const G4String &name, G4double halfZ, const std::vector< G4TwoVector > &vertices)
 
 ~G4GenericTrap () override
 
G4double GetZHalfLength () const
 
G4int GetNofVertices () const
 
G4TwoVector GetVertex (G4int index) const
 
const std::vector< G4TwoVector > & GetVertices () const
 
G4double GetTwistAngle (G4int index) const
 
G4bool IsTwisted () const
 
G4int GetVisSubdivisions () const
 
void SetVisSubdivisions (G4int subdiv)
 
EInside Inside (const G4ThreeVector &p) const override
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const override
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const override
 
G4double DistanceToIn (const G4ThreeVector &p) const override
 
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
 
G4double DistanceToOut (const G4ThreeVector &p) const override
 
void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const override
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const override
 
G4GeometryType GetEntityType () const override
 
G4VSolidClone () const override
 
std::ostream & StreamInfo (std::ostream &os) const override
 
G4ThreeVector GetPointOnSurface () const override
 
G4double GetCubicVolume () override
 
G4double GetSurfaceArea () override
 
G4PolyhedronGetPolyhedron () const override
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const override
 
G4VisExtent GetExtent () const override
 
G4PolyhedronCreatePolyhedron () const override
 
 G4GenericTrap (__void__ &)
 
 G4GenericTrap (const G4GenericTrap &rhs)
 
G4GenericTrapoperator= (const G4GenericTrap &rhs)
 
- Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
 
virtual ~G4VSolid ()
 
G4bool operator== (const G4VSolid &s) const
 
G4String GetName () const
 
void SetName (const G4String &name)
 
G4double GetTolerance () const
 
virtual void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
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
 

Protected Attributes

G4bool fRebuildPolyhedron = false
 
G4PolyhedronfpPolyhedron = nullptr
 
- Protected Attributes inherited from G4VSolid
G4double kCarTolerance
 

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
 

Detailed Description

Definition at line 79 of file G4GenericTrap.hh.

Constructor & Destructor Documentation

◆ G4GenericTrap() [1/3]

G4GenericTrap::G4GenericTrap ( const G4String & name,
G4double halfZ,
const std::vector< G4TwoVector > & vertices )

Definition at line 64 of file G4GenericTrap.cc.

66 : G4VSolid(name), fDz(halfZ), fVertices(),
67 fMinBBoxVector(G4ThreeVector(0,0,0)), fMaxBBoxVector(G4ThreeVector(0,0,0))
68{
69 // General constructor
70 const G4double min_length=5*1.e-6;
71 G4double length = 0.;
72 G4int k=0;
73 G4String errorDescription = "InvalidSetup in \" ";
74 errorDescription += name;
75 errorDescription += "\"";
76
77 halfCarTolerance = kCarTolerance*0.5;
78
79 // Check vertices size
80
81 if ( G4int(vertices.size()) != fgkNofVertices )
82 {
83 G4Exception("G4GenericTrap::G4GenericTrap()", "GeomSolids0002",
84 FatalErrorInArgument, "Number of vertices != 8");
85 }
86
87 // Check dZ
88 //
89 if (halfZ < kCarTolerance)
90 {
91 G4Exception("G4GenericTrap::G4GenericTrap()", "GeomSolids0002",
92 FatalErrorInArgument, "dZ is too small or negative");
93 }
94
95 // Check Ordering and Copy vertices
96 //
97 if(CheckOrder(vertices))
98 {
99 for (G4int i=0; i<fgkNofVertices; ++i) {fVertices.push_back(vertices[i]);}
100 }
101 else
102 {
103 for (auto i=0; i <4; ++i) {fVertices.push_back(vertices[3-i]);}
104 for (auto i=0; i <4; ++i) {fVertices.push_back(vertices[7-i]);}
105 }
106
107 // Check length of segments and Adjust
108 //
109 for (auto j=0; j<2; ++j)
110 {
111 for (auto i=1; i<4; ++i)
112 {
113 k = j*4+i;
114 length = (fVertices[k]-fVertices[k-1]).mag();
115 if ( ( length < min_length) && ( length > kCarTolerance ) )
116 {
117 std::ostringstream message;
118 message << "Length segment is too small." << G4endl
119 << "Distance between " << fVertices[k-1] << " and "
120 << fVertices[k] << " is only " << length << " mm !";
121 G4Exception("G4GenericTrap::G4GenericTrap()", "GeomSolids1001",
122 JustWarning, message, "Vertices will be collapsed.");
123 fVertices[k]=fVertices[k-1];
124 }
125 }
126 }
127
128 // Compute Twist
129 //
130 for(G4double & i : fTwist) { i=0.; }
131 fIsTwisted = ComputeIsTwisted();
132
133 // Compute Bounding Box
134 //
135 ComputeBBox();
136
137 // If not twisted - create tessellated solid
138 // (an alternative implementation for testing)
139 //
140#ifdef G4TESS_TEST
141 if ( !fIsTwisted ) { fTessellatedSolid = CreateTessellatedSolid(); }
142#endif
143}
@ JustWarning
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4VSolid(const G4String &name)
Definition G4VSolid.cc:57
G4double kCarTolerance
Definition G4VSolid.hh:299
const char * name(G4int ptype)

Referenced by Clone().

◆ ~G4GenericTrap()

G4GenericTrap::~G4GenericTrap ( )
override

Definition at line 157 of file G4GenericTrap.cc.

158{
159 delete fTessellatedSolid;
160}

◆ G4GenericTrap() [2/3]

G4GenericTrap::G4GenericTrap ( __void__ & a)

Definition at line 147 of file G4GenericTrap.cc.

148 : G4VSolid(a), halfCarTolerance(0.), fDz(0.), fVertices(),
149 fMinBBoxVector(G4ThreeVector(0,0,0)), fMaxBBoxVector(G4ThreeVector(0,0,0))
150{
151 // Fake default constructor - sets only member data and allocates memory
152 // for usage restricted to object persistency.
153}

◆ G4GenericTrap() [3/3]

G4GenericTrap::G4GenericTrap ( const G4GenericTrap & rhs)

Definition at line 164 of file G4GenericTrap.cc.

165 : G4VSolid(rhs),
166 halfCarTolerance(rhs.halfCarTolerance),
167 fDz(rhs.fDz), fVertices(rhs.fVertices),
168 fIsTwisted(rhs.fIsTwisted),
169 fMinBBoxVector(rhs.fMinBBoxVector), fMaxBBoxVector(rhs.fMaxBBoxVector),
170 fVisSubdivisions(rhs.fVisSubdivisions),
171 fSurfaceArea(rhs.fSurfaceArea), fCubicVolume(rhs.fCubicVolume)
172{
173 for (auto i=0; i<4; ++i) { fTwist[i] = rhs.fTwist[i]; }
174#ifdef G4TESS_TEST
175 if (rhs.fTessellatedSolid && !fIsTwisted )
176 { fTessellatedSolid = CreateTessellatedSolid(); }
177#endif
178}

Member Function Documentation

◆ BoundingLimits()

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

Reimplemented from G4VSolid.

Definition at line 1175 of file G4GenericTrap.cc.

1177{
1178 pMin = GetMinimumBBox();
1179 pMax = GetMaximumBBox();
1180
1181 // Check correctness of the bounding box
1182 //
1183 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
1184 {
1185 std::ostringstream message;
1186 message << "Bad bounding box (min >= max) for solid: "
1187 << GetName() << " !"
1188 << "\npMin = " << pMin
1189 << "\npMax = " << pMax;
1190 G4Exception("G4GenericTrap::BoundingLimits()", "GeomMgt0001",
1191 JustWarning, message);
1192 DumpInfo();
1193 }
1194}
double z() const
double x() const
double y() const
G4String GetName() const
void DumpInfo() const

Referenced by CalculateExtent().

◆ CalculateExtent()

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

Implements G4VSolid.

Definition at line 1199 of file G4GenericTrap.cc.

1203{
1204 G4ThreeVector bmin, bmax;
1205 G4bool exist;
1206
1207 // Check bounding box (bbox)
1208 //
1209 BoundingLimits(bmin,bmax);
1210 G4BoundingEnvelope bbox(bmin,bmax);
1211#ifdef G4BBOX_EXTENT
1212 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
1213#endif
1214 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
1215 {
1216 return exist = pMin < pMax;
1217 }
1218
1219 // Set bounding envelope (benv) and calculate extent
1220 //
1221 // To build the bounding envelope with plane faces each side face of
1222 // the trapezoid is subdivided in triangles. Subdivision is done by
1223 // duplication of vertices in the bases in a way that the envelope be
1224 // a convex polyhedron (some faces of the envelope can be degenerate)
1225 //
1226 G4double dz = GetZHalfLength();
1227 G4ThreeVectorList baseA(8), baseB(8);
1228 for (G4int i=0; i<4; ++i)
1229 {
1230 G4TwoVector va = GetVertex(i);
1231 G4TwoVector vb = GetVertex(i+4);
1232 baseA[2*i].set(va.x(),va.y(),-dz);
1233 baseB[2*i].set(vb.x(),vb.y(), dz);
1234 }
1235 for (G4int i=0; i<4; ++i)
1236 {
1237 G4int k1=2*i, k2=(2*i+2)%8;
1238 G4double ax = (baseA[k2].x()-baseA[k1].x());
1239 G4double ay = (baseA[k2].y()-baseA[k1].y());
1240 G4double bx = (baseB[k2].x()-baseB[k1].x());
1241 G4double by = (baseB[k2].y()-baseB[k1].y());
1242 G4double znorm = ax*by - ay*bx;
1243 baseA[k1+1] = (znorm < 0.0) ? baseA[k2] : baseA[k1];
1244 baseB[k1+1] = (znorm < 0.0) ? baseB[k1] : baseB[k2];
1245 }
1246
1247 std::vector<const G4ThreeVectorList *> polygons(2);
1248 polygons[0] = &baseA;
1249 polygons[1] = &baseB;
1250
1251 G4BoundingEnvelope benv(bmin,bmax,polygons);
1252 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
1253 return exist;
1254}
std::vector< G4ThreeVector > G4ThreeVectorList
bool G4bool
Definition G4Types.hh:86
double x() const
double y() const
G4double GetZHalfLength() const
G4TwoVector GetVertex(G4int index) const
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override

◆ Clone()

G4VSolid * G4GenericTrap::Clone ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 1265 of file G4GenericTrap.cc.

1266{
1267 return new G4GenericTrap(*this);
1268}
G4GenericTrap(const G4String &name, G4double halfZ, const std::vector< G4TwoVector > &vertices)

◆ CreatePolyhedron()

G4Polyhedron * G4GenericTrap::CreatePolyhedron ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 2012 of file G4GenericTrap.cc.

2013{
2014
2015#ifdef G4TESS_TEST
2016 if ( fTessellatedSolid != nullptr )
2017 {
2018 return fTessellatedSolid->CreatePolyhedron();
2019 }
2020#endif
2021
2022 // Approximation of Twisted Side
2023 // Construct extra Points, if Twisted Side
2024 //
2025 G4Polyhedron* polyhedron;
2026 G4int nVertices, nFacets;
2027
2028 G4int subdivisions = 0;
2029 if (fIsTwisted)
2030 {
2031 if (GetVisSubdivisions() != 0)
2032 {
2033 subdivisions = GetVisSubdivisions();
2034 }
2035 else
2036 {
2037 // Estimation of Number of Subdivisions for smooth visualisation
2038 //
2039 G4double maxTwist = 0.;
2040 for(G4int i = 0; i < 4; ++i)
2041 {
2042 if (GetTwistAngle(i) > maxTwist) { maxTwist = GetTwistAngle(i); }
2043 }
2044
2045 // Computes bounding vectors for the shape
2046 //
2047 G4double Dx, Dy;
2048 G4ThreeVector minVec = GetMinimumBBox();
2049 G4ThreeVector maxVec = GetMaximumBBox();
2050 Dx = 0.5*(maxVec.x() - minVec.y());
2051 Dy = 0.5*(maxVec.y() - minVec.y());
2052 if (Dy > Dx) { Dx = Dy; }
2053
2054 subdivisions = 8*G4int(maxTwist/(Dx*Dx*Dx)*fDz);
2055 if (subdivisions < 4) { subdivisions = 4; }
2056 if (subdivisions > 30) { subdivisions = 30; }
2057 }
2058 }
2059 G4int sub4 = 4*subdivisions;
2060 nVertices = 8 + subdivisions*4;
2061 nFacets = 6 + subdivisions*4;
2062 G4double cf = 1./(subdivisions + 1);
2063 polyhedron = new G4Polyhedron(nVertices, nFacets);
2064
2065 // Set vertices
2066 //
2067 G4int icur = 0;
2068 for (G4int i = 0; i < 4; ++i)
2069 {
2070 G4ThreeVector v(fVertices[i].x(),fVertices[i].y(),-fDz);
2071 polyhedron->SetVertex(++icur, v);
2072 }
2073 for (G4int i = 0; i < subdivisions; ++i)
2074 {
2075 for (G4int j = 0; j < 4; ++j)
2076 {
2077 G4TwoVector u = fVertices[j]+cf*(i+1)*(fVertices[j+4]-fVertices[j]);
2078 G4ThreeVector v(u.x(),u.y(),-fDz+cf*2*fDz*(i+1));
2079 polyhedron->SetVertex(++icur, v);
2080 }
2081 }
2082 for (G4int i = 4; i < 8; ++i)
2083 {
2084 G4ThreeVector v(fVertices[i].x(),fVertices[i].y(),fDz);
2085 polyhedron->SetVertex(++icur, v);
2086 }
2087
2088 // Set facets
2089 //
2090 icur = 0;
2091 polyhedron->SetFacet(++icur, 1, 4, 3, 2); // Z-plane
2092 for (G4int i = 0; i < subdivisions + 1; ++i)
2093 {
2094 G4int is = i*4;
2095 polyhedron->SetFacet(++icur, 5+is, 8+is, 4+is, 1+is);
2096 polyhedron->SetFacet(++icur, 8+is, 7+is, 3+is, 4+is);
2097 polyhedron->SetFacet(++icur, 7+is, 6+is, 2+is, 3+is);
2098 polyhedron->SetFacet(++icur, 6+is, 5+is, 1+is, 2+is);
2099 }
2100 polyhedron->SetFacet(++icur, 5+sub4, 6+sub4, 7+sub4, 8+sub4); // Z-plane
2101
2102 polyhedron->SetReferences();
2103 polyhedron->InvertFacets();
2104
2105 return polyhedron;
2106}
G4double GetTwistAngle(G4int index) const
G4int GetVisSubdivisions() const
G4Polyhedron * CreatePolyhedron() const override
void SetVertex(G4int index, const G4Point3D &v)
void SetFacet(G4int index, G4int iv1, G4int iv2, G4int iv3, G4int iv4=0)

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

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

Implements G4VSolid.

Definition at line 1977 of file G4GenericTrap.cc.

1978{
1979
1980#ifdef G4TESS_TEST
1981 if ( fTessellatedSolid )
1982 {
1983 return fTessellatedSolid->DescribeYourselfTo(scene);
1984 }
1985#endif
1986
1987 scene.AddSolid(*this);
1988}
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 799 of file G4GenericTrap.cc.

800{
801 // Computes the closest distance from given point to this shape
802
803#ifdef G4TESS_TEST
804 if ( fTessellatedSolid )
805 {
806 return fTessellatedSolid->DistanceToIn(p);
807 }
808#endif
809
810 G4double safz = std::fabs(p.z())-fDz;
811 if(safz<0) { safz=0; }
812
813 G4int iseg;
814 G4double safe = safz;
815 G4double safxy = safz;
816
817 for (iseg=0; iseg<4; ++iseg)
818 {
819 safxy = SafetyToFace(p,iseg);
820 if (safxy>safe) { safe=safxy; }
821 }
822
823 return safe;
824}
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override

◆ DistanceToIn() [2/2]

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

Implements G4VSolid.

Definition at line 731 of file G4GenericTrap.cc.

733{
734#ifdef G4TESS_TEST
735 if ( fTessellatedSolid )
736 {
737 return fTessellatedSolid->DistanceToIn(p, v);
738 }
739#endif
740
741 G4double dist[5];
743
744 // Check lateral faces
745 //
746 G4int i;
747 for (i=0; i<4; ++i)
748 {
749 dist[i]=DistToPlane(p, v, i);
750 }
751
752 // Check Z planes
753 //
754 dist[4]=kInfinity;
755 if (std::fabs(p.z())>fDz-halfCarTolerance)
756 {
757 if (v.z() != 0.0)
758 {
759 G4ThreeVector pt;
760 if (p.z()>0)
761 {
762 dist[4] = (fDz-p.z())/v.z();
763 }
764 else
765 {
766 dist[4] = (-fDz-p.z())/v.z();
767 }
768 if (dist[4]<-halfCarTolerance)
769 {
770 dist[4]=kInfinity;
771 }
772 else
773 {
774 if(dist[4]<halfCarTolerance)
775 {
776 if(p.z()>0) { n=G4ThreeVector(0,0,1); }
777 else { n=G4ThreeVector(0,0,-1); }
778 if (n.dot(v)<0) { dist[4]=0.; }
779 else { dist[4]=kInfinity; }
780 }
781 pt=p+dist[4]*v;
782 if (Inside(pt)==kOutside) { dist[4]=kInfinity; }
783 }
784 }
785 }
786 G4double distmin = dist[0];
787 for (i=1; i<5 ; ++i)
788 {
789 if (dist[i] < distmin) { distmin = dist[i]; }
790 }
791
792 if (distmin<halfCarTolerance) { distmin=0.; }
793
794 return distmin;
795}
EInside Inside(const G4ThreeVector &p) const override
@ kOutside
Definition geomdefs.hh:68

◆ DistanceToOut() [1/2]

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

Implements G4VSolid.

Definition at line 1148 of file G4GenericTrap.cc.

1149{
1150
1151#ifdef G4TESS_TEST
1152 if ( fTessellatedSolid )
1153 {
1154 return fTessellatedSolid->DistanceToOut(p);
1155 }
1156#endif
1157
1158 G4double safz = fDz-std::fabs(p.z());
1159 if (safz<0) { safz = 0; }
1160
1161 G4double safe = safz;
1162 G4double safxy = safz;
1163
1164 for (G4int iseg=0; iseg<4; ++iseg)
1165 {
1166 safxy = std::fabs(SafetyToFace(p,iseg));
1167 if (safxy < safe) { safe = safxy; }
1168 }
1169
1170 return safe;
1171}
G4double DistanceToOut(const G4ThreeVector &p) const override

◆ DistanceToOut() [2/2]

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

Implements G4VSolid.

Definition at line 903 of file G4GenericTrap.cc.

908{
909#ifdef G4TESS_TEST
910 if ( fTessellatedSolid )
911 {
912 return fTessellatedSolid->DistanceToOut(p, v, calcNorm, validNorm, n);
913 }
914#endif
915
916 G4double distmin;
917 G4bool lateral_cross = false;
918 ESide side = kUndef;
919
920 if (calcNorm) { *validNorm = true; } // All normals are valid
921
922 if (v.z() < 0)
923 {
924 distmin=(-fDz-p.z())/v.z();
925 if (calcNorm) { side=kMZ; *n=G4ThreeVector(0,0,-1); }
926 }
927 else
928 {
929 if (v.z() > 0)
930 {
931 distmin = (fDz-p.z())/v.z();
932 if (calcNorm) { side=kPZ; *n=G4ThreeVector(0,0,1); }
933 }
934 else { distmin = kInfinity; }
935 }
936
937 G4double dz2 =0.5/fDz;
938 G4double xa,xb,xc,xd;
939 G4double ya,yb,yc,yd;
940
941 for (G4int ipl=0; ipl<4; ++ipl)
942 {
943 G4int j = (ipl+1)%4;
944 xa=fVertices[ipl].x();
945 ya=fVertices[ipl].y();
946 xb=fVertices[ipl+4].x();
947 yb=fVertices[ipl+4].y();
948 xc=fVertices[j].x();
949 yc=fVertices[j].y();
950 xd=fVertices[4+j].x();
951 yd=fVertices[4+j].y();
952
953 if ( ((std::fabs(xb-xd)+std::fabs(yb-yd))<halfCarTolerance)
954 || ((std::fabs(xa-xc)+std::fabs(ya-yc))<halfCarTolerance) )
955 {
956 G4double q=DistToTriangle(p,v,ipl) ;
957 if ( (q>=0) && (q<distmin) )
958 {
959 distmin=q;
960 lateral_cross=true;
961 side=ESide(ipl+1);
962 }
963 continue;
964 }
965 G4double tx1 =dz2*(xb-xa);
966 G4double ty1 =dz2*(yb-ya);
967 G4double tx2 =dz2*(xd-xc);
968 G4double ty2 =dz2*(yd-yc);
969 G4double dzp =fDz+p.z();
970 G4double xs1 =xa+tx1*dzp;
971 G4double ys1 =ya+ty1*dzp;
972 G4double xs2 =xc+tx2*dzp;
973 G4double ys2 =yc+ty2*dzp;
974 G4double dxs =xs2-xs1;
975 G4double dys =ys2-ys1;
976 G4double dtx =tx2-tx1;
977 G4double dty =ty2-ty1;
978 G4double a = (dtx*v.y()-dty*v.x()+(tx1*ty2-tx2*ty1)*v.z())*v.z();
979 G4double b = dxs*v.y()-dys*v.x()+(dtx*p.y()-dty*p.x()+ty2*xs1-ty1*xs2
980 + tx1*ys2-tx2*ys1)*v.z();
981 G4double c=dxs*p.y()-dys*p.x()+xs1*ys2-xs2*ys1;
982 G4double q=kInfinity;
983
984 if (std::fabs(a) < kCarTolerance)
985 {
986 if (std::fabs(b) < kCarTolerance) { continue; }
987 q=-c/b;
988
989 // Check for Point on the Surface
990 //
991 if ((q > -halfCarTolerance) && (q < distmin))
992 {
993 if (q < halfCarTolerance)
994 {
995 if (NormalToPlane(p,ipl).dot(v)<0.) { continue; }
996 }
997 distmin =q;
998 lateral_cross=true;
999 side=ESide(ipl+1);
1000 }
1001 continue;
1002 }
1003 G4double d=b*b-4*a*c;
1004 if (d >= 0.)
1005 {
1006 if (a > 0) { q=0.5*(-b-std::sqrt(d))/a; }
1007 else { q=0.5*(-b+std::sqrt(d))/a; }
1008
1009 // Check for Point on the Surface
1010 //
1011 if (q > -halfCarTolerance )
1012 {
1013 if (q < distmin)
1014 {
1015 if(q < halfCarTolerance)
1016 {
1017 if (NormalToPlane(p,ipl).dot(v)<0.) // Check second root
1018 {
1019 if (a > 0) { q=0.5*(-b+std::sqrt(d))/a; }
1020 else { q=0.5*(-b-std::sqrt(d))/a; }
1021 if (( q > halfCarTolerance) && (q < distmin))
1022 {
1023 distmin=q;
1024 lateral_cross = true;
1025 side=ESide(ipl+1);
1026 }
1027 continue;
1028 }
1029 }
1030 distmin = q;
1031 lateral_cross = true;
1032 side=ESide(ipl+1);
1033 }
1034 }
1035 else
1036 {
1037 if (a > 0) { q=0.5*(-b+std::sqrt(d))/a; }
1038 else { q=0.5*(-b-std::sqrt(d))/a; }
1039
1040 // Check for Point on the Surface
1041 //
1042 if ((q > -halfCarTolerance) && (q < distmin))
1043 {
1044 if (q < halfCarTolerance)
1045 {
1046 if (NormalToPlane(p,ipl).dot(v)<0.) // Check second root
1047 {
1048 if (a > 0) { q=0.5*(-b-std::sqrt(d))/a; }
1049 else { q=0.5*(-b+std::sqrt(d))/a; }
1050 if ( ( q > halfCarTolerance) && (q < distmin) )
1051 {
1052 distmin=q;
1053 lateral_cross = true;
1054 side=ESide(ipl+1);
1055 }
1056 continue;
1057 }
1058 }
1059 distmin =q;
1060 lateral_cross = true;
1061 side=ESide(ipl+1);
1062 }
1063 }
1064 }
1065 }
1066 if (!lateral_cross) // Make sure that track crosses the top or bottom
1067 {
1068 if (distmin >= kInfinity) { distmin=kCarTolerance; }
1069 G4ThreeVector pt=p+distmin*v;
1070
1071 // Check if propagated point is in the polygon
1072 //
1073 G4int i=0;
1074 if (v.z()>0.) { i=4; }
1075 std::vector<G4TwoVector> xy;
1076 for ( G4int j=0; j<4; j++) { xy.push_back(fVertices[i+j]); }
1077
1078 // Check Inside
1079 //
1080 if (InsidePolygone(pt,xy)==kOutside)
1081 {
1082 if(calcNorm)
1083 {
1084 if (v.z()>0) {side= kPZ; *n = G4ThreeVector(0,0,1);}
1085 else { side=kMZ; *n = G4ThreeVector(0,0,-1);}
1086 }
1087 return 0.;
1088 }
1089 else
1090 {
1091 if(v.z()>0) {side=kPZ;}
1092 else {side=kMZ;}
1093 }
1094 }
1095
1096 if (calcNorm)
1097 {
1098 G4ThreeVector pt=p+v*distmin;
1099 switch (side)
1100 {
1101 case kXY0:
1102 *n=NormalToPlane(pt,0);
1103 break;
1104 case kXY1:
1105 *n=NormalToPlane(pt,1);
1106 break;
1107 case kXY2:
1108 *n=NormalToPlane(pt,2);
1109 break;
1110 case kXY3:
1111 *n=NormalToPlane(pt,3);
1112 break;
1113 case kPZ:
1114 *n=G4ThreeVector(0,0,1);
1115 break;
1116 case kMZ:
1117 *n=G4ThreeVector(0,0,-1);
1118 break;
1119 default:
1120 DumpInfo();
1121 std::ostringstream message;
1122 G4long oldprc = message.precision(16);
1123 message << "Undefined side for valid surface normal to solid." << G4endl
1124 << "Position:" << G4endl
1125 << " p.x() = " << p.x()/mm << " mm" << G4endl
1126 << " p.y() = " << p.y()/mm << " mm" << G4endl
1127 << " p.z() = " << p.z()/mm << " mm" << G4endl
1128 << "Direction:" << G4endl
1129 << " v.x() = " << v.x() << G4endl
1130 << " v.y() = " << v.y() << G4endl
1131 << " v.z() = " << v.z() << G4endl
1132 << "Proposed distance :" << G4endl
1133 << " distmin = " << distmin/mm << " mm";
1134 message.precision(oldprc);
1135 G4Exception("G4GenericTrap::DistanceToOut(p,v,..)",
1136 "GeomSolids1002", JustWarning, message);
1137 break;
1138 }
1139 }
1140
1141 if (distmin<halfCarTolerance) { distmin=0.; }
1142
1143 return distmin;
1144}
ESide
Definition G4Sphere.cc:61
long G4long
Definition G4Types.hh:87

◆ GetCubicVolume()

G4double G4GenericTrap::GetCubicVolume ( )
overridevirtual

Reimplemented from G4VSolid.

Definition at line 1382 of file G4GenericTrap.cc.

1383{
1384 if (fCubicVolume == 0.0)
1385 {
1386 // diagonals
1387 G4TwoVector A = fVertices[3] - fVertices[1];
1388 G4TwoVector B = fVertices[2] - fVertices[0];
1389 G4TwoVector C = fVertices[7] - fVertices[5];
1390 G4TwoVector D = fVertices[6] - fVertices[4];
1391
1392 // kross products
1393 G4double AB = A.x()*B.y() - A.y()*B.x();
1394 G4double CD = C.x()*D.y() - C.y()*D.x();
1395 G4double AD = A.x()*D.y() - A.y()*D.x();
1396 G4double CB = C.x()*B.y() - C.y()*B.x();
1397
1398 fCubicVolume = ((AB + CD)/3. + (AD + CB)/6.)*fDz;
1399 }
1400 return fCubicVolume;
1401}
G4double C(G4double temp)
G4double B(G4double temperature)
G4double D(G4double temp)
const G4double A[17]

◆ GetEntityType()

G4GeometryType G4GenericTrap::GetEntityType ( ) const
overridevirtual

Implements G4VSolid.

Definition at line 1258 of file G4GenericTrap.cc.

1259{
1260 return {"G4GenericTrap"};
1261}

Referenced by StreamInfo().

◆ GetExtent()

G4VisExtent G4GenericTrap::GetExtent ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 1992 of file G4GenericTrap.cc.

1993{
1994 // Computes bounding vectors for the shape
1995
1996#ifdef G4TESS_TEST
1997 if ( fTessellatedSolid != nullptr )
1998 {
1999 return fTessellatedSolid->GetExtent();
2000 }
2001#endif
2002
2003 G4ThreeVector minVec = GetMinimumBBox();
2004 G4ThreeVector maxVec = GetMaximumBBox();
2005 return { minVec.x(), maxVec.x(),
2006 minVec.y(), maxVec.y(),
2007 minVec.z(), maxVec.z() };
2008}
G4VisExtent GetExtent() const override

◆ GetNofVertices()

G4int G4GenericTrap::GetNofVertices ( ) const
inline

◆ GetPointOnSurface()

G4ThreeVector G4GenericTrap::GetPointOnSurface ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 1295 of file G4GenericTrap.cc.

1296{
1297
1298#ifdef G4TESS_TEST
1299 if ( fTessellatedSolid )
1300 {
1301 return fTessellatedSolid->GetPointOnSurface();
1302 }
1303#endif
1304
1305 G4ThreeVector point;
1306 G4TwoVector u,v,w;
1307 G4double rand,area,chose,cf,lambda0,lambda1,alfa,beta,zp;
1308 G4int ipl,j;
1309
1310 std::vector<G4ThreeVector> vertices;
1311 for (auto i=0; i<4; ++i)
1312 {
1313 vertices.emplace_back(fVertices[i].x(),fVertices[i].y(),-fDz);
1314 }
1315 for (auto i=4; i<8; ++i)
1316 {
1317 vertices.emplace_back(fVertices[i].x(),fVertices[i].y(),fDz);
1318 }
1319
1320 // Surface Area of Planes
1321 //
1322 G4TwoVector A = fVertices[3] - fVertices[1];
1323 G4TwoVector B = fVertices[2] - fVertices[0];
1324 G4TwoVector C = fVertices[7] - fVertices[5];
1325 G4TwoVector D = fVertices[6] - fVertices[4];
1326 G4double Surface0 = 0.5*(A.x()*B.y() - A.y()*B.x()); //-fDz plane
1327 G4double Surface1 = GetLateralFaceArea(0);
1328 G4double Surface2 = GetLateralFaceArea(1);
1329 G4double Surface3 = GetLateralFaceArea(2);
1330 G4double Surface4 = GetLateralFaceArea(3);
1331 G4double Surface5 = 0.5*(C.x()*D.y() - C.y()*D.x()); // fDz plane
1332
1333 rand = G4UniformRand();
1334 area = Surface0+Surface1+Surface2+Surface3+Surface4+Surface5;
1335 chose = rand*area;
1336
1337 if ( ( chose < Surface0)
1338 || ( chose > (Surface0+Surface1+Surface2+Surface3+Surface4)) )
1339 { // fDz or -fDz Plane
1340 ipl = G4int(G4UniformRand()*4);
1341 j = (ipl+1)%4;
1342 if(chose < Surface0)
1343 {
1344 zp = -fDz;
1345 u = fVertices[ipl]; v = fVertices[j];
1346 w = fVertices[(ipl+3)%4];
1347 }
1348 else
1349 {
1350 zp = fDz;
1351 u = fVertices[ipl+4]; v = fVertices[j+4];
1352 w = fVertices[(ipl+3)%4+4];
1353 }
1354 alfa = G4UniformRand();
1355 beta = G4UniformRand();
1356 lambda1=alfa*beta;
1357 lambda0=alfa-lambda1;
1358 v = v-u;
1359 w = w-u;
1360 v = u+lambda0*v+lambda1*w;
1361 }
1362 else // Lateral Plane Twisted or Not
1363 {
1364 if (chose < Surface0+Surface1) { ipl=0; }
1365 else if (chose < Surface0+Surface1+Surface2) { ipl=1; }
1366 else if (chose < Surface0+Surface1+Surface2+Surface3) { ipl=2; }
1367 else { ipl=3; }
1368 j = (ipl+1)%4;
1369 zp = -fDz+G4UniformRand()*2*fDz;
1370 cf = 0.5*(fDz-zp)/fDz;
1371 u = fVertices[ipl+4]+cf*( fVertices[ipl]-fVertices[ipl+4]);
1372 v = fVertices[j+4]+cf*(fVertices[j]-fVertices[j+4]);
1373 v = u+(v-u)*G4UniformRand();
1374 }
1375 point=G4ThreeVector(v.x(),v.y(),zp);
1376
1377 return point;
1378}
#define G4UniformRand()
Definition Randomize.hh:52
G4ThreeVector GetPointOnSurface() const override

◆ GetPolyhedron()

G4Polyhedron * G4GenericTrap::GetPolyhedron ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 1951 of file G4GenericTrap.cc.

1952{
1953
1954#ifdef G4TESS_TEST
1955 if ( fTessellatedSolid )
1956 {
1957 return fTessellatedSolid->GetPolyhedron();
1958 }
1959#endif
1960
1961 if ( (fpPolyhedron == nullptr)
1965 {
1966 G4AutoLock l(&polyhedronMutex);
1967 delete fpPolyhedron;
1969 fRebuildPolyhedron = false;
1970 l.unlock();
1971 }
1972 return fpPolyhedron;
1973}
G4bool fRebuildPolyhedron
G4Polyhedron * CreatePolyhedron() const override
G4Polyhedron * fpPolyhedron
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
G4Polyhedron * GetPolyhedron() const override
static G4int GetNumberOfRotationSteps()

◆ GetSurfaceArea()

G4double G4GenericTrap::GetSurfaceArea ( )
overridevirtual

Reimplemented from G4VSolid.

Definition at line 1468 of file G4GenericTrap.cc.

1469{
1470 if (fSurfaceArea == 0.0)
1471 {
1472 G4TwoVector A = fVertices[3] - fVertices[1];
1473 G4TwoVector B = fVertices[2] - fVertices[0];
1474 G4TwoVector C = fVertices[7] - fVertices[5];
1475 G4TwoVector D = fVertices[6] - fVertices[4];
1476 G4double S_bot = 0.5*(A.x()*B.y() - A.y()*B.x());
1477 G4double S_top = 0.5*(C.x()*D.y() - C.y()*D.x());
1478 fSurfaceArea = S_bot + S_top +
1479 GetLateralFaceArea(0) +
1480 GetLateralFaceArea(1) +
1481 GetLateralFaceArea(2) +
1482 GetLateralFaceArea(3);
1483 }
1484 return fSurfaceArea;
1485}

◆ GetTwistAngle()

G4double G4GenericTrap::GetTwistAngle ( G4int index) const
inline

Referenced by CreatePolyhedron(), and SurfaceNormal().

◆ GetVertex()

G4TwoVector G4GenericTrap::GetVertex ( G4int index) const
inline

Referenced by CalculateExtent().

◆ GetVertices()

const std::vector< G4TwoVector > & G4GenericTrap::GetVertices ( ) const
inline

◆ GetVisSubdivisions()

G4int G4GenericTrap::GetVisSubdivisions ( ) const
inline

Referenced by CreatePolyhedron().

◆ GetZHalfLength()

G4double G4GenericTrap::GetZHalfLength ( ) const
inline

◆ Inside()

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

Implements G4VSolid.

Definition at line 317 of file G4GenericTrap.cc.

318{
319 // Test if point is inside this shape
320
321#ifdef G4TESS_TEST
322 if ( fTessellatedSolid )
323 {
324 return fTessellatedSolid->Inside(p);
325 }
326#endif
327
328 EInside innew=kOutside;
329 std::vector<G4TwoVector> xy;
330
331 if (std::fabs(p.z()) <= fDz+halfCarTolerance) // First check Z range
332 {
333 // Compute intersection between Z plane containing point and the shape
334 //
335 G4double cf = 0.5*(fDz-p.z())/fDz;
336 for (auto i=0; i<4; ++i)
337 {
338 xy.push_back(fVertices[i+4]+cf*( fVertices[i]-fVertices[i+4]));
339 }
340
341 innew=InsidePolygone(p,xy);
342
343 if( (innew==kInside) || (innew==kSurface) )
344 {
345 if(std::fabs(p.z()) > fDz-halfCarTolerance) { innew=kSurface; }
346 }
347 }
348 return innew;
349}
EInside Inside(const G4ThreeVector &p) const override
EInside
Definition geomdefs.hh:67
@ kInside
Definition geomdefs.hh:70
@ kSurface
Definition geomdefs.hh:69

Referenced by DistanceToIn().

◆ IsTwisted()

G4bool G4GenericTrap::IsTwisted ( ) const
inline

◆ operator=()

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

Definition at line 182 of file G4GenericTrap.cc.

183{
184 // Check assignment to self
185 //
186 if (this == &rhs) { return *this; }
187
188 // Copy base class data
189 //
191
192 // Copy data
193 //
194 halfCarTolerance = rhs.halfCarTolerance;
195 fDz = rhs.fDz; fVertices = rhs.fVertices;
196 fIsTwisted = rhs.fIsTwisted; fTessellatedSolid = nullptr;
197 fMinBBoxVector = rhs.fMinBBoxVector; fMaxBBoxVector = rhs.fMaxBBoxVector;
198 fVisSubdivisions = rhs.fVisSubdivisions;
199 fSurfaceArea = rhs.fSurfaceArea; fCubicVolume = rhs.fCubicVolume;
200
201 for (auto i=0; i<4; ++i) { fTwist[i] = rhs.fTwist[i]; }
202#ifdef G4TESS_TEST
203 if (rhs.fTessellatedSolid && !fIsTwisted )
204 { delete fTessellatedSolid; fTessellatedSolid = CreateTessellatedSolid(); }
205#endif
206 fRebuildPolyhedron = false;
207 delete fpPolyhedron; fpPolyhedron = nullptr;
208
209 return *this;
210}
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:107

◆ SetVisSubdivisions()

void G4GenericTrap::SetVisSubdivisions ( G4int subdiv)
inline

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 1272 of file G4GenericTrap.cc.

1273{
1274 G4long oldprc = os.precision(16);
1275 os << "-----------------------------------------------------------\n"
1276 << " *** Dump for solid - " << GetName() << " *** \n"
1277 << " =================================================== \n"
1278 << " Solid geometry type: " << GetEntityType() << G4endl
1279 << " half length Z: " << fDz/mm << " mm \n"
1280 << " list of vertices:\n";
1281
1282 for ( G4int i=0; i<fgkNofVertices; ++i )
1283 {
1284 os << std::setw(5) << "#" << i
1285 << " vx = " << fVertices[i].x()/mm << " mm"
1286 << " vy = " << fVertices[i].y()/mm << " mm" << G4endl;
1287 }
1288 os.precision(oldprc);
1289
1290 return os;
1291}
G4GeometryType GetEntityType() const override

◆ SurfaceNormal()

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

Implements G4VSolid.

Definition at line 353 of file G4GenericTrap.cc.

354{
355 // Calculate side nearest to p, and return normal
356 // If two sides are equidistant, sum of the Normal is returned
357
358#ifdef G4TESS_TEST
359 if ( fTessellatedSolid )
360 {
361 return fTessellatedSolid->SurfaceNormal(p);
362 }
363#endif
364
365 G4ThreeVector lnorm, sumnorm(0.,0.,0.), apprnorm(0.,0.,1.),
366 p0, p1, p2, r1, r2, r3, r4;
367 G4int noSurfaces = 0;
368 G4double distxy,distz;
369 G4bool zPlusSide=false;
370
371 distz = fDz-std::fabs(p.z());
372 if (distz < halfCarTolerance)
373 {
374 if(p.z()>0)
375 {
376 zPlusSide=true;
377 sumnorm=G4ThreeVector(0,0,1);
378 }
379 else
380 {
381 sumnorm=G4ThreeVector(0,0,-1);
382 }
383 ++noSurfaces;
384 }
385
386 // Check lateral planes
387 //
388 std:: vector<G4TwoVector> vertices;
389 G4double cf = 0.5*(fDz-p.z())/fDz;
390 for (auto i=0; i<4; ++i)
391 {
392 vertices.push_back(fVertices[i+4]+cf*(fVertices[i]-fVertices[i+4]));
393 }
394
395 // Compute distance for lateral planes
396 //
397 for (G4int q=0; q<4; ++q)
398 {
399 p0=G4ThreeVector(vertices[q].x(),vertices[q].y(),p.z());
400 if(zPlusSide)
401 {
402 p1=G4ThreeVector(fVertices[q].x(),fVertices[q].y(),-fDz);
403 }
404 else
405 {
406 p1=G4ThreeVector(fVertices[q+4].x(),fVertices[q+4].y(),fDz);
407 }
408 p2=G4ThreeVector(vertices[(q+1)%4].x(),vertices[(q+1)%4].y(),p.z());
409
410 // Collapsed vertices
411 //
412 if ( (p2-p0).mag2() < kCarTolerance )
413 {
414 if ( std::fabs(p.z()+fDz) > kCarTolerance )
415 {
416 p2=G4ThreeVector(fVertices[(q+1)%4].x(),fVertices[(q+1)%4].y(),-fDz);
417 }
418 else
419 {
420 p2=G4ThreeVector(fVertices[(q+1)%4+4].x(),fVertices[(q+1)%4+4].y(),fDz);
421 }
422 }
423 lnorm = (p1-p0).cross(p2-p0);
424 lnorm = lnorm.unit();
425 if(zPlusSide) { lnorm=-lnorm; }
426
427 // Adjust Normal for Twisted Surface
428 //
429 if ( (fIsTwisted) && (GetTwistAngle(q)!=0) )
430 {
431 G4double normP=(p2-p0).mag();
432 if(normP != 0.0)
433 {
434 G4double proj=(p-p0).dot(p2-p0)/normP;
435 if(proj<0) { proj=0; }
436 if(proj>normP) { proj=normP; }
437 G4int j=(q+1)%4;
438 r1=G4ThreeVector(fVertices[q+4].x(),fVertices[q+4].y(),fDz);
439 r2=G4ThreeVector(fVertices[j+4].x(),fVertices[j+4].y(),fDz);
440 r3=G4ThreeVector(fVertices[q].x(),fVertices[q].y(),-fDz);
441 r4=G4ThreeVector(fVertices[j].x(),fVertices[j].y(),-fDz);
442 r1=r1+proj*(r2-r1)/normP;
443 r3=r3+proj*(r4-r3)/normP;
444 r2=r1-r3;
445 r4=r2.cross(p2-p0); r4=r4.unit();
446 lnorm=r4;
447 }
448 } // End if fIsTwisted
449
450 distxy=std::fabs((p0-p).dot(lnorm));
451 if ( distxy<halfCarTolerance )
452 {
453 ++noSurfaces;
454
455 // Negative sign for Normal is taken for Outside Normal
456 //
457 sumnorm=sumnorm+lnorm;
458 }
459
460 // For ApproxSurfaceNormal
461 //
462 if (distxy<distz)
463 {
464 distz=distxy;
465 apprnorm=lnorm;
466 }
467 } // End for loop
468
469 // Calculate final Normal, add Normal in the Corners and Touching Sides
470 //
471 if ( noSurfaces == 0 )
472 {
473#ifdef G4SPECSDEBUG
474 G4Exception("G4GenericTrap::SurfaceNormal(p)", "GeomSolids1002",
475 JustWarning, "Point p is not on surface !?" );
476#endif
477 sumnorm=apprnorm;
478 // Add Approximative Surface Normal Calculation?
479 }
480 else if ( noSurfaces == 1 ) { ; }
481 else { sumnorm = sumnorm.unit(); }
482
483 return sumnorm ;
484}
Hep3Vector unit() const
Hep3Vector cross(const Hep3Vector &) const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override

Member Data Documentation

◆ fpPolyhedron

G4Polyhedron* G4GenericTrap::fpPolyhedron = nullptr
mutableprotected

Definition at line 187 of file G4GenericTrap.hh.

Referenced by GetPolyhedron(), and operator=().

◆ fRebuildPolyhedron

G4bool G4GenericTrap::fRebuildPolyhedron = false
mutableprotected

Definition at line 186 of file G4GenericTrap.hh.

Referenced by GetPolyhedron(), and operator=().


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