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: G4EllipticalTube.cc,v 1.27 2006/10/20 13:45:21 gcosmo Exp $ |
---|
28 | // GEANT4 tag $Name: HEAD $ |
---|
29 | // |
---|
30 | // |
---|
31 | // -------------------------------------------------------------------- |
---|
32 | // GEANT 4 class source file |
---|
33 | // |
---|
34 | // |
---|
35 | // G4EllipticalTube.cc |
---|
36 | // |
---|
37 | // Implementation of a CSG volume representing a tube with elliptical cross |
---|
38 | // section (geant3 solid 'ELTU') |
---|
39 | // |
---|
40 | // -------------------------------------------------------------------- |
---|
41 | |
---|
42 | #include "G4EllipticalTube.hh" |
---|
43 | |
---|
44 | #include "G4ClippablePolygon.hh" |
---|
45 | #include "G4AffineTransform.hh" |
---|
46 | #include "G4SolidExtentList.hh" |
---|
47 | #include "G4VoxelLimits.hh" |
---|
48 | #include "meshdefs.hh" |
---|
49 | |
---|
50 | #include "Randomize.hh" |
---|
51 | |
---|
52 | #include "G4VGraphicsScene.hh" |
---|
53 | #include "G4Polyhedron.hh" |
---|
54 | #include "G4VisExtent.hh" |
---|
55 | |
---|
56 | using namespace CLHEP; |
---|
57 | |
---|
58 | // |
---|
59 | // Constructor |
---|
60 | // |
---|
61 | G4EllipticalTube::G4EllipticalTube( const G4String &name, |
---|
62 | G4double theDx, |
---|
63 | G4double theDy, |
---|
64 | G4double theDz ) |
---|
65 | : G4VSolid( name ), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0) |
---|
66 | { |
---|
67 | dx = theDx; |
---|
68 | dy = theDy; |
---|
69 | dz = theDz; |
---|
70 | } |
---|
71 | |
---|
72 | |
---|
73 | // |
---|
74 | // Fake default constructor - sets only member data and allocates memory |
---|
75 | // for usage restricted to object persistency. |
---|
76 | // |
---|
77 | G4EllipticalTube::G4EllipticalTube( __void__& a ) |
---|
78 | : G4VSolid(a), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0) |
---|
79 | { |
---|
80 | } |
---|
81 | |
---|
82 | |
---|
83 | // |
---|
84 | // Destructor |
---|
85 | // |
---|
86 | G4EllipticalTube::~G4EllipticalTube() |
---|
87 | { |
---|
88 | delete fpPolyhedron; |
---|
89 | } |
---|
90 | |
---|
91 | |
---|
92 | // |
---|
93 | // CalculateExtent |
---|
94 | // |
---|
95 | G4bool |
---|
96 | G4EllipticalTube::CalculateExtent( const EAxis axis, |
---|
97 | const G4VoxelLimits &voxelLimit, |
---|
98 | const G4AffineTransform &transform, |
---|
99 | G4double &min, G4double &max ) const |
---|
100 | { |
---|
101 | G4SolidExtentList extentList( axis, voxelLimit ); |
---|
102 | |
---|
103 | // |
---|
104 | // We are going to divide up our elliptical face into small |
---|
105 | // pieces |
---|
106 | // |
---|
107 | |
---|
108 | // |
---|
109 | // Choose phi size of our segment(s) based on constants as |
---|
110 | // defined in meshdefs.hh |
---|
111 | // |
---|
112 | G4int numPhi = kMaxMeshSections; |
---|
113 | G4double sigPhi = twopi/numPhi; |
---|
114 | |
---|
115 | // |
---|
116 | // We have to be careful to keep our segments completely outside |
---|
117 | // of the elliptical surface. To do so we imagine we have |
---|
118 | // a simple (unit radius) circular cross section (as in G4Tubs) |
---|
119 | // and then "stretch" the dimensions as necessary to fit the ellipse. |
---|
120 | // |
---|
121 | G4double rFudge = 1.0/std::cos(0.5*sigPhi); |
---|
122 | G4double dxFudge = dx*rFudge, |
---|
123 | dyFudge = dy*rFudge; |
---|
124 | |
---|
125 | // |
---|
126 | // As we work around the elliptical surface, we build |
---|
127 | // a "phi" segment on the way, and keep track of two |
---|
128 | // additional polygons for the two ends. |
---|
129 | // |
---|
130 | G4ClippablePolygon endPoly1, endPoly2, phiPoly; |
---|
131 | |
---|
132 | G4double phi = 0, |
---|
133 | cosPhi = std::cos(phi), |
---|
134 | sinPhi = std::sin(phi); |
---|
135 | G4ThreeVector v0( dxFudge*cosPhi, dyFudge*sinPhi, +dz ), |
---|
136 | v1( dxFudge*cosPhi, dyFudge*sinPhi, -dz ), |
---|
137 | w0, w1; |
---|
138 | transform.ApplyPointTransform( v0 ); |
---|
139 | transform.ApplyPointTransform( v1 ); |
---|
140 | do |
---|
141 | { |
---|
142 | phi += sigPhi; |
---|
143 | if (numPhi == 1) phi = 0; // Try to avoid roundoff |
---|
144 | cosPhi = std::cos(phi), |
---|
145 | sinPhi = std::sin(phi); |
---|
146 | |
---|
147 | w0 = G4ThreeVector( dxFudge*cosPhi, dyFudge*sinPhi, +dz ); |
---|
148 | w1 = G4ThreeVector( dxFudge*cosPhi, dyFudge*sinPhi, -dz ); |
---|
149 | transform.ApplyPointTransform( w0 ); |
---|
150 | transform.ApplyPointTransform( w1 ); |
---|
151 | |
---|
152 | // |
---|
153 | // Add a point to our z ends |
---|
154 | // |
---|
155 | endPoly1.AddVertexInOrder( v0 ); |
---|
156 | endPoly2.AddVertexInOrder( v1 ); |
---|
157 | |
---|
158 | // |
---|
159 | // Build phi polygon |
---|
160 | // |
---|
161 | phiPoly.ClearAllVertices(); |
---|
162 | |
---|
163 | phiPoly.AddVertexInOrder( v0 ); |
---|
164 | phiPoly.AddVertexInOrder( v1 ); |
---|
165 | phiPoly.AddVertexInOrder( w1 ); |
---|
166 | phiPoly.AddVertexInOrder( w0 ); |
---|
167 | |
---|
168 | if (phiPoly.PartialClip( voxelLimit, axis )) |
---|
169 | { |
---|
170 | // |
---|
171 | // Get unit normal |
---|
172 | // |
---|
173 | phiPoly.SetNormal( (v1-v0).cross(w0-v0).unit() ); |
---|
174 | |
---|
175 | extentList.AddSurface( phiPoly ); |
---|
176 | } |
---|
177 | |
---|
178 | // |
---|
179 | // Next vertex |
---|
180 | // |
---|
181 | v0 = w0; |
---|
182 | v1 = w1; |
---|
183 | } while( --numPhi > 0 ); |
---|
184 | |
---|
185 | // |
---|
186 | // Process the end pieces |
---|
187 | // |
---|
188 | if (endPoly1.PartialClip( voxelLimit, axis )) |
---|
189 | { |
---|
190 | static const G4ThreeVector normal(0,0,+1); |
---|
191 | endPoly1.SetNormal( transform.TransformAxis(normal) ); |
---|
192 | extentList.AddSurface( endPoly1 ); |
---|
193 | } |
---|
194 | |
---|
195 | if (endPoly2.PartialClip( voxelLimit, axis )) |
---|
196 | { |
---|
197 | static const G4ThreeVector normal(0,0,-1); |
---|
198 | endPoly2.SetNormal( transform.TransformAxis(normal) ); |
---|
199 | extentList.AddSurface( endPoly2 ); |
---|
200 | } |
---|
201 | |
---|
202 | // |
---|
203 | // Return min/max value |
---|
204 | // |
---|
205 | return extentList.GetExtent( min, max ); |
---|
206 | } |
---|
207 | |
---|
208 | |
---|
209 | // |
---|
210 | // Inside |
---|
211 | // |
---|
212 | // Note that for this solid, we've decided to define the tolerant |
---|
213 | // surface as that which is bounded by ellipses with axes |
---|
214 | // at +/- 0.5*kCarTolerance. |
---|
215 | // |
---|
216 | EInside G4EllipticalTube::Inside( const G4ThreeVector& p ) const |
---|
217 | { |
---|
218 | static const G4double halfTol = 0.5*kCarTolerance; |
---|
219 | |
---|
220 | // |
---|
221 | // Check z extents: are we outside? |
---|
222 | // |
---|
223 | G4double absZ = std::fabs(p.z()); |
---|
224 | if (absZ > dz+halfTol) return kOutside; |
---|
225 | |
---|
226 | // |
---|
227 | // Check x,y: are we outside? |
---|
228 | // |
---|
229 | // G4double x = p.x(), y = p.y(); |
---|
230 | |
---|
231 | if (CheckXY(p.x(), p.y(), +halfTol) > 1.0) return kOutside; |
---|
232 | |
---|
233 | // |
---|
234 | // We are either inside or on the surface: recheck z extents |
---|
235 | // |
---|
236 | if (absZ > dz-halfTol) return kSurface; |
---|
237 | |
---|
238 | // |
---|
239 | // Recheck x,y |
---|
240 | // |
---|
241 | if (CheckXY(p.x(), p.y(), -halfTol) > 1.0) return kSurface; |
---|
242 | |
---|
243 | return kInside; |
---|
244 | } |
---|
245 | |
---|
246 | |
---|
247 | // |
---|
248 | // SurfaceNormal |
---|
249 | // |
---|
250 | G4ThreeVector G4EllipticalTube::SurfaceNormal( const G4ThreeVector& p ) const |
---|
251 | { |
---|
252 | // |
---|
253 | // Which of the three surfaces are we closest to (approximately)? |
---|
254 | // |
---|
255 | G4double distZ = std::fabs(p.z()) - dz; |
---|
256 | |
---|
257 | G4double rxy = CheckXY( p.x(), p.y() ); |
---|
258 | G4double distR2 = (rxy < DBL_MIN) ? DBL_MAX : 1.0/rxy; |
---|
259 | |
---|
260 | // |
---|
261 | // Closer to z? |
---|
262 | // |
---|
263 | if (distZ*distZ < distR2) |
---|
264 | return G4ThreeVector( 0.0, 0.0, p.z() < 0 ? -1.0 : 1.0 ); |
---|
265 | |
---|
266 | // |
---|
267 | // Closer to x/y |
---|
268 | // |
---|
269 | return G4ThreeVector( p.x()*dy*dy, p.y()*dx*dx, 0.0 ).unit(); |
---|
270 | } |
---|
271 | |
---|
272 | |
---|
273 | // |
---|
274 | // DistanceToIn(p,v) |
---|
275 | // |
---|
276 | // Unlike DistanceToOut(p,v), it is possible for the trajectory |
---|
277 | // to miss. The geometric calculations here are quite simple. |
---|
278 | // More difficult is the logic required to prevent particles |
---|
279 | // from sneaking (or leaking) between the elliptical and end |
---|
280 | // surfaces. |
---|
281 | // |
---|
282 | // Keep in mind that the true distance is allowed to be |
---|
283 | // negative if the point is currently on the surface. For oblique |
---|
284 | // angles, it can be very negative. |
---|
285 | // |
---|
286 | G4double G4EllipticalTube::DistanceToIn( const G4ThreeVector& p, |
---|
287 | const G4ThreeVector& v ) const |
---|
288 | { |
---|
289 | static const G4double halfTol = 0.5*kCarTolerance; |
---|
290 | |
---|
291 | // |
---|
292 | // Check z = -dz planer surface |
---|
293 | // |
---|
294 | G4double sigz = p.z()+dz; |
---|
295 | |
---|
296 | if (sigz < halfTol) |
---|
297 | { |
---|
298 | // |
---|
299 | // We are "behind" the shape in z, and so can |
---|
300 | // potentially hit the rear face. Correct direction? |
---|
301 | // |
---|
302 | if (v.z() <= 0) |
---|
303 | { |
---|
304 | // |
---|
305 | // As long as we are far enough away, we know we |
---|
306 | // can't intersect |
---|
307 | // |
---|
308 | if (sigz < 0) return kInfinity; |
---|
309 | |
---|
310 | // |
---|
311 | // Otherwise, we don't intersect unless we are |
---|
312 | // on the surface of the ellipse |
---|
313 | // |
---|
314 | if (CheckXY(p.x(),p.y(),-halfTol) <= 1.0) return kInfinity; |
---|
315 | } |
---|
316 | else |
---|
317 | { |
---|
318 | // |
---|
319 | // How far? |
---|
320 | // |
---|
321 | G4double s = -sigz/v.z(); |
---|
322 | |
---|
323 | // |
---|
324 | // Where does that place us? |
---|
325 | // |
---|
326 | G4double xi = p.x() + s*v.x(), |
---|
327 | yi = p.y() + s*v.y(); |
---|
328 | |
---|
329 | // |
---|
330 | // Is this on the surface (within ellipse)? |
---|
331 | // |
---|
332 | if (CheckXY(xi,yi) <= 1.0) |
---|
333 | { |
---|
334 | // |
---|
335 | // Yup. Return s, unless we are on the surface |
---|
336 | // |
---|
337 | return (sigz < -halfTol) ? s : 0; |
---|
338 | } |
---|
339 | else if (xi*dy*dy*v.x() + yi*dx*dx*v.y() >= 0) |
---|
340 | { |
---|
341 | // |
---|
342 | // Else, if we are traveling outwards, we know |
---|
343 | // we must miss |
---|
344 | // |
---|
345 | return kInfinity; |
---|
346 | } |
---|
347 | } |
---|
348 | } |
---|
349 | |
---|
350 | // |
---|
351 | // Check z = +dz planer surface |
---|
352 | // |
---|
353 | sigz = p.z() - dz; |
---|
354 | |
---|
355 | if (sigz > -halfTol) |
---|
356 | { |
---|
357 | if (v.z() >= 0) |
---|
358 | { |
---|
359 | if (sigz > 0) return kInfinity; |
---|
360 | if (CheckXY(p.x(),p.y(),-halfTol) <= 1.0) return kInfinity; |
---|
361 | } |
---|
362 | else { |
---|
363 | G4double s = -sigz/v.z(); |
---|
364 | |
---|
365 | G4double xi = p.x() + s*v.x(), |
---|
366 | yi = p.y() + s*v.y(); |
---|
367 | |
---|
368 | if (CheckXY(xi,yi) <= 1.0) |
---|
369 | { |
---|
370 | return (sigz > -halfTol) ? s : 0; |
---|
371 | } |
---|
372 | else if (xi*dy*dy*v.x() + yi*dx*dx*v.y() >= 0) |
---|
373 | { |
---|
374 | return kInfinity; |
---|
375 | } |
---|
376 | } |
---|
377 | } |
---|
378 | |
---|
379 | // |
---|
380 | // Check intersection with the elliptical tube |
---|
381 | // |
---|
382 | G4double s[2]; |
---|
383 | G4int n = IntersectXY( p, v, s ); |
---|
384 | |
---|
385 | if (n==0) return kInfinity; |
---|
386 | |
---|
387 | // |
---|
388 | // Is the original point on the surface? |
---|
389 | // |
---|
390 | if (std::fabs(p.z()) < dz+halfTol) { |
---|
391 | if (CheckXY( p.x(), p.y(), halfTol ) < 1.0) |
---|
392 | { |
---|
393 | // |
---|
394 | // Well, yes, but are we traveling inwards at this point? |
---|
395 | // |
---|
396 | if (p.x()*dy*dy*v.x() + p.y()*dx*dx*v.y() < 0) return 0; |
---|
397 | } |
---|
398 | } |
---|
399 | |
---|
400 | // |
---|
401 | // We are now certain that point p is not on the surface of |
---|
402 | // the solid (and thus std::fabs(s[0]) > halfTol). |
---|
403 | // Return kInfinity if the intersection is "behind" the point. |
---|
404 | // |
---|
405 | if (s[0] < 0) return kInfinity; |
---|
406 | |
---|
407 | // |
---|
408 | // Check to see if we intersect the tube within |
---|
409 | // dz, but only when we know it might miss |
---|
410 | // |
---|
411 | G4double zi = p.z() + s[0]*v.z(); |
---|
412 | |
---|
413 | if (v.z() < 0) |
---|
414 | { |
---|
415 | if (zi < -dz) return kInfinity; |
---|
416 | } |
---|
417 | else if (v.z() > 0) |
---|
418 | { |
---|
419 | if (zi > +dz) return kInfinity; |
---|
420 | } |
---|
421 | |
---|
422 | return s[0]; |
---|
423 | } |
---|
424 | |
---|
425 | |
---|
426 | // |
---|
427 | // DistanceToIn(p) |
---|
428 | // |
---|
429 | // The distance from a point to an ellipse (in 2 dimensions) is a |
---|
430 | // surprisingly complicated quadric expression (this is easy to |
---|
431 | // appreciate once one understands that there may be up to |
---|
432 | // four lines normal to the ellipse intersecting any point). To |
---|
433 | // solve it exactly would be rather time consuming. This method, |
---|
434 | // however, is supposed to be a quick check, and is allowed to be an |
---|
435 | // underestimate. |
---|
436 | // |
---|
437 | // So, I will use the following underestimate of the distance |
---|
438 | // from an outside point to an ellipse. First: find the intersection "A" |
---|
439 | // of the line from the origin to the point with the ellipse. |
---|
440 | // Find the line passing through "A" and tangent to the ellipse |
---|
441 | // at A. The distance of the point p from the ellipse will be approximated |
---|
442 | // as the distance to this line. |
---|
443 | // |
---|
444 | G4double G4EllipticalTube::DistanceToIn( const G4ThreeVector& p ) const |
---|
445 | { |
---|
446 | static const G4double halfTol = 0.5*kCarTolerance; |
---|
447 | |
---|
448 | if (CheckXY( p.x(), p.y(), +halfTol ) < 1.0) |
---|
449 | { |
---|
450 | // |
---|
451 | // We are inside or on the surface of the |
---|
452 | // elliptical cross section in x/y. Check z |
---|
453 | // |
---|
454 | if (p.z() < -dz-halfTol) |
---|
455 | return -p.z()-dz; |
---|
456 | else if (p.z() > dz+halfTol) |
---|
457 | return p.z()-dz; |
---|
458 | else |
---|
459 | return 0; // On any surface here (or inside) |
---|
460 | } |
---|
461 | |
---|
462 | // |
---|
463 | // Find point on ellipse |
---|
464 | // |
---|
465 | G4double qnorm = CheckXY( p.x(), p.y() ); |
---|
466 | if (qnorm < DBL_MIN) return 0; // This should never happen |
---|
467 | |
---|
468 | G4double q = 1.0/std::sqrt(qnorm); |
---|
469 | |
---|
470 | G4double xe = q*p.x(), ye = q*p.y(); |
---|
471 | |
---|
472 | // |
---|
473 | // Get tangent to ellipse |
---|
474 | // |
---|
475 | G4double tx = -ye*dx*dx, ty = +xe*dy*dy; |
---|
476 | G4double tnorm = std::sqrt( tx*tx + ty*ty ); |
---|
477 | |
---|
478 | // |
---|
479 | // Calculate distance |
---|
480 | // |
---|
481 | G4double distR = ( (p.x()-xe)*ty - (p.y()-ye)*tx )/tnorm; |
---|
482 | |
---|
483 | // |
---|
484 | // Add the result in quadrature if we are, in addition, |
---|
485 | // outside the z bounds of the shape |
---|
486 | // |
---|
487 | // We could save some time by returning the maximum rather |
---|
488 | // than the quadrature sum |
---|
489 | // |
---|
490 | if (p.z() < -dz) |
---|
491 | return std::sqrt( (p.z()+dz)*(p.z()+dz) + distR*distR ); |
---|
492 | else if (p.z() > dz) |
---|
493 | return std::sqrt( (p.z()-dz)*(p.z()-dz) + distR*distR ); |
---|
494 | |
---|
495 | return distR; |
---|
496 | } |
---|
497 | |
---|
498 | |
---|
499 | // |
---|
500 | // DistanceToOut(p,v) |
---|
501 | // |
---|
502 | // This method can be somewhat complicated for a general shape. |
---|
503 | // For a convex one, like this, there are several simplifications, |
---|
504 | // the most important of which is that one can treat the surfaces |
---|
505 | // as infinite in extent when deciding if the p is on the surface. |
---|
506 | // |
---|
507 | G4double G4EllipticalTube::DistanceToOut( const G4ThreeVector& p, |
---|
508 | const G4ThreeVector& v, |
---|
509 | const G4bool calcNorm, |
---|
510 | G4bool *validNorm, |
---|
511 | G4ThreeVector *norm ) const |
---|
512 | { |
---|
513 | static const G4double halfTol = 0.5*kCarTolerance; |
---|
514 | |
---|
515 | // |
---|
516 | // Our normal is always valid |
---|
517 | // |
---|
518 | if (calcNorm) *validNorm = true; |
---|
519 | |
---|
520 | G4double sBest = kInfinity; |
---|
521 | const G4ThreeVector *nBest=0; |
---|
522 | |
---|
523 | // |
---|
524 | // Might we intersect the -dz surface? |
---|
525 | // |
---|
526 | if (v.z() < 0) |
---|
527 | { |
---|
528 | static const G4ThreeVector normHere(0.0,0.0,-1.0); |
---|
529 | // |
---|
530 | // Yup. What distance? |
---|
531 | // |
---|
532 | sBest = -(p.z()+dz)/v.z(); |
---|
533 | |
---|
534 | // |
---|
535 | // Are we on the surface? If so, return zero |
---|
536 | // |
---|
537 | if (p.z() < -dz+halfTol) { |
---|
538 | if (calcNorm) *norm = normHere; |
---|
539 | return 0; |
---|
540 | } |
---|
541 | else |
---|
542 | nBest = &normHere; |
---|
543 | } |
---|
544 | |
---|
545 | // |
---|
546 | // How about the +dz surface? |
---|
547 | // |
---|
548 | if (v.z() > 0) |
---|
549 | { |
---|
550 | static const G4ThreeVector normHere(0.0,0.0,+1.0); |
---|
551 | // |
---|
552 | // Yup. What distance? |
---|
553 | // |
---|
554 | G4double s = (dz-p.z())/v.z(); |
---|
555 | |
---|
556 | // |
---|
557 | // Are we on the surface? If so, return zero |
---|
558 | // |
---|
559 | if (p.z() > +dz-halfTol) { |
---|
560 | if (calcNorm) *norm = normHere; |
---|
561 | return 0; |
---|
562 | } |
---|
563 | |
---|
564 | // |
---|
565 | // Best so far? |
---|
566 | // |
---|
567 | if (s < sBest) { sBest = s; nBest = &normHere; } |
---|
568 | } |
---|
569 | |
---|
570 | // |
---|
571 | // Check furthest intersection with ellipse |
---|
572 | // |
---|
573 | G4double s[2]; |
---|
574 | G4int n = IntersectXY( p, v, s ); |
---|
575 | |
---|
576 | if (n == 0) |
---|
577 | { |
---|
578 | if (sBest == kInfinity) |
---|
579 | { |
---|
580 | G4cout.precision(16) ; |
---|
581 | G4cout << G4endl ; |
---|
582 | DumpInfo(); |
---|
583 | G4cout << "Position:" << G4endl << G4endl ; |
---|
584 | G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl ; |
---|
585 | G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl ; |
---|
586 | G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl ; |
---|
587 | G4cout << "Direction:" << G4endl << G4endl; |
---|
588 | G4cout << "v.x() = " << v.x() << G4endl; |
---|
589 | G4cout << "v.y() = " << v.y() << G4endl; |
---|
590 | G4cout << "v.z() = " << v.z() << G4endl << G4endl; |
---|
591 | G4cout << "Proposed distance :" << G4endl << G4endl; |
---|
592 | G4cout << "snxt = " << sBest/mm << " mm" << G4endl << G4endl; |
---|
593 | G4Exception( "G4EllipticalTube::DistanceToOut(p,v,...)", |
---|
594 | "Notification", JustWarning, "Point p is outside !?" ); |
---|
595 | } |
---|
596 | if (calcNorm) *norm = *nBest; |
---|
597 | return sBest; |
---|
598 | } |
---|
599 | else if (s[n-1] > sBest) |
---|
600 | { |
---|
601 | if (calcNorm) *norm = *nBest; |
---|
602 | return sBest; |
---|
603 | } |
---|
604 | sBest = s[n-1]; |
---|
605 | |
---|
606 | // |
---|
607 | // Intersection with ellipse. Get normal at intersection point. |
---|
608 | // |
---|
609 | if (calcNorm) |
---|
610 | { |
---|
611 | G4ThreeVector ip = p + sBest*v; |
---|
612 | *norm = G4ThreeVector( ip.x()*dy*dy, ip.y()*dx*dx, 0.0 ).unit(); |
---|
613 | } |
---|
614 | |
---|
615 | // |
---|
616 | // Do we start on the surface? |
---|
617 | // |
---|
618 | if (CheckXY( p.x(), p.y(), -halfTol ) > 1.0) |
---|
619 | { |
---|
620 | // |
---|
621 | // Well, yes, but are we traveling outwards at this point? |
---|
622 | // |
---|
623 | if (p.x()*dy*dy*v.x() + p.y()*dx*dx*v.y() > 0) return 0; |
---|
624 | } |
---|
625 | |
---|
626 | return sBest; |
---|
627 | } |
---|
628 | |
---|
629 | |
---|
630 | // |
---|
631 | // DistanceToOut(p) |
---|
632 | // |
---|
633 | // See DistanceToIn(p) for notes on the distance from a point |
---|
634 | // to an ellipse in two dimensions. |
---|
635 | // |
---|
636 | // The approximation used here for a point inside the ellipse |
---|
637 | // is to find the intersection with the ellipse of the lines |
---|
638 | // through the point and parallel to the x and y axes. The |
---|
639 | // distance of the point from the line connecting the two |
---|
640 | // intersecting points is then used. |
---|
641 | // |
---|
642 | G4double G4EllipticalTube::DistanceToOut( const G4ThreeVector& p ) const |
---|
643 | { |
---|
644 | static const G4double halfTol = 0.5*kCarTolerance; |
---|
645 | |
---|
646 | // |
---|
647 | // We need to calculate the distances to all surfaces, |
---|
648 | // and then return the smallest |
---|
649 | // |
---|
650 | // Check -dz and +dz surface |
---|
651 | // |
---|
652 | G4double sBest = dz - std::fabs(p.z()); |
---|
653 | if (sBest < halfTol) return 0; |
---|
654 | |
---|
655 | // |
---|
656 | // Check elliptical surface: find intersection of |
---|
657 | // line through p and parallel to x axis |
---|
658 | // |
---|
659 | G4double radical = 1.0 - p.y()*p.y()/dy/dy; |
---|
660 | if (radical < +DBL_MIN) return 0; |
---|
661 | |
---|
662 | G4double xi = dx*std::sqrt( radical ); |
---|
663 | if (p.x() < 0) xi = -xi; |
---|
664 | |
---|
665 | // |
---|
666 | // Do the same with y axis |
---|
667 | // |
---|
668 | radical = 1.0 - p.x()*p.x()/dx/dx; |
---|
669 | if (radical < +DBL_MIN) return 0; |
---|
670 | |
---|
671 | G4double yi = dy*std::sqrt( radical ); |
---|
672 | if (p.y() < 0) yi = -yi; |
---|
673 | |
---|
674 | // |
---|
675 | // Get distance from p to the line connecting |
---|
676 | // these two points |
---|
677 | // |
---|
678 | G4double xdi = p.x() - xi, |
---|
679 | ydi = yi - p.y(); |
---|
680 | |
---|
681 | G4double normi = std::sqrt( xdi*xdi + ydi*ydi ); |
---|
682 | if (normi < halfTol) return 0; |
---|
683 | xdi /= normi; |
---|
684 | ydi /= normi; |
---|
685 | |
---|
686 | G4double s = 0.5*(xdi*(p.y()-yi) - ydi*(p.x()-xi)); |
---|
687 | if (xi*yi < 0) s = -s; |
---|
688 | |
---|
689 | if (s < sBest) sBest = s; |
---|
690 | |
---|
691 | // |
---|
692 | // Return best answer |
---|
693 | // |
---|
694 | return sBest < halfTol ? 0 : sBest; |
---|
695 | } |
---|
696 | |
---|
697 | |
---|
698 | // |
---|
699 | // IntersectXY |
---|
700 | // |
---|
701 | // Decide if and where the x/y trajectory hits the elliptical cross |
---|
702 | // section. |
---|
703 | // |
---|
704 | // Arguments: |
---|
705 | // p - (in) Point on trajectory |
---|
706 | // v - (in) Vector along trajectory |
---|
707 | // s - (out) Up to two points of intersection, where the |
---|
708 | // intersection point is p + s*v, and if there are |
---|
709 | // two intersections, s[0] < s[1]. May be negative. |
---|
710 | // Returns: |
---|
711 | // The number of intersections. If 0, the trajectory misses. If 1, the |
---|
712 | // trajectory just grazes the surface. |
---|
713 | // |
---|
714 | // Solution: |
---|
715 | // One needs to solve: ( (p.x + s*v.x)/dx )**2 + ( (p.y + s*v.y)/dy )**2 = 1 |
---|
716 | // |
---|
717 | // The solution is quadratic: a*s**2 + b*s + c = 0 |
---|
718 | // |
---|
719 | // a = (v.x/dx)**2 + (v.y/dy)**2 |
---|
720 | // b = 2*p.x*v.x/dx**2 + 2*p.y*v.y/dy**2 |
---|
721 | // c = (p.x/dx)**2 + (p.y/dy)**2 - 1 |
---|
722 | // |
---|
723 | G4int G4EllipticalTube::IntersectXY( const G4ThreeVector &p, |
---|
724 | const G4ThreeVector &v, |
---|
725 | G4double s[2] ) const |
---|
726 | { |
---|
727 | G4double px = p.x(), py = p.y(); |
---|
728 | G4double vx = v.x(), vy = v.y(); |
---|
729 | |
---|
730 | G4double a = (vx/dx)*(vx/dx) + (vy/dy)*(vy/dy); |
---|
731 | G4double b = 2.0*( px*vx/dx/dx + py*vy/dy/dy ); |
---|
732 | G4double c = (px/dx)*(px/dx) + (py/dy)*(py/dy) - 1.0; |
---|
733 | |
---|
734 | if (a < DBL_MIN) return 0; // Trajectory parallel to z axis |
---|
735 | |
---|
736 | G4double radical = b*b - 4*a*c; |
---|
737 | |
---|
738 | if (radical < -DBL_MIN) return 0; // No solution |
---|
739 | |
---|
740 | if (radical < DBL_MIN) |
---|
741 | { |
---|
742 | // |
---|
743 | // Grazes surface |
---|
744 | // |
---|
745 | s[0] = -b/a/2.0; |
---|
746 | return 1; |
---|
747 | } |
---|
748 | |
---|
749 | radical = std::sqrt(radical); |
---|
750 | |
---|
751 | G4double q = -0.5*( b + (b < 0 ? -radical : +radical) ); |
---|
752 | G4double sa = q/a; |
---|
753 | G4double sb = c/q; |
---|
754 | if (sa < sb) { s[0] = sa; s[1] = sb; } else { s[0] = sb; s[1] = sa; } |
---|
755 | return 2; |
---|
756 | } |
---|
757 | |
---|
758 | |
---|
759 | // |
---|
760 | // GetEntityType |
---|
761 | // |
---|
762 | G4GeometryType G4EllipticalTube::GetEntityType() const |
---|
763 | { |
---|
764 | return G4String("G4EllipticalTube"); |
---|
765 | } |
---|
766 | |
---|
767 | |
---|
768 | // |
---|
769 | // GetCubicVolume |
---|
770 | // |
---|
771 | G4double G4EllipticalTube::GetCubicVolume() |
---|
772 | { |
---|
773 | if(fCubicVolume != 0.) {;} |
---|
774 | else { fCubicVolume = G4VSolid::GetCubicVolume(); } |
---|
775 | return fCubicVolume; |
---|
776 | } |
---|
777 | |
---|
778 | // |
---|
779 | // GetSurfaceArea |
---|
780 | // |
---|
781 | G4double G4EllipticalTube::GetSurfaceArea() |
---|
782 | { |
---|
783 | if(fSurfaceArea != 0.) {;} |
---|
784 | else { fSurfaceArea = G4VSolid::GetSurfaceArea(); } |
---|
785 | return fSurfaceArea; |
---|
786 | } |
---|
787 | |
---|
788 | // |
---|
789 | // Stream object contents to an output stream |
---|
790 | // |
---|
791 | std::ostream& G4EllipticalTube::StreamInfo(std::ostream& os) const |
---|
792 | { |
---|
793 | os << "-----------------------------------------------------------\n" |
---|
794 | << " *** Dump for solid - " << GetName() << " ***\n" |
---|
795 | << " ===================================================\n" |
---|
796 | << " Solid type: G4EllipticalTube\n" |
---|
797 | << " Parameters: \n" |
---|
798 | << " length Z: " << dz/mm << " mm \n" |
---|
799 | << " surface equation in X and Y: \n" |
---|
800 | << " (X / " << dx << ")^2 + (Y / " << dy << ")^2 = 1 \n" |
---|
801 | << "-----------------------------------------------------------\n"; |
---|
802 | |
---|
803 | return os; |
---|
804 | } |
---|
805 | |
---|
806 | |
---|
807 | // |
---|
808 | // GetPointOnSurface |
---|
809 | // |
---|
810 | // Randomly generates a point on the surface, |
---|
811 | // with ~ uniform distribution across surface. |
---|
812 | // |
---|
813 | G4ThreeVector G4EllipticalTube::GetPointOnSurface() const |
---|
814 | { |
---|
815 | G4double xRand, yRand, zRand, phi, cosphi, sinphi, zArea, cArea,p, chose; |
---|
816 | |
---|
817 | phi = RandFlat::shoot(0., 2.*pi); |
---|
818 | cosphi = std::cos(phi); |
---|
819 | sinphi = std::sin(phi); |
---|
820 | |
---|
821 | // the ellipse perimeter from: "http://mathworld.wolfram.com/Ellipse.html" |
---|
822 | // m = (dx - dy)/(dx + dy); |
---|
823 | // k = 1.+1./4.*m*m+1./64.*sqr(m)*sqr(m)+1./256.*sqr(m)*sqr(m)*sqr(m); |
---|
824 | // p = pi*(a+b)*k; |
---|
825 | |
---|
826 | // perimeter below from "http://www.efunda.com/math/areas/EllipseGen.cfm" |
---|
827 | |
---|
828 | p = 2.*pi*std::sqrt(0.5*(dx*dx+dy*dy)); |
---|
829 | |
---|
830 | cArea = 2.*dz*p; |
---|
831 | zArea = pi*dx*dy; |
---|
832 | |
---|
833 | xRand = dx*cosphi; |
---|
834 | yRand = dy*sinphi; |
---|
835 | zRand = RandFlat::shoot(dz, -1.*dz); |
---|
836 | |
---|
837 | chose = RandFlat::shoot(0.,2.*zArea+cArea); |
---|
838 | |
---|
839 | if( (chose>=0) && (chose < cArea) ) |
---|
840 | { |
---|
841 | return G4ThreeVector (xRand,yRand,zRand); |
---|
842 | } |
---|
843 | else if( (chose >= cArea) && (chose < cArea + zArea) ) |
---|
844 | { |
---|
845 | xRand = RandFlat::shoot(-1.*dx,dx); |
---|
846 | yRand = std::sqrt(1.-sqr(xRand/dx)); |
---|
847 | yRand = RandFlat::shoot(-1.*yRand, yRand); |
---|
848 | return G4ThreeVector (xRand,yRand,dz); |
---|
849 | } |
---|
850 | else |
---|
851 | { |
---|
852 | xRand = RandFlat::shoot(-1.*dx,dx); |
---|
853 | yRand = std::sqrt(1.-sqr(xRand/dx)); |
---|
854 | yRand = RandFlat::shoot(-1.*yRand, yRand); |
---|
855 | return G4ThreeVector (xRand,yRand,-1.*dz); |
---|
856 | } |
---|
857 | } |
---|
858 | |
---|
859 | |
---|
860 | // |
---|
861 | // CreatePolyhedron |
---|
862 | // |
---|
863 | G4Polyhedron* G4EllipticalTube::CreatePolyhedron() const |
---|
864 | { |
---|
865 | // create cylinder with radius=1... |
---|
866 | // |
---|
867 | G4Polyhedron* eTube = new G4PolyhedronTube(0.,1.,dz); |
---|
868 | |
---|
869 | // apply non-uniform scaling... |
---|
870 | // |
---|
871 | eTube->Transform(G4Scale3D(dx,dy,1.)); |
---|
872 | return eTube; |
---|
873 | } |
---|
874 | |
---|
875 | |
---|
876 | // |
---|
877 | // GetPolyhedron |
---|
878 | // |
---|
879 | G4Polyhedron* G4EllipticalTube::GetPolyhedron () const |
---|
880 | { |
---|
881 | if (!fpPolyhedron || |
---|
882 | fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() != |
---|
883 | fpPolyhedron->GetNumberOfRotationSteps()) |
---|
884 | { |
---|
885 | delete fpPolyhedron; |
---|
886 | fpPolyhedron = CreatePolyhedron(); |
---|
887 | } |
---|
888 | return fpPolyhedron; |
---|
889 | } |
---|
890 | |
---|
891 | |
---|
892 | // |
---|
893 | // DescribeYourselfTo |
---|
894 | // |
---|
895 | void G4EllipticalTube::DescribeYourselfTo( G4VGraphicsScene& scene ) const |
---|
896 | { |
---|
897 | scene.AddSolid (*this); |
---|
898 | } |
---|
899 | |
---|
900 | |
---|
901 | // |
---|
902 | // GetExtent |
---|
903 | // |
---|
904 | G4VisExtent G4EllipticalTube::GetExtent() const |
---|
905 | { |
---|
906 | return G4VisExtent( -dx, dx, -dy, dy, -dz, dz ); |
---|
907 | } |
---|