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

#include <G4Trap.hh>

+ Inheritance diagram for G4Trap:

Public Member Functions

 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, const G4ThreeVector pt[8])
 
 G4Trap (const G4String &pName, G4double pZ, G4double pY, G4double pX, G4double pLTX)
 
 G4Trap (const G4String &pName, G4double pDx1, G4double pDx2, G4double pDy1, G4double pDy2, G4double pDz)
 
 G4Trap (const G4String &pName, G4double pDx, G4double pDy, G4double pDz, G4double pAlpha, G4double pTheta, G4double pPhi)
 
 G4Trap (const G4String &pName)
 
virtual ~G4Trap ()
 
G4double GetZHalfLength () const
 
G4double GetYHalfLength1 () const
 
G4double GetXHalfLength1 () const
 
G4double GetXHalfLength2 () const
 
G4double GetTanAlpha1 () const
 
G4double GetYHalfLength2 () const
 
G4double GetXHalfLength3 () const
 
G4double GetXHalfLength4 () const
 
G4double GetTanAlpha2 () const
 
TrapSidePlane GetSidePlane (G4int n) const
 
G4ThreeVector GetSymAxis () 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)
 
G4double GetCubicVolume ()
 
G4double GetSurfaceArea ()
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
 
EInside Inside (const G4ThreeVector &p) const
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
G4double DistanceToIn (const G4ThreeVector &p) const
 
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
 
G4double DistanceToOut (const G4ThreeVector &p) const
 
G4GeometryType GetEntityType () const
 
G4ThreeVector GetPointOnSurface () const
 
G4VSolidClone () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
G4PolyhedronCreatePolyhedron () const
 
 G4Trap (__void__ &)
 
 G4Trap (const G4Trap &rhs)
 
G4Trapoperator= (const G4Trap &rhs)
 
- Public Member Functions inherited from G4CSGSolid
 G4CSGSolid (const G4String &pName)
 
virtual ~G4CSGSolid ()
 
virtual std::ostream & StreamInfo (std::ostream &os) const
 
virtual G4PolyhedronGetPolyhedron () const
 
 G4CSGSolid (__void__ &)
 
 G4CSGSolid (const G4CSGSolid &rhs)
 
G4CSGSolidoperator= (const G4CSGSolid &rhs)
 
- Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
 
virtual ~G4VSolid ()
 
G4bool operator== (const G4VSolid &s) const
 
G4String GetName () const
 
void SetName (const G4String &name)
 
G4double GetTolerance () const
 
virtual void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const
 
virtual G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
 
virtual EInside Inside (const G4ThreeVector &p) const =0
 
virtual G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p) const =0
 
virtual void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
virtual G4double GetCubicVolume ()
 
virtual G4double GetSurfaceArea ()
 
virtual G4GeometryType GetEntityType () const =0
 
virtual G4ThreeVector GetPointOnSurface () const
 
virtual G4VSolidClone () const
 
virtual std::ostream & StreamInfo (std::ostream &os) const =0
 
void DumpInfo () const
 
virtual void DescribeYourselfTo (G4VGraphicsScene &scene) const =0
 
virtual G4VisExtent GetExtent () const
 
virtual G4PolyhedronCreatePolyhedron () const
 
virtual G4PolyhedronGetPolyhedron () const
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
 G4VSolid (__void__ &)
 
 G4VSolid (const G4VSolid &rhs)
 
G4VSolidoperator= (const G4VSolid &rhs)
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 

Protected Member Functions

void MakePlanes ()
 
void MakePlanes (const G4ThreeVector pt[8])
 
G4bool MakePlane (const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3, const G4ThreeVector &p4, TrapSidePlane &plane)
 
void SetCachedValues ()
 
- Protected Member Functions inherited from G4CSGSolid
G4double GetRadiusInRing (G4double rmin, G4double rmax) const
 
- Protected Member Functions inherited from G4VSolid
void CalculateClippedPolygonExtent (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipCrossSection (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipBetweenSections (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipPolygon (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
 

Additional Inherited Members

- Protected Attributes inherited from G4CSGSolid
G4double fCubicVolume = 0.0
 
G4double fSurfaceArea = 0.0
 
G4bool fRebuildPolyhedron = false
 
G4PolyhedronfpPolyhedron = nullptr
 
- Protected Attributes inherited from G4VSolid
G4double kCarTolerance
 

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.

67 : G4CSGSolid(pName), halfCarTolerance(0.5*kCarTolerance)
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
76 CheckParameters();
77 MakePlanes();
78}
void MakePlanes()
Definition: G4Trap.cc:335
G4double kCarTolerance
Definition: G4VSolid.hh:302

◆ G4Trap() [2/8]

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

Definition at line 86 of file G4Trap.cc.

88 : G4CSGSolid(pName), halfCarTolerance(0.5*kCarTolerance)
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
137 CheckParameters();
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

◆ G4Trap() [3/8]

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

Definition at line 145 of file G4Trap.cc.

149 : G4CSGSolid(pName), halfCarTolerance(0.5*kCarTolerance)
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;
153 fDy2 = fDy1; fDx3 = fDx1; fDx4 = fDx2; fTalpha2 = fTalpha1;
154
155 CheckParameters();
156 MakePlanes();
157}

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

167 : G4CSGSolid(pName), halfCarTolerance(0.5*kCarTolerance), fTrapType(0)
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
173 CheckParameters();
174 MakePlanes();
175}

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

186 : G4CSGSolid(pName), halfCarTolerance(0.5*kCarTolerance)
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
195 CheckParameters();
196 MakePlanes();
197}

◆ G4Trap() [6/8]

G4Trap::G4Trap ( const G4String pName)

Definition at line 205 of file G4Trap.cc.

206 : G4CSGSolid (pName), halfCarTolerance(0.5*kCarTolerance),
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}

◆ ~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.

220 : G4CSGSolid(a), halfCarTolerance(0.5*kCarTolerance),
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}

◆ G4Trap() [8/8]

G4Trap::G4Trap ( const G4Trap rhs)

Definition at line 240 of file G4Trap.cc.

241 : G4CSGSolid(rhs), halfCarTolerance(rhs.halfCarTolerance),
242 fDz(rhs.fDz), fTthetaCphi(rhs.fTthetaCphi), fTthetaSphi(rhs.fTthetaSphi),
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

Member Function Documentation

◆ 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
double G4double
Definition: G4Types.hh:83
double z() const
double x() const
double y() const
void set(double x, double y, double z)
G4double GetZHalfLength() const
void DumpInfo() const

Referenced by CalculateExtent().

◆ 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

◆ Clone()

G4VSolid * G4Trap::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1116 of file G4Trap.cc.

1117{
1118 return new G4Trap(*this);
1119}

◆ 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

◆ CreatePolyhedron()

G4Polyhedron * G4Trap::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1225 of file G4Trap.cc.

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

◆ DescribeYourselfTo()

void G4Trap::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 1220 of file G4Trap.cc.

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

◆ 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}
G4double b
Definition: G4Trap.hh:92
G4double c
Definition: G4Trap.hh:92
G4double d
Definition: G4Trap.hh:92
G4double a
Definition: G4Trap.hh:92

◆ 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}
#define DBL_MAX
Definition: templates.hh:62

◆ 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}
#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

◆ 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}

◆ 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

◆ GetEntityType()

G4GeometryType G4Trap::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 1107 of file G4Trap.cc.

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

◆ GetPointOnSurface()

G4ThreeVector G4Trap::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1179 of file G4Trap.cc.

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

◆ GetSidePlane()

TrapSidePlane G4Trap::GetSidePlane ( G4int  n) const
inline

◆ 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

◆ GetSymAxis()

◆ GetTanAlpha1()

◆ GetTanAlpha2()

◆ GetXHalfLength1()

◆ GetXHalfLength2()

◆ GetXHalfLength3()

◆ GetXHalfLength4()

◆ GetYHalfLength1()

◆ GetYHalfLength2()

◆ GetZHalfLength()

◆ 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}
@ kInside
Definition: geomdefs.hh:70
@ kSurface
Definition: geomdefs.hh:69

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}
Hep3Vector unit() const
void setY(double)
double dot(const Hep3Vector &) const
void setZ(double)
void setX(double)
#define DBL_EPSILON
Definition: templates.hh:66

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}
CLHEP::Hep3Vector G4ThreeVector

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
375 G4ThreeVector normal(fPlanes[i].a,fPlanes[i].b,fPlanes[i].c);
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

◆ 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 //
267 halfCarTolerance = rhs.halfCarTolerance;
268 fDz = rhs.fDz; fTthetaCphi = rhs.fTthetaCphi; fTthetaSphi = rhs.fTthetaSphi;
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

◆ 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
307 CheckParameters();
308 MakePlanes();
309}
G4bool fRebuildPolyhedron
Definition: G4CSGSolid.hh:72

Referenced by G4ParameterisationTrdX::ComputeDimensions().

◆ 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}

Referenced by MakePlanes().

◆ 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 = std::atan2(fTthetaSphi,fTthetaCphi);
1128 G4double theta = std::atan(std::sqrt(fTthetaCphi*fTthetaCphi
1129 +fTthetaSphi*fTthetaSphi));
1130 G4double alpha1 = std::atan(fTalpha1);
1131 G4double alpha2 = std::atan(fTalpha2);
1132 G4String signDegree = "\u00B0";
1133
1134 G4int oldprc = os.precision(16);
1135 os << "-----------------------------------------------------------\n"
1136 << " *** Dump for solid: " << GetName() << " ***\n"
1137 << " ===================================================\n"
1138 << " Solid type: G4Trap\n"
1139 << " Parameters:\n"
1140 << " half length Z: " << fDz/mm << " mm\n"
1141 << " half length Y, face -Dz: " << fDy1/mm << " mm\n"
1142 << " half length X, face -Dz, side -Dy1: " << fDx1/mm << " mm\n"
1143 << " half length X, face -Dz, side +Dy1: " << fDx2/mm << " mm\n"
1144 << " half length Y, face +Dz: " << fDy2/mm << " mm\n"
1145 << " half length X, face +Dz, side -Dy2: " << fDx3/mm << " mm\n"
1146 << " half length X, face +Dz, side +Dy2: " << fDx4/mm << " mm\n"
1147 << " theta: " << theta/degree << signDegree << "\n"
1148 << " phi: " << phi/degree << signDegree << "\n"
1149 << " alpha, face -Dz: " << alpha1/degree << signDegree << "\n"
1150 << " alpha, face +Dz: " << alpha2/degree << signDegree << "\n"
1151 << "-----------------------------------------------------------\n";
1152 os.precision(oldprc);
1153
1154 return os;
1155}

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}

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