|
||||
File indexing completed on 2025-01-18 09:58:21
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // G4GeomTools 0027 // 0028 // Class description: 0029 // 0030 // A collection of utilities which can be helpfull for a wide range 0031 // of geometry-related tasks 0032 0033 // 10.10.2016, E.Tcherniaev: Initial version. 0034 // -------------------------------------------------------------------- 0035 #ifndef G4GEOMTOOLS_HH 0036 #define G4GEOMTOOLS_HH 0037 0038 #include <vector> 0039 #include "G4TwoVector.hh" 0040 #include "G4ThreeVector.hh" 0041 0042 using G4TwoVectorList = std::vector<G4TwoVector>; 0043 using G4ThreeVectorList = std::vector<G4ThreeVector>; 0044 0045 class G4GeomTools 0046 { 0047 public: 0048 0049 // ================================================================== 0050 // 2D Utilities 0051 // ------------------------------------------------------------------ 0052 0053 static G4double TriangleArea(G4double Ax, G4double Ay, 0054 G4double Bx, G4double By, 0055 G4double Cx, G4double Cy); 0056 0057 static G4double TriangleArea(const G4TwoVector& A, 0058 const G4TwoVector& B, 0059 const G4TwoVector& C); 0060 // Calculate area of 2D triangle, return value is positive if 0061 // vertices of the triangle are given in anticlockwise order, 0062 // otherwise it is negative 0063 0064 static G4double QuadArea(const G4TwoVector& A, 0065 const G4TwoVector& B, 0066 const G4TwoVector& C, 0067 const G4TwoVector& D); 0068 // Calculate area of 2D quadrilateral, return value is positive if 0069 // vertices of the quadrilateral are given in anticlockwise order, 0070 // otherwise it is negative 0071 0072 static G4double PolygonArea(const G4TwoVectorList& polygon); 0073 // Calculate area of 2D polygon, return value is positive if 0074 // vertices of the polygon are defined in anticlockwise order, 0075 // otherwise it is negative 0076 0077 static G4bool PointInTriangle(G4double Px, G4double Py, 0078 G4double Ax, G4double Ay, 0079 G4double Bx, G4double By, 0080 G4double Cx, G4double Cy); 0081 // Decide if point (Px,Py) is inside triangle (Ax,Ay)(Bx,By)(Cx,Cy) 0082 0083 static G4bool PointInTriangle(const G4TwoVector& P, 0084 const G4TwoVector& A, 0085 const G4TwoVector& B, 0086 const G4TwoVector& C); 0087 // Decide if point P is inside triangle ABC 0088 0089 static G4bool PointInPolygon(const G4TwoVector& P, 0090 const G4TwoVectorList& Polygon); 0091 // Decide if point P is inside Polygon 0092 0093 static G4bool IsConvex(const G4TwoVectorList& polygon); 0094 // Decide if 2D polygon is convex, i.e. all internal angles are 0095 // less than pi 0096 0097 static G4bool TriangulatePolygon(const G4TwoVectorList& polygon, 0098 G4TwoVectorList& result); 0099 0100 static G4bool TriangulatePolygon(const G4TwoVectorList& polygon, 0101 std::vector<G4int>& result); 0102 // Simple implementation of "ear clipping" algorithm for 0103 // triangulation of a simple contour/polygon, it places results 0104 // in a std::vector as triplets of vertices. If triangulation 0105 // is sucsessfull then the function returns true, otherwise false 0106 0107 static void RemoveRedundantVertices(G4TwoVectorList& polygon, 0108 std::vector<G4int>& iout, 0109 G4double tolerance = 0.0); 0110 // Remove collinear and coincident points from 2D polygon. 0111 // Indices of removed points are available in iout. 0112 0113 static G4bool DiskExtent(G4double rmin, G4double rmax, 0114 G4double startPhi, G4double delPhi, 0115 G4TwoVector& pmin, G4TwoVector& pmax); 0116 // Calculate bounding rectangle of a disk sector, 0117 // it returns false if input parameters do not meet the following: 0118 // rmin >= 0 0119 // rmax > rmin + kCarTolerance 0120 // delPhi > 0 + kCarTolerance 0121 0122 static void DiskExtent(G4double rmin, G4double rmax, 0123 G4double sinPhiStart, G4double cosPhiStart, 0124 G4double sinPhiEnd, G4double cosPhiEnd, 0125 G4TwoVector& pmin, G4TwoVector& pmax); 0126 // Calculate bounding rectangle of a disk sector, 0127 // faster version without check of parameters 0128 0129 static G4double EllipsePerimeter(G4double a, 0130 G4double b); 0131 // Compute the circumference (perimeter) of an ellipse 0132 0133 static G4double EllipticConeLateralArea(G4double a, 0134 G4double b, 0135 G4double h); 0136 // Compute the lateral surface area of an elliptic cone 0137 0138 // ================================================================== 0139 // 3D Utilities 0140 // ------------------------------------------------------------------ 0141 0142 static G4ThreeVector TriangleAreaNormal(const G4ThreeVector& A, 0143 const G4ThreeVector& B, 0144 const G4ThreeVector& C); 0145 // Find the normal to the plane of 3D triangle ABC, 0146 // length of the normal is equal to the area of the triangle 0147 0148 static G4ThreeVector QuadAreaNormal(const G4ThreeVector& A, 0149 const G4ThreeVector& B, 0150 const G4ThreeVector& C, 0151 const G4ThreeVector& D); 0152 // Find normal to the plane of 3D quadrilateral ABCD, 0153 // length of the normal is equal to the area of the quadrilateral 0154 0155 static G4ThreeVector PolygonAreaNormal(const G4ThreeVectorList& polygon); 0156 // Find normal to the plane of 3D polygon 0157 // length of the normal is equal to the area of the polygon 0158 0159 /* 0160 static G4bool IsPlanar(const G4ThreeVector& A, 0161 const G4ThreeVector& B, 0162 const G4ThreeVector& C, 0163 const G4ThreeVector& D); 0164 // Decide if 3D quadrilateral ABCD is planar 0165 0166 static G4bool IsPlanar(const G4ThreeVectorList& polygon 0167 const G4ThreeVector& normal); 0168 // Decide if 3D polygon is planar 0169 */ 0170 0171 static G4double DistancePointSegment(const G4ThreeVector& P, 0172 const G4ThreeVector& A, 0173 const G4ThreeVector& B); 0174 // Calculate distance between point P and line segment AB in 3D 0175 0176 static G4ThreeVector ClosestPointOnSegment(const G4ThreeVector& P, 0177 const G4ThreeVector& A, 0178 const G4ThreeVector& B); 0179 // Find point on 3D line segment AB closest to point P 0180 0181 static G4ThreeVector ClosestPointOnTriangle(const G4ThreeVector& P, 0182 const G4ThreeVector& A, 0183 const G4ThreeVector& B, 0184 const G4ThreeVector& C); 0185 // Find point on 3D triangle ABC closest to point P 0186 0187 static G4bool SphereExtent(G4double rmin, G4double rmax, 0188 G4double startTheta, G4double delTheta, 0189 G4double startPhi, G4double delPhi, 0190 G4ThreeVector& pmin, G4ThreeVector& pmax); 0191 // Calculate bounding box of a spherical sector, 0192 // it returns false if input parameters do not meet the following: 0193 // rmin >= 0 0194 // rmax > rmin + kCarTolerance 0195 // startTheta >= 0 && <= pi; 0196 // delTheta > 0 + kCarTolerance 0197 // delPhi > 0 + kCarTolerance 0198 0199 private: 0200 0201 static G4bool CheckSnip(const G4TwoVectorList& contour, 0202 G4int a, G4int b, G4int c, 0203 G4int n, const G4int* V); 0204 // Helper function for TriangulatePolygon() 0205 0206 static G4double comp_ellint_2(G4double e); 0207 // Complete Elliptic Integral of the Second Kind 0208 }; 0209 0210 #endif // G4GEOMTOOLS_HH
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |