Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/solids/CSG/src/G4Torus.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 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 // G4Torus implementation
 27 //
 28 // 30.10.96 V.Grichine: first implementation with G4Tubs elements in Fs
 29 // 26.05.00 V.Grichine: added new fuctions developed by O.Cremonesi
 30 // 31.08.00 E.Medernach: numerical computation of roots with bounding volume
 31 // 11.01.01 E.Medernach: Use G4PolynomialSolver to find roots
 32 // 03.05.05 V.Grichine: SurfaceNormal(p) according to J. Apostolakis proposal
 33 // 25.08.05 O.Link: new methods for DistanceToIn/Out using JTPolynomialSolver
 34 // 28.10.16 E.Tcherniaev: new CalculateExtent(); removed CreateRotatedVertices()
 35 // 16.12.16 H.Burkhardt: use radius differences and hypot to improve precision
 36 // --------------------------------------------------------------------
 37 
 38 #include "G4Torus.hh"
 39 
 40 #if !(defined(G4GEOM_USE_UTORUS) && defined(G4GEOM_USE_SYS_USOLIDS))
 41 
 42 #include "G4GeomTools.hh"
 43 #include "G4VoxelLimits.hh"
 44 #include "G4AffineTransform.hh"
 45 #include "G4BoundingEnvelope.hh"
 46 #include "G4GeometryTolerance.hh"
 47 #include "G4JTPolynomialSolver.hh"
 48 
 49 #include "G4VPVParameterisation.hh"
 50 
 51 #include "meshdefs.hh"
 52 
 53 #include "Randomize.hh"
 54 
 55 #include "G4VGraphicsScene.hh"
 56 #include "G4Polyhedron.hh"
 57 
 58 using namespace CLHEP;
 59 
 60 ///////////////////////////////////////////////////////////////
 61 //
 62 // Constructor - check parameters, convert angles so 0<sphi+dpshi<=2_PI
 63 //             - note if pdphi>2PI then reset to 2PI
 64 
 65 G4Torus::G4Torus( const G4String& pName,
 66                         G4double pRmin,
 67                         G4double pRmax,
 68                         G4double pRtor,
 69                         G4double pSPhi,
 70                         G4double pDPhi )
 71   : G4CSGSolid(pName)
 72 {
 73   SetAllParameters(pRmin, pRmax, pRtor, pSPhi, pDPhi);
 74 }
 75 
 76 ////////////////////////////////////////////////////////////////////////////
 77 //
 78 //
 79 
 80 void
 81 G4Torus::SetAllParameters( G4double pRmin,
 82                            G4double pRmax,
 83                            G4double pRtor,
 84                            G4double pSPhi,
 85                            G4double pDPhi )
 86 {
 87   const G4double fEpsilon = 4.e-11;  // relative tolerance of radii
 88 
 89   fCubicVolume = 0.;
 90   fSurfaceArea = 0.;
 91   fRebuildPolyhedron = true;
 92 
 93   kRadTolerance = G4GeometryTolerance::GetInstance()->GetRadialTolerance();
 94   kAngTolerance = G4GeometryTolerance::GetInstance()->GetAngularTolerance();
 95 
 96   halfCarTolerance = 0.5*kCarTolerance;
 97   halfAngTolerance = 0.5*kAngTolerance;
 98 
 99   if ( pRtor >= pRmax+1.e3*kCarTolerance )  // Check swept radius, as in G4Cons
100   {
101     fRtor = pRtor ;
102   }
103   else
104   {
105     std::ostringstream message;
106     message << "Invalid swept radius for Solid: " << GetName() << G4endl
107             << "        pRtor = " << pRtor << ", pRmax = " << pRmax;
108     G4Exception("G4Torus::SetAllParameters()",
109                 "GeomSolids0002", FatalException, message);
110   }
111 
112   // Check radii, as in G4Cons
113   //
114   if ( pRmin < pRmax - 1.e2*kCarTolerance && pRmin >= 0 )
115   {
116     if (pRmin >= 1.e2*kCarTolerance) { fRmin = pRmin ; }
117     else                             { fRmin = 0.0   ; }
118     fRmax = pRmax ;
119   }
120   else
121   {
122     std::ostringstream message;
123     message << "Invalid values of radii for Solid: " << GetName() << G4endl
124             << "        pRmin = " << pRmin << ", pRmax = " << pRmax;
125     G4Exception("G4Torus::SetAllParameters()",
126                 "GeomSolids0002", FatalException, message);
127   }
128 
129   // Relative tolerances
130   //
131   fRminTolerance = (fRmin) != 0.0
132                  ? 0.5*std::max( kRadTolerance, fEpsilon*(fRtor-fRmin )) : 0;
133   fRmaxTolerance = 0.5*std::max( kRadTolerance, fEpsilon*(fRtor+fRmax) );
134 
135   // Check angles
136   //
137   if ( pDPhi >= twopi )  { fDPhi = twopi ; }
138   else
139   {
140     if (pDPhi > 0)       { fDPhi = pDPhi ; }
141     else
142     {
143       std::ostringstream message;
144       message << "Invalid Z delta-Phi for Solid: " << GetName() << G4endl
145               << "        pDPhi = " << pDPhi;
146       G4Exception("G4Torus::SetAllParameters()",
147                   "GeomSolids0002", FatalException, message);
148     }
149   }
150   
151   // Ensure psphi in 0-2PI or -2PI-0 range if shape crosses 0
152   //
153   fSPhi = pSPhi;
154 
155   if (fSPhi < 0)  { fSPhi = twopi-std::fmod(std::fabs(fSPhi),twopi) ; }
156   else            { fSPhi = std::fmod(fSPhi,twopi) ; }
157 
158   if (fSPhi+fDPhi > twopi)  { fSPhi-=twopi ; }
159 }
160 
161 ///////////////////////////////////////////////////////////////////////
162 //
163 // Fake default constructor - sets only member data and allocates memory
164 //                            for usage restricted to object persistency.
165 //
166 G4Torus::G4Torus( __void__& a )
167   : G4CSGSolid(a)
168 {
169 }
170 
171 //////////////////////////////////////////////////////////////////////
172 //
173 // Destructor
174 
175 G4Torus::~G4Torus() = default;
176 
177 //////////////////////////////////////////////////////////////////////////
178 //
179 // Copy constructor
180 
181 G4Torus::G4Torus(const G4Torus&) = default;
182 
183 //////////////////////////////////////////////////////////////////////////
184 //
185 // Assignment operator
186 
187 G4Torus& G4Torus::operator = (const G4Torus& rhs) 
188 {
189    // Check assignment to self
190    //
191    if (this == &rhs)  { return *this; }
192 
193    // Copy base class data
194    //
195    G4CSGSolid::operator=(rhs);
196 
197    // Copy data
198    //
199    fRmin = rhs.fRmin; fRmax = rhs.fRmax;
200    fRtor = rhs.fRtor; fSPhi = rhs.fSPhi; fDPhi = rhs.fDPhi;
201    fRminTolerance = rhs.fRminTolerance; fRmaxTolerance = rhs.fRmaxTolerance;
202    kRadTolerance = rhs.kRadTolerance; kAngTolerance = rhs.kAngTolerance;
203    halfCarTolerance = rhs.halfCarTolerance;
204    halfAngTolerance = rhs.halfAngTolerance;
205 
206    return *this;
207 }
208 
209 //////////////////////////////////////////////////////////////////////
210 //
211 // Dispatch to parameterisation for replication mechanism dimension
212 // computation & modification.
213 
214 void G4Torus::ComputeDimensions(       G4VPVParameterisation* p,
215                                  const G4int n,
216                                  const G4VPhysicalVolume* pRep )
217 {
218   p->ComputeDimensions(*this,n,pRep);
219 }
220 
221 
222 
223 ////////////////////////////////////////////////////////////////////////////////
224 //
225 // Calculate the real roots to torus surface. 
226 // Returns negative solutions as well.
227 
228 void G4Torus::TorusRootsJT( const G4ThreeVector& p,
229                             const G4ThreeVector& v,
230                                   G4double r,
231                                   std::vector<G4double>& roots ) const
232 {
233 
234   G4int i, num ;
235   G4double c[5], srd[4], si[4] ;
236 
237   G4double Rtor2 = fRtor*fRtor, r2 = r*r  ;
238 
239   G4double pDotV = p.x()*v.x() + p.y()*v.y() + p.z()*v.z() ;
240   G4double pRad2 = p.x()*p.x() + p.y()*p.y() + p.z()*p.z() ;
241 
242   G4double d=pRad2 - Rtor2;
243   c[0] = 1.0 ;
244   c[1] = 4*pDotV ;
245   c[2] = 2*( (d + 2*pDotV*pDotV  - r2) + 2*Rtor2*v.z()*v.z());
246   c[3] = 4*(pDotV*(d - r2) + 2*Rtor2*p.z()*v.z()) ;
247   c[4] = (d-r2)*(d-r2) +4*Rtor2*(p.z()*p.z()-r2);
248 
249   G4JTPolynomialSolver  torusEq;
250 
251   num = torusEq.FindRoots( c, 4, srd, si );
252   
253   for ( i = 0; i < num; ++i ) 
254   {
255     if( si[i] == 0. )  { roots.push_back(srd[i]) ; }  // store real roots
256   }  
257 
258   std::sort(roots.begin() , roots.end() ) ;  // sorting  with <
259 }
260 
261 //////////////////////////////////////////////////////////////////////////////
262 //
263 // Interface for DistanceToIn and DistanceToOut.
264 // Calls TorusRootsJT and returns the smalles possible distance to 
265 // the surface.
266 // Attention: Difference in DistanceToIn/Out for points p on the surface.
267 
268 G4double G4Torus::SolveNumericJT( const G4ThreeVector& p,
269                                   const G4ThreeVector& v,
270                                         G4double r,
271                                         G4bool IsDistanceToIn ) const
272 {
273   G4double bigdist = 10*mm ;
274   G4double tmin = kInfinity ;
275   G4double t, scal ;
276 
277   // calculate the distances to the intersections with the Torus
278   // from a given point p and direction v.
279   //
280   std::vector<G4double> roots ;
281   std::vector<G4double> rootsrefined ;
282   TorusRootsJT(p,v,r,roots) ;
283 
284   G4ThreeVector ptmp ;
285 
286   // determine the smallest non-negative solution
287   //
288   for ( std::size_t k = 0 ; k<roots.size() ; ++k )
289   {
290     t = roots[k] ;
291 
292     if ( t < -halfCarTolerance )  { continue ; }  // skip negative roots
293 
294     if ( t > bigdist && t<kInfinity )    // problem with big distances
295     {
296       ptmp = p + t*v ;
297       TorusRootsJT(ptmp,v,r,rootsrefined) ;
298       if ( rootsrefined.size()==roots.size() )
299       {
300         t = t + rootsrefined[k] ;
301       }
302     }
303 
304     ptmp = p + t*v ;   // calculate the position of the proposed intersection
305 
306     G4double theta = std::atan2(ptmp.y(),ptmp.x());
307     
308     if ( fSPhi >= 0 )
309     {
310       if ( theta < - halfAngTolerance )  { theta += twopi; }
311       if ( (std::fabs(theta) < halfAngTolerance)
312         && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
313       { 
314         theta += twopi ; // 0 <= theta < 2pi
315       }
316     }
317     if ((fSPhi <= -pi )&&(theta>halfAngTolerance)) { theta = theta-twopi; }
318        
319     // We have to verify if this root is inside the region between
320     // fSPhi and fSPhi + fDPhi
321     //
322     if ( (theta - fSPhi >= - halfAngTolerance)
323       && (theta - (fSPhi + fDPhi) <=  halfAngTolerance) )
324     {
325       // check if P is on the surface, and called from DistanceToIn
326       // DistanceToIn has to return 0.0 if particle is going inside the solid
327 
328       if ( IsDistanceToIn )
329       {
330         if (std::fabs(t) < halfCarTolerance )
331         {
332           // compute scalar product at position p : v.n
333           // ( n taken from SurfaceNormal, not normalized )
334 
335           scal = v* G4ThreeVector( p.x()*(1-fRtor/std::hypot(p.x(),p.y())),
336                                    p.y()*(1-fRtor/std::hypot(p.x(),p.y())),
337                                    p.z() );
338 
339           // change sign in case of inner radius
340           //
341           if ( r == GetRmin() )  { scal = -scal ; }
342           if ( scal < 0 )  { return 0.0  ; }
343         }
344       }
345 
346       // check if P is on the surface, and called from DistanceToOut
347       // DistanceToIn has to return 0.0 if particle is leaving the solid
348 
349       if ( !IsDistanceToIn )
350       {
351         if (std::fabs(t) < halfCarTolerance )
352         {
353           // compute scalar product at position p : v.n   
354           //
355           scal = v* G4ThreeVector( p.x()*(1-fRtor/std::hypot(p.x(),p.y())),
356                                    p.y()*(1-fRtor/std::hypot(p.x(),p.y())),
357                                    p.z() );
358 
359           // change sign in case of inner radius
360           //
361           if ( r == GetRmin() )  { scal = -scal ; }
362           if ( scal > 0 )  { return 0.0  ; }
363         }
364       }
365 
366       // check if distance is larger than 1/2 kCarTolerance
367       //
368       if(  t > halfCarTolerance  )
369       {
370         tmin = t  ;
371         return tmin  ;
372       }
373     }
374   }
375 
376   return tmin;
377 }
378 
379 /////////////////////////////////////////////////////////////////////////////
380 //
381 // Get bounding box
382 
383 void G4Torus::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
384 {
385   G4double rmax = GetRmax();
386   G4double rtor = GetRtor();
387   G4double rint = rtor - rmax;
388   G4double rext = rtor + rmax;
389   G4double dz   = rmax;
390 
391   // Find bounding box
392   //
393   if (GetDPhi() >= twopi)
394   {
395     pMin.set(-rext,-rext,-dz);
396     pMax.set( rext, rext, dz);
397   }
398   else
399   {
400     G4TwoVector vmin,vmax;
401     G4GeomTools::DiskExtent(rint,rext,
402                             GetSinStartPhi(),GetCosStartPhi(),
403                             GetSinEndPhi(),GetCosEndPhi(),
404                             vmin,vmax);
405     pMin.set(vmin.x(),vmin.y(),-dz);
406     pMax.set(vmax.x(),vmax.y(), dz);
407   }
408 
409   // Check correctness of the bounding box
410   //
411   if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
412   {
413     std::ostringstream message;
414     message << "Bad bounding box (min >= max) for solid: "
415             << GetName() << " !"
416             << "\npMin = " << pMin
417             << "\npMax = " << pMax;
418     G4Exception("G4Torus::BoundingLimits()", "GeomMgt0001",
419                 JustWarning, message);
420     DumpInfo();
421   }
422 }
423 
424 /////////////////////////////////////////////////////////////////////////////
425 //
426 // Calculate extent under transform and specified limit
427 
428 G4bool G4Torus::CalculateExtent( const EAxis pAxis,
429                                  const G4VoxelLimits& pVoxelLimit,
430                                  const G4AffineTransform& pTransform,
431                                        G4double& pMin, G4double& pMax) const
432 {
433   G4ThreeVector bmin, bmax;
434   G4bool exist;
435 
436   // Get bounding box
437   BoundingLimits(bmin,bmax);
438 
439   // Check bounding box
440   G4BoundingEnvelope bbox(bmin,bmax);
441 #ifdef G4BBOX_EXTENT
442   return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
443 #endif
444   if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
445   {
446     return exist = pMin < pMax;
447   }
448 
449   // Get parameters of the solid
450   G4double rmin = GetRmin();
451   G4double rmax = GetRmax();
452   G4double rtor = GetRtor();
453   G4double dphi = GetDPhi();
454   G4double sinStart = GetSinStartPhi();
455   G4double cosStart = GetCosStartPhi();
456   G4double sinEnd   = GetSinEndPhi();
457   G4double cosEnd   = GetCosEndPhi();
458   G4double rint = rtor - rmax;
459   G4double rext = rtor + rmax;
460 
461   // Find bounding envelope and calculate extent
462   //
463   static const G4int NPHI  = 24; // number of steps for whole torus
464   static const G4int NDISK = 16; // number of steps for disk
465   static const G4double sinHalfDisk = std::sin(pi/NDISK);
466   static const G4double cosHalfDisk = std::cos(pi/NDISK);
467   static const G4double sinStepDisk = 2.*sinHalfDisk*cosHalfDisk;
468   static const G4double cosStepDisk = 1. - 2.*sinHalfDisk*sinHalfDisk;
469 
470   G4double astep = (360/NPHI)*deg; // max angle for one slice in phi
471   G4int    kphi  = (dphi <= astep) ? 1 : (G4int)((dphi-deg)/astep) + 1;
472   G4double ang   = dphi/kphi;
473 
474   G4double sinHalf = std::sin(0.5*ang);
475   G4double cosHalf = std::cos(0.5*ang);
476   G4double sinStep = 2.*sinHalf*cosHalf;
477   G4double cosStep = 1. - 2.*sinHalf*sinHalf;
478 
479   // define vectors for bounding envelope
480   G4ThreeVectorList pols[NDISK+1];
481   for (auto & pol : pols) pol.resize(4);
482 
483   std::vector<const G4ThreeVectorList *> polygons;
484   polygons.resize(NDISK+1);
485   for (G4int k=0; k<NDISK+1; ++k) polygons[k] = &pols[k];
486 
487   // set internal and external reference circles
488   G4TwoVector rzmin[NDISK];
489   G4TwoVector rzmax[NDISK];
490 
491   if ((rtor-rmin*sinHalfDisk)/cosHalf > (rtor+rmin*sinHalfDisk)) rmin = 0;
492   rmax /= cosHalfDisk;
493   G4double sinCurDisk = sinHalfDisk;
494   G4double cosCurDisk = cosHalfDisk;
495   for (G4int k=0; k<NDISK; ++k)
496   {
497     G4double rmincur = rtor + rmin*cosCurDisk;
498     if (cosCurDisk < 0 && rmin > 0) rmincur /= cosHalf;
499     rzmin[k].set(rmincur,rmin*sinCurDisk);
500 
501     G4double rmaxcur = rtor + rmax*cosCurDisk;
502     if (cosCurDisk > 0) rmaxcur /= cosHalf;
503     rzmax[k].set(rmaxcur,rmax*sinCurDisk);
504 
505     G4double sinTmpDisk = sinCurDisk;
506     sinCurDisk = sinCurDisk*cosStepDisk + cosCurDisk*sinStepDisk;
507     cosCurDisk = cosCurDisk*cosStepDisk - sinTmpDisk*sinStepDisk;
508   }
509 
510   // Loop along slices in Phi. The extent is calculated as cumulative
511   // extent of the slices
512   pMin =  kInfinity;
513   pMax = -kInfinity;
514   G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
515   G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
516   G4double sinCur1 = 0, cosCur1 = 0, sinCur2 = 0, cosCur2 = 0;
517   for (G4int i=0; i<kphi+1; ++i)
518   {
519     if (i == 0)
520     {
521       sinCur1 = sinStart;
522       cosCur1 = cosStart;
523       sinCur2 = sinCur1*cosHalf + cosCur1*sinHalf;
524       cosCur2 = cosCur1*cosHalf - sinCur1*sinHalf;
525     }
526     else
527     {
528       sinCur1 = sinCur2;
529       cosCur1 = cosCur2;
530       sinCur2 = (i == kphi) ? sinEnd : sinCur1*cosStep + cosCur1*sinStep;
531       cosCur2 = (i == kphi) ? cosEnd : cosCur1*cosStep - sinCur1*sinStep;
532     }
533     for (G4int k=0; k<NDISK; ++k)
534     {
535       G4double r1 = rzmin[k].x(), r2 = rzmax[k].x();
536       G4double z1 = rzmin[k].y(), z2 = rzmax[k].y();
537       pols[k][0].set(r1*cosCur1,r1*sinCur1,z1);
538       pols[k][1].set(r2*cosCur1,r2*sinCur1,z2);
539       pols[k][2].set(r2*cosCur2,r2*sinCur2,z2);
540       pols[k][3].set(r1*cosCur2,r1*sinCur2,z1);
541     }
542     pols[NDISK] = pols[0];
543 
544     // get bounding box of current slice
545     G4TwoVector vmin,vmax;
546     G4GeomTools::
547       DiskExtent(rint,rext,sinCur1,cosCur1,sinCur2,cosCur2,vmin,vmax);
548     bmin.setX(vmin.x()); bmin.setY(vmin.y());
549     bmax.setX(vmax.x()); bmax.setY(vmax.y());
550 
551     // set bounding envelope for current slice and adjust extent
552     G4double emin,emax;
553     G4BoundingEnvelope benv(bmin,bmax,polygons);
554     if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
555     if (emin < pMin) pMin = emin;
556     if (emax > pMax) pMax = emax;
557     if (eminlim > pMin && emaxlim < pMax) break; // max possible extent
558   }
559   return (pMin < pMax);
560 }
561 
562 //////////////////////////////////////////////////////////////////////////////
563 //
564 // Return whether point inside/outside/on surface
565 
566 EInside G4Torus::Inside( const G4ThreeVector& p ) const
567 {
568   G4double r, pt2, pPhi, tolRMin, tolRMax ;
569 
570   EInside in = kOutside ;
571 
572   // General precals
573   //
574   r   = std::hypot(p.x(),p.y());
575   pt2 = p.z()*p.z() + (r-fRtor)*(r-fRtor);
576 
577   if (fRmin != 0.0) tolRMin = fRmin + fRminTolerance ;
578   else       tolRMin = 0 ;
579 
580   tolRMax = fRmax - fRmaxTolerance;
581       
582   if (pt2 >= tolRMin*tolRMin && pt2 <= tolRMax*tolRMax )
583   {
584     if ( fDPhi == twopi || pt2 == 0 )  // on torus swept axis
585     {
586       in = kInside ;
587     }
588     else
589     {
590       // Try inner tolerant phi boundaries (=>inside)
591       // if not inside, try outer tolerant phi boundaries
592 
593       pPhi = std::atan2(p.y(),p.x()) ;
594 
595       if ( pPhi < -halfAngTolerance )  { pPhi += twopi ; }  // 0<=pPhi<2pi
596       if ( fSPhi >= 0 )
597       {
598         if ( (std::fabs(pPhi) < halfAngTolerance)
599             && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
600         { 
601             pPhi += twopi ; // 0 <= pPhi < 2pi
602         }
603         if ( (pPhi >= fSPhi + halfAngTolerance)
604             && (pPhi <= fSPhi + fDPhi - halfAngTolerance) )
605         {
606           in = kInside ;
607         }
608           else if ( (pPhi >= fSPhi - halfAngTolerance)
609                  && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
610         {
611           in = kSurface ;
612         }
613       }
614       else  // fSPhi < 0
615       {
616           if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
617             && (pPhi >= fSPhi + fDPhi  + halfAngTolerance) )  {;}
618           else
619           {
620             in = kSurface ;
621           }
622       }
623     }
624   }
625   else   // Try generous boundaries
626   {
627     tolRMin = fRmin - fRminTolerance ;
628     tolRMax = fRmax + fRmaxTolerance ;
629 
630     if (tolRMin < 0 )  { tolRMin = 0 ; }
631 
632     if ( (pt2 >= tolRMin*tolRMin) && (pt2 <= tolRMax*tolRMax) )
633     {
634       if ( (fDPhi == twopi) || (pt2 == 0) ) // Continuous in phi or on z-axis
635       {
636         in = kSurface ;
637       }
638       else // Try outer tolerant phi boundaries only
639       {
640         pPhi = std::atan2(p.y(),p.x()) ;
641 
642         if ( pPhi < -halfAngTolerance )  { pPhi += twopi ; }  // 0<=pPhi<2pi
643         if ( fSPhi >= 0 )
644         {
645           if ( (std::fabs(pPhi) < halfAngTolerance)
646             && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
647           { 
648             pPhi += twopi ; // 0 <= pPhi < 2pi
649           }
650           if ( (pPhi >= fSPhi - halfAngTolerance)
651             && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
652           {
653             in = kSurface;
654           }
655         }
656         else  // fSPhi < 0
657         {
658           if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
659             && (pPhi >= fSPhi + fDPhi  + halfAngTolerance) )  {;}
660           else
661           {
662             in = kSurface ;
663           }
664         }
665       }
666     }
667   }
668   return in ;
669 }
670 
671 /////////////////////////////////////////////////////////////////////////////
672 //
673 // Return unit normal of surface closest to p
674 // - note if point on z axis, ignore phi divided sides
675 // - unsafe if point close to z axis a rmin=0 - no explicit checks
676 
677 G4ThreeVector G4Torus::SurfaceNormal( const G4ThreeVector& p ) const
678 {
679   G4int noSurfaces = 0;  
680   G4double rho, pt, pPhi;
681   G4double distRMin = kInfinity;
682   G4double distSPhi = kInfinity, distEPhi = kInfinity;
683 
684   // To cope with precision loss
685   //
686   const G4double delta = std::max(10.0*kCarTolerance,
687                                   1.0e-8*(fRtor+fRmax));
688   const G4double dAngle = 10.0*kAngTolerance;
689 
690   G4ThreeVector nR, nPs, nPe;
691   G4ThreeVector norm, sumnorm(0.,0.,0.);
692 
693   rho = std::hypot(p.x(),p.y());
694   pt  = std::hypot(p.z(),rho-fRtor);
695 
696   G4double  distRMax = std::fabs(pt - fRmax);
697   if(fRmin != 0.0) distRMin = std::fabs(pt - fRmin);
698 
699   if( rho > delta && pt != 0.0 )
700   {
701     G4double redFactor= (rho-fRtor)/rho;
702     nR = G4ThreeVector( p.x()*redFactor,  // p.x()*(1.-fRtor/rho),
703                         p.y()*redFactor,  // p.y()*(1.-fRtor/rho),
704                         p.z()          );
705     nR *= 1.0/pt;
706   }
707 
708   if ( fDPhi < twopi ) // && rho ) // old limitation against (0,0,z)
709   {
710     if ( rho != 0.0 )
711     {
712       pPhi = std::atan2(p.y(),p.x());
713 
714       if(pPhi < fSPhi-delta)            { pPhi += twopi; }
715       else if(pPhi > fSPhi+fDPhi+delta) { pPhi -= twopi; }
716 
717       distSPhi = std::fabs( pPhi - fSPhi );
718       distEPhi = std::fabs(pPhi-fSPhi-fDPhi);
719     }
720     nPs = G4ThreeVector(std::sin(fSPhi),-std::cos(fSPhi),0);
721     nPe = G4ThreeVector(-std::sin(fSPhi+fDPhi),std::cos(fSPhi+fDPhi),0);
722   } 
723   if( distRMax <= delta )
724   {
725     ++noSurfaces;
726     sumnorm += nR;
727   }
728   else if( (fRmin != 0.0) && (distRMin <= delta) ) // Must not be on both Outer and Inner
729   {
730     ++noSurfaces;
731     sumnorm -= nR;
732   }
733 
734   //  To be on one of the 'phi' surfaces,
735   //  it must be within the 'tube' - with tolerance
736 
737   if( (fDPhi < twopi) && (fRmin-delta <= pt) && (pt <= (fRmax+delta)) )
738   {
739     if (distSPhi <= dAngle)
740     {
741       ++noSurfaces;
742       sumnorm += nPs;
743     }
744     if (distEPhi <= dAngle) 
745     {
746       ++noSurfaces;
747       sumnorm += nPe;
748     }
749   }
750   if ( noSurfaces == 0 )
751   {
752 #ifdef G4CSGDEBUG
753      G4ExceptionDescription ed;
754      ed.precision(16);
755 
756      EInside  inIt= Inside( p );
757      
758      if( inIt != kSurface )
759      {
760         ed << " ERROR>  Surface Normal was called for Torus,"
761            << " with point not on surface." << G4endl;
762      }
763      else
764      {
765         ed << " ERROR>  Surface Normal has not found a surface, "
766            << " despite the point being on the surface. " <<G4endl;
767      }
768 
769      if( inIt != kInside)
770      {
771          ed << " Safety (Dist To In)  = " << DistanceToIn(p) << G4endl;
772      }
773      if( inIt != kOutside)
774      {
775          ed << " Safety (Dist to Out) = " << DistanceToOut(p) << G4endl;
776      }
777      ed << " Coordinates of point : " << p << G4endl;
778      ed << " Parameters  of solid : " << G4endl << *this << G4endl;
779 
780      if( inIt == kSurface )
781      {
782         G4Exception("G4Torus::SurfaceNormal(p)", "GeomSolids1002",
783                     JustWarning, ed,
784                     "Failing to find normal, even though point is on surface!");
785      }
786      else
787      {
788         static const char* NameInside[3]= { "Inside", "Surface", "Outside" };
789         ed << "  The point is " << NameInside[inIt] << " the solid. "<< G4endl;
790         G4Exception("G4Torus::SurfaceNormal(p)", "GeomSolids1002",
791                     JustWarning, ed, "Point p is not on surface !?" );
792      }
793 #endif
794      norm = ApproxSurfaceNormal(p);
795   }
796   else if ( noSurfaces == 1 )  { norm = sumnorm; }
797   else                         { norm = sumnorm.unit(); }
798 
799   return norm ;
800 }
801 
802 //////////////////////////////////////////////////////////////////////////////
803 //
804 // Algorithm for SurfaceNormal() following the original specification
805 // for points not on the surface
806 
807 G4ThreeVector G4Torus::ApproxSurfaceNormal( const G4ThreeVector& p ) const
808 {
809   ENorm side ;
810   G4ThreeVector norm;
811   G4double rho,pt,phi;
812   G4double distRMin,distRMax,distSPhi,distEPhi,distMin;
813 
814   rho = std::hypot(p.x(),p.y());
815   pt  = std::hypot(p.z(),rho-fRtor);
816 
817 #ifdef G4CSGDEBUG
818   G4cout << " G4Torus::ApproximateSurfaceNormal called for point " << p
819          << G4endl;
820 #endif
821    
822   distRMax = std::fabs(pt - fRmax) ;
823 
824   if(fRmin != 0.0)  // First minimum radius
825   {
826     distRMin = std::fabs(pt - fRmin) ;
827 
828     if (distRMin < distRMax)
829     {
830       distMin = distRMin ;
831       side    = kNRMin ;
832     }
833     else
834     {
835       distMin = distRMax ;
836       side    = kNRMax ;
837     }
838   }
839   else
840   {
841     distMin = distRMax ;
842     side    = kNRMax ;
843   }    
844   if ( (fDPhi < twopi) && (rho != 0.0) )
845   {
846     phi = std::atan2(p.y(),p.x()) ; // Protected against (0,0,z) (above rho!=0)
847 
848     if (phi < 0)  { phi += twopi ; }
849 
850     if (fSPhi < 0 )  { distSPhi = std::fabs(phi-(fSPhi+twopi))*rho ; }
851     else             { distSPhi = std::fabs(phi-fSPhi)*rho ; }
852 
853     distEPhi = std::fabs(phi - fSPhi - fDPhi)*rho ;
854 
855     if (distSPhi < distEPhi) // Find new minimum
856     {
857       if (distSPhi<distMin) side = kNSPhi ;
858     }
859     else
860     {
861       if (distEPhi < distMin)  { side = kNEPhi ; }
862     }
863   }  
864   switch (side)
865   {
866     case kNRMin:      // Inner radius
867       norm = G4ThreeVector( -p.x()*(1-fRtor/rho)/pt,
868                             -p.y()*(1-fRtor/rho)/pt,
869                             -p.z()/pt                 ) ;
870       break ;
871     case kNRMax:      // Outer radius
872       norm = G4ThreeVector( p.x()*(1-fRtor/rho)/pt,
873                             p.y()*(1-fRtor/rho)/pt,
874                             p.z()/pt                  ) ;
875       break;
876     case kNSPhi:
877       norm = G4ThreeVector(std::sin(fSPhi),-std::cos(fSPhi),0) ;
878       break;
879     case kNEPhi:
880       norm = G4ThreeVector(-std::sin(fSPhi+fDPhi),std::cos(fSPhi+fDPhi),0) ;
881       break;
882     default:          // Should never reach this case ...
883       DumpInfo();
884       G4Exception("G4Torus::ApproxSurfaceNormal()",
885                   "GeomSolids1002", JustWarning,
886                   "Undefined side for valid surface normal to solid.");
887       break ;
888   } 
889   return norm ;
890 }
891 
892 ///////////////////////////////////////////////////////////////////////
893 //
894 // Calculate distance to shape from outside, along normalised vector
895 // - return kInfinity if no intersection, or intersection distance <= tolerance
896 //
897 // - Compute the intersection with the z planes 
898 //        - if at valid r, phi, return
899 //
900 // -> If point is outer outer radius, compute intersection with rmax
901 //        - if at valid phi,z return
902 //
903 // -> Compute intersection with inner radius, taking largest +ve root
904 //        - if valid (phi), save intersction
905 //
906 //    -> If phi segmented, compute intersections with phi half planes
907 //        - return smallest of valid phi intersections and
908 //          inner radius intersection
909 //
910 // NOTE:
911 // - Precalculations for phi trigonometry are Done `just in time'
912 // - `if valid' implies tolerant checking of intersection points
913 
914 G4double G4Torus::DistanceToIn( const G4ThreeVector& p,
915                                 const G4ThreeVector& v ) const
916 {
917   // Get bounding box of full torus
918   //
919   G4double boxDx  = fRtor + fRmax;
920   G4double boxDy  = boxDx;
921   G4double boxDz  = fRmax;
922   G4double boxMax = boxDx;
923   G4double boxMin = boxDz;
924 
925   // Check if point is traveling away
926   //
927   G4double distX = std::abs(p.x()) - boxDx;
928   G4double distY = std::abs(p.y()) - boxDy;
929   G4double distZ = std::abs(p.z()) - boxDz;
930   if (distX >= -halfCarTolerance && p.x()*v.x() >= 0) return kInfinity;
931   if (distY >= -halfCarTolerance && p.y()*v.y() >= 0) return kInfinity;
932   if (distZ >= -halfCarTolerance && p.z()*v.z() >= 0) return kInfinity;
933 
934   // Calculate safety distance to bounding box
935   // If point is too far, move it closer and calculate distance
936   //
937   G4double Dmax = 32*boxMax; 
938   G4double safe = std::max(std::max(distX,distY),distZ);
939   if (safe > Dmax)
940   {
941     G4double dist = safe - 1.e-8*safe - boxMin; // stay outside after the move
942     dist += DistanceToIn(p + dist*v, v);
943     return (dist >= kInfinity) ? kInfinity : dist;
944   }
945 
946   // Find intersection with torus
947   //
948   G4double snxt=kInfinity, sphi=kInfinity; // snxt = default return value
949 
950   G4double  sd[4] ;
951 
952   // Precalculated trig for phi intersections - used by r,z intersections to
953   //                                            check validity
954 
955   G4bool seg;        // true if segmented
956   G4double hDPhi;    // half dphi
957   G4double cPhi,sinCPhi=0.,cosCPhi=0.;  // central phi
958 
959   G4double tolORMin2;  // `generous' radii squared
960   G4double tolORMax2;
961 
962   G4double Dist,xi,yi,zi,rhoi,it2; // Intersection point variables
963 
964   G4double Comp;
965   G4double cosSPhi,sinSPhi;       // Trig for phi start intersect
966   G4double ePhi,cosEPhi,sinEPhi;  // for phi end intersect
967 
968   // Set phi divided flag and precalcs
969   //
970   if ( fDPhi < twopi )
971   {
972     seg        = true ;
973     hDPhi      = 0.5*fDPhi ;    // half delta phi
974     cPhi       = fSPhi + hDPhi ;
975     sinCPhi    = std::sin(cPhi) ;
976     cosCPhi    = std::cos(cPhi) ;
977   }
978   else
979   {
980     seg = false ;
981   }
982 
983   if (fRmin > fRminTolerance) // Calculate tolerant rmin and rmax
984   {
985     tolORMin2 = (fRmin - fRminTolerance)*(fRmin - fRminTolerance) ;
986   }
987   else
988   {
989     tolORMin2 = 0 ;
990   }
991   tolORMax2 = (fRmax + fRmaxTolerance)*(fRmax + fRmaxTolerance) ;
992 
993   // Intersection with Rmax (possible return) and Rmin (must also check phi)
994 
995   snxt = SolveNumericJT(p,v,fRmax,true);
996 
997   if (fRmin != 0.0)  // Possible Rmin intersection
998   {
999     sd[0] = SolveNumericJT(p,v,fRmin,true);
1000     if ( sd[0] < snxt )  { snxt = sd[0] ; }
1001   }
1002 
1003   //
1004   // Phi segment intersection
1005   //
1006   // o Tolerant of points inside phi planes by up to kCarTolerance*0.5
1007   //
1008   // o NOTE: Large duplication of code between sphi & ephi checks
1009   //         -> only diffs: sphi -> ephi, Comp -> -Comp and half-plane
1010   //            intersection check <=0 -> >=0
1011   //         -> use some form of loop Construct ?
1012 
1013   if (seg)
1014   {
1015     sinSPhi = std::sin(fSPhi) ; // First phi surface ('S'tarting phi)
1016     cosSPhi = std::cos(fSPhi) ;
1017     Comp    = v.x()*sinSPhi - v.y()*cosSPhi ;  // Component in outwards
1018                                                // normal direction
1019     if (Comp < 0 )
1020     {
1021       Dist = (p.y()*cosSPhi - p.x()*sinSPhi) ;
1022 
1023       if (Dist < halfCarTolerance)
1024       {
1025         sphi = Dist/Comp ;
1026         if (sphi < snxt)
1027         {
1028           if ( sphi < 0 )  { sphi = 0 ; }
1029 
1030           xi    = p.x() + sphi*v.x() ;
1031           yi    = p.y() + sphi*v.y() ;
1032           zi    = p.z() + sphi*v.z() ;
1033           rhoi = std::hypot(xi,yi);
1034           it2 = zi*zi + (rhoi-fRtor)*(rhoi-fRtor);
1035 
1036           if ( it2 >= tolORMin2 && it2 <= tolORMax2 )
1037           {
1038             // r intersection is good - check intersecting
1039             // with correct half-plane
1040             //
1041             if ((yi*cosCPhi-xi*sinCPhi)<=0)  { snxt=sphi; }
1042           }
1043         }
1044       }
1045     }
1046     ePhi=fSPhi+fDPhi;    // Second phi surface ('E'nding phi)
1047     sinEPhi=std::sin(ePhi);
1048     cosEPhi=std::cos(ePhi);
1049     Comp=-(v.x()*sinEPhi-v.y()*cosEPhi);
1050         
1051     if ( Comp < 0 )   // Component in outwards normal dirn
1052     {
1053       Dist = -(p.y()*cosEPhi - p.x()*sinEPhi) ;
1054 
1055       if (Dist < halfCarTolerance )
1056       {
1057         sphi = Dist/Comp ;
1058 
1059         if (sphi < snxt )
1060         {
1061           if (sphi < 0 )  { sphi = 0 ; }
1062        
1063           xi    = p.x() + sphi*v.x() ;
1064           yi    = p.y() + sphi*v.y() ;
1065           zi    = p.z() + sphi*v.z() ;
1066           rhoi = std::hypot(xi,yi);
1067           it2 = zi*zi + (rhoi-fRtor)*(rhoi-fRtor);
1068 
1069           if (it2 >= tolORMin2 && it2 <= tolORMax2)
1070           {
1071             // z and r intersections good - check intersecting
1072             // with correct half-plane
1073             //
1074             if ((yi*cosCPhi-xi*sinCPhi)>=0)  { snxt=sphi; }
1075           }    
1076         }
1077       }
1078     }
1079   }
1080   if(snxt < halfCarTolerance)  { snxt = 0.0 ; }
1081 
1082   return snxt ;
1083 }
1084 
1085 /////////////////////////////////////////////////////////////////////////////
1086 //
1087 // Calculate distance (<= actual) to closest surface of shape from outside
1088 // - Calculate distance to z, radial planes
1089 // - Only to phi planes if outside phi extent
1090 // - Return 0 if point inside
1091 
1092 G4double G4Torus::DistanceToIn( const G4ThreeVector& p ) const
1093 {
1094   G4double safe=0.0, safe1, safe2 ;
1095   G4double phiC, cosPhiC, sinPhiC, safePhi, ePhi, cosPsi ;
1096   G4double rho, pt ;
1097   
1098   rho = std::hypot(p.x(),p.y());
1099   pt  = std::hypot(p.z(),rho-fRtor);
1100   safe1 = fRmin - pt ;
1101   safe2 = pt - fRmax ;
1102 
1103   if (safe1 > safe2)  { safe = safe1; }
1104   else                { safe = safe2; }
1105 
1106   if ( fDPhi < twopi && (rho != 0.0) )
1107   {
1108     phiC    = fSPhi + fDPhi*0.5 ;
1109     cosPhiC = std::cos(phiC) ;
1110     sinPhiC = std::sin(phiC) ;
1111     cosPsi  = (p.x()*cosPhiC + p.y()*sinPhiC)/rho ;
1112 
1113     if (cosPsi < std::cos(fDPhi*0.5) ) // Psi=angle from central phi to point
1114     {                                  // Point lies outside phi range
1115       if ((p.y()*cosPhiC - p.x()*sinPhiC) <= 0 )
1116       {
1117         safePhi = std::fabs(p.x()*std::sin(fSPhi) - p.y()*std::cos(fSPhi)) ;
1118       }
1119       else
1120       {
1121         ePhi    = fSPhi + fDPhi ;
1122         safePhi = std::fabs(p.x()*std::sin(ePhi) - p.y()*std::cos(ePhi)) ;
1123       }
1124       if (safePhi > safe)  { safe = safePhi ; }
1125     }
1126   }
1127   if (safe < 0 )  { safe = 0 ; }
1128   return safe;
1129 }
1130 
1131 ///////////////////////////////////////////////////////////////////////////
1132 //
1133 // Calculate distance to surface of shape from `inside', allowing for tolerance
1134 // - Only Calc rmax intersection if no valid rmin intersection
1135 //
1136 
1137 G4double G4Torus::DistanceToOut( const G4ThreeVector& p,
1138                                  const G4ThreeVector& v,
1139                                  const G4bool calcNorm,
1140                                        G4bool* validNorm,
1141                                        G4ThreeVector* n ) const
1142 {
1143   ESide    side = kNull, sidephi = kNull ;
1144   G4double snxt = kInfinity, sphi, sd[4] ;
1145 
1146   // Vars for phi intersection
1147   //
1148   G4double sinSPhi, cosSPhi, ePhi, sinEPhi, cosEPhi;
1149   G4double cPhi, sinCPhi, cosCPhi ;
1150   G4double pDistS, compS, pDistE, compE, sphi2, xi, yi, zi, vphi ;
1151 
1152   // Radial Intersections Defenitions & General Precals
1153 
1154   //////////////////////// new calculation //////////////////////
1155 
1156 #if 1
1157 
1158   // This is the version with the calculation of CalcNorm = true 
1159   // To be done: Check the precision of this calculation.
1160   // If you want return always validNorm = false, then take the version below
1161   
1162   
1163   G4double rho = std::hypot(p.x(),p.y());
1164   G4double pt = hypot(p.z(),rho-fRtor);
1165 
1166   G4double pDotV = p.x()*v.x() + p.y()*v.y() + p.z()*v.z() ;
1167 
1168   G4double tolRMax = fRmax - fRmaxTolerance ;
1169    
1170   G4double vDotNmax   = pDotV - fRtor*(v.x()*p.x() + v.y()*p.y())/rho ;
1171   G4double pDotxyNmax = (1 - fRtor/rho) ;
1172 
1173   if( (pt*pt > tolRMax*tolRMax) && (vDotNmax >= 0) )
1174   {
1175     // On tolerant boundary & heading outwards (or perpendicular to) outer
1176     // radial surface -> leaving immediately with *n for really convex part
1177     // only
1178       
1179     if ( calcNorm && (pDotxyNmax >= -2.*fRmaxTolerance) ) 
1180     {
1181       *n = G4ThreeVector( p.x()*(1 - fRtor/rho)/pt,
1182                           p.y()*(1 - fRtor/rho)/pt,
1183                           p.z()/pt                  ) ;
1184       *validNorm = true ;
1185     }
1186      
1187     return snxt = 0 ; // Leaving by Rmax immediately
1188   }
1189   
1190   snxt = SolveNumericJT(p,v,fRmax,false);  
1191   side = kRMax ;
1192 
1193   // rmin
1194 
1195   if ( fRmin != 0.0 )
1196   {
1197     G4double tolRMin = fRmin + fRminTolerance ;
1198 
1199     if ( (pt*pt < tolRMin*tolRMin) && (vDotNmax < 0) )
1200     {
1201       if (calcNorm)  { *validNorm = false ; } // Concave surface of the torus
1202       return  snxt = 0 ;                      // Leaving by Rmin immediately
1203     }
1204     
1205     sd[0] = SolveNumericJT(p,v,fRmin,false);
1206     if ( sd[0] < snxt )
1207     {
1208       snxt = sd[0] ;
1209       side = kRMin ;
1210     }
1211   }
1212 
1213 #else
1214 
1215   // this is the "conservative" version which return always validnorm = false
1216   // NOTE: using this version the unit test testG4Torus will break
1217 
1218   snxt = SolveNumericJT(p,v,fRmax,false);  
1219   side = kRMax ;
1220 
1221   if ( fRmin )
1222   {
1223     sd[0] = SolveNumericJT(p,v,fRmin,false);
1224     if ( sd[0] < snxt )
1225     {
1226       snxt = sd[0] ;
1227       side = kRMin ;
1228     }
1229   }
1230 
1231   if ( calcNorm && (snxt == 0.0) )
1232   {
1233     *validNorm = false ;    // Leaving solid, but possible re-intersection
1234     return snxt  ;
1235   }
1236 
1237 #endif
1238   
1239   if (fDPhi < twopi)  // Phi Intersections
1240   {
1241     sinSPhi = std::sin(fSPhi) ;
1242     cosSPhi = std::cos(fSPhi) ;
1243     ePhi    = fSPhi + fDPhi ;
1244     sinEPhi = std::sin(ePhi) ;
1245     cosEPhi = std::cos(ePhi) ;
1246     cPhi    = fSPhi + fDPhi*0.5 ;
1247     sinCPhi = std::sin(cPhi) ;
1248     cosCPhi = std::cos(cPhi) ;
1249     
1250     // angle calculation with correction 
1251     // of difference in domain of atan2 and Sphi
1252     //
1253     vphi = std::atan2(v.y(),v.x()) ;
1254      
1255     if ( vphi < fSPhi - halfAngTolerance  )    { vphi += twopi; }
1256     else if ( vphi > ePhi + halfAngTolerance ) { vphi -= twopi; }
1257 
1258     if ( (p.x() != 0.0) || (p.y() != 0.0) ) // Check if on z axis (rho not needed later)
1259     {
1260       pDistS = p.x()*sinSPhi - p.y()*cosSPhi ; // pDist -ve when inside
1261       pDistE = -p.x()*sinEPhi + p.y()*cosEPhi ;
1262 
1263       // Comp -ve when in direction of outwards normal
1264       //
1265       compS   = -sinSPhi*v.x() + cosSPhi*v.y() ;
1266       compE   = sinEPhi*v.x() - cosEPhi*v.y() ;
1267       sidephi = kNull ;
1268      
1269       if( ( (fDPhi <= pi) && ( (pDistS <= halfCarTolerance)
1270                             && (pDistE <= halfCarTolerance) ) )
1271        || ( (fDPhi >  pi) && ((pDistS <=  halfCarTolerance)
1272                             || (pDistE <=  halfCarTolerance) ) )  )
1273       {
1274         // Inside both phi *full* planes
1275 
1276         if ( compS < 0 )
1277         {
1278           sphi = pDistS/compS ;
1279             
1280           if (sphi >= -halfCarTolerance)
1281           {
1282             xi = p.x() + sphi*v.x() ;
1283             yi = p.y() + sphi*v.y() ;
1284               
1285             // Check intersecting with correct half-plane
1286             // (if not -> no intersect)
1287             //
1288             if ( (std::fabs(xi)<=kCarTolerance)
1289               && (std::fabs(yi)<=kCarTolerance) )
1290             {
1291               sidephi = kSPhi;
1292               if ( ((fSPhi-halfAngTolerance)<=vphi)
1293                 && ((ePhi+halfAngTolerance)>=vphi) )
1294               {
1295                 sphi = kInfinity;
1296               }
1297             }
1298             else if ( yi*cosCPhi-xi*sinCPhi >=0 )
1299             {
1300               sphi = kInfinity ;
1301             }
1302             else
1303             {
1304               sidephi = kSPhi ;
1305             }       
1306           }
1307           else
1308           {
1309             sphi = kInfinity ;
1310           }
1311         }
1312         else
1313         {
1314           sphi = kInfinity ;
1315         }
1316 
1317         if ( compE < 0 )
1318         {
1319           sphi2 = pDistE/compE ;
1320             
1321           // Only check further if < starting phi intersection
1322           //
1323           if ( (sphi2 > -kCarTolerance) && (sphi2 < sphi) )
1324           {
1325             xi = p.x() + sphi2*v.x() ;
1326             yi = p.y() + sphi2*v.y() ;
1327               
1328             if ( (std::fabs(xi)<=kCarTolerance)
1329               && (std::fabs(yi)<=kCarTolerance) )
1330             {
1331               // Leaving via ending phi
1332               //
1333               if( (fSPhi-halfAngTolerance > vphi)
1334                   || (ePhi+halfAngTolerance < vphi) )
1335               {
1336                 sidephi = kEPhi ;
1337                 sphi = sphi2;
1338               }
1339             } 
1340             else    // Check intersecting with correct half-plane 
1341             {
1342               if ( (yi*cosCPhi-xi*sinCPhi) >= 0)
1343               {
1344                 // Leaving via ending phi
1345                 //
1346                 sidephi = kEPhi ;
1347                 sphi = sphi2;
1348                
1349               }
1350             }
1351           }
1352         }
1353       }
1354       else
1355       {
1356         sphi = kInfinity ;
1357       }
1358     } 
1359     else
1360     {
1361       // On z axis + travel not || to z axis -> if phi of vector direction
1362       // within phi of shape, Step limited by rmax, else Step =0
1363 
1364       vphi = std::atan2(v.y(),v.x());
1365  
1366       if ( ( fSPhi-halfAngTolerance <= vphi ) && 
1367            ( vphi <= ( ePhi+halfAngTolerance ) ) )
1368       {
1369         sphi = kInfinity;
1370       }
1371       else
1372       {
1373         sidephi = kSPhi ; // arbitrary 
1374         sphi=0;
1375       }
1376     }
1377 
1378     // Order intersections
1379 
1380     if (sphi<snxt)
1381     {
1382       snxt=sphi;
1383       side=sidephi;
1384     }     
1385   }
1386 
1387   G4double rhoi,it,iDotxyNmax ;
1388   // Note: by numerical computation we know where the ray hits the torus
1389   // So I propose to return the side where the ray hits
1390 
1391   if (calcNorm)
1392   {
1393     switch(side)
1394     {
1395       case kRMax:                     // n is unit vector 
1396         xi    = p.x() + snxt*v.x() ;
1397         yi    = p.y() + snxt*v.y() ;
1398         zi    = p.z() + snxt*v.z() ;
1399         rhoi = std::hypot(xi,yi);
1400         it = hypot(zi,rhoi-fRtor);
1401 
1402         iDotxyNmax = (1-fRtor/rhoi) ;
1403         if(iDotxyNmax >= -2.*fRmaxTolerance) // really convex part of Rmax
1404         {                       
1405           *n = G4ThreeVector( xi*(1-fRtor/rhoi)/it,
1406                               yi*(1-fRtor/rhoi)/it,
1407                               zi/it                 ) ;
1408           *validNorm = true ;
1409         }
1410         else
1411         {
1412           *validNorm = false ; // concave-convex part of Rmax
1413         }
1414         break ;
1415 
1416       case kRMin:
1417         *validNorm = false ;  // Rmin is concave or concave-convex
1418         break;
1419 
1420       case kSPhi:
1421         if (fDPhi <= pi )
1422         {
1423           *n=G4ThreeVector(std::sin(fSPhi),-std::cos(fSPhi),0);
1424           *validNorm=true;
1425         }
1426         else
1427         {
1428           *validNorm = false ;
1429         }
1430         break ;
1431 
1432       case kEPhi:
1433         if (fDPhi <= pi)
1434         {
1435           *n=G4ThreeVector(-std::sin(fSPhi+fDPhi),std::cos(fSPhi+fDPhi),0);
1436           *validNorm=true;
1437         }
1438         else
1439         {
1440           *validNorm = false ;
1441         }
1442         break;
1443 
1444       default:
1445 
1446         // It seems we go here from time to time ...
1447 
1448         G4cout << G4endl;
1449         DumpInfo();
1450         std::ostringstream message;
1451         G4long oldprc = message.precision(16);
1452         message << "Undefined side for valid surface normal to solid."
1453                 << G4endl
1454                 << "Position:"  << G4endl << G4endl
1455                 << "p.x() = "   << p.x()/mm << " mm" << G4endl
1456                 << "p.y() = "   << p.y()/mm << " mm" << G4endl
1457                 << "p.z() = "   << p.z()/mm << " mm" << G4endl << G4endl
1458                 << "Direction:" << G4endl << G4endl
1459                 << "v.x() = "   << v.x() << G4endl
1460                 << "v.y() = "   << v.y() << G4endl
1461                 << "v.z() = "   << v.z() << G4endl << G4endl
1462                 << "Proposed distance :" << G4endl << G4endl
1463                 << "snxt = " << snxt/mm << " mm" << G4endl;
1464         message.precision(oldprc);
1465         G4Exception("G4Torus::DistanceToOut(p,v,..)",
1466                     "GeomSolids1002",JustWarning, message);
1467         break;
1468     }
1469   }
1470   if ( snxt<halfCarTolerance )  { snxt=0 ; }
1471 
1472   return snxt;
1473 }
1474 
1475 /////////////////////////////////////////////////////////////////////////
1476 //
1477 // Calculate distance (<=actual) to closest surface of shape from inside
1478 
1479 G4double G4Torus::DistanceToOut( const G4ThreeVector& p ) const
1480 {
1481   G4double safe=0.0,safeR1,safeR2;
1482   G4double rho,pt ;
1483   G4double safePhi,phiC,cosPhiC,sinPhiC,ePhi;
1484   
1485   rho = std::hypot(p.x(),p.y());
1486   pt  = std::hypot(p.z(),rho-fRtor);
1487   
1488 #ifdef G4CSGDEBUG
1489   if( Inside(p) == kOutside )
1490   {
1491      G4long oldprc = G4cout.precision(16) ;
1492      G4cout << G4endl ;
1493      DumpInfo();
1494      G4cout << "Position:"  << G4endl << G4endl ;
1495      G4cout << "p.x() = "   << p.x()/mm << " mm" << G4endl ;
1496      G4cout << "p.y() = "   << p.y()/mm << " mm" << G4endl ;
1497      G4cout << "p.z() = "   << p.z()/mm << " mm" << G4endl << G4endl ;
1498      G4cout.precision(oldprc);
1499      G4Exception("G4Torus::DistanceToOut(p)", "GeomSolids1002",
1500                  JustWarning, "Point p is outside !?" );
1501   }
1502 #endif
1503 
1504   if (fRmin != 0.0)
1505   {
1506     safeR1 = pt - fRmin ;
1507     safeR2 = fRmax - pt ;
1508 
1509     if (safeR1 < safeR2)  { safe = safeR1 ; }
1510     else                  { safe = safeR2 ; }
1511   }
1512   else
1513   {
1514     safe = fRmax - pt ;
1515   }  
1516 
1517   // Check if phi divided, Calc distances closest phi plane
1518   //
1519   if (fDPhi < twopi) // Above/below central phi of Torus?
1520   {
1521     phiC    = fSPhi + fDPhi*0.5 ;
1522     cosPhiC = std::cos(phiC) ;
1523     sinPhiC = std::sin(phiC) ;
1524 
1525     if ((p.y()*cosPhiC-p.x()*sinPhiC)<=0)
1526     {
1527       safePhi = -(p.x()*std::sin(fSPhi) - p.y()*std::cos(fSPhi)) ;
1528     }
1529     else
1530     {
1531       ePhi    = fSPhi + fDPhi ;
1532       safePhi = (p.x()*std::sin(ePhi) - p.y()*std::cos(ePhi)) ;
1533     }
1534     if (safePhi < safe)  { safe = safePhi ; }
1535   }
1536   if (safe < 0)  { safe = 0 ; }
1537   return safe ;  
1538 }
1539 
1540 //////////////////////////////////////////////////////////////////////////
1541 //
1542 // Stream object contents to an output stream
1543 
1544 G4GeometryType G4Torus::GetEntityType() const
1545 {
1546   return {"G4Torus"};
1547 }
1548 
1549 //////////////////////////////////////////////////////////////////////////
1550 //
1551 // Make a clone of the object
1552 //
1553 G4VSolid* G4Torus::Clone() const
1554 {
1555   return new G4Torus(*this);
1556 }
1557 
1558 //////////////////////////////////////////////////////////////////////////
1559 //
1560 // Stream object contents to an output stream
1561 
1562 std::ostream& G4Torus::StreamInfo( std::ostream& os ) const
1563 {
1564   G4long oldprc = os.precision(16);
1565   os << "-----------------------------------------------------------\n"
1566      << "    *** Dump for solid - " << GetName() << " ***\n"
1567      << "    ===================================================\n"
1568      << " Solid type: G4Torus\n"
1569      << " Parameters: \n"
1570      << "    inner radius: " << fRmin/mm << " mm \n"
1571      << "    outer radius: " << fRmax/mm << " mm \n"
1572      << "    swept radius: " << fRtor/mm << " mm \n"
1573      << "    starting phi: " << fSPhi/degree << " degrees \n"
1574      << "    delta phi   : " << fDPhi/degree << " degrees \n"
1575      << "-----------------------------------------------------------\n";
1576   os.precision(oldprc);
1577 
1578   return os;
1579 }
1580 
1581 ////////////////////////////////////////////////////////////////////////////
1582 //
1583 // GetPointOnSurface
1584 
1585 G4ThreeVector G4Torus::GetPointOnSurface() const
1586 {
1587   G4double cosu, sinu,cosv, sinv, aOut, aIn, aSide, chose, phi, theta, rRand;
1588    
1589   phi   = G4RandFlat::shoot(fSPhi,fSPhi+fDPhi);
1590   theta = G4RandFlat::shoot(0.,twopi);
1591   
1592   cosu   = std::cos(phi);    sinu = std::sin(phi);
1593   cosv   = std::cos(theta);  sinv = std::sin(theta); 
1594 
1595   // compute the areas
1596 
1597   aOut   = (fDPhi)*twopi*fRtor*fRmax;
1598   aIn    = (fDPhi)*twopi*fRtor*fRmin;
1599   aSide  = pi*(fRmax*fRmax-fRmin*fRmin);
1600   
1601   if ((fSPhi == 0) && (fDPhi == twopi)){ aSide = 0; }
1602   chose = G4RandFlat::shoot(0.,aOut + aIn + 2.*aSide);
1603 
1604   if(chose < aOut)
1605   {
1606     return { (fRtor+fRmax*cosv)*cosu, (fRtor+fRmax*cosv)*sinu, fRmax*sinv };
1607   }
1608   else if( (chose >= aOut) && (chose < aOut + aIn) )
1609   {
1610     return { (fRtor+fRmin*cosv)*cosu, (fRtor+fRmin*cosv)*sinu, fRmin*sinv };
1611   }
1612   else if( (chose >= aOut + aIn) && (chose < aOut + aIn + aSide) )
1613   {
1614     rRand = GetRadiusInRing(fRmin,fRmax);
1615     return { (fRtor+rRand*cosv)*std::cos(fSPhi),
1616              (fRtor+rRand*cosv)*std::sin(fSPhi), rRand*sinv };
1617   }
1618   else
1619   {   
1620     rRand = GetRadiusInRing(fRmin,fRmax);
1621     return { (fRtor+rRand*cosv)*std::cos(fSPhi+fDPhi),
1622              (fRtor+rRand*cosv)*std::sin(fSPhi+fDPhi), rRand*sinv };
1623    }
1624 }
1625 
1626 ///////////////////////////////////////////////////////////////////////
1627 //
1628 // Visualisation Functions
1629 
1630 void G4Torus::DescribeYourselfTo ( G4VGraphicsScene& scene ) const 
1631 {
1632   scene.AddSolid (*this);
1633 }
1634 
1635 G4Polyhedron* G4Torus::CreatePolyhedron () const 
1636 {
1637   return new G4PolyhedronTorus (fRmin, fRmax, fRtor, fSPhi, fDPhi);
1638 }
1639 
1640 #endif // !defined(G4GEOM_USE_TORUS) || !defined(G4GEOM_USE_SYS_USOLIDS)
1641