Geant4  10.00.p02
G4Hype.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 //
27 // $Id: G4Hype.cc 81641 2014-06-04 09:11:38Z gcosmo $
28 // $Original: G4Hype.cc,v 1.0 1998/06/09 16:57:50 safai Exp $
29 //
30 //
31 // --------------------------------------------------------------------
32 // GEANT 4 class source file
33 //
34 //
35 // G4Hype.cc
36 //
37 // --------------------------------------------------------------------
38 //
39 // Authors:
40 // Ernesto Lamanna (Ernesto.Lamanna@roma1.infn.it) &
41 // Francesco Safai Tehrani (Francesco.SafaiTehrani@roma1.infn.it)
42 // Rome, INFN & University of Rome "La Sapienza", 9 June 1998.
43 //
44 // --------------------------------------------------------------------
45 
46 #include "G4Hype.hh"
47 
48 #include "G4VoxelLimits.hh"
49 #include "G4AffineTransform.hh"
50 #include "G4SolidExtentList.hh"
51 #include "G4ClippablePolygon.hh"
52 
53 #include "G4VPVParameterisation.hh"
54 
55 #include "meshdefs.hh"
56 
57 #include <cmath>
58 
59 #include "Randomize.hh"
60 
61 #include "G4VGraphicsScene.hh"
62 #include "G4VisExtent.hh"
63 
64 using namespace CLHEP;
65 
66 // Constructor - check parameters, and fills protected data members
67 G4Hype::G4Hype(const G4String& pName,
68  G4double newInnerRadius,
69  G4double newOuterRadius,
70  G4double newInnerStereo,
71  G4double newOuterStereo,
72  G4double newHalfLenZ)
73  : G4VSolid(pName), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
74 {
75  fHalfTol = 0.5*kCarTolerance;
76 
77  // Check z-len
78  //
79  if (newHalfLenZ<=0)
80  {
81  std::ostringstream message;
82  message << "Invalid Z half-length - " << GetName() << G4endl
83  << " Invalid Z half-length: "
84  << newHalfLenZ/mm << " mm";
85  G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
86  FatalErrorInArgument, message);
87  }
88  halfLenZ=newHalfLenZ;
89 
90  // Check radii
91  //
92  if (newInnerRadius<0 || newOuterRadius<0)
93  {
94  std::ostringstream message;
95  message << "Invalid radii - " << GetName() << G4endl
96  << " Invalid radii ! Inner radius: "
97  << newInnerRadius/mm << " mm" << G4endl
98  << " Outer radius: "
99  << newOuterRadius/mm << " mm";
100  G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
101  FatalErrorInArgument, message);
102  }
103  if (newInnerRadius >= newOuterRadius)
104  {
105  std::ostringstream message;
106  message << "Outer > inner radius - " << GetName() << G4endl
107  << " Invalid radii ! Inner radius: "
108  << newInnerRadius/mm << " mm" << G4endl
109  << " Outer radius: "
110  << newOuterRadius/mm << " mm";
111  G4Exception("G4Hype::G4Hype()", "GeomSolids0002",
112  FatalErrorInArgument, message);
113  }
114 
115  innerRadius=newInnerRadius;
116  outerRadius=newOuterRadius;
117 
120 
121  SetInnerStereo( newInnerStereo );
122  SetOuterStereo( newOuterStereo );
123 }
124 
125 
126 //
127 // Fake default constructor - sets only member data and allocates memory
128 // for usage restricted to object persistency.
129 //
130 G4Hype::G4Hype( __void__& a )
131  : G4VSolid(a), innerRadius(0.), outerRadius(0.), halfLenZ(0.), innerStereo(0.),
132  outerStereo(0.), tanInnerStereo(0.), tanOuterStereo(0.), tanInnerStereo2(0.),
133  tanOuterStereo2(0.), innerRadius2(0.), outerRadius2(0.), endInnerRadius2(0.),
134  endOuterRadius2(0.), endInnerRadius(0.), endOuterRadius(0.),
135  fCubicVolume(0.), fSurfaceArea(0.), fHalfTol(0.), fpPolyhedron(0)
136 {
137 }
138 
139 
140 //
141 // Destructor
142 //
144 {
145  delete fpPolyhedron;
146 }
147 
148 
149 //
150 // Copy constructor
151 //
153  : G4VSolid(rhs), innerRadius(rhs.innerRadius),
154  outerRadius(rhs.outerRadius), halfLenZ(rhs.halfLenZ),
155  innerStereo(rhs.innerStereo), outerStereo(rhs.outerStereo),
156  tanInnerStereo(rhs.tanInnerStereo), tanOuterStereo(rhs.tanOuterStereo),
157  tanInnerStereo2(rhs.tanInnerStereo2), tanOuterStereo2(rhs.tanOuterStereo2),
158  innerRadius2(rhs.innerRadius2), outerRadius2(rhs.outerRadius2),
159  endInnerRadius2(rhs.endInnerRadius2), endOuterRadius2(rhs.endOuterRadius2),
160  endInnerRadius(rhs.endInnerRadius), endOuterRadius(rhs.endOuterRadius),
161  fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
162  fHalfTol(rhs.fHalfTol), fpPolyhedron(0)
163 {
165 }
166 
167 
168 //
169 // Assignment operator
170 //
172 {
173  // Check assignment to self
174  //
175  if (this == &rhs) { return *this; }
176 
177  // Copy base class data
178  //
179  G4VSolid::operator=(rhs);
180 
181  // Copy data
182  //
184  halfLenZ = rhs.halfLenZ;
192  fHalfTol = rhs.fHalfTol;
194 
195  return *this;
196 }
197 
198 
199 //
200 // Dispatch to parameterisation for replication mechanism dimension
201 // computation & modification.
202 //
204  const G4int n,
205  const G4VPhysicalVolume* pRep)
206 {
207  p->ComputeDimensions(*this,n,pRep);
208 }
209 
210 
211 //
212 // CalculateExtent
213 //
215  const G4VoxelLimits &voxelLimit,
216  const G4AffineTransform &transform,
217  G4double &min, G4double &max ) const
218 {
219  G4SolidExtentList extentList( axis, voxelLimit );
220 
221  //
222  // Choose phi size of our segment(s) based on constants as
223  // defined in meshdefs.hh
224  //
225  G4int numPhi = kMaxMeshSections;
226  G4double sigPhi = twopi/numPhi;
227  G4double rFudge = 1.0/std::cos(0.5*sigPhi);
228 
229  //
230  // We work around in phi building polygons along the way.
231  // As a reasonable compromise between accuracy and
232  // complexity (=cpu time), the following facets are chosen:
233  //
234  // 1. If outerRadius/endOuterRadius > 0.95, approximate
235  // the outer surface as a cylinder, and use one
236  // rectangular polygon (0-1) to build its mesh.
237  //
238  // Otherwise, use two trapazoidal polygons that
239  // meet at z = 0 (0-4-1)
240  //
241  // 2. If there is no inner surface, then use one
242  // polygon for each entire endcap. (0) and (1)
243  //
244  // Otherwise, use a trapazoidal polygon for each
245  // phi segment of each endcap. (0-2) and (1-3)
246  //
247  // 3. For the inner surface, if innerRadius/endInnerRadius > 0.95,
248  // approximate the inner surface as a cylinder of
249  // radius innerRadius and use one rectangular polygon
250  // to build each phi segment of its mesh. (2-3)
251  //
252  // Otherwise, use one rectangular polygon centered
253  // at z = 0 (5-6) and two connecting trapazoidal polygons
254  // for each phi segment (2-5) and (3-6).
255  //
256 
257  G4bool splitOuter = (outerRadius/endOuterRadius < 0.95);
258  G4bool splitInner = 0;
259  if (InnerSurfaceExists())
260  {
261  splitInner = (innerRadius/endInnerRadius < 0.95);
262  }
263 
264  //
265  // Vertex assignments (v and w arrays)
266  // [0] and [1] are mandatory
267  // the rest are optional
268  //
269  // + -
270  // [0]------[4]------[1] <--- outer radius
271  // | |
272  // | |
273  // [2]---[5]---[6]---[3] <--- inner radius
274  //
275 
276 
277  G4ClippablePolygon endPoly1, endPoly2;
278 
279  G4double phi = 0,
280  cosPhi = std::cos(phi),
281  sinPhi = std::sin(phi);
282  G4ThreeVector v0( rFudge*endOuterRadius*cosPhi,
283  rFudge*endOuterRadius*sinPhi,
284  +halfLenZ ),
285  v1( rFudge*endOuterRadius*cosPhi,
286  rFudge*endOuterRadius*sinPhi,
287  -halfLenZ ),
288  v2, v3, v4, v5, v6,
289  w0, w1, w2, w3, w4, w5, w6;
290  transform.ApplyPointTransform( v0 );
291  transform.ApplyPointTransform( v1 );
292 
293  G4double zInnerSplit=0.;
294  if (InnerSurfaceExists())
295  {
296  if (splitInner)
297  {
298  v2 = transform.TransformPoint(
300  endInnerRadius*sinPhi, +halfLenZ ) );
301  v3 = transform.TransformPoint(
303  endInnerRadius*sinPhi, -halfLenZ ) );
304  //
305  // Find intersection of line normal to inner
306  // surface at z = halfLenZ and line r=innerRadius
307  //
310 
311  zInnerSplit = halfLenZ + (innerRadius - endInnerRadius)*zn/rn;
312 
313  //
314  // Build associated vertices
315  //
316  v5 = transform.TransformPoint(
317  G4ThreeVector( innerRadius*cosPhi,
318  innerRadius*sinPhi, +zInnerSplit ) );
319  v6 = transform.TransformPoint(
320  G4ThreeVector( innerRadius*cosPhi,
321  innerRadius*sinPhi, -zInnerSplit ) );
322  }
323  else
324  {
325  v2 = transform.TransformPoint(
326  G4ThreeVector( innerRadius*cosPhi,
327  innerRadius*sinPhi, +halfLenZ ) );
328  v3 = transform.TransformPoint(
329  G4ThreeVector( innerRadius*cosPhi,
330  innerRadius*sinPhi, -halfLenZ ) );
331  }
332  }
333 
334  if (splitOuter)
335  {
336  v4 = transform.TransformPoint(
337  G4ThreeVector( rFudge*outerRadius*cosPhi,
338  rFudge*outerRadius*sinPhi, 0 ) );
339  }
340 
341  //
342  // Loop over phi segments
343  //
344  do
345  {
346  phi += sigPhi;
347  if (numPhi == 1) phi = 0; // Try to avoid roundoff
348  cosPhi = std::cos(phi),
349  sinPhi = std::sin(phi);
350 
351  G4double r(rFudge*endOuterRadius);
352  w0 = G4ThreeVector( r*cosPhi, r*sinPhi, +halfLenZ );
353  w1 = G4ThreeVector( r*cosPhi, r*sinPhi, -halfLenZ );
354  transform.ApplyPointTransform( w0 );
355  transform.ApplyPointTransform( w1 );
356 
357  //
358  // Outer hyperbolic surface
359  //
360  if (splitOuter)
361  {
362  r = rFudge*outerRadius;
363  w4 = G4ThreeVector( r*cosPhi, r*sinPhi, 0 );
364  transform.ApplyPointTransform( w4 );
365 
366  AddPolyToExtent( v0, v4, w4, w0, voxelLimit, axis, extentList );
367  AddPolyToExtent( v4, v1, w1, w4, voxelLimit, axis, extentList );
368  }
369  else
370  {
371  AddPolyToExtent( v0, v1, w1, w0, voxelLimit, axis, extentList );
372  }
373 
374  if (InnerSurfaceExists())
375  {
376  //
377  // Inner hyperbolic surface
378  //
379  if (splitInner)
380  {
381  w2 = G4ThreeVector( endInnerRadius*cosPhi,
382  endInnerRadius*sinPhi, +halfLenZ );
383  w3 = G4ThreeVector( endInnerRadius*cosPhi,
384  endInnerRadius*sinPhi, -halfLenZ );
385  transform.ApplyPointTransform( w2 );
386  transform.ApplyPointTransform( w3 );
387 
388  w5 = G4ThreeVector( innerRadius*cosPhi,
389  innerRadius*sinPhi, +zInnerSplit );
390  w6 = G4ThreeVector( innerRadius*cosPhi,
391  innerRadius*sinPhi, -zInnerSplit );
392  transform.ApplyPointTransform( w5 );
393  transform.ApplyPointTransform( w6 );
394  AddPolyToExtent( v3, v6, w6, w3, voxelLimit, axis, extentList );
395  AddPolyToExtent( v6, v5, w5, w6, voxelLimit, axis, extentList );
396  AddPolyToExtent( v5, v2, w2, w5, voxelLimit, axis, extentList );
397  }
398  else
399  {
400  w2 = G4ThreeVector( innerRadius*cosPhi,
401  innerRadius*sinPhi, +halfLenZ );
402  w3 = G4ThreeVector( innerRadius*cosPhi,
403  innerRadius*sinPhi, -halfLenZ );
404  transform.ApplyPointTransform( w2 );
405  transform.ApplyPointTransform( w3 );
406 
407  AddPolyToExtent( v3, v2, w2, w3, voxelLimit, axis, extentList );
408  }
409 
410  //
411  // Endplate segments
412  //
413  AddPolyToExtent( v1, v3, w3, w1, voxelLimit, axis, extentList );
414  AddPolyToExtent( v2, v0, w0, w2, voxelLimit, axis, extentList );
415  }
416  else
417  {
418  //
419  // Continue building endplate polygons
420  //
421  endPoly1.AddVertexInOrder( v0 );
422  endPoly2.AddVertexInOrder( v1 );
423  }
424 
425  //
426  // Next phi segments
427  //
428  v0 = w0;
429  v1 = w1;
430  if (InnerSurfaceExists())
431  {
432  v2 = w2;
433  v3 = w3;
434  if (splitInner)
435  {
436  v5 = w5;
437  v6 = w6;
438  }
439  }
440  if (splitOuter) v4 = w4;
441 
442  } while( --numPhi > 0 );
443 
444 
445  //
446  // Don't forget about the endplate polygons, if
447  // we use them
448  //
449  if (!InnerSurfaceExists())
450  {
451  if (endPoly1.PartialClip( voxelLimit, axis ))
452  {
453  static const G4ThreeVector normal(0,0,+1);
454  endPoly1.SetNormal( transform.TransformAxis(normal) );
455  extentList.AddSurface( endPoly1 );
456  }
457 
458  if (endPoly2.PartialClip( voxelLimit, axis ))
459  {
460  static const G4ThreeVector normal(0,0,-1);
461  endPoly2.SetNormal( transform.TransformAxis(normal) );
462  extentList.AddSurface( endPoly2 );
463  }
464  }
465 
466  //
467  // Return min/max value
468  //
469  return extentList.GetExtent( min, max );
470 }
471 
472 
473 //
474 // AddPolyToExtent (static)
475 //
476 // Utility function for CalculateExtent
477 //
479  const G4ThreeVector &v1,
480  const G4ThreeVector &w1,
481  const G4ThreeVector &w0,
482  const G4VoxelLimits &voxelLimit,
483  const EAxis axis,
484  G4SolidExtentList &extentList )
485 {
486  G4ClippablePolygon phiPoly;
487 
488  phiPoly.AddVertexInOrder( v0 );
489  phiPoly.AddVertexInOrder( v1 );
490  phiPoly.AddVertexInOrder( w1 );
491  phiPoly.AddVertexInOrder( w0 );
492 
493  if (phiPoly.PartialClip( voxelLimit, axis ))
494  {
495  phiPoly.SetNormal( (v1-v0).cross(w0-v0).unit() );
496  extentList.AddSurface( phiPoly );
497  }
498 }
499 
500 
501 //
502 // Decides whether point is inside,outside or on the surface
503 //
505 {
506  //
507  // Check z extents: are we outside?
508  //
509  const G4double absZ(std::fabs(p.z()));
510  if (absZ > halfLenZ + fHalfTol) return kOutside;
511 
512  //
513  // Check outer radius
514  //
515  const G4double oRad2(HypeOuterRadius2(absZ));
516  const G4double xR2( p.x()*p.x()+p.y()*p.y() );
517 
518  if (xR2 > oRad2 + kCarTolerance*endOuterRadius) return kOutside;
519 
520  if (xR2 > oRad2 - kCarTolerance*endOuterRadius) return kSurface;
521 
522  if (InnerSurfaceExists())
523  {
524  //
525  // Check inner radius
526  //
527  const G4double iRad2(HypeInnerRadius2(absZ));
528 
529  if (xR2 < iRad2 - kCarTolerance*endInnerRadius) return kOutside;
530 
531  if (xR2 < iRad2 + kCarTolerance*endInnerRadius) return kSurface;
532  }
533 
534  //
535  // We are inside in radius, now check endplate surface
536  //
537  if (absZ > halfLenZ - fHalfTol) return kSurface;
538 
539  return kInside;
540 }
541 
542 
543 
544 //
545 // return the normal unit vector to the Hyperbolical Surface at a point
546 // p on (or nearly on) the surface
547 //
549 {
550  //
551  // Which of the three or four surfaces are we closest to?
552  //
553  const G4double absZ(std::fabs(p.z()));
554  const G4double distZ(absZ - halfLenZ);
555  const G4double dist2Z(distZ*distZ);
556 
557  const G4double xR2( p.x()*p.x()+p.y()*p.y() );
558  const G4double dist2Outer( std::fabs(xR2 - HypeOuterRadius2(absZ)) );
559 
560  if (InnerSurfaceExists())
561  {
562  //
563  // Has inner surface: is this closest?
564  //
565  const G4double dist2Inner( std::fabs(xR2 - HypeInnerRadius2(absZ)) );
566  if (dist2Inner < dist2Z && dist2Inner < dist2Outer)
567  return G4ThreeVector( -p.x(), -p.y(), p.z()*tanInnerStereo2 ).unit();
568  }
569 
570  //
571  // Do the "endcaps" win?
572  //
573  if (dist2Z < dist2Outer)
574  return G4ThreeVector( 0.0, 0.0, p.z() < 0 ? -1.0 : 1.0 );
575 
576 
577  //
578  // Outer surface wins
579  //
580  return G4ThreeVector( p.x(), p.y(), -p.z()*tanOuterStereo2 ).unit();
581 }
582 
583 
584 //
585 // Calculate distance to shape from outside, along normalised vector
586 // - return kInfinity if no intersection,
587 // or intersection distance <= tolerance
588 //
589 // Calculating the intersection of a line with the surfaces
590 // is fairly straight forward. The difficult problem is dealing
591 // with the intersections of the surfaces in a consistent manner,
592 // and this accounts for the complicated logic.
593 //
595  const G4ThreeVector& v ) const
596 {
597  //
598  // Quick test. Beware! This assumes v is a unit vector!
599  //
600  if (std::fabs(p.x()*v.y() - p.y()*v.x()) > endOuterRadius+kCarTolerance)
601  return kInfinity;
602 
603  //
604  // Take advantage of z symmetry, and reflect throught the
605  // z=0 plane so that pz is always positive
606  //
607  G4double pz(p.z()), vz(v.z());
608  if (pz < 0)
609  {
610  pz = -pz;
611  vz = -vz;
612  }
613 
614  //
615  // We must be very careful if we don't want to
616  // create subtle leaks at the edges where the
617  // hyperbolic surfaces connect to the endplate.
618  // The only reliable way to do so is to make sure
619  // that the decision as to when a track passes
620  // over the edge of one surface is exactly the
621  // same decision as to when a track passes into the
622  // other surface. By "exact", we don't mean algebraicly
623  // exact, but we mean the same machine instructions
624  // should be used.
625  //
626  G4bool couldMissOuter(true),
627  couldMissInner(true),
628  cantMissInnerCylinder(false);
629 
630  //
631  // Check endplate intersection
632  //
633  G4double sigz = pz-halfLenZ;
634 
635  if (sigz > -fHalfTol) // equivalent to: if (pz > halfLenZ - fHalfTol)
636  {
637  //
638  // We start in front of the endplate (within roundoff)
639  // Correct direction to intersect endplate?
640  //
641  if (vz >= 0)
642  {
643  //
644  // Nope. As long as we are far enough away, we
645  // can't intersect anything
646  //
647  if (sigz > 0) return kInfinity;
648 
649  //
650  // Otherwise, we may still hit a hyperbolic surface
651  // if the point is on the hyperbolic surface (within tolerance)
652  //
653  G4double pr2 = p.x()*p.x() + p.y()*p.y();
655  return kInfinity;
656 
657  if (InnerSurfaceExists())
658  {
660  return kInfinity;
661  if ( (pr2 < endOuterRadius2 - kCarTolerance*endOuterRadius)
662  && (pr2 > endInnerRadius2 + kCarTolerance*endInnerRadius) )
663  return kInfinity;
664  }
665  else
666  {
667  if (pr2 < endOuterRadius2 - kCarTolerance*endOuterRadius)
668  return kInfinity;
669  }
670  }
671  else
672  {
673  //
674  // Where do we intersect at z = halfLenZ?
675  //
676  G4double q = -sigz/vz;
677  G4double xi = p.x() + q*v.x(),
678  yi = p.y() + q*v.y();
679 
680  //
681  // Is this on the endplate? If so, return s, unless
682  // we are on the tolerant surface, in which case return 0
683  //
684  G4double pr2 = xi*xi + yi*yi;
685  if (pr2 <= endOuterRadius2)
686  {
687  if (InnerSurfaceExists())
688  {
689  if (pr2 >= endInnerRadius2) return (sigz < fHalfTol) ? 0 : q;
690  //
691  // This test is sufficient to ensure that the
692  // trajectory cannot miss the inner hyperbolic surface
693  // for z > 0, if the normal is correct.
694  //
695  G4double dot1 = (xi*v.x() + yi*v.y())*endInnerRadius/std::sqrt(pr2);
696  couldMissInner = (dot1 - halfLenZ*tanInnerStereo2*vz <= 0);
697 
698  if (pr2 > endInnerRadius2*(1 - 2*DBL_EPSILON) )
699  {
700  //
701  // There is a potential leak if the inner
702  // surface is a cylinder
703  //
704  if ( (innerStereo < DBL_MIN)
705  && ((std::fabs(v.x()) > DBL_MIN) || (std::fabs(v.y()) > DBL_MIN)) )
706  cantMissInnerCylinder = true;
707  }
708  }
709  else
710  {
711  return (sigz < fHalfTol) ? 0 : q;
712  }
713  }
714  else
715  {
716  G4double dotR( xi*v.x() + yi*v.y() );
717  if (dotR >= 0)
718  {
719  //
720  // Otherwise, if we are traveling outwards, we know
721  // we must miss the hyperbolic surfaces also, so
722  // we need not bother checking
723  //
724  return kInfinity;
725  }
726  else
727  {
728  //
729  // This test is sufficient to ensure that the
730  // trajectory cannot miss the outer hyperbolic surface
731  // for z > 0, if the normal is correct.
732  //
733  G4double dot1 = dotR*endOuterRadius/std::sqrt(pr2);
734  couldMissOuter = (dot1 - halfLenZ*tanOuterStereo2*vz>= 0);
735  }
736  }
737  }
738  }
739 
740  //
741  // Check intersection with outer hyperbolic surface, save
742  // distance to valid intersection into "best".
743  //
744  G4double best = kInfinity;
745 
746  G4double q[2];
748 
749  if (n > 0)
750  {
751  //
752  // Potential intersection: is p on this surface?
753  //
754  if (pz < halfLenZ+fHalfTol)
755  {
756  G4double dr2 = p.x()*p.x() + p.y()*p.y() - HypeOuterRadius2(pz);
757  if (std::fabs(dr2) < kCarTolerance*endOuterRadius)
758  {
759  //
760  // Sure, but make sure we're traveling inwards at
761  // this point
762  //
763  if (p.x()*v.x() + p.y()*v.y() - pz*tanOuterStereo2*vz < 0)
764  return 0;
765  }
766  }
767 
768  //
769  // We are now certain that p is not on the tolerant surface.
770  // Accept only position distance q
771  //
772  G4int i;
773  for( i=0; i<n; i++ )
774  {
775  if (q[i] >= 0)
776  {
777  //
778  // Check to make sure this intersection point is
779  // on the surface, but only do so if we haven't
780  // checked the endplate intersection already
781  //
782  G4double zi = pz + q[i]*vz;
783 
784  if (zi < -halfLenZ) continue;
785  if (zi > +halfLenZ && couldMissOuter) continue;
786 
787  //
788  // Check normal
789  //
790  G4double xi = p.x() + q[i]*v.x(),
791  yi = p.y() + q[i]*v.y();
792 
793  if (xi*v.x() + yi*v.y() - zi*tanOuterStereo2*vz > 0) continue;
794 
795  best = q[i];
796  break;
797  }
798  }
799  }
800 
801  if (!InnerSurfaceExists()) return best;
802 
803  //
804  // Check intersection with inner hyperbolic surface
805  //
806  n = IntersectHype( p, v, innerRadius2, tanInnerStereo2, q );
807  if (n == 0)
808  {
809  if (cantMissInnerCylinder) return (sigz < fHalfTol) ? 0 : -sigz/vz;
810 
811  return best;
812  }
813 
814  //
815  // P on this surface?
816  //
817  if (pz < halfLenZ+fHalfTol)
818  {
819  G4double dr2 = p.x()*p.x() + p.y()*p.y() - HypeInnerRadius2(pz);
820  if (std::fabs(dr2) < kCarTolerance*endInnerRadius)
821  {
822  //
823  // Sure, but make sure we're traveling outwards at
824  // this point
825  //
826  if (p.x()*v.x() + p.y()*v.y() - pz*tanInnerStereo2*vz > 0) return 0;
827  }
828  }
829 
830  //
831  // No, so only positive q is valid. Search for a valid intersection
832  // that is closer than the outer intersection (if it exists)
833  //
834  G4int i;
835  for( i=0; i<n; i++ )
836  {
837  if (q[i] > best) break;
838  if (q[i] >= 0)
839  {
840  //
841  // Check to make sure this intersection point is
842  // on the surface, but only do so if we haven't
843  // checked the endplate intersection already
844  //
845  G4double zi = pz + q[i]*vz;
846 
847  if (zi < -halfLenZ) continue;
848  if (zi > +halfLenZ && couldMissInner) continue;
849 
850  //
851  // Check normal
852  //
853  G4double xi = p.x() + q[i]*v.x(),
854  yi = p.y() + q[i]*v.y();
855 
856  if (xi*v.x() + yi*v.y() - zi*tanOuterStereo2*vz < 0) continue;
857 
858  best = q[i];
859  break;
860  }
861  }
862 
863  //
864  // Done
865  //
866  return best;
867 }
868 
869 
870 //
871 // Calculate distance to shape from outside, along perpendicular direction
872 // (if one exists). May be an underestimate.
873 //
874 // There are five (r,z) regions:
875 // 1. a point that is beyond the endcap but within the
876 // endcap radii
877 // 2. a point with r > outer endcap radius and with
878 // a z position that is beyond the cone formed by the
879 // normal of the outer hyperbolic surface at the
880 // edge at which it meets the endcap.
881 // 3. a point that is outside the outer surface and not in (1 or 2)
882 // 4. a point that is inside the inner surface and not in (5)
883 // 5. a point with radius < inner endcap radius and
884 // with a z position beyond the cone formed by the
885 // normal of the inner hyperbolic surface at the
886 // edge at which it meets the endcap.
887 // (regions 4 and 5 only exist if there is an inner surface)
888 //
890 {
891  G4double absZ(std::fabs(p.z()));
892 
893  //
894  // Check region
895  //
896  G4double r2 = p.x()*p.x() + p.y()*p.y();
897  G4double r = std::sqrt(r2);
898 
899  G4double sigz = absZ - halfLenZ;
900 
901  if (r < endOuterRadius)
902  {
903  if (sigz > -fHalfTol)
904  {
905  if (InnerSurfaceExists())
906  {
907  if (r > endInnerRadius)
908  return sigz < fHalfTol ? 0 : sigz; // Region 1
909 
910  G4double dr = endInnerRadius - r;
911  if (sigz > dr*tanInnerStereo2)
912  {
913  //
914  // In region 5
915  //
916  G4double answer = std::sqrt( dr*dr + sigz*sigz );
917  return answer < fHalfTol ? 0 : answer;
918  }
919  }
920  else
921  {
922  //
923  // In region 1 (no inner surface)
924  //
925  return sigz < fHalfTol ? 0 : sigz;
926  }
927  }
928  }
929  else
930  {
931  G4double dr = r - endOuterRadius;
932  if (sigz > -dr*tanOuterStereo2)
933  {
934  //
935  // In region 2
936  //
937  G4double answer = std::sqrt( dr*dr + sigz*sigz );
938  return answer < fHalfTol ? 0 : answer;
939  }
940  }
941 
942  if (InnerSurfaceExists())
943  {
945  {
946  //
947  // In region 4
948  //
950  return answer < fHalfTol ? 0 : answer;
951  }
952  }
953 
954  //
955  // We are left by elimination with region 3
956  //
958  return answer < fHalfTol ? 0 : answer;
959 }
960 
961 
962 //
963 // Calculate distance to surface of shape from `inside', allowing for tolerance
964 //
965 // The situation here is much simplier than DistanceToIn(p,v). For
966 // example, there is no need to even check whether an intersection
967 // point is inside the boundary of a surface, as long as all surfaces
968 // are checked and the smallest distance is used.
969 //
971  const G4bool calcNorm,
972  G4bool *validNorm, G4ThreeVector *norm ) const
973 {
974  static const G4ThreeVector normEnd1(0.0,0.0,+1.0);
975  static const G4ThreeVector normEnd2(0.0,0.0,-1.0);
976 
977  //
978  // Keep track of closest surface
979  //
980  G4double sBest; // distance to
981  const G4ThreeVector *nBest; // normal vector
982  G4bool vBest; // whether "valid"
983 
984  //
985  // Check endplate, taking advantage of symmetry.
986  // Note that the endcap is the only surface which
987  // has a "valid" normal, i.e. is a surface of which
988  // the entire solid is behind.
989  //
990  G4double pz(p.z()), vz(v.z());
991  if (vz < 0)
992  {
993  pz = -pz;
994  vz = -vz;
995  nBest = &normEnd2;
996  }
997  else
998  nBest = &normEnd1;
999 
1000  //
1001  // Possible intercept. Are we on the surface?
1002  //
1003  if (pz > halfLenZ-fHalfTol)
1004  {
1005  if (calcNorm) { *norm = *nBest; *validNorm = true; }
1006  return 0;
1007  }
1008 
1009  //
1010  // Nope. Get distance. Beware of zero vz.
1011  //
1012  sBest = (vz > DBL_MIN) ? (halfLenZ - pz)/vz : kInfinity;
1013  vBest = true;
1014 
1015  //
1016  // Check outer surface
1017  //
1018  G4double r2 = p.x()*p.x() + p.y()*p.y();
1019 
1020  G4double q[2];
1022 
1023  G4ThreeVector norm1, norm2;
1024 
1025  if (n > 0)
1026  {
1027  //
1028  // We hit somewhere. Are we on the surface?
1029  //
1030  G4double dr2 = r2 - HypeOuterRadius2(pz);
1031  if (std::fabs(dr2) < endOuterRadius*kCarTolerance)
1032  {
1033  G4ThreeVector normHere( p.x(), p.y(), -p.z()*tanOuterStereo2 );
1034  //
1035  // Sure. But are we going the right way?
1036  //
1037  if (normHere.dot(v) > 0)
1038  {
1039  if (calcNorm) { *norm = normHere.unit(); *validNorm = false; }
1040  return 0;
1041  }
1042  }
1043 
1044  //
1045  // Nope. Check closest positive intercept.
1046  //
1047  G4int i;
1048  for( i=0; i<n; i++ )
1049  {
1050  if (q[i] > sBest) break;
1051  if (q[i] > 0)
1052  {
1053  //
1054  // Make sure normal is correct (that this
1055  // solution is an outgoing solution)
1056  //
1057  G4ThreeVector pk(p+q[i]*v);
1058  norm1 = G4ThreeVector( pk.x(), pk.y(), -pk.z()*tanOuterStereo2 );
1059  if (norm1.dot(v) > 0)
1060  {
1061  sBest = q[i];
1062  nBest = &norm1;
1063  vBest = false;
1064  break;
1065  }
1066  }
1067  }
1068  }
1069 
1070  if (InnerSurfaceExists())
1071  {
1072  //
1073  // Check inner surface
1074  //
1075  n = IntersectHype( p, v, innerRadius2, tanInnerStereo2, q );
1076  if (n > 0)
1077  {
1078  //
1079  // On surface?
1080  //
1081  G4double dr2 = r2 - HypeInnerRadius2(pz);
1082  if (std::fabs(dr2) < endInnerRadius*kCarTolerance)
1083  {
1084  G4ThreeVector normHere( -p.x(), -p.y(), p.z()*tanInnerStereo2 );
1085  if (normHere.dot(v) > 0)
1086  {
1087  if (calcNorm)
1088  {
1089  *norm = normHere.unit();
1090  *validNorm = false;
1091  }
1092  return 0;
1093  }
1094  }
1095 
1096  //
1097  // Check closest positive
1098  //
1099  G4int i;
1100  for( i=0; i<n; i++ )
1101  {
1102  if (q[i] > sBest) break;
1103  if (q[i] > 0)
1104  {
1105  G4ThreeVector pk(p+q[i]*v);
1106  norm2 = G4ThreeVector( -pk.x(), -pk.y(), pk.z()*tanInnerStereo2 );
1107  if (norm2.dot(v) > 0)
1108  {
1109  sBest = q[i];
1110  nBest = &norm2;
1111  vBest = false;
1112  break;
1113  }
1114  }
1115  }
1116  }
1117  }
1118 
1119  //
1120  // Done!
1121  //
1122  if (calcNorm)
1123  {
1124  *validNorm = vBest;
1125 
1126  if (nBest == &norm1 || nBest == &norm2)
1127  *norm = nBest->unit();
1128  else
1129  *norm = *nBest;
1130  }
1131 
1132  return sBest;
1133 }
1134 
1135 
1136 //
1137 // Calculate distance (<=actual) to closest surface of shape from inside
1138 //
1139 // May be an underestimate
1140 //
1142 {
1143  //
1144  // Try each surface and remember the closest
1145  //
1146  G4double absZ(std::fabs(p.z()));
1147  G4double r(p.perp());
1148 
1149  G4double sBest = halfLenZ - absZ;
1150 
1151  G4double tryOuter = ApproxDistInside( r, absZ, outerRadius, tanOuterStereo2 );
1152  if (tryOuter < sBest)
1153  sBest = tryOuter;
1154 
1155  if (InnerSurfaceExists())
1156  {
1158  if (tryInner < sBest) sBest = tryInner;
1159  }
1160 
1161  return sBest < 0.5*kCarTolerance ? 0 : sBest;
1162 }
1163 
1164 
1165 //
1166 // IntersectHype (static)
1167 //
1168 // Decide if and where a line intersects with a hyperbolic
1169 // surface (of infinite extent)
1170 //
1171 // Arguments:
1172 // p - (in) Point on trajectory
1173 // v - (in) Vector along trajectory
1174 // r2 - (in) Square of radius at z = 0
1175 // tan2phi - (in) std::tan(phi)**2
1176 // q - (out) Up to two points of intersection, where the
1177 // intersection point is p + q*v, and if there are
1178 // two intersections, q[0] < q[1]. May be negative.
1179 // Returns:
1180 // The number of intersections. If 0, the trajectory misses.
1181 //
1182 //
1183 // Equation of a line:
1184 //
1185 // x = x0 + q*tx y = y0 + q*ty z = z0 + q*tz
1186 //
1187 // Equation of a hyperbolic surface:
1188 //
1189 // x**2 + y**2 = r**2 + (z*tanPhi)**2
1190 //
1191 // Solution is quadratic:
1192 //
1193 // a*q**2 + b*q + c = 0
1194 //
1195 // where:
1196 //
1197 // a = tx**2 + ty**2 - (tz*tanPhi)**2
1198 //
1199 // b = 2*( x0*tx + y0*ty - z0*tz*tanPhi**2 )
1200 //
1201 // c = x0**2 + y0**2 - r**2 - (z0*tanPhi)**2
1202 //
1203 //
1205  G4double r2, G4double tan2Phi, G4double ss[2] )
1206 {
1207  G4double x0 = p.x(), y0 = p.y(), z0 = p.z();
1208  G4double tx = v.x(), ty = v.y(), tz = v.z();
1209 
1210  G4double a = tx*tx + ty*ty - tz*tz*tan2Phi;
1211  G4double b = 2*( x0*tx + y0*ty - z0*tz*tan2Phi );
1212  G4double c = x0*x0 + y0*y0 - r2 - z0*z0*tan2Phi;
1213 
1214  if (std::fabs(a) < DBL_MIN)
1215  {
1216  //
1217  // The trajectory is parallel to the asympotic limit of
1218  // the surface: single solution
1219  //
1220  if (std::fabs(b) < DBL_MIN) return 0; // Unless we travel through exact center
1221 
1222  ss[0] = c/b;
1223  return 1;
1224  }
1225 
1226 
1227  G4double radical = b*b - 4*a*c;
1228 
1229  if (radical < -DBL_MIN) return 0; // No solution
1230 
1231  if (radical < DBL_MIN)
1232  {
1233  //
1234  // Grazes surface
1235  //
1236  ss[0] = -b/a/2.0;
1237  return 1;
1238  }
1239 
1240  radical = std::sqrt(radical);
1241 
1242  G4double q = -0.5*( b + (b < 0 ? -radical : +radical) );
1243  G4double sa = q/a;
1244  G4double sb = c/q;
1245  if (sa < sb) { ss[0] = sa; ss[1] = sb; } else { ss[0] = sb; ss[1] = sa; }
1246  return 2;
1247 }
1248 
1249 
1250 //
1251 // ApproxDistOutside (static)
1252 //
1253 // Find the approximate distance of a point outside
1254 // (greater radius) of a hyperbolic surface. The distance
1255 // must be an underestimate. It will also be nice (although
1256 // not necesary) that the estimate is always finite no
1257 // matter how close the point is.
1258 //
1259 // Our hyperbola approaches the asymptotic limit at z = +/- infinity
1260 // to the lines r = z*tanPhi. We call these lines the
1261 // asymptotic limit line.
1262 //
1263 // We need the distance of the 2d point p(r,z) to the
1264 // hyperbola r**2 = r0**2 + (z*tanPhi)**2. Find two
1265 // points that bracket the true normal and use the
1266 // distance to the line that connects these two points.
1267 // The first such point is z=p.z. The second point is
1268 // the z position on the asymptotic limit line that
1269 // contains the normal on the line through the point p.
1270 //
1272  G4double r0, G4double tanPhi )
1273 {
1274  if (tanPhi < DBL_MIN) return pr-r0;
1275 
1276  G4double tan2Phi = tanPhi*tanPhi;
1277 
1278  //
1279  // First point
1280  //
1281  G4double z1 = pz;
1282  G4double r1 = std::sqrt( r0*r0 + z1*z1*tan2Phi );
1283 
1284  //
1285  // Second point
1286  //
1287  G4double z2 = (pr*tanPhi + pz)/(1 + tan2Phi);
1288  G4double r2 = std::sqrt( r0*r0 + z2*z2*tan2Phi );
1289 
1290  //
1291  // Line between them
1292  //
1293  G4double dr = r2-r1;
1294  G4double dz = z2-z1;
1295 
1296  G4double len = std::sqrt(dr*dr + dz*dz);
1297  if (len < DBL_MIN)
1298  {
1299  //
1300  // The two points are the same?? I guess we
1301  // must have really bracketed the normal
1302  //
1303  dr = pr-r1;
1304  dz = pz-z1;
1305  return std::sqrt( dr*dr + dz*dz );
1306  }
1307 
1308  //
1309  // Distance
1310  //
1311  return std::fabs((pr-r1)*dz - (pz-z1)*dr)/len;
1312 }
1313 
1314 //
1315 // ApproxDistInside (static)
1316 //
1317 // Find the approximate distance of a point inside
1318 // of a hyperbolic surface. The distance
1319 // must be an underestimate. It will also be nice (although
1320 // not necesary) that the estimate is always finite no
1321 // matter how close the point is.
1322 //
1323 // This estimate uses the distance to a line tangent to
1324 // the hyperbolic function. The point of tangent is chosen
1325 // by the z position point
1326 //
1327 // Assumes pr and pz are positive
1328 //
1330  G4double r0, G4double tan2Phi )
1331 {
1332  if (tan2Phi < DBL_MIN) return r0 - pr;
1333 
1334  //
1335  // Corresponding position and normal on hyperbolic
1336  //
1337  G4double rh = std::sqrt( r0*r0 + pz*pz*tan2Phi );
1338 
1339  G4double dr = -rh;
1340  G4double dz = pz*tan2Phi;
1341  G4double len = std::sqrt(dr*dr + dz*dz);
1342 
1343  //
1344  // Answer
1345  //
1346  return std::fabs((pr-rh)*dr)/len;
1347 }
1348 
1349 
1350 //
1351 // GetEntityType
1352 //
1354 {
1355  return G4String("G4Hype");
1356 }
1357 
1358 
1359 //
1360 // Clone
1361 //
1363 {
1364  return new G4Hype(*this);
1365 }
1366 
1367 
1368 //
1369 // GetCubicVolume
1370 //
1372 {
1373  if(fCubicVolume != 0.) {;}
1375  return fCubicVolume;
1376 }
1377 
1378 
1379 //
1380 // GetSurfaceArea
1381 //
1383 {
1384  if(fSurfaceArea != 0.) {;}
1386  return fSurfaceArea;
1387 }
1388 
1389 
1390 //
1391 // Stream object contents to an output stream
1392 //
1393 std::ostream& G4Hype::StreamInfo(std::ostream& os) const
1394 {
1395  G4int oldprc = os.precision(16);
1396  os << "-----------------------------------------------------------\n"
1397  << " *** Dump for solid - " << GetName() << " ***\n"
1398  << " ===================================================\n"
1399  << " Solid type: G4Hype\n"
1400  << " Parameters: \n"
1401  << " half length Z: " << halfLenZ/mm << " mm \n"
1402  << " inner radius : " << innerRadius/mm << " mm \n"
1403  << " outer radius : " << outerRadius/mm << " mm \n"
1404  << " inner stereo angle : " << innerStereo/degree << " degrees \n"
1405  << " outer stereo angle : " << outerStereo/degree << " degrees \n"
1406  << "-----------------------------------------------------------\n";
1407  os.precision(oldprc);
1408 
1409  return os;
1410 }
1411 
1412 
1413 
1414 //
1415 // GetPointOnSurface
1416 //
1418 {
1419  G4double xRand, yRand, zRand, r2 , aOne, aTwo, aThree, chose, sinhu;
1420  G4double phi, cosphi, sinphi, rBar2Out, rBar2In, alpha, t, rOut, rIn2, rOut2;
1421 
1422  // we use the formula of the area of a surface of revolution to compute
1423  // the areas, using the equation of the hyperbola:
1424  // x^2 + y^2 = (z*tanphi)^2 + r^2
1425 
1426  rBar2Out = outerRadius2;
1427  alpha = 2.*pi*rBar2Out*std::cos(outerStereo)/tanOuterStereo;
1429  t = std::log(t+std::sqrt(sqr(t)+1));
1430  aOne = std::fabs(2.*alpha*(std::sinh(2.*t)/4.+t/2.));
1431 
1432 
1433  rBar2In = innerRadius2;
1434  alpha = 2.*pi*rBar2In*std::cos(innerStereo)/tanInnerStereo;
1436  t = std::log(t+std::sqrt(sqr(t)+1));
1437  aTwo = std::fabs(2.*alpha*(std::sinh(2.*t)/4.+t/2.));
1438 
1441 
1442  if(outerStereo == 0.) {aOne = std::fabs(2.*pi*outerRadius*2.*halfLenZ);}
1443  if(innerStereo == 0.) {aTwo = std::fabs(2.*pi*innerRadius*2.*halfLenZ);}
1444 
1445  phi = RandFlat::shoot(0.,2.*pi);
1446  cosphi = std::cos(phi);
1447  sinphi = std::sin(phi);
1448  sinhu = RandFlat::shoot(-1.*halfLenZ*tanOuterStereo/outerRadius,
1449  halfLenZ*tanOuterStereo/outerRadius);
1450 
1451  chose = RandFlat::shoot(0.,aOne+aTwo+2.*aThree);
1452  if(chose>=0. && chose < aOne)
1453  {
1454  if(outerStereo != 0.)
1455  {
1456  zRand = outerRadius*sinhu/tanOuterStereo;
1457  xRand = std::sqrt(sqr(sinhu)+1)*outerRadius*cosphi;
1458  yRand = std::sqrt(sqr(sinhu)+1)*outerRadius*sinphi;
1459  return G4ThreeVector (xRand, yRand, zRand);
1460  }
1461  else
1462  {
1463  return G4ThreeVector(outerRadius*cosphi,outerRadius*sinphi,
1465  }
1466  }
1467  else if(chose>=aOne && chose<aOne+aTwo)
1468  {
1469  if(innerStereo != 0.)
1470  {
1471  sinhu = RandFlat::shoot(-1.*halfLenZ*tanInnerStereo/innerRadius,
1472  halfLenZ*tanInnerStereo/innerRadius);
1473  zRand = innerRadius*sinhu/tanInnerStereo;
1474  xRand = std::sqrt(sqr(sinhu)+1)*innerRadius*cosphi;
1475  yRand = std::sqrt(sqr(sinhu)+1)*innerRadius*sinphi;
1476  return G4ThreeVector (xRand, yRand, zRand);
1477  }
1478  else
1479  {
1480  return G4ThreeVector(innerRadius*cosphi,innerRadius*sinphi,
1482  }
1483  }
1484  else if(chose>=aOne+aTwo && chose<aOne+aTwo+aThree)
1485  {
1487  rOut2 = outerRadius2+tanOuterStereo2*halfLenZ*halfLenZ;
1488  rOut = std::sqrt(rOut2) ;
1489 
1490  do {
1491  xRand = RandFlat::shoot(-rOut,rOut) ;
1492  yRand = RandFlat::shoot(-rOut,rOut) ;
1493  r2 = xRand*xRand + yRand*yRand ;
1494  } while ( ! ( r2 >= rIn2 && r2 <= rOut2 ) ) ;
1495 
1496  zRand = halfLenZ;
1497  return G4ThreeVector (xRand, yRand, zRand);
1498  }
1499  else
1500  {
1502  rOut2 = outerRadius2+tanOuterStereo2*halfLenZ*halfLenZ;
1503  rOut = std::sqrt(rOut2) ;
1504 
1505  do {
1506  xRand = RandFlat::shoot(-rOut,rOut) ;
1507  yRand = RandFlat::shoot(-rOut,rOut) ;
1508  r2 = xRand*xRand + yRand*yRand ;
1509  } while ( ! ( r2 >= rIn2 && r2 <= rOut2 ) ) ;
1510 
1511  zRand = -1.*halfLenZ;
1512  return G4ThreeVector (xRand, yRand, zRand);
1513  }
1514 }
1515 
1516 
1517 //
1518 // DescribeYourselfTo
1519 //
1521 {
1522  scene.AddSolid (*this);
1523 }
1524 
1525 
1526 //
1527 // GetExtent
1528 //
1530 {
1531  // Define the sides of the box into which the G4Tubs instance would fit.
1532  //
1535  -halfLenZ, halfLenZ );
1536 }
1537 
1538 
1539 //
1540 // CreatePolyhedron
1541 //
1543 {
1546 }
1547 
1548 
1549 //
1550 // GetPolyhedron
1551 //
1553 {
1554  if (!fpPolyhedron ||
1556  fpPolyhedron->GetNumberOfRotationSteps())
1557  {
1558  delete fpPolyhedron;
1560  }
1561  return fpPolyhedron;
1562 }
1563 
1564 
1565 //
1566 // asinh
1567 //
1569 {
1570  return std::log(arg+std::sqrt(sqr(arg)+1));
1571 }
G4double outerStereo
Definition: G4Hype.hh:172
G4String GetName() const
G4double GetCubicVolume()
Definition: G4Hype.cc:1371
G4double outerRadius2
Definition: G4Hype.hh:181
ThreeVector shoot(const G4int Ap, const G4int Af)
G4double fHalfTol
Definition: G4Hype.hh:200
G4Hype(const G4String &pName, G4double newInnerRadius, G4double newOuterRadius, G4double newInnerStereo, G4double newOuterStereo, G4double newHalfLenZ)
Definition: G4Hype.cc:67
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Hype.cc:594
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
G4double GetSurfaceArea()
Definition: G4Hype.cc:1382
G4double innerRadius2
Definition: G4Hype.hh:180
G4GeometryType GetEntityType() const
Definition: G4Hype.cc:1353
G4double innerStereo
Definition: G4Hype.hh:171
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Hype.cc:203
static G4double ApproxDistInside(G4double pr, G4double pz, G4double r0, G4double tan2Phi)
Definition: G4Hype.cc:1329
G4bool GetExtent(G4double &min, G4double &max) const
const G4double pi
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
void SetNormal(const G4ThreeVector &newNormal)
G4double fSurfaceArea
Definition: G4Hype.hh:198
virtual G4bool PartialClip(const G4VoxelLimits &voxelLimit, const EAxis IgnoreMe)
G4Polyhedron * CreatePolyhedron() const
Definition: G4Hype.cc:1542
virtual void AddVertexInOrder(const G4ThreeVector vertex)
G4double a
Definition: TRTMaterials.hh:39
G4double endOuterRadius2
Definition: G4Hype.hh:183
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Hype.cc:1520
virtual void AddSolid(const G4Box &)=0
G4double tanOuterStereo
Definition: G4Hype.hh:177
int G4int
Definition: G4Types.hh:78
G4bool InnerSurfaceExists() const
virtual ~G4Hype()
Definition: G4Hype.cc:143
G4double HypeOuterRadius2(G4double zVal) const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4double tanOuterStereo2
Definition: G4Hype.hh:179
G4VSolid * Clone() const
Definition: G4Hype.cc:1362
static G4double ApproxDistOutside(G4double pr, G4double pz, G4double r0, G4double tanPhi)
Definition: G4Hype.cc:1271
G4double tanInnerStereo2
Definition: G4Hype.hh:178
G4double outerRadius
Definition: G4Hype.hh:169
G4double HypeInnerRadius2(G4double zVal) const
Definition: G4Hype.hh:67
G4double endOuterRadius
Definition: G4Hype.hh:185
bool G4bool
Definition: G4Types.hh:79
void AddSurface(const G4ClippablePolygon &surface)
#define DBL_EPSILON
Definition: templates.hh:87
const G4int n
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Hype.cc:970
G4Hype & operator=(const G4Hype &rhs)
Definition: G4Hype.cc:171
G4double halfLenZ
Definition: G4Hype.hh:170
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ThreeVector GetPointOnSurface() const
Definition: G4Hype.cc:1417
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
static void AddPolyToExtent(const G4ThreeVector &v0, const G4ThreeVector &v1, const G4ThreeVector &w1, const G4ThreeVector &w0, const G4VoxelLimits &voxelLimit, const EAxis axis, G4SolidExtentList &extentList)
Definition: G4Hype.cc:478
G4VisExtent GetExtent() const
Definition: G4Hype.cc:1529
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
static G4int IntersectHype(const G4ThreeVector &p, const G4ThreeVector &v, G4double r2, G4double tan2Phi, G4double s[2])
Definition: G4Hype.cc:1204
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EInside
Definition: geomdefs.hh:58
G4double fCubicVolume
Definition: G4Hype.hh:197
EAxis
Definition: geomdefs.hh:54
G4double tanInnerStereo
Definition: G4Hype.hh:176
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
EInside Inside(const G4ThreeVector &p) const
Definition: G4Hype.cc:504
#define DBL_MIN
Definition: templates.hh:75
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static const double degree
Definition: G4SIunits.hh:125
G4Polyhedron * fpPolyhedron
Definition: G4Hype.hh:202
G4double asinh(G4double arg)
Definition: G4Hype.cc:1568
#define G4endl
Definition: G4ios.hh:61
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
G4double endInnerRadius2
Definition: G4Hype.hh:182
G4double kCarTolerance
Definition: G4VSolid.hh:305
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
void SetOuterStereo(G4double newOSte)
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
G4double innerRadius
Definition: G4Hype.hh:168
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Hype.cc:214
static const G4double alpha
G4double endInnerRadius
Definition: G4Hype.hh:184
static const double mm
Definition: G4SIunits.hh:102
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:250
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Hype.cc:1393
G4Polyhedron * GetPolyhedron() const
Definition: G4Hype.cc:1552
void SetInnerStereo(G4double newISte)
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Hype.cc:548
void ApplyPointTransform(G4ThreeVector &vec) const
const G4int kMaxMeshSections
Definition: meshdefs.hh:46