Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4MultiUnion.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// 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
51class G4Polyhedron;
52
53class G4MultiUnion : public G4VSolid
54{
55 friend class G4Voxelizer;
56
57 public:
58
60 G4MultiUnion(const G4String& name);
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);
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;
77
78 EInside InsideIterator(const G4ThreeVector& aPoint) const;
79
80 // Safety methods
81 G4double DistanceToIn(const G4ThreeVector& aPoint) const;
82 G4double DistanceToOut(const G4ThreeVector& aPoint) const;
83 inline void SetAccurateSafety(G4bool flag);
84
85 // Exact distance methods
87 const G4ThreeVector& aDirection) const;
89 const G4ThreeVector& aDirection,
90 const G4bool calcNorm = false,
91 G4bool* validNorm = nullptr,
92 G4ThreeVector* aNormalVector = nullptr) const;
93
95 const G4ThreeVector& aDirection) const;
97 const G4ThreeVector& aDirection,
98 G4ThreeVector* aNormalVector) const;
100 const G4ThreeVector& aDirection,
101 G4ThreeVector* aNormalVector,
102 G4bool& aConvex,
103 std::vector<G4int>& candidates) const;
105 const G4ThreeVector& aDirection,
106 G4ThreeVector* aNormalVector) const;
107
108 G4ThreeVector SurfaceNormal(const G4ThreeVector& aPoint) const;
109
110 void Extent(EAxis aAxis, G4double& aMin, G4double& aMax) const;
111 void BoundingLimits(G4ThreeVector& aMin, G4ThreeVector& aMax) const;
112 G4bool CalculateExtent(const EAxis pAxis,
113 const G4VoxelLimits& pVoxelLimit,
114 const G4AffineTransform& pTransform,
115 G4double& pMin, G4double& pMax) const;
118
119 G4VSolid* Clone() const ;
120
121 G4GeometryType GetEntityType() const { return "G4MultiUnion"; }
122
123 void Voxelize();
124 // Finalize and prepare for use. User MUST call it once before
125 // navigation use.
126
127 EInside InsideNoVoxels(const G4ThreeVector& aPoint) const;
128 inline G4Voxelizer& GetVoxels() const;
129
130 std::ostream& StreamInfo(std::ostream& os) const;
131
133
134 void DescribeYourselfTo ( G4VGraphicsScene& scene ) const ;
136 G4Polyhedron* GetPolyhedron () const;
137
138 G4MultiUnion(__void__&);
139 // Fake default constructor for usage restricted to direct object
140 // persistency for clients requiring preallocation of memory for
141 // persistifiable objects.
142
143 private:
144
145 EInside InsideWithExclusion(const G4ThreeVector& aPoint,
146 G4SurfBits* bits = 0) const;
147 G4int SafetyFromOutsideNumberNode(const G4ThreeVector& aPoint,
148 G4double& safety) const;
149 G4double DistanceToInCandidates(const G4ThreeVector& aPoint,
150 const G4ThreeVector& aDirection,
151 std::vector<G4int>& candidates,
152 G4SurfBits& bits) const;
153
154 // Conversion utilities
155 inline G4ThreeVector GetLocalPoint(const G4Transform3D& trans,
156 const G4ThreeVector& gpoint) const;
157 inline G4ThreeVector GetLocalVector(const G4Transform3D& trans,
158 const G4ThreeVector& gvec) const;
159 inline G4ThreeVector GetGlobalPoint(const G4Transform3D& trans,
160 const G4ThreeVector& lpoint) const;
161 inline G4ThreeVector GetGlobalVector(const G4Transform3D& trans,
162 const G4ThreeVector& lvec) const;
163 void TransformLimits(G4ThreeVector& min, G4ThreeVector& max,
164 const G4Transform3D& transformation) const;
165 private:
166
167 struct G4MultiUnionSurface
168 {
169 G4ThreeVector point;
170 G4VSolid* solid;
171 };
172
173 std::vector<G4VSolid*> fSolids;
174 std::vector<G4Transform3D> fTransformObjs;
175 G4Voxelizer fVoxels; // Vozelizer for the solid
176 G4double fCubicVolume = 0.0; // Cubic Volume
177 G4double fSurfaceArea = 0.0; // Surface Area
178 G4double kRadTolerance; // Cached radial tolerance
179 mutable G4bool fAccurate = false; // Accurate safety (off by default)
180
181 mutable G4bool fRebuildPolyhedron = false;
182 mutable G4Polyhedron* fpPolyhedron = nullptr;
183};
184
185//______________________________________________________________________________
187{
188 return (G4Voxelizer&)fVoxels;
189}
190
191//______________________________________________________________________________
193{
194 return fTransformObjs[index];
195}
196
197//______________________________________________________________________________
199{
200 return fSolids[index];
201}
202
203//______________________________________________________________________________
205{
206 return G4int(fSolids.size());
207}
208
209//______________________________________________________________________________
211{
212 fAccurate = flag;
213}
214
215//______________________________________________________________________________
216inline
217G4ThreeVector G4MultiUnion::GetLocalPoint(const G4Transform3D& trans,
218 const G4ThreeVector& global) const
219{
220 // Returns local point coordinates converted from the global frame defined
221 // by the transformation. This is defined by multiplying the inverse
222 // transformation with the global vector.
223
224 return trans.inverse()*G4Point3D(global);
225}
226
227//______________________________________________________________________________
228inline
229G4ThreeVector G4MultiUnion::GetLocalVector(const G4Transform3D& trans,
230 const G4ThreeVector& global) const
231{
232 // Returns local point coordinates converted from the global frame defined
233 // by the transformation. This is defined by multiplying the inverse
234 // transformation with the global vector.
235
236 G4Rotate3D rot;
237 G4Translate3D transl ;
238 G4Scale3D scale;
239
240 trans.getDecomposition(scale,rot,transl);
241 return rot.inverse()*G4Vector3D(global);
242}
243
244//______________________________________________________________________________
245inline
246G4ThreeVector G4MultiUnion::GetGlobalPoint(const G4Transform3D& trans,
247 const G4ThreeVector& local) const
248{
249 // Returns global point coordinates converted from the local frame defined
250 // by the transformation. This is defined by multiplying this transformation
251 // with the local vector.
252
253 return trans*G4Point3D(local);
254}
255
256//______________________________________________________________________________
257inline
258G4ThreeVector G4MultiUnion::GetGlobalVector(const G4Transform3D& trans,
259 const G4ThreeVector& local) const
260{
261 // Returns vector components converted from the local frame defined by the
262 // transformation to the global one. This is defined by multiplying this
263 // transformation with the local vector while ignoring the translation.
264
265 G4Rotate3D rot;
266 G4Translate3D transl ;
267 G4Scale3D scale;
268
269 trans.getDecomposition(scale,rot,transl);
270 return rot*G4Vector3D(local);
271}
272
273#endif
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:34
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:34
const G4Transform3D & GetTransformation(G4int index) const
G4int GetNumberOfSolids() const
EInside InsideIterator(const G4ThreeVector &aPoint) const
G4double DistanceToIn(const G4ThreeVector &aPoint) const
G4double GetSurfaceArea()
G4double DistanceToOutNoVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const
G4double GetCubicVolume()
EInside InsideNoVoxels(const G4ThreeVector &aPoint) const
G4double DistanceToInNoVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection) const
EInside Inside(const G4ThreeVector &aPoint) const
G4double DistanceToOut(const G4ThreeVector &aPoint) const
G4Polyhedron * CreatePolyhedron() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
void BoundingLimits(G4ThreeVector &aMin, G4ThreeVector &aMax) const
G4double DistanceToOutVoxelsCore(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector, G4bool &aConvex, std::vector< G4int > &candidates) const
G4Polyhedron * GetPolyhedron() const
G4GeometryType GetEntityType() const
G4double DistanceToOutVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const
std::ostream & StreamInfo(std::ostream &os) const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
void AddNode(G4VSolid &solid, const G4Transform3D &trans)
Definition: G4MultiUnion.cc:69
void Extent(EAxis aAxis, G4double &aMin, G4double &aMax) const
G4VSolid * Clone() const
Definition: G4MultiUnion.cc:83
G4ThreeVector GetPointOnSurface() const
G4VSolid * GetSolid(G4int index) const
void SetAccurateSafety(G4bool flag)
G4Voxelizer & GetVoxels() const
G4ThreeVector SurfaceNormal(const G4ThreeVector &aPoint) const
G4MultiUnion & operator=(const G4MultiUnion &rhs)
void getDecomposition(Scale3D &scale, Rotate3D &rotation, Translate3D &translation) const
Definition: Transform3D.cc:173
Transform3D inverse() const
Definition: Transform3D.cc:141
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67
#define local
Definition: gzguts.h:114