Geant4 11.1.1
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[])
 
virtual ~G4Polycone ()
 
EInside Inside (const G4ThreeVector &p) const
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
G4double DistanceToIn (const G4ThreeVector &p) const
 
void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
 
G4double GetCubicVolume ()
 
G4double GetSurfaceArea ()
 
G4ThreeVector GetPointOnSurface () const
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
G4GeometryType GetEntityType () const
 
G4VSolidClone () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4PolyhedronCreatePolyhedron () const
 
G4bool Reset ()
 
G4double GetStartPhi () const
 
G4double GetEndPhi () const
 
G4double GetSinStartPhi () const
 
G4double GetCosStartPhi () const
 
G4double GetSinEndPhi () const
 
G4double GetCosEndPhi () const
 
G4bool IsOpen () const
 
G4int GetNumRZCorner () const
 
G4PolyconeSideRZ GetCorner (G4int index) const
 
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)
 
virtual G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
 
virtual EInside Inside (const G4ThreeVector &p) const
 
virtual G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
virtual G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
virtual G4double DistanceToIn (const G4ThreeVector &p) const
 
virtual G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
 
virtual G4double DistanceToOut (const G4ThreeVector &p) const
 
virtual G4GeometryType GetEntityType () const
 
virtual std::ostream & StreamInfo (std::ostream &os) const
 
virtual G4PolyhedronCreatePolyhedron () const =0
 
virtual void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
virtual G4VisExtent GetExtent () const
 
virtual G4PolyhedronGetPolyhedron () const
 
G4int GetCubVolStatistics () const
 
G4double GetCubVolEpsilon () const
 
void SetCubVolStatistics (G4int st)
 
void SetCubVolEpsilon (G4double ep)
 
G4int GetAreaStatistics () const
 
G4double GetAreaAccuracy () const
 
void SetAreaStatistics (G4int st)
 
void SetAreaAccuracy (G4double ep)
 
virtual G4double GetCubicVolume ()
 
virtual G4double GetSurfaceArea ()
 
 G4VCSGfaceted (__void__ &)
 
- Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
 
virtual ~G4VSolid ()
 
G4bool operator== (const G4VSolid &s) const
 
G4String GetName () const
 
void SetName (const G4String &name)
 
G4double GetTolerance () const
 
virtual void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const
 
virtual G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
 
virtual EInside Inside (const G4ThreeVector &p) const =0
 
virtual G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p) const =0
 
virtual void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
virtual G4double GetCubicVolume ()
 
virtual G4double GetSurfaceArea ()
 
virtual G4GeometryType GetEntityType () const =0
 
virtual G4ThreeVector GetPointOnSurface () const
 
virtual G4VSolidClone () const
 
virtual std::ostream & StreamInfo (std::ostream &os) const =0
 
void DumpInfo () const
 
virtual void DescribeYourselfTo (G4VGraphicsScene &scene) const =0
 
virtual G4VisExtent GetExtent () const
 
virtual G4PolyhedronCreatePolyhedron () const
 
virtual G4PolyhedronGetPolyhedron () const
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
 G4VSolid (__void__ &)
 
 G4VSolid (const G4VSolid &rhs)
 
G4VSolidoperator= (const G4VSolid &rhs)
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 

Protected Member Functions

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 //
120 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)
Definition: G4Exception.cc:59
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
void Create(G4double phiStart, G4double phiTotal, G4ReduciblePolygon *rz)
Definition: G4Polycone.cc:172
G4PolyconeHistorical * original_parameters
Definition: G4Polycone.hh:178
void DumpInfo() const

◆ G4Polycone() [2/4]

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

Definition at line 132 of file G4Polycone.cc.

138 : G4VCSGfaceted( name )
139{
140
141 G4ReduciblePolygon* 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 ( )
virtual

Definition at line 342 of file G4Polycone.cc.

343{
344 delete [] corners;
345 delete original_parameters;
346 delete enclosingCylinder;
347 delete fElements;
348 delete fpPolyhedron;
349 corners = nullptr;
350 original_parameters = nullptr;
351 enclosingCylinder = nullptr;
352 fElements = nullptr;
353 fpPolyhedron = nullptr;
354}
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polycone.hh:180
G4PolyconeSideRZ * corners
Definition: G4Polycone.hh:177
std::vector< surface_element > * fElements
Definition: G4Polycone.hh:183
G4Polyhedron * fpPolyhedron

◆ G4Polycone() [3/4]

G4Polycone::G4Polycone ( __void__ &  a)

Definition at line 333 of file G4Polycone.cc.

334 : G4VCSGfaceted(a), startPhi(0.), endPhi(0.), numCorner(0)
335{
336}
G4double endPhi
Definition: G4Polycone.hh:174
G4int numCorner
Definition: G4Polycone.hh:176
G4double startPhi
Definition: G4Polycone.hh:173

◆ G4Polycone() [4/4]

G4Polycone::G4Polycone ( const G4Polycone source)

Definition at line 358 of file G4Polycone.cc.

359 : G4VCSGfaceted( source )
360{
361 CopyStuff( source );
362}
void CopyStuff(const G4Polycone &source)
Definition: G4Polycone.cc:384

Member Function Documentation

◆ BoundingLimits()

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

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
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)
Definition: G4GeomTools.cc:390
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
virtual

Reimplemented from G4VCSGfaceted.

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) ? true : false;
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.push_back(G4TwoVector(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) pols[k].resize(6);
636 for (G4int k=0; k<ksteps+2; ++k) polygons[k] = &pols[k];
637 G4double r0[6],z0[6]; // contour with original edges of triangle
638 G4double r1[6]; // shifted radii of external edges of triangle
639
640 // main loop along triangles
641 pMin = kInfinity;
642 pMax =-kInfinity;
643 G4int ntria = (G4int)triangles.size()/3;
644 for (G4int i=0; i<ntria; ++i)
645 {
646 G4int i3 = i*3;
647 for (G4int k=0; k<3; ++k)
648 {
649 G4int e0 = i3+k, e1 = (k<2) ? e0+1 : i3;
650 G4int k2 = k*2;
651 // set contour with original edges of triangle
652 r0[k2+0] = triangles[e0].x(); z0[k2+0] = triangles[e0].y();
653 r0[k2+1] = triangles[e1].x(); z0[k2+1] = triangles[e1].y();
654 // set shifted radii
655 r1[k2+0] = r0[k2+0];
656 r1[k2+1] = r0[k2+1];
657 if (z0[k2+1] - z0[k2+0] <= 0) continue;
658 r1[k2+0] /= cosHalf;
659 r1[k2+1] /= cosHalf;
660 }
661
662 // rotate countour, set sequence of 6-sided polygons
663 G4double sinCur = sinStart*cosHalf + cosStart*sinHalf;
664 G4double cosCur = cosStart*cosHalf - sinStart*sinHalf;
665 for (G4int j=0; j<6; ++j) pols[0][j].set(r0[j]*cosStart,r0[j]*sinStart,z0[j]);
666 for (G4int k=1; k<ksteps+1; ++k)
667 {
668 for (G4int j=0; j<6; ++j) pols[k][j].set(r1[j]*cosCur,r1[j]*sinCur,z0[j]);
669 G4double sinTmp = sinCur;
670 sinCur = sinCur*cosStep + cosCur*sinStep;
671 cosCur = cosCur*cosStep - sinTmp*sinStep;
672 }
673 for (G4int j=0; j<6; ++j) pols[ksteps+1][j].set(r0[j]*cosEnd,r0[j]*sinEnd,z0[j]);
674
675 // set sub-envelope and adjust extent
676 G4double emin,emax;
677 G4BoundingEnvelope benv(polygons);
678 if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
679 if (emin < pMin) pMin = emin;
680 if (emax > pMax) pMax = emax;
681 if (eminlim > pMin && emaxlim < pMax) return true; // max possible extent
682 }
683 return (pMin < pMax);
684}
std::vector< G4ThreeVector > G4ThreeVectorList
std::vector< G4TwoVector > G4TwoVectorList
Definition: G4GeomTools.hh:42
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:36
static G4bool TriangulatePolygon(const G4TwoVectorList &polygon, G4TwoVectorList &result)
Definition: G4GeomTools.cc:193
static void RemoveRedundantVertices(G4TwoVectorList &polygon, std::vector< G4int > &iout, G4double tolerance=0.0)
Definition: G4GeomTools.cc:305
static G4double PolygonArea(const G4TwoVectorList &polygon)
Definition: G4GeomTools.cc:76
G4double GetEndPhi() const
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Polycone.cc:511
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
virtual

Reimplemented from G4VSolid.

Definition at line 704 of file G4Polycone.cc.

705{
706 return new G4Polycone(*this);
707}

◆ ComputeDimensions()

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

Reimplemented from G4VSolid.

Definition at line 688 of file G4Polycone.cc.

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

◆ CopyStuff()

void G4Polycone::CopyStuff ( const G4Polycone source)
protected

Definition at line 384 of file G4Polycone.cc.

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

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

◆ Create()

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

Definition at line 172 of file G4Polycone.cc.

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

Referenced by G4Polycone(), and Reset().

◆ CreatePolyhedron()

G4Polyhedron * G4Polycone::CreatePolyhedron ( ) const
virtual

Implements G4VCSGfaceted.

Definition at line 951 of file G4Polycone.cc.

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

◆ DistanceToIn() [1/2]

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

Reimplemented from G4VCSGfaceted.

Definition at line 504 of file G4Polycone.cc.

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

◆ DistanceToIn() [2/2]

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

Reimplemented from G4VCSGfaceted.

Definition at line 487 of file G4Polycone.cc.

489{
490 //
491 // Quick test
492 //
494 return kInfinity;
495
496 //
497 // Long answer
498 //
499 return G4VCSGfaceted::DistanceToIn( p, v );
500}
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 ( )
virtual

Reimplemented from G4VCSGfaceted.

Definition at line 762 of file G4Polycone.cc.

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

◆ GetEndPhi()

◆ GetEntityType()

G4GeometryType G4Polycone::GetEntityType ( ) const
virtual

Reimplemented from G4VCSGfaceted.

Definition at line 697 of file G4Polycone.cc.

698{
699 return G4String("G4Polycone");
700}

◆ GetNumRZCorner()

◆ GetOriginalParameters()

◆ GetPointOnSurface()

G4ThreeVector G4Polycone::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 886 of file G4Polycone.cc.

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

◆ 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 ( )
virtual

Reimplemented from G4VCSGfaceted.

Definition at line 784 of file G4Polycone.cc.

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

◆ Inside()

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

Reimplemented from G4VCSGfaceted.

Definition at line 469 of file G4Polycone.cc.

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

◆ IsOpen()

G4bool G4Polycone::IsOpen ( ) const
inline

◆ operator=()

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

Definition at line 366 of file G4Polycone.cc.

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

◆ Reset()

G4bool G4Polycone::Reset ( )

Definition at line 436 of file G4Polycone.cc.

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

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

◆ SetOriginalParameters() [1/2]

◆ SetOriginalParameters() [2/2]

G4bool G4Polycone::SetOriginalParameters ( G4ReduciblePolygon rz)
protected

Definition at line 961 of file G4Polycone.cc.

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

◆ SetSurfaceElements()

void G4Polycone::SetSurfaceElements ( ) const
protected

Definition at line 823 of file G4Polycone.cc.

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

Referenced by GetPointOnSurface().

◆ StreamInfo()

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

Reimplemented from G4VCSGfaceted.

Definition at line 712 of file G4Polycone.cc.

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

Member Data Documentation

◆ corners

G4PolyconeSideRZ* G4Polycone::corners = nullptr
protected

◆ enclosingCylinder

G4EnclosingCylinder* G4Polycone::enclosingCylinder = nullptr
protected

Definition at line 180 of file G4Polycone.hh.

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

◆ endPhi

G4double G4Polycone::endPhi
protected

Definition at line 174 of file G4Polycone.hh.

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

◆ fElements

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

Definition at line 183 of file G4Polycone.hh.

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

◆ numCorner

G4int G4Polycone::numCorner
protected

Definition at line 176 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 175 of file G4Polycone.hh.

Referenced by CopyStuff(), and Create().

◆ startPhi

G4double G4Polycone::startPhi
protected

Definition at line 173 of file G4Polycone.hh.

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


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