Geant4 11.4.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// 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// Author: Evgueni Tcherniaev (CERN), 25.05.2016 - 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
58/**
59 * @brief G4BoundingEnvelope is a helper class to facilitate calculation of
60 * the extent of a solid within the limits defined by a G4VoxelLimits object.
61 * Calculation of the extent takes into account scaling and reflection, if any.
62 */
63
65{
66 public:
67
68 /**
69 * Constructor from an axis aligned bounding box (AABB).
70 * @param[in] pMin Lower boundary point.
71 * @param[in] pMax Higher boundary point.
72 */
74 const G4ThreeVector& pMax);
75
76 /**
77 * Constructor from a sequence of convex polygons, the polygons should have
78 * equal numbers of vertices except first and last polygons which may
79 * consist of a single vertex.
80 * @param[in] polygons The list of polygons.
81 */
82 G4BoundingEnvelope(const std::vector<const G4ThreeVectorList*>& polygons);
83
84 /**
85 * Constructor from AABB and a sequence of polygons.
86 * @param[in] pMin Lower boundary point.
87 * @param[in] pMax Higher boundary point.
88 * @param[in] polygons The list of polygons.
89 */
91 const G4ThreeVector& pMax,
92 const std::vector<const G4ThreeVectorList*>& polygons);
93
94 /**
95 * Default Destructor.
96 */
98
99 /**
100 * Analyses the position of the bounding box relative to the voxel.
101 * It returns "true" in the case where the value of the extent can be
102 * figured out directly from the dimensions of the bounding box, or
103 * it is clear that the bounding box and the voxel do not intersect.
104 * The reply "false" means that further calculations are needed.
105 * @param[in] pAxis The axis along which compute the extent.
106 * @param[in] pVoxelLimits The limiting space dictated by voxels.
107 * @param[in] pTransform3D Internal transformation applied to the envelope.
108 * @param[out] pMin The minimum extent value.
109 * @param[out] pMax The maximum extent value.
110 * @returns True if the envelope is intersected by the extent region.
111 */
113 const G4VoxelLimits& pVoxelLimits,
114 const G4Transform3D& pTransform3D,
115 G4double& pMin, G4double& pMax) const;
116
117 /**
118 * Calculates the minimum and maximum extent of the bounding envelope,
119 * when under the specified transform, and within the specified limits.
120 * @param[in] pAxis The axis along which compute the extent.
121 * @param[in] pVoxelLimits The limiting space dictated by voxels.
122 * @param[in] pTransform3D Internal transformation applied to the envelope.
123 * @param[out] pMin The minimum extent value.
124 * @param[out] pMax The maximum extent value.
125 * @returns True if the envelope is intersected by the extent region.
126 */
127 G4bool CalculateExtent(const EAxis pAxis,
128 const G4VoxelLimits& pVoxelLimits,
129 const G4Transform3D& pTransform3D,
130 G4double& pMin, G4double& pMax) const;
131
132 private:
133
134 /**
135 * Checks the correctness of the AABB (axis aligned bounding box).
136 */
137 void CheckBoundingBox();
138
139 /**
140 * Checks the correctness of the sequence of convex polygonal bases.
141 */
142 void CheckBoundingPolygons();
143
144 /**
145 * Finds the max scale factor of the transformation.
146 */
147 G4double FindScaleFactor(const G4Transform3D& pTransform3D) const;
148
149 /**
150 * Creates a list of transformed polygons.
151 */
152 void TransformVertices(const G4Transform3D& pTransform3D,
153 std::vector<G4Point3D>& pVertices,
154 std::vector<std::pair<G4int,G4int>>& pBases) const;
155
156 /**
157 * Finds the bounding box of a prism.
158 */
159 void GetPrismAABB(const G4Polygon3D& pBaseA,
160 const G4Polygon3D& pBaseB,
161 G4Segment3D& pAABB) const;
162
163 /**
164 * Creates a list of edges of a prism.
165 */
166 void CreateListOfEdges(const G4Polygon3D& baseA,
167 const G4Polygon3D& baseB,
168 std::vector<G4Segment3D>& pEdges) const;
169
170 /**
171 * Creates a list of planes bounding a prism.
172 */
173 void CreateListOfPlanes(const G4Polygon3D& baseA,
174 const G4Polygon3D& baseB,
175 std::vector<G4Plane3D>& pPlanes) const;
176
177 /**
178 * Clips a set of edges by G4VoxelLimits.
179 */
180 G4bool ClipEdgesByVoxel(const std::vector<G4Segment3D>& pEdges,
181 const G4VoxelLimits& pLimits,
182 G4Segment3D& pExtent) const;
183
184 /**
185 * Clips G4VoxelLimits by set of planes bounding a prism.
186 */
187 void ClipVoxelByPlanes(G4int pBits,
188 const G4VoxelLimits& pLimits,
189 const std::vector<G4Plane3D>& pPlanes,
190 const G4Segment3D& pAABB,
191 G4Segment3D& pExtent) const;
192
193 private:
194
195 /** Original bounding box. */
196 G4ThreeVector fMin, fMax;
197
198 /** Reference to the original sequence of polygonal bases. */
199 const std::vector<const G4ThreeVectorList*>* fPolygons = nullptr;
200};
201
202#endif
std::pair< G4Point3D, G4Point3D > G4Segment3D
std::vector< G4Point3D > G4Polygon3D
std::vector< G4ThreeVector > G4ThreeVectorList
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Transform3D G4Transform3D
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
~G4BoundingEnvelope()=default
G4BoundingEnvelope(const G4ThreeVector &pMin, const G4ThreeVector &pMax)
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
EAxis
Definition geomdefs.hh:54