Geant4 9.6.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//
27// $Id$
28//
29//
30// class G4Torus
31//
32// Implementation
33//
34// 05.04.12 M.Kelsey: Use sqrt(r) in GetPointOnSurface() for uniform points
35// 02.10.07 T.Nikitina: Bug fixed in SolveNumericJT(), b.969:segmentation fault.
36// rootsrefined is used only if the number of refined roots
37// is the same as for primary roots.
38// 02.10.07 T.Nikitina: Bug fixed in CalculateExtent() for case of non-rotated
39// full-phi torus:protect against negative value for sqrt,
40// correct formula for delta.
41// 20.11.05 V.Grichine: Bug fixed in Inside(p) for phi sections, b.810
42// 25.08.05 O.Link: new methods for DistanceToIn/Out using JTPolynomialSolver
43// 07.06.05 V.Grichine: SurfaceNormal(p) for rho=0, Constructor as G4Cons
44// 03.05.05 V.Grichine: SurfaceNormal(p) according to J. Apostolakis proposal
45// 18.03.04 V.Grichine: bug fixed in DistanceToIn(p)
46// 11.01.01 E.Medernach: Use G4PolynomialSolver to find roots
47// 03.10.00 E.Medernach: SafeNewton added
48// 31.08.00 E.Medernach: numerical computation of roots wuth bounding
49// volume technique
50// 26.05.00 V.Grichine: new fuctions developed by O.Cremonesi were added
51// 06.03.00 V.Grichine: modifications in Distance ToOut(p,v,...)
52// 19.11.99 V.Grichine: side = kNull in Distance ToOut(p,v,...)
53// 09.10.98 V.Grichine: modifications in Distance ToOut(p,v,...)
54// 30.10.96 V.Grichine: first implementation with G4Tubs elements in Fs
55//
56
57#include "G4Torus.hh"
58
59#include "G4VoxelLimits.hh"
60#include "G4AffineTransform.hh"
63
65
66#include "meshdefs.hh"
67
68#include "Randomize.hh"
69
70#include "G4VGraphicsScene.hh"
71#include "G4Polyhedron.hh"
72#include "G4NURBS.hh"
73#include "G4NURBStube.hh"
74#include "G4NURBScylinder.hh"
75#include "G4NURBStubesector.hh"
76
77using namespace CLHEP;
78
79///////////////////////////////////////////////////////////////
80//
81// Constructor - check parameters, convert angles so 0<sphi+dpshi<=2_PI
82// - note if pdphi>2PI then reset to 2PI
83
85 G4double pRmin,
86 G4double pRmax,
87 G4double pRtor,
88 G4double pSPhi,
89 G4double pDPhi)
90 : G4CSGSolid(pName)
91{
92 SetAllParameters(pRmin, pRmax, pRtor, pSPhi, pDPhi);
93}
94
95////////////////////////////////////////////////////////////////////////////
96//
97//
98
99void
101 G4double pRmax,
102 G4double pRtor,
103 G4double pSPhi,
104 G4double pDPhi )
105{
106 const G4double fEpsilon = 4.e-11; // relative tolerance of radii
107
108 fCubicVolume = 0.;
109 fSurfaceArea = 0.;
110 fpPolyhedron = 0;
111
114
115 if ( pRtor >= pRmax+1.e3*kCarTolerance ) // Check swept radius, as in G4Cons
116 {
117 fRtor = pRtor ;
118 }
119 else
120 {
121 std::ostringstream message;
122 message << "Invalid swept radius for Solid: " << GetName() << G4endl
123 << " pRtor = " << pRtor << ", pRmax = " << pRmax;
124 G4Exception("G4Torus::SetAllParameters()",
125 "GeomSolids0002", FatalException, message);
126 }
127
128 // Check radii, as in G4Cons
129 //
130 if ( pRmin < pRmax - 1.e2*kCarTolerance && pRmin >= 0 )
131 {
132 if (pRmin >= 1.e2*kCarTolerance) { fRmin = pRmin ; }
133 else { fRmin = 0.0 ; }
134 fRmax = pRmax ;
135 }
136 else
137 {
138 std::ostringstream message;
139 message << "Invalid values of radii for Solid: " << GetName() << G4endl
140 << " pRmin = " << pRmin << ", pRmax = " << pRmax;
141 G4Exception("G4Torus::SetAllParameters()",
142 "GeomSolids0002", FatalException, message);
143 }
144
145 // Relative tolerances
146 //
147 fRminTolerance = (fRmin)
148 ? 0.5*std::max( kRadTolerance, fEpsilon*(fRtor-fRmin )) : 0;
149 fRmaxTolerance = 0.5*std::max( kRadTolerance, fEpsilon*(fRtor+fRmax) );
150
151 // Check angles
152 //
153 if ( pDPhi >= twopi ) { fDPhi = twopi ; }
154 else
155 {
156 if (pDPhi > 0) { fDPhi = pDPhi ; }
157 else
158 {
159 std::ostringstream message;
160 message << "Invalid Z delta-Phi for Solid: " << GetName() << G4endl
161 << " pDPhi = " << pDPhi;
162 G4Exception("G4Torus::SetAllParameters()",
163 "GeomSolids0002", FatalException, message);
164 }
165 }
166
167 // Ensure psphi in 0-2PI or -2PI-0 range if shape crosses 0
168 //
169 fSPhi = pSPhi;
170
171 if (fSPhi < 0) { fSPhi = twopi-std::fmod(std::fabs(fSPhi),twopi) ; }
172 else { fSPhi = std::fmod(fSPhi,twopi) ; }
173
174 if (fSPhi+fDPhi > twopi) { fSPhi-=twopi ; }
175}
176
177///////////////////////////////////////////////////////////////////////
178//
179// Fake default constructor - sets only member data and allocates memory
180// for usage restricted to object persistency.
181//
182G4Torus::G4Torus( __void__& a )
183 : G4CSGSolid(a), fRmin(0.), fRmax(0.), fRtor(0.), fSPhi(0.),
184 fDPhi(0.), fRminTolerance(0.), fRmaxTolerance(0. ),
185 kRadTolerance(0.), kAngTolerance(0.)
186{
187}
188
189//////////////////////////////////////////////////////////////////////
190//
191// Destructor
192
194{}
195
196//////////////////////////////////////////////////////////////////////////
197//
198// Copy constructor
199
201 : G4CSGSolid(rhs), fRmin(rhs.fRmin),fRmax(rhs.fRmax),
202 fRtor(rhs.fRtor),fSPhi(rhs.fSPhi),fDPhi(rhs.fDPhi),
203 fRminTolerance(rhs.fRminTolerance), fRmaxTolerance(rhs.fRmaxTolerance),
204 kRadTolerance(rhs.kRadTolerance), kAngTolerance(rhs.kAngTolerance)
205{
206}
207
208//////////////////////////////////////////////////////////////////////////
209//
210// Assignment operator
211
213{
214 // Check assignment to self
215 //
216 if (this == &rhs) { return *this; }
217
218 // Copy base class data
219 //
221
222 // Copy data
223 //
224 fRmin = rhs.fRmin; fRmax = rhs.fRmax;
225 fRtor = rhs.fRtor; fSPhi = rhs.fSPhi; fDPhi = rhs.fDPhi;
226 fRminTolerance = rhs.fRminTolerance; fRmaxTolerance = rhs.fRmaxTolerance;
227 kRadTolerance = rhs.kRadTolerance; kAngTolerance = rhs.kAngTolerance;
228
229 return *this;
230}
231
232//////////////////////////////////////////////////////////////////////
233//
234// Dispatch to parameterisation for replication mechanism dimension
235// computation & modification.
236
238 const G4int n,
239 const G4VPhysicalVolume* pRep )
240{
241 p->ComputeDimensions(*this,n,pRep);
242}
243
244
245
246////////////////////////////////////////////////////////////////////////////////
247//
248// Calculate the real roots to torus surface.
249// Returns negative solutions as well.
250
251void G4Torus::TorusRootsJT( const G4ThreeVector& p,
252 const G4ThreeVector& v,
253 G4double r,
254 std::vector<G4double>& roots ) const
255{
256
257 G4int i, num ;
258 G4double c[5], srd[4], si[4] ;
259
260 G4double Rtor2 = fRtor*fRtor, r2 = r*r ;
261
262 G4double pDotV = p.x()*v.x() + p.y()*v.y() + p.z()*v.z() ;
263 G4double pRad2 = p.x()*p.x() + p.y()*p.y() + p.z()*p.z() ;
264
265 c[0] = 1.0 ;
266 c[1] = 4*pDotV ;
267 c[2] = 2*(pRad2 + 2*pDotV*pDotV - Rtor2 - r2 + 2*Rtor2*v.z()*v.z()) ;
268 c[3] = 4*(pDotV*(pRad2 - Rtor2 - r2) + 2*Rtor2*p.z()*v.z()) ;
269 c[4] = pRad2*pRad2 - 2*pRad2*(Rtor2+r2)
270 + 4*Rtor2*p.z()*p.z() + (Rtor2-r2)*(Rtor2-r2) ;
271
272 G4JTPolynomialSolver torusEq;
273
274 num = torusEq.FindRoots( c, 4, srd, si );
275
276 for ( i = 0; i < num; i++ )
277 {
278 if( si[i] == 0. ) { roots.push_back(srd[i]) ; } // store real roots
279 }
280
281 std::sort(roots.begin() , roots.end() ) ; // sorting with <
282}
283
284//////////////////////////////////////////////////////////////////////////////
285//
286// Interface for DistanceToIn and DistanceToOut.
287// Calls TorusRootsJT and returns the smalles possible distance to
288// the surface.
289// Attention: Difference in DistanceToIn/Out for points p on the surface.
290
291G4double G4Torus::SolveNumericJT( const G4ThreeVector& p,
292 const G4ThreeVector& v,
293 G4double r,
294 G4bool IsDistanceToIn ) const
295{
296 G4double bigdist = 10*mm ;
297 G4double tmin = kInfinity ;
298 G4double t, scal ;
299 static const G4double halfCarTolerance = 0.5*kCarTolerance;
300 static const G4double halfAngTolerance = 0.5*kAngTolerance;
301
302 // calculate the distances to the intersections with the Torus
303 // from a given point p and direction v.
304 //
305 std::vector<G4double> roots ;
306 std::vector<G4double> rootsrefined ;
307 TorusRootsJT(p,v,r,roots) ;
308
309 G4ThreeVector ptmp ;
310
311 // determine the smallest non-negative solution
312 //
313 for ( size_t k = 0 ; k<roots.size() ; k++ )
314 {
315 t = roots[k] ;
316
317 if ( t < -halfCarTolerance ) { continue ; } // skip negative roots
318
319 if ( t > bigdist && t<kInfinity ) // problem with big distances
320 {
321 ptmp = p + t*v ;
322 TorusRootsJT(ptmp,v,r,rootsrefined) ;
323 if ( rootsrefined.size()==roots.size() )
324 {
325 t = t + rootsrefined[k] ;
326 }
327 }
328
329 ptmp = p + t*v ; // calculate the position of the proposed intersection
330
331 G4double theta = std::atan2(ptmp.y(),ptmp.x());
332
333 if ( fSPhi >= 0 )
334 {
335 if ( theta < - halfAngTolerance ) { theta += twopi; }
336 if ( (std::fabs(theta) < halfAngTolerance)
337 && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
338 {
339 theta += twopi ; // 0 <= theta < 2pi
340 }
341 }
342 if ((fSPhi <= -pi )&&(theta>halfAngTolerance)) { theta = theta-twopi; }
343
344 // We have to verify if this root is inside the region between
345 // fSPhi and fSPhi + fDPhi
346 //
347 if ( (theta - fSPhi >= - halfAngTolerance)
348 && (theta - (fSPhi + fDPhi) <= halfAngTolerance) )
349 {
350 // check if P is on the surface, and called from DistanceToIn
351 // DistanceToIn has to return 0.0 if particle is going inside the solid
352
353 if ( IsDistanceToIn == true )
354 {
355 if (std::fabs(t) < halfCarTolerance )
356 {
357 // compute scalar product at position p : v.n
358 // ( n taken from SurfaceNormal, not normalized )
359
360 scal = v* G4ThreeVector( p.x()*(1-fRtor/std::sqrt(p.x()*p.x()
361 + p.y()*p.y())),
362 p.y()*(1-fRtor/std::sqrt(p.x()*p.x()
363 + p.y()*p.y())),
364 p.z() );
365
366 // change sign in case of inner radius
367 //
368 if ( r == GetRmin() ) { scal = -scal ; }
369 if ( scal < 0 ) { return 0.0 ; }
370 }
371 }
372
373 // check if P is on the surface, and called from DistanceToOut
374 // DistanceToIn has to return 0.0 if particle is leaving the solid
375
376 if ( IsDistanceToIn == false )
377 {
378 if (std::fabs(t) < halfCarTolerance )
379 {
380 // compute scalar product at position p : v.n
381 //
382 scal = v* G4ThreeVector( p.x()*(1-fRtor/std::sqrt(p.x()*p.x()
383 + p.y()*p.y())),
384 p.y()*(1-fRtor/std::sqrt(p.x()*p.x()
385 + p.y()*p.y())),
386 p.z() );
387
388 // change sign in case of inner radius
389 //
390 if ( r == GetRmin() ) { scal = -scal ; }
391 if ( scal > 0 ) { return 0.0 ; }
392 }
393 }
394
395 // check if distance is larger than 1/2 kCarTolerance
396 //
397 if( t > halfCarTolerance )
398 {
399 tmin = t ;
400 return tmin ;
401 }
402 }
403 }
404
405 return tmin;
406}
407
408/////////////////////////////////////////////////////////////////////////////
409//
410// Calculate extent under transform and specified limit
411
413 const G4VoxelLimits& pVoxelLimit,
414 const G4AffineTransform& pTransform,
415 G4double& pMin, G4double& pMax) const
416{
417 if ((!pTransform.IsRotated()) && (fDPhi==twopi) && (fRmin==0))
418 {
419 // Special case handling for unrotated solid torus
420 // Compute x/y/z mins and maxs for bounding box respecting limits,
421 // with early returns if outside limits. Then switch() on pAxis,
422 // and compute exact x and y limit for x/y case
423
424 G4double xoffset,xMin,xMax;
425 G4double yoffset,yMin,yMax;
426 G4double zoffset,zMin,zMax;
427
428 G4double RTorus,delta,diff1,diff2,maxDiff,newMin,newMax;
429 G4double xoff1,xoff2,yoff1,yoff2;
430
431 xoffset = pTransform.NetTranslation().x();
432 xMin = xoffset - fRmax - fRtor ;
433 xMax = xoffset + fRmax + fRtor ;
434
435 if (pVoxelLimit.IsXLimited())
436 {
437 if ( (xMin > pVoxelLimit.GetMaxXExtent()+kCarTolerance)
438 || (xMax < pVoxelLimit.GetMinXExtent()-kCarTolerance) )
439 return false ;
440 else
441 {
442 if (xMin < pVoxelLimit.GetMinXExtent())
443 {
444 xMin = pVoxelLimit.GetMinXExtent() ;
445 }
446 if (xMax > pVoxelLimit.GetMaxXExtent())
447 {
448 xMax = pVoxelLimit.GetMaxXExtent() ;
449 }
450 }
451 }
452 yoffset = pTransform.NetTranslation().y();
453 yMin = yoffset - fRmax - fRtor ;
454 yMax = yoffset + fRmax + fRtor ;
455
456 if (pVoxelLimit.IsYLimited())
457 {
458 if ( (yMin > pVoxelLimit.GetMaxYExtent()+kCarTolerance)
459 || (yMax < pVoxelLimit.GetMinYExtent()-kCarTolerance) )
460 {
461 return false ;
462 }
463 else
464 {
465 if (yMin < pVoxelLimit.GetMinYExtent() )
466 {
467 yMin = pVoxelLimit.GetMinYExtent() ;
468 }
469 if (yMax > pVoxelLimit.GetMaxYExtent() )
470 {
471 yMax = pVoxelLimit.GetMaxYExtent() ;
472 }
473 }
474 }
475 zoffset = pTransform.NetTranslation().z() ;
476 zMin = zoffset - fRmax ;
477 zMax = zoffset + fRmax ;
478
479 if (pVoxelLimit.IsZLimited())
480 {
481 if ( (zMin > pVoxelLimit.GetMaxZExtent()+kCarTolerance)
482 || (zMax < pVoxelLimit.GetMinZExtent()-kCarTolerance) )
483 {
484 return false ;
485 }
486 else
487 {
488 if (zMin < pVoxelLimit.GetMinZExtent() )
489 {
490 zMin = pVoxelLimit.GetMinZExtent() ;
491 }
492 if (zMax > pVoxelLimit.GetMaxZExtent() )
493 {
494 zMax = pVoxelLimit.GetMaxZExtent() ;
495 }
496 }
497 }
498
499 // Known to cut cylinder
500
501 switch (pAxis)
502 {
503 case kXAxis:
504 yoff1=yoffset-yMin;
505 yoff2=yMax-yoffset;
506 if ( yoff1 >= 0 && yoff2 >= 0 )
507 {
508 // Y limits cross max/min x => no change
509 //
510 pMin = xMin ;
511 pMax = xMax ;
512 }
513 else
514 {
515 // Y limits don't cross max/min x => compute max delta x,
516 // hence new mins/maxs
517 //
518
519 RTorus=fRmax+fRtor;
520 delta = RTorus*RTorus - yoff1*yoff1;
521 diff1 = (delta>0.) ? std::sqrt(delta) : 0.;
522 delta = RTorus*RTorus - yoff2*yoff2;
523 diff2 = (delta>0.) ? std::sqrt(delta) : 0.;
524 maxDiff = (diff1 > diff2) ? diff1:diff2 ;
525 newMin = xoffset - maxDiff ;
526 newMax = xoffset + maxDiff ;
527 pMin = (newMin < xMin) ? xMin : newMin ;
528 pMax = (newMax > xMax) ? xMax : newMax ;
529 }
530 break;
531
532 case kYAxis:
533 xoff1 = xoffset - xMin ;
534 xoff2 = xMax - xoffset ;
535 if (xoff1 >= 0 && xoff2 >= 0 )
536 {
537 // X limits cross max/min y => no change
538 //
539 pMin = yMin ;
540 pMax = yMax ;
541 }
542 else
543 {
544 // X limits don't cross max/min y => compute max delta y,
545 // hence new mins/maxs
546 //
547 RTorus=fRmax+fRtor;
548 delta = RTorus*RTorus - xoff1*xoff1;
549 diff1 = (delta>0.) ? std::sqrt(delta) : 0.;
550 delta = RTorus*RTorus - xoff2*xoff2;
551 diff2 = (delta>0.) ? std::sqrt(delta) : 0.;
552 maxDiff = (diff1 > diff2) ? diff1 : diff2 ;
553 newMin = yoffset - maxDiff ;
554 newMax = yoffset + maxDiff ;
555 pMin = (newMin < yMin) ? yMin : newMin ;
556 pMax = (newMax > yMax) ? yMax : newMax ;
557 }
558 break;
559
560 case kZAxis:
561 pMin=zMin;
562 pMax=zMax;
563 break;
564
565 default:
566 break;
567 }
568 pMin -= kCarTolerance ;
569 pMax += kCarTolerance ;
570
571 return true;
572 }
573 else
574 {
575 G4int i, noEntries, noBetweenSections4 ;
576 G4bool existsAfterClip = false ;
577
578 // Calculate rotated vertex coordinates
579
580 G4ThreeVectorList *vertices ;
581 G4int noPolygonVertices ; // will be 4
582 vertices = CreateRotatedVertices(pTransform,noPolygonVertices) ;
583
584 pMin = +kInfinity ;
585 pMax = -kInfinity ;
586
587 noEntries = vertices->size() ;
588 noBetweenSections4 = noEntries - noPolygonVertices ;
589
590 for (i=0;i<noEntries;i+=noPolygonVertices)
591 {
592 ClipCrossSection(vertices,i,pVoxelLimit,pAxis,pMin,pMax);
593 }
594 for (i=0;i<noBetweenSections4;i+=noPolygonVertices)
595 {
596 ClipBetweenSections(vertices,i,pVoxelLimit,pAxis,pMin,pMax);
597 }
598 if (pMin!=kInfinity||pMax!=-kInfinity)
599 {
600 existsAfterClip = true ; // Add 2*tolerance to avoid precision troubles
601 pMin -= kCarTolerance ;
602 pMax += kCarTolerance ;
603 }
604 else
605 {
606 // Check for case where completely enveloping clipping volume
607 // If point inside then we are confident that the solid completely
608 // envelopes the clipping volume. Hence set min/max extents according
609 // to clipping volume extents along the specified axis.
610
611 G4ThreeVector clipCentre(
612 (pVoxelLimit.GetMinXExtent()+pVoxelLimit.GetMaxXExtent())*0.5,
613 (pVoxelLimit.GetMinYExtent()+pVoxelLimit.GetMaxYExtent())*0.5,
614 (pVoxelLimit.GetMinZExtent()+pVoxelLimit.GetMaxZExtent())*0.5 ) ;
615
616 if (Inside(pTransform.Inverse().TransformPoint(clipCentre)) != kOutside )
617 {
618 existsAfterClip = true ;
619 pMin = pVoxelLimit.GetMinExtent(pAxis) ;
620 pMax = pVoxelLimit.GetMaxExtent(pAxis) ;
621 }
622 }
623 delete vertices;
624 return existsAfterClip;
625 }
626}
627
628//////////////////////////////////////////////////////////////////////////////
629//
630// Return whether point inside/outside/on surface
631
633{
634 G4double r2, pt2, pPhi, tolRMin, tolRMax ;
635
636 EInside in = kOutside ;
637 static const G4double halfAngTolerance = 0.5*kAngTolerance;
638
639 // General precals
640 r2 = p.x()*p.x() + p.y()*p.y() ;
641 pt2 = r2 + p.z()*p.z() + fRtor*fRtor - 2*fRtor*std::sqrt(r2) ;
642
643 if (fRmin) tolRMin = fRmin + fRminTolerance ;
644 else tolRMin = 0 ;
645
646 tolRMax = fRmax - fRmaxTolerance;
647
648 if (pt2 >= tolRMin*tolRMin && pt2 <= tolRMax*tolRMax )
649 {
650 if ( fDPhi == twopi || pt2 == 0 ) // on torus swept axis
651 {
652 in = kInside ;
653 }
654 else
655 {
656 // Try inner tolerant phi boundaries (=>inside)
657 // if not inside, try outer tolerant phi boundaries
658
659 pPhi = std::atan2(p.y(),p.x()) ;
660
661 if ( pPhi < -halfAngTolerance ) { pPhi += twopi ; } // 0<=pPhi<2pi
662 if ( fSPhi >= 0 )
663 {
664 if ( (std::fabs(pPhi) < halfAngTolerance)
665 && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
666 {
667 pPhi += twopi ; // 0 <= pPhi < 2pi
668 }
669 if ( (pPhi >= fSPhi + halfAngTolerance)
670 && (pPhi <= fSPhi + fDPhi - halfAngTolerance) )
671 {
672 in = kInside ;
673 }
674 else if ( (pPhi >= fSPhi - halfAngTolerance)
675 && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
676 {
677 in = kSurface ;
678 }
679 }
680 else // fSPhi < 0
681 {
682 if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
683 && (pPhi >= fSPhi + fDPhi + halfAngTolerance) ) {;}
684 else
685 {
686 in = kSurface ;
687 }
688 }
689 }
690 }
691 else // Try generous boundaries
692 {
693 tolRMin = fRmin - fRminTolerance ;
694 tolRMax = fRmax + fRmaxTolerance ;
695
696 if (tolRMin < 0 ) { tolRMin = 0 ; }
697
698 if ( (pt2 >= tolRMin*tolRMin) && (pt2 <= tolRMax*tolRMax) )
699 {
700 if ( (fDPhi == twopi) || (pt2 == 0) ) // Continuous in phi or on z-axis
701 {
702 in = kSurface ;
703 }
704 else // Try outer tolerant phi boundaries only
705 {
706 pPhi = std::atan2(p.y(),p.x()) ;
707
708 if ( pPhi < -halfAngTolerance ) { pPhi += twopi ; } // 0<=pPhi<2pi
709 if ( fSPhi >= 0 )
710 {
711 if ( (std::fabs(pPhi) < halfAngTolerance)
712 && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
713 {
714 pPhi += twopi ; // 0 <= pPhi < 2pi
715 }
716 if ( (pPhi >= fSPhi - halfAngTolerance)
717 && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
718 {
719 in = kSurface;
720 }
721 }
722 else // fSPhi < 0
723 {
724 if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
725 && (pPhi >= fSPhi + fDPhi + halfAngTolerance) ) {;}
726 else
727 {
728 in = kSurface ;
729 }
730 }
731 }
732 }
733 }
734 return in ;
735}
736
737/////////////////////////////////////////////////////////////////////////////
738//
739// Return unit normal of surface closest to p
740// - note if point on z axis, ignore phi divided sides
741// - unsafe if point close to z axis a rmin=0 - no explicit checks
742
744{
745 G4int noSurfaces = 0;
746 G4double rho2, rho, pt2, pt, pPhi;
747 G4double distRMin = kInfinity;
748 G4double distSPhi = kInfinity, distEPhi = kInfinity;
749
750 // To cope with precision loss
751 //
752 const G4double delta = std::max(10.0*kCarTolerance,
753 1.0e-8*(fRtor+fRmax));
754 const G4double dAngle = 10.0*kAngTolerance;
755
756 G4ThreeVector nR, nPs, nPe;
757 G4ThreeVector norm, sumnorm(0.,0.,0.);
758
759 rho2 = p.x()*p.x() + p.y()*p.y();
760 rho = std::sqrt(rho2);
761 pt2 = rho2+p.z()*p.z() +fRtor * (fRtor-2*rho);
762 pt2 = std::max(pt2, 0.0); // std::fabs(pt2);
763 pt = std::sqrt(pt2) ;
764
765 G4double distRMax = std::fabs(pt - fRmax);
766 if(fRmin) distRMin = std::fabs(pt - fRmin);
767
768 if( rho > delta && pt != 0.0 )
769 {
770 G4double redFactor= (rho-fRtor)/rho;
771 nR = G4ThreeVector( p.x()*redFactor, // p.x()*(1.-fRtor/rho),
772 p.y()*redFactor, // p.y()*(1.-fRtor/rho),
773 p.z() );
774 nR *= 1.0/pt;
775 }
776
777 if ( fDPhi < twopi ) // && rho ) // old limitation against (0,0,z)
778 {
779 if ( rho )
780 {
781 pPhi = std::atan2(p.y(),p.x());
782
783 if(pPhi < fSPhi-delta) { pPhi += twopi; }
784 else if(pPhi > fSPhi+fDPhi+delta) { pPhi -= twopi; }
785
786 distSPhi = std::fabs( pPhi - fSPhi );
787 distEPhi = std::fabs(pPhi-fSPhi-fDPhi);
788 }
789 nPs = G4ThreeVector(std::sin(fSPhi),-std::cos(fSPhi),0);
790 nPe = G4ThreeVector(-std::sin(fSPhi+fDPhi),std::cos(fSPhi+fDPhi),0);
791 }
792 if( distRMax <= delta )
793 {
794 noSurfaces ++;
795 sumnorm += nR;
796 }
797 else if( fRmin && (distRMin <= delta) ) // Must not be on both Outer and Inner
798 {
799 noSurfaces ++;
800 sumnorm -= nR;
801 }
802
803 // To be on one of the 'phi' surfaces,
804 // it must be within the 'tube' - with tolerance
805
806 if( (fDPhi < twopi) && (fRmin-delta <= pt) && (pt <= (fRmax+delta)) )
807 {
808 if (distSPhi <= dAngle)
809 {
810 noSurfaces ++;
811 sumnorm += nPs;
812 }
813 if (distEPhi <= dAngle)
814 {
815 noSurfaces ++;
816 sumnorm += nPe;
817 }
818 }
819 if ( noSurfaces == 0 )
820 {
821#ifdef G4CSGDEBUG
823 ed.precision(16);
824
825 EInside inIt= Inside( p );
826
827 if( inIt != kSurface )
828 {
829 ed << " ERROR> Surface Normal was called for Torus,"
830 << " with point not on surface." << G4endl;
831 }
832 else
833 {
834 ed << " ERROR> Surface Normal has not found a surface, "
835 << " despite the point being on the surface. " <<G4endl;
836 }
837
838 if( inIt != kInside)
839 {
840 ed << " Safety (Dist To In) = " << DistanceToIn(p) << G4endl;
841 }
842 if( inIt != kOutside)
843 {
844 ed << " Safety (Dist to Out) = " << DistanceToOut(p) << G4endl;
845 }
846 ed << " Coordinates of point : " << p << G4endl;
847 ed << " Parameters of solid : " << G4endl << *this << G4endl;
848
849 if( inIt == kSurface )
850 {
851 G4Exception("G4Torus::SurfaceNormal(p)", "GeomSolids1002",
852 JustWarning, ed,
853 "Failing to find normal, even though point is on surface!" );
854 }
855 else
856 {
857 static const char* NameInside[3]= { "Inside", "Surface", "Outside" };
858 ed << " The point is " << NameInside[inIt] << " the solid. "<< G4endl;
859 G4Exception("G4Torus::SurfaceNormal(p)", "GeomSolids1002",
860 JustWarning, ed, "Point p is not on surface !?" );
861 }
862#endif
863 norm = ApproxSurfaceNormal(p);
864 }
865 else if ( noSurfaces == 1 ) { norm = sumnorm; }
866 else { norm = sumnorm.unit(); }
867
868 // G4cout << "G4Torus::SurfaceNormal p= " << p << " returns norm= " << norm << G4endl;
869
870 return norm ;
871}
872
873//////////////////////////////////////////////////////////////////////////////
874//
875// Algorithm for SurfaceNormal() following the original specification
876// for points not on the surface
877
878G4ThreeVector G4Torus::ApproxSurfaceNormal( const G4ThreeVector& p ) const
879{
880 ENorm side ;
881 G4ThreeVector norm;
882 G4double rho2,rho,pt2,pt,phi;
883 G4double distRMin,distRMax,distSPhi,distEPhi,distMin;
884
885 rho2 = p.x()*p.x() + p.y()*p.y();
886 rho = std::sqrt(rho2) ;
887 pt2 = std::fabs(rho2+p.z()*p.z() +fRtor*fRtor - 2*fRtor*rho) ;
888 pt = std::sqrt(pt2) ;
889
890#ifdef G4CSGDEBUG
891 G4cout << " G4Torus::ApproximateSurfaceNormal called for point " << p
892 << G4endl;
893#endif
894
895 distRMax = std::fabs(pt - fRmax) ;
896
897 if(fRmin) // First minimum radius
898 {
899 distRMin = std::fabs(pt - fRmin) ;
900
901 if (distRMin < distRMax)
902 {
903 distMin = distRMin ;
904 side = kNRMin ;
905 }
906 else
907 {
908 distMin = distRMax ;
909 side = kNRMax ;
910 }
911 }
912 else
913 {
914 distMin = distRMax ;
915 side = kNRMax ;
916 }
917 if ( (fDPhi < twopi) && rho )
918 {
919 phi = std::atan2(p.y(),p.x()) ; // Protected against (0,0,z) (above rho!=0)
920
921 if (phi < 0) { phi += twopi ; }
922
923 if (fSPhi < 0 ) { distSPhi = std::fabs(phi-(fSPhi+twopi))*rho ; }
924 else { distSPhi = std::fabs(phi-fSPhi)*rho ; }
925
926 distEPhi = std::fabs(phi - fSPhi - fDPhi)*rho ;
927
928 if (distSPhi < distEPhi) // Find new minimum
929 {
930 if (distSPhi<distMin) side = kNSPhi ;
931 }
932 else
933 {
934 if (distEPhi < distMin) { side = kNEPhi ; }
935 }
936 }
937 switch (side)
938 {
939 case kNRMin: // Inner radius
940 norm = G4ThreeVector( -p.x()*(1-fRtor/rho)/pt,
941 -p.y()*(1-fRtor/rho)/pt,
942 -p.z()/pt ) ;
943 break ;
944 case kNRMax: // Outer radius
945 norm = G4ThreeVector( p.x()*(1-fRtor/rho)/pt,
946 p.y()*(1-fRtor/rho)/pt,
947 p.z()/pt ) ;
948 break;
949 case kNSPhi:
950 norm = G4ThreeVector(std::sin(fSPhi),-std::cos(fSPhi),0) ;
951 break;
952 case kNEPhi:
953 norm = G4ThreeVector(-std::sin(fSPhi+fDPhi),std::cos(fSPhi+fDPhi),0) ;
954 break;
955 default: // Should never reach this case ...
956 DumpInfo();
957 G4Exception("G4Torus::ApproxSurfaceNormal()",
958 "GeomSolids1002", JustWarning,
959 "Undefined side for valid surface normal to solid.");
960 break ;
961 }
962 return norm ;
963}
964
965///////////////////////////////////////////////////////////////////////
966//
967// Calculate distance to shape from outside, along normalised vector
968// - return kInfinity if no intersection, or intersection distance <= tolerance
969//
970// - Compute the intersection with the z planes
971// - if at valid r, phi, return
972//
973// -> If point is outer outer radius, compute intersection with rmax
974// - if at valid phi,z return
975//
976// -> Compute intersection with inner radius, taking largest +ve root
977// - if valid (phi), save intersction
978//
979// -> If phi segmented, compute intersections with phi half planes
980// - return smallest of valid phi intersections and
981// inner radius intersection
982//
983// NOTE:
984// - Precalculations for phi trigonometry are Done `just in time'
985// - `if valid' implies tolerant checking of intersection points
986
988 const G4ThreeVector& v ) const
989{
990
991 G4double snxt=kInfinity, sphi=kInfinity; // snxt = default return value
992
993 G4double sd[4] ;
994
995 // Precalculated trig for phi intersections - used by r,z intersections to
996 // check validity
997
998 G4bool seg; // true if segmented
999 G4double hDPhi; // half dphi
1000 G4double cPhi,sinCPhi=0.,cosCPhi=0.; // central phi
1001
1002 G4double tolORMin2; // `generous' radii squared
1003 G4double tolORMax2;
1004
1005 static const G4double halfCarTolerance = 0.5*kCarTolerance;
1006
1007 G4double Dist,xi,yi,zi,rhoi2,it2; // Intersection point variables
1008
1009 G4double Comp;
1010 G4double cosSPhi,sinSPhi; // Trig for phi start intersect
1011 G4double ePhi,cosEPhi,sinEPhi; // for phi end intersect
1012
1013 // Set phi divided flag and precalcs
1014 //
1015 if ( fDPhi < twopi )
1016 {
1017 seg = true ;
1018 hDPhi = 0.5*fDPhi ; // half delta phi
1019 cPhi = fSPhi + hDPhi ;
1020 sinCPhi = std::sin(cPhi) ;
1021 cosCPhi = std::cos(cPhi) ;
1022 }
1023 else
1024 {
1025 seg = false ;
1026 }
1027
1028 if (fRmin > fRminTolerance) // Calculate tolerant rmin and rmax
1029 {
1030 tolORMin2 = (fRmin - fRminTolerance)*(fRmin - fRminTolerance) ;
1031 }
1032 else
1033 {
1034 tolORMin2 = 0 ;
1035 }
1036 tolORMax2 = (fRmax + fRmaxTolerance)*(fRmax + fRmaxTolerance) ;
1037
1038 // Intersection with Rmax (possible return) and Rmin (must also check phi)
1039
1040 G4double Rtor2 = fRtor*fRtor ;
1041
1042 snxt = SolveNumericJT(p,v,fRmax,true);
1043
1044 if (fRmin) // Possible Rmin intersection
1045 {
1046 sd[0] = SolveNumericJT(p,v,fRmin,true);
1047 if ( sd[0] < snxt ) { snxt = sd[0] ; }
1048 }
1049
1050 //
1051 // Phi segment intersection
1052 //
1053 // o Tolerant of points inside phi planes by up to kCarTolerance*0.5
1054 //
1055 // o NOTE: Large duplication of code between sphi & ephi checks
1056 // -> only diffs: sphi -> ephi, Comp -> -Comp and half-plane
1057 // intersection check <=0 -> >=0
1058 // -> use some form of loop Construct ?
1059
1060 if (seg)
1061 {
1062 sinSPhi = std::sin(fSPhi) ; // First phi surface ('S'tarting phi)
1063 cosSPhi = std::cos(fSPhi) ;
1064 Comp = v.x()*sinSPhi - v.y()*cosSPhi ; // Component in outwards
1065 // normal direction
1066 if (Comp < 0 )
1067 {
1068 Dist = (p.y()*cosSPhi - p.x()*sinSPhi) ;
1069
1070 if (Dist < halfCarTolerance)
1071 {
1072 sphi = Dist/Comp ;
1073 if (sphi < snxt)
1074 {
1075 if ( sphi < 0 ) { sphi = 0 ; }
1076
1077 xi = p.x() + sphi*v.x() ;
1078 yi = p.y() + sphi*v.y() ;
1079 zi = p.z() + sphi*v.z() ;
1080 rhoi2 = xi*xi + yi*yi ;
1081 it2 = std::fabs(rhoi2 + zi*zi + Rtor2 - 2*fRtor*std::sqrt(rhoi2)) ;
1082
1083 if ( it2 >= tolORMin2 && it2 <= tolORMax2 )
1084 {
1085 // r intersection is good - check intersecting
1086 // with correct half-plane
1087 //
1088 if ((yi*cosCPhi-xi*sinCPhi)<=0) { snxt=sphi; }
1089 }
1090 }
1091 }
1092 }
1093 ePhi=fSPhi+fDPhi; // Second phi surface ('E'nding phi)
1094 sinEPhi=std::sin(ePhi);
1095 cosEPhi=std::cos(ePhi);
1096 Comp=-(v.x()*sinEPhi-v.y()*cosEPhi);
1097
1098 if ( Comp < 0 ) // Component in outwards normal dirn
1099 {
1100 Dist = -(p.y()*cosEPhi - p.x()*sinEPhi) ;
1101
1102 if (Dist < halfCarTolerance )
1103 {
1104 sphi = Dist/Comp ;
1105
1106 if (sphi < snxt )
1107 {
1108 if (sphi < 0 ) { sphi = 0 ; }
1109
1110 xi = p.x() + sphi*v.x() ;
1111 yi = p.y() + sphi*v.y() ;
1112 zi = p.z() + sphi*v.z() ;
1113 rhoi2 = xi*xi + yi*yi ;
1114 it2 = std::fabs(rhoi2 + zi*zi + Rtor2 - 2*fRtor*std::sqrt(rhoi2)) ;
1115
1116 if (it2 >= tolORMin2 && it2 <= tolORMax2)
1117 {
1118 // z and r intersections good - check intersecting
1119 // with correct half-plane
1120 //
1121 if ((yi*cosCPhi-xi*sinCPhi)>=0) { snxt=sphi; }
1122 }
1123 }
1124 }
1125 }
1126 }
1127 if(snxt < halfCarTolerance) { snxt = 0.0 ; }
1128
1129 return snxt ;
1130}
1131
1132/////////////////////////////////////////////////////////////////////////////
1133//
1134// Calculate distance (<= actual) to closest surface of shape from outside
1135// - Calculate distance to z, radial planes
1136// - Only to phi planes if outside phi extent
1137// - Return 0 if point inside
1138
1140{
1141 G4double safe=0.0, safe1, safe2 ;
1142 G4double phiC, cosPhiC, sinPhiC, safePhi, ePhi, cosPsi ;
1143 G4double rho2, rho, pt2, pt ;
1144
1145 rho2 = p.x()*p.x() + p.y()*p.y() ;
1146 rho = std::sqrt(rho2) ;
1147 pt2 = std::fabs(rho2 + p.z()*p.z() + fRtor*fRtor - 2*fRtor*rho) ;
1148 pt = std::sqrt(pt2) ;
1149
1150 safe1 = fRmin - pt ;
1151 safe2 = pt - fRmax ;
1152
1153 if (safe1 > safe2) { safe = safe1; }
1154 else { safe = safe2; }
1155
1156 if ( fDPhi < twopi && rho )
1157 {
1158 phiC = fSPhi + fDPhi*0.5 ;
1159 cosPhiC = std::cos(phiC) ;
1160 sinPhiC = std::sin(phiC) ;
1161 cosPsi = (p.x()*cosPhiC + p.y()*sinPhiC)/rho ;
1162
1163 if (cosPsi < std::cos(fDPhi*0.5) ) // Psi=angle from central phi to point
1164 { // Point lies outside phi range
1165 if ((p.y()*cosPhiC - p.x()*sinPhiC) <= 0 )
1166 {
1167 safePhi = std::fabs(p.x()*std::sin(fSPhi) - p.y()*std::cos(fSPhi)) ;
1168 }
1169 else
1170 {
1171 ePhi = fSPhi + fDPhi ;
1172 safePhi = std::fabs(p.x()*std::sin(ePhi) - p.y()*std::cos(ePhi)) ;
1173 }
1174 if (safePhi > safe) { safe = safePhi ; }
1175 }
1176 }
1177 if (safe < 0 ) { safe = 0 ; }
1178 return safe;
1179}
1180
1181///////////////////////////////////////////////////////////////////////////
1182//
1183// Calculate distance to surface of shape from `inside', allowing for tolerance
1184// - Only Calc rmax intersection if no valid rmin intersection
1185//
1186
1188 const G4ThreeVector& v,
1189 const G4bool calcNorm,
1190 G4bool *validNorm,
1191 G4ThreeVector *n ) const
1192{
1193 ESide side = kNull, sidephi = kNull ;
1194 G4double snxt = kInfinity, sphi, sd[4] ;
1195
1196 static const G4double halfCarTolerance = 0.5*kCarTolerance;
1197 static const G4double halfAngTolerance = 0.5*kAngTolerance;
1198
1199 // Vars for phi intersection
1200 //
1201 G4double sinSPhi, cosSPhi, ePhi, sinEPhi, cosEPhi;
1202 G4double cPhi, sinCPhi, cosCPhi ;
1203 G4double pDistS, compS, pDistE, compE, sphi2, xi, yi, zi, vphi ;
1204
1205 // Radial Intersections Defenitions & General Precals
1206
1207 //////////////////////// new calculation //////////////////////
1208
1209#if 1
1210
1211 // This is the version with the calculation of CalcNorm = true
1212 // To be done: Check the precision of this calculation.
1213 // If you want return always validNorm = false, then take the version below
1214
1215 G4double rho2 = p.x()*p.x()+p.y()*p.y();
1216 G4double rho = std::sqrt(rho2) ;
1217
1218 G4double pt2 = rho2 + p.z()*p.z() + fRtor * (fRtor - 2.0*rho);
1219 // Regroup for slightly better FP accuracy
1220
1221 if( pt2 < 0.0)
1222 {
1223 pt2= std::fabs( pt2 );
1224 }
1225
1226 G4double pt = std::sqrt(pt2) ;
1227
1228 G4double pDotV = p.x()*v.x() + p.y()*v.y() + p.z()*v.z() ;
1229
1230 G4double tolRMax = fRmax - fRmaxTolerance ;
1231
1232 G4double vDotNmax = pDotV - fRtor*(v.x()*p.x() + v.y()*p.y())/rho ;
1233 G4double pDotxyNmax = (1 - fRtor/rho) ;
1234
1235 if( (pt2 > tolRMax*tolRMax) && (vDotNmax >= 0) )
1236 {
1237 // On tolerant boundary & heading outwards (or perpendicular to) outer
1238 // radial surface -> leaving immediately with *n for really convex part
1239 // only
1240
1241 if ( calcNorm && (pDotxyNmax >= -2.*fRmaxTolerance) )
1242 {
1243 *n = G4ThreeVector( p.x()*(1 - fRtor/rho)/pt,
1244 p.y()*(1 - fRtor/rho)/pt,
1245 p.z()/pt ) ;
1246 *validNorm = true ;
1247 }
1248
1249 return snxt = 0 ; // Leaving by Rmax immediately
1250 }
1251
1252 snxt = SolveNumericJT(p,v,fRmax,false);
1253 side = kRMax ;
1254
1255 // rmin
1256
1257 if ( fRmin )
1258 {
1259 G4double tolRMin = fRmin + fRminTolerance ;
1260
1261 if ( (pt2 < tolRMin*tolRMin) && (vDotNmax < 0) )
1262 {
1263 if (calcNorm) { *validNorm = false ; } // Concave surface of the torus
1264 return snxt = 0 ; // Leaving by Rmin immediately
1265 }
1266
1267 sd[0] = SolveNumericJT(p,v,fRmin,false);
1268 if ( sd[0] < snxt )
1269 {
1270 snxt = sd[0] ;
1271 side = kRMin ;
1272 }
1273 }
1274
1275#else
1276
1277 // this is the "conservative" version which return always validnorm = false
1278 // NOTE: using this version the unit test testG4Torus will break
1279
1280 snxt = SolveNumericJT(p,v,fRmax,false);
1281 side = kRMax ;
1282
1283 if ( fRmin )
1284 {
1285 sd[0] = SolveNumericJT(p,v,fRmin,false);
1286 if ( sd[0] < snxt )
1287 {
1288 snxt = sd[0] ;
1289 side = kRMin ;
1290 }
1291 }
1292
1293 if ( calcNorm && (snxt == 0.0) )
1294 {
1295 *validNorm = false ; // Leaving solid, but possible re-intersection
1296 return snxt ;
1297 }
1298
1299#endif
1300
1301 if (fDPhi < twopi) // Phi Intersections
1302 {
1303 sinSPhi = std::sin(fSPhi) ;
1304 cosSPhi = std::cos(fSPhi) ;
1305 ePhi = fSPhi + fDPhi ;
1306 sinEPhi = std::sin(ePhi) ;
1307 cosEPhi = std::cos(ePhi) ;
1308 cPhi = fSPhi + fDPhi*0.5 ;
1309 sinCPhi = std::sin(cPhi) ;
1310 cosCPhi = std::cos(cPhi) ;
1311
1312 // angle calculation with correction
1313 // of difference in domain of atan2 and Sphi
1314 //
1315 vphi = std::atan2(v.y(),v.x()) ;
1316
1317 if ( vphi < fSPhi - halfAngTolerance ) { vphi += twopi; }
1318 else if ( vphi > ePhi + halfAngTolerance ) { vphi -= twopi; }
1319
1320 if ( p.x() || p.y() ) // Check if on z axis (rho not needed later)
1321 {
1322 pDistS = p.x()*sinSPhi - p.y()*cosSPhi ; // pDist -ve when inside
1323 pDistE = -p.x()*sinEPhi + p.y()*cosEPhi ;
1324
1325 // Comp -ve when in direction of outwards normal
1326 //
1327 compS = -sinSPhi*v.x() + cosSPhi*v.y() ;
1328 compE = sinEPhi*v.x() - cosEPhi*v.y() ;
1329 sidephi = kNull ;
1330
1331 if( ( (fDPhi <= pi) && ( (pDistS <= halfCarTolerance)
1332 && (pDistE <= halfCarTolerance) ) )
1333 || ( (fDPhi > pi) && !((pDistS > halfCarTolerance)
1334 && (pDistE > halfCarTolerance) ) ) )
1335 {
1336 // Inside both phi *full* planes
1337
1338 if ( compS < 0 )
1339 {
1340 sphi = pDistS/compS ;
1341
1342 if (sphi >= -halfCarTolerance)
1343 {
1344 xi = p.x() + sphi*v.x() ;
1345 yi = p.y() + sphi*v.y() ;
1346
1347 // Check intersecting with correct half-plane
1348 // (if not -> no intersect)
1349 //
1350 if ( (std::fabs(xi)<=kCarTolerance)
1351 && (std::fabs(yi)<=kCarTolerance) )
1352 {
1353 sidephi = kSPhi;
1354 if ( ((fSPhi-halfAngTolerance)<=vphi)
1355 && ((ePhi+halfAngTolerance)>=vphi) )
1356 {
1357 sphi = kInfinity;
1358 }
1359 }
1360 else if ( yi*cosCPhi-xi*sinCPhi >=0 )
1361 {
1362 sphi = kInfinity ;
1363 }
1364 else
1365 {
1366 sidephi = kSPhi ;
1367 }
1368 }
1369 else
1370 {
1371 sphi = kInfinity ;
1372 }
1373 }
1374 else
1375 {
1376 sphi = kInfinity ;
1377 }
1378
1379 if ( compE < 0 )
1380 {
1381 sphi2 = pDistE/compE ;
1382
1383 // Only check further if < starting phi intersection
1384 //
1385 if ( (sphi2 > -kCarTolerance) && (sphi2 < sphi) )
1386 {
1387 xi = p.x() + sphi2*v.x() ;
1388 yi = p.y() + sphi2*v.y() ;
1389
1390 if ( (std::fabs(xi)<=kCarTolerance)
1391 && (std::fabs(yi)<=kCarTolerance) )
1392 {
1393 // Leaving via ending phi
1394 //
1395 if( !( (fSPhi-halfAngTolerance <= vphi)
1396 && (ePhi+halfAngTolerance >= vphi) ) )
1397 {
1398 sidephi = kEPhi ;
1399 sphi = sphi2;
1400 }
1401 }
1402 else // Check intersecting with correct half-plane
1403 {
1404 if ( (yi*cosCPhi-xi*sinCPhi) >= 0)
1405 {
1406 // Leaving via ending phi
1407 //
1408 sidephi = kEPhi ;
1409 sphi = sphi2;
1410
1411 }
1412 }
1413 }
1414 }
1415 }
1416 else
1417 {
1418 sphi = kInfinity ;
1419 }
1420 }
1421 else
1422 {
1423 // On z axis + travel not || to z axis -> if phi of vector direction
1424 // within phi of shape, Step limited by rmax, else Step =0
1425
1426 vphi = std::atan2(v.y(),v.x());
1427
1428 if ( ( fSPhi-halfAngTolerance <= vphi ) &&
1429 ( vphi <= ( ePhi+halfAngTolerance ) ) )
1430 {
1431 sphi = kInfinity;
1432 }
1433 else
1434 {
1435 sidephi = kSPhi ; // arbitrary
1436 sphi=0;
1437 }
1438 }
1439
1440 // Order intersections
1441
1442 if (sphi<snxt)
1443 {
1444 snxt=sphi;
1445 side=sidephi;
1446 }
1447 }
1448
1449 G4double rhoi2,rhoi,it2,it,iDotxyNmax ;
1450 // Note: by numerical computation we know where the ray hits the torus
1451 // So I propose to return the side where the ray hits
1452
1453 if (calcNorm)
1454 {
1455 switch(side)
1456 {
1457 case kRMax: // n is unit vector
1458 xi = p.x() + snxt*v.x() ;
1459 yi =p.y() + snxt*v.y() ;
1460 zi = p.z() + snxt*v.z() ;
1461 rhoi2 = xi*xi + yi*yi ;
1462 rhoi = std::sqrt(rhoi2) ;
1463 it2 = std::fabs(rhoi2 + zi*zi + fRtor*fRtor - 2*fRtor*rhoi) ;
1464 it = std::sqrt(it2) ;
1465 iDotxyNmax = (1-fRtor/rhoi) ;
1466 if(iDotxyNmax >= -2.*fRmaxTolerance) // really convex part of Rmax
1467 {
1468 *n = G4ThreeVector( xi*(1-fRtor/rhoi)/it,
1469 yi*(1-fRtor/rhoi)/it,
1470 zi/it ) ;
1471 *validNorm = true ;
1472 }
1473 else
1474 {
1475 *validNorm = false ; // concave-convex part of Rmax
1476 }
1477 break ;
1478
1479 case kRMin:
1480 *validNorm = false ; // Rmin is concave or concave-convex
1481 break;
1482
1483 case kSPhi:
1484 if (fDPhi <= pi )
1485 {
1486 *n=G4ThreeVector(std::sin(fSPhi),-std::cos(fSPhi),0);
1487 *validNorm=true;
1488 }
1489 else
1490 {
1491 *validNorm = false ;
1492 }
1493 break ;
1494
1495 case kEPhi:
1496 if (fDPhi <= pi)
1497 {
1498 *n=G4ThreeVector(-std::sin(fSPhi+fDPhi),std::cos(fSPhi+fDPhi),0);
1499 *validNorm=true;
1500 }
1501 else
1502 {
1503 *validNorm = false ;
1504 }
1505 break;
1506
1507 default:
1508
1509 // It seems we go here from time to time ...
1510
1511 G4cout << G4endl;
1512 DumpInfo();
1513 std::ostringstream message;
1514 G4int oldprc = message.precision(16);
1515 message << "Undefined side for valid surface normal to solid."
1516 << G4endl
1517 << "Position:" << G4endl << G4endl
1518 << "p.x() = " << p.x()/mm << " mm" << G4endl
1519 << "p.y() = " << p.y()/mm << " mm" << G4endl
1520 << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl
1521 << "Direction:" << G4endl << G4endl
1522 << "v.x() = " << v.x() << G4endl
1523 << "v.y() = " << v.y() << G4endl
1524 << "v.z() = " << v.z() << G4endl << G4endl
1525 << "Proposed distance :" << G4endl << G4endl
1526 << "snxt = " << snxt/mm << " mm" << G4endl;
1527 message.precision(oldprc);
1528 G4Exception("G4Torus::DistanceToOut(p,v,..)",
1529 "GeomSolids1002",JustWarning, message);
1530 break;
1531 }
1532 }
1533 if ( snxt<halfCarTolerance ) { snxt=0 ; }
1534
1535 return snxt;
1536}
1537
1538/////////////////////////////////////////////////////////////////////////
1539//
1540// Calculate distance (<=actual) to closest surface of shape from inside
1541
1543{
1544 G4double safe=0.0,safeR1,safeR2;
1545 G4double rho2,rho,pt2,pt ;
1546 G4double safePhi,phiC,cosPhiC,sinPhiC,ePhi;
1547 rho2 = p.x()*p.x() + p.y()*p.y() ;
1548 rho = std::sqrt(rho2) ;
1549 pt2 = std::fabs(rho2 + p.z()*p.z() + fRtor*fRtor - 2*fRtor*rho) ;
1550 pt = std::sqrt(pt2) ;
1551
1552#ifdef G4CSGDEBUG
1553 if( Inside(p) == kOutside )
1554 {
1555 G4int oldprc = G4cout.precision(16) ;
1556 G4cout << G4endl ;
1557 DumpInfo();
1558 G4cout << "Position:" << G4endl << G4endl ;
1559 G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl ;
1560 G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl ;
1561 G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl ;
1562 G4cout.precision(oldprc);
1563 G4Exception("G4Torus::DistanceToOut(p)", "GeomSolids1002",
1564 JustWarning, "Point p is outside !?" );
1565 }
1566#endif
1567
1568 if (fRmin)
1569 {
1570 safeR1 = pt - fRmin ;
1571 safeR2 = fRmax - pt ;
1572
1573 if (safeR1 < safeR2) { safe = safeR1 ; }
1574 else { safe = safeR2 ; }
1575 }
1576 else
1577 {
1578 safe = fRmax - pt ;
1579 }
1580
1581 // Check if phi divided, Calc distances closest phi plane
1582 //
1583 if (fDPhi<twopi) // Above/below central phi of Torus?
1584 {
1585 phiC = fSPhi + fDPhi*0.5 ;
1586 cosPhiC = std::cos(phiC) ;
1587 sinPhiC = std::sin(phiC) ;
1588
1589 if ((p.y()*cosPhiC-p.x()*sinPhiC)<=0)
1590 {
1591 safePhi = -(p.x()*std::sin(fSPhi) - p.y()*std::cos(fSPhi)) ;
1592 }
1593 else
1594 {
1595 ePhi = fSPhi + fDPhi ;
1596 safePhi = (p.x()*std::sin(ePhi) - p.y()*std::cos(ePhi)) ;
1597 }
1598 if (safePhi < safe) { safe = safePhi ; }
1599 }
1600 if (safe < 0) { safe = 0 ; }
1601 return safe ;
1602}
1603
1604/////////////////////////////////////////////////////////////////////////////
1605//
1606// Create a List containing the transformed vertices
1607// Ordering [0-3] -fRtor cross section
1608// [4-7] +fRtor cross section such that [0] is below [4],
1609// [1] below [5] etc.
1610// Note:
1611// Caller has deletion resposibility
1612// Potential improvement: For last slice, use actual ending angle
1613// to avoid rounding error problems.
1614
1616G4Torus::CreateRotatedVertices( const G4AffineTransform& pTransform,
1617 G4int& noPolygonVertices ) const
1618{
1619 G4ThreeVectorList *vertices;
1620 G4ThreeVector vertex0,vertex1,vertex2,vertex3;
1621 G4double meshAngle,meshRMax,crossAngle,cosCrossAngle,sinCrossAngle,sAngle;
1622 G4double rMaxX,rMaxY,rMinX,rMinY;
1623 G4int crossSection,noCrossSections;
1624
1625 // Compute no of cross-sections necessary to mesh tube
1626 //
1627 noCrossSections = G4int (fDPhi/kMeshAngleDefault) + 1 ;
1628
1629 if (noCrossSections < kMinMeshSections)
1630 {
1631 noCrossSections = kMinMeshSections ;
1632 }
1633 else if (noCrossSections>kMaxMeshSections)
1634 {
1635 noCrossSections=kMaxMeshSections;
1636 }
1637 meshAngle = fDPhi/(noCrossSections - 1) ;
1638 meshRMax = (fRtor + fRmax)/std::cos(meshAngle*0.5) ;
1639
1640 // If complete in phi, set start angle such that mesh will be at fRmax
1641 // on the x axis. Will give better extent calculations when not rotated
1642
1643 if ( (fDPhi == twopi) && (fSPhi == 0) )
1644 {
1645 sAngle = -meshAngle*0.5 ;
1646 }
1647 else
1648 {
1649 sAngle = fSPhi ;
1650 }
1651 vertices = new G4ThreeVectorList();
1652
1653 if (vertices)
1654 {
1655 vertices->reserve(noCrossSections*4) ;
1656 for (crossSection=0;crossSection<noCrossSections;crossSection++)
1657 {
1658 // Compute coordinates of cross section at section crossSection
1659
1660 crossAngle=sAngle+crossSection*meshAngle;
1661 cosCrossAngle=std::cos(crossAngle);
1662 sinCrossAngle=std::sin(crossAngle);
1663
1664 rMaxX=meshRMax*cosCrossAngle;
1665 rMaxY=meshRMax*sinCrossAngle;
1666 rMinX=(fRtor-fRmax)*cosCrossAngle;
1667 rMinY=(fRtor-fRmax)*sinCrossAngle;
1668 vertex0=G4ThreeVector(rMinX,rMinY,-fRmax);
1669 vertex1=G4ThreeVector(rMaxX,rMaxY,-fRmax);
1670 vertex2=G4ThreeVector(rMaxX,rMaxY,+fRmax);
1671 vertex3=G4ThreeVector(rMinX,rMinY,+fRmax);
1672
1673 vertices->push_back(pTransform.TransformPoint(vertex0));
1674 vertices->push_back(pTransform.TransformPoint(vertex1));
1675 vertices->push_back(pTransform.TransformPoint(vertex2));
1676 vertices->push_back(pTransform.TransformPoint(vertex3));
1677 }
1678 noPolygonVertices = 4 ;
1679 }
1680 else
1681 {
1682 DumpInfo();
1683 G4Exception("G4Torus::CreateRotatedVertices()",
1684 "GeomSolids0003", FatalException,
1685 "Error in allocation of vertices. Out of memory !");
1686 }
1687 return vertices;
1688}
1689
1690//////////////////////////////////////////////////////////////////////////
1691//
1692// Stream object contents to an output stream
1693
1695{
1696 return G4String("G4Torus");
1697}
1698
1699//////////////////////////////////////////////////////////////////////////
1700//
1701// Make a clone of the object
1702//
1704{
1705 return new G4Torus(*this);
1706}
1707
1708//////////////////////////////////////////////////////////////////////////
1709//
1710// Stream object contents to an output stream
1711
1712std::ostream& G4Torus::StreamInfo( std::ostream& os ) const
1713{
1714 G4int oldprc = os.precision(16);
1715 os << "-----------------------------------------------------------\n"
1716 << " *** Dump for solid - " << GetName() << " ***\n"
1717 << " ===================================================\n"
1718 << " Solid type: G4Torus\n"
1719 << " Parameters: \n"
1720 << " inner radius: " << fRmin/mm << " mm \n"
1721 << " outer radius: " << fRmax/mm << " mm \n"
1722 << " swept radius: " << fRtor/mm << " mm \n"
1723 << " starting phi: " << fSPhi/degree << " degrees \n"
1724 << " delta phi : " << fDPhi/degree << " degrees \n"
1725 << "-----------------------------------------------------------\n";
1726 os.precision(oldprc);
1727
1728 return os;
1729}
1730
1731////////////////////////////////////////////////////////////////////////////
1732//
1733// GetPointOnSurface
1734
1736{
1737 G4double cosu, sinu,cosv, sinv, aOut, aIn, aSide, chose, phi, theta, rRand;
1738
1739 phi = RandFlat::shoot(fSPhi,fSPhi+fDPhi);
1740 theta = RandFlat::shoot(0.,twopi);
1741
1742 cosu = std::cos(phi); sinu = std::sin(phi);
1743 cosv = std::cos(theta); sinv = std::sin(theta);
1744
1745 // compute the areas
1746
1747 aOut = (fDPhi)*twopi*fRtor*fRmax;
1748 aIn = (fDPhi)*twopi*fRtor*fRmin;
1749 aSide = pi*(fRmax*fRmax-fRmin*fRmin);
1750
1751 if ((fSPhi == 0) && (fDPhi == twopi)){ aSide = 0; }
1752 chose = RandFlat::shoot(0.,aOut + aIn + 2.*aSide);
1753
1754 if(chose < aOut)
1755 {
1756 return G4ThreeVector ((fRtor+fRmax*cosv)*cosu,
1757 (fRtor+fRmax*cosv)*sinu, fRmax*sinv);
1758 }
1759 else if( (chose >= aOut) && (chose < aOut + aIn) )
1760 {
1761 return G4ThreeVector ((fRtor+fRmin*cosv)*cosu,
1762 (fRtor+fRmin*cosv)*sinu, fRmin*sinv);
1763 }
1764 else if( (chose >= aOut + aIn) && (chose < aOut + aIn + aSide) )
1765 {
1766 rRand = GetRadiusInRing(fRmin,fRmax);
1767 return G4ThreeVector ((fRtor+rRand*cosv)*std::cos(fSPhi),
1768 (fRtor+rRand*cosv)*std::sin(fSPhi), rRand*sinv);
1769 }
1770 else
1771 {
1772 rRand = GetRadiusInRing(fRmin,fRmax);
1773 return G4ThreeVector ((fRtor+rRand*cosv)*std::cos(fSPhi+fDPhi),
1774 (fRtor+rRand*cosv)*std::sin(fSPhi+fDPhi),
1775 rRand*sinv);
1776 }
1777}
1778
1779///////////////////////////////////////////////////////////////////////
1780//
1781// Visualisation Functions
1782
1784{
1785 scene.AddSolid (*this);
1786}
1787
1789{
1790 return new G4PolyhedronTorus (fRmin, fRmax, fRtor, fSPhi, fDPhi);
1791}
1792
1794{
1795 G4NURBS* pNURBS;
1796 if (fRmin != 0)
1797 {
1798 if (fDPhi >= twopi)
1799 {
1800 pNURBS = new G4NURBStube(fRmin, fRmax, fRtor);
1801 }
1802 else
1803 {
1804 pNURBS = new G4NURBStubesector(fRmin, fRmax, fRtor, fSPhi, fSPhi + fDPhi);
1805 }
1806 }
1807 else
1808 {
1809 if (fDPhi >= twopi)
1810 {
1811 pNURBS = new G4NURBScylinder (fRmax, fRtor);
1812 }
1813 else
1814 {
1815 const G4double epsilon = 1.e-4; // Cylinder sector not yet available!
1816 pNURBS = new G4NURBStubesector (epsilon, fRmax, fRtor,
1817 fSPhi, fSPhi + fDPhi);
1818 }
1819 }
1820 return pNURBS;
1821}
ENorm
Definition: G4Cons.cc:72
@ JustWarning
@ FatalException
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:85
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
double z() const
Hep3Vector unit() const
double x() const
double y() const
static double shoot()
Definition: RandFlat.cc:59
G4bool IsRotated() const
G4AffineTransform Inverse() const
G4ThreeVector NetTranslation() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4Polyhedron * fpPolyhedron
Definition: G4CSGSolid.hh:80
G4double fSurfaceArea
Definition: G4CSGSolid.hh:79
G4double fCubicVolume
Definition: G4CSGSolid.hh:78
G4double GetRadiusInRing(G4double rmin, G4double rmax) const
Definition: G4CSGSolid.cc:101
G4CSGSolid & operator=(const G4CSGSolid &rhs)
Definition: G4CSGSolid.cc:82
G4double GetRadialTolerance() const
static G4GeometryTolerance * GetInstance()
G4double GetAngularTolerance() const
G4int FindRoots(G4double *op, G4int degree, G4double *zeror, G4double *zeroi)
EInside Inside(const G4ThreeVector &p) const
Definition: G4Torus.cc:632
G4Torus & operator=(const G4Torus &rhs)
Definition: G4Torus.cc:212
G4GeometryType GetEntityType() const
Definition: G4Torus.cc:1694
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Torus.cc:1187
~G4Torus()
Definition: G4Torus.cc:193
G4VSolid * Clone() const
Definition: G4Torus.cc:1703
G4double GetRmin() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Torus.cc:412
void SetAllParameters(G4double pRmin, G4double pRmax, G4double pRtor, G4double pSPhi, G4double pDPhi)
Definition: G4Torus.cc:100
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Torus.cc:237
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Torus.cc:743
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Torus.cc:987
G4ThreeVector GetPointOnSurface() const
Definition: G4Torus.cc:1735
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Torus.cc:1783
G4NURBS * CreateNURBS() const
Definition: G4Torus.cc:1793
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Torus.cc:1712
G4Torus(const G4String &pName, G4double pRmin, G4double pRmax, G4double pRtor, G4double pSPhi, G4double pDPhi)
Definition: G4Torus.cc:84
G4Polyhedron * CreatePolyhedron() const
Definition: G4Torus.cc:1788
virtual void AddSolid(const G4Box &)=0
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4String GetName() const
void ClipBetweenSections(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:376
void DumpInfo() const
G4double kCarTolerance
Definition: G4VSolid.hh:307
void ClipCrossSection(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:345
G4double GetMinExtent(const EAxis pAxis) const
G4bool IsYLimited() const
G4double GetMinZExtent() const
G4bool IsXLimited() const
G4double GetMaxExtent(const EAxis pAxis) const
G4double GetMaxYExtent() const
G4double GetMaxZExtent() const
G4double GetMinYExtent() const
G4double GetMinXExtent() const
G4bool IsZLimited() const
G4double GetMaxXExtent() const
EAxis
Definition: geomdefs.hh:54
@ kYAxis
Definition: geomdefs.hh:54
@ kXAxis
Definition: geomdefs.hh:54
@ kZAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:58
@ kInside
Definition: geomdefs.hh:58
@ kOutside
Definition: geomdefs.hh:58
@ kSurface
Definition: geomdefs.hh:58
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4int kMaxMeshSections
Definition: meshdefs.hh:46
const G4int kMinMeshSections
Definition: meshdefs.hh:45
const G4double kMeshAngleDefault
Definition: meshdefs.hh:42
Definition: DoubConv.h:17