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