Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4UTessellatedSolid.cc
이 파일의 문서화 페이지로 가기
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 // Implementation of G4UTessellatedSolid wrapper class
31 // --------------------------------------------------------------------
32 
33 #include "G4TessellatedSolid.hh"
34 #include "G4UTessellatedSolid.hh"
35 
36 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
37 
38 #include "G4TriangularFacet.hh"
39 #include "G4QuadrangularFacet.hh"
40 
41 #include "G4GeomTools.hh"
42 #include "G4AffineTransform.hh"
43 #include "G4BoundingEnvelope.hh"
44 
45 #include "G4PolyhedronArbitrary.hh"
46 
48 //
49 // Constructors
50 //
51 G4UTessellatedSolid::G4UTessellatedSolid()
52  : Base_t("")
53 {
54 }
55 
56 
57 G4UTessellatedSolid::G4UTessellatedSolid(const G4String& name)
58  : Base_t(name)
59 {
60 }
61 
63 //
64 // Fake default constructor - sets only member data and allocates memory
65 // for usage restricted to object persistency.
66 //
67 G4UTessellatedSolid::G4UTessellatedSolid(__void__& a)
68  : Base_t(a)
69 {
70 }
71 
72 
74 //
75 // Destructor
76 //
77 G4UTessellatedSolid::~G4UTessellatedSolid()
78 {
79  G4int size = fFacets.size();
80  for (G4int i = 0; i < size; ++i) { delete fFacets[i]; }
81  fFacets.clear();
82 }
83 
84 
86 //
87 // Copy constructor
88 //
89 G4UTessellatedSolid::G4UTessellatedSolid(const G4UTessellatedSolid& source)
90  : Base_t(source)
91 {
92 }
93 
94 
96 //
97 // Assignment operator
98 //
99 G4UTessellatedSolid&
100 G4UTessellatedSolid::operator=(const G4UTessellatedSolid& source)
101 {
102  if (this == &source) return *this;
103 
104  Base_t::operator=( source );
105 
106  return *this;
107 }
108 
109 
111 //
112 // Accessors
113 
114 G4bool G4UTessellatedSolid::AddFacet(G4VFacet* aFacet)
115 {
116  // Add a facet to the structure, checking validity.
117  //
118  if (GetSolidClosed())
119  {
120  G4Exception("G4UTessellatedSolid::AddFacet()", "GeomSolids1002",
121  JustWarning, "Attempt to add facets when solid is closed.");
122  return false;
123  }
124  if (!aFacet->IsDefined())
125  {
126  G4Exception("G4UTessellatedSolid::AddFacet()", "GeomSolids1002",
127  JustWarning, "Attempt to add facet not properly defined.");
128  aFacet->StreamInfo(G4cout);
129  return false;
130  }
131  if (aFacet->GetNumberOfVertices() == 3)
132  {
133  G4TriangularFacet* a3Facet = dynamic_cast<G4TriangularFacet*>(aFacet);
134  return Base_t::AddTriangularFacet(UVector3(a3Facet->GetVertex(0).x(),
135  a3Facet->GetVertex(0).y(),
136  a3Facet->GetVertex(0).z()),
137  UVector3(a3Facet->GetVertex(1).x(),
138  a3Facet->GetVertex(1).y(),
139  a3Facet->GetVertex(1).z()),
140  UVector3(a3Facet->GetVertex(2).x(),
141  a3Facet->GetVertex(2).y(),
142  a3Facet->GetVertex(2).z()),
143  true);
144  }
145  else if (aFacet->GetNumberOfVertices() == 4)
146  {
147  G4QuadrangularFacet* a4Facet = dynamic_cast<G4QuadrangularFacet*>(aFacet);
148  return Base_t::AddQuadrilateralFacet(UVector3(a4Facet->GetVertex(0).x(),
149  a4Facet->GetVertex(0).y(),
150  a4Facet->GetVertex(0).z()),
151  UVector3(a4Facet->GetVertex(1).x(),
152  a4Facet->GetVertex(1).y(),
153  a4Facet->GetVertex(1).z()),
154  UVector3(a4Facet->GetVertex(2).x(),
155  a4Facet->GetVertex(2).y(),
156  a4Facet->GetVertex(2).z()),
157  UVector3(a4Facet->GetVertex(3).x(),
158  a4Facet->GetVertex(3).y(),
159  a4Facet->GetVertex(3).z()),
160  true);
161  }
162  else
163  {
164  G4Exception("G4UTessellatedSolid::AddFacet()", "GeomSolids1002",
165  JustWarning, "Attempt to add facet not properly defined.");
166  aFacet->StreamInfo(G4cout);
167  return false;
168  }
169 }
170 
171 G4VFacet* G4UTessellatedSolid::GetFacet(G4int i) const
172 {
173  return fFacets[i];
174 }
175 
176 G4int G4UTessellatedSolid::GetNumberOfFacets() const
177 {
178  return GetNFacets();
179 }
180 
181 void G4UTessellatedSolid::SetSolidClosed(const G4bool t)
182 {
183  if (t && !Base_t::IsClosed())
184  {
185  Base_t::Close();
186  G4int nVertices = fTessellated.fVertices.size();
187  G4int nFacets = fTessellated.fFacets.size();
188  for (G4int j = 0; j < nVertices; ++j)
189  {
190  UVector3 vt = fTessellated.fVertices[j];
191  fVertexList.push_back(G4ThreeVector(vt.x(), vt.y(), vt.z()));
192  }
193  for (G4int i = 0; i < nFacets; ++i)
194  {
195  vecgeom::TriangleFacet<G4double>* afacet = Base_t::GetFacet(i);
196  std::vector<G4ThreeVector> v;
197  for (G4int k=0; k<3; ++k)
198  {
199  v.push_back(G4ThreeVector(afacet->fVertices[k].x(),
200  afacet->fVertices[k].y(),
201  afacet->fVertices[k].z()));
202  }
203  G4VFacet* facet = new G4TriangularFacet(v[0], v[1], v[2],
205  facet->SetVertices(&fVertexList);
206  for (G4int k=0; k<3; ++k)
207  {
208  facet->SetVertexIndex(k, afacet->fIndices[k]);
209  }
210  fFacets.push_back(facet);
211  }
212  }
213 }
214 
215 G4bool G4UTessellatedSolid::GetSolidClosed() const
216 {
217  return Base_t::IsClosed();
218 }
219 
220 void G4UTessellatedSolid::SetMaxVoxels(G4int)
221 {
222  // Not yet implemented !
223 }
224 
225 G4double G4UTessellatedSolid::GetMinXExtent() const
226 {
227  UVector3 aMin, aMax;
228  Base_t::Extent(aMin, aMax);
229  return aMin.x();
230 }
231 G4double G4UTessellatedSolid::GetMaxXExtent() const
232 {
233  UVector3 aMin, aMax;
234  Base_t::Extent(aMin, aMax);
235  return aMax.x();
236 }
237 G4double G4UTessellatedSolid::GetMinYExtent() const
238 {
239  UVector3 aMin, aMax;
240  Base_t::Extent(aMin, aMax);
241  return aMin.y();
242 }
243 G4double G4UTessellatedSolid::GetMaxYExtent() const
244 {
245  UVector3 aMin, aMax;
246  Base_t::Extent(aMin, aMax);
247  return aMax.y();
248 }
249 G4double G4UTessellatedSolid::GetMinZExtent() const
250 {
251  UVector3 aMin, aMax;
252  Base_t::Extent(aMin, aMax);
253  return aMin.z();
254 }
255 G4double G4UTessellatedSolid::GetMaxZExtent() const
256 {
257  UVector3 aMin, aMax;
258  Base_t::Extent(aMin, aMax);
259  return aMax.z();
260 }
261 
262 G4int G4UTessellatedSolid::AllocatedMemoryWithoutVoxels()
263 {
264  G4int base = sizeof(*this);
265  base += fVertexList.capacity() * sizeof(G4ThreeVector);
266 
267  G4int limit = fFacets.size();
268  for (G4int i = 0; i < limit; i++)
269  {
270  G4VFacet &facet = *fFacets[i];
271  base += facet.AllocatedMemory();
272  }
273  return base;
274 }
275 G4int G4UTessellatedSolid::AllocatedMemory()
276 {
277  return AllocatedMemoryWithoutVoxels();
278 }
279 void G4UTessellatedSolid::DisplayAllocatedMemory()
280 {
281  G4int without = AllocatedMemoryWithoutVoxels();
282  // G4int with = AllocatedMemory();
283  // G4double ratio = (G4double) with / without;
284  // G4cout << "G4TessellatedSolid - Allocated memory without voxel overhead "
285  // << without << "; with " << with << "; ratio: " << ratio << G4endl;
286  G4cout << "G4TessellatedSolid - Allocated memory without voxel overhead "
287  << without << G4endl;
288 }
289 
290 
292 //
293 // Get bounding box
294 
295 void G4UTessellatedSolid::BoundingLimits(G4ThreeVector& pMin,
296  G4ThreeVector& pMax) const
297 {
298  UVector3 aMin, aMax;
299  Base_t::Extent(aMin, aMax);
300  pMin = G4ThreeVector(aMin.x(), aMin.y(), aMin.z());
301  pMax = G4ThreeVector(aMax.x(), aMax.y(), aMax.z());
302 
303  // Check correctness of the bounding box
304  //
305  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
306  {
307  std::ostringstream message;
308  message << "Bad bounding box (min >= max) for solid: "
309  << GetName() << " !"
310  << "\npMin = " << pMin
311  << "\npMax = " << pMax;
312  G4Exception("G4UTessellatedSolid::BoundingLimits()",
313  "GeomMgt0001", JustWarning, message);
314  StreamInfo(G4cout);
315  }
316 }
317 
318 
320 //
321 // Calculate extent under transform and specified limit
322 
323 G4bool
324 G4UTessellatedSolid::CalculateExtent(const EAxis pAxis,
325  const G4VoxelLimits& pVoxelLimit,
326  const G4AffineTransform& pTransform,
327  G4double& pMin, G4double& pMax) const
328 {
329  G4ThreeVector bmin, bmax;
330  G4bool exist;
331  G4double kCarToleranceHalf = 0.5*kCarTolerance;
332 
333  // Check bounding box (bbox)
334  //
335  BoundingLimits(bmin,bmax);
336  G4BoundingEnvelope bbox(bmin,bmax);
337 #ifdef G4BBOX_EXTENT
338  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
339 #endif
340  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
341  {
342  return exist = (pMin < pMax) ? true : false;
343  }
344 
345  // The extent is calculated as cumulative extent of the pyramids
346  // formed by facets and the center of the bounding box.
347  //
348  G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
349  G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
350 
352  G4ThreeVectorList apex(1);
353  std::vector<const G4ThreeVectorList *> pyramid(2);
354  pyramid[0] = &base;
355  pyramid[1] = &apex;
356  apex[0] = (bmin+bmax)*0.5;
357 
358  // main loop along facets
359  pMin = kInfinity;
360  pMax = -kInfinity;
361  for (G4int i=0; i<GetNumberOfFacets(); ++i)
362  {
363  G4VFacet* facet = GetFacet(i);
364  if (std::abs((facet->GetSurfaceNormal()).dot(facet->GetVertex(0)-apex[0]))
365  < kCarToleranceHalf) continue;
366 
367  base.resize(3);
368  for (G4int k=0; k<3; ++k) { base[k] = facet->GetVertex(k); }
369  G4double emin,emax;
370  G4BoundingEnvelope benv(pyramid);
371  if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
372  if (emin < pMin) pMin = emin;
373  if (emax > pMax) pMax = emax;
374  if (eminlim > pMin && emaxlim < pMax) break; // max possible extent
375  }
376  return (pMin < pMax);
377 }
378 
379 
381 //
382 // CreatePolyhedron()
383 //
384 G4Polyhedron* G4UTessellatedSolid::CreatePolyhedron () const
385 {
386  G4int nVertices = fVertexList.size();
387  G4int nFacets = fFacets.size();
388  G4PolyhedronArbitrary *polyhedron = new G4PolyhedronArbitrary (nVertices,
389  nFacets);
390  for (G4int j = 0; j < nVertices; ++j)
391  {
392  polyhedron->AddVertex(fVertexList[j]);
393  }
394 
395  for (G4int i = 0; i < nFacets; ++i)
396  {
397  G4int v[3]; // Only facets with 3 vertices are defined in VecGeom
398  G4VFacet* facet = GetFacet(i);
399  for (G4int j=0; j<3; ++j) // Retrieve indexing directly from VecGeom
400  {
401  v[j] = facet->GetVertexIndex(j) + 1;
402  }
403  polyhedron->AddFacet(v[0],v[1],v[2]);
404  }
405  polyhedron->SetReferences();
406 
407  return (G4Polyhedron*) polyhedron;
408 }
409 
410 #endif // G4GEOM_USE_USOLIDS
virtual G4int GetVertexIndex(G4int i) const =0
virtual G4ThreeVector GetSurfaceNormal() const =0
const XML_Char * name
Definition: expat.h:151
CLHEP::Hep3Vector G4ThreeVector
G4double GetMaxExtent(const EAxis pAxis) const
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
static const G4double kInfinity
Definition: geomdefs.hh:42
void AddFacet(const G4int iv1, const G4int iv2, const G4int iv3, const G4int iv4=0)
void AddVertex(const G4ThreeVector &v)
#define G4endl
Definition: G4ios.hh:61
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4VFacet.cc:98
void message(RunManager *runmanager)
Definition: ts_scorers.cc:72
double z() const
static const G4double emax
const G4double kCarTolerance
virtual void SetVertexIndex(G4int i, G4int j)=0
virtual G4ThreeVector GetVertex(G4int i) const =0
const XML_Char int const XML_Char int const XML_Char * base
Definition: expat.h:331
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
G4double GetMinExtent(const EAxis pAxis) const
virtual G4bool IsDefined() const =0
virtual G4int GetNumberOfVertices() const =0
c1 Close()
G4ThreeVector GetVertex(G4int i) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
int G4int
Definition: G4Types.hh:78
EAxis
Definition: geomdefs.hh:54
G4GLOB_DLL std::ostream G4cout
double x() const
G4ThreeVector GetVertex(G4int i) const
double y() const
virtual void SetVertices(std::vector< G4ThreeVector > *vertices)=0
virtual G4int AllocatedMemory()=0