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

#include <G4Polycone.hh>

+ Inheritance diagram for G4Polycone:

Classes

struct  surface_element
 

Public Member Functions

 G4Polycone (const G4String &name, G4double phiStart, G4double phiTotal, G4int numZPlanes, const G4double zPlane[], const G4double rInner[], const G4double rOuter[])
 
 G4Polycone (const G4String &name, G4double phiStart, G4double phiTotal, G4int numRZ, const G4double r[], const G4double z[])
 
 ~G4Polycone () override
 
EInside Inside (const G4ThreeVector &p) const override
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const override
 
G4double DistanceToIn (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
 
G4double GetCubicVolume () override
 
G4double GetSurfaceArea () override
 
G4ThreeVector GetPointOnSurface () const override
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
 
G4GeometryType GetEntityType () const override
 
G4VSolidClone () const override
 
std::ostream & StreamInfo (std::ostream &os) const override
 
G4PolyhedronCreatePolyhedron () const override
 
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
 
G4PolyconeHistoricalGetOriginalParameters () const
 
void SetOriginalParameters (G4PolyconeHistorical *pars)
 
 G4Polycone (__void__ &)
 
 G4Polycone (const G4Polycone &source)
 
G4Polyconeoperator= (const G4Polycone &source)
 
- Public Member Functions inherited from G4VCSGfaceted
 G4VCSGfaceted (const G4String &name)
 
virtual ~G4VCSGfaceted ()
 
 G4VCSGfaceted (const G4VCSGfaceted &source)
 
G4VCSGfacetedoperator= (const G4VCSGfaceted &source)
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const override
 
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
 
G4GeometryType GetEntityType () const override
 
std::ostream & StreamInfo (std::ostream &os) const override
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const override
 
G4VisExtent GetExtent () const override
 
G4PolyhedronGetPolyhedron () const override
 
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)
 
G4double GetCubicVolume () override
 
G4double GetSurfaceArea () override
 
 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
 
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 Member Functions

G4bool SetOriginalParameters (G4ReduciblePolygon *rz)
 
void Create (G4double phiStart, G4double phiTotal, G4ReduciblePolygon *rz)
 
void CopyStuff (const G4Polycone &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
 
G4PolyconeHistoricaloriginal_parameters = 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 75 of file G4Polycone.hh.

Constructor & Destructor Documentation

◆ G4Polycone() [1/4]

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

Definition at line 58 of file G4Polycone.cc.

65 : G4VCSGfaceted( name )
66{
67 //
68 // Some historical ugliness
69 //
71
75 original_parameters->Z_values = new G4double[numZPlanes];
76 original_parameters->Rmin = new G4double[numZPlanes];
77 original_parameters->Rmax = new G4double[numZPlanes];
78
79 for (G4int i=0; i<numZPlanes; ++i)
80 {
81 if(rInner[i]>rOuter[i])
82 {
83 DumpInfo();
84 std::ostringstream message;
85 message << "Cannot create a Polycone with rInner > rOuter for the same Z"
86 << G4endl
87 << " rInner > rOuter for the same Z !" << G4endl
88 << " rMin[" << i << "] = " << rInner[i]
89 << " -- rMax[" << i << "] = " << rOuter[i];
90 G4Exception("G4Polycone::G4Polycone()", "GeomSolids0002",
91 FatalErrorInArgument, message);
92 }
93 if (( i < numZPlanes-1) && ( zPlane[i] == zPlane[i+1] ))
94 {
95 if( (rInner[i] > rOuter[i+1])
96 ||(rInner[i+1] > rOuter[i]) )
97 {
98 DumpInfo();
99 std::ostringstream message;
100 message << "Cannot create a Polycone with no contiguous segments."
101 << G4endl
102 << " Segments are not contiguous !" << G4endl
103 << " rMin[" << i << "] = " << rInner[i]
104 << " -- rMax[" << i+1 << "] = " << rOuter[i+1] << G4endl
105 << " rMin[" << i+1 << "] = " << rInner[i+1]
106 << " -- rMax[" << i << "] = " << rOuter[i];
107 G4Exception("G4Polycone::G4Polycone()", "GeomSolids0002",
108 FatalErrorInArgument, message);
109 }
110 }
111 original_parameters->Z_values[i] = zPlane[i];
112 original_parameters->Rmin[i] = rInner[i];
113 original_parameters->Rmax[i] = rOuter[i];
114 }
115
116 //
117 // Build RZ polygon using special PCON/PGON GEANT3 constructor
118 //
119 auto rz = new G4ReduciblePolygon( rInner, rOuter, zPlane, numZPlanes );
120
121 //
122 // Do the real work
123 //
124 Create( phiStart, phiTotal, rz );
125
126 delete rz;
127}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
void Create(G4double phiStart, G4double phiTotal, G4ReduciblePolygon *rz)
G4PolyconeHistorical * original_parameters
G4VCSGfaceted(const G4String &name)
void DumpInfo() const

Referenced by Clone().

◆ G4Polycone() [2/4]

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

Definition at line 131 of file G4Polycone.cc.

137 : G4VCSGfaceted( name )
138{
139
140 auto rz = new G4ReduciblePolygon( r, z, numRZ );
141
142 Create( phiStart, phiTotal, rz );
143
144 // Set original_parameters struct for consistency
145 //
146
147 G4bool convertible = SetOriginalParameters(rz);
148
149 if(!convertible)
150 {
151 std::ostringstream message;
152 message << "Polycone " << GetName() << "cannot be converted" << G4endl
153 << "to Polycone with (Rmin,Rmaz,Z) parameters!";
154 G4Exception("G4Polycone::G4Polycone()", "GeomSolids0002",
155 FatalException, message, "Use G4GenericPolycone instead!");
156 }
157 else
158 {
159 G4cout << "INFO: Converting polycone " << GetName() << G4endl
160 << "to optimized polycone with (Rmin,Rmaz,Z) parameters !"
161 << G4endl;
162 }
163 delete rz;
164}
@ FatalException
bool G4bool
Definition G4Types.hh:86
G4GLOB_DLL std::ostream G4cout
void SetOriginalParameters(G4PolyconeHistorical *pars)
G4String GetName() const

◆ ~G4Polycone()

G4Polycone::~G4Polycone ( )
override

Definition at line 341 of file G4Polycone.cc.

342{
343 delete [] corners;
344 delete original_parameters;
345 delete enclosingCylinder;
346 delete fElements;
347 delete fpPolyhedron;
348 corners = nullptr;
349 original_parameters = nullptr;
350 enclosingCylinder = nullptr;
351 fElements = nullptr;
352 fpPolyhedron = nullptr;
353}
G4EnclosingCylinder * enclosingCylinder
G4PolyconeSideRZ * corners
std::vector< surface_element > * fElements
G4Polyhedron * fpPolyhedron

◆ G4Polycone() [3/4]

G4Polycone::G4Polycone ( __void__ & a)

Definition at line 332 of file G4Polycone.cc.

333 : G4VCSGfaceted(a), startPhi(0.), endPhi(0.), numCorner(0)
334{
335}
G4double endPhi
G4int numCorner
G4double startPhi

◆ G4Polycone() [4/4]

G4Polycone::G4Polycone ( const G4Polycone & source)

Definition at line 357 of file G4Polycone.cc.

358 : G4VCSGfaceted( source )
359{
360 CopyStuff( source );
361}
void CopyStuff(const G4Polycone &source)

Member Function Documentation

◆ BoundingLimits()

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

Reimplemented from G4VSolid.

Definition at line 509 of file G4Polycone.cc.

511{
512 G4double rmin = kInfinity, rmax = -kInfinity;
513 G4double zmin = kInfinity, zmax = -kInfinity;
514
515 for (G4int i=0; i<GetNumRZCorner(); ++i)
516 {
517 G4PolyconeSideRZ corner = GetCorner(i);
518 if (corner.r < rmin) rmin = corner.r;
519 if (corner.r > rmax) rmax = corner.r;
520 if (corner.z < zmin) zmin = corner.z;
521 if (corner.z > zmax) zmax = corner.z;
522 }
523
524 if (IsOpen())
525 {
526 G4TwoVector vmin,vmax;
527 G4GeomTools::DiskExtent(rmin,rmax,
530 vmin,vmax);
531 pMin.set(vmin.x(),vmin.y(),zmin);
532 pMax.set(vmax.x(),vmax.y(),zmax);
533 }
534 else
535 {
536 pMin.set(-rmax,-rmax, zmin);
537 pMax.set( rmax, rmax, zmax);
538 }
539
540 // Check correctness of the bounding box
541 //
542 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
543 {
544 std::ostringstream message;
545 message << "Bad bounding box (min >= max) for solid: "
546 << GetName() << " !"
547 << "\npMin = " << pMin
548 << "\npMax = " << pMax;
549 G4Exception("G4Polycone::BoundingLimits()", "GeomMgt0001",
550 JustWarning, message);
551 DumpInfo();
552 }
553}
@ JustWarning
double x() const
double y() const
double z() const
double x() const
double y() const
void set(double x, double y, double z)
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
G4double GetCosEndPhi() const
G4double GetSinEndPhi() const
G4double GetCosStartPhi() const
G4bool IsOpen() const
G4double GetSinStartPhi() const
G4int GetNumRZCorner() const
G4PolyconeSideRZ GetCorner(G4int index) const

Referenced by CalculateExtent().

◆ CalculateExtent()

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

Implements G4VSolid.

Definition at line 557 of file G4Polycone.cc.

561{
562 G4ThreeVector bmin, bmax;
563 G4bool exist;
564
565 // Check bounding box (bbox)
566 //
567 BoundingLimits(bmin,bmax);
568 G4BoundingEnvelope bbox(bmin,bmax);
569#ifdef G4BBOX_EXTENT
570 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
571#endif
572 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
573 {
574 return exist = pMin < pMax;
575 }
576
577 // To find the extent, RZ contour of the polycone is subdivided
578 // in triangles. The extent is calculated as cumulative extent of
579 // all sub-polycones formed by rotation of triangles around Z
580 //
581 G4TwoVectorList contourRZ;
582 G4TwoVectorList triangles;
583 std::vector<G4int> iout;
584 G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
585 G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
586
587 // get RZ contour, ensure anticlockwise order of corners
588 for (G4int i=0; i<GetNumRZCorner(); ++i)
589 {
590 G4PolyconeSideRZ corner = GetCorner(i);
591 contourRZ.emplace_back(corner.r,corner.z);
592 }
594 G4double area = G4GeomTools::PolygonArea(contourRZ);
595 if (area < 0.) std::reverse(contourRZ.begin(),contourRZ.end());
596
597 // triangulate RZ countour
598 if (!G4GeomTools::TriangulatePolygon(contourRZ,triangles))
599 {
600 std::ostringstream message;
601 message << "Triangulation of RZ contour has failed for solid: "
602 << GetName() << " !"
603 << "\nExtent has been calculated using boundary box";
604 G4Exception("G4Polycone::CalculateExtent()",
605 "GeomMgt1002", JustWarning, message);
606 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
607 }
608
609 // set trigonometric values
610 const G4int NSTEPS = 24; // number of steps for whole circle
611 G4double astep = twopi/NSTEPS; // max angle for one step
612
613 G4double sphi = GetStartPhi();
614 G4double ephi = GetEndPhi();
615 G4double dphi = IsOpen() ? ephi-sphi : twopi;
616 G4int ksteps = (dphi <= astep) ? 1 : (G4int)((dphi-deg)/astep) + 1;
617 G4double ang = dphi/ksteps;
618
619 G4double sinHalf = std::sin(0.5*ang);
620 G4double cosHalf = std::cos(0.5*ang);
621 G4double sinStep = 2.*sinHalf*cosHalf;
622 G4double cosStep = 1. - 2.*sinHalf*sinHalf;
623
624 G4double sinStart = GetSinStartPhi();
625 G4double cosStart = GetCosStartPhi();
626 G4double sinEnd = GetSinEndPhi();
627 G4double cosEnd = GetCosEndPhi();
628
629 // define vectors and arrays
630 std::vector<const G4ThreeVectorList *> polygons;
631 polygons.resize(ksteps+2);
632 G4ThreeVectorList pols[NSTEPS+2];
633 for (G4int k=0; k<ksteps+2; ++k) pols[k].resize(6);
634 for (G4int k=0; k<ksteps+2; ++k) polygons[k] = &pols[k];
635 G4double r0[6],z0[6]; // contour with original edges of triangle
636 G4double r1[6]; // shifted radii of external edges of triangle
637
638 // main loop along triangles
639 pMin = kInfinity;
640 pMax =-kInfinity;
641 G4int ntria = (G4int)triangles.size()/3;
642 for (G4int i=0; i<ntria; ++i)
643 {
644 G4int i3 = i*3;
645 for (G4int k=0; k<3; ++k)
646 {
647 G4int e0 = i3+k, e1 = (k<2) ? e0+1 : i3;
648 G4int k2 = k*2;
649 // set contour with original edges of triangle
650 r0[k2+0] = triangles[e0].x(); z0[k2+0] = triangles[e0].y();
651 r0[k2+1] = triangles[e1].x(); z0[k2+1] = triangles[e1].y();
652 // set shifted radii
653 r1[k2+0] = r0[k2+0];
654 r1[k2+1] = r0[k2+1];
655 if (z0[k2+1] - z0[k2+0] <= 0) continue;
656 r1[k2+0] /= cosHalf;
657 r1[k2+1] /= cosHalf;
658 }
659
660 // rotate countour, set sequence of 6-sided polygons
661 G4double sinCur = sinStart*cosHalf + cosStart*sinHalf;
662 G4double cosCur = cosStart*cosHalf - sinStart*sinHalf;
663 for (G4int j=0; j<6; ++j) pols[0][j].set(r0[j]*cosStart,r0[j]*sinStart,z0[j]);
664 for (G4int k=1; k<ksteps+1; ++k)
665 {
666 for (G4int j=0; j<6; ++j) pols[k][j].set(r1[j]*cosCur,r1[j]*sinCur,z0[j]);
667 G4double sinTmp = sinCur;
668 sinCur = sinCur*cosStep + cosCur*sinStep;
669 cosCur = cosCur*cosStep - sinTmp*sinStep;
670 }
671 for (G4int j=0; j<6; ++j) pols[ksteps+1][j].set(r0[j]*cosEnd,r0[j]*sinEnd,z0[j]);
672
673 // set sub-envelope and adjust extent
674 G4double emin,emax;
675 G4BoundingEnvelope benv(polygons);
676 if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
677 if (emin < pMin) pMin = emin;
678 if (emax > pMax) pMax = emax;
679 if (eminlim > pMin && emaxlim < pMax) return true; // max possible extent
680 }
681 return (pMin < pMax);
682}
std::vector< G4ThreeVector > G4ThreeVectorList
std::vector< G4TwoVector > G4TwoVectorList
static G4bool TriangulatePolygon(const G4TwoVectorList &polygon, G4TwoVectorList &result)
static void RemoveRedundantVertices(G4TwoVectorList &polygon, std::vector< G4int > &iout, G4double tolerance=0.0)
static G4double PolygonArea(const G4TwoVectorList &polygon)
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
G4double GetEndPhi() const
G4double GetStartPhi() const
G4double kCarTolerance
Definition G4VSolid.hh:299
G4double GetMinExtent(const EAxis pAxis) const
G4double GetMaxExtent(const EAxis pAxis) const

◆ Clone()

G4VSolid * G4Polycone::Clone ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 702 of file G4Polycone.cc.

703{
704 return new G4Polycone(*this);
705}
G4Polycone(const G4String &name, G4double phiStart, G4double phiTotal, G4int numZPlanes, const G4double zPlane[], const G4double rInner[], const G4double rOuter[])
Definition G4Polycone.cc:58

◆ ComputeDimensions()

void G4Polycone::ComputeDimensions ( G4VPVParameterisation * p,
const G4int n,
const G4VPhysicalVolume * pRep )
overridevirtual

Reimplemented from G4VSolid.

Definition at line 686 of file G4Polycone.cc.

689{
690 p->ComputeDimensions(*this,n,pRep);
691}
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const

◆ CopyStuff()

void G4Polycone::CopyStuff ( const G4Polycone & source)
protected

Definition at line 383 of file G4Polycone.cc.

384{
385 //
386 // Simple stuff
387 //
388 startPhi = source.startPhi;
389 endPhi = source.endPhi;
390 phiIsOpen = source.phiIsOpen;
391 numCorner = source.numCorner;
392
393 //
394 // The corner array
395 //
397
399 * sourceCorn = source.corners;
400 do // Loop checking, 13.08.2015, G.Cosmo
401 {
402 *corn = *sourceCorn;
403 } while( ++sourceCorn, ++corn < corners+numCorner );
404
405 //
406 // Original parameters
407 //
408 if (source.original_parameters != nullptr)
409 {
412 }
413
414 //
415 // Enclosing cylinder
416 //
418
419 //
420 // Surface elements
421 //
422 delete fElements;
423 fElements = nullptr;
424
425 //
426 // Polyhedron
427 //
428 fRebuildPolyhedron = false;
429 delete fpPolyhedron;
430 fpPolyhedron = nullptr;
431}
G4bool phiIsOpen
G4bool fRebuildPolyhedron

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

◆ Create()

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

Definition at line 171 of file G4Polycone.cc.

174{
175 //
176 // Perform checks of rz values
177 //
178 if (rz->Amin() < 0.0)
179 {
180 std::ostringstream message;
181 message << "Illegal input parameters - " << GetName() << G4endl
182 << " All R values must be >= 0 !";
183 G4Exception("G4Polycone::Create()", "GeomSolids0002",
184 FatalErrorInArgument, message);
185 }
186
187 G4double rzArea = rz->Area();
188 if (rzArea < -kCarTolerance)
189 {
190 rz->ReverseOrder();
191 }
192 else if (rzArea < kCarTolerance)
193 {
194 std::ostringstream message;
195 message << "Illegal input parameters - " << GetName() << G4endl
196 << " R/Z cross section is zero or near zero: " << rzArea;
197 G4Exception("G4Polycone::Create()", "GeomSolids0002",
198 FatalErrorInArgument, message);
199 }
200
203 {
204 std::ostringstream message;
205 message << "Illegal input parameters - " << GetName() << G4endl
206 << " Too few unique R/Z values !";
207 G4Exception("G4Polycone::Create()", "GeomSolids0002",
208 FatalErrorInArgument, message);
209 }
210
211 if (rz->CrossesItself(1/kInfinity))
212 {
213 std::ostringstream message;
214 message << "Illegal input parameters - " << GetName() << G4endl
215 << " R/Z segments cross !";
216 G4Exception("G4Polycone::Create()", "GeomSolids0002",
217 FatalErrorInArgument, message);
218 }
219
220 numCorner = rz->NumVertices();
221
222 startPhi = phiStart;
223 while( startPhi < 0. ) // Loop checking, 13.08.2015, G.Cosmo
224 startPhi += twopi;
225 //
226 // Phi opening? Account for some possible roundoff, and interpret
227 // nonsense value as representing no phi opening
228 //
229 if ( (phiTotal <= 0) || (phiTotal > twopi*(1-DBL_EPSILON)) )
230 {
231 phiIsOpen = false;
232 startPhi = 0.;
233 endPhi = twopi;
234 }
235 else
236 {
237 phiIsOpen = true;
238 endPhi = startPhi + phiTotal;
239 }
240
241 //
242 // Allocate corner array.
243 //
245
246 //
247 // Copy corners
248 //
250
252 iterRZ.Begin();
253 do // Loop checking, 13.08.2015, G.Cosmo
254 {
255 next->r = iterRZ.GetA();
256 next->z = iterRZ.GetB();
257 } while( ++next, iterRZ.Next() );
258
259 //
260 // Allocate face pointer array
261 //
263 faces = new G4VCSGface*[numFace];
264
265 //
266 // Construct conical faces
267 //
268 // But! Don't construct a face if both points are at zero radius!
269 //
270 G4PolyconeSideRZ* corner = corners,
271 * prev = corners + numCorner-1,
272 * nextNext;
273 G4VCSGface **face = faces;
274 do // Loop checking, 13.08.2015, G.Cosmo
275 {
276 next = corner+1;
277 if (next >= corners+numCorner) next = corners;
278 nextNext = next+1;
279 if (nextNext >= corners+numCorner) nextNext = corners;
280
281 if (corner->r < 1/kInfinity && next->r < 1/kInfinity) continue;
282
283 //
284 // We must decide here if we can dare declare one of our faces
285 // as having a "valid" normal (i.e. allBehind = true). This
286 // is never possible if the face faces "inward" in r.
287 //
288 G4bool allBehind;
289 if (corner->z > next->z)
290 {
291 allBehind = false;
292 }
293 else
294 {
295 //
296 // Otherwise, it is only true if the line passing
297 // through the two points of the segment do not
298 // split the r/z cross section
299 //
300 allBehind = !rz->BisectedBy( corner->r, corner->z,
301 next->r, next->z, kCarTolerance );
302 }
303
304 *face++ = new G4PolyconeSide( prev, corner, next, nextNext,
305 startPhi, endPhi-startPhi, phiIsOpen, allBehind );
306 } while( prev=corner, corner=next, corner > corners );
307
308 if (phiIsOpen)
309 {
310 //
311 // Construct phi open edges
312 //
313 *face++ = new G4PolyPhiFace( rz, startPhi, 0, endPhi );
314 *face++ = new G4PolyPhiFace( rz, endPhi, 0, startPhi );
315 }
316
317 //
318 // We might have dropped a face or two: recalculate numFace
319 //
320 numFace = (G4int)(face-faces);
321
322 //
323 // Make enclosingCylinder
324 //
326 new G4EnclosingCylinder( rz, phiIsOpen, phiStart, phiTotal );
327}
G4bool BisectedBy(G4double a1, G4double b1, G4double a2, G4double b2, G4double tolerance)
G4double Amin() const
G4bool RemoveDuplicateVertices(G4double tolerance)
G4bool RemoveRedundantVertices(G4double tolerance)
G4bool CrossesItself(G4double tolerance)
G4VCSGface ** faces
#define DBL_EPSILON
Definition templates.hh:66

Referenced by G4Polycone(), G4Polycone(), and Reset().

◆ CreatePolyhedron()

G4Polyhedron * G4Polycone::CreatePolyhedron ( ) const
overridevirtual

Implements G4VCSGfaceted.

Definition at line 949 of file G4Polycone.cc.

950{
951 std::vector<G4TwoVector> rz(numCorner);
952 for (G4int i = 0; i < numCorner; ++i)
953 rz[i].set(corners[i].r, corners[i].z);
954 return new G4PolyhedronPcon(startPhi, endPhi - startPhi, rz);
955}

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 502 of file G4Polycone.cc.

503{
505}
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override

◆ DistanceToIn() [2/2]

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

Implements G4VSolid.

Definition at line 485 of file G4Polycone.cc.

487{
488 //
489 // Quick test
490 //
492 return kInfinity;
493
494 //
495 // Long answer
496 //
497 return G4VCSGfaceted::DistanceToIn( p, v );
498}
G4bool ShouldMiss(const G4ThreeVector &p, const G4ThreeVector &v) const

◆ GetCorner()

◆ GetCosEndPhi()

G4double G4Polycone::GetCosEndPhi ( ) const
inline

Referenced by BoundingLimits(), and CalculateExtent().

◆ GetCosStartPhi()

G4double G4Polycone::GetCosStartPhi ( ) const
inline

Referenced by BoundingLimits(), and CalculateExtent().

◆ GetCubicVolume()

G4double G4Polycone::GetCubicVolume ( )
overridevirtual

Reimplemented from G4VSolid.

Definition at line 760 of file G4Polycone.cc.

761{
762 if (fCubicVolume == 0.)
763 {
764 G4double total = 0.;
765 G4int nrz = GetNumRZCorner();
766 G4PolyconeSideRZ a = GetCorner(nrz - 1);
767 for (G4int i=0; i<nrz; ++i)
768 {
770 total += (b.r*b.r + b.r*a.r + a.r*a.r)*(b.z - a.z);
771 a = b;
772 }
773 fCubicVolume = std::abs(total)*(GetEndPhi() - GetStartPhi())/6.;
774 }
775 return fCubicVolume;
776}
G4double fCubicVolume
G4double total(Particle const *const p1, Particle const *const p2)

◆ GetEndPhi()

◆ GetEntityType()

G4GeometryType G4Polycone::GetEntityType ( ) const
overridevirtual

Implements G4VSolid.

Definition at line 695 of file G4Polycone.cc.

696{
697 return {"G4Polycone"};
698}

◆ GetNumRZCorner()

◆ GetOriginalParameters()

◆ GetPointOnSurface()

G4ThreeVector G4Polycone::GetPointOnSurface ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 884 of file G4Polycone.cc.

885{
886 // Set surface elements
887 if (fElements == nullptr)
888 {
889 G4AutoLock l(&surface_elementsMutex);
891 l.unlock();
892 }
893
894 // Select surface element
896 selem = fElements->back();
897 G4double select = selem.area*G4QuickRand();
898 auto it = std::lower_bound(fElements->begin(), fElements->end(), select,
899 [](const G4Polycone::surface_element& x, G4double val)
900 -> G4bool { return x.area < val; });
901
902 // Generate random point
903 G4double r = 0, z = 0, phi = 0;
904 G4double u = G4QuickRand();
905 G4double v = G4QuickRand();
906 G4int i0 = (*it).i0;
907 G4int i1 = (*it).i1;
908 G4int i2 = (*it).i2;
909 if (i2 < 0) // lateral surface
910 {
913 if (p1.r < p0.r)
914 {
915 p0 = GetCorner(i1);
916 p1 = GetCorner(i0);
917 }
918 if (p1.r - p0.r < kCarTolerance) // cylindrical surface
919 {
920 r = (p1.r - p0.r)*u + p0.r;
921 z = (p1.z - p0.z)*u + p0.z;
922 }
923 else // conical surface
924 {
925 r = std::sqrt(p1.r*p1.r*u + p0.r*p0.r*(1. - u));
926 z = p0.z + (p1.z - p0.z)*(r - p0.r)/(p1.r - p0.r);
927 }
928 phi = (GetEndPhi() - GetStartPhi())*v + GetStartPhi();
929 }
930 else // phi cut
931 {
932 G4int nrz = GetNumRZCorner();
933 phi = (i0 < nrz) ? GetStartPhi() : GetEndPhi();
934 if (i0 >= nrz) { i0 -= nrz; }
938 if (u + v > 1.) { u = 1. - u; v = 1. - v; }
939 r = (p1.r - p0.r)*u + (p2.r - p0.r)*v + p0.r;
940 z = (p1.z - p0.z)*u + (p2.z - p0.z)*v + p0.z;
941 }
942 return { r*std::cos(phi), r*std::sin(phi), z };
943}
G4double G4QuickRand()
void SetSurfaceElements() const

◆ GetSinEndPhi()

G4double G4Polycone::GetSinEndPhi ( ) const
inline

Referenced by BoundingLimits(), and CalculateExtent().

◆ GetSinStartPhi()

G4double G4Polycone::GetSinStartPhi ( ) const
inline

Referenced by BoundingLimits(), and CalculateExtent().

◆ GetStartPhi()

◆ GetSurfaceArea()

G4double G4Polycone::GetSurfaceArea ( )
overridevirtual

Reimplemented from G4VSolid.

Definition at line 782 of file G4Polycone.cc.

783{
784 if (fSurfaceArea == 0.)
785 {
786 // phi cut area
787 G4int nrz = GetNumRZCorner();
788 G4double scut = 0.;
789 if (IsOpen())
790 {
791 G4PolyconeSideRZ a = GetCorner(nrz - 1);
792 for (G4int i=0; i<nrz; ++i)
793 {
795 scut += a.r*b.z - a.z*b.r;
796 a = b;
797 }
798 scut = std::abs(scut);
799 }
800 // lateral surface area
801 G4double slat = 0;
802 G4PolyconeSideRZ a = GetCorner(nrz - 1);
803 for (G4int i=0; i<nrz; ++i)
804 {
806 G4double h = std::sqrt((b.r - a.r)*(b.r - a.r) + (b.z - a.z)*(b.z - a.z));
807 slat += (b.r + a.r)*h;
808 a = b;
809 }
810 slat *= (GetEndPhi() - GetStartPhi())/2.;
811 fSurfaceArea = scut + slat;
812 }
813 return fSurfaceArea;
814}
G4double fSurfaceArea

◆ Inside()

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

Implements G4VSolid.

Definition at line 467 of file G4Polycone.cc.

468{
469 //
470 // Quick test
471 //
473
474 //
475 // Long answer
476 //
477 return G4VCSGfaceted::Inside(p);
478}
G4bool MustBeOutside(const G4ThreeVector &p) const
EInside Inside(const G4ThreeVector &p) const override
@ kOutside
Definition geomdefs.hh:68

◆ IsOpen()

G4bool G4Polycone::IsOpen ( ) const
inline

◆ operator=()

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

Definition at line 365 of file G4Polycone.cc.

366{
367 if (this == &source) return *this;
368
369 G4VCSGfaceted::operator=( source );
370
371 delete [] corners;
372 delete original_parameters;
373
374 delete enclosingCylinder;
375
376 CopyStuff( source );
377
378 return *this;
379}
G4VCSGfaceted & operator=(const G4VCSGfaceted &source)

◆ Reset()

G4bool G4Polycone::Reset ( )

Definition at line 435 of file G4Polycone.cc.

436{
437 //
438 // Clear old setup
439 //
441 delete [] corners;
442 delete enclosingCylinder;
443 delete fElements;
444 corners = nullptr;
445 fElements = nullptr;
446 enclosingCylinder = nullptr;
447
448 //
449 // Rebuild polycone
450 //
457 delete rz;
458
459 return false;
460}

Referenced by G4ParameterisationPolyconePhi::ComputeDimensions(), G4ParameterisationPolyconeRho::ComputeDimensions(), and G4ParameterisationPolyconeZ::ComputeDimensions().

◆ SetOriginalParameters() [1/2]

◆ SetOriginalParameters() [2/2]

G4bool G4Polycone::SetOriginalParameters ( G4ReduciblePolygon * rz)
protected

Definition at line 959 of file G4Polycone.cc.

960{
961 G4int numPlanes = numCorner;
962 G4bool isConvertible = true;
963 G4double Zmax=rz->Bmax();
964 rz->StartWithZMin();
965
966 // Prepare vectors for storage
967 //
968 std::vector<G4double> Z;
969 std::vector<G4double> Rmin;
970 std::vector<G4double> Rmax;
971
972 G4int countPlanes=1;
973 G4int icurr=0;
974 G4int icurl=0;
975
976 // first plane Z=Z[0]
977 //
978 Z.push_back(corners[0].z);
979 G4double Zprev=Z[0];
980 if (Zprev == corners[1].z)
981 {
982 Rmin.push_back(corners[0].r);
983 Rmax.push_back (corners[1].r);icurr=1;
984 }
985 else if (Zprev == corners[numPlanes-1].z)
986 {
987 Rmin.push_back(corners[numPlanes-1].r);
988 Rmax.push_back (corners[0].r);
989 icurl=numPlanes-1;
990 }
991 else
992 {
993 Rmin.push_back(corners[0].r);
994 Rmax.push_back (corners[0].r);
995 }
996
997 // next planes until last
998 //
999 G4int inextr=0, inextl=0;
1000 for (G4int i=0; i < numPlanes-2; ++i)
1001 {
1002 inextr=1+icurr;
1003 inextl=(icurl <= 0)? numPlanes-1 : icurl-1;
1004
1005 if((corners[inextr].z >= Zmax) & (corners[inextl].z >= Zmax)) { break; }
1006
1007 G4double Zleft = corners[inextl].z;
1008 G4double Zright = corners[inextr].z;
1009 if(Zright > Zleft) // Next plane will be Zleft
1010 {
1011 Z.push_back(Zleft);
1012 countPlanes++;
1013 G4double difZr=corners[inextr].z - corners[icurr].z;
1014 G4double difZl=corners[inextl].z - corners[icurl].z;
1015
1016 if(std::fabs(difZl) < kCarTolerance)
1017 {
1018 if(std::fabs(difZr) < kCarTolerance)
1019 {
1020 Rmin.push_back(corners[inextl].r);
1021 Rmax.push_back(corners[icurr].r);
1022 }
1023 else
1024 {
1025 Rmin.push_back(corners[inextl].r);
1026 Rmax.push_back(corners[icurr].r + (Zleft-corners[icurr].z)/difZr
1027 *(corners[inextr].r - corners[icurr].r));
1028 }
1029 }
1030 else if (difZl >= kCarTolerance)
1031 {
1032 if(std::fabs(difZr) < kCarTolerance)
1033 {
1034 Rmin.push_back(corners[icurl].r);
1035 Rmax.push_back(corners[icurr].r);
1036 }
1037 else
1038 {
1039 Rmin.push_back(corners[icurl].r);
1040 Rmax.push_back(corners[icurr].r + (Zleft-corners[icurr].z)/difZr
1041 *(corners[inextr].r - corners[icurr].r));
1042 }
1043 }
1044 else
1045 {
1046 isConvertible=false; break;
1047 }
1048 icurl=(icurl == 0)? numPlanes-1 : icurl-1;
1049 }
1050 else if(std::fabs(Zright-Zleft)<kCarTolerance) // Zright=Zleft
1051 {
1052 Z.push_back(Zleft);
1053 ++countPlanes;
1054 ++icurr;
1055
1056 icurl=(icurl == 0)? numPlanes-1 : icurl-1;
1057
1058 Rmin.push_back(corners[inextl].r);
1059 Rmax.push_back(corners[inextr].r);
1060 }
1061 else // Zright<Zleft
1062 {
1063 Z.push_back(Zright);
1064 ++countPlanes;
1065
1066 G4double difZr=corners[inextr].z - corners[icurr].z;
1067 G4double difZl=corners[inextl].z - corners[icurl].z;
1068 if(std::fabs(difZr) < kCarTolerance)
1069 {
1070 if(std::fabs(difZl) < kCarTolerance)
1071 {
1072 Rmax.push_back(corners[inextr].r);
1073 Rmin.push_back(corners[icurr].r);
1074 }
1075 else
1076 {
1077 Rmin.push_back(corners[icurl].r + (Zright-corners[icurl].z)/difZl
1078 *(corners[inextl].r - corners[icurl].r));
1079 Rmax.push_back(corners[inextr].r);
1080 }
1081 ++icurr;
1082 } // plate
1083 else if (difZr >= kCarTolerance)
1084 {
1085 if(std::fabs(difZl) < kCarTolerance)
1086 {
1087 Rmax.push_back(corners[inextr].r);
1088 Rmin.push_back (corners[icurr].r);
1089 }
1090 else
1091 {
1092 Rmax.push_back(corners[inextr].r);
1093 Rmin.push_back (corners[icurl].r+(Zright-corners[icurl].z)/difZl
1094 * (corners[inextl].r - corners[icurl].r));
1095 }
1096 ++icurr;
1097 }
1098 else
1099 {
1100 isConvertible=false; break;
1101 }
1102 }
1103 } // end for loop
1104
1105 // last plane Z=Zmax
1106 //
1107 Z.push_back(Zmax);
1108 ++countPlanes;
1109 inextr=1+icurr;
1110 inextl=(icurl <= 0)? numPlanes-1 : icurl-1;
1111
1112 Rmax.push_back(corners[inextr].r);
1113 Rmin.push_back(corners[inextl].r);
1114
1115 // Set original parameters Rmin,Rmax,Z
1116 //
1117 if(isConvertible)
1118 {
1120 original_parameters->Z_values = new G4double[countPlanes];
1121 original_parameters->Rmin = new G4double[countPlanes];
1122 original_parameters->Rmax = new G4double[countPlanes];
1123
1124 for(G4int j=0; j < countPlanes; ++j)
1125 {
1126 original_parameters->Z_values[j] = Z[j];
1127 original_parameters->Rmax[j] = Rmax[j];
1128 original_parameters->Rmin[j] = Rmin[j];
1129 }
1132 original_parameters->Num_z_planes = countPlanes;
1133
1134 }
1135 else // Set parameters(r,z) with Rmin==0 as convention
1136 {
1137#ifdef G4SPECSDEBUG
1138 std::ostringstream message;
1139 message << "Polycone " << GetName() << G4endl
1140 << "cannot be converted to Polycone with (Rmin,Rmaz,Z) parameters!";
1141 G4Exception("G4Polycone::SetOriginalParameters()", "GeomSolids0002",
1142 JustWarning, message);
1143#endif
1145 original_parameters->Z_values = new G4double[numPlanes];
1146 original_parameters->Rmin = new G4double[numPlanes];
1147 original_parameters->Rmax = new G4double[numPlanes];
1148
1149 for(G4int j=0; j < numPlanes; ++j)
1150 {
1153 original_parameters->Rmin[j] = 0.0;
1154 }
1157 original_parameters->Num_z_planes = numPlanes;
1158 }
1159 return isConvertible;
1160}
G4double Bmax() const

◆ SetSurfaceElements()

void G4Polycone::SetSurfaceElements ( ) const
protected

Definition at line 821 of file G4Polycone.cc.

822{
823 fElements = new std::vector<G4Polycone::surface_element>;
824 G4double total = 0.;
825 G4int nrz = GetNumRZCorner();
826
827 // set lateral surface elements
828 G4double dphi = GetEndPhi() - GetStartPhi();
829 G4int ia = nrz - 1;
830 for (G4int ib=0; ib<nrz; ++ib)
831 {
835 selem.i0 = ia;
836 selem.i1 = ib;
837 selem.i2 = -1;
838 ia = ib;
839 if (a.r == 0. && b.r == 0.) continue;
840 G4double h = std::sqrt((b.r - a.r)*(b.r - a.r) + (b.z - a.z)*(b.z - a.z));
841 total += 0.5*dphi*(b.r + a.r)*h;
842 selem.area = total;
843 fElements->push_back(selem);
844 }
845
846 // set elements for phi cuts
847 if (IsOpen())
848 {
849 G4TwoVectorList contourRZ;
850 std::vector<G4int> triangles;
851 for (G4int i=0; i<nrz; ++i)
852 {
853 G4PolyconeSideRZ corner = GetCorner(i);
854 contourRZ.emplace_back(corner.r, corner.z);
855 }
856 G4GeomTools::TriangulatePolygon(contourRZ, triangles);
857 auto ntria = (G4int)triangles.size();
858 for (G4int i=0; i<ntria; i+=3)
859 {
861 selem.i0 = triangles[i];
862 selem.i1 = triangles[i+1];
863 selem.i2 = triangles[i+2];
864 G4PolyconeSideRZ a = GetCorner(selem.i0);
865 G4PolyconeSideRZ b = GetCorner(selem.i1);
866 G4PolyconeSideRZ c = GetCorner(selem.i2);
867 G4double stria =
868 std::abs(G4GeomTools::TriangleArea(a.r, a.z, b.r, b.z, c.r, c.z));
869 total += stria;
870 selem.area = total;
871 fElements->push_back(selem); // start phi
872 total += stria;
873 selem.area = total;
874 selem.i0 += nrz;
875 fElements->push_back(selem); // end phi
876 }
877 }
878}
static G4double TriangleArea(G4double Ax, G4double Ay, G4double Bx, G4double By, G4double Cx, G4double Cy)

Referenced by GetPointOnSurface().

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 710 of file G4Polycone.cc.

711{
712 G4long oldprc = os.precision(16);
713 os << "-----------------------------------------------------------\n"
714 << " *** Dump for solid - " << GetName() << " ***\n"
715 << " ===================================================\n"
716 << " Solid type: G4Polycone\n"
717 << " Parameters: \n"
718 << " starting phi angle : " << startPhi/degree << " degrees \n"
719 << " ending phi angle : " << endPhi/degree << " degrees \n";
720 G4int i=0;
721
723 os << " number of Z planes: " << numPlanes << "\n"
724 << " Z values: \n";
725 for (i=0; i<numPlanes; ++i)
726 {
727 os << " Z plane " << i << ": "
728 << original_parameters->Z_values[i] << "\n";
729 }
730 os << " Tangent distances to inner surface (Rmin): \n";
731 for (i=0; i<numPlanes; ++i)
732 {
733 os << " Z plane " << i << ": "
734 << original_parameters->Rmin[i] << "\n";
735 }
736 os << " Tangent distances to outer surface (Rmax): \n";
737 for (i=0; i<numPlanes; ++i)
738 {
739 os << " Z plane " << i << ": "
740 << original_parameters->Rmax[i] << "\n";
741 }
742
743 os << " number of RZ points: " << numCorner << "\n"
744 << " RZ values (corners): \n";
745 for (i=0; i<numCorner; ++i)
746 {
747 os << " "
748 << corners[i].r << ", " << corners[i].z << "\n";
749 }
750 os << "-----------------------------------------------------------\n";
751 os.precision(oldprc);
752
753 return os;
754}
long G4long
Definition G4Types.hh:87

Member Data Documentation

◆ corners

G4PolyconeSideRZ* G4Polycone::corners = nullptr
protected

◆ enclosingCylinder

G4EnclosingCylinder* G4Polycone::enclosingCylinder = nullptr
protected

Definition at line 176 of file G4Polycone.hh.

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

◆ endPhi

G4double G4Polycone::endPhi
protected

Definition at line 170 of file G4Polycone.hh.

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

◆ fElements

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

Definition at line 179 of file G4Polycone.hh.

Referenced by CopyStuff(), GetPointOnSurface(), Reset(), SetSurfaceElements(), and ~G4Polycone().

◆ numCorner

G4int G4Polycone::numCorner
protected

Definition at line 172 of file G4Polycone.hh.

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

◆ original_parameters

G4PolyconeHistorical* G4Polycone::original_parameters = nullptr
protected

◆ phiIsOpen

G4bool G4Polycone::phiIsOpen = false
protected

Definition at line 171 of file G4Polycone.hh.

Referenced by CopyStuff(), and Create().

◆ startPhi

G4double G4Polycone::startPhi
protected

Definition at line 169 of file G4Polycone.hh.

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


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