Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BoundingEnvelope.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// class G4BoundingEnvelope
27//
28// Class description:
29//
30// Helper class to facilitate calculation of the extent of a solid
31// within the limits defined by a G4VoxelLimits object.
32//
33// The function CalculateExtent() of a particular solid can create
34// a G4BoundingEnvelope object that bounds the solid and then call
35// CalculateExtent() of the G4BoundingEnvelope object.
36//
37// Calculation of extent uses G4Transform3D, thus takes into account
38// scaling and reflection, if any.
39
40// 2016.05.25, E.Tcherniaev - initial version
41// --------------------------------------------------------------------
42#ifndef G4BOUNDINGENVELOPE_HH
43#define G4BOUNDINGENVELOPE_HH
44
45#include <vector>
46#include "geomdefs.hh"
47
48#include "G4ThreeVector.hh"
49#include "G4VoxelLimits.hh"
50#include "G4Transform3D.hh"
51#include "G4Point3D.hh"
52#include "G4Plane3D.hh"
53
54using G4ThreeVectorList = std::vector<G4ThreeVector>;
55using G4Polygon3D = std::vector<G4Point3D>;
56using G4Segment3D = std::pair<G4Point3D,G4Point3D>;
57
59{
60 public:
61
63 const G4ThreeVector& pMax);
64 // Constructor from an axis aligned bounding box (AABB)
65
66 G4BoundingEnvelope(const std::vector<const G4ThreeVectorList*>& polygons);
67 // Constructor from a sequence of convex polygons, the polygons
68 // should have equal numbers of vertices except first and last
69 // polygons which may consist of a single vertex
70
72 const G4ThreeVector& pMax,
73 const std::vector<const G4ThreeVectorList*>& polygons);
74 // Constructor from AABB and a sequence of polygons
75
77 // Destructor
78
80 const G4VoxelLimits& pVoxelLimits,
81 const G4Transform3D& pTransform3D,
82 G4double& pMin, G4double& pMax) const;
83 // Analyse the position of the bounding box relative to the voxel.
84 // It returns "true" in the case where the value of the extent can be
85 // figured out directly from the dimensions of the bounding box, or
86 // it is clear that the bounding box and the voxel do not intersect.
87 // The reply "false" means that further calculations are needed.
88
89 G4bool CalculateExtent(const EAxis pAxis,
90 const G4VoxelLimits& pVoxelLimits,
91 const G4Transform3D& pTransform3D,
92 G4double& pMin, G4double& pMax) const;
93 // Calculate extent of the bounding envelope
94
95 private:
96
97 void CheckBoundingBox();
98 // Check correctness of the AABB (axis aligned bounding box)
99
100 void CheckBoundingPolygons();
101 // Check correctness of the sequence of convex polygonal bases
102
103 G4double FindScaleFactor(const G4Transform3D& pTransform3D) const;
104 // Find max scale factor of the transformation
105
106 void TransformVertices(const G4Transform3D& pTransform3D,
107 std::vector<G4Polygon3D*>& pBases) const;
108 // Create list of transformed polygons
109
110 void GetPrismAABB(const G4Polygon3D& pBaseA,
111 const G4Polygon3D& pBaseB,
112 G4Segment3D& pAABB) const;
113 // Find bounding box of a prism
114
115 void CreateListOfEdges(const G4Polygon3D& baseA,
116 const G4Polygon3D& baseB,
117 std::vector<G4Segment3D>& pEdges) const;
118 // Create list of edges of a prism
119
120 void CreateListOfPlanes(const G4Polygon3D& baseA,
121 const G4Polygon3D& baseB,
122 std::vector<G4Plane3D>& pPlanes) const;
123 // Create list of planes bounding a prism
124
125 G4bool ClipEdgesByVoxel(const std::vector<G4Segment3D>& pEdges,
126 const G4VoxelLimits& pLimits,
127 G4Segment3D& pExtent) const;
128 // Clip set of edges by G4VoxelLimits
129
130 void ClipVoxelByPlanes(G4int pBits,
131 const G4VoxelLimits& pLimits,
132 const std::vector<G4Plane3D>& pPlanes,
133 const G4Segment3D& pAABB,
134 G4Segment3D& pExtent) const;
135 // Clip G4VoxelLimits by set of planes bounding a prism
136
137 private:
138
139 G4ThreeVector fMin, fMax;
140 // original bounding box
141
142 const std::vector<const G4ThreeVectorList*>* fPolygons = nullptr;
143 // ref to original sequence of polygonal bases
144};
145
146#endif // G4BOUNDINGENVELOPE_HH
std::pair< G4Point3D, G4Point3D > G4Segment3D
std::vector< G4Point3D > G4Polygon3D
std::vector< G4ThreeVector > G4ThreeVectorList
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4bool BoundingBoxVsVoxelLimits(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
EAxis
Definition: geomdefs.hh:54