Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
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.

48 {
49  const G4double halfCarTolerance
51 
52  //
53  // What type of cone are we?
54  //
55  type1 = (std::fabs(z[1]-z[0]) > std::fabs(r[1]-r[0]));
56 
57  if (type1)
58  {
59  B = (r[1]-r[0])/(z[1]-z[0]); // tube like
60  A = 0.5*( r[1]+r[0] - B*(z[1]+z[0]) );
61  }
62  else
63  {
64  B = (z[1]-z[0])/(r[1]-r[0]); // disk like
65  A = 0.5*( z[1]+z[0] - B*(r[1]+r[0]) );
66  }
67  //
68  // Calculate extent
69  //
70  if (r[0] < r[1])
71  {
72  rLo = r[0]-halfCarTolerance; rHi = r[1]+halfCarTolerance;
73  }
74  else
75  {
76  rLo = r[1]-halfCarTolerance; rHi = r[0]+halfCarTolerance;
77  }
78 
79  if (z[0] < z[1])
80  {
81  zLo = z[0]-halfCarTolerance; zHi = z[1]+halfCarTolerance;
82  }
83  else
84  {
85  zLo = z[1]-halfCarTolerance; zHi = z[0]+halfCarTolerance;
86  }
87 }
G4double z
Definition: TRTMaterials.hh:39
G4double GetSurfaceTolerance() const
double G4double
Definition: G4Types.hh:76
static G4GeometryTolerance * GetInstance()
G4IntersectingCone::~G4IntersectingCone ( )
virtual

Definition at line 103 of file G4IntersectingCone.cc.

104 {
105 }
G4IntersectingCone::G4IntersectingCone ( __void__ &  )

Definition at line 94 of file G4IntersectingCone.cc.

95  : zLo(0.), zHi(0.), rLo(0.), rHi(0.), type1(false), A(0.), B(0.)
96 {
97 }

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().

116 {
117  //
118  // Be careful! The inequalities cannot be "<=" and ">=" here without
119  // punching a tiny hole in our shape!
120  //
121  if (type1)
122  {
123  if (z < zLo || z > zHi) return false;
124  }
125  else
126  {
127  if (r < rLo || r > rHi) return false;
128  }
129 
130  return true;
131 }
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().

143 {
144  if (type1)
145  {
146  return LineHitsCone1( p, v, s1, s2 );
147  }
148  else
149  {
150  return LineHitsCone2( p, v, s1, s2 );
151  }
152 }
G4int LineHitsCone1(const G4ThreeVector &p, const G4ThreeVector &v, G4double *s1, G4double *s2)
G4int LineHitsCone2(const G4ThreeVector &p, const G4ThreeVector &v, G4double *s1, G4double *s2)
G4int G4IntersectingCone::LineHitsCone1 ( const G4ThreeVector p,
const G4ThreeVector v,
G4double s1,
G4double s2 
)
protected

Definition at line 214 of file G4IntersectingCone.cc.

References test::a, A, test::b, B, test::c, sqr(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), CLHEP::Hep3Vector::z(), and G4InuclParticleNames::z0.

Referenced by LineHitsCone().

217 {
218  G4double x0 = p.x(), y0 = p.y(), z0 = p.z();
219  G4double tx = v.x(), ty = v.y(), tz = v.z();
220 
221  G4double a = tx*tx + ty*ty - sqr(B*tz);
222  G4double b = 2*( x0*tx + y0*ty - (A*B + B*B*z0)*tz);
223  G4double c = x0*x0 + y0*y0 - sqr(A + B*z0);
224 
225  G4double radical = b*b - 4*a*c;
226 
227  if (radical < -1E-6*std::fabs(b)) { return 0; } // No solution
228 
229  if (radical < 1E-6*std::fabs(b))
230  {
231  //
232  // The radical is roughly zero: check for special, very rare, cases
233  //
234  if (std::fabs(a) > 1/kInfinity)
235  {
236  if(B==0.) { return 0; }
237  if ( std::fabs(x0*ty - y0*tx) < std::fabs(1E-6/B) )
238  {
239  *s1 = -0.5*b/a;
240  return 1;
241  }
242  return 0;
243  }
244  }
245  else
246  {
247  radical = std::sqrt(radical);
248  }
249 
250  if (a > 1/kInfinity)
251  {
252  G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) );
253  sa = q/a;
254  sb = c/q;
255  if (sa < sb) { *s1 = sa; *s2 = sb; } else { *s1 = sb; *s2 = sa; }
256  if (A + B*(z0+(*s1)*tz) < 0) { return 0; }
257  return 2;
258  }
259  else if (a < -1/kInfinity)
260  {
261  G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) );
262  sa = q/a;
263  sb = c/q;
264  *s1 = (B*tz > 0)^(sa > sb) ? sb : sa;
265  return 1;
266  }
267  else if (std::fabs(b) < 1/kInfinity)
268  {
269  return 0;
270  }
271  else
272  {
273  *s1 = -c/b;
274  if (A + B*(z0+(*s1)*tz) < 0) { return 0; }
275  return 1;
276  }
277 }
double x() const
double z() const
double y() const
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
G4int G4IntersectingCone::LineHitsCone2 ( const G4ThreeVector p,
const G4ThreeVector v,
G4double s1,
G4double s2 
)
protected

Definition at line 304 of file G4IntersectingCone.cc.

References test::a, A, test::b, B, test::c, sqr(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), CLHEP::Hep3Vector::z(), and G4InuclParticleNames::z0.

Referenced by LineHitsCone().

307 {
308  G4double x0 = p.x(), y0 = p.y(), z0 = p.z();
309  G4double tx = v.x(), ty = v.y(), tz = v.z();
310 
311 
312  // Special case which might not be so rare: B = 0 (precisely)
313  //
314  if (B==0)
315  {
316  if (std::fabs(tz) < 1/kInfinity) { return 0; }
317 
318  *s1 = (A-z0)/tz;
319  return 1;
320  }
321 
322  G4double B2 = B*B;
323 
324  G4double a = tz*tz - B2*(tx*tx + ty*ty);
325  G4double b = 2*( (z0-A)*tz - B2*(x0*tx + y0*ty) );
326  G4double c = sqr(z0-A) - B2*( x0*x0 + y0*y0 );
327 
328  G4double radical = b*b - 4*a*c;
329 
330  if (radical < -1E-6*std::fabs(b)) { return 0; } // No solution
331 
332  if (radical < 1E-6*std::fabs(b))
333  {
334  //
335  // The radical is roughly zero: check for special, very rare, cases
336  //
337  if (std::fabs(a) > 1/kInfinity)
338  {
339  if ( std::fabs(x0*ty - y0*tx) < std::fabs(1E-6/B) )
340  {
341  *s1 = -0.5*b/a;
342  return 1;
343  }
344  return 0;
345  }
346  }
347  else
348  {
349  radical = std::sqrt(radical);
350  }
351 
352  if (a < -1/kInfinity)
353  {
354  G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) );
355  sa = q/a;
356  sb = c/q;
357  if (sa < sb) { *s1 = sa; *s2 = sb; } else { *s1 = sb; *s2 = sa; }
358  if ((z0 + (*s1)*tz - A)/B < 0) { return 0; }
359  return 2;
360  }
361  else if (a > 1/kInfinity)
362  {
363  G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) );
364  sa = q/a;
365  sb = c/q;
366  *s1 = (tz*B > 0)^(sa > sb) ? sb : sa;
367  return 1;
368  }
369  else if (std::fabs(b) < 1/kInfinity)
370  {
371  return 0;
372  }
373  else
374  {
375  *s1 = -c/b;
376  if ((z0 + (*s1)*tz - A)/B < 0) { return 0; }
377  return 1;
378  }
379 }
double x() const
double z() const
double y() const
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
G4double G4IntersectingCone::RHi ( ) const
inline

Definition at line 64 of file G4IntersectingCone.hh.

References rHi.

64 { return rHi; }
G4double G4IntersectingCone::RLo ( ) const
inline

Definition at line 63 of file G4IntersectingCone.hh.

References rLo.

63 { return rLo; }
G4double G4IntersectingCone::ZHi ( ) const
inline

Definition at line 66 of file G4IntersectingCone.hh.

References zHi.

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

66 { return zHi; }
G4double G4IntersectingCone::ZLo ( ) const
inline

Definition at line 65 of file G4IntersectingCone.hh.

References zLo.

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

65 { 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: