|
||||
File indexing completed on 2025-01-18 09:57:59
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 // class G4BoundingEnvelope 0027 // 0028 // Class description: 0029 // 0030 // Helper class to facilitate calculation of the extent of a solid 0031 // within the limits defined by a G4VoxelLimits object. 0032 // 0033 // The function CalculateExtent() of a particular solid can create 0034 // a G4BoundingEnvelope object that bounds the solid and then call 0035 // CalculateExtent() of the G4BoundingEnvelope object. 0036 // 0037 // Calculation of extent uses G4Transform3D, thus takes into account 0038 // scaling and reflection, if any. 0039 0040 // 2016.05.25, E.Tcherniaev - initial version 0041 // -------------------------------------------------------------------- 0042 #ifndef G4BOUNDINGENVELOPE_HH 0043 #define G4BOUNDINGENVELOPE_HH 0044 0045 #include <vector> 0046 #include "geomdefs.hh" 0047 0048 #include "G4ThreeVector.hh" 0049 #include "G4VoxelLimits.hh" 0050 #include "G4Transform3D.hh" 0051 #include "G4Point3D.hh" 0052 #include "G4Plane3D.hh" 0053 0054 using G4ThreeVectorList = std::vector<G4ThreeVector>; 0055 using G4Polygon3D = std::vector<G4Point3D>; 0056 using G4Segment3D = std::pair<G4Point3D,G4Point3D>; 0057 0058 class G4BoundingEnvelope 0059 { 0060 public: 0061 0062 G4BoundingEnvelope(const G4ThreeVector& pMin, 0063 const G4ThreeVector& pMax); 0064 // Constructor from an axis aligned bounding box (AABB) 0065 0066 G4BoundingEnvelope(const std::vector<const G4ThreeVectorList*>& polygons); 0067 // Constructor from a sequence of convex polygons, the polygons 0068 // should have equal numbers of vertices except first and last 0069 // polygons which may consist of a single vertex 0070 0071 G4BoundingEnvelope(const G4ThreeVector& pMin, 0072 const G4ThreeVector& pMax, 0073 const std::vector<const G4ThreeVectorList*>& polygons); 0074 // Constructor from AABB and a sequence of polygons 0075 0076 ~G4BoundingEnvelope() = default; 0077 // Destructor 0078 0079 G4bool BoundingBoxVsVoxelLimits(const EAxis pAxis, 0080 const G4VoxelLimits& pVoxelLimits, 0081 const G4Transform3D& pTransform3D, 0082 G4double& pMin, G4double& pMax) const; 0083 // Analyse the position of the bounding box relative to the voxel. 0084 // It returns "true" in the case where the value of the extent can be 0085 // figured out directly from the dimensions of the bounding box, or 0086 // it is clear that the bounding box and the voxel do not intersect. 0087 // The reply "false" means that further calculations are needed. 0088 0089 G4bool CalculateExtent(const EAxis pAxis, 0090 const G4VoxelLimits& pVoxelLimits, 0091 const G4Transform3D& pTransform3D, 0092 G4double& pMin, G4double& pMax) const; 0093 // Calculate extent of the bounding envelope 0094 0095 private: 0096 0097 void CheckBoundingBox(); 0098 // Check correctness of the AABB (axis aligned bounding box) 0099 0100 void CheckBoundingPolygons(); 0101 // Check correctness of the sequence of convex polygonal bases 0102 0103 G4double FindScaleFactor(const G4Transform3D& pTransform3D) const; 0104 // Find max scale factor of the transformation 0105 0106 void TransformVertices(const G4Transform3D& pTransform3D, 0107 std::vector<G4Point3D>& pVertices, 0108 std::vector<std::pair<G4int,G4int>>& pBases) const; 0109 // Create list of transformed polygons 0110 0111 void GetPrismAABB(const G4Polygon3D& pBaseA, 0112 const G4Polygon3D& pBaseB, 0113 G4Segment3D& pAABB) const; 0114 // Find bounding box of a prism 0115 0116 void CreateListOfEdges(const G4Polygon3D& baseA, 0117 const G4Polygon3D& baseB, 0118 std::vector<G4Segment3D>& pEdges) const; 0119 // Create list of edges of a prism 0120 0121 void CreateListOfPlanes(const G4Polygon3D& baseA, 0122 const G4Polygon3D& baseB, 0123 std::vector<G4Plane3D>& pPlanes) const; 0124 // Create list of planes bounding a prism 0125 0126 G4bool ClipEdgesByVoxel(const std::vector<G4Segment3D>& pEdges, 0127 const G4VoxelLimits& pLimits, 0128 G4Segment3D& pExtent) const; 0129 // Clip set of edges by G4VoxelLimits 0130 0131 void ClipVoxelByPlanes(G4int pBits, 0132 const G4VoxelLimits& pLimits, 0133 const std::vector<G4Plane3D>& pPlanes, 0134 const G4Segment3D& pAABB, 0135 G4Segment3D& pExtent) const; 0136 // Clip G4VoxelLimits by set of planes bounding a prism 0137 0138 private: 0139 0140 G4ThreeVector fMin, fMax; 0141 // original bounding box 0142 0143 const std::vector<const G4ThreeVectorList*>* fPolygons = nullptr; 0144 // ref to original sequence of polygonal bases 0145 }; 0146 0147 #endif // G4BOUNDINGENVELOPE_HH
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |