G4IntersectingCone Class Reference

#include <G4IntersectingCone.hh>


Public Member Functions

 G4IntersectingCone (const G4double r[2], const G4double z[2])
virtual ~G4IntersectingCone ()
G4int LineHitsCone (const G4ThreeVector &p, const G4ThreeVector &v, G4double *s1, G4double *s2)
G4bool HitOn (const G4double r, const G4double z)
G4double RLo () const
G4double RHi () const
G4double ZLo () const
G4double ZHi () const
 G4IntersectingCone (__void__ &)

Protected Member Functions

G4int LineHitsCone1 (const G4ThreeVector &p, const G4ThreeVector &v, G4double *s1, G4double *s2)
G4int LineHitsCone2 (const G4ThreeVector &p, const G4ThreeVector &v, G4double *s1, G4double *s2)

Protected Attributes

G4double zLo
G4double zHi
G4double rLo
G4double rHi
G4bool type1
G4double A
G4double B


Detailed Description

Definition at line 51 of file G4IntersectingCone.hh.


Constructor & Destructor Documentation

G4IntersectingCone::G4IntersectingCone ( const G4double  r[2],
const G4double  z[2] 
)

Definition at line 46 of file G4IntersectingCone.cc.

References A, B, G4GeometryTolerance::GetInstance(), G4GeometryTolerance::GetSurfaceTolerance(), rHi, rLo, type1, zHi, and zLo.

00048 { 
00049   static const G4double halfCarTolerance
00050     = 0.5 * G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
00051 
00052   //
00053   // What type of cone are we?
00054   //
00055   type1 = (std::fabs(z[1]-z[0]) > std::fabs(r[1]-r[0]));
00056   
00057   if (type1)
00058   {
00059     B = (r[1]-r[0])/(z[1]-z[0]);      // tube like
00060     A = 0.5*( r[1]+r[0] - B*(z[1]+z[0]) );
00061   }
00062   else
00063   {
00064     B = (z[1]-z[0])/(r[1]-r[0]);      // disk like
00065     A = 0.5*( z[1]+z[0] - B*(r[1]+r[0]) );
00066   }
00067   //
00068   // Calculate extent
00069   //
00070   if (r[0] < r[1])
00071   {
00072     rLo = r[0]-halfCarTolerance; rHi = r[1]+halfCarTolerance;
00073   }
00074   else
00075   {
00076     rLo = r[1]-halfCarTolerance; rHi = r[0]+halfCarTolerance;
00077   }
00078   
00079   if (z[0] < z[1])
00080   {
00081     zLo = z[0]-halfCarTolerance; zHi = z[1]+halfCarTolerance;
00082   }
00083   else
00084   {
00085     zLo = z[1]-halfCarTolerance; zHi = z[0]+halfCarTolerance;
00086   }
00087 }

G4IntersectingCone::~G4IntersectingCone (  )  [virtual]

Definition at line 103 of file G4IntersectingCone.cc.

00104 {
00105 }

G4IntersectingCone::G4IntersectingCone ( __void__ &   ) 

Definition at line 94 of file G4IntersectingCone.cc.

00095   : zLo(0.), zHi(0.), rLo(0.), rHi(0.), type1(false), A(0.), B(0.)
00096 {
00097 }


Member Function Documentation

G4bool G4IntersectingCone::HitOn ( const G4double  r,
const G4double  z 
)

Definition at line 114 of file G4IntersectingCone.cc.

References rHi, type1, and zHi.

Referenced by G4PolyconeSide::PointOnCone().

00116 {
00117   //
00118   // Be careful! The inequalities cannot be "<=" and ">=" here without
00119   // punching a tiny hole in our shape!
00120   //
00121   if (type1)
00122   {
00123     if (z < zLo || z > zHi) return false;
00124   }
00125   else
00126   {
00127     if (r < rLo || r > rHi) return false;
00128   }
00129 
00130   return true;
00131 }

G4int G4IntersectingCone::LineHitsCone ( const G4ThreeVector p,
const G4ThreeVector v,
G4double s1,
G4double s2 
)

Definition at line 140 of file G4IntersectingCone.cc.

References LineHitsCone1(), LineHitsCone2(), and type1.

Referenced by G4PolyconeSide::Intersect(), and G4PolyhedraSide::LineHitsSegments().

00143 {
00144   if (type1)
00145   {
00146     return LineHitsCone1( p, v, s1, s2 );
00147   }
00148   else
00149   {
00150     return LineHitsCone2( p, v, s1, s2 );
00151   }
00152 }

G4int G4IntersectingCone::LineHitsCone1 ( const G4ThreeVector p,
const G4ThreeVector v,
G4double s1,
G4double s2 
) [protected]

Definition at line 214 of file G4IntersectingCone.cc.

References A, B, and sqr().

Referenced by LineHitsCone().

00217 {
00218   G4double x0 = p.x(), y0 = p.y(), z0 = p.z();
00219   G4double tx = v.x(), ty = v.y(), tz = v.z();
00220 
00221   G4double a = tx*tx + ty*ty - sqr(B*tz);
00222   G4double b = 2*( x0*tx + y0*ty - (A*B + B*B*z0)*tz);
00223   G4double c = x0*x0 + y0*y0 - sqr(A + B*z0);
00224   
00225   G4double radical = b*b - 4*a*c;
00226  
00227   if (radical < -1E-6*std::fabs(b))  { return 0; }    // No solution
00228   
00229   if (radical < 1E-6*std::fabs(b))
00230   {
00231     //
00232     // The radical is roughly zero: check for special, very rare, cases
00233     //
00234     if (std::fabs(a) > 1/kInfinity)
00235       {
00236       if(B==0.) { return 0; }
00237       if ( std::fabs(x0*ty - y0*tx) < std::fabs(1E-6/B) )
00238       {
00239          *s1 = -0.5*b/a;
00240          return 1;
00241       }
00242       return 0;
00243     }
00244   }
00245   else
00246   {
00247     radical = std::sqrt(radical);
00248   }
00249   
00250   if (a > 1/kInfinity)
00251   {
00252     G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) );
00253     sa = q/a;
00254     sb = c/q;
00255     if (sa < sb) { *s1 = sa; *s2 = sb; } else { *s1 = sb; *s2 = sa; }
00256     if (A + B*(z0+(*s1)*tz) < 0)  { return 0; }
00257     return 2;
00258   }
00259   else if (a < -1/kInfinity)
00260   {
00261     G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) );
00262     sa = q/a;
00263     sb = c/q;
00264     *s1 = (B*tz > 0)^(sa > sb) ? sb : sa;
00265     return 1;
00266   }
00267   else if (std::fabs(b) < 1/kInfinity)
00268   {
00269     return 0;
00270   }
00271   else
00272   {
00273     *s1 = -c/b;
00274     if (A + B*(z0+(*s1)*tz) < 0)  { return 0; }
00275     return 1;
00276   }
00277 }

G4int G4IntersectingCone::LineHitsCone2 ( const G4ThreeVector p,
const G4ThreeVector v,
G4double s1,
G4double s2 
) [protected]

Definition at line 304 of file G4IntersectingCone.cc.

References A, B, and sqr().

Referenced by LineHitsCone().

00307 {
00308   G4double x0 = p.x(), y0 = p.y(), z0 = p.z();
00309   G4double tx = v.x(), ty = v.y(), tz = v.z();
00310   
00311   
00312   // Special case which might not be so rare: B = 0 (precisely)
00313   //
00314   if (B==0)
00315   {
00316     if (std::fabs(tz) < 1/kInfinity)  { return 0; }
00317     
00318     *s1 = (A-z0)/tz;
00319     return 1;
00320   }
00321 
00322   G4double B2 = B*B;
00323 
00324   G4double a = tz*tz - B2*(tx*tx + ty*ty);
00325   G4double b = 2*( (z0-A)*tz - B2*(x0*tx + y0*ty) );
00326   G4double c = sqr(z0-A) - B2*( x0*x0 + y0*y0 );
00327   
00328   G4double radical = b*b - 4*a*c;
00329  
00330   if (radical < -1E-6*std::fabs(b)) { return 0; }   // No solution
00331   
00332   if (radical < 1E-6*std::fabs(b))
00333   {
00334     //
00335     // The radical is roughly zero: check for special, very rare, cases
00336     //
00337     if (std::fabs(a) > 1/kInfinity)
00338     {
00339       if ( std::fabs(x0*ty - y0*tx) < std::fabs(1E-6/B) )
00340       {
00341         *s1 = -0.5*b/a;
00342         return 1;
00343       }
00344       return 0;
00345     }
00346   }
00347   else
00348   {
00349     radical = std::sqrt(radical);
00350   }
00351   
00352   if (a < -1/kInfinity)
00353   {
00354     G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) );
00355     sa = q/a;
00356     sb = c/q;
00357     if (sa < sb) { *s1 = sa; *s2 = sb; } else { *s1 = sb; *s2 = sa; }
00358     if ((z0 + (*s1)*tz  - A)/B < 0)  { return 0; }
00359     return 2;
00360   }
00361   else if (a > 1/kInfinity)
00362   {
00363     G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) );
00364     sa = q/a;
00365     sb = c/q;
00366     *s1 = (tz*B > 0)^(sa > sb) ? sb : sa;
00367     return 1;
00368   }
00369   else if (std::fabs(b) < 1/kInfinity)
00370   {
00371     return 0;
00372   }
00373   else
00374   {
00375     *s1 = -c/b;
00376     if ((z0 + (*s1)*tz  - A)/B < 0)  { return 0; }
00377     return 1;
00378   }
00379 }

G4double G4IntersectingCone::RHi (  )  const [inline]

Definition at line 64 of file G4IntersectingCone.hh.

References rHi.

00064 { return rHi; }

G4double G4IntersectingCone::RLo (  )  const [inline]

Definition at line 63 of file G4IntersectingCone.hh.

References rLo.

00063 { return rLo; }

G4double G4IntersectingCone::ZHi (  )  const [inline]

Definition at line 66 of file G4IntersectingCone.hh.

References zHi.

Referenced by G4PolyhedraSide::Extent(), and G4PolyconeSide::Extent().

00066 { return zHi; }

G4double G4IntersectingCone::ZLo (  )  const [inline]

Definition at line 65 of file G4IntersectingCone.hh.

References zLo.

Referenced by G4PolyhedraSide::Extent(), and G4PolyconeSide::Extent().

00065 { return zLo; }


Field Documentation

G4double G4IntersectingCone::A [protected]

Definition at line 83 of file G4IntersectingCone.hh.

Referenced by G4IntersectingCone(), LineHitsCone1(), and LineHitsCone2().

G4double G4IntersectingCone::B [protected]

Definition at line 83 of file G4IntersectingCone.hh.

Referenced by G4IntersectingCone(), LineHitsCone1(), and LineHitsCone2().

G4double G4IntersectingCone::rHi [protected]

Definition at line 78 of file G4IntersectingCone.hh.

Referenced by G4IntersectingCone(), HitOn(), and RHi().

G4double G4IntersectingCone::rLo [protected]

Definition at line 78 of file G4IntersectingCone.hh.

Referenced by G4IntersectingCone(), and RLo().

G4bool G4IntersectingCone::type1 [protected]

Definition at line 81 of file G4IntersectingCone.hh.

Referenced by G4IntersectingCone(), HitOn(), and LineHitsCone().

G4double G4IntersectingCone::zHi [protected]

Definition at line 78 of file G4IntersectingCone.hh.

Referenced by G4IntersectingCone(), HitOn(), and ZHi().

G4double G4IntersectingCone::zLo [protected]

Definition at line 78 of file G4IntersectingCone.hh.

Referenced by G4IntersectingCone(), and ZLo().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:16 2013 for Geant4 by  doxygen 1.4.7