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

G4EllipticalTube is a tube with elliptical cross section. The equation of the lateral surface is: (x/dx)^2 + (y/dy)^2 = 1. More...

#include <G4EllipticalTube.hh>

Inheritance diagram for G4EllipticalTube:

Public Member Functions

 G4EllipticalTube (const G4String &name, G4double Dx, G4double Dy, G4double Dz)
 ~G4EllipticalTube () override
G4double GetDx () const
G4double GetDy () const
G4double GetDz () const
void SetDx (G4double Dx)
void SetDy (G4double Dy)
void SetDz (G4double Dz)
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
G4PolyhedronCreatePolyhedron () const override
G4PolyhedronGetPolyhedron () const override
void DescribeYourselfTo (G4VGraphicsScene &scene) const override
G4VisExtent GetExtent () const override
 G4EllipticalTube (__void__ &)
 G4EllipticalTube (const G4EllipticalTube &rhs)
G4EllipticalTubeoperator= (const G4EllipticalTube &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 void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
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

G4EllipticalTube is a tube with elliptical cross section. The equation of the lateral surface is: (x/dx)^2 + (y/dy)^2 = 1.

Definition at line 64 of file G4EllipticalTube.hh.

Constructor & Destructor Documentation

◆ G4EllipticalTube() [1/3]

G4EllipticalTube::G4EllipticalTube ( const G4String & name,
G4double Dx,
G4double Dy,
G4double Dz )

Constructs an elliptical tube, given its parameters.

Parameters
[in]nameThe solid name.
[in]DxHalf length of axis along X.
[in]DyHalf length of axis along Y.
[in]DzHalf length in Z.

Definition at line 59 of file G4EllipticalTube.cc.

63 : G4VSolid(name), fDx(Dx), fDy(Dy), fDz(Dz)
64{
65 CheckParameters();
66 fSurfaceArea = GetCachedSurfaceArea();
67}
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59

Referenced by Clone(), G4EllipticalTube(), operator=(), and SetDz().

◆ ~G4EllipticalTube()

G4EllipticalTube::~G4EllipticalTube ( )
override

Destructor.

Definition at line 85 of file G4EllipticalTube.cc.

86{
87 delete fpPolyhedron; fpPolyhedron = nullptr;
88}

◆ G4EllipticalTube() [2/3]

G4EllipticalTube::G4EllipticalTube ( __void__ & a)

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

Definition at line 74 of file G4EllipticalTube.cc.

75 : G4VSolid(a), halfTolerance(0.), fDx(0.), fDy(0.), fDz(0.),
76 fRsph(0.), fDDx(0.), fDDy(0.), fSx(0.), fSy(0.), fR(0.),
77 fQ1(0.), fQ2(0.), fScratch(0.)
78{
79}

◆ G4EllipticalTube() [3/3]

G4EllipticalTube::G4EllipticalTube ( const G4EllipticalTube & rhs)

Copy constructor and assignment operator.

Definition at line 94 of file G4EllipticalTube.cc.

95 : G4VSolid(rhs), halfTolerance(rhs.halfTolerance),
96 fDx(rhs.fDx), fDy(rhs.fDy), fDz(rhs.fDz),
97 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
98 fRsph(rhs.fRsph), fDDx(rhs.fDDx), fDDy(rhs.fDDy),
99 fSx(rhs.fSx), fSy(rhs.fSy), fR(rhs.fR),
100 fQ1(rhs.fQ1), fQ2(rhs.fQ2), fScratch(rhs.fScratch)
101{
102}

Member Function Documentation

◆ BoundingLimits()

void G4EllipticalTube::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 186 of file G4EllipticalTube.cc.

188{
189 pMin.set(-fDx,-fDy,-fDz);
190 pMax.set( fDx, fDy, fDz);
191}
void set(double x, double y, double z)

Referenced by CalculateExtent().

◆ CalculateExtent()

G4bool G4EllipticalTube::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 198 of file G4EllipticalTube.cc.

202{
203 G4ThreeVector bmin, bmax;
204 G4bool exist;
205
206 // Check bounding box (bbox)
207 //
208 BoundingLimits(bmin,bmax);
209 G4BoundingEnvelope bbox(bmin,bmax);
210#ifdef G4BBOX_EXTENT
211 return bbox.CalculateExtent(pAxis,pVoxelLimit, pTransform, pMin, pMax);
212#endif
213 if (bbox.BoundingBoxVsVoxelLimits(pAxis, pVoxelLimit, pTransform, pMin, pMax))
214 {
215 return exist = pMin < pMax;
216 }
217
218 G4double dx = fDx;
219 G4double dy = fDy;
220 G4double dz = fDz;
221
222 // Set bounding envelope (benv) and calculate extent
223 //
224 const G4int NSTEPS = 24; // number of steps for whole circle
225 G4double ang = twopi/NSTEPS;
226
227 G4double sinHalf = std::sin(0.5*ang);
228 G4double cosHalf = std::cos(0.5*ang);
229 G4double sinStep = 2.*sinHalf*cosHalf;
230 G4double cosStep = 1. - 2.*sinHalf*sinHalf;
231 G4double sx = dx/cosHalf;
232 G4double sy = dy/cosHalf;
233
234 G4double sinCur = sinHalf;
235 G4double cosCur = cosHalf;
236 G4ThreeVectorList baseA(NSTEPS),baseB(NSTEPS);
237 for (G4int k=0; k<NSTEPS; ++k)
238 {
239 baseA[k].set(sx*cosCur,sy*sinCur,-dz);
240 baseB[k].set(sx*cosCur,sy*sinCur, dz);
241
242 G4double sinTmp = sinCur;
243 sinCur = sinCur*cosStep + cosCur*sinStep;
244 cosCur = cosCur*cosStep - sinTmp*sinStep;
245 }
246
247 std::vector<const G4ThreeVectorList *> polygons(2);
248 polygons[0] = &baseA;
249 polygons[1] = &baseB;
250 G4BoundingEnvelope benv(bmin, bmax, polygons);
251 exist = benv.CalculateExtent(pAxis, pVoxelLimit, pTransform, pMin, pMax);
252 return exist;
253}
std::vector< G4ThreeVector > G4ThreeVectorList
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override

◆ Clone()

G4VSolid * G4EllipticalTube::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 641 of file G4EllipticalTube.cc.

642{
643 return new G4EllipticalTube(*this);
644}
G4EllipticalTube(const G4String &name, G4double Dx, G4double Dy, G4double Dz)

◆ CreatePolyhedron()

G4Polyhedron * G4EllipticalTube::CreatePolyhedron ( ) const
overridevirtual

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

Reimplemented from G4VSolid.

Definition at line 758 of file G4EllipticalTube.cc.

759{
760 // create cylinder with radius=1...
761 //
762 G4Polyhedron* eTube = new G4PolyhedronTube(0., 1., fDz);
763
764 // apply non-uniform scaling...
765 //
766 eTube->Transform(G4Scale3D(fDx, fDy, 1.));
767 return eTube;
768}
HepGeom::Scale3D G4Scale3D
HepPolyhedron & Transform(const G4Transform3D &t)

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

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

A "double dispatch" function which identifies the solid to the graphics scene for visualization.

Implements G4VSolid.

Definition at line 794 of file G4EllipticalTube.cc.

795{
796 scene.AddSolid (*this);
797}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

G4double G4EllipticalTube::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 441 of file G4EllipticalTube.cc.

442{
443 // safety distance to bounding box
444 G4double distX = std::abs(p.x()) - fDx;
445 G4double distY = std::abs(p.y()) - fDy;
446 G4double distZ = std::abs(p.z()) - fDz;
447 G4double distB = std::max(std::max(distX, distY), distZ);
448 // return (distB < 0) ? 0 : distB;
449
450 // safety distance to lateral surface
451 G4double x = p.x() * fSx;
452 G4double y = p.y() * fSy;
453 G4double distR = std::sqrt(x * x + y * y) - fR;
454
455 // return SafetyToIn
456 G4double dist = std::max(distB, distR);
457 return (dist < 0) ? 0 : dist;
458}
double z() const
double x() const
double y() const

◆ DistanceToIn() [2/2]

G4double G4EllipticalTube::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 353 of file G4EllipticalTube.cc.

355{
356 G4double offset = 0.;
357 G4ThreeVector pcur = p;
358
359 // Check if point is flying away
360 //
361 G4double safex = std::abs(pcur.x()) - fDx;
362 G4double safey = std::abs(pcur.y()) - fDy;
363 G4double safez = std::abs(pcur.z()) - fDz;
364
365 if (safez >= -halfTolerance && pcur.z() * v.z() >= 0.) { return kInfinity; }
366 if (safey >= -halfTolerance && pcur.y() * v.y() >= 0.) { return kInfinity; }
367 if (safex >= -halfTolerance && pcur.x() * v.x() >= 0.) { return kInfinity; }
368
369 // Relocate point, if required
370 //
371 G4double Dmax = 32. * fRsph;
372 if (std::max(std::max(safex, safey), safez) > Dmax)
373 {
374 offset = (1. - 1.e-08) * pcur.mag() - 2. * fRsph;
375 pcur += offset * v;
376 G4double dist = DistanceToIn(pcur, v);
377 return (dist == kInfinity) ? kInfinity : dist + offset;
378 }
379
380 // Scale elliptical tube to cylinder
381 //
382 G4double px = pcur.x() * fSx;
383 G4double py = pcur.y() * fSy;
384 G4double pz = pcur.z();
385 G4double vx = v.x() * fSx;
386 G4double vy = v.y() * fSy;
387 G4double vz = v.z();
388
389 // Set coefficients of quadratic equation: A t^2 + 2B t + C = 0
390 //
391 G4double rr = px * px + py * py;
392 G4double A = vx * vx + vy * vy;
393 G4double B = px * vx + py * vy;
394 G4double C = rr - fR * fR;
395 G4double D = B * B - A * C;
396
397 // Check if point is flying away relative to lateral surface
398 //
399 G4double distR = fQ1 * rr - fQ2;
400 G4bool parallelToZ = (A < DBL_EPSILON || std::abs(vz) >= 1.);
401 if (distR >= -halfTolerance && (B >= 0. || parallelToZ)) { return kInfinity; }
402
403 // Find intersection with Z planes
404 //
405 G4double invz = (vz == 0) ? DBL_MAX : -1./vz;
406 G4double dz = std::copysign(fDz, invz);
407 G4double tzmin = (pz - dz) * invz;
408 G4double tzmax = (pz + dz) * invz;
409
410 // Solve qudratic equation. There are two cases special where D <= 0:
411 // 1) trajectory parallel to Z axis (A = 0, B = 0, C - any, D = 0)
412 // 2) touch (D = 0) or no intersection (D < 0) with lateral surface
413 //
414 if (parallelToZ)
415 {
416 return (tzmin<halfTolerance) ? offset : tzmin + offset; // 1)
417 }
418 if (D <= A * A * fScratch) { return kInfinity; } // 2)
419
420 // Find roots of quadratic equation
421 G4double tmp = -B - std::copysign(std::sqrt(D), B);
422 G4double t1 = tmp / A;
423 G4double t2 = C / tmp;
424 G4double trmin = std::min(t1, t2);
425 G4double trmax = std::max(t1, t2);
426
427 // Return distance
428 G4double tin = std::max(tzmin, trmin);
429 G4double tout = std::min(tzmax, trmax);
430
431 if (tout <= tin + halfTolerance) { return kInfinity; } // touch or no hit
432
433 return (tin<halfTolerance) ? offset : tin + offset;
434}
G4double C(G4double temp)
G4double B(G4double temperature)
G4double D(G4double temp)
G4ThreadLocal T * G4GeomSplitter< T >::offset
const G4double A[17]
double mag() const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
#define DBL_EPSILON
Definition templates.hh:66
#define DBL_MAX
Definition templates.hh:62

Referenced by DistanceToIn().

◆ DistanceToOut() [1/2]

G4double G4EllipticalTube::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 597 of file G4EllipticalTube.cc.

598{
599#ifdef G4SPECDEBUG
600 if( Inside(p) == kOutside )
601 {
602 std::ostringstream message;
603 G4long oldprc = message.precision(16);
604 message << "Point p is outside (!?) of solid: " << GetName() << "\n"
605 << "Position:\n"
606 << " p.x() = " << p.x()/mm << " mm\n"
607 << " p.y() = " << p.y()/mm << " mm\n"
608 << " p.z() = " << p.z()/mm << " mm";
609 message.precision(oldprc) ;
610 G4Exception("G4ElliptocalTube::DistanceToOut(p)", "GeomSolids1002",
611 JustWarning, message);
612 DumpInfo();
613 }
614#endif
615 // safety distance to Z-bases
616 G4double distZ = fDz - std::abs(p.z());
617
618 // safety distance lateral surface
619 G4double x = p.x() * fSx;
620 G4double y = p.y() * fSy;
621 G4double distR = fR - std::sqrt(x * x + y * y);
622
623 // return SafetyToOut
624 G4double dist = std::min(distZ, distR);
625 return (dist < 0) ? 0 : dist;
626}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
long G4long
Definition G4Types.hh:87
EInside Inside(const G4ThreeVector &p) const override
G4String GetName() const
void DumpInfo() const
@ kOutside
Definition geomdefs.hh:68

◆ DistanceToOut() [2/2]

G4double G4EllipticalTube::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 466 of file G4EllipticalTube.cc.

471{
472 // Check if point flying away relative to Z planes
473 //
474 G4double pz = p.z();
475 G4double vz = v.z();
476 G4double distZ = std::abs(pz) - fDz;
477 if (distZ >= -halfTolerance && pz * vz > 0)
478 {
479 if (calcNorm)
480 {
481 *validNorm = true;
482 n->set(0, 0, (pz < 0) ? -1. : 1.);
483 }
484 return 0.;
485 }
486 G4double tzmax = (vz == 0) ? DBL_MAX : (std::copysign(fDz, vz) - pz) / vz;
487
488 // Scale elliptical tube to cylinder
489 //
490 G4double px = p.x() * fSx;
491 G4double py = p.y() * fSy;
492 G4double vx = v.x() * fSx;
493 G4double vy = v.y() * fSy;
494
495 // Check if point is flying away relative to lateral surface
496 //
497 G4double rr = px * px + py * py;
498 G4double B = px * vx + py * vy;
499 G4double distR = fQ1 * rr - fQ2;
500 if (distR >= -halfTolerance && B > 0.)
501 {
502 if (calcNorm)
503 {
504 *validNorm = true;
505 *n = G4ThreeVector(px * fDDy, py * fDDx, 0.).unit();
506 }
507 return 0.;
508 }
509
510 // Just in case check if point is outside, normally it should never be
511 //
512 if (std::max(distZ, distR) > halfTolerance)
513 {
514#ifdef G4SPECDEBUG
515 std::ostringstream message;
516 G4long oldprc = message.precision(16);
517 message << "Point p is outside (!?) of solid: "
518 << GetName() << G4endl;
519 message << "Position: " << p << G4endl;;
520 message << "Direction: " << v;
521 G4cout.precision(oldprc);
522 G4Exception("G4EllipticalTube::DistanceToOut(p,v)", "GeomSolids1002",
523 JustWarning, message );
524 DumpInfo();
525#endif
526 if (calcNorm)
527 {
528 *validNorm = true;
529 *n = ApproxSurfaceNormal(p);
530 }
531 return 0.;
532 }
533
534 // Set coefficients of quadratic equation: A t^2 + 2B t + C = 0
535 //
536 G4double A = vx * vx + vy * vy;
537 G4double C = rr - fR * fR;
538 G4double D = B * B - A * C;
539
540 // Solve qudratic equation. There are two special cases where D <= 0:
541 // 1) trajectory parallel to Z axis (A = 0, B = 0, C - any, D = 0)
542 // 2) touch (D = 0) or no intersection (D < 0) with lateral surface
543 //
544 G4bool parallelToZ = (A < DBL_EPSILON || std::abs(vz) >= 1.);
545 if (parallelToZ) // 1)
546 {
547 if (calcNorm)
548 {
549 *validNorm = true;
550 n->set(0, 0, (vz < 0) ? -1. : 1.);
551 }
552 return tzmax;
553 }
554 if (D <= A * A * fScratch) // 2)
555 {
556 if (calcNorm)
557 {
558 *validNorm = true;
559 *n = G4ThreeVector(px * fDDy, py * fDDx, 0.).unit();
560 }
561 return 0.;
562 }
563
564 // Find roots of quadratic equation
565 G4double tmp = -B - std::copysign(std::sqrt(D), B);
566 G4double t1 = tmp / A;
567 G4double t2 = C / tmp;
568 G4double trmax = std::max(t1, t2);
569
570 // Return distance
571 G4double tmax = std::min(tzmax, trmax);
572
573 // Set normal, if required, and return distance
574 //
575 if (calcNorm)
576 {
577 *validNorm = true;
578 G4ThreeVector pnew = p + tmax * v;
579 if (tmax == tzmax)
580 {
581 n->set(0, 0, (pnew.z() < 0) ? -1. : 1.);
582 }
583 else
584 {
585 *n = G4ThreeVector(pnew.x() * fDDy, pnew.y() * fDDx, 0.).unit();
586 }
587 }
588 return tmax;
589}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const

◆ GetCubicVolume()

G4double G4EllipticalTube::GetCubicVolume ( )
overridevirtual

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

Reimplemented from G4VSolid.

Definition at line 670 of file G4EllipticalTube.cc.

671{
672 if (fCubicVolume == 0)
673 {
674 G4AutoLock l(&eltubeMutex);
675 fCubicVolume = twopi * fDx * fDy * fDz;
676 l.unlock();
677 }
678 return fCubicVolume;
679}
G4TemplateAutoLock< G4Mutex > G4AutoLock

◆ GetDx()

G4double G4EllipticalTube::GetDx ( ) const
inline

◆ GetDy()

G4double G4EllipticalTube::GetDy ( ) const
inline

◆ GetDz()

G4double G4EllipticalTube::GetDz ( ) const
inline

◆ GetEntityType()

G4GeometryType G4EllipticalTube::GetEntityType ( ) const
overridevirtual

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

Implements G4VSolid.

Definition at line 632 of file G4EllipticalTube.cc.

633{
634 return {"G4EllipticalTube"};
635}

◆ GetExtent()

G4VisExtent G4EllipticalTube::GetExtent ( ) const
overridevirtual

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

Reimplemented from G4VSolid.

Definition at line 803 of file G4EllipticalTube.cc.

804{
805 return { -fDx, fDx, -fDy, fDy, -fDz, fDz };
806}

◆ GetPointOnSurface()

G4ThreeVector G4EllipticalTube::GetPointOnSurface ( ) const
overridevirtual

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

Reimplemented from G4VSolid.

Definition at line 722 of file G4EllipticalTube.cc.

723{
724 // Pick random point on selected surface
725 //
726 G4double sbase = pi * fDx * fDy; // base area
727 G4double select = fSurfaceArea * G4QuickRand();
728 G4double x, y, z;
729 if (select < 2. * sbase)
730 { // base (ellipse)
731 G4double phi = CLHEP::twopi * G4QuickRand();
732 G4double rho = std::sqrt(G4QuickRand());
733 x = rho * std::cos(phi);
734 y = rho * std::sin(phi);
735 z = (select < sbase) ? fDz : -fDz;
736 }
737 else
738 { // lateral surface (rejection sampling)
739 G4double s_max = std::max(fDx, fDy);
740 for (auto i = 0; i < 10000; ++i)
741 {
742 G4double phi = CLHEP::twopi * G4QuickRand();
743 x = std::cos(phi);
744 y = std::sin(phi);
745 G4double ss = sqr(fDy * x) + sqr(fDx * y);
746 if (sqr(s_max*G4QuickRand()) <= ss) { break; }
747 }
748 z = (2. * G4QuickRand() - 1.) * fDz;
749 }
750 return { fDx * x, fDy * y, z };
751}
G4double G4QuickRand(uint32_t seed=0)
const G4double pi
T sqr(const T &x)
Definition templates.hh:128

◆ GetPolyhedron()

G4Polyhedron * G4EllipticalTube::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 774 of file G4EllipticalTube.cc.

775{
776 if (fpPolyhedron == nullptr ||
777 fRebuildPolyhedron ||
778 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
779 fpPolyhedron->GetNumberOfRotationSteps())
780 {
781 G4AutoLock l(&polyhedronMutex);
782 delete fpPolyhedron;
783 fpPolyhedron = CreatePolyhedron();
784 fRebuildPolyhedron = false;
785 l.unlock();
786 }
787 return fpPolyhedron;
788}
G4Polyhedron * CreatePolyhedron() const override

◆ GetSurfaceArea()

G4double G4EllipticalTube::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 685 of file G4EllipticalTube.cc.

686{
687 if(fSurfaceArea == 0)
688 {
689 G4AutoLock l(&eltubeMutex);
690 fSurfaceArea = GetCachedSurfaceArea();
691 l.unlock();
692 }
693 return fSurfaceArea;
694}

◆ Inside()

EInside G4EllipticalTube::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 260 of file G4EllipticalTube.cc.

261{
262 G4double x = p.x() * fSx;
263 G4double y = p.y() * fSy;
264 G4double distR = fQ1 * (x * x + y * y) - fQ2;
265 G4double distZ = std::abs(p.z()) - fDz;
266 G4double dist = std::max(distR, distZ);
267
268 if (dist > halfTolerance) { return kOutside; }
269 return (dist > -halfTolerance) ? kSurface : kInside;
270}
@ kInside
Definition geomdefs.hh:70
@ kSurface
Definition geomdefs.hh:69

Referenced by DistanceToOut().

◆ operator=()

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

Definition at line 108 of file G4EllipticalTube.cc.

109{
110 // Check assignment to self
111 //
112 if (this == &rhs) { return *this; }
113
114 // Copy base class data
115 //
117
118 // Copy data
119 //
120 halfTolerance = rhs.halfTolerance;
121 fDx = rhs.fDx;
122 fDy = rhs.fDy;
123 fDz = rhs.fDz;
124 fCubicVolume = rhs.fCubicVolume;
125 fSurfaceArea = rhs.fSurfaceArea;
126
127 fRsph = rhs.fRsph;
128 fDDx = rhs.fDDx;
129 fDDy = rhs.fDDy;
130 fSx = rhs.fSx;
131 fSy = rhs.fSy;
132 fR = rhs.fR;
133 fQ1 = rhs.fQ1;
134 fQ2 = rhs.fQ2;
135 fScratch = rhs.fScratch;
136
137 fRebuildPolyhedron = false;
138 delete fpPolyhedron; fpPolyhedron = nullptr;
139
140 return *this;
141}
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:108

◆ SetDx()

void G4EllipticalTube::SetDx ( G4double Dx)
inline

Modifiers.

◆ SetDy()

void G4EllipticalTube::SetDy ( G4double Dy)
inline

◆ SetDz()

void G4EllipticalTube::SetDz ( G4double Dz)
inline

◆ StreamInfo()

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

Streams the object contents to an output stream.

Implements G4VSolid.

Definition at line 700 of file G4EllipticalTube.cc.

701{
702 G4long oldprc = os.precision(16);
703 os << "-----------------------------------------------------------\n"
704 << " *** Dump for solid - " << GetName() << " ***\n"
705 << " ===================================================\n"
706 << " Solid type: G4EllipticalTube\n"
707 << " Parameters: \n"
708 << " length Z: " << fDz/mm << " mm \n"
709 << " lateral surface equation: \n"
710 << " (X / " << fDx << ")^2 + (Y / " << fDy << ")^2 = 1 \n"
711 << "-----------------------------------------------------------\n";
712 os.precision(oldprc);
713
714 return os;
715}

◆ SurfaceNormal()

G4ThreeVector G4EllipticalTube::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 276 of file G4EllipticalTube.cc.

277{
278 G4ThreeVector norm(0, 0, 0);
279 G4int nsurf = 0;
280
281 // check lateral surface
282 G4double x = p.x() * fSx;
283 G4double y = p.y() * fSy;
284 G4double distR = fQ1 * (x * x + y * y) - fQ2;
285 if (std::abs(distR) <= halfTolerance)
286 {
287 norm = G4ThreeVector(p.x() * fDDy, p.y() * fDDx, 0.).unit();
288 ++nsurf;
289 }
290
291 // check lateral bases
292 G4double distZ = std::abs(p.z()) - fDz;
293 if (std::abs(distZ) <= halfTolerance)
294 {
295 norm.setZ(p.z() < 0 ? -1. : 1.);
296 ++nsurf;
297 }
298
299 // return normal
300 if (nsurf == 1)
301 {
302 return norm;
303 }
304 if (nsurf > 1)
305 {
306 return norm.unit(); // edge
307 }
308
309 // Point is not on the surface
310 //
311#ifdef G4SPECDEBUG
312 std::ostringstream message;
313 G4long oldprc = message.precision(16);
314 message << "Point p is not on surface (!?) of solid: "
315 << GetName() << G4endl;
316 message << "Position:\n";
317 message << " p.x() = " << p.x()/mm << " mm\n";
318 message << " p.y() = " << p.y()/mm << " mm\n";
319 message << " p.z() = " << p.z()/mm << " mm";
320 G4cout.precision(oldprc);
321 G4Exception("G4EllipticalTube::SurfaceNormal(p)", "GeomSolids1002",
322 JustWarning, message );
323 DumpInfo();
324#endif
325 return ApproxSurfaceNormal(p);
326}

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