Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Ellipsoid.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 // $Id: G4Ellipsoid.cc 67011 2013-01-29 16:17:41Z gcosmo $
27 //
28 // class G4Ellipsoid
29 //
30 // Implementation for G4Ellipsoid class
31 //
32 // History:
33 //
34 // 10.11.99 G.Horton-Smith -- first writing, based on G4Sphere class
35 // 25.02.05 G.Guerrieri -- Modified for future Geant4 release
36 //
37 // --------------------------------------------------------------------
38 
39 #include "globals.hh"
40 
41 #include "G4Ellipsoid.hh"
42 
43 #include "G4VoxelLimits.hh"
44 #include "G4AffineTransform.hh"
45 #include "G4GeometryTolerance.hh"
46 
47 #include "meshdefs.hh"
48 
49 #include "Randomize.hh"
50 
51 #include "G4VGraphicsScene.hh"
52 #include "G4Polyhedron.hh"
53 #include "G4NURBS.hh"
54 #include "G4NURBSbox.hh"
55 #include "G4VisExtent.hh"
56 
57 using namespace CLHEP;
58 
60 //
61 // constructor - check parameters, convert angles so 0<sphi+dpshi<=2_PI
62 // - note if pDPhi>2PI then reset to 2PI
63 
65  G4double pxSemiAxis,
66  G4double pySemiAxis,
67  G4double pzSemiAxis,
68  G4double pzBottomCut,
69  G4double pzTopCut)
70  : G4VSolid(pName), fpPolyhedron(0), fCubicVolume(0.), fSurfaceArea(0.),
71  zBottomCut(0.), zTopCut(0.)
72 {
73  // note: for users that want to use the full ellipsoid it is useful
74  // to include a default for the cuts
75 
77 
78  // Check Semi-Axis
79  if ( (pxSemiAxis<=0.) || (pySemiAxis<=0.) || (pzSemiAxis<=0.) )
80  {
81  std::ostringstream message;
82  message << "Invalid semi-axis - " << GetName();
83  G4Exception("G4Ellipsoid::G4Ellipsoid()", "GeomSolids0002",
84  FatalErrorInArgument, message);
85  }
86  SetSemiAxis(pxSemiAxis, pySemiAxis, pzSemiAxis);
87 
88  if ( pzBottomCut == 0 && pzTopCut == 0 )
89  {
90  SetZCuts(-pzSemiAxis, pzSemiAxis);
91  }
92  else if ( (pzBottomCut < pzSemiAxis) && (pzTopCut > -pzSemiAxis)
93  && (pzBottomCut < pzTopCut) )
94  {
95  SetZCuts(pzBottomCut, pzTopCut);
96  }
97  else
98  {
99  std::ostringstream message;
100  message << "Invalid z-coordinate for cutting plane - " << GetName();
101  G4Exception("G4Ellipsoid::G4Ellipsoid()", "GeomSolids0002",
102  FatalErrorInArgument, message);
103  }
104 }
105 
107 //
108 // Fake default constructor - sets only member data and allocates memory
109 // for usage restricted to object persistency.
110 //
112  : G4VSolid(a), fpPolyhedron(0), kRadTolerance(0.), fCubicVolume(0.),
113  fSurfaceArea(0.), xSemiAxis(0.), ySemiAxis(0.), zSemiAxis(0.),
114  semiAxisMax(0.), zBottomCut(0.), zTopCut(0.)
115 {
116 }
117 
119 //
120 // Destructor
121 
123 {
124 }
125 
127 //
128 // Copy constructor
129 
131  : G4VSolid(rhs),
132  fpPolyhedron(0), kRadTolerance(rhs.kRadTolerance),
133  fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
134  xSemiAxis(rhs.xSemiAxis), ySemiAxis(rhs.ySemiAxis),
135  zSemiAxis(rhs.zSemiAxis), semiAxisMax(rhs.semiAxisMax),
136  zBottomCut(rhs.zBottomCut), zTopCut(rhs.zTopCut)
137 {
138 }
139 
141 //
142 // Assignment operator
143 
145 {
146  // Check assignment to self
147  //
148  if (this == &rhs) { return *this; }
149 
150  // Copy base class data
151  //
152  G4VSolid::operator=(rhs);
153 
154  // Copy data
155  //
156  fpPolyhedron = 0; kRadTolerance = rhs.kRadTolerance;
157  fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
158  xSemiAxis = rhs.xSemiAxis; ySemiAxis = rhs.ySemiAxis;
159  zSemiAxis = rhs.zSemiAxis; semiAxisMax = rhs.semiAxisMax;
160  zBottomCut = rhs.zBottomCut; zTopCut = rhs.zTopCut;
161 
162  return *this;
163 }
164 
166 //
167 // Calculate extent under transform and specified limit
168 
169 G4bool
171  const G4VoxelLimits& pVoxelLimit,
172  const G4AffineTransform& pTransform,
173  G4double& pMin, G4double& pMax) const
174 {
175  if (!pTransform.IsRotated())
176  {
177  // Special case handling for unrotated solid ellipsoid
178  // Compute x/y/z mins and maxs for bounding box respecting limits,
179  // with early returns if outside limits. Then switch() on pAxis,
180  // and compute exact x and y limit for x/y case
181 
182  G4double xoffset,xMin,xMax;
183  G4double yoffset,yMin,yMax;
184  G4double zoffset,zMin,zMax;
185 
186  G4double maxDiff,newMin,newMax;
187  G4double xoff,yoff;
188 
189  xoffset=pTransform.NetTranslation().x();
190  xMin=xoffset - xSemiAxis;
191  xMax=xoffset + xSemiAxis;
192  if (pVoxelLimit.IsXLimited())
193  {
194  if ( (xMin>pVoxelLimit.GetMaxXExtent()+kCarTolerance)
195  || (xMax<pVoxelLimit.GetMinXExtent()-kCarTolerance) )
196  {
197  return false;
198  }
199  else
200  {
201  if (xMin<pVoxelLimit.GetMinXExtent())
202  {
203  xMin=pVoxelLimit.GetMinXExtent();
204  }
205  if (xMax>pVoxelLimit.GetMaxXExtent())
206  {
207  xMax=pVoxelLimit.GetMaxXExtent();
208  }
209  }
210  }
211 
212  yoffset=pTransform.NetTranslation().y();
213  yMin=yoffset - ySemiAxis;
214  yMax=yoffset + ySemiAxis;
215  if (pVoxelLimit.IsYLimited())
216  {
217  if ( (yMin>pVoxelLimit.GetMaxYExtent()+kCarTolerance)
218  || (yMax<pVoxelLimit.GetMinYExtent()-kCarTolerance) )
219  {
220  return false;
221  }
222  else
223  {
224  if (yMin<pVoxelLimit.GetMinYExtent())
225  {
226  yMin=pVoxelLimit.GetMinYExtent();
227  }
228  if (yMax>pVoxelLimit.GetMaxYExtent())
229  {
230  yMax=pVoxelLimit.GetMaxYExtent();
231  }
232  }
233  }
234 
235  zoffset=pTransform.NetTranslation().z();
236  zMin=zoffset + (-zSemiAxis > zBottomCut ? -zSemiAxis : zBottomCut);
237  zMax=zoffset + ( zSemiAxis < zTopCut ? zSemiAxis : zTopCut);
238  if (pVoxelLimit.IsZLimited())
239  {
240  if ( (zMin>pVoxelLimit.GetMaxZExtent()+kCarTolerance)
241  || (zMax<pVoxelLimit.GetMinZExtent()-kCarTolerance) )
242  {
243  return false;
244  }
245  else
246  {
247  if (zMin<pVoxelLimit.GetMinZExtent())
248  {
249  zMin=pVoxelLimit.GetMinZExtent();
250  }
251  if (zMax>pVoxelLimit.GetMaxZExtent())
252  {
253  zMax=pVoxelLimit.GetMaxZExtent();
254  }
255  }
256  }
257 
258  // if here, then known to cut bounding box around ellipsoid
259  //
260  xoff = (xoffset < xMin) ? (xMin-xoffset)
261  : (xoffset > xMax) ? (xoffset-xMax) : 0.0;
262  yoff = (yoffset < yMin) ? (yMin-yoffset)
263  : (yoffset > yMax) ? (yoffset-yMax) : 0.0;
264 
265  // detailed calculations
266  // NOTE: does not use X or Y offsets to adjust Z range,
267  // and does not use Z offset to adjust X or Y range,
268  // which is consistent with G4Sphere::CalculateExtent behavior
269  //
270  switch (pAxis)
271  {
272  case kXAxis:
273  if (yoff==0.)
274  {
275  // YZ limits cross max/min x => no change
276  //
277  pMin=xMin;
278  pMax=xMax;
279  }
280  else
281  {
282  // YZ limits don't cross max/min x => compute max delta x,
283  // hence new mins/maxs
284  //
285  maxDiff= 1.0-sqr(yoff/ySemiAxis);
286  if (maxDiff < 0.0) { return false; }
287  maxDiff= xSemiAxis * std::sqrt(maxDiff);
288  newMin=xoffset-maxDiff;
289  newMax=xoffset+maxDiff;
290  pMin=(newMin<xMin) ? xMin : newMin;
291  pMax=(newMax>xMax) ? xMax : newMax;
292  }
293  break;
294  case kYAxis:
295  if (xoff==0.)
296  {
297  // XZ limits cross max/min y => no change
298  //
299  pMin=yMin;
300  pMax=yMax;
301  }
302  else
303  {
304  // XZ limits don't cross max/min y => compute max delta y,
305  // hence new mins/maxs
306  //
307  maxDiff= 1.0-sqr(xoff/xSemiAxis);
308  if (maxDiff < 0.0) { return false; }
309  maxDiff= ySemiAxis * std::sqrt(maxDiff);
310  newMin=yoffset-maxDiff;
311  newMax=yoffset+maxDiff;
312  pMin=(newMin<yMin) ? yMin : newMin;
313  pMax=(newMax>yMax) ? yMax : newMax;
314  }
315  break;
316  case kZAxis:
317  pMin=zMin;
318  pMax=zMax;
319  break;
320  default:
321  break;
322  }
323 
324  pMin-=kCarTolerance;
325  pMax+=kCarTolerance;
326  return true;
327  }
328  else // not rotated
329  {
330  G4int i,j,noEntries,noBetweenSections;
331  G4bool existsAfterClip=false;
332 
333  // Calculate rotated vertex coordinates
334 
335  G4int noPolygonVertices=0;
336  G4ThreeVectorList* vertices =
337  CreateRotatedVertices(pTransform,noPolygonVertices);
338 
339  pMin=+kInfinity;
340  pMax=-kInfinity;
341 
342  noEntries=vertices->size(); // noPolygonVertices*noPhiCrossSections
343  noBetweenSections=noEntries-noPolygonVertices;
344 
345  G4ThreeVectorList ThetaPolygon;
346  for (i=0;i<noEntries;i+=noPolygonVertices)
347  {
348  for(j=0;j<(noPolygonVertices/2)-1;j++)
349  {
350  ThetaPolygon.push_back((*vertices)[i+j]);
351  ThetaPolygon.push_back((*vertices)[i+j+1]);
352  ThetaPolygon.push_back((*vertices)[i+noPolygonVertices-2-j]);
353  ThetaPolygon.push_back((*vertices)[i+noPolygonVertices-1-j]);
354  CalculateClippedPolygonExtent(ThetaPolygon,pVoxelLimit,pAxis,pMin,pMax);
355  ThetaPolygon.clear();
356  }
357  }
358  for (i=0;i<noBetweenSections;i+=noPolygonVertices)
359  {
360  for(j=0;j<noPolygonVertices-1;j++)
361  {
362  ThetaPolygon.push_back((*vertices)[i+j]);
363  ThetaPolygon.push_back((*vertices)[i+j+1]);
364  ThetaPolygon.push_back((*vertices)[i+noPolygonVertices+j+1]);
365  ThetaPolygon.push_back((*vertices)[i+noPolygonVertices+j]);
366  CalculateClippedPolygonExtent(ThetaPolygon,pVoxelLimit,pAxis,pMin,pMax);
367  ThetaPolygon.clear();
368  }
369  ThetaPolygon.push_back((*vertices)[i+noPolygonVertices-1]);
370  ThetaPolygon.push_back((*vertices)[i]);
371  ThetaPolygon.push_back((*vertices)[i+noPolygonVertices]);
372  ThetaPolygon.push_back((*vertices)[i+2*noPolygonVertices-1]);
373  CalculateClippedPolygonExtent(ThetaPolygon,pVoxelLimit,pAxis,pMin,pMax);
374  ThetaPolygon.clear();
375  }
376  if ( (pMin!=kInfinity) || (pMax!=-kInfinity) )
377  {
378  existsAfterClip=true;
379 
380  // Add 2*tolerance to avoid precision troubles
381  //
382  pMin-=kCarTolerance;
383  pMax+=kCarTolerance;
384 
385  }
386  else
387  {
388  // Check for case where completely enveloping clipping volume
389  // If point inside then we are confident that the solid completely
390  // envelopes the clipping volume. Hence set min/max extents according
391  // to clipping volume extents along the specified axis.
392  //
394  clipCentre((pVoxelLimit.GetMinXExtent()+pVoxelLimit.GetMaxXExtent())*0.5,
395  (pVoxelLimit.GetMinYExtent()+pVoxelLimit.GetMaxYExtent())*0.5,
396  (pVoxelLimit.GetMinZExtent()+pVoxelLimit.GetMaxZExtent())*0.5);
397 
398  if (Inside(pTransform.Inverse().TransformPoint(clipCentre))!=kOutside)
399  {
400  existsAfterClip=true;
401  pMin=pVoxelLimit.GetMinExtent(pAxis);
402  pMax=pVoxelLimit.GetMaxExtent(pAxis);
403  }
404  }
405  delete vertices;
406  return existsAfterClip;
407  }
408 }
409 
411 //
412 // Return whether point inside/outside/on surface
413 // Split into radius, phi, theta checks
414 // Each check modifies `in', or returns as approprate
415 
417 {
418  G4double rad2oo, // outside surface outer tolerance
419  rad2oi; // outside surface inner tolerance
420  EInside in;
421 
422  static const G4double halfRadTolerance=kRadTolerance*0.5;
423 
424  // check this side of z cut first, because that's fast
425  //
426  if (p.z() < zBottomCut-halfRadTolerance) { return in=kOutside; }
427  if (p.z() > zTopCut+halfRadTolerance) { return in=kOutside; }
428 
429  rad2oo= sqr(p.x()/(xSemiAxis+halfRadTolerance))
430  + sqr(p.y()/(ySemiAxis+halfRadTolerance))
431  + sqr(p.z()/(zSemiAxis+halfRadTolerance));
432 
433  if (rad2oo > 1.0) { return in=kOutside; }
434 
435  rad2oi= sqr(p.x()*(1.0+halfRadTolerance/xSemiAxis)/xSemiAxis)
436  + sqr(p.y()*(1.0+halfRadTolerance/ySemiAxis)/ySemiAxis)
437  + sqr(p.z()*(1.0+halfRadTolerance/zSemiAxis)/zSemiAxis);
438 
439  // Check radial surfaces
440  // sets `in' (already checked for rad2oo > 1.0)
441  //
442  if (rad2oi < 1.0)
443  {
444  in = ( (p.z() < zBottomCut+halfRadTolerance)
445  || (p.z() > zTopCut-halfRadTolerance) ) ? kSurface : kInside;
446  if ( rad2oi > 1.0-halfRadTolerance ) { in=kSurface; }
447  }
448  else
449  {
450  in = kSurface;
451  }
452  return in;
453 
454 }
455 
457 //
458 // Return unit normal of surface closest to p not protected against p=0
459 
461 {
462  G4double distR, distZBottom, distZTop;
463 
464  // normal vector with special magnitude: parallel to normal, units 1/length
465  // norm*p == 1.0 if on surface, >1.0 if outside, <1.0 if inside
466  //
467  G4ThreeVector norm(p.x()/(xSemiAxis*xSemiAxis),
468  p.y()/(ySemiAxis*ySemiAxis),
469  p.z()/(zSemiAxis*zSemiAxis));
470  G4double radius = 1.0/norm.mag();
471 
472  // approximate distance to curved surface
473  //
474  distR = std::fabs( (p*norm - 1.0) * radius ) / 2.0;
475 
476  // Distance to z-cut plane
477  //
478  distZBottom = std::fabs( p.z() - zBottomCut );
479  distZTop = std::fabs( p.z() - zTopCut );
480 
481  if ( (distZBottom < distR) || (distZTop < distR) )
482  {
483  return G4ThreeVector(0.,0.,(distZBottom < distZTop) ? -1.0 : 1.0);
484  }
485  return ( norm *= radius );
486 }
487 
489 //
490 // Calculate distance to shape from outside, along normalised vector
491 // - return kInfinity if no intersection, or intersection distance <= tolerance
492 //
493 
495  const G4ThreeVector& v ) const
496 {
497  static const G4double halfCarTolerance=kCarTolerance*0.5;
498  static const G4double halfRadTolerance=kRadTolerance*0.5;
499 
500  G4double distMin = std::min(xSemiAxis,ySemiAxis);
501  const G4double dRmax = 100.*std::min(distMin,zSemiAxis);
502  distMin= kInfinity;
503 
504  // check to see if Z plane is relevant
505  if (p.z() <= zBottomCut+halfCarTolerance)
506  {
507  if (v.z() <= 0.0) { return distMin; }
508  G4double distZ = (zBottomCut - p.z()) / v.z();
509 
510  if ( (distZ > -halfRadTolerance) && (Inside(p+distZ*v) != kOutside) )
511  {
512  // early exit since can't intercept curved surface if we reach here
513  if ( std::fabs(distZ) < halfRadTolerance ) { distZ=0.; }
514  return distMin= distZ;
515  }
516  }
517  if (p.z() >= zTopCut-halfCarTolerance)
518  {
519  if (v.z() >= 0.0) { return distMin;}
520  G4double distZ = (zTopCut - p.z()) / v.z();
521  if ( (distZ > -halfRadTolerance) && (Inside(p+distZ*v) != kOutside) )
522  {
523  // early exit since can't intercept curved surface if we reach here
524  if ( std::fabs(distZ) < halfRadTolerance ) { distZ=0.; }
525  return distMin= distZ;
526  }
527  }
528  // if fZCut1 <= p.z() <= fZCut2, then must hit curved surface
529 
530  // now check curved surface intercept
531  G4double A,B,C;
532 
533  A= sqr(v.x()/xSemiAxis) + sqr(v.y()/ySemiAxis) + sqr(v.z()/zSemiAxis);
534  C= sqr(p.x()/xSemiAxis) + sqr(p.y()/ySemiAxis) + sqr(p.z()/zSemiAxis) - 1.0;
535  B= 2.0 * ( p.x()*v.x()/(xSemiAxis*xSemiAxis)
536  + p.y()*v.y()/(ySemiAxis*ySemiAxis)
537  + p.z()*v.z()/(zSemiAxis*zSemiAxis) );
538 
539  C= B*B - 4.0*A*C;
540  if (C > 0.0)
541  {
542  G4double distR= (-B - std::sqrt(C)) / (2.0*A);
543  G4double intZ = p.z()+distR*v.z();
544  if ( (distR > halfRadTolerance)
545  && (intZ >= zBottomCut-halfRadTolerance)
546  && (intZ <= zTopCut+halfRadTolerance) )
547  {
548  distMin = distR;
549  }
550  else if( (distR >- halfRadTolerance)
551  && (intZ >= zBottomCut-halfRadTolerance)
552  && (intZ <= zTopCut+halfRadTolerance) )
553  {
554  // p is on the curved surface, DistanceToIn returns 0 or kInfinity:
555  // DistanceToIn returns 0, if second root is positive (means going inside)
556  // If second root is negative, DistanceToIn returns kInfinity (outside)
557  //
558  distR = (-B + std::sqrt(C) ) / (2.0*A);
559  if(distR>0.) { distMin=0.; }
560  }
561  else
562  {
563  distR= (-B + std::sqrt(C)) / (2.0*A);
564  intZ = p.z()+distR*v.z();
565  if ( (distR > halfRadTolerance)
566  && (intZ >= zBottomCut-halfRadTolerance)
567  && (intZ <= zTopCut+halfRadTolerance) )
568  {
570  if (norm.dot(v)<0.) { distMin = distR; }
571  }
572  }
573  if ( (distMin!=kInfinity) && (distMin>dRmax) )
574  { // Avoid rounding errors due to precision issues on
575  // 64 bits systems. Split long distances and recompute
576  G4double fTerm = distMin-std::fmod(distMin,dRmax);
577  distMin = fTerm + DistanceToIn(p+fTerm*v,v);
578  }
579  }
580 
581  if (std::fabs(distMin)<halfRadTolerance) { distMin=0.; }
582  return distMin;
583 }
584 
586 //
587 // Calculate distance (<= actual) to closest surface of shape from outside
588 // - Return 0 if point inside
589 
591 {
592  G4double distR, distZ;
593 
594  // normal vector: parallel to normal, magnitude 1/(characteristic radius)
595  //
596  G4ThreeVector norm(p.x()/(xSemiAxis*xSemiAxis),
597  p.y()/(ySemiAxis*ySemiAxis),
598  p.z()/(zSemiAxis*zSemiAxis));
599  G4double radius= 1.0/norm.mag();
600 
601  // approximate distance to curved surface ( <= actual distance )
602  //
603  distR= (p*norm - 1.0) * radius / 2.0;
604 
605  // Distance to z-cut plane
606  //
607  distZ= zBottomCut - p.z();
608  if (distZ < 0.0)
609  {
610  distZ = p.z() - zTopCut;
611  }
612 
613  // Distance to closest surface from outside
614  //
615  if (distZ < 0.0)
616  {
617  return (distR < 0.0) ? 0.0 : distR;
618  }
619  else if (distR < 0.0)
620  {
621  return distZ;
622  }
623  else
624  {
625  return (distZ < distR) ? distZ : distR;
626  }
627 }
628 
630 //
631 // Calculate distance to surface of shape from `inside', allowing for tolerance
632 
634  const G4ThreeVector& v,
635  const G4bool calcNorm,
636  G4bool *validNorm,
637  G4ThreeVector *n ) const
638 {
639  G4double distMin;
640  enum surface_e {kPlaneSurf, kCurvedSurf, kNoSurf} surface;
641 
642  distMin= kInfinity;
643  surface= kNoSurf;
644 
645  // check to see if Z plane is relevant
646  //
647  if (v.z() < 0.0)
648  {
649  G4double distZ = (zBottomCut - p.z()) / v.z();
650  if (distZ < 0.0)
651  {
652  distZ= 0.0;
653  if (!calcNorm) {return 0.0;}
654  }
655  distMin= distZ;
656  surface= kPlaneSurf;
657  }
658  if (v.z() > 0.0)
659  {
660  G4double distZ = (zTopCut - p.z()) / v.z();
661  if (distZ < 0.0)
662  {
663  distZ= 0.0;
664  if (!calcNorm) {return 0.0;}
665  }
666  distMin= distZ;
667  surface= kPlaneSurf;
668  }
669 
670  // normal vector: parallel to normal, magnitude 1/(characteristic radius)
671  //
672  G4ThreeVector nearnorm(p.x()/(xSemiAxis*xSemiAxis),
673  p.y()/(ySemiAxis*ySemiAxis),
674  p.z()/(zSemiAxis*zSemiAxis));
675 
676  // now check curved surface intercept
677  //
678  G4double A,B,C;
679 
680  A= sqr(v.x()/xSemiAxis) + sqr(v.y()/ySemiAxis) + sqr(v.z()/zSemiAxis);
681  C= (p * nearnorm) - 1.0;
682  B= 2.0 * (v * nearnorm);
683 
684  C= B*B - 4.0*A*C;
685  if (C > 0.0)
686  {
687  G4double distR= (-B + std::sqrt(C) ) / (2.0*A);
688  if (distR < 0.0)
689  {
690  distR= 0.0;
691  if (!calcNorm) {return 0.0;}
692  }
693  if (distR < distMin)
694  {
695  distMin= distR;
696  surface= kCurvedSurf;
697  }
698  }
699 
700  // set normal if requested
701  //
702  if (calcNorm)
703  {
704  if (surface == kNoSurf)
705  {
706  *validNorm = false;
707  }
708  else
709  {
710  *validNorm = true;
711  switch (surface)
712  {
713  case kPlaneSurf:
714  *n= G4ThreeVector(0.,0.,(v.z() > 0.0 ? 1. : -1.));
715  break;
716  case kCurvedSurf:
717  {
718  G4ThreeVector pexit= p + distMin*v;
719  G4ThreeVector truenorm(pexit.x()/(xSemiAxis*xSemiAxis),
720  pexit.y()/(ySemiAxis*ySemiAxis),
721  pexit.z()/(zSemiAxis*zSemiAxis));
722  truenorm *= 1.0/truenorm.mag();
723  *n= truenorm;
724  } break;
725  default: // Should never reach this case ...
726  DumpInfo();
727  std::ostringstream message;
728  G4int oldprc = message.precision(16);
729  message << "Undefined side for valid surface normal to solid."
730  << G4endl
731  << "Position:" << G4endl
732  << " p.x() = " << p.x()/mm << " mm" << G4endl
733  << " p.y() = " << p.y()/mm << " mm" << G4endl
734  << " p.z() = " << p.z()/mm << " mm" << G4endl
735  << "Direction:" << G4endl << G4endl
736  << " v.x() = " << v.x() << G4endl
737  << " v.y() = " << v.y() << G4endl
738  << " v.z() = " << v.z() << G4endl
739  << "Proposed distance :" << G4endl
740  << " distMin = " << distMin/mm << " mm";
741  message.precision(oldprc);
742  G4Exception("G4Ellipsoid::DistanceToOut(p,v,..)",
743  "GeomSolids1002", JustWarning, message);
744  break;
745  }
746  }
747  }
748 
749  return distMin;
750 }
751 
753 //
754 // Calculate distance (<=actual) to closest surface of shape from inside
755 
757 {
758  G4double distR, distZ;
759 
760 #ifdef G4SPECSDEBUG
761  if( Inside(p) == kOutside )
762  {
763  DumpInfo();
764  std::ostringstream message;
765  G4int oldprc = message.precision(16);
766  message << "Point p is outside !?" << G4endl
767  << "Position:" << G4endl
768  << " p.x() = " << p.x()/mm << " mm" << G4endl
769  << " p.y() = " << p.y()/mm << " mm" << G4endl
770  << " p.z() = " << p.z()/mm << " mm";
771  message.precision(oldprc) ;
772  G4Exception("G4Ellipsoid::DistanceToOut(p)", "GeomSolids1002",
773  JustWarning, message);
774  }
775 #endif
776 
777  // Normal vector: parallel to normal, magnitude 1/(characteristic radius)
778  //
779  G4ThreeVector norm(p.x()/(xSemiAxis*xSemiAxis),
780  p.y()/(ySemiAxis*ySemiAxis),
781  p.z()/(zSemiAxis*zSemiAxis));
782 
783  // the following is a safe inlined "radius= min(1.0/norm.mag(),p.mag())
784  //
785  G4double radius= p.mag();
786  G4double tmp= norm.mag();
787  if ( (tmp > 0.0) && (1.0 < radius*tmp) ) {radius = 1.0/tmp;}
788 
789  // Approximate distance to curved surface ( <= actual distance )
790  //
791  distR = (1.0 - p*norm) * radius / 2.0;
792 
793  // Distance to z-cut plane
794  //
795  distZ = p.z() - zBottomCut;
796  if (distZ < 0.0) {distZ= zTopCut - p.z();}
797 
798  // Distance to closest surface from inside
799  //
800  if ( (distZ < 0.0) || (distR < 0.0) )
801  {
802  return 0.0;
803  }
804  else
805  {
806  return (distZ < distR) ? distZ : distR;
807  }
808 }
809 
811 //
812 // Create a List containing the transformed vertices
813 // Ordering [0-3] -fDz cross section
814 // [4-7] +fDz cross section such that [0] is below [4],
815 // [1] below [5] etc.
816 // Note:
817 // Caller has deletion resposibility
818 // Potential improvement: For last slice, use actual ending angle
819 // to avoid rounding error problems.
820 
823  G4int& noPolygonVertices) const
824 {
825  G4ThreeVectorList *vertices;
826  G4ThreeVector vertex;
827  G4double meshAnglePhi, meshRMaxFactor,
828  crossAnglePhi, coscrossAnglePhi, sincrossAnglePhi, sAnglePhi;
829  G4double meshTheta, crossTheta, startTheta;
830  G4double rMaxX, rMaxY, rMaxZ, rMaxMax, rx, ry, rz;
831  G4int crossSectionPhi, noPhiCrossSections, crossSectionTheta, noThetaSections;
832 
833  // Phi cross sections
834  //
835  noPhiCrossSections=G4int (twopi/kMeshAngleDefault)+1; // = 9!
836 
837 /*
838  if (noPhiCrossSections<kMinMeshSections) // <3
839  {
840  noPhiCrossSections=kMinMeshSections;
841  }
842  else if (noPhiCrossSections>kMaxMeshSections) // >37
843  {
844  noPhiCrossSections=kMaxMeshSections;
845  }
846 */
847  meshAnglePhi=twopi/(noPhiCrossSections-1);
848 
849  // Set start angle such that mesh will be at fRMax
850  // on the x axis. Will give better extent calculations when not rotated.
851 
852  sAnglePhi = -meshAnglePhi*0.5;
853 
854  // Theta cross sections
855 
856  noThetaSections = G4int(pi/kMeshAngleDefault)+3; // = 7!
857 
858 /*
859  if (noThetaSections<kMinMeshSections) // <3
860  {
861  noThetaSections=kMinMeshSections;
862  }
863  else if (noThetaSections>kMaxMeshSections) // >37
864  {
865  noThetaSections=kMaxMeshSections;
866  }
867 */
868  meshTheta= pi/(noThetaSections-2);
869 
870  // Set start angle such that mesh will be at fRMax
871  // on the z axis. Will give better extent calculations when not rotated.
872 
873  startTheta = -meshTheta*0.5;
874 
875  meshRMaxFactor = 1.0/std::cos(0.5*
876  std::sqrt(meshAnglePhi*meshAnglePhi+meshTheta*meshTheta));
877  rMaxMax= (xSemiAxis > ySemiAxis ? xSemiAxis : ySemiAxis);
878  if (zSemiAxis > rMaxMax) rMaxMax= zSemiAxis;
879  rMaxX= xSemiAxis + rMaxMax*(meshRMaxFactor-1.0);
880  rMaxY= ySemiAxis + rMaxMax*(meshRMaxFactor-1.0);
881  rMaxZ= zSemiAxis + rMaxMax*(meshRMaxFactor-1.0);
882  G4double* cosCrossTheta = new G4double[noThetaSections];
883  G4double* sinCrossTheta = new G4double[noThetaSections];
884  vertices=new G4ThreeVectorList(noPhiCrossSections*noThetaSections);
885  if (vertices && cosCrossTheta && sinCrossTheta)
886  {
887  for (crossSectionTheta=0; crossSectionTheta<noThetaSections;
888  crossSectionTheta++)
889  {
890  // Compute sine and cosine table (for historical reasons)
891  //
892  crossTheta=startTheta+crossSectionTheta*meshTheta;
893  cosCrossTheta[crossSectionTheta]=std::cos(crossTheta);
894  sinCrossTheta[crossSectionTheta]=std::sin(crossTheta);
895  }
896  for (crossSectionPhi=0; crossSectionPhi<noPhiCrossSections;
897  crossSectionPhi++)
898  {
899  crossAnglePhi=sAnglePhi+crossSectionPhi*meshAnglePhi;
900  coscrossAnglePhi=std::cos(crossAnglePhi);
901  sincrossAnglePhi=std::sin(crossAnglePhi);
902  for (crossSectionTheta=0; crossSectionTheta<noThetaSections;
903  crossSectionTheta++)
904  {
905  // Compute coordinates of cross section at section crossSectionPhi
906  //
907  rx= sinCrossTheta[crossSectionTheta]*coscrossAnglePhi*rMaxX;
908  ry= sinCrossTheta[crossSectionTheta]*sincrossAnglePhi*rMaxY;
909  rz= cosCrossTheta[crossSectionTheta]*rMaxZ;
910  if (rz < zBottomCut)
911  { rz= zBottomCut; }
912  if (rz > zTopCut)
913  { rz= zTopCut; }
914  vertex= G4ThreeVector(rx,ry,rz);
915  vertices->push_back(pTransform.TransformPoint(vertex));
916  } // Theta forward
917  } // Phi
918  noPolygonVertices = noThetaSections ;
919  }
920  else
921  {
922  DumpInfo();
923  G4Exception("G4Ellipsoid::CreateRotatedVertices()",
924  "GeomSolids0003", FatalException,
925  "Error in allocation of vertices. Out of memory !");
926  }
927 
928  delete[] cosCrossTheta;
929  delete[] sinCrossTheta;
930 
931  return vertices;
932 }
933 
935 //
936 // G4EntityType
937 
939 {
940  return G4String("G4Ellipsoid");
941 }
942 
944 //
945 // Make a clone of the object
946 
948 {
949  return new G4Ellipsoid(*this);
950 }
951 
953 //
954 // Stream object contents to an output stream
955 
956 std::ostream& G4Ellipsoid::StreamInfo( std::ostream& os ) const
957 {
958  G4int oldprc = os.precision(16);
959  os << "-----------------------------------------------------------\n"
960  << " *** Dump for solid - " << GetName() << " ***\n"
961  << " ===================================================\n"
962  << " Solid type: G4Ellipsoid\n"
963  << " Parameters: \n"
964 
965  << " semi-axis x: " << xSemiAxis/mm << " mm \n"
966  << " semi-axis y: " << ySemiAxis/mm << " mm \n"
967  << " semi-axis z: " << zSemiAxis/mm << " mm \n"
968  << " max semi-axis: " << semiAxisMax/mm << " mm \n"
969  << " lower cut plane level z: " << zBottomCut/mm << " mm \n"
970  << " upper cut plane level z: " << zTopCut/mm << " mm \n"
971  << "-----------------------------------------------------------\n";
972  os.precision(oldprc);
973 
974  return os;
975 }
976 
978 //
979 // GetPointOnSurface
980 
982 {
983  G4double aTop, aBottom, aCurved, chose, xRand, yRand, zRand, phi;
984  G4double cosphi, sinphi, costheta, sintheta, alpha, beta, max1, max2, max3;
985 
986  max1 = xSemiAxis > ySemiAxis ? xSemiAxis : ySemiAxis;
987  max1 = max1 > zSemiAxis ? max1 : zSemiAxis;
988  if (max1 == xSemiAxis) { max2 = ySemiAxis; max3 = zSemiAxis; }
989  else if (max1 == ySemiAxis) { max2 = xSemiAxis; max3 = zSemiAxis; }
990  else { max2 = xSemiAxis; max3 = ySemiAxis; }
991 
992  phi = RandFlat::shoot(0.,twopi);
993 
994  cosphi = std::cos(phi); sinphi = std::sin(phi);
995  costheta = RandFlat::shoot(zBottomCut,zTopCut)/zSemiAxis;
996  sintheta = std::sqrt(1.-sqr(costheta));
997 
998  alpha = 1.-sqr(max2/max1); beta = 1.-sqr(max3/max1);
999 
1000  aTop = pi*xSemiAxis*ySemiAxis*(1 - sqr(zTopCut/zSemiAxis));
1001  aBottom = pi*xSemiAxis*ySemiAxis*(1 - sqr(zBottomCut/zSemiAxis));
1002 
1003  // approximation
1004  // from:" http://www.citr.auckland.ac.nz/techreports/2004/CITR-TR-139.pdf"
1005  aCurved = 4.*pi*max1*max2*(1.-1./6.*(alpha+beta)-
1006  1./120.*(3.*sqr(alpha)+2.*alpha*beta+3.*sqr(beta)));
1007 
1008  aCurved *= 0.5*(1.2*zTopCut/zSemiAxis - 1.2*zBottomCut/zSemiAxis);
1009 
1010  if( ( zTopCut >= zSemiAxis && zBottomCut <= -1.*zSemiAxis )
1011  || ( zTopCut == 0 && zBottomCut ==0 ) )
1012  {
1013  aTop = 0; aBottom = 0;
1014  }
1015 
1016  chose = RandFlat::shoot(0.,aTop + aBottom + aCurved);
1017 
1018  if(chose < aCurved)
1019  {
1020  xRand = xSemiAxis*sintheta*cosphi;
1021  yRand = ySemiAxis*sintheta*sinphi;
1022  zRand = zSemiAxis*costheta;
1023  return G4ThreeVector (xRand,yRand,zRand);
1024  }
1025  else if(chose >= aCurved && chose < aCurved + aTop)
1026  {
1027  xRand = RandFlat::shoot(-1.,1.)*xSemiAxis
1028  * std::sqrt(1-sqr(zTopCut/zSemiAxis));
1029  yRand = RandFlat::shoot(-1.,1.)*ySemiAxis
1030  * std::sqrt(1.-sqr(zTopCut/zSemiAxis)-sqr(xRand/xSemiAxis));
1031  zRand = zTopCut;
1032  return G4ThreeVector (xRand,yRand,zRand);
1033  }
1034  else
1035  {
1036  xRand = RandFlat::shoot(-1.,1.)*xSemiAxis
1037  * std::sqrt(1-sqr(zBottomCut/zSemiAxis));
1038  yRand = RandFlat::shoot(-1.,1.)*ySemiAxis
1039  * std::sqrt(1.-sqr(zBottomCut/zSemiAxis)-sqr(xRand/xSemiAxis));
1040  zRand = zBottomCut;
1041  return G4ThreeVector (xRand,yRand,zRand);
1042  }
1043 }
1044 
1046 //
1047 // Methods for visualisation
1048 
1050 {
1051  scene.AddSolid(*this);
1052 }
1053 
1055 {
1056  // Define the sides of the box into which the G4Ellipsoid instance would fit.
1057  //
1058  return G4VisExtent (-semiAxisMax, semiAxisMax,
1059  -semiAxisMax, semiAxisMax,
1060  -semiAxisMax, semiAxisMax);
1061 }
1062 
1064 {
1065  // Box for now!!!
1066  //
1067  return new G4NURBSbox(semiAxisMax, semiAxisMax, semiAxisMax);
1068 }
1069 
1071 {
1072  return new G4PolyhedronEllipsoid(xSemiAxis, ySemiAxis, zSemiAxis,
1073  zBottomCut, zTopCut);
1074 }
1075 
1077 {
1078  if (!fpPolyhedron ||
1081  {
1082  delete fpPolyhedron;
1084  }
1085  return fpPolyhedron;
1086 }