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