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