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

#include <G4EllipticalTube.hh>

+ Inheritance diagram for G4EllipticalTube:

Public Member Functions

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

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

Definition at line 60 of file G4EllipticalTube.hh.

Constructor & Destructor Documentation

◆ G4EllipticalTube() [1/3]

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

Definition at line 61 of file G4EllipticalTube.cc.

65 : G4VSolid(name), fDx(Dx), fDy(Dy), fDz(Dz)
66{
67 CheckParameters();
68}
G4VSolid(const G4String &name)
Definition G4VSolid.cc:57

Referenced by Clone().

◆ ~G4EllipticalTube()

G4EllipticalTube::~G4EllipticalTube ( )
override

Definition at line 86 of file G4EllipticalTube.cc.

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

◆ G4EllipticalTube() [2/3]

G4EllipticalTube::G4EllipticalTube ( __void__ & a)

Definition at line 75 of file G4EllipticalTube.cc.

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

◆ G4EllipticalTube() [3/3]

G4EllipticalTube::G4EllipticalTube ( const G4EllipticalTube & rhs)

Definition at line 95 of file G4EllipticalTube.cc.

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

Member Function Documentation

◆ BoundingLimits()

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

Reimplemented from G4VSolid.

Definition at line 187 of file G4EllipticalTube.cc.

189{
190 pMin.set(-fDx,-fDy,-fDz);
191 pMax.set( fDx, fDy, fDz);
192}
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

Implements G4VSolid.

Definition at line 199 of file G4EllipticalTube.cc.

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

Reimplemented from G4VSolid.

Definition at line 629 of file G4EllipticalTube.cc.

630{
631 return new G4EllipticalTube(*this);
632}
G4EllipticalTube(const G4String &name, G4double Dx, G4double Dy, G4double Dz)

◆ CreatePolyhedron()

G4Polyhedron * G4EllipticalTube::CreatePolyhedron ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 749 of file G4EllipticalTube.cc.

750{
751 // create cylinder with radius=1...
752 //
753 G4Polyhedron* eTube = new G4PolyhedronTube(0., 1., fDz);
754
755 // apply non-uniform scaling...
756 //
757 eTube->Transform(G4Scale3D(fDx, fDy, 1.));
758 return eTube;
759}
HepGeom::Scale3D G4Scale3D
HepPolyhedron & Transform(const G4Transform3D &t)

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

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

Implements G4VSolid.

Definition at line 785 of file G4EllipticalTube.cc.

786{
787 scene.AddSolid (*this);
788}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 433 of file G4EllipticalTube.cc.

434{
435 // safety distance to bounding box
436 G4double distX = std::abs(p.x()) - fDx;
437 G4double distY = std::abs(p.y()) - fDy;
438 G4double distZ = std::abs(p.z()) - fDz;
439 G4double distB = std::max(std::max(distX, distY), distZ);
440 // return (distB < 0) ? 0 : distB;
441
442 // safety distance to lateral surface
443 G4double x = p.x() * fSx;
444 G4double y = p.y() * fSy;
445 G4double distR = std::sqrt(x * x + y * y) - fR;
446
447 // return SafetyToIn
448 G4double dist = std::max(distB, distR);
449 return (dist < 0) ? 0 : dist;
450}
double z() const
double x() const
double y() const

◆ DistanceToIn() [2/2]

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

Implements G4VSolid.

Definition at line 349 of file G4EllipticalTube.cc.

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

Implements G4VSolid.

Definition at line 585 of file G4EllipticalTube.cc.

586{
587#ifdef G4SPECDEBUG
588 if( Inside(p) == kOutside )
589 {
590 std::ostringstream message;
591 G4long oldprc = message.precision(16);
592 message << "Point p is outside (!?) of solid: " << GetName() << "\n"
593 << "Position:\n"
594 << " p.x() = " << p.x()/mm << " mm\n"
595 << " p.y() = " << p.y()/mm << " mm\n"
596 << " p.z() = " << p.z()/mm << " mm";
597 message.precision(oldprc) ;
598 G4Exception("G4ElliptocalTube::DistanceToOut(p)", "GeomSolids1002",
599 JustWarning, message);
600 DumpInfo();
601 }
602#endif
603 // safety distance to Z-bases
604 G4double distZ = fDz - std::abs(p.z());
605
606 // safety distance lateral surface
607 G4double x = p.x() * fSx;
608 G4double y = p.y() * fSy;
609 G4double distR = fR - std::sqrt(x * x + y * y);
610
611 // return SafetyToOut
612 G4double dist = std::min(distZ, distR);
613 return (dist < 0) ? 0 : dist;
614}
@ 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

Implements G4VSolid.

Definition at line 458 of file G4EllipticalTube.cc.

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

◆ GetCubicVolume()

G4double G4EllipticalTube::GetCubicVolume ( )
overridevirtual

Reimplemented from G4VSolid.

Definition at line 638 of file G4EllipticalTube.cc.

639{
640 if (fCubicVolume == 0.)
641 {
642 fCubicVolume = twopi * fDx * fDy * fDz;
643 }
644 return fCubicVolume;
645}

◆ GetDx()

G4double G4EllipticalTube::GetDx ( ) const
inline

◆ GetDy()

G4double G4EllipticalTube::GetDy ( ) const
inline

◆ GetDz()

G4double G4EllipticalTube::GetDz ( ) const
inline

◆ GetEntityType()

G4GeometryType G4EllipticalTube::GetEntityType ( ) const
overridevirtual

Implements G4VSolid.

Definition at line 620 of file G4EllipticalTube.cc.

621{
622 return {"G4EllipticalTube"};
623}

◆ GetExtent()

G4VisExtent G4EllipticalTube::GetExtent ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 794 of file G4EllipticalTube.cc.

795{
796 return { -fDx, fDx, -fDy, fDy, -fDz, fDz };
797}

◆ GetPointOnSurface()

G4ThreeVector G4EllipticalTube::GetPointOnSurface ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 706 of file G4EllipticalTube.cc.

707{
708 // Select surface (0 - base at -Z, 1 - base at +Z, 2 - lateral surface)
709 //
710 G4double sbase = pi * fDx * fDy;
711 G4double ssurf = GetCachedSurfaceArea();
712 G4double select = ssurf * G4UniformRand();
713
714 G4int k = 0;
715 if (select > sbase) k = 1;
716 if (select > 2. * sbase) k = 2;
717
718 // Pick random point on selected surface (rejection sampling)
719 //
721 switch (k) {
722 case 0: // base at -Z
723 {
724 G4TwoVector rho = G4RandomPointInEllipse(fDx, fDy);
725 p.set(rho.x(), rho.y(), -fDz);
726 break;
727 }
728 case 1: // base at +Z
729 {
730 G4TwoVector rho = G4RandomPointInEllipse(fDx, fDy);
731 p.set(rho.x(), rho.y(), fDz);
732 break;
733 }
734 case 2: // lateral surface
735 {
736 G4TwoVector rho = G4RandomPointOnEllipse(fDx, fDy);
737 p.set(rho.x(), rho.y(), (2. * G4UniformRand() - 1.) * fDz);
738 break;
739 }
740 }
741 return p;
742}
G4TwoVector G4RandomPointInEllipse(G4double a, G4double b)
G4TwoVector G4RandomPointOnEllipse(G4double a, G4double b)
#define G4UniformRand()
Definition Randomize.hh:52
double x() const
double y() const

◆ GetPolyhedron()

G4Polyhedron * G4EllipticalTube::GetPolyhedron ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 765 of file G4EllipticalTube.cc.

766{
767 if (fpPolyhedron == nullptr ||
768 fRebuildPolyhedron ||
770 fpPolyhedron->GetNumberOfRotationSteps())
771 {
772 G4AutoLock l(&polyhedronMutex);
773 delete fpPolyhedron;
774 fpPolyhedron = CreatePolyhedron();
775 fRebuildPolyhedron = false;
776 l.unlock();
777 }
778 return fpPolyhedron;
779}
G4Polyhedron * CreatePolyhedron() const override
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
static G4int GetNumberOfRotationSteps()

◆ GetSurfaceArea()

G4double G4EllipticalTube::GetSurfaceArea ( )
overridevirtual

Reimplemented from G4VSolid.

Definition at line 671 of file G4EllipticalTube.cc.

672{
673 if(fSurfaceArea == 0.)
674 {
675 fSurfaceArea = GetCachedSurfaceArea();
676 }
677 return fSurfaceArea;
678}

◆ Inside()

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

Implements G4VSolid.

Definition at line 261 of file G4EllipticalTube.cc.

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

Referenced by DistanceToOut().

◆ operator=()

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

Definition at line 109 of file G4EllipticalTube.cc.

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

◆ SetDx()

void G4EllipticalTube::SetDx ( G4double Dx)
inline

◆ SetDy()

void G4EllipticalTube::SetDy ( G4double Dy)
inline

◆ SetDz()

void G4EllipticalTube::SetDz ( G4double Dz)
inline

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 684 of file G4EllipticalTube.cc.

685{
686 G4long oldprc = os.precision(16);
687 os << "-----------------------------------------------------------\n"
688 << " *** Dump for solid - " << GetName() << " ***\n"
689 << " ===================================================\n"
690 << " Solid type: G4EllipticalTube\n"
691 << " Parameters: \n"
692 << " length Z: " << fDz/mm << " mm \n"
693 << " lateral surface equation: \n"
694 << " (X / " << fDx << ")^2 + (Y / " << fDy << ")^2 = 1 \n"
695 << "-----------------------------------------------------------\n";
696 os.precision(oldprc);
697
698 return os;
699}

◆ SurfaceNormal()

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

Implements G4VSolid.

Definition at line 277 of file G4EllipticalTube.cc.

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

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