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

#include <G4Hype.hh>

+ Inheritance diagram for G4Hype:

Public Member Functions

 G4Hype (const G4String &pName, G4double newInnerRadius, G4double newOuterRadius, G4double newInnerStereo, G4double newOuterStereo, G4double newHalfLenZ)
 
virtual ~G4Hype ()
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
 
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)
 
EInside Inside (const G4ThreeVector &p) const
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
G4double DistanceToIn (const G4ThreeVector &p) const
 
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
 
G4double DistanceToOut (const G4ThreeVector &p) const
 
G4GeometryType GetEntityType () const
 
G4VSolidClone () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4double GetCubicVolume ()
 
G4double GetSurfaceArea ()
 
G4ThreeVector GetPointOnSurface () const
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
G4VisExtent GetExtent () const
 
G4PolyhedronCreatePolyhedron () const
 
G4PolyhedronGetPolyhedron () const
 
 G4Hype (__void__ &)
 
 G4Hype (const G4Hype &rhs)
 
G4Hypeoperator= (const G4Hype &rhs)
 
- 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 Types

enum  ESide { outerFace , innerFace , leftCap , rightCap }
 

Protected Member Functions

G4bool InnerSurfaceExists () const
 
G4double HypeInnerRadius2 (G4double zVal) const
 
G4double HypeOuterRadius2 (G4double zVal) const
 
- 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
 

Static Protected Member Functions

static G4double ApproxDistOutside (G4double pr, G4double pz, G4double r0, G4double tanPhi)
 
static G4double ApproxDistInside (G4double pr, G4double pz, G4double r0, G4double tan2Phi)
 
static G4int IntersectHype (const G4ThreeVector &p, const G4ThreeVector &v, G4double r2, G4double tan2Phi, G4double s[2])
 

Protected Attributes

G4double innerRadius
 
G4double outerRadius
 
G4double halfLenZ
 
G4double innerStereo
 
G4double outerStereo
 
G4double tanInnerStereo
 
G4double tanOuterStereo
 
G4double tanInnerStereo2
 
G4double tanOuterStereo2
 
G4double innerRadius2
 
G4double outerRadius2
 
G4double endInnerRadius2
 
G4double endOuterRadius2
 
G4double endInnerRadius
 
G4double endOuterRadius
 
- Protected Attributes inherited from G4VSolid
G4double kCarTolerance
 

Detailed Description

Definition at line 68 of file G4Hype.hh.

Member Enumeration Documentation

◆ ESide

enum G4Hype::ESide
protected
Enumerator
outerFace 
innerFace 
leftCap 
rightCap 

Definition at line 185 of file G4Hype.hh.

@ innerFace
Definition: G4Hype.hh:185
@ leftCap
Definition: G4Hype.hh:185
@ rightCap
Definition: G4Hype.hh:185
@ outerFace
Definition: G4Hype.hh:185

Constructor & Destructor Documentation

◆ G4Hype() [1/3]

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

Definition at line 65 of file G4Hype.cc.

71 : G4VSolid(pName)
72{
73 fHalfTol = 0.5*kCarTolerance;
74
75 // Check z-len
76 //
77 if (newHalfLenZ<=0)
78 {
79 std::ostringstream message;
80 message << "Invalid Z half-length - " << GetName() << G4endl
81 << " Invalid Z half-length: "
82 << newHalfLenZ/mm << " mm";
83 G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
84 FatalErrorInArgument, message);
85 }
86 halfLenZ=newHalfLenZ;
87
88 // Check radii
89 //
90 if (newInnerRadius<0 || newOuterRadius<0)
91 {
92 std::ostringstream message;
93 message << "Invalid radii - " << GetName() << G4endl
94 << " Invalid radii ! Inner radius: "
95 << newInnerRadius/mm << " mm" << G4endl
96 << " Outer radius: "
97 << newOuterRadius/mm << " mm";
98 G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
99 FatalErrorInArgument, message);
100 }
101 if (newInnerRadius >= newOuterRadius)
102 {
103 std::ostringstream message;
104 message << "Outer > inner radius - " << GetName() << G4endl
105 << " Invalid radii ! Inner radius: "
106 << newInnerRadius/mm << " mm" << G4endl
107 << " Outer radius: "
108 << newOuterRadius/mm << " mm";
109 G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
110 FatalErrorInArgument, message);
111 }
112
113 innerRadius=newInnerRadius;
114 outerRadius=newOuterRadius;
115
118
119 SetInnerStereo( newInnerStereo );
120 SetOuterStereo( newOuterStereo );
121}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
#define G4endl
Definition: G4ios.hh:57
G4double halfLenZ
Definition: G4Hype.hh:166
G4double innerRadius2
Definition: G4Hype.hh:176
void SetOuterStereo(G4double newOSte)
G4double outerRadius
Definition: G4Hype.hh:165
G4double outerRadius2
Definition: G4Hype.hh:177
G4double innerRadius
Definition: G4Hype.hh:164
void SetInnerStereo(G4double newISte)
G4String GetName() const
G4double kCarTolerance
Definition: G4VSolid.hh:299

◆ ~G4Hype()

G4Hype::~G4Hype ( )
virtual

Definition at line 136 of file G4Hype.cc.

137{
138 delete fpPolyhedron; fpPolyhedron = 0;
139}

◆ G4Hype() [2/3]

G4Hype::G4Hype ( __void__ &  a)

Definition at line 126 of file G4Hype.cc.

127 : G4VSolid(a), innerRadius(0.), outerRadius(0.), halfLenZ(0.), innerStereo(0.),
130 endOuterRadius2(0.), endInnerRadius(0.), endOuterRadius(0.), fHalfTol(0.)
131{
132}
G4double endInnerRadius
Definition: G4Hype.hh:180
G4double tanOuterStereo
Definition: G4Hype.hh:173
G4double tanOuterStereo2
Definition: G4Hype.hh:175
G4double tanInnerStereo2
Definition: G4Hype.hh:174
G4double endInnerRadius2
Definition: G4Hype.hh:178
G4double endOuterRadius
Definition: G4Hype.hh:181
G4double outerStereo
Definition: G4Hype.hh:168
G4double tanInnerStereo
Definition: G4Hype.hh:172
G4double endOuterRadius2
Definition: G4Hype.hh:179
G4double innerStereo
Definition: G4Hype.hh:167

◆ G4Hype() [3/3]

G4Hype::G4Hype ( const G4Hype rhs)

Definition at line 143 of file G4Hype.cc.

144 : G4VSolid(rhs), innerRadius(rhs.innerRadius),
152 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
153 fHalfTol(rhs.fHalfTol)
154{
155}

Member Function Documentation

◆ ApproxDistInside()

G4double G4Hype::ApproxDistInside ( G4double  pr,
G4double  pz,
G4double  r0,
G4double  tan2Phi 
)
staticprotected

Definition at line 1043 of file G4Hype.cc.

1045{
1046 if (tan2Phi < DBL_MIN) return r0 - pr;
1047
1048 //
1049 // Corresponding position and normal on hyperbolic
1050 //
1051 G4double rh = std::sqrt( r0*r0 + pz*pz*tan2Phi );
1052
1053 G4double dr = -rh;
1054 G4double dz = pz*tan2Phi;
1055 G4double len = std::sqrt(dr*dr + dz*dz);
1056
1057 //
1058 // Answer
1059 //
1060 return std::fabs((pr-rh)*dr)/len;
1061}
double G4double
Definition: G4Types.hh:83
#define DBL_MIN
Definition: templates.hh:54

Referenced by DistanceToIn(), and DistanceToOut().

◆ ApproxDistOutside()

G4double G4Hype::ApproxDistOutside ( G4double  pr,
G4double  pz,
G4double  r0,
G4double  tanPhi 
)
staticprotected

Definition at line 986 of file G4Hype.cc.

988{
989 if (tanPhi < DBL_MIN) return pr-r0;
990
991 G4double tan2Phi = tanPhi*tanPhi;
992
993 //
994 // First point
995 //
996 G4double z1 = pz;
997 G4double r1 = std::sqrt( r0*r0 + z1*z1*tan2Phi );
998
999 //
1000 // Second point
1001 //
1002 G4double z2 = (pr*tanPhi + pz)/(1 + tan2Phi);
1003 G4double r2 = std::sqrt( r0*r0 + z2*z2*tan2Phi );
1004
1005 //
1006 // Line between them
1007 //
1008 G4double dr = r2-r1;
1009 G4double dz = z2-z1;
1010
1011 G4double len = std::sqrt(dr*dr + dz*dz);
1012 if (len < DBL_MIN)
1013 {
1014 //
1015 // The two points are the same?? I guess we
1016 // must have really bracketed the normal
1017 //
1018 dr = pr-r1;
1019 dz = pz-z1;
1020 return std::sqrt( dr*dr + dz*dz );
1021 }
1022
1023 //
1024 // Distance
1025 //
1026 return std::fabs((pr-r1)*dz - (pz-z1)*dr)/len;
1027}

Referenced by DistanceToIn(), and DistanceToOut().

◆ BoundingLimits()

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

Reimplemented from G4VSolid.

Definition at line 199 of file G4Hype.cc.

200{
203
204 // Check correctness of the bounding box
205 //
206 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
207 {
208 std::ostringstream message;
209 message << "Bad bounding box (min >= max) for solid: "
210 << GetName() << " !"
211 << "\npMin = " << pMin
212 << "\npMax = " << pMax;
213 G4Exception("G4Hype::BoundingLimits()", "GeomMgt0001",
214 JustWarning, message);
215 DumpInfo();
216 }
217}
@ 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
virtual

Implements G4VSolid.

Definition at line 221 of file G4Hype.cc.

225{
226 G4ThreeVector bmin, bmax;
227
228 // Get bounding box
229 BoundingLimits(bmin,bmax);
230
231 // Find extent
232 G4BoundingEnvelope bbox(bmin,bmax);
233 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
234}
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Hype.cc:199

◆ Clone()

G4VSolid * G4Hype::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1072 of file G4Hype.cc.

1073{
1074 return new G4Hype(*this);
1075}
Definition: G4Hype.hh:69

◆ ComputeDimensions()

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

Reimplemented from G4VSolid.

Definition at line 190 of file G4Hype.cc.

193{
194 p->ComputeDimensions(*this,n,pRep);
195}
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const

◆ CreatePolyhedron()

G4Polyhedron * G4Hype::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1269 of file G4Hype.cc.

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

void G4Hype::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 1251 of file G4Hype.cc.

1252{
1253 scene.AddSolid (*this);
1254}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 614 of file G4Hype.cc.

615{
616 G4double absZ(std::fabs(p.z()));
617
618 //
619 // Check region
620 //
621 G4double r2 = p.x()*p.x() + p.y()*p.y();
622 G4double r = std::sqrt(r2);
623
624 G4double sigz = absZ - halfLenZ;
625
626 if (r < endOuterRadius)
627 {
628 if (sigz > -fHalfTol)
629 {
630 if (InnerSurfaceExists())
631 {
632 if (r > endInnerRadius)
633 return sigz < fHalfTol ? 0 : sigz; // Region 1
634
635 G4double dr = endInnerRadius - r;
636 if (sigz > dr*tanInnerStereo2)
637 {
638 //
639 // In region 5
640 //
641 G4double answer = std::sqrt( dr*dr + sigz*sigz );
642 return answer < fHalfTol ? 0 : answer;
643 }
644 }
645 else
646 {
647 //
648 // In region 1 (no inner surface)
649 //
650 return sigz < fHalfTol ? 0 : sigz;
651 }
652 }
653 }
654 else
655 {
656 G4double dr = r - endOuterRadius;
657 if (sigz > -dr*tanOuterStereo2)
658 {
659 //
660 // In region 2
661 //
662 G4double answer = std::sqrt( dr*dr + sigz*sigz );
663 return answer < fHalfTol ? 0 : answer;
664 }
665 }
666
667 if (InnerSurfaceExists())
668 {
670 {
671 //
672 // In region 4
673 //
675 return answer < fHalfTol ? 0 : answer;
676 }
677 }
678
679 //
680 // We are left by elimination with region 3
681 //
683 return answer < fHalfTol ? 0 : answer;
684}
static G4double ApproxDistOutside(G4double pr, G4double pz, G4double r0, G4double tanPhi)
Definition: G4Hype.cc:986
G4bool InnerSurfaceExists() const
static G4double ApproxDistInside(G4double pr, G4double pz, G4double r0, G4double tan2Phi)
Definition: G4Hype.cc:1043
G4double HypeInnerRadius2(G4double zVal) const

◆ DistanceToIn() [2/2]

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

Implements G4VSolid.

Definition at line 323 of file G4Hype.cc.

325{
326 //
327 // Quick test. Beware! This assumes v is a unit vector!
328 //
329 if (std::fabs(p.x()*v.y() - p.y()*v.x()) > endOuterRadius+kCarTolerance)
330 return kInfinity;
331
332 //
333 // Take advantage of z symmetry, and reflect throught the
334 // z=0 plane so that pz is always positive
335 //
336 G4double pz(p.z()), vz(v.z());
337 if (pz < 0)
338 {
339 pz = -pz;
340 vz = -vz;
341 }
342
343 //
344 // We must be very careful if we don't want to
345 // create subtle leaks at the edges where the
346 // hyperbolic surfaces connect to the endplate.
347 // The only reliable way to do so is to make sure
348 // that the decision as to when a track passes
349 // over the edge of one surface is exactly the
350 // same decision as to when a track passes into the
351 // other surface. By "exact", we don't mean algebraicly
352 // exact, but we mean the same machine instructions
353 // should be used.
354 //
355 G4bool couldMissOuter(true),
356 couldMissInner(true),
357 cantMissInnerCylinder(false);
358
359 //
360 // Check endplate intersection
361 //
362 G4double sigz = pz-halfLenZ;
363
364 if (sigz > -fHalfTol) // equivalent to: if (pz > halfLenZ - fHalfTol)
365 {
366 //
367 // We start in front of the endplate (within roundoff)
368 // Correct direction to intersect endplate?
369 //
370 if (vz >= 0)
371 {
372 //
373 // Nope. As long as we are far enough away, we
374 // can't intersect anything
375 //
376 if (sigz > 0) return kInfinity;
377
378 //
379 // Otherwise, we may still hit a hyperbolic surface
380 // if the point is on the hyperbolic surface (within tolerance)
381 //
382 G4double pr2 = p.x()*p.x() + p.y()*p.y();
384 return kInfinity;
385
386 if (InnerSurfaceExists())
387 {
389 return kInfinity;
392 return kInfinity;
393 }
394 else
395 {
397 return kInfinity;
398 }
399 }
400 else
401 {
402 //
403 // Where do we intersect at z = halfLenZ?
404 //
405 G4double q = -sigz/vz;
406 G4double xi = p.x() + q*v.x(),
407 yi = p.y() + q*v.y();
408
409 //
410 // Is this on the endplate? If so, return s, unless
411 // we are on the tolerant surface, in which case return 0
412 //
413 G4double pr2 = xi*xi + yi*yi;
414 if (pr2 <= endOuterRadius2)
415 {
416 if (InnerSurfaceExists())
417 {
418 if (pr2 >= endInnerRadius2) return (sigz < fHalfTol) ? 0 : q;
419 //
420 // This test is sufficient to ensure that the
421 // trajectory cannot miss the inner hyperbolic surface
422 // for z > 0, if the normal is correct.
423 //
424 G4double dot1 = (xi*v.x() + yi*v.y())*endInnerRadius/std::sqrt(pr2);
425 couldMissInner = (dot1 - halfLenZ*tanInnerStereo2*vz <= 0);
426
427 if (pr2 > endInnerRadius2*(1 - 2*DBL_EPSILON) )
428 {
429 //
430 // There is a potential leak if the inner
431 // surface is a cylinder
432 //
433 if ( (innerStereo < DBL_MIN)
434 && ((std::fabs(v.x()) > DBL_MIN) || (std::fabs(v.y()) > DBL_MIN)))
435 cantMissInnerCylinder = true;
436 }
437 }
438 else
439 {
440 return (sigz < fHalfTol) ? 0 : q;
441 }
442 }
443 else
444 {
445 G4double dotR( xi*v.x() + yi*v.y() );
446 if (dotR >= 0)
447 {
448 //
449 // Otherwise, if we are traveling outwards, we know
450 // we must miss the hyperbolic surfaces also, so
451 // we need not bother checking
452 //
453 return kInfinity;
454 }
455 else
456 {
457 //
458 // This test is sufficient to ensure that the
459 // trajectory cannot miss the outer hyperbolic surface
460 // for z > 0, if the normal is correct.
461 //
462 G4double dot1 = dotR*endOuterRadius/std::sqrt(pr2);
463 couldMissOuter = (dot1 - halfLenZ*tanOuterStereo2*vz>= 0);
464 }
465 }
466 }
467 }
468
469 //
470 // Check intersection with outer hyperbolic surface, save
471 // distance to valid intersection into "best".
472 //
473 G4double best = kInfinity;
474
475 G4double q[2];
477
478 if (n > 0)
479 {
480 //
481 // Potential intersection: is p on this surface?
482 //
483 if (pz < halfLenZ+fHalfTol)
484 {
485 G4double dr2 = p.x()*p.x() + p.y()*p.y() - HypeOuterRadius2(pz);
486 if (std::fabs(dr2) < kCarTolerance*endOuterRadius)
487 {
488 //
489 // Sure, but make sure we're traveling inwards at
490 // this point
491 //
492 if (p.x()*v.x() + p.y()*v.y() - pz*tanOuterStereo2*vz < 0)
493 return 0;
494 }
495 }
496
497 //
498 // We are now certain that p is not on the tolerant surface.
499 // Accept only position distance q
500 //
501 for( G4int i=0; i<n; ++i )
502 {
503 if (q[i] >= 0)
504 {
505 //
506 // Check to make sure this intersection point is
507 // on the surface, but only do so if we haven't
508 // checked the endplate intersection already
509 //
510 G4double zi = pz + q[i]*vz;
511
512 if (zi < -halfLenZ) continue;
513 if (zi > +halfLenZ && couldMissOuter) continue;
514
515 //
516 // Check normal
517 //
518 G4double xi = p.x() + q[i]*v.x(),
519 yi = p.y() + q[i]*v.y();
520
521 if (xi*v.x() + yi*v.y() - zi*tanOuterStereo2*vz > 0) continue;
522
523 best = q[i];
524 break;
525 }
526 }
527 }
528
529 if (!InnerSurfaceExists()) return best;
530
531 //
532 // Check intersection with inner hyperbolic surface
533 //
535 if (n == 0)
536 {
537 if (cantMissInnerCylinder) return (sigz < fHalfTol) ? 0 : -sigz/vz;
538
539 return best;
540 }
541
542 //
543 // P on this surface?
544 //
545 if (pz < halfLenZ+fHalfTol)
546 {
547 G4double dr2 = p.x()*p.x() + p.y()*p.y() - HypeInnerRadius2(pz);
548 if (std::fabs(dr2) < kCarTolerance*endInnerRadius)
549 {
550 //
551 // Sure, but make sure we're traveling outwards at
552 // this point
553 //
554 if (p.x()*v.x() + p.y()*v.y() - pz*tanInnerStereo2*vz > 0) return 0;
555 }
556 }
557
558 //
559 // No, so only positive q is valid. Search for a valid intersection
560 // that is closer than the outer intersection (if it exists)
561 //
562 for( G4int i=0; i<n; ++i )
563 {
564 if (q[i] > best) break;
565 if (q[i] >= 0)
566 {
567 //
568 // Check to make sure this intersection point is
569 // on the surface, but only do so if we haven't
570 // checked the endplate intersection already
571 //
572 G4double zi = pz + q[i]*vz;
573
574 if (zi < -halfLenZ) continue;
575 if (zi > +halfLenZ && couldMissInner) continue;
576
577 //
578 // Check normal
579 //
580 G4double xi = p.x() + q[i]*v.x(),
581 yi = p.y() + q[i]*v.y();
582
583 if (xi*v.x() + yi*v.y() - zi*tanOuterStereo2*vz < 0) continue;
584
585 best = q[i];
586 break;
587 }
588 }
589
590 //
591 // Done
592 //
593 return best;
594}
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
static G4int IntersectHype(const G4ThreeVector &p, const G4ThreeVector &v, G4double r2, G4double tan2Phi, G4double s[2])
Definition: G4Hype.cc:921
G4double HypeOuterRadius2(G4double zVal) const
#define DBL_EPSILON
Definition: templates.hh:66

◆ DistanceToOut() [1/2]

G4double G4Hype::DistanceToOut ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 860 of file G4Hype.cc.

861{
862 //
863 // Try each surface and remember the closest
864 //
865 G4double absZ(std::fabs(p.z()));
866 G4double r(p.perp());
867
868 G4double sBest = halfLenZ - absZ;
869
871 if (tryOuter < sBest)
872 sBest = tryOuter;
873
874 if (InnerSurfaceExists())
875 {
877 if (tryInner < sBest) sBest = tryInner;
878 }
879
880 return sBest < 0.5*kCarTolerance ? 0 : sBest;
881}
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
virtual

Implements G4VSolid.

Definition at line 693 of file G4Hype.cc.

696{
697 static const G4ThreeVector normEnd1(0.0,0.0,+1.0);
698 static const G4ThreeVector normEnd2(0.0,0.0,-1.0);
699
700 //
701 // Keep track of closest surface
702 //
703 G4double sBest; // distance to
704 const G4ThreeVector* nBest; // normal vector
705 G4bool vBest; // whether "valid"
706
707 //
708 // Check endplate, taking advantage of symmetry.
709 // Note that the endcap is the only surface which
710 // has a "valid" normal, i.e. is a surface of which
711 // the entire solid is behind.
712 //
713 G4double pz(p.z()), vz(v.z());
714 if (vz < 0)
715 {
716 pz = -pz;
717 vz = -vz;
718 nBest = &normEnd2;
719 }
720 else
721 nBest = &normEnd1;
722
723 //
724 // Possible intercept. Are we on the surface?
725 //
726 if (pz > halfLenZ-fHalfTol)
727 {
728 if (calcNorm) { *norm = *nBest; *validNorm = true; }
729 return 0;
730 }
731
732 //
733 // Nope. Get distance. Beware of zero vz.
734 //
735 sBest = (vz > DBL_MIN) ? (halfLenZ - pz)/vz : kInfinity;
736 vBest = true;
737
738 //
739 // Check outer surface
740 //
741 G4double r2 = p.x()*p.x() + p.y()*p.y();
742
743 G4double q[2];
745
746 G4ThreeVector norm1, norm2;
747
748 if (n > 0)
749 {
750 //
751 // We hit somewhere. Are we on the surface?
752 //
753 G4double dr2 = r2 - HypeOuterRadius2(pz);
754 if (std::fabs(dr2) < endOuterRadius*kCarTolerance)
755 {
756 G4ThreeVector normHere( p.x(), p.y(), -p.z()*tanOuterStereo2 );
757 //
758 // Sure. But are we going the right way?
759 //
760 if (normHere.dot(v) > 0)
761 {
762 if (calcNorm) { *norm = normHere.unit(); *validNorm = false; }
763 return 0;
764 }
765 }
766
767 //
768 // Nope. Check closest positive intercept.
769 //
770 for( G4int i=0; i<n; ++i )
771 {
772 if (q[i] > sBest) break;
773 if (q[i] > 0)
774 {
775 //
776 // Make sure normal is correct (that this
777 // solution is an outgoing solution)
778 //
779 G4ThreeVector pk(p+q[i]*v);
780 norm1 = G4ThreeVector( pk.x(), pk.y(), -pk.z()*tanOuterStereo2 );
781 if (norm1.dot(v) > 0)
782 {
783 sBest = q[i];
784 nBest = &norm1;
785 vBest = false;
786 break;
787 }
788 }
789 }
790 }
791
792 if (InnerSurfaceExists())
793 {
794 //
795 // Check inner surface
796 //
798 if (n > 0)
799 {
800 //
801 // On surface?
802 //
803 G4double dr2 = r2 - HypeInnerRadius2(pz);
804 if (std::fabs(dr2) < endInnerRadius*kCarTolerance)
805 {
806 G4ThreeVector normHere( -p.x(), -p.y(), p.z()*tanInnerStereo2 );
807 if (normHere.dot(v) > 0)
808 {
809 if (calcNorm)
810 {
811 *norm = normHere.unit();
812 *validNorm = false;
813 }
814 return 0;
815 }
816 }
817
818 //
819 // Check closest positive
820 //
821 for( G4int i=0; i<n; ++i )
822 {
823 if (q[i] > sBest) break;
824 if (q[i] > 0)
825 {
826 G4ThreeVector pk(p+q[i]*v);
827 norm2 = G4ThreeVector( -pk.x(), -pk.y(), pk.z()*tanInnerStereo2 );
828 if (norm2.dot(v) > 0)
829 {
830 sBest = q[i];
831 nBest = &norm2;
832 vBest = false;
833 break;
834 }
835 }
836 }
837 }
838 }
839
840 //
841 // Done!
842 //
843 if (calcNorm)
844 {
845 *validNorm = vBest;
846
847 if (nBest == &norm1 || nBest == &norm2)
848 *norm = nBest->unit();
849 else
850 *norm = *nBest;
851 }
852
853 return sBest;
854}
CLHEP::Hep3Vector G4ThreeVector
Hep3Vector unit() const
double dot(const Hep3Vector &) const

◆ GetCubicVolume()

G4double G4Hype::GetCubicVolume ( )
virtual

Reimplemented from G4VSolid.

Definition at line 1081 of file G4Hype.cc.

1082{
1083 if (fCubicVolume == 0.)
1084 {
1085 fCubicVolume = CLHEP::twopi*halfLenZ*
1087 }
1088 return fCubicVolume;
1089}

◆ GetEntityType()

G4GeometryType G4Hype::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 1065 of file G4Hype.cc.

1066{
1067 return G4String("G4Hype");
1068}

◆ GetExtent()

G4VisExtent G4Hype::GetExtent ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1258 of file G4Hype.cc.

1259{
1260 // Define the sides of the box into which the G4Tubs instance would fit.
1261 //
1264 -halfLenZ, halfLenZ );
1265}

◆ GetInnerRadius()

◆ GetInnerStereo()

◆ GetOuterRadius()

◆ GetOuterStereo()

◆ GetPointOnSurface()

G4ThreeVector G4Hype::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1148 of file G4Hype.cc.

1149{
1150 G4double xRand, yRand, zRand, r2 , aOne, aTwo, aThree, chose, sinhu;
1151 G4double phi, cosphi, sinphi, rBar2Out, rBar2In, alpha, t, rOut, rIn2, rOut2;
1152
1153 // we use the formula of the area of a surface of revolution to compute
1154 // the areas, using the equation of the hyperbola:
1155 // x^2 + y^2 = (z*tanphi)^2 + r^2
1156
1157 rBar2Out = outerRadius2;
1158 alpha = 2.*pi*rBar2Out*std::cos(outerStereo)/tanOuterStereo;
1160 t = std::log(t+std::sqrt(sqr(t)+1));
1161 aOne = std::fabs(2.*alpha*(std::sinh(2.*t)/4.+t/2.));
1162
1163
1164 rBar2In = innerRadius2;
1165 alpha = 2.*pi*rBar2In*std::cos(innerStereo)/tanInnerStereo;
1167 t = std::log(t+std::sqrt(sqr(t)+1));
1168 aTwo = std::fabs(2.*alpha*(std::sinh(2.*t)/4.+t/2.));
1169
1172
1173 if(outerStereo == 0.) {aOne = std::fabs(2.*pi*outerRadius*2.*halfLenZ);}
1174 if(innerStereo == 0.) {aTwo = std::fabs(2.*pi*innerRadius*2.*halfLenZ);}
1175
1176 phi = G4RandFlat::shoot(0.,2.*pi);
1177 cosphi = std::cos(phi);
1178 sinphi = std::sin(phi);
1179 sinhu = G4RandFlat::shoot(-1.*halfLenZ*tanOuterStereo/outerRadius,
1181
1182 chose = G4RandFlat::shoot(0.,aOne+aTwo+2.*aThree);
1183 if(chose>=0. && chose < aOne)
1184 {
1185 if(outerStereo != 0.)
1186 {
1187 zRand = outerRadius*sinhu/tanOuterStereo;
1188 xRand = std::sqrt(sqr(sinhu)+1)*outerRadius*cosphi;
1189 yRand = std::sqrt(sqr(sinhu)+1)*outerRadius*sinphi;
1190 return G4ThreeVector (xRand, yRand, zRand);
1191 }
1192 else
1193 {
1194 return G4ThreeVector(outerRadius*cosphi,outerRadius*sinphi,
1195 G4RandFlat::shoot(-halfLenZ,halfLenZ));
1196 }
1197 }
1198 else if(chose>=aOne && chose<aOne+aTwo)
1199 {
1200 if(innerStereo != 0.)
1201 {
1202 sinhu = G4RandFlat::shoot(-1.*halfLenZ*tanInnerStereo/innerRadius,
1204 zRand = innerRadius*sinhu/tanInnerStereo;
1205 xRand = std::sqrt(sqr(sinhu)+1)*innerRadius*cosphi;
1206 yRand = std::sqrt(sqr(sinhu)+1)*innerRadius*sinphi;
1207 return G4ThreeVector (xRand, yRand, zRand);
1208 }
1209 else
1210 {
1211 return G4ThreeVector(innerRadius*cosphi,innerRadius*sinphi,
1212 G4RandFlat::shoot(-1.*halfLenZ,halfLenZ));
1213 }
1214 }
1215 else if(chose>=aOne+aTwo && chose<aOne+aTwo+aThree)
1216 {
1219 rOut = std::sqrt(rOut2) ;
1220
1221 do // Loop checking, 13.08.2015, G.Cosmo
1222 {
1223 xRand = G4RandFlat::shoot(-rOut,rOut) ;
1224 yRand = G4RandFlat::shoot(-rOut,rOut) ;
1225 r2 = xRand*xRand + yRand*yRand ;
1226 } while ( ! ( r2 >= rIn2 && r2 <= rOut2 ) ) ;
1227
1228 zRand = halfLenZ;
1229 return G4ThreeVector (xRand, yRand, zRand);
1230 }
1231 else
1232 {
1235 rOut = std::sqrt(rOut2) ;
1236
1237 do // Loop checking, 13.08.2015, G.Cosmo
1238 {
1239 xRand = G4RandFlat::shoot(-rOut,rOut) ;
1240 yRand = G4RandFlat::shoot(-rOut,rOut) ;
1241 r2 = xRand*xRand + yRand*yRand ;
1242 } while ( ! ( r2 >= rIn2 && r2 <= rOut2 ) ) ;
1243
1244 zRand = -1.*halfLenZ;
1245 return G4ThreeVector (xRand, yRand, zRand);
1246 }
1247}
const G4double pi
T sqr(const T &x)
Definition: templates.hh:128

◆ GetPolyhedron()

G4Polyhedron * G4Hype::GetPolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1277 of file G4Hype.cc.

1278{
1279 if (fpPolyhedron == nullptr ||
1280 fRebuildPolyhedron ||
1282 fpPolyhedron->GetNumberOfRotationSteps())
1283 {
1284 G4AutoLock l(&polyhedronMutex);
1285 delete fpPolyhedron;
1286 fpPolyhedron = CreatePolyhedron();
1287 fRebuildPolyhedron = false;
1288 l.unlock();
1289 }
1290 return fpPolyhedron;
1291}
G4Polyhedron * CreatePolyhedron() const
Definition: G4Hype.cc:1269
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
static G4int GetNumberOfRotationSteps()

◆ GetSurfaceArea()

G4double G4Hype::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 1093 of file G4Hype.cc.

1094{
1095 if (fSurfaceArea == 0.)
1096 {
1097 G4double h = halfLenZ;
1098 G4double innS = 2.*h*innerRadius;
1099 if (std::abs(endInnerRadius - innerRadius) > kCarTolerance)
1100 {
1104 G4double CC = AA*h*h/(RR - AA);
1105 G4double K = std::sqrt(AA + CC)/CC;
1106 G4double Kh = K*h;
1107 innS = A*(h*std::sqrt(1. + Kh*Kh) + std::asinh(Kh)/K);
1108 }
1109 G4double outS = 2.*h*outerRadius;
1110 if (std::abs(endOuterRadius - outerRadius) > kCarTolerance)
1111 {
1115 G4double CC = AA*h*h/(RR - AA);
1116 G4double K = std::sqrt(AA + CC)/CC;
1117 G4double Kh = K*h;
1118 outS = A*(h*std::sqrt(1. + Kh*Kh) + std::asinh(Kh)/K);
1119 }
1120 fSurfaceArea = CLHEP::twopi*(endOuterRadius2 - endInnerRadius2 + innS + outS);
1121 }
1122 return fSurfaceArea;
1123}
const G4double A[17]

◆ GetZHalfLength()

◆ HypeInnerRadius2()

G4double G4Hype::HypeInnerRadius2 ( G4double  zVal) const
inlineprotected

◆ HypeOuterRadius2()

G4double G4Hype::HypeOuterRadius2 ( G4double  zVal) const
inlineprotected

◆ InnerSurfaceExists()

G4bool G4Hype::InnerSurfaceExists ( ) const
inlineprotected

◆ Inside()

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

Implements G4VSolid.

Definition at line 238 of file G4Hype.cc.

239{
240 //
241 // Check z extents: are we outside?
242 //
243 const G4double absZ(std::fabs(p.z()));
244 if (absZ > halfLenZ + fHalfTol) return kOutside;
245
246 //
247 // Check outer radius
248 //
249 const G4double oRad2(HypeOuterRadius2(absZ));
250 const G4double xR2( p.x()*p.x()+p.y()*p.y() );
251
252 if (xR2 > oRad2 + kCarTolerance*endOuterRadius) return kOutside;
253
254 if (xR2 > oRad2 - kCarTolerance*endOuterRadius) return kSurface;
255
256 if (InnerSurfaceExists())
257 {
258 //
259 // Check inner radius
260 //
261 const G4double iRad2(HypeInnerRadius2(absZ));
262
263 if (xR2 < iRad2 - kCarTolerance*endInnerRadius) return kOutside;
264
265 if (xR2 < iRad2 + kCarTolerance*endInnerRadius) return kSurface;
266 }
267
268 //
269 // We are inside in radius, now check endplate surface
270 //
271 if (absZ > halfLenZ - fHalfTol) return kSurface;
272
273 return kInside;
274}
@ kInside
Definition: geomdefs.hh:70
@ kOutside
Definition: geomdefs.hh:68
@ kSurface
Definition: geomdefs.hh:69

◆ IntersectHype()

G4int G4Hype::IntersectHype ( const G4ThreeVector p,
const G4ThreeVector v,
G4double  r2,
G4double  tan2Phi,
G4double  s[2] 
)
staticprotected

Definition at line 921 of file G4Hype.cc.

923{
924 G4double x0 = p.x(), y0 = p.y(), z0 = p.z();
925 G4double tx = v.x(), ty = v.y(), tz = v.z();
926
927 G4double a = tx*tx + ty*ty - tz*tz*tan2Phi;
928 G4double b = 2*( x0*tx + y0*ty - z0*tz*tan2Phi );
929 G4double c = x0*x0 + y0*y0 - r2 - z0*z0*tan2Phi;
930
931 if (std::fabs(a) < DBL_MIN)
932 {
933 //
934 // The trajectory is parallel to the asympotic limit of
935 // the surface: single solution
936 //
937 if (std::fabs(b) < DBL_MIN) return 0;
938 // Unless we travel through exact center
939
940 ss[0] = c/b;
941 return 1;
942 }
943
944 G4double radical = b*b - 4*a*c;
945
946 if (radical < -DBL_MIN) return 0; // No solution
947
948 if (radical < DBL_MIN)
949 {
950 //
951 // Grazes surface
952 //
953 ss[0] = -b/a/2.0;
954 return 1;
955 }
956
957 radical = std::sqrt(radical);
958
959 G4double q = -0.5*( b + (b < 0 ? -radical : +radical) );
960 G4double sa = q/a;
961 G4double sb = c/q;
962 if (sa < sb) { ss[0] = sa; ss[1] = sb; } else { ss[0] = sb; ss[1] = sa; }
963 return 2;
964}

Referenced by DistanceToIn(), and DistanceToOut().

◆ operator=()

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

Definition at line 159 of file G4Hype.cc.

160{
161 // Check assignment to self
162 //
163 if (this == &rhs) { return *this; }
164
165 // Copy base class data
166 //
168
169 // Copy data
170 //
172 halfLenZ = rhs.halfLenZ;
179 fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
180 fHalfTol = rhs.fHalfTol;
181 fRebuildPolyhedron = false;
182 delete fpPolyhedron; fpPolyhedron = nullptr;
183
184 return *this;
185}
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:107

◆ SetInnerRadius()

void G4Hype::SetInnerRadius ( G4double  newIRad)
inline

◆ SetInnerStereo()

void G4Hype::SetInnerStereo ( G4double  newISte)
inline

Referenced by G4Hype().

◆ SetOuterRadius()

void G4Hype::SetOuterRadius ( G4double  newORad)
inline

◆ SetOuterStereo()

void G4Hype::SetOuterStereo ( G4double  newOSte)
inline

Referenced by G4Hype().

◆ SetZHalfLength()

void G4Hype::SetZHalfLength ( G4double  newHLZ)
inline

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 1127 of file G4Hype.cc.

1128{
1129 G4long oldprc = os.precision(16);
1130 os << "-----------------------------------------------------------\n"
1131 << " *** Dump for solid - " << GetName() << " ***\n"
1132 << " ===================================================\n"
1133 << " Solid type: G4Hype\n"
1134 << " Parameters: \n"
1135 << " half length Z: " << halfLenZ/mm << " mm \n"
1136 << " inner radius : " << innerRadius/mm << " mm \n"
1137 << " outer radius : " << outerRadius/mm << " mm \n"
1138 << " inner stereo angle : " << innerStereo/degree << " degrees \n"
1139 << " outer stereo angle : " << outerStereo/degree << " degrees \n"
1140 << "-----------------------------------------------------------\n";
1141 os.precision(oldprc);
1142
1143 return os;
1144}
long G4long
Definition: G4Types.hh:87

◆ SurfaceNormal()

G4ThreeVector G4Hype::SurfaceNormal ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 279 of file G4Hype.cc.

280{
281 //
282 // Which of the three or four surfaces are we closest to?
283 //
284 const G4double absZ(std::fabs(p.z()));
285 const G4double distZ(absZ - halfLenZ);
286 const G4double dist2Z(distZ*distZ);
287
288 const G4double xR2( p.x()*p.x()+p.y()*p.y() );
289 const G4double dist2Outer( std::fabs(xR2 - HypeOuterRadius2(absZ)) );
290
291 if (InnerSurfaceExists())
292 {
293 //
294 // Has inner surface: is this closest?
295 //
296 const G4double dist2Inner( std::fabs(xR2 - HypeInnerRadius2(absZ)) );
297 if (dist2Inner < dist2Z && dist2Inner < dist2Outer)
298 return G4ThreeVector( -p.x(), -p.y(), p.z()*tanInnerStereo2 ).unit();
299 }
300
301 //
302 // Do the "endcaps" win?
303 //
304 if (dist2Z < dist2Outer)
305 return G4ThreeVector( 0.0, 0.0, p.z() < 0 ? -1.0 : 1.0 );
306
307
308 //
309 // Outer surface wins
310 //
311 return G4ThreeVector( p.x(), p.y(), -p.z()*tanOuterStereo2 ).unit();
312}

Member Data Documentation

◆ endInnerRadius

G4double G4Hype::endInnerRadius
protected

Definition at line 180 of file G4Hype.hh.

Referenced by DistanceToIn(), DistanceToOut(), GetSurfaceArea(), Inside(), and operator=().

◆ endInnerRadius2

G4double G4Hype::endInnerRadius2
protected

Definition at line 178 of file G4Hype.hh.

Referenced by DistanceToIn(), GetCubicVolume(), GetSurfaceArea(), and operator=().

◆ endOuterRadius

G4double G4Hype::endOuterRadius
protected

◆ endOuterRadius2

G4double G4Hype::endOuterRadius2
protected

Definition at line 179 of file G4Hype.hh.

Referenced by DistanceToIn(), GetCubicVolume(), GetSurfaceArea(), and operator=().

◆ halfLenZ

◆ innerRadius

G4double G4Hype::innerRadius
protected

◆ innerRadius2

G4double G4Hype::innerRadius2
protected

◆ innerStereo

G4double G4Hype::innerStereo
protected

Definition at line 167 of file G4Hype.hh.

Referenced by DistanceToIn(), GetPointOnSurface(), operator=(), and StreamInfo().

◆ outerRadius

G4double G4Hype::outerRadius
protected

◆ outerRadius2

G4double G4Hype::outerRadius2
protected

◆ outerStereo

G4double G4Hype::outerStereo
protected

Definition at line 168 of file G4Hype.hh.

Referenced by GetPointOnSurface(), operator=(), and StreamInfo().

◆ tanInnerStereo

G4double G4Hype::tanInnerStereo
protected

Definition at line 172 of file G4Hype.hh.

Referenced by DistanceToOut(), GetPointOnSurface(), and operator=().

◆ tanInnerStereo2

G4double G4Hype::tanInnerStereo2
protected

◆ tanOuterStereo

G4double G4Hype::tanOuterStereo
protected

Definition at line 173 of file G4Hype.hh.

Referenced by DistanceToIn(), GetPointOnSurface(), and operator=().

◆ tanOuterStereo2

G4double G4Hype::tanOuterStereo2
protected

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