Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolyhedraSide.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4PolyhedraSide
27//
28// Class description:
29//
30// Class implementing a face that represents one segmented side
31// of a polyhedra:
32//
33// G4PolyhedraSide( const G4PolyhedraSideRZ* prevRZ,
34// const G4PolyhedraSideRZ* tail,
35// const G4PolyhedraSideRZ* head,
36// const G4PolyhedraSideRZ* nextRZ,
37// G4int numSide,
38// G4double phiStart, G4double phiTotal,
39// G4bool phiIsOpen, G4bool isAllBehind = false )
40//
41// Values for r1,z1 and r2,z2 should be specified in clockwise
42// order in (r,z).
43
44// Author: David C. Williams ([email protected])
45// --------------------------------------------------------------------
46
47#ifndef G4PolyhedraSide_hh
48#define G4PolyhedraSide_hh
49
50#include "G4VCSGface.hh"
51
53
55{
56 G4double r, z; // start of vector
57};
58
59// ----------------------------------------------------------------------------
60// MT-specific utility code
61
62#include "G4GeomSplitter.hh"
63
64// The class G4PhSideData is introduced to encapsulate the
65// fields of the class G4PolyhedraSide that may not be read-only.
66//
68{
69 public:
70
72 {
73 fPhix = 0.; fPhiy = 0.; fPhiz = 0.; fPhik = 0.;
74 }
75
76 G4double fPhix=0., fPhiy=0., fPhiz=0., fPhik=0.; // Cached values for phi
77};
78
79// The type G4PhSideManager is introduced to encapsulate the methods used
80// by both the master thread and worker threads to allocate memory space
81// for the fields encapsulated by the class G4PhSideData.
82//
84
85//
86// ----------------------------------------------------------------------------
87
89{
90
91 public:
92
94 const G4PolyhedraSideRZ* tail,
95 const G4PolyhedraSideRZ* head,
96 const G4PolyhedraSideRZ* nextRZ,
98 G4double phiStart, G4double phiTotal,
99 G4bool phiIsOpen, G4bool isAllBehind = false );
100 ~G4PolyhedraSide() override;
101
102 G4PolyhedraSide( const G4PolyhedraSide& source );
103 G4PolyhedraSide& operator=( const G4PolyhedraSide& source );
104
105 G4bool Intersect( const G4ThreeVector& p, const G4ThreeVector& v,
106 G4bool outgoing, G4double surfTolerance,
107 G4double& distance, G4double& distFromSurface,
108 G4ThreeVector& normal, G4bool& allBehind ) override;
109
110 G4double Distance( const G4ThreeVector& p, G4bool outgoing ) override;
111
112 EInside Inside( const G4ThreeVector &p, G4double tolerance,
113 G4double *bestDistance ) override;
114
116 G4double* bestDistance ) override;
117
118 G4double Extent( const G4ThreeVector axis ) override;
119
120 void CalculateExtent( const EAxis axis,
121 const G4VoxelLimits &voxelLimit,
122 const G4AffineTransform& tranform,
123 G4SolidExtentList& extentList ) override;
124
125 G4VCSGface* Clone() override { return new G4PolyhedraSide( *this ); }
126
127 // Methods used for GetPointOnSurface()
128
130 const G4ThreeVector& p2,
131 const G4ThreeVector& p3,
132 G4ThreeVector* p4 );
134 const G4ThreeVector& p2, const G4ThreeVector& p3,
135 G4double* Area );
136 G4double SurfaceArea() override;
137 G4ThreeVector GetPointOnFace() override;
138
139 G4PolyhedraSide(__void__&);
140 // Fake default constructor for usage restricted to direct object
141 // persistency for clients requiring preallocation of memory for
142 // persistifiable objects.
143
144 inline G4int GetInstanceID() const { return instanceID; }
145 // Returns the instance ID.
146
148 // Returns the private data instance manager.
149
150 //
151 // A couple internal data structures
152 //
153 struct sG4PolyhedraSideVec; // Secret recipe for allowing
154 friend struct sG4PolyhedraSideVec; // protected nested structures
155
156 using G4PolyhedraSideEdge = struct sG4PolyhedraSideEdge
157 {
158 G4ThreeVector normal; // Unit normal to this edge
159 G4ThreeVector corner[2]; // The two corners of this phi edge
160 G4ThreeVector cornNorm[2]; // The normals of these corners
161 };
162
164 {
165 G4ThreeVector normal, // Normal (point out of the shape)
166 center, // Point in center of side
167 surfPhi, // Unit vector on surface pointing along phi
168 surfRZ; // Unit vector on surface pointing along R/Z
169 G4PolyhedraSideEdge* edges[2]; // The phi boundary edges to this side
170 // [0]=low phi [1]=high phi
171 G4ThreeVector edgeNorm[2]; // RZ edge normals [i] at {r[i],z[i]}
172 };
173
174 protected:
175
177 const G4PolyhedraSideVec& vec,
178 G4double normSign,
179 G4double surfTolerance,
180 G4double &distance,
181 G4double &distFromSurface );
182
184 const G4ThreeVector& v,
185 G4int* i1, G4int* i2 );
186
188
190
191 G4double GetPhi( const G4ThreeVector& p );
192
194 const G4PolyhedraSideVec& vec,
195 G4double* normDist );
196
198 const G4PolyhedraSideVec& vec,
199 G4double* normDist );
200
201 void CopyStuff( const G4PolyhedraSide& source );
202
203 protected:
204
205 G4int numSide = 0; // Number sides
206 G4double r[2], z[2]; // r, z parameters, in specified order
207 G4double startPhi, // Start phi (0 to 2pi), if phiIsOpen
208 deltaPhi, // Delta phi (0 to 2pi), if phiIsOpen
209 endPhi; // End phi (>startPhi), if phiIsOpen
210 G4bool phiIsOpen = false; // True if there is a phi slice
211 G4bool allBehind = false; // True if the entire solid is "behind" this face
212
213 G4IntersectingCone* cone = nullptr; // Our intersecting cone
214
215 G4PolyhedraSideVec* vecs = nullptr; // Vector set for each facet of our face
216 G4PolyhedraSideEdge* edges = nullptr; // The edges belong to vecs
217 G4double lenRZ, // RZ length of each side
218 lenPhi[2]; // Phi dimensions of each side
219 G4double edgeNorm; // Normal in RZ/Phi space to each side
220
221 private:
222
223 G4double kCarTolerance; // Geometrical surface thickness
224 G4double fSurfaceArea = 0.0; // Surface Area
225
226 G4int instanceID;
227 // This field is used as instance ID.
228 G4GEOM_DLL static G4PhSideManager subInstanceManager;
229 // This field helps to use the class G4PhSideManager introduced above.
230};
231
232#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4PolyhedraSideVec * vecs
G4PolyhedraSideEdge * edges
struct sG4PolyhedraSideVec { G4ThreeVector normal, center, surfPhi, surfRZ; G4PolyhedraSideEdge *edges[2]; G4ThreeVector edgeNorm[2]; } G4PolyhedraSideVec
G4ThreeVector Normal(const G4ThreeVector &p, G4double *bestDistance) override
G4ThreeVector GetPointOnPlane(const G4ThreeVector &p0, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3, G4double *Area)
static const G4PhSideManager & GetSubInstanceManager()
G4PolyhedraSide & operator=(const G4PolyhedraSide &source)
G4double GetPhi(const G4ThreeVector &p)
G4bool IntersectSidePlane(const G4ThreeVector &p, const G4ThreeVector &v, const G4PolyhedraSideVec &vec, G4double normSign, G4double surfTolerance, G4double &distance, G4double &distFromSurface)
G4double SurfaceTriangle(const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3, G4ThreeVector *p4)
G4PolyhedraSide(const G4PolyhedraSideRZ *prevRZ, const G4PolyhedraSideRZ *tail, const G4PolyhedraSideRZ *head, const G4PolyhedraSideRZ *nextRZ, G4int numSide, G4double phiStart, G4double phiTotal, G4bool phiIsOpen, G4bool isAllBehind=false)
G4int GetInstanceID() const
G4int ClosestPhiSegment(G4double phi)
struct sG4PolyhedraSideEdge { G4ThreeVector normal; G4ThreeVector corner[2]; G4ThreeVector cornNorm[2]; } G4PolyhedraSideEdge
G4int PhiSegment(G4double phi)
G4int LineHitsSegments(const G4ThreeVector &p, const G4ThreeVector &v, G4int *i1, G4int *i2)
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, G4bool outgoing, G4double surfTolerance, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal, G4bool &allBehind) override
G4double SurfaceArea() override
EInside Inside(const G4ThreeVector &p, G4double tolerance, G4double *bestDistance) override
G4double Extent(const G4ThreeVector axis) override
G4double Distance(const G4ThreeVector &p, G4bool outgoing) override
G4double DistanceAway(const G4ThreeVector &p, const G4PolyhedraSideVec &vec, G4double *normDist)
void CalculateExtent(const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList) override
G4double DistanceToOneSide(const G4ThreeVector &p, const G4PolyhedraSideVec &vec, G4double *normDist)
~G4PolyhedraSide() override
friend struct sG4PolyhedraSideVec
G4IntersectingCone * cone
G4VCSGface * Clone() override
void CopyStuff(const G4PolyhedraSide &source)
G4ThreeVector GetPointOnFace() override
EAxis
Definition geomdefs.hh:54
EInside
Definition geomdefs.hh:67
#define G4GEOM_DLL
Definition geomwdefs.hh:44