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

G4Trap is a general trapezoid: the faces perpendicular to the Z planes are trapezia, and their centres are not necessarily on a line parallel to the Z axis. A check for planarity is made in the calculation of the equation for each plane. If the planes are not parallel, a call to G4Exception is made. More...

#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=default
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
G4bool IsFaceted () 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)
 ~G4CSGSolid () override
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 ()
 G4VSolid (const G4VSolid &rhs)
G4VSolidoperator= (const G4VSolid &rhs)
G4bool operator== (const G4VSolid &s) const
G4String GetName () const
void SetName (const G4String &name)
G4double GetTolerance () const
virtual G4int GetNumOfConstituents () 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__ &)
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
G4double EstimateSurfaceArea (G4int nStat, G4double epsilon) 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

G4Trap is a general trapezoid: the faces perpendicular to the Z planes are trapezia, and their centres are not necessarily on a line parallel to the Z axis. A check for planarity is made in the calculation of the equation for each plane. If the planes are not parallel, a call to G4Exception is made.

Definition at line 115 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 )

The most general constructor for G4Trap which prepares plane equations and corner coordinates from parameters.

Parameters
[in]pNameThe name of the solid.
[in]pDzHalf-length along the Z-axis.
[in]pThetaPolar angle of the line joining the centres of the faces at -/+pDz.
[in]pPhiAzimuthal angle of the line joining the centre of the face at -pDz to the centre of the face at +pDz.
[in]pDy1Half-length along Y of the face at -pDz.
[in]pDx1Half-length along X of the side at y=-pDy1 of the face at -pDz.
[in]pDx2Half-length along X of the side at y=+pDy1 of the face at -pDz.
[in]pAlp1Angle with respect to the Y axis from the centre of the side at y=-pDy1 to the centre at y=+pDy1 of the face at -pDz.
[in]pDy2Half-length along Y of the face at +pDz.
[in]pDx3Half-length along X of the side at y=-pDy2 of the face at +pDz.
[in]pDx4Half-length along X of the side at y=+pDy2 of the face at +pDz.
[in]pAlp2Angle with respect to the Y axis from the centre of the side at y=-pDy2 to the centre at y=+pDy2 of the face at +pDz.

Definition at line 66 of file G4Trap.cc.

73 : G4CSGSolid(pName), halfCarTolerance(0.5*kCarTolerance)
74{
75 fDz = pDz;
76 fTthetaCphi = std::tan(pTheta)*std::cos(pPhi);
77 fTthetaSphi = std::tan(pTheta)*std::sin(pPhi);
78
79 fDy1 = pDy1; fDx1 = pDx1; fDx2 = pDx2; fTalpha1 = std::tan(pAlp1);
80 fDy2 = pDy2; fDx3 = pDx3; fDx4 = pDx4; fTalpha2 = std::tan(pAlp2);
81
82 CheckParameters();
83 MakePlanes();
84}
G4CSGSolid(const G4String &pName)
Definition G4CSGSolid.cc:49
void MakePlanes()
Definition G4Trap.cc:333
G4double kCarTolerance
Definition G4VSolid.hh:418

Referenced by Clone(), G4Trap(), GetAlpha2(), and operator=().

◆ G4Trap() [2/8]

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

Prepares plane equations and parameters from corner coordinates.

Parameters
[in]pNameThe name of the solid.
[in]ptPoints of the 8 vertices.

Definition at line 92 of file G4Trap.cc.

94 : G4CSGSolid(pName), halfCarTolerance(0.5*kCarTolerance)
95{
96 // Start with check of centering - the center of gravity trap line
97 // should cross the origin of frame
98 //
99 if ( pt[0].z() >= 0
100 || pt[0].z() != pt[1].z()
101 || pt[0].z() != pt[2].z()
102 || pt[0].z() != pt[3].z()
103
104 || pt[4].z() <= 0
105 || pt[4].z() != pt[5].z()
106 || pt[4].z() != pt[6].z()
107 || pt[4].z() != pt[7].z()
108
109 || std::fabs( pt[0].z() + pt[4].z() ) >= kCarTolerance
110
111 || pt[0].y() != pt[1].y()
112 || pt[2].y() != pt[3].y()
113 || pt[4].y() != pt[5].y()
114 || pt[6].y() != pt[7].y()
115
116 || std::fabs(pt[0].y()+pt[2].y()+pt[4].y()+pt[6].y()) >= kCarTolerance
117 || std::fabs(pt[0].x()+pt[1].x()+pt[4].x()+pt[5].x() +
118 pt[2].x()+pt[3].x()+pt[6].x()+pt[7].x()) >= kCarTolerance )
119 {
120 std::ostringstream message;
121 message << "Invalid vertice coordinates for Solid: " << GetName();
122 G4Exception("G4Trap::G4Trap()", "GeomSolids0002",
123 FatalException, message);
124 }
125
126 // Set parameters
127 //
128 fDz = (pt[7]).z();
129
130 fDy1 = ((pt[2]).y()-(pt[1]).y())*0.5;
131 fDx1 = ((pt[1]).x()-(pt[0]).x())*0.5;
132 fDx2 = ((pt[3]).x()-(pt[2]).x())*0.5;
133 fTalpha1 = ((pt[2]).x()+(pt[3]).x()-(pt[1]).x()-(pt[0]).x())*0.25/fDy1;
134
135 fDy2 = ((pt[6]).y()-(pt[5]).y())*0.5;
136 fDx3 = ((pt[5]).x()-(pt[4]).x())*0.5;
137 fDx4 = ((pt[7]).x()-(pt[6]).x())*0.5;
138 fTalpha2 = ((pt[6]).x()+(pt[7]).x()-(pt[5]).x()-(pt[4]).x())*0.25/fDy2;
139
140 fTthetaCphi = ((pt[4]).x()+fDy2*fTalpha2+fDx3)/fDz;
141 fTthetaSphi = ((pt[4]).y()+fDy2)/fDz;
142
143 CheckParameters();
144 MakePlanes(pt);
145}
@ 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 )

Constructor for Right Angular Wedge from STEP (assumes pLTX<=pX).

Parameters
[in]pNameThe name of the solid.
[in]pZLength along Z.
[in]pYLength along Y.
[in]pXLength along X at the wider side.
[in]pLTXLength along X at the narrower side (plTX<=pX).

Definition at line 151 of file G4Trap.cc.

155 : G4CSGSolid(pName), halfCarTolerance(0.5*kCarTolerance)
156{
157 fDz = 0.5*pZ; fTthetaCphi = 0; fTthetaSphi = 0;
158 fDy1 = 0.5*pY; fDx1 = 0.5*pX; fDx2 = 0.5*pLTX; fTalpha1 = 0.5*(pLTX - pX)/pY;
159 fDy2 = fDy1; fDx3 = fDx1; fDx4 = fDx2; fTalpha2 = fTalpha1;
160
161 CheckParameters();
162 MakePlanes();
163}

◆ G4Trap() [4/8]

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

Constructor for G4Trd.

Parameters
[in]pNameThe name of the solid.
[in]pDx1Half-length along X at the surface positioned at -dz.
[in]pDx2Half-length along X at the surface positioned at +dz.
[in]pDy1Half-length along Y at the surface positioned at -dz.
[in]pDy2Half-length along Y at the surface positioned at +dz.
[in]pDzHalf-length along Z axis.

Definition at line 169 of file G4Trap.cc.

173 : G4CSGSolid(pName), halfCarTolerance(0.5*kCarTolerance), fTrapType(0)
174{
175 fDz = pDz; fTthetaCphi = 0; fTthetaSphi = 0;
176 fDy1 = pDy1; fDx1 = pDx1; fDx2 = pDx1; fTalpha1 = 0;
177 fDy2 = pDy2; fDx3 = pDx2; fDx4 = pDx2; fTalpha2 = 0;
178
179 CheckParameters();
180 MakePlanes();
181}

◆ G4Trap() [5/8]

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

Constructor for G4Para.

Parameters
[in]pNameThe name of the solid.
[in]pDxHalf-length in X.
[in]pDyHalf-length in Y.
[in]pDzHalf-length in Z.
[in]pAlphaAngle formed by the Y axis and the plane joining the centre of the faces parallel to the Z-X plane at -dy and +dy.
[in]pThetaPolar angle of the line joining the centres of the faces at -dz and +dz in Z.
[in]pPhiAzimuthal angle of the line joining the centres of the faces at -dz and +dz in Z.

Definition at line 187 of file G4Trap.cc.

192 : G4CSGSolid(pName), halfCarTolerance(0.5*kCarTolerance)
193{
194 fDz = pDz;
195 fTthetaCphi = std::tan(pTheta)*std::cos(pPhi);
196 fTthetaSphi = std::tan(pTheta)*std::sin(pPhi);
197
198 fDy1 = pDy; fDx1 = pDx; fDx2 = pDx; fTalpha1 = std::tan(pAlpha);
199 fDy2 = pDy; fDx3 = pDx; fDx4 = pDx; fTalpha2 = fTalpha1;
200
201 CheckParameters();
202 MakePlanes();
203}

◆ G4Trap() [6/8]

G4Trap::G4Trap ( const G4String & pName)

Constructor for "nominal" G4Trap whose parameters are to be set by a G4VPVParamaterisation later on.

Parameters
[in]pNameThe name of the solid.

Definition at line 211 of file G4Trap.cc.

212 : G4CSGSolid (pName), halfCarTolerance(0.5*kCarTolerance),
213 fDz(1.), fTthetaCphi(0.), fTthetaSphi(0.),
214 fDy1(1.), fDx1(1.), fDx2(1.), fTalpha1(0.),
215 fDy2(1.), fDx3(1.), fDx4(1.), fTalpha2(0.)
216{
217 MakePlanes();
218}

◆ ~G4Trap()

G4Trap::~G4Trap ( )
overridedefault

Default destructor.

◆ G4Trap() [7/8]

G4Trap::G4Trap ( __void__ & a)

Fake default constructor for usage restricted to direct object persistency for clients requiring preallocation of memory for persistifiable objects.

Definition at line 225 of file G4Trap.cc.

226 : G4CSGSolid(a), halfCarTolerance(0.5*kCarTolerance),
227 fDz(1.), fTthetaCphi(0.), fTthetaSphi(0.),
228 fDy1(1.), fDx1(1.), fDx2(1.), fTalpha1(0.),
229 fDy2(1.), fDx3(1.), fDx4(1.), fTalpha2(0.)
230{
231 MakePlanes();
232}

◆ G4Trap() [8/8]

G4Trap::G4Trap ( const G4Trap & rhs)

Copy constructor and assignment operator.

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

Computes the bounding limits of the solid.

Parameters
[out]pMinThe minimum bounding limit point.
[out]pMaxThe maximum bounding limit point.

Reimplemented from G4VSolid.

Definition at line 498 of file G4Trap.cc.

499{
500 G4ThreeVector pt[8];
501 GetVertices(pt);
502
503 G4double xmin = kInfinity, xmax = -kInfinity;
504 G4double ymin = kInfinity, ymax = -kInfinity;
505 for (const auto & i : pt)
506 {
507 G4double x = i.x();
508 if (x < xmin) { xmin = x; }
509 if (x > xmax) { xmax = x; }
510 G4double y = i.y();
511 if (y < ymin) { ymin = y; }
512 if (y > ymax) { ymax = y; }
513 }
514
516 pMin.set(xmin,ymin,-dz);
517 pMax.set(xmax,ymax, dz);
518
519 // Check correctness of the bounding box
520 //
521 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
522 {
523 std::ostringstream message;
524 message << "Bad bounding box (min >= max) for solid: "
525 << GetName() << " !"
526 << "\npMin = " << pMin
527 << "\npMax = " << pMax;
528 G4Exception("G4Trap::BoundingLimits()", "GeomMgt0001",
529 JustWarning, message);
530 DumpInfo();
531 }
532}
@ JustWarning
CLHEP::Hep3Vector G4ThreeVector
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

Calculates the minimum and maximum extent of the solid, when under the specified transform, and within the specified limits.

Parameters
[in]pAxisThe axis along which compute the extent.
[in]pVoxelLimitThe limiting space dictated by voxels.
[in]pTransformThe internal transformation applied to the solid.
[out]pMinThe minimum extent value.
[out]pMaxThe maximum extent value.
Returns
True if the solid is intersected by the extent region.

Implements G4VSolid.

Definition at line 538 of file G4Trap.cc.

542{
543 G4ThreeVector bmin, bmax;
544 G4bool exist;
545
546 // Check bounding box (bbox)
547 //
548 BoundingLimits(bmin,bmax);
549 G4BoundingEnvelope bbox(bmin,bmax);
550#ifdef G4BBOX_EXTENT
551 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
552#endif
553 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
554 {
555 return exist = pMin < pMax;
556 }
557
558 // Set bounding envelope (benv) and calculate extent
559 //
560 G4ThreeVector pt[8];
561 GetVertices(pt);
562
563 G4ThreeVectorList baseA(4), baseB(4);
564 baseA[0] = pt[0];
565 baseA[1] = pt[1];
566 baseA[2] = pt[3];
567 baseA[3] = pt[2];
568
569 baseB[0] = pt[4];
570 baseB[1] = pt[5];
571 baseB[2] = pt[7];
572 baseB[3] = pt[6];
573
574 std::vector<const G4ThreeVectorList *> polygons(2);
575 polygons[0] = &baseA;
576 polygons[1] = &baseB;
577
578 G4BoundingEnvelope benv(bmin,bmax,polygons);
579 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
580 return exist;
581}
std::vector< G4ThreeVector > G4ThreeVectorList
bool G4bool
Definition G4Types.hh:86
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
Definition G4Trap.cc:498

◆ Clone()

G4VSolid * G4Trap::Clone ( ) const
overridevirtual

Makes a clone of the object for use in multi-treading.

Returns
A pointer to the new cloned allocated solid.

Reimplemented from G4VSolid.

Definition at line 1088 of file G4Trap.cc.

1089{
1090 return new G4Trap(*this);
1091}
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:66

◆ ComputeDimensions()

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

Dispatch method for parameterisation replication mechanism and dimension computation.

Reimplemented from G4VSolid.

Definition at line 487 of file G4Trap.cc.

490{
491 p->ComputeDimensions(*this,n,pRep);
492}
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const

◆ CreatePolyhedron()

G4Polyhedron * G4Trap::CreatePolyhedron ( ) const
overridevirtual

Creates a Polyhedron used for Visualisation. It is the caller's responsibility to delete it. A null pointer means "not created".

Reimplemented from G4VSolid.

Definition at line 1248 of file G4Trap.cc.

1249{
1250 G4double phi = std::atan2(fTthetaSphi, fTthetaCphi);
1251 G4double alpha1 = std::atan(fTalpha1);
1252 G4double alpha2 = std::atan(fTalpha2);
1253 G4double theta = std::atan(std::sqrt(fTthetaCphi*fTthetaCphi
1254 +fTthetaSphi*fTthetaSphi));
1255
1256 return new G4PolyhedronTrap(fDz, theta, phi,
1257 fDy1, fDx1, fDx2, alpha1,
1258 fDy2, fDx3, fDx4, alpha2);
1259}
const G4double alpha2

◆ DescribeYourselfTo()

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

Methods for creating graphical representations (i.e. for visualisation).

Implements G4VSolid.

Definition at line 1243 of file G4Trap.cc.

1244{
1245 scene.AddSolid (*this);
1246}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

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

Calculates the distance to the nearest surface of a shape from an outside point. The distance can be an underestimate.

Parameters
[in]pThe point at offset p.
Returns
The safety distance to enter the shape.

Implements G4VSolid.

Definition at line 855 of file G4Trap.cc.

856{
857 switch (fTrapType)
858 {
859 case 0: // General case
860 {
861 G4double dz = std::abs(p.z())-fDz;
862 G4double dy1 = fPlanes[0].b*p.y()+fPlanes[0].c*p.z()+fPlanes[0].d;
863 G4double dy2 = fPlanes[1].b*p.y()+fPlanes[1].c*p.z()+fPlanes[1].d;
864 G4double dy = std::max(dz,std::max(dy1,dy2));
865
866 G4double dx1 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()
867 + fPlanes[2].c*p.z()+fPlanes[2].d;
868 G4double dx2 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
869 + fPlanes[3].c*p.z()+fPlanes[3].d;
870 G4double dist = std::max(dy,std::max(dx1,dx2));
871 return (dist > 0) ? dist : 0.;
872 }
873 case 1: // YZ section is a rectangle
874 {
875 G4double dz = std::abs(p.z())-fDz;
876 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
877 G4double dx1 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()
878 + fPlanes[2].c*p.z()+fPlanes[2].d;
879 G4double dx2 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
880 + fPlanes[3].c*p.z()+fPlanes[3].d;
881 G4double dist = std::max(dy,std::max(dx1,dx2));
882 return (dist > 0) ? dist : 0.;
883 }
884 case 2: // YZ section is a rectangle and
885 { // XZ section is an isosceles trapezoid
886 G4double dz = std::abs(p.z())-fDz;
887 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
888 G4double dx = fPlanes[3].a*std::abs(p.x())
889 + fPlanes[3].c*p.z()+fPlanes[3].d;
890 G4double dist = std::max(dy,dx);
891 return (dist > 0) ? dist : 0.;
892 }
893 case 3: // YZ section is a rectangle and
894 { // XY section is an isosceles trapezoid
895 G4double dz = std::abs(p.z())-fDz;
896 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
897 G4double dx = fPlanes[3].a*std::abs(p.x())
898 + fPlanes[3].b*p.y()+fPlanes[3].d;
899 G4double dist = std::max(dy,dx);
900 return (dist > 0) ? dist : 0.;
901 }
902 }
903 return 0.;
904}

◆ DistanceToIn() [2/2]

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

Returns the distance along the normalised vector 'v' to the shape, from the point at offset 'p'. If there is no intersection, returns kInfinity. The first intersection resulting from 'leaving' a surface/volume is discarded. Hence, it is tolerant of points on the surface of the shape.

Parameters
[in]pThe point at offset p.
[in]vThe normalised direction vector.
Returns
The distance to enter the shape.

Implements G4VSolid.

Definition at line 783 of file G4Trap.cc.

785{
786 // Z intersections
787 //
788 if ((std::abs(p.z()) - fDz) >= -halfCarTolerance && p.z()*v.z() >= 0)
789 {
790 return kInfinity;
791 }
792 G4double invz = (-v.z() == 0) ? DBL_MAX : -1./v.z();
793 G4double dz = (invz < 0) ? fDz : -fDz;
794 G4double tzmin = (p.z() + dz)*invz;
795 G4double tzmax = (p.z() - dz)*invz;
796
797 // Y intersections
798 //
799 G4double tymin = 0, tymax = DBL_MAX;
800 G4int i = 0;
801 for ( ; i<2; ++i)
802 {
803 G4double cosa = fPlanes[i].b*v.y() + fPlanes[i].c*v.z();
804 G4double dist = fPlanes[i].b*p.y() + fPlanes[i].c*p.z() + fPlanes[i].d;
805 if (dist >= -halfCarTolerance)
806 {
807 if (cosa >= 0) { return kInfinity; }
808 G4double tmp = -dist/cosa;
809 if (tymin < tmp) { tymin = tmp; }
810 }
811 else if (cosa > 0)
812 {
813 G4double tmp = -dist/cosa;
814 if (tymax > tmp) { tymax = tmp; }
815 }
816 }
817
818 // Z intersections
819 //
820 G4double txmin = 0, txmax = DBL_MAX;
821 for ( ; i<4; ++i)
822 {
823 G4double cosa = fPlanes[i].a*v.x()+fPlanes[i].b*v.y()+fPlanes[i].c*v.z();
824 G4double dist = fPlanes[i].a*p.x()+fPlanes[i].b*p.y()+fPlanes[i].c*p.z() +
825 fPlanes[i].d;
826 if (dist >= -halfCarTolerance)
827 {
828 if (cosa >= 0) { return kInfinity; }
829 G4double tmp = -dist/cosa;
830 if (txmin < tmp) { txmin = tmp; }
831 }
832 else if (cosa > 0)
833 {
834 G4double tmp = -dist/cosa;
835 if (txmax > tmp) { txmax = tmp; }
836 }
837 }
838
839 // Find distance
840 //
841 G4double tmin = std::max(std::max(txmin,tymin),tzmin);
842 G4double tmax = std::min(std::min(txmax,tymax),tzmax);
843
844 if (tmax <= tmin + halfCarTolerance) { return kInfinity; } // touch or no hit
845
846 return (tmin < halfCarTolerance ) ? 0. : tmin;
847}
#define DBL_MAX
Definition templates.hh:62

◆ DistanceToOut() [1/2]

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

Calculates the distance to the nearest surface of a shape from an inside point 'p'. The distance can be an underestimate.

Parameters
[in]pThe point at offset p.
Returns
The safety distance to exit the shape.

Implements G4VSolid.

Definition at line 999 of file G4Trap.cc.

1000{
1001#ifdef G4CSGDEBUG
1002 if( Inside(p) == kOutside )
1003 {
1004 std::ostringstream message;
1005 G4long oldprc = message.precision(16);
1006 message << "Point p is outside (!?) of solid: " << GetName() << G4endl;
1007 message << "Position:\n";
1008 message << " p.x() = " << p.x()/mm << " mm\n";
1009 message << " p.y() = " << p.y()/mm << " mm\n";
1010 message << " p.z() = " << p.z()/mm << " mm";
1011 G4cout.precision(oldprc);
1012 G4Exception("G4Trap::DistanceToOut(p)", "GeomSolids1002",
1013 JustWarning, message );
1014 DumpInfo();
1015 }
1016#endif
1017 switch (fTrapType)
1018 {
1019 case 0: // General case
1020 {
1021 G4double dz = std::abs(p.z())-fDz;
1022 G4double dy1 = fPlanes[0].b*p.y()+fPlanes[0].c*p.z()+fPlanes[0].d;
1023 G4double dy2 = fPlanes[1].b*p.y()+fPlanes[1].c*p.z()+fPlanes[1].d;
1024 G4double dy = std::max(dz,std::max(dy1,dy2));
1025
1026 G4double dx1 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()
1027 + fPlanes[2].c*p.z()+fPlanes[2].d;
1028 G4double dx2 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
1029 + fPlanes[3].c*p.z()+fPlanes[3].d;
1030 G4double dist = std::max(dy,std::max(dx1,dx2));
1031 return (dist < 0) ? -dist : 0.;
1032 }
1033 case 1: // YZ section is a rectangle
1034 {
1035 G4double dz = std::abs(p.z())-fDz;
1036 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
1037 G4double dx1 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()
1038 + fPlanes[2].c*p.z()+fPlanes[2].d;
1039 G4double dx2 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
1040 + fPlanes[3].c*p.z()+fPlanes[3].d;
1041 G4double dist = std::max(dy,std::max(dx1,dx2));
1042 return (dist < 0) ? -dist : 0.;
1043 }
1044 case 2: // YZ section is a rectangle and
1045 { // XZ section is an isosceles trapezoid
1046 G4double dz = std::abs(p.z())-fDz;
1047 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
1048 G4double dx = fPlanes[3].a*std::abs(p.x())
1049 + fPlanes[3].c*p.z()+fPlanes[3].d;
1050 G4double dist = std::max(dy,dx);
1051 return (dist < 0) ? -dist : 0.;
1052 }
1053 case 3: // YZ section is a rectangle and
1054 { // XY section is an isosceles trapezoid
1055 G4double dz = std::abs(p.z())-fDz;
1056 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
1057 G4double dx = fPlanes[3].a*std::abs(p.x())
1058 + fPlanes[3].b*p.y()+fPlanes[3].d;
1059 G4double dist = std::max(dy,dx);
1060 return (dist < 0) ? -dist : 0.;
1061 }
1062 }
1063 return 0.;
1064}
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:587
@ 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

Returns the distance along the normalised vector 'v' to the shape, from a point at an offset 'p' inside or on the surface of the shape. Intersections with surfaces, when the point is less than Tolerance/2 from a surface must be ignored.

Parameters
[in]pThe point at offset p.
[in]vThe normalised direction vector.
[in]calcNormFlag to indicate if to calculate the normal or not.
[out]validNormFlag set to true if the solid lies entirely behind or on the exiting surface. It is set false if the solid does not lie entirely behind or on the exiting surface. 'calcNorm' must be true, otherwise it is unused.
[out]nThe exiting outwards normal vector (undefined Magnitude). 'calcNorm' must be true, otherwise it is unused.
Returns
The distance to exit the shape.

Implements G4VSolid.

Definition at line 912 of file G4Trap.cc.

915{
916 // Z intersections
917 //
918 if ((std::abs(p.z()) - fDz) >= -halfCarTolerance && p.z()*v.z() > 0)
919 {
920 if (calcNorm)
921 {
922 *validNorm = true;
923 n->set(0, 0, (p.z() < 0) ? -1 : 1);
924 }
925 return 0;
926 }
927 G4double vz = v.z();
928 G4double tmax = (vz == 0) ? DBL_MAX : (std::copysign(fDz,vz) - p.z())/vz;
929 G4int iside = (vz < 0) ? -4 : -2; // little trick: (-4+3)=-1, (-2+3)=+1
930
931 // Y intersections
932 //
933 G4int i = 0;
934 for ( ; i<2; ++i)
935 {
936 G4double cosa = fPlanes[i].b*v.y() + fPlanes[i].c*v.z();
937 if (cosa > 0)
938 {
939 G4double dist = fPlanes[i].b*p.y() + fPlanes[i].c*p.z() + fPlanes[i].d;
940 if (dist >= -halfCarTolerance)
941 {
942 if (calcNorm)
943 {
944 *validNorm = true;
945 n->set(0, fPlanes[i].b, fPlanes[i].c);
946 }
947 return 0;
948 }
949 G4double tmp = -dist/cosa;
950 if (tmax > tmp) { tmax = tmp; iside = i; }
951 }
952 }
953
954 // X intersections
955 //
956 for ( ; i<4; ++i)
957 {
958 G4double cosa = fPlanes[i].a*v.x()+fPlanes[i].b*v.y()+fPlanes[i].c*v.z();
959 if (cosa > 0)
960 {
961 G4double dist = fPlanes[i].a*p.x() +
962 fPlanes[i].b*p.y() + fPlanes[i].c*p.z() + fPlanes[i].d;
963 if (dist >= -halfCarTolerance)
964 {
965 if (calcNorm)
966 {
967 *validNorm = true;
968 n->set(fPlanes[i].a, fPlanes[i].b, fPlanes[i].c);
969 }
970 return 0;
971 }
972 G4double tmp = -dist/cosa;
973 if (tmax > tmp) { tmax = tmp; iside = i; }
974 }
975 }
976
977 // Set normal, if required, and return distance
978 //
979 if (calcNorm)
980 {
981 *validNorm = true;
982 if (iside < 0)
983 {
984 n->set(0, 0, iside + 3); // (-4+3)=-1, (-2+3)=+1
985 }
986 else
987 {
988 n->set(fPlanes[iside].a, fPlanes[iside].b, fPlanes[iside].c);
989 }
990 }
991 return tmax;
992}

◆ GetAlpha1()

G4double G4Trap::GetAlpha1 ( ) const
inline

Referenced by StreamInfo().

◆ GetAlpha2()

G4double G4Trap::GetAlpha2 ( ) const
inline

Referenced by StreamInfo().

◆ GetCubicVolume()

G4double G4Trap::GetCubicVolume ( )
overridevirtual

Returning an estimation of the solid volume (capacity) and surface area, in internal units.

Reimplemented from G4VSolid.

Definition at line 1190 of file G4Trap.cc.

1191{
1192 if (fCubicVolume == 0)
1193 {
1194 G4AutoLock l(&trapMutex);
1195 G4ThreeVector pt[8];
1196 GetVertices(pt);
1197
1198 G4double dz = pt[4].z() - pt[0].z();
1199 G4double dy1 = pt[2].y() - pt[0].y();
1200 G4double dx1 = pt[1].x() - pt[0].x();
1201 G4double dx2 = pt[3].x() - pt[2].x();
1202 G4double dy2 = pt[6].y() - pt[4].y();
1203 G4double dx3 = pt[5].x() - pt[4].x();
1204 G4double dx4 = pt[7].x() - pt[6].x();
1205
1206 fCubicVolume = ((dx1 + dx2 + dx3 + dx4)*(dy1 + dy2) +
1207 (dx4 + dx3 - dx2 - dx1)*(dy2 - dy1)/3)*dz*0.125;
1208 l.unlock();
1209 }
1210 return fCubicVolume;
1211}
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4double fCubicVolume
Definition G4CSGSolid.hh:92

◆ GetEntityType()

G4GeometryType G4Trap::GetEntityType ( ) const
overridevirtual

Returns the type ID, "G4Trap" of the solid.

Implements G4VSolid.

Definition at line 1070 of file G4Trap.cc.

1071{
1072 return {"G4Trap"};
1073}

◆ GetPhi()

G4double G4Trap::GetPhi ( ) const
inline

Accessors obtaining (re)computed values of the original parameters.

Referenced by StreamInfo().

◆ GetPointOnSurface()

G4ThreeVector G4Trap::GetPointOnSurface ( ) const
overridevirtual

Returns a random point located and uniformly distributed on the surface of the solid.

Reimplemented from G4VSolid.

Definition at line 1149 of file G4Trap.cc.

1150{
1151 // Set indeces
1152 constexpr G4int iface [6][4] =
1153 { {0,1,3,2}, {0,4,5,1}, {2,3,7,6}, {0,2,6,4}, {1,5,7,3}, {4,6,7,5} };
1154
1155 // Set vertices
1156 G4ThreeVector pt[8];
1157 GetVertices(pt);
1158
1159 // Select face
1160 //
1161 G4double select = fAreas[5]*G4QuickRand();
1162 G4int k = 5;
1163 k -= (G4int)(select <= fAreas[4]);
1164 k -= (G4int)(select <= fAreas[3]);
1165 k -= (G4int)(select <= fAreas[2]);
1166 k -= (G4int)(select <= fAreas[1]);
1167 k -= (G4int)(select <= fAreas[0]);
1168
1169 // Select sub-triangle
1170 //
1171 G4int i0 = iface[k][0];
1172 G4int i1 = iface[k][1];
1173 G4int i2 = iface[k][2];
1174 G4int i3 = iface[k][3];
1175 G4double s2 = G4GeomTools::TriangleAreaNormal(pt[i2],pt[i1],pt[i3]).mag();
1176 if (select > fAreas[k] - s2) { i0 = i2; }
1177
1178 // Generate point
1179 //
1180 G4double u = G4QuickRand();
1181 G4double v = G4QuickRand();
1182 if (u + v > 1.) { u = 1. - u; v = 1. - v; }
1183 return (1.-u-v)*pt[i0] + u*pt[i1] + v*pt[i3];
1184}
G4double G4QuickRand(uint32_t seed=0)
double mag() const
static G4ThreeVector TriangleAreaNormal(const G4ThreeVector &A, const G4ThreeVector &B, const G4ThreeVector &C)

◆ GetSidePlane()

TrapSidePlane G4Trap::GetSidePlane ( G4int n) const
inline

More accessors.

◆ GetSurfaceArea()

G4double G4Trap::GetSurfaceArea ( )
overridevirtual

Returns an estimation of the solid surface area in internal units. This method may be overloaded by derived classes to compute the exact geometrical quantity for solids where this is possible, or anyway to cache the computed value. Note: the computed value is NOT cached.

Reimplemented from G4VSolid.

Definition at line 1217 of file G4Trap.cc.

1218{
1219 if (fSurfaceArea == 0)
1220 {
1221 G4AutoLock l(&trapMutex);
1222 G4ThreeVector pt[8];
1223 G4int iface [6][4] =
1224 { {0,1,3,2}, {0,4,5,1}, {2,3,7,6}, {0,2,6,4}, {1,5,7,3}, {4,6,7,5} };
1225
1226 GetVertices(pt);
1227 for (const auto & i : iface)
1228 {
1230 pt[i[1]],
1231 pt[i[2]],
1232 pt[i[3]]).mag();
1233 }
1234 l.unlock();
1235 }
1236 return fSurfaceArea;
1237}
G4double fSurfaceArea
Definition G4CSGSolid.hh:93
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()

G4double G4Trap::GetZHalfLength ( ) const
inline

Accessors. Returning the coordinates of a unit vector along a straight line joining centers of -/+fDz planes.

Referenced by BoundingLimits(), G4tgbGeometryDumper::GetSolidParams(), G4GDMLWriteParamvol::Trap_dimensionsWrite(), and G4GDMLWriteSolids::TrapWrite().

◆ Inside()

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

Concrete implementations of the expected query interfaces for solids, as defined in the base class G4VSolid.

Implements G4VSolid.

Definition at line 587 of file G4Trap.cc.

588{
589 switch (fTrapType)
590 {
591 case 0: // General case
592 {
593 G4double dz = std::abs(p.z())-fDz;
594 G4double dy1 = fPlanes[0].b*p.y()+fPlanes[0].c*p.z()+fPlanes[0].d;
595 G4double dy2 = fPlanes[1].b*p.y()+fPlanes[1].c*p.z()+fPlanes[1].d;
596 G4double dy = std::max(dz,std::max(dy1,dy2));
597
598 G4double dx1 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()
599 + fPlanes[2].c*p.z()+fPlanes[2].d;
600 G4double dx2 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
601 + fPlanes[3].c*p.z()+fPlanes[3].d;
602 G4double dist = std::max(dy,std::max(dx1,dx2));
603
604 return (dist > halfCarTolerance) ? kOutside :
605 ((dist > -halfCarTolerance) ? kSurface : kInside);
606 }
607 case 1: // YZ section is a rectangle
608 {
609 G4double dz = std::abs(p.z())-fDz;
610 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
611 G4double dx1 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()
612 + fPlanes[2].c*p.z()+fPlanes[2].d;
613 G4double dx2 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
614 + fPlanes[3].c*p.z()+fPlanes[3].d;
615 G4double dist = std::max(dy,std::max(dx1,dx2));
616
617 return (dist > halfCarTolerance) ? kOutside :
618 ((dist > -halfCarTolerance) ? kSurface : kInside);
619 }
620 case 2: // YZ section is a rectangle and
621 { // XZ section is an isosceles trapezoid
622 G4double dz = std::abs(p.z())-fDz;
623 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
624 G4double dx = fPlanes[3].a*std::abs(p.x())
625 + fPlanes[3].c*p.z()+fPlanes[3].d;
626 G4double dist = std::max(dy,dx);
627
628 return (dist > halfCarTolerance) ? kOutside :
629 ((dist > -halfCarTolerance) ? kSurface : kInside);
630 }
631 case 3: // YZ section is a rectangle and
632 { // XY section is an isosceles trapezoid
633 G4double dz = std::abs(p.z())-fDz;
634 G4double dy = std::max(dz,std::abs(p.y())+fPlanes[1].d);
635 G4double dx = fPlanes[3].a*std::abs(p.x())
636 + fPlanes[3].b*p.y()+fPlanes[3].d;
637 G4double dist = std::max(dy,dx);
638
639 return (dist > halfCarTolerance) ? kOutside :
640 ((dist > -halfCarTolerance) ? kSurface : kInside);
641 }
642 }
643 return kOutside;
644}
@ kInside
Definition geomdefs.hh:70
@ kSurface
Definition geomdefs.hh:69

Referenced by DistanceToOut().

◆ IsFaceted()

G4bool G4Trap::IsFaceted ( ) const
overridevirtual

Returns true as the solid has only planar faces.

Reimplemented from G4VSolid.

Definition at line 1079 of file G4Trap.cc.

1080{
1081 return true;
1082}

◆ MakePlane()

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

Calculates the coefficents of the plane p1->p2->p3->p4->p1 where the ThreeVectors 1-4 are in anti-clockwise order when viewed from infront of the plane (i.e. from normal direction).

Returns
true if the points are co-planar, false otherwise.

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)
G4double b
Definition G4Trap.hh:90
G4double c
Definition G4Trap.hh:90
G4double d
Definition G4Trap.hh:90
G4double a
Definition G4Trap.hh:90
#define DBL_EPSILON
Definition templates.hh:66

Referenced by MakePlanes().

◆ MakePlanes() [1/2]

void G4Trap::MakePlanes ( )
protected

Internal methods for checking and building planes. Computing the vertices and setting side planes, checking for planarity.

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}

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:1097
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 )

Sets all parameters, as for constructor. Checks and sets half-widths as well as angles. Makes a final check of co-planarity.

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:94

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

◆ SetCachedValues()

void G4Trap::SetCachedValues ( )
protected

Recomputes parameters using planes.

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

Streams the object contents to an output stream.

Reimplemented from G4CSGSolid.

Definition at line 1097 of file G4Trap.cc.

1098{
1099 G4double phi = GetPhi();
1100 G4double theta = GetTheta();
1101 G4double alpha1 = GetAlpha1();
1103
1104 G4long oldprc = os.precision(16);
1105 os << "-----------------------------------------------------------\n"
1106 << " *** Dump for solid: " << GetName() << " ***\n"
1107 << " ===================================================\n"
1108 << " Solid type: G4Trap\n"
1109 << " Parameters:\n"
1110 << " half length Z: " << fDz/mm << " mm\n"
1111 << " half length Y, face -Dz: " << fDy1/mm << " mm\n"
1112 << " half length X, face -Dz, side -Dy1: " << fDx1/mm << " mm\n"
1113 << " half length X, face -Dz, side +Dy1: " << fDx2/mm << " mm\n"
1114 << " half length Y, face +Dz: " << fDy2/mm << " mm\n"
1115 << " half length X, face +Dz, side -Dy2: " << fDx3/mm << " mm\n"
1116 << " half length X, face +Dz, side +Dy2: " << fDx4/mm << " mm\n"
1117 << " theta: " << theta/degree << " degrees\n"
1118 << " phi: " << phi/degree << " degrees\n"
1119 << " alpha, face -Dz: " << alpha1/degree << " degrees\n"
1120 << " alpha, face +Dz: " << alpha2/degree << " degrees\n"
1121 << "-----------------------------------------------------------\n";
1122 os.precision(oldprc);
1123
1124 return os;
1125}
G4double GetAlpha2() const
G4double GetAlpha1() const
G4double GetTheta() const
G4double GetPhi() const

Referenced by MakePlanes().

◆ SurfaceNormal()

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

Returns the outwards pointing unit normal of the shape for the surface closest to the point at offset 'p'.

Parameters
[in]pThe point at offset p.
Returns
The outwards pointing unit normal.

Implements G4VSolid.

Definition at line 650 of file G4Trap.cc.

651{
652 G4double nx = 0, ny = 0, nz = 0;
653 G4double dz = std::abs(p.z()) - fDz;
654 nz = std::copysign(G4double(std::abs(dz) <= halfCarTolerance), p.z());
655
656 switch (fTrapType)
657 {
658 case 0: // General case
659 {
660 for (G4int i=0; i<2; ++i)
661 {
662 G4double dy = fPlanes[i].b*p.y() + fPlanes[i].c*p.z() + fPlanes[i].d;
663 if (std::abs(dy) > halfCarTolerance) { continue; }
664 ny = fPlanes[i].b;
665 nz += fPlanes[i].c;
666 break;
667 }
668 for (G4int i=2; i<4; ++i)
669 {
670 G4double dx = fPlanes[i].a*p.x() +
671 fPlanes[i].b*p.y() + fPlanes[i].c*p.z() + fPlanes[i].d;
672 if (std::abs(dx) > halfCarTolerance) { continue; }
673 nx = fPlanes[i].a;
674 ny += fPlanes[i].b;
675 nz += fPlanes[i].c;
676 break;
677 }
678 break;
679 }
680 case 1: // YZ section - rectangle
681 {
682 G4double dy = std::abs(p.y()) + fPlanes[1].d;
683 ny = std::copysign(G4double(std::abs(dy) <= halfCarTolerance), p.y());
684 for (G4int i=2; i<4; ++i)
685 {
686 G4double dx = fPlanes[i].a*p.x() +
687 fPlanes[i].b*p.y() + fPlanes[i].c*p.z() + fPlanes[i].d;
688 if (std::abs(dx) > halfCarTolerance) { continue; }
689 nx = fPlanes[i].a;
690 ny += fPlanes[i].b;
691 nz += fPlanes[i].c;
692 break;
693 }
694 break;
695 }
696 case 2: // YZ section - rectangle, XZ section - isosceles trapezoid
697 {
698 G4double dy = std::abs(p.y()) + fPlanes[1].d;
699 ny = std::copysign(G4double(std::abs(dy) <= halfCarTolerance), p.y());
700 G4double dx = fPlanes[3].a*std::abs(p.x()) +
701 fPlanes[3].c*p.z() + fPlanes[3].d;
702 G4double k = static_cast<G4double>(std::abs(dx) <= halfCarTolerance);
703 nx = std::copysign(k, p.x())*fPlanes[3].a;
704 nz += k*fPlanes[3].c;
705 break;
706 }
707 case 3: // YZ section - rectangle, XY section - isosceles trapezoid
708 {
709 G4double dy = std::abs(p.y()) + fPlanes[1].d;
710 ny = std::copysign(G4double(std::abs(dy) <= halfCarTolerance), p.y());
711 G4double dx = fPlanes[3].a*std::abs(p.x()) +
712 fPlanes[3].b*p.y() + fPlanes[3].d;
713 G4double k = static_cast<G4double>(std::abs(dx) <= halfCarTolerance);
714 nx = std::copysign(k, p.x())*fPlanes[3].a;
715 ny += k*fPlanes[3].b;
716 break;
717 }
718 }
719
720 // Return normal
721 //
722 G4double mag2 = nx*nx + ny*ny + nz*nz;
723 if (mag2 == 1)
724 {
725 return { nx,ny,nz };
726 }
727 if (mag2 != 0)
728 {
729 return G4ThreeVector(nx,ny,nz).unit(); // edge or corner
730 }
731
732 // Point is not on the surface
733 //
734#ifdef G4CSGDEBUG
735 std::ostringstream message;
736 G4long oldprc = message.precision(16);
737 message << "Point p is not on surface (!?) of solid: "
738 << GetName() << G4endl;
739 message << "Position:\n";
740 message << " p.x() = " << p.x()/mm << " mm\n";
741 message << " p.y() = " << p.y()/mm << " mm\n";
742 message << " p.z() = " << p.z()/mm << " mm";
743 G4cout.precision(oldprc) ;
744 G4Exception("G4Trap::SurfaceNormal(p)", "GeomSolids1002",
745 JustWarning, message );
746 DumpInfo();
747#endif
748
749 return ApproxSurfaceNormal(p);
750}

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