File indexing completed on 2025-09-18 09:15:07
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
0036
0037
0038 #ifndef G4MULTIUNION_HH
0039 #define G4MULTIUNION_HH
0040
0041 #include <vector>
0042
0043 #include "G4VSolid.hh"
0044 #include "G4ThreeVector.hh"
0045 #include "G4Transform3D.hh"
0046 #include "G4Point3D.hh"
0047 #include "G4Vector3D.hh"
0048 #include "G4SurfBits.hh"
0049 #include "G4Voxelizer.hh"
0050
0051 class G4Polyhedron;
0052
0053 class G4MultiUnion : public G4VSolid
0054 {
0055 friend class G4Voxelizer;
0056
0057 public:
0058
0059 G4MultiUnion() : G4VSolid("") {}
0060 G4MultiUnion(const G4String& name);
0061 ~G4MultiUnion() override;
0062
0063
0064 void AddNode(G4VSolid& solid, const G4Transform3D& trans);
0065 void AddNode(G4VSolid* solid, const G4Transform3D& trans);
0066
0067 G4MultiUnion(const G4MultiUnion& rhs);
0068 G4MultiUnion& operator=(const G4MultiUnion& rhs);
0069
0070
0071 inline const G4Transform3D& GetTransformation(G4int index) const;
0072 inline G4VSolid* GetSolid(G4int index) const;
0073 inline G4int GetNumberOfSolids()const;
0074
0075
0076 EInside Inside(const G4ThreeVector& aPoint) const override;
0077
0078 EInside InsideIterator(const G4ThreeVector& aPoint) const;
0079
0080
0081 G4double DistanceToIn(const G4ThreeVector& aPoint) const override;
0082 G4double DistanceToOut(const G4ThreeVector& aPoint) const override;
0083 inline void SetAccurateSafety(G4bool flag);
0084
0085
0086 G4double DistanceToIn(const G4ThreeVector& aPoint,
0087 const G4ThreeVector& aDirection) const override;
0088 G4double DistanceToOut(const G4ThreeVector& aPoint,
0089 const G4ThreeVector& aDirection,
0090 const G4bool calcNorm = false,
0091 G4bool* validNorm = nullptr,
0092 G4ThreeVector* aNormalVector = nullptr) const override;
0093
0094 G4double DistanceToInNoVoxels(const G4ThreeVector& aPoint,
0095 const G4ThreeVector& aDirection) const;
0096 G4double DistanceToOutVoxels(const G4ThreeVector& aPoint,
0097 const G4ThreeVector& aDirection,
0098 G4ThreeVector* aNormalVector) const;
0099 G4double DistanceToOutVoxelsCore(const G4ThreeVector& aPoint,
0100 const G4ThreeVector& aDirection,
0101 G4ThreeVector* aNormalVector,
0102 G4bool& aConvex,
0103 std::vector<G4int>& candidates) const;
0104 G4double DistanceToOutNoVoxels(const G4ThreeVector& aPoint,
0105 const G4ThreeVector& aDirection,
0106 G4ThreeVector* aNormalVector) const;
0107
0108 G4ThreeVector SurfaceNormal(const G4ThreeVector& aPoint) const override;
0109
0110 void Extent(EAxis aAxis, G4double& aMin, G4double& aMax) const;
0111 void BoundingLimits(G4ThreeVector& aMin, G4ThreeVector& aMax) const override;
0112 G4bool CalculateExtent(const EAxis pAxis,
0113 const G4VoxelLimits& pVoxelLimit,
0114 const G4AffineTransform& pTransform,
0115 G4double& pMin, G4double& pMax) const override;
0116 G4double GetCubicVolume() override;
0117 G4double GetSurfaceArea() override;
0118
0119 G4int GetNumOfConstituents() const override;
0120 G4bool IsFaceted() const override;
0121
0122 G4VSolid* Clone() const override ;
0123
0124 G4GeometryType GetEntityType() const override { return "G4MultiUnion"; }
0125
0126 void Voxelize();
0127
0128
0129
0130 EInside InsideNoVoxels(const G4ThreeVector& aPoint) const;
0131 inline G4Voxelizer& GetVoxels() const;
0132
0133 std::ostream& StreamInfo(std::ostream& os) const override;
0134
0135 G4ThreeVector GetPointOnSurface() const override;
0136
0137 void DescribeYourselfTo ( G4VGraphicsScene& scene ) const override ;
0138 G4Polyhedron* CreatePolyhedron () const override ;
0139 G4Polyhedron* GetPolyhedron () const override;
0140
0141 G4MultiUnion(__void__&);
0142
0143
0144
0145
0146 private:
0147
0148 EInside InsideWithExclusion(const G4ThreeVector& aPoint,
0149 G4SurfBits* bits = nullptr) const;
0150 G4int SafetyFromOutsideNumberNode(const G4ThreeVector& aPoint,
0151 G4double& safety) const;
0152 G4double DistanceToInCandidates(const G4ThreeVector& aPoint,
0153 const G4ThreeVector& aDirection,
0154 std::vector<G4int>& candidates,
0155 G4SurfBits& bits) const;
0156
0157
0158 inline G4ThreeVector GetLocalPoint(const G4Transform3D& trans,
0159 const G4ThreeVector& gpoint) const;
0160 inline G4ThreeVector GetLocalVector(const G4Transform3D& trans,
0161 const G4ThreeVector& gvec) const;
0162 inline G4ThreeVector GetGlobalPoint(const G4Transform3D& trans,
0163 const G4ThreeVector& lpoint) const;
0164 inline G4ThreeVector GetGlobalVector(const G4Transform3D& trans,
0165 const G4ThreeVector& lvec) const;
0166 void TransformLimits(G4ThreeVector& min, G4ThreeVector& max,
0167 const G4Transform3D& transformation) const;
0168 private:
0169
0170 struct G4MultiUnionSurface
0171 {
0172 G4ThreeVector point;
0173 G4VSolid* solid;
0174 };
0175
0176 std::vector<G4VSolid*> fSolids;
0177 std::vector<G4Transform3D> fTransformObjs;
0178 G4Voxelizer fVoxels;
0179 G4double fCubicVolume = 0.0;
0180 G4double fSurfaceArea = 0.0;
0181 G4double kRadTolerance;
0182 mutable G4bool fAccurate = false;
0183
0184 mutable G4bool fRebuildPolyhedron = false;
0185 mutable G4Polyhedron* fpPolyhedron = nullptr;
0186 };
0187
0188
0189 inline G4Voxelizer& G4MultiUnion::GetVoxels() const
0190 {
0191 return (G4Voxelizer&)fVoxels;
0192 }
0193
0194
0195 inline const G4Transform3D& G4MultiUnion::GetTransformation(G4int index) const
0196 {
0197 return fTransformObjs[index];
0198 }
0199
0200
0201 inline G4VSolid* G4MultiUnion::GetSolid(G4int index) const
0202 {
0203 return fSolids[index];
0204 }
0205
0206
0207 inline G4int G4MultiUnion::GetNumberOfSolids() const
0208 {
0209 return G4int(fSolids.size());
0210 }
0211
0212
0213 inline void G4MultiUnion::SetAccurateSafety(G4bool flag)
0214 {
0215 fAccurate = flag;
0216 }
0217
0218
0219 inline
0220 G4ThreeVector G4MultiUnion::GetLocalPoint(const G4Transform3D& trans,
0221 const G4ThreeVector& global) const
0222 {
0223
0224
0225
0226
0227 return trans.inverse()*G4Point3D(global);
0228 }
0229
0230
0231 inline
0232 G4ThreeVector G4MultiUnion::GetLocalVector(const G4Transform3D& trans,
0233 const G4ThreeVector& global) const
0234 {
0235
0236
0237
0238
0239 G4Rotate3D rot;
0240 G4Translate3D transl ;
0241 G4Scale3D scale;
0242
0243 trans.getDecomposition(scale,rot,transl);
0244 return rot.inverse()*G4Vector3D(global);
0245 }
0246
0247
0248 inline
0249 G4ThreeVector G4MultiUnion::GetGlobalPoint(const G4Transform3D& trans,
0250 const G4ThreeVector& local) const
0251 {
0252
0253
0254
0255
0256 return trans*G4Point3D(local);
0257 }
0258
0259
0260 inline
0261 G4ThreeVector G4MultiUnion::GetGlobalVector(const G4Transform3D& trans,
0262 const G4ThreeVector& local) const
0263 {
0264
0265
0266
0267
0268 G4Rotate3D rot;
0269 G4Translate3D transl ;
0270 G4Scale3D scale;
0271
0272 trans.getDecomposition(scale,rot,transl);
0273 return rot*G4Vector3D(local);
0274 }
0275
0276 #endif