Geant4 10.7.0
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:35
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 349 of file G4Polycone.cc.

350{
351 delete [] corners;
352 delete original_parameters;
353 delete enclosingCylinder;
354 delete fElements;
355 delete fpPolyhedron;
356 corners = nullptr;
357 original_parameters = nullptr;
358 enclosingCylinder = nullptr;
359 fElements = nullptr;
360 fpPolyhedron = nullptr;
361}
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 340 of file G4Polycone.cc.

341 : G4VCSGfaceted(a), startPhi(0.), endPhi(0.), numCorner(0)
342{
343}
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 365 of file G4Polycone.cc.

366 : G4VCSGfaceted( source )
367{
368 CopyStuff( source );
369}
void CopyStuff(const G4Polycone &source)
Definition: G4Polycone.cc:391

Member Function Documentation

◆ BoundingLimits()

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

Reimplemented from G4VSolid.

Definition at line 518 of file G4Polycone.cc.

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

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

◆ Clone()

G4VSolid * G4Polycone::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 711 of file G4Polycone.cc.

712{
713 return new G4Polycone(*this);
714}

◆ ComputeDimensions()

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

Reimplemented from G4VSolid.

Definition at line 695 of file G4Polycone.cc.

698{
699 p->ComputeDimensions(*this,n,pRep);
700}
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const

◆ CopyStuff()

void G4Polycone::CopyStuff ( const G4Polycone source)
protected

Definition at line 391 of file G4Polycone.cc.

392{
393 //
394 // Simple stuff
395 //
396 startPhi = source.startPhi;
397 endPhi = source.endPhi;
398 phiIsOpen = source.phiIsOpen;
399 numCorner = source.numCorner;
400
401 //
402 // The corner array
403 //
405
407 * sourceCorn = source.corners;
408 do // Loop checking, 13.08.2015, G.Cosmo
409 {
410 *corn = *sourceCorn;
411 } while( ++sourceCorn, ++corn < corners+numCorner );
412
413 //
414 // Original parameters
415 //
416 if (source.original_parameters)
417 {
420 }
421
422 //
423 // Enclosing cylinder
424 //
426
427 //
428 // Surface elements
429 //
430 delete fElements;
431 fElements = nullptr;
432
433 //
434 // Polyhedron
435 //
436 fRebuildPolyhedron = false;
437 delete fpPolyhedron;
438 fpPolyhedron = nullptr;
439}
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 //
224 // Phi opening? Account for some possible roundoff, and interpret
225 // nonsense value as representing no phi opening
226 //
227 if (phiTotal <= 0 || phiTotal > twopi-1E-10)
228 {
229 phiIsOpen = false;
230 startPhi = 0.;
231 endPhi = twopi;
232 }
233 else
234 {
235 phiIsOpen = true;
236
237 //
238 // Convert phi into our convention
239 //
240 startPhi = phiStart;
241 while( startPhi < 0. ) // Loop checking, 13.08.2015, G.Cosmo
242 startPhi += twopi;
243
244 endPhi = phiStart+phiTotal;
245 while( endPhi < startPhi ) // Loop checking, 13.08.2015, G.Cosmo
246 endPhi += twopi;
247 }
248
249 //
250 // Allocate corner array.
251 //
253
254 //
255 // Copy corners
256 //
258
260 iterRZ.Begin();
261 do // Loop checking, 13.08.2015, G.Cosmo
262 {
263 next->r = iterRZ.GetA();
264 next->z = iterRZ.GetB();
265 } while( ++next, iterRZ.Next() );
266
267 //
268 // Allocate face pointer array
269 //
271 faces = new G4VCSGface*[numFace];
272
273 //
274 // Construct conical faces
275 //
276 // But! Don't construct a face if both points are at zero radius!
277 //
278 G4PolyconeSideRZ* corner = corners,
279 * prev = corners + numCorner-1,
280 * nextNext;
281 G4VCSGface **face = faces;
282 do // Loop checking, 13.08.2015, G.Cosmo
283 {
284 next = corner+1;
285 if (next >= corners+numCorner) next = corners;
286 nextNext = next+1;
287 if (nextNext >= corners+numCorner) nextNext = corners;
288
289 if (corner->r < 1/kInfinity && next->r < 1/kInfinity) continue;
290
291 //
292 // We must decide here if we can dare declare one of our faces
293 // as having a "valid" normal (i.e. allBehind = true). This
294 // is never possible if the face faces "inward" in r.
295 //
296 G4bool allBehind;
297 if (corner->z > next->z)
298 {
299 allBehind = false;
300 }
301 else
302 {
303 //
304 // Otherwise, it is only true if the line passing
305 // through the two points of the segment do not
306 // split the r/z cross section
307 //
308 allBehind = !rz->BisectedBy( corner->r, corner->z,
309 next->r, next->z, kCarTolerance );
310 }
311
312 *face++ = new G4PolyconeSide( prev, corner, next, nextNext,
313 startPhi, endPhi-startPhi, phiIsOpen, allBehind );
314 } while( prev=corner, corner=next, corner > corners );
315
316 if (phiIsOpen)
317 {
318 //
319 // Construct phi open edges
320 //
321 *face++ = new G4PolyPhiFace( rz, startPhi, 0, endPhi );
322 *face++ = new G4PolyPhiFace( rz, endPhi, 0, startPhi );
323 }
324
325 //
326 // We might have dropped a face or two: recalculate numFace
327 //
328 numFace = face-faces;
329
330 //
331 // Make enclosingCylinder
332 //
334 new G4EnclosingCylinder( rz, phiIsOpen, phiStart, phiTotal );
335}
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

Referenced by G4Polycone(), and Reset().

◆ CreatePolyhedron()

G4Polyhedron * G4Polycone::CreatePolyhedron ( ) const
virtual

Implements G4VCSGfaceted.

Definition at line 958 of file G4Polycone.cc.

959{
960 //
961 // This has to be fixed in visualization. Fake it for the moment.
962 //
963
970}

◆ DistanceToIn() [1/2]

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

Reimplemented from G4VCSGfaceted.

Definition at line 511 of file G4Polycone.cc.

512{
514}
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 494 of file G4Polycone.cc.

496{
497 //
498 // Quick test
499 //
501 return kInfinity;
502
503 //
504 // Long answer
505 //
506 return G4VCSGfaceted::DistanceToIn( p, v );
507}
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 769 of file G4Polycone.cc.

770{
771 if (fCubicVolume == 0.)
772 {
773 G4double total = 0.;
774 G4int nrz = GetNumRZCorner();
775 G4PolyconeSideRZ a = GetCorner(nrz - 1);
776 for (G4int i=0; i<nrz; ++i)
777 {
779 total += (b.r*b.r + b.r*a.r + a.r*a.r)*(b.z - a.z);
780 a = b;
781 }
782 fCubicVolume = std::abs(total)*(GetEndPhi() - GetStartPhi())/6.;
783 }
784 return fCubicVolume;
785}
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 704 of file G4Polycone.cc.

705{
706 return G4String("G4Polycone");
707}

◆ GetNumRZCorner()

◆ GetOriginalParameters()

◆ GetPointOnSurface()

G4ThreeVector G4Polycone::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 893 of file G4Polycone.cc.

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

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

792{
793 if (fSurfaceArea == 0.)
794 {
795 // phi cut area
796 G4int nrz = GetNumRZCorner();
797 G4double scut = 0.;
798 if (IsOpen())
799 {
800 G4PolyconeSideRZ a = GetCorner(nrz - 1);
801 for (G4int i=0; i<nrz; ++i)
802 {
804 scut += a.r*b.z - a.z*b.r;
805 a = b;
806 }
807 scut = std::abs(scut);
808 }
809 // lateral surface area
810 G4double slat = 0;
811 G4PolyconeSideRZ a = GetCorner(nrz - 1);
812 for (G4int i=0; i<nrz; ++i)
813 {
815 G4double h = std::sqrt((b.r - a.r)*(b.r - a.r) + (b.z - a.z)*(b.z - a.z));
816 slat += (b.r + a.r)*h;
817 a = b;
818 }
819 slat *= (GetEndPhi() - GetStartPhi())/2.;
820 fSurfaceArea = scut + slat;
821 }
822 return fSurfaceArea;
823}
G4double fSurfaceArea

◆ Inside()

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

Reimplemented from G4VCSGfaceted.

Definition at line 476 of file G4Polycone.cc.

477{
478 //
479 // Quick test
480 //
482
483 //
484 // Long answer
485 //
486 return G4VCSGfaceted::Inside(p);
487}
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 373 of file G4Polycone.cc.

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

◆ Reset()

G4bool G4Polycone::Reset ( )

Definition at line 443 of file G4Polycone.cc.

444{
445 //
446 // Clear old setup
447 //
449 delete [] corners;
450 delete enclosingCylinder;
451 delete fElements;
452 corners = nullptr;
453 fElements = nullptr;
454 enclosingCylinder = nullptr;
455
456 //
457 // Rebuild polycone
458 //
466 delete rz;
467
468 return false;
469}

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

◆ SetOriginalParameters() [1/2]

◆ SetOriginalParameters() [2/2]

G4bool G4Polycone::SetOriginalParameters ( G4ReduciblePolygon rz)
protected

Definition at line 974 of file G4Polycone.cc.

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

◆ SetSurfaceElements()

void G4Polycone::SetSurfaceElements ( ) const
protected

Definition at line 830 of file G4Polycone.cc.

831{
832 fElements = new std::vector<G4Polycone::surface_element>;
833 G4double total = 0.;
834 G4int nrz = GetNumRZCorner();
835
836 // set lateral surface elements
837 G4double dphi = GetEndPhi() - GetStartPhi();
838 G4int ia = nrz - 1;
839 for (G4int ib=0; ib<nrz; ++ib)
840 {
844 selem.i0 = ia;
845 selem.i1 = ib;
846 selem.i2 = -1;
847 ia = ib;
848 if (a.r == 0. && b.r == 0.) continue;
849 G4double h = std::sqrt((b.r - a.r)*(b.r - a.r) + (b.z - a.z)*(b.z - a.z));
850 total += 0.5*dphi*(b.r + a.r)*h;
851 selem.area = total;
852 fElements->push_back(selem);
853 }
854
855 // set elements for phi cuts
856 if (IsOpen())
857 {
858 G4TwoVectorList contourRZ;
859 std::vector<G4int> triangles;
860 for (G4int i=0; i<nrz; ++i)
861 {
862 G4PolyconeSideRZ corner = GetCorner(i);
863 contourRZ.push_back(G4TwoVector(corner.r, corner.z));
864 }
865 G4GeomTools::TriangulatePolygon(contourRZ, triangles);
866 G4int ntria = triangles.size();
867 for (G4int i=0; i<ntria; i+=3)
868 {
870 selem.i0 = triangles[i];
871 selem.i1 = triangles[i+1];
872 selem.i2 = triangles[i+2];
873 G4PolyconeSideRZ a = GetCorner(selem.i0);
874 G4PolyconeSideRZ b = GetCorner(selem.i1);
875 G4PolyconeSideRZ c = GetCorner(selem.i2);
876 G4double stria =
877 std::abs(G4GeomTools::TriangleArea(a.r, a.z, b.r, b.z, c.r, c.z));
878 total += stria;
879 selem.area = total;
880 fElements->push_back(selem); // start phi
881 total += stria;
882 selem.area = total;
883 selem.i0 += nrz;
884 fElements->push_back(selem); // end phi
885 }
886 }
887}
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 719 of file G4Polycone.cc.

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

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(), 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(), 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(), SetOriginalParameters(), and StreamInfo().


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