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
G4PolyPhiFace.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// G4PolyPhiFace
27//
28// Class description:
29//
30// Definition of a face that bounds a polycone or polyhedra when
31// it has a phi opening:
32//
33// G4PolyPhiFace( const G4ReduciblePolygon* rz,
34// G4double phi,
35// G4double deltaPhi,
36// G4double phiOther )
37//
38// Specifically: a face that lies on a plane that passes through
39// the z axis. It has boundaries that are straight lines of arbitrary
40// length and direction, but with corners aways on the same side of
41// the z axis.
42
43// Author: David C. Williams (davidw@scipp.ucsc.edu)
44// --------------------------------------------------------------------
45#ifndef G4POLYPHIFACE_HH
46#define G4POLYPHIFACE_HH
47
48#include "G4VCSGface.hh"
49#include "G4TwoVector.hh"
50
52
54{
55 G4double x, y, r, z; // position
57 zNorm; // r/z normal
58 G4ThreeVector norm3D; // 3D normal
59
60 // Needed for Triangulation Algorithm
61 //
64};
65
67{
68 G4PolyPhiFaceEdge(): v0(0), v1(0), tr(.0), tz(0.), length(0.) {}
69 G4PolyPhiFaceVertex *v0, *v1; // Corners
70 G4double tr, tz, // Unit vector along edge
71 length; // Length of edge
72 G4ThreeVector norm3D; // 3D edge normal vector
73};
74
76{
77
78 public: // with description
79
81 G4double phi, G4double deltaPhi, G4double phiOther );
82 // Constructor.
83 // Points r,z should be supplied in clockwise order in r,z.
84 // For example:
85 // [1]---------[2] ^ R
86 // | | |
87 // | | +--> z
88 // [0]---------[3]
89
90 virtual ~G4PolyPhiFace();
91 // Destructor. Removes edges and corners.
92
93 G4PolyPhiFace( const G4PolyPhiFace &source );
94 G4PolyPhiFace& operator=( const G4PolyPhiFace &source );
95 // Copy constructor and assgnment operator.
96
97 G4bool Intersect( const G4ThreeVector& p, const G4ThreeVector& v,
98 G4bool outgoing, G4double surfTolerance,
99 G4double& distance, G4double& distFromSurface,
101
102 G4double Distance( const G4ThreeVector& p, G4bool outgoing );
103
104 EInside Inside( const G4ThreeVector& p, G4double tolerance,
105 G4double* bestDistance );
106
107 G4ThreeVector Normal( const G4ThreeVector& p, G4double* bestDistance );
108
109 G4double Extent( const G4ThreeVector axis );
110
111 void CalculateExtent( const EAxis axis,
112 const G4VoxelLimits &voxelLimit,
113 const G4AffineTransform& tranform,
114 G4SolidExtentList& extentList );
115
116 inline G4VCSGface* Clone();
117 // Allocates on the heap a clone of this face.
118
123 // Auxiliary methods for determination of points on surface.
124
125 public: // without description
126
127 G4PolyPhiFace(__void__&);
128 // Fake default constructor for usage restricted to direct object
129 // persistency for clients requiring preallocation of memory for
130 // persistifiable objects.
131
132 void Diagnose( G4VSolid* solid );
133 // Throw an exception if something is found inconsistent with
134 // the solid. For debugging purposes only
135
136 protected:
137
139 const G4ThreeVector& p, const G4ThreeVector& v );
140 // Decide if the point in r,z is inside the edges of our face,
141 // **but** do so consistently with other faces.
142
145 G4PolyPhiFaceVertex** base3Dnorm = nullptr,
146 G4ThreeVector** head3Dnorm = nullptr );
147 // Decide if the point in r,z is inside the edges of our face.
148
150 G4double qx, G4double qy, G4double qz,
151 const G4ThreeVector& v,
152 G4double normSign,
153 const G4PolyPhiFaceVertex* vert ) const;
154 // Decide precisely whether a trajectory passes to the left, right,
155 // or exactly passes through the z position of a vertex point in face.
156
157 void CopyStuff( const G4PolyPhiFace& source );
158
159 protected:
160
161 // Functions used for Triangulation in Case of generic Polygone.
162 // The triangulation is used for GetPointOnFace()
163
165 // Calculation of 2*Area of Triangle with Sign
166
170 // Boolean functions for sign of Surface
171
174 // Boolean function for finding proper intersection of two
175 // line segments (a,b) and (c,d).
176
178 // Boolean function for determining if point c is between a and b
179 // where the three points (a,b,c) are on the same line.
180
183 // Boolean function for finding proper intersection or not
184 // of two line segments (a,b) and (c,d).
185
187 // Boolean Diagonalie help to determine if diagonal s
188 // of segment (a,b) is convex or reflex.
189
191 // Boolean function for determining if b is inside the cone (a0,a,a1)
192 // where a is the center of the cone.
193
195 // Boolean function for determining if Diagonal is possible
196 // inside Polycone or PolyHedra.
197
198 void EarInit();
199 // Initialisation for Triangulisation by ear tips.
200 // For details see "Computational Geometry in C" by Joseph O'Rourke.
201
202 void Triangulate();
203 // Triangularisation by ear tips for Polycone or Polyhedra.
204 // For details see "Computational Geometry in C" by Joseph O'Rourke.
205 // NOTE: a copy of the shape is made and this copy is reordered in
206 // order to have a list of triangles. This list is used by the
207 // method GetPointOnFace().
208
209 protected:
210
211 G4int numEdges; // Number of edges
212 G4PolyPhiFaceEdge* edges = nullptr; // The edges of the face
213 G4PolyPhiFaceVertex* corners = nullptr; // And the corners
214 G4ThreeVector normal; // Normal unit vector
215 G4ThreeVector radial; // Unit vector along radial direction
216 G4ThreeVector surface; // Point on surface
217 G4ThreeVector surface_point; // Auxiliary point on surface used for
218 // method GetPointOnFace()
219 G4double rMin, rMax, // Extent in r
220 zMin, zMax; // Extent in z
221 G4bool allBehind = false; // True if the polycone/polyhedra
222 // is behind the place of this face
223 G4double kCarTolerance; // Surface thickness
224 G4double fSurfaceArea = 0.0; // Surface Area of PolyPhiFace
226 // Auxiliary pointer to 'corners' used for triangulation.
227 // Copy structure, changing the structure of 'corners' (ear removal)
228};
229
230#include "G4PolyPhiFace.icc"
231
232#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double SurfaceArea()
G4bool Left(G4TwoVector a, G4TwoVector b, G4TwoVector c)
G4bool Collinear(G4TwoVector a, G4TwoVector b, G4TwoVector c)
void Diagnose(G4VSolid *solid)
G4PolyPhiFace & operator=(const G4PolyPhiFace &source)
virtual ~G4PolyPhiFace()
G4bool Diagonal(G4PolyPhiFaceVertex *a, G4PolyPhiFaceVertex *b)
G4ThreeVector surface_point
G4bool InsideEdges(G4double r, G4double z)
G4bool Diagonalie(G4PolyPhiFaceVertex *a, G4PolyPhiFaceVertex *b)
G4VCSGface * Clone()
G4bool IntersectProp(G4TwoVector a, G4TwoVector b, G4TwoVector c, G4TwoVector d)
G4ThreeVector surface
G4bool LeftOn(G4TwoVector a, G4TwoVector b, G4TwoVector c)
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, G4bool outgoing, G4double surfTolerance, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal, G4bool &allBehind)
G4ThreeVector Normal(const G4ThreeVector &p, G4double *bestDistance)
G4double fSurfaceArea
void CopyStuff(const G4PolyPhiFace &source)
void CalculateExtent(const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList)
G4double Area2(G4TwoVector a, G4TwoVector b, G4TwoVector c)
G4double ExactZOrder(G4double z, G4double qx, G4double qy, G4double qz, const G4ThreeVector &v, G4double normSign, const G4PolyPhiFaceVertex *vert) const
G4double Distance(const G4ThreeVector &p, G4bool outgoing)
G4double Extent(const G4ThreeVector axis)
G4PolyPhiFaceEdge * edges
G4ThreeVector radial
G4bool InCone(G4PolyPhiFaceVertex *a, G4PolyPhiFaceVertex *b)
G4bool Between(G4TwoVector a, G4TwoVector b, G4TwoVector c)
G4PolyPhiFaceVertex * corners
G4PolyPhiFaceVertex * triangles
G4ThreeVector normal
EInside Inside(const G4ThreeVector &p, G4double tolerance, G4double *bestDistance)
G4double kCarTolerance
G4ThreeVector GetPointOnFace()
G4double SurfaceTriangle(G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector *p4)
G4bool InsideEdgesExact(G4double r, G4double z, G4double normSign, const G4ThreeVector &p, const G4ThreeVector &v)
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67
G4PolyPhiFaceVertex * v1
G4ThreeVector norm3D
G4PolyPhiFaceVertex * v0
G4ThreeVector norm3D
G4PolyPhiFaceVertex * next
G4PolyPhiFaceVertex * prev