Geant4 11.2.2
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 helpfull for a wide range
31// of geometry-related tasks
32
33// 10.10.2016, E.Tcherniaev: Initial version.
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
46{
47 public:
48
49 // ==================================================================
50 // 2D Utilities
51 // ------------------------------------------------------------------
52
54 G4double Bx, G4double By,
55 G4double Cx, G4double Cy);
56
57 static G4double TriangleArea(const G4TwoVector& A,
58 const G4TwoVector& B,
59 const G4TwoVector& C);
60 // Calculate area of 2D triangle, return value is positive if
61 // vertices of the triangle are given in anticlockwise order,
62 // otherwise it is negative
63
64 static G4double QuadArea(const G4TwoVector& A,
65 const G4TwoVector& B,
66 const G4TwoVector& C,
67 const G4TwoVector& D);
68 // Calculate area of 2D quadrilateral, return value is positive if
69 // vertices of the quadrilateral are given in anticlockwise order,
70 // otherwise it is negative
71
72 static G4double PolygonArea(const G4TwoVectorList& polygon);
73 // Calculate area of 2D polygon, return value is positive if
74 // vertices of the polygon are defined in anticlockwise order,
75 // otherwise it is negative
76
78 G4double Ax, G4double Ay,
79 G4double Bx, G4double By,
80 G4double Cx, G4double Cy);
81 // Decide if point (Px,Py) is inside triangle (Ax,Ay)(Bx,By)(Cx,Cy)
82
83 static G4bool PointInTriangle(const G4TwoVector& P,
84 const G4TwoVector& A,
85 const G4TwoVector& B,
86 const G4TwoVector& C);
87 // Decide if point P is inside triangle ABC
88
89 static G4bool PointInPolygon(const G4TwoVector& P,
90 const G4TwoVectorList& Polygon);
91 // Decide if point P is inside Polygon
92
93 static G4bool IsConvex(const G4TwoVectorList& polygon);
94 // Decide if 2D polygon is convex, i.e. all internal angles are
95 // less than pi
96
97 static G4bool TriangulatePolygon(const G4TwoVectorList& polygon,
98 G4TwoVectorList& result);
99
100 static G4bool TriangulatePolygon(const G4TwoVectorList& polygon,
101 std::vector<G4int>& result);
102 // Simple implementation of "ear clipping" algorithm for
103 // triangulation of a simple contour/polygon, it places results
104 // in a std::vector as triplets of vertices. If triangulation
105 // is sucsessfull then the function returns true, otherwise false
106
107 static void RemoveRedundantVertices(G4TwoVectorList& polygon,
108 std::vector<G4int>& iout,
109 G4double tolerance = 0.0);
110 // Remove collinear and coincident points from 2D polygon.
111 // Indices of removed points are available in iout.
112
113 static G4bool DiskExtent(G4double rmin, G4double rmax,
114 G4double startPhi, G4double delPhi,
115 G4TwoVector& pmin, G4TwoVector& pmax);
116 // Calculate bounding rectangle of a disk sector,
117 // it returns false if input parameters do not meet the following:
118 // rmin >= 0
119 // rmax > rmin + kCarTolerance
120 // delPhi > 0 + kCarTolerance
121
122 static void DiskExtent(G4double rmin, G4double rmax,
123 G4double sinPhiStart, G4double cosPhiStart,
124 G4double sinPhiEnd, G4double cosPhiEnd,
125 G4TwoVector& pmin, G4TwoVector& pmax);
126 // Calculate bounding rectangle of a disk sector,
127 // faster version without check of parameters
128
130 G4double b);
131 // Compute the circumference (perimeter) of an ellipse
132
134 G4double b,
135 G4double h);
136 // Compute the lateral surface area of an elliptic cone
137
138 // ==================================================================
139 // 3D Utilities
140 // ------------------------------------------------------------------
141
143 const G4ThreeVector& B,
144 const G4ThreeVector& C);
145 // Find the normal to the plane of 3D triangle ABC,
146 // length of the normal is equal to the area of the triangle
147
149 const G4ThreeVector& B,
150 const G4ThreeVector& C,
151 const G4ThreeVector& D);
152 // Find normal to the plane of 3D quadrilateral ABCD,
153 // length of the normal is equal to the area of the quadrilateral
154
156 // Find normal to the plane of 3D polygon
157 // length of the normal is equal to the area of the polygon
158
159 /*
160 static G4bool IsPlanar(const G4ThreeVector& A,
161 const G4ThreeVector& B,
162 const G4ThreeVector& C,
163 const G4ThreeVector& D);
164 // Decide if 3D quadrilateral ABCD is planar
165
166 static G4bool IsPlanar(const G4ThreeVectorList& polygon
167 const G4ThreeVector& normal);
168 // Decide if 3D polygon is planar
169 */
170
172 const G4ThreeVector& A,
173 const G4ThreeVector& B);
174 // Calculate distance between point P and line segment AB in 3D
175
177 const G4ThreeVector& A,
178 const G4ThreeVector& B);
179 // Find point on 3D line segment AB closest to point P
180
182 const G4ThreeVector& A,
183 const G4ThreeVector& B,
184 const G4ThreeVector& C);
185 // Find point on 3D triangle ABC closest to point P
186
187 static G4bool SphereExtent(G4double rmin, G4double rmax,
188 G4double startTheta, G4double delTheta,
189 G4double startPhi, G4double delPhi,
190 G4ThreeVector& pmin, G4ThreeVector& pmax);
191 // Calculate bounding box of a spherical sector,
192 // it returns false if input parameters do not meet the following:
193 // rmin >= 0
194 // rmax > rmin + kCarTolerance
195 // startTheta >= 0 && <= pi;
196 // delTheta > 0 + kCarTolerance
197 // delPhi > 0 + kCarTolerance
198
199 private:
200
201 static G4bool CheckSnip(const G4TwoVectorList& contour,
202 G4int a, G4int b, G4int c,
203 G4int n, const G4int* V);
204 // Helper function for TriangulatePolygon()
205
206 static G4double comp_ellint_2(G4double e);
207 // Complete Elliptic Integral of the Second Kind
208};
209
210#endif // G4GEOMTOOLS_HH
std::vector< G4ThreeVector > G4ThreeVectorList
G4double C(G4double temp)
G4double B(G4double temperature)
G4double D(G4double temp)
std::vector< G4ThreeVector > G4ThreeVectorList
std::vector< G4TwoVector > G4TwoVectorList
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
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 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 G4bool TriangulatePolygon(const G4TwoVectorList &polygon, G4TwoVectorList &result)
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 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)