File indexing completed on 2025-01-18 09:59:25
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035 #ifndef G4VOXELIZER_HH
0036 #define G4VOXELIZER_HH
0037
0038 #include <vector>
0039 #include <string>
0040 #include <map>
0041
0042 #include "G4Transform3D.hh"
0043 #include "G4RotationMatrix.hh"
0044 #include "G4SurfBits.hh"
0045 #include "G4Box.hh"
0046 #include "G4VFacet.hh"
0047 #include "G4VSolid.hh"
0048
0049 struct G4VoxelBox
0050 {
0051 G4ThreeVector hlen;
0052 G4ThreeVector pos;
0053 };
0054
0055 struct G4VoxelInfo
0056 {
0057 G4int count;
0058 G4int previous;
0059 G4int next;
0060 };
0061
0062 class G4Voxelizer
0063 {
0064 public:
0065
0066 template <typename T>
0067 static inline G4int BinarySearch(const std::vector<T>& vec, T value);
0068
0069 void Voxelize(std::vector<G4VSolid*>& solids,
0070 std::vector<G4Transform3D>& transforms);
0071 void Voxelize(std::vector<G4VFacet*>& facets);
0072
0073 void DisplayVoxelLimits() const;
0074 void DisplayBoundaries();
0075 void DisplayListNodes() const;
0076
0077 G4Voxelizer();
0078 ~G4Voxelizer();
0079
0080 void GetCandidatesVoxel(std::vector<G4int>& voxels);
0081
0082
0083
0084 G4int GetCandidatesVoxelArray(const G4ThreeVector& point,
0085 std::vector<G4int>& list,
0086 G4SurfBits* crossed = nullptr) const;
0087
0088
0089 G4int GetCandidatesVoxelArray(const std::vector<G4int>& voxels,
0090 const G4SurfBits bitmasks[],
0091 std::vector<G4int>& list,
0092 G4SurfBits* crossed = nullptr) const;
0093 G4int GetCandidatesVoxelArray(const std::vector<G4int>& voxels,
0094 std::vector<G4int>& list,
0095 G4SurfBits* crossed = nullptr)const;
0096
0097 inline const std::vector<G4VoxelBox>& GetBoxes() const;
0098
0099
0100
0101 inline const std::vector<G4double>& GetBoundary(G4int index) const;
0102
0103 G4bool UpdateCurrentVoxel(const G4ThreeVector& point,
0104 const G4ThreeVector& direction,
0105 std::vector<G4int>& curVoxel) const;
0106
0107 inline void GetVoxel(std::vector<G4int>& curVoxel,
0108 const G4ThreeVector& point) const;
0109
0110 inline G4int GetBitsPerSlice () const;
0111
0112 G4bool Contains(const G4ThreeVector& point) const;
0113
0114 G4double DistanceToNext(const G4ThreeVector& point,
0115 const G4ThreeVector& direction,
0116 std::vector<G4int>& curVoxel) const;
0117
0118 G4double DistanceToFirst(const G4ThreeVector& point,
0119 const G4ThreeVector& direction) const;
0120
0121 G4double DistanceToBoundingBox(const G4ThreeVector& point) const;
0122
0123 inline G4int GetVoxelsIndex(G4int x, G4int y, G4int z) const;
0124 inline G4int GetVoxelsIndex(const std::vector<G4int>& voxels) const;
0125 inline G4bool GetPointVoxel(const G4ThreeVector& p,
0126 std::vector<G4int>& voxels) const;
0127 inline G4int GetPointIndex(const G4ThreeVector& p) const;
0128
0129 inline const G4SurfBits& Empty() const;
0130 inline G4bool IsEmpty(G4int index) const;
0131
0132 void SetMaxVoxels(G4int max);
0133 void SetMaxVoxels(const G4ThreeVector& reductionRatio);
0134
0135 inline G4int GetMaxVoxels(G4ThreeVector& ratioOfReduction);
0136
0137 G4int AllocatedMemory();
0138
0139 inline long long GetCountOfVoxels() const;
0140
0141 inline long long CountVoxels(std::vector<G4double> boundaries[]) const;
0142
0143 inline const std::vector<G4int>&
0144 GetCandidates(std::vector<G4int>& curVoxel) const;
0145
0146 inline G4int GetVoxelBoxesSize() const;
0147
0148 inline const G4VoxelBox &GetVoxelBox(G4int i) const;
0149
0150 inline const std::vector<G4int>& GetVoxelBoxCandidates(G4int i) const;
0151
0152 inline G4int GetTotalCandidates() const;
0153
0154 static G4double MinDistanceToBox (const G4ThreeVector& aPoint,
0155 const G4ThreeVector& f);
0156
0157 static void SetDefaultVoxelsCount(G4int count);
0158
0159 static G4int GetDefaultVoxelsCount();
0160
0161 private:
0162
0163 class G4VoxelComparator
0164 {
0165 public:
0166
0167 std::vector<G4VoxelInfo>& fVoxels;
0168
0169 G4VoxelComparator(std::vector<G4VoxelInfo>& voxels) : fVoxels(voxels) {}
0170
0171 G4bool operator()(const G4int& l, const G4int& r) const
0172 {
0173 G4VoxelInfo &lv = fVoxels[l], &rv = fVoxels[r];
0174 G4int left = lv.count + fVoxels[lv.next].count;
0175 G4int right = rv.count + fVoxels[rv.next].count;
0176 return (left == right) ? l < r : left < right;
0177 }
0178 };
0179
0180 void BuildEmpty ();
0181
0182 G4String GetCandidatesAsString(const G4SurfBits& bits) const;
0183
0184 void CreateSortedBoundary(std::vector<G4double>& boundaryRaw, G4int axis);
0185
0186 void BuildBoundaries();
0187
0188 void BuildReduceVoxels(std::vector<G4double> fBoundaries[],
0189 G4ThreeVector reductionRatio);
0190 void BuildReduceVoxels2(std::vector<G4double> fBoundaries[],
0191 G4ThreeVector reductionRatio);
0192
0193 void BuildVoxelLimits(std::vector<G4VSolid*>& solids,
0194 std::vector<G4Transform3D>& transforms);
0195 void BuildVoxelLimits(std::vector<G4VFacet*>& facets);
0196
0197 void DisplayBoundaries(std::vector<G4double>& fBoundaries);
0198
0199 void BuildBitmasks(std::vector<G4double> fBoundaries[],
0200 G4SurfBits bitmasks[], G4bool countsOnly = false);
0201
0202 void BuildBoundingBox();
0203 void BuildBoundingBox(G4ThreeVector& amin, G4ThreeVector& amax,
0204 G4double tolerance = 0.0);
0205
0206 void SetReductionRatio(G4int maxVoxels, G4ThreeVector& reductionRatio);
0207
0208 void CreateMiniVoxels(std::vector<G4double> fBoundaries[],
0209 G4SurfBits bitmasks[]);
0210 static void FindComponentsFastest(unsigned int mask,
0211 std::vector<G4int>& list, G4int i);
0212
0213 inline G4ThreeVector GetGlobalPoint(const G4Transform3D& trans,
0214 const G4ThreeVector& lpoint) const;
0215 void TransformLimits(G4ThreeVector& min, G4ThreeVector& max,
0216 const G4Transform3D& transformation) const;
0217
0218 private:
0219
0220 static G4ThreadLocal G4int fDefaultVoxelsCount;
0221
0222 std::vector<G4VoxelBox> fVoxelBoxes;
0223 std::vector<std::vector<G4int> > fVoxelBoxesCandidates;
0224 mutable std::map<G4int, std::vector<G4int> > fCandidates;
0225
0226 const std::vector<G4int> fNoCandidates;
0227
0228 long long fCountOfVoxels;
0229
0230 G4int fNPerSlice;
0231
0232 std::vector<G4VoxelBox> fBoxes;
0233
0234
0235 std::vector<G4double> fBoundaries[3];
0236
0237
0238 std::vector<G4int> fCandidatesCounts[3];
0239
0240 G4int fTotalCandidates;
0241
0242 G4SurfBits fBitmasks[3];
0243
0244 G4ThreeVector fBoundingBoxCenter;
0245
0246 G4Box fBoundingBox;
0247
0248 G4ThreeVector fBoundingBoxSize;
0249
0250 G4ThreeVector fReductionRatio;
0251
0252 G4int fMaxVoxels;
0253
0254 G4double fTolerance;
0255
0256 G4SurfBits fEmpty;
0257 };
0258
0259 #include "G4Voxelizer.icc"
0260
0261 #endif