Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EllipticalTube Class Reference

#include <G4EllipticalTube.hh>

Inheritance diagram for G4EllipticalTube:
Collaboration diagram for G4EllipticalTube:

Public Member Functions

 G4EllipticalTube (const G4String &name, G4double theDx, G4double theDy, G4double theDz)
 
virtual ~G4EllipticalTube ()
 
void Extent (G4ThreeVector &pMin, G4ThreeVector &pMax) const
 
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 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
 

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
 

Protected Attributes

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

Detailed Description

Definition at line 56 of file G4EllipticalTube.hh.

Constructor & Destructor Documentation

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

Definition at line 67 of file G4EllipticalTube.cc.

71  : G4VSolid( name ), fCubicVolume(0.), fSurfaceArea(0.),
72  fRebuildPolyhedron(false), fpPolyhedron(0)
73 {
74  halfTol = 0.5*kCarTolerance;
75 
76  dx = theDx;
77  dy = theDy;
78  dz = theDz;
79 }
G4VSolid(const G4String &name)
Definition: G4VSolid.cc:61
G4double kCarTolerance
Definition: G4VSolid.hh:307

Here is the caller graph for this function:

G4EllipticalTube::~G4EllipticalTube ( )
virtual

Definition at line 97 of file G4EllipticalTube.cc.

98 {
99  delete fpPolyhedron; fpPolyhedron = 0;
100 }
G4EllipticalTube::G4EllipticalTube ( __void__ &  a)

Definition at line 86 of file G4EllipticalTube.cc.

87  : G4VSolid(a), dx(0.), dy(0.), dz(0.), halfTol(0.),
88  fCubicVolume(0.), fSurfaceArea(0.),
89  fRebuildPolyhedron(false), fpPolyhedron(0)
90 {
91 }
G4VSolid(const G4String &name)
Definition: G4VSolid.cc:61
G4EllipticalTube::G4EllipticalTube ( const G4EllipticalTube rhs)

Definition at line 106 of file G4EllipticalTube.cc.

107  : G4VSolid(rhs), dx(rhs.dx), dy(rhs.dy), dz(rhs.dz), halfTol(rhs.halfTol),
108  fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
109  fRebuildPolyhedron(false), fpPolyhedron(0)
110 {
111 }
G4VSolid(const G4String &name)
Definition: G4VSolid.cc:61

Member Function Documentation

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

Implements G4VSolid.

Definition at line 154 of file G4EllipticalTube.cc.

158 {
159  G4ThreeVector bmin, bmax;
160  G4bool exist;
161 
162  // Check bounding box (bbox)
163  //
164  Extent(bmin,bmax);
165  G4BoundingEnvelope bbox(bmin,bmax);
166 #ifdef G4BBOX_EXTENT
167  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
168 #endif
169  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
170  {
171  return exist = (pMin < pMax) ? true : false;
172  }
173 
174  // Set bounding envelope (benv) and calculate extent
175  //
176  const G4int NSTEPS = 48; // number of steps for whole circle
177  G4double ang = twopi/NSTEPS;
178 
179  G4double sinHalf = std::sin(0.5*ang);
180  G4double cosHalf = std::cos(0.5*ang);
181  G4double sinStep = 2.*sinHalf*cosHalf;
182  G4double cosStep = 1. - 2.*sinHalf*sinHalf;
183  G4double sx = dx/cosHalf;
184  G4double sy = dy/cosHalf;
185 
186  G4double sinCur = sinHalf;
187  G4double cosCur = cosHalf;
188  G4ThreeVectorList baseA(NSTEPS),baseB(NSTEPS);
189  for (G4int k=0; k<NSTEPS; ++k)
190  {
191  baseA[k].set(sx*cosCur,sy*sinCur,-dz);
192  baseB[k].set(sx*cosCur,sy*sinCur, dz);
193 
194  G4double sinTmp = sinCur;
195  sinCur = sinCur*cosStep + cosCur*sinStep;
196  cosCur = cosCur*cosStep - sinTmp*sinStep;
197  }
198 
199  std::vector<const G4ThreeVectorList *> polygons(2);
200  polygons[0] = &baseA;
201  polygons[1] = &baseB;
202  G4BoundingEnvelope benv(bmin,bmax,polygons);
203  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
204  return exist;
205 }
int G4int
Definition: G4Types.hh:78
static constexpr double twopi
Definition: G4SIunits.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
void Extent(G4ThreeVector &pMin, G4ThreeVector &pMax) const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Here is the caller graph for this function:

G4double G4EllipticalTube::CheckXY ( const G4double  x,
const G4double  y 
) const
inlineprotected
G4VSolid * G4EllipticalTube::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 811 of file G4EllipticalTube.cc.

812 {
813  return new G4EllipticalTube(*this);
814 }
G4EllipticalTube(const G4String &name, G4double theDx, G4double theDy, G4double theDz)

Here is the call graph for this function:

G4Polyhedron * G4EllipticalTube::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 914 of file G4EllipticalTube.cc.

915 {
916  // create cylinder with radius=1...
917  //
918  G4Polyhedron* eTube = new G4PolyhedronTube(0.,1.,dz);
919 
920  // apply non-uniform scaling...
921  //
922  eTube->Transform(G4Scale3D(dx,dy,1.));
923  return eTube;
924 }
HepPolyhedron & Transform(const G4Transform3D &t)
HepGeom::Scale3D G4Scale3D

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EllipticalTube::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 950 of file G4EllipticalTube.cc.

951 {
952  scene.AddSolid (*this);
953 }
virtual void AddSolid(const G4Box &)=0

Here is the call graph for this function:

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

Implements G4VSolid.

Definition at line 328 of file G4EllipticalTube.cc.

330 {
331  //
332  // Check z = -dz planer surface
333  //
334  G4double sigz = p.z()+dz;
335 
336  if (sigz < halfTol)
337  {
338  //
339  // We are "behind" the shape in z, and so can
340  // potentially hit the rear face. Correct direction?
341  //
342  if (v.z() <= 0)
343  {
344  //
345  // As long as we are far enough away, we know we
346  // can't intersect
347  //
348  if (sigz < 0) return kInfinity;
349 
350  //
351  // Otherwise, we don't intersect unless we are
352  // on the surface of the ellipse
353  //
354  if (CheckXY(p.x(),p.y(),-halfTol) <= 1.0) return kInfinity;
355  }
356  else
357  {
358  //
359  // How far?
360  //
361  G4double q = -sigz/v.z();
362 
363  //
364  // Where does that place us?
365  //
366  G4double xi = p.x() + q*v.x(),
367  yi = p.y() + q*v.y();
368 
369  //
370  // Is this on the surface (within ellipse)?
371  //
372  if (CheckXY(xi,yi) <= 1.0)
373  {
374  //
375  // Yup. Return q, unless we are on the surface
376  //
377  return (sigz < -halfTol) ? q : 0;
378  }
379  else if (xi*dy*dy*v.x() + yi*dx*dx*v.y() >= 0)
380  {
381  //
382  // Else, if we are traveling outwards, we know
383  // we must miss
384  //
385  return kInfinity;
386  }
387  }
388  }
389 
390  //
391  // Check z = +dz planer surface
392  //
393  sigz = p.z() - dz;
394 
395  if (sigz > -halfTol)
396  {
397  if (v.z() >= 0)
398  {
399  if (sigz > 0) return kInfinity;
400  if (CheckXY(p.x(),p.y(),-halfTol) <= 1.0) return kInfinity;
401  }
402  else {
403  G4double q = -sigz/v.z();
404 
405  G4double xi = p.x() + q*v.x(),
406  yi = p.y() + q*v.y();
407 
408  if (CheckXY(xi,yi) <= 1.0)
409  {
410  return (sigz > -halfTol) ? q : 0;
411  }
412  else if (xi*dy*dy*v.x() + yi*dx*dx*v.y() >= 0)
413  {
414  return kInfinity;
415  }
416  }
417  }
418 
419  //
420  // Check intersection with the elliptical tube
421  //
422  G4double q[2];
423  G4int n = IntersectXY( p, v, q );
424 
425  if (n==0) return kInfinity;
426 
427  //
428  // Is the original point on the surface?
429  //
430  if (std::fabs(p.z()) < dz+halfTol) {
431  if (CheckXY( p.x(), p.y(), halfTol ) < 1.0)
432  {
433  //
434  // Well, yes, but are we traveling inwards at this point?
435  //
436  if (p.x()*dy*dy*v.x() + p.y()*dx*dx*v.y() < 0) return 0;
437  }
438  }
439 
440  //
441  // We are now certain that point p is not on the surface of
442  // the solid (and thus std::fabs(q[0]) > halfTol).
443  // Return kInfinity if the intersection is "behind" the point.
444  //
445  if (q[0] < 0) return kInfinity;
446 
447  //
448  // Check to see if we intersect the tube within
449  // dz, but only when we know it might miss
450  //
451  G4double zi = p.z() + q[0]*v.z();
452 
453  if (v.z() < 0)
454  {
455  if (zi < -dz) return kInfinity;
456  }
457  else if (v.z() > 0)
458  {
459  if (zi > +dz) return kInfinity;
460  }
461 
462  return q[0];
463 }
static const G4double kInfinity
Definition: geomdefs.hh:42
double x() const
int G4int
Definition: G4Types.hh:78
double z() const
G4double CheckXY(const G4double x, const G4double y, const G4double toler) const
G4int IntersectXY(const G4ThreeVector &p, const G4ThreeVector &v, G4double s[2]) const
double y() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VSolid.

Definition at line 484 of file G4EllipticalTube.cc.

485 {
486  if (CheckXY( p.x(), p.y(), +halfTol ) < 1.0)
487  {
488  //
489  // We are inside or on the surface of the
490  // elliptical cross section in x/y. Check z
491  //
492  if (p.z() < -dz-halfTol)
493  return -p.z()-dz;
494  else if (p.z() > dz+halfTol)
495  return p.z()-dz;
496  else
497  return 0; // On any surface here (or inside)
498  }
499 
500  //
501  // Find point on ellipse
502  //
503  G4double qnorm = CheckXY( p.x(), p.y() );
504  if (qnorm < DBL_MIN) return 0; // This should never happen
505 
506  G4double q = 1.0/std::sqrt(qnorm);
507 
508  G4double xe = q*p.x(), ye = q*p.y();
509 
510  //
511  // Get tangent to ellipse
512  //
513  G4double tx = -ye*dx*dx, ty = +xe*dy*dy;
514  G4double tnorm = std::sqrt( tx*tx + ty*ty );
515 
516  //
517  // Calculate distance
518  //
519  G4double distR = ( (p.x()-xe)*ty - (p.y()-ye)*tx )/tnorm;
520 
521  //
522  // Add the result in quadrature if we are, in addition,
523  // outside the z bounds of the shape
524  //
525  // We could save some time by returning the maximum rather
526  // than the quadrature sum
527  //
528  if (p.z() < -dz)
529  return std::sqrt( (p.z()+dz)*(p.z()+dz) + distR*distR );
530  else if (p.z() > dz)
531  return std::sqrt( (p.z()-dz)*(p.z()-dz) + distR*distR );
532 
533  return distR;
534 }
double x() const
double z() const
G4double CheckXY(const G4double x, const G4double y, const G4double toler) const
double y() const
#define DBL_MIN
Definition: templates.hh:75
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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 545 of file G4EllipticalTube.cc.

550 {
551  //
552  // Our normal is always valid
553  //
554  if (calcNorm) { *validNorm = true; }
555 
556  G4double sBest = kInfinity;
557  G4ThreeVector nBest(0,0,0);
558 
559  //
560  // Might we intersect the -dz surface?
561  //
562  if (v.z() < 0)
563  {
564  static const G4ThreeVector normHere(0.0,0.0,-1.0);
565  //
566  // Yup. What distance?
567  //
568  sBest = -(p.z()+dz)/v.z();
569 
570  //
571  // Are we on the surface? If so, return zero
572  //
573  if (p.z() < -dz+halfTol)
574  {
575  if (calcNorm) { *norm = normHere; }
576  return 0;
577  }
578  else
579  {
580  nBest = normHere;
581  }
582  }
583 
584  //
585  // How about the +dz surface?
586  //
587  if (v.z() > 0)
588  {
589  static const G4ThreeVector normHere(0.0,0.0,+1.0);
590  //
591  // Yup. What distance?
592  //
593  G4double q = (dz-p.z())/v.z();
594 
595  //
596  // Are we on the surface? If so, return zero
597  //
598  if (p.z() > +dz-halfTol)
599  {
600  if (calcNorm) { *norm = normHere; }
601  return 0;
602  }
603 
604  //
605  // Best so far?
606  //
607  if (q < sBest) { sBest = q; nBest = normHere; }
608  }
609 
610  //
611  // Check furthest intersection with ellipse
612  //
613  G4double q[2];
614  G4int n = IntersectXY( p, v, q );
615 
616  if (n == 0)
617  {
618  if (sBest == kInfinity)
619  {
620  DumpInfo();
621  std::ostringstream message;
622  G4int oldprc = message.precision(16) ;
623  message << "Point p is outside !?" << G4endl
624  << "Position:" << G4endl
625  << " p.x() = " << p.x()/mm << " mm" << G4endl
626  << " p.y() = " << p.y()/mm << " mm" << G4endl
627  << " p.z() = " << p.z()/mm << " mm" << G4endl
628  << "Direction:" << G4endl << G4endl
629  << " v.x() = " << v.x() << G4endl
630  << " v.y() = " << v.y() << G4endl
631  << " v.z() = " << v.z() << G4endl
632  << "Proposed distance :" << G4endl
633  << " snxt = " << sBest/mm << " mm";
634  message.precision(oldprc) ;
635  G4Exception( "G4EllipticalTube::DistanceToOut(p,v,...)",
636  "GeomSolids1002", JustWarning, message);
637  }
638  if (calcNorm) { *norm = nBest; }
639  return sBest;
640  }
641  else if (q[n-1] > sBest)
642  {
643  if (calcNorm) { *norm = nBest; }
644  return sBest;
645  }
646  sBest = q[n-1];
647 
648  //
649  // Intersection with ellipse. Get normal at intersection point.
650  //
651  if (calcNorm)
652  {
653  G4ThreeVector ip = p + sBest*v;
654  *norm = G4ThreeVector( ip.x()*dy*dy, ip.y()*dx*dx, 0.0 ).unit();
655  }
656 
657  //
658  // Do we start on the surface?
659  //
660  if (CheckXY( p.x(), p.y(), -halfTol ) > 1.0)
661  {
662  //
663  // Well, yes, but are we traveling outwards at this point?
664  //
665  if (p.x()*dy*dy*v.x() + p.y()*dx*dx*v.y() > 0) return 0;
666  }
667 
668  return sBest;
669 }
static constexpr double mm
Definition: G4SIunits.hh:115
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
double x() const
int G4int
Definition: G4Types.hh:78
double z() const
void DumpInfo() const
G4double CheckXY(const G4double x, const G4double y, const G4double toler) const
G4int IntersectXY(const G4ThreeVector &p, const G4ThreeVector &v, G4double s[2]) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double y() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VSolid.

Definition at line 684 of file G4EllipticalTube.cc.

685 {
686  //
687  // We need to calculate the distances to all surfaces,
688  // and then return the smallest
689  //
690  // Check -dz and +dz surface
691  //
692  G4double sBest = dz - std::fabs(p.z());
693  if (sBest < halfTol) return 0;
694 
695  //
696  // Check elliptical surface: find intersection of
697  // line through p and parallel to x axis
698  //
699  G4double radical = 1.0 - p.y()*p.y()/dy/dy;
700  if (radical < +DBL_MIN) return 0;
701 
702  G4double xi = dx*std::sqrt( radical );
703  if (p.x() < 0) xi = -xi;
704 
705  //
706  // Do the same with y axis
707  //
708  radical = 1.0 - p.x()*p.x()/dx/dx;
709  if (radical < +DBL_MIN) return 0;
710 
711  G4double yi = dy*std::sqrt( radical );
712  if (p.y() < 0) yi = -yi;
713 
714  //
715  // Get distance from p to the line connecting
716  // these two points
717  //
718  G4double xdi = p.x() - xi,
719  ydi = yi - p.y();
720 
721  G4double normi = std::sqrt( xdi*xdi + ydi*ydi );
722  if (normi < halfTol) return 0;
723  xdi /= normi;
724  ydi /= normi;
725 
726  G4double q = 0.5*(xdi*(p.y()-yi) - ydi*(p.x()-xi));
727  if (xi*yi < 0) q = -q;
728 
729  if (q < sBest) sBest = q;
730 
731  //
732  // Return best answer
733  //
734  return sBest < halfTol ? 0 : sBest;
735 }
double x() const
double z() const
double y() const
#define DBL_MIN
Definition: templates.hh:75
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4EllipticalTube::Extent ( G4ThreeVector pMin,
G4ThreeVector pMax 
) const
virtual

Reimplemented from G4VSolid.

Definition at line 142 of file G4EllipticalTube.cc.

144 {
145  pMin.set(-dx,-dy,-dz);
146  pMax.set( dx, dy, dz);
147 }
void set(double x, double y, double z)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EllipticalTube::GetCubicVolume ( )
virtual

Reimplemented from G4VSolid.

Definition at line 820 of file G4EllipticalTube.cc.

821 {
822  if(fCubicVolume != 0.) {;}
823  else { fCubicVolume = G4VSolid::GetCubicVolume(); }
824  return fCubicVolume;
825 }
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:189

Here is the call graph for this function:

G4double G4EllipticalTube::GetDx ( ) const
inline

Here is the caller graph for this function:

G4double G4EllipticalTube::GetDy ( ) const
inline

Here is the caller graph for this function:

G4double G4EllipticalTube::GetDz ( ) const
inline

Here is the caller graph for this function:

G4GeometryType G4EllipticalTube::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 802 of file G4EllipticalTube.cc.

803 {
804  return G4String("G4EllipticalTube");
805 }
G4VisExtent G4EllipticalTube::GetExtent ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 959 of file G4EllipticalTube.cc.

960 {
961  return G4VisExtent( -dx, dx, -dy, dy, -dz, dz );
962 }
G4ThreeVector G4EllipticalTube::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 864 of file G4EllipticalTube.cc.

865 {
866  G4double xRand, yRand, zRand, phi, cosphi, sinphi, zArea, cArea,p, chose;
867 
868  phi = G4RandFlat::shoot(0., 2.*pi);
869  cosphi = std::cos(phi);
870  sinphi = std::sin(phi);
871 
872  // the ellipse perimeter from: "http://mathworld.wolfram.com/Ellipse.html"
873  // m = (dx - dy)/(dx + dy);
874  // k = 1.+1./4.*m*m+1./64.*sqr(m)*sqr(m)+1./256.*sqr(m)*sqr(m)*sqr(m);
875  // p = pi*(a+b)*k;
876 
877  // perimeter below from "http://www.efunda.com/math/areas/EllipseGen.cfm"
878 
879  p = 2.*pi*std::sqrt(0.5*(dx*dx+dy*dy));
880 
881  cArea = 2.*dz*p;
882  zArea = pi*dx*dy;
883 
884  xRand = dx*cosphi;
885  yRand = dy*sinphi;
886  zRand = G4RandFlat::shoot(dz, -1.*dz);
887 
888  chose = G4RandFlat::shoot(0.,2.*zArea+cArea);
889 
890  if( (chose>=0) && (chose < cArea) )
891  {
892  return G4ThreeVector (xRand,yRand,zRand);
893  }
894  else if( (chose >= cArea) && (chose < cArea + zArea) )
895  {
896  xRand = G4RandFlat::shoot(-1.*dx,dx);
897  yRand = std::sqrt(1.-sqr(xRand/dx));
898  yRand = G4RandFlat::shoot(-1.*yRand, yRand);
899  return G4ThreeVector (xRand,yRand,dz);
900  }
901  else
902  {
903  xRand = G4RandFlat::shoot(-1.*dx,dx);
904  yRand = std::sqrt(1.-sqr(xRand/dx));
905  yRand = G4RandFlat::shoot(-1.*yRand, yRand);
906  return G4ThreeVector (xRand,yRand,-1.*dz);
907  }
908 }
ThreeVector shoot(const G4int Ap, const G4int Af)
CLHEP::Hep3Vector G4ThreeVector
const char * p
Definition: xmltok.h:285
static constexpr double pi
Definition: G4SIunits.hh:75
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4Polyhedron * G4EllipticalTube::GetPolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 930 of file G4EllipticalTube.cc.

931 {
932  if (!fpPolyhedron ||
933  fRebuildPolyhedron ||
935  fpPolyhedron->GetNumberOfRotationSteps())
936  {
937  G4AutoLock l(&polyhedronMutex);
938  delete fpPolyhedron;
939  fpPolyhedron = CreatePolyhedron();
940  fRebuildPolyhedron = false;
941  l.unlock();
942  }
943  return fpPolyhedron;
944 }
G4Polyhedron * CreatePolyhedron() const
static G4int GetNumberOfRotationSteps()
G4int GetNumberOfRotationStepsAtTimeOfCreation() const

Here is the call graph for this function:

G4double G4EllipticalTube::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 830 of file G4EllipticalTube.cc.

831 {
832  if(fSurfaceArea != 0.) {;}
833  else { fSurfaceArea = G4VSolid::GetSurfaceArea(); }
834  return fSurfaceArea;
835 }
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:251

Here is the call graph for this function:

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

Implements G4VSolid.

Definition at line 214 of file G4EllipticalTube.cc.

215 {
216  //
217  // Check z extents: are we outside?
218  //
219  G4double absZ = std::fabs(p.z());
220  if (absZ > dz+halfTol) return kOutside;
221 
222  //
223  // Check x,y: are we outside?
224  //
225  // G4double x = p.x(), y = p.y();
226 
227  if (CheckXY(p.x(), p.y(), +halfTol) > 1.0) return kOutside;
228 
229  //
230  // We are either inside or on the surface: recheck z extents
231  //
232  if (absZ > dz-halfTol) return kSurface;
233 
234  //
235  // Recheck x,y
236  //
237  if (CheckXY(p.x(), p.y(), -halfTol) > 1.0) return kSurface;
238 
239  return kInside;
240 }
double x() const
double z() const
G4double CheckXY(const G4double x, const G4double y, const G4double toler) const
double y() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Definition at line 763 of file G4EllipticalTube.cc.

766 {
767  G4double px = p.x(), py = p.y();
768  G4double vx = v.x(), vy = v.y();
769 
770  G4double a = (vx/dx)*(vx/dx) + (vy/dy)*(vy/dy);
771  G4double b = 2.0*( px*vx/dx/dx + py*vy/dy/dy );
772  G4double c = (px/dx)*(px/dx) + (py/dy)*(py/dy) - 1.0;
773 
774  if (a < DBL_MIN) return 0; // Trajectory parallel to z axis
775 
776  G4double radical = b*b - 4*a*c;
777 
778  if (radical < -DBL_MIN) return 0; // No solution
779 
780  if (radical < DBL_MIN)
781  {
782  //
783  // Grazes surface
784  //
785  ss[0] = -b/a/2.0;
786  return 1;
787  }
788 
789  radical = std::sqrt(radical);
790 
791  G4double q = -0.5*( b + (b < 0 ? -radical : +radical) );
792  G4double sa = q/a;
793  G4double sb = c/q;
794  if (sa < sb) { ss[0] = sa; ss[1] = sb; } else { ss[0] = sb; ss[1] = sa; }
795  return 2;
796 }
double x() const
double y() const
#define DBL_MIN
Definition: templates.hh:75
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 117 of file G4EllipticalTube.cc.

118 {
119  // Check assignment to self
120  //
121  if (this == &rhs) { return *this; }
122 
123  // Copy base class data
124  //
125  G4VSolid::operator=(rhs);
126 
127  // Copy data
128  //
129  dx = rhs.dx; dy = rhs.dy; dz = rhs.dz;
130  halfTol = rhs.halfTol;
131  fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
132  fRebuildPolyhedron = false;
133  delete fpPolyhedron; fpPolyhedron = 0;
134 
135  return *this;
136 }
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:111

Here is the call graph for this function:

void G4EllipticalTube::SetDx ( const G4double  newDx)
inline
void G4EllipticalTube::SetDy ( const G4double  newDy)
inline
void G4EllipticalTube::SetDz ( const G4double  newDz)
inline
std::ostream & G4EllipticalTube::StreamInfo ( std::ostream &  os) const
virtual

Implements G4VSolid.

Definition at line 840 of file G4EllipticalTube.cc.

841 {
842  G4int oldprc = os.precision(16);
843  os << "-----------------------------------------------------------\n"
844  << " *** Dump for solid - " << GetName() << " ***\n"
845  << " ===================================================\n"
846  << " Solid type: G4EllipticalTube\n"
847  << " Parameters: \n"
848  << " length Z: " << dz/mm << " mm \n"
849  << " surface equation in X and Y: \n"
850  << " (X / " << dx << ")^2 + (Y / " << dy << ")^2 = 1 \n"
851  << "-----------------------------------------------------------\n";
852  os.precision(oldprc);
853 
854  return os;
855 }
G4String GetName() const
static constexpr double mm
Definition: G4SIunits.hh:115
int G4int
Definition: G4Types.hh:78

Here is the call graph for this function:

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

Implements G4VSolid.

Definition at line 246 of file G4EllipticalTube.cc.

247 {
248  //
249  // SurfaceNormal for the point On the Surface, sum the normals on the Corners
250  //
251 
252  G4int noSurfaces=0;
253  G4ThreeVector norm, sumnorm(0.,0.,0.);
254 
255  G4double distZ = std::fabs(std::fabs(p.z()) - dz);
256 
257  G4double distR1 = CheckXY( p.x(), p.y(),+ halfTol );
258  G4double distR2 = CheckXY( p.x(), p.y(),- halfTol );
259 
260  if ( (distZ < halfTol ) && ( distR1 <= 1 ) )
261  {
262  noSurfaces++;
263  sumnorm=G4ThreeVector( 0.0, 0.0, p.z() < 0 ? -1.0 : 1.0 );
264  }
265  if( (distR1 <= 1 ) && ( distR2 >= 1 ) )
266  {
267  noSurfaces++;
268  norm= G4ThreeVector( p.x()*dy*dy, p.y()*dx*dx, 0.0 ).unit();
269  sumnorm+=norm;
270  }
271  if ( noSurfaces == 0 )
272  {
273 #ifdef G4SPECSDEBUG
274  G4Exception("G4EllipticalTube::SurfaceNormal(p)", "GeomSolids1002",
275  JustWarning, "Point p is not on surface !?" );
276 #endif
277  norm = ApproxSurfaceNormal(p);
278  }
279  else if ( noSurfaces == 1 ) { norm = sumnorm; }
280  else { norm = sumnorm.unit(); }
281 
282  return norm;
283 }
CLHEP::Hep3Vector G4ThreeVector
double x() const
int G4int
Definition: G4Types.hh:78
double z() const
G4double CheckXY(const G4double x, const G4double y, const G4double toler) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
Hep3Vector unit() const
double y() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Member Data Documentation

G4double G4EllipticalTube::dx
protected

Definition at line 131 of file G4EllipticalTube.hh.

G4double G4EllipticalTube::dy
protected

Definition at line 131 of file G4EllipticalTube.hh.

G4double G4EllipticalTube::dz
protected

Definition at line 131 of file G4EllipticalTube.hh.


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