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

G4Polyhedra represents a composed closed polyhedra (PGON) made of planar sizes along the Z axis, with or without cut in Phi. More...

#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[])
 ~G4Polyhedra () override
G4int GetNumSide () const
G4double GetStartPhi () const
G4double GetEndPhi () const
G4double GetSinStartPhi () const
G4double GetCosStartPhi () const
G4double GetSinEndPhi () const
G4double GetCosEndPhi () const
G4bool IsOpen () const
G4bool IsGeneric () const
G4int GetNumRZCorner () const
G4PolyhedraSideRZ GetCorner (const G4int index) const
G4PolyhedraHistoricalGetOriginalParameters () const
void SetOriginalParameters (G4PolyhedraHistorical *pars)
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
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
G4GeometryType GetEntityType () const override
G4bool IsFaceted () const override
G4VSolidClone () const override
G4double GetCubicVolume () override
G4double GetSurfaceArea () override
G4ThreeVector GetPointOnSurface () const override
std::ostream & StreamInfo (std::ostream &os) const override
G4PolyhedronCreatePolyhedron () const override
G4bool Reset ()
 G4Polyhedra (__void__ &)
 G4Polyhedra (const G4Polyhedra &source)
G4Polyhedraoperator= (const G4Polyhedra &source)
Public Member Functions inherited from G4VCSGfaceted
 G4VCSGfaceted (const G4String &name)
 ~G4VCSGfaceted () override
 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 ()
 G4VSolid (const G4VSolid &rhs)
G4VSolidoperator= (const G4VSolid &rhs)
G4bool operator== (const G4VSolid &s) const
G4String GetName () const
void SetName (const G4String &name)
G4double GetTolerance () const
virtual G4int GetNumOfConstituents () const
void DumpInfo () const
virtual const G4VSolidGetConstituentSolid (G4int no) const
virtual G4VSolidGetConstituentSolid (G4int no)
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 G4VSolid (__void__ &)
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
G4double EstimateSurfaceArea (G4int nStat, G4double epsilon) const

Additional Inherited Members

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 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

G4Polyhedra represents a composed closed polyhedra (PGON) made of planar sizes along the Z axis, with or without cut in Phi.

Definition at line 79 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[] )

Constructs a polyhedra, given its parameters.

Parameters
[in]nameThe solid name.
[in]phiStartInitial Phi starting angle.
[in]phiTotalTotal Phi angle.
[in]numSideNumber of sides.
[in]numZPlanesNumber of Z planes.
[in]zPlanePosition of Z planes.
[in]rInnerTangent distance to inner surface.
[in]rOuterTangent distance to outer surface.

Definition at line 72 of file G4Polyhedra.cc.

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

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

◆ G4Polyhedra() [2/4]

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

Alternative constructor of a polyhedra, given corners coordinates.

Parameters
[in]nameThe solid name.
[in]phiStartInitial Phi starting angle.
[in]phiTotalTotal Phi angle.
[in]numSideNumber of sides.
[in]numRZNumber of corners in r,Z space.
[in]rr coordinates of corners.
[in]zZ coordinates of corners.

Definition at line 154 of file G4Polyhedra.cc.

161 : G4VCSGfaceted( name ), genericPgon(true)
162{
163 if (theNumSide <= 0)
164 {
165 std::ostringstream message;
166 message << "Solid must have at least one side - " << GetName() << G4endl
167 << " No sides specified !";
168 G4Exception("G4Polyhedra::G4Polyhedra()", "GeomSolids0002",
169 FatalErrorInArgument, message);
170 }
171
172 auto rz = new G4ReduciblePolygon( r, z, numRZ );
173
174 Create( phiStart, phiTotal, theNumSide, rz );
175
176 // Set original_parameters struct for consistency
177 //
179
180 delete rz;
181}
void SetOriginalParameters(G4PolyhedraHistorical *pars)

◆ ~G4Polyhedra()

G4Polyhedra::~G4Polyhedra ( )
override

Destructor.

Definition at line 347 of file G4Polyhedra.cc.

348{
349 delete [] corners;
350 delete original_parameters;
351 delete enclosingCylinder;
352 delete fElements;
353 delete fpPolyhedron;
354 corners = nullptr;
355 original_parameters = nullptr;
356 enclosingCylinder = nullptr;
357 fElements = nullptr;
358 fpPolyhedron = nullptr;
359}
G4Polyhedron * fpPolyhedron

◆ G4Polyhedra() [3/4]

G4Polyhedra::G4Polyhedra ( __void__ & a)

Fake default constructor for usage restricted to direct object persistency for clients requiring preallocation of memory for persistifiable objects.

Definition at line 340 of file G4Polyhedra.cc.

341 : G4VCSGfaceted(a), startPhi(0.), endPhi(0.)
342{
343}

◆ G4Polyhedra() [4/4]

G4Polyhedra::G4Polyhedra ( const G4Polyhedra & source)

Copy constructor and assignment operator.

Definition at line 363 of file G4Polyhedra.cc.

364 : G4VCSGfaceted( source )
365{
366 CopyStuff( source );
367}

Member Function Documentation

◆ BoundingLimits()

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

Computes the bounding limits of the solid.

Parameters
[out]pMinThe minimum bounding limit point.
[out]pMaxThe maximum bounding limit point.

Reimplemented from G4VSolid.

Definition at line 529 of file G4Polyhedra.cc.

531{
532 G4double rmin = kInfinity, rmax = -kInfinity;
533 G4double zmin = kInfinity, zmax = -kInfinity;
534 for (G4int i=0; i<GetNumRZCorner(); ++i)
535 {
536 G4PolyhedraSideRZ corner = GetCorner(i);
537 if (corner.r < rmin) { rmin = corner.r; }
538 if (corner.r > rmax) { rmax = corner.r; }
539 if (corner.z < zmin) { zmin = corner.z; }
540 if (corner.z > zmax) { zmax = corner.z; }
541 }
542
543 G4double sphi = GetStartPhi();
544 G4double ephi = GetEndPhi();
545 G4double dphi = IsOpen() ? ephi-sphi : twopi;
546 G4int ksteps = GetNumSide();
547 G4double astep = dphi/ksteps;
548 G4double sinStep = std::sin(astep);
549 G4double cosStep = std::cos(astep);
550
551 G4double sinCur = GetSinStartPhi();
552 G4double cosCur = GetCosStartPhi();
553 if (!IsOpen()) { rmin = 0.; }
554 G4double xmin = rmin*cosCur, xmax = xmin;
555 G4double ymin = rmin*sinCur, ymax = ymin;
556 for (G4int k=0; k<ksteps+1; ++k)
557 {
558 G4double x = rmax*cosCur;
559 if (x < xmin) { xmin = x; }
560 if (x > xmax) { xmax = x; }
561 G4double y = rmax*sinCur;
562 if (y < ymin) { ymin = y; }
563 if (y > ymax) { ymax = y; }
564 if (rmin > 0)
565 {
566 G4double xx = rmin*cosCur;
567 if (xx < xmin) { xmin = xx; }
568 if (xx > xmax) { xmax = xx; }
569 G4double yy = rmin*sinCur;
570 if (yy < ymin) { ymin = yy; }
571 if (yy > ymax) { ymax = yy; }
572 }
573 G4double sinTmp = sinCur;
574 sinCur = sinCur*cosStep + cosCur*sinStep;
575 cosCur = cosCur*cosStep - sinTmp*sinStep;
576 }
577 pMin.set(xmin,ymin,zmin);
578 pMax.set(xmax,ymax,zmax);
579
580 // Check correctness of the bounding box
581 //
582 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
583 {
584 std::ostringstream message;
585 message << "Bad bounding box (min >= max) for solid: "
586 << GetName() << " !"
587 << "\npMin = " << pMin
588 << "\npMax = " << pMax;
589 G4Exception("G4Polyhedra::BoundingLimits()", "GeomMgt0001",
590 JustWarning, message);
591 DumpInfo();
592 }
593}
@ JustWarning
double z() const
double x() const
double y() const
void set(double x, double y, double z)
G4int GetNumRZCorner() const
G4double GetSinStartPhi() const
G4bool IsOpen() const
G4double GetEndPhi() const
G4int GetNumSide() const
G4PolyhedraSideRZ GetCorner(const G4int index) const
G4double GetStartPhi() const
G4double GetCosStartPhi() const

Referenced by CalculateExtent().

◆ CalculateExtent()

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

Calculates the minimum and maximum extent of the solid, when under the specified transform, and within the specified limits.

Parameters
[in]pAxisThe axis along which compute the extent.
[in]pVoxelLimitThe limiting space dictated by voxels.
[in]pTransformThe internal transformation applied to the solid.
[out]pMinThe minimum extent value.
[out]pMaxThe maximum extent value.
Returns
True if the solid is intersected by the extent region.

Implements G4VSolid.

Definition at line 597 of file G4Polyhedra.cc.

601{
602 G4ThreeVector bmin, bmax;
603 G4bool exist;
604
605 // Check bounding box (bbox)
606 //
607 BoundingLimits(bmin,bmax);
608 G4BoundingEnvelope bbox(bmin,bmax);
609#ifdef G4BBOX_EXTENT
610 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
611#endif
612 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
613 {
614 return exist = pMin < pMax;
615 }
616
617 // To find the extent, RZ contour of the polycone is subdivided
618 // in triangles. The extent is calculated as cumulative extent of
619 // all sub-polycones formed by rotation of triangles around Z
620 //
621 G4TwoVectorList contourRZ;
622 G4TwoVectorList triangles;
623 std::vector<G4int> iout;
624 G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
625 G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
626
627 // get RZ contour, ensure anticlockwise order of corners
628 for (G4int i=0; i<GetNumRZCorner(); ++i)
629 {
630 G4PolyhedraSideRZ corner = GetCorner(i);
631 contourRZ.emplace_back(corner.r,corner.z);
632 }
634 G4double area = G4GeomTools::PolygonArea(contourRZ);
635 if (area < 0.) { std::reverse(contourRZ.begin(),contourRZ.end()); }
636
637 // triangulate RZ countour
638 if (!G4GeomTools::TriangulatePolygon(contourRZ,triangles))
639 {
640 std::ostringstream message;
641 message << "Triangulation of RZ contour has failed for solid: "
642 << GetName() << " !"
643 << "\nExtent has been calculated using boundary box";
644 G4Exception("G4Polyhedra::CalculateExtent()",
645 "GeomMgt1002",JustWarning,message);
646 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
647 }
648
649 // set trigonometric values
650 G4double sphi = GetStartPhi();
651 G4double ephi = GetEndPhi();
652 G4double dphi = IsOpen() ? ephi-sphi : twopi;
653 G4int ksteps = GetNumSide();
654 G4double astep = dphi/ksteps;
655 G4double sinStep = std::sin(astep);
656 G4double cosStep = std::cos(astep);
657 G4double sinStart = GetSinStartPhi();
658 G4double cosStart = GetCosStartPhi();
659
660 // allocate vector lists
661 std::vector<const G4ThreeVectorList *> polygons;
662 polygons.resize(ksteps+1);
663 for (G4int k=0; k<ksteps+1; ++k)
664 {
665 polygons[k] = new G4ThreeVectorList(3);
666 }
667
668 // main loop along triangles
669 pMin = kInfinity;
670 pMax = -kInfinity;
671 G4int ntria = (G4int)triangles.size()/3;
672 for (G4int i=0; i<ntria; ++i)
673 {
674 G4double sinCur = sinStart;
675 G4double cosCur = cosStart;
676 G4int i3 = i*3;
677 for (G4int k=0; k<ksteps+1; ++k) // rotate triangle
678 {
679 auto ptr = const_cast<G4ThreeVectorList*>(polygons[k]);
680 auto iter = ptr->begin();
681 iter->set(triangles[i3+0].x()*cosCur,
682 triangles[i3+0].x()*sinCur,
683 triangles[i3+0].y());
684 iter++;
685 iter->set(triangles[i3+1].x()*cosCur,
686 triangles[i3+1].x()*sinCur,
687 triangles[i3+1].y());
688 iter++;
689 iter->set(triangles[i3+2].x()*cosCur,
690 triangles[i3+2].x()*sinCur,
691 triangles[i3+2].y());
692
693 G4double sinTmp = sinCur;
694 sinCur = sinCur*cosStep + cosCur*sinStep;
695 cosCur = cosCur*cosStep - sinTmp*sinStep;
696 }
697
698 // set sub-envelope and adjust extent
699 G4double emin,emax;
700 G4BoundingEnvelope benv(polygons);
701 if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax))
702 {
703 continue;
704 }
705 if (emin < pMin) { pMin = emin; }
706 if (emax > pMax) { pMax = emax; }
707 if (eminlim > pMin && emaxlim < pMax) { break; } // max possible extent
708 }
709 // free memory
710 for (G4int k=0; k<ksteps+1; ++k)
711 {
712 delete polygons[k];
713 polygons[k]=nullptr;
714 }
715 return (pMin < pMax);
716}
std::vector< G4ThreeVector > G4ThreeVectorList
std::vector< G4TwoVector > G4TwoVectorList
CLHEP::Hep3Vector G4ThreeVector
bool G4bool
Definition G4Types.hh:86
static void RemoveRedundantVertices(G4TwoVectorList &polygon, std::vector< G4int > &iout, G4double tolerance=0.0)
static G4double PolygonArea(const G4TwoVectorList &polygon)
static G4bool TriangulatePolygon(const G4TwoVectorList &polygon, std::vector< G4int > &result)
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
G4double kCarTolerance
Definition G4VSolid.hh:418
G4double GetMinExtent(const EAxis pAxis) const
G4double GetMaxExtent(const EAxis pAxis) const

◆ Clone()

G4VSolid * G4Polyhedra::Clone ( ) const
overridevirtual

Makes a clone of the object for use in multi-treading.

Returns
A pointer to the new cloned allocated solid.

Reimplemented from G4VSolid.

Definition at line 743 of file G4Polyhedra.cc.

744{
745 return new G4Polyhedra(*this);
746}
G4Polyhedra(const G4String &name, G4double phiStart, G4double phiTotal, G4int numSide, G4int numZPlanes, const G4double zPlane[], const G4double rInner[], const G4double rOuter[])

◆ ComputeDimensions()

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

Dispatch method for parameterisation replication mechanism and dimension computation.

Reimplemented from G4VSolid.

Definition at line 720 of file G4Polyhedra.cc.

723{
724 p->ComputeDimensions(*this,n,pRep);
725}
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const

◆ CreatePolyhedron()

G4Polyhedron * G4Polyhedra::CreatePolyhedron ( ) const
overridevirtual

Returns a pointer to a generated polyhedron used for visualisation.

Implements G4VCSGfaceted.

Definition at line 1029 of file G4Polyhedra.cc.

1030{
1031 std::vector<G4TwoVector> rz(numCorner);
1032 for (G4int i = 0; i < numCorner; ++i)
1033 {
1034 rz[i].set(corners[i].r, corners[i].z);
1035 }
1036
1037 // Check the validity of the delta phi
1038 G4double wrDelta = endPhi - startPhi;
1039 if (wrDelta <= 0. || wrDelta >= twopi*(1-DBL_EPSILON))
1040 {
1041 wrDelta = twopi;
1042 }
1043 return new G4PolyhedronPgon(startPhi, wrDelta, numSide, rz);
1044}

◆ DistanceToIn() [1/2]

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

Calculates the distance to the nearest surface of a shape from an outside point. The distance can be an underestimate.

Parameters
[in]pThe point at offset p.
Returns
The safety distance to enter the shape.

Implements G4VSolid.

Definition at line 522 of file G4Polyhedra.cc.

523{
525}
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override

◆ DistanceToIn() [2/2]

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

Returns the distance along the normalised vector 'v' to the shape, from the point at offset 'p'. If there is no intersection, returns kInfinity. The first intersection resulting from 'leaving' a surface/volume is discarded. Hence, it is tolerant of points on the surface of the shape.

Parameters
[in]pThe point at offset p.
[in]vThe normalised direction vector.
Returns
The distance to enter the shape.

Implements G4VSolid.

Definition at line 506 of file G4Polyhedra.cc.

508{
509 //
510 // Quick test
511 //
512 if (enclosingCylinder->ShouldMiss(p,v)) { return kInfinity; }
513
514 //
515 // Long answer
516 //
517 return G4VCSGfaceted::DistanceToIn( p, v );
518}

◆ GetCorner()

◆ GetCosEndPhi()

G4double G4Polyhedra::GetCosEndPhi ( ) const
inline

◆ GetCosStartPhi()

G4double G4Polyhedra::GetCosStartPhi ( ) const
inline

Referenced by BoundingLimits(), and CalculateExtent().

◆ GetCubicVolume()

G4double G4Polyhedra::GetCubicVolume ( )
overridevirtual

Returning an estimation of the solid volume (capacity) and surface area, in internal units.

Reimplemented from G4VSolid.

Definition at line 802 of file G4Polyhedra.cc.

803{
804 if (fCubicVolume == 0)
805 {
806 G4AutoLock l(&polyhedraMutex);
807 G4double total = 0.;
808 G4int nrz = GetNumRZCorner();
809 G4PolyhedraSideRZ a = GetCorner(nrz - 1);
810 for (G4int i=0; i<nrz; ++i)
811 {
812 G4PolyhedraSideRZ b = GetCorner(i);
813 total += (b.r*b.r + b.r*a.r + a.r*a.r)*(b.z - a.z);
814 a = b;
815 }
816 fCubicVolume = std::abs(total)*
817 std::sin((GetEndPhi() - GetStartPhi())/GetNumSide())*GetNumSide()/6.;
818 l.unlock();
819 }
820 return fCubicVolume;
821}
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4double fCubicVolume
G4double total(Particle const *const p1, Particle const *const p2)

◆ GetEndPhi()

◆ GetEntityType()

G4GeometryType G4Polyhedra::GetEntityType ( ) const
overridevirtual

Returns the type ID, "G4Polyhedra" of the solid.

Implements G4VSolid.

Definition at line 729 of file G4Polyhedra.cc.

730{
731 return {"G4Polyhedra"};
732}

◆ GetNumRZCorner()

◆ GetNumSide()

◆ GetOriginalParameters()

G4PolyhedraHistorical * G4Polyhedra::GetOriginalParameters ( ) const
inline

◆ GetPointOnSurface()

G4ThreeVector G4Polyhedra::GetPointOnSurface ( ) const
overridevirtual

Returns a random point located and uniformly distributed on the surface of the solid.

Reimplemented from G4VSolid.

Definition at line 948 of file G4Polyhedra.cc.

949{
950 // Set surface elements
951 if (fElements == nullptr)
952 {
953 G4AutoLock l(&surface_elementsMutex);
954 SetSurfaceElements();
955 l.unlock();
956 }
957
958 // Select surface element
959 G4Polyhedra::surface_element selem;
960 selem = fElements->back();
961 G4double select = selem.area*G4QuickRand();
962 auto it = std::lower_bound(fElements->begin(), fElements->end(), select,
963 [](const G4Polyhedra::surface_element& x, G4double val)
964 -> G4bool { return x.area < val; });
965
966 // Generate random point
967 G4double x = 0, y = 0, z = 0;
968 G4double u = G4QuickRand();
969 G4double v = G4QuickRand();
970 if (u + v > 1.) { u = 1. - u; v = 1. - v; }
971 G4int i0 = (*it).i0;
972 G4int i1 = (*it).i1;
973 G4int i2 = (*it).i2;
974 if (i2 < 0) // lateral surface
975 {
976 // sample point
977 G4int nside = GetNumSide();
978 G4double dphi = (GetEndPhi() - GetStartPhi())/nside;
979 G4double cosa = std::cos(dphi);
980 G4double sina = std::sin(dphi);
981 G4PolyhedraSideRZ a = GetCorner(i0);
982 G4PolyhedraSideRZ b = GetCorner(i1);
983 G4ThreeVector p0(a.r, 0, a.z);
984 G4ThreeVector p1(b.r, 0, b.z);
985 G4ThreeVector p2(b.r*cosa, b.r*sina, b.z);
986 if (i2 == -1) { p1.set(a.r*cosa, a.r*sina, a.z); }
987 p0 += (p1 - p0)*u + (p2 - p0)*v;
988 // find selected side and rotate point
989 G4double scurr = (*it).area;
990 G4double sprev = (it == fElements->begin()) ? 0. : (*(--it)).area;
991 G4int iside = nside*(select - sprev)/(scurr - sprev);
992 if (iside == 0 && GetStartPhi() == 0.)
993 {
994 x = p0.x();
995 y = p0.y();
996 z = p0.z();
997 }
998 else
999 {
1000 if (iside == nside) { --iside; } // iside must be less then nside
1001 G4double phi = iside*dphi + GetStartPhi();
1002 G4double cosphi = std::cos(phi);
1003 G4double sinphi = std::sin(phi);
1004 x = p0.x()*cosphi - p0.y()*sinphi;
1005 y = p0.x()*sinphi + p0.y()*cosphi;
1006 z = p0.z();
1007 }
1008 }
1009 else // phi cut
1010 {
1011 G4int nrz = GetNumRZCorner();
1012 G4double phi = (i0 < nrz) ? GetStartPhi() : GetEndPhi();
1013 if (i0 >= nrz) { i0 -= nrz; }
1014 G4PolyhedraSideRZ p0 = GetCorner(i0);
1015 G4PolyhedraSideRZ p1 = GetCorner(i1);
1016 G4PolyhedraSideRZ p2 = GetCorner(i2);
1017 G4double r = (p1.r - p0.r)*u + (p2.r - p0.r)*v + p0.r;
1018 x = r*std::cos(phi);
1019 y = r*std::sin(phi);
1020 z = (p1.z - p0.z)*u + (p2.z - p0.z)*v + p0.z;
1021 }
1022 return {x, y, z};
1023}
G4double G4QuickRand(uint32_t seed=0)

◆ GetSinEndPhi()

G4double G4Polyhedra::GetSinEndPhi ( ) const
inline

◆ GetSinStartPhi()

G4double G4Polyhedra::GetSinStartPhi ( ) const
inline

Referenced by BoundingLimits(), and CalculateExtent().

◆ GetStartPhi()

◆ GetSurfaceArea()

G4double G4Polyhedra::GetSurfaceArea ( )
overridevirtual

Returns an estimation of the solid surface area in internal units. This method may be overloaded by derived classes to compute the exact geometrical quantity for solids where this is possible, or anyway to cache the computed value. Note: the computed value is NOT cached.

Reimplemented from G4VSolid.

Definition at line 827 of file G4Polyhedra.cc.

828{
829 if (fSurfaceArea == 0)
830 {
831 G4AutoLock l(&polyhedraMutex);
832 G4double total = 0.;
833 G4int nrz = GetNumRZCorner();
834 if (IsOpen())
835 {
836 G4PolyhedraSideRZ a = GetCorner(nrz - 1);
837 for (G4int i=0; i<nrz; ++i)
838 {
839 G4PolyhedraSideRZ b = GetCorner(i);
840 total += a.r*b.z - a.z*b.r;
841 a = b;
842 }
843 total = std::abs(total);
844 }
846 G4double cosa = std::cos(alp);
847 G4double sina = std::sin(alp);
848 G4PolyhedraSideRZ a = GetCorner(nrz - 1);
849 for (G4int i=0; i<nrz; ++i)
850 {
851 G4PolyhedraSideRZ b = GetCorner(i);
852 G4ThreeVector p1(a.r, 0, a.z);
853 G4ThreeVector p2(a.r*cosa, a.r*sina, a.z);
854 G4ThreeVector p3(b.r*cosa, b.r*sina, b.z);
855 G4ThreeVector p4(b.r, 0, b.z);
856 total += GetNumSide()*(G4GeomTools::QuadAreaNormal(p1, p2, p3, p4)).mag();
857 a = b;
858 }
860 l.unlock();
861 }
862 return fSurfaceArea;
863}
static G4ThreeVector QuadAreaNormal(const G4ThreeVector &A, const G4ThreeVector &B, const G4ThreeVector &C, const G4ThreeVector &D)
G4double fSurfaceArea

◆ Inside()

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

Concrete implementations of the expected query interfaces for solids, as defined in G4VSolid. Remaining functions are concretely defined in the base class G4VCSGfaceted.

Implements G4VSolid.

Definition at line 488 of file G4Polyhedra.cc.

489{
490 //
491 // Quick test
492 //
493 if (enclosingCylinder->MustBeOutside(p)) { return kOutside; }
494
495 //
496 // Long answer
497 //
498 return G4VCSGfaceted::Inside(p);
499}
EInside Inside(const G4ThreeVector &p) const override
@ kOutside
Definition geomdefs.hh:68

◆ IsFaceted()

G4bool G4Polyhedra::IsFaceted ( ) const
overridevirtual

Returns true as the solid has only planar faces.

Reimplemented from G4VSolid.

Definition at line 736 of file G4Polyhedra.cc.

737{
738 return true;
739}

◆ IsGeneric()

G4bool G4Polyhedra::IsGeneric ( ) const
inline

◆ IsOpen()

G4bool G4Polyhedra::IsOpen ( ) const
inline

◆ operator=()

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

Definition at line 371 of file G4Polyhedra.cc.

372{
373 if (this == &source) { return *this; }
374
375 G4VCSGfaceted::operator=( source );
376
377 delete [] corners;
378 delete original_parameters;
379 delete enclosingCylinder;
380
381 CopyStuff( source );
382
383 return *this;
384}
G4VCSGfaceted & operator=(const G4VCSGfaceted &source)

◆ Reset()

G4bool G4Polyhedra::Reset ( )

Clears all parameters and rebuild the shape, for use in divisions.

Definition at line 445 of file G4Polyhedra.cc.

446{
447 if (genericPgon)
448 {
449 std::ostringstream message;
450 message << "Solid " << GetName() << " built using generic construct."
451 << G4endl << "Not applicable to the generic construct !";
452 G4Exception("G4Polyhedra::Reset()", "GeomSolids1001",
453 JustWarning, message, "Parameters NOT resetted.");
454 return true;
455 }
456
457 //
458 // Clear old setup
459 //
461 delete [] corners;
462 delete enclosingCylinder;
463 delete fElements;
464 corners = nullptr;
465 fElements = nullptr;
466 enclosingCylinder = nullptr;
467
468 //
469 // Rebuild polyhedra
470 //
471 auto rz = new G4ReduciblePolygon( original_parameters->Rmin,
472 original_parameters->Rmax,
473 original_parameters->Z_values,
474 original_parameters->Num_z_planes );
475 Create( original_parameters->Start_angle,
476 original_parameters->Opening_angle,
477 original_parameters->numSide, rz );
478 delete rz;
479
480 return false;
481}

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

◆ SetOriginalParameters()

void G4Polyhedra::SetOriginalParameters ( G4PolyhedraHistorical * pars)
inline

Sets internal parameters. Parameters 'Rmin' and 'Rmax' in input must be scaled first by a factor computed as 'cos(0.5*phiTotal/theNumSide)', if not already scaled.

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

◆ StreamInfo()

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

Streams the object contents to an output stream.

Implements G4VSolid.

Definition at line 750 of file G4Polyhedra.cc.

751{
752 G4long oldprc = os.precision(16);
753 os << "-----------------------------------------------------------\n"
754 << " *** Dump for solid - " << GetName() << " ***\n"
755 << " ===================================================\n"
756 << " Solid type: G4Polyhedra\n"
757 << " Parameters: \n"
758 << " starting phi angle : " << startPhi/degree << " degrees \n"
759 << " ending phi angle : " << endPhi/degree << " degrees \n"
760 << " number of sides : " << numSide << " \n";
761 G4int i=0;
762 if (!genericPgon)
763 {
764 G4int numPlanes = original_parameters->Num_z_planes;
765 os << " number of Z planes: " << numPlanes << "\n"
766 << " Z values: \n";
767 for (i=0; i<numPlanes; ++i)
768 {
769 os << " Z plane " << i << ": "
770 << original_parameters->Z_values[i] << "\n";
771 }
772 os << " Tangent distances to inner surface (Rmin): \n";
773 for (i=0; i<numPlanes; ++i)
774 {
775 os << " Z plane " << i << ": "
776 << original_parameters->Rmin[i] << "\n";
777 }
778 os << " Tangent distances to outer surface (Rmax): \n";
779 for (i=0; i<numPlanes; ++i)
780 {
781 os << " Z plane " << i << ": "
782 << original_parameters->Rmax[i] << "\n";
783 }
784 }
785 os << " number of RZ points: " << numCorner << "\n"
786 << " RZ values (corners): \n";
787 for (i=0; i<numCorner; ++i)
788 {
789 os << " "
790 << corners[i].r << ", " << corners[i].z << "\n";
791 }
792 os << "-----------------------------------------------------------\n";
793 os.precision(oldprc);
794
795 return os;
796}
long G4long
Definition G4Types.hh:87

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