Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4Trd Class Reference

#include <G4Trd.hh>

Inheritance diagram for G4Trd:
G4CSGSolid G4VSolid

Public Member Functions

void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
 
G4VSolidClone () const
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
G4PolyhedronCreatePolyhedron () const
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
G4double DistanceToIn (const G4ThreeVector &p) const
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
G4double DistanceToOut (const G4ThreeVector &p) const
 
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
 
void DumpInfo () const
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 
 G4Trd (__void__ &)
 
 G4Trd (const G4String &pName, G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
 
 G4Trd (const G4Trd &rhs)
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
G4double GetCubicVolume ()
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
G4GeometryType GetEntityType () const
 
virtual G4VisExtent GetExtent () const
 
G4String GetName () const
 
G4ThreeVector GetPointOnSurface () const
 
virtual G4PolyhedronGetPolyhedron () const
 
G4double GetSurfaceArea ()
 
G4double GetTolerance () const
 
G4double GetXHalfLength1 () const
 
G4double GetXHalfLength2 () const
 
G4double GetYHalfLength1 () const
 
G4double GetYHalfLength2 () const
 
G4double GetZHalfLength () const
 
EInside Inside (const G4ThreeVector &p) const
 
G4Trdoperator= (const G4Trd &rhs)
 
G4bool operator== (const G4VSolid &s) const
 
void SetAllParameters (G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
 
void SetName (const G4String &name)
 
void SetXHalfLength1 (G4double val)
 
void SetXHalfLength2 (G4double val)
 
void SetYHalfLength1 (G4double val)
 
void SetYHalfLength2 (G4double val)
 
void SetZHalfLength (G4double val)
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
 ~G4Trd ()
 

Protected Member Functions

void CalculateClippedPolygonExtent (G4ThreeVectorList &pPolygon, 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 ClipCrossSection (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 GetRadiusInRing (G4double rmin, G4double rmax) const
 

Protected Attributes

G4double fCubicVolume = 0.0
 
G4PolyhedronfpPolyhedron = nullptr
 
G4bool fRebuildPolyhedron = false
 
G4double fSurfaceArea = 0.0
 
G4double kCarTolerance
 

Private Member Functions

G4ThreeVector ApproxSurfaceNormal (const G4ThreeVector &p) const
 
void CheckParameters ()
 
void ClipPolygonToSimpleLimits (G4ThreeVectorList &pPolygon, G4ThreeVectorList &outputPolygon, const G4VoxelLimits &pVoxelLimit) const
 
void MakePlanes ()
 

Private Attributes

G4double fDx1
 
G4double fDx2
 
G4double fDy1
 
G4double fDy2
 
G4double fDz
 
G4double fHx
 
G4double fHy
 
struct {
   G4double   a
 
   G4double   b
 
   G4double   c
 
   G4double   d
 
fPlanes [4]
 
G4String fshapeName
 
G4double halfCarTolerance
 

Detailed Description

Definition at line 62 of file G4Trd.hh.

Constructor & Destructor Documentation

◆ G4Trd() [1/3]

G4Trd::G4Trd ( const G4String pName,
G4double  pdx1,
G4double  pdx2,
G4double  pdy1,
G4double  pdy2,
G4double  pdz 
)

Definition at line 54 of file G4Trd.cc.

59 fDx1(pdx1), fDx2(pdx2), fDy1(pdy1), fDy2(pdy2), fDz(pdz)
60{
62 MakePlanes();
63}
G4CSGSolid(const G4String &pName)
Definition: G4CSGSolid.cc:49
void CheckParameters()
Definition: G4Trd.cc:149
G4double fDy1
Definition: G4Trd.hh:169
G4double fDx1
Definition: G4Trd.hh:169
G4double halfCarTolerance
Definition: G4Trd.hh:168
G4double fDz
Definition: G4Trd.hh:169
G4double fDy2
Definition: G4Trd.hh:169
G4double fDx2
Definition: G4Trd.hh:169
void MakePlanes()
Definition: G4Trd.cc:171
G4double kCarTolerance
Definition: G4VSolid.hh:299

References CheckParameters(), and MakePlanes().

Referenced by Clone().

◆ ~G4Trd()

G4Trd::~G4Trd ( )

Definition at line 81 of file G4Trd.cc.

82{
83}

◆ G4Trd() [2/3]

G4Trd::G4Trd ( __void__ &  a)

Definition at line 70 of file G4Trd.cc.

72 fDx1(1.), fDx2(1.), fDy1(1.), fDy2(1.), fDz(1.)
73{
74 MakePlanes();
75}
G4double a
Definition: G4Trd.hh:170

References MakePlanes().

◆ G4Trd() [3/3]

G4Trd::G4Trd ( const G4Trd rhs)

Definition at line 89 of file G4Trd.cc.

91 fDx1(rhs.fDx1), fDx2(rhs.fDx2),
92 fDy1(rhs.fDy1), fDy2(rhs.fDy2), fDz(rhs.fDz),
93 fHx(rhs.fHx), fHy(rhs.fHy)
94{
95 for (G4int i=0; i<4; ++i) { fPlanes[i] = rhs.fPlanes[i]; }
96}
int G4int
Definition: G4Types.hh:85
G4double fHx
Definition: G4Trd.hh:169
G4double fHy
Definition: G4Trd.hh:169
struct G4Trd::@23 fPlanes[4]

References fPlanes.

Member Function Documentation

◆ ApproxSurfaceNormal()

G4ThreeVector G4Trd::ApproxSurfaceNormal ( const G4ThreeVector p) const
private

Definition at line 428 of file G4Trd.cc.

429{
430 G4double dist = -DBL_MAX;
431 G4int iside = 0;
432 for (G4int i=0; i<4; ++i)
433 {
434 G4double d = fPlanes[i].a*p.x() +
435 fPlanes[i].b*p.y() +
436 fPlanes[i].c*p.z() + fPlanes[i].d;
437 if (d > dist) { dist = d; iside = i; }
438 }
439
440 G4double distz = std::abs(p.z()) - fDz;
441 if (dist > distz)
442 return G4ThreeVector(fPlanes[iside].a, fPlanes[iside].b, fPlanes[iside].c);
443 else
444 return G4ThreeVector(0, 0, (p.z() < 0) ? -1 : 1);
445}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
double z() const
double x() const
double y() const
G4double b
Definition: G4Trd.hh:170
G4double c
Definition: G4Trd.hh:170
G4double d
Definition: G4Trd.hh:170
#define DBL_MAX
Definition: templates.hh:62

References a, b, c, d, DBL_MAX, fDz, fPlanes, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by SurfaceNormal().

◆ BoundingLimits()

void G4Trd::BoundingLimits ( G4ThreeVector pMin,
G4ThreeVector pMax 
) const
virtual

Reimplemented from G4VSolid.

Definition at line 248 of file G4Trd.cc.

249{
255
256 G4double xmax = std::max(dx1,dx2);
257 G4double ymax = std::max(dy1,dy2);
258 pMin.set(-xmax,-ymax,-dz);
259 pMax.set( xmax, ymax, dz);
260
261 // Check correctness of the bounding box
262 //
263 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
264 {
265 std::ostringstream message;
266 message << "Bad bounding box (min >= max) for solid: "
267 << GetName() << " !"
268 << "\npMin = " << pMin
269 << "\npMax = " << pMax;
270 G4Exception("G4Trd::BoundingLimits()", "GeomMgt0001", JustWarning, message);
271 DumpInfo();
272 }
273}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static const G4double pMax
static const G4double pMin
G4double GetXHalfLength2() const
G4double GetYHalfLength2() const
G4double GetXHalfLength1() const
G4double GetYHalfLength1() const
G4double GetZHalfLength() const
G4String GetName() const
void DumpInfo() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References G4VSolid::DumpInfo(), G4Exception(), G4VSolid::GetName(), GetXHalfLength1(), GetXHalfLength2(), GetYHalfLength1(), GetYHalfLength2(), GetZHalfLength(), JustWarning, G4INCL::Math::max(), pMax, and pMin.

Referenced by CalculateExtent().

◆ CalculateClippedPolygonExtent()

void G4VSolid::CalculateClippedPolygonExtent ( G4ThreeVectorList pPolygon,
const G4VoxelLimits pVoxelLimit,
const EAxis  pAxis,
G4double pMin,
G4double pMax 
) const
protectedinherited

Definition at line 489 of file G4VSolid.cc.

494{
495 G4int noLeft,i;
496 G4double component;
497
498 ClipPolygon(pPolygon,pVoxelLimit,pAxis);
499 noLeft = pPolygon.size();
500
501 if ( noLeft )
502 {
503 for (i=0; i<noLeft; ++i)
504 {
505 component = pPolygon[i].operator()(pAxis);
506
507 if (component < pMin)
508 {
509 pMin = component;
510 }
511 if (component > pMax)
512 {
513 pMax = component;
514 }
515 }
516 }
517}
void ClipPolygon(G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
Definition: G4VSolid.cc:539

References G4VSolid::ClipPolygon(), pMax, and pMin.

Referenced by G4VSolid::ClipBetweenSections(), and G4VSolid::ClipCrossSection().

◆ CalculateExtent()

G4bool G4Trd::CalculateExtent ( const EAxis  pAxis,
const G4VoxelLimits pVoxelLimit,
const G4AffineTransform pTransform,
G4double pMin,
G4double pMax 
) const
virtual

Implements G4VSolid.

Definition at line 279 of file G4Trd.cc.

283{
284 G4ThreeVector bmin, bmax;
285 G4bool exist;
286
287 // Check bounding box (bbox)
288 //
289 BoundingLimits(bmin,bmax);
290 G4BoundingEnvelope bbox(bmin,bmax);
291#ifdef G4BBOX_EXTENT
292 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
293#endif
294 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
295 {
296 return exist = (pMin < pMax) ? true : false;
297 }
298
299 // Set bounding envelope (benv) and calculate extent
300 //
306
307 G4ThreeVectorList baseA(4), baseB(4);
308 baseA[0].set(-dx1,-dy1,-dz);
309 baseA[1].set( dx1,-dy1,-dz);
310 baseA[2].set( dx1, dy1,-dz);
311 baseA[3].set(-dx1, dy1,-dz);
312 baseB[0].set(-dx2,-dy2, dz);
313 baseB[1].set( dx2,-dy2, dz);
314 baseB[2].set( dx2, dy2, dz);
315 baseB[3].set(-dx2, dy2, dz);
316
317 std::vector<const G4ThreeVectorList *> polygons(2);
318 polygons[0] = &baseA;
319 polygons[1] = &baseB;
320
321 G4BoundingEnvelope benv(bmin,bmax,polygons);
322 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
323 return exist;
324}
std::vector< G4ThreeVector > G4ThreeVectorList
bool G4bool
Definition: G4Types.hh:86
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Trd.cc:248

References G4BoundingEnvelope::BoundingBoxVsVoxelLimits(), BoundingLimits(), G4BoundingEnvelope::CalculateExtent(), GetXHalfLength1(), GetXHalfLength2(), GetYHalfLength1(), GetYHalfLength2(), GetZHalfLength(), pMax, and pMin.

◆ CheckParameters()

void G4Trd::CheckParameters ( )
private

Definition at line 149 of file G4Trd.cc.

150{
151 G4double dmin = 2*kCarTolerance;
152 if ((fDx1 < 0 || fDx2 < 0 || fDy1 < 0 || fDy2 < 0 || fDz < dmin) ||
153 (fDx1 < dmin && fDx2 < dmin) ||
154 (fDy1 < dmin && fDy2 < dmin))
155 {
156 std::ostringstream message;
157 message << "Invalid (too small or negative) dimensions for Solid: "
158 << GetName()
159 << "\n X - " << fDx1 << ", " << fDx2
160 << "\n Y - " << fDy1 << ", " << fDy2
161 << "\n Z - " << fDz;
162 G4Exception("G4Trd::CheckParameters()", "GeomSolids0002",
163 FatalException, message);
164 }
165}
@ FatalException

References FatalException, fDx1, fDx2, fDy1, fDy2, fDz, G4Exception(), G4VSolid::GetName(), and G4VSolid::kCarTolerance.

Referenced by G4Trd(), and SetAllParameters().

◆ ClipBetweenSections()

void G4VSolid::ClipBetweenSections ( G4ThreeVectorList pVertices,
const G4int  pSectionIndex,
const G4VoxelLimits pVoxelLimit,
const EAxis  pAxis,
G4double pMin,
G4double pMax 
) const
protectedinherited

Definition at line 444 of file G4VSolid.cc.

449{
450 G4ThreeVectorList polygon;
451 polygon.reserve(4);
452 polygon.push_back((*pVertices)[pSectionIndex]);
453 polygon.push_back((*pVertices)[pSectionIndex+4]);
454 polygon.push_back((*pVertices)[pSectionIndex+5]);
455 polygon.push_back((*pVertices)[pSectionIndex+1]);
456 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
457 polygon.clear();
458
459 polygon.push_back((*pVertices)[pSectionIndex+1]);
460 polygon.push_back((*pVertices)[pSectionIndex+5]);
461 polygon.push_back((*pVertices)[pSectionIndex+6]);
462 polygon.push_back((*pVertices)[pSectionIndex+2]);
463 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
464 polygon.clear();
465
466 polygon.push_back((*pVertices)[pSectionIndex+2]);
467 polygon.push_back((*pVertices)[pSectionIndex+6]);
468 polygon.push_back((*pVertices)[pSectionIndex+7]);
469 polygon.push_back((*pVertices)[pSectionIndex+3]);
470 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
471 polygon.clear();
472
473 polygon.push_back((*pVertices)[pSectionIndex+3]);
474 polygon.push_back((*pVertices)[pSectionIndex+7]);
475 polygon.push_back((*pVertices)[pSectionIndex+4]);
476 polygon.push_back((*pVertices)[pSectionIndex]);
477 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
478 return;
479}
void CalculateClippedPolygonExtent(G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:489

References G4VSolid::CalculateClippedPolygonExtent(), pMax, and pMin.

◆ ClipCrossSection()

void G4VSolid::ClipCrossSection ( G4ThreeVectorList pVertices,
const G4int  pSectionIndex,
const G4VoxelLimits pVoxelLimit,
const EAxis  pAxis,
G4double pMin,
G4double pMax 
) const
protectedinherited

Definition at line 414 of file G4VSolid.cc.

419{
420
421 G4ThreeVectorList polygon;
422 polygon.reserve(4);
423 polygon.push_back((*pVertices)[pSectionIndex]);
424 polygon.push_back((*pVertices)[pSectionIndex+1]);
425 polygon.push_back((*pVertices)[pSectionIndex+2]);
426 polygon.push_back((*pVertices)[pSectionIndex+3]);
427 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
428 return;
429}

References G4VSolid::CalculateClippedPolygonExtent(), pMax, and pMin.

◆ ClipPolygon()

void G4VSolid::ClipPolygon ( G4ThreeVectorList pPolygon,
const G4VoxelLimits pVoxelLimit,
const EAxis  pAxis 
) const
protectedinherited

Definition at line 539 of file G4VSolid.cc.

542{
543 G4ThreeVectorList outputPolygon;
544
545 if ( pVoxelLimit.IsLimited() )
546 {
547 if (pVoxelLimit.IsXLimited() ) // && pAxis != kXAxis)
548 {
549 G4VoxelLimits simpleLimit1;
550 simpleLimit1.AddLimit(kXAxis,pVoxelLimit.GetMinXExtent(),kInfinity);
551 ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
552
553 pPolygon.clear();
554
555 if ( !outputPolygon.size() ) return;
556
557 G4VoxelLimits simpleLimit2;
558 simpleLimit2.AddLimit(kXAxis,-kInfinity,pVoxelLimit.GetMaxXExtent());
559 ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
560
561 if ( !pPolygon.size() ) return;
562 else outputPolygon.clear();
563 }
564 if ( pVoxelLimit.IsYLimited() ) // && pAxis != kYAxis)
565 {
566 G4VoxelLimits simpleLimit1;
567 simpleLimit1.AddLimit(kYAxis,pVoxelLimit.GetMinYExtent(),kInfinity);
568 ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
569
570 // Must always clear pPolygon - for clip to simpleLimit2 and in case of
571 // early exit
572
573 pPolygon.clear();
574
575 if ( !outputPolygon.size() ) return;
576
577 G4VoxelLimits simpleLimit2;
578 simpleLimit2.AddLimit(kYAxis,-kInfinity,pVoxelLimit.GetMaxYExtent());
579 ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
580
581 if ( !pPolygon.size() ) return;
582 else outputPolygon.clear();
583 }
584 if ( pVoxelLimit.IsZLimited() ) // && pAxis != kZAxis)
585 {
586 G4VoxelLimits simpleLimit1;
587 simpleLimit1.AddLimit(kZAxis,pVoxelLimit.GetMinZExtent(),kInfinity);
588 ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
589
590 // Must always clear pPolygon - for clip to simpleLimit2 and in case of
591 // early exit
592
593 pPolygon.clear();
594
595 if ( !outputPolygon.size() ) return;
596
597 G4VoxelLimits simpleLimit2;
598 simpleLimit2.AddLimit(kZAxis,-kInfinity,pVoxelLimit.GetMaxZExtent());
599 ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
600
601 // Return after final clip - no cleanup
602 }
603 }
604}
void ClipPolygonToSimpleLimits(G4ThreeVectorList &pPolygon, G4ThreeVectorList &outputPolygon, const G4VoxelLimits &pVoxelLimit) const
Definition: G4VSolid.cc:612
G4bool IsYLimited() const
G4double GetMinZExtent() const
void AddLimit(const EAxis pAxis, const G4double pMin, const G4double pMax)
G4bool IsXLimited() const
G4double GetMaxYExtent() const
G4double GetMaxZExtent() const
G4double GetMinYExtent() const
G4double GetMinXExtent() const
G4bool IsZLimited() const
G4bool IsLimited() const
G4double GetMaxXExtent() const
@ kYAxis
Definition: geomdefs.hh:56
@ kXAxis
Definition: geomdefs.hh:55
@ kZAxis
Definition: geomdefs.hh:57
static const G4double kInfinity
Definition: geomdefs.hh:41

References G4VoxelLimits::AddLimit(), G4VSolid::ClipPolygonToSimpleLimits(), G4VoxelLimits::GetMaxXExtent(), G4VoxelLimits::GetMaxYExtent(), G4VoxelLimits::GetMaxZExtent(), G4VoxelLimits::GetMinXExtent(), G4VoxelLimits::GetMinYExtent(), G4VoxelLimits::GetMinZExtent(), G4VoxelLimits::IsLimited(), G4VoxelLimits::IsXLimited(), G4VoxelLimits::IsYLimited(), G4VoxelLimits::IsZLimited(), kInfinity, kXAxis, kYAxis, and kZAxis.

Referenced by G4VSolid::CalculateClippedPolygonExtent().

◆ ClipPolygonToSimpleLimits()

void G4VSolid::ClipPolygonToSimpleLimits ( G4ThreeVectorList pPolygon,
G4ThreeVectorList outputPolygon,
const G4VoxelLimits pVoxelLimit 
) const
privateinherited

Definition at line 612 of file G4VSolid.cc.

615{
616 G4int i;
617 G4int noVertices=pPolygon.size();
618 G4ThreeVector vEnd,vStart;
619
620 for (i = 0 ; i < noVertices ; ++i )
621 {
622 vStart = pPolygon[i];
623 if ( i == noVertices-1 ) vEnd = pPolygon[0];
624 else vEnd = pPolygon[i+1];
625
626 if ( pVoxelLimit.Inside(vStart) )
627 {
628 if (pVoxelLimit.Inside(vEnd))
629 {
630 // vStart and vEnd inside -> output end point
631 //
632 outputPolygon.push_back(vEnd);
633 }
634 else
635 {
636 // vStart inside, vEnd outside -> output crossing point
637 //
638 pVoxelLimit.ClipToLimits(vStart,vEnd);
639 outputPolygon.push_back(vEnd);
640 }
641 }
642 else
643 {
644 if (pVoxelLimit.Inside(vEnd))
645 {
646 // vStart outside, vEnd inside -> output inside section
647 //
648 pVoxelLimit.ClipToLimits(vStart,vEnd);
649 outputPolygon.push_back(vStart);
650 outputPolygon.push_back(vEnd);
651 }
652 else // Both point outside -> no output
653 {
654 // outputPolygon.push_back(vStart);
655 // outputPolygon.push_back(vEnd);
656 }
657 }
658 }
659}
G4bool ClipToLimits(G4ThreeVector &pStart, G4ThreeVector &pEnd) const
G4bool Inside(const G4ThreeVector &pVec) const

References G4VoxelLimits::ClipToLimits(), and G4VoxelLimits::Inside().

Referenced by G4VSolid::ClipPolygon().

◆ Clone()

G4VSolid * G4Trd::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 686 of file G4Trd.cc.

687{
688 return new G4Trd(*this);
689}
G4Trd(const G4String &pName, G4double pdx1, G4double pdx2, G4double pdy1, G4double pdy2, G4double pdz)
Definition: G4Trd.cc:54

References G4Trd().

◆ ComputeDimensions()

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

Reimplemented from G4VSolid.

Definition at line 237 of file G4Trd.cc.

240{
241 p->ComputeDimensions(*this,n,pRep);
242}
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const

References G4VPVParameterisation::ComputeDimensions(), and CLHEP::detail::n.

◆ CreatePolyhedron()

G4Polyhedron * G4Trd::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 804 of file G4Trd.cc.

805{
806 return new G4PolyhedronTrd2 (fDx1, fDx2, fDy1, fDy2, fDz);
807}

References fDx1, fDx2, fDy1, fDy2, and fDz.

◆ DescribeYourselfTo()

void G4Trd::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 799 of file G4Trd.cc.

800{
801 scene.AddSolid (*this);
802}
virtual void AddSolid(const G4Box &)=0

References G4VGraphicsScene::AddSolid().

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 545 of file G4Trd.cc.

546{
547 G4double dx = fPlanes[3].a*std::abs(p.x())+fPlanes[3].c*p.z()+fPlanes[3].d;
548 G4double dy = fPlanes[1].b*std::abs(p.y())+fPlanes[1].c*p.z()+fPlanes[1].d;
549 G4double dxy = std::max(dx,dy);
550
551 G4double dz = std::abs(p.z())-fDz;
552 G4double dist = std::max(dz,dxy);
553
554 return (dist > 0) ? dist : 0.;
555}

References fDz, fPlanes, G4INCL::Math::max(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ DistanceToIn() [2/2]

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

Implements G4VSolid.

Definition at line 452 of file G4Trd.cc.

454{
455 // Z intersections
456 //
457 if ((std::abs(p.z()) - fDz) >= -halfCarTolerance && p.z()*v.z() >= 0)
458 return kInfinity;
459 G4double invz = (-v.z() == 0) ? DBL_MAX : -1./v.z();
460 G4double dz = (invz < 0) ? fDz : -fDz;
461 G4double tzmin = (p.z() + dz)*invz;
462 G4double tzmax = (p.z() - dz)*invz;
463
464 // Y intersections
465 //
466 G4double tmin0 = tzmin, tmax0 = tzmax;
467 G4double ya = fPlanes[0].b*v.y(), yb = fPlanes[0].c*v.z();
468 G4double yc = fPlanes[0].b*p.y(), yd = fPlanes[0].c*p.z()+fPlanes[0].d;
469 G4double cos0 = yb + ya;
470 G4double dis0 = yd + yc;
471 if (dis0 >= -halfCarTolerance)
472 {
473 if (cos0 >= 0) return kInfinity;
474 G4double tmp = -dis0/cos0;
475 if (tmin0 < tmp) tmin0 = tmp;
476 }
477 else if (cos0 > 0)
478 {
479 G4double tmp = -dis0/cos0;
480 if (tmax0 > tmp) tmax0 = tmp;
481 }
482
483 G4double tmin1 = tmin0, tmax1 = tmax0;
484 G4double cos1 = yb - ya;
485 G4double dis1 = yd - yc;
486 if (dis1 >= -halfCarTolerance)
487 {
488 if (cos1 >= 0) return kInfinity;
489 G4double tmp = -dis1/cos1;
490 if (tmin1 < tmp) tmin1 = tmp;
491 }
492 else if (cos1 > 0)
493 {
494 G4double tmp = -dis1/cos1;
495 if (tmax1 > tmp) tmax1 = tmp;
496 }
497
498 // X intersections
499 //
500 G4double tmin2 = tmin1, tmax2 = tmax1;
501 G4double xa = fPlanes[2].a*v.x(), xb = fPlanes[2].c*v.z();
502 G4double xc = fPlanes[2].a*p.x(), xd = fPlanes[2].c*p.z()+fPlanes[2].d;
503 G4double cos2 = xb + xa;
504 G4double dis2 = xd + xc;
505 if (dis2 >= -halfCarTolerance)
506 {
507 if (cos2 >= 0) return kInfinity;
508 G4double tmp = -dis2/cos2;
509 if (tmin2 < tmp) tmin2 = tmp;
510 }
511 else if (cos2 > 0)
512 {
513 G4double tmp = -dis2/cos2;
514 if (tmax2 > tmp) tmax2 = tmp;
515 }
516
517 G4double tmin3 = tmin2, tmax3 = tmax2;
518 G4double cos3 = xb - xa;
519 G4double dis3 = xd - xc;
520 if (dis3 >= -halfCarTolerance)
521 {
522 if (cos3 >= 0) return kInfinity;
523 G4double tmp = -dis3/cos3;
524 if (tmin3 < tmp) tmin3 = tmp;
525 }
526 else if (cos3 > 0)
527 {
528 G4double tmp = -dis3/cos3;
529 if (tmax3 > tmp) tmax3 = tmp;
530 }
531
532 // Find distance
533 //
534 G4double tmin = tmin3, tmax = tmax3;
535 if (tmax <= tmin + halfCarTolerance) return kInfinity; // touch or no hit
536 return (tmin < halfCarTolerance ) ? 0. : tmin;
537}

References DBL_MAX, fDz, fPlanes, halfCarTolerance, kInfinity, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ DistanceToOut() [1/2]

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

Implements G4VSolid.

Definition at line 645 of file G4Trd.cc.

646{
647#ifdef G4CSGDEBUG
648 if( Inside(p) == kOutside )
649 {
650 std::ostringstream message;
651 G4int oldprc = message.precision(16);
652 message << "Point p is outside (!?) of solid: " << GetName() << G4endl;
653 message << "Position:\n";
654 message << " p.x() = " << p.x()/mm << " mm\n";
655 message << " p.y() = " << p.y()/mm << " mm\n";
656 message << " p.z() = " << p.z()/mm << " mm";
657 G4cout.precision(oldprc);
658 G4Exception("G4Trd::DistanceToOut(p)", "GeomSolids1002",
659 JustWarning, message );
660 DumpInfo();
661 }
662#endif
663 G4double dx = fPlanes[3].a*std::abs(p.x())+fPlanes[3].c*p.z()+fPlanes[3].d;
664 G4double dy = fPlanes[1].b*std::abs(p.y())+fPlanes[1].c*p.z()+fPlanes[1].d;
665 G4double dxy = std::max(dx,dy);
666
667 G4double dz = std::abs(p.z())-fDz;
668 G4double dist = std::max(dz,dxy);
669
670 return (dist < 0) ? -dist : 0.;
671}
static constexpr double mm
Definition: G4SIunits.hh:95
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
EInside Inside(const G4ThreeVector &p) const
Definition: G4Trd.cc:330
@ kOutside
Definition: geomdefs.hh:68

References G4VSolid::DumpInfo(), fDz, fPlanes, G4cout, G4endl, G4Exception(), G4VSolid::GetName(), Inside(), JustWarning, kOutside, G4INCL::Math::max(), mm, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ DistanceToOut() [2/2]

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

Implements G4VSolid.

Definition at line 563 of file G4Trd.cc.

566{
567 // Z intersections
568 //
569 if ((std::abs(p.z()) - fDz) >= -halfCarTolerance && p.z()*v.z() > 0)
570 {
571 if (calcNorm)
572 {
573 *validNorm = true;
574 n->set(0, 0, (p.z() < 0) ? -1 : 1);
575 }
576 return 0;
577 }
578 G4double vz = v.z();
579 G4double tmax = (vz == 0) ? DBL_MAX : (std::copysign(fDz,vz) - p.z())/vz;
580 G4int iside = (vz < 0) ? -4 : -2; // little trick: (-4+3)=-1, (-2+3)=+1
581
582 // Y intersections
583 //
584 G4int i = 0;
585 for ( ; i<2; ++i)
586 {
587 G4double cosa = fPlanes[i].b*v.y() + fPlanes[i].c*v.z();
588 if (cosa > 0)
589 {
590 G4double dist = fPlanes[i].b*p.y()+fPlanes[i].c*p.z()+fPlanes[i].d;
591 if (dist >= -halfCarTolerance)
592 {
593 if (calcNorm)
594 {
595 *validNorm = true;
596 n->set(0, fPlanes[i].b, fPlanes[i].c);
597 }
598 return 0;
599 }
600 G4double tmp = -dist/cosa;
601 if (tmax > tmp) { tmax = tmp; iside = i; }
602 }
603 }
604
605 // X intersections
606 //
607 for ( ; i<4; ++i)
608 {
609 G4double cosa = fPlanes[i].a*v.x()+fPlanes[i].c*v.z();
610 if (cosa > 0)
611 {
612 G4double dist = fPlanes[i].a*p.x()+fPlanes[i].c*p.z()+fPlanes[i].d;
613 if (dist >= -halfCarTolerance)
614 {
615 if (calcNorm)
616 {
617 *validNorm = true;
618 n->set(fPlanes[i].a, fPlanes[i].b, fPlanes[i].c);
619 }
620 return 0;
621 }
622 G4double tmp = -dist/cosa;
623 if (tmax > tmp) { tmax = tmp; iside = i; }
624 }
625 }
626
627 // Set normal, if required, and return distance
628 //
629 if (calcNorm)
630 {
631 *validNorm = true;
632 if (iside < 0)
633 n->set(0, 0, iside + 3); // (-4+3)=-1, (-2+3)=+1
634 else
635 n->set(fPlanes[iside].a, fPlanes[iside].b, fPlanes[iside].c);
636 }
637 return tmax;
638}

References a, b, c, DBL_MAX, fDz, fPlanes, halfCarTolerance, CLHEP::detail::n, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ DumpInfo()

void G4VSolid::DumpInfo ( ) const
inlineinherited

Referenced by G4Cons::ApproxSurfaceNormal(), G4CutTubs::ApproxSurfaceNormal(), G4Sphere::ApproxSurfaceNormal(), G4Torus::ApproxSurfaceNormal(), G4Tubs::ApproxSurfaceNormal(), G4ReflectedSolid::BoundingLimits(), G4DisplacedSolid::BoundingLimits(), G4IntersectionSolid::BoundingLimits(), G4ScaledSolid::BoundingLimits(), G4SubtractionSolid::BoundingLimits(), G4UnionSolid::BoundingLimits(), G4Box::BoundingLimits(), G4Cons::BoundingLimits(), G4CutTubs::BoundingLimits(), G4Orb::BoundingLimits(), G4Para::BoundingLimits(), G4Sphere::BoundingLimits(), G4Torus::BoundingLimits(), G4Trap::BoundingLimits(), BoundingLimits(), G4Tubs::BoundingLimits(), G4EllipticalCone::BoundingLimits(), G4ExtrudedSolid::BoundingLimits(), G4GenericPolycone::BoundingLimits(), G4GenericTrap::BoundingLimits(), G4Hype::BoundingLimits(), G4Paraboloid::BoundingLimits(), G4Polycone::BoundingLimits(), G4Polyhedra::BoundingLimits(), G4TessellatedSolid::BoundingLimits(), G4TwistedTubs::BoundingLimits(), G4ParameterisationBoxX::ComputeDimensions(), G4ParameterisationBoxY::ComputeDimensions(), G4ParameterisationBoxZ::ComputeDimensions(), G4ParameterisationConsRho::ComputeDimensions(), G4ParameterisationConsPhi::ComputeDimensions(), G4ParameterisationConsZ::ComputeDimensions(), G4ParameterisationParaX::ComputeDimensions(), G4ParameterisationParaY::ComputeDimensions(), G4ParameterisationParaZ::ComputeDimensions(), G4ParameterisationPolyconeRho::ComputeDimensions(), G4ParameterisationPolyconePhi::ComputeDimensions(), G4ParameterisationPolyconeZ::ComputeDimensions(), G4ParameterisationPolyhedraRho::ComputeDimensions(), G4ParameterisationPolyhedraPhi::ComputeDimensions(), G4ParameterisationPolyhedraZ::ComputeDimensions(), G4ParameterisationTrdX::ComputeDimensions(), G4ParameterisationTrdY::ComputeDimensions(), G4ParameterisationTrdZ::ComputeDimensions(), G4ParameterisationTubsRho::ComputeDimensions(), G4ParameterisationTubsPhi::ComputeDimensions(), G4ParameterisationTubsZ::ComputeDimensions(), G4ReflectedSolid::ComputeDimensions(), G4DisplacedSolid::ComputeDimensions(), G4ScaledSolid::ComputeDimensions(), G4ParameterisedNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4DisplacedSolid::CreatePolyhedron(), G4ScaledSolid::CreatePolyhedron(), G4SubtractionSolid::DistanceToIn(), G4Box::DistanceToOut(), G4Orb::DistanceToOut(), G4Para::DistanceToOut(), G4Trap::DistanceToOut(), DistanceToOut(), G4Paraboloid::DistanceToOut(), G4VTwistedFaceted::DistanceToOut(), G4Cons::DistanceToOut(), G4CutTubs::DistanceToOut(), G4Sphere::DistanceToOut(), G4Torus::DistanceToOut(), G4Tubs::DistanceToOut(), G4Ellipsoid::DistanceToOut(), G4EllipticalCone::DistanceToOut(), G4EllipticalTube::DistanceToOut(), G4GenericTrap::DistanceToOut(), export_G4VSolid(), G4Polycone::G4Polycone(), G4Polyhedra::G4Polyhedra(), G4BooleanSolid::GetConstituentSolid(), G4NavigationLogger::PostComputeStepLog(), G4Box::SurfaceNormal(), G4Para::SurfaceNormal(), G4Trap::SurfaceNormal(), SurfaceNormal(), G4Ellipsoid::SurfaceNormal(), G4EllipticalCone::SurfaceNormal(), G4EllipticalTube::SurfaceNormal(), G4ExtrudedSolid::SurfaceNormal(), and G4Tet::SurfaceNormal().

◆ EstimateCubicVolume()

G4double G4VSolid::EstimateCubicVolume ( G4int  nStat,
G4double  epsilon 
) const
inherited

Definition at line 203 of file G4VSolid.cc.

204{
205 G4int iInside=0;
206 G4double px,py,pz,minX,maxX,minY,maxY,minZ,maxZ,volume,halfepsilon;
208 EInside in;
209
210 // values needed for CalculateExtent signature
211
212 G4VoxelLimits limit; // Unlimited
213 G4AffineTransform origin;
214
215 // min max extents of pSolid along X,Y,Z
216
217 CalculateExtent(kXAxis,limit,origin,minX,maxX);
218 CalculateExtent(kYAxis,limit,origin,minY,maxY);
219 CalculateExtent(kZAxis,limit,origin,minZ,maxZ);
220
221 // limits
222
223 if(nStat < 100) nStat = 100;
224 if(epsilon > 0.01) epsilon = 0.01;
225 halfepsilon = 0.5*epsilon;
226
227 for(auto i = 0; i < nStat; ++i )
228 {
229 px = minX-halfepsilon+(maxX-minX+epsilon)*G4QuickRand();
230 py = minY-halfepsilon+(maxY-minY+epsilon)*G4QuickRand();
231 pz = minZ-halfepsilon+(maxZ-minZ+epsilon)*G4QuickRand();
232 p = G4ThreeVector(px,py,pz);
233 in = Inside(p);
234 if(in != kOutside) ++iInside;
235 }
236 volume = (maxX-minX+epsilon)*(maxY-minY+epsilon)
237 * (maxZ-minZ+epsilon)*iInside/nStat;
238 return volume;
239}
G4double epsilon(G4double density, G4double temperature)
static const G4int maxZ
G4double G4QuickRand()
Definition: G4QuickRand.hh:34
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
EInside
Definition: geomdefs.hh:67

References G4VSolid::CalculateExtent(), epsilon(), G4QuickRand(), G4VSolid::Inside(), kOutside, kXAxis, kYAxis, kZAxis, and maxZ.

Referenced by G4VSolid::GetCubicVolume(), G4BooleanSolid::GetCubicVolume(), and G4VCSGfaceted::GetCubicVolume().

◆ EstimateSurfaceArea()

G4double G4VSolid::EstimateSurfaceArea ( G4int  nStat,
G4double  ell 
) const
inherited

Definition at line 265 of file G4VSolid.cc.

266{
267 static const G4double s2 = 1./std::sqrt(2.);
268 static const G4double s3 = 1./std::sqrt(3.);
269 static const G4ThreeVector directions[64] =
270 {
271 G4ThreeVector( 0, 0, 0), G4ThreeVector( -1, 0, 0), // ( , , ) ( -, , )
272 G4ThreeVector( 1, 0, 0), G4ThreeVector( -1, 0, 0), // ( +, , ) (-+, , )
273 G4ThreeVector( 0, -1, 0), G4ThreeVector(-s2,-s2, 0), // ( , -, ) ( -, -, )
274 G4ThreeVector( s2, -s2, 0), G4ThreeVector( 0, -1, 0), // ( +, -, ) (-+, -, )
275
276 G4ThreeVector( 0, 1, 0), G4ThreeVector( -s2, s2, 0), // ( , +, ) ( -, +, )
277 G4ThreeVector( s2, s2, 0), G4ThreeVector( 0, 1, 0), // ( +, +, ) (-+, +, )
278 G4ThreeVector( 0, -1, 0), G4ThreeVector( -1, 0, 0), // ( ,-+, ) ( -,-+, )
279 G4ThreeVector( 1, 0, 0), G4ThreeVector( -1, 0, 0), // ( +,-+, ) (-+,-+, )
280
281 G4ThreeVector( 0, 0, -1), G4ThreeVector(-s2, 0,-s2), // ( , , -) ( -, , -)
282 G4ThreeVector( s2, 0,-s2), G4ThreeVector( 0, 0, -1), // ( +, , -) (-+, , -)
283 G4ThreeVector( 0,-s2,-s2), G4ThreeVector(-s3,-s3,-s3), // ( , -, -) ( -, -, -)
284 G4ThreeVector( s3,-s3,-s3), G4ThreeVector( 0,-s2,-s2), // ( +, -, -) (-+, -, -)
285
286 G4ThreeVector( 0, s2,-s2), G4ThreeVector(-s3, s3,-s3), // ( , +, -) ( -, +, -)
287 G4ThreeVector( s3, s3,-s3), G4ThreeVector( 0, s2,-s2), // ( +, +, -) (-+, +, -)
288 G4ThreeVector( 0, 0, -1), G4ThreeVector(-s2, 0,-s2), // ( ,-+, -) ( -,-+, -)
289 G4ThreeVector( s2, 0,-s2), G4ThreeVector( 0, 0, -1), // ( +,-+, -) (-+,-+, -)
290
291 G4ThreeVector( 0, 0, 1), G4ThreeVector(-s2, 0, s2), // ( , , +) ( -, , +)
292 G4ThreeVector( s2, 0, s2), G4ThreeVector( 0, 0, 1), // ( +, , +) (-+, , +)
293 G4ThreeVector( 0,-s2, s2), G4ThreeVector(-s3,-s3, s3), // ( , -, +) ( -, -, +)
294 G4ThreeVector( s3,-s3, s3), G4ThreeVector( 0,-s2, s2), // ( +, -, +) (-+, -, +)
295
296 G4ThreeVector( 0, s2, s2), G4ThreeVector(-s3, s3, s3), // ( , +, +) ( -, +, +)
297 G4ThreeVector( s3, s3, s3), G4ThreeVector( 0, s2, s2), // ( +, +, +) (-+, +, +)
298 G4ThreeVector( 0, 0, 1), G4ThreeVector(-s2, 0, s2), // ( ,-+, +) ( -,-+, +)
299 G4ThreeVector( s2, 0, s2), G4ThreeVector( 0, 0, 1), // ( +,-+, +) (-+,-+, +)
300
301 G4ThreeVector( 0, 0, -1), G4ThreeVector( -1, 0, 0), // ( , ,-+) ( -, ,-+)
302 G4ThreeVector( 1, 0, 0), G4ThreeVector( -1, 0, 0), // ( +, ,-+) (-+, ,-+)
303 G4ThreeVector( 0, -1, 0), G4ThreeVector(-s2,-s2, 0), // ( , -,-+) ( -, -,-+)
304 G4ThreeVector( s2, -s2, 0), G4ThreeVector( 0, -1, 0), // ( +, -,-+) (-+, -,-+)
305
306 G4ThreeVector( 0, 1, 0), G4ThreeVector( -s2, s2, 0), // ( , +,-+) ( -, +,-+)
307 G4ThreeVector( s2, s2, 0), G4ThreeVector( 0, 1, 0), // ( +, +,-+) (-+, +,-+)
308 G4ThreeVector( 0, -1, 0), G4ThreeVector( -1, 0, 0), // ( ,-+,-+) ( -,-+,-+)
309 G4ThreeVector( 1, 0, 0), G4ThreeVector( -1, 0, 0), // ( +,-+,-+) (-+,-+,-+)
310 };
311
312 G4ThreeVector bmin, bmax;
313 BoundingLimits(bmin, bmax);
314
315 G4double dX = bmax.x() - bmin.x();
316 G4double dY = bmax.y() - bmin.y();
317 G4double dZ = bmax.z() - bmin.z();
318
319 // Define statistics and shell thickness
320 //
321 G4int npoints = (nstat < 1000) ? 1000 : nstat;
322 G4double coeff = 0.5 / std::cbrt(G4double(npoints));
323 G4double eps = (ell > 0) ? ell : coeff * std::min(std::min(dX, dY), dZ);
324 G4double del = 1.8 * eps; // shold be more than sqrt(3.)
325
326 G4double minX = bmin.x() - eps;
327 G4double minY = bmin.y() - eps;
328 G4double minZ = bmin.z() - eps;
329
330 G4double dd = 2. * eps;
331 dX += dd;
332 dY += dd;
333 dZ += dd;
334
335 // Calculate surface area
336 //
337 G4int icount = 0;
338 for(auto i = 0; i < npoints; ++i)
339 {
340 G4double px = minX + dX*G4QuickRand();
341 G4double py = minY + dY*G4QuickRand();
342 G4double pz = minZ + dZ*G4QuickRand();
343 G4ThreeVector p = G4ThreeVector(px, py, pz);
344 EInside in = Inside(p);
345 G4double dist = 0;
346 if (in == kInside)
347 {
348 if (DistanceToOut(p) >= eps) continue;
349 G4int icase = 0;
350 if (Inside(G4ThreeVector(px-del, py, pz)) != kInside) icase += 1;
351 if (Inside(G4ThreeVector(px+del, py, pz)) != kInside) icase += 2;
352 if (Inside(G4ThreeVector(px, py-del, pz)) != kInside) icase += 4;
353 if (Inside(G4ThreeVector(px, py+del, pz)) != kInside) icase += 8;
354 if (Inside(G4ThreeVector(px, py, pz-del)) != kInside) icase += 16;
355 if (Inside(G4ThreeVector(px, py, pz+del)) != kInside) icase += 32;
356 if (icase == 0) continue;
357 G4ThreeVector v = directions[icase];
358 dist = DistanceToOut(p, v);
359 G4ThreeVector n = SurfaceNormal(p + v*dist);
360 dist *= v.dot(n);
361 }
362 else if (in == kOutside)
363 {
364 if (DistanceToIn(p) >= eps) continue;
365 G4int icase = 0;
366 if (Inside(G4ThreeVector(px-del, py, pz)) != kOutside) icase += 1;
367 if (Inside(G4ThreeVector(px+del, py, pz)) != kOutside) icase += 2;
368 if (Inside(G4ThreeVector(px, py-del, pz)) != kOutside) icase += 4;
369 if (Inside(G4ThreeVector(px, py+del, pz)) != kOutside) icase += 8;
370 if (Inside(G4ThreeVector(px, py, pz-del)) != kOutside) icase += 16;
371 if (Inside(G4ThreeVector(px, py, pz+del)) != kOutside) icase += 32;
372 if (icase == 0) continue;
373 G4ThreeVector v = directions[icase];
374 dist = DistanceToIn(p, v);
375 if (dist == kInfinity) continue;
376 G4ThreeVector n = SurfaceNormal(p + v*dist);
377 dist *= -(v.dot(n));
378 }
379 if (dist < eps) ++icount;
380 }
381 return dX*dY*dZ*icount/npoints/dd;
382}
static const G4double eps
double dot(const Hep3Vector &) const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4VSolid.cc:665
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
@ kInside
Definition: geomdefs.hh:70
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References G4VSolid::BoundingLimits(), G4VSolid::DistanceToIn(), G4VSolid::DistanceToOut(), CLHEP::Hep3Vector::dot(), eps, G4QuickRand(), G4VSolid::Inside(), kInfinity, kInside, kOutside, G4INCL::Math::min(), CLHEP::detail::n, G4VSolid::SurfaceNormal(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by G4VSolid::GetSurfaceArea(), G4MultiUnion::GetSurfaceArea(), and G4VCSGfaceted::GetSurfaceArea().

◆ GetConstituentSolid() [1/2]

G4VSolid * G4VSolid::GetConstituentSolid ( G4int  no)
virtualinherited

Reimplemented in G4BooleanSolid.

Definition at line 170 of file G4VSolid.cc.

171{ return nullptr; }

◆ GetConstituentSolid() [2/2]

const G4VSolid * G4VSolid::GetConstituentSolid ( G4int  no) const
virtualinherited

Reimplemented in G4BooleanSolid.

Definition at line 167 of file G4VSolid.cc.

168{ return nullptr; }

Referenced by G4BooleanSolid::StackPolyhedron().

◆ GetCubicVolume()

G4double G4Trd::GetCubicVolume ( )
virtual

Reimplemented from G4VSolid.

Definition at line 208 of file G4Trd.cc.

209{
210 if (fCubicVolume == 0.)
211 {
212 fCubicVolume = 2*fDz*( (fDx1+fDx2)*(fDy1+fDy2) +
213 (fDx2-fDx1)*(fDy2-fDy1)/3 );
214 }
215 return fCubicVolume;
216}
G4double fCubicVolume
Definition: G4CSGSolid.hh:70

References G4CSGSolid::fCubicVolume, fDx1, fDx2, fDy1, fDy2, and fDz.

◆ GetDisplacedSolidPtr() [1/2]

G4DisplacedSolid * G4VSolid::GetDisplacedSolidPtr ( )
virtualinherited

Reimplemented in G4DisplacedSolid.

Definition at line 176 of file G4VSolid.cc.

177{ return nullptr; }

◆ GetDisplacedSolidPtr() [2/2]

const G4DisplacedSolid * G4VSolid::GetDisplacedSolidPtr ( ) const
virtualinherited

Reimplemented in G4DisplacedSolid.

Definition at line 173 of file G4VSolid.cc.

174{ return nullptr; }

◆ GetEntityType()

G4GeometryType G4Trd::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 677 of file G4Trd.cc.

678{
679 return G4String("G4Trd");
680}

◆ GetExtent()

G4VisExtent G4VSolid::GetExtent ( ) const
virtualinherited

Reimplemented in G4Box, G4Orb, G4Sphere, G4Ellipsoid, G4EllipticalCone, G4EllipticalTube, G4GenericTrap, G4Hype, G4TessellatedSolid, G4Tet, G4TwistedTubs, G4VCSGfaceted, and G4VTwistedFaceted.

Definition at line 682 of file G4VSolid.cc.

683{
684 G4VisExtent extent;
685 G4VoxelLimits voxelLimits; // Defaults to "infinite" limits.
686 G4AffineTransform affineTransform;
687 G4double vmin, vmax;
688 CalculateExtent(kXAxis,voxelLimits,affineTransform,vmin,vmax);
689 extent.SetXmin (vmin);
690 extent.SetXmax (vmax);
691 CalculateExtent(kYAxis,voxelLimits,affineTransform,vmin,vmax);
692 extent.SetYmin (vmin);
693 extent.SetYmax (vmax);
694 CalculateExtent(kZAxis,voxelLimits,affineTransform,vmin,vmax);
695 extent.SetZmin (vmin);
696 extent.SetZmax (vmax);
697 return extent;
698}
void SetYmin(G4double ymin)
Definition: G4VisExtent.hh:114
void SetYmax(G4double ymax)
Definition: G4VisExtent.hh:116
void SetXmax(G4double xmax)
Definition: G4VisExtent.hh:112
void SetXmin(G4double xmin)
Definition: G4VisExtent.hh:110
void SetZmax(G4double zmax)
Definition: G4VisExtent.hh:120
void SetZmin(G4double zmin)
Definition: G4VisExtent.hh:118

References G4VSolid::CalculateExtent(), kXAxis, kYAxis, kZAxis, G4VisExtent::SetXmax(), G4VisExtent::SetXmin(), G4VisExtent::SetYmax(), G4VisExtent::SetYmin(), G4VisExtent::SetZmax(), and G4VisExtent::SetZmin().

Referenced by G4tgbVolume::BuildSolidForDivision(), G4BoundingExtentScene::ProcessVolume(), G4BoundingSphereScene::ProcessVolume(), and G4VisCommandsTouchable::SetNewValue().

◆ GetName()

G4String G4VSolid::GetName ( ) const
inlineinherited

Referenced by G4GMocrenFileSceneHandler::AddDetector(), G4HepRepFileSceneHandler::AddHepRepInstance(), G4GMocrenFileSceneHandler::AddPrimitive(), G4HepRepFileSceneHandler::AddSolid(), G4GMocrenFileSceneHandler::AddSolid(), G4VtkSceneHandler::AddSolid(), G4GDMLWriteSolids::AddSolid(), G4NavigationLogger::AlongComputeStepLog(), G4GDMLWriteSolids::BooleanWrite(), G4ReflectedSolid::BoundingLimits(), G4DisplacedSolid::BoundingLimits(), G4IntersectionSolid::BoundingLimits(), G4ScaledSolid::BoundingLimits(), G4SubtractionSolid::BoundingLimits(), G4UnionSolid::BoundingLimits(), G4Box::BoundingLimits(), G4Cons::BoundingLimits(), G4CutTubs::BoundingLimits(), G4Orb::BoundingLimits(), G4Para::BoundingLimits(), G4Sphere::BoundingLimits(), G4Torus::BoundingLimits(), G4Trap::BoundingLimits(), BoundingLimits(), G4Tubs::BoundingLimits(), G4EllipticalCone::BoundingLimits(), G4ExtrudedSolid::BoundingLimits(), G4GenericPolycone::BoundingLimits(), G4GenericTrap::BoundingLimits(), G4Hype::BoundingLimits(), G4Paraboloid::BoundingLimits(), G4Polycone::BoundingLimits(), G4Polyhedra::BoundingLimits(), G4TessellatedSolid::BoundingLimits(), G4TwistedTubs::BoundingLimits(), G4GDMLWriteSolids::BoxWrite(), G4ExtrudedSolid::CalculateExtent(), G4GenericPolycone::CalculateExtent(), G4Polycone::CalculateExtent(), G4Polyhedra::CalculateExtent(), G4NavigationLogger::CheckDaughterEntryPoint(), G4VDivisionParameterisation::CheckNDivAndWidth(), G4VDivisionParameterisation::CheckOffset(), G4GenericTrap::CheckOrder(), G4Para::CheckParameters(), G4Trap::CheckParameters(), CheckParameters(), G4Ellipsoid::CheckParameters(), G4EllipticalTube::CheckParameters(), G4ParameterisationPolyconeRho::CheckParametersValidity(), G4ParameterisationPolyconeZ::CheckParametersValidity(), G4ParameterisationPolyhedraRho::CheckParametersValidity(), G4ParameterisationPolyhedraPhi::CheckParametersValidity(), G4ParameterisationPolyhedraZ::CheckParametersValidity(), G4PhantomParameterisation::CheckVoxelsFillContainer(), G4GenericTrap::ComputeIsTwisted(), G4VoxelNavigation::ComputeSafety(), G4VoxelSafety::ComputeSafety(), G4NavigationLogger::ComputeSafetyLog(), G4ParameterisedNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4GDMLWriteSolids::ConeWrite(), G4Polyhedra::Create(), G4GenericPolycone::Create(), G4Polycone::Create(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G4ReflectedSolid::CreatePolyhedron(), G4ReflectionFactory::CreateReflectedLV(), G4GenericTrap::CreateTessellatedSolid(), G4GDMLWriteSolids::CutTubeWrite(), G4SolidStore::DeRegister(), G4PhysicalVolumeModel::DescribeSolid(), G4SubtractionSolid::DistanceToIn(), G4Paraboloid::DistanceToIn(), G4TessellatedSolid::DistanceToIn(), G4Box::DistanceToOut(), G4Orb::DistanceToOut(), G4Para::DistanceToOut(), G4Trap::DistanceToOut(), DistanceToOut(), G4EllipticalCone::DistanceToOut(), G4TessellatedSolid::DistanceToOut(), G4Ellipsoid::DistanceToOut(), G4EllipticalTube::DistanceToOut(), G4tgbGeometryDumper::DumpMultiUnionVolume(), G4tgbGeometryDumper::DumpScaledVolume(), G4tgbGeometryDumper::DumpSolid(), G4GDMLWriteSolids::ElconeWrite(), G4GDMLWriteSolids::EllipsoidWrite(), G4GDMLWriteSolids::EltubeWrite(), G4PVDivision::ErrorInAxis(), G4ReplicatedSlice::ErrorInAxis(), export_G4VSolid(), G4Box::G4Box(), G4Cons::G4Cons(), G4CutTubs::G4CutTubs(), G4EllipticalCone::G4EllipticalCone(), G4Hype::G4Hype(), G4Para::G4Para(), G4Paraboloid::G4Paraboloid(), G4Polycone::G4Polycone(), G4Polyhedra::G4Polyhedra(), G4Sphere::G4Sphere(), G4Tet::G4Tet(), G4Trap::G4Trap(), G4Tubs::G4Tubs(), G4VParameterisationCons::G4VParameterisationCons(), G4VParameterisationPara::G4VParameterisationPara(), G4VParameterisationPolycone::G4VParameterisationPolycone(), G4VParameterisationPolyhedra::G4VParameterisationPolyhedra(), G4VParameterisationTrd::G4VParameterisationTrd(), G4VTwistedFaceted::G4VTwistedFaceted(), G4GDMLWriteSolids::GenericPolyconeWrite(), G4GDMLWriteSolids::GenTrapWrite(), G4Navigator::GetGlobalExitNormal(), G4Navigator::GetLocalExitNormal(), G4ITNavigator1::GetLocalExitNormal(), G4ITNavigator2::GetLocalExitNormal(), G4BooleanSolid::GetPointOnSurface(), G4PhantomParameterisation::GetReplicaNo(), G4GDMLWriteSolids::HypeWrite(), G4TessellatedSolid::InsideNoVoxels(), G4TessellatedSolid::InsideVoxels(), G4ITNavigator1::LocateGlobalPointAndSetup(), G4ITNavigator2::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointAndSetup(), G4GenericTrap::MakeDownFacet(), G4Trap::MakePlanes(), G4GenericTrap::MakeUpFacet(), G4GDMLWriteSolids::MultiUnionWrite(), G4GDMLWriteSolids::OrbWrite(), G4GDMLWriteSolids::ParaboloidWrite(), G4GDMLWriteParamvol::ParametersWrite(), G4GDMLWriteSolids::ParaWrite(), G4GDMLWriteSolids::PolyconeWrite(), G4GDMLWriteSolids::PolyhedraWrite(), G4NavigationLogger::PostComputeStepLog(), G4NavigationLogger::PreComputeStepLog(), G4NavigationLogger::PrintDaughterLog(), G4PseudoScene::ProcessVolume(), G4SolidStore::Register(), G4tgbVolumeMgr::RegisterMe(), G4NavigationLogger::ReportOutsideMother(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4VSceneHandler::RequestPrimitives(), G4GenericPolycone::Reset(), G4Polyhedra::Reset(), G4VoxelSafety::SafetyForVoxelNode(), G4GDMLWriteSolids::ScaledWrite(), G4Torus::SetAllParameters(), G4Tet::SetBoundingLimits(), G4Polycone::SetOriginalParameters(), G4Polyhedra::SetOriginalParameters(), G4TessellatedSolid::SetSolidClosed(), G4Tet::SetVertices(), G4Box::SetXHalfLength(), G4Box::SetYHalfLength(), G4Box::SetZHalfLength(), G4GDMLWriteSolids::SphereWrite(), G4BooleanSolid::StackPolyhedron(), G4ReflectedSolid::StreamInfo(), G4BooleanSolid::StreamInfo(), G4DisplacedSolid::StreamInfo(), G4MultiUnion::StreamInfo(), G4ScaledSolid::StreamInfo(), G4Box::StreamInfo(), G4Cons::StreamInfo(), G4CSGSolid::StreamInfo(), G4CutTubs::StreamInfo(), G4Orb::StreamInfo(), G4Para::StreamInfo(), G4Sphere::StreamInfo(), G4Torus::StreamInfo(), G4Trap::StreamInfo(), StreamInfo(), G4Tubs::StreamInfo(), G4Ellipsoid::StreamInfo(), G4EllipticalCone::StreamInfo(), G4EllipticalTube::StreamInfo(), G4ExtrudedSolid::StreamInfo(), G4GenericPolycone::StreamInfo(), G4GenericTrap::StreamInfo(), G4Hype::StreamInfo(), G4Paraboloid::StreamInfo(), G4Polycone::StreamInfo(), G4Polyhedra::StreamInfo(), G4TessellatedSolid::StreamInfo(), G4Tet::StreamInfo(), G4TwistedBox::StreamInfo(), G4TwistedTrap::StreamInfo(), G4TwistedTrd::StreamInfo(), G4TwistedTubs::StreamInfo(), G4VCSGfaceted::StreamInfo(), G4VTwistedFaceted::StreamInfo(), G4GDMLRead::StripNames(), SubstractSolids(), G4UnionSolid::SurfaceNormal(), G4Box::SurfaceNormal(), G4Para::SurfaceNormal(), G4Trap::SurfaceNormal(), SurfaceNormal(), G4Ellipsoid::SurfaceNormal(), G4EllipticalCone::SurfaceNormal(), G4EllipticalTube::SurfaceNormal(), G4ExtrudedSolid::SurfaceNormal(), G4Tet::SurfaceNormal(), G4GDMLWriteSolids::TessellatedWrite(), G4GDMLWriteSolids::TetWrite(), G4GDMLWriteSolids::TorusWrite(), G4GDMLWriteSolids::TrapWrite(), G4GDMLWriteStructure::TraverseVolumeTree(), G4GDMLWriteSolids::TrdWrite(), G4GDMLWriteSolids::TubeWrite(), G4GDMLWriteSolids::TwistedboxWrite(), G4GDMLWriteSolids::TwistedtrapWrite(), G4GDMLWriteSolids::TwistedtrdWrite(), G4GDMLWriteSolids::TwistedtubsWrite(), G4PhysicalVolumeModel::VisitGeometryAndGetVisReps(), and G4GDMLWriteSolids::XtruWrite().

◆ GetPointOnSurface()

G4ThreeVector G4Trd::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 718 of file G4Trd.cc.

719{
720 // Set areas
721 //
722 G4double sxz = (fDx1 + fDx2)*fHx;
723 G4double syz = (fDy1 + fDy2)*fHy;
724 G4double ssurf[6] = { 4.*fDx1*fDy1, sxz, sxz, syz, syz, 4.*fDx2*fDy2 };
725 ssurf[1] += ssurf[0];
726 ssurf[2] += ssurf[1];
727 ssurf[3] += ssurf[2];
728 ssurf[4] += ssurf[3];
729 ssurf[5] += ssurf[4];
730
731 // Select face
732 //
733 G4double select = ssurf[5]*G4QuickRand();
734 G4int k = 5;
735 k -= (select <= ssurf[4]);
736 k -= (select <= ssurf[3]);
737 k -= (select <= ssurf[2]);
738 k -= (select <= ssurf[1]);
739 k -= (select <= ssurf[0]);
740
741 // Generate point on selected surface
742 //
743 G4double u = G4QuickRand();
744 G4double v = G4QuickRand();
745 switch(k)
746 {
747 case 0: // base at -Z
748 {
749 return G4ThreeVector((2.*u - 1.)*fDx1, (2.*v - 1.)*fDy1, -fDz);
750 }
751 case 1: // X face at -Y
752 {
753 if (u + v > 1.) { u = 1. - u; v = 1. - v; }
754 G4ThreeVector p0(-fDx1,-fDy1,-fDz);
755 G4ThreeVector p1( fDx2,-fDy2, fDz);
756 return (select <= ssurf[0] + fDx1*fHx) ?
757 (1. - u - v)*p0 + u*p1 + v*G4ThreeVector( fDx1,-fDy1,-fDz) :
758 (1. - u - v)*p0 + u*p1 + v*G4ThreeVector(-fDx2,-fDy2, fDz);
759 }
760 case 2: // X face at +Y
761 {
762 if (u + v > 1.) { u = 1. - u; v = 1. - v; }
763 G4ThreeVector p0( fDx1, fDy1,-fDz);
764 G4ThreeVector p1(-fDx2, fDy2, fDz);
765 return (select <= ssurf[1] + fDx1*fHx) ?
766 (1. - u - v)*p0 + u*p1 + v*G4ThreeVector(-fDx1, fDy1,-fDz) :
767 (1. - u - v)*p0 + u*p1 + v*G4ThreeVector( fDx2, fDy2, fDz);
768 }
769 case 3: // Y face at -X
770 {
771 if (u + v > 1.) { u = 1. - u; v = 1. - v; }
772 G4ThreeVector p0(-fDx1, fDy1,-fDz);
773 G4ThreeVector p1(-fDx2,-fDy2, fDz);
774 return (select <= ssurf[2] + fDy1*fHy) ?
775 (1. - u - v)*p0 + u*p1 + v*G4ThreeVector(-fDx1,-fDy1,-fDz) :
776 (1. - u - v)*p0 + u*p1 + v*G4ThreeVector(-fDx2, fDy2, fDz);
777 }
778 case 4: // Y face at +X
779 {
780 if (u + v > 1.) { u = 1. - u; v = 1. - v; }
781 G4ThreeVector p0( fDx1,-fDy1,-fDz);
782 G4ThreeVector p1( fDx2, fDy2, fDz);
783 return (select <= ssurf[3] + fDy1*fHy) ?
784 (1. - u - v)*p0 + u*p1 + v*G4ThreeVector( fDx1, fDy1,-fDz) :
785 (1. - u - v)*p0 + u*p1 + v*G4ThreeVector( fDx2,-fDy2, fDz);
786 }
787 case 5: // base at +Z
788 {
789 return G4ThreeVector((2.*u - 1.)*fDx2, (2.*v - 1.)*fDy2, fDz);
790 }
791 }
792 return G4ThreeVector(0., 0., 0.);
793}

References fDx1, fDx2, fDy1, fDy2, fDz, fHx, fHy, and G4QuickRand().

◆ GetPolyhedron()

G4Polyhedron * G4CSGSolid::GetPolyhedron ( ) const
virtualinherited

Reimplemented from G4VSolid.

Definition at line 129 of file G4CSGSolid.cc.

130{
131 if (fpPolyhedron == nullptr ||
135 {
137 delete fpPolyhedron;
139 fRebuildPolyhedron = false;
140 l.unlock();
141 }
142 return fpPolyhedron;
143}
G4Polyhedron * fpPolyhedron
Definition: G4CSGSolid.hh:73
G4bool fRebuildPolyhedron
Definition: G4CSGSolid.hh:72
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
virtual G4Polyhedron * CreatePolyhedron() const
Definition: G4VSolid.cc:700
static G4int GetNumberOfRotationSteps()

References G4VSolid::CreatePolyhedron(), G4CSGSolid::fpPolyhedron, G4CSGSolid::fRebuildPolyhedron, HepPolyhedron::GetNumberOfRotationSteps(), G4Polyhedron::GetNumberOfRotationStepsAtTimeOfCreation(), anonymous_namespace{G4CSGSolid.cc}::polyhedronMutex, and G4TemplateAutoLock< _Mutex_t >::unlock().

Referenced by G4ScoringBox::Draw(), G4ScoringCylinder::Draw(), G4ScoringBox::DrawColumn(), and G4ScoringCylinder::DrawColumn().

◆ GetRadiusInRing()

G4double G4CSGSolid::GetRadiusInRing ( G4double  rmin,
G4double  rmax 
) const
protectedinherited

Definition at line 109 of file G4CSGSolid.cc.

110{
111 G4double k = G4QuickRand();
112 return (rmin <= 0) ? rmax*std::sqrt(k)
113 : std::sqrt(k*rmax*rmax + (1. - k)*rmin*rmin);
114}

References G4QuickRand().

Referenced by G4Cons::GetPointOnSurface(), and G4Torus::GetPointOnSurface().

◆ GetSurfaceArea()

G4double G4Trd::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 222 of file G4Trd.cc.

223{
224 if (fSurfaceArea == 0.)
225 {
227 4*(fDx1*fDy1 + fDx2*fDy2) + 2*(fDx1+fDx2)*fHx + 2*(fDy1+fDy2)*fHy;
228 }
229 return fSurfaceArea;
230}
G4double fSurfaceArea
Definition: G4CSGSolid.hh:71

References fDx1, fDx2, fDy1, fDy2, fHx, fHy, and G4CSGSolid::fSurfaceArea.

◆ GetTolerance()

G4double G4VSolid::GetTolerance ( ) const
inlineinherited

◆ GetXHalfLength1()

G4double G4Trd::GetXHalfLength1 ( ) const
inline

◆ GetXHalfLength2()

G4double G4Trd::GetXHalfLength2 ( ) const
inline

◆ GetYHalfLength1()

G4double G4Trd::GetYHalfLength1 ( ) const
inline

◆ GetYHalfLength2()

G4double G4Trd::GetYHalfLength2 ( ) const
inline

◆ GetZHalfLength()

G4double G4Trd::GetZHalfLength ( ) const
inline

◆ Inside()

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

Implements G4VSolid.

Definition at line 330 of file G4Trd.cc.

331{
332 G4double dx = fPlanes[3].a*std::abs(p.x())+fPlanes[3].c*p.z()+fPlanes[3].d;
333 G4double dy = fPlanes[1].b*std::abs(p.y())+fPlanes[1].c*p.z()+fPlanes[1].d;
334 G4double dxy = std::max(dx,dy);
335
336 G4double dz = std::abs(p.z())-fDz;
337 G4double dist = std::max(dz,dxy);
338
339 return (dist > halfCarTolerance) ? kOutside :
340 ((dist > -halfCarTolerance) ? kSurface : kInside);
341}
@ kSurface
Definition: geomdefs.hh:69

References fDz, fPlanes, halfCarTolerance, kInside, kOutside, kSurface, G4INCL::Math::max(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by DistanceToOut().

◆ MakePlanes()

void G4Trd::MakePlanes ( )
private

Definition at line 171 of file G4Trd.cc.

172{
173 G4double dx = fDx1 - fDx2;
174 G4double dy = fDy1 - fDy2;
175 G4double dz = 2*fDz;
176 fHx = std::sqrt(dy*dy + dz*dz);
177 fHy = std::sqrt(dx*dx + dz*dz);
178
179 // Set X planes at -Y & +Y
180 //
181 fPlanes[0].a = 0.;
182 fPlanes[0].b = -dz/fHx;
183 fPlanes[0].c = dy/fHx;
184 fPlanes[0].d = fPlanes[0].b*fDy1 + fPlanes[0].c*fDz;
185
186 fPlanes[1].a = fPlanes[0].a;
187 fPlanes[1].b = -fPlanes[0].b;
188 fPlanes[1].c = fPlanes[0].c;
189 fPlanes[1].d = fPlanes[0].d;
190
191 // Set Y planes at -X & +X
192 //
193 fPlanes[2].a = -dz/fHy;
194 fPlanes[2].b = 0.;
195 fPlanes[2].c = dx/fHy;
196 fPlanes[2].d = fPlanes[2].a*fDx1 + fPlanes[2].c*fDz;
197
198 fPlanes[3].a = -fPlanes[2].a;
199 fPlanes[3].b = fPlanes[2].b;
200 fPlanes[3].c = fPlanes[2].c;
201 fPlanes[3].d = fPlanes[2].d;
202}

References fDx1, fDx2, fDy1, fDy2, fDz, fHx, fHy, and fPlanes.

Referenced by G4Trd(), and SetAllParameters().

◆ operator=()

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

Definition at line 102 of file G4Trd.cc.

103{
104 // Check assignment to self
105 //
106 if (this == &rhs) { return *this; }
107
108 // Copy base class data
109 //
111
112 // Copy data
113 //
115 fDx1 = rhs.fDx1; fDx2 = rhs.fDx2;
116 fDy1 = rhs.fDy1; fDy2 = rhs.fDy2;
117 fDz = rhs.fDz;
118 fHx = rhs.fHx; fHy = rhs.fHy;
119 for (G4int i=0; i<4; ++i) { fPlanes[i] = rhs.fPlanes[i]; }
120
121 return *this;
122}
G4CSGSolid & operator=(const G4CSGSolid &rhs)
Definition: G4CSGSolid.cc:89

References fDx1, fDx2, fDy1, fDy2, fDz, fHx, fHy, fPlanes, halfCarTolerance, and G4CSGSolid::operator=().

◆ operator==()

G4bool G4VSolid::operator== ( const G4VSolid s) const
inlineinherited

◆ SetAllParameters()

void G4Trd::SetAllParameters ( G4double  pdx1,
G4double  pdx2,
G4double  pdy1,
G4double  pdy2,
G4double  pdz 
)

Definition at line 128 of file G4Trd.cc.

130{
131 // Reset data of the base class
132 fCubicVolume = 0.;
133 fSurfaceArea = 0.;
134 fRebuildPolyhedron = true;
135
136 // Set parameters
137 fDx1 = pdx1; fDx2 = pdx2;
138 fDy1 = pdy1; fDy2 = pdy2;
139 fDz = pdz;
140
142 MakePlanes();
143}

References CheckParameters(), G4CSGSolid::fCubicVolume, fDx1, fDx2, fDy1, fDy2, fDz, G4CSGSolid::fRebuildPolyhedron, G4CSGSolid::fSurfaceArea, and MakePlanes().

Referenced by G4ParameterisationTrdX::ComputeDimensions(), G4ParameterisationTrdY::ComputeDimensions(), and G4ParameterisationTrdZ::ComputeDimensions().

◆ SetName()

void G4VSolid::SetName ( const G4String name)
inherited

Definition at line 127 of file G4VSolid.cc.

128{
131}
void SetMapValid(G4bool val)
Definition: G4SolidStore.hh:76
static G4SolidStore * GetInstance()
G4String fshapeName
Definition: G4VSolid.hh:312
const char * name(G4int ptype)

References G4VSolid::fshapeName, G4SolidStore::GetInstance(), G4InuclParticleNames::name(), and G4SolidStore::SetMapValid().

Referenced by export_G4VSolid(), G4MultiUnion::G4MultiUnion(), and G4GDMLRead::StripNames().

◆ SetXHalfLength1()

void G4Trd::SetXHalfLength1 ( G4double  val)
inline

◆ SetXHalfLength2()

void G4Trd::SetXHalfLength2 ( G4double  val)
inline

◆ SetYHalfLength1()

void G4Trd::SetYHalfLength1 ( G4double  val)
inline

◆ SetYHalfLength2()

void G4Trd::SetYHalfLength2 ( G4double  val)
inline

◆ SetZHalfLength()

void G4Trd::SetZHalfLength ( G4double  val)
inline

◆ StreamInfo()

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

Reimplemented from G4CSGSolid.

Definition at line 695 of file G4Trd.cc.

696{
697 G4int oldprc = os.precision(16);
698 os << "-----------------------------------------------------------\n"
699 << " *** Dump for solid - " << GetName() << " ***\n"
700 << " ===================================================\n"
701 << " Solid type: G4Trd\n"
702 << " Parameters: \n"
703 << " half length X, surface -dZ: " << fDx1/mm << " mm \n"
704 << " half length X, surface +dZ: " << fDx2/mm << " mm \n"
705 << " half length Y, surface -dZ: " << fDy1/mm << " mm \n"
706 << " half length Y, surface +dZ: " << fDy2/mm << " mm \n"
707 << " half length Z : " << fDz/mm << " mm \n"
708 << "-----------------------------------------------------------\n";
709 os.precision(oldprc);
710
711 return os;
712}

References fDx1, fDx2, fDy1, fDy2, fDz, G4VSolid::GetName(), and mm.

◆ SurfaceNormal()

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

Implements G4VSolid.

Definition at line 347 of file G4Trd.cc.

348{
349 G4int nsurf = 0; // number of surfaces where p is placed
350
351 // Check Z faces
352 //
353 G4double nz = 0;
354 G4double dz = std::abs(p.z()) - fDz;
355 if (std::abs(dz) <= halfCarTolerance)
356 {
357 nz = (p.z() < 0) ? -1 : 1;
358 ++nsurf;
359 }
360
361 // Check Y faces
362 //
363 G4double ny = 0;
364 G4double dy1 = fPlanes[0].b*p.y();
365 G4double dy2 = fPlanes[0].c*p.z() + fPlanes[0].d;
366 if (std::abs(dy2 + dy1) <= halfCarTolerance)
367 {
368 ny += fPlanes[0].b;
369 nz += fPlanes[0].c;
370 ++nsurf;
371 }
372 if (std::abs(dy2 - dy1) <= halfCarTolerance)
373 {
374 ny += fPlanes[1].b;
375 nz += fPlanes[1].c;
376 ++nsurf;
377 }
378
379 // Check X faces
380 //
381 G4double nx = 0;
382 G4double dx1 = fPlanes[2].a*p.x();
383 G4double dx2 = fPlanes[2].c*p.z() + fPlanes[2].d;
384 if (std::abs(dx2 + dx1) <= halfCarTolerance)
385 {
386 nx += fPlanes[2].a;
387 nz += fPlanes[2].c;
388 ++nsurf;
389 }
390 if (std::abs(dx2 - dx1) <= halfCarTolerance)
391 {
392 nx += fPlanes[3].a;
393 nz += fPlanes[3].c;
394 ++nsurf;
395 }
396
397 // Return normal
398 //
399 if (nsurf == 1) return G4ThreeVector(nx,ny,nz);
400 else if (nsurf != 0) return G4ThreeVector(nx,ny,nz).unit(); // edge or corner
401 else
402 {
403 // Point is not on the surface
404 //
405#ifdef G4CSGDEBUG
406 std::ostringstream message;
407 G4int oldprc = message.precision(16);
408 message << "Point p is not on surface (!?) of solid: "
409 << GetName() << G4endl;
410 message << "Position:\n";
411 message << " p.x() = " << p.x()/mm << " mm\n";
412 message << " p.y() = " << p.y()/mm << " mm\n";
413 message << " p.z() = " << p.z()/mm << " mm";
414 G4cout.precision(oldprc) ;
415 G4Exception("G4Trd::SurfaceNormal(p)", "GeomSolids1002",
416 JustWarning, message );
417 DumpInfo();
418#endif
419 return ApproxSurfaceNormal(p);
420 }
421}
Hep3Vector unit() const
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Trd.cc:428

References ApproxSurfaceNormal(), G4VSolid::DumpInfo(), fDz, fPlanes, G4cout, G4endl, G4Exception(), G4VSolid::GetName(), halfCarTolerance, JustWarning, mm, CLHEP::Hep3Vector::unit(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Field Documentation

◆ a

G4double G4Trd::a

Definition at line 170 of file G4Trd.hh.

Referenced by ApproxSurfaceNormal(), and DistanceToOut().

◆ b

G4double G4Trd::b

Definition at line 170 of file G4Trd.hh.

Referenced by ApproxSurfaceNormal(), and DistanceToOut().

◆ c

G4double G4Trd::c

Definition at line 170 of file G4Trd.hh.

Referenced by ApproxSurfaceNormal(), and DistanceToOut().

◆ d

G4double G4Trd::d

Definition at line 170 of file G4Trd.hh.

Referenced by ApproxSurfaceNormal().

◆ fCubicVolume

G4double G4CSGSolid::fCubicVolume = 0.0
protectedinherited

◆ fDx1

G4double G4Trd::fDx1
private

◆ fDx2

G4double G4Trd::fDx2
private

◆ fDy1

G4double G4Trd::fDy1
private

◆ fDy2

G4double G4Trd::fDy2
private

◆ fDz

G4double G4Trd::fDz
private

◆ fHx

G4double G4Trd::fHx
private

Definition at line 169 of file G4Trd.hh.

Referenced by GetPointOnSurface(), GetSurfaceArea(), MakePlanes(), and operator=().

◆ fHy

G4double G4Trd::fHy
private

Definition at line 169 of file G4Trd.hh.

Referenced by GetPointOnSurface(), GetSurfaceArea(), MakePlanes(), and operator=().

◆ 

struct { ... } G4Trd::fPlanes[4]

◆ fpPolyhedron

G4Polyhedron* G4CSGSolid::fpPolyhedron = nullptr
mutableprotectedinherited

◆ fRebuildPolyhedron

G4bool G4CSGSolid::fRebuildPolyhedron = false
mutableprotectedinherited

◆ fshapeName

G4String G4VSolid::fshapeName
privateinherited

Definition at line 312 of file G4VSolid.hh.

Referenced by G4VSolid::operator=(), and G4VSolid::SetName().

◆ fSurfaceArea

G4double G4CSGSolid::fSurfaceArea = 0.0
protectedinherited

◆ halfCarTolerance

G4double G4Trd::halfCarTolerance
private

Definition at line 168 of file G4Trd.hh.

Referenced by DistanceToIn(), DistanceToOut(), Inside(), operator=(), and SurfaceNormal().

◆ kCarTolerance

G4double G4VSolid::kCarTolerance
protectedinherited

Definition at line 299 of file G4VSolid.hh.

Referenced by G4TessellatedSolid::AddFacet(), G4Polycone::CalculateExtent(), G4Polyhedra::CalculateExtent(), G4Tet::CheckDegeneracy(), G4Para::CheckParameters(), CheckParameters(), G4Ellipsoid::CheckParameters(), G4EllipticalTube::CheckParameters(), G4GenericTrap::ComputeIsTwisted(), G4Polyhedra::Create(), G4GenericPolycone::Create(), G4Polycone::Create(), G4CutTubs::CreatePolyhedron(), G4TessellatedSolid::CreateVertexList(), G4VCSGfaceted::DistanceTo(), G4Sphere::DistanceToIn(), G4Ellipsoid::DistanceToIn(), G4Hype::DistanceToIn(), G4Paraboloid::DistanceToIn(), G4VCSGfaceted::DistanceToIn(), G4TessellatedSolid::DistanceToInCore(), G4Cons::DistanceToOut(), G4CutTubs::DistanceToOut(), G4Sphere::DistanceToOut(), G4Torus::DistanceToOut(), G4Tubs::DistanceToOut(), G4GenericTrap::DistanceToOut(), G4Hype::DistanceToOut(), G4Paraboloid::DistanceToOut(), G4VCSGfaceted::DistanceToOut(), G4TessellatedSolid::DistanceToOutCandidates(), G4TessellatedSolid::DistanceToOutCore(), G4TessellatedSolid::DistanceToOutNoVoxels(), G4GenericTrap::DistToPlane(), G4GenericTrap::DistToTriangle(), G4Box::G4Box(), G4Cons::G4Cons(), G4CutTubs::G4CutTubs(), G4EllipticalCone::G4EllipticalCone(), G4ExtrudedSolid::G4ExtrudedSolid(), G4GenericTrap::G4GenericTrap(), G4Hype::G4Hype(), G4Para::G4Para(), G4Sphere::G4Sphere(), G4Tet::G4Tet(), G4Trap::G4Trap(), G4Tubs::G4Tubs(), G4UnionSolid::G4UnionSolid(), G4VSolid::G4VSolid(), G4VTwistedFaceted::G4VTwistedFaceted(), G4GenericPolycone::GetPointOnSurface(), G4Polycone::GetPointOnSurface(), G4UnionSolid::Init(), G4Orb::Initialize(), G4TessellatedSolid::Initialize(), G4SubtractionSolid::Inside(), G4Hype::Inside(), G4Paraboloid::Inside(), G4VCSGfaceted::Inside(), G4VTwistedFaceted::Inside(), G4TessellatedSolid::InsideNoVoxels(), G4GenericTrap::InsidePolygone(), G4TessellatedSolid::InsideVoxels(), G4CutTubs::IsCrossingCutPlanes(), G4GenericTrap::IsSegCrossingZ(), G4Trap::MakePlane(), G4GenericTrap::NormalToPlane(), G4VSolid::operator=(), G4TessellatedSolid::SafetyFromInside(), G4TessellatedSolid::SafetyFromOutside(), G4Torus::SetAllParameters(), G4Polycone::SetOriginalParameters(), G4Polyhedra::SetOriginalParameters(), G4Box::SetXHalfLength(), G4Box::SetYHalfLength(), G4Box::SetZHalfLength(), G4Torus::SurfaceNormal(), G4GenericTrap::SurfaceNormal(), and G4Paraboloid::SurfaceNormal().


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