Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GeomTools.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// G4GeomTools
27//
28// Class description:
29//
30// A collection of utilities which can be helpful for a wide range
31// of geometry-related tasks
32
33// Author: Evgueni Tcherniaev (CERN), 10.10.2016
34// --------------------------------------------------------------------
35#ifndef G4GEOMTOOLS_HH
36#define G4GEOMTOOLS_HH
37
38#include <vector>
39#include "G4TwoVector.hh"
40#include "G4ThreeVector.hh"
41
42using G4TwoVectorList = std::vector<G4TwoVector>;
43using G4ThreeVectorList = std::vector<G4ThreeVector>;
44
45/**
46 * @brief G4GeomTools is a collecting utilities which can be helpful
47 * for a wide range of geometry-related tasks.
48 */
49
51{
52 public:
53
54 // ==================================================================
55 // 2D Utilities
56 // ------------------------------------------------------------------
57
58 /**
59 * Functions to calculate the area of 2D triangle, returned value is
60 * positive if the vertices of the triangle are given in anticlockwise
61 * order, otherwise it is negative.
62 */
64 G4double Bx, G4double By,
65 G4double Cx, G4double Cy);
66 static G4double TriangleArea(const G4TwoVector& A,
67 const G4TwoVector& B,
68 const G4TwoVector& C);
69
70 /**
71 * Calculates the area of a 2D quadrilateral, returned value is positive if
72 * the vertices of the quadrilateral are given in anticlockwise order,
73 * otherwise it is negative.
74 */
75 static G4double QuadArea(const G4TwoVector& A,
76 const G4TwoVector& B,
77 const G4TwoVector& C,
78 const G4TwoVector& D);
79
80 /**
81 * Calculates the area of a 2D polygon, returned value is positive if
82 * the vertices of the polygon are defined in anticlockwise order,
83 * otherwise it is negative.
84 */
85 static G4double PolygonArea(const G4TwoVectorList& polygon);
86
87 /**
88 * Decides if a point (Px,Py) is inside the triangle (Ax,Ay)(Bx,By)(Cx,Cy).
89 */
91 G4double Ax, G4double Ay,
92 G4double Bx, G4double By,
93 G4double Cx, G4double Cy);
94
95 /**
96 * Decides if a point P is inside the triangle ABC.
97 */
98 static G4bool PointInTriangle(const G4TwoVector& P,
99 const G4TwoVector& A,
100 const G4TwoVector& B,
101 const G4TwoVector& C);
102
103 /**
104 * Decides if a point P is inside the 'Polygon'.
105 */
106 static G4bool PointInPolygon(const G4TwoVector& P,
107 const G4TwoVectorList& Polygon);
108
109 /**
110 * Decides if a 2D 'polygon' is convex, i.e. if all internal angles are
111 * less than pi.
112 */
113 static G4bool IsConvex(const G4TwoVectorList& polygon);
114
115 /**
116 * Simple implementation of "ear clipping" algorithm for triangulation
117 * of a simple contour/polygon, it places results in a std::vector as
118 * triplets of vertices. If triangulation is successful the function
119 * returns true, otherwise false.
120 */
121 static G4bool TriangulatePolygon(const G4TwoVectorList& polygon,
122 std::vector<G4int>& result);
123
124 /**
125 * Same using the function above and returning as 'result' a list
126 * of triangles.
127 */
128 static G4bool TriangulatePolygon(const G4TwoVectorList& polygon,
129 G4TwoVectorList& result);
130
131 /**
132 * Removes collinear and coincident points from a 2D 'polygon'.
133 * Indices of removed points are available in 'iout'.
134 * Allows to specify a 'tolerance'.
135 */
136 static void RemoveRedundantVertices(G4TwoVectorList& polygon,
137 std::vector<G4int>& iout,
138 G4double tolerance = 0.0);
139
140 /**
141 * Calculates the bounding rectangle of a disk sector. It returns false
142 * if the input parameters do not meet the following criteria:
143 * rmin >= 0
144 * rmax > rmin + kCarTolerance
145 * delPhi > 0 + kCarTolerance.
146 */
147 static G4bool DiskExtent(G4double rmin, G4double rmax,
148 G4double startPhi, G4double delPhi,
149 G4TwoVector& pmin, G4TwoVector& pmax);
150
151 /**
152 * Calculates the bounding rectangle of a disk sector.
153 * Faster version without check of parameters.
154 */
155 static void DiskExtent(G4double rmin, G4double rmax,
156 G4double sinPhiStart, G4double cosPhiStart,
157 G4double sinPhiEnd, G4double cosPhiEnd,
158 G4TwoVector& pmin, G4TwoVector& pmax);
159
160 /**
161 * Computes the circumference (perimeter) of an ellipse.
162 */
164 G4double b);
165
166 /**
167 * Computes the lateral surface area of an elliptic cone.
168 */
170 G4double b,
171 G4double h);
172
173 // ==================================================================
174 // 3D Utilities
175 // ------------------------------------------------------------------
176
177 /**
178 * Finds the normal to the plane of a 3D triangle ABC;
179 * the length of the normal is equal to the area of the triangle.
180 */
182 const G4ThreeVector& B,
183 const G4ThreeVector& C);
184
185 /**
186 * Finds the normal to the plane of a 3D quadrilateral ABCD;
187 * the length of the normal is equal to the area of the quadrilateral.
188 */
190 const G4ThreeVector& B,
191 const G4ThreeVector& C,
192 const G4ThreeVector& D);
193
194 /**
195 * Finds the normal to the plane of a 3D polygon; the length of the
196 * normal is equal to the area of the polygon.
197 */
199
200 /**
201 * Calculates the distance between a point 'P' and line segment AB in 3D.
202 */
204 const G4ThreeVector& A,
205 const G4ThreeVector& B);
206
207 /**
208 * Finds a point on a 3D line segment AB closest to point 'P'.
209 */
211 const G4ThreeVector& A,
212 const G4ThreeVector& B);
213
214 /**
215 * Finds a point on a 3D triangle ABC closest to point 'P'.
216 */
218 const G4ThreeVector& A,
219 const G4ThreeVector& B,
220 const G4ThreeVector& C);
221
222 /**
223 * Calculates the bounding box of a spherical sector,
224 * @returns false if input parameters do not meet the following criteria:
225 * rmin >= 0
226 * rmax > rmin + kCarTolerance
227 * startTheta >= 0 && <= pi;
228 * delTheta > 0 + kCarTolerance
229 * delPhi > 0 + kCarTolerance.
230 */
231 static G4bool SphereExtent(G4double rmin, G4double rmax,
232 G4double startTheta, G4double delTheta,
233 G4double startPhi, G4double delPhi,
234 G4ThreeVector& pmin, G4ThreeVector& pmax);
235
236 /**
237 * Calculates the hyperbolic surface stereo. Stereo is a half angle at the
238 * intersection point of the two lines in the tangent plane cross-section.
239 */
240 static G4double HypeStereo(G4double r0, // radius at z = 0
241 G4double r, // radius at z = h
242 G4double h);
243
244 /**
245 * Finds the XY-coordinates of the corners of a generic trap that bounds
246 * specified twisted tube.
247 * @param[in] twistAng The twist angle.
248 * @param[in] endInnerRad The inner radius at z = halfZ.
249 * @param[in] endOuterRad The outer radius at z = halfZ.
250 * @param[in] dPhi Delta phi.
251 * @param[in] vertices The corners of the generic trap.
252 */
253 static void TwistedTubeBoundingTrap(G4double twistAng,
254 G4double endInnerRad,
255 G4double endOuterRad,
256 G4double dPhi,
257 G4TwoVectorList& vertices);
258
259 /**
260 * Calculate surface area of the hyperboloid between 'zmin' and 'zmax'.
261 * @param[in] dphi Delta phi.
262 * @param[in] r0 The radius at z = 0.
263 * @param[in] tanstereo The tangent of the stereo angle.
264 * @param[in] zmin Minimum Z.
265 * @param[in] zmax Maximum Z.
266 */
268 G4double r0,
269 G4double tanstereo,
270 G4double zmin,
271 G4double zmax);
272
273 private:
274
275 /**
276 * Helper function for use by TriangulatePolygon().
277 */
278 static G4bool CheckSnip(const G4TwoVectorList& contour,
279 G4int a, G4int b, G4int c,
280 G4int n, const G4int* V);
281
282 /**
283 * Complete Elliptic Integral of the Second Kind.
284 */
285 static G4double comp_ellint_2(G4double e);
286};
287
288#endif
std::vector< G4ThreeVector > G4ThreeVectorList
G4double C(G4double temp)
G4double B(G4double temperature)
G4double D(G4double temp)
std::vector< G4ThreeVector > G4ThreeVectorList
std::vector< G4TwoVector > G4TwoVectorList
CLHEP::Hep3Vector G4ThreeVector
CLHEP::Hep2Vector G4TwoVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4GeomTools is a collecting utilities which can be helpful for a wide range of geometry-related tasks...
static G4bool PointInPolygon(const G4TwoVector &P, const G4TwoVectorList &Polygon)
static G4ThreeVector QuadAreaNormal(const G4ThreeVector &A, const G4ThreeVector &B, const G4ThreeVector &C, const G4ThreeVector &D)
static G4bool SphereExtent(G4double rmin, G4double rmax, G4double startTheta, G4double delTheta, G4double startPhi, G4double delPhi, G4ThreeVector &pmin, G4ThreeVector &pmax)
static G4double HypeStereo(G4double r0, G4double r, G4double h)
static void TwistedTubeBoundingTrap(G4double twistAng, G4double endInnerRad, G4double endOuterRad, G4double dPhi, G4TwoVectorList &vertices)
static G4bool PointInTriangle(G4double Px, G4double Py, G4double Ax, G4double Ay, G4double Bx, G4double By, G4double Cx, G4double Cy)
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
static G4double HyperboloidSurfaceArea(G4double dphi, G4double r0, G4double tanstereo, G4double zmin, G4double zmax)
static G4double EllipticConeLateralArea(G4double a, G4double b, G4double h)
static G4double DistancePointSegment(const G4ThreeVector &P, const G4ThreeVector &A, const G4ThreeVector &B)
static G4ThreeVector ClosestPointOnTriangle(const G4ThreeVector &P, const G4ThreeVector &A, const G4ThreeVector &B, const G4ThreeVector &C)
static void RemoveRedundantVertices(G4TwoVectorList &polygon, std::vector< G4int > &iout, G4double tolerance=0.0)
static G4double TriangleArea(G4double Ax, G4double Ay, G4double Bx, G4double By, G4double Cx, G4double Cy)
static G4double QuadArea(const G4TwoVector &A, const G4TwoVector &B, const G4TwoVector &C, const G4TwoVector &D)
static G4double EllipsePerimeter(G4double a, G4double b)
static G4double PolygonArea(const G4TwoVectorList &polygon)
static G4bool TriangulatePolygon(const G4TwoVectorList &polygon, std::vector< G4int > &result)
static G4ThreeVector ClosestPointOnSegment(const G4ThreeVector &P, const G4ThreeVector &A, const G4ThreeVector &B)
static G4bool IsConvex(const G4TwoVectorList &polygon)
static G4ThreeVector PolygonAreaNormal(const G4ThreeVectorList &polygon)
static G4ThreeVector TriangleAreaNormal(const G4ThreeVector &A, const G4ThreeVector &B, const G4ThreeVector &C)