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