Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/solids/specific/src/G4EllipticalCone.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

Diff markup

Differences between /geometry/solids/specific/src/G4EllipticalCone.cc (Version 11.3.0) and /geometry/solids/specific/src/G4EllipticalCone.cc (Version 9.6.p4)


  1 //                                                  1 //
  2 // *******************************************      2 // ********************************************************************
  3 // * License and Disclaimer                         3 // * License and Disclaimer                                           *
  4 // *                                                4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of th      5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided      6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License      7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/      8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.           9 // * include a list of copyright holders.                             *
 10 // *                                               10 // *                                                                  *
 11 // * Neither the authors of this software syst     11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing fin     12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warran     13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assum     14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file      15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitatio     16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                               17 // *                                                                  *
 18 // * This  code  implementation is the result      18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboratio     19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distri     20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  ag     21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publicati     22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Sof     23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // *******************************************     24 // ********************************************************************
 25 //                                                 25 //
                                                   >>  26 // $Id: G4EllipticalCone.cc 67011 2013-01-29 16:17:41Z gcosmo $
                                                   >>  27 //
 26 // Implementation of G4EllipticalCone class        28 // Implementation of G4EllipticalCone class
 27 //                                                 29 //
 28 // This code implements an Elliptical Cone giv     30 // This code implements an Elliptical Cone given explicitly by the
 29 // equation:                                       31 // equation:
 30 //   x^2/a^2 + y^2/b^2 = (z-h)^2                   32 //   x^2/a^2 + y^2/b^2 = (z-h)^2
 31 // and specified by the parameters (a,b,h) and     33 // and specified by the parameters (a,b,h) and a cut parallel to the
 32 // xy plane above z = 0.                           34 // xy plane above z = 0.
 33 //                                                 35 //
 34 // Author: Dionysios Anninos                       36 // Author: Dionysios Anninos
 35 // Revised: Evgueni Tcherniaev                 <<  37 //
 36 // -------------------------------------------     38 // --------------------------------------------------------------------
 37                                                    39 
 38 #if !(defined(G4GEOM_USE_UELLIPTICALCONE) && d << 
 39                                                << 
 40 #include "globals.hh"                              40 #include "globals.hh"
 41                                                    41 
 42 #include "G4EllipticalCone.hh"                     42 #include "G4EllipticalCone.hh"
 43                                                    43 
 44 #include "G4RandomTools.hh"                    << 
 45 #include "G4GeomTools.hh"                      << 
 46 #include "G4ClippablePolygon.hh"                   44 #include "G4ClippablePolygon.hh"
                                                   >>  45 #include "G4SolidExtentList.hh"
 47 #include "G4VoxelLimits.hh"                        46 #include "G4VoxelLimits.hh"
 48 #include "G4AffineTransform.hh"                    47 #include "G4AffineTransform.hh"
 49 #include "G4BoundingEnvelope.hh"               << 
 50 #include "G4GeometryTolerance.hh"                  48 #include "G4GeometryTolerance.hh"
 51                                                    49 
 52 #include "meshdefs.hh"                             50 #include "meshdefs.hh"
 53                                                    51 
 54 #include "Randomize.hh"                            52 #include "Randomize.hh"
 55                                                    53 
 56 #include "G4VGraphicsScene.hh"                     54 #include "G4VGraphicsScene.hh"
                                                   >>  55 #include "G4Polyhedron.hh"
                                                   >>  56 #include "G4NURBS.hh"
                                                   >>  57 #include "G4NURBSbox.hh"
 57 #include "G4VisExtent.hh"                          58 #include "G4VisExtent.hh"
 58                                                    59 
 59 #include "G4AutoLock.hh"                       <<  60 //#define G4SPECSDEBUG 1    
 60                                                << 
 61 namespace                                      << 
 62 {                                              << 
 63   G4Mutex polyhedronMutex = G4MUTEX_INITIALIZE << 
 64 }                                              << 
 65                                                    61 
 66 using namespace CLHEP;                             62 using namespace CLHEP;
 67                                                    63 
 68 ////////////////////////////////////////////// <<  64 //////////////////////////////////////////////////////////////////////
 69 //                                                 65 //
 70 // Constructor - check parameters                  66 // Constructor - check parameters
 71                                                <<  67 //
 72 G4EllipticalCone::G4EllipticalCone(const G4Str     68 G4EllipticalCone::G4EllipticalCone(const G4String& pName,
 73                                          G4dou     69                                          G4double  pxSemiAxis,
 74                                          G4dou     70                                          G4double  pySemiAxis,
 75                                          G4dou     71                                          G4double  pzMax,
 76                                          G4dou     72                                          G4double  pzTopCut)
 77   : G4VSolid(pName), zTopCut(0.)               <<  73   : G4VSolid(pName), fpPolyhedron(0), fCubicVolume(0.), fSurfaceArea(0.),
                                                   >>  74     zTopCut(0.)
 78 {                                                  75 {
 79   halfCarTol = 0.5*kCarTolerance;              <<  76 
                                                   >>  77   kRadTolerance = G4GeometryTolerance::GetInstance()->GetRadialTolerance();
 80                                                    78 
 81   // Check Semi-Axis & Z-cut                       79   // Check Semi-Axis & Z-cut
 82   //                                               80   //
 83   if ( (pxSemiAxis <= 0.) || (pySemiAxis <= 0.     81   if ( (pxSemiAxis <= 0.) || (pySemiAxis <= 0.) || (pzMax <= 0.) )
 84   {                                                82   {
 85      std::ostringstream message;                   83      std::ostringstream message;
 86      message << "Invalid semi-axis or height f <<  84      message << "Invalid semi-axis or height - " << GetName();
 87              << "\n   X semi-axis, Y semi-axis << 
 88              << pxSemiAxis << ", " << pySemiAx << 
 89      G4Exception("G4EllipticalCone::G4Elliptic     85      G4Exception("G4EllipticalCone::G4EllipticalCone()", "GeomSolids0002",
 90                  FatalErrorInArgument, message     86                  FatalErrorInArgument, message);
 91    }                                           <<  87   }
 92                                                << 
 93   if ( pzTopCut <= 0 )                             88   if ( pzTopCut <= 0 )
 94   {                                                89   {
 95      std::ostringstream message;                   90      std::ostringstream message;
 96      message << "Invalid z-coordinate for cutt <<  91      message << "Invalid z-coordinate for cutting plane - " << GetName();
 97              << "\n   Z top cut = " << pzTopCu <<  92      G4Exception("G4EllipticalCone::G4EllipticalCone()", "InvalidSetup",
 98      G4Exception("G4EllipticalCone::G4Elliptic << 
 99                  FatalErrorInArgument, message     93                  FatalErrorInArgument, message);
100   }                                                94   }
101                                                    95 
102   SetSemiAxis( pxSemiAxis, pySemiAxis, pzMax )     96   SetSemiAxis( pxSemiAxis, pySemiAxis, pzMax );
103   SetZCut(pzTopCut);                               97   SetZCut(pzTopCut);
104 }                                                  98 }
105                                                    99 
106 ////////////////////////////////////////////// << 100 ///////////////////////////////////////////////////////////////////////////////
107 //                                                101 //
108 // Fake default constructor - sets only member    102 // Fake default constructor - sets only member data and allocates memory
109 //                            for usage restri    103 //                            for usage restricted to object persistency.
110                                                << 104 //
111 G4EllipticalCone::G4EllipticalCone( __void__&     105 G4EllipticalCone::G4EllipticalCone( __void__& a )
112   : G4VSolid(a), halfCarTol(0.),               << 106   : G4VSolid(a), fpPolyhedron(0), kRadTolerance(0.), fCubicVolume(0.),
113     xSemiAxis(0.), ySemiAxis(0.), zheight(0.), << 107     fSurfaceArea(0.), xSemiAxis(0.), ySemiAxis(0.), zheight(0.),
114     cosAxisMin(0.), invXX(0.), invYY(0.)       << 108     semiAxisMax(0.), zTopCut(0.)
115 {                                                 109 {
116 }                                                 110 }
117                                                   111 
118 ////////////////////////////////////////////// << 112 ///////////////////////////////////////////////////////////////////////////////
119 //                                                113 //
120 // Destructor                                     114 // Destructor
121                                                << 115 //
122 G4EllipticalCone::~G4EllipticalCone()             116 G4EllipticalCone::~G4EllipticalCone()
123 {                                                 117 {
124   delete fpPolyhedron; fpPolyhedron = nullptr; << 
125 }                                                 118 }
126                                                   119 
127 ////////////////////////////////////////////// << 120 ///////////////////////////////////////////////////////////////////////////////
128 //                                                121 //
129 // Copy constructor                               122 // Copy constructor
130                                                << 123 //
131 G4EllipticalCone::G4EllipticalCone(const G4Ell    124 G4EllipticalCone::G4EllipticalCone(const G4EllipticalCone& rhs)
132   : G4VSolid(rhs), halfCarTol(rhs.halfCarTol), << 125   : G4VSolid(rhs),
                                                   >> 126     fpPolyhedron(0), kRadTolerance(rhs.kRadTolerance),
133     fCubicVolume(rhs.fCubicVolume), fSurfaceAr    127     fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
134     xSemiAxis(rhs.xSemiAxis), ySemiAxis(rhs.yS << 128     xSemiAxis(rhs.xSemiAxis), ySemiAxis(rhs.ySemiAxis), zheight(rhs.zheight),
135     zheight(rhs.zheight), zTopCut(rhs.zTopCut) << 129     semiAxisMax(rhs.semiAxisMax), zTopCut(rhs.zTopCut)
136     cosAxisMin(rhs.cosAxisMin), invXX(rhs.invX << 
137 {                                                 130 {
138 }                                                 131 }
139                                                   132 
140 ////////////////////////////////////////////// << 133 ///////////////////////////////////////////////////////////////////////////////
141 //                                                134 //
142 // Assignment operator                            135 // Assignment operator
143                                                << 136 //
144 G4EllipticalCone& G4EllipticalCone::operator =    137 G4EllipticalCone& G4EllipticalCone::operator = (const G4EllipticalCone& rhs) 
145 {                                                 138 {
146    // Check assignment to self                    139    // Check assignment to self
147    //                                             140    //
148    if (this == &rhs)  { return *this; }           141    if (this == &rhs)  { return *this; }
149                                                   142 
150    // Copy base class data                        143    // Copy base class data
151    //                                             144    //
152    G4VSolid::operator=(rhs);                      145    G4VSolid::operator=(rhs);
153                                                   146 
154    // Copy data                                   147    // Copy data
155    //                                             148    //
156    halfCarTol = rhs.halfCarTol;                << 149    fpPolyhedron = 0; kRadTolerance = rhs.kRadTolerance;
157    fCubicVolume = rhs.fCubicVolume; fSurfaceAr    150    fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
158    xSemiAxis = rhs.xSemiAxis; ySemiAxis = rhs.    151    xSemiAxis = rhs.xSemiAxis; ySemiAxis = rhs.ySemiAxis;
159    zheight = rhs.zheight; zTopCut = rhs.zTopCu << 152    zheight = rhs.zheight; semiAxisMax = rhs.semiAxisMax; zTopCut = rhs.zTopCut;
160    cosAxisMin = rhs.cosAxisMin; invXX = rhs.in << 
161                                                << 
162    fRebuildPolyhedron = false;                 << 
163    delete fpPolyhedron; fpPolyhedron = nullptr << 
164                                                   153 
165    return *this;                                  154    return *this;
166 }                                                 155 }
167                                                   156 
168 ////////////////////////////////////////////// << 157 ///////////////////////////////////////////////////////////////////////////////
169 //                                                158 //
170 // Get bounding box                            << 159 // Calculate extent under transform and specified limit
171                                                << 160 //
172 void G4EllipticalCone::BoundingLimits(G4ThreeV << 161 G4bool
173                                       G4ThreeV << 162 G4EllipticalCone::CalculateExtent( const EAxis axis,
                                                   >> 163                                    const G4VoxelLimits &voxelLimit,
                                                   >> 164                                    const G4AffineTransform &transform,
                                                   >> 165                                          G4double &min, G4double &max ) const
174 {                                                 166 {
175   G4double zcut   = GetZTopCut();              << 167   G4SolidExtentList  extentList( axis, voxelLimit );
176   G4double height = GetZMax();                 << 168   
177   G4double xmax   = GetSemiAxisX()*(height+zcu << 169   //
178   G4double ymax   = GetSemiAxisY()*(height+zcu << 170   // We are going to divide up our elliptical face into small pieces
179   pMin.set(-xmax,-ymax,-zcut);                 << 171   //
180   pMax.set( xmax, ymax, zcut);                 << 172   
181                                                << 173   //
182   // Check correctness of the bounding box     << 174   // Choose phi size of our segment(s) based on constants as
183   //                                           << 175   // defined in meshdefs.hh
184   if (pMin.x() >= pMax.x() || pMin.y() >= pMax << 176   //
185   {                                            << 177   G4int numPhi = kMaxMeshSections;
186     std::ostringstream message;                << 178   G4double sigPhi = twopi/numPhi;
187     message << "Bad bounding box (min >= max)  << 179   
188             << GetName() << " !"               << 180   //
189             << "\npMin = " << pMin             << 181   // We have to be careful to keep our segments completely outside
190             << "\npMax = " << pMax;            << 182   // of the elliptical surface. To do so we imagine we have
191     G4Exception("G4EllipticalCone::BoundingLim << 183   // a simple (unit radius) circular cross section (as in G4Tubs) 
192                 JustWarning, message);         << 184   // and then "stretch" the dimensions as necessary to fit the ellipse.
193     DumpInfo();                                << 185   //
194   }                                            << 186   G4double rFudge = 1.0/std::cos(0.5*sigPhi);
195 }                                              << 187   G4double dxFudgeBot = xSemiAxis*2.*zheight*rFudge,
                                                   >> 188            dyFudgeBot = ySemiAxis*2.*zheight*rFudge;
                                                   >> 189   G4double dxFudgeTop = xSemiAxis*(zheight-zTopCut)*rFudge,
                                                   >> 190            dyFudgeTop = ySemiAxis*(zheight-zTopCut)*rFudge;
                                                   >> 191   
                                                   >> 192   //
                                                   >> 193   // As we work around the elliptical surface, we build
                                                   >> 194   // a "phi" segment on the way, and keep track of two
                                                   >> 195   // additional polygons for the two ends.
                                                   >> 196   //
                                                   >> 197   G4ClippablePolygon endPoly1, endPoly2, phiPoly;
                                                   >> 198   
                                                   >> 199   G4double phi = 0, 
                                                   >> 200            cosPhi = std::cos(phi),
                                                   >> 201            sinPhi = std::sin(phi);
                                                   >> 202   G4ThreeVector v0( dxFudgeTop*cosPhi, dyFudgeTop*sinPhi, +zTopCut ),
                                                   >> 203                 v1( dxFudgeBot*cosPhi, dyFudgeBot*sinPhi, -zTopCut ),
                                                   >> 204                 w0, w1;
                                                   >> 205   transform.ApplyPointTransform( v0 );
                                                   >> 206   transform.ApplyPointTransform( v1 );
                                                   >> 207   do
                                                   >> 208   {
                                                   >> 209     phi += sigPhi;
                                                   >> 210     if (numPhi == 1) phi = 0;  // Try to avoid roundoff
                                                   >> 211     cosPhi = std::cos(phi), 
                                                   >> 212     sinPhi = std::sin(phi);
                                                   >> 213     
                                                   >> 214     w0 = G4ThreeVector( dxFudgeTop*cosPhi, dyFudgeTop*sinPhi, +zTopCut );
                                                   >> 215     w1 = G4ThreeVector( dxFudgeBot*cosPhi, dyFudgeBot*sinPhi, -zTopCut );
                                                   >> 216     transform.ApplyPointTransform( w0 );
                                                   >> 217     transform.ApplyPointTransform( w1 );
                                                   >> 218     
                                                   >> 219     //
                                                   >> 220     // Add a point to our z ends
                                                   >> 221     //
                                                   >> 222     endPoly1.AddVertexInOrder( v0 );
                                                   >> 223     endPoly2.AddVertexInOrder( v1 );
                                                   >> 224     
                                                   >> 225     //
                                                   >> 226     // Build phi polygon
                                                   >> 227     //
                                                   >> 228     phiPoly.ClearAllVertices();
                                                   >> 229     
                                                   >> 230     phiPoly.AddVertexInOrder( v0 );
                                                   >> 231     phiPoly.AddVertexInOrder( v1 );
                                                   >> 232     phiPoly.AddVertexInOrder( w1 );
                                                   >> 233     phiPoly.AddVertexInOrder( w0 );
                                                   >> 234     
                                                   >> 235     if (phiPoly.PartialClip( voxelLimit, axis ))
                                                   >> 236     {
                                                   >> 237       //
                                                   >> 238       // Get unit normal
                                                   >> 239       //
                                                   >> 240       phiPoly.SetNormal( (v1-v0).cross(w0-v0).unit() );
                                                   >> 241       
                                                   >> 242       extentList.AddSurface( phiPoly );
                                                   >> 243     }
196                                                   244 
197 ////////////////////////////////////////////// << 245     //
198 //                                             << 246     // Next vertex
199 // Calculate extent under transform and specif << 247     //    
                                                   >> 248     v0 = w0;
                                                   >> 249     v1 = w1;
                                                   >> 250   } while( --numPhi > 0 );
200                                                   251 
201 G4bool                                         << 252   //
202 G4EllipticalCone::CalculateExtent(const EAxis  << 253   // Process the end pieces
203                                   const G4Voxe << 254   //
204                                   const G4Affi << 255   if (endPoly1.PartialClip( voxelLimit, axis ))
205                                         G4doub << 
206 {                                              << 
207   G4ThreeVector bmin,bmax;                     << 
208   G4bool exist;                                << 
209                                                << 
210   // Check bounding box (bbox)                 << 
211   //                                           << 
212   BoundingLimits(bmin,bmax);                   << 
213   G4BoundingEnvelope bbox(bmin,bmax);          << 
214 #ifdef G4BBOX_EXTENT                           << 
215   return bbox.CalculateExtent(pAxis,pVoxelLimi << 
216 #endif                                         << 
217   if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVox << 
218   {                                               256   {
219     return exist = pMin < pMax;                << 257     static const G4ThreeVector normal(0,0,+1);
                                                   >> 258     endPoly1.SetNormal( transform.TransformAxis(normal) );
                                                   >> 259     extentList.AddSurface( endPoly1 );
220   }                                               260   }
221                                                << 261   
222   // Set bounding envelope (benv) and calculat << 262   if (endPoly2.PartialClip( voxelLimit, axis ))
                                                   >> 263   {
                                                   >> 264     static const G4ThreeVector normal(0,0,-1);
                                                   >> 265     endPoly2.SetNormal( transform.TransformAxis(normal) );
                                                   >> 266     extentList.AddSurface( endPoly2 );
                                                   >> 267   }
                                                   >> 268   
223   //                                              269   //
224   static const G4int NSTEPS = 48; // number of << 270   // Return min/max value
225   static const G4double ang = twopi/NSTEPS;    << 271   //
226   static const G4double sinHalf = std::sin(0.5 << 272   return extentList.GetExtent( min, max );
227   static const G4double cosHalf = std::cos(0.5 << 
228   static const G4double sinStep = 2.*sinHalf*c << 
229   static const G4double cosStep = 1. - 2.*sinH << 
230   G4double zcut   = bmax.z();                  << 
231   G4double height = GetZMax();                 << 
232   G4double sxmin  = GetSemiAxisX()*(height-zcu << 
233   G4double symin  = GetSemiAxisY()*(height-zcu << 
234   G4double sxmax  = bmax.x()/cosHalf;          << 
235   G4double symax  = bmax.y()/cosHalf;          << 
236                                                << 
237   G4double sinCur = sinHalf;                   << 
238   G4double cosCur = cosHalf;                   << 
239   G4ThreeVectorList baseA(NSTEPS),baseB(NSTEPS << 
240   for (G4int k=0; k<NSTEPS; ++k)               << 
241   {                                            << 
242     baseA[k].set(sxmax*cosCur,symax*sinCur,-zc << 
243     baseB[k].set(sxmin*cosCur,symin*sinCur, zc << 
244                                                << 
245     G4double sinTmp = sinCur;                  << 
246     sinCur = sinCur*cosStep + cosCur*sinStep;  << 
247     cosCur = cosCur*cosStep - sinTmp*sinStep;  << 
248   }                                            << 
249                                                << 
250   std::vector<const G4ThreeVectorList *> polyg << 
251   polygons[0] = &baseA;                        << 
252   polygons[1] = &baseB;                        << 
253   G4BoundingEnvelope benv(bmin,bmax,polygons); << 
254   exist = benv.CalculateExtent(pAxis,pVoxelLim << 
255   return exist;                                << 
256 }                                                 273 }
257                                                   274 
258 ////////////////////////////////////////////// << 275 ////////////////////////////////////////////////////////////////////////
                                                   >> 276 //
                                                   >> 277 // Return whether point inside/outside/on surface
                                                   >> 278 // Split into radius, phi, theta checks
                                                   >> 279 // Each check modifies `in', or returns as approprate
259 //                                                280 //
260 // Determine where is point: inside, outside o << 
261                                                << 
262 EInside G4EllipticalCone::Inside(const G4Three    281 EInside G4EllipticalCone::Inside(const G4ThreeVector& p) const
263 {                                                 282 {
264   G4double hp = std::sqrt(p.x()*p.x()*invXX +  << 283   G4double rad2oo,  // outside surface outer tolerance
265   G4double ds = (hp - zheight)*cosAxisMin;     << 284            rad2oi;  // outside surface inner tolerance
266   G4double dz = std::abs(p.z()) - zTopCut;     << 285   
267   G4double dist = std::max(ds,dz);             << 286   EInside in;
                                                   >> 287 
                                                   >> 288   static const G4double halfRadTol = 0.5*kRadTolerance;
                                                   >> 289   static const G4double halfCarTol = 0.5*kCarTolerance;
                                                   >> 290 
                                                   >> 291   // check this side of z cut first, because that's fast
                                                   >> 292   //
268                                                   293 
269   if (dist > halfCarTol) return kOutside;      << 294   if ( (p.z() < -zTopCut - halfCarTol)
270   return (dist > -halfCarTol) ? kSurface : kIn << 295     || (p.z() > zTopCut + halfCarTol ) )
                                                   >> 296   {
                                                   >> 297     return in = kOutside; 
                                                   >> 298   }
                                                   >> 299 
                                                   >> 300   rad2oo= sqr(p.x()/( xSemiAxis + halfRadTol ))
                                                   >> 301         + sqr(p.y()/( ySemiAxis + halfRadTol ));
                                                   >> 302 
                                                   >> 303   if ( rad2oo > sqr( zheight-p.z() ) )
                                                   >> 304   {
                                                   >> 305     return in = kOutside; 
                                                   >> 306   }
                                                   >> 307 
                                                   >> 308   //  rad2oi= sqr( p.x()*(1.0 + 0.5*kRadTolerance/(xSemiAxis*xSemiAxis)) )
                                                   >> 309   //      + sqr( p.y()*(1.0 + 0.5*kRadTolerance/(ySemiAxis*ySemiAxis)) );
                                                   >> 310   rad2oi = sqr(p.x()/( xSemiAxis - halfRadTol ))
                                                   >> 311         + sqr(p.y()/( ySemiAxis - halfRadTol ));
                                                   >> 312      
                                                   >> 313   if (rad2oi < sqr( zheight-p.z() ) )
                                                   >> 314   {
                                                   >> 315     in = ( ( p.z() < -zTopCut + halfRadTol )
                                                   >> 316         || ( p.z() >  zTopCut - halfRadTol ) ) ? kSurface : kInside;
                                                   >> 317   }
                                                   >> 318   else 
                                                   >> 319   {
                                                   >> 320     in = kSurface;
                                                   >> 321   }
                                                   >> 322 
                                                   >> 323   return in;
271 }                                                 324 }
272                                                   325 
273 //////////////////////////////////////////////    326 /////////////////////////////////////////////////////////////////////////
274 //                                                327 //
275 // Return unit normal at surface closest to p  << 328 // Return unit normal of surface closest to p not protected against p=0
276                                                << 329 //
277 G4ThreeVector G4EllipticalCone::SurfaceNormal(    330 G4ThreeVector G4EllipticalCone::SurfaceNormal( const G4ThreeVector& p) const
278 {                                                 331 {
279   G4ThreeVector norm(0,0,0);                   << 
280   G4int nsurf = 0;  // number of surfaces wher << 
281                                                   332 
282   G4double hp = std::sqrt(p.x()*p.x()*invXX +  << 333   G4double rx = sqr(p.x()/xSemiAxis), 
283   G4double ds = (hp - zheight)*cosAxisMin;     << 334            ry = sqr(p.y()/ySemiAxis);
284   if (std::abs(ds) <= halfCarTol)              << 335 
                                                   >> 336   G4double rds = std::sqrt(rx + ry); 
                                                   >> 337 
                                                   >> 338   G4ThreeVector norm;
                                                   >> 339 
                                                   >> 340   if( (p.z() < -zTopCut) && ((rx+ry) < sqr(zTopCut + zheight)) )
285   {                                               341   {
286     norm = G4ThreeVector(p.x()*invXX, p.y()*in << 342     return G4ThreeVector( 0., 0., -1. ); 
287     G4double mag = norm.mag();                 << 
288     if (mag == 0) return {0,0,1}; // apex      << 
289     norm *= (1/mag);                           << 
290     ++nsurf;                                   << 
291   }                                               343   }
292   G4double dz = std::abs(p.z()) - zTopCut;     << 344 
293   if (std::abs(dz) <= halfCarTol)              << 345   if( (p.z() > (zheight > zTopCut ? zheight : zTopCut)) &&
                                                   >> 346       ((rx+ry) < sqr(zheight-zTopCut)) )
294   {                                               347   {
295     norm += G4ThreeVector(0., 0.,(p.z() < 0) ? << 348     return G4ThreeVector( 0., 0., 1. );
296     ++nsurf;                                   << 
297   }                                               349   }
298                                                   350 
299   if      (nsurf == 1) return norm;            << 351   if( p.z() > rds + 2.*zTopCut - zheight ) 
300   else if (nsurf >  1) return norm.unit(); //  << 
301   else                                         << 
302   {                                               352   {
303     // Point is not on the surface             << 353     if ( p.z() > zTopCut )
304     //                                         << 354     {
305 #ifdef G4CSGDEBUG                              << 355       if( p.x() == 0. ) 
306     std::ostringstream message;                << 356       {
307     G4long oldprc = message.precision(16);     << 357         norm = G4ThreeVector( 0., p.y() < 0. ? -1. : 1., 1. ); 
308     message << "Point p is not on surface (!?) << 358         return norm /= norm.mag();
309             << GetName() << G4endl;            << 359       } 
310     message << "Position:\n";                  << 360       if( p.y() == 0. )
311     message << "   p.x() = " << p.x()/mm << "  << 361       {
312     message << "   p.y() = " << p.y()/mm << "  << 362         norm = G4ThreeVector( p.x() < 0. ? -1. : 1., 0., 1. ); 
313     message << "   p.z() = " << p.z()/mm << "  << 363         return norm /= norm.mag();
314     G4cout.precision(oldprc);                  << 364       } 
315     G4Exception("G4EllipticalCone::SurfaceNorm << 365       
316                 JustWarning, message );        << 366       G4double k =  std::fabs(p.x()/p.y());
317     DumpInfo();                                << 367       G4double c2 = sqr(zheight-zTopCut)/(1./sqr(xSemiAxis)+sqr(k/ySemiAxis));
318 #endif                                         << 368       G4double x  = std::sqrt(c2);
319     return ApproxSurfaceNormal(p);             << 369       G4double y  = k*x;
                                                   >> 370         
                                                   >> 371       x /= sqr(xSemiAxis);
                                                   >> 372       y /= sqr(ySemiAxis);
                                                   >> 373       
                                                   >> 374       norm = G4ThreeVector( p.x() < 0. ? -x : x, 
                                                   >> 375                             p.y() < 0. ? -y : y,
                                                   >> 376                             - ( zheight - zTopCut ) );
                                                   >> 377       norm /= norm.mag();
                                                   >> 378       norm += G4ThreeVector( 0., 0., 1. );
                                                   >> 379       return norm /= norm.mag();      
                                                   >> 380     }
                                                   >> 381     
                                                   >> 382     return G4ThreeVector( 0., 0., 1. );    
320   }                                               383   }
321 }                                              << 384   
                                                   >> 385   if( p.z() < rds - 2.*zTopCut - zheight )
                                                   >> 386   {
                                                   >> 387     if( p.x() == 0. ) 
                                                   >> 388     {
                                                   >> 389       norm = G4ThreeVector( 0., p.y() < 0. ? -1. : 1., -1. ); 
                                                   >> 390       return norm /= norm.mag();
                                                   >> 391     } 
                                                   >> 392     if( p.y() == 0. )
                                                   >> 393     {
                                                   >> 394       norm = G4ThreeVector( p.x() < 0. ? -1. : 1., 0., -1. ); 
                                                   >> 395       return norm /= norm.mag();
                                                   >> 396     } 
                                                   >> 397     
                                                   >> 398     G4double k =  std::fabs(p.x()/p.y());
                                                   >> 399     G4double c2 = sqr(zheight+zTopCut)/(1./sqr(xSemiAxis)+sqr(k/ySemiAxis));
                                                   >> 400     G4double x  = std::sqrt(c2);
                                                   >> 401     G4double y  = k*x;
                                                   >> 402     
                                                   >> 403     x /= sqr(xSemiAxis);
                                                   >> 404     y /= sqr(ySemiAxis);
                                                   >> 405     
                                                   >> 406     norm = G4ThreeVector( p.x() < 0. ? -x : x, 
                                                   >> 407                           p.y() < 0. ? -y : y,
                                                   >> 408                           - ( zheight - zTopCut ) );
                                                   >> 409     norm /= norm.mag();
                                                   >> 410     norm += G4ThreeVector( 0., 0., -1. );
                                                   >> 411     return norm /= norm.mag();      
                                                   >> 412   }
                                                   >> 413     
                                                   >> 414   norm  = G4ThreeVector(p.x()/sqr(xSemiAxis), p.y()/sqr(ySemiAxis), rds);
                                                   >> 415    
                                                   >> 416   G4double k = std::tan(pi/8.);
                                                   >> 417   G4double c = -zTopCut - k*(zTopCut + zheight);
322                                                   418 
323 ////////////////////////////////////////////// << 419   if( p.z() < -k*rds + c )
324 //                                             << 420     return G4ThreeVector (0.,0.,-1.);
325 // Find surface nearest to point and return co << 421 
326 // The algorithm is similar to the algorithm u << 422   return norm /= norm.mag();
327 // This method normally should not be called.  << 
328                                                << 
329 G4ThreeVector                                  << 
330 G4EllipticalCone::ApproxSurfaceNormal(const G4 << 
331 {                                              << 
332   G4double hp = std::sqrt(p.x()*p.x()*invXX +  << 
333   G4double ds = (hp - zheight)*cosAxisMin;     << 
334   G4double dz = std::abs(p.z()) - zTopCut;     << 
335   if (ds > dz && std::abs(hp - p.z()) > halfCa << 
336     return G4ThreeVector(p.x()*invXX, p.y()*in << 
337   else                                         << 
338     return { 0., 0., (G4double)((p.z() < 0) ?  << 
339 }                                                 423 }
340                                                   424 
341 ////////////////////////////////////////////// << 425 //////////////////////////////////////////////////////////////////////////
342 //                                                426 //
343 // Calculate distance to shape from outside, a    427 // Calculate distance to shape from outside, along normalised vector
344 // return kInfinity if no intersection, or int    428 // return kInfinity if no intersection, or intersection distance <= tolerance
345                                                << 429 //
346 G4double G4EllipticalCone::DistanceToIn( const    430 G4double G4EllipticalCone::DistanceToIn( const G4ThreeVector& p,
347                                          const    431                                          const G4ThreeVector& v  ) const
348 {                                                 432 {
                                                   >> 433   static const G4double halfTol = 0.5*kCarTolerance;
                                                   >> 434 
349   G4double distMin = kInfinity;                   435   G4double distMin = kInfinity;
350                                                   436 
351   // code from EllipticalTube                     437   // code from EllipticalTube
352                                                   438 
353   G4double sigz = p.z()+zTopCut;                  439   G4double sigz = p.z()+zTopCut;
354                                                   440 
355   //                                              441   //
356   // Check z = -dz planer surface                 442   // Check z = -dz planer surface
357   //                                              443   //
358                                                   444 
359   if (sigz < halfCarTol)                       << 445   if (sigz < halfTol)
360   {                                               446   {
361     //                                            447     //
362     // We are "behind" the shape in z, and so     448     // We are "behind" the shape in z, and so can
363     // potentially hit the rear face. Correct     449     // potentially hit the rear face. Correct direction?
364     //                                            450     //
365     if (v.z() <= 0)                               451     if (v.z() <= 0)
366     {                                             452     {
367       //                                          453       //
368       // As long as we are far enough away, we    454       // As long as we are far enough away, we know we
369       // can't intersect                          455       // can't intersect
370       //                                          456       //
371       if (sigz < 0) return kInfinity;             457       if (sigz < 0) return kInfinity;
372                                                   458       
373       //                                          459       //
374       // Otherwise, we don't intersect unless     460       // Otherwise, we don't intersect unless we are
375       // on the surface of the ellipse            461       // on the surface of the ellipse
376       //                                          462       //
377                                                   463 
378       if ( sqr(p.x()/( xSemiAxis - halfCarTol  << 464       if ( sqr(p.x()/( xSemiAxis - halfTol ))
379          + sqr(p.y()/( ySemiAxis - halfCarTol  << 465          + sqr(p.y()/( ySemiAxis - halfTol )) <= sqr( zheight+zTopCut ) )
380         return kInfinity;                         466         return kInfinity;
381                                                   467 
382     }                                             468     }
383     else                                          469     else
384     {                                             470     {
385       //                                          471       //
386       // How far?                                 472       // How far?
387       //                                          473       //
388       G4double q = -sigz/v.z();                   474       G4double q = -sigz/v.z();
389                                                   475       
390       //                                          476       //
391       // Where does that place us?                477       // Where does that place us?
392       //                                          478       //
393       G4double xi = p.x() + q*v.x(),              479       G4double xi = p.x() + q*v.x(),
394                yi = p.y() + q*v.y();              480                yi = p.y() + q*v.y();
395                                                   481       
396       //                                          482       //
397       // Is this on the surface (within ellips    483       // Is this on the surface (within ellipse)?
398       //                                          484       //
399       if ( sqr(xi/xSemiAxis) + sqr(yi/ySemiAxi    485       if ( sqr(xi/xSemiAxis) + sqr(yi/ySemiAxis) <= sqr( zheight + zTopCut ) )
400       {                                           486       {
401         //                                        487         //
402         // Yup. Return q, unless we are on the    488         // Yup. Return q, unless we are on the surface
403         //                                        489         //
404         return (sigz < -halfCarTol) ? q : 0;   << 490         return (sigz < -halfTol) ? q : 0;
405       }                                           491       }
406       else if (xi/(xSemiAxis*xSemiAxis)*v.x()     492       else if (xi/(xSemiAxis*xSemiAxis)*v.x()
407              + yi/(ySemiAxis*ySemiAxis)*v.y()     493              + yi/(ySemiAxis*ySemiAxis)*v.y() >= 0)
408       {                                           494       {
409         //                                        495         //
410         // Else, if we are traveling outwards,    496         // Else, if we are traveling outwards, we know
411         // we must miss                           497         // we must miss
412         //                                        498         //
413         //        return kInfinity;               499         //        return kInfinity;
414       }                                           500       }
415     }                                             501     }
416   }                                               502   }
417                                                   503 
418   //                                              504   //
419   // Check z = +dz planer surface                 505   // Check z = +dz planer surface
420   //                                              506   //
421   sigz = p.z() - zTopCut;                         507   sigz = p.z() - zTopCut;
422                                                   508   
423   if (sigz > -halfCarTol)                      << 509   if (sigz > -halfTol)
424   {                                               510   {
425     if (v.z() >= 0)                               511     if (v.z() >= 0)
426     {                                             512     {
427                                                   513 
428       if (sigz > 0) return kInfinity;             514       if (sigz > 0) return kInfinity;
429                                                   515 
430       if ( sqr(p.x()/( xSemiAxis - halfCarTol  << 516       if ( sqr(p.x()/( xSemiAxis - halfTol ))
431          + sqr(p.y()/( ySemiAxis - halfCarTol  << 517          + sqr(p.y()/( ySemiAxis - halfTol )) <= sqr( zheight-zTopCut ) )
432         return kInfinity;                         518         return kInfinity;
433                                                   519 
434     }                                             520     }
435     else {                                        521     else {
436       G4double q = -sigz/v.z();                   522       G4double q = -sigz/v.z();
437                                                   523 
438       G4double xi = p.x() + q*v.x(),              524       G4double xi = p.x() + q*v.x(),
439                yi = p.y() + q*v.y();              525                yi = p.y() + q*v.y();
440                                                   526 
441       if ( sqr(xi/xSemiAxis) + sqr(yi/ySemiAxi    527       if ( sqr(xi/xSemiAxis) + sqr(yi/ySemiAxis) <= sqr( zheight - zTopCut ) )
442       {                                           528       {
443         return (sigz > -halfCarTol) ? q : 0;   << 529         return (sigz > -halfTol) ? q : 0;
444       }                                           530       }
445       else if (xi/(xSemiAxis*xSemiAxis)*v.x()     531       else if (xi/(xSemiAxis*xSemiAxis)*v.x()
446              + yi/(ySemiAxis*ySemiAxis)*v.y()     532              + yi/(ySemiAxis*ySemiAxis)*v.y() >= 0)
447       {                                           533       {
448         //        return kInfinity;               534         //        return kInfinity;
449       }                                           535       }
450     }                                             536     }
451   }                                               537   }
452                                                   538 
453                                                   539 
454 #if 0                                             540 #if 0
455                                                   541 
456   // check to see if Z plane is relevant          542   // check to see if Z plane is relevant
457   //                                              543   //
458   if (p.z() < -zTopCut - halfCarTol)           << 544   if (p.z() < -zTopCut - 0.5*kCarTolerance)
459   {                                               545   {
460     if (v.z() <= 0.0)                             546     if (v.z() <= 0.0)
461       return distMin;                             547       return distMin; 
462                                                   548 
463     G4double lambda = (-zTopCut - p.z())/v.z()    549     G4double lambda = (-zTopCut - p.z())/v.z();
464                                                   550     
465     if ( sqr((lambda*v.x()+p.x())/xSemiAxis) +    551     if ( sqr((lambda*v.x()+p.x())/xSemiAxis) + 
466          sqr((lambda*v.y()+p.y())/ySemiAxis) <    552          sqr((lambda*v.y()+p.y())/ySemiAxis) <=
467          sqr(zTopCut + zheight + halfCarTol) ) << 553          sqr(zTopCut + zheight + 0.5*kRadTolerance) ) 
468     {                                             554     { 
469       return distMin = std::fabs(lambda);         555       return distMin = std::fabs(lambda);    
470     }                                             556     }
471   }                                               557   }
472                                                   558 
473   if (p.z() > zTopCut + halfCarTol)            << 559   if (p.z() > zTopCut+0.5*kCarTolerance) 
474   {                                               560   {
475     if (v.z() >= 0.0)                             561     if (v.z() >= 0.0)
476       { return distMin; }                         562       { return distMin; }
477                                                   563 
478     G4double lambda  = (zTopCut - p.z()) / v.z    564     G4double lambda  = (zTopCut - p.z()) / v.z();
479                                                   565 
480     if ( sqr((lambda*v.x() + p.x())/xSemiAxis)    566     if ( sqr((lambda*v.x() + p.x())/xSemiAxis) + 
481          sqr((lambda*v.y() + p.y())/ySemiAxis)    567          sqr((lambda*v.y() + p.y())/ySemiAxis) <=
482          sqr(zheight - zTopCut + halfCarTol) ) << 568          sqr(zheight - zTopCut + 0.5*kRadTolerance) )
483       {                                           569       {
484         return distMin = std::fabs(lambda);       570         return distMin = std::fabs(lambda);
485       }                                           571       }
486   }                                               572   }
487                                                   573   
488   if (p.z() > zTopCut - halfCarTol             << 574   if (p.z() > zTopCut - halfTol
489    && p.z() < zTopCut + halfCarTol )           << 575    && p.z() < zTopCut + halfTol )
490   {                                               576   {
491     if (v.z() > 0.)                               577     if (v.z() > 0.) 
492       { return kInfinity; }                       578       { return kInfinity; }
493                                                   579 
494     return distMin = 0.;                          580     return distMin = 0.;
495   }                                               581   }
496                                                   582   
497   if (p.z() < -zTopCut + halfCarTol            << 583   if (p.z() < -zTopCut + halfTol
498    && p.z() > -zTopCut - halfCarTol)           << 584    && p.z() > -zTopCut - halfTol)
499   {                                               585   {
500     if (v.z() < 0.)                               586     if (v.z() < 0.)
501       { return distMin = kInfinity; }             587       { return distMin = kInfinity; }
502                                                   588     
503     return distMin = 0.;                          589     return distMin = 0.;
504   }                                               590   }
505                                                   591   
506 #endif                                            592 #endif
507                                                   593 
508   // if we are here then it either intersects     594   // if we are here then it either intersects or grazes the curved surface 
509   // or it does not intersect at all              595   // or it does not intersect at all
510   //                                              596   //
511   G4double A = sqr(v.x()/xSemiAxis) + sqr(v.y(    597   G4double A = sqr(v.x()/xSemiAxis) + sqr(v.y()/ySemiAxis) - sqr(v.z());
512   G4double B = 2*(v.x()*p.x()/sqr(xSemiAxis) +    598   G4double B = 2*(v.x()*p.x()/sqr(xSemiAxis) + 
513                   v.y()*p.y()/sqr(ySemiAxis) +    599                   v.y()*p.y()/sqr(ySemiAxis) + v.z()*(zheight-p.z()));
514   G4double C = sqr(p.x()/xSemiAxis) + sqr(p.y(    600   G4double C = sqr(p.x()/xSemiAxis) + sqr(p.y()/ySemiAxis) - 
515                sqr(zheight - p.z());              601                sqr(zheight - p.z());
516                                                   602  
517   G4double discr = B*B - 4.*A*C;                  603   G4double discr = B*B - 4.*A*C;
518                                                   604    
519   // if the discriminant is negative it never     605   // if the discriminant is negative it never hits the curved object
520   //                                              606   //
521   if ( discr < -halfCarTol )                   << 607   if ( discr < -halfTol )
522     { return distMin; }                           608     { return distMin; }
523                                                   609   
524   // case below is when it hits or grazes the     610   // case below is when it hits or grazes the surface
525   //                                              611   //
526   if ( (discr >= -halfCarTol ) && (discr < hal << 612   if ( (discr >= - halfTol ) && (discr < halfTol ) )
527   {                                               613   {
528     return distMin = std::fabs(-B/(2.*A));        614     return distMin = std::fabs(-B/(2.*A)); 
529   }                                               615   }
530                                                   616   
531   G4double plus  = (-B+std::sqrt(discr))/(2.*A    617   G4double plus  = (-B+std::sqrt(discr))/(2.*A);
532   G4double minus = (-B-std::sqrt(discr))/(2.*A    618   G4double minus = (-B-std::sqrt(discr))/(2.*A);
533                                                   619  
534   // Special case::Point on Surface, Check nor    620   // Special case::Point on Surface, Check norm.dot(v)
535                                                   621 
536   if ( ( std::fabs(plus) < halfCarTol )||( std << 622   if ( ( std::fabs(plus) < halfTol )||( std::fabs(minus) < halfTol ) )
537   {                                               623   {
538     G4ThreeVector truenorm(p.x()/(xSemiAxis*xS    624     G4ThreeVector truenorm(p.x()/(xSemiAxis*xSemiAxis),
539                            p.y()/(ySemiAxis*yS    625                            p.y()/(ySemiAxis*ySemiAxis),
540                            -( p.z() - zheight     626                            -( p.z() - zheight ));
541     if ( truenorm*v >= 0)  //  going outside t    627     if ( truenorm*v >= 0)  //  going outside the solid from surface
542     {                                             628     {
543       return kInfinity;                           629       return kInfinity;
544     }                                             630     }
545     else                                          631     else
546     {                                             632     {
547       return 0;                                   633       return 0;
548     }                                             634     }
549   }                                               635   }
550                                                   636 
551   // G4double lambda = std::fabs(plus) < std::    637   // G4double lambda = std::fabs(plus) < std::fabs(minus) ? plus : minus;  
552   G4double lambda = 0;                            638   G4double lambda = 0;
553                                                   639 
554   if ( minus > halfCarTol && minus < distMin ) << 640   if ( minus > halfTol && minus < distMin ) 
555   {                                               641   {
556     lambda = minus ;                              642     lambda = minus ;
557     // check normal vector   n * v < 0            643     // check normal vector   n * v < 0
558     G4ThreeVector pin = p + lambda*v;             644     G4ThreeVector pin = p + lambda*v;
559     if(std::fabs(pin.z())< zTopCut + halfCarTo << 645     if(std::fabs(pin.z())<zTopCut+0.5*kCarTolerance)
560     {                                             646     {
561       G4ThreeVector truenorm(pin.x()/(xSemiAxi    647       G4ThreeVector truenorm(pin.x()/(xSemiAxis*xSemiAxis),
562                              pin.y()/(ySemiAxi    648                              pin.y()/(ySemiAxis*ySemiAxis),
563                              - ( pin.z() - zhe    649                              - ( pin.z() - zheight ));
564       if ( truenorm*v < 0)                        650       if ( truenorm*v < 0)
565       {   // yes, going inside the solid          651       {   // yes, going inside the solid
566         distMin = lambda;                         652         distMin = lambda;
567       }                                           653       }
568     }                                             654     }
569   }                                               655   }
570   if ( plus > halfCarTol  && plus < distMin )  << 656   if ( plus > halfTol  && plus < distMin )
571   {                                               657   {
572     lambda = plus ;                               658     lambda = plus ;
573     // check normal vector   n * v < 0            659     // check normal vector   n * v < 0
574     G4ThreeVector pin = p + lambda*v;             660     G4ThreeVector pin = p + lambda*v;
575     if(std::fabs(pin.z()) < zTopCut + halfCarT << 661     if(std::fabs(pin.z())<zTopCut+0.5*kCarTolerance)
576     {                                             662     {
577       G4ThreeVector truenorm(pin.x()/(xSemiAxi    663       G4ThreeVector truenorm(pin.x()/(xSemiAxis*xSemiAxis),
578                              pin.y()/(ySemiAxi    664                              pin.y()/(ySemiAxis*ySemiAxis),
579                              - ( pin.z() - zhe    665                              - ( pin.z() - zheight ) );
580       if ( truenorm*v < 0)                        666       if ( truenorm*v < 0)
581       {   // yes, going inside the solid          667       {   // yes, going inside the solid
582         distMin = lambda;                         668         distMin = lambda;
583       }                                           669       }
584     }                                             670     }
585   }                                               671   }
586   if (distMin < halfCarTol) distMin=0.;        << 672   if (distMin < halfTol) distMin=0.;
587   return distMin ;                                673   return distMin ;
588 }                                                 674 }
589                                                   675 
590 ////////////////////////////////////////////// << 676 //////////////////////////////////////////////////////////////////////////
591 //                                                677 //
592 // Calculate distance (<= actual) to closest s    678 // Calculate distance (<= actual) to closest surface of shape from outside
593 // Return 0 if point inside                       679 // Return 0 if point inside
594                                                << 680 //
595 G4double G4EllipticalCone::DistanceToIn(const     681 G4double G4EllipticalCone::DistanceToIn(const G4ThreeVector& p) const
596 {                                                 682 {
597   G4double hp = std::sqrt(p.x()*p.x()*invXX +  << 683   G4double distR, distR2, distZ, maxDim;
598   G4double ds = (hp - zheight)*cosAxisMin;     << 684   G4double distRad;  
599   G4double dz = std::abs(p.z()) - zTopCut;     << 685 
600   G4double dist = std::max(ds,dz);             << 686   // check if the point lies either below z=-zTopCut in bottom elliptical
601   return (dist > 0) ? dist : 0.;               << 687   // region or on top within cut elliptical region
                                                   >> 688   //
                                                   >> 689   if( (p.z() <= -zTopCut) && (sqr(p.x()/xSemiAxis) + sqr(p.y()/ySemiAxis)
                                                   >> 690                            <= sqr(zTopCut + zheight + 0.5*kCarTolerance )) )
                                                   >> 691   {  
                                                   >> 692     //return distZ = std::fabs(zTopCut - p.z());
                                                   >> 693      return distZ = std::fabs(zTopCut + p.z());
                                                   >> 694   } 
                                                   >> 695   
                                                   >> 696   if( (p.z() >= zTopCut) && (sqr(p.x()/xSemiAxis)+sqr(p.y()/ySemiAxis)
                                                   >> 697                           <= sqr(zheight - zTopCut + kCarTolerance/2.0 )) )
                                                   >> 698   {
                                                   >> 699     return distZ = std::fabs(p.z() - zTopCut);
                                                   >> 700   } 
                                                   >> 701   
                                                   >> 702   // below we use the following approximation: we take the largest of the
                                                   >> 703   // axes and find the shortest distance to the circular (cut) cone of that
                                                   >> 704   // radius.  
                                                   >> 705   //
                                                   >> 706   maxDim = xSemiAxis >= ySemiAxis ? xSemiAxis:ySemiAxis;
                                                   >> 707   distRad = std::sqrt(p.x()*p.x()+p.y()*p.y());
                                                   >> 708 
                                                   >> 709   if( p.z() > maxDim*distRad + zTopCut*(1.+maxDim)-sqr(maxDim)*zheight )
                                                   >> 710   {
                                                   >> 711     distR2 = sqr(p.z() - zTopCut) + sqr(distRad - maxDim*(zheight - zTopCut));
                                                   >> 712     return std::sqrt( distR2 );
                                                   >> 713   } 
                                                   >> 714 
                                                   >> 715   if( distRad > maxDim*( zheight - p.z() ) )
                                                   >> 716   {
                                                   >> 717     if( p.z() > maxDim*distRad - (zTopCut*(1.+maxDim)+sqr(maxDim)*zheight) )
                                                   >> 718     {
                                                   >> 719       G4double zVal = (p.z()-maxDim*(distRad-maxDim*zheight))/(1.+sqr(maxDim));
                                                   >> 720       G4double rVal = maxDim*(zheight - zVal);
                                                   >> 721       return distR  = std::sqrt(sqr(p.z() - zVal) + sqr(distRad - rVal));
                                                   >> 722     }
                                                   >> 723   }
                                                   >> 724 
                                                   >> 725   if( distRad <= maxDim*(zheight - p.z()) )
                                                   >> 726   {
                                                   >> 727     distR2 = sqr(distRad - maxDim*(zheight + zTopCut)) + sqr(p.z() + zTopCut);
                                                   >> 728     return std::sqrt( distR2 );    
                                                   >> 729   }   
                                                   >> 730   
                                                   >> 731   return distR = 0;
602 }                                                 732 }
603                                                   733 
604 ////////////////////////////////////////////// << 734 /////////////////////////////////////////////////////////////////////////
605 //                                                735 //
606 // Calculate distance to surface of shape from    736 // Calculate distance to surface of shape from `inside',
607 // allowing for tolerance                         737 // allowing for tolerance
608                                                << 738 //
609 G4double G4EllipticalCone::DistanceToOut(const    739 G4double G4EllipticalCone::DistanceToOut(const G4ThreeVector& p,
610                                          const    740                                          const G4ThreeVector& v,
611                                          const    741                                          const G4bool calcNorm,
612                                                << 742                                                G4bool *validNorm,
613                                                << 743                                                G4ThreeVector *n  ) const
614 {                                                 744 {
615   G4double distMin, lambda;                       745   G4double distMin, lambda;
616   enum surface_e {kPlaneSurf, kCurvedSurf, kNo    746   enum surface_e {kPlaneSurf, kCurvedSurf, kNoSurf} surface;
617                                                   747   
618   distMin = kInfinity;                            748   distMin = kInfinity;
619   surface = kNoSurf;                              749   surface = kNoSurf;
620                                                   750 
621   if (v.z() < 0.0)                                751   if (v.z() < 0.0)
622   {                                               752   {
623     lambda = (-p.z() - zTopCut)/v.z();            753     lambda = (-p.z() - zTopCut)/v.z();
624                                                   754 
625     if ( (sqr((p.x() + lambda*v.x())/xSemiAxis    755     if ( (sqr((p.x() + lambda*v.x())/xSemiAxis) + 
626           sqr((p.y() + lambda*v.y())/ySemiAxis    756           sqr((p.y() + lambda*v.y())/ySemiAxis)) < 
627           sqr(zheight + zTopCut + halfCarTol)  << 757           sqr(zheight + zTopCut + 0.5*kCarTolerance) )
628     {                                             758     {
629       distMin = std::fabs(lambda);                759       distMin = std::fabs(lambda);
630                                                   760 
631       if (!calcNorm) { return distMin; }          761       if (!calcNorm) { return distMin; }
632     }                                             762     } 
633     distMin = std::fabs(lambda);                  763     distMin = std::fabs(lambda);
634     surface = kPlaneSurf;                         764     surface = kPlaneSurf;
635   }                                               765   }
636                                                   766 
637   if (v.z() > 0.0)                                767   if (v.z() > 0.0)
638   {                                               768   {
639     lambda = (zTopCut - p.z()) / v.z();           769     lambda = (zTopCut - p.z()) / v.z();
640                                                   770 
641     if ( (sqr((p.x() + lambda*v.x())/xSemiAxis    771     if ( (sqr((p.x() + lambda*v.x())/xSemiAxis)
642         + sqr((p.y() + lambda*v.y())/ySemiAxis    772         + sqr((p.y() + lambda*v.y())/ySemiAxis) )
643        < (sqr(zheight - zTopCut + halfCarTol)) << 773        < (sqr(zheight - zTopCut + 0.5*kCarTolerance)) )
644     {                                             774     {
645       distMin = std::fabs(lambda);                775       distMin = std::fabs(lambda);
646       if (!calcNorm) { return distMin; }          776       if (!calcNorm) { return distMin; }
647     }                                             777     }
648     distMin = std::fabs(lambda);                  778     distMin = std::fabs(lambda);
649     surface = kPlaneSurf;                         779     surface = kPlaneSurf;
650   }                                               780   }
651                                                   781   
652   // if we are here then it either intersects     782   // if we are here then it either intersects or grazes the 
653   // curved surface...                            783   // curved surface...
654   //                                              784   //
655   G4double A = sqr(v.x()/xSemiAxis) + sqr(v.y(    785   G4double A = sqr(v.x()/xSemiAxis) + sqr(v.y()/ySemiAxis) - sqr(v.z());
656   G4double B = 2.*(v.x()*p.x()/sqr(xSemiAxis)     786   G4double B = 2.*(v.x()*p.x()/sqr(xSemiAxis) +  
657                    v.y()*p.y()/sqr(ySemiAxis)     787                    v.y()*p.y()/sqr(ySemiAxis) + v.z()*(zheight-p.z()));
658   G4double C = sqr(p.x()/xSemiAxis) + sqr(p.y(    788   G4double C = sqr(p.x()/xSemiAxis) + sqr(p.y()/ySemiAxis)
659              - sqr(zheight - p.z());              789              - sqr(zheight - p.z());
660                                                   790  
661   G4double discr = B*B - 4.*A*C;                  791   G4double discr = B*B - 4.*A*C;
662                                                   792   
663   if ( discr >= - halfCarTol && discr < halfCa << 793   if ( discr >= - 0.5*kCarTolerance && discr < 0.5*kCarTolerance )
664   {                                               794   { 
665     if(!calcNorm) { return distMin = std::fabs    795     if(!calcNorm) { return distMin = std::fabs(-B/(2.*A)); }
666   }                                               796   }
667                                                   797 
668   else if ( discr > halfCarTol )               << 798   else if ( discr > 0.5*kCarTolerance )
669   {                                               799   {
670     G4double plus  = (-B+std::sqrt(discr))/(2.    800     G4double plus  = (-B+std::sqrt(discr))/(2.*A);
671     G4double minus = (-B-std::sqrt(discr))/(2.    801     G4double minus = (-B-std::sqrt(discr))/(2.*A);
672                                                   802 
673     if ( plus > halfCarTol && minus > halfCarT << 803     if ( plus > 0.5*kCarTolerance && minus > 0.5*kCarTolerance )
674     {                                             804     {
675       // take the shorter distance                805       // take the shorter distance
676       //                                          806       //
677       lambda   = std::fabs(plus) < std::fabs(m    807       lambda   = std::fabs(plus) < std::fabs(minus) ? plus : minus;
678     }                                             808     }
679     else                                          809     else
680     {                                             810     {
681       // at least one solution is close to zer    811       // at least one solution is close to zero or negative
682       // so, take small positive solution or z    812       // so, take small positive solution or zero 
683       //                                          813       //
684       lambda   = plus > -halfCarTol ? plus : 0 << 814       lambda   = plus > -0.5*kCarTolerance ? plus : 0;
685     }                                             815     }
686                                                   816 
687     if ( std::fabs(lambda) < distMin )            817     if ( std::fabs(lambda) < distMin )
688     {                                             818     {
689       if( std::fabs(lambda) > halfCarTol)      << 819       if( std::fabs(lambda) > 0.5*kCarTolerance)
690       {                                           820       {
691         distMin  = std::fabs(lambda);             821         distMin  = std::fabs(lambda);
692         surface  = kCurvedSurf;                   822         surface  = kCurvedSurf;
693       }                                           823       }
694       else  // Point is On the Surface, Check     824       else  // Point is On the Surface, Check Normal
695       {                                           825       {
696         G4ThreeVector truenorm(p.x()/(xSemiAxi    826         G4ThreeVector truenorm(p.x()/(xSemiAxis*xSemiAxis),
697                                p.y()/(ySemiAxi    827                                p.y()/(ySemiAxis*ySemiAxis),
698                                -( p.z() - zhei    828                                -( p.z() - zheight ));
699         if( truenorm.dot(v) > 0 )                 829         if( truenorm.dot(v) > 0 )
700         {                                         830         {
701           distMin  = 0.0;                         831           distMin  = 0.0;
702           surface  = kCurvedSurf;                 832           surface  = kCurvedSurf;
703         }                                         833         }
704       }                                           834       } 
705     }                                             835     }
706   }                                               836   }
707                                                   837 
708   // set normal if requested                      838   // set normal if requested
709   //                                              839   //
710   if (calcNorm)                                   840   if (calcNorm)
711   {                                               841   {
712     if (surface == kNoSurf)                       842     if (surface == kNoSurf)
713     {                                             843     {
714       *validNorm = false;                         844       *validNorm = false;
715     }                                             845     }
716     else                                          846     else
717     {                                             847     {
718       *validNorm = true;                          848       *validNorm = true;
719       switch (surface)                            849       switch (surface)
720       {                                           850       {
721         case kPlaneSurf:                          851         case kPlaneSurf:
722         {                                         852         {
723           *n = G4ThreeVector(0.,0.,(v.z() > 0.    853           *n = G4ThreeVector(0.,0.,(v.z() > 0.0 ? 1. : -1.));
724         }                                         854         }
725         break;                                    855         break;
726                                                   856 
727         case kCurvedSurf:                         857         case kCurvedSurf:
728         {                                         858         {
729           G4ThreeVector pexit = p + distMin*v;    859           G4ThreeVector pexit = p + distMin*v;
730           G4ThreeVector truenorm( pexit.x()/(x    860           G4ThreeVector truenorm( pexit.x()/(xSemiAxis*xSemiAxis),
731                                   pexit.y()/(y    861                                   pexit.y()/(ySemiAxis*ySemiAxis),
732                                   -( pexit.z()    862                                   -( pexit.z() - zheight ) );
733           truenorm /= truenorm.mag();             863           truenorm /= truenorm.mag();
734           *n= truenorm;                           864           *n= truenorm;
735         }                                         865         } 
736         break;                                    866         break;
737                                                   867 
738         default:            // Should never re    868         default:            // Should never reach this case ...
739           DumpInfo();                             869           DumpInfo();
740           std::ostringstream message;             870           std::ostringstream message;
741           G4long oldprc = message.precision(16 << 871           G4int oldprc = message.precision(16);
742           message << "Undefined side for valid    872           message << "Undefined side for valid surface normal to solid."
743                   << G4endl                       873                   << G4endl
744                   << "Position:"  << G4endl       874                   << "Position:"  << G4endl
745                   << "   p.x() = "   << p.x()/    875                   << "   p.x() = "   << p.x()/mm << " mm" << G4endl
746                   << "   p.y() = "   << p.y()/    876                   << "   p.y() = "   << p.y()/mm << " mm" << G4endl
747                   << "   p.z() = "   << p.z()/    877                   << "   p.z() = "   << p.z()/mm << " mm" << G4endl
748                   << "Direction:" << G4endl       878                   << "Direction:" << G4endl
749                   << "   v.x() = "   << v.x()     879                   << "   v.x() = "   << v.x() << G4endl
750                   << "   v.y() = "   << v.y()     880                   << "   v.y() = "   << v.y() << G4endl
751                   << "   v.z() = "   << v.z()     881                   << "   v.z() = "   << v.z() << G4endl
752                   << "Proposed distance :" <<     882                   << "Proposed distance :" << G4endl
753                   << "   distMin = "    << dis    883                   << "   distMin = "    << distMin/mm << " mm";
754           message.precision(oldprc);              884           message.precision(oldprc);
755           G4Exception("G4EllipticalCone::Dista    885           G4Exception("G4EllipticalCone::DistanceToOut(p,v,..)",
756                       "GeomSolids1002", JustWa    886                       "GeomSolids1002", JustWarning, message);
757           break;                                  887           break;
758       }                                           888       }
759     }                                             889     }
760   }                                               890   }
761                                                   891 
762   if (distMin < halfCarTol) { distMin=0; }     << 892   if (distMin<0.5*kCarTolerance) { distMin=0; }
763                                                   893 
764   return distMin;                                 894   return distMin;
765 }                                                 895 }
766                                                   896 
767 //////////////////////////////////////////////    897 /////////////////////////////////////////////////////////////////////////
768 //                                                898 //
769 // Calculate distance (<=actual) to closest su    899 // Calculate distance (<=actual) to closest surface of shape from inside
770                                                << 900 //
771 G4double G4EllipticalCone::DistanceToOut(const    901 G4double G4EllipticalCone::DistanceToOut(const G4ThreeVector& p) const
772 {                                                 902 {
                                                   >> 903   G4double rds,roo,roo1, distR, distZ, distMin=0.;
                                                   >> 904   G4double minAxis = xSemiAxis < ySemiAxis ? xSemiAxis : ySemiAxis;
                                                   >> 905 
773 #ifdef G4SPECSDEBUG                               906 #ifdef G4SPECSDEBUG
774   if( Inside(p) == kOutside )                     907   if( Inside(p) == kOutside )
775   {                                               908   {
                                                   >> 909      DumpInfo();
776      std::ostringstream message;                  910      std::ostringstream message;
777      G4long oldprc = message.precision(16);    << 911      G4int oldprc = message.precision(16);
778      message << "Point p is outside (!?) of so << 912      message << "Point p is outside !?" << G4endl
779              << "Position:\n"                  << 913              << "Position:"  << G4endl
780              << "   p.x() = "  << p.x()/mm <<  << 914              << "   p.x() = "   << p.x()/mm << " mm" << G4endl
781              << "   p.y() = "  << p.y()/mm <<  << 915              << "   p.y() = "   << p.y()/mm << " mm" << G4endl
782              << "   p.z() = "  << p.z()/mm <<  << 916              << "   p.z() = "   << p.z()/mm << " mm";
783      message.precision(oldprc) ;                  917      message.precision(oldprc) ;
784      G4Exception("G4Ellipsoid::DistanceToOut(p    918      G4Exception("G4Ellipsoid::DistanceToOut(p)", "GeomSolids1002",
785                  JustWarning, message);           919                  JustWarning, message);
786      DumpInfo();                               << 
787   }                                               920   }
788 #endif                                            921 #endif
789   G4double hp = std::sqrt(p.x()*p.x()*invXX +  << 922     
790   G4double ds = (zheight - hp)*cosAxisMin;     << 923   // since we have made the above warning, below we are working assuming p
791   G4double dz = zTopCut - std::abs(p.z());     << 924   // is inside check how close it is to the circular cone with radius equal
792   G4double dist = std::min(ds,dz);             << 925   // to the smaller of the axes
793   return (dist > 0) ? dist : 0.;               << 926   //
                                                   >> 927   if( sqr(p.x()/minAxis)+sqr(p.y()/minAxis) < sqr(zheight - p.z()) )
                                                   >> 928   {
                                                   >> 929     rds     = std::sqrt(sqr(p.x()) + sqr(p.y()));
                                                   >> 930     roo     = minAxis*(zheight-p.z()); // radius of cone at z= p.z()
                                                   >> 931     roo1    = minAxis*(zheight-zTopCut); // radius of cone at z=+zTopCut
                                                   >> 932 
                                                   >> 933     distZ=zTopCut - std::fabs(p.z()) ;
                                                   >> 934     distR=(roo-rds)/(std::sqrt(1+sqr(minAxis)));
                                                   >> 935 
                                                   >> 936     if(rds>roo1)
                                                   >> 937     {
                                                   >> 938       distMin=(zTopCut-p.z())*(roo-rds)/(roo-roo1);
                                                   >> 939       distMin=std::min(distMin,distR);
                                                   >> 940     }      
                                                   >> 941     distMin=std::min(distR,distZ);
                                                   >> 942   }
                                                   >> 943 
                                                   >> 944   return distMin;
794 }                                                 945 }
795                                                   946 
796 ////////////////////////////////////////////// << 947 //////////////////////////////////////////////////////////////////////////
797 //                                                948 //
798 // GetEntityType                                  949 // GetEntityType
799                                                << 950 //
800 G4GeometryType G4EllipticalCone::GetEntityType    951 G4GeometryType G4EllipticalCone::GetEntityType() const
801 {                                                 952 {
802   return {"G4EllipticalCone"};                 << 953   return G4String("G4EllipticalCone");
803 }                                                 954 }
804                                                   955 
805 ////////////////////////////////////////////// << 956 //////////////////////////////////////////////////////////////////////////
806 //                                                957 //
807 // Make a clone of the object                     958 // Make a clone of the object
808                                                << 959 //
809 G4VSolid* G4EllipticalCone::Clone() const         960 G4VSolid* G4EllipticalCone::Clone() const
810 {                                                 961 {
811   return new G4EllipticalCone(*this);             962   return new G4EllipticalCone(*this);
812 }                                                 963 }
813                                                   964 
814 ////////////////////////////////////////////// << 965 //////////////////////////////////////////////////////////////////////////
815 //                                                966 //
816 // Stream object contents to an output stream     967 // Stream object contents to an output stream
817                                                << 968 //
818 std::ostream& G4EllipticalCone::StreamInfo( st    969 std::ostream& G4EllipticalCone::StreamInfo( std::ostream& os ) const
819 {                                                 970 {
820   G4long oldprc = os.precision(16);            << 971   G4int oldprc = os.precision(16);
821   os << "-------------------------------------    972   os << "-----------------------------------------------------------\n"
822      << "    *** Dump for solid - " << GetName    973      << "    *** Dump for solid - " << GetName() << " ***\n"
823      << "    =================================    974      << "    ===================================================\n"
824      << " Solid type: G4EllipticalCone\n"         975      << " Solid type: G4EllipticalCone\n"
825      << " Parameters: \n"                         976      << " Parameters: \n"
826                                                   977 
827      << "    semi-axis x: " << xSemiAxis/mm <<    978      << "    semi-axis x: " << xSemiAxis/mm << " mm \n"
828      << "    semi-axis y: " << ySemiAxis/mm <<    979      << "    semi-axis y: " << ySemiAxis/mm << " mm \n"
829      << "    height    z: " << zheight/mm << "    980      << "    height    z: " << zheight/mm << " mm \n"
830      << "    half length in  z: " << zTopCut/m    981      << "    half length in  z: " << zTopCut/mm << " mm \n"
831      << "-------------------------------------    982      << "-----------------------------------------------------------\n";
832   os.precision(oldprc);                           983   os.precision(oldprc);
833                                                   984 
834   return os;                                      985   return os;
835 }                                                 986 }
836                                                   987 
837 //////////////////////////////////////////////    988 /////////////////////////////////////////////////////////////////////////
838 //                                                989 //
839 // Return random point on the surface of the s << 990 // GetPointOnSurface
840                                                << 991 //
                                                   >> 992 // returns quasi-uniformly distributed point on surface of elliptical cone
                                                   >> 993 //
841 G4ThreeVector G4EllipticalCone::GetPointOnSurf    994 G4ThreeVector G4EllipticalCone::GetPointOnSurface() const
842 {                                                 995 {
843   G4double x0 = xSemiAxis*zheight; // x semi a << 
844   G4double y0 = ySemiAxis*zheight; // y semi a << 
845   G4double s0 = G4GeomTools::EllipticConeLater << 
846   G4double kmin = (zTopCut >= zheight ) ? 0. : << 
847   G4double kmax = (zTopCut >= zheight ) ? 2. : << 
848                                                << 
849   // Set areas (base at -Z, side surface, base << 
850   //                                           << 
851   G4double szmin =  pi*x0*y0*kmax*kmax;        << 
852   G4double szmax =  pi*x0*y0*kmin*kmin;        << 
853   G4double sside =  s0*(kmax*kmax - kmin*kmin) << 
854   G4double ssurf[3] = { szmin, sside, szmax }; << 
855   for (auto i=1; i<3; ++i) { ssurf[i] += ssurf << 
856                                                   996 
857   // Select surface                            << 997   G4double phi, sinphi, cosphi, aOne, aTwo, aThree,
858   //                                           << 998            chose, zRand, rRand1, rRand2;
859   G4double select = ssurf[2]*G4UniformRand();  << 999   
860   G4int k = 2;                                 << 1000   G4double rOne = std::sqrt(sqr(xSemiAxis)
861   if (select <= ssurf[1]) k = 1;               << 1001                 + sqr(ySemiAxis))*(zheight - zTopCut);
862   if (select <= ssurf[0]) k = 0;               << 1002   G4double rTwo = std::sqrt(sqr(xSemiAxis)
                                                   >> 1003                 + sqr(ySemiAxis))*(zheight + zTopCut);
                                                   >> 1004   
                                                   >> 1005   aOne   = pi*(rOne + rTwo)*std::sqrt(sqr(rOne - rTwo)+sqr(2.*zTopCut));
                                                   >> 1006   aTwo   = pi*xSemiAxis*ySemiAxis*sqr(zheight+zTopCut);
                                                   >> 1007   aThree = pi*xSemiAxis*ySemiAxis*sqr(zheight-zTopCut);  
                                                   >> 1008 
                                                   >> 1009   phi = RandFlat::shoot(0.,twopi);
                                                   >> 1010   cosphi = std::cos(phi);
                                                   >> 1011   sinphi = std::sin(phi);
                                                   >> 1012   
                                                   >> 1013   if(zTopCut >= zheight) aThree = 0.;
863                                                   1014 
864   // Pick random point on selected surface     << 1015   chose = RandFlat::shoot(0.,aOne+aTwo+aThree);
865   //                                           << 1016   if((chose>=0.) && (chose<aOne))
866   G4ThreeVector p;                             << 
867   switch(k)                                    << 
868   {                                               1017   {
869     case 0: // base at -Z, uniform distributio << 1018     zRand = RandFlat::shoot(-zTopCut,zTopCut);
870     {                                          << 1019     return G4ThreeVector(xSemiAxis*(zheight-zRand)*cosphi,
871       G4double zh = zheight + zTopCut;         << 1020                          ySemiAxis*(zheight-zRand)*sinphi,zRand);    
872       G4TwoVector rho = G4RandomPointInEllipse << 1021   }
873       p.set(rho.x(),rho.y(),-zTopCut);         << 1022   else if((chose>=aOne) && (chose<aOne+aTwo))
874       break;                                   << 1023   {
875     }                                          << 1024     do
876     case 1: // side surface, uniform distribut << 
877     {                                             1025     {
878       G4double zh = G4RandomRadiusInRing(zheig << 1026       rRand1 = RandFlat::shoot(0.,1.) ;
879       G4double a = x0;                         << 1027       rRand2 = RandFlat::shoot(0.,1.) ;
880       G4double b = y0;                         << 1028     } while ( rRand2 >= rRand1  ) ;
881                                                   1029 
882       G4double hh = zheight*zheight;           << 1030     //    rRand2 = RandFlat::shoot(0.,std::sqrt(1.-sqr(rRand1)));
883       G4double aa = a*a;                       << 1031     return G4ThreeVector(rRand1*xSemiAxis*(zheight+zTopCut)*cosphi,
884       G4double bb = b*b;                       << 1032                          rRand1*ySemiAxis*(zheight+zTopCut)*sinphi, -zTopCut);
885       G4double R  = std::max(a,b);             << 
886       G4double mu_max = R*std::sqrt(hh + R*R); << 
887                                                   1033 
888       G4double x,y;                            << 
889       for (auto i=0; i<1000; ++i)              << 
890       {                                        << 
891   G4double phi = CLHEP::twopi*G4UniformRand(); << 
892         x = std::cos(phi);                     << 
893         y = std::sin(phi);                     << 
894         G4double xx = x*x;                     << 
895         G4double yy = y*y;                     << 
896         G4double E = hh + aa*xx + bb*yy;       << 
897         G4double F = (aa-bb)*x*y;              << 
898         G4double G = aa*yy + bb*xx;            << 
899         G4double mu = std::sqrt(E*G - F*F);    << 
900         if (mu_max*G4UniformRand() <= mu) brea << 
901       }                                        << 
902       p.set(zh*xSemiAxis*x,zh*ySemiAxis*y,zhei << 
903       break;                                   << 
904     }                                          << 
905     case 2: // base at +Z, uniform distributio << 
906     {                                          << 
907       G4double zh = zheight - zTopCut;         << 
908       G4TwoVector rho = G4RandomPointInEllipse << 
909       p.set(rho.x(),rho.y(),zTopCut);          << 
910       break;                                   << 
911     }                                          << 
912   }                                               1034   }
913   return p;                                    << 1035   // else
914 }                                              << 1036   //
915                                                << 
916 ////////////////////////////////////////////// << 
917 //                                             << 
918 // Get cubic volume                            << 
919                                                   1037 
920 G4double G4EllipticalCone::GetCubicVolume()    << 1038   do
921 {                                              << 
922   if (fCubicVolume == 0.0)                     << 
923   {                                               1039   {
924     G4double x0 = xSemiAxis*zheight; // x semi << 1040     rRand1 = RandFlat::shoot(0.,1.) ;
925     G4double y0 = ySemiAxis*zheight; // y semi << 1041     rRand2 = RandFlat::shoot(0.,1.) ;
926     G4double v0 = CLHEP::pi*x0*y0*zheight/3.;  << 1042   } while ( rRand2 >= rRand1  ) ;
927     G4double kmin = (zTopCut >= zheight ) ? 0. << 
928     G4double kmax = (zTopCut >= zheight ) ? 2. << 
929     fCubicVolume = (kmax - kmin)*(kmax*kmax +  << 
930   }                                            << 
931   return fCubicVolume;                         << 
932 }                                              << 
933                                                   1043 
934 ////////////////////////////////////////////// << 1044   return G4ThreeVector(rRand1*xSemiAxis*(zheight-zTopCut)*cosphi,
935 //                                             << 1045                        rRand1*ySemiAxis*(zheight-zTopCut)*sinphi, zTopCut);
936 // Get surface area                            << 
937                                                << 
938 G4double G4EllipticalCone::GetSurfaceArea()    << 
939 {                                              << 
940   if (fSurfaceArea == 0.0)                     << 
941   {                                            << 
942     G4double x0 = xSemiAxis*zheight; // x semi << 
943     G4double y0 = ySemiAxis*zheight; // y semi << 
944     G4double s0 = G4GeomTools::EllipticConeLat << 
945     G4double kmin = (zTopCut >= zheight ) ? 0. << 
946     G4double kmax = (zTopCut >= zheight ) ? 2. << 
947     fSurfaceArea = (kmax - kmin)*(kmax + kmin) << 
948                  + CLHEP::pi*x0*y0*(kmin*kmin  << 
949   }                                            << 
950   return fSurfaceArea;                         << 
951 }                                                 1046 }
952                                                   1047 
953 ////////////////////////////////////////////// << 
954 //                                                1048 //
955 // Methods for visualisation                      1049 // Methods for visualisation
                                                   >> 1050 //
956                                                   1051 
957 void G4EllipticalCone::DescribeYourselfTo (G4V    1052 void G4EllipticalCone::DescribeYourselfTo (G4VGraphicsScene& scene) const
958 {                                                 1053 {
959   scene.AddSolid(*this);                          1054   scene.AddSolid(*this);
960 }                                                 1055 }
961                                                   1056 
962 G4VisExtent G4EllipticalCone::GetExtent() cons    1057 G4VisExtent G4EllipticalCone::GetExtent() const
963 {                                                 1058 {
964   // Define the sides of the box into which th    1059   // Define the sides of the box into which the solid instance would fit.
965   //                                              1060   //
966   G4ThreeVector pmin,pmax;                     << 1061   G4double maxDim;
967   BoundingLimits(pmin,pmax);                   << 1062   maxDim = xSemiAxis > ySemiAxis ? xSemiAxis : ySemiAxis;
968   return { pmin.x(), pmax.x(), pmin.y(), pmax. << 1063   maxDim = maxDim > zTopCut ? maxDim : zTopCut;
                                                   >> 1064   
                                                   >> 1065   return G4VisExtent (-maxDim, maxDim,
                                                   >> 1066                       -maxDim, maxDim,
                                                   >> 1067                       -maxDim, maxDim);
                                                   >> 1068 }
                                                   >> 1069 
                                                   >> 1070 G4NURBS* G4EllipticalCone::CreateNURBS () const
                                                   >> 1071 {
                                                   >> 1072   // Box for now!!!
                                                   >> 1073   //
                                                   >> 1074   return new G4NURBSbox(xSemiAxis, ySemiAxis,zheight);
969 }                                                 1075 }
970                                                   1076 
971 G4Polyhedron* G4EllipticalCone::CreatePolyhedr    1077 G4Polyhedron* G4EllipticalCone::CreatePolyhedron () const
972 {                                                 1078 {
973   return new G4PolyhedronEllipticalCone(xSemiA    1079   return new G4PolyhedronEllipticalCone(xSemiAxis, ySemiAxis, zheight, zTopCut);
974 }                                                 1080 }
975                                                   1081 
976 G4Polyhedron* G4EllipticalCone::GetPolyhedron     1082 G4Polyhedron* G4EllipticalCone::GetPolyhedron () const
977 {                                                 1083 {
978   if ( (fpPolyhedron == nullptr)               << 1084   if ( (!fpPolyhedron)
979     || fRebuildPolyhedron                      << 
980     || (fpPolyhedron->GetNumberOfRotationSteps    1085     || (fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
981         fpPolyhedron->GetNumberOfRotationSteps    1086         fpPolyhedron->GetNumberOfRotationSteps()) )
982     {                                             1087     {
983       G4AutoLock l(&polyhedronMutex);          << 
984       delete fpPolyhedron;                        1088       delete fpPolyhedron;
985       fpPolyhedron = CreatePolyhedron();          1089       fpPolyhedron = CreatePolyhedron();
986       fRebuildPolyhedron = false;              << 
987       l.unlock();                              << 
988     }                                             1090     }
989   return fpPolyhedron;                            1091   return fpPolyhedron;
990 }                                                 1092 }
991                                                << 
992 #endif // !defined(G4GEOM_USE_UELLIPTICALCONE) << 
993                                                   1093