Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 (davidw@scipp.ucsc.edu)
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: // with description
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 virtual ~G4PolyhedraSide();
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 );
109
110 G4double Distance( const G4ThreeVector& p, G4bool outgoing );
111
112 EInside Inside( const G4ThreeVector &p, G4double tolerance,
113 G4double *bestDistance );
114
115 G4ThreeVector Normal( const G4ThreeVector& p, G4double* bestDistance );
116
117 G4double Extent( const G4ThreeVector axis );
118
119 void CalculateExtent( const EAxis axis,
120 const G4VoxelLimits &voxelLimit,
121 const G4AffineTransform& tranform,
122 G4SolidExtentList& extentList );
123
124 G4VCSGface* Clone() { return new G4PolyhedraSide( *this ); }
125
126 public: // without description
127
128 // Methods used for GetPointOnSurface()
129
131 G4ThreeVector p2,
132 G4ThreeVector p3,
133 G4ThreeVector* p4 );
136 G4double* Area );
139
140 public: // without description
141
142 G4PolyhedraSide(__void__&);
143 // Fake default constructor for usage restricted to direct object
144 // persistency for clients requiring preallocation of memory for
145 // persistifiable objects.
146
147 inline G4int GetInstanceID() const { return instanceID; }
148 // Returns the instance ID.
149
151 // Returns the private data instance manager.
152
153 //
154 // A couple internal data structures
155 //
156 struct sG4PolyhedraSideVec; // Secret recipe for allowing
157 friend struct sG4PolyhedraSideVec; // protected nested structures
158
159 typedef struct sG4PolyhedraSideEdge
160 {
161 G4ThreeVector normal; // Unit normal to this edge
162 G4ThreeVector corner[2]; // The two corners of this phi edge
163 G4ThreeVector cornNorm[2]; // The normals of these corners
165
166 typedef struct sG4PolyhedraSideVec
167 {
168 G4ThreeVector normal, // Normal (point out of the shape)
169 center, // Point in center of side
170 surfPhi, // Unit vector on surface pointing along phi
171 surfRZ; // Unit vector on surface pointing along R/Z
172 G4PolyhedraSideEdge* edges[2]; // The phi boundary edges to this side
173 // [0]=low phi [1]=high phi
174 G4ThreeVector edgeNorm[2]; // RZ edge normals [i] at {r[i],z[i]}
176
177 protected:
178
180 const G4PolyhedraSideVec& vec,
181 G4double normSign,
182 G4double surfTolerance,
183 G4double &distance,
184 G4double &distFromSurface );
185
187 const G4ThreeVector& v,
188 G4int* i1, G4int* i2 );
189
191
193
194 G4double GetPhi( const G4ThreeVector& p );
195
197 const G4PolyhedraSideVec& vec,
198 G4double* normDist );
199
201 const G4PolyhedraSideVec& vec,
202 G4double* normDist );
203
204 void CopyStuff( const G4PolyhedraSide& source );
205
206 protected:
207
208 G4int numSide = 0; // Number sides
209 G4double r[2], z[2]; // r, z parameters, in specified order
210 G4double startPhi, // Start phi (0 to 2pi), if phiIsOpen
211 deltaPhi, // Delta phi (0 to 2pi), if phiIsOpen
212 endPhi; // End phi (>startPhi), if phiIsOpen
213 G4bool phiIsOpen = false; // True if there is a phi slice
214 G4bool allBehind = false; // True if the entire solid is "behind" this face
215
216 G4IntersectingCone* cone = nullptr; // Our intersecting cone
217
218 G4PolyhedraSideVec* vecs = nullptr; // Vector set for each facet of our face
219 G4PolyhedraSideEdge* edges = nullptr; // The edges belong to vecs
220 G4double lenRZ, // RZ length of each side
221 lenPhi[2]; // Phi dimensions of each side
222 G4double edgeNorm; // Normal in RZ/Phi space to each side
223
224 private:
225
226 G4double kCarTolerance; // Geometrical surface thickness
227 G4double fSurfaceArea = 0.0; // Surface Area
228
229 G4int instanceID;
230 // This field is used as instance ID.
231 G4GEOM_DLL static G4PhSideManager subInstanceManager;
232 // This field helps to use the class G4PhSideManager introduced above.
233};
234
235#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
EInside Inside(const G4ThreeVector &p, G4double tolerance, G4double *bestDistance)
G4PolyhedraSideVec * vecs
G4double Distance(const G4ThreeVector &p, G4bool outgoing)
G4PolyhedraSideEdge * edges
G4ThreeVector GetPointOnPlane(G4ThreeVector p0, G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4double *Area)
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, G4bool outgoing, G4double surfTolerance, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal, G4bool &allBehind)
static const G4PhSideManager & GetSubInstanceManager()
G4VCSGface * Clone()
G4PolyhedraSide & operator=(const G4PolyhedraSide &source)
G4double SurfaceTriangle(G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector *p4)
struct G4PolyhedraSide::sG4PolyhedraSideVec G4PolyhedraSideVec
G4double GetPhi(const G4ThreeVector &p)
G4bool IntersectSidePlane(const G4ThreeVector &p, const G4ThreeVector &v, const G4PolyhedraSideVec &vec, G4double normSign, G4double surfTolerance, G4double &distance, G4double &distFromSurface)
G4int GetInstanceID() const
G4ThreeVector GetPointOnFace()
G4double Extent(const G4ThreeVector axis)
G4int ClosestPhiSegment(G4double phi)
G4double lenPhi[2]
struct G4PolyhedraSide::sG4PolyhedraSideEdge G4PolyhedraSideEdge
G4int PhiSegment(G4double phi)
G4int LineHitsSegments(const G4ThreeVector &p, const G4ThreeVector &v, G4int *i1, G4int *i2)
G4ThreeVector Normal(const G4ThreeVector &p, G4double *bestDistance)
G4double DistanceAway(const G4ThreeVector &p, const G4PolyhedraSideVec &vec, G4double *normDist)
G4double DistanceToOneSide(const G4ThreeVector &p, const G4PolyhedraSideVec &vec, G4double *normDist)
G4double SurfaceArea()
virtual ~G4PolyhedraSide()
G4IntersectingCone * cone
void CalculateExtent(const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList)
void CopyStuff(const G4PolyhedraSide &source)
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67
#define G4GEOM_DLL
Definition: geomwdefs.hh:44