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

#include <G4GenericPolycone.hh>

+ Inheritance diagram for G4GenericPolycone:

Classes

struct  surface_element
 

Public Member Functions

 G4GenericPolycone (const G4String &name, G4double phiStart, G4double phiTotal, G4int numRZ, const G4double r[], const G4double z[])
 
virtual ~G4GenericPolycone ()
 
EInside Inside (const G4ThreeVector &p) const
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
G4double DistanceToIn (const G4ThreeVector &p) const
 
void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
 
G4double GetCubicVolume ()
 
G4double GetSurfaceArea ()
 
G4ThreeVector GetPointOnSurface () const
 
G4GeometryType GetEntityType () const
 
G4VSolidClone () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4PolyhedronCreatePolyhedron () const
 
G4bool Reset ()
 
G4double GetStartPhi () const
 
G4double GetEndPhi () const
 
G4double GetSinStartPhi () const
 
G4double GetCosStartPhi () const
 
G4double GetSinEndPhi () const
 
G4double GetCosEndPhi () const
 
G4bool IsOpen () const
 
G4int GetNumRZCorner () const
 
G4PolyconeSideRZ GetCorner (G4int index) const
 
 G4GenericPolycone (__void__ &)
 
 G4GenericPolycone (const G4GenericPolycone &source)
 
G4GenericPolyconeoperator= (const G4GenericPolycone &source)
 
- Public Member Functions inherited from G4VCSGfaceted
 G4VCSGfaceted (const G4String &name)
 
virtual ~G4VCSGfaceted ()
 
 G4VCSGfaceted (const G4VCSGfaceted &source)
 
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=nullptr, G4ThreeVector *n=nullptr) 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 void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) 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=nullptr, G4ThreeVector *n=nullptr) 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 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)
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 

Protected Member Functions

void Create (G4double phiStart, G4double phiTotal, G4ReduciblePolygon *rz)
 
void CopyStuff (const G4GenericPolycone &source)
 
void SetSurfaceElements () 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
 

Protected Attributes

G4double startPhi
 
G4double endPhi
 
G4bool phiIsOpen = false
 
G4int numCorner
 
G4PolyconeSideRZcorners = nullptr
 
G4EnclosingCylinderenclosingCylinder = nullptr
 
std::vector< surface_element > * fElements = nullptr
 
- Protected Attributes inherited from G4VCSGfaceted
G4int numFace = 0
 
G4VCSGface ** faces = nullptr
 
G4double fCubicVolume = 0.0
 
G4double fSurfaceArea = 0.0
 
G4bool fRebuildPolyhedron = false
 
G4PolyhedronfpPolyhedron = nullptr
 
- Protected Attributes inherited from G4VSolid
G4double kCarTolerance
 

Detailed Description

Definition at line 65 of file G4GenericPolycone.hh.

Constructor & Destructor Documentation

◆ G4GenericPolycone() [1/3]

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

Definition at line 59 of file G4GenericPolycone.cc.

65 : G4VCSGfaceted( name )
66{
67
68 G4ReduciblePolygon *rz = new G4ReduciblePolygon( r, z, numRZ );
69
70 Create( phiStart, phiTotal, rz );
71
72 // Set original_parameters struct for consistency
73 //
74 //SetOriginalParameters(rz);
75
76 delete rz;
77}
void Create(G4double phiStart, G4double phiTotal, G4ReduciblePolygon *rz)

◆ ~G4GenericPolycone()

G4GenericPolycone::~G4GenericPolycone ( )
virtual

Definition at line 259 of file G4GenericPolycone.cc.

260{
261 delete [] corners;
262 delete enclosingCylinder;
263 delete fElements;
264 delete fpPolyhedron;
265 corners = nullptr;
266 enclosingCylinder = nullptr;
267 fElements = nullptr;
268 fpPolyhedron = nullptr;
269}
G4PolyconeSideRZ * corners
G4EnclosingCylinder * enclosingCylinder
std::vector< surface_element > * fElements
G4Polyhedron * fpPolyhedron

◆ G4GenericPolycone() [2/3]

G4GenericPolycone::G4GenericPolycone ( __void__ &  a)

Definition at line 252 of file G4GenericPolycone.cc.

◆ G4GenericPolycone() [3/3]

G4GenericPolycone::G4GenericPolycone ( const G4GenericPolycone source)

Definition at line 273 of file G4GenericPolycone.cc.

274 : G4VCSGfaceted( source )
275{
276 CopyStuff( source );
277}
void CopyStuff(const G4GenericPolycone &source)

Member Function Documentation

◆ BoundingLimits()

void G4GenericPolycone::BoundingLimits ( G4ThreeVector pMin,
G4ThreeVector pMax 
) const
virtual

Reimplemented from G4VSolid.

Definition at line 402 of file G4GenericPolycone.cc.

404{
405 G4double rmin = kInfinity, rmax = -kInfinity;
406 G4double zmin = kInfinity, zmax = -kInfinity;
407
408 for (G4int i=0; i<GetNumRZCorner(); ++i)
409 {
410 G4PolyconeSideRZ corner = GetCorner(i);
411 if (corner.r < rmin) rmin = corner.r;
412 if (corner.r > rmax) rmax = corner.r;
413 if (corner.z < zmin) zmin = corner.z;
414 if (corner.z > zmax) zmax = corner.z;
415 }
416
417 if (IsOpen())
418 {
419 G4TwoVector vmin,vmax;
420 G4GeomTools::DiskExtent(rmin,rmax,
423 vmin,vmax);
424 pMin.set(vmin.x(),vmin.y(),zmin);
425 pMax.set(vmax.x(),vmax.y(),zmax);
426 }
427 else
428 {
429 pMin.set(-rmax,-rmax, zmin);
430 pMax.set( rmax, rmax, zmax);
431 }
432
433 // Check correctness of the bounding box
434 //
435 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
436 {
437 std::ostringstream message;
438 message << "Bad bounding box (min >= max) for solid: "
439 << GetName() << " !"
440 << "\npMin = " << pMin
441 << "\npMax = " << pMax;
442 G4Exception("GenericG4Polycone::BoundingLimits()", "GeomMgt0001",
443 JustWarning, message);
444 DumpInfo();
445 }
446}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
double x() const
double y() const
double z() const
double x() const
double y() const
void set(double x, double y, double z)
G4bool IsOpen() const
G4int GetNumRZCorner() const
G4double GetCosEndPhi() const
G4double GetSinStartPhi() const
G4PolyconeSideRZ GetCorner(G4int index) const
G4double GetSinEndPhi() const
G4double GetCosStartPhi() const
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
Definition: G4GeomTools.cc:390
G4String GetName() const
void DumpInfo() const

Referenced by CalculateExtent().

◆ CalculateExtent()

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

Reimplemented from G4VCSGfaceted.

Definition at line 453 of file G4GenericPolycone.cc.

457{
458 G4ThreeVector bmin, bmax;
459 G4bool exist;
460
461 // Check bounding box (bbox)
462 //
463 BoundingLimits(bmin,bmax);
464 G4BoundingEnvelope bbox(bmin,bmax);
465#ifdef G4BBOX_EXTENT
466 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
467#endif
468 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
469 {
470 return exist = (pMin < pMax) ? true : false;
471 }
472
473 // To find the extent, RZ contour of the polycone is subdivided
474 // in triangles. The extent is calculated as cumulative extent of
475 // all sub-polycones formed by rotation of triangles around Z
476 //
477 G4TwoVectorList contourRZ;
478 G4TwoVectorList triangles;
479 G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
480 G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
481
482 // get RZ contour, ensure anticlockwise order of corners
483 for (G4int i=0; i<GetNumRZCorner(); ++i)
484 {
485 G4PolyconeSideRZ corner = GetCorner(i);
486 contourRZ.push_back(G4TwoVector(corner.r,corner.z));
487 }
488 G4double area = G4GeomTools::PolygonArea(contourRZ);
489 if (area < 0.) std::reverse(contourRZ.begin(),contourRZ.end());
490
491 // triangulate RZ countour
492 if (!G4GeomTools::TriangulatePolygon(contourRZ,triangles))
493 {
494 std::ostringstream message;
495 message << "Triangulation of RZ contour has failed for solid: "
496 << GetName() << " !"
497 << "\nExtent has been calculated using boundary box";
498 G4Exception("G4GenericPolycone::CalculateExtent()",
499 "GeomMgt1002", JustWarning, message);
500 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
501 }
502
503 // set trigonometric values
504 const G4int NSTEPS = 24; // number of steps for whole circle
505 G4double astep = twopi/NSTEPS; // max angle for one step
506
507 G4double sphi = GetStartPhi();
508 G4double ephi = GetEndPhi();
509 G4double dphi = IsOpen() ? ephi-sphi : twopi;
510 G4int ksteps = (dphi <= astep) ? 1 : (G4int)((dphi-deg)/astep) + 1;
511 G4double ang = dphi/ksteps;
512
513 G4double sinHalf = std::sin(0.5*ang);
514 G4double cosHalf = std::cos(0.5*ang);
515 G4double sinStep = 2.*sinHalf*cosHalf;
516 G4double cosStep = 1. - 2.*sinHalf*sinHalf;
517
518 G4double sinStart = GetSinStartPhi();
519 G4double cosStart = GetCosStartPhi();
520 G4double sinEnd = GetSinEndPhi();
521 G4double cosEnd = GetCosEndPhi();
522
523 // define vectors and arrays
524 std::vector<const G4ThreeVectorList *> polygons;
525 polygons.resize(ksteps+2);
526 G4ThreeVectorList pols[NSTEPS+2];
527 for (G4int k=0; k<ksteps+2; ++k) pols[k].resize(6);
528 for (G4int k=0; k<ksteps+2; ++k) polygons[k] = &pols[k];
529 G4double r0[6],z0[6]; // contour with original edges of triangle
530 G4double r1[6]; // shifted radii of external edges of triangle
531
532 // main loop along triangles
533 pMin = kInfinity;
534 pMax =-kInfinity;
535 G4int ntria = triangles.size()/3;
536 for (G4int i=0; i<ntria; ++i)
537 {
538 G4int i3 = i*3;
539 for (G4int k=0; k<3; ++k)
540 {
541 G4int e0 = i3+k, e1 = (k<2) ? e0+1 : i3;
542 G4int k2 = k*2;
543 // set contour with original edges of triangle
544 r0[k2+0] = triangles[e0].x(); z0[k2+0] = triangles[e0].y();
545 r0[k2+1] = triangles[e1].x(); z0[k2+1] = triangles[e1].y();
546 // set shifted radii
547 r1[k2+0] = r0[k2+0];
548 r1[k2+1] = r0[k2+1];
549 if (z0[k2+1] - z0[k2+0] <= 0) continue;
550 r1[k2+0] /= cosHalf;
551 r1[k2+1] /= cosHalf;
552 }
553
554 // rotate countour, set sequence of 6-sided polygons
555 G4double sinCur = sinStart*cosHalf + cosStart*sinHalf;
556 G4double cosCur = cosStart*cosHalf - sinStart*sinHalf;
557 for (G4int j=0; j<6; ++j)
558 {
559 pols[0][j].set(r0[j]*cosStart,r0[j]*sinStart,z0[j]);
560 }
561 for (G4int k=1; k<ksteps+1; ++k)
562 {
563 for (G4int j=0; j<6; ++j)
564 {
565 pols[k][j].set(r1[j]*cosCur,r1[j]*sinCur,z0[j]);
566 }
567 G4double sinTmp = sinCur;
568 sinCur = sinCur*cosStep + cosCur*sinStep;
569 cosCur = cosCur*cosStep - sinTmp*sinStep;
570 }
571 for (G4int j=0; j<6; ++j)
572 {
573 pols[ksteps+1][j].set(r0[j]*cosEnd,r0[j]*sinEnd,z0[j]);
574 }
575
576 // set sub-envelope and adjust extent
577 G4double emin,emax;
578 G4BoundingEnvelope benv(polygons);
579 if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
580 if (emin < pMin) pMin = emin;
581 if (emax > pMax) pMax = emax;
582 if (eminlim > pMin && emaxlim < pMax) return true; // max possible extent
583 }
584 return (pMin < pMax);
585}
std::vector< G4ThreeVector > G4ThreeVectorList
std::vector< G4TwoVector > G4TwoVectorList
Definition: G4GeomTools.hh:42
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:36
bool G4bool
Definition: G4Types.hh:86
G4double GetStartPhi() const
G4double GetEndPhi() const
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
static G4bool TriangulatePolygon(const G4TwoVectorList &polygon, G4TwoVectorList &result)
Definition: G4GeomTools.cc:193
static G4double PolygonArea(const G4TwoVectorList &polygon)
Definition: G4GeomTools.cc:76
G4double GetMinExtent(const EAxis pAxis) const
G4double GetMaxExtent(const EAxis pAxis) const

◆ Clone()

G4VSolid * G4GenericPolycone::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 598 of file G4GenericPolycone.cc.

599{
600 return new G4GenericPolycone(*this);
601}

◆ CopyStuff()

void G4GenericPolycone::CopyStuff ( const G4GenericPolycone source)
protected

Definition at line 300 of file G4GenericPolycone.cc.

301{
302 //
303 // Simple stuff
304 //
305 startPhi = source.startPhi;
306 endPhi = source.endPhi;
307 phiIsOpen = source.phiIsOpen;
308 numCorner = source.numCorner;
309
310 //
311 // The corner array
312 //
314
316 *sourceCorn = source.corners;
317 do // Loop checking, 13.08.2015, G.Cosmo
318 {
319 *corn = *sourceCorn;
320 } while( ++sourceCorn, ++corn < corners+numCorner );
321
322 //
323 // Enclosing cylinder
324 //
326
327 //
328 // Surface elements
329 //
330 delete fElements;
331 fElements = nullptr;
332
333 // Polyhedron
334 //
335 fRebuildPolyhedron = false;
336 delete fpPolyhedron;
337 fpPolyhedron = nullptr;
338}
G4bool fRebuildPolyhedron

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

◆ Create()

void G4GenericPolycone::Create ( G4double  phiStart,
G4double  phiTotal,
G4ReduciblePolygon rz 
)
protected

Definition at line 84 of file G4GenericPolycone.cc.

87{
88 //
89 // Perform checks of rz values
90 //
91 if (rz->Amin() < 0.0)
92 {
93 std::ostringstream message;
94 message << "Illegal input parameters - " << GetName() << G4endl
95 << " All R values must be >= 0 !";
96 G4Exception("G4GenericPolycone::Create()", "GeomSolids0002",
97 FatalErrorInArgument, message);
98 }
99
100 G4double rzArea = rz->Area();
101 if (rzArea < -kCarTolerance)
102 {
103 rz->ReverseOrder();
104 }
105 else if (rzArea < kCarTolerance)
106 {
107 std::ostringstream message;
108 message << "Illegal input parameters - " << GetName() << G4endl
109 << " R/Z cross section is zero or near zero: " << rzArea;
110 G4Exception("G4GenericPolycone::Create()", "GeomSolids0002",
111 FatalErrorInArgument, message);
112 }
113
116 {
117 std::ostringstream message;
118 message << "Illegal input parameters - " << GetName() << G4endl
119 << " Too few unique R/Z values !";
120 G4Exception("G4GenericPolycone::Create()", "GeomSolids0002",
121 FatalErrorInArgument, message);
122 }
123
124 if (rz->CrossesItself(1/kInfinity))
125 {
126 std::ostringstream message;
127 message << "Illegal input parameters - " << GetName() << G4endl
128 << " R/Z segments cross !";
129 G4Exception("G4GenericPolycone::Create()", "GeomSolids0002",
130 FatalErrorInArgument, message);
131 }
132
133 numCorner = rz->NumVertices();
134
135 //
136 // Phi opening? Account for some possible roundoff, and interpret
137 // nonsense value as representing no phi opening
138 //
139 if (phiTotal <= 0 || phiTotal > twopi-1E-10)
140 {
141 phiIsOpen = false;
142 startPhi = 0;
143 endPhi = twopi;
144 }
145 else
146 {
147 phiIsOpen = true;
148
149 //
150 // Convert phi into our convention
151 //
152 startPhi = phiStart;
153 while( startPhi < 0 ) // Loop checking, 13.08.2015, G.Cosmo
154 startPhi += twopi;
155
156 endPhi = phiStart+phiTotal;
157 while( endPhi < startPhi ) // Loop checking, 13.08.2015, G.Cosmo
158 endPhi += twopi;
159 }
160
161 //
162 // Allocate corner array.
163 //
165
166 //
167 // Copy corners
168 //
170
172 iterRZ.Begin();
173 do // Loop checking, 13.08.2015, G.Cosmo
174 {
175 next->r = iterRZ.GetA();
176 next->z = iterRZ.GetB();
177 } while( ++next, iterRZ.Next() );
178
179 //
180 // Allocate face pointer array
181 //
183 faces = new G4VCSGface*[numFace];
184
185 //
186 // Construct conical faces
187 //
188 // But! Don't construct a face if both points are at zero radius!
189 //
190 G4PolyconeSideRZ *corner = corners,
191 *prev = corners + numCorner-1,
192 *nextNext;
193 G4VCSGface** face = faces;
194 do // Loop checking, 13.08.2015, G.Cosmo
195 {
196 next = corner+1;
197 if (next >= corners+numCorner) next = corners;
198 nextNext = next+1;
199 if (nextNext >= corners+numCorner) nextNext = corners;
200
201 if (corner->r < 1/kInfinity && next->r < 1/kInfinity) continue;
202
203 //
204 // We must decide here if we can dare declare one of our faces
205 // as having a "valid" normal (i.e. allBehind = true). This
206 // is never possible if the face faces "inward" in r.
207 //
208 G4bool allBehind;
209 if (corner->z > next->z)
210 {
211 allBehind = false;
212 }
213 else
214 {
215 //
216 // Otherwise, it is only true if the line passing
217 // through the two points of the segment do not
218 // split the r/z cross section
219 //
220 allBehind = !rz->BisectedBy( corner->r, corner->z,
221 next->r, next->z, kCarTolerance );
222 }
223
224 *face++ = new G4PolyconeSide( prev, corner, next, nextNext,
225 startPhi, endPhi-startPhi, phiIsOpen, allBehind );
226 } while( prev=corner, corner=next, corner > corners );
227
228 if (phiIsOpen)
229 {
230 //
231 // Construct phi open edges
232 //
233 *face++ = new G4PolyPhiFace( rz, startPhi, 0, endPhi );
234 *face++ = new G4PolyPhiFace( rz, endPhi, 0, startPhi );
235 }
236
237 //
238 // We might have dropped a face or two: recalculate numFace
239 //
240 numFace = face-faces;
241
242 //
243 // Make enclosingCylinder
244 //
246 new G4EnclosingCylinder( rz, phiIsOpen, phiStart, phiTotal );
247}
@ FatalErrorInArgument
#define G4endl
Definition: G4ios.hh:57
G4bool BisectedBy(G4double a1, G4double b1, G4double a2, G4double b2, G4double tolerance)
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:302

Referenced by G4GenericPolycone().

◆ CreatePolyhedron()

G4Polyhedron * G4GenericPolycone::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 825 of file G4GenericPolycone.cc.

826{
827 // The following code prepares for:
828 // HepPolyhedron::createPolyhedron(int Nnodes, int Nfaces,
829 // const double xyz[][3],
830 // const int faces_vec[][4])
831 // Here is an extract from the header file HepPolyhedron.h:
832 /**
833 * Creates user defined polyhedron.
834 * This function allows to the user to define arbitrary polyhedron.
835 * The faces of the polyhedron should be either triangles or planar
836 * quadrilateral. Nodes of a face are defined by indexes pointing to
837 * the elements in the xyz array. Numeration of the elements in the
838 * array starts from 1 (like in fortran). The indexes can be positive
839 * or negative. Negative sign means that the corresponding edge is
840 * invisible. The normal of the face should be directed to exterior
841 * of the polyhedron.
842 *
843 * @param Nnodes number of nodes
844 * @param Nfaces number of faces
845 * @param xyz nodes
846 * @param faces_vec faces (quadrilaterals or triangles)
847 * @return status of the operation - is non-zero in case of problem
848 */
849 const G4int numSide =
851 * (endPhi - startPhi) / twopi) + 1;
852 G4int nNodes;
853 G4int nFaces;
854 typedef G4double double3[3];
855 double3* xyz;
856 typedef G4int int4[4];
857 int4* faces_vec;
858 if (phiIsOpen)
859 {
860 // Triangulate open ends. Simple ear-chopping algorithm...
861 // I'm not sure how robust this algorithm is (J.Allison).
862 //
863 std::vector<G4bool> chopped(numCorner, false);
864 std::vector<G4int*> triQuads;
865 G4int remaining = numCorner;
866 G4int iStarter = 0;
867 while (remaining >= 3) // Loop checking, 13.08.2015, G.Cosmo
868 {
869 // Find unchopped corners...
870 //
871 G4int A = -1, B = -1, C = -1;
872 G4int iStepper = iStarter;
873 do // Loop checking, 13.08.2015, G.Cosmo
874 {
875 if (A < 0) { A = iStepper; }
876 else if (B < 0) { B = iStepper; }
877 else if (C < 0) { C = iStepper; }
878 do // Loop checking, 13.08.2015, G.Cosmo
879 {
880 if (++iStepper >= numCorner) { iStepper = 0; }
881 }
882 while (chopped[iStepper]);
883 }
884 while (C < 0 && iStepper != iStarter);
885
886 // Check triangle at B is pointing outward (an "ear").
887 // Sign of z cross product determines...
888 //
889 G4double BAr = corners[A].r - corners[B].r;
890 G4double BAz = corners[A].z - corners[B].z;
891 G4double BCr = corners[C].r - corners[B].r;
892 G4double BCz = corners[C].z - corners[B].z;
893 if (BAr * BCz - BAz * BCr < kCarTolerance)
894 {
895 G4int* tq = new G4int[3];
896 tq[0] = A + 1;
897 tq[1] = B + 1;
898 tq[2] = C + 1;
899 triQuads.push_back(tq);
900 chopped[B] = true;
901 --remaining;
902 }
903 else
904 {
905 do // Loop checking, 13.08.2015, G.Cosmo
906 {
907 if (++iStarter >= numCorner) { iStarter = 0; }
908 }
909 while (chopped[iStarter]);
910 }
911 }
912 // Transfer to faces...
913 //
914 nNodes = (numSide + 1) * numCorner;
915 nFaces = numSide * numCorner + 2 * triQuads.size();
916 faces_vec = new int4[nFaces];
917 G4int iface = 0;
918 G4int addition = numCorner * numSide;
919 G4int d = numCorner - 1;
920 for (G4int iEnd = 0; iEnd < 2; ++iEnd)
921 {
922 for (size_t i = 0; i < triQuads.size(); ++i)
923 {
924 // Negative for soft/auxiliary/normally invisible edges...
925 //
926 G4int a, b, c;
927 if (iEnd == 0)
928 {
929 a = triQuads[i][0];
930 b = triQuads[i][1];
931 c = triQuads[i][2];
932 }
933 else
934 {
935 a = triQuads[i][0] + addition;
936 b = triQuads[i][2] + addition;
937 c = triQuads[i][1] + addition;
938 }
939 G4int ab = std::abs(b - a);
940 G4int bc = std::abs(c - b);
941 G4int ca = std::abs(a - c);
942 faces_vec[iface][0] = (ab == 1 || ab == d)? a: -a;
943 faces_vec[iface][1] = (bc == 1 || bc == d)? b: -b;
944 faces_vec[iface][2] = (ca == 1 || ca == d)? c: -c;
945 faces_vec[iface][3] = 0;
946 ++iface;
947 }
948 }
949
950 // Continue with sides...
951
952 xyz = new double3[nNodes];
953 const G4double dPhi = (endPhi - startPhi) / numSide;
954 G4double phi = startPhi;
955 G4int ixyz = 0;
956 for (G4int iSide = 0; iSide < numSide; ++iSide)
957 {
958 for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
959 {
960 xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
961 xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
962 xyz[ixyz][2] = corners[iCorner].z;
963 if (iSide == 0) // startPhi
964 {
965 if (iCorner < numCorner - 1)
966 {
967 faces_vec[iface][0] = ixyz + 1;
968 faces_vec[iface][1] = -(ixyz + numCorner + 1);
969 faces_vec[iface][2] = ixyz + numCorner + 2;
970 faces_vec[iface][3] = ixyz + 2;
971 }
972 else
973 {
974 faces_vec[iface][0] = ixyz + 1;
975 faces_vec[iface][1] = -(ixyz + numCorner + 1);
976 faces_vec[iface][2] = ixyz + 2;
977 faces_vec[iface][3] = ixyz - numCorner + 2;
978 }
979 }
980 else if (iSide == numSide - 1) // endPhi
981 {
982 if (iCorner < numCorner - 1)
983 {
984 faces_vec[iface][0] = ixyz + 1;
985 faces_vec[iface][1] = ixyz + numCorner + 1;
986 faces_vec[iface][2] = ixyz + numCorner + 2;
987 faces_vec[iface][3] = -(ixyz + 2);
988 }
989 else
990 {
991 faces_vec[iface][0] = ixyz + 1;
992 faces_vec[iface][1] = ixyz + numCorner + 1;
993 faces_vec[iface][2] = ixyz + 2;
994 faces_vec[iface][3] = -(ixyz - numCorner + 2);
995 }
996 }
997 else
998 {
999 if (iCorner < numCorner - 1)
1000 {
1001 faces_vec[iface][0] = ixyz + 1;
1002 faces_vec[iface][1] = -(ixyz + numCorner + 1);
1003 faces_vec[iface][2] = ixyz + numCorner + 2;
1004 faces_vec[iface][3] = -(ixyz + 2);
1005 }
1006 else
1007 {
1008 faces_vec[iface][0] = ixyz + 1;
1009 faces_vec[iface][1] = -(ixyz + numCorner + 1);
1010 faces_vec[iface][2] = ixyz + 2;
1011 faces_vec[iface][3] = -(ixyz - numCorner + 2);
1012 }
1013 }
1014 ++iface;
1015 ++ixyz;
1016 }
1017 phi += dPhi;
1018 }
1019
1020 // Last corners...
1021
1022 for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1023 {
1024 xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1025 xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1026 xyz[ixyz][2] = corners[iCorner].z;
1027 ++ixyz;
1028 }
1029 }
1030 else // !phiIsOpen - i.e., a complete 360 degrees.
1031 {
1032 nNodes = numSide * numCorner;
1033 nFaces = numSide * numCorner;;
1034 xyz = new double3[nNodes];
1035 faces_vec = new int4[nFaces];
1036 const G4double dPhi = (endPhi - startPhi) / numSide;
1037 G4double phi = startPhi;
1038 G4int ixyz = 0, iface = 0;
1039 for (G4int iSide = 0; iSide < numSide; ++iSide)
1040 {
1041 for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1042 {
1043 xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1044 xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1045 xyz[ixyz][2] = corners[iCorner].z;
1046
1047 if (iSide < numSide - 1)
1048 {
1049 if (iCorner < numCorner - 1)
1050 {
1051 faces_vec[iface][0] = ixyz + 1;
1052 faces_vec[iface][1] = -(ixyz + numCorner + 1);
1053 faces_vec[iface][2] = ixyz + numCorner + 2;
1054 faces_vec[iface][3] = -(ixyz + 2);
1055 }
1056 else
1057 {
1058 faces_vec[iface][0] = ixyz + 1;
1059 faces_vec[iface][1] = -(ixyz + numCorner + 1);
1060 faces_vec[iface][2] = ixyz + 2;
1061 faces_vec[iface][3] = -(ixyz - numCorner + 2);
1062 }
1063 }
1064 else // Last side joins ends...
1065 {
1066 if (iCorner < numCorner - 1)
1067 {
1068 faces_vec[iface][0] = ixyz + 1;
1069 faces_vec[iface][1] = -(ixyz + numCorner - nFaces + 1);
1070 faces_vec[iface][2] = ixyz + numCorner - nFaces + 2;
1071 faces_vec[iface][3] = -(ixyz + 2);
1072 }
1073 else
1074 {
1075 faces_vec[iface][0] = ixyz + 1;
1076 faces_vec[iface][1] = -(ixyz - nFaces + numCorner + 1);
1077 faces_vec[iface][2] = ixyz - nFaces + 2;
1078 faces_vec[iface][3] = -(ixyz - numCorner + 2);
1079 }
1080 }
1081 ++ixyz;
1082 ++iface;
1083 }
1084 phi += dPhi;
1085 }
1086 }
1087 G4Polyhedron* polyhedron = new G4Polyhedron;
1088 G4int prob = polyhedron->createPolyhedron(nNodes, nFaces, xyz, faces_vec);
1089 delete [] faces_vec;
1090 delete [] xyz;
1091 if (prob)
1092 {
1093 std::ostringstream message;
1094 message << "Problem creating G4Polyhedron for: " << GetName();
1095 G4Exception("G4GenericPolycone::CreatePolyhedron()", "GeomSolids1002",
1096 JustWarning, message);
1097 delete polyhedron;
1098 return nullptr;
1099 }
1100 else
1101 {
1102 return polyhedron;
1103 }
1104}
double B(double temperature)
double C(double temp)
double A(double temperature)
static G4int GetNumberOfRotationSteps()
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])

◆ DistanceToIn() [1/2]

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

Reimplemented from G4VCSGfaceted.

Definition at line 392 of file G4GenericPolycone.cc.

393{
395}
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const

◆ DistanceToIn() [2/2]

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

Reimplemented from G4VCSGfaceted.

Definition at line 375 of file G4GenericPolycone.cc.

377{
378 //
379 // Quick test
380 //
382 return kInfinity;
383
384 //
385 // Long answer
386 //
387 return G4VCSGfaceted::DistanceToIn( p, v );
388}
G4bool ShouldMiss(const G4ThreeVector &p, const G4ThreeVector &v) const

◆ GetCorner()

◆ GetCosEndPhi()

G4double G4GenericPolycone::GetCosEndPhi ( ) const
inline

Referenced by BoundingLimits(), and CalculateExtent().

◆ GetCosStartPhi()

G4double G4GenericPolycone::GetCosStartPhi ( ) const
inline

Referenced by BoundingLimits(), and CalculateExtent().

◆ GetCubicVolume()

G4double G4GenericPolycone::GetCubicVolume ( )
virtual

Reimplemented from G4VCSGfaceted.

Definition at line 636 of file G4GenericPolycone.cc.

637{
638 if (fCubicVolume == 0.)
639 {
640 G4double total = 0.;
641 G4int nrz = GetNumRZCorner();
642 G4PolyconeSideRZ a = GetCorner(nrz - 1);
643 for (G4int i=0; i<nrz; ++i)
644 {
646 total += (b.r*b.r + b.r*a.r + a.r*a.r)*(b.z - a.z);
647 a = b;
648 }
649 fCubicVolume = std::abs(total)*(GetEndPhi() - GetStartPhi())/6.;
650 }
651 return fCubicVolume;
652}
G4double fCubicVolume
G4double total(Particle const *const p1, Particle const *const p2)

◆ GetEndPhi()

◆ GetEntityType()

G4GeometryType G4GenericPolycone::GetEntityType ( ) const
virtual

Reimplemented from G4VCSGfaceted.

Definition at line 589 of file G4GenericPolycone.cc.

590{
591 return G4String("G4GenericPolycone");
592}

◆ GetNumRZCorner()

◆ GetPointOnSurface()

G4ThreeVector G4GenericPolycone::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 760 of file G4GenericPolycone.cc.

761{
762 // Set surface elements
763 if (!fElements)
764 {
765 G4AutoLock l(&surface_elementsMutex);
767 l.unlock();
768 }
769
770 // Select surface element
772 selem = fElements->back();
773 G4double select = selem.area*G4QuickRand();
774 auto it = std::lower_bound(fElements->begin(), fElements->end(), select,
776 -> G4bool { return x.area < val; });
777
778 // Generate random point
779 G4double r = 0, z = 0, phi = 0;
780 G4double u = G4QuickRand();
781 G4double v = G4QuickRand();
782 G4int i0 = (*it).i0;
783 G4int i1 = (*it).i1;
784 G4int i2 = (*it).i2;
785 if (i2 < 0) // lateral surface
786 {
789 if (p1.r < p0.r)
790 {
791 p0 = GetCorner(i1);
792 p1 = GetCorner(i0);
793 }
794 if (p1.r - p0.r < kCarTolerance) // cylindrical surface
795 {
796 r = (p1.r - p0.r)*u + p0.r;
797 z = (p1.z - p0.z)*u + p0.z;
798 }
799 else // conical surface
800 {
801 r = std::sqrt(p1.r*p1.r*u + p0.r*p0.r*(1. - u));
802 z = p0.z + (p1.z - p0.z)*(r - p0.r)/(p1.r - p0.r);
803 }
804 phi = (GetEndPhi() - GetStartPhi())*v + GetStartPhi();
805 }
806 else // phi cut
807 {
808 G4int nrz = GetNumRZCorner();
809 phi = (i0 < nrz) ? GetStartPhi() : GetEndPhi();
810 if (i0 >= nrz) { i0 -= nrz; }
814 if (u + v > 1.) { u = 1. - u; v = 1. - v; }
815 r = (p1.r - p0.r)*u + (p2.r - p0.r)*v + p0.r;
816 z = (p1.z - p0.z)*u + (p2.z - p0.z)*v + p0.z;
817 }
818 return G4ThreeVector(r*std::cos(phi), r*std::sin(phi), z);
819}
G4double G4QuickRand()
Definition: G4QuickRand.hh:34
CLHEP::Hep3Vector G4ThreeVector
void SetSurfaceElements() const

◆ GetSinEndPhi()

G4double G4GenericPolycone::GetSinEndPhi ( ) const
inline

Referenced by BoundingLimits(), and CalculateExtent().

◆ GetSinStartPhi()

G4double G4GenericPolycone::GetSinStartPhi ( ) const
inline

Referenced by BoundingLimits(), and CalculateExtent().

◆ GetStartPhi()

◆ GetSurfaceArea()

G4double G4GenericPolycone::GetSurfaceArea ( )
virtual

Reimplemented from G4VCSGfaceted.

Definition at line 658 of file G4GenericPolycone.cc.

659{
660 if (fSurfaceArea == 0.)
661 {
662 // phi cut area
663 G4int nrz = GetNumRZCorner();
664 G4double scut = 0.;
665 if (IsOpen())
666 {
667 G4PolyconeSideRZ a = GetCorner(nrz - 1);
668 for (G4int i=0; i<nrz; ++i)
669 {
671 scut += a.r*b.z - a.z*b.r;
672 a = b;
673 }
674 scut = std::abs(scut);
675 }
676 // lateral surface area
677 G4double slat = 0;
678 G4PolyconeSideRZ a = GetCorner(nrz - 1);
679 for (G4int i=0; i<nrz; ++i)
680 {
682 G4double h = std::sqrt((b.r - a.r)*(b.r - a.r) + (b.z - a.z)*(b.z - a.z));
683 slat += (b.r + a.r)*h;
684 a = b;
685 }
686 slat *= (GetEndPhi() - GetStartPhi())/2.;
687 fSurfaceArea = scut + slat;
688 }
689 return fSurfaceArea;
690}
G4double fSurfaceArea

◆ Inside()

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

Reimplemented from G4VCSGfaceted.

Definition at line 357 of file G4GenericPolycone.cc.

358{
359 //
360 // Quick test
361 //
363
364 //
365 // Long answer
366 //
367 return G4VCSGfaceted::Inside(p);
368}
G4bool MustBeOutside(const G4ThreeVector &p) const
virtual EInside Inside(const G4ThreeVector &p) const
@ kOutside
Definition: geomdefs.hh:68

◆ IsOpen()

G4bool G4GenericPolycone::IsOpen ( ) const
inline

◆ operator=()

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

Definition at line 282 of file G4GenericPolycone.cc.

283{
284 if (this == &source) return *this;
285
286 G4VCSGfaceted::operator=( source );
287
288 delete [] corners;
289 // if (original_parameters) delete original_parameters;
290
291 delete enclosingCylinder;
292
293 CopyStuff( source );
294
295 return *this;
296}
G4VCSGfaceted & operator=(const G4VCSGfaceted &source)

◆ Reset()

G4bool G4GenericPolycone::Reset ( )

Definition at line 342 of file G4GenericPolycone.cc.

343{
344 std::ostringstream message;
345 message << "Solid " << GetName() << " built using generic construct."
346 << G4endl << "Not applicable to the generic construct !";
347 G4Exception("G4GenericPolycone::Reset()", "GeomSolids1001",
348 JustWarning, message, "Parameters NOT resetted.");
349 return true;
350}

◆ SetSurfaceElements()

void G4GenericPolycone::SetSurfaceElements ( ) const
protected

Definition at line 697 of file G4GenericPolycone.cc.

698{
699 fElements = new std::vector<G4GenericPolycone::surface_element>;
700 G4double sarea = 0.;
701 G4int nrz = GetNumRZCorner();
702
703 // set lateral surface elements
704 G4double dphi = GetEndPhi() - GetStartPhi();
705 G4int ia = nrz - 1;
706 for (G4int ib=0; ib<nrz; ++ib)
707 {
711 selem.i0 = ia;
712 selem.i1 = ib;
713 selem.i2 = -1;
714 ia = ib;
715 if (a.r == 0. && b.r == 0.) continue;
716 G4double h = std::sqrt((b.r - a.r)*(b.r - a.r) + (b.z - a.z)*(b.z - a.z));
717 sarea += 0.5*dphi*(b.r + a.r)*h;
718 selem.area = sarea;
719 fElements->push_back(selem);
720 }
721
722 // set elements for phi cuts
723 if (IsOpen())
724 {
725 G4TwoVectorList contourRZ;
726 std::vector<G4int> triangles;
727 for (G4int i=0; i<nrz; ++i)
728 {
729 G4PolyconeSideRZ corner = GetCorner(i);
730 contourRZ.push_back(G4TwoVector(corner.r, corner.z));
731 }
732 G4GeomTools::TriangulatePolygon(contourRZ, triangles);
733 G4int ntria = triangles.size();
734 for (G4int i=0; i<ntria; i+=3)
735 {
737 selem.i0 = triangles[i];
738 selem.i1 = triangles[i+1];
739 selem.i2 = triangles[i+2];
740 G4PolyconeSideRZ a = GetCorner(selem.i0);
741 G4PolyconeSideRZ b = GetCorner(selem.i1);
742 G4PolyconeSideRZ c = GetCorner(selem.i2);
743 G4double stria =
744 std::abs(G4GeomTools::TriangleArea(a.r, a.z, b.r, b.z, c.r, c.z));
745 sarea += stria;
746 selem.area = sarea;
747 fElements->push_back(selem); // start phi
748 sarea += stria;
749 selem.area = sarea;
750 selem.i0 += nrz;
751 fElements->push_back(selem); // end phi
752 }
753 }
754}
static G4double TriangleArea(G4double Ax, G4double Ay, G4double Bx, G4double By, G4double Cx, G4double Cy)
Definition: G4GeomTools.cc:41

Referenced by GetPointOnSurface().

◆ StreamInfo()

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

Reimplemented from G4VCSGfaceted.

Definition at line 607 of file G4GenericPolycone.cc.

608{
609 G4int oldprc = os.precision(16);
610 os << "-----------------------------------------------------------\n"
611 << " *** Dump for solid - " << GetName() << " ***\n"
612 << " ===================================================\n"
613 << " Solid type: G4GenericPolycone\n"
614 << " Parameters: \n"
615 << " starting phi angle : " << startPhi/degree << " degrees \n"
616 << " ending phi angle : " << endPhi/degree << " degrees \n";
617 G4int i=0;
618
619 os << " number of RZ points: " << numCorner << "\n"
620 << " RZ values (corners): \n";
621 for (i=0; i<numCorner; i++)
622 {
623 os << " "
624 << corners[i].r << ", " << corners[i].z << "\n";
625 }
626 os << "-----------------------------------------------------------\n";
627 os.precision(oldprc);
628
629 return os;
630}

Member Data Documentation

◆ corners

G4PolyconeSideRZ* G4GenericPolycone::corners = nullptr
protected

◆ enclosingCylinder

G4EnclosingCylinder* G4GenericPolycone::enclosingCylinder = nullptr
protected

◆ endPhi

G4double G4GenericPolycone::endPhi
protected

Definition at line 148 of file G4GenericPolycone.hh.

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

◆ fElements

std::vector<surface_element>* G4GenericPolycone::fElements = nullptr
mutableprotected

◆ numCorner

G4int G4GenericPolycone::numCorner
protected

Definition at line 150 of file G4GenericPolycone.hh.

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

◆ phiIsOpen

G4bool G4GenericPolycone::phiIsOpen = false
protected

Definition at line 149 of file G4GenericPolycone.hh.

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

◆ startPhi

G4double G4GenericPolycone::startPhi
protected

Definition at line 147 of file G4GenericPolycone.hh.

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


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