Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Tubs Class Reference

#include <G4Tubs.hh>

+ Inheritance diagram for G4Tubs:

Public Member Functions

 G4Tubs (const G4String &pName, G4double pRMin, G4double pRMax, G4double pDz, G4double pSPhi, G4double pDPhi)
 
virtual ~G4Tubs ()
 
G4double GetInnerRadius () const
 
G4double GetOuterRadius () const
 
G4double GetZHalfLength () const
 
G4double GetStartPhiAngle () const
 
G4double GetDeltaPhiAngle () const
 
void SetInnerRadius (G4double newRMin)
 
void SetOuterRadius (G4double newRMax)
 
void SetZHalfLength (G4double newDz)
 
void SetStartPhiAngle (G4double newSPhi, G4bool trig=true)
 
void SetDeltaPhiAngle (G4double newDPhi)
 
G4double GetCubicVolume ()
 
G4double GetSurfaceArea ()
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
 
EInside Inside (const G4ThreeVector &p) const
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
G4double DistanceToIn (const G4ThreeVector &p) const
 
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
 
G4double DistanceToOut (const G4ThreeVector &p) const
 
G4GeometryType GetEntityType () const
 
G4ThreeVector GetPointOnSurface () const
 
G4VSolidClone () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
G4PolyhedronCreatePolyhedron () const
 
G4NURBSCreateNURBS () const
 
 G4Tubs (__void__ &)
 
 G4Tubs (const G4Tubs &rhs)
 
G4Tubsoperator= (const G4Tubs &rhs)
 
G4double GetRMin () const
 
G4double GetRMax () const
 
G4double GetDz () const
 
G4double GetSPhi () const
 
G4double GetDPhi () const
 
- Public Member Functions inherited from G4CSGSolid
 G4CSGSolid (const G4String &pName)
 
virtual ~G4CSGSolid ()
 
virtual std::ostream & StreamInfo (std::ostream &os) const
 
virtual G4PolyhedronGetPolyhedron () const
 
 G4CSGSolid (__void__ &)
 
 G4CSGSolid (const G4CSGSolid &rhs)
 
G4CSGSolidoperator= (const G4CSGSolid &rhs)
 
- Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
 
virtual ~G4VSolid ()
 
G4bool operator== (const G4VSolid &s) const
 
G4String GetName () const
 
void SetName (const G4String &name)
 
G4double GetTolerance () const
 
virtual G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
 
virtual EInside Inside (const G4ThreeVector &p) const =0
 
virtual G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p) const =0
 
virtual void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
virtual G4double GetCubicVolume ()
 
virtual G4double GetSurfaceArea ()
 
virtual G4GeometryType GetEntityType () const =0
 
virtual G4ThreeVector GetPointOnSurface () const
 
virtual G4VSolidClone () const
 
virtual std::ostream & StreamInfo (std::ostream &os) const =0
 
void DumpInfo () const
 
virtual void DescribeYourselfTo (G4VGraphicsScene &scene) const =0
 
virtual G4VisExtent GetExtent () const
 
virtual G4PolyhedronCreatePolyhedron () const
 
virtual G4NURBSCreateNURBS () const
 
virtual G4PolyhedronGetPolyhedron () const
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
 G4VSolid (__void__ &)
 
 G4VSolid (const G4VSolid &rhs)
 
G4VSolidoperator= (const G4VSolid &rhs)
 

Protected Types

enum  ESide {
  kNull , kRMin , kRMax , kSPhi ,
  kEPhi , kPZ , kMZ
}
 
enum  ENorm {
  kNRMin , kNRMax , kNSPhi , kNEPhi ,
  kNZ
}
 

Protected Member Functions

G4ThreeVectorListCreateRotatedVertices (const G4AffineTransform &pTransform) const
 
void Initialize ()
 
void CheckSPhiAngle (G4double sPhi)
 
void CheckDPhiAngle (G4double dPhi)
 
void CheckPhiAngles (G4double sPhi, G4double dPhi)
 
void InitializeTrigonometry ()
 
virtual G4ThreeVector ApproxSurfaceNormal (const G4ThreeVector &p) const
 
- Protected Member Functions inherited from G4CSGSolid
G4double GetRadiusInRing (G4double rmin, G4double rmax) const
 
- Protected Member Functions inherited from G4VSolid
void CalculateClippedPolygonExtent (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipCrossSection (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipBetweenSections (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipPolygon (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 

Protected Attributes

G4double kRadTolerance
 
G4double kAngTolerance
 
G4double fRMin
 
G4double fRMax
 
G4double fDz
 
G4double fSPhi
 
G4double fDPhi
 
G4double sinCPhi
 
G4double cosCPhi
 
G4double cosHDPhiOT
 
G4double cosHDPhiIT
 
G4double sinSPhi
 
G4double cosSPhi
 
G4double sinEPhi
 
G4double cosEPhi
 
G4bool fPhiFullTube
 
- Protected Attributes inherited from G4CSGSolid
G4double fCubicVolume
 
G4double fSurfaceArea
 
G4PolyhedronfpPolyhedron
 
- Protected Attributes inherited from G4VSolid
G4double kCarTolerance
 

Detailed Description

Definition at line 76 of file G4Tubs.hh.

Member Enumeration Documentation

◆ ENorm

enum G4Tubs::ENorm
protected
Enumerator
kNRMin 
kNRMax 
kNSPhi 
kNEPhi 
kNZ 

Definition at line 203 of file G4Tubs.hh.

@ kNRMax
Definition: G4Tubs.hh:203
@ kNZ
Definition: G4Tubs.hh:203
@ kNRMin
Definition: G4Tubs.hh:203
@ kNEPhi
Definition: G4Tubs.hh:203
@ kNSPhi
Definition: G4Tubs.hh:203

◆ ESide

enum G4Tubs::ESide
protected
Enumerator
kNull 
kRMin 
kRMax 
kSPhi 
kEPhi 
kPZ 
kMZ 

Definition at line 199 of file G4Tubs.hh.

@ kEPhi
Definition: G4Tubs.hh:199
@ kRMax
Definition: G4Tubs.hh:199
@ kPZ
Definition: G4Tubs.hh:199
@ kMZ
Definition: G4Tubs.hh:199
@ kRMin
Definition: G4Tubs.hh:199
@ kSPhi
Definition: G4Tubs.hh:199
@ kNull
Definition: G4Tubs.hh:199

Constructor & Destructor Documentation

◆ G4Tubs() [1/3]

G4Tubs::G4Tubs ( const G4String pName,
G4double  pRMin,
G4double  pRMax,
G4double  pDz,
G4double  pSPhi,
G4double  pDPhi 
)

Definition at line 88 of file G4Tubs.cc.

92 : G4CSGSolid(pName), fRMin(pRMin), fRMax(pRMax), fDz(pDz), fSPhi(0), fDPhi(0)
93{
94
97
98 if (pDz<=0) // Check z-len
99 {
100 std::ostringstream message;
101 message << "Negative Z half-length (" << pDz << ") in solid: " << GetName();
102 G4Exception("G4Tubs::G4Tubs()", "GeomSolids0002", FatalException, message);
103 }
104 if ( (pRMin >= pRMax) || (pRMin < 0) ) // Check radii
105 {
106 std::ostringstream message;
107 message << "Invalid values for radii in solid: " << GetName()
108 << G4endl
109 << " pRMin = " << pRMin << ", pRMax = " << pRMax;
110 G4Exception("G4Tubs::G4Tubs()", "GeomSolids0002", FatalException, message);
111 }
112
113 // Check angles
114 //
115 CheckPhiAngles(pSPhi, pDPhi);
116}
@ FatalException
#define G4endl
Definition: G4ios.hh:52
G4double GetRadialTolerance() const
static G4GeometryTolerance * GetInstance()
G4double GetAngularTolerance() const
void CheckPhiAngles(G4double sPhi, G4double dPhi)
G4double fRMin
Definition: G4Tubs.hh:209
G4double kAngTolerance
Definition: G4Tubs.hh:205
G4double fDPhi
Definition: G4Tubs.hh:209
G4double fRMax
Definition: G4Tubs.hh:209
G4double kRadTolerance
Definition: G4Tubs.hh:205
G4double fDz
Definition: G4Tubs.hh:209
G4double fSPhi
Definition: G4Tubs.hh:209
G4String GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ ~G4Tubs()

G4Tubs::~G4Tubs ( )
virtual

Definition at line 136 of file G4Tubs.cc.

137{
138}

◆ G4Tubs() [2/3]

G4Tubs::G4Tubs ( __void__ &  a)

Definition at line 123 of file G4Tubs.cc.

125 fRMin(0.), fRMax(0.), fDz(0.), fSPhi(0.), fDPhi(0.),
126 sinCPhi(0.), cosCPhi(0.), cosHDPhiOT(0.), cosHDPhiIT(0.),
127 sinSPhi(0.), cosSPhi(0.), sinEPhi(0.), cosEPhi(0.),
128 fPhiFullTube(false)
129{
130}
G4double cosHDPhiIT
Definition: G4Tubs.hh:213
G4double sinCPhi
Definition: G4Tubs.hh:213
G4double cosEPhi
Definition: G4Tubs.hh:214
G4double sinSPhi
Definition: G4Tubs.hh:214
G4double cosCPhi
Definition: G4Tubs.hh:213
G4double cosSPhi
Definition: G4Tubs.hh:214
G4double sinEPhi
Definition: G4Tubs.hh:214
G4bool fPhiFullTube
Definition: G4Tubs.hh:218
G4double cosHDPhiOT
Definition: G4Tubs.hh:213

◆ G4Tubs() [3/3]

G4Tubs::G4Tubs ( const G4Tubs rhs)

Definition at line 144 of file G4Tubs.cc.

145 : G4CSGSolid(rhs),
147 fRMin(rhs.fRMin), fRMax(rhs.fRMax), fDz(rhs.fDz),
148 fSPhi(rhs.fSPhi), fDPhi(rhs.fDPhi),
149 sinCPhi(rhs.sinCPhi), cosCPhi(rhs.sinCPhi),
151 sinSPhi(rhs.sinSPhi), cosSPhi(rhs.cosSPhi),
153{
154}

Member Function Documentation

◆ ApproxSurfaceNormal()

G4ThreeVector G4Tubs::ApproxSurfaceNormal ( const G4ThreeVector p) const
protectedvirtual

Reimplemented in G4CutTubs.

Definition at line 678 of file G4Tubs.cc.

679{
680 ENorm side ;
681 G4ThreeVector norm ;
682 G4double rho, phi ;
683 G4double distZ, distRMin, distRMax, distSPhi, distEPhi, distMin ;
684
685 rho = std::sqrt(p.x()*p.x() + p.y()*p.y()) ;
686
687 distRMin = std::fabs(rho - fRMin) ;
688 distRMax = std::fabs(rho - fRMax) ;
689 distZ = std::fabs(std::fabs(p.z()) - fDz) ;
690
691 if (distRMin < distRMax) // First minimum
692 {
693 if ( distZ < distRMin )
694 {
695 distMin = distZ ;
696 side = kNZ ;
697 }
698 else
699 {
700 distMin = distRMin ;
701 side = kNRMin ;
702 }
703 }
704 else
705 {
706 if ( distZ < distRMax )
707 {
708 distMin = distZ ;
709 side = kNZ ;
710 }
711 else
712 {
713 distMin = distRMax ;
714 side = kNRMax ;
715 }
716 }
717 if (!fPhiFullTube && rho ) // Protected against (0,0,z)
718 {
719 phi = std::atan2(p.y(),p.x()) ;
720
721 if ( phi < 0 ) { phi += twopi; }
722
723 if ( fSPhi < 0 )
724 {
725 distSPhi = std::fabs(phi - (fSPhi + twopi))*rho ;
726 }
727 else
728 {
729 distSPhi = std::fabs(phi - fSPhi)*rho ;
730 }
731 distEPhi = std::fabs(phi - fSPhi - fDPhi)*rho ;
732
733 if (distSPhi < distEPhi) // Find new minimum
734 {
735 if ( distSPhi < distMin )
736 {
737 side = kNSPhi ;
738 }
739 }
740 else
741 {
742 if ( distEPhi < distMin )
743 {
744 side = kNEPhi ;
745 }
746 }
747 }
748 switch ( side )
749 {
750 case kNRMin : // Inner radius
751 {
752 norm = G4ThreeVector(-p.x()/rho, -p.y()/rho, 0) ;
753 break ;
754 }
755 case kNRMax : // Outer radius
756 {
757 norm = G4ThreeVector(p.x()/rho, p.y()/rho, 0) ;
758 break ;
759 }
760 case kNZ : // + or - dz
761 {
762 if ( p.z() > 0 ) { norm = G4ThreeVector(0,0,1) ; }
763 else { norm = G4ThreeVector(0,0,-1); }
764 break ;
765 }
766 case kNSPhi:
767 {
768 norm = G4ThreeVector(std::sin(fSPhi), -std::cos(fSPhi), 0) ;
769 break ;
770 }
771 case kNEPhi:
772 {
773 norm = G4ThreeVector(-std::sin(fSPhi+fDPhi), std::cos(fSPhi+fDPhi), 0) ;
774 break;
775 }
776 default: // Should never reach this case ...
777 {
778 DumpInfo();
779 G4Exception("G4Tubs::ApproxSurfaceNormal()",
780 "GeomSolids1002", JustWarning,
781 "Undefined side for valid surface normal to solid.");
782 break ;
783 }
784 }
785 return norm;
786}
ENorm
Definition: G4Cons.cc:72
@ JustWarning
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
double z() const
double x() const
double y() const
void DumpInfo() const

Referenced by SurfaceNormal().

◆ CalculateExtent()

G4bool G4Tubs::CalculateExtent ( const EAxis  pAxis,
const G4VoxelLimits pVoxelLimit,
const G4AffineTransform pTransform,
G4double pmin,
G4double pmax 
) const
virtual

Implements G4VSolid.

Definition at line 200 of file G4Tubs.cc.

205{
206
207 if ( (!pTransform.IsRotated()) && (fDPhi == twopi) && (fRMin == 0) )
208 {
209 // Special case handling for unrotated solid tubes
210 // Compute x/y/z mins and maxs fro bounding box respecting limits,
211 // with early returns if outside limits. Then switch() on pAxis,
212 // and compute exact x and y limit for x/y case
213
214 G4double xoffset, xMin, xMax;
215 G4double yoffset, yMin, yMax;
216 G4double zoffset, zMin, zMax;
217
218 G4double diff1, diff2, maxDiff, newMin, newMax;
219 G4double xoff1, xoff2, yoff1, yoff2, delta;
220
221 xoffset = pTransform.NetTranslation().x();
222 xMin = xoffset - fRMax;
223 xMax = xoffset + fRMax;
224
225 if (pVoxelLimit.IsXLimited())
226 {
227 if ( (xMin > pVoxelLimit.GetMaxXExtent())
228 || (xMax < pVoxelLimit.GetMinXExtent()) )
229 {
230 return false;
231 }
232 else
233 {
234 if (xMin < pVoxelLimit.GetMinXExtent())
235 {
236 xMin = pVoxelLimit.GetMinXExtent();
237 }
238 if (xMax > pVoxelLimit.GetMaxXExtent())
239 {
240 xMax = pVoxelLimit.GetMaxXExtent();
241 }
242 }
243 }
244 yoffset = pTransform.NetTranslation().y();
245 yMin = yoffset - fRMax;
246 yMax = yoffset + fRMax;
247
248 if ( pVoxelLimit.IsYLimited() )
249 {
250 if ( (yMin > pVoxelLimit.GetMaxYExtent())
251 || (yMax < pVoxelLimit.GetMinYExtent()) )
252 {
253 return false;
254 }
255 else
256 {
257 if (yMin < pVoxelLimit.GetMinYExtent())
258 {
259 yMin = pVoxelLimit.GetMinYExtent();
260 }
261 if (yMax > pVoxelLimit.GetMaxYExtent())
262 {
263 yMax=pVoxelLimit.GetMaxYExtent();
264 }
265 }
266 }
267 zoffset = pTransform.NetTranslation().z();
268 zMin = zoffset - fDz;
269 zMax = zoffset + fDz;
270
271 if ( pVoxelLimit.IsZLimited() )
272 {
273 if ( (zMin > pVoxelLimit.GetMaxZExtent())
274 || (zMax < pVoxelLimit.GetMinZExtent()) )
275 {
276 return false;
277 }
278 else
279 {
280 if (zMin < pVoxelLimit.GetMinZExtent())
281 {
282 zMin = pVoxelLimit.GetMinZExtent();
283 }
284 if (zMax > pVoxelLimit.GetMaxZExtent())
285 {
286 zMax = pVoxelLimit.GetMaxZExtent();
287 }
288 }
289 }
290 switch ( pAxis ) // Known to cut cylinder
291 {
292 case kXAxis :
293 {
294 yoff1 = yoffset - yMin;
295 yoff2 = yMax - yoffset;
296
297 if ( (yoff1 >= 0) && (yoff2 >= 0) ) // Y limits cross max/min x
298 { // => no change
299 pMin = xMin;
300 pMax = xMax;
301 }
302 else
303 {
304 // Y limits don't cross max/min x => compute max delta x,
305 // hence new mins/maxs
306
307 delta = fRMax*fRMax - yoff1*yoff1;
308 diff1 = (delta>0.) ? std::sqrt(delta) : 0.;
309 delta = fRMax*fRMax - yoff2*yoff2;
310 diff2 = (delta>0.) ? std::sqrt(delta) : 0.;
311 maxDiff = (diff1 > diff2) ? diff1:diff2;
312 newMin = xoffset - maxDiff;
313 newMax = xoffset + maxDiff;
314 pMin = (newMin < xMin) ? xMin : newMin;
315 pMax = (newMax > xMax) ? xMax : newMax;
316 }
317 break;
318 }
319 case kYAxis :
320 {
321 xoff1 = xoffset - xMin;
322 xoff2 = xMax - xoffset;
323
324 if ( (xoff1 >= 0) && (xoff2 >= 0) ) // X limits cross max/min y
325 { // => no change
326 pMin = yMin;
327 pMax = yMax;
328 }
329 else
330 {
331 // X limits don't cross max/min y => compute max delta y,
332 // hence new mins/maxs
333
334 delta = fRMax*fRMax - xoff1*xoff1;
335 diff1 = (delta>0.) ? std::sqrt(delta) : 0.;
336 delta = fRMax*fRMax - xoff2*xoff2;
337 diff2 = (delta>0.) ? std::sqrt(delta) : 0.;
338 maxDiff = (diff1 > diff2) ? diff1 : diff2;
339 newMin = yoffset - maxDiff;
340 newMax = yoffset + maxDiff;
341 pMin = (newMin < yMin) ? yMin : newMin;
342 pMax = (newMax > yMax) ? yMax : newMax;
343 }
344 break;
345 }
346 case kZAxis:
347 {
348 pMin = zMin;
349 pMax = zMax;
350 break;
351 }
352 default:
353 break;
354 }
355 pMin -= kCarTolerance;
356 pMax += kCarTolerance;
357 return true;
358 }
359 else // Calculate rotated vertex coordinates
360 {
361 G4int i, noEntries, noBetweenSections4;
362 G4bool existsAfterClip = false;
363 G4ThreeVectorList* vertices = CreateRotatedVertices(pTransform);
364
365 pMin = kInfinity;
366 pMax = -kInfinity;
367
368 noEntries = vertices->size();
369 noBetweenSections4 = noEntries - 4;
370
371 for ( i = 0 ; i < noEntries ; i += 4 )
372 {
373 ClipCrossSection(vertices, i, pVoxelLimit, pAxis, pMin, pMax);
374 }
375 for ( i = 0 ; i < noBetweenSections4 ; i += 4 )
376 {
377 ClipBetweenSections(vertices, i, pVoxelLimit, pAxis, pMin, pMax);
378 }
379 if ( (pMin != kInfinity) || (pMax != -kInfinity) )
380 {
381 existsAfterClip = true;
382 pMin -= kCarTolerance; // Add 2*tolerance to avoid precision troubles
383 pMax += kCarTolerance;
384 }
385 else
386 {
387 // Check for case where completely enveloping clipping volume
388 // If point inside then we are confident that the solid completely
389 // envelopes the clipping volume. Hence set min/max extents according
390 // to clipping volume extents along the specified axis.
391
392 G4ThreeVector clipCentre(
393 (pVoxelLimit.GetMinXExtent()+pVoxelLimit.GetMaxXExtent())*0.5,
394 (pVoxelLimit.GetMinYExtent()+pVoxelLimit.GetMaxYExtent())*0.5,
395 (pVoxelLimit.GetMinZExtent()+pVoxelLimit.GetMaxZExtent())*0.5 );
396
397 if ( Inside(pTransform.Inverse().TransformPoint(clipCentre)) != kOutside )
398 {
399 existsAfterClip = true;
400 pMin = pVoxelLimit.GetMinExtent(pAxis);
401 pMax = pVoxelLimit.GetMaxExtent(pAxis);
402 }
403 }
404 delete vertices;
405 return existsAfterClip;
406 }
407}
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:85
G4bool IsRotated() const
G4AffineTransform Inverse() const
G4ThreeVector NetTranslation() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
EInside Inside(const G4ThreeVector &p) const
Definition: G4Tubs.cc:414
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform) const
Definition: G4Tubs.cc:1721
void ClipBetweenSections(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:376
G4double kCarTolerance
Definition: G4VSolid.hh:307
void ClipCrossSection(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:345
G4double GetMinExtent(const EAxis pAxis) const
G4bool IsYLimited() const
G4double GetMinZExtent() const
G4bool IsXLimited() const
G4double GetMaxExtent(const EAxis pAxis) const
G4double GetMaxYExtent() const
G4double GetMaxZExtent() const
G4double GetMinYExtent() const
G4double GetMinXExtent() const
G4bool IsZLimited() const
G4double GetMaxXExtent() const
@ kYAxis
Definition: geomdefs.hh:54
@ kXAxis
Definition: geomdefs.hh:54
@ kZAxis
Definition: geomdefs.hh:54
@ kOutside
Definition: geomdefs.hh:58

◆ CheckDPhiAngle()

void G4Tubs::CheckDPhiAngle ( G4double  dPhi)
inlineprotected

◆ CheckPhiAngles()

void G4Tubs::CheckPhiAngles ( G4double  sPhi,
G4double  dPhi 
)
inlineprotected

Referenced by G4Tubs().

◆ CheckSPhiAngle()

void G4Tubs::CheckSPhiAngle ( G4double  sPhi)
inlineprotected

◆ Clone()

G4VSolid * G4Tubs::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1816 of file G4Tubs.cc.

1817{
1818 return new G4Tubs(*this);
1819}
Definition: G4Tubs.hh:77

◆ ComputeDimensions()

void G4Tubs::ComputeDimensions ( G4VPVParameterisation p,
const G4int  n,
const G4VPhysicalVolume pRep 
)
virtual

Reimplemented from G4VSolid.

Definition at line 189 of file G4Tubs.cc.

192{
193 p->ComputeDimensions(*this,n,pRep) ;
194}
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const

◆ CreateNURBS()

G4NURBS * G4Tubs::CreateNURBS ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1928 of file G4Tubs.cc.

1929{
1930 G4NURBS* pNURBS ;
1931 if (fRMin != 0)
1932 {
1933 if (fPhiFullTube)
1934 {
1935 pNURBS = new G4NURBStube (fRMin,fRMax,fDz) ;
1936 }
1937 else
1938 {
1940 }
1941 }
1942 else
1943 {
1944 if (fPhiFullTube)
1945 {
1946 pNURBS = new G4NURBScylinder (fRMax,fDz) ;
1947 }
1948 else
1949 {
1950 const G4double epsilon = 1.e-4 ; // Cylinder sector not yet available!
1951 pNURBS = new G4NURBStubesector (epsilon,fRMax,fDz,fSPhi,fSPhi+fDPhi) ;
1952 }
1953 }
1954 return pNURBS ;
1955}

◆ CreatePolyhedron()

G4Polyhedron * G4Tubs::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1923 of file G4Tubs.cc.

1924{
1925 return new G4PolyhedronTubs (fRMin, fRMax, fDz, fSPhi, fDPhi) ;
1926}

Referenced by G4GMocrenFileSceneHandler::AddSolid(), G4CutTubs::CreatePolyhedron(), and G4ArrowModel::G4ArrowModel().

◆ CreateRotatedVertices()

G4ThreeVectorList * G4Tubs::CreateRotatedVertices ( const G4AffineTransform pTransform) const
protected

Definition at line 1721 of file G4Tubs.cc.

1722{
1723 G4ThreeVectorList* vertices ;
1724 G4ThreeVector vertex0, vertex1, vertex2, vertex3 ;
1725 G4double meshAngle, meshRMax, crossAngle,
1726 cosCrossAngle, sinCrossAngle, sAngle;
1727 G4double rMaxX, rMaxY, rMinX, rMinY, meshRMin ;
1728 G4int crossSection, noCrossSections;
1729
1730 // Compute no of cross-sections necessary to mesh tube
1731 //
1732 noCrossSections = G4int(fDPhi/kMeshAngleDefault) + 1 ;
1733
1734 if ( noCrossSections < kMinMeshSections )
1735 {
1736 noCrossSections = kMinMeshSections ;
1737 }
1738 else if (noCrossSections>kMaxMeshSections)
1739 {
1740 noCrossSections = kMaxMeshSections ;
1741 }
1742 // noCrossSections = 4 ;
1743
1744 meshAngle = fDPhi/(noCrossSections - 1) ;
1745 // meshAngle = fDPhi/(noCrossSections) ;
1746
1747 meshRMax = (fRMax+100*kCarTolerance)/std::cos(meshAngle*0.5) ;
1748 meshRMin = fRMin - 100*kCarTolerance ;
1749
1750 // If complete in phi, set start angle such that mesh will be at fRMax
1751 // on the x axis. Will give better extent calculations when not rotated.
1752
1753 if (fPhiFullTube && (fSPhi == 0) ) { sAngle = -meshAngle*0.5 ; }
1754 else { sAngle = fSPhi ; }
1755
1756 vertices = new G4ThreeVectorList();
1757
1758 if ( vertices )
1759 {
1760 vertices->reserve(noCrossSections*4);
1761 for (crossSection = 0 ; crossSection < noCrossSections ; crossSection++ )
1762 {
1763 // Compute coordinates of cross section at section crossSection
1764
1765 crossAngle = sAngle + crossSection*meshAngle ;
1766 cosCrossAngle = std::cos(crossAngle) ;
1767 sinCrossAngle = std::sin(crossAngle) ;
1768
1769 rMaxX = meshRMax*cosCrossAngle ;
1770 rMaxY = meshRMax*sinCrossAngle ;
1771
1772 if(meshRMin <= 0.0)
1773 {
1774 rMinX = 0.0 ;
1775 rMinY = 0.0 ;
1776 }
1777 else
1778 {
1779 rMinX = meshRMin*cosCrossAngle ;
1780 rMinY = meshRMin*sinCrossAngle ;
1781 }
1782 vertex0 = G4ThreeVector(rMinX,rMinY,-fDz) ;
1783 vertex1 = G4ThreeVector(rMaxX,rMaxY,-fDz) ;
1784 vertex2 = G4ThreeVector(rMaxX,rMaxY,+fDz) ;
1785 vertex3 = G4ThreeVector(rMinX,rMinY,+fDz) ;
1786
1787 vertices->push_back(pTransform.TransformPoint(vertex0)) ;
1788 vertices->push_back(pTransform.TransformPoint(vertex1)) ;
1789 vertices->push_back(pTransform.TransformPoint(vertex2)) ;
1790 vertices->push_back(pTransform.TransformPoint(vertex3)) ;
1791 }
1792 }
1793 else
1794 {
1795 DumpInfo();
1796 G4Exception("G4Tubs::CreateRotatedVertices()",
1797 "GeomSolids0003", FatalException,
1798 "Error in allocation of vertices. Out of memory !");
1799 }
1800 return vertices ;
1801}
const G4int kMaxMeshSections
Definition: meshdefs.hh:46
const G4int kMinMeshSections
Definition: meshdefs.hh:45
const G4double kMeshAngleDefault
Definition: meshdefs.hh:42

Referenced by CalculateExtent().

◆ DescribeYourselfTo()

void G4Tubs::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 1918 of file G4Tubs.cc.

1919{
1920 scene.AddSolid (*this) ;
1921}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

G4double G4Tubs::DistanceToIn ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 1189 of file G4Tubs.cc.

1190{
1191 G4double safe=0.0, rho, safe1, safe2, safe3 ;
1192 G4double safePhi, cosPsi ;
1193
1194 rho = std::sqrt(p.x()*p.x() + p.y()*p.y()) ;
1195 safe1 = fRMin - rho ;
1196 safe2 = rho - fRMax ;
1197 safe3 = std::fabs(p.z()) - fDz ;
1198
1199 if ( safe1 > safe2 ) { safe = safe1; }
1200 else { safe = safe2; }
1201 if ( safe3 > safe ) { safe = safe3; }
1202
1203 if ( (!fPhiFullTube) && (rho) )
1204 {
1205 // Psi=angle from central phi to point
1206 //
1207 cosPsi = (p.x()*cosCPhi + p.y()*sinCPhi)/rho ;
1208
1209 if ( cosPsi < std::cos(fDPhi*0.5) )
1210 {
1211 // Point lies outside phi range
1212
1213 if ( (p.y()*cosCPhi - p.x()*sinCPhi) <= 0 )
1214 {
1215 safePhi = std::fabs(p.x()*sinSPhi - p.y()*cosSPhi) ;
1216 }
1217 else
1218 {
1219 safePhi = std::fabs(p.x()*sinEPhi - p.y()*cosEPhi) ;
1220 }
1221 if ( safePhi > safe ) { safe = safePhi; }
1222 }
1223 }
1224 if ( safe < 0 ) { safe = 0; }
1225 return safe ;
1226}

◆ DistanceToIn() [2/2]

G4double G4Tubs::DistanceToIn ( const G4ThreeVector p,
const G4ThreeVector v 
) const
virtual

Implements G4VSolid.

Definition at line 810 of file G4Tubs.cc.

812{
813 G4double snxt = kInfinity ; // snxt = default return value
814 G4double tolORMin2, tolIRMax2 ; // 'generous' radii squared
815 G4double tolORMax2, tolIRMin2, tolODz, tolIDz ;
816 const G4double dRmax = 100.*fRMax;
817
818 static const G4double halfCarTolerance = 0.5*kCarTolerance;
819 static const G4double halfRadTolerance = 0.5*kRadTolerance;
820
821 // Intersection point variables
822 //
823 G4double Dist, sd, xi, yi, zi, rho2, inum, iden, cosPsi, Comp ;
824 G4double t1, t2, t3, b, c, d ; // Quadratic solver variables
825
826 // Calculate tolerant rmin and rmax
827
828 if (fRMin > kRadTolerance)
829 {
830 tolORMin2 = (fRMin - halfRadTolerance)*(fRMin - halfRadTolerance) ;
831 tolIRMin2 = (fRMin + halfRadTolerance)*(fRMin + halfRadTolerance) ;
832 }
833 else
834 {
835 tolORMin2 = 0.0 ;
836 tolIRMin2 = 0.0 ;
837 }
838 tolORMax2 = (fRMax + halfRadTolerance)*(fRMax + halfRadTolerance) ;
839 tolIRMax2 = (fRMax - halfRadTolerance)*(fRMax - halfRadTolerance) ;
840
841 // Intersection with Z surfaces
842
843 tolIDz = fDz - halfCarTolerance ;
844 tolODz = fDz + halfCarTolerance ;
845
846 if (std::fabs(p.z()) >= tolIDz)
847 {
848 if ( p.z()*v.z() < 0 ) // at +Z going in -Z or visa versa
849 {
850 sd = (std::fabs(p.z()) - fDz)/std::fabs(v.z()) ; // Z intersect distance
851
852 if(sd < 0.0) { sd = 0.0; }
853
854 xi = p.x() + sd*v.x() ; // Intersection coords
855 yi = p.y() + sd*v.y() ;
856 rho2 = xi*xi + yi*yi ;
857
858 // Check validity of intersection
859
860 if ((tolIRMin2 <= rho2) && (rho2 <= tolIRMax2))
861 {
862 if (!fPhiFullTube && rho2)
863 {
864 // Psi = angle made with central (average) phi of shape
865 //
866 inum = xi*cosCPhi + yi*sinCPhi ;
867 iden = std::sqrt(rho2) ;
868 cosPsi = inum/iden ;
869 if (cosPsi >= cosHDPhiIT) { return sd ; }
870 }
871 else
872 {
873 return sd ;
874 }
875 }
876 }
877 else
878 {
879 if ( snxt<halfCarTolerance ) { snxt=0; }
880 return snxt ; // On/outside extent, and heading away
881 // -> cannot intersect
882 }
883 }
884
885 // -> Can not intersect z surfaces
886 //
887 // Intersection with rmax (possible return) and rmin (must also check phi)
888 //
889 // Intersection point (xi,yi,zi) on line x=p.x+t*v.x etc.
890 //
891 // Intersects with x^2+y^2=R^2
892 //
893 // 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
894 // t1 t2 t3
895
896 t1 = 1.0 - v.z()*v.z() ;
897 t2 = p.x()*v.x() + p.y()*v.y() ;
898 t3 = p.x()*p.x() + p.y()*p.y() ;
899
900 if ( t1 > 0 ) // Check not || to z axis
901 {
902 b = t2/t1 ;
903 c = t3 - fRMax*fRMax ;
904 if ((t3 >= tolORMax2) && (t2<0)) // This also handles the tangent case
905 {
906 // Try outer cylinder intersection
907 // c=(t3-fRMax*fRMax)/t1;
908
909 c /= t1 ;
910 d = b*b - c ;
911
912 if (d >= 0) // If real root
913 {
914 sd = c/(-b+std::sqrt(d));
915 if (sd >= 0) // If 'forwards'
916 {
917 if ( sd>dRmax ) // Avoid rounding errors due to precision issues on
918 { // 64 bits systems. Split long distances and recompute
919 G4double fTerm = sd-std::fmod(sd,dRmax);
920 sd = fTerm + DistanceToIn(p+fTerm*v,v);
921 }
922 // Check z intersection
923 //
924 zi = p.z() + sd*v.z() ;
925 if (std::fabs(zi)<=tolODz)
926 {
927 // Z ok. Check phi intersection if reqd
928 //
929 if (fPhiFullTube)
930 {
931 return sd ;
932 }
933 else
934 {
935 xi = p.x() + sd*v.x() ;
936 yi = p.y() + sd*v.y() ;
937 cosPsi = (xi*cosCPhi + yi*sinCPhi)/fRMax ;
938 if (cosPsi >= cosHDPhiIT) { return sd ; }
939 }
940 } // end if std::fabs(zi)
941 } // end if (sd>=0)
942 } // end if (d>=0)
943 } // end if (r>=fRMax)
944 else
945 {
946 // Inside outer radius :
947 // check not inside, and heading through tubs (-> 0 to in)
948
949 if ((t3 > tolIRMin2) && (t2 < 0) && (std::fabs(p.z()) <= tolIDz))
950 {
951 // Inside both radii, delta r -ve, inside z extent
952
953 if (!fPhiFullTube)
954 {
955 inum = p.x()*cosCPhi + p.y()*sinCPhi ;
956 iden = std::sqrt(t3) ;
957 cosPsi = inum/iden ;
958 if (cosPsi >= cosHDPhiIT)
959 {
960 // In the old version, the small negative tangent for the point
961 // on surface was not taken in account, and returning 0.0 ...
962 // New version: check the tangent for the point on surface and
963 // if no intersection, return kInfinity, if intersection instead
964 // return sd.
965 //
966 c = t3-fRMax*fRMax;
967 if ( c<=0.0 )
968 {
969 return 0.0;
970 }
971 else
972 {
973 c = c/t1 ;
974 d = b*b-c;
975 if ( d>=0.0 )
976 {
977 snxt = c/(-b+std::sqrt(d)); // using safe solution
978 // for quadratic equation
979 if ( snxt < halfCarTolerance ) { snxt=0; }
980 return snxt ;
981 }
982 else
983 {
984 return kInfinity;
985 }
986 }
987 }
988 }
989 else
990 {
991 // In the old version, the small negative tangent for the point
992 // on surface was not taken in account, and returning 0.0 ...
993 // New version: check the tangent for the point on surface and
994 // if no intersection, return kInfinity, if intersection instead
995 // return sd.
996 //
997 c = t3 - fRMax*fRMax;
998 if ( c<=0.0 )
999 {
1000 return 0.0;
1001 }
1002 else
1003 {
1004 c = c/t1 ;
1005 d = b*b-c;
1006 if ( d>=0.0 )
1007 {
1008 snxt= c/(-b+std::sqrt(d)); // using safe solution
1009 // for quadratic equation
1010 if ( snxt < halfCarTolerance ) { snxt=0; }
1011 return snxt ;
1012 }
1013 else
1014 {
1015 return kInfinity;
1016 }
1017 }
1018 } // end if (!fPhiFullTube)
1019 } // end if (t3>tolIRMin2)
1020 } // end if (Inside Outer Radius)
1021 if ( fRMin ) // Try inner cylinder intersection
1022 {
1023 c = (t3 - fRMin*fRMin)/t1 ;
1024 d = b*b - c ;
1025 if ( d >= 0.0 ) // If real root
1026 {
1027 // Always want 2nd root - we are outside and know rmax Hit was bad
1028 // - If on surface of rmin also need farthest root
1029
1030 sd =( b > 0. )? c/(-b - std::sqrt(d)) : (-b + std::sqrt(d));
1031 if (sd >= -halfCarTolerance) // check forwards
1032 {
1033 // Check z intersection
1034 //
1035 if(sd < 0.0) { sd = 0.0; }
1036 if ( sd>dRmax ) // Avoid rounding errors due to precision issues seen
1037 { // 64 bits systems. Split long distances and recompute
1038 G4double fTerm = sd-std::fmod(sd,dRmax);
1039 sd = fTerm + DistanceToIn(p+fTerm*v,v);
1040 }
1041 zi = p.z() + sd*v.z() ;
1042 if (std::fabs(zi) <= tolODz)
1043 {
1044 // Z ok. Check phi
1045 //
1046 if ( fPhiFullTube )
1047 {
1048 return sd ;
1049 }
1050 else
1051 {
1052 xi = p.x() + sd*v.x() ;
1053 yi = p.y() + sd*v.y() ;
1054 cosPsi = (xi*cosCPhi + yi*sinCPhi)/fRMin ;
1055 if (cosPsi >= cosHDPhiIT)
1056 {
1057 // Good inner radius isect
1058 // - but earlier phi isect still possible
1059
1060 snxt = sd ;
1061 }
1062 }
1063 } // end if std::fabs(zi)
1064 } // end if (sd>=0)
1065 } // end if (d>=0)
1066 } // end if (fRMin)
1067 }
1068
1069 // Phi segment intersection
1070 //
1071 // o Tolerant of points inside phi planes by up to kCarTolerance*0.5
1072 //
1073 // o NOTE: Large duplication of code between sphi & ephi checks
1074 // -> only diffs: sphi -> ephi, Comp -> -Comp and half-plane
1075 // intersection check <=0 -> >=0
1076 // -> use some form of loop Construct ?
1077 //
1078 if ( !fPhiFullTube )
1079 {
1080 // First phi surface (Starting phi)
1081 //
1082 Comp = v.x()*sinSPhi - v.y()*cosSPhi ;
1083
1084 if ( Comp < 0 ) // Component in outwards normal dirn
1085 {
1086 Dist = (p.y()*cosSPhi - p.x()*sinSPhi) ;
1087
1088 if ( Dist < halfCarTolerance )
1089 {
1090 sd = Dist/Comp ;
1091
1092 if (sd < snxt)
1093 {
1094 if ( sd < 0 ) { sd = 0.0; }
1095 zi = p.z() + sd*v.z() ;
1096 if ( std::fabs(zi) <= tolODz )
1097 {
1098 xi = p.x() + sd*v.x() ;
1099 yi = p.y() + sd*v.y() ;
1100 rho2 = xi*xi + yi*yi ;
1101
1102 if ( ( (rho2 >= tolIRMin2) && (rho2 <= tolIRMax2) )
1103 || ( (rho2 > tolORMin2) && (rho2 < tolIRMin2)
1104 && ( v.y()*cosSPhi - v.x()*sinSPhi > 0 )
1105 && ( v.x()*cosSPhi + v.y()*sinSPhi >= 0 ) )
1106 || ( (rho2 > tolIRMax2) && (rho2 < tolORMax2)
1107 && (v.y()*cosSPhi - v.x()*sinSPhi > 0)
1108 && (v.x()*cosSPhi + v.y()*sinSPhi < 0) ) )
1109 {
1110 // z and r intersections good
1111 // - check intersecting with correct half-plane
1112 //
1113 if ((yi*cosCPhi-xi*sinCPhi) <= halfCarTolerance) { snxt = sd; }
1114 }
1115 }
1116 }
1117 }
1118 }
1119
1120 // Second phi surface (Ending phi)
1121
1122 Comp = -(v.x()*sinEPhi - v.y()*cosEPhi) ;
1123
1124 if (Comp < 0 ) // Component in outwards normal dirn
1125 {
1126 Dist = -(p.y()*cosEPhi - p.x()*sinEPhi) ;
1127
1128 if ( Dist < halfCarTolerance )
1129 {
1130 sd = Dist/Comp ;
1131
1132 if (sd < snxt)
1133 {
1134 if ( sd < 0 ) { sd = 0; }
1135 zi = p.z() + sd*v.z() ;
1136 if ( std::fabs(zi) <= tolODz )
1137 {
1138 xi = p.x() + sd*v.x() ;
1139 yi = p.y() + sd*v.y() ;
1140 rho2 = xi*xi + yi*yi ;
1141 if ( ( (rho2 >= tolIRMin2) && (rho2 <= tolIRMax2) )
1142 || ( (rho2 > tolORMin2) && (rho2 < tolIRMin2)
1143 && (v.x()*sinEPhi - v.y()*cosEPhi > 0)
1144 && (v.x()*cosEPhi + v.y()*sinEPhi >= 0) )
1145 || ( (rho2 > tolIRMax2) && (rho2 < tolORMax2)
1146 && (v.x()*sinEPhi - v.y()*cosEPhi > 0)
1147 && (v.x()*cosEPhi + v.y()*sinEPhi < 0) ) )
1148 {
1149 // z and r intersections good
1150 // - check intersecting with correct half-plane
1151 //
1152 if ( (yi*cosCPhi-xi*sinCPhi) >= 0 ) { snxt = sd; }
1153 } //?? >=-halfCarTolerance
1154 }
1155 }
1156 }
1157 } // Comp < 0
1158 } // !fPhiFullTube
1159 if ( snxt<halfCarTolerance ) { snxt=0; }
1160 return snxt ;
1161}
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Tubs.cc:810

Referenced by DistanceToIn().

◆ DistanceToOut() [1/2]

G4double G4Tubs::DistanceToOut ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 1653 of file G4Tubs.cc.

1654{
1655 G4double safe=0.0, rho, safeR1, safeR2, safeZ, safePhi ;
1656 rho = std::sqrt(p.x()*p.x() + p.y()*p.y()) ;
1657
1658#ifdef G4CSGDEBUG
1659 if( Inside(p) == kOutside )
1660 {
1661 G4int oldprc = G4cout.precision(16) ;
1662 G4cout << G4endl ;
1663 DumpInfo();
1664 G4cout << "Position:" << G4endl << G4endl ;
1665 G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl ;
1666 G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl ;
1667 G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl ;
1668 G4cout.precision(oldprc) ;
1669 G4Exception("G4Tubs::DistanceToOut(p)", "GeomSolids1002",
1670 JustWarning, "Point p is outside !?");
1671 }
1672#endif
1673
1674 if ( fRMin )
1675 {
1676 safeR1 = rho - fRMin ;
1677 safeR2 = fRMax - rho ;
1678
1679 if ( safeR1 < safeR2 ) { safe = safeR1 ; }
1680 else { safe = safeR2 ; }
1681 }
1682 else
1683 {
1684 safe = fRMax - rho ;
1685 }
1686 safeZ = fDz - std::fabs(p.z()) ;
1687
1688 if ( safeZ < safe ) { safe = safeZ ; }
1689
1690 // Check if phi divided, Calc distances closest phi plane
1691 //
1692 if ( !fPhiFullTube )
1693 {
1694 if ( p.y()*cosCPhi-p.x()*sinCPhi <= 0 )
1695 {
1696 safePhi = -(p.x()*sinSPhi - p.y()*cosSPhi) ;
1697 }
1698 else
1699 {
1700 safePhi = (p.x()*sinEPhi - p.y()*cosEPhi) ;
1701 }
1702 if (safePhi < safe) { safe = safePhi ; }
1703 }
1704 if ( safe < 0 ) { safe = 0 ; }
1705
1706 return safe ;
1707}
G4DLLIMPORT std::ostream G4cout

◆ DistanceToOut() [2/2]

G4double G4Tubs::DistanceToOut ( const G4ThreeVector p,
const G4ThreeVector v,
const G4bool  calcNorm = G4bool(false),
G4bool validNorm = 0,
G4ThreeVector n = 0 
) const
virtual

Implements G4VSolid.

Definition at line 1233 of file G4Tubs.cc.

1238{
1239 ESide side=kNull , sider=kNull, sidephi=kNull ;
1240 G4double snxt, srd=kInfinity, sphi=kInfinity, pdist ;
1241 G4double deltaR, t1, t2, t3, b, c, d2, roMin2 ;
1242
1243 static const G4double halfCarTolerance = kCarTolerance*0.5;
1244 static const G4double halfAngTolerance = kAngTolerance*0.5;
1245
1246 // Vars for phi intersection:
1247
1248 G4double pDistS, compS, pDistE, compE, sphi2, xi, yi, vphi, roi2 ;
1249
1250 // Z plane intersection
1251
1252 if (v.z() > 0 )
1253 {
1254 pdist = fDz - p.z() ;
1255 if ( pdist > halfCarTolerance )
1256 {
1257 snxt = pdist/v.z() ;
1258 side = kPZ ;
1259 }
1260 else
1261 {
1262 if (calcNorm)
1263 {
1264 *n = G4ThreeVector(0,0,1) ;
1265 *validNorm = true ;
1266 }
1267 return snxt = 0 ;
1268 }
1269 }
1270 else if ( v.z() < 0 )
1271 {
1272 pdist = fDz + p.z() ;
1273
1274 if ( pdist > halfCarTolerance )
1275 {
1276 snxt = -pdist/v.z() ;
1277 side = kMZ ;
1278 }
1279 else
1280 {
1281 if (calcNorm)
1282 {
1283 *n = G4ThreeVector(0,0,-1) ;
1284 *validNorm = true ;
1285 }
1286 return snxt = 0.0 ;
1287 }
1288 }
1289 else
1290 {
1291 snxt = kInfinity ; // Travel perpendicular to z axis
1292 side = kNull;
1293 }
1294
1295 // Radial Intersections
1296 //
1297 // Find intersection with cylinders at rmax/rmin
1298 // Intersection point (xi,yi,zi) on line x=p.x+t*v.x etc.
1299 //
1300 // Intersects with x^2+y^2=R^2
1301 //
1302 // 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
1303 //
1304 // t1 t2 t3
1305
1306 t1 = 1.0 - v.z()*v.z() ; // since v normalised
1307 t2 = p.x()*v.x() + p.y()*v.y() ;
1308 t3 = p.x()*p.x() + p.y()*p.y() ;
1309
1310 if ( snxt > 10*(fDz+fRMax) ) { roi2 = 2*fRMax*fRMax; }
1311 else { roi2 = snxt*snxt*t1 + 2*snxt*t2 + t3; } // radius^2 on +-fDz
1312
1313 if ( t1 > 0 ) // Check not parallel
1314 {
1315 // Calculate srd, r exit distance
1316
1317 if ( (t2 >= 0.0) && (roi2 > fRMax*(fRMax + kRadTolerance)) )
1318 {
1319 // Delta r not negative => leaving via rmax
1320
1321 deltaR = t3 - fRMax*fRMax ;
1322
1323 // NOTE: Should use rho-fRMax<-kRadTolerance*0.5
1324 // - avoid sqrt for efficiency
1325
1326 if ( deltaR < -kRadTolerance*fRMax )
1327 {
1328 b = t2/t1 ;
1329 c = deltaR/t1 ;
1330 d2 = b*b-c;
1331 if( d2 >= 0 ) { srd = c/( -b - std::sqrt(d2)); }
1332 else { srd = 0.; }
1333 sider = kRMax ;
1334 }
1335 else
1336 {
1337 // On tolerant boundary & heading outwards (or perpendicular to)
1338 // outer radial surface -> leaving immediately
1339
1340 if ( calcNorm )
1341 {
1342 *n = G4ThreeVector(p.x()/fRMax,p.y()/fRMax,0) ;
1343 *validNorm = true ;
1344 }
1345 return snxt = 0 ; // Leaving by rmax immediately
1346 }
1347 }
1348 else if ( t2 < 0. ) // i.e. t2 < 0; Possible rmin intersection
1349 {
1350 roMin2 = t3 - t2*t2/t1 ; // min ro2 of the plane of movement
1351
1352 if ( fRMin && (roMin2 < fRMin*(fRMin - kRadTolerance)) )
1353 {
1354 deltaR = t3 - fRMin*fRMin ;
1355 b = t2/t1 ;
1356 c = deltaR/t1 ;
1357 d2 = b*b - c ;
1358
1359 if ( d2 >= 0 ) // Leaving via rmin
1360 {
1361 // NOTE: SHould use rho-rmin>kRadTolerance*0.5
1362 // - avoid sqrt for efficiency
1363
1364 if (deltaR > kRadTolerance*fRMin)
1365 {
1366 srd = c/(-b+std::sqrt(d2));
1367 sider = kRMin ;
1368 }
1369 else
1370 {
1371 if ( calcNorm ) { *validNorm = false; } // Concave side
1372 return snxt = 0.0;
1373 }
1374 }
1375 else // No rmin intersect -> must be rmax intersect
1376 {
1377 deltaR = t3 - fRMax*fRMax ;
1378 c = deltaR/t1 ;
1379 d2 = b*b-c;
1380 if( d2 >=0. )
1381 {
1382 srd = -b + std::sqrt(d2) ;
1383 sider = kRMax ;
1384 }
1385 else // Case: On the border+t2<kRadTolerance
1386 // (v is perpendicular to the surface)
1387 {
1388 if (calcNorm)
1389 {
1390 *n = G4ThreeVector(p.x()/fRMax,p.y()/fRMax,0) ;
1391 *validNorm = true ;
1392 }
1393 return snxt = 0.0;
1394 }
1395 }
1396 }
1397 else if ( roi2 > fRMax*(fRMax + kRadTolerance) )
1398 // No rmin intersect -> must be rmax intersect
1399 {
1400 deltaR = t3 - fRMax*fRMax ;
1401 b = t2/t1 ;
1402 c = deltaR/t1;
1403 d2 = b*b-c;
1404 if( d2 >= 0 )
1405 {
1406 srd = -b + std::sqrt(d2) ;
1407 sider = kRMax ;
1408 }
1409 else // Case: On the border+t2<kRadTolerance
1410 // (v is perpendicular to the surface)
1411 {
1412 if (calcNorm)
1413 {
1414 *n = G4ThreeVector(p.x()/fRMax,p.y()/fRMax,0) ;
1415 *validNorm = true ;
1416 }
1417 return snxt = 0.0;
1418 }
1419 }
1420 }
1421
1422 // Phi Intersection
1423
1424 if ( !fPhiFullTube )
1425 {
1426 // add angle calculation with correction
1427 // of the difference in domain of atan2 and Sphi
1428 //
1429 vphi = std::atan2(v.y(),v.x()) ;
1430
1431 if ( vphi < fSPhi - halfAngTolerance ) { vphi += twopi; }
1432 else if ( vphi > fSPhi + fDPhi + halfAngTolerance ) { vphi -= twopi; }
1433
1434
1435 if ( p.x() || p.y() ) // Check if on z axis (rho not needed later)
1436 {
1437 // pDist -ve when inside
1438
1439 pDistS = p.x()*sinSPhi - p.y()*cosSPhi ;
1440 pDistE = -p.x()*sinEPhi + p.y()*cosEPhi ;
1441
1442 // Comp -ve when in direction of outwards normal
1443
1444 compS = -sinSPhi*v.x() + cosSPhi*v.y() ;
1445 compE = sinEPhi*v.x() - cosEPhi*v.y() ;
1446
1447 sidephi = kNull;
1448
1449 if( ( (fDPhi <= pi) && ( (pDistS <= halfCarTolerance)
1450 && (pDistE <= halfCarTolerance) ) )
1451 || ( (fDPhi > pi) && !((pDistS > halfCarTolerance)
1452 && (pDistE > halfCarTolerance) ) ) )
1453 {
1454 // Inside both phi *full* planes
1455
1456 if ( compS < 0 )
1457 {
1458 sphi = pDistS/compS ;
1459
1460 if (sphi >= -halfCarTolerance)
1461 {
1462 xi = p.x() + sphi*v.x() ;
1463 yi = p.y() + sphi*v.y() ;
1464
1465 // Check intersecting with correct half-plane
1466 // (if not -> no intersect)
1467 //
1468 if( (std::fabs(xi)<=kCarTolerance)&&(std::fabs(yi)<=kCarTolerance) )
1469 {
1470 sidephi = kSPhi;
1471 if (((fSPhi-halfAngTolerance)<=vphi)
1472 &&((fSPhi+fDPhi+halfAngTolerance)>=vphi))
1473 {
1474 sphi = kInfinity;
1475 }
1476 }
1477 else if ( yi*cosCPhi-xi*sinCPhi >=0 )
1478 {
1479 sphi = kInfinity ;
1480 }
1481 else
1482 {
1483 sidephi = kSPhi ;
1484 if ( pDistS > -halfCarTolerance )
1485 {
1486 sphi = 0.0 ; // Leave by sphi immediately
1487 }
1488 }
1489 }
1490 else
1491 {
1492 sphi = kInfinity ;
1493 }
1494 }
1495 else
1496 {
1497 sphi = kInfinity ;
1498 }
1499
1500 if ( compE < 0 )
1501 {
1502 sphi2 = pDistE/compE ;
1503
1504 // Only check further if < starting phi intersection
1505 //
1506 if ( (sphi2 > -halfCarTolerance) && (sphi2 < sphi) )
1507 {
1508 xi = p.x() + sphi2*v.x() ;
1509 yi = p.y() + sphi2*v.y() ;
1510
1511 if ((std::fabs(xi)<=kCarTolerance)&&(std::fabs(yi)<=kCarTolerance))
1512 {
1513 // Leaving via ending phi
1514 //
1515 if( !((fSPhi-halfAngTolerance <= vphi)
1516 &&(fSPhi+fDPhi+halfAngTolerance >= vphi)) )
1517 {
1518 sidephi = kEPhi ;
1519 if ( pDistE <= -halfCarTolerance ) { sphi = sphi2 ; }
1520 else { sphi = 0.0 ; }
1521 }
1522 }
1523 else // Check intersecting with correct half-plane
1524
1525 if ( (yi*cosCPhi-xi*sinCPhi) >= 0)
1526 {
1527 // Leaving via ending phi
1528 //
1529 sidephi = kEPhi ;
1530 if ( pDistE <= -halfCarTolerance ) { sphi = sphi2 ; }
1531 else { sphi = 0.0 ; }
1532 }
1533 }
1534 }
1535 }
1536 else
1537 {
1538 sphi = kInfinity ;
1539 }
1540 }
1541 else
1542 {
1543 // On z axis + travel not || to z axis -> if phi of vector direction
1544 // within phi of shape, Step limited by rmax, else Step =0
1545
1546 if ( (fSPhi - halfAngTolerance <= vphi)
1547 && (vphi <= fSPhi + fDPhi + halfAngTolerance ) )
1548 {
1549 sphi = kInfinity ;
1550 }
1551 else
1552 {
1553 sidephi = kSPhi ; // arbitrary
1554 sphi = 0.0 ;
1555 }
1556 }
1557 if (sphi < snxt) // Order intersecttions
1558 {
1559 snxt = sphi ;
1560 side = sidephi ;
1561 }
1562 }
1563 if (srd < snxt) // Order intersections
1564 {
1565 snxt = srd ;
1566 side = sider ;
1567 }
1568 }
1569 if (calcNorm)
1570 {
1571 switch(side)
1572 {
1573 case kRMax:
1574 // Note: returned vector not normalised
1575 // (divide by fRMax for unit vector)
1576 //
1577 xi = p.x() + snxt*v.x() ;
1578 yi = p.y() + snxt*v.y() ;
1579 *n = G4ThreeVector(xi/fRMax,yi/fRMax,0) ;
1580 *validNorm = true ;
1581 break ;
1582
1583 case kRMin:
1584 *validNorm = false ; // Rmin is inconvex
1585 break ;
1586
1587 case kSPhi:
1588 if ( fDPhi <= pi )
1589 {
1591 *validNorm = true ;
1592 }
1593 else
1594 {
1595 *validNorm = false ;
1596 }
1597 break ;
1598
1599 case kEPhi:
1600 if (fDPhi <= pi)
1601 {
1603 *validNorm = true ;
1604 }
1605 else
1606 {
1607 *validNorm = false ;
1608 }
1609 break ;
1610
1611 case kPZ:
1612 *n = G4ThreeVector(0,0,1) ;
1613 *validNorm = true ;
1614 break ;
1615
1616 case kMZ:
1617 *n = G4ThreeVector(0,0,-1) ;
1618 *validNorm = true ;
1619 break ;
1620
1621 default:
1622 G4cout << G4endl ;
1623 DumpInfo();
1624 std::ostringstream message;
1625 G4int oldprc = message.precision(16);
1626 message << "Undefined side for valid surface normal to solid."
1627 << G4endl
1628 << "Position:" << G4endl << G4endl
1629 << "p.x() = " << p.x()/mm << " mm" << G4endl
1630 << "p.y() = " << p.y()/mm << " mm" << G4endl
1631 << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl
1632 << "Direction:" << G4endl << G4endl
1633 << "v.x() = " << v.x() << G4endl
1634 << "v.y() = " << v.y() << G4endl
1635 << "v.z() = " << v.z() << G4endl << G4endl
1636 << "Proposed distance :" << G4endl << G4endl
1637 << "snxt = " << snxt/mm << " mm" << G4endl ;
1638 message.precision(oldprc) ;
1639 G4Exception("G4Tubs::DistanceToOut(p,v,..)", "GeomSolids1002",
1640 JustWarning, message);
1641 break ;
1642 }
1643 }
1644 if ( snxt<halfCarTolerance ) { snxt=0 ; }
1645
1646 return snxt ;
1647}
ESide
Definition: G4Cons.cc:68

◆ GetCubicVolume()

G4double G4Tubs::GetCubicVolume ( )
inlinevirtual

Reimplemented from G4VSolid.

◆ GetDeltaPhiAngle()

◆ GetDPhi()

G4double G4Tubs::GetDPhi ( ) const
inline

◆ GetDz()

G4double G4Tubs::GetDz ( ) const
inline

◆ GetEntityType()

G4GeometryType G4Tubs::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 1807 of file G4Tubs.cc.

1808{
1809 return G4String("G4Tubs");
1810}

◆ GetInnerRadius()

◆ GetOuterRadius()

◆ GetPointOnSurface()

G4ThreeVector G4Tubs::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1848 of file G4Tubs.cc.

1849{
1850 G4double xRand, yRand, zRand, phi, cosphi, sinphi, chose,
1851 aOne, aTwo, aThr, aFou;
1852 G4double rRand;
1853
1854 aOne = 2.*fDz*fDPhi*fRMax;
1855 aTwo = 2.*fDz*fDPhi*fRMin;
1856 aThr = 0.5*fDPhi*(fRMax*fRMax-fRMin*fRMin);
1857 aFou = 2.*fDz*(fRMax-fRMin);
1858
1860 cosphi = std::cos(phi);
1861 sinphi = std::sin(phi);
1862
1863 rRand = GetRadiusInRing(fRMin,fRMax);
1864
1865 if( (fSPhi == 0) && (fDPhi == twopi) ) { aFou = 0; }
1866
1867 chose = RandFlat::shoot(0.,aOne+aTwo+2.*aThr+2.*aFou);
1868
1869 if( (chose >=0) && (chose < aOne) )
1870 {
1871 xRand = fRMax*cosphi;
1872 yRand = fRMax*sinphi;
1873 zRand = RandFlat::shoot(-1.*fDz,fDz);
1874 return G4ThreeVector (xRand, yRand, zRand);
1875 }
1876 else if( (chose >= aOne) && (chose < aOne + aTwo) )
1877 {
1878 xRand = fRMin*cosphi;
1879 yRand = fRMin*sinphi;
1880 zRand = RandFlat::shoot(-1.*fDz,fDz);
1881 return G4ThreeVector (xRand, yRand, zRand);
1882 }
1883 else if( (chose >= aOne + aTwo) && (chose < aOne + aTwo + aThr) )
1884 {
1885 xRand = rRand*cosphi;
1886 yRand = rRand*sinphi;
1887 zRand = fDz;
1888 return G4ThreeVector (xRand, yRand, zRand);
1889 }
1890 else if( (chose >= aOne + aTwo + aThr) && (chose < aOne + aTwo + 2.*aThr) )
1891 {
1892 xRand = rRand*cosphi;
1893 yRand = rRand*sinphi;
1894 zRand = -1.*fDz;
1895 return G4ThreeVector (xRand, yRand, zRand);
1896 }
1897 else if( (chose >= aOne + aTwo + 2.*aThr)
1898 && (chose < aOne + aTwo + 2.*aThr + aFou) )
1899 {
1900 xRand = rRand*std::cos(fSPhi);
1901 yRand = rRand*std::sin(fSPhi);
1902 zRand = RandFlat::shoot(-1.*fDz,fDz);
1903 return G4ThreeVector (xRand, yRand, zRand);
1904 }
1905 else
1906 {
1907 xRand = rRand*std::cos(fSPhi+fDPhi);
1908 yRand = rRand*std::sin(fSPhi+fDPhi);
1909 zRand = RandFlat::shoot(-1.*fDz,fDz);
1910 return G4ThreeVector (xRand, yRand, zRand);
1911 }
1912}
static double shoot()
Definition: RandFlat.cc:59
G4double GetRadiusInRing(G4double rmin, G4double rmax) const
Definition: G4CSGSolid.cc:101

◆ GetRMax()

G4double G4Tubs::GetRMax ( ) const
inline

◆ GetRMin()

G4double G4Tubs::GetRMin ( ) const
inline

◆ GetSPhi()

G4double G4Tubs::GetSPhi ( ) const
inline

◆ GetStartPhiAngle()

◆ GetSurfaceArea()

G4double G4Tubs::GetSurfaceArea ( )
inlinevirtual

Reimplemented from G4VSolid.

◆ GetZHalfLength()

◆ Initialize()

void G4Tubs::Initialize ( )
inlineprotected

◆ InitializeTrigonometry()

void G4Tubs::InitializeTrigonometry ( )
inlineprotected

◆ Inside()

EInside G4Tubs::Inside ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 414 of file G4Tubs.cc.

415{
416 G4double r2,pPhi,tolRMin,tolRMax;
417 EInside in = kOutside ;
418 static const G4double halfCarTolerance=kCarTolerance*0.5;
419 static const G4double halfRadTolerance=kRadTolerance*0.5;
420 static const G4double halfAngTolerance=kAngTolerance*0.5;
421
422 if (std::fabs(p.z()) <= fDz - halfCarTolerance)
423 {
424 r2 = p.x()*p.x() + p.y()*p.y() ;
425
426 if (fRMin) { tolRMin = fRMin + halfRadTolerance ; }
427 else { tolRMin = 0 ; }
428
429 tolRMax = fRMax - halfRadTolerance ;
430
431 if ((r2 >= tolRMin*tolRMin) && (r2 <= tolRMax*tolRMax))
432 {
433 if ( fPhiFullTube )
434 {
435 in = kInside ;
436 }
437 else
438 {
439 // Try inner tolerant phi boundaries (=>inside)
440 // if not inside, try outer tolerant phi boundaries
441
442 if ( (tolRMin==0) && (std::fabs(p.x())<=halfCarTolerance)
443 && (std::fabs(p.y())<=halfCarTolerance) )
444 {
445 in=kSurface;
446 }
447 else
448 {
449 pPhi = std::atan2(p.y(),p.x()) ;
450 if ( pPhi < -halfAngTolerance ) { pPhi += twopi; } // 0<=pPhi<2pi
451
452 if ( fSPhi >= 0 )
453 {
454 if ( (std::fabs(pPhi) < halfAngTolerance)
455 && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
456 {
457 pPhi += twopi ; // 0 <= pPhi < 2pi
458 }
459 if ( (pPhi >= fSPhi + halfAngTolerance)
460 && (pPhi <= fSPhi + fDPhi - halfAngTolerance) )
461 {
462 in = kInside ;
463 }
464 else if ( (pPhi >= fSPhi - halfAngTolerance)
465 && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
466 {
467 in = kSurface ;
468 }
469 }
470 else // fSPhi < 0
471 {
472 if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
473 && (pPhi >= fSPhi + fDPhi + halfAngTolerance) ) {;} //kOutside
474 else if ( (pPhi <= fSPhi + twopi + halfAngTolerance)
475 && (pPhi >= fSPhi + fDPhi - halfAngTolerance) )
476 {
477 in = kSurface ;
478 }
479 else
480 {
481 in = kInside ;
482 }
483 }
484 }
485 }
486 }
487 else // Try generous boundaries
488 {
489 tolRMin = fRMin - halfRadTolerance ;
490 tolRMax = fRMax + halfRadTolerance ;
491
492 if ( tolRMin < 0 ) { tolRMin = 0; }
493
494 if ( (r2 >= tolRMin*tolRMin) && (r2 <= tolRMax*tolRMax) )
495 {
496 if (fPhiFullTube || (r2 <=halfRadTolerance*halfRadTolerance) )
497 { // Continuous in phi or on z-axis
498 in = kSurface ;
499 }
500 else // Try outer tolerant phi boundaries only
501 {
502 pPhi = std::atan2(p.y(),p.x()) ;
503
504 if ( pPhi < -halfAngTolerance) { pPhi += twopi; } // 0<=pPhi<2pi
505 if ( fSPhi >= 0 )
506 {
507 if ( (std::fabs(pPhi) < halfAngTolerance)
508 && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
509 {
510 pPhi += twopi ; // 0 <= pPhi < 2pi
511 }
512 if ( (pPhi >= fSPhi - halfAngTolerance)
513 && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
514 {
515 in = kSurface ;
516 }
517 }
518 else // fSPhi < 0
519 {
520 if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
521 && (pPhi >= fSPhi + fDPhi + halfAngTolerance) ) {;} // kOutside
522 else
523 {
524 in = kSurface ;
525 }
526 }
527 }
528 }
529 }
530 }
531 else if (std::fabs(p.z()) <= fDz + halfCarTolerance)
532 { // Check within tolerant r limits
533 r2 = p.x()*p.x() + p.y()*p.y() ;
534 tolRMin = fRMin - halfRadTolerance ;
535 tolRMax = fRMax + halfRadTolerance ;
536
537 if ( tolRMin < 0 ) { tolRMin = 0; }
538
539 if ( (r2 >= tolRMin*tolRMin) && (r2 <= tolRMax*tolRMax) )
540 {
541 if (fPhiFullTube || (r2 <=halfRadTolerance*halfRadTolerance))
542 { // Continuous in phi or on z-axis
543 in = kSurface ;
544 }
545 else // Try outer tolerant phi boundaries
546 {
547 pPhi = std::atan2(p.y(),p.x()) ;
548
549 if ( pPhi < -halfAngTolerance ) { pPhi += twopi; } // 0<=pPhi<2pi
550 if ( fSPhi >= 0 )
551 {
552 if ( (std::fabs(pPhi) < halfAngTolerance)
553 && (std::fabs(fSPhi + fDPhi - twopi) < halfAngTolerance) )
554 {
555 pPhi += twopi ; // 0 <= pPhi < 2pi
556 }
557 if ( (pPhi >= fSPhi - halfAngTolerance)
558 && (pPhi <= fSPhi + fDPhi + halfAngTolerance) )
559 {
560 in = kSurface;
561 }
562 }
563 else // fSPhi < 0
564 {
565 if ( (pPhi <= fSPhi + twopi - halfAngTolerance)
566 && (pPhi >= fSPhi + fDPhi + halfAngTolerance) ) {;}
567 else
568 {
569 in = kSurface ;
570 }
571 }
572 }
573 }
574 }
575 return in;
576}
EInside
Definition: geomdefs.hh:58
@ kInside
Definition: geomdefs.hh:58
@ kSurface
Definition: geomdefs.hh:58

Referenced by CalculateExtent(), and DistanceToOut().

◆ operator=()

G4Tubs & G4Tubs::operator= ( const G4Tubs rhs)

Definition at line 160 of file G4Tubs.cc.

161{
162 // Check assignment to self
163 //
164 if (this == &rhs) { return *this; }
165
166 // Copy base class data
167 //
169
170 // Copy data
171 //
173 fRMin = rhs.fRMin; fRMax = rhs.fRMax; fDz = rhs.fDz;
174 fSPhi = rhs.fSPhi; fDPhi = rhs.fDPhi;
175 sinCPhi = rhs.sinCPhi; cosCPhi = rhs.sinCPhi;
177 sinSPhi = rhs.sinSPhi; cosSPhi = rhs.cosSPhi;
178 sinEPhi = rhs.sinEPhi; cosEPhi = rhs.cosEPhi;
180
181 return *this;
182}
G4CSGSolid & operator=(const G4CSGSolid &rhs)
Definition: G4CSGSolid.cc:82

Referenced by G4CutTubs::operator=().

◆ SetDeltaPhiAngle()

◆ SetInnerRadius()

◆ SetOuterRadius()

◆ SetStartPhiAngle()

◆ SetZHalfLength()

◆ StreamInfo()

std::ostream & G4Tubs::StreamInfo ( std::ostream &  os) const
virtual

Reimplemented from G4CSGSolid.

Definition at line 1825 of file G4Tubs.cc.

1826{
1827 G4int oldprc = os.precision(16);
1828 os << "-----------------------------------------------------------\n"
1829 << " *** Dump for solid - " << GetName() << " ***\n"
1830 << " ===================================================\n"
1831 << " Solid type: G4Tubs\n"
1832 << " Parameters: \n"
1833 << " inner radius : " << fRMin/mm << " mm \n"
1834 << " outer radius : " << fRMax/mm << " mm \n"
1835 << " half length Z: " << fDz/mm << " mm \n"
1836 << " starting phi : " << fSPhi/degree << " degrees \n"
1837 << " delta phi : " << fDPhi/degree << " degrees \n"
1838 << "-----------------------------------------------------------\n";
1839 os.precision(oldprc);
1840
1841 return os;
1842}

◆ SurfaceNormal()

G4ThreeVector G4Tubs::SurfaceNormal ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 584 of file G4Tubs.cc.

585{
586 G4int noSurfaces = 0;
587 G4double rho, pPhi;
588 G4double distZ, distRMin, distRMax;
589 G4double distSPhi = kInfinity, distEPhi = kInfinity;
590
591 static const G4double halfCarTolerance = 0.5*kCarTolerance;
592 static const G4double halfAngTolerance = 0.5*kAngTolerance;
593
594 G4ThreeVector norm, sumnorm(0.,0.,0.);
595 G4ThreeVector nZ = G4ThreeVector(0, 0, 1.0);
596 G4ThreeVector nR, nPs, nPe;
597
598 rho = std::sqrt(p.x()*p.x() + p.y()*p.y());
599
600 distRMin = std::fabs(rho - fRMin);
601 distRMax = std::fabs(rho - fRMax);
602 distZ = std::fabs(std::fabs(p.z()) - fDz);
603
604 if (!fPhiFullTube) // Protected against (0,0,z)
605 {
606 if ( rho > halfCarTolerance )
607 {
608 pPhi = std::atan2(p.y(),p.x());
609
610 if(pPhi < fSPhi- halfCarTolerance) { pPhi += twopi; }
611 else if(pPhi > fSPhi+fDPhi+ halfCarTolerance) { pPhi -= twopi; }
612
613 distSPhi = std::fabs(pPhi - fSPhi);
614 distEPhi = std::fabs(pPhi - fSPhi - fDPhi);
615 }
616 else if( !fRMin )
617 {
618 distSPhi = 0.;
619 distEPhi = 0.;
620 }
621 nPs = G4ThreeVector(std::sin(fSPhi),-std::cos(fSPhi),0);
622 nPe = G4ThreeVector(-std::sin(fSPhi+fDPhi),std::cos(fSPhi+fDPhi),0);
623 }
624 if ( rho > halfCarTolerance ) { nR = G4ThreeVector(p.x()/rho,p.y()/rho,0); }
625
626 if( distRMax <= halfCarTolerance )
627 {
628 noSurfaces ++;
629 sumnorm += nR;
630 }
631 if( fRMin && (distRMin <= halfCarTolerance) )
632 {
633 noSurfaces ++;
634 sumnorm -= nR;
635 }
636 if( fDPhi < twopi )
637 {
638 if (distSPhi <= halfAngTolerance)
639 {
640 noSurfaces ++;
641 sumnorm += nPs;
642 }
643 if (distEPhi <= halfAngTolerance)
644 {
645 noSurfaces ++;
646 sumnorm += nPe;
647 }
648 }
649 if (distZ <= halfCarTolerance)
650 {
651 noSurfaces ++;
652 if ( p.z() >= 0.) { sumnorm += nZ; }
653 else { sumnorm -= nZ; }
654 }
655 if ( noSurfaces == 0 )
656 {
657#ifdef G4CSGDEBUG
658 G4Exception("G4Tubs::SurfaceNormal(p)", "GeomSolids1002",
659 JustWarning, "Point p is not on surface !?" );
660 G4int oldprc = G4cout.precision(20);
661 G4cout<< "G4Tubs::SN ( "<<p.x()<<", "<<p.y()<<", "<<p.z()<<" ); "
662 << G4endl << G4endl;
663 G4cout.precision(oldprc) ;
664#endif
665 norm = ApproxSurfaceNormal(p);
666 }
667 else if ( noSurfaces == 1 ) { norm = sumnorm; }
668 else { norm = sumnorm.unit(); }
669
670 return norm;
671}
Hep3Vector unit() const
virtual G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tubs.cc:678

Member Data Documentation

◆ cosCPhi

G4double G4Tubs::cosCPhi
protected

◆ cosEPhi

G4double G4Tubs::cosEPhi
protected

◆ cosHDPhiIT

G4double G4Tubs::cosHDPhiIT
protected

Definition at line 213 of file G4Tubs.hh.

Referenced by G4CutTubs::DistanceToIn(), DistanceToIn(), and operator=().

◆ cosHDPhiOT

G4double G4Tubs::cosHDPhiOT
protected

Definition at line 213 of file G4Tubs.hh.

Referenced by operator=().

◆ cosSPhi

G4double G4Tubs::cosSPhi
protected

◆ fDPhi

◆ fDz

◆ fPhiFullTube

G4bool G4Tubs::fPhiFullTube
protected

◆ fRMax

◆ fRMin

◆ fSPhi

◆ kAngTolerance

◆ kRadTolerance

◆ sinCPhi

G4double G4Tubs::sinCPhi
protected

◆ sinEPhi

G4double G4Tubs::sinEPhi
protected

◆ sinSPhi

G4double G4Tubs::sinSPhi
protected

The documentation for this class was generated from the following files: