Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Tubs.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// G4Tubs implementation
27//
28// 1994-95 P.Kent: first implementation
29// 08.08.00 V.Grichine: more stable roots of 2-equation in DistanceToOut(p,v,..)
30// 07.12.00 V.Grichine: phi-section algorithm was changed in Inside(p)
31// 03.05.05 V.Grichine: SurfaceNormal(p) according to J.Apostolakis proposal
32// 24.08.16 E.Tcherniaev: reimplemented CalculateExtent().
33// --------------------------------------------------------------------
34
35#include "G4Tubs.hh"
36
37#if !defined(G4GEOM_USE_UTUBS)
38
39#include "G4GeomTools.hh"
40#include "G4VoxelLimits.hh"
41#include "G4AffineTransform.hh"
43#include "G4BoundingEnvelope.hh"
44
46#include "G4QuickRand.hh"
47
48#include "G4VGraphicsScene.hh"
49#include "G4Polyhedron.hh"
50#include "G4AutoLock.hh"
51
52namespace
53{
54 G4Mutex tubsMutex = G4MUTEX_INITIALIZER;
55}
56
57using namespace CLHEP;
58
59// Private enums: Not for external use
60namespace {
61// Used by distanceToOut
62enum ESide {kNull,kRMin,kRMax,kSPhi,kEPhi,kPZ,kMZ};
63
64// Used by normal
65enum ENorm {kNRMin,kNRMax,kNSPhi,kNEPhi,kNZ};
66}
67
68/////////////////////////////////////////////////////////////////////////
69//
70// Constructor - check parameters, convert angles so 0<sphi+dpshi<=2_PI
71// - note if pdphi>2PI then reset to 2PI
72
74 G4double pRMin, G4double pRMax,
75 G4double pDz,
76 G4double pSPhi, G4double pDPhi )
77 : G4CSGSolid(pName), fRMin(pRMin), fRMax(pRMax), fDz(pDz),
78 fSPhi(0), fDPhi(0),
79 fInvRmax( pRMax > 0.0 ? 1.0/pRMax : 0.0 ),
80 fInvRmin( pRMin > 0.0 ? 1.0/pRMin : 0.0 )
81{
84
88
89 if (pDz<=0) // Check z-len
90 {
91 std::ostringstream message;
92 message << "Negative Z half-length (" << pDz << ") in solid: " << GetName();
93 G4Exception("G4Tubs::G4Tubs()", "GeomSolids0002", FatalException, message);
94 }
95 if ( (pRMin >= pRMax) || (pRMin < 0) ) // Check radii
96 {
97 std::ostringstream message;
98 message << "Invalid values for radii in solid: " << GetName()
99 << G4endl
100 << " pRMin = " << pRMin << ", pRMax = " << pRMax;
101 G4Exception("G4Tubs::G4Tubs()", "GeomSolids0002", FatalException, message);
102 }
103
104 // Check angles
105 //
106 CheckPhiAngles(pSPhi, pDPhi);
107}
108
109///////////////////////////////////////////////////////////////////////
110//
111// Fake default constructor - sets only member data and allocates memory
112// for usage restricted to object persistency.
113//
114G4Tubs::G4Tubs( __void__& a )
115 : G4CSGSolid(a)
116{
117}
118
119//////////////////////////////////////////////////////////////////////////
120//
121// Assignment operator
122
124{
125 // Check assignment to self
126 //
127 if (this == &rhs) { return *this; }
128
129 // Copy base class data
130 //
132
133 // Copy data
134 //
136 fRMin = rhs.fRMin; fRMax = rhs.fRMax; fDz = rhs.fDz;
137 fSPhi = rhs.fSPhi; fDPhi = rhs.fDPhi;
138 sinCPhi = rhs.sinCPhi; cosCPhi = rhs.cosCPhi; cosHDPhi = rhs.cosHDPhi;
140 sinSPhi = rhs.sinSPhi; cosSPhi = rhs.cosSPhi;
141 sinEPhi = rhs.sinEPhi; cosEPhi = rhs.cosEPhi;
143 fInvRmax = rhs.fInvRmax;
144 fInvRmin = rhs.fInvRmin;
148
149 return *this;
150}
151
152/////////////////////////////////////////////////////////////////////////
153//
154// Dispatch to parameterisation for replication mechanism dimension
155// computation & modification.
156
158 const G4int n,
159 const G4VPhysicalVolume* pRep )
160{
161 p->ComputeDimensions(*this,n,pRep) ;
162}
163
164/////////////////////////////////////////////////////////////////////////
165//
166// Get bounding box
167
169{
170 G4double rmin = GetInnerRadius();
171 G4double rmax = GetOuterRadius();
173
174 // Find bounding box
175 //
176 if (GetDeltaPhiAngle() < twopi)
177 {
178 G4TwoVector vmin,vmax;
179 G4GeomTools::DiskExtent(rmin,rmax,
182 vmin,vmax);
183 pMin.set(vmin.x(),vmin.y(),-dz);
184 pMax.set(vmax.x(),vmax.y(), dz);
185 }
186 else
187 {
188 pMin.set(-rmax,-rmax,-dz);
189 pMax.set( rmax, rmax, dz);
190 }
191
192 // Check correctness of the bounding box
193 //
194 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
195 {
196 std::ostringstream message;
197 message << "Bad bounding box (min >= max) for solid: "
198 << GetName() << " !"
199 << "\npMin = " << pMin
200 << "\npMax = " << pMax;
201 G4Exception("G4Tubs::BoundingLimits()", "GeomMgt0001",
202 JustWarning, message);
203 DumpInfo();
204 }
205}
206
207/////////////////////////////////////////////////////////////////////////
208//
209// Calculate extent under transform and specified limit
210
212 const G4VoxelLimits& pVoxelLimit,
213 const G4AffineTransform& pTransform,
214 G4double& pMin,
215 G4double& pMax ) const
216{
217 G4ThreeVector bmin, bmax;
218 G4bool exist;
219
220 // Get bounding box
221 BoundingLimits(bmin,bmax);
222
223 // Check bounding box
224 G4BoundingEnvelope bbox(bmin,bmax);
225#ifdef G4BBOX_EXTENT
226 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
227#endif
228 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
229 {
230 return exist = pMin < pMax;
231 }
232
233 // Get parameters of the solid
234 G4double rmin = GetInnerRadius();
235 G4double rmax = GetOuterRadius();
237 G4double dphi = GetDeltaPhiAngle();
238
239 // Find bounding envelope and calculate extent
240 //
241 const G4int NSTEPS = 24; // number of steps for whole circle
242 G4double astep = twopi/NSTEPS; // max angle for one step
243 G4int ksteps = (dphi <= astep) ? 1 : (G4int)((dphi-deg)/astep) + 1;
244 G4double ang = dphi/ksteps;
245
246 G4double sinHalf = std::sin(0.5*ang);
247 G4double cosHalf = std::cos(0.5*ang);
248 G4double sinStep = 2.*sinHalf*cosHalf;
249 G4double cosStep = 1. - 2.*sinHalf*sinHalf;
250 G4double rext = rmax/cosHalf;
251
252 // bounding envelope for full cylinder consists of two polygons,
253 // in other cases it is a sequence of quadrilaterals
254 if (rmin == 0 && dphi == twopi)
255 {
256 G4double sinCur = sinHalf;
257 G4double cosCur = cosHalf;
258
259 G4ThreeVectorList baseA(NSTEPS),baseB(NSTEPS);
260 for (G4int k=0; k<NSTEPS; ++k)
261 {
262 baseA[k].set(rext*cosCur,rext*sinCur,-dz);
263 baseB[k].set(rext*cosCur,rext*sinCur, dz);
264
265 G4double sinTmp = sinCur;
266 sinCur = sinCur*cosStep + cosCur*sinStep;
267 cosCur = cosCur*cosStep - sinTmp*sinStep;
268 }
269 std::vector<const G4ThreeVectorList *> polygons(2);
270 polygons[0] = &baseA;
271 polygons[1] = &baseB;
272 G4BoundingEnvelope benv(bmin,bmax,polygons);
273 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
274 }
275 else
276 {
277 G4double sinStart = GetSinStartPhi();
278 G4double cosStart = GetCosStartPhi();
279 G4double sinEnd = GetSinEndPhi();
280 G4double cosEnd = GetCosEndPhi();
281 G4double sinCur = sinStart*cosHalf + cosStart*sinHalf;
282 G4double cosCur = cosStart*cosHalf - sinStart*sinHalf;
283
284 // set quadrilaterals
285 G4ThreeVectorList pols[NSTEPS+2];
286 for (G4int k=0; k<ksteps+2; ++k)
287 {
288 pols[k].resize(4);
289 }
290 pols[0][0].set(rmin*cosStart,rmin*sinStart, dz);
291 pols[0][1].set(rmin*cosStart,rmin*sinStart,-dz);
292 pols[0][2].set(rmax*cosStart,rmax*sinStart,-dz);
293 pols[0][3].set(rmax*cosStart,rmax*sinStart, dz);
294 for (G4int k=1; k<ksteps+1; ++k)
295 {
296 pols[k][0].set(rmin*cosCur,rmin*sinCur, dz);
297 pols[k][1].set(rmin*cosCur,rmin*sinCur,-dz);
298 pols[k][2].set(rext*cosCur,rext*sinCur,-dz);
299 pols[k][3].set(rext*cosCur,rext*sinCur, dz);
300
301 G4double sinTmp = sinCur;
302 sinCur = sinCur*cosStep + cosCur*sinStep;
303 cosCur = cosCur*cosStep - sinTmp*sinStep;
304 }
305 pols[ksteps+1][0].set(rmin*cosEnd,rmin*sinEnd, dz);
306 pols[ksteps+1][1].set(rmin*cosEnd,rmin*sinEnd,-dz);
307 pols[ksteps+1][2].set(rmax*cosEnd,rmax*sinEnd,-dz);
308 pols[ksteps+1][3].set(rmax*cosEnd,rmax*sinEnd, dz);
309
310 // set envelope and calculate extent
311 std::vector<const G4ThreeVectorList *> polygons;
312 polygons.resize(ksteps+2);
313 for (G4int k=0; k<ksteps+2; ++k)
314 {
315 polygons[k] = &pols[k];
316 }
317 G4BoundingEnvelope benv(bmin,bmax,polygons);
318 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
319 }
320 return exist;
321}
322
323///////////////////////////////////////////////////////////////////////////
324//
325// Return whether point inside/outside/on surface
326
328{
329 G4double r2,pPhi,tolRMin,tolRMax;
330 EInside in = kOutside ;
331
332 if (std::fabs(p.z()) <= fDz - halfCarTolerance)
333 {
334 r2 = p.x()*p.x() + p.y()*p.y() ;
335
336 if (fRMin != 0.0) { tolRMin = fRMin + halfRadTolerance ; }
337 else { tolRMin = 0 ; }
338
339 tolRMax = fRMax - halfRadTolerance ;
340
341 if ((r2 >= tolRMin*tolRMin) && (r2 <= tolRMax*tolRMax))
342 {
343 if ( fPhiFullTube )
344 {
345 in = kInside ;
346 }
347 else
348 {
349 // Try inner tolerant phi boundaries (=>inside)
350 // if not inside, try outer tolerant phi boundaries
351
352 if ( (tolRMin==0) && (std::fabs(p.x())<=halfCarTolerance)
353 && (std::fabs(p.y())<=halfCarTolerance) )
354 {
355 in=kSurface;
356 }
357 else
358 {
359 pPhi = std::atan2(p.y(),p.x()) ;
360 if ( pPhi < -halfAngTolerance ) { pPhi += twopi; } // 0<=pPhi<2pi
361
362 if ( fSPhi >= 0 )
363 {
364 if ( (std::fabs(pPhi) < halfAngTolerance)
365 && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
366 {
367 pPhi += twopi ; // 0 <= pPhi < 2pi
368 }
369 if ( (pPhi >= fSPhi + halfAngTolerance)
370 && (pPhi <= fSPhi + fDPhi - halfAngTolerance) )
371 {
372 in = kInside ;
373 }
374 else if ( (pPhi >= fSPhi - halfAngTolerance)
375 && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
376 {
377 in = kSurface ;
378 }
379 }
380 else // fSPhi < 0
381 {
382 if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
383 && (pPhi >= fSPhi + fDPhi + halfAngTolerance) ) {;} //kOutside
384 else if ( (pPhi <= fSPhi + twopi + halfAngTolerance)
385 && (pPhi >= fSPhi + fDPhi - halfAngTolerance) )
386 {
387 in = kSurface ;
388 }
389 else
390 {
391 in = kInside ;
392 }
393 }
394 }
395 }
396 }
397 else // Try generous boundaries
398 {
399 tolRMin = fRMin - halfRadTolerance ;
400 tolRMax = fRMax + halfRadTolerance ;
401
402 if ( tolRMin < 0 ) { tolRMin = 0; }
403
404 if ( (r2 >= tolRMin*tolRMin) && (r2 <= tolRMax*tolRMax) )
405 {
407 { // Continuous in phi or on z-axis
408 in = kSurface ;
409 }
410 else // Try outer tolerant phi boundaries only
411 {
412 pPhi = std::atan2(p.y(),p.x()) ;
413
414 if ( pPhi < -halfAngTolerance) { pPhi += twopi; } // 0<=pPhi<2pi
415 if ( fSPhi >= 0 )
416 {
417 if ( (std::fabs(pPhi) < halfAngTolerance)
418 && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
419 {
420 pPhi += twopi ; // 0 <= pPhi < 2pi
421 }
422 if ( (pPhi >= fSPhi - halfAngTolerance)
423 && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
424 {
425 in = kSurface ;
426 }
427 }
428 else // fSPhi < 0
429 {
430 if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
431 && (pPhi >= fSPhi + fDPhi + halfAngTolerance) ) {;} // kOutside
432 else
433 {
434 in = kSurface ;
435 }
436 }
437 }
438 }
439 }
440 }
441 else if (std::fabs(p.z()) <= fDz + halfCarTolerance)
442 { // Check within tolerant r limits
443 r2 = p.x()*p.x() + p.y()*p.y() ;
444 tolRMin = fRMin - halfRadTolerance ;
445 tolRMax = fRMax + halfRadTolerance ;
446
447 if ( tolRMin < 0 ) { tolRMin = 0; }
448
449 if ( (r2 >= tolRMin*tolRMin) && (r2 <= tolRMax*tolRMax) )
450 {
452 { // Continuous in phi or on z-axis
453 in = kSurface ;
454 }
455 else // Try outer tolerant phi boundaries
456 {
457 pPhi = std::atan2(p.y(),p.x()) ;
458
459 if ( pPhi < -halfAngTolerance ) { pPhi += twopi; } // 0<=pPhi<2pi
460 if ( fSPhi >= 0 )
461 {
462 if ( (std::fabs(pPhi) < halfAngTolerance)
463 && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
464 {
465 pPhi += twopi ; // 0 <= pPhi < 2pi
466 }
467 if ( (pPhi >= fSPhi - halfAngTolerance)
468 && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
469 {
470 in = kSurface;
471 }
472 }
473 else // fSPhi < 0
474 {
475 if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
476 && (pPhi >= fSPhi + fDPhi + halfAngTolerance) ) {;}
477 else
478 {
479 in = kSurface ;
480 }
481 }
482 }
483 }
484 }
485 return in;
486}
487
488///////////////////////////////////////////////////////////////////////////
489//
490// Return unit normal of surface closest to p
491// - note if point on z axis, ignore phi divided sides
492// - unsafe if point close to z axis a rmin=0 - no explicit checks
493
495{
496 G4int noSurfaces = 0;
497 G4double rho, pPhi;
498 G4double distZ, distRMin, distRMax;
499 G4double distSPhi = kInfinity, distEPhi = kInfinity;
500
501 G4ThreeVector norm, sumnorm(0.,0.,0.);
502 G4ThreeVector nZ = G4ThreeVector(0, 0, 1.0);
503 G4ThreeVector nR, nPs, nPe;
504
505 rho = std::sqrt(p.x()*p.x() + p.y()*p.y());
506
507 distRMin = std::fabs(rho - fRMin);
508 distRMax = std::fabs(rho - fRMax);
509 distZ = std::fabs(std::fabs(p.z()) - fDz);
510
511 if (!fPhiFullTube) // Protected against (0,0,z)
512 {
513 if ( rho > halfCarTolerance )
514 {
515 pPhi = std::atan2(p.y(),p.x());
516
517 if (pPhi < fSPhi-halfCarTolerance) { pPhi += twopi; }
518 else if (pPhi > fSPhi+fDPhi+halfCarTolerance) { pPhi -= twopi; }
519
520 distSPhi = std::fabs( pPhi - fSPhi );
521 distEPhi = std::fabs( pPhi - fSPhi - fDPhi );
522 }
523 else if ( fRMin == 0.0 )
524 {
525 distSPhi = 0.;
526 distEPhi = 0.;
527 }
528 nPs = G4ThreeVector( sinSPhi, -cosSPhi, 0 );
529 nPe = G4ThreeVector( -sinEPhi, cosEPhi, 0 );
530 }
531 if ( rho > halfCarTolerance ) { nR = G4ThreeVector(p.x()/rho,p.y()/rho,0); }
532
533 if( distRMax <= halfCarTolerance )
534 {
535 ++noSurfaces;
536 sumnorm += nR;
537 }
538 if( (fRMin != 0.0) && (distRMin <= halfCarTolerance) )
539 {
540 ++noSurfaces;
541 sumnorm -= nR;
542 }
543 if( fDPhi < twopi )
544 {
545 if (distSPhi <= halfAngTolerance)
546 {
547 ++noSurfaces;
548 sumnorm += nPs;
549 }
550 if (distEPhi <= halfAngTolerance)
551 {
552 ++noSurfaces;
553 sumnorm += nPe;
554 }
555 }
556 if (distZ <= halfCarTolerance)
557 {
558 ++noSurfaces;
559 if ( p.z() >= 0.) { sumnorm += nZ; }
560 else { sumnorm -= nZ; }
561 }
562 if ( noSurfaces == 0 )
563 {
564#ifdef G4CSGDEBUG
565 G4Exception("G4Tubs::SurfaceNormal(p)", "GeomSolids1002",
566 JustWarning, "Point p is not on surface !?" );
567 G4long oldprc = G4cout.precision(20);
568 G4cout<< "G4Tubs::SN ( "<<p.x()<<", "<<p.y()<<", "<<p.z()<<" ); "
569 << G4endl << G4endl;
570 G4cout.precision(oldprc) ;
571#endif
572 norm = ApproxSurfaceNormal(p);
573 }
574 else if ( noSurfaces == 1 ) { norm = sumnorm; }
575 else { norm = sumnorm.unit(); }
576
577 return norm;
578}
579
580/////////////////////////////////////////////////////////////////////////////
581//
582// Algorithm for SurfaceNormal() following the original specification
583// for points not on the surface
584
586{
587 ENorm side ;
588 G4ThreeVector norm ;
589 G4double rho, phi ;
590 G4double distZ, distRMin, distRMax, distSPhi, distEPhi, distMin ;
591
592 rho = std::sqrt(p.x()*p.x() + p.y()*p.y()) ;
593
594 distRMin = std::fabs(rho - fRMin) ;
595 distRMax = std::fabs(rho - fRMax) ;
596 distZ = std::fabs(std::fabs(p.z()) - fDz) ;
597
598 if (distRMin < distRMax) // First minimum
599 {
600 if ( distZ < distRMin )
601 {
602 distMin = distZ ;
603 side = kNZ ;
604 }
605 else
606 {
607 distMin = distRMin ;
608 side = kNRMin ;
609 }
610 }
611 else
612 {
613 if ( distZ < distRMax )
614 {
615 distMin = distZ ;
616 side = kNZ ;
617 }
618 else
619 {
620 distMin = distRMax ;
621 side = kNRMax ;
622 }
623 }
624 if (!fPhiFullTube && (rho != 0.0) ) // Protected against (0,0,z)
625 {
626 phi = std::atan2(p.y(),p.x()) ;
627
628 if ( phi < 0 ) { phi += twopi; }
629
630 if ( fSPhi < 0 )
631 {
632 distSPhi = std::fabs(phi - (fSPhi + twopi))*rho ;
633 }
634 else
635 {
636 distSPhi = std::fabs(phi - fSPhi)*rho ;
637 }
638 distEPhi = std::fabs(phi - fSPhi - fDPhi)*rho ;
639
640 if (distSPhi < distEPhi) // Find new minimum
641 {
642 if ( distSPhi < distMin )
643 {
644 side = kNSPhi ;
645 }
646 }
647 else
648 {
649 if ( distEPhi < distMin )
650 {
651 side = kNEPhi ;
652 }
653 }
654 }
655 switch ( side )
656 {
657 case kNRMin : // Inner radius
658 {
659 norm = G4ThreeVector(-p.x()/rho, -p.y()/rho, 0) ;
660 break ;
661 }
662 case kNRMax : // Outer radius
663 {
664 norm = G4ThreeVector(p.x()/rho, p.y()/rho, 0) ;
665 break ;
666 }
667 case kNZ : // + or - dz
668 {
669 if ( p.z() > 0 ) { norm = G4ThreeVector(0,0,1) ; }
670 else { norm = G4ThreeVector(0,0,-1); }
671 break ;
672 }
673 case kNSPhi:
674 {
675 norm = G4ThreeVector(sinSPhi, -cosSPhi, 0) ;
676 break ;
677 }
678 case kNEPhi:
679 {
680 norm = G4ThreeVector(-sinEPhi, cosEPhi, 0) ;
681 break;
682 }
683 default: // Should never reach this case ...
684 {
685 DumpInfo();
686 G4Exception("G4Tubs::ApproxSurfaceNormal()",
687 "GeomSolids1002", JustWarning,
688 "Undefined side for valid surface normal to solid.");
689 break ;
690 }
691 }
692 return norm;
693}
694
695////////////////////////////////////////////////////////////////////
696//
697//
698// Calculate distance to shape from outside, along normalised vector
699// - return kInfinity if no intersection, or intersection distance <= tolerance
700//
701// - Compute the intersection with the z planes
702// - if at valid r, phi, return
703//
704// -> If point is outer outer radius, compute intersection with rmax
705// - if at valid phi,z return
706//
707// -> Compute intersection with inner radius, taking largest +ve root
708// - if valid (in z,phi), save intersction
709//
710// -> If phi segmented, compute intersections with phi half planes
711// - return smallest of valid phi intersections and
712// inner radius intersection
713//
714// NOTE:
715// - 'if valid' implies tolerant checking of intersection points
716
718 const G4ThreeVector& v ) const
719{
720 G4double snxt = kInfinity ; // snxt = default return value
721 G4double tolORMin2, tolIRMax2 ; // 'generous' radii squared
722 G4double tolORMax2, tolIRMin2, tolODz, tolIDz ;
723 const G4double dRmax = 100.*fRMax;
724
725 // Intersection point variables
726 //
727 G4double Dist, sd, xi, yi, zi, rho2, inum, iden, cosPsi, Comp ;
728 G4double t1, t2, t3, b, c, d ; // Quadratic solver variables
729
730 // Calculate tolerant rmin and rmax
731
732 if (fRMin > kRadTolerance)
733 {
734 tolORMin2 = (fRMin - halfRadTolerance)*(fRMin - halfRadTolerance) ;
735 tolIRMin2 = (fRMin + halfRadTolerance)*(fRMin + halfRadTolerance) ;
736 }
737 else
738 {
739 tolORMin2 = 0.0 ;
740 tolIRMin2 = 0.0 ;
741 }
742 tolORMax2 = (fRMax + halfRadTolerance)*(fRMax + halfRadTolerance) ;
743 tolIRMax2 = (fRMax - halfRadTolerance)*(fRMax - halfRadTolerance) ;
744
745 // Intersection with Z surfaces
746
747 tolIDz = fDz - halfCarTolerance ;
748 tolODz = fDz + halfCarTolerance ;
749
750 if (std::fabs(p.z()) >= tolIDz)
751 {
752 if ( p.z()*v.z() < 0 ) // at +Z going in -Z or visa versa
753 {
754 sd = (std::fabs(p.z()) - fDz)/std::fabs(v.z()) ; // Z intersect distance
755
756 if(sd < 0.0) { sd = 0.0; }
757
758 xi = p.x() + sd*v.x() ; // Intersection coords
759 yi = p.y() + sd*v.y() ;
760 rho2 = xi*xi + yi*yi ;
761
762 // Check validity of intersection
763
764 if ((tolIRMin2 <= rho2) && (rho2 <= tolIRMax2))
765 {
766 if (!fPhiFullTube && (rho2 != 0.0))
767 {
768 // Psi = angle made with central (average) phi of shape
769 //
770 inum = xi*cosCPhi + yi*sinCPhi ;
771 iden = std::sqrt(rho2) ;
772 cosPsi = inum/iden ;
773 if (cosPsi >= cosHDPhiIT) { return sd ; }
774 }
775 else
776 {
777 return sd ;
778 }
779 }
780 }
781 else
782 {
783 if ( snxt<halfCarTolerance ) { snxt=0; }
784 return snxt ; // On/outside extent, and heading away
785 // -> cannot intersect
786 }
787 }
788
789 // -> Can not intersect z surfaces
790 //
791 // Intersection with rmax (possible return) and rmin (must also check phi)
792 //
793 // Intersection point (xi,yi,zi) on line x=p.x+t*v.x etc.
794 //
795 // Intersects with x^2+y^2=R^2
796 //
797 // Hence (v.x^2+v.y^2)t^2+ 2t(p.x*v.x+p.y*v.y)+p.x^2+p.y^2-R^2=0
798 // t1 t2 t3
799
800 t1 = 1.0 - v.z()*v.z() ;
801 t2 = p.x()*v.x() + p.y()*v.y() ;
802 t3 = p.x()*p.x() + p.y()*p.y() ;
803
804 if ( t1 > 0 ) // Check not || to z axis
805 {
806 b = t2/t1 ;
807 c = t3 - fRMax*fRMax ;
808 if ((t3 >= tolORMax2) && (t2<0)) // This also handles the tangent case
809 {
810 // Try outer cylinder intersection
811 // c=(t3-fRMax*fRMax)/t1;
812
813 c /= t1 ;
814 d = b*b - c ;
815
816 if (d >= 0) // If real root
817 {
818 sd = c/(-b+std::sqrt(d));
819 if (sd >= 0) // If 'forwards'
820 {
821 if ( sd>dRmax ) // Avoid rounding errors due to precision issues on
822 { // 64 bits systems. Split long distances and recompute
823 G4double fTerm = sd-std::fmod(sd,dRmax);
824 sd = fTerm + DistanceToIn(p+fTerm*v,v);
825 }
826 // Check z intersection
827 //
828 zi = p.z() + sd*v.z() ;
829 if (std::fabs(zi)<=tolODz)
830 {
831 // Z ok. Check phi intersection if reqd
832 //
833 if (fPhiFullTube)
834 {
835 return sd ;
836 }
837 xi = p.x() + sd*v.x() ;
838 yi = p.y() + sd*v.y() ;
839 cosPsi = (xi*cosCPhi + yi*sinCPhi)/fRMax ;
840 if (cosPsi >= cosHDPhiIT)
841 {
842 return sd ;
843 }
844 } // end if std::fabs(zi)
845 } // end if (sd>=0)
846 } // end if (d>=0)
847 } // end if (r>=fRMax)
848 else
849 {
850 // Inside outer radius :
851 // check not inside, and heading through tubs (-> 0 to in)
852
853 if ((t3 > tolIRMin2) && (t2 < 0) && (std::fabs(p.z()) <= tolIDz))
854 {
855 // Inside both radii, delta r -ve, inside z extent
856
857 if (!fPhiFullTube)
858 {
859 inum = p.x()*cosCPhi + p.y()*sinCPhi ;
860 iden = std::sqrt(t3) ;
861 cosPsi = inum/iden ;
862 if (cosPsi >= cosHDPhiIT)
863 {
864 // In the old version, the small negative tangent for the point
865 // on surface was not taken in account, and returning 0.0 ...
866 // New version: check the tangent for the point on surface and
867 // if no intersection, return kInfinity, if intersection instead
868 // return sd.
869 //
870 c = t3-fRMax*fRMax;
871 if ( c<=0.0 )
872 {
873 return 0.0;
874 }
875
876 c = c/t1 ;
877 d = b*b-c;
878 if ( d>=0.0 )
879 {
880 snxt = c/(-b+std::sqrt(d)); // using safe solution
881 // for quadratic equation
882 if ( snxt < halfCarTolerance ) { snxt=0; }
883 return snxt ;
884 }
885
886 return kInfinity;
887 }
888 }
889 else
890 {
891 // In the old version, the small negative tangent for the point
892 // on surface was not taken in account, and returning 0.0 ...
893 // New version: check the tangent for the point on surface and
894 // if no intersection, return kInfinity, if intersection instead
895 // return sd.
896 //
897 c = t3 - fRMax*fRMax;
898 if ( c<=0.0 )
899 {
900 return 0.0;
901 }
902
903 c = c/t1 ;
904 d = b*b-c;
905 if ( d>=0.0 )
906 {
907 snxt= c/(-b+std::sqrt(d)); // using safe solution
908 // for quadratic equation
909 if ( snxt < halfCarTolerance ) { snxt=0; }
910 return snxt ;
911 }
912
913 return kInfinity;
914 } // end if (!fPhiFullTube)
915 } // end if (t3>tolIRMin2)
916 } // end if (Inside Outer Radius)
917 if ( fRMin != 0.0 ) // Try inner cylinder intersection
918 {
919 c = (t3 - fRMin*fRMin)/t1 ;
920 d = b*b - c ;
921 if ( d >= 0.0 ) // If real root
922 {
923 // Always want 2nd root - we are outside and know rmax Hit was bad
924 // - If on surface of rmin also need farthest root
925
926 sd =( b > 0. )? c/(-b - std::sqrt(d)) : (-b + std::sqrt(d));
927 if (sd >= -halfCarTolerance) // check forwards
928 {
929 // Check z intersection
930 //
931 if(sd < 0.0) { sd = 0.0; }
932 if ( sd>dRmax ) // Avoid rounding errors due to precision issues seen
933 { // 64 bits systems. Split long distances and recompute
934 G4double fTerm = sd-std::fmod(sd,dRmax);
935 sd = fTerm + DistanceToIn(p+fTerm*v,v);
936 }
937 zi = p.z() + sd*v.z() ;
938 if (std::fabs(zi) <= tolODz)
939 {
940 // Z ok. Check phi
941 //
942 if ( fPhiFullTube )
943 {
944 return sd ;
945 }
946
947 xi = p.x() + sd*v.x() ;
948 yi = p.y() + sd*v.y() ;
949 cosPsi = (xi*cosCPhi + yi*sinCPhi)*fInvRmin;
950 if (cosPsi >= cosHDPhiIT)
951 {
952 // Good inner radius isect
953 // - but earlier phi isect still possible
954
955 snxt = sd ;
956 }
957 } // end if std::fabs(zi)
958 } // end if (sd>=0)
959 } // end if (d>=0)
960 } // end if (fRMin)
961 }
962
963 // Phi segment intersection
964 //
965 // o Tolerant of points inside phi planes by up to kCarTolerance*0.5
966 //
967 // o NOTE: Large duplication of code between sphi & ephi checks
968 // -> only diffs: sphi -> ephi, Comp -> -Comp and half-plane
969 // intersection check <=0 -> >=0
970 // -> use some form of loop Construct ?
971 //
972 if ( !fPhiFullTube )
973 {
974 // First phi surface (Starting phi)
975 //
976 Comp = v.x()*sinSPhi - v.y()*cosSPhi ;
977
978 if ( Comp < 0 ) // Component in outwards normal dirn
979 {
980 Dist = (p.y()*cosSPhi - p.x()*sinSPhi) ;
981
982 if ( Dist < halfCarTolerance )
983 {
984 sd = Dist/Comp ;
985
986 if (sd < snxt)
987 {
988 if ( sd < 0 ) { sd = 0.0; }
989 zi = p.z() + sd*v.z() ;
990 if ( std::fabs(zi) <= tolODz )
991 {
992 xi = p.x() + sd*v.x() ;
993 yi = p.y() + sd*v.y() ;
994 rho2 = xi*xi + yi*yi ;
995
996 if ( ( (rho2 >= tolIRMin2) && (rho2 <= tolIRMax2) )
997 || ( (rho2 > tolORMin2) && (rho2 < tolIRMin2)
998 && ( v.y()*cosSPhi - v.x()*sinSPhi > 0 )
999 && ( v.x()*cosSPhi + v.y()*sinSPhi >= 0 ) )
1000 || ( (rho2 > tolIRMax2) && (rho2 < tolORMax2)
1001 && (v.y()*cosSPhi - v.x()*sinSPhi > 0)
1002 && (v.x()*cosSPhi + v.y()*sinSPhi < 0) ) )
1003 {
1004 // z and r intersections good
1005 // - check intersecting with correct half-plane
1006 //
1007 if ((yi*cosCPhi-xi*sinCPhi) <= halfCarTolerance) { snxt = sd; }
1008 }
1009 }
1010 }
1011 }
1012 }
1013
1014 // Second phi surface (Ending phi)
1015
1016 Comp = -(v.x()*sinEPhi - v.y()*cosEPhi) ;
1017
1018 if (Comp < 0 ) // Component in outwards normal dirn
1019 {
1020 Dist = -(p.y()*cosEPhi - p.x()*sinEPhi) ;
1021
1022 if ( Dist < halfCarTolerance )
1023 {
1024 sd = Dist/Comp ;
1025
1026 if (sd < snxt)
1027 {
1028 if ( sd < 0 ) { sd = 0; }
1029 zi = p.z() + sd*v.z() ;
1030 if ( std::fabs(zi) <= tolODz )
1031 {
1032 xi = p.x() + sd*v.x() ;
1033 yi = p.y() + sd*v.y() ;
1034 rho2 = xi*xi + yi*yi ;
1035 if ( ( (rho2 >= tolIRMin2) && (rho2 <= tolIRMax2) )
1036 || ( (rho2 > tolORMin2) && (rho2 < tolIRMin2)
1037 && (v.x()*sinEPhi - v.y()*cosEPhi > 0)
1038 && (v.x()*cosEPhi + v.y()*sinEPhi >= 0) )
1039 || ( (rho2 > tolIRMax2) && (rho2 < tolORMax2)
1040 && (v.x()*sinEPhi - v.y()*cosEPhi > 0)
1041 && (v.x()*cosEPhi + v.y()*sinEPhi < 0) ) )
1042 {
1043 // z and r intersections good
1044 // - check intersecting with correct half-plane
1045 //
1046 if ( (yi*cosCPhi-xi*sinCPhi) >= 0 ) { snxt = sd; }
1047 } //?? >=-halfCarTolerance
1048 }
1049 }
1050 }
1051 } // Comp < 0
1052 } // !fPhiFullTube
1053 if ( snxt<halfCarTolerance ) { snxt=0; }
1054 return snxt ;
1055}
1056
1057//////////////////////////////////////////////////////////////////
1058//
1059// Calculate distance to shape from outside, along normalised vector
1060// - return kInfinity if no intersection, or intersection distance <= tolerance
1061//
1062// - Compute the intersection with the z planes
1063// - if at valid r, phi, return
1064//
1065// -> If point is outer outer radius, compute intersection with rmax
1066// - if at valid phi,z return
1067//
1068// -> Compute intersection with inner radius, taking largest +ve root
1069// - if valid (in z,phi), save intersction
1070//
1071// -> If phi segmented, compute intersections with phi half planes
1072// - return smallest of valid phi intersections and
1073// inner radius intersection
1074//
1075// NOTE:
1076// - Precalculations for phi trigonometry are Done `just in time'
1077// - `if valid' implies tolerant checking of intersection points
1078// Calculate distance (<= actual) to closest surface of shape from outside
1079// - Calculate distance to z, radial planes
1080// - Only to phi planes if outside phi extent
1081// - Return 0 if point inside
1082
1084{
1085 G4double safe=0.0, rho, safe1, safe2, safe3 ;
1086 G4double safePhi, cosPsi ;
1087
1088 rho = std::sqrt(p.x()*p.x() + p.y()*p.y()) ;
1089 safe1 = fRMin - rho ;
1090 safe2 = rho - fRMax ;
1091 safe3 = std::fabs(p.z()) - fDz ;
1092
1093 if ( safe1 > safe2 ) { safe = safe1; }
1094 else { safe = safe2; }
1095 if ( safe3 > safe ) { safe = safe3; }
1096
1097 if ( (!fPhiFullTube) && ((rho) != 0.0) )
1098 {
1099 // Psi=angle from central phi to point
1100 //
1101 cosPsi = (p.x()*cosCPhi + p.y()*sinCPhi)/rho ;
1102
1103 if ( cosPsi < cosHDPhi )
1104 {
1105 // Point lies outside phi range
1106
1107 if ( (p.y()*cosCPhi - p.x()*sinCPhi) <= 0 )
1108 {
1109 safePhi = std::fabs(p.x()*sinSPhi - p.y()*cosSPhi) ;
1110 }
1111 else
1112 {
1113 safePhi = std::fabs(p.x()*sinEPhi - p.y()*cosEPhi) ;
1114 }
1115 if ( safePhi > safe ) { safe = safePhi; }
1116 }
1117 }
1118 if ( safe < 0 ) { safe = 0; }
1119 return safe ;
1120}
1121
1122//////////////////////////////////////////////////////////////////////////////
1123//
1124// Calculate distance to surface of shape from `inside', allowing for tolerance
1125// - Only Calc rmax intersection if no valid rmin intersection
1126
1128 const G4ThreeVector& v,
1129 const G4bool calcNorm,
1130 G4bool* validNorm,
1131 G4ThreeVector* n ) const
1132{
1133 ESide side=kNull , sider=kNull, sidephi=kNull ;
1134 G4double snxt, srd=kInfinity, sphi=kInfinity, pdist ;
1135 G4double deltaR, t1, t2, t3, b, c, d2, roMin2 ;
1136
1137 // Vars for phi intersection:
1138
1139 G4double pDistS, compS, pDistE, compE, sphi2, xi, yi, vphi, roi2 ;
1140
1141 // Z plane intersection
1142
1143 if (v.z() > 0 )
1144 {
1145 pdist = fDz - p.z() ;
1146 if ( pdist > halfCarTolerance )
1147 {
1148 snxt = pdist/v.z() ;
1149 side = kPZ ;
1150 }
1151 else
1152 {
1153 if (calcNorm)
1154 {
1155 *n = G4ThreeVector(0,0,1) ;
1156 *validNorm = true ;
1157 }
1158 return snxt = 0 ;
1159 }
1160 }
1161 else if ( v.z() < 0 )
1162 {
1163 pdist = fDz + p.z() ;
1164
1165 if ( pdist > halfCarTolerance )
1166 {
1167 snxt = -pdist/v.z() ;
1168 side = kMZ ;
1169 }
1170 else
1171 {
1172 if (calcNorm)
1173 {
1174 *n = G4ThreeVector(0,0,-1) ;
1175 *validNorm = true ;
1176 }
1177 return snxt = 0.0 ;
1178 }
1179 }
1180 else
1181 {
1182 snxt = kInfinity ; // Travel perpendicular to z axis
1183 side = kNull;
1184 }
1185
1186 // Radial Intersections
1187 //
1188 // Find intersection with cylinders at rmax/rmin
1189 // Intersection point (xi,yi,zi) on line x=p.x+t*v.x etc.
1190 //
1191 // Intersects with x^2+y^2=R^2
1192 //
1193 // Hence (v.x^2+v.y^2)t^2+ 2t(p.x*v.x+p.y*v.y)+p.x^2+p.y^2-R^2=0
1194 //
1195 // t1 t2 t3
1196
1197 t1 = 1.0 - v.z()*v.z() ; // since v normalised
1198 t2 = p.x()*v.x() + p.y()*v.y() ;
1199 t3 = p.x()*p.x() + p.y()*p.y() ;
1200
1201 if ( snxt > 10*(fDz+fRMax) ) { roi2 = 2*fRMax*fRMax; }
1202 else { roi2 = snxt*snxt*t1 + 2*snxt*t2 + t3; } // radius^2 on +-fDz
1203
1204 if ( t1 > 0 ) // Check not parallel
1205 {
1206 // Calculate srd, r exit distance
1207
1208 if ( (t2 >= 0.0) && (roi2 > fRMax*(fRMax + kRadTolerance)) )
1209 {
1210 // Delta r not negative => leaving via rmax
1211
1212 deltaR = t3 - fRMax*fRMax ;
1213
1214 // NOTE: Should use rho-fRMax<-kRadTolerance*0.5
1215 // - avoid sqrt for efficiency
1216
1217 if ( deltaR < -kRadTolerance*fRMax )
1218 {
1219 b = t2/t1 ;
1220 c = deltaR/t1 ;
1221 d2 = b*b-c;
1222 if( d2 >= 0 ) { srd = c/( -b - std::sqrt(d2)); }
1223 else { srd = 0.; }
1224 sider = kRMax ;
1225 }
1226 else
1227 {
1228 // On tolerant boundary & heading outwards (or perpendicular to)
1229 // outer radial surface -> leaving immediately
1230
1231 if ( calcNorm )
1232 {
1234 *n = G4ThreeVector(p.x()*invRho,p.y()*invRho,0) ;
1235 *validNorm = true ;
1236 }
1237 return snxt = 0 ; // Leaving by rmax immediately
1238 }
1239 }
1240 else if ( t2 < 0. ) // i.e. t2 < 0; Possible rmin intersection
1241 {
1242 roMin2 = t3 - t2*t2/t1 ; // min ro2 of the plane of movement
1243
1244 if ( (fRMin != 0.0) && (roMin2 < fRMin*(fRMin - kRadTolerance)) )
1245 {
1246 deltaR = t3 - fRMin*fRMin ;
1247 b = t2/t1 ;
1248 c = deltaR/t1 ;
1249 d2 = b*b - c ;
1250
1251 if ( d2 >= 0 ) // Leaving via rmin
1252 {
1253 // NOTE: SHould use rho-rmin>kRadTolerance*0.5
1254 // - avoid sqrt for efficiency
1255
1256 if (deltaR > kRadTolerance*fRMin)
1257 {
1258 srd = c/(-b+std::sqrt(d2));
1259 sider = kRMin ;
1260 }
1261 else
1262 {
1263 if ( calcNorm ) {
1264 *validNorm = false;
1265 } // Concave side
1266 return snxt = 0.0;
1267 }
1268 }
1269 else // No rmin intersect -> must be rmax intersect
1270 {
1271 deltaR = t3 - fRMax*fRMax ;
1272 c = deltaR/t1 ;
1273 d2 = b*b-c;
1274 if( d2 >=0. )
1275 {
1276 srd = -b + std::sqrt(d2) ;
1277 sider = kRMax ;
1278 }
1279 else // Case: On the border+t2<kRadTolerance
1280 // (v is perpendicular to the surface)
1281 {
1282 if (calcNorm)
1283 {
1285 *n = G4ThreeVector(p.x()*invRho,p.y()*invRho,0) ;
1286 *validNorm = true ;
1287 }
1288 return snxt = 0.0;
1289 }
1290 }
1291 }
1292 else if ( roi2 > fRMax*(fRMax + kRadTolerance) )
1293 // No rmin intersect -> must be rmax intersect
1294 {
1295 deltaR = t3 - fRMax*fRMax ;
1296 b = t2/t1 ;
1297 c = deltaR/t1;
1298 d2 = b*b-c;
1299 if( d2 >= 0 )
1300 {
1301 srd = -b + std::sqrt(d2) ;
1302 sider = kRMax ;
1303 }
1304 else // Case: On the border+t2<kRadTolerance
1305 // (v is perpendicular to the surface)
1306 {
1307 if (calcNorm)
1308 {
1310 *n = G4ThreeVector(p.x()*invRho,p.y()*invRho,0) ;
1311 *validNorm = true ;
1312 }
1313 return snxt = 0.0;
1314 }
1315 }
1316 }
1317
1318 // Phi Intersection
1319
1320 if ( !fPhiFullTube )
1321 {
1322 // add angle calculation with correction
1323 // of the difference in domain of atan2 and Sphi
1324 //
1325 vphi = std::atan2(v.y(),v.x()) ;
1326
1327 if ( vphi < fSPhi - halfAngTolerance ) { vphi += twopi; }
1328 else if ( vphi > fSPhi + fDPhi + halfAngTolerance ) { vphi -= twopi; }
1329
1330
1331 if ( (p.x() != 0.0) || (p.y() != 0.0) ) // Check if on z axis (rho not needed later)
1332 {
1333 // pDist -ve when inside
1334
1335 pDistS = p.x()*sinSPhi - p.y()*cosSPhi ;
1336 pDistE = -p.x()*sinEPhi + p.y()*cosEPhi ;
1337
1338 // Comp -ve when in direction of outwards normal
1339
1340 compS = -sinSPhi*v.x() + cosSPhi*v.y() ;
1341 compE = sinEPhi*v.x() - cosEPhi*v.y() ;
1342
1343 sidephi = kNull;
1344
1345 if( ( (fDPhi <= pi) && ( (pDistS <= halfCarTolerance)
1346 && (pDistE <= halfCarTolerance) ) )
1347 || ( (fDPhi > pi) && ((pDistS <= halfCarTolerance)
1348 || (pDistE <= halfCarTolerance) ) ) )
1349 {
1350 // Inside both phi *full* planes
1351
1352 if ( compS < 0 )
1353 {
1354 sphi = pDistS/compS ;
1355
1356 if (sphi >= -halfCarTolerance)
1357 {
1358 xi = p.x() + sphi*v.x() ;
1359 yi = p.y() + sphi*v.y() ;
1360
1361 // Check intersecting with correct half-plane
1362 // (if not -> no intersect)
1363 //
1364 if((std::fabs(xi)<=kCarTolerance)&&(std::fabs(yi)<=kCarTolerance))
1365 {
1366 sidephi = kSPhi;
1367 if (((fSPhi-halfAngTolerance)<=vphi)
1368 &&((fSPhi+fDPhi+halfAngTolerance)>=vphi))
1369 {
1370 sphi = kInfinity;
1371 }
1372 }
1373 else if ( yi*cosCPhi-xi*sinCPhi >=0 )
1374 {
1375 sphi = kInfinity ;
1376 }
1377 else
1378 {
1379 sidephi = kSPhi ;
1380 if ( pDistS > -halfCarTolerance )
1381 {
1382 sphi = 0.0 ; // Leave by sphi immediately
1383 }
1384 }
1385 }
1386 else
1387 {
1388 sphi = kInfinity ;
1389 }
1390 }
1391 else
1392 {
1393 sphi = kInfinity ;
1394 }
1395
1396 if ( compE < 0 )
1397 {
1398 sphi2 = pDistE/compE ;
1399
1400 // Only check further if < starting phi intersection
1401 //
1402 if ( (sphi2 > -halfCarTolerance) && (sphi2 < sphi) )
1403 {
1404 xi = p.x() + sphi2*v.x() ;
1405 yi = p.y() + sphi2*v.y() ;
1406
1407 if((std::fabs(xi)<=kCarTolerance)&&(std::fabs(yi)<=kCarTolerance))
1408 {
1409 // Leaving via ending phi
1410 //
1411 if( (fSPhi-halfAngTolerance > vphi)
1412 ||(fSPhi+fDPhi+halfAngTolerance < vphi) )
1413 {
1414 sidephi = kEPhi ;
1415 if ( pDistE <= -halfCarTolerance ) { sphi = sphi2 ; }
1416 else { sphi = 0.0 ; }
1417 }
1418 }
1419 else // Check intersecting with correct half-plane
1420
1421 if ( (yi*cosCPhi-xi*sinCPhi) >= 0)
1422 {
1423 // Leaving via ending phi
1424 //
1425 sidephi = kEPhi ;
1426 if ( pDistE <= -halfCarTolerance ) { sphi = sphi2 ; }
1427 else { sphi = 0.0 ; }
1428 }
1429 }
1430 }
1431 }
1432 else
1433 {
1434 sphi = kInfinity ;
1435 }
1436 }
1437 else
1438 {
1439 // On z axis + travel not || to z axis -> if phi of vector direction
1440 // within phi of shape, Step limited by rmax, else Step =0
1441
1442 if ( (fSPhi - halfAngTolerance <= vphi)
1443 && (vphi <= fSPhi + fDPhi + halfAngTolerance ) )
1444 {
1445 sphi = kInfinity ;
1446 }
1447 else
1448 {
1449 sidephi = kSPhi ; // arbitrary
1450 sphi = 0.0 ;
1451 }
1452 }
1453 if (sphi < snxt) // Order intersecttions
1454 {
1455 snxt = sphi ;
1456 side = sidephi ;
1457 }
1458 }
1459 if (srd < snxt) // Order intersections
1460 {
1461 snxt = srd ;
1462 side = sider ;
1463 }
1464 }
1465 if (calcNorm)
1466 {
1467 switch(side)
1468 {
1469 case kRMax:
1470 // Note: returned vector not normalised
1471 // (divide by fRMax for unit vector)
1472 //
1473 xi = p.x() + snxt*v.x() ;
1474 yi = p.y() + snxt*v.y() ;
1475 *n = G4ThreeVector(xi/fRMax,yi/fRMax,0) ;
1476 *validNorm = true ;
1477 break ;
1478
1479 case kRMin:
1480 *validNorm = false ; // Rmin is inconvex
1481 break ;
1482
1483 case kSPhi:
1484 if ( fDPhi <= pi )
1485 {
1486 *n = G4ThreeVector(sinSPhi,-cosSPhi,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(-sinEPhi,cosEPhi,0) ;
1499 *validNorm = true ;
1500 }
1501 else
1502 {
1503 *validNorm = false ;
1504 }
1505 break ;
1506
1507 case kPZ:
1508 *n = G4ThreeVector(0,0,1) ;
1509 *validNorm = true ;
1510 break ;
1511
1512 case kMZ:
1513 *n = G4ThreeVector(0,0,-1) ;
1514 *validNorm = true ;
1515 break ;
1516
1517 default:
1518 G4cout << G4endl ;
1519 DumpInfo();
1520 std::ostringstream message;
1521 G4long oldprc = message.precision(16);
1522 message << "Undefined side for valid surface normal to solid."
1523 << G4endl
1524 << "Position:" << G4endl << G4endl
1525 << "p.x() = " << p.x()/mm << " mm" << G4endl
1526 << "p.y() = " << p.y()/mm << " mm" << G4endl
1527 << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl
1528 << "Direction:" << G4endl << G4endl
1529 << "v.x() = " << v.x() << G4endl
1530 << "v.y() = " << v.y() << G4endl
1531 << "v.z() = " << v.z() << G4endl << G4endl
1532 << "Proposed distance :" << G4endl << G4endl
1533 << "snxt = " << snxt/mm << " mm" << G4endl ;
1534 message.precision(oldprc) ;
1535 G4Exception("G4Tubs::DistanceToOut(p,v,..)", "GeomSolids1002",
1536 JustWarning, message);
1537 break ;
1538 }
1539 }
1540 if ( snxt<halfCarTolerance ) { snxt=0 ; }
1541
1542 return snxt ;
1543}
1544
1545//////////////////////////////////////////////////////////////////////////
1546//
1547// Calculate distance (<=actual) to closest surface of shape from inside
1548
1550{
1551 G4double safe=0.0, rho, safeR1, safeR2, safeZ, safePhi ;
1552 rho = std::sqrt(p.x()*p.x() + p.y()*p.y()) ;
1553
1554#ifdef G4CSGDEBUG
1555 if( Inside(p) == kOutside )
1556 {
1557 G4long oldprc = G4cout.precision(16) ;
1558 G4cout << G4endl ;
1559 DumpInfo();
1560 G4cout << "Position:" << G4endl << G4endl ;
1561 G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl ;
1562 G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl ;
1563 G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl ;
1564 G4cout.precision(oldprc) ;
1565 G4Exception("G4Tubs::DistanceToOut(p)", "GeomSolids1002",
1566 JustWarning, "Point p is outside !?");
1567 }
1568#endif
1569
1570 if ( fRMin != 0.0 )
1571 {
1572 safeR1 = rho - fRMin ;
1573 safeR2 = fRMax - rho ;
1574
1575 if ( safeR1 < safeR2 ) { safe = safeR1 ; }
1576 else { safe = safeR2 ; }
1577 }
1578 else
1579 {
1580 safe = fRMax - rho ;
1581 }
1582 safeZ = fDz - std::fabs(p.z()) ;
1583
1584 if ( safeZ < safe ) { safe = safeZ ; }
1585
1586 // Check if phi divided, Calc distances closest phi plane
1587 //
1588 if ( !fPhiFullTube )
1589 {
1590 if ( p.y()*cosCPhi-p.x()*sinCPhi <= 0 )
1591 {
1592 safePhi = -(p.x()*sinSPhi - p.y()*cosSPhi) ;
1593 }
1594 else
1595 {
1596 safePhi = (p.x()*sinEPhi - p.y()*cosEPhi) ;
1597 }
1598 if (safePhi < safe) { safe = safePhi ; }
1599 }
1600 if ( safe < 0 ) { safe = 0 ; }
1601
1602 return safe ;
1603}
1604
1605//////////////////////////////////////////////////////////////////////////
1606//
1607// Stream object contents to an output stream
1608
1610{
1611 return {"G4Tubs"};
1612}
1613
1614//////////////////////////////////////////////////////////////////////////
1615//
1616// Make a clone of the object
1617//
1619{
1620 return new G4Tubs(*this);
1621}
1622
1623//////////////////////////////////////////////////////////////////////////
1624//
1625// Stream object contents to an output stream
1626
1627std::ostream& G4Tubs::StreamInfo( std::ostream& os ) const
1628{
1629 G4long oldprc = os.precision(16);
1630 os << "-----------------------------------------------------------\n"
1631 << " *** Dump for solid - " << GetName() << " ***\n"
1632 << " ===================================================\n"
1633 << " Solid type: G4Tubs\n"
1634 << " Parameters: \n"
1635 << " inner radius : " << fRMin/mm << " mm \n"
1636 << " outer radius : " << fRMax/mm << " mm \n"
1637 << " half length Z: " << fDz/mm << " mm \n"
1638 << " starting phi : " << fSPhi/degree << " degrees \n"
1639 << " delta phi : " << fDPhi/degree << " degrees \n"
1640 << "-----------------------------------------------------------\n";
1641 os.precision(oldprc);
1642
1643 return os;
1644}
1645
1646/////////////////////////////////////////////////////////////////////////
1647//
1648// GetPointOnSurface
1649
1651{
1652 G4double Rmax = fRMax;
1653 G4double Rmin = fRMin;
1654 G4double hz = 2.*fDz; // height
1655 G4double lext = fDPhi*Rmax; // length of external circular arc
1656 G4double lint = fDPhi*Rmin; // length of internal circular arc
1657
1658 // Set array of surface areas
1659 //
1660 G4double RRmax = Rmax * Rmax;
1661 G4double RRmin = Rmin * Rmin;
1662 G4double sbase = 0.5*fDPhi*(RRmax - RRmin);
1663 G4double scut = (fDPhi == twopi) ? 0. : hz*(Rmax - Rmin);
1664 G4double ssurf[6] = { scut, scut, sbase, sbase, hz*lext, hz*lint };
1665 ssurf[1] += ssurf[0];
1666 ssurf[2] += ssurf[1];
1667 ssurf[3] += ssurf[2];
1668 ssurf[4] += ssurf[3];
1669 ssurf[5] += ssurf[4];
1670
1671 // Select surface
1672 //
1673 G4double select = ssurf[5]*G4QuickRand();
1674 G4int k = 5;
1675 k -= (G4int)(select <= ssurf[4]);
1676 k -= (G4int)(select <= ssurf[3]);
1677 k -= (G4int)(select <= ssurf[2]);
1678 k -= (G4int)(select <= ssurf[1]);
1679 k -= (G4int)(select <= ssurf[0]);
1680
1681 // Generate point on selected surface
1682 //
1683 switch(k)
1684 {
1685 case 0: // start phi cut
1686 {
1687 G4double r = Rmin + (Rmax - Rmin)*G4QuickRand();
1688 return { r*cosSPhi, r*sinSPhi, hz*G4QuickRand() - fDz };
1689 }
1690 case 1: // end phi cut
1691 {
1692 G4double r = Rmin + (Rmax - Rmin)*G4QuickRand();
1693 return { r*cosEPhi, r*sinEPhi, hz*G4QuickRand() - fDz };
1694 }
1695 case 2: // base at -dz
1696 {
1697 G4double r = std::sqrt(RRmin + (RRmax - RRmin)*G4QuickRand());
1698 G4double phi = fSPhi + fDPhi*G4QuickRand();
1699 return { r*std::cos(phi), r*std::sin(phi), -fDz };
1700 }
1701 case 3: // base at +dz
1702 {
1703 G4double r = std::sqrt(RRmin + (RRmax - RRmin)*G4QuickRand());
1704 G4double phi = fSPhi + fDPhi*G4QuickRand();
1705 return { r*std::cos(phi), r*std::sin(phi), fDz };
1706 }
1707 case 4: // external lateral surface
1708 {
1709 G4double phi = fSPhi + fDPhi*G4QuickRand();
1710 G4double z = hz*G4QuickRand() - fDz;
1711 G4double x = Rmax*std::cos(phi);
1712 G4double y = Rmax*std::sin(phi);
1713 return { x,y,z };
1714 }
1715 case 5: // internal lateral surface
1716 {
1717 G4double phi = fSPhi + fDPhi*G4QuickRand();
1718 G4double z = hz*G4QuickRand() - fDz;
1719 G4double x = Rmin*std::cos(phi);
1720 G4double y = Rmin*std::sin(phi);
1721 return { x,y,z };
1722 }
1723 }
1724 return {0., 0., 0.};
1725}
1726
1727/////////////////////////////////////////////////////////////////////////
1728//
1729// GetCubicVolume
1730
1732{
1733 if (fCubicVolume == 0)
1734 {
1735 G4AutoLock l(&tubsMutex);
1737 l.unlock();
1738 }
1739 return fCubicVolume;
1740}
1741
1742/////////////////////////////////////////////////////////////////////////
1743//
1744// GetSurfaceArea
1745
1747{
1748 if (fSurfaceArea == 0)
1749 {
1750 G4AutoLock l(&tubsMutex);
1752 if (!fPhiFullTube)
1753 {
1755 }
1756 l.unlock();
1757 }
1758 return fSurfaceArea;
1759}
1760
1761///////////////////////////////////////////////////////////////////////////
1762//
1763// Methods for visualisation
1764
1766{
1767 scene.AddSolid (*this) ;
1768}
1769
1771{
1772 return new G4PolyhedronTubs (fRMin, fRMax, fDz, fSPhi, fDPhi) ;
1773}
1774
1775#endif
G4TemplateAutoLock< G4Mutex > G4AutoLock
std::vector< G4ThreeVector > G4ThreeVectorList
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
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
double z() const
Hep3Vector unit() const
double x() const
double y() const
void set(double x, double y, double z)
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
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
void CheckPhiAngles(G4double sPhi, G4double dPhi)
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition G4Tubs.cc:585
G4Tubs(const G4String &pName, G4double pRMin, G4double pRMax, G4double pDz, G4double pSPhi, G4double pDPhi)
Definition G4Tubs.cc:73
G4double GetZHalfLength() const
G4double cosHDPhiIT
Definition G4Tubs.hh:271
G4double sinCPhi
Definition G4Tubs.hh:271
G4double cosEPhi
Definition G4Tubs.hh:272
G4ThreeVector GetPointOnSurface() const override
Definition G4Tubs.cc:1650
G4double GetSurfaceArea() override
Definition G4Tubs.cc:1746
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
Definition G4Tubs.cc:168
G4double fRMin
Definition G4Tubs.hh:268
G4double kAngTolerance
Definition G4Tubs.hh:262
G4double fDPhi
Definition G4Tubs.hh:268
G4double GetCosStartPhi() const
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
Definition G4Tubs.cc:1765
G4double GetCosEndPhi() const
G4double fRMax
Definition G4Tubs.hh:268
G4double GetInnerRadius() const
G4double fInvRmin
Definition G4Tubs.hh:278
G4double GetCubicVolume() override
Definition G4Tubs.cc:1731
G4Tubs & operator=(const G4Tubs &rhs)
Definition G4Tubs.cc:123
G4double GetOuterRadius() const
G4double sinSPhi
Definition G4Tubs.hh:272
G4double fInvRmax
Definition G4Tubs.hh:278
static constexpr G4double kNormTolerance
Definition G4Tubs.hh:265
G4double cosCPhi
Definition G4Tubs.hh:271
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
Definition G4Tubs.cc:494
EInside Inside(const G4ThreeVector &p) const override
Definition G4Tubs.cc:327
G4double cosHDPhi
Definition G4Tubs.hh:271
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const override
Definition G4Tubs.cc:211
G4Polyhedron * CreatePolyhedron() const override
Definition G4Tubs.cc:1770
G4double cosSPhi
Definition G4Tubs.hh:272
G4double GetSinEndPhi() const
G4double sinEPhi
Definition G4Tubs.hh:272
std::ostream & StreamInfo(std::ostream &os) const override
Definition G4Tubs.cc:1627
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
Definition G4Tubs.cc:1127
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
Definition G4Tubs.cc:157
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
Definition G4Tubs.cc:717
G4double kRadTolerance
Definition G4Tubs.hh:262
G4double FastInverseRxy(const G4ThreeVector &pos, G4double invRad, G4double normalTolerance) const
G4double fDz
Definition G4Tubs.hh:268
G4bool fPhiFullTube
Definition G4Tubs.hh:275
G4GeometryType GetEntityType() const override
Definition G4Tubs.cc:1609
G4double halfRadTolerance
Definition G4Tubs.hh:281
G4double halfAngTolerance
Definition G4Tubs.hh:281
G4double halfCarTolerance
Definition G4Tubs.hh:281
G4double GetDeltaPhiAngle() const
G4double GetSinStartPhi() const
G4VSolid * Clone() const override
Definition G4Tubs.cc:1618
G4double cosHDPhiOT
Definition G4Tubs.hh:271
G4double fSPhi
Definition G4Tubs.hh:268
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...
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