Geant4 Cross Reference |
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 // G4MultiUnion 27 // 28 // Class description: 29 // 30 // An instance of "G4MultiUnion" constitutes a grouping of several solids. 31 // The constituent solids are stored with their respective location in an 32 // instance of "G4Node". An instance of "G4MultiUnion" is subsequently 33 // composed of one or several nodes. 34 35 // 19.10.12 M.Gayer - Original implementation from USolids module 36 // 06.04.17 G.Cosmo - Adapted implementation in Geant4 for VecGeom migration 37 // -------------------------------------------------------------------- 38 #ifndef G4MULTIUNION_HH 39 #define G4MULTIUNION_HH 40 41 #include <vector> 42 43 #include "G4VSolid.hh" 44 #include "G4ThreeVector.hh" 45 #include "G4Transform3D.hh" 46 #include "G4Point3D.hh" 47 #include "G4Vector3D.hh" 48 #include "G4SurfBits.hh" 49 #include "G4Voxelizer.hh" 50 51 class G4Polyhedron; 52 53 class G4MultiUnion : public G4VSolid 54 { 55 friend class G4Voxelizer; 56 57 public: 58 59 G4MultiUnion() : G4VSolid("") {} 60 G4MultiUnion(const G4String& name); 61 ~G4MultiUnion() override; 62 63 // Build the multiple union by adding nodes 64 void AddNode(G4VSolid& solid, const G4Transform3D& trans); 65 void AddNode(G4VSolid* solid, const G4Transform3D& trans); 66 67 G4MultiUnion(const G4MultiUnion& rhs); 68 G4MultiUnion& operator=(const G4MultiUnion& rhs); 69 70 // Accessors 71 inline const G4Transform3D& GetTransformation(G4int index) const; 72 inline G4VSolid* GetSolid(G4int index) const; 73 inline G4int GetNumberOfSolids()const; 74 75 // Navigation methods 76 EInside Inside(const G4ThreeVector& aPoint) const override; 77 78 EInside InsideIterator(const G4ThreeVector& aPoint) const; 79 80 // Safety methods 81 G4double DistanceToIn(const G4ThreeVector& aPoint) const override; 82 G4double DistanceToOut(const G4ThreeVector& aPoint) const override; 83 inline void SetAccurateSafety(G4bool flag); 84 85 // Exact distance methods 86 G4double DistanceToIn(const G4ThreeVector& aPoint, 87 const G4ThreeVector& aDirection) const override; 88 G4double DistanceToOut(const G4ThreeVector& aPoint, 89 const G4ThreeVector& aDirection, 90 const G4bool calcNorm = false, 91 G4bool* validNorm = nullptr, 92 G4ThreeVector* aNormalVector = nullptr) const override; 93 94 G4double DistanceToInNoVoxels(const G4ThreeVector& aPoint, 95 const G4ThreeVector& aDirection) const; 96 G4double DistanceToOutVoxels(const G4ThreeVector& aPoint, 97 const G4ThreeVector& aDirection, 98 G4ThreeVector* aNormalVector) const; 99 G4double DistanceToOutVoxelsCore(const G4ThreeVector& aPoint, 100 const G4ThreeVector& aDirection, 101 G4ThreeVector* aNormalVector, 102 G4bool& aConvex, 103 std::vector<G4int>& candidates) const; 104 G4double DistanceToOutNoVoxels(const G4ThreeVector& aPoint, 105 const G4ThreeVector& aDirection, 106 G4ThreeVector* aNormalVector) const; 107 108 G4ThreeVector SurfaceNormal(const G4ThreeVector& aPoint) const override; 109 110 void Extent(EAxis aAxis, G4double& aMin, G4double& aMax) const; 111 void BoundingLimits(G4ThreeVector& aMin, G4ThreeVector& aMax) const override; 112 G4bool CalculateExtent(const EAxis pAxis, 113 const G4VoxelLimits& pVoxelLimit, 114 const G4AffineTransform& pTransform, 115 G4double& pMin, G4double& pMax) const override; 116 G4double GetCubicVolume() override; 117 G4double GetSurfaceArea() override; 118 119 G4int GetNumOfConstituents() const override; 120 G4bool IsFaceted() const override; 121 122 G4VSolid* Clone() const override ; 123 124 G4GeometryType GetEntityType() const override { return "G4MultiUnion"; } 125 126 void Voxelize(); 127 // Finalize and prepare for use. User MUST call it once before 128 // navigation use. 129 130 EInside InsideNoVoxels(const G4ThreeVector& aPoint) const; 131 inline G4Voxelizer& GetVoxels() const; 132 133 std::ostream& StreamInfo(std::ostream& os) const override; 134 135 G4ThreeVector GetPointOnSurface() const override; 136 137 void DescribeYourselfTo ( G4VGraphicsScene& scene ) const override ; 138 G4Polyhedron* CreatePolyhedron () const override ; 139 G4Polyhedron* GetPolyhedron () const override; 140 141 G4MultiUnion(__void__&); 142 // Fake default constructor for usage restricted to direct object 143 // persistency for clients requiring preallocation of memory for 144 // persistifiable objects. 145 146 private: 147 148 EInside InsideWithExclusion(const G4ThreeVector& aPoint, 149 G4SurfBits* bits = nullptr) const; 150 G4int SafetyFromOutsideNumberNode(const G4ThreeVector& aPoint, 151 G4double& safety) const; 152 G4double DistanceToInCandidates(const G4ThreeVector& aPoint, 153 const G4ThreeVector& aDirection, 154 std::vector<G4int>& candidates, 155 G4SurfBits& bits) const; 156 157 // Conversion utilities 158 inline G4ThreeVector GetLocalPoint(const G4Transform3D& trans, 159 const G4ThreeVector& gpoint) const; 160 inline G4ThreeVector GetLocalVector(const G4Transform3D& trans, 161 const G4ThreeVector& gvec) const; 162 inline G4ThreeVector GetGlobalPoint(const G4Transform3D& trans, 163 const G4ThreeVector& lpoint) const; 164 inline G4ThreeVector GetGlobalVector(const G4Transform3D& trans, 165 const G4ThreeVector& lvec) const; 166 void TransformLimits(G4ThreeVector& min, G4ThreeVector& max, 167 const G4Transform3D& transformation) const; 168 private: 169 170 struct G4MultiUnionSurface 171 { 172 G4ThreeVector point; 173 G4VSolid* solid; 174 }; 175 176 std::vector<G4VSolid*> fSolids; 177 std::vector<G4Transform3D> fTransformObjs; 178 G4Voxelizer fVoxels; // Vozelizer for the solid 179 G4double fCubicVolume = 0.0; // Cubic Volume 180 G4double fSurfaceArea = 0.0; // Surface Area 181 G4double kRadTolerance; // Cached radial tolerance 182 mutable G4bool fAccurate = false; // Accurate safety (off by default) 183 184 mutable G4bool fRebuildPolyhedron = false; 185 mutable G4Polyhedron* fpPolyhedron = nullptr; 186 }; 187 188 //______________________________________________________________________________ 189 inline G4Voxelizer& G4MultiUnion::GetVoxels() const 190 { 191 return (G4Voxelizer&)fVoxels; 192 } 193 194 //______________________________________________________________________________ 195 inline const G4Transform3D& G4MultiUnion::GetTransformation(G4int index) const 196 { 197 return fTransformObjs[index]; 198 } 199 200 //______________________________________________________________________________ 201 inline G4VSolid* G4MultiUnion::GetSolid(G4int index) const 202 { 203 return fSolids[index]; 204 } 205 206 //______________________________________________________________________________ 207 inline G4int G4MultiUnion::GetNumberOfSolids() const 208 { 209 return G4int(fSolids.size()); 210 } 211 212 //______________________________________________________________________________ 213 inline void G4MultiUnion::SetAccurateSafety(G4bool flag) 214 { 215 fAccurate = flag; 216 } 217 218 //______________________________________________________________________________ 219 inline 220 G4ThreeVector G4MultiUnion::GetLocalPoint(const G4Transform3D& trans, 221 const G4ThreeVector& global) const 222 { 223 // Returns local point coordinates converted from the global frame defined 224 // by the transformation. This is defined by multiplying the inverse 225 // transformation with the global vector. 226 227 return trans.inverse()*G4Point3D(global); 228 } 229 230 //______________________________________________________________________________ 231 inline 232 G4ThreeVector G4MultiUnion::GetLocalVector(const G4Transform3D& trans, 233 const G4ThreeVector& global) const 234 { 235 // Returns local point coordinates converted from the global frame defined 236 // by the transformation. This is defined by multiplying the inverse 237 // transformation with the global vector. 238 239 G4Rotate3D rot; 240 G4Translate3D transl ; 241 G4Scale3D scale; 242 243 trans.getDecomposition(scale,rot,transl); 244 return rot.inverse()*G4Vector3D(global); 245 } 246 247 //______________________________________________________________________________ 248 inline 249 G4ThreeVector G4MultiUnion::GetGlobalPoint(const G4Transform3D& trans, 250 const G4ThreeVector& local) const 251 { 252 // Returns global point coordinates converted from the local frame defined 253 // by the transformation. This is defined by multiplying this transformation 254 // with the local vector. 255 256 return trans*G4Point3D(local); 257 } 258 259 //______________________________________________________________________________ 260 inline 261 G4ThreeVector G4MultiUnion::GetGlobalVector(const G4Transform3D& trans, 262 const G4ThreeVector& local) const 263 { 264 // Returns vector components converted from the local frame defined by the 265 // transformation to the global one. This is defined by multiplying this 266 // transformation with the local vector while ignoring the translation. 267 268 G4Rotate3D rot; 269 G4Translate3D transl ; 270 G4Scale3D scale; 271 272 trans.getDecomposition(scale,rot,transl); 273 return rot*G4Vector3D(local); 274 } 275 276 #endif 277