Geant4 9.6.0
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 theDx, G4double theDy, G4double theDz)
 
virtual ~G4EllipticalTube ()
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
 
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=0, G4ThreeVector *n=0) 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
 
G4PolyhedronCreatePolyhedron () const
 
G4PolyhedronGetPolyhedron () const
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
G4VisExtent GetExtent () const
 
G4double GetDx () const
 
G4double GetDy () const
 
G4double GetDz () const
 
void SetDx (const G4double newDx)
 
void SetDy (const G4double newDy)
 
void SetDz (const G4double newDz)
 
 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 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=0, G4ThreeVector *n=0) 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 G4NURBSCreateNURBS () 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)
 

Protected Member Functions

G4double CheckXY (const G4double x, const G4double y, const G4double toler) const
 
G4double CheckXY (const G4double x, const G4double y) const
 
G4int IntersectXY (const G4ThreeVector &p, const G4ThreeVector &v, G4double s[2]) 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
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 

Protected Attributes

G4double dx
 
G4double dy
 
G4double dz
 
- Protected Attributes inherited from G4VSolid
G4double kCarTolerance
 

Detailed Description

Definition at line 55 of file G4EllipticalTube.hh.

Constructor & Destructor Documentation

◆ G4EllipticalTube() [1/3]

G4EllipticalTube::G4EllipticalTube ( const G4String name,
G4double  theDx,
G4double  theDy,
G4double  theDz 
)

Definition at line 60 of file G4EllipticalTube.cc.

64 : G4VSolid( name ), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
65{
66 dx = theDx;
67 dy = theDy;
68 dz = theDz;
69}

◆ ~G4EllipticalTube()

G4EllipticalTube::~G4EllipticalTube ( )
virtual

Definition at line 86 of file G4EllipticalTube.cc.

87{
88 delete fpPolyhedron;
89}

◆ G4EllipticalTube() [2/3]

G4EllipticalTube::G4EllipticalTube ( __void__ &  a)

Definition at line 76 of file G4EllipticalTube.cc.

77 : G4VSolid(a), dx(0.), dy(0.), dz(0.),
78 fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
79{
80}

◆ G4EllipticalTube() [3/3]

G4EllipticalTube::G4EllipticalTube ( const G4EllipticalTube rhs)

Definition at line 95 of file G4EllipticalTube.cc.

96 : G4VSolid(rhs), dx(rhs.dx), dy(rhs.dy), dz(rhs.dz),
97 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
98 fpPolyhedron(0)
99{
100}

Member Function Documentation

◆ CalculateExtent()

G4bool G4EllipticalTube::CalculateExtent ( const EAxis  pAxis,
const G4VoxelLimits pVoxelLimit,
const G4AffineTransform pTransform,
G4double pmin,
G4double pmax 
) const
virtual

Implements G4VSolid.

Definition at line 130 of file G4EllipticalTube.cc.

134{
135 G4SolidExtentList extentList( axis, voxelLimit );
136
137 //
138 // We are going to divide up our elliptical face into small
139 // pieces
140 //
141
142 //
143 // Choose phi size of our segment(s) based on constants as
144 // defined in meshdefs.hh
145 //
146 G4int numPhi = kMaxMeshSections;
147 G4double sigPhi = twopi/numPhi;
148
149 //
150 // We have to be careful to keep our segments completely outside
151 // of the elliptical surface. To do so we imagine we have
152 // a simple (unit radius) circular cross section (as in G4Tubs)
153 // and then "stretch" the dimensions as necessary to fit the ellipse.
154 //
155 G4double rFudge = 1.0/std::cos(0.5*sigPhi);
156 G4double dxFudge = dx*rFudge,
157 dyFudge = dy*rFudge;
158
159 //
160 // As we work around the elliptical surface, we build
161 // a "phi" segment on the way, and keep track of two
162 // additional polygons for the two ends.
163 //
164 G4ClippablePolygon endPoly1, endPoly2, phiPoly;
165
166 G4double phi = 0,
167 cosPhi = std::cos(phi),
168 sinPhi = std::sin(phi);
169 G4ThreeVector v0( dxFudge*cosPhi, dyFudge*sinPhi, +dz ),
170 v1( dxFudge*cosPhi, dyFudge*sinPhi, -dz ),
171 w0, w1;
172 transform.ApplyPointTransform( v0 );
173 transform.ApplyPointTransform( v1 );
174 do
175 {
176 phi += sigPhi;
177 if (numPhi == 1) phi = 0; // Try to avoid roundoff
178 cosPhi = std::cos(phi),
179 sinPhi = std::sin(phi);
180
181 w0 = G4ThreeVector( dxFudge*cosPhi, dyFudge*sinPhi, +dz );
182 w1 = G4ThreeVector( dxFudge*cosPhi, dyFudge*sinPhi, -dz );
183 transform.ApplyPointTransform( w0 );
184 transform.ApplyPointTransform( w1 );
185
186 //
187 // Add a point to our z ends
188 //
189 endPoly1.AddVertexInOrder( v0 );
190 endPoly2.AddVertexInOrder( v1 );
191
192 //
193 // Build phi polygon
194 //
195 phiPoly.ClearAllVertices();
196
197 phiPoly.AddVertexInOrder( v0 );
198 phiPoly.AddVertexInOrder( v1 );
199 phiPoly.AddVertexInOrder( w1 );
200 phiPoly.AddVertexInOrder( w0 );
201
202 if (phiPoly.PartialClip( voxelLimit, axis ))
203 {
204 //
205 // Get unit normal
206 //
207 phiPoly.SetNormal( (v1-v0).cross(w0-v0).unit() );
208
209 extentList.AddSurface( phiPoly );
210 }
211
212 //
213 // Next vertex
214 //
215 v0 = w0;
216 v1 = w1;
217 } while( --numPhi > 0 );
218
219 //
220 // Process the end pieces
221 //
222 if (endPoly1.PartialClip( voxelLimit, axis ))
223 {
224 static const G4ThreeVector normal(0,0,+1);
225 endPoly1.SetNormal( transform.TransformAxis(normal) );
226 extentList.AddSurface( endPoly1 );
227 }
228
229 if (endPoly2.PartialClip( voxelLimit, axis ))
230 {
231 static const G4ThreeVector normal(0,0,-1);
232 endPoly2.SetNormal( transform.TransformAxis(normal) );
233 extentList.AddSurface( endPoly2 );
234 }
235
236 //
237 // Return min/max value
238 //
239 return extentList.GetExtent( min, max );
240}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
virtual G4bool PartialClip(const G4VoxelLimits &voxelLimit, const EAxis IgnoreMe)
virtual void AddVertexInOrder(const G4ThreeVector vertex)
virtual void ClearAllVertices()
void SetNormal(const G4ThreeVector &newNormal)
const G4int kMaxMeshSections
Definition: meshdefs.hh:46

◆ CheckXY() [1/2]

G4double G4EllipticalTube::CheckXY ( const G4double  x,
const G4double  y 
) const
inlineprotected

◆ CheckXY() [2/2]

G4double G4EllipticalTube::CheckXY ( const G4double  x,
const G4double  y,
const G4double  toler 
) const
inlineprotected

◆ Clone()

G4VSolid * G4EllipticalTube::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 858 of file G4EllipticalTube.cc.

859{
860 return new G4EllipticalTube(*this);
861}

◆ CreatePolyhedron()

G4Polyhedron * G4EllipticalTube::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 961 of file G4EllipticalTube.cc.

962{
963 // create cylinder with radius=1...
964 //
965 G4Polyhedron* eTube = new G4PolyhedronTube(0.,1.,dz);
966
967 // apply non-uniform scaling...
968 //
969 eTube->Transform(G4Scale3D(dx,dy,1.));
970 return eTube;
971}
HepGeom::Scale3D G4Scale3D
HepPolyhedron & Transform(const G4Transform3D &t)

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

void G4EllipticalTube::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 993 of file G4EllipticalTube.cc.

994{
995 scene.AddSolid (*this);
996}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 525 of file G4EllipticalTube.cc.

526{
527 static const G4double halfTol = 0.5*kCarTolerance;
528
529 if (CheckXY( p.x(), p.y(), +halfTol ) < 1.0)
530 {
531 //
532 // We are inside or on the surface of the
533 // elliptical cross section in x/y. Check z
534 //
535 if (p.z() < -dz-halfTol)
536 return -p.z()-dz;
537 else if (p.z() > dz+halfTol)
538 return p.z()-dz;
539 else
540 return 0; // On any surface here (or inside)
541 }
542
543 //
544 // Find point on ellipse
545 //
546 G4double qnorm = CheckXY( p.x(), p.y() );
547 if (qnorm < DBL_MIN) return 0; // This should never happen
548
549 G4double q = 1.0/std::sqrt(qnorm);
550
551 G4double xe = q*p.x(), ye = q*p.y();
552
553 //
554 // Get tangent to ellipse
555 //
556 G4double tx = -ye*dx*dx, ty = +xe*dy*dy;
557 G4double tnorm = std::sqrt( tx*tx + ty*ty );
558
559 //
560 // Calculate distance
561 //
562 G4double distR = ( (p.x()-xe)*ty - (p.y()-ye)*tx )/tnorm;
563
564 //
565 // Add the result in quadrature if we are, in addition,
566 // outside the z bounds of the shape
567 //
568 // We could save some time by returning the maximum rather
569 // than the quadrature sum
570 //
571 if (p.z() < -dz)
572 return std::sqrt( (p.z()+dz)*(p.z()+dz) + distR*distR );
573 else if (p.z() > dz)
574 return std::sqrt( (p.z()-dz)*(p.z()-dz) + distR*distR );
575
576 return distR;
577}
double z() const
double x() const
double y() const
G4double CheckXY(const G4double x, const G4double y, const G4double toler) const
G4double kCarTolerance
Definition: G4VSolid.hh:307
#define DBL_MIN
Definition: templates.hh:75

◆ DistanceToIn() [2/2]

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

Implements G4VSolid.

Definition at line 367 of file G4EllipticalTube.cc.

369{
370 static const G4double halfTol = 0.5*kCarTolerance;
371
372 //
373 // Check z = -dz planer surface
374 //
375 G4double sigz = p.z()+dz;
376
377 if (sigz < halfTol)
378 {
379 //
380 // We are "behind" the shape in z, and so can
381 // potentially hit the rear face. Correct direction?
382 //
383 if (v.z() <= 0)
384 {
385 //
386 // As long as we are far enough away, we know we
387 // can't intersect
388 //
389 if (sigz < 0) return kInfinity;
390
391 //
392 // Otherwise, we don't intersect unless we are
393 // on the surface of the ellipse
394 //
395 if (CheckXY(p.x(),p.y(),-halfTol) <= 1.0) return kInfinity;
396 }
397 else
398 {
399 //
400 // How far?
401 //
402 G4double q = -sigz/v.z();
403
404 //
405 // Where does that place us?
406 //
407 G4double xi = p.x() + q*v.x(),
408 yi = p.y() + q*v.y();
409
410 //
411 // Is this on the surface (within ellipse)?
412 //
413 if (CheckXY(xi,yi) <= 1.0)
414 {
415 //
416 // Yup. Return q, unless we are on the surface
417 //
418 return (sigz < -halfTol) ? q : 0;
419 }
420 else if (xi*dy*dy*v.x() + yi*dx*dx*v.y() >= 0)
421 {
422 //
423 // Else, if we are traveling outwards, we know
424 // we must miss
425 //
426 return kInfinity;
427 }
428 }
429 }
430
431 //
432 // Check z = +dz planer surface
433 //
434 sigz = p.z() - dz;
435
436 if (sigz > -halfTol)
437 {
438 if (v.z() >= 0)
439 {
440 if (sigz > 0) return kInfinity;
441 if (CheckXY(p.x(),p.y(),-halfTol) <= 1.0) return kInfinity;
442 }
443 else {
444 G4double q = -sigz/v.z();
445
446 G4double xi = p.x() + q*v.x(),
447 yi = p.y() + q*v.y();
448
449 if (CheckXY(xi,yi) <= 1.0)
450 {
451 return (sigz > -halfTol) ? q : 0;
452 }
453 else if (xi*dy*dy*v.x() + yi*dx*dx*v.y() >= 0)
454 {
455 return kInfinity;
456 }
457 }
458 }
459
460 //
461 // Check intersection with the elliptical tube
462 //
463 G4double q[2];
464 G4int n = IntersectXY( p, v, q );
465
466 if (n==0) return kInfinity;
467
468 //
469 // Is the original point on the surface?
470 //
471 if (std::fabs(p.z()) < dz+halfTol) {
472 if (CheckXY( p.x(), p.y(), halfTol ) < 1.0)
473 {
474 //
475 // Well, yes, but are we traveling inwards at this point?
476 //
477 if (p.x()*dy*dy*v.x() + p.y()*dx*dx*v.y() < 0) return 0;
478 }
479 }
480
481 //
482 // We are now certain that point p is not on the surface of
483 // the solid (and thus std::fabs(q[0]) > halfTol).
484 // Return kInfinity if the intersection is "behind" the point.
485 //
486 if (q[0] < 0) return kInfinity;
487
488 //
489 // Check to see if we intersect the tube within
490 // dz, but only when we know it might miss
491 //
492 G4double zi = p.z() + q[0]*v.z();
493
494 if (v.z() < 0)
495 {
496 if (zi < -dz) return kInfinity;
497 }
498 else if (v.z() > 0)
499 {
500 if (zi > +dz) return kInfinity;
501 }
502
503 return q[0];
504}
G4int IntersectXY(const G4ThreeVector &p, const G4ThreeVector &v, G4double s[2]) const

◆ DistanceToOut() [1/2]

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

Implements G4VSolid.

Definition at line 729 of file G4EllipticalTube.cc.

730{
731 static const G4double halfTol = 0.5*kCarTolerance;
732
733 //
734 // We need to calculate the distances to all surfaces,
735 // and then return the smallest
736 //
737 // Check -dz and +dz surface
738 //
739 G4double sBest = dz - std::fabs(p.z());
740 if (sBest < halfTol) return 0;
741
742 //
743 // Check elliptical surface: find intersection of
744 // line through p and parallel to x axis
745 //
746 G4double radical = 1.0 - p.y()*p.y()/dy/dy;
747 if (radical < +DBL_MIN) return 0;
748
749 G4double xi = dx*std::sqrt( radical );
750 if (p.x() < 0) xi = -xi;
751
752 //
753 // Do the same with y axis
754 //
755 radical = 1.0 - p.x()*p.x()/dx/dx;
756 if (radical < +DBL_MIN) return 0;
757
758 G4double yi = dy*std::sqrt( radical );
759 if (p.y() < 0) yi = -yi;
760
761 //
762 // Get distance from p to the line connecting
763 // these two points
764 //
765 G4double xdi = p.x() - xi,
766 ydi = yi - p.y();
767
768 G4double normi = std::sqrt( xdi*xdi + ydi*ydi );
769 if (normi < halfTol) return 0;
770 xdi /= normi;
771 ydi /= normi;
772
773 G4double q = 0.5*(xdi*(p.y()-yi) - ydi*(p.x()-xi));
774 if (xi*yi < 0) q = -q;
775
776 if (q < sBest) sBest = q;
777
778 //
779 // Return best answer
780 //
781 return sBest < halfTol ? 0 : sBest;
782}

◆ DistanceToOut() [2/2]

G4double G4EllipticalTube::DistanceToOut ( const G4ThreeVector p,
const G4ThreeVector v,
const G4bool  calcNorm = false,
G4bool validNorm = 0,
G4ThreeVector n = 0 
) const
virtual

Implements G4VSolid.

Definition at line 588 of file G4EllipticalTube.cc.

593{
594 static const G4double halfTol = 0.5*kCarTolerance;
595
596 //
597 // Our normal is always valid
598 //
599 if (calcNorm) { *validNorm = true; }
600
601 G4double sBest = kInfinity;
602 G4ThreeVector nBest(0,0,0);
603
604 //
605 // Might we intersect the -dz surface?
606 //
607 if (v.z() < 0)
608 {
609 static const G4ThreeVector normHere(0.0,0.0,-1.0);
610 //
611 // Yup. What distance?
612 //
613 sBest = -(p.z()+dz)/v.z();
614
615 //
616 // Are we on the surface? If so, return zero
617 //
618 if (p.z() < -dz+halfTol)
619 {
620 if (calcNorm) { *norm = normHere; }
621 return 0;
622 }
623 else
624 {
625 nBest = normHere;
626 }
627 }
628
629 //
630 // How about the +dz surface?
631 //
632 if (v.z() > 0)
633 {
634 static const G4ThreeVector normHere(0.0,0.0,+1.0);
635 //
636 // Yup. What distance?
637 //
638 G4double q = (dz-p.z())/v.z();
639
640 //
641 // Are we on the surface? If so, return zero
642 //
643 if (p.z() > +dz-halfTol)
644 {
645 if (calcNorm) { *norm = normHere; }
646 return 0;
647 }
648
649 //
650 // Best so far?
651 //
652 if (q < sBest) { sBest = q; nBest = normHere; }
653 }
654
655 //
656 // Check furthest intersection with ellipse
657 //
658 G4double q[2];
659 G4int n = IntersectXY( p, v, q );
660
661 if (n == 0)
662 {
663 if (sBest == kInfinity)
664 {
665 DumpInfo();
666 std::ostringstream message;
667 G4int oldprc = message.precision(16) ;
668 message << "Point p is outside !?" << G4endl
669 << "Position:" << G4endl
670 << " p.x() = " << p.x()/mm << " mm" << G4endl
671 << " p.y() = " << p.y()/mm << " mm" << G4endl
672 << " p.z() = " << p.z()/mm << " mm" << G4endl
673 << "Direction:" << G4endl << G4endl
674 << " v.x() = " << v.x() << G4endl
675 << " v.y() = " << v.y() << G4endl
676 << " v.z() = " << v.z() << G4endl
677 << "Proposed distance :" << G4endl
678 << " snxt = " << sBest/mm << " mm";
679 message.precision(oldprc) ;
680 G4Exception( "G4EllipticalTube::DistanceToOut(p,v,...)",
681 "GeomSolids1002", JustWarning, message);
682 }
683 if (calcNorm) { *norm = nBest; }
684 return sBest;
685 }
686 else if (q[n-1] > sBest)
687 {
688 if (calcNorm) { *norm = nBest; }
689 return sBest;
690 }
691 sBest = q[n-1];
692
693 //
694 // Intersection with ellipse. Get normal at intersection point.
695 //
696 if (calcNorm)
697 {
698 G4ThreeVector ip = p + sBest*v;
699 *norm = G4ThreeVector( ip.x()*dy*dy, ip.y()*dx*dx, 0.0 ).unit();
700 }
701
702 //
703 // Do we start on the surface?
704 //
705 if (CheckXY( p.x(), p.y(), -halfTol ) > 1.0)
706 {
707 //
708 // Well, yes, but are we traveling outwards at this point?
709 //
710 if (p.x()*dy*dy*v.x() + p.y()*dx*dx*v.y() > 0) return 0;
711 }
712
713 return sBest;
714}
@ JustWarning
#define G4endl
Definition: G4ios.hh:52
Hep3Vector unit() const
void DumpInfo() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ GetCubicVolume()

G4double G4EllipticalTube::GetCubicVolume ( )
virtual

Reimplemented from G4VSolid.

Definition at line 867 of file G4EllipticalTube.cc.

868{
869 if(fCubicVolume != 0.) {;}
870 else { fCubicVolume = G4VSolid::GetCubicVolume(); }
871 return fCubicVolume;
872}
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188

◆ GetDx()

G4double G4EllipticalTube::GetDx ( ) const
inline

◆ GetDy()

G4double G4EllipticalTube::GetDy ( ) const
inline

◆ GetDz()

G4double G4EllipticalTube::GetDz ( ) const
inline

◆ GetEntityType()

G4GeometryType G4EllipticalTube::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 849 of file G4EllipticalTube.cc.

850{
851 return G4String("G4EllipticalTube");
852}

◆ GetExtent()

G4VisExtent G4EllipticalTube::GetExtent ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1002 of file G4EllipticalTube.cc.

1003{
1004 return G4VisExtent( -dx, dx, -dy, dy, -dz, dz );
1005}

◆ GetPointOnSurface()

G4ThreeVector G4EllipticalTube::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 911 of file G4EllipticalTube.cc.

912{
913 G4double xRand, yRand, zRand, phi, cosphi, sinphi, zArea, cArea,p, chose;
914
915 phi = RandFlat::shoot(0., 2.*pi);
916 cosphi = std::cos(phi);
917 sinphi = std::sin(phi);
918
919 // the ellipse perimeter from: "http://mathworld.wolfram.com/Ellipse.html"
920 // m = (dx - dy)/(dx + dy);
921 // k = 1.+1./4.*m*m+1./64.*sqr(m)*sqr(m)+1./256.*sqr(m)*sqr(m)*sqr(m);
922 // p = pi*(a+b)*k;
923
924 // perimeter below from "http://www.efunda.com/math/areas/EllipseGen.cfm"
925
926 p = 2.*pi*std::sqrt(0.5*(dx*dx+dy*dy));
927
928 cArea = 2.*dz*p;
929 zArea = pi*dx*dy;
930
931 xRand = dx*cosphi;
932 yRand = dy*sinphi;
933 zRand = RandFlat::shoot(dz, -1.*dz);
934
935 chose = RandFlat::shoot(0.,2.*zArea+cArea);
936
937 if( (chose>=0) && (chose < cArea) )
938 {
939 return G4ThreeVector (xRand,yRand,zRand);
940 }
941 else if( (chose >= cArea) && (chose < cArea + zArea) )
942 {
943 xRand = RandFlat::shoot(-1.*dx,dx);
944 yRand = std::sqrt(1.-sqr(xRand/dx));
945 yRand = RandFlat::shoot(-1.*yRand, yRand);
946 return G4ThreeVector (xRand,yRand,dz);
947 }
948 else
949 {
950 xRand = RandFlat::shoot(-1.*dx,dx);
951 yRand = std::sqrt(1.-sqr(xRand/dx));
952 yRand = RandFlat::shoot(-1.*yRand, yRand);
953 return G4ThreeVector (xRand,yRand,-1.*dz);
954 }
955}
static double shoot()
Definition: RandFlat.cc:59
const G4double pi
T sqr(const T &x)
Definition: templates.hh:145

◆ GetPolyhedron()

G4Polyhedron * G4EllipticalTube::GetPolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 977 of file G4EllipticalTube.cc.

978{
979 if (!fpPolyhedron ||
981 fpPolyhedron->GetNumberOfRotationSteps())
982 {
983 delete fpPolyhedron;
984 fpPolyhedron = CreatePolyhedron();
985 }
986 return fpPolyhedron;
987}
G4Polyhedron * CreatePolyhedron() const
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
static G4int GetNumberOfRotationSteps()

◆ GetSurfaceArea()

G4double G4EllipticalTube::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 877 of file G4EllipticalTube.cc.

878{
879 if(fSurfaceArea != 0.) {;}
880 else { fSurfaceArea = G4VSolid::GetSurfaceArea(); }
881 return fSurfaceArea;
882}
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:248

◆ Inside()

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

Implements G4VSolid.

Definition at line 250 of file G4EllipticalTube.cc.

251{
252 static const G4double halfTol = 0.5*kCarTolerance;
253
254 //
255 // Check z extents: are we outside?
256 //
257 G4double absZ = std::fabs(p.z());
258 if (absZ > dz+halfTol) return kOutside;
259
260 //
261 // Check x,y: are we outside?
262 //
263 // G4double x = p.x(), y = p.y();
264
265 if (CheckXY(p.x(), p.y(), +halfTol) > 1.0) return kOutside;
266
267 //
268 // We are either inside or on the surface: recheck z extents
269 //
270 if (absZ > dz-halfTol) return kSurface;
271
272 //
273 // Recheck x,y
274 //
275 if (CheckXY(p.x(), p.y(), -halfTol) > 1.0) return kSurface;
276
277 return kInside;
278}
@ kInside
Definition: geomdefs.hh:58
@ kOutside
Definition: geomdefs.hh:58
@ kSurface
Definition: geomdefs.hh:58

◆ IntersectXY()

G4int G4EllipticalTube::IntersectXY ( const G4ThreeVector p,
const G4ThreeVector v,
G4double  s[2] 
) const
protected

Definition at line 810 of file G4EllipticalTube.cc.

813{
814 G4double px = p.x(), py = p.y();
815 G4double vx = v.x(), vy = v.y();
816
817 G4double a = (vx/dx)*(vx/dx) + (vy/dy)*(vy/dy);
818 G4double b = 2.0*( px*vx/dx/dx + py*vy/dy/dy );
819 G4double c = (px/dx)*(px/dx) + (py/dy)*(py/dy) - 1.0;
820
821 if (a < DBL_MIN) return 0; // Trajectory parallel to z axis
822
823 G4double radical = b*b - 4*a*c;
824
825 if (radical < -DBL_MIN) return 0; // No solution
826
827 if (radical < DBL_MIN)
828 {
829 //
830 // Grazes surface
831 //
832 ss[0] = -b/a/2.0;
833 return 1;
834 }
835
836 radical = std::sqrt(radical);
837
838 G4double q = -0.5*( b + (b < 0 ? -radical : +radical) );
839 G4double sa = q/a;
840 G4double sb = c/q;
841 if (sa < sb) { ss[0] = sa; ss[1] = sb; } else { ss[0] = sb; ss[1] = sa; }
842 return 2;
843}

Referenced by DistanceToIn(), and DistanceToOut().

◆ operator=()

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

Definition at line 106 of file G4EllipticalTube.cc.

107{
108 // Check assignment to self
109 //
110 if (this == &rhs) { return *this; }
111
112 // Copy base class data
113 //
115
116 // Copy data
117 //
118 dx = rhs.dx; dy = rhs.dy; dz = rhs.dz;
119 fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
120 fpPolyhedron = 0;
121
122 return *this;
123}
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110

◆ SetDx()

void G4EllipticalTube::SetDx ( const G4double  newDx)
inline

◆ SetDy()

void G4EllipticalTube::SetDy ( const G4double  newDy)
inline

◆ SetDz()

void G4EllipticalTube::SetDz ( const G4double  newDz)
inline

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 887 of file G4EllipticalTube.cc.

888{
889 G4int oldprc = os.precision(16);
890 os << "-----------------------------------------------------------\n"
891 << " *** Dump for solid - " << GetName() << " ***\n"
892 << " ===================================================\n"
893 << " Solid type: G4EllipticalTube\n"
894 << " Parameters: \n"
895 << " length Z: " << dz/mm << " mm \n"
896 << " surface equation in X and Y: \n"
897 << " (X / " << dx << ")^2 + (Y / " << dy << ")^2 = 1 \n"
898 << "-----------------------------------------------------------\n";
899 os.precision(oldprc);
900
901 return os;
902}
G4String GetName() const

◆ SurfaceNormal()

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

Implements G4VSolid.

Definition at line 284 of file G4EllipticalTube.cc.

285{
286 //
287 // SurfaceNormal for the point On the Surface, sum the normals on the Corners
288 //
289 static const G4double halfTol = 0.5*kCarTolerance;
290
291 G4int noSurfaces=0;
292 G4ThreeVector norm, sumnorm(0.,0.,0.);
293
294 G4double distZ = std::fabs(std::fabs(p.z()) - dz);
295
296 G4double distR1 = CheckXY( p.x(), p.y(),+ halfTol );
297 G4double distR2 = CheckXY( p.x(), p.y(),- halfTol );
298
299 if ( (distZ < halfTol ) && ( distR1 <= 1 ) )
300 {
301 noSurfaces++;
302 sumnorm=G4ThreeVector( 0.0, 0.0, p.z() < 0 ? -1.0 : 1.0 );
303 }
304 if( (distR1 <= 1 ) && ( distR2 >= 1 ) )
305 {
306 noSurfaces++;
307 norm= G4ThreeVector( p.x()*dy*dy, p.y()*dx*dx, 0.0 ).unit();
308 sumnorm+=norm;
309 }
310 if ( noSurfaces == 0 )
311 {
312#ifdef G4SPECSDEBUG
313 G4Exception("G4EllipticalTube::SurfaceNormal(p)", "GeomSolids1002",
314 JustWarning, "Point p is not on surface !?" );
315#endif
316 norm = ApproxSurfaceNormal(p);
317 }
318 else if ( noSurfaces == 1 ) { norm = sumnorm; }
319 else { norm = sumnorm.unit(); }
320
321 return norm;
322}

Member Data Documentation

◆ dx

◆ dy

◆ dz


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