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

G4Hype is a tube with hyperbolic profile; it describes an hyperbolic volume with curved sides parallel to the Z axis. The solid has a specified half-length along the Z axis, about which it is centered, and a given minimum and maximum radii. More...

#include <G4Hype.hh>

Inheritance diagram for G4Hype:

Public Member Functions

 G4Hype (const G4String &pName, G4double newInnerRadius, G4double newOuterRadius, G4double newInnerStereo, G4double newOuterStereo, G4double newHalfLenZ)
 ~G4Hype () override
G4double GetInnerRadius () const
G4double GetOuterRadius () const
G4double GetZHalfLength () const
G4double GetInnerStereo () const
G4double GetOuterStereo () const
void SetInnerRadius (G4double newIRad)
void SetOuterRadius (G4double newORad)
void SetZHalfLength (G4double newHLZ)
void SetInnerStereo (G4double newISte)
void SetOuterStereo (G4double newOSte)
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const override
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const override
EInside Inside (const G4ThreeVector &p) const override
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const override
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const override
G4double DistanceToIn (const G4ThreeVector &p) const override
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
G4double DistanceToOut (const G4ThreeVector &p) const override
G4GeometryType GetEntityType () const override
G4VSolidClone () const override
std::ostream & StreamInfo (std::ostream &os) const override
G4double GetCubicVolume () override
G4double GetSurfaceArea () override
G4ThreeVector GetPointOnSurface () const override
void DescribeYourselfTo (G4VGraphicsScene &scene) const override
G4VisExtent GetExtent () const override
G4PolyhedronCreatePolyhedron () const override
G4PolyhedronGetPolyhedron () const override
 G4Hype (__void__ &)
 G4Hype (const G4Hype &rhs)
G4Hypeoperator= (const G4Hype &rhs)
Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
virtual ~G4VSolid ()
 G4VSolid (const G4VSolid &rhs)
G4VSolidoperator= (const G4VSolid &rhs)
G4bool operator== (const G4VSolid &s) const
G4String GetName () const
void SetName (const G4String &name)
G4double GetTolerance () const
virtual G4int GetNumOfConstituents () const
virtual G4bool IsFaceted () const
void DumpInfo () const
virtual const G4VSolidGetConstituentSolid (G4int no) const
virtual G4VSolidGetConstituentSolid (G4int no)
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 G4VSolid (__void__ &)
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
G4double EstimateSurfaceArea (G4int nStat, G4double epsilon) const

Additional Inherited Members

Protected Member Functions inherited from G4VSolid
void CalculateClippedPolygonExtent (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
void ClipCrossSection (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
void ClipBetweenSections (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
void ClipPolygon (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
Protected Attributes inherited from G4VSolid
G4double kCarTolerance

Detailed Description

G4Hype is a tube with hyperbolic profile; it describes an hyperbolic volume with curved sides parallel to the Z axis. The solid has a specified half-length along the Z axis, about which it is centered, and a given minimum and maximum radii.

Definition at line 71 of file G4Hype.hh.

Constructor & Destructor Documentation

◆ G4Hype() [1/3]

G4Hype::G4Hype ( const G4String & pName,
G4double newInnerRadius,
G4double newOuterRadius,
G4double newInnerStereo,
G4double newOuterStereo,
G4double newHalfLenZ )

Constructs a hyperbolic tube, given its parameters.

Parameters
[in]pNameThe solid name.
[in]newInnerRadiusInner radius.
[in]newOuterRadiusOuter radius.
[in]newInnerStereoInner stereo angle in radians.
[in]newOuterStereoOuter stereo angle in radians.
[in]newHalfLenZHalf length in Z.

Definition at line 61 of file G4Hype.cc.

67 : G4VSolid(pName)
68{
69 fHalfTol = 0.5*kCarTolerance;
70
71 // Check z-len
72 //
73 if (newHalfLenZ <= 0)
74 {
75 std::ostringstream message;
76 message << "Invalid Z half-length in solid: " << GetName() << " !"
77 << "\nZ half-length: " << newHalfLenZ/mm << " mm";
78 G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
79 FatalErrorInArgument, message);
80 }
81 halfLenZ = newHalfLenZ;
82
83 // Check radii
84 //
85 if (newInnerRadius < 0 || newOuterRadius < 0 || newInnerRadius >= newOuterRadius)
86 {
87 std::ostringstream message;
88 message << "Invalid radii in solid: " << GetName() << " !"
89 << "\nInner radius: " << newInnerRadius/mm << " mm"
90 << "\nOuter radius: " << newOuterRadius/mm << " mm";
91 G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
92 FatalErrorInArgument, message);
93 }
94
95 innerRadius = newInnerRadius;
96 outerRadius = newOuterRadius;
97
98 innerRadius2 = innerRadius*innerRadius;
99 outerRadius2 = outerRadius*outerRadius;
100
101 SetInnerStereo(newInnerStereo);
102 SetOuterStereo(newOuterStereo);
103
104 // Check end radii
105 //
106 if (endInnerRadius > endOuterRadius)
107 {
108 std::ostringstream message;
109 message << "Inconsistent stereo angles in solid: " << GetName() << " !"
110 << "\nEnd inner radius: " << endInnerRadius/mm << " mm"
111 << "\nEnd outer radius: " << endOuterRadius/mm << " mm";
112 G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
113 FatalErrorInArgument, message);
114 }
115}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
void SetOuterStereo(G4double newOSte)
Definition G4Hype.cc:249
void SetInnerStereo(G4double newISte)
Definition G4Hype.cc:233
G4String GetName() const
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59
G4double kCarTolerance
Definition G4VSolid.hh:418

Referenced by Clone(), G4Hype(), GetOuterStereo(), and operator=().

◆ ~G4Hype()

G4Hype::~G4Hype ( )
override

Destructor.

Definition at line 130 of file G4Hype.cc.

131{
132 delete fpPolyhedron; fpPolyhedron = nullptr;
133}

◆ G4Hype() [2/3]

G4Hype::G4Hype ( __void__ & a)

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

Definition at line 120 of file G4Hype.cc.

121 : G4VSolid(a), innerRadius(0.), outerRadius(0.), halfLenZ(0.), innerStereo(0.),
122 outerStereo(0.), tanInnerStereo(0.), tanOuterStereo(0.), tanInnerStereo2(0.),
123 tanOuterStereo2(0.), innerRadius2(0.), outerRadius2(0.), endInnerRadius2(0.),
124 endOuterRadius2(0.), endInnerRadius(0.), endOuterRadius(0.), fHalfTol(0.)
125{
126}

◆ G4Hype() [3/3]

G4Hype::G4Hype ( const G4Hype & rhs)

Copy constructor and assignment operator.

Definition at line 137 of file G4Hype.cc.

138 : G4VSolid(rhs), innerRadius(rhs.innerRadius),
139 outerRadius(rhs.outerRadius), halfLenZ(rhs.halfLenZ),
140 innerStereo(rhs.innerStereo), outerStereo(rhs.outerStereo),
141 tanInnerStereo(rhs.tanInnerStereo), tanOuterStereo(rhs.tanOuterStereo),
142 tanInnerStereo2(rhs.tanInnerStereo2), tanOuterStereo2(rhs.tanOuterStereo2),
143 innerRadius2(rhs.innerRadius2), outerRadius2(rhs.outerRadius2),
144 endInnerRadius2(rhs.endInnerRadius2), endOuterRadius2(rhs.endOuterRadius2),
145 endInnerRadius(rhs.endInnerRadius), endOuterRadius(rhs.endOuterRadius),
146 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
147 fInnerSurfaceArea(rhs.fInnerSurfaceArea), fOuterSurfaceArea(rhs.fOuterSurfaceArea),
148 fHalfTol(rhs.fHalfTol)
149{
150}

Member Function Documentation

◆ BoundingLimits()

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

Computes the bounding limits of the solid.

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

Reimplemented from G4VSolid.

Definition at line 275 of file G4Hype.cc.

276{
277 pMin.set(-endOuterRadius,-endOuterRadius,-halfLenZ);
278 pMax.set( endOuterRadius, endOuterRadius, halfLenZ);
279
280 // Check correctness of the bounding box
281 //
282 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
283 {
284 std::ostringstream message;
285 message << "Bad bounding box (min >= max) for solid: " << GetName() << " !"
286 << "\npMin = " << pMin
287 << "\npMax = " << pMax;
288 G4Exception("G4Hype::BoundingLimits()", "GeomMgt0001",
289 JustWarning, message);
290 DumpInfo();
291 }
292}
@ JustWarning
double z() const
double x() const
double y() const
void set(double x, double y, double z)
void DumpInfo() const

Referenced by CalculateExtent().

◆ CalculateExtent()

G4bool G4Hype::CalculateExtent ( const EAxis pAxis,
const G4VoxelLimits & pVoxelLimit,
const G4AffineTransform & pTransform,
G4double & pMin,
G4double & pMax ) const
overridevirtual

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

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

Implements G4VSolid.

Definition at line 296 of file G4Hype.cc.

300{
301 G4ThreeVector bmin, bmax;
302
303 // Get bounding box
304 BoundingLimits(bmin,bmax);
305
306 // Find extent
307 G4BoundingEnvelope bbox(bmin,bmax);
308 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
309}
CLHEP::Hep3Vector G4ThreeVector
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
Definition G4Hype.cc:275

◆ Clone()

G4VSolid * G4Hype::Clone ( ) const
overridevirtual

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

Returns
A pointer to the new cloned allocated solid.

Reimplemented from G4VSolid.

Definition at line 1176 of file G4Hype.cc.

1177{
1178 return new G4Hype(*this);
1179}
G4Hype(const G4String &pName, G4double newInnerRadius, G4double newOuterRadius, G4double newInnerStereo, G4double newOuterStereo, G4double newHalfLenZ)
Definition G4Hype.cc:61

◆ ComputeDimensions()

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

Dispatch method for parameterisation replication mechanism and dimension computation.

Reimplemented from G4VSolid.

Definition at line 266 of file G4Hype.cc.

269{
270 p->ComputeDimensions(*this,n,pRep);
271}
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const

◆ CreatePolyhedron()

G4Polyhedron * G4Hype::CreatePolyhedron ( ) const
overridevirtual

Creates a Polyhedron used for Visualisation. It is the caller's responsibility to delete it. A null pointer means "not created".

Reimplemented from G4VSolid.

Definition at line 1290 of file G4Hype.cc.

1291{
1292 return new G4PolyhedronHype(innerRadius, outerRadius,
1293 tanInnerStereo2, tanOuterStereo2, halfLenZ);
1294}

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

void G4Hype::DescribeYourselfTo ( G4VGraphicsScene & scene) const
overridevirtual

Methods for creating graphical representations (i.e. for visualisation).

Implements G4VSolid.

Definition at line 1272 of file G4Hype.cc.

1273{
1274 scene.AddSolid (*this);
1275}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

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

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

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

Implements G4VSolid.

Definition at line 704 of file G4Hype.cc.

705{
706 G4double absZ(std::fabs(p.z()));
707
708 //
709 // Check region
710 //
711 G4double r2 = p.x()*p.x() + p.y()*p.y();
712 G4double r = std::sqrt(r2);
713
714 G4double sigz = absZ - halfLenZ;
715
716 if (r < endOuterRadius)
717 {
718 if (sigz > -fHalfTol)
719 {
720 if (InnerSurfaceExists())
721 {
722 if (r > endInnerRadius)
723 {
724 return sigz < fHalfTol ? 0 : sigz; // Region 1
725 }
726
727 G4double dr = endInnerRadius - r;
728 if (sigz > dr*tanInnerStereo2)
729 {
730 //
731 // In region 5
732 //
733 G4double answer = std::sqrt( dr*dr + sigz*sigz );
734 return answer < fHalfTol ? 0 : answer;
735 }
736 }
737 else
738 {
739 //
740 // In region 1 (no inner surface)
741 //
742 return sigz < fHalfTol ? 0 : sigz;
743 }
744 }
745 }
746 else
747 {
748 G4double dr = r - endOuterRadius;
749 if (sigz > -dr*tanOuterStereo2)
750 {
751 //
752 // In region 2
753 //
754 G4double answer = std::sqrt( dr*dr + sigz*sigz );
755 return answer < fHalfTol ? 0 : answer;
756 }
757 }
758
759 if (InnerSurfaceExists())
760 {
761 if (r2 < HypeInnerRadius2(absZ)+kCarTolerance*endInnerRadius)
762 {
763 //
764 // In region 4
765 //
766 G4double answer = ApproxDistInside( r,absZ,innerRadius,tanInnerStereo2 );
767 return answer < fHalfTol ? 0 : answer;
768 }
769 }
770
771 //
772 // We are left by elimination with region 3
773 //
774 G4double answer = ApproxDistOutside( r, absZ, outerRadius, tanOuterStereo );
775 return answer < fHalfTol ? 0 : answer;
776}
double G4double
Definition G4Types.hh:83

◆ DistanceToIn() [2/2]

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

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

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

Implements G4VSolid.

Definition at line 401 of file G4Hype.cc.

403{
404 //
405 // Quick test. Beware! This assumes v is a unit vector!
406 //
407 if (std::fabs(p.x()*v.y() - p.y()*v.x()) > endOuterRadius+kCarTolerance)
408 {
409 return kInfinity;
410 }
411
412 //
413 // Take advantage of z symmetry, and reflect throught the
414 // z=0 plane so that pz is always positive
415 //
416 G4double pz(p.z()), vz(v.z());
417 if (pz < 0)
418 {
419 pz = -pz;
420 vz = -vz;
421 }
422
423 //
424 // We must be very careful if we don't want to
425 // create subtle leaks at the edges where the
426 // hyperbolic surfaces connect to the endplate.
427 // The only reliable way to do so is to make sure
428 // that the decision as to when a track passes
429 // over the edge of one surface is exactly the
430 // same decision as to when a track passes into the
431 // other surface. By "exact", we don't mean algebraicly
432 // exact, but we mean the same machine instructions
433 // should be used.
434 //
435 G4bool couldMissOuter(true),
436 couldMissInner(true),
437 cantMissInnerCylinder(false);
438
439 //
440 // Check endplate intersection
441 //
442 G4double sigz = pz-halfLenZ;
443
444 if (sigz > -fHalfTol) // equivalent to: if (pz > halfLenZ - fHalfTol)
445 {
446 //
447 // We start in front of the endplate (within roundoff)
448 // Correct direction to intersect endplate?
449 //
450 if (vz >= 0)
451 {
452 //
453 // Nope. As long as we are far enough away, we
454 // can't intersect anything
455 //
456 if (sigz > 0) { return kInfinity; }
457
458 //
459 // Otherwise, we may still hit a hyperbolic surface
460 // if the point is on the hyperbolic surface (within tolerance)
461 //
462 G4double pr2 = p.x()*p.x() + p.y()*p.y();
463 if (pr2 > endOuterRadius2 + kCarTolerance*endOuterRadius)
464 {
465 return kInfinity;
466 }
467
468 if (InnerSurfaceExists())
469 {
470 if (pr2 < endInnerRadius2 - kCarTolerance*endInnerRadius)
471 {
472 return kInfinity;
473 }
474 if ( (pr2 < endOuterRadius2 - kCarTolerance*endOuterRadius)
475 && (pr2 > endInnerRadius2 + kCarTolerance*endInnerRadius) )
476 {
477 return kInfinity;
478 }
479 }
480 else
481 {
482 if (pr2 < endOuterRadius2 - kCarTolerance*endOuterRadius)
483 {
484 return kInfinity;
485 }
486 }
487 }
488 else
489 {
490 //
491 // Where do we intersect at z = halfLenZ?
492 //
493 G4double q = -sigz/vz;
494 G4double xi = p.x() + q*v.x(),
495 yi = p.y() + q*v.y();
496
497 //
498 // Is this on the endplate? If so, return s, unless
499 // we are on the tolerant surface, in which case return 0
500 //
501 G4double pr2 = xi*xi + yi*yi;
502 if (pr2 <= endOuterRadius2)
503 {
504 if (InnerSurfaceExists())
505 {
506 if (pr2 >= endInnerRadius2) { return (sigz < fHalfTol) ? 0 : q; }
507 //
508 // This test is sufficient to ensure that the
509 // trajectory cannot miss the inner hyperbolic surface
510 // for z > 0, if the normal is correct.
511 //
512 G4double dot1 = (xi*v.x() + yi*v.y())*endInnerRadius/std::sqrt(pr2);
513 couldMissInner = (dot1 - halfLenZ*tanInnerStereo2*vz <= 0);
514
515 if (pr2 > endInnerRadius2*(1 - 2*DBL_EPSILON) )
516 {
517 //
518 // There is a potential leak if the inner
519 // surface is a cylinder
520 //
521 if ( (innerStereo < DBL_MIN)
522 && ((std::fabs(v.x()) > DBL_MIN) || (std::fabs(v.y()) > DBL_MIN)))
523 {
524 cantMissInnerCylinder = true;
525 }
526 }
527 }
528 else
529 {
530 return (sigz < fHalfTol) ? 0 : q;
531 }
532 }
533 else
534 {
535 G4double dotR( xi*v.x() + yi*v.y() );
536 if (dotR >= 0)
537 {
538 //
539 // Otherwise, if we are traveling outwards, we know
540 // we must miss the hyperbolic surfaces also, so
541 // we need not bother checking
542 //
543 return kInfinity;
544 }
545
546 //
547 // This test is sufficient to ensure that the
548 // trajectory cannot miss the outer hyperbolic surface
549 // for z > 0, if the normal is correct.
550 //
551 G4double dot1 = dotR*endOuterRadius/std::sqrt(pr2);
552 couldMissOuter = (dot1 - halfLenZ*tanOuterStereo2*vz>= 0);
553 }
554 }
555 }
556
557 //
558 // Check intersection with outer hyperbolic surface, save
559 // distance to valid intersection into "best".
560 //
561 G4double best = kInfinity;
562
563 G4double q[2];
564 G4int n = IntersectHype( p, v, outerRadius2, tanOuterStereo2, q );
565
566 if (n > 0)
567 {
568 //
569 // Potential intersection: is p on this surface?
570 //
571 if (pz < halfLenZ+fHalfTol)
572 {
573 G4double dr2 = p.x()*p.x() + p.y()*p.y() - HypeOuterRadius2(pz);
574 if (std::fabs(dr2) < kCarTolerance*endOuterRadius)
575 {
576 //
577 // Sure, but make sure we're traveling inwards at
578 // this point
579 //
580 if (p.x()*v.x() + p.y()*v.y() - pz*tanOuterStereo2*vz < 0) { return 0; }
581 }
582 }
583
584 //
585 // We are now certain that p is not on the tolerant surface.
586 // Accept only position distance q
587 //
588 for( G4int i=0; i<n; ++i )
589 {
590 if (q[i] >= 0)
591 {
592 //
593 // Check to make sure this intersection point is
594 // on the surface, but only do so if we haven't
595 // checked the endplate intersection already
596 //
597 G4double zi = pz + q[i]*vz;
598
599 if (zi < -halfLenZ) { continue; }
600 if (zi > +halfLenZ && couldMissOuter) { continue; }
601
602 //
603 // Check normal
604 //
605 G4double xi = p.x() + q[i]*v.x(),
606 yi = p.y() + q[i]*v.y();
607
608 if (xi*v.x() + yi*v.y() - zi*tanOuterStereo2*vz > 0) { continue; }
609
610 best = q[i];
611 break;
612 }
613 }
614 }
615
616 if (!InnerSurfaceExists()) { return best; }
617
618 //
619 // Check intersection with inner hyperbolic surface
620 //
621 n = IntersectHype( p, v, innerRadius2, tanInnerStereo2, q );
622 if (n == 0)
623 {
624 if (cantMissInnerCylinder)
625 {
626 return (sigz < fHalfTol) ? 0 : -sigz/vz;
627 }
628
629 return best;
630 }
631
632 //
633 // P on this surface?
634 //
635 if (pz < halfLenZ+fHalfTol)
636 {
637 G4double dr2 = p.x()*p.x() + p.y()*p.y() - HypeInnerRadius2(pz);
638 if (std::fabs(dr2) < kCarTolerance*endInnerRadius)
639 {
640 //
641 // Sure, but make sure we're traveling outwards at
642 // this point
643 //
644 if (p.x()*v.x() + p.y()*v.y() - pz*tanInnerStereo2*vz > 0) { return 0; }
645 }
646 }
647
648 //
649 // No, so only positive q is valid. Search for a valid intersection
650 // that is closer than the outer intersection (if it exists)
651 //
652 for( G4int i=0; i<n; ++i )
653 {
654 if (q[i] > best) { break; }
655 if (q[i] >= 0)
656 {
657 //
658 // Check to make sure this intersection point is
659 // on the surface, but only do so if we haven't
660 // checked the endplate intersection already
661 //
662 G4double zi = pz + q[i]*vz;
663
664 if (zi < -halfLenZ) { continue; }
665 if (zi > +halfLenZ && couldMissInner) { continue; }
666
667 //
668 // Check normal
669 //
670 G4double xi = p.x() + q[i]*v.x(),
671 yi = p.y() + q[i]*v.y();
672
673 if (xi*v.x() + yi*v.y() - zi*tanOuterStereo2*vz < 0) { continue; }
674
675 best = q[i];
676 break;
677 }
678 }
679
680 //
681 // Done
682 //
683 return best;
684}
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define DBL_MIN
Definition templates.hh:54
#define DBL_EPSILON
Definition templates.hh:66

◆ DistanceToOut() [1/2]

G4double G4Hype::DistanceToOut ( const G4ThreeVector & p) const
overridevirtual

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

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

Implements G4VSolid.

Definition at line 958 of file G4Hype.cc.

959{
960 //
961 // Try each surface and remember the closest
962 //
963 G4double absZ(std::fabs(p.z()));
964 G4double r(p.perp());
965
966 G4double sBest = halfLenZ - absZ;
967
968 G4double tryOuter = ApproxDistInside( r, absZ, outerRadius, tanOuterStereo2 );
969 if (tryOuter < sBest) { sBest = tryOuter; }
970
971 if (InnerSurfaceExists())
972 {
973 G4double tryInner = ApproxDistOutside( r,absZ,innerRadius,tanInnerStereo );
974 if (tryInner < sBest) { sBest = tryInner; }
975 }
976
977 return sBest < 0.5*kCarTolerance ? 0 : sBest;
978}
double perp() const

◆ DistanceToOut() [2/2]

G4double G4Hype::DistanceToOut ( const G4ThreeVector & p,
const G4ThreeVector & v,
const G4bool calcNorm = false,
G4bool * validNorm = nullptr,
G4ThreeVector * n = nullptr ) const
overridevirtual

Returns the distance along the normalised vector 'v' to the shape, from a point at an offset 'p' inside or on the surface of the shape. Intersections with surfaces, when the point is less than Tolerance/2 from a surface must be ignored.

Parameters
[in]pThe point at offset p.
[in]vThe normalised direction vector.
[in]calcNormFlag to indicate if to calculate the normal or not.
[out]validNormFlag set to true if the solid lies entirely behind or on the exiting surface. It is set false if the solid does not lie entirely behind or on the exiting surface. 'calcNorm' must be true, otherwise it is unused.
[out]nThe exiting outwards normal vector (undefined Magnitude). 'calcNorm' must be true, otherwise it is unused.
Returns
The distance to exit the shape.

Implements G4VSolid.

Definition at line 785 of file G4Hype.cc.

788{
789 static const G4ThreeVector normEnd1(0.0,0.0,+1.0);
790 static const G4ThreeVector normEnd2(0.0,0.0,-1.0);
791
792 //
793 // Keep track of closest surface
794 //
795 G4double sBest; // distance to
796 const G4ThreeVector* nBest; // normal vector
797 G4bool vBest; // whether "valid"
798
799 //
800 // Check endplate, taking advantage of symmetry.
801 // Note that the endcap is the only surface which
802 // has a "valid" normal, i.e. is a surface of which
803 // the entire solid is behind.
804 //
805 G4double pz(p.z()), vz(v.z());
806 if (vz < 0)
807 {
808 pz = -pz;
809 vz = -vz;
810 nBest = &normEnd2;
811 }
812 else
813 {
814 nBest = &normEnd1;
815 }
816
817 //
818 // Possible intercept. Are we on the surface?
819 //
820 if (pz > halfLenZ-fHalfTol)
821 {
822 if (calcNorm) { *norm = *nBest; *validNorm = true; }
823 return 0;
824 }
825
826 //
827 // Nope. Get distance. Beware of zero vz.
828 //
829 sBest = (vz > DBL_MIN) ? (halfLenZ - pz)/vz : kInfinity;
830 vBest = true;
831
832 //
833 // Check outer surface
834 //
835 G4double r2 = p.x()*p.x() + p.y()*p.y();
836
837 G4double q[2];
838 G4int n = IntersectHype( p, v, outerRadius2, tanOuterStereo2, q );
839
840 G4ThreeVector norm1, norm2;
841
842 if (n > 0)
843 {
844 //
845 // We hit somewhere. Are we on the surface?
846 //
847 G4double dr2 = r2 - HypeOuterRadius2(pz);
848 if (std::fabs(dr2) < endOuterRadius*kCarTolerance)
849 {
850 G4ThreeVector normHere( p.x(), p.y(), -p.z()*tanOuterStereo2 );
851 //
852 // Sure. But are we going the right way?
853 //
854 if (normHere.dot(v) > 0)
855 {
856 if (calcNorm) { *norm = normHere.unit(); *validNorm = false; }
857 return 0;
858 }
859 }
860
861 //
862 // Nope. Check closest positive intercept.
863 //
864 for( G4int i=0; i<n; ++i )
865 {
866 if (q[i] > sBest) { break; }
867 if (q[i] > 0)
868 {
869 //
870 // Make sure normal is correct (that this
871 // solution is an outgoing solution)
872 //
873 G4ThreeVector pk(p+q[i]*v);
874 norm1 = G4ThreeVector( pk.x(), pk.y(), -pk.z()*tanOuterStereo2 );
875 if (norm1.dot(v) > 0)
876 {
877 sBest = q[i];
878 nBest = &norm1;
879 vBest = false;
880 break;
881 }
882 }
883 }
884 }
885
886 if (InnerSurfaceExists())
887 {
888 //
889 // Check inner surface
890 //
891 n = IntersectHype( p, v, innerRadius2, tanInnerStereo2, q );
892 if (n > 0)
893 {
894 //
895 // On surface?
896 //
897 G4double dr2 = r2 - HypeInnerRadius2(pz);
898 if (std::fabs(dr2) < endInnerRadius*kCarTolerance)
899 {
900 G4ThreeVector normHere( -p.x(), -p.y(), p.z()*tanInnerStereo2 );
901 if (normHere.dot(v) > 0)
902 {
903 if (calcNorm)
904 {
905 *norm = normHere.unit();
906 *validNorm = false;
907 }
908 return 0;
909 }
910 }
911
912 //
913 // Check closest positive
914 //
915 for( G4int i=0; i<n; ++i )
916 {
917 if (q[i] > sBest) { break; }
918 if (q[i] > 0)
919 {
920 G4ThreeVector pk(p+q[i]*v);
921 norm2 = G4ThreeVector( -pk.x(), -pk.y(), pk.z()*tanInnerStereo2 );
922 if (norm2.dot(v) > 0)
923 {
924 sBest = q[i];
925 nBest = &norm2;
926 vBest = false;
927 break;
928 }
929 }
930 }
931 }
932 }
933
934 //
935 // Done!
936 //
937 if (calcNorm)
938 {
939 *validNorm = vBest;
940
941 if (nBest == &norm1 || nBest == &norm2)
942 {
943 *norm = nBest->unit();
944 }
945 else
946 {
947 *norm = *nBest;
948 }
949 }
950
951 return sBest;
952}
Hep3Vector unit() const
double dot(const Hep3Vector &) const

◆ GetCubicVolume()

G4double G4Hype::GetCubicVolume ( )
overridevirtual

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

Reimplemented from G4VSolid.

Definition at line 1245 of file G4Hype.cc.

1246{
1247 if (fCubicVolume == 0)
1248 {
1249 G4AutoLock l(&hypeMutex);
1250 fCubicVolume =
1251 twopi*halfLenZ*(2.*(outerRadius2 - innerRadius2) + endOuterRadius2 - endInnerRadius2)/3.;
1252 l.unlock();
1253 }
1254 return fCubicVolume;
1255}
G4TemplateAutoLock< G4Mutex > G4AutoLock

◆ GetEntityType()

G4GeometryType G4Hype::GetEntityType ( ) const
overridevirtual

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

Implements G4VSolid.

Definition at line 1169 of file G4Hype.cc.

1170{
1171 return {"G4Hype"};
1172}

◆ GetExtent()

G4VisExtent G4Hype::GetExtent ( ) const
overridevirtual

Provides extent (bounding box) as possible hint to the graphics view.

Reimplemented from G4VSolid.

Definition at line 1279 of file G4Hype.cc.

1280{
1281 // Define the sides of the box into which the G4Tubs instance would fit.
1282 //
1283 return { -endOuterRadius, endOuterRadius,
1284 -endOuterRadius, endOuterRadius,
1285 -halfLenZ, halfLenZ };
1286}

◆ GetInnerRadius()

◆ GetInnerStereo()

◆ GetOuterRadius()

◆ GetOuterStereo()

◆ GetPointOnSurface()

G4ThreeVector G4Hype::GetPointOnSurface ( ) const
overridevirtual

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

Reimplemented from G4VSolid.

Definition at line 1204 of file G4Hype.cc.

1205{
1206 G4double sbases = twopi*(endOuterRadius2 - endInnerRadius2);
1207 G4double stotal = fInnerSurfaceArea + fOuterSurfaceArea + sbases;
1208 G4double select = stotal*G4QuickRand();
1209
1210 G4double h = halfLenZ;
1211 G4double phi = twopi*G4QuickRand();
1212 G4double cosphi = std::cos(phi);
1213 G4double sinphi = std::sin(phi);
1214 if (select < sbases)
1215 {
1216 // circular bases
1217 G4double u = G4QuickRand();
1218 G4double rho = std::sqrt(u*endOuterRadius2 + (1. - u)*endInnerRadius2);
1219 G4double z = (select < 0.5*sbases) ? -h : h;
1220 return { rho*cosphi, rho*sinphi, z };
1221 }
1222
1223 // lateral surfaces (rejection sampling)
1224 G4double hh = h*h;
1225 G4double rr = (select < stotal - fInnerSurfaceArea) ? outerRadius2 : innerRadius2;
1226 G4double endrr = (select < stotal - fInnerSurfaceArea) ? endOuterRadius2 : endInnerRadius2;
1227 G4double delrr = endrr - rr;
1228 // max surface element: std::sqrt(aa*hh + (RR - aa)*(RR - aa + hh)*1.0)
1229 G4double ss_max = rr*hh + delrr*(delrr + hh);
1230 G4double u = 0.;
1231 for (auto i = 0; i < 10000; ++i)
1232 {
1233 u = 2.*G4QuickRand() - 1.;
1234 // surface element: std::sqrt(aa*hh + (RR - aa)*(RR - aa + hh)*uu)
1235 G4double ss = rr*hh + delrr*(delrr + hh)*u*u;
1236 if (ss_max*sqr(G4QuickRand()) <= ss) { break; }
1237 }
1238 G4double z = h*u;
1239 G4double rho = std::sqrt(rr + delrr*u*u);
1240 return { rho*cosphi, rho*sinphi, z };
1241}
G4double G4QuickRand(uint32_t seed=0)
T sqr(const T &x)
Definition templates.hh:128

◆ GetPolyhedron()

G4Polyhedron * G4Hype::GetPolyhedron ( ) const
overridevirtual

Smart access function - creates on request and stores for future access. A null pointer means "not available".

Reimplemented from G4VSolid.

Definition at line 1298 of file G4Hype.cc.

1299{
1300 if (fpPolyhedron == nullptr ||
1301 fRebuildPolyhedron ||
1302 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
1303 fpPolyhedron->GetNumberOfRotationSteps())
1304 {
1305 G4AutoLock l(&polyhedronMutex);
1306 delete fpPolyhedron;
1307 fpPolyhedron = CreatePolyhedron();
1308 fRebuildPolyhedron = false;
1309 l.unlock();
1310 }
1311 return fpPolyhedron;
1312}
G4Polyhedron * CreatePolyhedron() const override
Definition G4Hype.cc:1290

◆ GetSurfaceArea()

G4double G4Hype::GetSurfaceArea ( )
overridevirtual

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

Reimplemented from G4VSolid.

Definition at line 1259 of file G4Hype.cc.

1260{
1261 if (fSurfaceArea == 0)
1262 {
1263 G4AutoLock l(&hypeMutex);
1264 fSurfaceArea = fInnerSurfaceArea + fOuterSurfaceArea + twopi*(endOuterRadius2 - endInnerRadius2);
1265 l.unlock();
1266 }
1267 return fSurfaceArea;
1268}

◆ GetZHalfLength()

◆ Inside()

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

Concrete implementations of the expected query interfaces for solids, as defined in the base class G4VSolid.

Implements G4VSolid.

Definition at line 313 of file G4Hype.cc.

314{
315 //
316 // Check z extents: are we outside?
317 //
318 const G4double absZ(std::fabs(p.z()));
319 if (absZ > halfLenZ + fHalfTol) { return kOutside; }
320
321 //
322 // Check outer radius
323 //
324 const G4double oRad2(HypeOuterRadius2(absZ));
325 const G4double xR2( p.x()*p.x()+p.y()*p.y() );
326
327 if (xR2 > oRad2 + kCarTolerance*endOuterRadius) { return kOutside; }
328
329 if (xR2 > oRad2 - kCarTolerance*endOuterRadius) { return kSurface; }
330
331 if (InnerSurfaceExists())
332 {
333 //
334 // Check inner radius
335 //
336 const G4double iRad2(HypeInnerRadius2(absZ));
337
338 if (xR2 < iRad2 - kCarTolerance*endInnerRadius) { return kOutside; }
339
340 if (xR2 < iRad2 + kCarTolerance*endInnerRadius) { return kSurface; }
341 }
342
343 //
344 // We are inside in radius, now check endplate surface
345 //
346 if (absZ > halfLenZ - fHalfTol) { return kSurface; }
347
348 return kInside;
349}
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68
@ kSurface
Definition geomdefs.hh:69

◆ operator=()

G4Hype & G4Hype::operator= ( const G4Hype & rhs)

Definition at line 154 of file G4Hype.cc.

155{
156 // Check assignment to self
157 //
158 if (this == &rhs) { return *this; }
159
160 // Copy base class data
161 //
163
164 // Copy data
165 //
166 innerRadius = rhs.innerRadius; outerRadius = rhs.outerRadius;
167 halfLenZ = rhs.halfLenZ;
168 innerStereo = rhs.innerStereo; outerStereo = rhs.outerStereo;
169 tanInnerStereo = rhs.tanInnerStereo; tanOuterStereo = rhs.tanOuterStereo;
170 tanInnerStereo2 = rhs.tanInnerStereo2; tanOuterStereo2 = rhs.tanOuterStereo2;
171 innerRadius2 = rhs.innerRadius2; outerRadius2 = rhs.outerRadius2;
172 endInnerRadius2 = rhs.endInnerRadius2; endOuterRadius2 = rhs.endOuterRadius2;
173 endInnerRadius = rhs.endInnerRadius; endOuterRadius = rhs.endOuterRadius;
174 fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
175 fInnerSurfaceArea = rhs.fInnerSurfaceArea; fOuterSurfaceArea = rhs.fOuterSurfaceArea;
176 fHalfTol = rhs.fHalfTol;
177 fRebuildPolyhedron = false;
178 delete fpPolyhedron; fpPolyhedron = nullptr;
179
180 return *this;
181}
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:108

◆ SetInnerRadius()

void G4Hype::SetInnerRadius ( G4double newIRad)

Modifiers.

Definition at line 185 of file G4Hype.cc.

186{
187 innerRadius = newIRad;
188 innerRadius2 = newIRad*newIRad;
189 endInnerRadius2 = HypeInnerRadius2(halfLenZ);
190 endInnerRadius = std::sqrt(endInnerRadius2);
191 fCubicVolume = 0.;
192 fSurfaceArea = 0.;
193 fInnerSurfaceArea =
194 G4GeomTools::HyperboloidSurfaceArea(twopi, innerRadius, tanInnerStereo, -halfLenZ, halfLenZ);
195 fRebuildPolyhedron = true;
196}
static G4double HyperboloidSurfaceArea(G4double dphi, G4double r0, G4double tanstereo, G4double zmin, G4double zmax)

◆ SetInnerStereo()

void G4Hype::SetInnerStereo ( G4double newISte)

Definition at line 233 of file G4Hype.cc.

234{
235 innerStereo = std::abs(newISte);
236 tanInnerStereo = std::tan(innerStereo);
237 tanInnerStereo2 = tanInnerStereo*tanInnerStereo;
238 endInnerRadius2 = HypeInnerRadius2(halfLenZ);
239 endInnerRadius = std::sqrt(endInnerRadius2);
240 fCubicVolume = 0.;
241 fSurfaceArea = 0.;
242 fInnerSurfaceArea =
243 G4GeomTools::HyperboloidSurfaceArea(twopi, innerRadius, tanInnerStereo, -halfLenZ, halfLenZ);
244 fRebuildPolyhedron = true;
245}

Referenced by G4Hype().

◆ SetOuterRadius()

void G4Hype::SetOuterRadius ( G4double newORad)

Definition at line 200 of file G4Hype.cc.

201{
202 outerRadius = newORad;
203 outerRadius2 = newORad*newORad;
204 endOuterRadius2 = HypeOuterRadius2(halfLenZ);
205 endOuterRadius = std::sqrt(endOuterRadius2);
206 fCubicVolume = 0.;
207 fSurfaceArea = 0.;
208 fOuterSurfaceArea =
209 G4GeomTools::HyperboloidSurfaceArea(twopi, outerRadius, tanOuterStereo, -halfLenZ, halfLenZ);
210 fRebuildPolyhedron = true;
211}

◆ SetOuterStereo()

void G4Hype::SetOuterStereo ( G4double newOSte)

Definition at line 249 of file G4Hype.cc.

250{
251 outerStereo = std::abs(newOSte);
252 tanOuterStereo = std::tan(outerStereo);
253 tanOuterStereo2 = tanOuterStereo*tanOuterStereo;
254 endOuterRadius2 = HypeOuterRadius2(halfLenZ);
255 endOuterRadius = std::sqrt(endOuterRadius2);
256 fCubicVolume = 0.;
257 fSurfaceArea = 0.;
258 fOuterSurfaceArea =
259 G4GeomTools::HyperboloidSurfaceArea(twopi, outerRadius, tanOuterStereo, -halfLenZ, halfLenZ);
260 fRebuildPolyhedron = true;
261}

Referenced by G4Hype().

◆ SetZHalfLength()

void G4Hype::SetZHalfLength ( G4double newHLZ)

Definition at line 215 of file G4Hype.cc.

216{
217 halfLenZ = newHLZ;
218 endInnerRadius2 = HypeInnerRadius2(halfLenZ);
219 endInnerRadius = std::sqrt(endInnerRadius2);
220 endOuterRadius2 = HypeOuterRadius2(halfLenZ);
221 endOuterRadius = std::sqrt(endOuterRadius2);
222 fCubicVolume = 0.;
223 fSurfaceArea = 0.;
224 fInnerSurfaceArea =
225 G4GeomTools::HyperboloidSurfaceArea(twopi, innerRadius, tanInnerStereo, -halfLenZ, halfLenZ);
226 fOuterSurfaceArea =
227 G4GeomTools::HyperboloidSurfaceArea(twopi, outerRadius, tanOuterStereo, -halfLenZ, halfLenZ);
228 fRebuildPolyhedron = true;
229}

◆ StreamInfo()

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

Streams the object contents to an output stream.

Implements G4VSolid.

Definition at line 1183 of file G4Hype.cc.

1184{
1185 G4long oldprc = os.precision(16);
1186 os << "-----------------------------------------------------------\n"
1187 << " *** Dump for solid - " << GetName() << " ***\n"
1188 << " ===================================================\n"
1189 << " Solid type: G4Hype\n"
1190 << " Parameters: \n"
1191 << " half length Z: " << halfLenZ/mm << " mm \n"
1192 << " inner radius : " << innerRadius/mm << " mm \n"
1193 << " outer radius : " << outerRadius/mm << " mm \n"
1194 << " inner stereo angle : " << innerStereo/degree << " degrees \n"
1195 << " outer stereo angle : " << outerStereo/degree << " degrees \n"
1196 << "-----------------------------------------------------------\n";
1197 os.precision(oldprc);
1198
1199 return os;
1200}
long G4long
Definition G4Types.hh:87

◆ SurfaceNormal()

G4ThreeVector G4Hype::SurfaceNormal ( const G4ThreeVector & p) const
overridevirtual

Returns the outwards pointing unit normal of the shape for the surface closest to the point at offset 'p'.

Parameters
[in]pThe point at offset p.
Returns
The outwards pointing unit normal.

Implements G4VSolid.

Definition at line 354 of file G4Hype.cc.

355{
356 //
357 // Which of the three or four surfaces are we closest to?
358 //
359 const G4double absZ(std::fabs(p.z()));
360 const G4double distZ(absZ - halfLenZ);
361 const G4double dist2Z(distZ*distZ);
362
363 const G4double xR2( p.x()*p.x()+p.y()*p.y() );
364 const G4double dist2Outer( std::fabs(xR2 - HypeOuterRadius2(absZ)) );
365
366 if (InnerSurfaceExists())
367 {
368 //
369 // Has inner surface: is this closest?
370 //
371 const G4double dist2Inner( std::fabs(xR2 - HypeInnerRadius2(absZ)) );
372 if (dist2Inner < dist2Z && dist2Inner < dist2Outer)
373 {
374 return G4ThreeVector( -p.x(), -p.y(), p.z()*tanInnerStereo2 ).unit();
375 }
376 }
377
378 //
379 // Do the "endcaps" win?
380 //
381 if (dist2Z < dist2Outer)
382 {
383 return { 0.0, 0.0, (G4double)(p.z() < 0 ? -1.0 : 1.0) };
384 }
385
386 //
387 // Outer surface wins
388 //
389 return G4ThreeVector( p.x(), p.y(), -p.z()*tanOuterStereo2 ).unit();
390}

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