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

G4Polycone represents a composed closed shape (PCON) made of cones and cylinders, along the Z axis with increasing Z, with or without cut in Phi. More...

#include <G4Polycone.hh>

Inheritance diagram for G4Polycone:

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
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)
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 ()
 G4Polycone (__void__ &)
 G4Polycone (const G4Polycone &source)
G4Polyconeoperator= (const G4Polycone &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
virtual G4bool IsFaceted () 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

G4Polycone represents a composed closed shape (PCON) made of cones and cylinders, along the Z axis with increasing Z, with or without cut in Phi.

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

Constructs a polycone shape, given its parameters.

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

Definition at line 59 of file G4Polycone.cc.

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

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

◆ G4Polycone() [2/4]

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

Alternative constructor of a polycone shape, given corners coordinates.

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

Definition at line 132 of file G4Polycone.cc.

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

◆ ~G4Polycone()

G4Polycone::~G4Polycone ( )
override

Destructor.

Definition at line 344 of file G4Polycone.cc.

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

◆ G4Polycone() [3/4]

G4Polycone::G4Polycone ( __void__ & a)

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

Definition at line 335 of file G4Polycone.cc.

336 : G4VCSGfaceted(a), startPhi(0.), endPhi(0.), numCorner(0)
337{
338}

◆ G4Polycone() [4/4]

G4Polycone::G4Polycone ( const G4Polycone & source)

Copy constructor and assignment operator.

Definition at line 360 of file G4Polycone.cc.

361 : G4VCSGfaceted( source )
362{
363 CopyStuff( source );
364}

Member Function Documentation

◆ BoundingLimits()

void G4Polycone::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 511 of file G4Polycone.cc.

513{
514 G4double rmin = kInfinity, rmax = -kInfinity;
515 G4double zmin = kInfinity, zmax = -kInfinity;
516
517 for (G4int i=0; i<GetNumRZCorner(); ++i)
518 {
519 G4PolyconeSideRZ corner = GetCorner(i);
520 if (corner.r < rmin) { rmin = corner.r; }
521 if (corner.r > rmax) { rmax = corner.r; }
522 if (corner.z < zmin) { zmin = corner.z; }
523 if (corner.z > zmax) { zmax = corner.z; }
524 }
525
526 if (IsOpen())
527 {
528 G4TwoVector vmin,vmax;
529 G4GeomTools::DiskExtent(rmin,rmax,
532 vmin,vmax);
533 pMin.set(vmin.x(),vmin.y(),zmin);
534 pMax.set(vmax.x(),vmax.y(),zmax);
535 }
536 else
537 {
538 pMin.set(-rmax,-rmax, zmin);
539 pMax.set( rmax, rmax, zmax);
540 }
541
542 // Check correctness of the bounding box
543 //
544 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
545 {
546 std::ostringstream message;
547 message << "Bad bounding box (min >= max) for solid: "
548 << GetName() << " !"
549 << "\npMin = " << pMin
550 << "\npMax = " << pMax;
551 G4Exception("G4Polycone::BoundingLimits()", "GeomMgt0001",
552 JustWarning, message);
553 DumpInfo();
554 }
555}
@ JustWarning
CLHEP::Hep2Vector G4TwoVector
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

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 559 of file G4Polycone.cc.

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

◆ Clone()

G4VSolid * G4Polycone::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 722 of file G4Polycone.cc.

723{
724 return new G4Polycone(*this);
725}
G4Polycone(const G4String &name, G4double phiStart, G4double phiTotal, G4int numZPlanes, const G4double zPlane[], const G4double rInner[], const G4double rOuter[])
Definition G4Polycone.cc:59

◆ ComputeDimensions()

void G4Polycone::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 706 of file G4Polycone.cc.

709{
710 p->ComputeDimensions(*this,n,pRep);
711}
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const

◆ CreatePolyhedron()

G4Polyhedron * G4Polycone::CreatePolyhedron ( ) const
overridevirtual

Returns a pointer to a generated polyhedron used for visualisation.

Implements G4VCSGfaceted.

Definition at line 973 of file G4Polycone.cc.

974{
975 std::vector<G4TwoVector> rz(numCorner);
976 for (G4int i = 0; i < numCorner; ++i)
977 {
978 rz[i].set(corners[i].r, corners[i].z);
979 }
980 return new G4PolyhedronPcon(startPhi, endPhi - startPhi, rz);
981}

◆ DistanceToIn() [1/2]

G4double G4Polycone::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 504 of file G4Polycone.cc.

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

◆ DistanceToIn() [2/2]

G4double G4Polycone::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 488 of file G4Polycone.cc.

490{
491 //
492 // Quick test
493 //
494 if (enclosingCylinder->ShouldMiss(p,v)) { return kInfinity; }
495
496 //
497 // Long answer
498 //
499 return G4VCSGfaceted::DistanceToIn( p, v );
500}

◆ 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

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

Reimplemented from G4VSolid.

Definition at line 780 of file G4Polycone.cc.

781{
782 if (fCubicVolume == 0)
783 {
784 G4AutoLock l(&polyconeMutex);
785 G4double total = 0.;
786 G4int nrz = GetNumRZCorner();
787 G4PolyconeSideRZ a = GetCorner(nrz - 1);
788 for (G4int i=0; i<nrz; ++i)
789 {
790 G4PolyconeSideRZ b = GetCorner(i);
791 total += (b.r*b.r + b.r*a.r + a.r*a.r)*(b.z - a.z);
792 a = b;
793 }
794 fCubicVolume = std::abs(total)*(GetEndPhi() - GetStartPhi())/6.;
795 l.unlock();
796 }
797 return fCubicVolume;
798}
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4double fCubicVolume
G4double total(Particle const *const p1, Particle const *const p2)

◆ GetEndPhi()

◆ GetEntityType()

G4GeometryType G4Polycone::GetEntityType ( ) const
overridevirtual

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

Implements G4VSolid.

Definition at line 715 of file G4Polycone.cc.

716{
717 return {"G4Polycone"};
718}

◆ GetNumRZCorner()

◆ GetOriginalParameters()

G4PolyconeHistorical * G4Polycone::GetOriginalParameters ( ) const
inline

◆ GetPointOnSurface()

G4ThreeVector G4Polycone::GetPointOnSurface ( ) const
overridevirtual

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

Reimplemented from G4VSolid.

Definition at line 908 of file G4Polycone.cc.

909{
910 // Set surface elements
911 if (fElements == nullptr)
912 {
913 G4AutoLock l(&surface_elementsMutex);
914 SetSurfaceElements();
915 l.unlock();
916 }
917
918 // Select surface element
919 G4Polycone::surface_element selem;
920 selem = fElements->back();
921 G4double select = selem.area*G4QuickRand();
922 auto it = std::lower_bound(fElements->begin(), fElements->end(), select,
923 [](const G4Polycone::surface_element& x, G4double val)
924 -> G4bool { return x.area < val; });
925
926 // Generate random point
927 G4double r = 0, z = 0, phi = 0;
928 G4double u = G4QuickRand();
929 G4double v = G4QuickRand();
930 G4int i0 = (*it).i0;
931 G4int i1 = (*it).i1;
932 G4int i2 = (*it).i2;
933 if (i2 < 0) // lateral surface
934 {
935 G4PolyconeSideRZ p0 = GetCorner(i0);
936 G4PolyconeSideRZ p1 = GetCorner(i1);
937 if (p1.r < p0.r)
938 {
939 p0 = GetCorner(i1);
940 p1 = GetCorner(i0);
941 }
942 if (p1.r - p0.r < kCarTolerance) // cylindrical surface
943 {
944 r = (p1.r - p0.r)*u + p0.r;
945 z = (p1.z - p0.z)*u + p0.z;
946 }
947 else // conical surface
948 {
949 r = std::sqrt(p1.r*p1.r*u + p0.r*p0.r*(1. - u));
950 z = p0.z + (p1.z - p0.z)*(r - p0.r)/(p1.r - p0.r);
951 }
952 phi = (GetEndPhi() - GetStartPhi())*v + GetStartPhi();
953 }
954 else // phi cut
955 {
956 G4int nrz = GetNumRZCorner();
957 phi = (i0 < nrz) ? GetStartPhi() : GetEndPhi();
958 if (i0 >= nrz) { i0 -= nrz; }
959 G4PolyconeSideRZ p0 = GetCorner(i0);
960 G4PolyconeSideRZ p1 = GetCorner(i1);
961 G4PolyconeSideRZ p2 = GetCorner(i2);
962 if (u + v > 1.) { u = 1. - u; v = 1. - v; }
963 r = (p1.r - p0.r)*u + (p2.r - p0.r)*v + p0.r;
964 z = (p1.z - p0.z)*u + (p2.z - p0.z)*v + p0.z;
965 }
966 return { r*std::cos(phi), r*std::sin(phi), z };
967}
G4double G4QuickRand(uint32_t seed=0)

◆ 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

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 804 of file G4Polycone.cc.

805{
806 if (fSurfaceArea == 0)
807 {
808 G4AutoLock l(&polyconeMutex);
809 // phi cut area
810 G4int nrz = GetNumRZCorner();
811 G4double scut = 0.;
812 if (IsOpen())
813 {
814 G4PolyconeSideRZ a = GetCorner(nrz - 1);
815 for (G4int i=0; i<nrz; ++i)
816 {
817 G4PolyconeSideRZ b = GetCorner(i);
818 scut += a.r*b.z - a.z*b.r;
819 a = b;
820 }
821 scut = std::abs(scut);
822 }
823 // lateral surface area
824 G4double slat = 0;
825 G4PolyconeSideRZ a = GetCorner(nrz - 1);
826 for (G4int i=0; i<nrz; ++i)
827 {
828 G4PolyconeSideRZ b = GetCorner(i);
829 G4double h = std::sqrt((b.r - a.r)*(b.r - a.r) + (b.z - a.z)*(b.z - a.z));
830 slat += (b.r + a.r)*h;
831 a = b;
832 }
833 slat *= (GetEndPhi() - GetStartPhi())/2.;
834 fSurfaceArea = scut + slat;
835 l.unlock();
836 }
837 return fSurfaceArea;
838}
G4double fSurfaceArea

◆ Inside()

EInside G4Polycone::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 470 of file G4Polycone.cc.

471{
472 //
473 // Quick test
474 //
475 if (enclosingCylinder->MustBeOutside(p)) { return kOutside; }
476
477 //
478 // Long answer
479 //
480 return G4VCSGfaceted::Inside(p);
481}
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 368 of file G4Polycone.cc.

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

◆ Reset()

G4bool G4Polycone::Reset ( )

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

Definition at line 438 of file G4Polycone.cc.

439{
440 //
441 // Clear old setup
442 //
444 delete [] corners;
445 delete enclosingCylinder;
446 delete fElements;
447 corners = nullptr;
448 fElements = nullptr;
449 enclosingCylinder = nullptr;
450
451 //
452 // Rebuild polycone
453 //
454 auto rz = new G4ReduciblePolygon( original_parameters->Rmin,
455 original_parameters->Rmax,
456 original_parameters->Z_values,
457 original_parameters->Num_z_planes );
458 Create( original_parameters->Start_angle,
459 original_parameters->Opening_angle, rz );
460 delete rz;
461
462 return false;
463}

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

◆ SetOriginalParameters()

◆ StreamInfo()

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

Streams the object contents to an output stream.

Implements G4VSolid.

Definition at line 730 of file G4Polycone.cc.

731{
732 G4long oldprc = os.precision(16);
733 os << "-----------------------------------------------------------\n"
734 << " *** Dump for solid - " << GetName() << " ***\n"
735 << " ===================================================\n"
736 << " Solid type: G4Polycone\n"
737 << " Parameters: \n"
738 << " starting phi angle : " << startPhi/degree << " degrees \n"
739 << " ending phi angle : " << endPhi/degree << " degrees \n";
740 G4int i=0;
741
742 G4int numPlanes = original_parameters->Num_z_planes;
743 os << " number of Z planes: " << numPlanes << "\n"
744 << " Z values: \n";
745 for (i=0; i<numPlanes; ++i)
746 {
747 os << " Z plane " << i << ": "
748 << original_parameters->Z_values[i] << "\n";
749 }
750 os << " Tangent distances to inner surface (Rmin): \n";
751 for (i=0; i<numPlanes; ++i)
752 {
753 os << " Z plane " << i << ": "
754 << original_parameters->Rmin[i] << "\n";
755 }
756 os << " Tangent distances to outer surface (Rmax): \n";
757 for (i=0; i<numPlanes; ++i)
758 {
759 os << " Z plane " << i << ": "
760 << original_parameters->Rmax[i] << "\n";
761 }
762
763 os << " number of RZ points: " << numCorner << "\n"
764 << " RZ values (corners): \n";
765 for (i=0; i<numCorner; ++i)
766 {
767 os << " "
768 << corners[i].r << ", " << corners[i].z << "\n";
769 }
770 os << "-----------------------------------------------------------\n";
771 os.precision(oldprc);
772
773 return os;
774}
long G4long
Definition G4Types.hh:87

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