Geant4  10.00.p01
G4EllipticalCone.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: G4EllipticalCone.cc 72937 2013-08-14 13:20:38Z gcosmo $
27 //
28 // Implementation of G4EllipticalCone class
29 //
30 // This code implements an Elliptical Cone given explicitly by the
31 // equation:
32 // x^2/a^2 + y^2/b^2 = (z-h)^2
33 // and specified by the parameters (a,b,h) and a cut parallel to the
34 // xy plane above z = 0.
35 //
36 // Author: Dionysios Anninos
37 //
38 // --------------------------------------------------------------------
39 
40 #include "globals.hh"
41 
42 #include "G4EllipticalCone.hh"
43 
44 #include "G4ClippablePolygon.hh"
45 #include "G4SolidExtentList.hh"
46 #include "G4VoxelLimits.hh"
47 #include "G4AffineTransform.hh"
48 #include "G4GeometryTolerance.hh"
49 
50 #include "meshdefs.hh"
51 
52 #include "Randomize.hh"
53 
54 #include "G4VGraphicsScene.hh"
55 #include "G4Polyhedron.hh"
56 #include "G4VisExtent.hh"
57 
58 //#define G4SPECSDEBUG 1
59 
60 using namespace CLHEP;
61 
63 //
64 // Constructor - check parameters
65 //
67  G4double pxSemiAxis,
68  G4double pySemiAxis,
69  G4double pzMax,
70  G4double pzTopCut)
71  : G4VSolid(pName), fpPolyhedron(0), fCubicVolume(0.), fSurfaceArea(0.),
72  zTopCut(0.)
73 {
74 
76 
79 
80  // Check Semi-Axis & Z-cut
81  //
82  if ( (pxSemiAxis <= 0.) || (pySemiAxis <= 0.) || (pzMax <= 0.) )
83  {
84  std::ostringstream message;
85  message << "Invalid semi-axis or height - " << GetName();
86  G4Exception("G4EllipticalCone::G4EllipticalCone()", "GeomSolids0002",
87  FatalErrorInArgument, message);
88  }
89  if ( pzTopCut <= 0 )
90  {
91  std::ostringstream message;
92  message << "Invalid z-coordinate for cutting plane - " << GetName();
93  G4Exception("G4EllipticalCone::G4EllipticalCone()", "InvalidSetup",
94  FatalErrorInArgument, message);
95  }
96 
97  SetSemiAxis( pxSemiAxis, pySemiAxis, pzMax );
98  SetZCut(pzTopCut);
99 }
100 
102 //
103 // Fake default constructor - sets only member data and allocates memory
104 // for usage restricted to object persistency.
105 //
107  : G4VSolid(a), fpPolyhedron(0), kRadTolerance(0.),
108  halfRadTol(0.), halfCarTol(0.), fCubicVolume(0.),
109  fSurfaceArea(0.), xSemiAxis(0.), ySemiAxis(0.), zheight(0.),
110  semiAxisMax(0.), zTopCut(0.)
111 {
112 }
113 
115 //
116 // Destructor
117 //
119 {
120 }
121 
123 //
124 // Copy constructor
125 //
127  : G4VSolid(rhs),
128  fpPolyhedron(0), kRadTolerance(rhs.kRadTolerance),
129  halfRadTol(rhs.halfRadTol), halfCarTol(rhs.halfCarTol),
130  fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
131  xSemiAxis(rhs.xSemiAxis), ySemiAxis(rhs.ySemiAxis), zheight(rhs.zheight),
132  semiAxisMax(rhs.semiAxisMax), zTopCut(rhs.zTopCut)
133 {
134 }
135 
137 //
138 // Assignment operator
139 //
141 {
142  // Check assignment to self
143  //
144  if (this == &rhs) { return *this; }
145 
146  // Copy base class data
147  //
148  G4VSolid::operator=(rhs);
149 
150  // Copy data
151  //
155  xSemiAxis = rhs.xSemiAxis; ySemiAxis = rhs.ySemiAxis;
157 
158  return *this;
159 }
160 
162 //
163 // Calculate extent under transform and specified limit
164 //
165 G4bool
167  const G4VoxelLimits &voxelLimit,
168  const G4AffineTransform &transform,
169  G4double &min, G4double &max ) const
170 {
171  G4SolidExtentList extentList( axis, voxelLimit );
172 
173  //
174  // We are going to divide up our elliptical face into small pieces
175  //
176 
177  //
178  // Choose phi size of our segment(s) based on constants as
179  // defined in meshdefs.hh
180  //
181  G4int numPhi = kMaxMeshSections;
182  G4double sigPhi = twopi/numPhi;
183 
184  //
185  // We have to be careful to keep our segments completely outside
186  // of the elliptical surface. To do so we imagine we have
187  // a simple (unit radius) circular cross section (as in G4Tubs)
188  // and then "stretch" the dimensions as necessary to fit the ellipse.
189  //
190  G4double rFudge = 1.0/std::cos(0.5*sigPhi);
191  G4double dxFudgeBot = xSemiAxis*2.*zheight*rFudge,
192  dyFudgeBot = ySemiAxis*2.*zheight*rFudge;
193  G4double dxFudgeTop = xSemiAxis*(zheight-zTopCut)*rFudge,
194  dyFudgeTop = ySemiAxis*(zheight-zTopCut)*rFudge;
195 
196  //
197  // As we work around the elliptical surface, we build
198  // a "phi" segment on the way, and keep track of two
199  // additional polygons for the two ends.
200  //
201  G4ClippablePolygon endPoly1, endPoly2, phiPoly;
202 
203  G4double phi = 0,
204  cosPhi = std::cos(phi),
205  sinPhi = std::sin(phi);
206  G4ThreeVector v0( dxFudgeTop*cosPhi, dyFudgeTop*sinPhi, +zTopCut ),
207  v1( dxFudgeBot*cosPhi, dyFudgeBot*sinPhi, -zTopCut ),
208  w0, w1;
209  transform.ApplyPointTransform( v0 );
210  transform.ApplyPointTransform( v1 );
211  do
212  {
213  phi += sigPhi;
214  if (numPhi == 1) phi = 0; // Try to avoid roundoff
215  cosPhi = std::cos(phi),
216  sinPhi = std::sin(phi);
217 
218  w0 = G4ThreeVector( dxFudgeTop*cosPhi, dyFudgeTop*sinPhi, +zTopCut );
219  w1 = G4ThreeVector( dxFudgeBot*cosPhi, dyFudgeBot*sinPhi, -zTopCut );
220  transform.ApplyPointTransform( w0 );
221  transform.ApplyPointTransform( w1 );
222 
223  //
224  // Add a point to our z ends
225  //
226  endPoly1.AddVertexInOrder( v0 );
227  endPoly2.AddVertexInOrder( v1 );
228 
229  //
230  // Build phi polygon
231  //
232  phiPoly.ClearAllVertices();
233 
234  phiPoly.AddVertexInOrder( v0 );
235  phiPoly.AddVertexInOrder( v1 );
236  phiPoly.AddVertexInOrder( w1 );
237  phiPoly.AddVertexInOrder( w0 );
238 
239  if (phiPoly.PartialClip( voxelLimit, axis ))
240  {
241  //
242  // Get unit normal
243  //
244  phiPoly.SetNormal( (v1-v0).cross(w0-v0).unit() );
245 
246  extentList.AddSurface( phiPoly );
247  }
248 
249  //
250  // Next vertex
251  //
252  v0 = w0;
253  v1 = w1;
254  } while( --numPhi > 0 );
255 
256  //
257  // Process the end pieces
258  //
259  if (endPoly1.PartialClip( voxelLimit, axis ))
260  {
261  static const G4ThreeVector normal(0,0,+1);
262  endPoly1.SetNormal( transform.TransformAxis(normal) );
263  extentList.AddSurface( endPoly1 );
264  }
265 
266  if (endPoly2.PartialClip( voxelLimit, axis ))
267  {
268  static const G4ThreeVector normal(0,0,-1);
269  endPoly2.SetNormal( transform.TransformAxis(normal) );
270  extentList.AddSurface( endPoly2 );
271  }
272 
273  //
274  // Return min/max value
275  //
276  return extentList.GetExtent( min, max );
277 }
278 
280 //
281 // Return whether point inside/outside/on surface
282 // Split into radius, phi, theta checks
283 // Each check modifies `in', or returns as approprate
284 //
286 {
287  G4double rad2oo, // outside surface outer tolerance
288  rad2oi; // outside surface inner tolerance
289 
290  EInside in;
291 
292  // check this side of z cut first, because that's fast
293  //
294 
295  if ( (p.z() < -zTopCut - halfCarTol)
296  || (p.z() > zTopCut + halfCarTol ) )
297  {
298  return in = kOutside;
299  }
300 
301  rad2oo= sqr(p.x()/( xSemiAxis + halfRadTol ))
302  + sqr(p.y()/( ySemiAxis + halfRadTol ));
303 
304  if ( rad2oo > sqr( zheight-p.z() ) )
305  {
306  return in = kOutside;
307  }
308 
309  // rad2oi= sqr( p.x()*(1.0 + 0.5*kRadTolerance/(xSemiAxis*xSemiAxis)) )
310  // + sqr( p.y()*(1.0 + 0.5*kRadTolerance/(ySemiAxis*ySemiAxis)) );
311  rad2oi = sqr(p.x()/( xSemiAxis - halfRadTol ))
312  + sqr(p.y()/( ySemiAxis - halfRadTol ));
313 
314  if (rad2oi < sqr( zheight-p.z() ) )
315  {
316  in = ( ( p.z() < -zTopCut + halfRadTol )
317  || ( p.z() > zTopCut - halfRadTol ) ) ? kSurface : kInside;
318  }
319  else
320  {
321  in = kSurface;
322  }
323 
324  return in;
325 }
326 
328 //
329 // Return unit normal of surface closest to p not protected against p=0
330 //
332 {
333 
334  G4double rx = sqr(p.x()/xSemiAxis),
335  ry = sqr(p.y()/ySemiAxis);
336 
337  G4double rds = std::sqrt(rx + ry);
338 
339  G4ThreeVector norm;
340 
341  if( (p.z() < -zTopCut) && ((rx+ry) < sqr(zTopCut + zheight)) )
342  {
343  return G4ThreeVector( 0., 0., -1. );
344  }
345 
346  if( (p.z() > (zheight > zTopCut ? zheight : zTopCut)) &&
347  ((rx+ry) < sqr(zheight-zTopCut)) )
348  {
349  return G4ThreeVector( 0., 0., 1. );
350  }
351 
352  if( p.z() > rds + 2.*zTopCut - zheight )
353  {
354  if ( p.z() > zTopCut )
355  {
356  if( p.x() == 0. )
357  {
358  norm = G4ThreeVector( 0., p.y() < 0. ? -1. : 1., 1. );
359  return norm /= norm.mag();
360  }
361  if( p.y() == 0. )
362  {
363  norm = G4ThreeVector( p.x() < 0. ? -1. : 1., 0., 1. );
364  return norm /= norm.mag();
365  }
366 
367  G4double k = std::fabs(p.x()/p.y());
369  G4double x = std::sqrt(c2);
370  G4double y = k*x;
371 
372  x /= sqr(xSemiAxis);
373  y /= sqr(ySemiAxis);
374 
375  norm = G4ThreeVector( p.x() < 0. ? -x : x,
376  p.y() < 0. ? -y : y,
377  - ( zheight - zTopCut ) );
378  norm /= norm.mag();
379  norm += G4ThreeVector( 0., 0., 1. );
380  return norm /= norm.mag();
381  }
382 
383  return G4ThreeVector( 0., 0., 1. );
384  }
385 
386  if( p.z() < rds - 2.*zTopCut - zheight )
387  {
388  if( p.x() == 0. )
389  {
390  norm = G4ThreeVector( 0., p.y() < 0. ? -1. : 1., -1. );
391  return norm /= norm.mag();
392  }
393  if( p.y() == 0. )
394  {
395  norm = G4ThreeVector( p.x() < 0. ? -1. : 1., 0., -1. );
396  return norm /= norm.mag();
397  }
398 
399  G4double k = std::fabs(p.x()/p.y());
401  G4double x = std::sqrt(c2);
402  G4double y = k*x;
403 
404  x /= sqr(xSemiAxis);
405  y /= sqr(ySemiAxis);
406 
407  norm = G4ThreeVector( p.x() < 0. ? -x : x,
408  p.y() < 0. ? -y : y,
409  - ( zheight - zTopCut ) );
410  norm /= norm.mag();
411  norm += G4ThreeVector( 0., 0., -1. );
412  return norm /= norm.mag();
413  }
414 
415  norm = G4ThreeVector(p.x()/sqr(xSemiAxis), p.y()/sqr(ySemiAxis), rds);
416 
417  G4double k = std::tan(pi/8.);
418  G4double c = -zTopCut - k*(zTopCut + zheight);
419 
420  if( p.z() < -k*rds + c )
421  return G4ThreeVector (0.,0.,-1.);
422 
423  return norm /= norm.mag();
424 }
425 
427 //
428 // Calculate distance to shape from outside, along normalised vector
429 // return kInfinity if no intersection, or intersection distance <= tolerance
430 //
432  const G4ThreeVector& v ) const
433 {
434  G4double distMin = kInfinity;
435 
436  // code from EllipticalTube
437 
438  G4double sigz = p.z()+zTopCut;
439 
440  //
441  // Check z = -dz planer surface
442  //
443 
444  if (sigz < halfCarTol)
445  {
446  //
447  // We are "behind" the shape in z, and so can
448  // potentially hit the rear face. Correct direction?
449  //
450  if (v.z() <= 0)
451  {
452  //
453  // As long as we are far enough away, we know we
454  // can't intersect
455  //
456  if (sigz < 0) return kInfinity;
457 
458  //
459  // Otherwise, we don't intersect unless we are
460  // on the surface of the ellipse
461  //
462 
463  if ( sqr(p.x()/( xSemiAxis - halfCarTol ))
464  + sqr(p.y()/( ySemiAxis - halfCarTol )) <= sqr( zheight+zTopCut ) )
465  return kInfinity;
466 
467  }
468  else
469  {
470  //
471  // How far?
472  //
473  G4double q = -sigz/v.z();
474 
475  //
476  // Where does that place us?
477  //
478  G4double xi = p.x() + q*v.x(),
479  yi = p.y() + q*v.y();
480 
481  //
482  // Is this on the surface (within ellipse)?
483  //
484  if ( sqr(xi/xSemiAxis) + sqr(yi/ySemiAxis) <= sqr( zheight + zTopCut ) )
485  {
486  //
487  // Yup. Return q, unless we are on the surface
488  //
489  return (sigz < -halfCarTol) ? q : 0;
490  }
491  else if (xi/(xSemiAxis*xSemiAxis)*v.x()
492  + yi/(ySemiAxis*ySemiAxis)*v.y() >= 0)
493  {
494  //
495  // Else, if we are traveling outwards, we know
496  // we must miss
497  //
498  // return kInfinity;
499  }
500  }
501  }
502 
503  //
504  // Check z = +dz planer surface
505  //
506  sigz = p.z() - zTopCut;
507 
508  if (sigz > -halfCarTol)
509  {
510  if (v.z() >= 0)
511  {
512 
513  if (sigz > 0) return kInfinity;
514 
515  if ( sqr(p.x()/( xSemiAxis - halfCarTol ))
516  + sqr(p.y()/( ySemiAxis - halfCarTol )) <= sqr( zheight-zTopCut ) )
517  return kInfinity;
518 
519  }
520  else {
521  G4double q = -sigz/v.z();
522 
523  G4double xi = p.x() + q*v.x(),
524  yi = p.y() + q*v.y();
525 
526  if ( sqr(xi/xSemiAxis) + sqr(yi/ySemiAxis) <= sqr( zheight - zTopCut ) )
527  {
528  return (sigz > -halfCarTol) ? q : 0;
529  }
530  else if (xi/(xSemiAxis*xSemiAxis)*v.x()
531  + yi/(ySemiAxis*ySemiAxis)*v.y() >= 0)
532  {
533  // return kInfinity;
534  }
535  }
536  }
537 
538 
539 #if 0
540 
541  // check to see if Z plane is relevant
542  //
543  if (p.z() < -zTopCut - 0.5*kCarTolerance)
544  {
545  if (v.z() <= 0.0)
546  return distMin;
547 
548  G4double lambda = (-zTopCut - p.z())/v.z();
549 
550  if ( sqr((lambda*v.x()+p.x())/xSemiAxis) +
551  sqr((lambda*v.y()+p.y())/ySemiAxis) <=
552  sqr(zTopCut + zheight + 0.5*kRadTolerance) )
553  {
554  return distMin = std::fabs(lambda);
555  }
556  }
557 
558  if (p.z() > zTopCut+0.5*kCarTolerance)
559  {
560  if (v.z() >= 0.0)
561  { return distMin; }
562 
563  G4double lambda = (zTopCut - p.z()) / v.z();
564 
565  if ( sqr((lambda*v.x() + p.x())/xSemiAxis) +
566  sqr((lambda*v.y() + p.y())/ySemiAxis) <=
567  sqr(zheight - zTopCut + 0.5*kRadTolerance) )
568  {
569  return distMin = std::fabs(lambda);
570  }
571  }
572 
573  if (p.z() > zTopCut - halfCarTol
574  && p.z() < zTopCut + halfCarTol )
575  {
576  if (v.z() > 0.)
577  { return kInfinity; }
578 
579  return distMin = 0.;
580  }
581 
582  if (p.z() < -zTopCut + halfCarTol
583  && p.z() > -zTopCut - halfCarTol)
584  {
585  if (v.z() < 0.)
586  { return distMin = kInfinity; }
587 
588  return distMin = 0.;
589  }
590 
591 #endif
592 
593  // if we are here then it either intersects or grazes the curved surface
594  // or it does not intersect at all
595  //
596  G4double A = sqr(v.x()/xSemiAxis) + sqr(v.y()/ySemiAxis) - sqr(v.z());
597  G4double B = 2*(v.x()*p.x()/sqr(xSemiAxis) +
598  v.y()*p.y()/sqr(ySemiAxis) + v.z()*(zheight-p.z()));
599  G4double C = sqr(p.x()/xSemiAxis) + sqr(p.y()/ySemiAxis) -
600  sqr(zheight - p.z());
601 
602  G4double discr = B*B - 4.*A*C;
603 
604  // if the discriminant is negative it never hits the curved object
605  //
606  if ( discr < -halfCarTol )
607  { return distMin; }
608 
609  // case below is when it hits or grazes the surface
610  //
611  if ( (discr >= - halfCarTol ) && (discr < halfCarTol ) )
612  {
613  return distMin = std::fabs(-B/(2.*A));
614  }
615 
616  G4double plus = (-B+std::sqrt(discr))/(2.*A);
617  G4double minus = (-B-std::sqrt(discr))/(2.*A);
618 
619  // Special case::Point on Surface, Check norm.dot(v)
620 
621  if ( ( std::fabs(plus) < halfCarTol )||( std::fabs(minus) < halfCarTol ) )
622  {
623  G4ThreeVector truenorm(p.x()/(xSemiAxis*xSemiAxis),
624  p.y()/(ySemiAxis*ySemiAxis),
625  -( p.z() - zheight ));
626  if ( truenorm*v >= 0) // going outside the solid from surface
627  {
628  return kInfinity;
629  }
630  else
631  {
632  return 0;
633  }
634  }
635 
636  // G4double lambda = std::fabs(plus) < std::fabs(minus) ? plus : minus;
637  G4double lambda = 0;
638 
639  if ( minus > halfCarTol && minus < distMin )
640  {
641  lambda = minus ;
642  // check normal vector n * v < 0
643  G4ThreeVector pin = p + lambda*v;
644  if(std::fabs(pin.z())<zTopCut+0.5*kCarTolerance)
645  {
646  G4ThreeVector truenorm(pin.x()/(xSemiAxis*xSemiAxis),
647  pin.y()/(ySemiAxis*ySemiAxis),
648  - ( pin.z() - zheight ));
649  if ( truenorm*v < 0)
650  { // yes, going inside the solid
651  distMin = lambda;
652  }
653  }
654  }
655  if ( plus > halfCarTol && plus < distMin )
656  {
657  lambda = plus ;
658  // check normal vector n * v < 0
659  G4ThreeVector pin = p + lambda*v;
660  if(std::fabs(pin.z())<zTopCut+0.5*kCarTolerance)
661  {
662  G4ThreeVector truenorm(pin.x()/(xSemiAxis*xSemiAxis),
663  pin.y()/(ySemiAxis*ySemiAxis),
664  - ( pin.z() - zheight ) );
665  if ( truenorm*v < 0)
666  { // yes, going inside the solid
667  distMin = lambda;
668  }
669  }
670  }
671  if (distMin < halfCarTol) distMin=0.;
672  return distMin ;
673 }
674 
676 //
677 // Calculate distance (<= actual) to closest surface of shape from outside
678 // Return 0 if point inside
679 //
681 {
682  G4double distR, distR2, distZ, maxDim;
683  G4double distRad;
684 
685  // check if the point lies either below z=-zTopCut in bottom elliptical
686  // region or on top within cut elliptical region
687  //
688  if( (p.z() <= -zTopCut) && (sqr(p.x()/xSemiAxis) + sqr(p.y()/ySemiAxis)
689  <= sqr(zTopCut + zheight + 0.5*kCarTolerance )) )
690  {
691  //return distZ = std::fabs(zTopCut - p.z());
692  return distZ = std::fabs(zTopCut + p.z());
693  }
694 
695  if( (p.z() >= zTopCut) && (sqr(p.x()/xSemiAxis)+sqr(p.y()/ySemiAxis)
696  <= sqr(zheight - zTopCut + kCarTolerance/2.0 )) )
697  {
698  return distZ = std::fabs(p.z() - zTopCut);
699  }
700 
701  // below we use the following approximation: we take the largest of the
702  // axes and find the shortest distance to the circular (cut) cone of that
703  // radius.
704  //
705  maxDim = xSemiAxis >= ySemiAxis ? xSemiAxis:ySemiAxis;
706  distRad = std::sqrt(p.x()*p.x()+p.y()*p.y());
707 
708  if( p.z() > maxDim*distRad + zTopCut*(1.+maxDim)-sqr(maxDim)*zheight )
709  {
710  distR2 = sqr(p.z() - zTopCut) + sqr(distRad - maxDim*(zheight - zTopCut));
711  return std::sqrt( distR2 );
712  }
713 
714  if( distRad > maxDim*( zheight - p.z() ) )
715  {
716  if( p.z() > maxDim*distRad - (zTopCut*(1.+maxDim)+sqr(maxDim)*zheight) )
717  {
718  G4double zVal = (p.z()-maxDim*(distRad-maxDim*zheight))/(1.+sqr(maxDim));
719  G4double rVal = maxDim*(zheight - zVal);
720  return distR = std::sqrt(sqr(p.z() - zVal) + sqr(distRad - rVal));
721  }
722  }
723 
724  if( distRad <= maxDim*(zheight - p.z()) )
725  {
726  distR2 = sqr(distRad - maxDim*(zheight + zTopCut)) + sqr(p.z() + zTopCut);
727  return std::sqrt( distR2 );
728  }
729 
730  return distR = 0;
731 }
732 
734 //
735 // Calculate distance to surface of shape from `inside',
736 // allowing for tolerance
737 //
739  const G4ThreeVector& v,
740  const G4bool calcNorm,
741  G4bool *validNorm,
742  G4ThreeVector *n ) const
743 {
744  G4double distMin, lambda;
745  enum surface_e {kPlaneSurf, kCurvedSurf, kNoSurf} surface;
746 
747  distMin = kInfinity;
748  surface = kNoSurf;
749 
750  if (v.z() < 0.0)
751  {
752  lambda = (-p.z() - zTopCut)/v.z();
753 
754  if ( (sqr((p.x() + lambda*v.x())/xSemiAxis) +
755  sqr((p.y() + lambda*v.y())/ySemiAxis)) <
756  sqr(zheight + zTopCut + 0.5*kCarTolerance) )
757  {
758  distMin = std::fabs(lambda);
759 
760  if (!calcNorm) { return distMin; }
761  }
762  distMin = std::fabs(lambda);
763  surface = kPlaneSurf;
764  }
765 
766  if (v.z() > 0.0)
767  {
768  lambda = (zTopCut - p.z()) / v.z();
769 
770  if ( (sqr((p.x() + lambda*v.x())/xSemiAxis)
771  + sqr((p.y() + lambda*v.y())/ySemiAxis) )
772  < (sqr(zheight - zTopCut + 0.5*kCarTolerance)) )
773  {
774  distMin = std::fabs(lambda);
775  if (!calcNorm) { return distMin; }
776  }
777  distMin = std::fabs(lambda);
778  surface = kPlaneSurf;
779  }
780 
781  // if we are here then it either intersects or grazes the
782  // curved surface...
783  //
784  G4double A = sqr(v.x()/xSemiAxis) + sqr(v.y()/ySemiAxis) - sqr(v.z());
785  G4double B = 2.*(v.x()*p.x()/sqr(xSemiAxis) +
786  v.y()*p.y()/sqr(ySemiAxis) + v.z()*(zheight-p.z()));
787  G4double C = sqr(p.x()/xSemiAxis) + sqr(p.y()/ySemiAxis)
788  - sqr(zheight - p.z());
789 
790  G4double discr = B*B - 4.*A*C;
791 
792  if ( discr >= - 0.5*kCarTolerance && discr < 0.5*kCarTolerance )
793  {
794  if(!calcNorm) { return distMin = std::fabs(-B/(2.*A)); }
795  }
796 
797  else if ( discr > 0.5*kCarTolerance )
798  {
799  G4double plus = (-B+std::sqrt(discr))/(2.*A);
800  G4double minus = (-B-std::sqrt(discr))/(2.*A);
801 
802  if ( plus > 0.5*kCarTolerance && minus > 0.5*kCarTolerance )
803  {
804  // take the shorter distance
805  //
806  lambda = std::fabs(plus) < std::fabs(minus) ? plus : minus;
807  }
808  else
809  {
810  // at least one solution is close to zero or negative
811  // so, take small positive solution or zero
812  //
813  lambda = plus > -0.5*kCarTolerance ? plus : 0;
814  }
815 
816  if ( std::fabs(lambda) < distMin )
817  {
818  if( std::fabs(lambda) > 0.5*kCarTolerance)
819  {
820  distMin = std::fabs(lambda);
821  surface = kCurvedSurf;
822  }
823  else // Point is On the Surface, Check Normal
824  {
825  G4ThreeVector truenorm(p.x()/(xSemiAxis*xSemiAxis),
826  p.y()/(ySemiAxis*ySemiAxis),
827  -( p.z() - zheight ));
828  if( truenorm.dot(v) > 0 )
829  {
830  distMin = 0.0;
831  surface = kCurvedSurf;
832  }
833  }
834  }
835  }
836 
837  // set normal if requested
838  //
839  if (calcNorm)
840  {
841  if (surface == kNoSurf)
842  {
843  *validNorm = false;
844  }
845  else
846  {
847  *validNorm = true;
848  switch (surface)
849  {
850  case kPlaneSurf:
851  {
852  *n = G4ThreeVector(0.,0.,(v.z() > 0.0 ? 1. : -1.));
853  }
854  break;
855 
856  case kCurvedSurf:
857  {
858  G4ThreeVector pexit = p + distMin*v;
859  G4ThreeVector truenorm( pexit.x()/(xSemiAxis*xSemiAxis),
860  pexit.y()/(ySemiAxis*ySemiAxis),
861  -( pexit.z() - zheight ) );
862  truenorm /= truenorm.mag();
863  *n= truenorm;
864  }
865  break;
866 
867  default: // Should never reach this case ...
868  DumpInfo();
869  std::ostringstream message;
870  G4int oldprc = message.precision(16);
871  message << "Undefined side for valid surface normal to solid."
872  << G4endl
873  << "Position:" << G4endl
874  << " p.x() = " << p.x()/mm << " mm" << G4endl
875  << " p.y() = " << p.y()/mm << " mm" << G4endl
876  << " p.z() = " << p.z()/mm << " mm" << G4endl
877  << "Direction:" << G4endl
878  << " v.x() = " << v.x() << G4endl
879  << " v.y() = " << v.y() << G4endl
880  << " v.z() = " << v.z() << G4endl
881  << "Proposed distance :" << G4endl
882  << " distMin = " << distMin/mm << " mm";
883  message.precision(oldprc);
884  G4Exception("G4EllipticalCone::DistanceToOut(p,v,..)",
885  "GeomSolids1002", JustWarning, message);
886  break;
887  }
888  }
889  }
890 
891  if (distMin<0.5*kCarTolerance) { distMin=0; }
892 
893  return distMin;
894 }
895 
897 //
898 // Calculate distance (<=actual) to closest surface of shape from inside
899 //
901 {
902  G4double rds,roo,roo1, distR, distZ, distMin=0.;
904 
905 #ifdef G4SPECSDEBUG
906  if( Inside(p) == kOutside )
907  {
908  DumpInfo();
909  std::ostringstream message;
910  G4int oldprc = message.precision(16);
911  message << "Point p is outside !?" << G4endl
912  << "Position:" << G4endl
913  << " p.x() = " << p.x()/mm << " mm" << G4endl
914  << " p.y() = " << p.y()/mm << " mm" << G4endl
915  << " p.z() = " << p.z()/mm << " mm";
916  message.precision(oldprc) ;
917  G4Exception("G4Ellipsoid::DistanceToOut(p)", "GeomSolids1002",
918  JustWarning, message);
919  }
920 #endif
921 
922  // since we have made the above warning, below we are working assuming p
923  // is inside check how close it is to the circular cone with radius equal
924  // to the smaller of the axes
925  //
926  if( sqr(p.x()/minAxis)+sqr(p.y()/minAxis) < sqr(zheight - p.z()) )
927  {
928  rds = std::sqrt(sqr(p.x()) + sqr(p.y()));
929  roo = minAxis*(zheight-p.z()); // radius of cone at z= p.z()
930  roo1 = minAxis*(zheight-zTopCut); // radius of cone at z=+zTopCut
931 
932  distZ=zTopCut - std::fabs(p.z()) ;
933  distR=(roo-rds)/(std::sqrt(1+sqr(minAxis)));
934 
935  if(rds>roo1)
936  {
937  distMin=(zTopCut-p.z())*(roo-rds)/(roo-roo1);
938  distMin=std::min(distMin,distR);
939  }
940  distMin=std::min(distR,distZ);
941  }
942 
943  return distMin;
944 }
945 
947 //
948 // GetEntityType
949 //
951 {
952  return G4String("G4EllipticalCone");
953 }
954 
956 //
957 // Make a clone of the object
958 //
960 {
961  return new G4EllipticalCone(*this);
962 }
963 
965 //
966 // Stream object contents to an output stream
967 //
968 std::ostream& G4EllipticalCone::StreamInfo( std::ostream& os ) const
969 {
970  G4int oldprc = os.precision(16);
971  os << "-----------------------------------------------------------\n"
972  << " *** Dump for solid - " << GetName() << " ***\n"
973  << " ===================================================\n"
974  << " Solid type: G4EllipticalCone\n"
975  << " Parameters: \n"
976 
977  << " semi-axis x: " << xSemiAxis/mm << " mm \n"
978  << " semi-axis y: " << ySemiAxis/mm << " mm \n"
979  << " height z: " << zheight/mm << " mm \n"
980  << " half length in z: " << zTopCut/mm << " mm \n"
981  << "-----------------------------------------------------------\n";
982  os.precision(oldprc);
983 
984  return os;
985 }
986 
988 //
989 // GetPointOnSurface
990 //
991 // returns quasi-uniformly distributed point on surface of elliptical cone
992 //
994 {
995 
996  G4double phi, sinphi, cosphi, aOne, aTwo, aThree,
997  chose, zRand, rRand1, rRand2;
998 
999  G4double rOne = std::sqrt(sqr(xSemiAxis)
1000  + sqr(ySemiAxis))*(zheight - zTopCut);
1001  G4double rTwo = std::sqrt(sqr(xSemiAxis)
1002  + sqr(ySemiAxis))*(zheight + zTopCut);
1003 
1004  aOne = pi*(rOne + rTwo)*std::sqrt(sqr(rOne - rTwo)+sqr(2.*zTopCut));
1006  aThree = pi*xSemiAxis*ySemiAxis*sqr(zheight-zTopCut);
1007 
1008  phi = RandFlat::shoot(0.,twopi);
1009  cosphi = std::cos(phi);
1010  sinphi = std::sin(phi);
1011 
1012  if(zTopCut >= zheight) aThree = 0.;
1013 
1014  chose = RandFlat::shoot(0.,aOne+aTwo+aThree);
1015  if((chose>=0.) && (chose<aOne))
1016  {
1017  zRand = RandFlat::shoot(-zTopCut,zTopCut);
1018  return G4ThreeVector(xSemiAxis*(zheight-zRand)*cosphi,
1019  ySemiAxis*(zheight-zRand)*sinphi,zRand);
1020  }
1021  else if((chose>=aOne) && (chose<aOne+aTwo))
1022  {
1023  do
1024  {
1025  rRand1 = RandFlat::shoot(0.,1.) ;
1026  rRand2 = RandFlat::shoot(0.,1.) ;
1027  } while ( rRand2 >= rRand1 ) ;
1028 
1029  // rRand2 = RandFlat::shoot(0.,std::sqrt(1.-sqr(rRand1)));
1030  return G4ThreeVector(rRand1*xSemiAxis*(zheight+zTopCut)*cosphi,
1031  rRand1*ySemiAxis*(zheight+zTopCut)*sinphi, -zTopCut);
1032 
1033  }
1034  // else
1035  //
1036 
1037  do
1038  {
1039  rRand1 = RandFlat::shoot(0.,1.) ;
1040  rRand2 = RandFlat::shoot(0.,1.) ;
1041  } while ( rRand2 >= rRand1 ) ;
1042 
1043  return G4ThreeVector(rRand1*xSemiAxis*(zheight-zTopCut)*cosphi,
1044  rRand1*ySemiAxis*(zheight-zTopCut)*sinphi, zTopCut);
1045 }
1046 
1047 //
1048 // Methods for visualisation
1049 //
1050 
1052 {
1053  scene.AddSolid(*this);
1054 }
1055 
1057 {
1058  // Define the sides of the box into which the solid instance would fit.
1059  //
1060  G4double maxDim;
1061  maxDim = xSemiAxis > ySemiAxis ? xSemiAxis : ySemiAxis;
1062  maxDim = maxDim > zTopCut ? maxDim : zTopCut;
1063 
1064  return G4VisExtent (-maxDim, maxDim,
1065  -maxDim, maxDim,
1066  -maxDim, maxDim);
1067 }
1068 
1070 {
1072 }
1073 
1075 {
1076  if ( (!fpPolyhedron)
1078  fpPolyhedron->GetNumberOfRotationSteps()) )
1079  {
1080  delete fpPolyhedron;
1082  }
1083  return fpPolyhedron;
1084 }
G4String GetName() const
G4VSolid * Clone() const
ThreeVector shoot(const G4int Ap, const G4int Af)
static c2_factory< G4double > c2
G4EllipticalCone & operator=(const G4EllipticalCone &rhs)
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
G4bool GetExtent(G4double &min, G4double &max) const
const G4double pi
EInside Inside(const G4ThreeVector &p) const
void SetNormal(const G4ThreeVector &newNormal)
virtual G4bool PartialClip(const G4VoxelLimits &voxelLimit, const EAxis IgnoreMe)
virtual void AddVertexInOrder(const G4ThreeVector vertex)
G4double a
Definition: TRTMaterials.hh:39
virtual void ClearAllVertices()
virtual void AddSolid(const G4Box &)=0
G4Polyhedron * CreatePolyhedron() const
int G4int
Definition: G4Types.hh:78
void DumpInfo() const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
std::ostream & StreamInfo(std::ostream &os) const
G4ThreeVector GetPointOnSurface() const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
G4Polyhedron * fpPolyhedron
G4double GetRadialTolerance() const
bool G4bool
Definition: G4Types.hh:79
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4EllipticalCone(const G4String &pName, G4double pxSemiAxis, G4double pySemiAxis, G4double zMax, G4double pzTopCut)
void AddSurface(const G4ClippablePolygon &surface)
const G4int n
static const G4double A[nN]
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void DescribeYourselfTo(G4VGraphicsScene &scene) const
void SetSemiAxis(G4double x, G4double y, G4double z)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
void SetZCut(G4double newzTopCut)
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
#define G4endl
Definition: G4ios.hh:61
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
double G4double
Definition: G4Types.hh:76
G4Polyhedron * GetPolyhedron() const
G4GeometryType GetEntityType() const
G4VisExtent GetExtent() const
virtual ~G4EllipticalCone()
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
static const double mm
Definition: G4SIunits.hh:102
static G4GeometryTolerance * GetInstance()
void ApplyPointTransform(G4ThreeVector &vec) const
const G4int kMaxMeshSections
Definition: meshdefs.hh:46