Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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, G4Transform3D& trans);
65
66 G4MultiUnion(const G4MultiUnion& rhs);
68
69 // Accessors
70 inline const G4Transform3D& GetTransformation(G4int index) const;
71 inline G4VSolid* GetSolid(G4int index) const;
72 inline G4int GetNumberOfSolids()const;
73
74 // Navigation methods
75 EInside Inside(const G4ThreeVector& aPoint) const;
76
77 EInside InsideIterator(const G4ThreeVector& aPoint) const;
78
79 // Safety methods
80 G4double DistanceToIn(const G4ThreeVector& aPoint) const;
81 G4double DistanceToOut(const G4ThreeVector& aPoint) const;
82 inline void SetAccurateSafety(G4bool flag);
83
84 // Exact distance methods
86 const G4ThreeVector& aDirection) const;
88 const G4ThreeVector& aDirection,
89 const G4bool calcNorm = false,
90 G4bool* validNorm = nullptr,
91 G4ThreeVector* aNormalVector = nullptr) const;
92
94 const G4ThreeVector& aDirection) const;
96 const G4ThreeVector& aDirection,
97 G4ThreeVector* aNormalVector) const;
99 const G4ThreeVector& aDirection,
100 G4ThreeVector* aNormalVector,
101 G4bool& aConvex,
102 std::vector<G4int>& candidates) const;
104 const G4ThreeVector& aDirection,
105 G4ThreeVector* aNormalVector) const;
106
107 G4ThreeVector SurfaceNormal(const G4ThreeVector& aPoint) const;
108
109 void Extent(EAxis aAxis, G4double& aMin, G4double& aMax) const;
110 void BoundingLimits(G4ThreeVector& aMin, G4ThreeVector& aMax) const;
111 G4bool CalculateExtent(const EAxis pAxis,
112 const G4VoxelLimits& pVoxelLimit,
113 const G4AffineTransform& pTransform,
114 G4double& pMin, G4double& pMax) const;
117
118 G4VSolid* Clone() const ;
119
120 G4GeometryType GetEntityType() const { return "G4MultiUnion"; }
121
122 void Voxelize();
123 // Finalize and prepare for use. User MUST call it once before
124 // navigation use.
125
126 EInside InsideNoVoxels(const G4ThreeVector& aPoint) const;
127 inline G4Voxelizer& GetVoxels() const;
128
129 std::ostream& StreamInfo(std::ostream& os) const;
130
132
133 void DescribeYourselfTo ( G4VGraphicsScene& scene ) const ;
135 G4Polyhedron* GetPolyhedron () const;
136
137 G4MultiUnion(__void__&);
138 // Fake default constructor for usage restricted to direct object
139 // persistency for clients requiring preallocation of memory for
140 // persistifiable objects.
141
142 private:
143
144 EInside InsideWithExclusion(const G4ThreeVector& aPoint,
145 G4SurfBits* bits = 0) const;
146 G4int SafetyFromOutsideNumberNode(const G4ThreeVector& aPoint,
147 G4double& safety) const;
148 G4double DistanceToInCandidates(const G4ThreeVector& aPoint,
149 const G4ThreeVector& aDirection,
150 std::vector<G4int>& candidates,
151 G4SurfBits& bits) const;
152
153 // Conversion utilities
154 inline G4ThreeVector GetLocalPoint(const G4Transform3D& trans,
155 const G4ThreeVector& gpoint) const;
156 inline G4ThreeVector GetLocalVector(const G4Transform3D& trans,
157 const G4ThreeVector& gvec) const;
158 inline G4ThreeVector GetGlobalPoint(const G4Transform3D& trans,
159 const G4ThreeVector& lpoint) const;
160 inline G4ThreeVector GetGlobalVector(const G4Transform3D& trans,
161 const G4ThreeVector& lvec) const;
162 void TransformLimits(G4ThreeVector& min, G4ThreeVector& max,
163 const G4Transform3D& transformation) const;
164 private:
165
166 struct G4MultiUnionSurface
167 {
168 G4ThreeVector point;
169 G4VSolid* solid;
170 };
171
172 std::vector<G4VSolid*> fSolids;
173 std::vector<G4Transform3D> fTransformObjs;
174 G4Voxelizer fVoxels; // Vozelizer for the solid
175 G4double fCubicVolume = 0.0; // Cubic Volume
176 G4double fSurfaceArea = 0.0; // Surface Area
177 G4double kRadTolerance; // Cached radial tolerance
178 mutable G4bool fAccurate = false; // Accurate safety (off by default)
179
180 mutable G4bool fRebuildPolyhedron = false;
181 mutable G4Polyhedron* fpPolyhedron = nullptr;
182};
183
184//______________________________________________________________________________
186{
187 return (G4Voxelizer&)fVoxels;
188}
189
190//______________________________________________________________________________
192{
193 return fTransformObjs[index];
194}
195
196//______________________________________________________________________________
198{
199 return fSolids[index];
200}
201
202//______________________________________________________________________________
204{
205 return G4int(fSolids.size());
206}
207
208//______________________________________________________________________________
210{
211 fAccurate = flag;
212}
213
214//______________________________________________________________________________
215inline
216G4ThreeVector G4MultiUnion::GetLocalPoint(const G4Transform3D& trans,
217 const G4ThreeVector& global) const
218{
219 // Returns local point coordinates converted from the global frame defined
220 // by the transformation. This is defined by multiplying the inverse
221 // transformation with the global vector.
222
223 return trans.inverse()*G4Point3D(global);
224}
225
226//______________________________________________________________________________
227inline
228G4ThreeVector G4MultiUnion::GetLocalVector(const G4Transform3D& trans,
229 const G4ThreeVector& global) const
230{
231 // Returns local point coordinates converted from the global frame defined
232 // by the transformation. This is defined by multiplying the inverse
233 // transformation with the global vector.
234
235 G4Rotate3D rot;
236 G4Translate3D transl ;
237 G4Scale3D scale;
238
239 trans.getDecomposition(scale,rot,transl);
240 return rot.inverse()*G4Vector3D(global);
241}
242
243//______________________________________________________________________________
244inline
245G4ThreeVector G4MultiUnion::GetGlobalPoint(const G4Transform3D& trans,
246 const G4ThreeVector& local) const
247{
248 // Returns global point coordinates converted from the local frame defined
249 // by the transformation. This is defined by multiplying this transformation
250 // with the local vector.
251
252 return trans*G4Point3D(local);
253}
254
255//______________________________________________________________________________
256inline
257G4ThreeVector G4MultiUnion::GetGlobalVector(const G4Transform3D& trans,
258 const G4ThreeVector& local) const
259{
260 // Returns vector components converted from the local frame defined by the
261 // transformation to the global one. This is defined by multiplying this
262 // transformation with the local vector while ignoring the translation.
263
264 G4Rotate3D rot;
265 G4Translate3D transl ;
266 G4Scale3D scale;
267
268 trans.getDecomposition(scale,rot,transl);
269 return rot*G4Vector3D(local);
270}
271
272#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 AddNode(G4VSolid &solid, G4Transform3D &trans)
Definition: G4MultiUnion.cc:69
void DescribeYourselfTo(G4VGraphicsScene &scene) const
void Extent(EAxis aAxis, G4double &aMin, G4double &aMax) const
G4VSolid * Clone() const
Definition: G4MultiUnion.cc:76
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)
Definition: G4MultiUnion.cc:99
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