Geant4 11.4.0
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 of a polyhedra:
31//
32// G4PolyhedraSide( const G4PolyhedraSideRZ* prevRZ,
33// const G4PolyhedraSideRZ* tail,
34// const G4PolyhedraSideRZ* head,
35// const G4PolyhedraSideRZ* nextRZ,
36// G4int numSide,
37// G4double phiStart, G4double phiTotal,
38// G4bool phiIsOpen, G4bool isAllBehind=false )
39//
40// Values for r1,z1 and r2,z2 should be specified in clockwise order in (r,z).
41
42// Author: David C. Williams (UCSC), 1998 - Created
43// --------------------------------------------------------------------
44#ifndef G4POLYHEDRASIDE_HH
45#define G4POLYHEDRASIDE_HH
46
47#include "G4VCSGface.hh"
48
50
52{
53 G4double r, z; // start of vector
54};
55
56// ----------------------------------------------------------------------------
57// MT-specific utility code
58
59#include "G4GeomSplitter.hh"
60
61// The class G4PhSideData is introduced to encapsulate the
62// fields of the class G4PolyhedraSide that may not be read-only.
63//
65{
66 public:
67
69 {
70 fPhix = 0.; fPhiy = 0.; fPhiz = 0.; fPhik = 0.;
71 }
72
73 G4double fPhix=0., fPhiy=0., fPhiz=0., fPhik=0.; // Cached values for phi
74};
75
76// The type G4PhSideManager is introduced to encapsulate the methods used
77// by both the master thread and worker threads to allocate memory space
78// for the fields encapsulated by the class G4PhSideData.
79//
81
82//
83// ----------------------------------------------------------------------------
84
85/**
86 * @brief G4PolyhedraSide is a utility class implementing a face that
87 * represents one segmented side of a polyhedra.
88 */
89
91{
92
93 public:
94
95 /**
96 * Constructor for the segmented side of a polyhedra.
97 * @param[in] prevRZ Pointer to previous r,Z section.
98 * @param[in] tail Pointer to r,Z tail of section.
99 * @param[in] head Pointer to r,Z head of section.
100 * @param[in] nextRZ Pointer to next r,Z section.
101 * @param[in] numSide The number od sides.
102 * @param[in] phiStart Initial Phi starting angle.
103 * @param[in] phiTotal Total Phi angle.
104 * @param[in] phiIsOpen Flag indicating if it is a Phi section.
105 * @param[in] isAllBehind Indicating if entire surface is behind normal.
106 */
107 G4PolyhedraSide( const G4PolyhedraSideRZ* prevRZ,
108 const G4PolyhedraSideRZ* tail,
109 const G4PolyhedraSideRZ* head,
110 const G4PolyhedraSideRZ* nextRZ,
111 G4int numSide,
112 G4double phiStart, G4double phiTotal,
113 G4bool phiIsOpen, G4bool isAllBehind = false );
114
115 /**
116 * Destructor.
117 */
118 ~G4PolyhedraSide() override;
119
120 /**
121 * Copy constructor and assignment operator.
122 */
123 G4PolyhedraSide( const G4PolyhedraSide& source );
124 G4PolyhedraSide& operator=( const G4PolyhedraSide& source );
125
126 /**
127 * Determines the distance along a line to the face.
128 * @param[in] p Position.
129 * @param[in] v Direction (assumed to be a unit vector).
130 * @param[in] outgoing Flag true, to consider only inside surfaces;
131 * false, to consider only outside surfaces.
132 * @param[in] surfTolerance Minimum distance from the surface.
133 * @param[out] distance Distance to intersection.
134 * @param[out] distFromSurface Distance from surface (along surface normal),
135 * < 0 if the point is in front of the surface.
136 * @param[out] normal Normal of surface at intersection point.
137 * @param[out] allBehind Flag, true, if entire surface is behind normal.
138 * @returns true if there is an intersection, false otherwise.
139 */
140 G4bool Intersect( const G4ThreeVector& p, const G4ThreeVector& v,
141 G4bool outgoing, G4double surfTolerance,
142 G4double& distance, G4double& distFromSurface,
143 G4ThreeVector& normal, G4bool& allBehind ) override;
144
145 /**
146 * Determines the distance of a point from either the inside or outside
147 * surfaces of the face.
148 * @param[in] p Position.
149 * @param[in] outgoing Flag, true, to consider only inside surfaces
150 * or false, to consider only outside surfaces.
151 * @returns The distance to the closest surface satisfying requirements
152 * or kInfinity if no such surface exists.
153 */
154 G4double Distance( const G4ThreeVector& p, G4bool outgoing ) override;
155
156 /**
157 * Determines whether a point is inside, outside, or on the surface of
158 * the face.
159 * @param[in] p Position.
160 * @param[in] tolerance Tolerance defining the bounds of the "kSurface",
161 * nominally equal to kCarTolerance/2.
162 * @param[out] bestDistance Distance to the closest surface (in or out).
163 * @returns kInside if the point is closest to the inside surface;
164 * kOutside if the point is closest to the outside surface;
165 * kSurface if the point is withing tolerance of the surface.
166 */
167 EInside Inside( const G4ThreeVector &p, G4double tolerance,
168 G4double *bestDistance ) override;
169
170 /**
171 * Returns the normal of surface closest to the point.
172 * @param[in] p Position.
173 * @param[out] bestDistance Distance to the closest surface (in or out).
174 * @returns The normal of the surface nearest the point.
175 */
177 G4double* bestDistance ) override;
178
179 /**
180 * Returns the face extent along the axis.
181 * @param[in] axis Unit vector defining the direction.
182 * @returns The largest point along the given axis of the face's extent.
183 */
184 G4double Extent( const G4ThreeVector axis ) override;
185
186 /**
187 * Calculates the extent of the face for the voxel navigator.
188 * @param[in] axis The axis in which to check the shapes 3D extent against.
189 * @param[in] voxelLimit Limits along x, y, and/or z axes.
190 * @param[in] tranform A coordinate transformation on which to apply to
191 * the shape before testing.
192 * @param[out] extentList The list of (voxel) extents along the axis.
193 */
194 void CalculateExtent( const EAxis axis,
195 const G4VoxelLimits &voxelLimit,
196 const G4AffineTransform& tranform,
197 G4SolidExtentList& extentList ) override;
198
199 /**
200 * Method invoked by the copy constructor or the assignment operator.
201 * Its purpose is to return a pointer to a duplicate copy of the face.
202 */
203 inline G4VCSGface* Clone() override { return new G4PolyhedraSide( *this ); }
204
205 /**
206 * Returning an estimation of the face surface area, in internal units.
207 */
208 G4double SurfaceArea() override;
209
210 /**
211 * Fake default constructor for usage restricted to direct object
212 * persistency for clients requiring preallocation of memory for
213 * persistifiable objects.
214 */
215 G4PolyhedraSide(__void__&);
216
217 /**
218 * Returns the instance ID.
219 */
220 inline G4int GetInstanceID() const { return instanceID; }
221
222 /**
223 * Returns the private data instance manager.
224 */
226
227 private:
228
229 /**
230 * Internal data structures.
231 */
232 struct sG4PolyhedraSideVec; // Secret recipe for allowing
233 friend struct sG4PolyhedraSideVec; // protected nested structures
234
235 using G4PolyhedraSideEdge = struct sG4PolyhedraSideEdge
236 {
237 G4ThreeVector normal; // Unit normal to this edge
238 G4ThreeVector corner[2]; // The two corners of this phi edge
239 G4ThreeVector cornNorm[2]; // The normals of these corners
240 };
241
242 using G4PolyhedraSideVec = struct sG4PolyhedraSideVec
243 {
244 G4ThreeVector normal, // Normal (point out of the shape)
245 center, // Point in center of side
246 surfPhi, // Unit vector on surface pointing along phi
247 surfRZ; // Unit vector on surface pointing along R/Z
248 G4PolyhedraSideEdge* edges[2]; // The phi boundary edges to this side
249 // [0]=low phi [1]=high phi
250 G4ThreeVector edgeNorm[2]; // RZ edge normals [i] at {r[i],z[i]}
251 };
252
253 /**
254 * Calculates the surface area of a triangle.
255 * At the same time a random point in the triangle is given.
256 */
257 G4double SurfaceTriangle( const G4ThreeVector& p1,
258 const G4ThreeVector& p2,
259 const G4ThreeVector& p3,
260 G4ThreeVector* p4 );
261
262 /**
263 * Returns a random point located and uniformly distributed on the face.
264 */
265 G4ThreeVector GetPointOnFace() override;
266
267 /**
268 * Auxiliary method for GetPointOnSurface().
269 */
270 G4ThreeVector GetPointOnPlane( const G4ThreeVector& p0, const G4ThreeVector& p1,
271 const G4ThreeVector& p2, const G4ThreeVector& p3,
272 G4double* Area );
273
274 /**
275 * Decides if a line correctly intersects one side plane of a segment.
276 * It is assumed that the correct side has been chosen, and thus only
277 * the Z bounds (of the entire segment) are checked.
278 * @param[in] p The point to check.
279 * @param[in] v The direction.
280 * @param[in] vec Description record of the side plane.
281 * @param[in] normSign Sign (+/- 1) to apply to normal.
282 * @param[in] surfTolerance Surface tolerance (generally > 0, see below).
283 * @param[out] distance Distance along v to intersection.
284 * @param[out] distFromSurface Distance from surface normal.
285 * @returns true is the line is intersecting, false otherwise.
286 */
287 G4bool IntersectSidePlane( const G4ThreeVector& p, const G4ThreeVector& v,
288 const G4PolyhedraSideVec& vec,
289 G4double normSign,
290 G4double surfTolerance,
291 G4double& distance,
292 G4double& distFromSurface );
293
294 /**
295 * Calculates which phi segments a line intersects in three dimensions.
296 * No check is made as to whether the intersections are within the Z bounds
297 * of the segment.
298 */
299 G4int LineHitsSegments( const G4ThreeVector& p,
300 const G4ThreeVector& v,
301 G4int* i1, G4int* i2 );
302
303 /**
304 * Decides which phi segment an angle belongs to, counting from zero.
305 * A returned value of -1 indicates that the phi value is outside the
306 * shape (only possible if phiTotal < 360 degrees).
307 */
308 G4int PhiSegment( G4double phi );
309
310 /**
311 * Decides which phi segment is closest in phi to the point.
312 * The result is the same as PhiSegment() if there is no phi opening.
313 */
314 G4int ClosestPhiSegment( G4double phi );
315
316 /**
317 * Calculates Phi for a given 3-vector (point 'p'), if not already cached
318 * for the same point, in the attempt to avoid consecutive computation of
319 * the same quantity.
320 */
321 G4double GetPhi( const G4ThreeVector& p );
322
323 /**
324 * Computes the total distance from the side.
325 * @param[in] p The point to check.
326 * @param[in] vec The vector set of this side.
327 * @param[out] normDist The returned distance normal to the side or edge,
328 * as appropriate, signed.
329 * @returns The total distance from the side.
330 */
331 G4double DistanceToOneSide( const G4ThreeVector& p,
332 const G4PolyhedraSideVec& vec,
333 G4double* normDist );
334
335 /**
336 * Calculates the distance of a point from the segmented surface, including
337 * the effect of any phi segmentation.
338 * Adds distance from side edges, if necessary, to the total distance,
339 * and updates 'normDist' appropriate depending on edge normals.
340 * @param[in] p The point to check.
341 * @param[in] opposite If true, check the opposite hemisphere (see below).
342 * @param[out] normDist The returned normal distance.
343 * @returns The distance from the segmented plane.
344 */
345 G4double DistanceAway( const G4ThreeVector& p,
346 const G4PolyhedraSideVec& vec,
347 G4double* normDist );
348
349 /**
350 * Copies parameters from other object; used in copy constructor and
351 * assignment operator.
352 */
353 void CopyStuff( const G4PolyhedraSide& source );
354
355 private:
356
357 G4int numSide = 0; // Number sides
358 G4double r[2], z[2]; // r, z parameters, in specified order
359 G4double startPhi, // Start phi (0 to 2pi), if phiIsOpen
360 deltaPhi, // Delta phi (0 to 2pi), if phiIsOpen
361 endPhi; // End phi (>startPhi), if phiIsOpen
362 G4bool phiIsOpen = false; // True if there is a phi slice
363 G4bool allBehind = false; // True if the entire solid is "behind" this face
364
365 G4IntersectingCone* cone = nullptr; // Our intersecting cone
366
367 G4PolyhedraSideVec* vecs = nullptr; // Vector set for each facet of our face
368 G4PolyhedraSideEdge* edges = nullptr; // The edges belong to vecs
369 G4double lenRZ, // RZ length of each side
370 lenPhi[2]; // Phi dimensions of each side
371 G4double edgeNorm; // Normal in RZ/Phi space to each side
372
373 G4double kCarTolerance; // Geometrical surface thickness
374 G4double fSurfaceArea = 0.0; // Surface Area
375
376 G4int instanceID;
377 // This field is used as instance ID.
378 G4GEOM_DLL static G4PhSideManager subInstanceManager;
379 // This field helps to use the class G4PhSideManager introduced above.
380};
381
382#endif
G4GeomSplitter< G4PhSideData > G4PhSideManager
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4GeomSplitter is an utility class for splitting of R/W data for thread-safety from geometry classes....
G4IntersectingCone is a utility class used to calculate the intersection of an arbitrary line with a ...
G4ThreeVector Normal(const G4ThreeVector &p, G4double *bestDistance) override
static const G4PhSideManager & GetSubInstanceManager()
G4PolyhedraSide & operator=(const G4PolyhedraSide &source)
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
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
void CalculateExtent(const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList) override
~G4PolyhedraSide() override
friend struct sG4PolyhedraSideVec
G4VCSGface * Clone() override
G4SolidExtentList is utility class designed for calculating the extent of a CSG-like solid for a voxe...
virtual G4ThreeVector GetPointOnFace()=0
G4VCSGface()=default
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
EAxis
Definition geomdefs.hh:54
EInside
Definition geomdefs.hh:67
#define G4GEOM_DLL
Definition geomwdefs.hh:45
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668