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