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

G4PolyconeSide is a utility class implementing a face that represents one conical side of a polycone. More...

#include <G4PolyconeSide.hh>

Inheritance diagram for G4PolyconeSide:

Public Member Functions

 G4PolyconeSide (const G4PolyconeSideRZ *prevRZ, const G4PolyconeSideRZ *tail, const G4PolyconeSideRZ *head, const G4PolyconeSideRZ *nextRZ, G4double phiStart, G4double deltaPhi, G4bool phiIsOpen, G4bool isAllBehind=false)
 ~G4PolyconeSide () override
 G4PolyconeSide (const G4PolyconeSide &source)
G4PolyconeSideoperator= (const G4PolyconeSide &source)
G4bool Intersect (const G4ThreeVector &p, const G4ThreeVector &v, G4bool outgoing, G4double surfTolerance, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal, G4bool &isAllBehind) override
G4double Distance (const G4ThreeVector &p, G4bool outgoing) override
EInside Inside (const G4ThreeVector &p, G4double tolerance, G4double *bestDistance) override
G4ThreeVector Normal (const G4ThreeVector &p, G4double *bestDistance) override
G4double Extent (const G4ThreeVector axis) override
void CalculateExtent (const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList) override
G4VCSGfaceClone () override
G4double SurfaceArea () override
G4ThreeVector GetPointOnFace () override
 G4PolyconeSide (__void__ &)
G4int GetInstanceID () const
Public Member Functions inherited from G4VCSGface
 G4VCSGface ()=default
virtual ~G4VCSGface ()=default

Static Public Member Functions

static const G4PlSideManagerGetSubInstanceManager ()

Detailed Description

G4PolyconeSide is a utility class implementing a face that represents one conical side of a polycone.

Definition at line 90 of file G4PolyconeSide.hh.

Constructor & Destructor Documentation

◆ G4PolyconeSide() [1/3]

G4PolyconeSide::G4PolyconeSide ( const G4PolyconeSideRZ * prevRZ,
const G4PolyconeSideRZ * tail,
const G4PolyconeSideRZ * head,
const G4PolyconeSideRZ * nextRZ,
G4double phiStart,
G4double deltaPhi,
G4bool phiIsOpen,
G4bool isAllBehind = false )

Constructor for the conical side of a polycone.

Parameters
[in]prevRZPointer to previous r,Z section.
[in]tailPointer to r,Z tail of section.
[in]headPointer to r,Z head of section.
[in]nextRZPointer to next r,Z section.
[in]phiStartInitial Phi starting angle.
[in]deltaPhiTotal Phi angle.
[in]phiIsOpenFlag indicating if it is a Phi section.
[in]isAllBehindIndicating if entire surface is behind normal.

Definition at line 67 of file G4PolyconeSide.cc.

75{
76 instanceID = subInstanceManager.CreateSubInstance();
77
79 G4MT_pcphix = 0.0; G4MT_pcphiy = 0.0; G4MT_pcphiz = 0.0; G4MT_pcphik = 0.0;
80
81 //
82 // Record values
83 //
84 r[0] = tail->r; z[0] = tail->z;
85 r[1] = head->r; z[1] = head->z;
86
87 phiIsOpen = thePhiIsOpen;
88 if (phiIsOpen)
89 {
90 deltaPhi = theDeltaPhi;
91 startPhi = thePhiStart;
92
93 //
94 // Set phi values to our conventions
95 //
96 while (deltaPhi < 0.0) // Loop checking, 13.08.2015, G.Cosmo
97 {
98 deltaPhi += twopi;
99 }
100 while (startPhi < 0.0) // Loop checking, 13.08.2015, G.Cosmo
101 {
102 startPhi += twopi;
103 }
104
105 //
106 // Calculate corner coordinates
107 //
108 ncorners = 4;
109 corners = new G4ThreeVector[ncorners];
110
111 corners[0] = G4ThreeVector( tail->r*std::cos(startPhi),
112 tail->r*std::sin(startPhi), tail->z );
113 corners[1] = G4ThreeVector( head->r*std::cos(startPhi),
114 head->r*std::sin(startPhi), head->z );
115 corners[2] = G4ThreeVector( tail->r*std::cos(startPhi+deltaPhi),
116 tail->r*std::sin(startPhi+deltaPhi), tail->z );
117 corners[3] = G4ThreeVector( head->r*std::cos(startPhi+deltaPhi),
118 head->r*std::sin(startPhi+deltaPhi), head->z );
119 }
120 else
121 {
122 deltaPhi = twopi;
123 startPhi = 0.0;
124 }
125
126 allBehind = isAllBehind;
127
128 //
129 // Make our intersecting cone
130 //
131 cone = new G4IntersectingCone( r, z );
132
133 //
134 // Calculate vectors in r,z space
135 //
136 rS = r[1]-r[0]; zS = z[1]-z[0];
137 length = std::sqrt( rS*rS + zS*zS);
138 rS /= length; zS /= length;
139
140 rNorm = +zS;
141 zNorm = -rS;
142
143 G4double lAdj;
144
145 prevRS = r[0]-prevRZ->r;
146 prevZS = z[0]-prevRZ->z;
147 lAdj = std::sqrt( prevRS*prevRS + prevZS*prevZS );
148 prevRS /= lAdj;
149 prevZS /= lAdj;
150
151 rNormEdge[0] = rNorm + prevZS;
152 zNormEdge[0] = zNorm - prevRS;
153 lAdj = std::sqrt( rNormEdge[0]*rNormEdge[0] + zNormEdge[0]*zNormEdge[0] );
154 rNormEdge[0] /= lAdj;
155 zNormEdge[0] /= lAdj;
156
157 nextRS = nextRZ->r-r[1];
158 nextZS = nextRZ->z-z[1];
159 lAdj = std::sqrt( nextRS*nextRS + nextZS*nextZS );
160 nextRS /= lAdj;
161 nextZS /= lAdj;
162
163 rNormEdge[1] = rNorm + nextZS;
164 zNormEdge[1] = zNorm - nextRS;
165 lAdj = std::sqrt( rNormEdge[1]*rNormEdge[1] + zNormEdge[1]*zNormEdge[1] );
166 rNormEdge[1] /= lAdj;
167 zNormEdge[1] /= lAdj;
168}
#define G4MT_pcphiz
#define G4MT_pcphik
#define G4MT_pcphix
#define G4MT_pcphiy
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()

Referenced by Clone(), G4PolyconeSide(), and operator=().

◆ ~G4PolyconeSide()

G4PolyconeSide::~G4PolyconeSide ( )
override

Destructor.

Definition at line 187 of file G4PolyconeSide.cc.

188{
189 delete cone;
190 if (phiIsOpen) { delete [] corners; }
191}

◆ G4PolyconeSide() [2/3]

G4PolyconeSide::G4PolyconeSide ( const G4PolyconeSide & source)

Copy constructor and assignment operator.

Definition at line 195 of file G4PolyconeSide.cc.

196{
197 instanceID = subInstanceManager.CreateSubInstance();
198
199 CopyStuff( source );
200}

◆ G4PolyconeSide() [3/3]

G4PolyconeSide::G4PolyconeSide ( __void__ & )

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

Definition at line 173 of file G4PolyconeSide.cc.

174 : startPhi(0.), deltaPhi(0.),
175 rNorm(0.), zNorm(0.), rS(0.), zS(0.), length(0.),
176 prevRS(0.), prevZS(0.), nextRS(0.), nextZS(0.),
177 kCarTolerance(0.), instanceID(0)
178{
179 r[0] = r[1] = 0.;
180 z[0] = z[1] = 0.;
181 rNormEdge[0]= rNormEdge[1] = 0.;
182 zNormEdge[0]= zNormEdge[1] = 0.;
183}
const G4double kCarTolerance

Member Function Documentation

◆ CalculateExtent()

void G4PolyconeSide::CalculateExtent ( const EAxis axis,
const G4VoxelLimits & voxelLimit,
const G4AffineTransform & tranform,
G4SolidExtentList & extentList )
overridevirtual

Calculates the extent of the face for the voxel navigator.

Parameters
[in]axisThe axis in which to check the shapes 3D extent against.
[in]voxelLimitLimits along x, y, and/or z axes.
[in]tranformA coordinate transformation on which to apply to the shape before testing.
[out]extentListThe list of (voxel) extents along the axis.

Implements G4VCSGface.

Definition at line 549 of file G4PolyconeSide.cc.

553{
554 G4ClippablePolygon polygon;
555
556 //
557 // Here we will approximate (ala G4Cons) and divide our conical section
558 // into segments, like G4Polyhedra. When doing so, the radius
559 // is extented far enough such that the segments always lie
560 // just outside the surface of the conical section we are
561 // approximating.
562 //
563
564 //
565 // Choose phi size of our segment(s) based on constants as
566 // defined in meshdefs.hh
567 //
568 G4int numPhi = (G4int)(deltaPhi/kMeshAngleDefault) + 1;
569 if (numPhi < kMinMeshSections)
570 {
571 numPhi = kMinMeshSections;
572 }
573 else if (numPhi > kMaxMeshSections)
574 {
575 numPhi = kMaxMeshSections;
576 }
577
578 G4double sigPhi = deltaPhi/numPhi;
579
580 //
581 // Determine radius factor to keep segments outside
582 //
583 G4double rFudge = 1.0/std::cos(0.5*sigPhi);
584
585 //
586 // Decide which radius to use on each end of the side,
587 // and whether a transition mesh is required
588 //
589 // {r0,z0} - Beginning of this side
590 // {r1,z1} - Ending of this side
591 // {r2,z0} - Beginning of transition piece connecting previous
592 // side (and ends at beginning of this side)
593 //
594 // So, order is 2 --> 0 --> 1.
595 // -------
596 //
597 // r2 < 0 indicates that no transition piece is required
598 //
599 G4double r0, r1, r2, z0, z1;
600
601 r2 = -1; // By default: no transition piece
602
603 if (rNorm < -DBL_MIN)
604 {
605 //
606 // This side faces *inward*, and so our mesh has
607 // the same radius
608 //
609 r1 = r[1];
610 z1 = z[1];
611 z0 = z[0];
612 r0 = r[0];
613
614 r2 = -1;
615
616 if (prevZS > DBL_MIN)
617 {
618 //
619 // The previous side is facing outwards
620 //
621 if ( prevRS*zS - prevZS*rS > 0 )
622 {
623 //
624 // Transition was convex: build transition piece
625 //
626 if (r[0] > DBL_MIN) { r2 = r[0]*rFudge; }
627 }
628 else
629 {
630 //
631 // Transition was concave: short this side
632 //
633 FindLineIntersect( z0, r0, zS, rS,
634 z0, r0*rFudge, prevZS, prevRS*rFudge, z0, r0 );
635 }
636 }
637
638 if ( nextZS > DBL_MIN && (rS*nextZS - zS*nextRS < 0) )
639 {
640 //
641 // The next side is facing outwards, forming a
642 // concave transition: short this side
643 //
644 FindLineIntersect( z1, r1, zS, rS,
645 z1, r1*rFudge, nextZS, nextRS*rFudge, z1, r1 );
646 }
647 }
648 else if (rNorm > DBL_MIN)
649 {
650 //
651 // This side faces *outward* and is given a boost to
652 // it radius
653 //
654 r0 = r[0]*rFudge;
655 z0 = z[0];
656 r1 = r[1]*rFudge;
657 z1 = z[1];
658
659 if (prevZS < -DBL_MIN)
660 {
661 //
662 // The previous side is facing inwards
663 //
664 if ( prevRS*zS - prevZS*rS > 0 )
665 {
666 //
667 // Transition was convex: build transition piece
668 //
669 if (r[0] > DBL_MIN) { r2 = r[0]; }
670 }
671 else
672 {
673 //
674 // Transition was concave: short this side
675 //
676 FindLineIntersect( z0, r0, zS, rS*rFudge,
677 z0, r[0], prevZS, prevRS, z0, r0 );
678 }
679 }
680
681 if ( nextZS < -DBL_MIN && (rS*nextZS - zS*nextRS < 0) )
682 {
683 //
684 // The next side is facing inwards, forming a
685 // concave transition: short this side
686 //
687 FindLineIntersect( z1, r1, zS, rS*rFudge,
688 z1, r[1], nextZS, nextRS, z1, r1 );
689 }
690 }
691 else
692 {
693 //
694 // This side is perpendicular to the z axis (is a disk)
695 //
696 // Whether or not r0 needs a rFudge factor depends
697 // on the normal of the previous edge. Similar with r1
698 // and the next edge. No transition piece is required.
699 //
700 r0 = r[0];
701 r1 = r[1];
702 z0 = z[0];
703 z1 = z[1];
704
705 if (prevZS > DBL_MIN) { r0 *= rFudge; }
706 if (nextZS > DBL_MIN) { r1 *= rFudge; }
707 }
708
709 //
710 // Loop
711 //
712 G4double phi = startPhi,
713 cosPhi = std::cos(phi),
714 sinPhi = std::sin(phi);
715
716 G4ThreeVector v0( r0*cosPhi, r0*sinPhi, z0 ),
717 v1( r1*cosPhi, r1*sinPhi, z1 ),
718 v2, w0, w1, w2;
719 transform.ApplyPointTransform( v0 );
720 transform.ApplyPointTransform( v1 );
721
722 if (r2 >= 0)
723 {
724 v2 = G4ThreeVector( r2*cosPhi, r2*sinPhi, z0 );
725 transform.ApplyPointTransform( v2 );
726 }
727
728 do // Loop checking, 13.08.2015, G.Cosmo
729 {
730 phi += sigPhi;
731 if (numPhi == 1) { phi = startPhi+deltaPhi; } // Try to avoid roundoff
732 cosPhi = std::cos(phi),
733 sinPhi = std::sin(phi);
734
735 w0 = G4ThreeVector( r0*cosPhi, r0*sinPhi, z0 );
736 w1 = G4ThreeVector( r1*cosPhi, r1*sinPhi, z1 );
737 transform.ApplyPointTransform( w0 );
738 transform.ApplyPointTransform( w1 );
739
740 G4ThreeVector deltaV = r0 > r1 ? w0-v0 : w1-v1;
741
742 //
743 // Build polygon, taking special care to keep the vertices
744 // in order
745 //
746 polygon.ClearAllVertices();
747
748 polygon.AddVertexInOrder( v0 );
749 polygon.AddVertexInOrder( v1 );
750 polygon.AddVertexInOrder( w1 );
751 polygon.AddVertexInOrder( w0 );
752
753 //
754 // Get extent
755 //
756 if (polygon.PartialClip( voxelLimit, axis ))
757 {
758 //
759 // Get dot product of normal with target axis
760 //
761 polygon.SetNormal( deltaV.cross(v1-v0).unit() );
762
763 extentList.AddSurface( polygon );
764 }
765
766 if (r2 >= 0)
767 {
768 //
769 // Repeat, for transition piece
770 //
771 w2 = G4ThreeVector( r2*cosPhi, r2*sinPhi, z0 );
772 transform.ApplyPointTransform( w2 );
773
774 polygon.ClearAllVertices();
775
776 polygon.AddVertexInOrder( v2 );
777 polygon.AddVertexInOrder( v0 );
778 polygon.AddVertexInOrder( w0 );
779 polygon.AddVertexInOrder( w2 );
780
781 if (polygon.PartialClip( voxelLimit, axis ))
782 {
783 polygon.SetNormal( deltaV.cross(v0-v2).unit() );
784
785 extentList.AddSurface( polygon );
786 }
787
788 v2 = w2;
789 }
790
791 //
792 // Next vertex
793 //
794 v0 = w0;
795 v1 = w1;
796 } while( --numPhi > 0 );
797
798 //
799 // We are almost done. But, it is important that we leave no
800 // gaps in the surface of our solid. By using rFudge, however,
801 // we've done exactly that, if we have a phi segment.
802 // Add two additional faces if necessary
803 //
804 if (phiIsOpen && rNorm > DBL_MIN)
805 {
806 cosPhi = std::cos(startPhi);
807 sinPhi = std::sin(startPhi);
808
809 G4ThreeVector a0( r[0]*cosPhi, r[0]*sinPhi, z[0] ),
810 a1( r[1]*cosPhi, r[1]*sinPhi, z[1] ),
811 b0( r0*cosPhi, r0*sinPhi, z[0] ),
812 b1( r1*cosPhi, r1*sinPhi, z[1] );
813
814 transform.ApplyPointTransform( a0 );
815 transform.ApplyPointTransform( a1 );
816 transform.ApplyPointTransform( b0 );
817 transform.ApplyPointTransform( b1 );
818
819 polygon.ClearAllVertices();
820
821 polygon.AddVertexInOrder( a0 );
822 polygon.AddVertexInOrder( a1 );
823 polygon.AddVertexInOrder( b0 );
824 polygon.AddVertexInOrder( b1 );
825
826 if (polygon.PartialClip( voxelLimit , axis))
827 {
828 G4ThreeVector normal( sinPhi, -cosPhi, 0 );
829 polygon.SetNormal( transform.TransformAxis( normal ) );
830
831 extentList.AddSurface( polygon );
832 }
833
834 cosPhi = std::cos(startPhi+deltaPhi);
835 sinPhi = std::sin(startPhi+deltaPhi);
836
837 a0 = G4ThreeVector( r[0]*cosPhi, r[0]*sinPhi, z[0] ),
838 a1 = G4ThreeVector( r[1]*cosPhi, r[1]*sinPhi, z[1] ),
839 b0 = G4ThreeVector( r0*cosPhi, r0*sinPhi, z[0] ),
840 b1 = G4ThreeVector( r1*cosPhi, r1*sinPhi, z[1] );
841 transform.ApplyPointTransform( a0 );
842 transform.ApplyPointTransform( a1 );
843 transform.ApplyPointTransform( b0 );
844 transform.ApplyPointTransform( b1 );
845
846 polygon.ClearAllVertices();
847
848 polygon.AddVertexInOrder( a0 );
849 polygon.AddVertexInOrder( a1 );
850 polygon.AddVertexInOrder( b0 );
851 polygon.AddVertexInOrder( b1 );
852
853 if (polygon.PartialClip( voxelLimit, axis ))
854 {
855 G4ThreeVector normal( -sinPhi, cosPhi, 0 );
856 polygon.SetNormal( transform.TransformAxis( normal ) );
857
858 extentList.AddSurface( polygon );
859 }
860 }
861
862 return;
863}
const G4double a0
int G4int
Definition G4Types.hh:85
Hep3Vector unit() const
Hep3Vector cross(const Hep3Vector &) const
G4bool PartialClip(const G4VoxelLimits &voxelLimit, const EAxis IgnoreMe)
void AddVertexInOrder(const G4ThreeVector &vertex)
void SetNormal(const G4ThreeVector &newNormal)
void AddSurface(const G4ClippablePolygon &surface)
const G4int kMaxMeshSections
Definition meshdefs.hh:44
const G4int kMinMeshSections
Definition meshdefs.hh:41
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668
#define DBL_MIN
Definition templates.hh:54

◆ Clone()

G4VCSGface * G4PolyconeSide::Clone ( )
inlineoverridevirtual

Method invoked by the copy constructor or the assignment operator. Its purpose is to return a pointer to a duplicate copy of the face.

Implements G4VCSGface.

Definition at line 200 of file G4PolyconeSide.hh.

200{ return new G4PolyconeSide( *this ); }
G4PolyconeSide(const G4PolyconeSideRZ *prevRZ, const G4PolyconeSideRZ *tail, const G4PolyconeSideRZ *head, const G4PolyconeSideRZ *nextRZ, G4double phiStart, G4double deltaPhi, G4bool phiIsOpen, G4bool isAllBehind=false)

◆ Distance()

G4double G4PolyconeSide::Distance ( const G4ThreeVector & p,
G4bool outgoing )
overridevirtual

Determines the distance of a point from either the inside or outside surfaces of the face.

Parameters
[in]pPosition.
[in]outgoingFlag, true, to consider only inside surfaces or false, to consider only outside surfaces.
Returns
The distance to the closest surface satisfying requirements or kInfinity if no such surface exists.

Implements G4VCSGface.

Definition at line 399 of file G4PolyconeSide.cc.

400{
401 G4double normSign = outgoing ? -1 : +1;
402 G4double distFrom, distOut2;
403
404 //
405 // We have two tries for each hemisphere. Try the closest first.
406 //
407 distFrom = normSign*DistanceAway( p, false, distOut2 );
408 if (distFrom > -0.5*kCarTolerance )
409 {
410 //
411 // Good answer
412 //
413 if (distOut2 > 0)
414 {
415 return std::sqrt( distFrom*distFrom + distOut2 );
416 }
417 return std::fabs(distFrom);
418 }
419
420 //
421 // Try second side.
422 //
423 distFrom = normSign*DistanceAway( p, true, distOut2 );
424 if (distFrom > -0.5*kCarTolerance)
425 {
426 if (distOut2 > 0)
427 {
428 return std::sqrt( distFrom*distFrom + distOut2 );
429 }
430 return std::fabs(distFrom);
431 }
432
433 return kInfinity;
434}

◆ Extent()

G4double G4PolyconeSide::Extent ( const G4ThreeVector axis)
overridevirtual

Returns the face extent along the axis.

Parameters
[in]axisUnit vector defining the direction.
Returns
The largest point along the given axis of the face's extent.

Implements G4VCSGface.

Definition at line 485 of file G4PolyconeSide.cc.

486{
487 if (axis.perp2() < DBL_MIN)
488 {
489 //
490 // Special case
491 //
492 return axis.z() < 0 ? -cone->ZLo() : cone->ZHi();
493 }
494
495 //
496 // Is the axis pointing inside our phi gap?
497 //
498 if (phiIsOpen)
499 {
500 G4double phi = GetPhi(axis);
501 while( phi < startPhi ) // Loop checking, 13.08.2015, G.Cosmo
502 {
503 phi += twopi;
504 }
505
506 if (phi > deltaPhi+startPhi)
507 {
508 //
509 // Yeah, looks so. Make four three vectors defining the phi
510 // opening
511 //
512 G4double cosP = std::cos(startPhi), sinP = std::sin(startPhi);
513 G4ThreeVector a( r[0]*cosP, r[0]*sinP, z[0] );
514 G4ThreeVector b( r[1]*cosP, r[1]*sinP, z[1] );
515 cosP = std::cos(startPhi+deltaPhi); sinP = std::sin(startPhi+deltaPhi);
516 G4ThreeVector c( r[0]*cosP, r[0]*sinP, z[0] );
517 G4ThreeVector d( r[1]*cosP, r[1]*sinP, z[1] );
518
519 G4double ad = axis.dot(a),
520 bd = axis.dot(b),
521 cd = axis.dot(c),
522 dd = axis.dot(d);
523
524 if (bd > ad) { ad = bd; }
525 if (cd > ad) { ad = cd; }
526 if (dd > ad) { ad = dd; }
527
528 return ad;
529 }
530 }
531
532 //
533 // Check either end
534 //
535 G4double aPerp = axis.perp();
536
537 G4double a = aPerp*r[0] + axis.z()*z[0];
538 G4double b = aPerp*r[1] + axis.z()*z[1];
539
540 if (b > a) { a = b; }
541
542 return a;
543}

◆ GetInstanceID()

G4int G4PolyconeSide::GetInstanceID ( ) const
inline

Returns the instance ID.

Definition at line 222 of file G4PolyconeSide.hh.

222{ return instanceID; }

◆ GetPointOnFace()

G4ThreeVector G4PolyconeSide::GetPointOnFace ( )
overridevirtual

Returns a random point located and uniformly distributed on the face.

Implements G4VCSGface.

Definition at line 1205 of file G4PolyconeSide.cc.

1206{
1207 G4double x,y,zz;
1208 G4double rr,phi,dz,dr;
1209 dr=r[1]-r[0];dz=z[1]-z[0];
1210 phi=startPhi+deltaPhi*G4QuickRand();
1211 rr=r[0]+dr*G4QuickRand();
1212
1213 x=rr*std::cos(phi);
1214 y=rr*std::sin(phi);
1215
1216 // PolyconeSide has a Ring Form
1217 //
1218 if (dz==0.)
1219 {
1220 zz=z[0];
1221 }
1222 else
1223 {
1224 if(dr==0.) // PolyconeSide has a Tube Form
1225 {
1226 zz = z[0]+dz*G4QuickRand();
1227 }
1228 else
1229 {
1230 zz = z[0]+(rr-r[0])*dz/dr;
1231 }
1232 }
1233
1234 return {x,y,zz};
1235}
G4double G4QuickRand(uint32_t seed=0)

◆ GetSubInstanceManager()

const G4PlSideManager & G4PolyconeSide::GetSubInstanceManager ( )
static

Returns the private data instance manager.

Definition at line 57 of file G4PolyconeSide.cc.

58{
59 return subInstanceManager;
60}

Referenced by G4SolidsWorkspace::G4SolidsWorkspace().

◆ Inside()

EInside G4PolyconeSide::Inside ( const G4ThreeVector & p,
G4double tolerance,
G4double * bestDistance )
overridevirtual

Determines whether a point is inside, outside, or on the surface of the face.

Parameters
[in]pPosition.
[in]toleranceTolerance defining the bounds of the "kSurface", nominally equal to kCarTolerance/2.
[out]bestDistanceDistance to the closest surface (in or out).
Returns
kInside if the point is closest to the inside surface; kOutside if the point is closest to the outside surface; kSurface if the point is withing tolerance of the surface.

Implements G4VCSGface.

Definition at line 438 of file G4PolyconeSide.cc.

441{
442 G4double distFrom, distOut2, dist2;
443 G4double edgeRZnorm;
444
445 distFrom = DistanceAway( p, distOut2, &edgeRZnorm );
446 dist2 = distFrom*distFrom + distOut2;
447
448 *bestDistance = std::sqrt( dist2);
449
450 // Okay then, inside or out?
451 //
452 if ( (std::fabs(edgeRZnorm) < tolerance)
453 && (distOut2< tolerance*tolerance) )
454 {
455 return kSurface;
456 }
457 if (edgeRZnorm < 0)
458 {
459 return kInside;
460 }
461
462 return kOutside;
463}
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68
@ kSurface
Definition geomdefs.hh:69

◆ Intersect()

G4bool G4PolyconeSide::Intersect ( const G4ThreeVector & p,
const G4ThreeVector & v,
G4bool outgoing,
G4double surfTolerance,
G4double & distance,
G4double & distFromSurface,
G4ThreeVector & normal,
G4bool & isAllBehind )
overridevirtual

Determines the distance along a line to the face.

Parameters
[in]pPosition.
[in]vDirection (assumed to be a unit vector).
[in]outgoingFlag true, to consider only inside surfaces; false, to consider only outside surfaces.
[in]surfToleranceMinimum distance from the surface.
[out]distanceDistance to intersection.
[out]distFromSurfaceDistance from surface (along surface normal), < 0 if the point is in front of the surface.
[out]normalNormal of surface at intersection point.
[out]allBehindFlag, true, if entire surface is behind normal.
Returns
true if there is an intersection, false otherwise.

Implements G4VCSGface.

Definition at line 264 of file G4PolyconeSide.cc.

272{
273 G4double s1=0., s2=0.;
274 G4double normSign = outgoing ? +1 : -1;
275
276 isAllBehind = allBehind;
277
278 //
279 // Check for two possible intersections
280 //
281 G4int nside = cone->LineHitsCone( p, v, &s1, &s2 );
282 if (nside == 0) { return false; }
283
284 //
285 // Check the first side first, since it is (supposed to be) closest
286 //
287 G4ThreeVector hit = p + s1*v;
288
289 if (PointOnCone( hit, normSign, p, v, normal ))
290 {
291 //
292 // Good intersection! What about the normal?
293 //
294 if (normSign*v.dot(normal) > 0)
295 {
296 //
297 // We have a valid intersection, but it could very easily
298 // be behind the point. To decide if we tolerate this,
299 // we have to see if the point p is on the surface near
300 // the intersecting point.
301 //
302 // What does it mean exactly for the point p to be "near"
303 // the intersection? It means that if we draw a line from
304 // p to the hit, the line remains entirely within the
305 // tolerance bounds of the cone. To test this, we can
306 // ask if the normal is correct near p.
307 //
308 G4double pr = p.perp();
309 if (pr < DBL_MIN) { pr = DBL_MIN; }
310 G4ThreeVector pNormal( rNorm*p.x()/pr, rNorm*p.y()/pr, zNorm );
311 if (normSign*v.dot(pNormal) > 0)
312 {
313 //
314 // p and intersection in same hemisphere
315 //
316 G4double distOutside2;
317 distFromSurface = -normSign*DistanceAway( p, false, distOutside2 );
318 if (distOutside2 < surfTolerance*surfTolerance)
319 {
320 if (distFromSurface > -surfTolerance)
321 {
322 //
323 // We are just inside or away from the
324 // surface. Accept *any* value of distance.
325 //
326 distance = s1;
327 return true;
328 }
329 }
330 }
331 else
332 {
333 distFromSurface = s1;
334 }
335
336 //
337 // Accept positive distances
338 //
339 if (s1 > 0)
340 {
341 distance = s1;
342 return true;
343 }
344 }
345 }
346
347 if (nside==1) { return false;
348}
349
350 //
351 // Well, try the second hit
352 //
353 hit = p + s2*v;
354
355 if (PointOnCone( hit, normSign, p, v, normal ))
356 {
357 //
358 // Good intersection! What about the normal?
359 //
360 if (normSign*v.dot(normal) > 0)
361 {
362 G4double pr = p.perp();
363 if (pr < DBL_MIN) { pr = DBL_MIN; }
364 G4ThreeVector pNormal( rNorm*p.x()/pr, rNorm*p.y()/pr, zNorm );
365 if (normSign*v.dot(pNormal) > 0)
366 {
367 G4double distOutside2;
368 distFromSurface = -normSign*DistanceAway( p, false, distOutside2 );
369 if (distOutside2 < surfTolerance*surfTolerance)
370 {
371 if (distFromSurface > -surfTolerance)
372 {
373 distance = s2;
374 return true;
375 }
376 }
377 }
378 else
379 {
380 distFromSurface = s2;
381 }
382
383 if (s2 > 0)
384 {
385 distance = s2;
386 return true;
387 }
388 }
389 }
390
391 //
392 // Better luck next time
393 //
394 return false;
395}
double x() const
double y() const
double perp() const

◆ Normal()

G4ThreeVector G4PolyconeSide::Normal ( const G4ThreeVector & p,
G4double * bestDistance )
overridevirtual

Returns the normal of surface closest to the point.

Parameters
[in]pPosition.
[out]bestDistanceDistance to the closest surface (in or out).
Returns
The normal of the surface nearest the point.

Implements G4VCSGface.

Definition at line 467 of file G4PolyconeSide.cc.

469{
470 if (p == G4ThreeVector(0.,0.,0.)) { return p; }
471
472 G4double dFrom, dOut2;
473
474 dFrom = DistanceAway( p, false, dOut2 );
475
476 *bestDistance = std::sqrt( dFrom*dFrom + dOut2 );
477
478 G4double rds = p.perp();
479 if (rds!=0.) { return {rNorm*p.x()/rds,rNorm*p.y()/rds,zNorm}; }
480 return G4ThreeVector( 0.,0., zNorm ).unit();
481}

◆ operator=()

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

Definition at line 204 of file G4PolyconeSide.cc.

205{
206 if (this == &source) { return *this; }
207
208 delete cone;
209 if (phiIsOpen) { delete [] corners; }
210
211 CopyStuff( source );
212
213 return *this;
214}

◆ SurfaceArea()

G4double G4PolyconeSide::SurfaceArea ( )
overridevirtual

Returning an estimation of the face surface area, in internal units.

Implements G4VCSGface.

Definition at line 1193 of file G4PolyconeSide.cc.

1194{
1195 if(fSurfaceArea==0.)
1196 {
1197 fSurfaceArea = (r[0]+r[1])* std::sqrt(sqr(r[0]-r[1])+sqr(z[0]-z[1]));
1198 fSurfaceArea *= 0.5*(deltaPhi);
1199 }
1200 return fSurfaceArea;
1201}
T sqr(const T &x)
Definition templates.hh:128

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