Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4SurfaceVoxelizer.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 // $Id: G4SurfaceVoxelizer.hh 67011 2013-01-29 16:17:41Z gcosmo $
27 //
28 // --------------------------------------------------------------------
29 // GEANT 4 class header file
30 //
31 // G4SurfaceVoxelizer
32 //
33 // Class description:
34 //
35 // Voxelizer for tessellated surfaces, used in G4TessellatedSolid.
36 
37 // History:
38 // 19.10.12 Marek Gayer, created
39 // --------------------------------------------------------------------
40 
41 #ifndef G4SurfaceVoxelizer_HH
42 #define G4SurfaceVoxelizer_HH
43 
44 #include <vector>
45 #include <string>
46 #include <map>
47 
48 #include "G4SurfBits.hh"
49 #include "G4Box.hh"
50 #include "G4VFacet.hh"
51 
52 #include "Randomize.hh"
53 
54 struct G4VoxelBox
55 {
56  G4ThreeVector hlen; // half length of the box
57  G4ThreeVector pos; // position of the box
58 };
59 
61 {
65 };
66 
68 {
70 
71  public:
72 
73  template <typename T>
74  static inline G4int BinarySearch(const std::vector<T> &vec, T value);
75 
76  void Voxelize(std::vector<G4VFacet *> &facets);
77 
78  void DisplayVoxelLimits();
79  void DisplayBoundaries();
80  void DisplayListNodes();
81 
84 
85  void GetCandidatesVoxel(std::vector<G4int> &voxels);
86  // Method displaying the nodes located in a voxel characterized
87  // by its three indexes.
88 
90  std::vector<G4int> &list,
91  G4SurfBits *crossed=0) const;
92  // Method returning in a vector container the nodes located in a voxel
93  // characterized by its three indexes.
94  G4int GetCandidatesVoxelArray(const std::vector<G4int> &voxels,
95  const G4SurfBits bitmasks[],
96  std::vector<G4int> &list,
97  G4SurfBits *crossed=0) const;
98  G4int GetCandidatesVoxelArray(const std::vector<G4int> &voxels,
99  std::vector<G4int> &list,
100  G4SurfBits *crossed=0)const;
101 
102  inline const std::vector<G4VoxelBox> &GetBoxes() const;
103  // Method returning the pointer to the array containing the
104  // characteristics of each box.
105 
106  inline const std::vector<G4double> &GetBoundary(G4int index) const;
107 
109  const G4ThreeVector &direction,
110  std::vector<G4int> &curVoxel) const;
111 
112  inline void GetVoxel(std::vector<G4int> &curVoxel,
113  const G4ThreeVector &point) const;
114 
115  inline G4int GetBitsPerSlice () const;
116 
117  G4bool Contains(const G4ThreeVector &point) const;
118 
120  const G4ThreeVector &direction,
121  const std::vector<G4int> &curVoxel) const;
122 
124  const G4ThreeVector &direction) const;
125 
126  G4double DistanceToBoundingBox(const G4ThreeVector &point) const;
127 
128  inline G4int GetVoxelsIndex(G4int x, G4int y, G4int z) const;
129  inline G4int GetVoxelsIndex(const std::vector<G4int> &voxels) const;
130 
131  inline G4int GetPointIndex(const G4ThreeVector &p) const;
132 
133  inline const G4SurfBits &Empty() const;
134  inline G4bool IsEmpty(G4int index) const;
135 
136  void SetMaxVoxels(G4int max);
137  void SetMaxVoxels(const G4ThreeVector &reductionRatio);
138 
139  inline G4int GetMaxVoxels(G4ThreeVector &ratioOfReduction);
140 
142 
143  inline long long GetCountOfVoxels() const;
144 
145  inline long long CountVoxels(std::vector<G4double> boundaries[]) const;
146 
147  inline G4int GetCandidates(std::vector<G4int> &curVoxel,
148  std::vector<G4int> *&candidates,
149  std::vector<G4int> &space) const;
150  inline const std::vector<G4int> &
151  GetCandidates(std::vector<G4int> &curVoxel) const;
152 
153  inline G4int GetVoxelBoxesSize() const;
154 
155  inline const G4VoxelBox &GetVoxelBox(G4int i) const;
156 
157  inline const std::vector<G4int> &GetVoxelBoxCandidates(G4int i) const;
158 
159  static G4double MinDistanceToBox (const G4ThreeVector &aPoint,
160  const G4ThreeVector &f);
161 
162  static G4int SetDefaultVoxelsCount(G4int count);
163 
164  static G4int GetDefaultVoxelsCount();
165 
166  private:
167 
168  class G4VoxelComparator
169  {
170  public:
171 
172  std::vector<G4VoxelInfo> &fVoxels;
173 
174  G4VoxelComparator(std::vector<G4VoxelInfo> &voxels) : fVoxels(voxels) {}
175 
176  G4bool operator()(const G4int& l, const G4int& r) const
177  {
178  G4VoxelInfo &lv = fVoxels[l], &rv = fVoxels[r];
179  G4int left = lv.count + fVoxels[lv.next].count;
180  G4int right = rv.count + fVoxels[rv.next].count;
181  return (left == right) ? l < r : left < right;
182  }
183  };
184 
185  private:
186 
187  static int fDefaultVoxelsCount;
188 
189  void BuildEmpty ();
190 
191  G4String GetCandidatesAsString(const G4SurfBits &bits);
192 
193  void CreateSortedBoundary(std::vector<G4double> &boundaryRaw, G4int axis);
194 
195  void BuildBoundaries();
196 
197  void BuildReduceVoxels(std::vector<G4double> fBoundaries[],
198  G4ThreeVector reductionRatio);
199  void BuildReduceVoxels2(std::vector<G4double> fBoundaries[],
200  G4ThreeVector reductionRatio);
201 
202  void BuildVoxelLimits(std::vector<G4VFacet *> &facets);
203 
204  void DisplayBoundaries(std::vector<G4double> &fBoundaries);
205 
206  void BuildBitmasks(std::vector<G4double> fBoundaries[],
207  G4SurfBits bitmasks[], G4bool countsOnly=false);
208 
209  void BuildBoundingBox();
210 
211  void SetReductionRatio(G4int maxVoxels, G4ThreeVector &reductionRatio);
212 
213  void CreateMiniVoxels(std::vector<G4double> fBoundaries[],
214  G4SurfBits bitmasks[]);
215  static void FindComponentsFastest(unsigned int mask,
216  std::vector<G4int> &list, G4int i);
217 
218  private:
219 
220  std::vector<G4VoxelBox> fVoxelBoxes;
221  std::vector<std::vector<G4int> > fVoxelBoxesCandidates;
222  mutable std::map<G4int, std::vector<G4int> > fCandidates;
223 
224  const std::vector<G4int> fNoCandidates;
225 
226  long long fCountOfVoxels;
227 
228  G4int fNPerSlice;
229 
230  std::vector<G4VoxelBox> fBoxes;
231  // Array of box limits on the 3 cartesian axis
232 
233  std::vector<G4double> fBoundaries[3];
234  // Sorted and if need skimmed fBoundaries along X,Y,Z axis
235 
236  std::vector<G4int> fCandidatesCounts[3];
237 
238  G4int fTotalCandidates;
239 
240  G4SurfBits fBitmasks[3];
241 
242  G4ThreeVector fBoundingBoxCenter;
243 
244  G4Box fBoundingBox;
245 
246  G4ThreeVector fBoundingBoxSize;
247 
248  G4ThreeVector fReductionRatio;
249 
250  G4int fMaxVoxels;
251 
252  G4double fTolerance;
253 
254  G4SurfBits fEmpty;
255 };
256 
257 #include "G4SurfaceVoxelizer.icc"
258 
259 #endif