source: trunk/source/geometry/solids/specific/src/G4Hype.cc @ 1202

Last change on this file since 1202 was 1058, checked in by garnier, 15 years ago

file release beta

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