Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4IntersectingCone.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// Implementation of G4IntersectingCone, a utility class which calculates
27// the intersection of an arbitrary line with a fixed cone
28//
29// Author: David C. Williams (UCSC), 1998
30// --------------------------------------------------------------------
31
32#include "G4IntersectingCone.hh"
34
35// Constructor
36//
38 const G4double z[2] )
39{
40 const G4double halfCarTolerance
42
43 // What type of cone are we?
44 //
45 type1 = (std::abs(z[1]-z[0]) > std::abs(r[1]-r[0]));
46
47 if (type1) // tube like
48 {
49 B = (r[1] - r[0]) / (z[1] - z[0]);
50 A = (r[0]*z[1] - r[1]*z[0]) / (z[1] -z[0]);
51 }
52 else // disk like
53 {
54 B = (z[1] - z[0]) / (r[1] - r[0]);
55 A = (z[0]*r[1] - z[1]*r[0]) / (r[1] - r[0]);
56 }
57
58 // Calculate extent
59 //
60 rLo = std::min(r[0], r[1]) - halfCarTolerance;
61 rHi = std::max(r[0], r[1]) + halfCarTolerance;
62 zLo = std::min(z[0], z[1]) - halfCarTolerance;
63 zHi = std::max(z[0], z[1]) + halfCarTolerance;
64}
65
66// Fake default constructor - sets only member data and allocates memory
67// for usage restricted to object persistency.
68//
70 : zLo(0.), zHi(0.), rLo(0.), rHi(0.), A(0.), B(0.)
71{
72}
73
74// HitOn
75//
76// Check r or z extent, as appropriate, to see if the point is possibly
77// on the cone.
78//
80 const G4double z )
81{
82 //
83 // Be careful! The inequalities cannot be "<=" and ">=" here without
84 // punching a tiny hole in our shape!
85 //
86 if (type1)
87 {
88 if (z < zLo || z > zHi) { return false; }
89 }
90 else
91 {
92 if (r < rLo || r > rHi) { return false; }
93 }
94
95 return true;
96}
97
98// LineHitsCone
99//
100// Calculate the intersection of a line with our conical surface, ignoring
101// any phi division
102//
104 const G4ThreeVector& v,
105 G4double* s1, G4double* s2 )
106{
107 if (type1)
108 {
109 return LineHitsCone1( p, v, s1, s2 );
110 }
111 return LineHitsCone2( p, v, s1, s2 );
112}
113
114// LineHitsCone1
115//
116// Calculate the intersections of a line with a conical surface. Only
117// suitable if zPlane[0] != zPlane[1].
118//
119// Equation of a line:
120//
121// x = x0 + s*tx y = y0 + s*ty z = z0 + s*tz
122//
123// Equation of a conical surface:
124//
125// x**2 + y**2 = (A + B*z)**2
126//
127// Solution is quadratic:
128//
129// a*s**2 + b*s + c = 0
130//
131// where:
132//
133// a = tx**2 + ty**2 - (B*tz)**2
134//
135// b = 2*( px*vx + py*vy - B*(A + B*pz)*vz )
136//
137// c = x0**2 + y0**2 - (A + B*z0)**2
138//
139// Notice, that if a < 0, this indicates that the two solutions (assuming
140// they exist) are in opposite cones (that is, given z0 = -A/B, one z < z0
141// and the other z > z0). For our shapes, the invalid solution is one
142// which produces A + Bz < 0, or the one where Bz is smallest (most negative).
143// Since Bz = B*s*tz, if B*tz > 0, we want the largest s, otherwise,
144// the smaller.
145//
146// If there are two solutions on one side of the cone, we want to make
147// sure that they are on the "correct" side, that is A + B*z0 + s*B*tz >= 0.
148//
149// If a = 0, we have a linear problem: s = c/b, which again gives one solution.
150// This should be rare.
151//
152// For b*b - 4*a*c = 0, we also have one solution, which is almost always
153// a line just grazing the surface of a the cone, which we want to ignore.
154// However, there are two other, very rare, possibilities:
155// a line intersecting the z axis and either:
156// 1. At the same angle std::atan(B) to just miss one side of the cone, or
157// 2. Intersecting the cone apex (0,0,-A/B)
158// We *don't* want to miss these! How do we identify them? Well, since
159// this case is rare, we can at least swallow a little more CPU than we would
160// normally be comfortable with. Intersection with the z axis means
161// x0*ty - y0*tx = 0. Case (1) means a==0, and we've already dealt with that
162// above. Case (2) means a < 0.
163//
164// Now: x0*tx + y0*ty = 0 in terms of roundoff error. We can write:
165// Delta = x0*tx + y0*ty
166// b = 2*( Delta - B*(A + B*z0)*tz )
167// For:
168// b*b - 4*a*c = epsilon
169// where epsilon is small, then:
170// Delta = epsilon/2/B
171//
172G4int G4IntersectingCone::LineHitsCone1( const G4ThreeVector& p,
173 const G4ThreeVector& v,
174 G4double* s1, G4double* s2 )
175{
176 static const G4double EPS = DBL_EPSILON; // Precision constant,
177 // originally it was 1E-6
178 G4double x0 = p.x(), y0 = p.y(), z0 = p.z();
179 G4double tx = v.x(), ty = v.y(), tz = v.z();
180
181 // Value of radical can be inaccurate due to loss of precision
182 // if to calculate the coefficiets a,b,c like the following:
183 // G4double a = tx*tx + ty*ty - sqr(B*tz);
184 // G4double b = 2*( x0*tx + y0*ty - B*(A + B*z0)*tz);
185 // G4double c = x0*x0 + y0*y0 - sqr(A + B*z0);
186 //
187 // For more accurate calculation of radical the coefficients
188 // are splitted in two components, radial and along z-axis
189 //
190 G4double ar = tx*tx + ty*ty;
191 G4double az = sqr(B*tz);
192 G4double br = 2*(x0*tx + y0*ty);
193 G4double bz = 2*B*(A + B*z0)*tz;
194 G4double cr = x0*x0 + y0*y0;
195 G4double cz = sqr(A + B*z0);
196
197 // Instead radical = b*b - 4*a*c
198 G4double arcz = 4*ar*cz;
199 G4double azcr = 4*az*cr;
200 G4double radical = (br*br - 4*ar*cr) + ((std::max(arcz,azcr) - 2*bz*br) + std::min(arcz,azcr));
201
202 // Find the coefficients
203 G4double a = ar - az;
204 G4double b = br - bz;
205 G4double c = cr - cz;
206
207 if (radical < -EPS*std::fabs(b)) { return 0; } // No solution
208
209 if (radical < EPS*std::fabs(b))
210 {
211 //
212 // The radical is roughly zero: check for special, very rare, cases
213 //
214 if (std::fabs(a) > 1/kInfinity)
215 {
216 if(B==0.) { return 0; }
217 if ( std::fabs(x0*ty - y0*tx) < std::fabs(EPS/B) )
218 {
219 *s1 = -0.5*b/a;
220 return 1;
221 }
222 return 0;
223 }
224 }
225 else
226 {
227 radical = std::sqrt(radical);
228 }
229
230 if (a > 1/kInfinity)
231 {
232 G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) );
233 sa = q/a;
234 sb = c/q;
235 if (sa < sb) { *s1 = sa; *s2 = sb; } else { *s1 = sb; *s2 = sa; }
236 if (A + B*(z0+(*s1)*tz) < 0) { return 0; }
237 return 2;
238 }
239 if (a < -1/kInfinity)
240 {
241 G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) );
242 sa = q/a;
243 sb = c/q;
244 *s1 = (B*tz > 0)^(sa > sb) ? sb : sa;
245 return 1;
246 }
247
248 if (std::fabs(b) < 1/kInfinity) { return 0; }
249
250 *s1 = -c/b;
251 if (A + B*(z0+(*s1)*tz) < 0) { return 0; }
252 return 1;
253}
254
255// LineHitsCone2
256//
257// See comments under LineHitsCone1. In this routine, case2, we have:
258//
259// Z = A + B*R
260//
261// The solution is still quadratic:
262//
263// a = tz**2 - B*B*(tx**2 + ty**2)
264//
265// b = 2*( (z0-A)*tz - B*B*(x0*tx+y0*ty) )
266//
267// c = ( (z0-A)**2 - B*B*(x0**2 + y0**2) )
268//
269// The rest is much the same, except some details.
270//
271// a > 0 now means we intersect only once in the correct hemisphere.
272//
273// a > 0 ? We only want solution which produces R > 0.
274// since R = (z0+s*tz-A)/B, for tz/B > 0, this is the largest s
275// for tz/B < 0, this is the smallest s
276// thus, same as in case 1 ( since sign(tz/B) = sign(tz*B) )
277//
278G4int G4IntersectingCone::LineHitsCone2( const G4ThreeVector& p,
279 const G4ThreeVector& v,
280 G4double* s1, G4double* s2 )
281{
282 static const G4double EPS = DBL_EPSILON; // Precision constant,
283 // originally it was 1E-6
284 G4double x0 = p.x(), y0 = p.y(), z0 = p.z();
285 G4double tx = v.x(), ty = v.y(), tz = v.z();
286
287 // Special case which might not be so rare: B = 0 (precisely)
288 //
289 if (B==0)
290 {
291 if (std::fabs(tz) < 1/kInfinity) { return 0; }
292
293 *s1 = (A-z0)/tz;
294 return 1;
295 }
296
297 // Value of radical can be inaccurate due to loss of precision
298 // if to calculate the coefficiets a,b,c like the following:
299 // G4double a = tz*tz - B2*(tx*tx + ty*ty);
300 // G4double b = 2*( (z0-A)*tz - B2*(x0*tx + y0*ty) );
301 // G4double c = sqr(z0-A) - B2*( x0*x0 + y0*y0 );
302 //
303 // For more accurate calculation of radical the coefficients
304 // are splitted in two components, radial and along z-axis
305 //
306 G4double B2 = B*B;
307
308 G4double az = tz*tz;
309 G4double ar = B2*(tx*tx + ty*ty);
310 G4double bz = 2*(z0-A)*tz;
311 G4double br = 2*B2*(x0*tx + y0*ty);
312 G4double cz = sqr(z0-A);
313 G4double cr = B2*(x0*x0 + y0*y0);
314
315 // Instead radical = b*b - 4*a*c
316 G4double arcz = 4*ar*cz;
317 G4double azcr = 4*az*cr;
318 G4double radical = (br*br - 4*ar*cr) + ((std::max(arcz,azcr) - 2*bz*br) + std::min(arcz,azcr));
319
320 // Find the coefficients
321 G4double a = az - ar;
322 G4double b = bz - br;
323 G4double c = cz - cr;
324
325 if (radical < -EPS*std::fabs(b)) { return 0; } // No solution
326
327 if (radical < EPS*std::fabs(b))
328 {
329 //
330 // The radical is roughly zero: check for special, very rare, cases
331 //
332 if (std::fabs(a) > 1/kInfinity)
333 {
334 if ( std::fabs(x0*ty - y0*tx) < std::fabs(EPS/B) )
335 {
336 *s1 = -0.5*b/a;
337 return 1;
338 }
339 return 0;
340 }
341 }
342 else
343 {
344 radical = std::sqrt(radical);
345 }
346
347 if (a < -1/kInfinity)
348 {
349 G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) );
350 sa = q/a;
351 sb = c/q;
352 if (sa < sb) { *s1 = sa; *s2 = sb; } else { *s1 = sb; *s2 = sa; }
353 if ((z0 + (*s1)*tz - A)/B < 0) { return 0; }
354 return 2;
355 }
356 if (a > 1/kInfinity)
357 {
358 G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) );
359 sa = q/a;
360 sb = c/q;
361 *s1 = (tz*B > 0)^(sa > sb) ? sb : sa;
362 return 1;
363 }
364
365 if (std::fabs(b) < 1/kInfinity) { return 0; }
366
367 *s1 = -c/b;
368 if ((z0 + (*s1)*tz - A)/B < 0) { return 0; }
369 return 1;
370}
G4double B(G4double temperature)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
double z() const
double x() const
double y() const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4bool HitOn(const G4double r, const G4double z)
G4int LineHitsCone(const G4ThreeVector &p, const G4ThreeVector &v, G4double *s1, G4double *s2)
G4IntersectingCone(const G4double r[2], const G4double z[2])
T sqr(const T &x)
Definition templates.hh:128
#define DBL_EPSILON
Definition templates.hh:66