Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4MultiUnion.hh
이 파일의 문서화 페이지로 가기
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 //
27 // $Id:$
28 //
29 // --------------------------------------------------------------------
30 // GEANT 4 class header file
31 //
32 //
33 // G4MultiUnion
34 //
35 // Class description:
36 //
37 // An instance of "G4MultiUnion" constitutes a grouping of several solids.
38 // The constituent solids are stored with their respective location in an
39 // instance of "G4Node". An instance of "G4MultiUnion" is subsequently
40 // composed of one or several nodes.
41 
42 // History:
43 // 06.04.17 G.Cosmo - Imported implementation in Geant4 for VecGeom migration
44 // 19.10.12 M.Gayer - Original implementation from USolids module
45 // --------------------------------------------------------------------
46 #ifndef G4MULTIUNION_HH
47 #define G4MULTIUNION_HH
48 
49 #include <vector>
50 
51 #include "G4VSolid.hh"
52 #include "G4ThreeVector.hh"
53 #include "G4Transform3D.hh"
54 #include "G4Point3D.hh"
55 #include "G4Vector3D.hh"
56 #include "G4SurfBits.hh"
57 #include "G4Voxelizer.hh"
58 
59 class G4Polyhedron;
60 
61 class G4MultiUnion : public G4VSolid
62 {
63  friend class G4Voxelizer;
64 
65  public:
66 
68  G4MultiUnion(const G4String& name);
69  ~G4MultiUnion();
70 
71  // Build the multiple union by adding nodes
72  void AddNode(G4VSolid& solid, G4Transform3D& trans);
73 
74  G4MultiUnion(const G4MultiUnion& rhs);
75  G4MultiUnion& operator=(const G4MultiUnion& rhs);
76 
77  // Accessors
78  inline const G4Transform3D& GetTransformation(G4int index) const;
79  inline G4VSolid* GetSolid(G4int index) const;
80  inline G4int GetNumberOfSolids()const;
81 
82  // Navigation methods
83  EInside Inside(const G4ThreeVector& aPoint) const;
84 
85  EInside InsideIterator(const G4ThreeVector& aPoint) const;
86 
87  // Safety methods
88  G4double DistanceToIn(const G4ThreeVector& aPoint) const;
89  G4double DistanceToOut(const G4ThreeVector& aPoint) const;
90  inline void SetAccurateSafety(G4bool flag);
91 
92  // Exact distance methods
93  G4double DistanceToIn(const G4ThreeVector& aPoint,
94  const G4ThreeVector& aDirection) const;
95  G4double DistanceToOut(const G4ThreeVector& aPoint,
96  const G4ThreeVector& aDirection,
97  const G4bool calcNorm=false,
98  G4bool *validNorm=0,
99  G4ThreeVector *aNormalVector=0) const;
100 
102  const G4ThreeVector& aDirection) const;
104  const G4ThreeVector& aDirection,
105  G4ThreeVector* aNormalVector) const;
107  const G4ThreeVector& aDirection,
108  G4ThreeVector* aNormalVector,
109  G4bool& aConvex,
110  std::vector<G4int>& candidates) const;
112  const G4ThreeVector& aDirection,
113  G4ThreeVector* aNormalVector) const;
114 
115  G4ThreeVector SurfaceNormal(const G4ThreeVector& aPoint) const;
116 
117  void Extent(EAxis aAxis, G4double& aMin, G4double& aMax) const;
118  void BoundingLimits(G4ThreeVector& aMin, G4ThreeVector& aMax) const;
119  G4bool CalculateExtent(const EAxis pAxis,
120  const G4VoxelLimits& pVoxelLimit,
121  const G4AffineTransform& pTransform,
122  G4double& pMin, G4double& pMax) const;
125 
126  G4VSolid* Clone() const ;
127 
128  G4GeometryType GetEntityType() const { return "G4MultiUnion"; }
129 
130  void Voxelize();
131  // Finalize and prepare for use. User MUST call it once before
132  // navigation use.
133 
134  EInside InsideNoVoxels(const G4ThreeVector& aPoint) const;
135  inline G4Voxelizer& GetVoxels() const;
136 
137  std::ostream& StreamInfo(std::ostream& os) const;
138 
140 
141  void DescribeYourselfTo ( G4VGraphicsScene& scene ) const ;
142  G4Polyhedron* CreatePolyhedron () const ;
143  G4Polyhedron* GetPolyhedron () const;
144 
145  G4MultiUnion(__void__&);
146  // Fake default constructor for usage restricted to direct object
147  // persistency for clients requiring preallocation of memory for
148  // persistifiable objects.
149 
150  private:
151 
153  G4SurfBits* bits = 0) const;
155  G4double& safety) const;
157  const G4ThreeVector& aDirection,
158  std::vector<G4int>& candidates,
159  G4SurfBits& bits) const;
160 
161  // Conversion utilities
162  inline G4ThreeVector GetLocalPoint(const G4Transform3D& trans,
163  const G4ThreeVector& gpoint) const;
164  inline G4ThreeVector GetLocalVector(const G4Transform3D& trans,
165  const G4ThreeVector& gvec) const;
166  inline G4ThreeVector GetGlobalPoint(const G4Transform3D& trans,
167  const G4ThreeVector& lpoint) const;
168  inline G4ThreeVector GetGlobalVector(const G4Transform3D& trans,
169  const G4ThreeVector& lvec) const;
171  const G4Transform3D& transformation) const;
172  private:
173 
175  {
178  };
179 
180  std::vector<G4VSolid*> fSolids;
181  std::vector<G4Transform3D> fTransformObjs;
182  G4Voxelizer fVoxels; // Pointer to the vozelized solid
183  G4double fCubicVolume; // Cubic Volume
184  G4double fSurfaceArea; // Surface Area
185  G4double kRadTolerance; // Cached radial tolerance
186  mutable G4bool fAccurate; // Accurate safety (off by default)
187 
190 };
191 
192 //______________________________________________________________________________
194 {
195  return (G4Voxelizer&)fVoxels;
196 }
197 
198 //______________________________________________________________________________
200 {
201  return fTransformObjs[index];
202 }
203 
204 //______________________________________________________________________________
206 {
207  return fSolids[index];
208 }
209 
210 //______________________________________________________________________________
212 {
213  return fSolids.size();
214 }
215 
216 //______________________________________________________________________________
218 {
219  fAccurate = flag;
220 }
221 
222 //______________________________________________________________________________
223 inline
225  const G4ThreeVector& global) const
226 {
227  // Returns local point coordinates converted from the global frame defined
228  // by the transformation. This is defined by multiplying the inverse
229  // transformation with the global vector.
230 
231  return trans.inverse()*G4Point3D(global);
232 }
233 
234 //______________________________________________________________________________
235 inline
237  const G4ThreeVector& global) const
238 {
239  // Returns local point coordinates converted from the global frame defined
240  // by the transformation. This is defined by multiplying the inverse
241  // transformation with the global vector.
242 
243  G4Rotate3D rot;
244  G4Translate3D transl ;
246 
247  trans.getDecomposition(scale,rot,transl);
248  return rot.inverse()*G4Vector3D(global);
249 }
250 
251 //______________________________________________________________________________
252 inline
254  const G4ThreeVector& local) const
255 {
256  // Returns global point coordinates converted from the local frame defined
257  // by the transformation. This is defined by multiplying this transformation
258  // with the local vector.
259 
260  return trans*G4Point3D(local);
261 }
262 
263 //______________________________________________________________________________
264 inline
266  const G4ThreeVector& local) const
267 {
268  // Returns vector components converted from the local frame defined by the
269  // transformation to the global one. This is defined by multiplying this
270  // transformation with the local vector while ignoring the translation.
271 
272  G4Rotate3D rot;
273  G4Translate3D transl ;
275 
276  trans.getDecomposition(scale,rot,transl);
277  return rot*G4Vector3D(local);
278 }
279 
280 #endif
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const XML_Char * name
Definition: expat.h:151
void SetAccurateSafety(G4bool flag)
G4VSolid * GetSolid(G4int index) const
G4GeometryType GetEntityType() const
G4Polyhedron * CreatePolyhedron() const
std::vector< G4Transform3D > fTransformObjs
G4Voxelizer & GetVoxels() const
G4ThreeVector GetGlobalPoint(const G4Transform3D &trans, const G4ThreeVector &lpoint) const
void AddNode(G4VSolid &solid, G4Transform3D &trans)
Definition: G4MultiUnion.cc:77
G4bool fAccurate
G4ThreeVector SurfaceNormal(const G4ThreeVector &aPoint) const
std::ostream & StreamInfo(std::ostream &os) const
Transform3D inverse() const
Definition: Transform3D.cc:142
G4ThreeVector GetPointOnSurface() const
G4VSolid * Clone() const
Definition: G4MultiUnion.cc:84
G4double DistanceToInNoVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection) const
G4Polyhedron * GetPolyhedron() const
G4double fSurfaceArea
G4double DistanceToIn(const G4ThreeVector &aPoint) const
EInside InsideNoVoxels(const G4ThreeVector &aPoint) const
G4double GetSurfaceArea()
G4ThreeVector GetLocalPoint(const G4Transform3D &trans, const G4ThreeVector &gpoint) const
G4ThreeVector GetGlobalVector(const G4Transform3D &trans, const G4ThreeVector &lvec) const
G4double kRadTolerance
G4double DistanceToOutVoxelsCore(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector, G4bool &aConvex, std::vector< G4int > &candidates) const
EInside InsideIterator(const G4ThreeVector &aPoint) const
G4double DistanceToOutNoVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
Double_t scale
void BoundingLimits(G4ThreeVector &aMin, G4ThreeVector &aMax) const
void TransformLimits(G4ThreeVector &min, G4ThreeVector &max, const G4Transform3D &transformation) const
#define local
Definition: gzguts.h:114
const G4Transform3D & GetTransformation(G4int index) const
G4double GetCubicVolume()
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35
G4bool fRebuildPolyhedron
EInside InsideWithExclusion(const G4ThreeVector &aPoint, G4SurfBits *bits=0) const
G4ThreeVector GetLocalVector(const G4Transform3D &trans, const G4ThreeVector &gvec) const
G4Polyhedron * fpPolyhedron
G4double fCubicVolume
int G4int
Definition: G4Types.hh:78
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
G4double DistanceToOutVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const
std::vector< G4VSolid * > fSolids
G4MultiUnion & operator=(const G4MultiUnion &rhs)
void getDecomposition(Scale3D &scale, Rotate3D &rotation, Translate3D &translation) const
Definition: Transform3D.cc:174
void DescribeYourselfTo(G4VGraphicsScene &scene) const
void Extent(EAxis aAxis, G4double &aMin, G4double &aMax) const
EInside Inside(const G4ThreeVector &aPoint) const
G4Voxelizer fVoxels
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
G4double DistanceToOut(const G4ThreeVector &aPoint) const
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:35
G4int GetNumberOfSolids() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double DistanceToInCandidates(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, std::vector< G4int > &candidates, G4SurfBits &bits) const
G4int SafetyFromOutsideNumberNode(const G4ThreeVector &aPoint, G4double &safety) const