Geant4 11.2.2
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)
 
 ~G4Trap () override
 
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
 
G4double GetPhi () const
 
G4double GetTheta () const
 
G4double GetAlpha1 () const
 
G4double GetAlpha2 () 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 () override
 
G4double GetSurfaceArea () override
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
 
void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const override
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const override
 
EInside Inside (const G4ThreeVector &p) const override
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const override
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const override
 
G4double DistanceToIn (const G4ThreeVector &p) const override
 
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
 
G4double DistanceToOut (const G4ThreeVector &p) const override
 
G4GeometryType GetEntityType () const override
 
G4ThreeVector GetPointOnSurface () const override
 
G4VSolidClone () const override
 
std::ostream & StreamInfo (std::ostream &os) const override
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const override
 
G4PolyhedronCreatePolyhedron () const override
 
 G4Trap (__void__ &)
 
 G4Trap (const G4Trap &rhs)
 
G4Trapoperator= (const G4Trap &rhs)
 
- Public Member Functions inherited from G4CSGSolid
 G4CSGSolid (const G4String &pName)
 
virtual ~G4CSGSolid ()
 
G4PolyhedronGetPolyhedron () const override
 
 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
 
void DumpInfo () const
 
virtual G4VisExtent GetExtent () 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}
G4CSGSolid(const G4String &pName)
Definition G4CSGSolid.cc:49
void MakePlanes()
Definition G4Trap.cc:333
G4double kCarTolerance
Definition G4VSolid.hh:299

Referenced by Clone().

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

◆ 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 238 of file G4Trap.cc.

239 : G4CSGSolid(rhs), halfCarTolerance(rhs.halfCarTolerance),
240 fDz(rhs.fDz), fTthetaCphi(rhs.fTthetaCphi), fTthetaSphi(rhs.fTthetaSphi),
241 fDy1(rhs.fDy1), fDx1(rhs.fDx1), fDx2(rhs.fDx2), fTalpha1(rhs.fTalpha1),
242 fDy2(rhs.fDy2), fDx3(rhs.fDx3), fDx4(rhs.fDx4), fTalpha2(rhs.fTalpha2)
243{
244 for (G4int i=0; i<4; ++i) { fPlanes[i] = rhs.fPlanes[i]; }
245 for (G4int i=0; i<6; ++i) { fAreas[i] = rhs.fAreas[i]; }
246 fTrapType = rhs.fTrapType;
247}
int G4int
Definition G4Types.hh:85

Member Function Documentation

◆ BoundingLimits()

void G4Trap::BoundingLimits ( G4ThreeVector & pMin,
G4ThreeVector & pMax ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 547 of file G4Trap.cc.

548{
549 G4ThreeVector pt[8];
550 GetVertices(pt);
551
552 G4double xmin = kInfinity, xmax = -kInfinity;
553 G4double ymin = kInfinity, ymax = -kInfinity;
554 for (const auto & i : pt)
555 {
556 G4double x = i.x();
557 if (x < xmin) xmin = x;
558 if (x > xmax) xmax = x;
559 G4double y = i.y();
560 if (y < ymin) ymin = y;
561 if (y > ymax) ymax = y;
562 }
563
565 pMin.set(xmin,ymin,-dz);
566 pMax.set(xmax,ymax, dz);
567
568 // Check correctness of the bounding box
569 //
570 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
571 {
572 std::ostringstream message;
573 message << "Bad bounding box (min >= max) for solid: "
574 << GetName() << " !"
575 << "\npMin = " << pMin
576 << "\npMax = " << pMax;
577 G4Exception("G4Trap::BoundingLimits()", "GeomMgt0001",
578 JustWarning, message);
579 DumpInfo();
580 }
581}
@ 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
overridevirtual

Implements G4VSolid.

Definition at line 587 of file G4Trap.cc.

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

◆ Clone()

G4VSolid * G4Trap::Clone ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 1114 of file G4Trap.cc.

1115{
1116 return new G4Trap(*this);
1117}
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

◆ ComputeDimensions()

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

Reimplemented from G4VSolid.

Definition at line 536 of file G4Trap.cc.

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

◆ CreatePolyhedron()

G4Polyhedron * G4Trap::CreatePolyhedron ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 1221 of file G4Trap.cc.

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

◆ DescribeYourselfTo()

void G4Trap::DescribeYourselfTo ( G4VGraphicsScene & scene) const
overridevirtual

Implements G4VSolid.

Definition at line 1216 of file G4Trap.cc.

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

◆ DistanceToIn() [1/2]

G4double G4Trap::DistanceToIn ( const G4ThreeVector & p) const
overridevirtual

Implements G4VSolid.

Definition at line 894 of file G4Trap.cc.

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

Implements G4VSolid.

Definition at line 825 of file G4Trap.cc.

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

◆ DistanceToOut() [1/2]

G4double G4Trap::DistanceToOut ( const G4ThreeVector & p) const
overridevirtual

Implements G4VSolid.

Definition at line 1034 of file G4Trap.cc.

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

Implements G4VSolid.

Definition at line 951 of file G4Trap.cc.

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

◆ GetAlpha1()

G4double G4Trap::GetAlpha1 ( ) const
inline

Referenced by StreamInfo().

◆ GetAlpha2()

G4double G4Trap::GetAlpha2 ( ) const
inline

Referenced by StreamInfo().

◆ GetCubicVolume()

G4double G4Trap::GetCubicVolume ( )
overridevirtual

Reimplemented from G4VSolid.

Definition at line 486 of file G4Trap.cc.

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

◆ GetEntityType()

G4GeometryType G4Trap::GetEntityType ( ) const
overridevirtual

Implements G4VSolid.

Definition at line 1105 of file G4Trap.cc.

1106{
1107 return {"G4Trap"};
1108}

◆ GetPhi()

G4double G4Trap::GetPhi ( ) const
inline

Referenced by StreamInfo().

◆ GetPointOnSurface()

G4ThreeVector G4Trap::GetPointOnSurface ( ) const
overridevirtual

Reimplemented from G4VSolid.

Definition at line 1175 of file G4Trap.cc.

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

◆ GetSidePlane()

TrapSidePlane G4Trap::GetSidePlane ( G4int n) const
inline

◆ GetSurfaceArea()

G4double G4Trap::GetSurfaceArea ( )
overridevirtual

Reimplemented from G4VSolid.

Definition at line 511 of file G4Trap.cc.

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

◆ GetSymAxis()

◆ GetTanAlpha1()

◆ GetTanAlpha2()

◆ GetTheta()

G4double G4Trap::GetTheta ( ) const
inline

Referenced by StreamInfo().

◆ GetXHalfLength1()

◆ GetXHalfLength2()

◆ GetXHalfLength3()

◆ GetXHalfLength4()

◆ GetYHalfLength1()

◆ GetYHalfLength2()

◆ GetZHalfLength()

◆ Inside()

EInside G4Trap::Inside ( const G4ThreeVector & p) const
overridevirtual

Implements G4VSolid.

Definition at line 636 of file G4Trap.cc.

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

405{
406 G4ThreeVector normal = ((p4 - p2).cross(p3 - p1)).unit();
407 if (std::abs(normal.x()) < DBL_EPSILON) normal.setX(0);
408 if (std::abs(normal.y()) < DBL_EPSILON) normal.setY(0);
409 if (std::abs(normal.z()) < DBL_EPSILON) normal.setZ(0);
410 normal = normal.unit();
411
412 G4ThreeVector centre = (p1 + p2 + p3 + p4)*0.25;
413 plane.a = normal.x();
414 plane.b = normal.y();
415 plane.c = normal.z();
416 plane.d = -normal.dot(centre);
417
418 // compute distances and check planarity
419 G4double d1 = std::abs(normal.dot(p1) + plane.d);
420 G4double d2 = std::abs(normal.dot(p2) + plane.d);
421 G4double d3 = std::abs(normal.dot(p3) + plane.d);
422 G4double d4 = std::abs(normal.dot(p4) + plane.d);
423 G4double dmax = std::max(std::max(std::max(d1,d2),d3),d4);
424
425 return dmax <= 1000 * kCarTolerance;
426}
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 333 of file G4Trap.cc.

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

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

◆ MakePlanes() [2/2]

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

Definition at line 359 of file G4Trap.cc.

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

◆ operator=()

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

Definition at line 253 of file G4Trap.cc.

254{
255 // Check assignment to self
256 //
257 if (this == &rhs) { return *this; }
258
259 // Copy base class data
260 //
262
263 // Copy data
264 //
265 halfCarTolerance = rhs.halfCarTolerance;
266 fDz = rhs.fDz; fTthetaCphi = rhs.fTthetaCphi; fTthetaSphi = rhs.fTthetaSphi;
267 fDy1 = rhs.fDy1; fDx1 = rhs.fDx1; fDx2 = rhs.fDx2; fTalpha1 = rhs.fTalpha1;
268 fDy2 = rhs.fDy2; fDx3 = rhs.fDx3; fDx4 = rhs.fDx4; fTalpha2 = rhs.fTalpha2;
269 for (G4int i=0; i<4; ++i) { fPlanes[i] = rhs.fPlanes[i]; }
270 for (G4int i=0; i<6; ++i) { fAreas[i] = rhs.fAreas[i]; }
271 fTrapType = rhs.fTrapType;
272 return *this;
273}
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 280 of file G4Trap.cc.

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

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

◆ SetCachedValues()

void G4Trap::SetCachedValues ( )
protected

Definition at line 432 of file G4Trap.cc.

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

Referenced by MakePlanes().

◆ StreamInfo()

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

Reimplemented from G4CSGSolid.

Definition at line 1123 of file G4Trap.cc.

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

Referenced by MakePlanes().

◆ SurfaceNormal()

G4ThreeVector G4Trap::SurfaceNormal ( const G4ThreeVector & p) const
overridevirtual

Implements G4VSolid.

Definition at line 699 of file G4Trap.cc.

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

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