Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Paraboloid.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: G4Paraboloid.cc 67011 2013-01-29 16:17:41Z gcosmo $
27 //
28 // class G4Paraboloid
29 //
30 // Implementation for G4Paraboloid class
31 //
32 // Author : Lukas Lindroos (CERN), July 2007
33 // Revised: Tatiana Nikitina (CERN)
34 // --------------------------------------------------------------------
35 
36 #include "globals.hh"
37 
38 #include "G4Paraboloid.hh"
39 
40 #include "G4VoxelLimits.hh"
41 #include "G4AffineTransform.hh"
42 
43 #include "meshdefs.hh"
44 
45 #include "Randomize.hh"
46 
47 #include "G4VGraphicsScene.hh"
48 #include "G4Polyhedron.hh"
49 #include "G4NURBS.hh"
50 #include "G4NURBSbox.hh"
51 #include "G4VisExtent.hh"
52 
53 using namespace CLHEP;
54 
56 //
57 // constructor - check parameters
58 
60  G4double pDz,
61  G4double pR1,
62  G4double pR2)
63  : G4VSolid(pName), fpPolyhedron(0), fSurfaceArea(0.), fCubicVolume(0.)
64 
65 {
66  if( (pDz <= 0.) || (pR2 <= pR1) || (pR1 < 0.) )
67  {
68  std::ostringstream message;
69  message << "Invalid dimensions. Negative Input Values or R1>=R2 - "
70  << GetName();
71  G4Exception("G4Paraboloid::G4Paraboloid()", "GeomSolids0002",
72  FatalErrorInArgument, message,
73  "Z half-length must be larger than zero or R1>=R2.");
74  }
75 
76  r1 = pR1;
77  r2 = pR2;
78  dz = pDz;
79 
80  // r1^2 = k1 * (-dz) + k2
81  // r2^2 = k1 * ( dz) + k2
82  // => r1^2 + r2^2 = k2 + k2 => k2 = (r2^2 + r1^2) / 2
83  // and r2^2 - r1^2 = k1 * dz - k1 * (-dz) => k1 = (r2^2 - r1^2) / 2 / dz
84 
85  k1 = (r2 * r2 - r1 * r1) / 2 / dz;
86  k2 = (r2 * r2 + r1 * r1) / 2;
87 }
88 
90 //
91 // Fake default constructor - sets only member data and allocates memory
92 // for usage restricted to object persistency.
93 //
95  : G4VSolid(a), fpPolyhedron(0), fSurfaceArea(0.), fCubicVolume(0.),
96  dz(0.), r1(0.), r2(0.), k1(0.), k2(0.)
97 {
98 }
99 
101 //
102 // Destructor
103 
105 {
106 }
107 
109 //
110 // Copy constructor
111 
113  : G4VSolid(rhs), fpPolyhedron(0),
114  fSurfaceArea(rhs.fSurfaceArea), fCubicVolume(rhs.fCubicVolume),
115  dz(rhs.dz), r1(rhs.r1), r2(rhs.r2), k1(rhs.k1), k2(rhs.k2)
116 {
117 }
118 
119 
121 //
122 // Assignment operator
123 
125 {
126  // Check assignment to self
127  //
128  if (this == &rhs) { return *this; }
129 
130  // Copy base class data
131  //
132  G4VSolid::operator=(rhs);
133 
134  // Copy data
135  //
136  fpPolyhedron = 0;
137  fSurfaceArea = rhs.fSurfaceArea; fCubicVolume = rhs.fCubicVolume;
138  dz = rhs.dz; r1 = rhs.r1; r2 = rhs.r2; k1 = rhs.k1; k2 = rhs.k2;
139 
140  return *this;
141 }
142 
144 //
145 // Dispatch to parameterisation for replication mechanism dimension
146 // computation & modification.
147 
148 //void ComputeDimensions( G4VPVParamerisation p,
149 // const G4Int n,
150 // const G4VPhysicalVolume* pRep )
151 //{
152 // p->ComputeDimensions(*this,n,pRep) ;
153 //}
154 
155 
157 //
158 // Calculate extent under transform and specified limit
159 
160 G4bool
162  const G4VoxelLimits& pVoxelLimit,
163  const G4AffineTransform& pTransform,
164  G4double& pMin, G4double& pMax) const
165 {
166  G4double xMin = -r2 + pTransform.NetTranslation().x(),
167  xMax = r2 + pTransform.NetTranslation().x(),
168  yMin = -r2 + pTransform.NetTranslation().y(),
169  yMax = r2 + pTransform.NetTranslation().y(),
170  zMin = -dz + pTransform.NetTranslation().z(),
171  zMax = dz + pTransform.NetTranslation().z();
172 
173  if(!pTransform.IsRotated()
174  || pTransform.NetRotation()(G4ThreeVector(0, 0, 1)) == G4ThreeVector(0, 0, 1))
175  {
176  if(pVoxelLimit.IsXLimited())
177  {
178  if(pVoxelLimit.GetMaxXExtent() < xMin - 0.5 * kCarTolerance
179  || pVoxelLimit.GetMinXExtent() > xMax + 0.5 * kCarTolerance)
180  {
181  return false;
182  }
183  else
184  {
185  if(pVoxelLimit.GetMinXExtent() > xMin)
186  {
187  xMin = pVoxelLimit.GetMinXExtent();
188  }
189  if(pVoxelLimit.GetMaxXExtent() < xMax)
190  {
191  xMax = pVoxelLimit.GetMaxXExtent();
192  }
193  }
194  }
195  if(pVoxelLimit.IsYLimited())
196  {
197  if(pVoxelLimit.GetMaxYExtent() < yMin - 0.5 * kCarTolerance
198  || pVoxelLimit.GetMinYExtent() > yMax + 0.5 * kCarTolerance)
199  {
200  return false;
201  }
202  else
203  {
204  if(pVoxelLimit.GetMinYExtent() > yMin)
205  {
206  yMin = pVoxelLimit.GetMinYExtent();
207  }
208  if(pVoxelLimit.GetMaxYExtent() < yMax)
209  {
210  yMax = pVoxelLimit.GetMaxYExtent();
211  }
212  }
213  }
214  if(pVoxelLimit.IsZLimited())
215  {
216  if(pVoxelLimit.GetMaxZExtent() < zMin - 0.5 * kCarTolerance
217  || pVoxelLimit.GetMinZExtent() > zMax + 0.5 * kCarTolerance)
218  {
219  return false;
220  }
221  else
222  {
223  if(pVoxelLimit.GetMinZExtent() > zMin)
224  {
225  zMin = pVoxelLimit.GetMinZExtent();
226  }
227  if(pVoxelLimit.GetMaxZExtent() < zMax)
228  {
229  zMax = pVoxelLimit.GetMaxZExtent();
230  }
231  }
232  }
233  switch(pAxis)
234  {
235  case kXAxis:
236  pMin = xMin;
237  pMax = xMax;
238  break;
239  case kYAxis:
240  pMin = yMin;
241  pMax = yMax;
242  break;
243  case kZAxis:
244  pMin = zMin;
245  pMax = zMax;
246  break;
247  default:
248  pMin = 0;
249  pMax = 0;
250  return false;
251  }
252  }
253  else
254  {
255  G4bool existsAfterClip=true;
256 
257  // Calculate rotated vertex coordinates
258 
259  G4int noPolygonVertices=0;
260  G4ThreeVectorList* vertices
261  = CreateRotatedVertices(pTransform,noPolygonVertices);
262 
263  if(pAxis == kXAxis || pAxis == kYAxis || pAxis == kZAxis)
264  {
265 
266  pMin = kInfinity;
267  pMax = -kInfinity;
268 
269  for(G4ThreeVectorList::iterator it = vertices->begin();
270  it < vertices->end(); it++)
271  {
272  if(pMin > (*it)[pAxis]) pMin = (*it)[pAxis];
273  if((*it)[pAxis] < pVoxelLimit.GetMinExtent(pAxis))
274  {
275  pMin = pVoxelLimit.GetMinExtent(pAxis);
276  }
277  if(pMax < (*it)[pAxis])
278  {
279  pMax = (*it)[pAxis];
280  }
281  if((*it)[pAxis] > pVoxelLimit.GetMaxExtent(pAxis))
282  {
283  pMax = pVoxelLimit.GetMaxExtent(pAxis);
284  }
285  }
286 
287  if(pMin > pVoxelLimit.GetMaxExtent(pAxis)
288  || pMax < pVoxelLimit.GetMinExtent(pAxis)) { existsAfterClip = false; }
289  }
290  else
291  {
292  pMin = 0;
293  pMax = 0;
294  existsAfterClip = false;
295  }
296  delete vertices;
297  return existsAfterClip;
298  }
299  return true;
300 }
301 
303 //
304 // Return whether point inside/outside/on surface
305 
307 {
308  // First check is the point is above or below the solid.
309  //
310  if(std::fabs(p.z()) > dz + 0.5 * kCarTolerance) { return kOutside; }
311 
312  G4double rho2 = p.perp2(),
313  rhoSurfTimesTol2 = (k1 * p.z() + k2) * sqr(kCarTolerance),
314  A = rho2 - ((k1 *p.z() + k2) + 0.25 * kCarTolerance * kCarTolerance);
315 
316  if(A < 0 && sqr(A) > rhoSurfTimesTol2)
317  {
318  // Actually checking rho < radius of paraboloid at z = p.z().
319  // We're either inside or in lower/upper cutoff area.
320 
321  if(std::fabs(p.z()) > dz - 0.5 * kCarTolerance)
322  {
323  // We're in the upper/lower cutoff area, sides have a paraboloid shape
324  // maybe further checks should be made to make these nicer
325 
326  return kSurface;
327  }
328  else
329  {
330  return kInside;
331  }
332  }
333  else if(A <= 0 || sqr(A) < rhoSurfTimesTol2)
334  {
335  // We're in the parabolic surface.
336 
337  return kSurface;
338  }
339  else
340  {
341  return kOutside;
342  }
343 }
344 
346 //
347 
349 {
350  G4ThreeVector n(0, 0, 0);
351  if(std::fabs(p.z()) > dz + 0.5 * kCarTolerance)
352  {
353  // If above or below just return normal vector for the cutoff plane.
354 
355  n = G4ThreeVector(0, 0, p.z()/std::fabs(p.z()));
356  }
357  else if(std::fabs(p.z()) > dz - 0.5 * kCarTolerance)
358  {
359  // This means we're somewhere in the plane z = dz or z = -dz.
360  // (As far as the program is concerned anyway.
361 
362  if(p.z() < 0) // Are we in upper or lower plane?
363  {
364  if(p.perp2() > sqr(r1 + 0.5 * kCarTolerance))
365  {
366  n = G4ThreeVector(p.x(), p.y(), -k1 / 2).unit();
367  }
368  else if(r1 < 0.5 * kCarTolerance
369  || p.perp2() > sqr(r1 - 0.5 * kCarTolerance))
370  {
371  n = G4ThreeVector(p.x(), p.y(), 0.).unit()
372  + G4ThreeVector(0., 0., -1.).unit();
373  n = n.unit();
374  }
375  else
376  {
377  n = G4ThreeVector(0., 0., -1.);
378  }
379  }
380  else
381  {
382  if(p.perp2() > sqr(r2 + 0.5 * kCarTolerance))
383  {
384  n = G4ThreeVector(p.x(), p.y(), 0.).unit();
385  }
386  else if(r2 < 0.5 * kCarTolerance
387  || p.perp2() > sqr(r2 - 0.5 * kCarTolerance))
388  {
389  n = G4ThreeVector(p.x(), p.y(), 0.).unit()
390  + G4ThreeVector(0., 0., 1.).unit();
391  n = n.unit();
392  }
393  else
394  {
395  n = G4ThreeVector(0., 0., 1.);
396  }
397  }
398  }
399  else
400  {
401  G4double rho2 = p.perp2();
402  G4double rhoSurfTimesTol2 = (k1 * p.z() + k2) * sqr(kCarTolerance);
403  G4double A = rho2 - ((k1 *p.z() + k2)
404  + 0.25 * kCarTolerance * kCarTolerance);
405 
406  if(A < 0 && sqr(A) > rhoSurfTimesTol2)
407  {
408  // Actually checking rho < radius of paraboloid at z = p.z().
409  // We're inside.
410 
411  if(p.mag2() != 0) { n = p.unit(); }
412  }
413  else if(A <= 0 || sqr(A) < rhoSurfTimesTol2)
414  {
415  // We're in the parabolic surface.
416 
417  n = G4ThreeVector(p.x(), p.y(), - k1 / 2).unit();
418  }
419  else
420  {
421  n = G4ThreeVector(p.x(), p.y(), - k1 / 2).unit();
422  }
423  }
424 
425  if(n.mag2() == 0)
426  {
427  std::ostringstream message;
428  message << "No normal defined for this point p." << G4endl
429  << " p = " << 1 / mm * p << " mm";
430  G4Exception("G4Paraboloid::SurfaceNormal(p)", "GeomSolids1002",
431  JustWarning, message);
432  }
433  return n;
434 }
435 
437 //
438 // Calculate distance to shape from outside, along normalised vector
439 // - return kInfinity if no intersection
440 //
441 
443  const G4ThreeVector& v ) const
444 {
445  G4double rho2 = p.perp2(), paraRho2 = std::fabs(k1 * p.z() + k2);
447  G4double tolh = 0.5*kCarTolerance;
448 
449  if(r2 && p.z() > - tolh + dz)
450  {
451  // If the points is above check for intersection with upper edge.
452 
453  if(v.z() < 0)
454  {
455  G4double intersection = (dz - p.z()) / v.z(); // With plane z = dz.
456  if(sqr(p.x() + v.x()*intersection)
457  + sqr(p.y() + v.y()*intersection) < sqr(r2 + 0.5 * kCarTolerance))
458  {
459  if(p.z() < tolh + dz)
460  { return 0; }
461  else
462  { return intersection; }
463  }
464  }
465  else // Direction away, no possibility of intersection
466  {
467  return kInfinity;
468  }
469  }
470  else if(r1 && p.z() < tolh - dz)
471  {
472  // If the points is belove check for intersection with lower edge.
473 
474  if(v.z() > 0)
475  {
476  G4double intersection = (-dz - p.z()) / v.z(); // With plane z = -dz.
477  if(sqr(p.x() + v.x()*intersection)
478  + sqr(p.y() + v.y()*intersection) < sqr(r1 + 0.5 * kCarTolerance))
479  {
480  if(p.z() > -tolh - dz)
481  {
482  return 0;
483  }
484  else
485  {
486  return intersection;
487  }
488  }
489  }
490  else // Direction away, no possibility of intersection
491  {
492  return kInfinity;
493  }
494  }
495 
496  G4double A = k1 / 2 * v.z() - p.x() * v.x() - p.y() * v.y(),
497  vRho2 = v.perp2(), intersection,
498  B = (k1 * p.z() + k2 - rho2) * vRho2;
499 
500  if ( ( (rho2 > paraRho2) && (sqr(rho2-paraRho2-0.25*tol2) > tol2*paraRho2) )
501  || (p.z() < - dz+kCarTolerance)
502  || (p.z() > dz-kCarTolerance) ) // Make sure it's safely outside.
503  {
504  // Is there a problem with squaring rho twice?
505 
506  if(vRho2<tol2) // Needs to be treated seperately.
507  {
508  intersection = ((rho2 - k2)/k1 - p.z())/v.z();
509  if(intersection < 0) { return kInfinity; }
510  else if(std::fabs(p.z() + v.z() * intersection) <= dz)
511  {
512  return intersection;
513  }
514  else
515  {
516  return kInfinity;
517  }
518  }
519  else if(A*A + B < 0) // No real intersections.
520  {
521  return kInfinity;
522  }
523  else
524  {
525  intersection = (A - std::sqrt(B + sqr(A))) / vRho2;
526  if(intersection < 0)
527  {
528  return kInfinity;
529  }
530  else if(std::fabs(p.z() + intersection * v.z()) < dz + tolh)
531  {
532  return intersection;
533  }
534  else
535  {
536  return kInfinity;
537  }
538  }
539  }
540  else if(sqr(rho2 - paraRho2 - .25 * tol2) <= tol2 * paraRho2)
541  {
542  // If this is true we're somewhere in the border.
543 
544  G4ThreeVector normal(p.x(), p.y(), -k1/2);
545  if(normal.dot(v) <= 0)
546  { return 0; }
547  else
548  { return kInfinity; }
549  }
550  else
551  {
552  std::ostringstream message;
553  if(Inside(p) == kInside)
554  {
555  message << "Point p is inside! - " << GetName() << G4endl;
556  }
557  else
558  {
559  message << "Likely a problem in this function, for solid: " << GetName()
560  << G4endl;
561  }
562  message << " p = " << p * (1/mm) << " mm" << G4endl
563  << " v = " << v * (1/mm) << " mm";
564  G4Exception("G4Paraboloid::DistanceToIn(p,v)", "GeomSolids1002",
565  JustWarning, message);
566  return 0;
567  }
568 }
569 
571 //
572 // Calculate distance (<= actual) to closest surface of shape from outside
573 // - Return 0 if point inside
574 
576 {
577  G4double safz = -dz+std::fabs(p.z());
578  if(safz<0) { safz=0; }
579  G4double safr = kInfinity;
580 
581  G4double rho = p.x()*p.x()+p.y()*p.y();
582  G4double paraRho = (p.z()-k2)/k1;
583  G4double sqrho = std::sqrt(rho);
584 
585  if(paraRho<0)
586  {
587  safr=sqrho-r2;
588  if(safr>safz) { safz=safr; }
589  return safz;
590  }
591 
592  G4double sqprho = std::sqrt(paraRho);
593  G4double dRho = sqrho-sqprho;
594  if(dRho<0) { return safz; }
595 
596  G4double talf = -2.*k1*sqprho;
597  G4double tmp = 1+talf*talf;
598  if(tmp<0.) { return safz; }
599 
600  G4double salf = talf/std::sqrt(tmp);
601  safr = std::fabs(dRho*salf);
602  if(safr>safz) { safz=safr; }
603 
604  return safz;
605 }
606 
608 //
609 // Calculate distance to surface of shape from 'inside'
610 
612  const G4ThreeVector& v,
613  const G4bool calcNorm,
614  G4bool *validNorm,
615  G4ThreeVector *n ) const
616 {
617  G4double rho2 = p.perp2(), paraRho2 = std::fabs(k1 * p.z() + k2);
618  G4double vRho2 = v.perp2(), intersection;
620  G4double tolh = 0.5*kCarTolerance;
621 
622  if(calcNorm) { *validNorm = false; }
623 
624  // We have that the particle p follows the line x = p + s * v
625  // meaning x = p.x() + s * v.x(), y = p.y() + s * v.y() and
626  // z = p.z() + s * v.z()
627  // The equation for all points on the surface (surface expanded for
628  // to include all z) x^2 + y^2 = k1 * z + k2 => .. =>
629  // => s = (A +- std::sqrt(A^2 + B)) / vRho2
630  // where:
631  //
632  G4double A = k1 / 2 * v.z() - p.x() * v.x() - p.y() * v.y();
633  //
634  // and:
635  //
636  G4double B = (-rho2 + paraRho2) * vRho2;
637 
638  if ( rho2 < paraRho2 && sqr(rho2 - paraRho2 - 0.25 * tol2) > tol2 * paraRho2
639  && std::fabs(p.z()) < dz - kCarTolerance)
640  {
641  // Make sure it's safely inside.
642 
643  if(v.z() > 0)
644  {
645  // It's heading upwards, check where it colides with the plane z = dz.
646  // When it does, is that in the surface of the paraboloid.
647  // z = p.z() + variable * v.z() for all points the particle can go.
648  // => variable = (z - p.z()) / v.z() so intersection must be:
649 
650  intersection = (dz - p.z()) / v.z();
651  G4ThreeVector ip = p + intersection * v; // Point of intersection.
652 
653  if(ip.perp2() < sqr(r2 + kCarTolerance))
654  {
655  if(calcNorm)
656  {
657  *n = G4ThreeVector(0, 0, 1);
658  if(r2 < tolh || ip.perp2() > sqr(r2 - tolh))
659  {
660  *n += G4ThreeVector(ip.x(), ip.y(), - k1 / 2).unit();
661  *n = n->unit();
662  }
663  *validNorm = true;
664  }
665  return intersection;
666  }
667  }
668  else if(v.z() < 0)
669  {
670  // It's heading downwards, check were it colides with the plane z = -dz.
671  // When it does, is that in the surface of the paraboloid.
672  // z = p.z() + variable * v.z() for all points the particle can go.
673  // => variable = (z - p.z()) / v.z() so intersection must be:
674 
675  intersection = (-dz - p.z()) / v.z();
676  G4ThreeVector ip = p + intersection * v; // Point of intersection.
677 
678  if(ip.perp2() < sqr(r1 + tolh))
679  {
680  if(calcNorm)
681  {
682  *n = G4ThreeVector(0, 0, -1);
683  if(r1 < tolh || ip.perp2() > sqr(r1 - tolh))
684  {
685  *n += G4ThreeVector(ip.x(), ip.y(), - k1 / 2).unit();
686  *n = n->unit();
687  }
688  *validNorm = true;
689  }
690  return intersection;
691  }
692  }
693 
694  // Now check for collisions with paraboloid surface.
695 
696  if(vRho2 == 0) // Needs to be treated seperately.
697  {
698  intersection = ((rho2 - k2)/k1 - p.z())/v.z();
699  if(calcNorm)
700  {
701  G4ThreeVector intersectionP = p + v * intersection;
702  *n = G4ThreeVector(intersectionP.x(), intersectionP.y(), -k1/2);
703  *n = n->unit();
704 
705  *validNorm = true;
706  }
707  return intersection;
708  }
709  else if( ((A <= 0) && (B >= sqr(A) * (sqr(vRho2) - 1))) || (A >= 0))
710  {
711  // intersection = (A + std::sqrt(B + sqr(A))) / vRho2;
712  // The above calculation has a precision problem:
713  // known problem of solving quadratic equation with small A
714 
715  A = A/vRho2;
716  B = (k1 * p.z() + k2 - rho2)/vRho2;
717  intersection = B/(-A + std::sqrt(B + sqr(A)));
718  if(calcNorm)
719  {
720  G4ThreeVector intersectionP = p + v * intersection;
721  *n = G4ThreeVector(intersectionP.x(), intersectionP.y(), -k1/2);
722  *n = n->unit();
723  *validNorm = true;
724  }
725  return intersection;
726  }
727  std::ostringstream message;
728  message << "There is no intersection between given line and solid!"
729  << G4endl
730  << " p = " << p << G4endl
731  << " v = " << v;
732  G4Exception("G4Paraboloid::DistanceToOut(p,v,...)", "GeomSolids1002",
733  JustWarning, message);
734 
735  return kInfinity;
736  }
737  else if ( (rho2 < paraRho2 + kCarTolerance
738  || sqr(rho2 - paraRho2 - 0.25 * tol2) < tol2 * paraRho2 )
739  && std::fabs(p.z()) < dz + tolh)
740  {
741  // If this is true we're somewhere in the border.
742 
743  G4ThreeVector normal = G4ThreeVector (p.x(), p.y(), -k1/2);
744 
745  if(std::fabs(p.z()) > dz - tolh)
746  {
747  // We're in the lower or upper edge
748  //
749  if( ((v.z() > 0) && (p.z() > 0)) || ((v.z() < 0) && (p.z() < 0)) )
750  { // If we're heading out of the object that is treated here
751  if(calcNorm)
752  {
753  *validNorm = true;
754  if(p.z() > 0)
755  { *n = G4ThreeVector(0, 0, 1); }
756  else
757  { *n = G4ThreeVector(0, 0, -1); }
758  }
759  return 0;
760  }
761 
762  if(v.z() == 0)
763  {
764  // Case where we're moving inside the surface needs to be
765  // treated separately.
766  // Distance until it goes out through a side is returned.
767 
768  G4double r = (p.z() > 0)? r2 : r1;
769  G4double pDotV = p.dot(v);
770  A = vRho2 * ( sqr(r) - sqr(p.x()) - sqr(p.y()));
771  intersection = (-pDotV + std::sqrt(A + sqr(pDotV))) / vRho2;
772 
773  if(calcNorm)
774  {
775  *validNorm = true;
776 
777  *n = (G4ThreeVector(0, 0, p.z()/std::fabs(p.z()))
778  + G4ThreeVector(p.x() + v.x() * intersection, p.y() + v.y()
779  * intersection, -k1/2).unit()).unit();
780  }
781  return intersection;
782  }
783  }
784  //
785  // Problem in the Logic :: Following condition for point on upper surface
786  // and Vz<0 will return 0 (Problem #1015), but
787  // it has to return intersection with parabolic
788  // surface or with lower plane surface (z = -dz)
789  // The logic has to be :: If not found intersection until now,
790  // do not exit but continue to search for possible intersection.
791  // Only for point situated on both borders (Z and parabolic)
792  // this condition has to be taken into account and done later
793  //
794  //
795  // else if(normal.dot(v) >= 0)
796  // {
797  // if(calcNorm)
798  // {
799  // *validNorm = true;
800  // *n = normal.unit();
801  // }
802  // return 0;
803  // }
804 
805  if(v.z() > 0)
806  {
807  // Check for collision with upper edge.
808 
809  intersection = (dz - p.z()) / v.z();
810  G4ThreeVector ip = p + intersection * v;
811 
812  if(ip.perp2() < sqr(r2 - tolh))
813  {
814  if(calcNorm)
815  {
816  *validNorm = true;
817  *n = G4ThreeVector(0, 0, 1);
818  }
819  return intersection;
820  }
821  else if(ip.perp2() < sqr(r2 + tolh))
822  {
823  if(calcNorm)
824  {
825  *validNorm = true;
826  *n = G4ThreeVector(0, 0, 1)
827  + G4ThreeVector(ip.x(), ip.y(), - k1 / 2).unit();
828  *n = n->unit();
829  }
830  return intersection;
831  }
832  }
833  if( v.z() < 0)
834  {
835  // Check for collision with lower edge.
836 
837  intersection = (-dz - p.z()) / v.z();
838  G4ThreeVector ip = p + intersection * v;
839 
840  if(ip.perp2() < sqr(r1 - tolh))
841  {
842  if(calcNorm)
843  {
844  *validNorm = true;
845  *n = G4ThreeVector(0, 0, -1);
846  }
847  return intersection;
848  }
849  else if(ip.perp2() < sqr(r1 + tolh))
850  {
851  if(calcNorm)
852  {
853  *validNorm = true;
854  *n = G4ThreeVector(0, 0, -1)
855  + G4ThreeVector(ip.x(), ip.y(), - k1 / 2).unit();
856  *n = n->unit();
857  }
858  return intersection;
859  }
860  }
861 
862  // Note: comparison with zero below would not be correct !
863  //
864  if(std::fabs(vRho2) > tol2) // precision error in the calculation of
865  { // intersection = (A+std::sqrt(B+sqr(A)))/vRho2
866  A = A/vRho2;
867  B = (k1 * p.z() + k2 - rho2);
868  if(std::fabs(B)>kCarTolerance)
869  {
870  B = (B)/vRho2;
871  intersection = B/(-A + std::sqrt(B + sqr(A)));
872  }
873  else // Point is On both borders: Z and parabolic
874  { // solution depends on normal.dot(v) sign
875  if(normal.dot(v) >= 0)
876  {
877  if(calcNorm)
878  {
879  *validNorm = true;
880  *n = normal.unit();
881  }
882  return 0;
883  }
884  intersection = 2.*A;
885  }
886  }
887  else
888  {
889  intersection = ((rho2 - k2) / k1 - p.z()) / v.z();
890  }
891 
892  if(calcNorm)
893  {
894  *validNorm = true;
895  *n = G4ThreeVector(p.x() + intersection * v.x(), p.y()
896  + intersection * v.y(), - k1 / 2);
897  *n = n->unit();
898  }
899  return intersection;
900  }
901  else
902  {
903 #ifdef G4SPECSDEBUG
904  if(kOutside == Inside(p))
905  {
906  G4Exception("G4Paraboloid::DistanceToOut(p,v,...)", "GeomSolids1002",
907  JustWarning, "Point p is outside!");
908  }
909  else
910  G4Exception("G4Paraboloid::DistanceToOut(p,v,...)", "GeomSolids1002",
911  JustWarning, "There's an error in this functions code.");
912 #endif
913  return kInfinity;
914  }
915  return 0;
916 }
917 
919 //
920 // Calculate distance (<=actual) to closest surface of shape from inside
921 
923 {
924  G4double safe=0.0,rho,safeR,safeZ ;
925  G4double tanRMax,secRMax,pRMax ;
926 
927 #ifdef G4SPECSDEBUG
928  if( Inside(p) == kOutside )
929  {
930  G4cout << G4endl ;
931  DumpInfo();
932  std::ostringstream message;
933  G4int oldprc = message.precision(16);
934  message << "Point p is outside !?" << G4endl
935  << "Position:" << G4endl
936  << " p.x() = " << p.x()/mm << " mm" << G4endl
937  << " p.y() = " << p.y()/mm << " mm" << G4endl
938  << " p.z() = " << p.z()/mm << " mm";
939  message.precision(oldprc) ;
940  G4Exception("G4Paraboloid::DistanceToOut(p)", "GeomSolids1002",
941  JustWarning, message);
942  }
943 #endif
944 
945  rho = p.perp();
946  safeZ = dz - std::fabs(p.z()) ;
947 
948  tanRMax = (r2 - r1)*0.5/dz ;
949  secRMax = std::sqrt(1.0 + tanRMax*tanRMax) ;
950  pRMax = tanRMax*p.z() + (r1+r2)*0.5 ;
951  safeR = (pRMax - rho)/secRMax ;
952 
953  if (safeZ < safeR) { safe = safeZ; }
954  else { safe = safeR; }
955  if ( safe < 0.5 * kCarTolerance ) { safe = 0; }
956  return safe ;
957 }
958 
960 //
961 // G4EntityType
962 
964 {
965  return G4String("G4Paraboloid");
966 }
967 
969 //
970 // Make a clone of the object
971 
973 {
974  return new G4Paraboloid(*this);
975 }
976 
978 //
979 // Stream object contents to an output stream
980 
981 std::ostream& G4Paraboloid::StreamInfo( std::ostream& os ) const
982 {
983  G4int oldprc = os.precision(16);
984  os << "-----------------------------------------------------------\n"
985  << " *** Dump for solid - " << GetName() << " ***\n"
986  << " ===================================================\n"
987  << " Solid type: G4Paraboloid\n"
988  << " Parameters: \n"
989  << " z half-axis: " << dz/mm << " mm \n"
990  << " radius at -dz: " << r1/mm << " mm \n"
991  << " radius at dz: " << r2/mm << " mm \n"
992  << "-----------------------------------------------------------\n";
993  os.precision(oldprc);
994 
995  return os;
996 }
997 
999 //
1000 // GetPointOnSurface
1001 
1003 {
1004  G4double A = (fSurfaceArea == 0)? CalculateSurfaceArea(): fSurfaceArea;
1005  G4double z = RandFlat::shoot(0.,1.);
1006  G4double phi = RandFlat::shoot(0., twopi);
1007  if(pi*(sqr(r1) + sqr(r2))/A >= z)
1008  {
1009  G4double rho;
1010  if(pi * sqr(r1) / A > z)
1011  {
1012  rho = r1 * std::sqrt(RandFlat::shoot(0., 1.));
1013  return G4ThreeVector(rho * std::cos(phi), rho * std::sin(phi), -dz);
1014  }
1015  else
1016  {
1017  rho = r2 * std::sqrt(RandFlat::shoot(0., 1));
1018  return G4ThreeVector(rho * std::cos(phi), rho * std::sin(phi), dz);
1019  }
1020  }
1021  else
1022  {
1023  z = RandFlat::shoot(0., 1.)*2*dz - dz;
1024  return G4ThreeVector(std::sqrt(z*k1 + k2)*std::cos(phi),
1025  std::sqrt(z*k1 + k2)*std::sin(phi), z);
1026  }
1027 }
1028 
1030 G4Paraboloid::CreateRotatedVertices(const G4AffineTransform& pTransform,
1031  G4int& noPolygonVertices) const
1032 {
1033  G4ThreeVectorList *vertices;
1034  G4ThreeVector vertex;
1035  G4double meshAnglePhi, cosMeshAnglePhiPer2,
1036  crossAnglePhi, coscrossAnglePhi, sincrossAnglePhi, sAnglePhi,
1037  sRho, dRho, rho, lastRho = 0., swapRho;
1038  G4double rx, ry, rz, k3, k4, zm;
1039  G4int crossSectionPhi, noPhiCrossSections, noRhoSections;
1040 
1041  // Phi cross sections
1042  //
1043  noPhiCrossSections = G4int(twopi/kMeshAngleDefault)+1; // =9!
1044 /*
1045  if (noPhiCrossSections<kMinMeshSections) // <3
1046  {
1047  noPhiCrossSections=kMinMeshSections;
1048  }
1049  else if (noPhiCrossSections>kMaxMeshSections) // >37
1050  {
1051  noPhiCrossSections=kMaxMeshSections;
1052  }
1053 */
1054  meshAnglePhi=twopi/(noPhiCrossSections-1);
1055 
1056  sAnglePhi = -meshAnglePhi*0.5*0;
1057  cosMeshAnglePhiPer2 = std::cos(meshAnglePhi / 2.);
1058 
1059  noRhoSections = G4int(pi/2/kMeshAngleDefault) + 1;
1060 
1061  // There is no obvious value for noRhoSections, at the moment the parabola is
1062  // viewed as a quarter circle mean this formula for it.
1063 
1064  // An alternetive would be to calculate max deviation from parabola and
1065  // keep adding new vertices there until it was under a decided constant.
1066 
1067  // maxDeviation on a line between points (rho1, z1) and (rho2, z2) is given
1068  // by rhoMax = sqrt(k1 * z + k2) - z * (rho2 - rho1)
1069  // / (z2 - z1) - (rho1 * z2 - rho2 * z1) / (z2 - z1)
1070  // where z is k1 / 2 * (rho1 + rho2) - k2 / k1
1071 
1072  sRho = r1;
1073  dRho = (r2 - r1) / double(noRhoSections - 1);
1074 
1075  vertices=new G4ThreeVectorList();
1076 
1077  if (vertices)
1078  {
1079  for (crossSectionPhi=0; crossSectionPhi<noPhiCrossSections;
1080  crossSectionPhi++)
1081  {
1082  crossAnglePhi=sAnglePhi+crossSectionPhi*meshAnglePhi;
1083  coscrossAnglePhi=std::cos(crossAnglePhi);
1084  sincrossAnglePhi=std::sin(crossAnglePhi);
1085  lastRho = 0;
1086  for (int iRho=0; iRho < noRhoSections;
1087  iRho++)
1088  {
1089  // Compute coordinates of cross section at section crossSectionPhi
1090  //
1091  if(iRho == noRhoSections - 1)
1092  {
1093  rho = r2;
1094  }
1095  else
1096  {
1097  rho = iRho * dRho + sRho;
1098 
1099  // This part is to ensure that the vertices
1100  // will form a volume larger than the paraboloid
1101 
1102  k3 = k1 / (2*rho + dRho);
1103  k4 = rho - k3 * (sqr(rho) - k2) / k1;
1104  zm = (sqr(k1 / (2 * k3)) - k2) / k1;
1105  rho += std::sqrt(k1 * zm + k2) - zm * k3 - k4;
1106  }
1107 
1108  rho += (1 / cosMeshAnglePhiPer2 - 1) * (iRho * dRho + sRho);
1109 
1110  if(rho < lastRho)
1111  {
1112  swapRho = lastRho;
1113  lastRho = rho + dRho;
1114  rho = swapRho;
1115  }
1116  else
1117  {
1118  lastRho = rho + dRho;
1119  }
1120 
1121  rx = coscrossAnglePhi*rho;
1122  ry = sincrossAnglePhi*rho;
1123  rz = (sqr(iRho * dRho + sRho) - k2) / k1;
1124  vertex = G4ThreeVector(rx,ry,rz);
1125  vertices->push_back(pTransform.TransformPoint(vertex));
1126  }
1127  } // Phi
1128  noPolygonVertices = noRhoSections ;
1129  }
1130  else
1131  {
1132  DumpInfo();
1133  G4Exception("G4Paraboloid::CreateRotatedVertices()",
1134  "GeomSolids0003", FatalException,
1135  "Error in allocation of vertices. Out of memory !");
1136  }
1137  return vertices;
1138 }
1139 
1141 //
1142 // Methods for visualisation
1143 
1145 {
1146  scene.AddSolid(*this);
1147 }
1148 
1150 {
1151  // Box for now!!!
1152  //
1153  return new G4NURBSbox(r1, r1, dz);
1154 }
1155 
1157 {
1158  return new G4PolyhedronParaboloid(r1, r2, dz, 0., twopi);
1159 }
1160 
1161 
1163 {
1164  if (!fpPolyhedron ||
1167  {
1168  delete fpPolyhedron;
1170  }
1171  return fpPolyhedron;
1172 }