Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Voxelizer.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// G4Voxelizer
27//
28// Class description:
29//
30// Voxelizer for tessellated surfaces and solids positioning in 3D space,
31// used in G4TessellatedSolid and G4MultiUnion.
32
33// Author: Marek Gayer (CERN), 19.10.2012 - Created
34// --------------------------------------------------------------------
35#ifndef G4VOXELIZER_HH
36#define G4VOXELIZER_HH
37
38#include <vector>
39#include <string>
40#include <map>
41
42#include "G4Transform3D.hh"
43#include "G4RotationMatrix.hh"
44#include "G4SurfBits.hh"
45#include "G4Box.hh"
46#include "G4VFacet.hh"
47#include "G4VSolid.hh"
48
50{
51 G4ThreeVector hlen; // half length of the box
52 G4ThreeVector pos; // position of the box
53};
54
61
62/**
63 * @brief G4Voxelizer is a tool for generating the optimisation structure
64 * of tessellated surfaces and solids positioned in 3D space; it is used in
65 * G4TessellatedSolid and G4MultiUnion.
66 */
67
69{
70 public:
71
72 /**
73 * Constructor and default Destructor.
74 */
76 ~G4Voxelizer() = default;
77
78 /**
79 * Builds the voxelisation structure for solids positioned in space.
80 * @param[in] solids The list of solids.
81 * @param[in] transforms The associated transformation in space.
82 */
83 void Voxelize(std::vector<G4VSolid*>& solids,
84 std::vector<G4Transform3D>& transforms);
85
86 /**
87 * Builds the voxelisation structure for facets forming a shape.
88 * @param[in] facets The list of facets.
89 */
90 void Voxelize(std::vector<G4VFacet*>& facets);
91
92 /**
93 * Displays the dX, dY, dZ, pX, pY and pZ for each node.
94 */
95 void DisplayVoxelLimits() const;
96
97 /**
98 * Prints the positions of the boundaries of the slices on the three axes.
99 */
100 void DisplayBoundaries();
101
102 /**
103 * Prints which solids are present in the slices previously elaborated.
104 */
105 void DisplayListNodes() const;
106
107 /**
108 * Displays the nodes located in a voxel characterised by its three indexes.
109 */
110 void GetCandidatesVoxel(std::vector<G4int>& voxels);
111
112 /**
113 * Methods returning in a vector container the nodes located in a voxel
114 * characterised by its three indexes.
115 * @returns The total candidates number.
116 */
118 std::vector<G4int>& list,
119 G4SurfBits* crossed = nullptr) const;
120 G4int GetCandidatesVoxelArray(const std::vector<G4int>& voxels,
121 const G4SurfBits bitmasks[],
122 std::vector<G4int>& list,
123 G4SurfBits* crossed = nullptr) const;
124 G4int GetCandidatesVoxelArray(const std::vector<G4int>& voxels,
125 std::vector<G4int>& list,
126 G4SurfBits* crossed = nullptr)const;
127
128 /**
129 * Returns the pointer to the array containing the characteristics
130 * of each box.
131 */
132 inline const std::vector<G4VoxelBox>& GetBoxes() const;
133
134 /**
135 * Returns the boundary vector, given an 'index'.
136 */
137 inline const std::vector<G4double>& GetBoundary(G4int index) const;
138
139 /**
140 * Utility method for checking/updating current voxel given in input.
141 */
143 const G4ThreeVector& direction,
144 std::vector<G4int>& curVoxel) const;
145
146 /**
147 * Updates current voxel based on provided 'point'.
148 */
149 inline void GetVoxel(std::vector<G4int>& curVoxel,
150 const G4ThreeVector& point) const;
151
152 /**
153 * Returns memory size of a slice.
154 */
155 inline G4int GetBitsPerSlice () const;
156
157 /**
158 * Returns true if 'point' is contained within boundaries.
159 */
160 G4bool Contains(const G4ThreeVector& point) const;
161
162 /**
163 * Returns the distance to next boundary, given 'point' and 'direction'
164 * and updates current voxel 'curVoxel', using the index corresponding
165 * to the closest voxel boundary on the ray.
166 */
168 const G4ThreeVector& direction,
169 std::vector<G4int>& curVoxel) const;
170
171 /**
172 * Returns the distance to first bounding box, given 'point' and 'direction'.
173 */
175 const G4ThreeVector& direction) const;
176
177 /**
178 * Returns the minimum distance of 'point' to the bounding box.
179 */
180 G4double DistanceToBoundingBox(const G4ThreeVector& point) const;
181
182 /**
183 * Utility for estimating the isotropic safety from a point 'p' outside
184 * the current solid to any of its surfaces. The algorithm may be accurate
185 * or should provide a fast underestimate, based on safety point 'f'.
186 */
188 const G4ThreeVector& f) const;
189
190 /**
191 * Accessors for voxels and points.
192 */
193 inline G4int GetVoxelsIndex(G4int x, G4int y, G4int z) const;
194 inline G4int GetVoxelsIndex(const std::vector<G4int>& voxels) const;
196 std::vector<G4int>& voxels) const;
197 inline G4int GetPointIndex(const G4ThreeVector& p) const;
198
199 /**
200 * Returns the empty bits container.
201 */
202 inline const G4SurfBits& Empty() const;
203
204 /**
205 * Returns true if empty bit in container, given an 'index'.
206 */
207 inline G4bool IsEmpty(G4int index) const;
208
209 /**
210 * Setters/getter for the maximum number of voxels.
211 */
212 void SetMaxVoxels(G4int max);
213 void SetMaxVoxels(const G4ThreeVector& reductionRatio);
214 inline G4int GetMaxVoxels(G4ThreeVector& ratioOfReduction);
215
216 /**
217 * Logger returning the size of total allocated memory.
218 */
220
221 /**
222 * Utility accessors/functions for voxels.
223 */
224 inline long long GetCountOfVoxels() const;
225 inline long long CountVoxels(std::vector<G4double> boundaries[]) const;
226 inline const std::vector<G4int>&
227 GetCandidates(std::vector<G4int>& curVoxel) const;
228 inline G4int GetVoxelBoxesSize() const;
229 inline const G4VoxelBox& GetVoxelBox(G4int i) const;
230 inline const std::vector<G4int>& GetVoxelBoxCandidates(G4int i) const;
231 inline G4int GetTotalCandidates() const;
232 static void SetDefaultVoxelsCount(G4int count);
234
235 private:
236
237 /**
238 * Binary search function for retrieving a value in a vector.
239 */
240 template <typename T>
241 inline G4int BinarySearch(const std::vector<T>& vec, T value) const;
242
243 /**
244 * Utilities.
245 */
246 G4String GetCandidatesAsString(const G4SurfBits& bits) const;
247 void CreateSortedBoundary(std::vector<G4double>& boundaryRaw, G4int axis);
248 void DisplayBoundaries(std::vector<G4double>& fBoundaries);
249 void FindComponentsFastest(unsigned int mask,
250 std::vector<G4int>& list, G4int i) const;
251 inline G4ThreeVector GetGlobalPoint(const G4Transform3D& trans,
252 const G4ThreeVector& lpoint) const;
253 void TransformLimits(G4ThreeVector& min, G4ThreeVector& max,
254 const G4Transform3D& transformation) const;
255
256 /**
257 * Build utilities.
258 */
259 void BuildEmpty ();
260 void BuildBoundaries();
261 void BuildReduceVoxels(std::vector<G4double> fBoundaries[],
262 G4ThreeVector reductionRatio);
263 void BuildReduceVoxels2(std::vector<G4double> fBoundaries[],
264 G4ThreeVector reductionRatio);
265 void BuildVoxelLimits(std::vector<G4VSolid*>& solids,
266 std::vector<G4Transform3D>& transforms);
267 void BuildVoxelLimits(std::vector<G4VFacet*>& facets);
268 void CreateMiniVoxels(std::vector<G4double> fBoundaries[],
269 G4SurfBits bitmasks[]);
270 void BuildBitmasks(std::vector<G4double> fBoundaries[],
271 G4SurfBits bitmasks[], G4bool countsOnly = false);
272 void BuildBoundingBox();
273 void BuildBoundingBox(G4ThreeVector& amin, G4ThreeVector& amax,
274 G4double tolerance = 0.0);
275 void SetReductionRatio(G4int maxVoxels, G4ThreeVector& reductionRatio);
276
277
278 private:
279
280 /**
281 * @brief G4VoxelComparator is utility class used for comparing voxels.
282 */
283
284 class G4VoxelComparator
285 {
286 public:
287
288 std::vector<G4VoxelInfo>& fVoxels;
289
290 G4VoxelComparator(std::vector<G4VoxelInfo>& voxels) : fVoxels(voxels) {}
291
292 G4bool operator()(const G4int& l, const G4int& r) const
293 {
294 G4VoxelInfo &lv = fVoxels[l], &rv = fVoxels[r];
295 G4int left = lv.count + fVoxels[lv.next].count;
296 G4int right = rv.count + fVoxels[rv.next].count;
297 return (left == right) ? l < r : left < right;
298 }
299 };
300
301 static G4int fDefaultVoxelsCount;
302
303 std::vector<G4VoxelBox> fVoxelBoxes;
304 std::vector<std::vector<G4int> > fVoxelBoxesCandidates;
305 mutable std::map<G4int, std::vector<G4int> > fCandidates;
306
307 const std::vector<G4int> fNoCandidates;
308
309 long long fCountOfVoxels;
310
311 G4int fNPerSlice;
312
313 std::vector<G4VoxelBox> fBoxes;
314 // Array of box limits on the 3 cartesian axis
315
316 std::vector<G4double> fBoundaries[3];
317 // Sorted and if need skimmed fBoundaries along X,Y,Z axis
318
319 std::vector<G4int> fCandidatesCounts[3];
320
321 G4int fTotalCandidates;
322
323 G4SurfBits fBitmasks[3];
324
325 G4ThreeVector fBoundingBoxCenter;
326
327 G4Box fBoundingBox;
328
329 G4ThreeVector fBoundingBoxSize;
330
331 G4ThreeVector fReductionRatio;
332
333 G4int fMaxVoxels;
334
335 G4double fTolerance;
336
337 G4SurfBits fEmpty;
338};
339
340#include "G4Voxelizer.icc"
341
342#endif
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Transform3D G4Transform3D
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4SurfBits provides a simple container of bits, to be used for optimization of tessellated surfaces....
Definition G4SurfBits.hh:57
const G4SurfBits & Empty() const
long long CountVoxels(std::vector< G4double > boundaries[]) const
G4double DistanceToBoundingBox(const G4ThreeVector &point) const
G4double MinDistanceToBox(const G4ThreeVector &p, const G4ThreeVector &f) const
G4bool GetPointVoxel(const G4ThreeVector &p, std::vector< G4int > &voxels) const
long long GetCountOfVoxels() const
const std::vector< G4double > & GetBoundary(G4int index) const
const std::vector< G4VoxelBox > & GetBoxes() const
G4bool IsEmpty(G4int index) const
G4bool UpdateCurrentVoxel(const G4ThreeVector &point, const G4ThreeVector &direction, std::vector< G4int > &curVoxel) const
void DisplayListNodes() const
void GetVoxel(std::vector< G4int > &curVoxel, const G4ThreeVector &point) const
void DisplayVoxelLimits() const
G4int GetMaxVoxels(G4ThreeVector &ratioOfReduction)
static void SetDefaultVoxelsCount(G4int count)
G4int GetBitsPerSlice() const
G4double DistanceToFirst(const G4ThreeVector &point, const G4ThreeVector &direction) const
G4int GetCandidatesVoxelArray(const G4ThreeVector &point, std::vector< G4int > &list, G4SurfBits *crossed=nullptr) const
G4int GetVoxelBoxesSize() const
void SetMaxVoxels(G4int max)
const G4VoxelBox & GetVoxelBox(G4int i) const
void GetCandidatesVoxel(std::vector< G4int > &voxels)
G4int AllocatedMemory()
G4int GetPointIndex(const G4ThreeVector &p) const
G4int GetVoxelsIndex(const std::vector< G4int > &voxels) const
~G4Voxelizer()=default
static G4int GetDefaultVoxelsCount()
void DisplayBoundaries()
G4double DistanceToNext(const G4ThreeVector &point, const G4ThreeVector &direction, std::vector< G4int > &curVoxel) const
void Voxelize(std::vector< G4VSolid * > &solids, std::vector< G4Transform3D > &transforms)
const std::vector< G4int > & GetVoxelBoxCandidates(G4int i) const
G4int GetVoxelsIndex(G4int x, G4int y, G4int z) const
G4bool Contains(const G4ThreeVector &point) const
const std::vector< G4int > & GetCandidates(std::vector< G4int > &curVoxel) const
G4int GetTotalCandidates() const
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668
G4ThreeVector hlen
G4ThreeVector pos