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

#include <G4Trap.hh>

Inheritance diagram for G4Trap:
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
 
 G4Trap (__void__ &)
 
 G4Trap (const G4String &pName)
 
 G4Trap (const G4String &pName, const G4ThreeVector pt[8])
 
 G4Trap (const G4String &pName, G4double pDx, G4double pDy, G4double pDz, G4double pAlpha, G4double pTheta, G4double pPhi)
 
 G4Trap (const G4String &pName, G4double pDx1, G4double pDx2, G4double pDy1, G4double pDy2, G4double pDz)
 
 G4Trap (const G4String &pName, G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pAlp1, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlp2)
 
 G4Trap (const G4String &pName, G4double pZ, G4double pY, G4double pX, G4double pLTX)
 
 G4Trap (const G4Trap &rhs)
 
G4double GetAlpha1 () const
 
G4double GetAlpha2 () const
 
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
 
G4double GetPhi () const
 
G4ThreeVector GetPointOnSurface () const
 
virtual G4PolyhedronGetPolyhedron () const
 
TrapSidePlane GetSidePlane (G4int n) const
 
G4double GetSurfaceArea ()
 
G4ThreeVector GetSymAxis () const
 
G4double GetTanAlpha1 () const
 
G4double GetTanAlpha2 () const
 
G4double GetTheta () const
 
G4double GetTolerance () const
 
G4double GetXHalfLength1 () const
 
G4double GetXHalfLength2 () const
 
G4double GetXHalfLength3 () const
 
G4double GetXHalfLength4 () const
 
G4double GetYHalfLength1 () const
 
G4double GetYHalfLength2 () const
 
G4double GetZHalfLength () const
 
EInside Inside (const G4ThreeVector &p) const
 
G4Trapoperator= (const G4Trap &rhs)
 
G4bool operator== (const G4VSolid &s) const
 
void SetAllParameters (G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pAlp1, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlp2)
 
void SetName (const G4String &name)
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
virtual ~G4Trap ()
 

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
 
G4bool MakePlane (const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3, const G4ThreeVector &p4, TrapSidePlane &plane)
 
void MakePlanes ()
 
void MakePlanes (const G4ThreeVector pt[8])
 
void SetCachedValues ()
 

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 GetVertices (G4ThreeVector pt[8]) const
 

Private Attributes

G4double fAreas [6]
 
G4double fDx1
 
G4double fDx2
 
G4double fDx3
 
G4double fDx4
 
G4double fDy1
 
G4double fDy2
 
G4double fDz
 
TrapSidePlane fPlanes [4]
 
G4String fshapeName
 
G4double fTalpha1
 
G4double fTalpha2
 
G4int fTrapType
 
G4double fTthetaCphi
 
G4double fTthetaSphi
 
G4double halfCarTolerance
 

Detailed Description

Definition at line 109 of file G4Trap.hh.

Constructor & Destructor Documentation

◆ G4Trap() [1/8]

G4Trap::G4Trap ( const G4String pName,
G4double  pDz,
G4double  pTheta,
G4double  pPhi,
G4double  pDy1,
G4double  pDx1,
G4double  pDx2,
G4double  pAlp1,
G4double  pDy2,
G4double  pDx3,
G4double  pDx4,
G4double  pAlp2 
)

Definition at line 60 of file G4Trap.cc.

68{
69 fDz = pDz;
70 fTthetaCphi = std::tan(pTheta)*std::cos(pPhi);
71 fTthetaSphi = std::tan(pTheta)*std::sin(pPhi);
72
73 fDy1 = pDy1; fDx1 = pDx1; fDx2 = pDx2; fTalpha1 = std::tan(pAlp1);
74 fDy2 = pDy2; fDx3 = pDx3; fDx4 = pDx4; fTalpha2 = std::tan(pAlp2);
75
77 MakePlanes();
78}
G4CSGSolid(const G4String &pName)
Definition: G4CSGSolid.cc:49
G4double halfCarTolerance
Definition: G4Trap.hh:277
G4double fDy1
Definition: G4Trap.hh:279
G4double fDx3
Definition: G4Trap.hh:280
G4double fDy2
Definition: G4Trap.hh:280
G4double fTthetaCphi
Definition: G4Trap.hh:278
G4double fDz
Definition: G4Trap.hh:278
G4double fDx1
Definition: G4Trap.hh:279
void CheckParameters()
Definition: G4Trap.cc:315
G4double fDx4
Definition: G4Trap.hh:280
void MakePlanes()
Definition: G4Trap.cc:335
G4double fTalpha1
Definition: G4Trap.hh:279
G4double fTthetaSphi
Definition: G4Trap.hh:278
G4double fTalpha2
Definition: G4Trap.hh:280
G4double fDx2
Definition: G4Trap.hh:279
G4double kCarTolerance
Definition: G4VSolid.hh:299

References CheckParameters(), fDx1, fDx2, fDx3, fDx4, fDy1, fDy2, fDz, fTalpha1, fTalpha2, fTthetaCphi, fTthetaSphi, and MakePlanes().

Referenced by Clone().

◆ G4Trap() [2/8]

G4Trap::G4Trap ( const G4String pName,
const G4ThreeVector  pt[8] 
)

Definition at line 86 of file G4Trap.cc.

89{
90 // Start with check of centering - the center of gravity trap line
91 // should cross the origin of frame
92 //
93 if (!( pt[0].z() < 0
94 && pt[0].z() == pt[1].z()
95 && pt[0].z() == pt[2].z()
96 && pt[0].z() == pt[3].z()
97
98 && pt[4].z() > 0
99 && pt[4].z() == pt[5].z()
100 && pt[4].z() == pt[6].z()
101 && pt[4].z() == pt[7].z()
102
103 && std::fabs( pt[0].z() + pt[4].z() ) < kCarTolerance
104
105 && pt[0].y() == pt[1].y()
106 && pt[2].y() == pt[3].y()
107 && pt[4].y() == pt[5].y()
108 && pt[6].y() == pt[7].y()
109
110 && std::fabs(pt[0].y()+pt[2].y()+pt[4].y()+pt[6].y()) < kCarTolerance
111 && std::fabs(pt[0].x()+pt[1].x()+pt[4].x()+pt[5].x() +
112 pt[2].x()+pt[3].x()+pt[6].x()+pt[7].x()) < kCarTolerance ))
113 {
114 std::ostringstream message;
115 message << "Invalid vertice coordinates for Solid: " << GetName();
116 G4Exception("G4Trap::G4Trap()", "GeomSolids0002",
117 FatalException, message);
118 }
119
120 // Set parameters
121 //
122 fDz = (pt[7]).z();
123
124 fDy1 = ((pt[2]).y()-(pt[1]).y())*0.5;
125 fDx1 = ((pt[1]).x()-(pt[0]).x())*0.5;
126 fDx2 = ((pt[3]).x()-(pt[2]).x())*0.5;
127 fTalpha1 = ((pt[2]).x()+(pt[3]).x()-(pt[1]).x()-(pt[0]).x())*0.25/fDy1;
128
129 fDy2 = ((pt[6]).y()-(pt[5]).y())*0.5;
130 fDx3 = ((pt[5]).x()-(pt[4]).x())*0.5;
131 fDx4 = ((pt[7]).x()-(pt[6]).x())*0.5;
132 fTalpha2 = ((pt[6]).x()+(pt[7]).x()-(pt[5]).x()-(pt[4]).x())*0.25/fDy2;
133
134 fTthetaCphi = ((pt[4]).x()+fDy2*fTalpha2+fDx3)/fDz;
135 fTthetaSphi = ((pt[4]).y()+fDy2)/fDz;
136
138 MakePlanes(pt);
139}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4String GetName() const

References CheckParameters(), FatalException, fDx1, fDx2, fDx3, fDx4, fDy1, fDy2, fDz, fTalpha1, fTalpha2, fTthetaCphi, fTthetaSphi, G4Exception(), G4VSolid::GetName(), G4VSolid::kCarTolerance, and MakePlanes().

◆ G4Trap() [3/8]

G4Trap::G4Trap ( const G4String pName,
G4double  pZ,
G4double  pY,
G4double  pX,
G4double  pLTX 
)

Definition at line 145 of file G4Trap.cc.

150{
151 fDz = 0.5*pZ; fTthetaCphi = 0; fTthetaSphi = 0;
152 fDy1 = 0.5*pY; fDx1 = 0.5*pX; fDx2 = 0.5*pLTX; fTalpha1 = 0.5*(pLTX - pX)/pY;
154
156 MakePlanes();
157}

References CheckParameters(), fDx1, fDx2, fDx3, fDx4, fDy1, fDy2, fDz, fTalpha1, fTalpha2, fTthetaCphi, fTthetaSphi, and MakePlanes().

◆ G4Trap() [4/8]

G4Trap::G4Trap ( const G4String pName,
G4double  pDx1,
G4double  pDx2,
G4double  pDy1,
G4double  pDy2,
G4double  pDz 
)

Definition at line 163 of file G4Trap.cc.

168{
169 fDz = pDz; fTthetaCphi = 0; fTthetaSphi = 0;
170 fDy1 = pDy1; fDx1 = pDx1; fDx2 = pDx1; fTalpha1 = 0;
171 fDy2 = pDy2; fDx3 = pDx2; fDx4 = pDx2; fTalpha2 = 0;
172
174 MakePlanes();
175}
G4int fTrapType
Definition: G4Trap.hh:283

References CheckParameters(), fDx1, fDx2, fDx3, fDx4, fDy1, fDy2, fDz, fTalpha1, fTalpha2, fTthetaCphi, fTthetaSphi, and MakePlanes().

◆ G4Trap() [5/8]

G4Trap::G4Trap ( const G4String pName,
G4double  pDx,
G4double  pDy,
G4double  pDz,
G4double  pAlpha,
G4double  pTheta,
G4double  pPhi 
)

Definition at line 181 of file G4Trap.cc.

187{
188 fDz = pDz;
189 fTthetaCphi = std::tan(pTheta)*std::cos(pPhi);
190 fTthetaSphi = std::tan(pTheta)*std::sin(pPhi);
191
192 fDy1 = pDy; fDx1 = pDx; fDx2 = pDx; fTalpha1 = std::tan(pAlpha);
193 fDy2 = pDy; fDx3 = pDx; fDx4 = pDx; fTalpha2 = fTalpha1;
194
196 MakePlanes();
197}

References CheckParameters(), fDx1, fDx2, fDx3, fDx4, fDy1, fDy2, fDz, fTalpha1, fTalpha2, fTthetaCphi, fTthetaSphi, and MakePlanes().

◆ G4Trap() [6/8]

G4Trap::G4Trap ( const G4String pName)

Definition at line 205 of file G4Trap.cc.

207 fDz(1.), fTthetaCphi(0.), fTthetaSphi(0.),
208 fDy1(1.), fDx1(1.), fDx2(1.), fTalpha1(0.),
209 fDy2(1.), fDx3(1.), fDx4(1.), fTalpha2(0.)
210{
211 MakePlanes();
212}

References MakePlanes().

◆ ~G4Trap()

G4Trap::~G4Trap ( )
virtual

Definition at line 232 of file G4Trap.cc.

233{
234}

◆ G4Trap() [7/8]

G4Trap::G4Trap ( __void__ &  a)

Definition at line 219 of file G4Trap.cc.

221 fDz(1.), fTthetaCphi(0.), fTthetaSphi(0.),
222 fDy1(1.), fDx1(1.), fDx2(1.), fTalpha1(0.),
223 fDy2(1.), fDx3(1.), fDx4(1.), fTalpha2(0.)
224{
225 MakePlanes();
226}

References MakePlanes().

◆ G4Trap() [8/8]

G4Trap::G4Trap ( const G4Trap rhs)

Definition at line 240 of file G4Trap.cc.

243 fDy1(rhs.fDy1), fDx1(rhs.fDx1), fDx2(rhs.fDx2), fTalpha1(rhs.fTalpha1),
244 fDy2(rhs.fDy2), fDx3(rhs.fDx3), fDx4(rhs.fDx4), fTalpha2(rhs.fTalpha2)
245{
246 for (G4int i=0; i<4; ++i) { fPlanes[i] = rhs.fPlanes[i]; }
247 for (G4int i=0; i<6; ++i) { fAreas[i] = rhs.fAreas[i]; }
248 fTrapType = rhs.fTrapType;
249}
int G4int
Definition: G4Types.hh:85
G4double fAreas[6]
Definition: G4Trap.hh:282
TrapSidePlane fPlanes[4]
Definition: G4Trap.hh:281

References fAreas, fPlanes, and fTrapType.

Member Function Documentation

◆ ApproxSurfaceNormal()

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

Definition at line 803 of file G4Trap.cc.

804{
805 G4double dist = -DBL_MAX;
806 G4int iside = 0;
807 for (G4int i=0; i<4; ++i)
808 {
809 G4double d = fPlanes[i].a*p.x() +
810 fPlanes[i].b*p.y() +
811 fPlanes[i].c*p.z() + fPlanes[i].d;
812 if (d > dist) { dist = d; iside = i; }
813 }
814
815 G4double distz = std::abs(p.z()) - fDz;
816 if (dist > distz)
817 return G4ThreeVector(fPlanes[iside].a, fPlanes[iside].b, fPlanes[iside].c);
818 else
819 return G4ThreeVector(0, 0, (p.z() < 0) ? -1 : 1);
820}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
double z() const
double x() const
double y() const
G4double b
Definition: G4Trap.hh:92
G4double c
Definition: G4Trap.hh:92
G4double d
Definition: G4Trap.hh:92
G4double a
Definition: G4Trap.hh:92
#define DBL_MAX
Definition: templates.hh:62

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

Referenced by SurfaceNormal().

◆ BoundingLimits()

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

Reimplemented from G4VSolid.

Definition at line 549 of file G4Trap.cc.

550{
551 G4ThreeVector pt[8];
552 GetVertices(pt);
553
554 G4double xmin = kInfinity, xmax = -kInfinity;
555 G4double ymin = kInfinity, ymax = -kInfinity;
556 for (G4int i=0; i<8; ++i)
557 {
558 G4double x = pt[i].x();
559 if (x < xmin) xmin = x;
560 if (x > xmax) xmax = x;
561 G4double y = pt[i].y();
562 if (y < ymin) ymin = y;
563 if (y > ymax) ymax = y;
564 }
565
567 pMin.set(xmin,ymin,-dz);
568 pMax.set(xmax,ymax, dz);
569
570 // Check correctness of the bounding box
571 //
572 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
573 {
574 std::ostringstream message;
575 message << "Bad bounding box (min >= max) for solid: "
576 << GetName() << " !"
577 << "\npMin = " << pMin
578 << "\npMax = " << pMax;
579 G4Exception("G4Trap::BoundingLimits()", "GeomMgt0001",
580 JustWarning, message);
581 DumpInfo();
582 }
583}
@ JustWarning
static const G4double pMax
static const G4double pMin
G4double GetZHalfLength() const
void GetVertices(G4ThreeVector pt[8]) const
Definition: G4Trap.cc:1159
void DumpInfo() const
static const G4double kInfinity
Definition: geomdefs.hh:41

References G4VSolid::DumpInfo(), G4Exception(), G4VSolid::GetName(), GetVertices(), GetZHalfLength(), JustWarning, kInfinity, pMax, pMin, CLHEP::Hep3Vector::x(), and CLHEP::Hep3Vector::y().

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 G4Trap::CalculateExtent ( const EAxis  pAxis,
const G4VoxelLimits pVoxelLimit,
const G4AffineTransform pTransform,
G4double pMin,
G4double pMax 
) const
virtual

Implements G4VSolid.

Definition at line 589 of file G4Trap.cc.

593{
594 G4ThreeVector bmin, bmax;
595 G4bool exist;
596
597 // Check bounding box (bbox)
598 //
599 BoundingLimits(bmin,bmax);
600 G4BoundingEnvelope bbox(bmin,bmax);
601#ifdef G4BBOX_EXTENT
602 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
603#endif
604 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
605 {
606 return exist = (pMin < pMax) ? true : false;
607 }
608
609 // Set bounding envelope (benv) and calculate extent
610 //
611 G4ThreeVector pt[8];
612 GetVertices(pt);
613
614 G4ThreeVectorList baseA(4), baseB(4);
615 baseA[0] = pt[0];
616 baseA[1] = pt[1];
617 baseA[2] = pt[3];
618 baseA[3] = pt[2];
619
620 baseB[0] = pt[4];
621 baseB[1] = pt[5];
622 baseB[2] = pt[7];
623 baseB[3] = pt[6];
624
625 std::vector<const G4ThreeVectorList *> polygons(2);
626 polygons[0] = &baseA;
627 polygons[1] = &baseB;
628
629 G4BoundingEnvelope benv(bmin,bmax,polygons);
630 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
631 return exist;
632}
std::vector< G4ThreeVector > G4ThreeVectorList
bool G4bool
Definition: G4Types.hh:86
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Trap.cc:549

References G4BoundingEnvelope::BoundingBoxVsVoxelLimits(), BoundingLimits(), G4BoundingEnvelope::CalculateExtent(), GetVertices(), pMax, and pMin.

◆ CheckParameters()

void G4Trap::CheckParameters ( )
private

Definition at line 315 of file G4Trap.cc.

316{
317 if (fDz<=0 ||
318 fDy1<=0 || fDx1<=0 || fDx2<=0 ||
319 fDy2<=0 || fDx3<=0 || fDx4<=0)
320 {
321 std::ostringstream message;
322 message << "Invalid Length Parameters for Solid: " << GetName()
323 << "\n X - " <<fDx1<<", "<<fDx2<<", "<<fDx3<<", "<<fDx4
324 << "\n Y - " <<fDy1<<", "<<fDy2
325 << "\n Z - " <<fDz;
326 G4Exception("G4Trap::CheckParameters()", "GeomSolids0002",
327 FatalException, message);
328 }
329}

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

Referenced by G4Trap(), 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

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 * G4Trap::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1116 of file G4Trap.cc.

1117{
1118 return new G4Trap(*this);
1119}
G4Trap(const G4String &pName, G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pAlp1, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlp2)
Definition: G4Trap.cc:60

References G4Trap().

◆ ComputeDimensions()

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

Reimplemented from G4VSolid.

Definition at line 538 of file G4Trap.cc.

541{
542 p->ComputeDimensions(*this,n,pRep);
543}
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const

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

◆ CreatePolyhedron()

G4Polyhedron * G4Trap::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1223 of file G4Trap.cc.

1224{
1225 G4double phi = std::atan2(fTthetaSphi, fTthetaCphi);
1226 G4double alpha1 = std::atan(fTalpha1);
1227 G4double alpha2 = std::atan(fTalpha2);
1228 G4double theta = std::atan(std::sqrt(fTthetaCphi*fTthetaCphi
1230
1231 return new G4PolyhedronTrap(fDz, theta, phi,
1232 fDy1, fDx1, fDx2, alpha1,
1233 fDy2, fDx3, fDx4, alpha2);
1234}
const G4double alpha2

References alpha2, fDx1, fDx2, fDx3, fDx4, fDy1, fDy2, fDz, fTalpha1, fTalpha2, fTthetaCphi, and fTthetaSphi.

◆ DescribeYourselfTo()

void G4Trap::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 1218 of file G4Trap.cc.

1219{
1220 scene.AddSolid (*this);
1221}
virtual void AddSolid(const G4Box &)=0

References G4VGraphicsScene::AddSolid().

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 896 of file G4Trap.cc.

897{
898 switch (fTrapType)
899 {
900 case 0: // General case
901 {
902 G4double dz = std::abs(p.z())-fDz;
903 G4double dy1 = fPlanes[0].b*p.y()+fPlanes[0].c*p.z()+fPlanes[0].d;
904 G4double dy2 = fPlanes[1].b*p.y()+fPlanes[1].c*p.z()+fPlanes[1].d;
905 G4double dy = std::max(dz,std::max(dy1,dy2));
906
907 G4double dx1 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()
908 + fPlanes[2].c*p.z()+fPlanes[2].d;
909 G4double dx2 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
910 + fPlanes[3].c*p.z()+fPlanes[3].d;
911 G4double dist = std::max(dy,std::max(dx1,dx2));
912 return (dist > 0) ? dist : 0.;
913 }
914 case 1: // YZ section is a rectangle
915 {
916 G4double dz = std::abs(p.z())-fDz;
917 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
918 G4double dx1 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()
919 + fPlanes[2].c*p.z()+fPlanes[2].d;
920 G4double dx2 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
921 + fPlanes[3].c*p.z()+fPlanes[3].d;
922 G4double dist = std::max(dy,std::max(dx1,dx2));
923 return (dist > 0) ? dist : 0.;
924 }
925 case 2: // YZ section is a rectangle and
926 { // XZ section is an isosceles trapezoid
927 G4double dz = std::abs(p.z())-fDz;
928 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
929 G4double dx = fPlanes[3].a*std::abs(p.x())
930 + fPlanes[3].c*p.z()+fPlanes[3].d;
931 G4double dist = std::max(dy,dx);
932 return (dist > 0) ? dist : 0.;
933 }
934 case 3: // YZ section is a rectangle and
935 { // XY section is an isosceles trapezoid
936 G4double dz = std::abs(p.z())-fDz;
937 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
938 G4double dx = fPlanes[3].a*std::abs(p.x())
939 + fPlanes[3].b*p.y()+fPlanes[3].d;
940 G4double dist = std::max(dy,dx);
941 return (dist > 0) ? dist : 0.;
942 }
943 }
944 return 0.;
945}
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References TrapSidePlane::a, TrapSidePlane::b, TrapSidePlane::c, TrapSidePlane::d, fDz, fPlanes, fTrapType, G4INCL::Math::max(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ DistanceToIn() [2/2]

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

Implements G4VSolid.

Definition at line 827 of file G4Trap.cc.

829{
830 // Z intersections
831 //
832 if ((std::abs(p.z()) - fDz) >= -halfCarTolerance && p.z()*v.z() >= 0)
833 return kInfinity;
834 G4double invz = (-v.z() == 0) ? DBL_MAX : -1./v.z();
835 G4double dz = (invz < 0) ? fDz : -fDz;
836 G4double tzmin = (p.z() + dz)*invz;
837 G4double tzmax = (p.z() - dz)*invz;
838
839 // Y intersections
840 //
841 G4double tymin = 0, tymax = DBL_MAX;
842 G4int i = 0;
843 for ( ; i<2; ++i)
844 {
845 G4double cosa = fPlanes[i].b*v.y() + fPlanes[i].c*v.z();
846 G4double dist = fPlanes[i].b*p.y() + fPlanes[i].c*p.z() + fPlanes[i].d;
847 if (dist >= -halfCarTolerance)
848 {
849 if (cosa >= 0) return kInfinity;
850 G4double tmp = -dist/cosa;
851 if (tymin < tmp) tymin = tmp;
852 }
853 else if (cosa > 0)
854 {
855 G4double tmp = -dist/cosa;
856 if (tymax > tmp) tymax = tmp;
857 }
858 }
859
860 // Z intersections
861 //
862 G4double txmin = 0, txmax = DBL_MAX;
863 for ( ; i<4; ++i)
864 {
865 G4double cosa = fPlanes[i].a*v.x()+fPlanes[i].b*v.y()+fPlanes[i].c*v.z();
866 G4double dist = fPlanes[i].a*p.x()+fPlanes[i].b*p.y()+fPlanes[i].c*p.z() +
867 fPlanes[i].d;
868 if (dist >= -halfCarTolerance)
869 {
870 if (cosa >= 0) return kInfinity;
871 G4double tmp = -dist/cosa;
872 if (txmin < tmp) txmin = tmp;
873 }
874 else if (cosa > 0)
875 {
876 G4double tmp = -dist/cosa;
877 if (txmax > tmp) txmax = tmp;
878 }
879 }
880
881 // Find distance
882 //
883 G4double tmin = std::max(std::max(txmin,tymin),tzmin);
884 G4double tmax = std::min(std::min(txmax,tymax),tzmax);
885
886 if (tmax <= tmin + halfCarTolerance) return kInfinity; // touch or no hit
887 return (tmin < halfCarTolerance ) ? 0. : tmin;
888}
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References TrapSidePlane::a, TrapSidePlane::b, TrapSidePlane::c, TrapSidePlane::d, DBL_MAX, fDz, fPlanes, halfCarTolerance, kInfinity, G4INCL::Math::max(), G4INCL::Math::min(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ DistanceToOut() [1/2]

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

Implements G4VSolid.

Definition at line 1036 of file G4Trap.cc.

1037{
1038#ifdef G4CSGDEBUG
1039 if( Inside(p) == kOutside )
1040 {
1041 std::ostringstream message;
1042 G4int oldprc = message.precision(16);
1043 message << "Point p is outside (!?) of solid: " << GetName() << G4endl;
1044 message << "Position:\n";
1045 message << " p.x() = " << p.x()/mm << " mm\n";
1046 message << " p.y() = " << p.y()/mm << " mm\n";
1047 message << " p.z() = " << p.z()/mm << " mm";
1048 G4cout.precision(oldprc);
1049 G4Exception("G4Trap::DistanceToOut(p)", "GeomSolids1002",
1050 JustWarning, message );
1051 DumpInfo();
1052 }
1053#endif
1054 switch (fTrapType)
1055 {
1056 case 0: // General case
1057 {
1058 G4double dz = std::abs(p.z())-fDz;
1059 G4double dy1 = fPlanes[0].b*p.y()+fPlanes[0].c*p.z()+fPlanes[0].d;
1060 G4double dy2 = fPlanes[1].b*p.y()+fPlanes[1].c*p.z()+fPlanes[1].d;
1061 G4double dy = std::max(dz,std::max(dy1,dy2));
1062
1063 G4double dx1 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()
1064 + fPlanes[2].c*p.z()+fPlanes[2].d;
1065 G4double dx2 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
1066 + fPlanes[3].c*p.z()+fPlanes[3].d;
1067 G4double dist = std::max(dy,std::max(dx1,dx2));
1068 return (dist < 0) ? -dist : 0.;
1069 }
1070 case 1: // YZ section is a rectangle
1071 {
1072 G4double dz = std::abs(p.z())-fDz;
1073 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
1074 G4double dx1 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()
1075 + fPlanes[2].c*p.z()+fPlanes[2].d;
1076 G4double dx2 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
1077 + fPlanes[3].c*p.z()+fPlanes[3].d;
1078 G4double dist = std::max(dy,std::max(dx1,dx2));
1079 return (dist < 0) ? -dist : 0.;
1080 }
1081 case 2: // YZ section is a rectangle and
1082 { // XZ section is an isosceles trapezoid
1083 G4double dz = std::abs(p.z())-fDz;
1084 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
1085 G4double dx = fPlanes[3].a*std::abs(p.x())
1086 + fPlanes[3].c*p.z()+fPlanes[3].d;
1087 G4double dist = std::max(dy,dx);
1088 return (dist < 0) ? -dist : 0.;
1089 }
1090 case 3: // YZ section is a rectangle and
1091 { // XY section is an isosceles trapezoid
1092 G4double dz = std::abs(p.z())-fDz;
1093 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
1094 G4double dx = fPlanes[3].a*std::abs(p.x())
1095 + fPlanes[3].b*p.y()+fPlanes[3].d;
1096 G4double dist = std::max(dy,dx);
1097 return (dist < 0) ? -dist : 0.;
1098 }
1099 }
1100 return 0.;
1101}
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: G4Trap.cc:638
@ kOutside
Definition: geomdefs.hh:68

References TrapSidePlane::a, TrapSidePlane::b, TrapSidePlane::c, TrapSidePlane::d, G4VSolid::DumpInfo(), fDz, fPlanes, fTrapType, 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 G4Trap::DistanceToOut ( const G4ThreeVector p,
const G4ThreeVector v,
const G4bool  calcNorm = false,
G4bool validNorm = nullptr,
G4ThreeVector n = nullptr 
) const
virtual

Implements G4VSolid.

Definition at line 953 of file G4Trap.cc.

956{
957 // Z intersections
958 //
959 if ((std::abs(p.z()) - fDz) >= -halfCarTolerance && p.z()*v.z() > 0)
960 {
961 if (calcNorm)
962 {
963 *validNorm = true;
964 n->set(0, 0, (p.z() < 0) ? -1 : 1);
965 }
966 return 0;
967 }
968 G4double vz = v.z();
969 G4double tmax = (vz == 0) ? DBL_MAX : (std::copysign(fDz,vz) - p.z())/vz;
970 G4int iside = (vz < 0) ? -4 : -2; // little trick: (-4+3)=-1, (-2+3)=+1
971
972 // Y intersections
973 //
974 G4int i = 0;
975 for ( ; i<2; ++i)
976 {
977 G4double cosa = fPlanes[i].b*v.y() + fPlanes[i].c*v.z();
978 if (cosa > 0)
979 {
980 G4double dist = fPlanes[i].b*p.y() + fPlanes[i].c*p.z() + fPlanes[i].d;
981 if (dist >= -halfCarTolerance)
982 {
983 if (calcNorm)
984 {
985 *validNorm = true;
986 n->set(0, fPlanes[i].b, fPlanes[i].c);
987 }
988 return 0;
989 }
990 G4double tmp = -dist/cosa;
991 if (tmax > tmp) { tmax = tmp; iside = i; }
992 }
993 }
994
995 // X intersections
996 //
997 for ( ; i<4; ++i)
998 {
999 G4double cosa = fPlanes[i].a*v.x()+fPlanes[i].b*v.y()+fPlanes[i].c*v.z();
1000 if (cosa > 0)
1001 {
1002 G4double dist = fPlanes[i].a*p.x() +
1003 fPlanes[i].b*p.y() + fPlanes[i].c*p.z() + fPlanes[i].d;
1004 if (dist >= -halfCarTolerance)
1005 {
1006 if (calcNorm)
1007 {
1008 *validNorm = true;
1009 n->set(fPlanes[i].a, fPlanes[i].b, fPlanes[i].c);
1010 }
1011 return 0;
1012 }
1013 G4double tmp = -dist/cosa;
1014 if (tmax > tmp) { tmax = tmp; iside = i; }
1015 }
1016 }
1017
1018 // Set normal, if required, and return distance
1019 //
1020 if (calcNorm)
1021 {
1022 *validNorm = true;
1023 if (iside < 0)
1024 n->set(0, 0, iside + 3); // (-4+3)=-1, (-2+3)=+1
1025 else
1026 n->set(fPlanes[iside].a, fPlanes[iside].b, fPlanes[iside].c);
1027 }
1028 return tmax;
1029}

References TrapSidePlane::a, TrapSidePlane::b, TrapSidePlane::c, TrapSidePlane::d, 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(), BoundingLimits(), G4Trd::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(), DistanceToOut(), G4Trd::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(), SurfaceNormal(), G4Trd::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

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().

◆ GetAlpha1()

G4double G4Trap::GetAlpha1 ( ) const
inline

Referenced by StreamInfo().

◆ GetAlpha2()

G4double G4Trap::GetAlpha2 ( ) const
inline

Referenced by StreamInfo().

◆ 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 G4Trap::GetCubicVolume ( )
virtual

Reimplemented from G4VSolid.

Definition at line 488 of file G4Trap.cc.

489{
490 if (fCubicVolume == 0)
491 {
492 G4ThreeVector pt[8];
493 GetVertices(pt);
494
495 G4double dz = pt[4].z() - pt[0].z();
496 G4double dy1 = pt[2].y() - pt[0].y();
497 G4double dx1 = pt[1].x() - pt[0].x();
498 G4double dx2 = pt[3].x() - pt[2].x();
499 G4double dy2 = pt[6].y() - pt[4].y();
500 G4double dx3 = pt[5].x() - pt[4].x();
501 G4double dx4 = pt[7].x() - pt[6].x();
502
503 fCubicVolume = ((dx1 + dx2 + dx3 + dx4)*(dy1 + dy2) +
504 (dx4 + dx3 - dx2 - dx1)*(dy2 - dy1)/3)*dz*0.125;
505 }
506 return fCubicVolume;
507}
G4double fCubicVolume
Definition: G4CSGSolid.hh:70

References G4CSGSolid::fCubicVolume, GetVertices(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ 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 G4Trap::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 1107 of file G4Trap.cc.

1108{
1109 return G4String("G4Trap");
1110}

◆ 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(), BoundingLimits(), G4Trd::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(), CheckParameters(), G4Trd::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(), DistanceToOut(), G4Trd::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(), 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(), 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(), StreamInfo(), G4Trd::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(), SurfaceNormal(), G4Trd::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().

◆ GetPhi()

G4double G4Trap::GetPhi ( ) const
inline

Referenced by StreamInfo().

◆ GetPointOnSurface()

G4ThreeVector G4Trap::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1177 of file G4Trap.cc.

1178{
1179 // Set indeces
1180 constexpr G4int iface [6][4] =
1181 { {0,1,3,2}, {0,4,5,1}, {2,3,7,6}, {0,2,6,4}, {1,5,7,3}, {4,6,7,5} };
1182
1183 // Set vertices
1184 G4ThreeVector pt[8];
1185 GetVertices(pt);
1186
1187 // Select face
1188 //
1189 G4double select = fAreas[5]*G4QuickRand();
1190 G4int k = 5;
1191 k -= (select <= fAreas[4]);
1192 k -= (select <= fAreas[3]);
1193 k -= (select <= fAreas[2]);
1194 k -= (select <= fAreas[1]);
1195 k -= (select <= fAreas[0]);
1196
1197 // Select sub-triangle
1198 //
1199 G4int i0 = iface[k][0];
1200 G4int i1 = iface[k][1];
1201 G4int i2 = iface[k][2];
1202 G4int i3 = iface[k][3];
1203 G4double s2 = G4GeomTools::TriangleAreaNormal(pt[i2],pt[i1],pt[i3]).mag();
1204 if (select > fAreas[k] - s2) i0 = i2;
1205
1206 // Generate point
1207 //
1208 G4double u = G4QuickRand();
1209 G4double v = G4QuickRand();
1210 if (u + v > 1.) { u = 1. - u; v = 1. - v; }
1211 return (1.-u-v)*pt[i0] + u*pt[i1] + v*pt[i3];
1212}
double mag() const
static G4ThreeVector TriangleAreaNormal(const G4ThreeVector &A, const G4ThreeVector &B, const G4ThreeVector &C)
Definition: G4GeomTools.cc:598

References fAreas, G4QuickRand(), GetVertices(), CLHEP::Hep3Vector::mag(), and G4GeomTools::TriangleAreaNormal().

◆ 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().

◆ GetSidePlane()

TrapSidePlane G4Trap::GetSidePlane ( G4int  n) const
inline

Referenced by export_G4Trap().

◆ GetSurfaceArea()

G4double G4Trap::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 513 of file G4Trap.cc.

514{
515 if (fSurfaceArea == 0)
516 {
517 G4ThreeVector pt[8];
518 G4int iface [6][4] =
519 { {0,1,3,2}, {0,4,5,1}, {2,3,7,6}, {0,2,6,4}, {1,5,7,3}, {4,6,7,5} };
520
521 GetVertices(pt);
522 for (G4int i=0; i<6; ++i)
523 {
524 fSurfaceArea += G4GeomTools::QuadAreaNormal(pt[iface[i][0]],
525 pt[iface[i][1]],
526 pt[iface[i][2]],
527 pt[iface[i][3]]).mag();
528 }
529 }
530 return fSurfaceArea;
531}
G4double fSurfaceArea
Definition: G4CSGSolid.hh:71
static G4ThreeVector QuadAreaNormal(const G4ThreeVector &A, const G4ThreeVector &B, const G4ThreeVector &C, const G4ThreeVector &D)
Definition: G4GeomTools.cc:609

References G4CSGSolid::fSurfaceArea, GetVertices(), CLHEP::Hep3Vector::mag(), and G4GeomTools::QuadAreaNormal().

◆ GetSymAxis()

G4ThreeVector G4Trap::GetSymAxis ( ) const
inline

◆ GetTanAlpha1()

G4double G4Trap::GetTanAlpha1 ( ) const
inline

◆ GetTanAlpha2()

G4double G4Trap::GetTanAlpha2 ( ) const
inline

◆ GetTheta()

G4double G4Trap::GetTheta ( ) const
inline

Referenced by StreamInfo().

◆ GetTolerance()

G4double G4VSolid::GetTolerance ( ) const
inlineinherited

◆ GetVertices()

void G4Trap::GetVertices ( G4ThreeVector  pt[8]) const
private

Definition at line 1159 of file G4Trap.cc.

1160{
1161 for (G4int i=0; i<8; ++i)
1162 {
1163 G4int iy = (i==0 || i==1 || i==4 || i==5) ? 0 : 1;
1164 G4int ix = (i==0 || i==2 || i==4 || i==6) ? 2 : 3;
1165 G4double z = (i < 4) ? -fDz : fDz;
1166 G4double y = -(fPlanes[iy].c*z + fPlanes[iy].d)/fPlanes[iy].b;
1167 G4double x = -(fPlanes[ix].b*y + fPlanes[ix].c*z
1168 + fPlanes[ix].d)/fPlanes[ix].a;
1169 pt[i].set(x,y,z);
1170 }
1171}
void set(double x, double y, double z)

References TrapSidePlane::b, TrapSidePlane::c, TrapSidePlane::d, fDz, fPlanes, and CLHEP::Hep3Vector::set().

Referenced by BoundingLimits(), CalculateExtent(), GetCubicVolume(), GetPointOnSurface(), GetSurfaceArea(), and SetCachedValues().

◆ GetXHalfLength1()

G4double G4Trap::GetXHalfLength1 ( ) const
inline

◆ GetXHalfLength2()

G4double G4Trap::GetXHalfLength2 ( ) const
inline

◆ GetXHalfLength3()

G4double G4Trap::GetXHalfLength3 ( ) const
inline

◆ GetXHalfLength4()

G4double G4Trap::GetXHalfLength4 ( ) const
inline

◆ GetYHalfLength1()

G4double G4Trap::GetYHalfLength1 ( ) const
inline

◆ GetYHalfLength2()

G4double G4Trap::GetYHalfLength2 ( ) const
inline

◆ GetZHalfLength()

G4double G4Trap::GetZHalfLength ( ) const
inline

◆ Inside()

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

Implements G4VSolid.

Definition at line 638 of file G4Trap.cc.

639{
640 switch (fTrapType)
641 {
642 case 0: // General case
643 {
644 G4double dz = std::abs(p.z())-fDz;
645 G4double dy1 = fPlanes[0].b*p.y()+fPlanes[0].c*p.z()+fPlanes[0].d;
646 G4double dy2 = fPlanes[1].b*p.y()+fPlanes[1].c*p.z()+fPlanes[1].d;
647 G4double dy = std::max(dz,std::max(dy1,dy2));
648
649 G4double dx1 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()
650 + fPlanes[2].c*p.z()+fPlanes[2].d;
651 G4double dx2 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
652 + fPlanes[3].c*p.z()+fPlanes[3].d;
653 G4double dist = std::max(dy,std::max(dx1,dx2));
654
655 return (dist > halfCarTolerance) ? kOutside :
656 ((dist > -halfCarTolerance) ? kSurface : kInside);
657 }
658 case 1: // YZ section is a rectangle
659 {
660 G4double dz = std::abs(p.z())-fDz;
661 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
662 G4double dx1 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()
663 + fPlanes[2].c*p.z()+fPlanes[2].d;
664 G4double dx2 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
665 + fPlanes[3].c*p.z()+fPlanes[3].d;
666 G4double dist = std::max(dy,std::max(dx1,dx2));
667
668 return (dist > halfCarTolerance) ? kOutside :
669 ((dist > -halfCarTolerance) ? kSurface : kInside);
670 }
671 case 2: // YZ section is a rectangle and
672 { // XZ section is an isosceles trapezoid
673 G4double dz = std::abs(p.z())-fDz;
674 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
675 G4double dx = fPlanes[3].a*std::abs(p.x())
676 + fPlanes[3].c*p.z()+fPlanes[3].d;
677 G4double dist = std::max(dy,dx);
678
679 return (dist > halfCarTolerance) ? kOutside :
680 ((dist > -halfCarTolerance) ? kSurface : kInside);
681 }
682 case 3: // YZ section is a rectangle and
683 { // XY section is an isosceles trapezoid
684 G4double dz = std::abs(p.z())-fDz;
685 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
686 G4double dx = fPlanes[3].a*std::abs(p.x())
687 + fPlanes[3].b*p.y()+fPlanes[3].d;
688 G4double dist = std::max(dy,dx);
689
690 return (dist > halfCarTolerance) ? kOutside :
691 ((dist > -halfCarTolerance) ? kSurface : kInside);
692 }
693 }
694 return kOutside;
695}
@ kSurface
Definition: geomdefs.hh:69

References TrapSidePlane::a, TrapSidePlane::b, TrapSidePlane::c, TrapSidePlane::d, fDz, fPlanes, fTrapType, halfCarTolerance, kInside, kOutside, kSurface, G4INCL::Math::max(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by DistanceToOut().

◆ MakePlane()

G4bool G4Trap::MakePlane ( const G4ThreeVector p1,
const G4ThreeVector p2,
const G4ThreeVector p3,
const G4ThreeVector p4,
TrapSidePlane plane 
)
protected

Definition at line 402 of file G4Trap.cc.

407{
408 G4ThreeVector normal = ((p4 - p2).cross(p3 - p1)).unit();
409 if (std::abs(normal.x()) < DBL_EPSILON) normal.setX(0);
410 if (std::abs(normal.y()) < DBL_EPSILON) normal.setY(0);
411 if (std::abs(normal.z()) < DBL_EPSILON) normal.setZ(0);
412 normal = normal.unit();
413
414 G4ThreeVector centre = (p1 + p2 + p3 + p4)*0.25;
415 plane.a = normal.x();
416 plane.b = normal.y();
417 plane.c = normal.z();
418 plane.d = -normal.dot(centre);
419
420 // compute distances and check planarity
421 G4double d1 = std::abs(normal.dot(p1) + plane.d);
422 G4double d2 = std::abs(normal.dot(p2) + plane.d);
423 G4double d3 = std::abs(normal.dot(p3) + plane.d);
424 G4double d4 = std::abs(normal.dot(p4) + plane.d);
425 G4double dmax = std::max(std::max(std::max(d1,d2),d3),d4);
426
427 return (dmax > 1000 * kCarTolerance) ? false : true;
428}
static const G4double d1
static const G4double d2
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:79
#define DBL_EPSILON
Definition: templates.hh:66

References TrapSidePlane::a, TrapSidePlane::b, TrapSidePlane::c, TrapSidePlane::d, d1, d2, DBL_EPSILON, G4VSolid::kCarTolerance, G4INCL::Math::max(), and CLHEP::normal().

Referenced by MakePlanes().

◆ MakePlanes() [1/2]

void G4Trap::MakePlanes ( )
protected

Definition at line 335 of file G4Trap.cc.

336{
337 G4double DzTthetaCphi = fDz*fTthetaCphi;
338 G4double DzTthetaSphi = fDz*fTthetaSphi;
339 G4double Dy1Talpha1 = fDy1*fTalpha1;
340 G4double Dy2Talpha2 = fDy2*fTalpha2;
341
342 G4ThreeVector pt[8] =
343 {
344 G4ThreeVector(-DzTthetaCphi-Dy1Talpha1-fDx1,-DzTthetaSphi-fDy1,-fDz),
345 G4ThreeVector(-DzTthetaCphi-Dy1Talpha1+fDx1,-DzTthetaSphi-fDy1,-fDz),
346 G4ThreeVector(-DzTthetaCphi+Dy1Talpha1-fDx2,-DzTthetaSphi+fDy1,-fDz),
347 G4ThreeVector(-DzTthetaCphi+Dy1Talpha1+fDx2,-DzTthetaSphi+fDy1,-fDz),
348 G4ThreeVector( DzTthetaCphi-Dy2Talpha2-fDx3, DzTthetaSphi-fDy2, fDz),
349 G4ThreeVector( DzTthetaCphi-Dy2Talpha2+fDx3, DzTthetaSphi-fDy2, fDz),
350 G4ThreeVector( DzTthetaCphi+Dy2Talpha2-fDx4, DzTthetaSphi+fDy2, fDz),
351 G4ThreeVector( DzTthetaCphi+Dy2Talpha2+fDx4, DzTthetaSphi+fDy2, fDz)
352 };
353
354 MakePlanes(pt);
355}

References fDx1, fDx2, fDx3, fDx4, fDy1, fDy2, fDz, fTalpha1, fTalpha2, fTthetaCphi, fTthetaSphi, and MakePlanes().

Referenced by G4Trap(), MakePlanes(), and SetAllParameters().

◆ MakePlanes() [2/2]

void G4Trap::MakePlanes ( const G4ThreeVector  pt[8])
protected

Definition at line 361 of file G4Trap.cc.

362{
363 constexpr G4int iface[4][4] = { {0,4,5,1}, {2,3,7,6}, {0,2,6,4}, {1,5,7,3} };
364 const static G4String side[4] = { "~-Y", "~+Y", "~-X", "~+X" };
365
366 for (G4int i=0; i<4; ++i)
367 {
368 if (MakePlane(pt[iface[i][0]],
369 pt[iface[i][1]],
370 pt[iface[i][2]],
371 pt[iface[i][3]],
372 fPlanes[i])) continue;
373
374 // Non planar side face
376 G4double dmax = 0;
377 for (G4int k=0; k<4; ++k)
378 {
379 G4double dist = normal.dot(pt[iface[i][k]]) + fPlanes[i].d;
380 if (std::abs(dist) > std::abs(dmax)) dmax = dist;
381 }
382 std::ostringstream message;
383 message << "Side face " << side[i] << " is not planar for solid: "
384 << GetName() << "\nDiscrepancy: " << dmax/mm << " mm\n";
385 StreamInfo(message);
386 G4Exception("G4Trap::MakePlanes()", "GeomSolids0002",
387 FatalException, message);
388 }
389
390 // Re-compute parameters
392}
void SetCachedValues()
Definition: G4Trap.cc:434
G4bool MakePlane(const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3, const G4ThreeVector &p4, TrapSidePlane &plane)
Definition: G4Trap.cc:402
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Trap.cc:1125

References TrapSidePlane::d, FatalException, fPlanes, G4Exception(), G4VSolid::GetName(), MakePlane(), mm, CLHEP::normal(), SetCachedValues(), and StreamInfo().

◆ operator=()

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

Definition at line 255 of file G4Trap.cc.

256{
257 // Check assignment to self
258 //
259 if (this == &rhs) { return *this; }
260
261 // Copy base class data
262 //
264
265 // Copy data
266 //
269 fDy1 = rhs.fDy1; fDx1 = rhs.fDx1; fDx2 = rhs.fDx2; fTalpha1 = rhs.fTalpha1;
270 fDy2 = rhs.fDy2; fDx3 = rhs.fDx3; fDx4 = rhs.fDx4; fTalpha2 = rhs.fTalpha2;
271 for (G4int i=0; i<4; ++i) { fPlanes[i] = rhs.fPlanes[i]; }
272 for (G4int i=0; i<6; ++i) { fAreas[i] = rhs.fAreas[i]; }
273 fTrapType = rhs.fTrapType;
274 return *this;
275}
G4CSGSolid & operator=(const G4CSGSolid &rhs)
Definition: G4CSGSolid.cc:89

References fAreas, fDx1, fDx2, fDx3, fDx4, fDy1, fDy2, fDz, fPlanes, fTalpha1, fTalpha2, fTrapType, fTthetaCphi, fTthetaSphi, halfCarTolerance, and G4CSGSolid::operator=().

◆ operator==()

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

◆ SetAllParameters()

void G4Trap::SetAllParameters ( G4double  pDz,
G4double  pTheta,
G4double  pPhi,
G4double  pDy1,
G4double  pDx1,
G4double  pDx2,
G4double  pAlp1,
G4double  pDy2,
G4double  pDx3,
G4double  pDx4,
G4double  pAlp2 
)

Definition at line 282 of file G4Trap.cc.

293{
294 // Reset data of the base class
295 fCubicVolume = 0;
296 fSurfaceArea = 0;
297 fRebuildPolyhedron = true;
298
299 // Set parameters
300 fDz = pDz;
301 fTthetaCphi = std::tan(pTheta)*std::cos(pPhi);
302 fTthetaSphi = std::tan(pTheta)*std::sin(pPhi);
303
304 fDy1 = pDy1; fDx1 = pDx1; fDx2 = pDx2; fTalpha1 = std::tan(pAlp1);
305 fDy2 = pDy2; fDx3 = pDx3; fDx4 = pDx4; fTalpha2 = std::tan(pAlp2);
306
308 MakePlanes();
309}

References CheckParameters(), G4CSGSolid::fCubicVolume, fDx1, fDx2, fDx3, fDx4, fDy1, fDy2, fDz, G4CSGSolid::fRebuildPolyhedron, G4CSGSolid::fSurfaceArea, fTalpha1, fTalpha2, fTthetaCphi, fTthetaSphi, and MakePlanes().

Referenced by G4ParameterisationTrdY::ComputeDimensions(), G4GDMLParameterisation::ComputeDimensions(), G4ParameterisationTrdX::ComputeDimensions(), and export_G4Trap().

◆ SetCachedValues()

void G4Trap::SetCachedValues ( )
protected

Definition at line 434 of file G4Trap.cc.

435{
436 // Set indeces
437 constexpr G4int iface[6][4] =
438 { {0,1,3,2}, {0,4,5,1}, {2,3,7,6}, {0,2,6,4}, {1,5,7,3}, {4,6,7,5} };
439
440 // Get vertices
441 G4ThreeVector pt[8];
442 GetVertices(pt);
443
444 // Set face areas
445 for (G4int i=0; i<6; ++i)
446 {
447 fAreas[i] = G4GeomTools::QuadAreaNormal(pt[iface[i][0]],
448 pt[iface[i][1]],
449 pt[iface[i][2]],
450 pt[iface[i][3]]).mag();
451 }
452 for (G4int i=1; i<6; ++i) { fAreas[i] += fAreas[i - 1]; }
453
454 // Define type of trapezoid
455 fTrapType = 0;
456 if (fPlanes[0].b == -1 && fPlanes[1].b == 1 &&
457 std::abs(fPlanes[0].a) < DBL_EPSILON &&
458 std::abs(fPlanes[0].c) < DBL_EPSILON &&
459 std::abs(fPlanes[1].a) < DBL_EPSILON &&
460 std::abs(fPlanes[1].c) < DBL_EPSILON)
461 {
462 fTrapType = 1; // YZ section is a rectangle ...
463 if (std::abs(fPlanes[2].a + fPlanes[3].a) < DBL_EPSILON &&
464 std::abs(fPlanes[2].c - fPlanes[3].c) < DBL_EPSILON &&
465 fPlanes[2].b == 0 &&
466 fPlanes[3].b == 0)
467 {
468 fTrapType = 2; // ... and XZ section is a isosceles trapezoid
469 fPlanes[2].a = -fPlanes[3].a;
470 fPlanes[2].c = fPlanes[3].c;
471 }
472 if (std::abs(fPlanes[2].a + fPlanes[3].a) < DBL_EPSILON &&
473 std::abs(fPlanes[2].b - fPlanes[3].b) < DBL_EPSILON &&
474 fPlanes[2].c == 0 &&
475 fPlanes[3].c == 0)
476 {
477 fTrapType = 3; // ... and XY section is a isosceles trapezoid
478 fPlanes[2].a = -fPlanes[3].a;
479 fPlanes[2].b = fPlanes[3].b;
480 }
481 }
482}

References TrapSidePlane::a, TrapSidePlane::b, TrapSidePlane::c, DBL_EPSILON, fAreas, fPlanes, fTrapType, GetVertices(), CLHEP::Hep3Vector::mag(), and G4GeomTools::QuadAreaNormal().

Referenced by MakePlanes().

◆ 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().

◆ StreamInfo()

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

Reimplemented from G4CSGSolid.

Definition at line 1125 of file G4Trap.cc.

1126{
1127 G4double phi = GetPhi();
1128 G4double theta = GetTheta();
1129 G4double alpha1 = GetAlpha1();
1131
1132 G4int oldprc = os.precision(16);
1133 os << "-----------------------------------------------------------\n"
1134 << " *** Dump for solid: " << GetName() << " ***\n"
1135 << " ===================================================\n"
1136 << " Solid type: G4Trap\n"
1137 << " Parameters:\n"
1138 << " half length Z: " << fDz/mm << " mm\n"
1139 << " half length Y, face -Dz: " << fDy1/mm << " mm\n"
1140 << " half length X, face -Dz, side -Dy1: " << fDx1/mm << " mm\n"
1141 << " half length X, face -Dz, side +Dy1: " << fDx2/mm << " mm\n"
1142 << " half length Y, face +Dz: " << fDy2/mm << " mm\n"
1143 << " half length X, face +Dz, side -Dy2: " << fDx3/mm << " mm\n"
1144 << " half length X, face +Dz, side +Dy2: " << fDx4/mm << " mm\n"
1145 << " theta: " << theta/degree << " degrees\n"
1146 << " phi: " << phi/degree << " degrees\n"
1147 << " alpha, face -Dz: " << alpha1/degree << " degrees\n"
1148 << " alpha, face +Dz: " << alpha2/degree << " degrees\n"
1149 << "-----------------------------------------------------------\n";
1150 os.precision(oldprc);
1151
1152 return os;
1153}
static constexpr double degree
Definition: G4SIunits.hh:124
G4double GetAlpha2() const
G4double GetAlpha1() const
G4double GetTheta() const
G4double GetPhi() const

References alpha2, degree, fDx1, fDx2, fDx3, fDx4, fDy1, fDy2, fDz, GetAlpha1(), GetAlpha2(), G4VSolid::GetName(), GetPhi(), GetTheta(), and mm.

Referenced by MakePlanes().

◆ SurfaceNormal()

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

Implements G4VSolid.

Definition at line 701 of file G4Trap.cc.

702{
703 G4double nx = 0, ny = 0, nz = 0;
704 G4double dz = std::abs(p.z()) - fDz;
705 nz = std::copysign(G4double(std::abs(dz) <= halfCarTolerance), p.z());
706
707 switch (fTrapType)
708 {
709 case 0: // General case
710 {
711 for (G4int i=0; i<2; ++i)
712 {
713 G4double dy = fPlanes[i].b*p.y() + fPlanes[i].c*p.z() + fPlanes[i].d;
714 if (std::abs(dy) > halfCarTolerance) continue;
715 ny = fPlanes[i].b;
716 nz += fPlanes[i].c;
717 break;
718 }
719 for (G4int i=2; i<4; ++i)
720 {
721 G4double dx = fPlanes[i].a*p.x() +
722 fPlanes[i].b*p.y() + fPlanes[i].c*p.z() + fPlanes[i].d;
723 if (std::abs(dx) > halfCarTolerance) continue;
724 nx = fPlanes[i].a;
725 ny += fPlanes[i].b;
726 nz += fPlanes[i].c;
727 break;
728 }
729 break;
730 }
731 case 1: // YZ section - rectangle
732 {
733 G4double dy = std::abs(p.y()) + fPlanes[1].d;
734 ny = std::copysign(G4double(std::abs(dy) <= halfCarTolerance), p.y());
735 for (G4int i=2; i<4; ++i)
736 {
737 G4double dx = fPlanes[i].a*p.x() +
738 fPlanes[i].b*p.y() + fPlanes[i].c*p.z() + fPlanes[i].d;
739 if (std::abs(dx) > halfCarTolerance) continue;
740 nx = fPlanes[i].a;
741 ny += fPlanes[i].b;
742 nz += fPlanes[i].c;
743 break;
744 }
745 break;
746 }
747 case 2: // YZ section - rectangle, XZ section - isosceles trapezoid
748 {
749 G4double dy = std::abs(p.y()) + fPlanes[1].d;
750 ny = std::copysign(G4double(std::abs(dy) <= halfCarTolerance), p.y());
751 G4double dx = fPlanes[3].a*std::abs(p.x()) +
752 fPlanes[3].c*p.z() + fPlanes[3].d;
753 G4double k = std::abs(dx) <= halfCarTolerance;
754 nx = std::copysign(k, p.x())*fPlanes[3].a;
755 nz += k*fPlanes[3].c;
756 break;
757 }
758 case 3: // YZ section - rectangle, XY section - isosceles trapezoid
759 {
760 G4double dy = std::abs(p.y()) + fPlanes[1].d;
761 ny = std::copysign(G4double(std::abs(dy) <= halfCarTolerance), p.y());
762 G4double dx = fPlanes[3].a*std::abs(p.x()) +
763 fPlanes[3].b*p.y() + fPlanes[3].d;
764 G4double k = std::abs(dx) <= halfCarTolerance;
765 nx = std::copysign(k, p.x())*fPlanes[3].a;
766 ny += k*fPlanes[3].b;
767 break;
768 }
769 }
770
771 // Return normal
772 //
773 G4double mag2 = nx*nx + ny*ny + nz*nz;
774 if (mag2 == 1) return G4ThreeVector(nx,ny,nz);
775 else if (mag2 != 0) return G4ThreeVector(nx,ny,nz).unit(); // edge or corner
776 else
777 {
778 // Point is not on the surface
779 //
780#ifdef G4CSGDEBUG
781 std::ostringstream message;
782 G4int oldprc = message.precision(16);
783 message << "Point p is not on surface (!?) of solid: "
784 << GetName() << G4endl;
785 message << "Position:\n";
786 message << " p.x() = " << p.x()/mm << " mm\n";
787 message << " p.y() = " << p.y()/mm << " mm\n";
788 message << " p.z() = " << p.z()/mm << " mm";
789 G4cout.precision(oldprc) ;
790 G4Exception("G4Trap::SurfaceNormal(p)", "GeomSolids1002",
791 JustWarning, message );
792 DumpInfo();
793#endif
794 return ApproxSurfaceNormal(p);
795 }
796}
Hep3Vector unit() const
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Trap.cc:803

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

Field Documentation

◆ fAreas

G4double G4Trap::fAreas[6]
private

Definition at line 282 of file G4Trap.hh.

Referenced by G4Trap(), GetPointOnSurface(), operator=(), and SetCachedValues().

◆ fCubicVolume

G4double G4CSGSolid::fCubicVolume = 0.0
protectedinherited

◆ fDx1

G4double G4Trap::fDx1
private

◆ fDx2

G4double G4Trap::fDx2
private

◆ fDx3

G4double G4Trap::fDx3
private

◆ fDx4

G4double G4Trap::fDx4
private

◆ fDy1

G4double G4Trap::fDy1
private

◆ fDy2

G4double G4Trap::fDy2
private

◆ fDz

G4double G4Trap::fDz
private

◆ fPlanes

TrapSidePlane G4Trap::fPlanes[4]
private

◆ 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

◆ fTalpha1

G4double G4Trap::fTalpha1
private

Definition at line 279 of file G4Trap.hh.

Referenced by CreatePolyhedron(), G4Trap(), MakePlanes(), operator=(), and SetAllParameters().

◆ fTalpha2

G4double G4Trap::fTalpha2
private

Definition at line 280 of file G4Trap.hh.

Referenced by CreatePolyhedron(), G4Trap(), MakePlanes(), operator=(), and SetAllParameters().

◆ fTrapType

G4int G4Trap::fTrapType
private

◆ fTthetaCphi

G4double G4Trap::fTthetaCphi
private

Definition at line 278 of file G4Trap.hh.

Referenced by CreatePolyhedron(), G4Trap(), MakePlanes(), operator=(), and SetAllParameters().

◆ fTthetaSphi

G4double G4Trap::fTthetaSphi
private

Definition at line 278 of file G4Trap.hh.

Referenced by CreatePolyhedron(), G4Trap(), MakePlanes(), operator=(), and SetAllParameters().

◆ halfCarTolerance

G4double G4Trap::halfCarTolerance
private

Definition at line 277 of file G4Trap.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(), G4Trd::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(), 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(), 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: