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 | // |
---|
27 | // $Id: G4IntersectingCone.cc,v 1.12 2008/04/28 08:59:47 gcosmo Exp $ |
---|
28 | // GEANT4 tag $Name: geant4-09-02-ref-02 $ |
---|
29 | // |
---|
30 | // |
---|
31 | // -------------------------------------------------------------------- |
---|
32 | // GEANT 4 class source file |
---|
33 | // |
---|
34 | // |
---|
35 | // G4IntersectingCone.cc |
---|
36 | // |
---|
37 | // Implementation of a utility class which calculates the intersection |
---|
38 | // of an arbitrary line with a fixed cone |
---|
39 | // -------------------------------------------------------------------- |
---|
40 | |
---|
41 | #include "G4IntersectingCone.hh" |
---|
42 | #include "G4GeometryTolerance.hh" |
---|
43 | |
---|
44 | // |
---|
45 | // Constructor |
---|
46 | // |
---|
47 | G4IntersectingCone::G4IntersectingCone( const G4double r[2], |
---|
48 | const G4double z[2] ) |
---|
49 | { |
---|
50 | |
---|
51 | |
---|
52 | half_kCarTolerance = 0.5 * G4GeometryTolerance::GetInstance() |
---|
53 | ->GetSurfaceTolerance(); |
---|
54 | // |
---|
55 | // What type of cone are we? |
---|
56 | // |
---|
57 | type1 = (std::fabs(z[1]-z[0]) > std::fabs(r[1]-r[0])); |
---|
58 | |
---|
59 | if (type1) |
---|
60 | { |
---|
61 | B = (r[1]-r[0])/(z[1]-z[0]); // tube like |
---|
62 | A = 0.5*( r[1]+r[0] - B*(z[1]+z[0]) ); |
---|
63 | } |
---|
64 | else |
---|
65 | { |
---|
66 | B = (z[1]-z[0])/(r[1]-r[0]); // disk like |
---|
67 | A = 0.5*( z[1]+z[0] - B*(r[1]+r[0]) ); |
---|
68 | } |
---|
69 | // |
---|
70 | // Calculate extent |
---|
71 | // |
---|
72 | if (r[0] < r[1]) |
---|
73 | { |
---|
74 | rLo = r[0]-half_kCarTolerance; rHi = r[1]+half_kCarTolerance; |
---|
75 | } |
---|
76 | else |
---|
77 | { |
---|
78 | rLo = r[1]-half_kCarTolerance; rHi = r[0]+half_kCarTolerance; |
---|
79 | } |
---|
80 | |
---|
81 | if (z[0] < z[1]) |
---|
82 | { |
---|
83 | zLo = z[0]-half_kCarTolerance; zHi = z[1]+half_kCarTolerance; |
---|
84 | } |
---|
85 | else |
---|
86 | { |
---|
87 | zLo = z[1]-half_kCarTolerance; zHi = z[0]+half_kCarTolerance; |
---|
88 | } |
---|
89 | } |
---|
90 | |
---|
91 | |
---|
92 | // |
---|
93 | // Fake default constructor - sets only member data and allocates memory |
---|
94 | // for usage restricted to object persistency. |
---|
95 | // |
---|
96 | G4IntersectingCone::G4IntersectingCone( __void__& ) |
---|
97 | { |
---|
98 | } |
---|
99 | |
---|
100 | |
---|
101 | // |
---|
102 | // Destructor |
---|
103 | // |
---|
104 | G4IntersectingCone::~G4IntersectingCone() |
---|
105 | { |
---|
106 | } |
---|
107 | |
---|
108 | |
---|
109 | // |
---|
110 | // HitOn |
---|
111 | // |
---|
112 | // Check r or z extent, as appropriate, to see if the point is possibly |
---|
113 | // on the cone. |
---|
114 | // |
---|
115 | G4bool G4IntersectingCone::HitOn( const G4double r, |
---|
116 | const G4double z ) |
---|
117 | { |
---|
118 | // |
---|
119 | // Be careful! The inequalities cannot be "<=" and ">=" here without |
---|
120 | // punching a tiny hole in our shape! |
---|
121 | // |
---|
122 | if (type1) |
---|
123 | { |
---|
124 | if (z < zLo || z > zHi) return false; |
---|
125 | } |
---|
126 | else |
---|
127 | { |
---|
128 | if (r < rLo || r > rHi) return false; |
---|
129 | } |
---|
130 | |
---|
131 | return true; |
---|
132 | } |
---|
133 | |
---|
134 | |
---|
135 | // |
---|
136 | // LineHitsCone |
---|
137 | // |
---|
138 | // Calculate the intersection of a line with our conical surface, ignoring |
---|
139 | // any phi division |
---|
140 | // |
---|
141 | G4int G4IntersectingCone::LineHitsCone( const G4ThreeVector &p, |
---|
142 | const G4ThreeVector &v, |
---|
143 | G4double *s1, G4double *s2 ) |
---|
144 | { |
---|
145 | if (type1) |
---|
146 | { |
---|
147 | return LineHitsCone1( p, v, s1, s2 ); |
---|
148 | } |
---|
149 | else |
---|
150 | { |
---|
151 | return LineHitsCone2( p, v, s1, s2 ); |
---|
152 | } |
---|
153 | } |
---|
154 | |
---|
155 | |
---|
156 | // |
---|
157 | // LineHitsCone1 |
---|
158 | // |
---|
159 | // Calculate the intersections of a line with a conical surface. Only |
---|
160 | // suitable if zPlane[0] != zPlane[1]. |
---|
161 | // |
---|
162 | // Equation of a line: |
---|
163 | // |
---|
164 | // x = x0 + s*tx y = y0 + s*ty z = z0 + s*tz |
---|
165 | // |
---|
166 | // Equation of a conical surface: |
---|
167 | // |
---|
168 | // x**2 + y**2 = (A + B*z)**2 |
---|
169 | // |
---|
170 | // Solution is quadratic: |
---|
171 | // |
---|
172 | // a*s**2 + b*s + c = 0 |
---|
173 | // |
---|
174 | // where: |
---|
175 | // |
---|
176 | // a = x0**2 + y0**2 - (A + B*z0)**2 |
---|
177 | // |
---|
178 | // b = 2*( x0*tx + y0*ty - (A*B - B*B*z0)*tz) |
---|
179 | // |
---|
180 | // c = tx**2 + ty**2 - (B*tz)**2 |
---|
181 | // |
---|
182 | // Notice, that if a < 0, this indicates that the two solutions (assuming |
---|
183 | // they exist) are in opposite cones (that is, given z0 = -A/B, one z < z0 |
---|
184 | // and the other z > z0). For our shapes, the invalid solution is one |
---|
185 | // which produces A + Bz < 0, or the one where Bz is smallest (most negative). |
---|
186 | // Since Bz = B*s*tz, if B*tz > 0, we want the largest s, otherwise, |
---|
187 | // the smaller. |
---|
188 | // |
---|
189 | // If there are two solutions on one side of the cone, we want to make |
---|
190 | // sure that they are on the "correct" side, that is A + B*z0 + s*B*tz >= 0. |
---|
191 | // |
---|
192 | // If a = 0, we have a linear problem: s = c/b, which again gives one solution. |
---|
193 | // This should be rare. |
---|
194 | // |
---|
195 | // For b*b - 4*a*c = 0, we also have one solution, which is almost always |
---|
196 | // a line just grazing the surface of a the cone, which we want to ignore. |
---|
197 | // However, there are two other, very rare, possibilities: |
---|
198 | // a line intersecting the z axis and either: |
---|
199 | // 1. At the same angle std::atan(B) to just miss one side of the cone, or |
---|
200 | // 2. Intersecting the cone apex (0,0,-A/B) |
---|
201 | // We *don't* want to miss these! How do we identify them? Well, since |
---|
202 | // this case is rare, we can at least swallow a little more CPU than we would |
---|
203 | // normally be comfortable with. Intersection with the z axis means |
---|
204 | // x0*ty - y0*tx = 0. Case (1) means a==0, and we've already dealt with that |
---|
205 | // above. Case (2) means a < 0. |
---|
206 | // |
---|
207 | // Now: x0*tx + y0*ty = 0 in terms of roundoff error. We can write: |
---|
208 | // Delta = x0*tx + y0*ty |
---|
209 | // b = 2*( Delta - (A*B + B*B*z0)*tz ) |
---|
210 | // For: |
---|
211 | // b*b - 4*a*c = epsilon |
---|
212 | // where epsilon is small, then: |
---|
213 | // Delta = epsilon/2/B |
---|
214 | // |
---|
215 | G4int G4IntersectingCone::LineHitsCone1( const G4ThreeVector &p, |
---|
216 | const G4ThreeVector &v, |
---|
217 | G4double *s1, G4double *s2 ) |
---|
218 | { |
---|
219 | G4double x0 = p.x(), y0 = p.y(), z0 = p.z(); |
---|
220 | G4double tx = v.x(), ty = v.y(), tz = v.z(); |
---|
221 | |
---|
222 | G4double a = tx*tx + ty*ty - sqr(B*tz); |
---|
223 | G4double b = 2*( x0*tx + y0*ty - (A*B + B*B*z0)*tz); |
---|
224 | G4double c = x0*x0 + y0*y0 - sqr(A + B*z0); |
---|
225 | |
---|
226 | G4double radical = b*b - 4*a*c; |
---|
227 | |
---|
228 | if (radical < -1E-6*std::fabs(b)) { return 0; } // No solution |
---|
229 | |
---|
230 | if (radical < 1E-6*std::fabs(b)) |
---|
231 | { |
---|
232 | // |
---|
233 | // The radical is roughly zero: check for special, very rare, cases |
---|
234 | // |
---|
235 | if (std::fabs(a) > 1/kInfinity) |
---|
236 | { |
---|
237 | if(B==0.) { return 0; } |
---|
238 | if ( std::fabs(x0*ty - y0*tx) < std::fabs(1E-6/B) ) |
---|
239 | { |
---|
240 | *s1 = -0.5*b/a; |
---|
241 | return 1; |
---|
242 | } |
---|
243 | return 0; |
---|
244 | } |
---|
245 | } |
---|
246 | else |
---|
247 | { |
---|
248 | radical = std::sqrt(radical); |
---|
249 | } |
---|
250 | |
---|
251 | if (a > 1/kInfinity) |
---|
252 | { |
---|
253 | G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) ); |
---|
254 | sa = q/a; |
---|
255 | sb = c/q; |
---|
256 | if (sa < sb) { *s1 = sa; *s2 = sb; } else { *s1 = sb; *s2 = sa; } |
---|
257 | if (A + B*(z0+(*s1)*tz) < 0) { return 0; } |
---|
258 | return 2; |
---|
259 | } |
---|
260 | else if (a < -1/kInfinity) |
---|
261 | { |
---|
262 | G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) ); |
---|
263 | sa = q/a; |
---|
264 | sb = c/q; |
---|
265 | *s1 = (B*tz > 0)^(sa > sb) ? sb : sa; |
---|
266 | return 1; |
---|
267 | } |
---|
268 | else if (std::fabs(b) < 1/kInfinity) |
---|
269 | { |
---|
270 | return 0; |
---|
271 | } |
---|
272 | else |
---|
273 | { |
---|
274 | *s1 = -c/b; |
---|
275 | if (A + B*(z0+(*s1)*tz) < 0) { return 0; } |
---|
276 | return 1; |
---|
277 | } |
---|
278 | } |
---|
279 | |
---|
280 | |
---|
281 | // |
---|
282 | // LineHitsCone2 |
---|
283 | // |
---|
284 | // See comments under LineHitsCone1. In this routine, case2, we have: |
---|
285 | // |
---|
286 | // Z = A + B*R |
---|
287 | // |
---|
288 | // The solution is still quadratic: |
---|
289 | // |
---|
290 | // a = tz**2 - B*B*(tx**2 + ty**2) |
---|
291 | // |
---|
292 | // b = 2*( (z0-A)*tz - B*B*(x0*tx+y0*ty) ) |
---|
293 | // |
---|
294 | // c = ( (z0-A)**2 - B*B*(x0**2 + y0**2) ) |
---|
295 | // |
---|
296 | // The rest is much the same, except some details. |
---|
297 | // |
---|
298 | // a > 0 now means we intersect only once in the correct hemisphere. |
---|
299 | // |
---|
300 | // a > 0 ? We only want solution which produces R > 0. |
---|
301 | // since R = (z0+s*tz-A)/B, for tz/B > 0, this is the largest s |
---|
302 | // for tz/B < 0, this is the smallest s |
---|
303 | // thus, same as in case 1 ( since sign(tz/B) = sign(tz*B) ) |
---|
304 | // |
---|
305 | G4int G4IntersectingCone::LineHitsCone2( const G4ThreeVector &p, |
---|
306 | const G4ThreeVector &v, |
---|
307 | G4double *s1, G4double *s2 ) |
---|
308 | { |
---|
309 | G4double x0 = p.x(), y0 = p.y(), z0 = p.z(); |
---|
310 | G4double tx = v.x(), ty = v.y(), tz = v.z(); |
---|
311 | |
---|
312 | |
---|
313 | // Special case which might not be so rare: B = 0 (precisely) |
---|
314 | // |
---|
315 | if (B==0) |
---|
316 | { |
---|
317 | if (std::fabs(tz) < 1/kInfinity) { return 0; } |
---|
318 | |
---|
319 | *s1 = (A-z0)/tz; |
---|
320 | return 1; |
---|
321 | } |
---|
322 | |
---|
323 | G4double B2 = B*B; |
---|
324 | |
---|
325 | G4double a = tz*tz - B2*(tx*tx + ty*ty); |
---|
326 | G4double b = 2*( (z0-A)*tz - B2*(x0*tx + y0*ty) ); |
---|
327 | G4double c = sqr(z0-A) - B2*( x0*x0 + y0*y0 ); |
---|
328 | |
---|
329 | G4double radical = b*b - 4*a*c; |
---|
330 | |
---|
331 | if (radical < -1E-6*std::fabs(b)) { return 0; } // No solution |
---|
332 | |
---|
333 | if (radical < 1E-6*std::fabs(b)) |
---|
334 | { |
---|
335 | // |
---|
336 | // The radical is roughly zero: check for special, very rare, cases |
---|
337 | // |
---|
338 | if (std::fabs(a) > 1/kInfinity) |
---|
339 | { |
---|
340 | if ( std::fabs(x0*ty - y0*tx) < std::fabs(1E-6/B) ) |
---|
341 | { |
---|
342 | *s1 = -0.5*b/a; |
---|
343 | return 1; |
---|
344 | } |
---|
345 | return 0; |
---|
346 | } |
---|
347 | } |
---|
348 | else |
---|
349 | { |
---|
350 | radical = std::sqrt(radical); |
---|
351 | } |
---|
352 | |
---|
353 | if (a < -1/kInfinity) |
---|
354 | { |
---|
355 | G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) ); |
---|
356 | sa = q/a; |
---|
357 | sb = c/q; |
---|
358 | if (sa < sb) { *s1 = sa; *s2 = sb; } else { *s1 = sb; *s2 = sa; } |
---|
359 | if ((z0 + (*s1)*tz - A)/B < 0) { return 0; } |
---|
360 | return 2; |
---|
361 | } |
---|
362 | else if (a > 1/kInfinity) |
---|
363 | { |
---|
364 | G4double sa, sb, q = -0.5*( b + (b < 0 ? -radical : +radical) ); |
---|
365 | sa = q/a; |
---|
366 | sb = c/q; |
---|
367 | *s1 = (tz*B > 0)^(sa > sb) ? sb : sa; |
---|
368 | return 1; |
---|
369 | } |
---|
370 | else if (std::fabs(b) < 1/kInfinity) |
---|
371 | { |
---|
372 | return 0; |
---|
373 | } |
---|
374 | else |
---|
375 | { |
---|
376 | *s1 = -c/b; |
---|
377 | if ((z0 + (*s1)*tz - A)/B < 0) { return 0; } |
---|
378 | return 1; |
---|
379 | } |
---|
380 | } |
---|