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

#include <G4Polyhedra.hh>

+ Inheritance diagram for G4Polyhedra:

Public Member Functions

 G4Polyhedra (const G4String &name, G4double phiStart, G4double phiTotal, G4int numSide, G4int numZPlanes, const G4double zPlane[], const G4double rInner[], const G4double rOuter[])
 
 G4Polyhedra (const G4String &name, G4double phiStart, G4double phiTotal, G4int numSide, G4int numRZ, const G4double r[], const G4double z[])
 
virtual ~G4Polyhedra ()
 
EInside Inside (const G4ThreeVector &p) const
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
G4double DistanceToIn (const G4ThreeVector &p) const
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
G4GeometryType GetEntityType () const
 
G4VSolidClone () const
 
G4ThreeVector GetPointOnSurface () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4PolyhedronCreatePolyhedron () const
 
G4NURBSCreateNURBS () const
 
G4bool Reset ()
 
G4int GetNumSide () const
 
G4double GetStartPhi () const
 
G4double GetEndPhi () const
 
G4bool IsOpen () const
 
G4bool IsGeneric () const
 
G4int GetNumRZCorner () const
 
G4PolyhedraSideRZ GetCorner (const G4int index) const
 
G4PolyhedraHistoricalGetOriginalParameters () const
 
void SetOriginalParameters (G4PolyhedraHistorical *pars)
 
 G4Polyhedra (__void__ &)
 
 G4Polyhedra (const G4Polyhedra &source)
 
const G4Polyhedraoperator= (const G4Polyhedra &source)
 
- Public Member Functions inherited from G4VCSGfaceted
 G4VCSGfaceted (const G4String &name)
 
virtual ~G4VCSGfaceted ()
 
 G4VCSGfaceted (const G4VCSGfaceted &source)
 
const G4VCSGfacetedoperator= (const G4VCSGfaceted &source)
 
virtual G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
 
virtual EInside Inside (const G4ThreeVector &p) const
 
virtual G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
virtual G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
virtual G4double DistanceToIn (const G4ThreeVector &p) const
 
virtual G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
 
virtual G4double DistanceToOut (const G4ThreeVector &p) const
 
virtual G4GeometryType GetEntityType () const
 
virtual std::ostream & StreamInfo (std::ostream &os) const
 
virtual G4PolyhedronCreatePolyhedron () const =0
 
virtual void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
virtual G4VisExtent GetExtent () const
 
virtual G4PolyhedronGetPolyhedron () const
 
G4int GetCubVolStatistics () const
 
G4double GetCubVolEpsilon () const
 
void SetCubVolStatistics (G4int st)
 
void SetCubVolEpsilon (G4double ep)
 
G4int GetAreaStatistics () const
 
G4double GetAreaAccuracy () const
 
void SetAreaStatistics (G4int st)
 
void SetAreaAccuracy (G4double ep)
 
virtual G4double GetCubicVolume ()
 
virtual G4double GetSurfaceArea ()
 
 G4VCSGfaceted (__void__ &)
 
- Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
 
virtual ~G4VSolid ()
 
G4bool operator== (const G4VSolid &s) const
 
G4String GetName () const
 
void SetName (const G4String &name)
 
G4double GetTolerance () const
 
virtual G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
 
virtual EInside Inside (const G4ThreeVector &p) const =0
 
virtual G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p) const =0
 
virtual void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
virtual G4double GetCubicVolume ()
 
virtual G4double GetSurfaceArea ()
 
virtual G4GeometryType GetEntityType () const =0
 
virtual G4ThreeVector GetPointOnSurface () const
 
virtual G4VSolidClone () const
 
virtual std::ostream & StreamInfo (std::ostream &os) const =0
 
void DumpInfo () const
 
virtual void DescribeYourselfTo (G4VGraphicsScene &scene) const =0
 
virtual G4VisExtent GetExtent () const
 
virtual G4PolyhedronCreatePolyhedron () const
 
virtual G4NURBSCreateNURBS () const
 
virtual G4PolyhedronGetPolyhedron () const
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
 G4VSolid (__void__ &)
 
 G4VSolid (const G4VSolid &rhs)
 
G4VSolidoperator= (const G4VSolid &rhs)
 

Protected Member Functions

void SetOriginalParameters ()
 
void Create (G4double phiStart, G4double phiTotal, G4int numSide, G4ReduciblePolygon *rz)
 
void CopyStuff (const G4Polyhedra &source)
 
void DeleteStuff ()
 
G4ThreeVector GetPointOnPlane (G4ThreeVector p0, G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3) const
 
G4ThreeVector GetPointOnTriangle (G4ThreeVector p0, G4ThreeVector p1, G4ThreeVector p2) const
 
G4ThreeVector GetPointOnSurfaceCorners () const
 
- Protected Member Functions inherited from G4VCSGfaceted
virtual G4double DistanceTo (const G4ThreeVector &p, const G4bool outgoing) const
 
G4ThreeVector GetPointOnSurfaceGeneric () const
 
void CopyStuff (const G4VCSGfaceted &source)
 
void DeleteStuff ()
 
- Protected Member Functions inherited from G4VSolid
void CalculateClippedPolygonExtent (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipCrossSection (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipBetweenSections (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipPolygon (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 

Protected Attributes

G4int numSide
 
G4double startPhi
 
G4double endPhi
 
G4bool phiIsOpen
 
G4bool genericPgon
 
G4int numCorner
 
G4PolyhedraSideRZcorners
 
G4PolyhedraHistoricaloriginal_parameters
 
G4EnclosingCylinderenclosingCylinder
 
- Protected Attributes inherited from G4VCSGfaceted
G4int numFace
 
G4VCSGface ** faces
 
G4double fCubicVolume
 
G4double fSurfaceArea
 
G4PolyhedronfpPolyhedron
 
- Protected Attributes inherited from G4VSolid
G4double kCarTolerance
 

Detailed Description

Definition at line 88 of file G4Polyhedra.hh.

Constructor & Destructor Documentation

◆ G4Polyhedra() [1/4]

G4Polyhedra::G4Polyhedra ( const G4String name,
G4double  phiStart,
G4double  phiTotal,
G4int  numSide,
G4int  numZPlanes,
const G4double  zPlane[],
const G4double  rInner[],
const G4double  rOuter[] 
)

Definition at line 78 of file G4Polyhedra.cc.

86 : G4VCSGfaceted( name ), genericPgon(false)
87{
88 if (theNumSide <= 0)
89 {
90 std::ostringstream message;
91 message << "Solid must have at least one side - " << GetName() << G4endl
92 << " No sides specified !";
93 G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
94 FatalErrorInArgument, message);
95 }
96
97 //
98 // Calculate conversion factor from G3 radius to G4 radius
99 //
100 G4double phiTotal = thePhiTotal;
101 if ( (phiTotal <=0) || (phiTotal >= twopi*(1-DBL_EPSILON)) )
102 { phiTotal = twopi; }
103 G4double convertRad = std::cos(0.5*phiTotal/theNumSide);
104
105 //
106 // Some historical stuff
107 //
109
110 original_parameters->numSide = theNumSide;
113 original_parameters->Num_z_planes = numZPlanes;
114 original_parameters->Z_values = new G4double[numZPlanes];
115 original_parameters->Rmin = new G4double[numZPlanes];
116 original_parameters->Rmax = new G4double[numZPlanes];
117
118 G4int i;
119 for (i=0; i<numZPlanes; i++)
120 {
121 if (( i < numZPlanes-1) && ( zPlane[i] == zPlane[i+1] ))
122 {
123 if( (rInner[i] > rOuter[i+1])
124 ||(rInner[i+1] > rOuter[i]) )
125 {
126 DumpInfo();
127 std::ostringstream message;
128 message << "Cannot create a Polyhedra with no contiguous segments."
129 << G4endl
130 << " Segments are not contiguous !" << G4endl
131 << " rMin[" << i << "] = " << rInner[i]
132 << " -- rMax[" << i+1 << "] = " << rOuter[i+1] << G4endl
133 << " rMin[" << i+1 << "] = " << rInner[i+1]
134 << " -- rMax[" << i << "] = " << rOuter[i];
135 G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
136 FatalErrorInArgument, message);
137 }
138 }
139 original_parameters->Z_values[i] = zPlane[i];
140 original_parameters->Rmin[i] = rInner[i]/convertRad;
141 original_parameters->Rmax[i] = rOuter[i]/convertRad;
142 }
143
144
145 //
146 // Build RZ polygon using special PCON/PGON GEANT3 constructor
147 //
149 new G4ReduciblePolygon( rInner, rOuter, zPlane, numZPlanes );
150 rz->ScaleA( 1/convertRad );
151
152 //
153 // Do the real work
154 //
155 Create( phiStart, phiTotal, theNumSide, rz );
156
157 delete rz;
158}
@ FatalErrorInArgument
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4PolyhedraHistorical * original_parameters
Definition: G4Polyhedra.hh:195
G4bool genericPgon
Definition: G4Polyhedra.hh:192
void Create(G4double phiStart, G4double phiTotal, G4int numSide, G4ReduciblePolygon *rz)
Definition: G4Polyhedra.cc:191
void ScaleA(G4double scale)
G4String GetName() const
void DumpInfo() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define DBL_EPSILON
Definition: templates.hh:87

◆ G4Polyhedra() [2/4]

G4Polyhedra::G4Polyhedra ( const G4String name,
G4double  phiStart,
G4double  phiTotal,
G4int  numSide,
G4int  numRZ,
const G4double  r[],
const G4double  z[] 
)

Definition at line 164 of file G4Polyhedra.cc.

171 : G4VCSGfaceted( name ), genericPgon(true)
172{
173 G4ReduciblePolygon *rz = new G4ReduciblePolygon( r, z, numRZ );
174
175 Create( phiStart, phiTotal, theNumSide, rz );
176
177 // Set original_parameters struct for consistency
178 //
180
181 delete rz;
182}
void SetOriginalParameters()

◆ ~G4Polyhedra()

G4Polyhedra::~G4Polyhedra ( )
virtual

Definition at line 377 of file G4Polyhedra.cc.

378{
379 delete [] corners;
381
382 delete enclosingCylinder;
383}
G4PolyhedraSideRZ * corners
Definition: G4Polyhedra.hh:194
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polyhedra.hh:197

◆ G4Polyhedra() [3/4]

G4Polyhedra::G4Polyhedra ( __void__ &  a)

Definition at line 366 of file G4Polyhedra.cc.

367 : G4VCSGfaceted(a), numSide(0), startPhi(0.), endPhi(0.),
368 phiIsOpen(false), genericPgon(false), numCorner(0), corners(0),
370{
371}
G4double endPhi
Definition: G4Polyhedra.hh:190
G4bool phiIsOpen
Definition: G4Polyhedra.hh:191
G4double startPhi
Definition: G4Polyhedra.hh:189
G4int numCorner
Definition: G4Polyhedra.hh:193

◆ G4Polyhedra() [4/4]

G4Polyhedra::G4Polyhedra ( const G4Polyhedra source)

Definition at line 389 of file G4Polyhedra.cc.

390 : G4VCSGfaceted( source )
391{
392 CopyStuff( source );
393}
void CopyStuff(const G4Polyhedra &source)
Definition: G4Polyhedra.cc:419

Member Function Documentation

◆ Clone()

G4VSolid * G4Polyhedra::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 575 of file G4Polyhedra.cc.

576{
577 return new G4Polyhedra(*this);
578}

◆ ComputeDimensions()

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

Reimplemented from G4VSolid.

Definition at line 555 of file G4Polyhedra.cc.

558{
559 p->ComputeDimensions(*this,n,pRep);
560}
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const

◆ CopyStuff()

void G4Polyhedra::CopyStuff ( const G4Polyhedra source)
protected

Definition at line 419 of file G4Polyhedra.cc.

420{
421 //
422 // Simple stuff
423 //
424 numSide = source.numSide;
425 startPhi = source.startPhi;
426 endPhi = source.endPhi;
427 phiIsOpen = source.phiIsOpen;
428 numCorner = source.numCorner;
429 genericPgon= source.genericPgon;
430
431 //
432 // The corner array
433 //
435
437 *sourceCorn = source.corners;
438 do
439 {
440 *corn = *sourceCorn;
441 } while( ++sourceCorn, ++corn < corners+numCorner );
442
443 //
444 // Original parameters
445 //
446 if (source.original_parameters)
447 {
450 }
451
452 //
453 // Enclosing cylinder
454 //
456}

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

◆ Create()

void G4Polyhedra::Create ( G4double  phiStart,
G4double  phiTotal,
G4int  numSide,
G4ReduciblePolygon rz 
)
protected

Definition at line 191 of file G4Polyhedra.cc.

195{
196 //
197 // Perform checks of rz values
198 //
199 if (rz->Amin() < 0.0)
200 {
201 std::ostringstream message;
202 message << "Illegal input parameters - " << GetName() << G4endl
203 << " All R values must be >= 0 !";
204 G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
205 FatalErrorInArgument, message);
206 }
207
208 G4double rzArea = rz->Area();
209 if (rzArea < -kCarTolerance)
210 rz->ReverseOrder();
211
212 else if (rzArea < -kCarTolerance)
213 {
214 std::ostringstream message;
215 message << "Illegal input parameters - " << GetName() << G4endl
216 << " R/Z cross section is zero or near zero: " << rzArea;
217 G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
218 FatalErrorInArgument, message);
219 }
220
223 {
224 std::ostringstream message;
225 message << "Illegal input parameters - " << GetName() << G4endl
226 << " Too few unique R/Z values !";
227 G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
228 FatalErrorInArgument, message);
229 }
230
231 if (rz->CrossesItself( 1/kInfinity ))
232 {
233 std::ostringstream message;
234 message << "Illegal input parameters - " << GetName() << G4endl
235 << " R/Z segments cross !";
236 G4Exception("G4Polyhedra::Create()", "GeomSolids0002",
237 FatalErrorInArgument, message);
238 }
239
240 numCorner = rz->NumVertices();
241
242
243 startPhi = phiStart;
244 while( startPhi < 0 ) startPhi += twopi;
245 //
246 // Phi opening? Account for some possible roundoff, and interpret
247 // nonsense value as representing no phi opening
248 //
249 if ( (phiTotal <= 0) || (phiTotal > twopi*(1-DBL_EPSILON)) )
250 {
251 phiIsOpen = false;
252 endPhi = phiStart+twopi;
253 }
254 else
255 {
256 phiIsOpen = true;
257
258 //
259 // Convert phi into our convention
260 //
261 endPhi = phiStart+phiTotal;
262 while( endPhi < startPhi ) endPhi += twopi;
263 }
264
265 //
266 // Save number sides
267 //
268 numSide = theNumSide;
269
270 //
271 // Allocate corner array.
272 //
274
275 //
276 // Copy corners
277 //
279
281 iterRZ.Begin();
282 do
283 {
284 next->r = iterRZ.GetA();
285 next->z = iterRZ.GetB();
286 } while( ++next, iterRZ.Next() );
287
288 //
289 // Allocate face pointer array
290 //
292 faces = new G4VCSGface*[numFace];
293
294 //
295 // Construct side faces
296 //
297 // To do so properly, we need to keep track of four successive RZ
298 // corners.
299 //
300 // But! Don't construct a face if both points are at zero radius!
301 //
302 G4PolyhedraSideRZ *corner = corners,
303 *prev = corners + numCorner-1,
304 *nextNext;
305 G4VCSGface **face = faces;
306 do
307 {
308 next = corner+1;
309 if (next >= corners+numCorner) next = corners;
310 nextNext = next+1;
311 if (nextNext >= corners+numCorner) nextNext = corners;
312
313 if (corner->r < 1/kInfinity && next->r < 1/kInfinity) continue;
314/*
315 // We must decide here if we can dare declare one of our faces
316 // as having a "valid" normal (i.e. allBehind = true). This
317 // is never possible if the face faces "inward" in r *unless*
318 // we have only one side
319 //
320 G4bool allBehind;
321 if ((corner->z > next->z) && (numSide > 1))
322 {
323 allBehind = false;
324 }
325 else
326 {
327 //
328 // Otherwise, it is only true if the line passing
329 // through the two points of the segment do not
330 // split the r/z cross section
331 //
332 allBehind = !rz->BisectedBy( corner->r, corner->z,
333 next->r, next->z, kCarTolerance );
334 }
335*/
336 *face++ = new G4PolyhedraSide( prev, corner, next, nextNext,
338 } while( prev=corner, corner=next, corner > corners );
339
340 if (phiIsOpen)
341 {
342 //
343 // Construct phi open edges
344 //
345 *face++ = new G4PolyPhiFace( rz, startPhi, phiTotal/numSide, endPhi );
346 *face++ = new G4PolyPhiFace( rz, endPhi, phiTotal/numSide, startPhi );
347 }
348
349 //
350 // We might have dropped a face or two: recalculate numFace
351 //
352 numFace = face-faces;
353
354 //
355 // Make enclosingCylinder
356 //
358 new G4EnclosingCylinder( rz, phiIsOpen, phiStart, phiTotal );
359}
G4double Amin() const
G4bool RemoveDuplicateVertices(G4double tolerance)
G4int NumVertices() const
G4bool RemoveRedundantVertices(G4double tolerance)
G4bool CrossesItself(G4double tolerance)
G4VCSGface ** faces
G4double kCarTolerance
Definition: G4VSolid.hh:307

Referenced by G4Polyhedra(), and Reset().

◆ CreateNURBS()

G4NURBS * G4Polyhedra::CreateNURBS ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1161 of file G4Polyhedra.cc.

1162{
1163 return 0;
1164}

◆ CreatePolyhedron()

G4Polyhedron * G4Polyhedra::CreatePolyhedron ( ) const
virtual

Creates user defined polyhedron. This function allows to the user to define arbitrary polyhedron. The faces of the polyhedron should be either triangles or planar quadrilateral. Nodes of a face are defined by indexes pointing to the elements in the xyz array. Numeration of the elements in the array starts from 1 (like in fortran). The indexes can be positive or negative. Negative sign means that the corresponding edge is invisible. The normal of the face should be directed to exterior of the polyhedron.

Parameters
Nnodesnumber of nodes
Nfacesnumber of faces
xyznodes
faces_vecfaces (quadrilaterals or triangles)
Returns
status of the operation - is non-zero in case of problem

Implements G4VCSGfaceted.

Definition at line 903 of file G4Polyhedra.cc.

904{
905 if (!genericPgon)
906 {
914 }
915 else
916 {
917 // The following code prepares for:
918 // HepPolyhedron::createPolyhedron(int Nnodes, int Nfaces,
919 // const double xyz[][3],
920 // const int faces_vec[][4])
921 // Here is an extract from the header file HepPolyhedron.h:
922 /**
923 * Creates user defined polyhedron.
924 * This function allows to the user to define arbitrary polyhedron.
925 * The faces of the polyhedron should be either triangles or planar
926 * quadrilateral. Nodes of a face are defined by indexes pointing to
927 * the elements in the xyz array. Numeration of the elements in the
928 * array starts from 1 (like in fortran). The indexes can be positive
929 * or negative. Negative sign means that the corresponding edge is
930 * invisible. The normal of the face should be directed to exterior
931 * of the polyhedron.
932 *
933 * @param Nnodes number of nodes
934 * @param Nfaces number of faces
935 * @param xyz nodes
936 * @param faces_vec faces (quadrilaterals or triangles)
937 * @return status of the operation - is non-zero in case of problem
938 */
939 G4int nNodes;
940 G4int nFaces;
941 typedef G4double double3[3];
942 double3* xyz;
943 typedef G4int int4[4];
944 int4* faces_vec;
945 if (phiIsOpen)
946 {
947 // Triangulate open ends. Simple ear-chopping algorithm...
948 // I'm not sure how robust this algorithm is (J.Allison).
949 //
950 std::vector<G4bool> chopped(numCorner, false);
951 std::vector<G4int*> triQuads;
952 G4int remaining = numCorner;
953 G4int iStarter = 0;
954 while (remaining >= 3)
955 {
956 // Find unchopped corners...
957 //
958 G4int A = -1, B = -1, C = -1;
959 G4int iStepper = iStarter;
960 do
961 {
962 if (A < 0) { A = iStepper; }
963 else if (B < 0) { B = iStepper; }
964 else if (C < 0) { C = iStepper; }
965 do
966 {
967 if (++iStepper >= numCorner) iStepper = 0;
968 }
969 while (chopped[iStepper]);
970 }
971 while (C < 0 && iStepper != iStarter);
972
973 // Check triangle at B is pointing outward (an "ear").
974 // Sign of z cross product determines...
975
976 G4double BAr = corners[A].r - corners[B].r;
977 G4double BAz = corners[A].z - corners[B].z;
978 G4double BCr = corners[C].r - corners[B].r;
979 G4double BCz = corners[C].z - corners[B].z;
980 if (BAr * BCz - BAz * BCr < kCarTolerance)
981 {
982 G4int* tq = new G4int[3];
983 tq[0] = A + 1;
984 tq[1] = B + 1;
985 tq[2] = C + 1;
986 triQuads.push_back(tq);
987 chopped[B] = true;
988 --remaining;
989 }
990 else
991 {
992 do
993 {
994 if (++iStarter >= numCorner) { iStarter = 0; }
995 }
996 while (chopped[iStarter]);
997 }
998 }
999
1000 // Transfer to faces...
1001
1002 nNodes = (numSide + 1) * numCorner;
1003 nFaces = numSide * numCorner + 2 * triQuads.size();
1004 faces_vec = new int4[nFaces];
1005 G4int iface = 0;
1006 G4int addition = numCorner * numSide;
1007 G4int d = numCorner - 1;
1008 for (G4int iEnd = 0; iEnd < 2; ++iEnd)
1009 {
1010 for (size_t i = 0; i < triQuads.size(); ++i)
1011 {
1012 // Negative for soft/auxiliary/normally invisible edges...
1013 //
1014 G4int a, b, c;
1015 if (iEnd == 0)
1016 {
1017 a = triQuads[i][0];
1018 b = triQuads[i][1];
1019 c = triQuads[i][2];
1020 }
1021 else
1022 {
1023 a = triQuads[i][0] + addition;
1024 b = triQuads[i][2] + addition;
1025 c = triQuads[i][1] + addition;
1026 }
1027 G4int ab = std::abs(b - a);
1028 G4int bc = std::abs(c - b);
1029 G4int ca = std::abs(a - c);
1030 faces_vec[iface][0] = (ab == 1 || ab == d)? a: -a;
1031 faces_vec[iface][1] = (bc == 1 || bc == d)? b: -b;
1032 faces_vec[iface][2] = (ca == 1 || ca == d)? c: -c;
1033 faces_vec[iface][3] = 0;
1034 ++iface;
1035 }
1036 }
1037
1038 // Continue with sides...
1039
1040 xyz = new double3[nNodes];
1041 const G4double dPhi = (endPhi - startPhi) / numSide;
1042 G4double phi = startPhi;
1043 G4int ixyz = 0;
1044 for (G4int iSide = 0; iSide < numSide; ++iSide)
1045 {
1046 for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1047 {
1048 xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1049 xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1050 xyz[ixyz][2] = corners[iCorner].z;
1051 if (iCorner < numCorner - 1)
1052 {
1053 faces_vec[iface][0] = ixyz + 1;
1054 faces_vec[iface][1] = ixyz + numCorner + 1;
1055 faces_vec[iface][2] = ixyz + numCorner + 2;
1056 faces_vec[iface][3] = ixyz + 2;
1057 }
1058 else
1059 {
1060 faces_vec[iface][0] = ixyz + 1;
1061 faces_vec[iface][1] = ixyz + numCorner + 1;
1062 faces_vec[iface][2] = ixyz + 2;
1063 faces_vec[iface][3] = ixyz - numCorner + 2;
1064 }
1065 ++iface;
1066 ++ixyz;
1067 }
1068 phi += dPhi;
1069 }
1070
1071 // Last corners...
1072
1073 for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1074 {
1075 xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1076 xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1077 xyz[ixyz][2] = corners[iCorner].z;
1078 ++ixyz;
1079 }
1080 }
1081 else // !phiIsOpen - i.e., a complete 360 degrees.
1082 {
1083 nNodes = numSide * numCorner;
1084 nFaces = numSide * numCorner;;
1085 xyz = new double3[nNodes];
1086 faces_vec = new int4[nFaces];
1087 // const G4double dPhi = (endPhi - startPhi) / numSide;
1088 const G4double dPhi = twopi / numSide; // !phiIsOpen endPhi-startPhi = 360 degrees.
1089 G4double phi = startPhi;
1090 G4int ixyz = 0, iface = 0;
1091 for (G4int iSide = 0; iSide < numSide; ++iSide)
1092 {
1093 for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1094 {
1095 xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1096 xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1097 xyz[ixyz][2] = corners[iCorner].z;
1098 if (iSide < numSide - 1)
1099 {
1100 if (iCorner < numCorner - 1)
1101 {
1102 faces_vec[iface][0] = ixyz + 1;
1103 faces_vec[iface][1] = ixyz + numCorner + 1;
1104 faces_vec[iface][2] = ixyz + numCorner + 2;
1105 faces_vec[iface][3] = ixyz + 2;
1106 }
1107 else
1108 {
1109 faces_vec[iface][0] = ixyz + 1;
1110 faces_vec[iface][1] = ixyz + numCorner + 1;
1111 faces_vec[iface][2] = ixyz + 2;
1112 faces_vec[iface][3] = ixyz - numCorner + 2;
1113 }
1114 }
1115 else // Last side joins ends...
1116 {
1117 if (iCorner < numCorner - 1)
1118 {
1119 faces_vec[iface][0] = ixyz + 1;
1120 faces_vec[iface][1] = ixyz + numCorner - nFaces + 1;
1121 faces_vec[iface][2] = ixyz + numCorner - nFaces + 2;
1122 faces_vec[iface][3] = ixyz + 2;
1123 }
1124 else
1125 {
1126 faces_vec[iface][0] = ixyz + 1;
1127 faces_vec[iface][1] = ixyz - nFaces + numCorner + 1;
1128 faces_vec[iface][2] = ixyz - nFaces + 2;
1129 faces_vec[iface][3] = ixyz - numCorner + 2;
1130 }
1131 }
1132 ++ixyz;
1133 ++iface;
1134 }
1135 phi += dPhi;
1136 }
1137 }
1138 G4Polyhedron* polyhedron = new G4Polyhedron;
1139 G4int problem = polyhedron->createPolyhedron(nNodes, nFaces, xyz, faces_vec);
1140 delete [] faces_vec;
1141 delete [] xyz;
1142 if (problem)
1143 {
1144 std::ostringstream message;
1145 message << "Problem creating G4Polyhedron for: " << GetName();
1146 G4Exception("G4Polyhedra::CreatePolyhedron()", "GeomSolids1002",
1147 JustWarning, message);
1148 delete polyhedron;
1149 return 0;
1150 }
1151 else
1152 {
1153 return polyhedron;
1154 }
1155 }
1156}
@ JustWarning
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])

◆ DeleteStuff()

void G4Polyhedra::DeleteStuff ( )
protected

◆ DistanceToIn() [1/2]

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

Reimplemented from G4VCSGfaceted.

Definition at line 546 of file G4Polyhedra.cc.

547{
549}
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const

◆ DistanceToIn() [2/2]

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

Reimplemented from G4VCSGfaceted.

Definition at line 527 of file G4Polyhedra.cc.

529{
530 //
531 // Quick test
532 //
534 return kInfinity;
535
536 //
537 // Long answer
538 //
539 return G4VCSGfaceted::DistanceToIn( p, v );
540}
G4bool ShouldMiss(const G4ThreeVector &p, const G4ThreeVector &v) const

◆ GetCorner()

G4PolyhedraSideRZ G4Polyhedra::GetCorner ( const G4int  index) const
inline

◆ GetEndPhi()

◆ GetEntityType()

G4GeometryType G4Polyhedra::GetEntityType ( ) const
virtual

Reimplemented from G4VCSGfaceted.

Definition at line 566 of file G4Polyhedra.cc.

567{
568 return G4String("G4Polyhedra");
569}

◆ GetNumRZCorner()

G4int G4Polyhedra::GetNumRZCorner ( ) const
inline

◆ GetNumSide()

◆ GetOriginalParameters()

◆ GetPointOnPlane()

G4ThreeVector G4Polyhedra::GetPointOnPlane ( G4ThreeVector  p0,
G4ThreeVector  p1,
G4ThreeVector  p2,
G4ThreeVector  p3 
) const
protected

Definition at line 637 of file G4Polyhedra.cc.

639{
640 G4double lambda1, lambda2, chose,aOne,aTwo;
641 G4ThreeVector t, u, v, w, Area, normal;
642 aOne = 1.;
643 aTwo = 1.;
644
645 t = p1 - p0;
646 u = p2 - p1;
647 v = p3 - p2;
648 w = p0 - p3;
649
650 chose = RandFlat::shoot(0.,aOne+aTwo);
651 if( (chose>=0.) && (chose < aOne) )
652 {
653 lambda1 = RandFlat::shoot(0.,1.);
654 lambda2 = RandFlat::shoot(0.,lambda1);
655 return (p2+lambda1*v+lambda2*w);
656 }
657
658 lambda1 = RandFlat::shoot(0.,1.);
659 lambda2 = RandFlat::shoot(0.,lambda1);
660 return (p0+lambda1*t+lambda2*u);
661}
static double shoot()
Definition: RandFlat.cc:59

Referenced by GetPointOnSurface().

◆ GetPointOnSurface()

G4ThreeVector G4Polyhedra::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 686 of file G4Polyhedra.cc.

687{
688 if( !genericPgon ) // Polyhedra by faces
689 {
690 G4int j, numPlanes = original_parameters->Num_z_planes, Flag=0;
691 G4double chose, totArea=0., Achose1, Achose2,
692 rad1, rad2, sinphi1, sinphi2, cosphi1, cosphi2;
693 G4double a, b, l2, rang, totalPhi, ksi,
694 area, aTop=0., aBottom=0., zVal=0.;
695
696 G4ThreeVector p0, p1, p2, p3;
697 std::vector<G4double> aVector1;
698 std::vector<G4double> aVector2;
699 std::vector<G4double> aVector3;
700
701 totalPhi= (phiIsOpen) ? (endPhi-startPhi) : twopi;
702 ksi = totalPhi/numSide;
703 G4double cosksi = std::cos(ksi/2.);
704
705 // Below we generate the areas relevant to our solid
706 //
707 for(j=0; j<numPlanes-1; j++)
708 {
709 a = original_parameters->Rmax[j+1];
712 -original_parameters->Z_values[j+1]) + sqr(b-a);
713 area = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
714 aVector1.push_back(area);
715 }
716
717 for(j=0; j<numPlanes-1; j++)
718 {
719 a = original_parameters->Rmin[j+1];//*cosksi;
720 b = original_parameters->Rmin[j];//*cosksi;
722 -original_parameters->Z_values[j+1]) + sqr(b-a);
723 area = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
724 aVector2.push_back(area);
725 }
726
727 for(j=0; j<numPlanes-1; j++)
728 {
729 if(phiIsOpen == true)
730 {
731 aVector3.push_back(0.5*(original_parameters->Rmax[j]
735 *std::fabs(original_parameters->Z_values[j+1]
737 }
738 else { aVector3.push_back(0.); }
739 }
740
741 for(j=0; j<numPlanes-1; j++)
742 {
743 totArea += numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
744 }
745
746 // Must include top and bottom areas
747 //
748 if(original_parameters->Rmax[numPlanes-1] != 0.)
749 {
750 a = original_parameters->Rmax[numPlanes-1];
751 b = original_parameters->Rmin[numPlanes-1];
752 l2 = sqr(a-b);
753 aTop = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
754 }
755
756 if(original_parameters->Rmax[0] != 0.)
757 {
760 l2 = sqr(a-b);
761 aBottom = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi;
762 }
763
764 Achose1 = 0.;
765 Achose2 = numSide*(aVector1[0]+aVector2[0])+2.*aVector3[0];
766
767 chose = RandFlat::shoot(0.,totArea+aTop+aBottom);
768 if( (chose >= 0.) && (chose < aTop + aBottom) )
769 {
770 chose = RandFlat::shoot(startPhi,startPhi+totalPhi);
771 rang = std::floor((chose-startPhi)/ksi-0.01);
772 if(rang<0) { rang=0; }
773 rang = std::fabs(rang);
774 sinphi1 = std::sin(startPhi+rang*ksi);
775 sinphi2 = std::sin(startPhi+(rang+1)*ksi);
776 cosphi1 = std::cos(startPhi+rang*ksi);
777 cosphi2 = std::cos(startPhi+(rang+1)*ksi);
778 chose = RandFlat::shoot(0., aTop + aBottom);
779 if(chose>=0. && chose<aTop)
780 {
781 rad1 = original_parameters->Rmin[numPlanes-1];
782 rad2 = original_parameters->Rmax[numPlanes-1];
783 zVal = original_parameters->Z_values[numPlanes-1];
784 }
785 else
786 {
787 rad1 = original_parameters->Rmin[0];
788 rad2 = original_parameters->Rmax[0];
789 zVal = original_parameters->Z_values[0];
790 }
791 p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
792 p1 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
793 p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
794 p3 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
795 return GetPointOnPlane(p0,p1,p2,p3);
796 }
797 else
798 {
799 for (j=0; j<numPlanes-1; j++)
800 {
801 if( ((chose >= Achose1) && (chose < Achose2)) || (j == numPlanes-2) )
802 {
803 Flag = j; break;
804 }
805 Achose1 += numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
806 Achose2 = Achose1 + numSide*(aVector1[j+1]+aVector2[j+1])
807 + 2.*aVector3[j+1];
808 }
809 }
810
811 // At this point we have chosen a subsection
812 // between to adjacent plane cuts...
813
814 j = Flag;
815
816 totArea = numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j];
817 chose = RandFlat::shoot(0.,totArea);
818
819 if( (chose>=0.) && (chose<numSide*aVector1[j]) )
820 {
821 chose = RandFlat::shoot(startPhi,startPhi+totalPhi);
822 rang = std::floor((chose-startPhi)/ksi-0.01);
823 if(rang<0) { rang=0; }
824 rang = std::fabs(rang);
825 rad1 = original_parameters->Rmax[j];
826 rad2 = original_parameters->Rmax[j+1];
827 sinphi1 = std::sin(startPhi+rang*ksi);
828 sinphi2 = std::sin(startPhi+(rang+1)*ksi);
829 cosphi1 = std::cos(startPhi+rang*ksi);
830 cosphi2 = std::cos(startPhi+(rang+1)*ksi);
831 zVal = original_parameters->Z_values[j];
832
833 p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
834 p1 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
835
836 zVal = original_parameters->Z_values[j+1];
837
838 p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
839 p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
840 return GetPointOnPlane(p0,p1,p2,p3);
841 }
842 else if ( (chose >= numSide*aVector1[j])
843 && (chose <= numSide*(aVector1[j]+aVector2[j])) )
844 {
845 chose = RandFlat::shoot(startPhi,startPhi+totalPhi);
846 rang = std::floor((chose-startPhi)/ksi-0.01);
847 if(rang<0) { rang=0; }
848 rang = std::fabs(rang);
849 rad1 = original_parameters->Rmin[j];
850 rad2 = original_parameters->Rmin[j+1];
851 sinphi1 = std::sin(startPhi+rang*ksi);
852 sinphi2 = std::sin(startPhi+(rang+1)*ksi);
853 cosphi1 = std::cos(startPhi+rang*ksi);
854 cosphi2 = std::cos(startPhi+(rang+1)*ksi);
855 zVal = original_parameters->Z_values[j];
856
857 p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal);
858 p1 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal);
859
860 zVal = original_parameters->Z_values[j+1];
861
862 p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal);
863 p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal);
864 return GetPointOnPlane(p0,p1,p2,p3);
865 }
866
867 chose = RandFlat::shoot(0.,2.2);
868 if( (chose>=0.) && (chose < 1.) )
869 {
870 rang = startPhi;
871 }
872 else
873 {
874 rang = endPhi;
875 }
876
877 cosphi1 = std::cos(rang); rad1 = original_parameters->Rmin[j];
878 sinphi1 = std::sin(rang); rad2 = original_parameters->Rmax[j];
879
880 p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,
882 p1 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,
884
885 rad1 = original_parameters->Rmax[j+1];
886 rad2 = original_parameters->Rmin[j+1];
887
888 p2 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,
890 p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,
892 return GetPointOnPlane(p0,p1,p2,p3);
893 }
894 else // Generic polyhedra
895 {
896 return GetPointOnSurfaceGeneric();
897 }
898}
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector GetPointOnPlane(G4ThreeVector p0, G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3) const
Definition: G4Polyhedra.cc:637
G4ThreeVector GetPointOnSurfaceGeneric() const
T sqr(const T &x)
Definition: templates.hh:145

◆ GetPointOnSurfaceCorners()

G4ThreeVector G4Polyhedra::GetPointOnSurfaceCorners ( ) const
protected

◆ GetPointOnTriangle()

G4ThreeVector G4Polyhedra::GetPointOnTriangle ( G4ThreeVector  p0,
G4ThreeVector  p1,
G4ThreeVector  p2 
) const
protected

Definition at line 669 of file G4Polyhedra.cc.

672{
673 G4double lambda1,lambda2;
674 G4ThreeVector v=p3-p1, w=p1-p2;
675
676 lambda1 = RandFlat::shoot(0.,1.);
677 lambda2 = RandFlat::shoot(0.,lambda1);
678
679 return (p2 + lambda1*w + lambda2*v);
680}

◆ GetStartPhi()

◆ Inside()

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

Reimplemented from G4VCSGfaceted.

Definition at line 507 of file G4Polyhedra.cc.

508{
509 //
510 // Quick test
511 //
513
514 //
515 // Long answer
516 //
517 return G4VCSGfaceted::Inside(p);
518}
G4bool MustBeOutside(const G4ThreeVector &p) const
virtual EInside Inside(const G4ThreeVector &p) const
@ kOutside
Definition: geomdefs.hh:58

◆ IsGeneric()

G4bool G4Polyhedra::IsGeneric ( ) const
inline

◆ IsOpen()

G4bool G4Polyhedra::IsOpen ( ) const
inline

◆ operator=()

const G4Polyhedra & G4Polyhedra::operator= ( const G4Polyhedra source)

Definition at line 399 of file G4Polyhedra.cc.

400{
401 if (this == &source) return *this;
402
403 G4VCSGfaceted::operator=( source );
404
405 delete [] corners;
407
408 delete enclosingCylinder;
409
410 CopyStuff( source );
411
412 return *this;
413}
const G4VCSGfaceted & operator=(const G4VCSGfaceted &source)

◆ Reset()

G4bool G4Polyhedra::Reset ( )

Definition at line 465 of file G4Polyhedra.cc.

466{
467 if (genericPgon)
468 {
469 std::ostringstream message;
470 message << "Solid " << GetName() << " built using generic construct."
471 << G4endl << "Not applicable to the generic construct !";
472 G4Exception("G4Polyhedra::Reset()", "GeomSolids1001",
473 JustWarning, message, "Parameters NOT resetted.");
474 return 1;
475 }
476
477 //
478 // Clear old setup
479 //
481 delete [] corners;
482 delete enclosingCylinder;
483
484 //
485 // Rebuild polyhedra
486 //
495 delete rz;
496
497 return 0;
498}

Referenced by G4ParameterisationPolyhedraRho::ComputeDimensions(), G4ParameterisationPolyhedraPhi::ComputeDimensions(), and G4ParameterisationPolyhedraZ::ComputeDimensions().

◆ SetOriginalParameters() [1/2]

void G4Polyhedra::SetOriginalParameters ( )
inlineprotected

Referenced by G4Polyhedra().

◆ SetOriginalParameters() [2/2]

◆ StreamInfo()

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

Reimplemented from G4VCSGfaceted.

Definition at line 584 of file G4Polyhedra.cc.

585{
586 G4int oldprc = os.precision(16);
587 os << "-----------------------------------------------------------\n"
588 << " *** Dump for solid - " << GetName() << " ***\n"
589 << " ===================================================\n"
590 << " Solid type: G4Polyhedra\n"
591 << " Parameters: \n"
592 << " starting phi angle : " << startPhi/degree << " degrees \n"
593 << " ending phi angle : " << endPhi/degree << " degrees \n";
594 G4int i=0;
595 if (!genericPgon)
596 {
598 os << " number of Z planes: " << numPlanes << "\n"
599 << " Z values: \n";
600 for (i=0; i<numPlanes; i++)
601 {
602 os << " Z plane " << i << ": "
603 << original_parameters->Z_values[i] << "\n";
604 }
605 os << " Tangent distances to inner surface (Rmin): \n";
606 for (i=0; i<numPlanes; i++)
607 {
608 os << " Z plane " << i << ": "
609 << original_parameters->Rmin[i] << "\n";
610 }
611 os << " Tangent distances to outer surface (Rmax): \n";
612 for (i=0; i<numPlanes; i++)
613 {
614 os << " Z plane " << i << ": "
615 << original_parameters->Rmax[i] << "\n";
616 }
617 }
618 os << " number of RZ points: " << numCorner << "\n"
619 << " RZ values (corners): \n";
620 for (i=0; i<numCorner; i++)
621 {
622 os << " "
623 << corners[i].r << ", " << corners[i].z << "\n";
624 }
625 os << "-----------------------------------------------------------\n";
626 os.precision(oldprc);
627
628 return os;
629}

Member Data Documentation

◆ corners

G4PolyhedraSideRZ* G4Polyhedra::corners
protected

◆ enclosingCylinder

G4EnclosingCylinder* G4Polyhedra::enclosingCylinder
protected

Definition at line 197 of file G4Polyhedra.hh.

Referenced by CopyStuff(), Create(), DistanceToIn(), Inside(), operator=(), Reset(), and ~G4Polyhedra().

◆ endPhi

G4double G4Polyhedra::endPhi
protected

Definition at line 190 of file G4Polyhedra.hh.

Referenced by CopyStuff(), Create(), CreatePolyhedron(), GetPointOnSurface(), and StreamInfo().

◆ genericPgon

G4bool G4Polyhedra::genericPgon
protected

Definition at line 192 of file G4Polyhedra.hh.

Referenced by CopyStuff(), CreatePolyhedron(), GetPointOnSurface(), Reset(), and StreamInfo().

◆ numCorner

G4int G4Polyhedra::numCorner
protected

Definition at line 193 of file G4Polyhedra.hh.

Referenced by CopyStuff(), Create(), CreatePolyhedron(), and StreamInfo().

◆ numSide

G4int G4Polyhedra::numSide
protected

Definition at line 188 of file G4Polyhedra.hh.

Referenced by CopyStuff(), Create(), CreatePolyhedron(), and GetPointOnSurface().

◆ original_parameters

G4PolyhedraHistorical* G4Polyhedra::original_parameters
protected

◆ phiIsOpen

G4bool G4Polyhedra::phiIsOpen
protected

Definition at line 191 of file G4Polyhedra.hh.

Referenced by CopyStuff(), Create(), CreatePolyhedron(), and GetPointOnSurface().

◆ startPhi

G4double G4Polyhedra::startPhi
protected

Definition at line 189 of file G4Polyhedra.hh.

Referenced by CopyStuff(), Create(), CreatePolyhedron(), GetPointOnSurface(), and StreamInfo().


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