Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4Voxelizer.hh
이 파일의 문서화 페이지로 가기
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:$
27 //
28 // --------------------------------------------------------------------
29 // GEANT 4 class header file
30 //
31 // G4Voxelizer
32 //
33 // Class description:
34 //
35 // Voxelizer for tessellated surfaces and solids positioning in 3D space,
36 // used in G4TessellatedSolid and G4MultiUnion.
37 
38 // History:
39 // 19.10.12 Marek Gayer, created
40 // --------------------------------------------------------------------
41 
42 #ifndef G4Voxelizer_HH
43 #define G4Voxelizer_HH
44 
45 #include <vector>
46 #include <string>
47 #include <map>
48 
49 #include "G4Transform3D.hh"
50 #include "G4RotationMatrix.hh"
51 #include "G4SurfBits.hh"
52 #include "G4Box.hh"
53 #include "G4VFacet.hh"
54 #include "G4VSolid.hh"
55 
56 struct G4VoxelBox
57 {
58  G4ThreeVector hlen; // half length of the box
59  G4ThreeVector pos; // position of the box
60 };
61 
63 {
67 };
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<G4VSolid *> &solids,
77  std::vector<G4Transform3D>& transforms);
78  void Voxelize(std::vector<G4VFacet *> &facets);
79 
80  void DisplayVoxelLimits() const;
81  void DisplayBoundaries();
82  void DisplayListNodes() const;
83 
84  G4Voxelizer();
85  ~G4Voxelizer();
86 
87  void GetCandidatesVoxel(std::vector<G4int> &voxels);
88  // Method displaying the nodes located in a voxel characterized
89  // by its three indexes.
90 
92  std::vector<G4int> &list,
93  G4SurfBits *crossed=0) const;
94  // Method returning in a vector container the nodes located in a voxel
95  // characterized by its three indexes.
96  G4int GetCandidatesVoxelArray(const std::vector<G4int> &voxels,
97  const G4SurfBits bitmasks[],
98  std::vector<G4int> &list,
99  G4SurfBits *crossed=0) const;
100  G4int GetCandidatesVoxelArray(const std::vector<G4int> &voxels,
101  std::vector<G4int> &list,
102  G4SurfBits *crossed=0)const;
103 
104  inline const std::vector<G4VoxelBox> &GetBoxes() const;
105  // Method returning the pointer to the array containing the
106  // characteristics of each box.
107 
108  inline const std::vector<G4double> &GetBoundary(G4int index) const;
109 
111  const G4ThreeVector &direction,
112  std::vector<G4int> &curVoxel) const;
113 
114  inline void GetVoxel(std::vector<G4int> &curVoxel,
115  const G4ThreeVector &point) const;
116 
117  inline G4int GetBitsPerSlice () const;
118 
119  G4bool Contains(const G4ThreeVector &point) const;
120 
122  const G4ThreeVector &direction,
123  std::vector<G4int> &curVoxel) const;
124 
126  const G4ThreeVector &direction) const;
127 
128  G4double DistanceToBoundingBox(const G4ThreeVector &point) const;
129 
130  inline G4int GetVoxelsIndex(G4int x, G4int y, G4int z) const;
131  inline G4int GetVoxelsIndex(const std::vector<G4int> &voxels) const;
132  inline G4bool GetPointVoxel(const G4ThreeVector& p,
133  std::vector<G4int>& voxels) const;
134  inline G4int GetPointIndex(const G4ThreeVector &p) const;
135 
136  inline const G4SurfBits &Empty() const;
137  inline G4bool IsEmpty(G4int index) const;
138 
139  void SetMaxVoxels(G4int max);
140  void SetMaxVoxels(const G4ThreeVector &reductionRatio);
141 
142  inline G4int GetMaxVoxels(G4ThreeVector &ratioOfReduction);
143 
145 
146  inline long long GetCountOfVoxels() const;
147 
148  inline long long CountVoxels(std::vector<G4double> boundaries[]) const;
149 
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  inline G4int GetTotalCandidates() const;
160 
161  static G4double MinDistanceToBox (const G4ThreeVector &aPoint,
162  const G4ThreeVector &f);
163 
164  static void SetDefaultVoxelsCount(G4int count);
165 
166  static G4int GetDefaultVoxelsCount();
167 
168  private:
169 
171  {
172  public:
173 
174  std::vector<G4VoxelInfo> &fVoxels;
175 
176  G4VoxelComparator(std::vector<G4VoxelInfo> &voxels) : fVoxels(voxels) {}
177 
178  G4bool operator()(const G4int& l, const G4int& r) const
179  {
180  G4VoxelInfo &lv = fVoxels[l], &rv = fVoxels[r];
181  G4int left = lv.count + fVoxels[lv.next].count;
182  G4int right = rv.count + fVoxels[rv.next].count;
183  return (left == right) ? l < r : left < right;
184  }
185  };
186 
187  private:
188 
189  void BuildEmpty ();
190 
191  G4String GetCandidatesAsString(const G4SurfBits &bits) const;
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<G4VSolid*>& solids,
203  std::vector<G4Transform3D>& transforms);
204  void BuildVoxelLimits(std::vector<G4VFacet *> &facets);
205 
206  void DisplayBoundaries(std::vector<G4double> &fBoundaries);
207 
208  void BuildBitmasks(std::vector<G4double> fBoundaries[],
209  G4SurfBits bitmasks[], G4bool countsOnly=false);
210 
211  void BuildBoundingBox();
212  void BuildBoundingBox(G4ThreeVector& amin, G4ThreeVector& amax,
213  G4double tolerance = 0);
214 
215  void SetReductionRatio(G4int maxVoxels, G4ThreeVector &reductionRatio);
216 
217  void CreateMiniVoxels(std::vector<G4double> fBoundaries[],
218  G4SurfBits bitmasks[]);
219  static void FindComponentsFastest(unsigned int mask,
220  std::vector<G4int> &list, G4int i);
221 
222  inline G4ThreeVector GetGlobalPoint(const G4Transform3D& trans,
223  const G4ThreeVector& lpoint) const;
225  const G4Transform3D& transformation) const;
226 
227  private:
228 
230 
231  std::vector<G4VoxelBox> fVoxelBoxes;
232  std::vector<std::vector<G4int> > fVoxelBoxesCandidates;
233  mutable std::map<G4int, std::vector<G4int> > fCandidates;
234 
235  const std::vector<G4int> fNoCandidates;
236 
237  long long fCountOfVoxels;
238 
240 
241  std::vector<G4VoxelBox> fBoxes;
242  // Array of box limits on the 3 cartesian axis
243 
244  std::vector<G4double> fBoundaries[3];
245  // Sorted and if need skimmed fBoundaries along X,Y,Z axis
246 
247  std::vector<G4int> fCandidatesCounts[3];
248 
250 
252 
254 
256 
258 
260 
262 
264 
266 };
267 
268 #include "G4Voxelizer.icc"
269 
270 #endif
Float_t x
Definition: compare.C:6
G4int GetPointIndex(const G4ThreeVector &p) const
std::vector< G4VoxelInfo > & fVoxels
Definition: G4Voxelizer.hh:174
T max(const T t1, const T t2)
brief Return the largest of the two arguments
long long fCountOfVoxels
Definition: G4Voxelizer.hh:237
const G4SurfBits & Empty() const
const std::vector< G4int > & GetVoxelBoxCandidates(G4int i) const
G4ThreeVector fBoundingBoxSize
Definition: G4Voxelizer.hh:257
std::vector< G4double > fBoundaries[3]
Definition: G4Voxelizer.hh:244
void Voxelize(std::vector< G4VSolid * > &solids, std::vector< G4Transform3D > &transforms)
Definition: G4Voxelizer.cc:724
G4VoxelComparator(std::vector< G4VoxelInfo > &voxels)
Definition: G4Voxelizer.hh:176
G4bool Contains(const G4ThreeVector &point) const
G4SurfBits fEmpty
Definition: G4Voxelizer.hh:265
void BuildReduceVoxels(std::vector< G4double > fBoundaries[], G4ThreeVector reductionRatio)
Definition: G4Voxelizer.cc:526
G4bool IsEmpty(G4int index) const
G4ThreeVector GetGlobalPoint(const G4Transform3D &trans, const G4ThreeVector &lpoint) const
G4ThreeVector fReductionRatio
Definition: G4Voxelizer.hh:259
G4int GetCandidatesVoxelArray(const G4ThreeVector &point, std::vector< G4int > &list, G4SurfBits *crossed=0) const
Definition: G4Voxelizer.cc:984
static ush mask[]
Definition: csz_inflate.cc:317
Float_t y
Definition: compare.C:6
const char * p
Definition: xmltok.h:285
std::map< G4int, std::vector< G4int > > fCandidates
Definition: G4Voxelizer.hh:233
void BuildVoxelLimits(std::vector< G4VSolid * > &solids, std::vector< G4Transform3D > &transforms)
Definition: G4Voxelizer.cc:122
void DisplayVoxelLimits() const
Definition: G4Voxelizer.cc:208
Double_t z
const std::vector< G4int > & GetCandidates(std::vector< G4int > &curVoxel) const
const std::vector< G4double > & GetBoundary(G4int index) const
G4double fTolerance
Definition: G4Voxelizer.hh:263
void DisplayBoundaries()
Definition: G4Voxelizer.cc:335
long long GetCountOfVoxels() const
#define G4ThreadLocal
Definition: tls.hh:69
G4int previous
Definition: G4Voxelizer.hh:65
G4ThreeVector hlen
Definition: G4Voxelizer.hh:58
void BuildBitmasks(std::vector< G4double > fBoundaries[], G4SurfBits bitmasks[], G4bool countsOnly=false)
Definition: G4Voxelizer.cc:362
G4int GetTotalCandidates() const
std::vector< G4int > fCandidatesCounts[3]
Definition: G4Voxelizer.hh:247
void SetMaxVoxels(G4int max)
void DisplayListNodes() const
Definition: G4Voxelizer.cc:442
void SetReductionRatio(G4int maxVoxels, G4ThreeVector &reductionRatio)
Definition: G4Voxelizer.cc:510
G4double DistanceToFirst(const G4ThreeVector &point, const G4ThreeVector &direction) const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< std::vector< G4int > > fVoxelBoxesCandidates
Definition: G4Voxelizer.hh:232
const std::vector< G4int > fNoCandidates
Definition: G4Voxelizer.hh:235
G4bool GetPointVoxel(const G4ThreeVector &p, std::vector< G4int > &voxels) const
G4int GetMaxVoxels(G4ThreeVector &ratioOfReduction)
const XML_Char int const XML_Char * value
Definition: expat.h:331
G4String GetCandidatesAsString(const G4SurfBits &bits) const
Definition: G4Voxelizer.cc:427
std::vector< G4VoxelBox > fVoxelBoxes
Definition: G4Voxelizer.hh:231
G4int fNPerSlice
Definition: G4Voxelizer.hh:239
G4SurfBits fBitmasks[3]
Definition: G4Voxelizer.hh:251
Definition: G4Box.hh:64
void GetCandidatesVoxel(std::vector< G4int > &voxels)
Definition: G4Voxelizer.cc:907
G4int GetVoxelBoxesSize() const
static G4double MinDistanceToBox(const G4ThreeVector &aPoint, const G4ThreeVector &f)
static void FindComponentsFastest(unsigned int mask, std::vector< G4int > &list, G4int i)
Definition: G4Voxelizer.cc:922
const std::vector< G4VoxelBox > & GetBoxes() const
G4ThreeVector pos
Definition: G4Voxelizer.hh:59
void CreateSortedBoundary(std::vector< G4double > &boundaryRaw, G4int axis)
Definition: G4Voxelizer.cc:225
G4int fTotalCandidates
Definition: G4Voxelizer.hh:249
G4int AllocatedMemory()
long long CountVoxels(std::vector< G4double > boundaries[]) const
G4double DistanceToNext(const G4ThreeVector &point, const G4ThreeVector &direction, std::vector< G4int > &curVoxel) const
void BuildEmpty()
Definition: G4Voxelizer.cc:79
G4Box fBoundingBox
Definition: G4Voxelizer.hh:255
int G4int
Definition: G4Types.hh:78
std::vector< G4VoxelBox > fBoxes
Definition: G4Voxelizer.hh:241
G4int fMaxVoxels
Definition: G4Voxelizer.hh:261
G4int GetBitsPerSlice() const
G4ThreeVector fBoundingBoxCenter
Definition: G4Voxelizer.hh:253
G4bool UpdateCurrentVoxel(const G4ThreeVector &point, const G4ThreeVector &direction, std::vector< G4int > &curVoxel) const
static G4ThreadLocal G4int fDefaultVoxelsCount
Definition: G4Voxelizer.hh:229
void BuildBoundingBox()
Definition: G4Voxelizer.cc:467
G4double DistanceToBoundingBox(const G4ThreeVector &point) const
void BuildReduceVoxels2(std::vector< G4double > fBoundaries[], G4ThreeVector reductionRatio)
Definition: G4Voxelizer.cc:681
G4bool operator()(const G4int &l, const G4int &r) const
Definition: G4Voxelizer.hh:178
void TransformLimits(G4ThreeVector &min, G4ThreeVector &max, const G4Transform3D &transformation) const
Definition: G4Voxelizer.cc:941
void GetVoxel(std::vector< G4int > &curVoxel, const G4ThreeVector &point) const
G4int GetVoxelsIndex(G4int x, G4int y, G4int z) const
void CreateMiniVoxels(std::vector< G4double > fBoundaries[], G4SurfBits bitmasks[])
Definition: G4Voxelizer.cc:742
static G4int BinarySearch(const std::vector< T > &vec, T value)
void BuildBoundaries()
Definition: G4Voxelizer.cc:255
const G4VoxelBox & GetVoxelBox(G4int i) const
static G4int GetDefaultVoxelsCount()
static void SetDefaultVoxelsCount(G4int count)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments