Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TessellatedSolid.cc
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 and of QinetiQ Ltd, *
20// * subject to DEFCON 705 IPR conditions. *
21// * By using, copying, modifying or distributing the software (or *
22// * any work based on the software) you agree to acknowledge its *
23// * use in resulting scientific publications, and indicate your *
24// * acceptance of all terms of the Geant4 Software license. *
25// ********************************************************************
26//
27// G4TessellatedSolid implementation
28//
29// 31.10.2004, P R Truscott, QinetiQ Ltd, UK - Created.
30// 17.09.2007, P R Truscott, QinetiQ Ltd & Richard Holmberg
31// Updated extensively prior to this date to deal with
32// concaved tessellated surfaces, based on the algorithm
33// of Richard Holmberg. This had been slightly modified
34// to determine with inside the geometry by projecting
35// random rays from the point provided. Now random rays
36// are predefined rather than making use of random
37// number generator at run-time.
38// 12.10.2012, M Gayer, CERN, complete rewrite reducing memory
39// requirements more than 50% and speedup by a factor of
40// tens or more depending on the number of facets, thanks
41// to voxelization of surface and improvements.
42// Speedup factor of thousands for solids with number of
43// facets in hundreds of thousands.
44// 23.10.2016, E Tcherniaev, reimplemented CalculateExtent() to make
45// use of G4BoundingEnvelope.
46// --------------------------------------------------------------------
47
48#include "G4TessellatedSolid.hh"
49
50#if !defined(G4GEOM_USE_UTESSELLATEDSOLID)
51
52#include <algorithm>
53#include <fstream>
54#include <iomanip>
55#include <iostream>
56#include <list>
57#include <random>
58#include <stack>
59
60#include "geomdefs.hh"
61#include "G4SystemOfUnits.hh"
64#include "G4VoxelLimits.hh"
65#include "G4AffineTransform.hh"
66#include "G4BoundingEnvelope.hh"
67#include "G4QuickRand.hh"
68
69#include "G4VGraphicsScene.hh"
70#include "G4VisExtent.hh"
71
72#include "G4AutoLock.hh"
73
74namespace
75{
76 G4Mutex polyhedronMutex = G4MUTEX_INITIALIZER;
77 G4Mutex tessMutex = G4MUTEX_INITIALIZER;
78}
79
80using namespace std;
81
82///////////////////////////////////////////////////////////////////////////////
83//
84// Standard contructor has blank name and defines no fFacets.
85//
87{
88 Initialize();
89}
90
91///////////////////////////////////////////////////////////////////////////////
92//
93// Alternative constructor. Simple define name and geometry type - no fFacets
94// to detine.
95//
97 : G4VSolid(name)
98{
99 Initialize();
100}
101
102///////////////////////////////////////////////////////////////////////////////
103//
104// Fake default constructor - sets only member data and allocates memory
105// for usage restricted to object persistency.
106//
108{
109 Initialize();
110 fMinExtent.set(0,0,0);
111 fMaxExtent.set(0,0,0);
112}
113
114
115///////////////////////////////////////////////////////////////////////////////
117{
118 DeleteObjects();
119}
120
121///////////////////////////////////////////////////////////////////////////////
122//
123// Copy constructor.
124//
126 : G4VSolid(ts)
127{
128 Initialize();
129
130 CopyObjects(ts);
131}
132
133///////////////////////////////////////////////////////////////////////////////
134//
135// Assignment operator.
136//
139{
140 if (&ts == this) { return *this; }
141
142 // Copy base class data
144
145 DeleteObjects ();
146
147 Initialize();
148
149 CopyObjects (ts);
150
151 return *this;
152}
153
154///////////////////////////////////////////////////////////////////////////////
155//
156void G4TessellatedSolid::Initialize()
157{
159
160 fRebuildPolyhedron = false; fpPolyhedron = nullptr;
161 fCubicVolume = 0.; fSurfaceArea = 0.;
162
163 fGeometryType = "G4TessellatedSolid";
164 fSolidClosed = false;
165
166 fMinExtent.set(kInfinity,kInfinity,kInfinity);
167 fMaxExtent.set(-kInfinity,-kInfinity,-kInfinity);
168
169 SetRandomVectors();
170}
171
172///////////////////////////////////////////////////////////////////////////////
173//
174void G4TessellatedSolid::DeleteObjects()
175{
176 std::size_t size = fFacets.size();
177 for (std::size_t i = 0; i < size; ++i) { delete fFacets[i]; }
178 fFacets.clear();
179 delete fpPolyhedron; fpPolyhedron = nullptr;
180}
181
182///////////////////////////////////////////////////////////////////////////////
183//
184void G4TessellatedSolid::CopyObjects (const G4TessellatedSolid &ts)
185{
186 G4ThreeVector reductionRatio;
187 G4int fmaxVoxels = fVoxels.GetMaxVoxels(reductionRatio);
188 if (fmaxVoxels < 0)
189 {
190 fVoxels.SetMaxVoxels(reductionRatio);
191 }
192 else
193 {
194 fVoxels.SetMaxVoxels(fmaxVoxels);
195 }
196
198 for (G4int i = 0; i < n; ++i)
199 {
200 G4VFacet *facetClone = (ts.GetFacet(i))->GetClone();
201 AddFacet(facetClone);
202 }
203 if (ts.GetSolidClosed()) { SetSolidClosed(true); }
204}
205
206///////////////////////////////////////////////////////////////////////////////
207//
208// Add a facet to the facet list.
209// Note that you can add, but you cannot delete.
210//
212{
213 // Add the facet to the vector.
214 //
215 if (fSolidClosed)
216 {
217 G4Exception("G4TessellatedSolid::AddFacet()", "GeomSolids1002",
218 JustWarning, "Attempt to add facets when solid is closed.");
219 return false;
220 }
221 if (aFacet->IsDefined())
222 {
223 set<G4VertexInfo,G4VertexComparator>::iterator begin
224 = fFacetList.begin(), end = fFacetList.end(), pos, it;
225 G4ThreeVector p = aFacet->GetCircumcentre();
226 G4VertexInfo value;
227 value.id = (G4int)fFacetList.size();
228 value.mag2 = p.x() + p.y() + p.z();
229
230 G4bool found = false;
231 if (!OutsideOfExtent(p, kCarTolerance))
232 {
233 G4double kCarTolerance3 = 3 * kCarTolerance;
234 pos = fFacetList.lower_bound(value);
235
236 it = pos;
237 while (!found && it != end) // Loop checking, 13.08.2015, G.Cosmo
238 {
239 G4int id = (*it).id;
240 G4VFacet *facet = fFacets[id];
241 G4ThreeVector q = facet->GetCircumcentre();
242 if ((found = (facet == aFacet))) { break; }
243 G4double dif = q.x() + q.y() + q.z() - value.mag2;
244 if (dif > kCarTolerance3) { break; }
245 it++;
246 }
247
248 if (fFacets.size() > 1)
249 {
250 it = pos;
251 while (!found && it != begin) // Loop checking, 13.08.2015, G.Cosmo
252 {
253 --it;
254 G4int id = (*it).id;
255 G4VFacet *facet = fFacets[id];
256 G4ThreeVector q = facet->GetCircumcentre();
257 found = (facet == aFacet);
258 if (found) { break; }
259 G4double dif = value.mag2 - (q.x() + q.y() + q.z());
260 if (dif > kCarTolerance3) { break; }
261 }
262 }
263 }
264
265 if (!found)
266 {
267 fFacets.push_back(aFacet);
268 fFacetList.insert(value);
269 }
270 return true;
271 }
272
273 G4Exception("G4TessellatedSolid::AddFacet()", "GeomSolids1002",
274 JustWarning, "Attempt to add facet not properly defined.");
275 aFacet->StreamInfo(G4cout);
276 return false;
277}
278
279///////////////////////////////////////////////////////////////////////////////
280//
281G4int G4TessellatedSolid::SetAllUsingStack(const std::vector<G4int>& voxel,
282 const std::vector<G4int>& max,
283 G4bool status, G4SurfBits& checked)
284{
285 vector<G4int> xyz = voxel;
286 stack<vector<G4int> > pos;
287 pos.push(xyz);
288 G4int filled = 0;
289
290 vector<G4int> candidates;
291
292 while (!pos.empty()) // Loop checking, 13.08.2015, G.Cosmo
293 {
294 xyz = pos.top();
295 pos.pop();
296 G4int index = fVoxels.GetVoxelsIndex(xyz);
297 if (!checked[index])
298 {
299 checked.SetBitNumber(index, true);
300
301 if (fVoxels.IsEmpty(index))
302 {
303 ++filled;
304
305 fInsides.SetBitNumber(index, status);
306
307 for (auto i = 0; i <= 2; ++i)
308 {
309 if (xyz[i] < max[i] - 1)
310 {
311 xyz[i]++;
312 pos.push(xyz);
313 xyz[i]--;
314 }
315
316 if (xyz[i] > 0)
317 {
318 xyz[i]--;
319 pos.push(xyz);
320 xyz[i]++;
321 }
322 }
323 }
324 }
325 }
326 return filled;
327}
328
329///////////////////////////////////////////////////////////////////////////////
330//
331void G4TessellatedSolid::PrecalculateInsides()
332{
333 vector<G4int> voxel(3), maxVoxels(3);
334 for (auto i = 0; i <= 2; ++i)
335 {
336 maxVoxels[i] = (G4int)fVoxels.GetBoundary(i).size();
337 }
338 G4int size = maxVoxels[0] * maxVoxels[1] * maxVoxels[2];
339
340 G4SurfBits checked(size-1);
341 fInsides.Clear();
342 fInsides.ResetBitNumber(size-1);
343
344 G4ThreeVector point;
345 for (voxel[2] = 0; voxel[2] < maxVoxels[2] - 1; ++voxel[2])
346 {
347 for (voxel[1] = 0; voxel[1] < maxVoxels[1] - 1; ++voxel[1])
348 {
349 for (voxel[0] = 0; voxel[0] < maxVoxels[0] - 1; ++voxel[0])
350 {
351 G4int index = fVoxels.GetVoxelsIndex(voxel);
352 if (!checked[index] && fVoxels.IsEmpty(index))
353 {
354 for (auto i = 0; i <= 2; ++i)
355 {
356 point[i] = fVoxels.GetBoundary(i)[voxel[i]];
357 }
358 auto inside = (G4bool) (InsideNoVoxels(point) == kInside);
359 SetAllUsingStack(voxel, maxVoxels, inside, checked);
360 }
361 else
362 {
363 checked.SetBitNumber(index);
364 }
365 }
366 }
367 }
368}
369
370///////////////////////////////////////////////////////////////////////////////
371//
372void G4TessellatedSolid::Voxelize ()
373{
374#ifdef G4SPECSDEBUG
375 G4cout << "Voxelizing..." << G4endl;
376#endif
377 fVoxels.Voxelize(fFacets);
378
379 if (fVoxels.Empty().GetNbits() != 0u)
380 {
381#ifdef G4SPECSDEBUG
382 G4cout << "Precalculating Insides..." << G4endl;
383#endif
384 PrecalculateInsides();
385 }
386}
387
388///////////////////////////////////////////////////////////////////////////////
389//
390// Compute extremeFacets, i.e. find those facets that have surface
391// planes that bound the volume.
392// Note that this is going to reject concaved surfaces as being extreme. Also
393// note that if the vertex is on the facet, displacement is zero, so IsInside
394// returns true. So will this work?? Need non-equality
395// "G4bool inside = displacement < 0.0;"
396// or
397// "G4bool inside = displacement <= -0.5*kCarTolerance"
398// (Notes from PT 13/08/2007).
399//
400void G4TessellatedSolid::SetExtremeFacets()
401{
402 // Copy vertices to local array
403 std::size_t vsize = fVertexList.size();
404 std::vector<G4ThreeVector> vertices(vsize);
405 for (std::size_t i = 0; i < vsize; ++i) { vertices[i] = fVertexList[i]; }
406
407 // Shuffle vertices
408 std::mt19937 gen(12345678);
409 std::shuffle(vertices.begin(), vertices.end(), gen);
410
411 // Select six extreme vertices in different directions
412 G4ThreeVector points[6];
413 for (auto & point : points) { point = vertices[0]; }
414 for (std::size_t i=1; i < vsize; ++i)
415 {
416 if (vertices[i].x() < points[0].x()) { points[0] = vertices[i]; }
417 if (vertices[i].x() > points[1].x()) { points[1] = vertices[i]; }
418 if (vertices[i].y() < points[2].y()) { points[2] = vertices[i]; }
419 if (vertices[i].y() > points[3].y()) { points[3] = vertices[i]; }
420 if (vertices[i].z() < points[4].z()) { points[4] = vertices[i]; }
421 if (vertices[i].z() > points[5].z()) { points[5] = vertices[i]; }
422 }
423
424 // Find extreme facets
425 std::size_t size = fFacets.size();
426 for (std::size_t j = 0; j < size; ++j)
427 {
428 G4VFacet &facet = *fFacets[j];
429
430 // Check extreme vertices
431 if (!facet.IsInside(points[0])) { continue; }
432 if (!facet.IsInside(points[1])) { continue; }
433 if (!facet.IsInside(points[2])) { continue; }
434 if (!facet.IsInside(points[3])) { continue; }
435 if (!facet.IsInside(points[4])) { continue; }
436 if (!facet.IsInside(points[5])) { continue; }
437
438 // Check vertices
439 G4bool isExtreme = true;
440 for (std::size_t i=0; i < vsize; ++i)
441 {
442 if (!facet.IsInside(vertices[i]))
443 {
444 isExtreme = false;
445 break;
446 }
447 }
448 if (isExtreme) { fExtremeFacets.insert(&facet); }
449 }
450}
451
452///////////////////////////////////////////////////////////////////////////////
453//
454void G4TessellatedSolid::CreateVertexList()
455{
456 // The algorithm:
457 // we will have additional vertexListSorted, where all the items will be
458 // sorted by magnitude of vertice vector.
459 // New candidate for fVertexList - we will determine the position fo first
460 // item which would be within its magnitude - 0.5*kCarTolerance.
461 // We will go trough until we will reach > +0.5 kCarTolerance.
462 // Comparison (q-p).mag() < 0.5*kCarTolerance will be made.
463 // They can be just stored in std::vector, with custom insertion based
464 // on binary search.
465
466 set<G4VertexInfo,G4VertexComparator> vertexListSorted;
467 set<G4VertexInfo,G4VertexComparator>::iterator begin
468 = vertexListSorted.begin(), end = vertexListSorted.end(), pos, it;
470 G4VertexInfo value;
471
472 fVertexList.clear();
473 std::size_t size = fFacets.size();
474
475 G4double kCarTolerance24 = kCarTolerance * kCarTolerance / 4.0;
476 G4double kCarTolerance3 = 3 * kCarTolerance;
477 vector<G4int> newIndex(100);
478
479 for (std::size_t k = 0; k < size; ++k)
480 {
481 G4VFacet &facet = *fFacets[k];
482 G4int max = facet.GetNumberOfVertices();
483
484 for (G4int i = 0; i < max; ++i)
485 {
486 p = facet.GetVertex(i);
487 value.id = (G4int)fVertexList.size();
488 value.mag2 = p.x() + p.y() + p.z();
489
490 G4bool found = false;
491 G4int id = 0;
492 if (!OutsideOfExtent(p, kCarTolerance))
493 {
494 pos = vertexListSorted.lower_bound(value);
495 it = pos;
496 while (it != end) // Loop checking, 13.08.2015, G.Cosmo
497 {
498 id = (*it).id;
499 G4ThreeVector q = fVertexList[id];
500 G4double dif = (q-p).mag2();
501 found = (dif < kCarTolerance24);
502 if (found) { break; }
503 dif = q.x() + q.y() + q.z() - value.mag2;
504 if (dif > kCarTolerance3) { break; }
505 ++it;
506 }
507
508 if (!found && (fVertexList.size() > 1))
509 {
510 it = pos;
511 while (it != begin) // Loop checking, 13.08.2015, G.Cosmo
512 {
513 --it;
514 id = (*it).id;
515 G4ThreeVector q = fVertexList[id];
516 G4double dif = (q-p).mag2();
517 found = (dif < kCarTolerance24);
518 if (found) { break; }
519 dif = value.mag2 - (q.x() + q.y() + q.z());
520 if (dif > kCarTolerance3) { break; }
521 }
522 }
523 }
524
525 if (!found)
526 {
527#ifdef G4SPECSDEBUG
528 G4cout << p.x() << ":" << p.y() << ":" << p.z() << G4endl;
529 G4cout << "Adding new vertex #" << i << " of facet " << k
530 << " id " << value.id << G4endl;
531 G4cout << "===" << G4endl;
532#endif
533 fVertexList.push_back(p);
534 vertexListSorted.insert(value);
535 begin = vertexListSorted.begin();
536 end = vertexListSorted.end();
537 newIndex[i] = value.id;
538 //
539 // Now update the maximum x, y and z limits of the volume.
540 //
541 if (value.id == 0)
542 {
543 fMinExtent = fMaxExtent = p;
544 }
545 else
546 {
547 if (p.x() > fMaxExtent.x())
548 {
549 fMaxExtent.setX(p.x());
550 }
551 else if (p.x() < fMinExtent.x())
552 {
553 fMinExtent.setX(p.x());
554 }
555 if (p.y() > fMaxExtent.y())
556 {
557 fMaxExtent.setY(p.y());
558 }
559 else if (p.y() < fMinExtent.y())
560 {
561 fMinExtent.setY(p.y());
562 }
563 if (p.z() > fMaxExtent.z())
564 {
565 fMaxExtent.setZ(p.z());
566 }
567 else if (p.z() < fMinExtent.z())
568 {
569 fMinExtent.setZ(p.z());
570 }
571 }
572 }
573 else
574 {
575#ifdef G4SPECSDEBUG
576 G4cout << p.x() << ":" << p.y() << ":" << p.z() << G4endl;
577 G4cout << "Vertex #" << i << " of facet " << k
578 << " found, redirecting to " << id << G4endl;
579 G4cout << "===" << G4endl;
580#endif
581 newIndex[i] = id;
582 }
583 }
584 // only now it is possible to change vertices pointer
585 //
586 facet.SetVertices(&fVertexList);
587 for (G4int i = 0; i < max; ++i)
588 {
589 facet.SetVertexIndex(i,newIndex[i]);
590 }
591 }
592 vector<G4ThreeVector>(fVertexList).swap(fVertexList);
593
594#ifdef G4SPECSDEBUG
595 G4double previousValue = 0.;
596 for (auto res=vertexListSorted.cbegin(); res!=vertexListSorted.cend(); ++res)
597 {
598 G4int id = (*res).id;
599 G4ThreeVector vec = fVertexList[id];
600 G4double mvalue = vec.x() + vec.y() + vec.z();
601 if (previousValue && (previousValue - 1e-9 > mvalue))
602 G4cout << "Error in CreateVertexList: previousValue " << previousValue
603 << " is smaller than mvalue " << mvalue << G4endl;
604 previousValue = mvalue;
605 }
606#endif
607}
608
609///////////////////////////////////////////////////////////////////////////////
610//
612{
614 G4int with = AllocatedMemory();
615 G4double ratio = (G4double) with / without;
616 G4cout << "G4TessellatedSolid - Allocated memory without voxel overhead "
617 << without << "; with " << with << "; ratio: " << ratio << G4endl;
618}
619
620///////////////////////////////////////////////////////////////////////////////
621//
623{
624 if (t)
625 {
626#ifdef G4SPECSDEBUG
627 G4cout << "Creating vertex list..." << G4endl;
628#endif
629 CreateVertexList();
630
631#ifdef G4SPECSDEBUG
632 G4cout << "Setting extreme facets..." << G4endl;
633#endif
634 SetExtremeFacets();
635
636#ifdef G4SPECSDEBUG
637 G4cout << "Voxelizing..." << G4endl;
638#endif
639 Voxelize();
640
641#ifdef G4SPECSDEBUG
643#endif
644
645#ifdef G4SPECSDEBUG
646 G4cout << "Checking Structure..." << G4endl;
647#endif
648 G4int irep = CheckStructure();
649 if (irep != 0)
650 {
651 if ((irep & 1) != 0)
652 {
653 std::ostringstream message;
654 message << "Defects in solid: " << GetName()
655 << " - negative cubic volume, please check orientation of facets!";
656 G4Exception("G4TessellatedSolid::SetSolidClosed()",
657 "GeomSolids1001", JustWarning, message);
658 }
659 if ((irep & 2) != 0)
660 {
661 std::ostringstream message;
662 message << "Defects in solid: " << GetName()
663 << " - some facets have wrong orientation!";
664 G4Exception("G4TessellatedSolid::SetSolidClosed()",
665 "GeomSolids1001", JustWarning, message);
666 }
667 if ((irep & 4) != 0)
668 {
669 std::ostringstream message;
670 message << "Defects in solid: " << GetName()
671 << " - there are holes in the surface!";
672 G4Exception("G4TessellatedSolid::SetSolidClosed()",
673 "GeomSolids1001", JustWarning, message);
674 }
675 }
676 }
677 fSolidClosed = t;
678}
679
680///////////////////////////////////////////////////////////////////////////////
681//
682// GetSolidClosed
683//
684// Used to determine whether the solid is closed to adding further fFacets.
685//
687{
688 return fSolidClosed;
689}
690
691///////////////////////////////////////////////////////////////////////////////
692//
693// CheckStructure
694//
695// Checks structure of the solid. Return value is a sum of the following
696// defect indicators, if any (0 means no defects):
697// 1 - cubic volume is negative, wrong orientation of facets
698// 2 - some facets have wrong orientation
699// 4 - holes in the surface
700//
702{
703 G4int nedge = 0;
704 std::size_t nface = fFacets.size();
705
706 // Calculate volume
707 //
708 G4double volume = 0.;
709 for (std::size_t i = 0; i < nface; ++i)
710 {
711 G4VFacet& facet = *fFacets[i];
712 nedge += facet.GetNumberOfVertices();
713 volume += facet.GetArea()*(facet.GetVertex(0).dot(facet.GetSurfaceNormal()));
714 }
715 auto ivolume = static_cast<G4int>(volume <= 0.);
716
717 // Create sorted vector of edges
718 //
719 std::vector<int64_t> iedge(nedge);
720 G4int kk = 0;
721 for (std::size_t i = 0; i < nface; ++i)
722 {
723 G4VFacet& facet = *fFacets[i];
724 G4int nnode = facet.GetNumberOfVertices();
725 for (G4int k = 0; k < nnode; ++k)
726 {
727 int64_t i1 = facet.GetVertexIndex((k == 0) ? nnode - 1 : k - 1);
728 int64_t i2 = facet.GetVertexIndex(k);
729 auto inverse = static_cast<int64_t>(i2 > i1);
730 if (inverse != 0) { std::swap(i1, i2); }
731 iedge[kk++] = i1*1000000000 + i2*2 + inverse;
732 }
733 }
734 std::sort(iedge.begin(), iedge.end());
735
736 // Check edges, correct structure should consist of paired edges
737 // with different orientation
738 //
739 G4int iorder = 0;
740 G4int ihole = 0;
741 G4int i = 0;
742 while (i < nedge - 1)
743 {
744 if (iedge[i + 1] - iedge[i] == 1) // paired edges with different orientation
745 {
746 i += 2;
747 }
748 else if (iedge[i + 1] == iedge[i]) // paired edges with the same orientation
749 {
750 iorder = 2;
751 i += 2;
752 }
753 else // unpaired edge
754 {
755 ihole = 4;
756 i++;
757 }
758 }
759 return ivolume + iorder + ihole;
760}
761
762///////////////////////////////////////////////////////////////////////////////
763//
764// operator+=
765//
766// This operator allows the user to add two tessellated solids together, so
767// that the solid on the left then includes all of the facets in the solid
768// on the right. Note that copies of the facets are generated, rather than
769// using the original facet set of the solid on the right.
770//
773{
774 G4int size = right.GetNumberOfFacets();
775 for (G4int i = 0; i < size; ++i) {
776 AddFacet(right.GetFacet(i)->GetClone());
777}
778
779 return *this;
780}
781
782///////////////////////////////////////////////////////////////////////////////
783//
784// GetNumberOfFacets
785//
787{
788 return (G4int)fFacets.size();
789}
790
791///////////////////////////////////////////////////////////////////////////////
792//
793EInside G4TessellatedSolid::InsideVoxels(const G4ThreeVector& p) const
794{
795 //
796 // First the simple test - check if we're outside of the X-Y-Z extremes
797 // of the tessellated solid.
798 //
799 if (OutsideOfExtent(p, kCarTolerance)) { return kOutside; }
800
801 vector<G4int> startingVoxel(3);
802 fVoxels.GetVoxel(startingVoxel, p);
803
804 const G4double dirTolerance = 1.0E-14;
805
806 const vector<G4int> &startingCandidates =
807 fVoxels.GetCandidates(startingVoxel);
808 std::size_t limit = startingCandidates.size();
809 if (limit == 0 && (fInsides.GetNbits() != 0u))
810 {
811 G4int index = fVoxels.GetPointIndex(p);
812 EInside location = fInsides[index] ? kInside : kOutside;
813 return location;
814 }
815
816 G4double minDist = kInfinity;
817
818 for(std::size_t i = 0; i < limit; ++i)
819 {
820 G4int candidate = startingCandidates[i];
821 G4VFacet &facet = *fFacets[candidate];
822 G4double dist = facet.Distance(p,minDist);
823 if (dist < minDist) { minDist = dist; }
824 if (dist <= kCarToleranceHalf) { return kSurface; }
825 }
826
827 // The following is something of an adaptation of the method implemented by
828 // Rickard Holmberg augmented with information from Schneider & Eberly,
829 // "Geometric Tools for Computer Graphics," pp700-701, 2003. In essence,
830 // we're trying to determine whether we're inside the volume by projecting
831 // a few rays and determining if the first surface crossed is has a normal
832 // vector between 0 to pi/2 (out-going) or pi/2 to pi (in-going).
833 // We should also avoid rays which are nearly within the plane of the
834 // tessellated surface, and therefore produce rays randomly.
835 // For the moment, this is a bit over-engineered (belt-braces-and-ducttape).
836 //
837 G4double distOut = kInfinity;
838 G4double distIn = kInfinity;
839 G4double distO = 0.0;
840 G4double distI = 0.0;
841 G4double distFromSurfaceO = 0.0;
842 G4double distFromSurfaceI = 0.0;
843 G4ThreeVector normalO, normalI;
844 G4bool crossingO = false;
845 G4bool crossingI = false;
846 EInside location = kOutside;
847 G4int sm = 0;
848
849 G4bool nearParallel = false;
850 do // Loop checking, 13.08.2015, G.Cosmo
851 {
852 // We loop until we find direction where the vector is not nearly parallel
853 // to the surface of any facet since this causes ambiguities. The usual
854 // case is that the angles should be sufficiently different, but there
855 // are 20 random directions to select from - hopefully sufficient.
856 //
857 distOut = distIn = kInfinity;
858 const G4ThreeVector& v = fRandir[sm];
859 ++sm;
860 //
861 // This code could be voxelized by the same algorithm, which is used for
862 // DistanceToOut(). We will traverse through fVoxels. we will call
863 // intersect only for those, which would be candidates and was not
864 // checked before.
865 //
866 G4ThreeVector currentPoint = p;
867 G4ThreeVector direction = v.unit();
868 // G4SurfBits exclusion(fVoxels.GetBitsPerSlice());
869 vector<G4int> curVoxel(3);
870 curVoxel = startingVoxel;
871 G4double shiftBonus = kCarTolerance;
872
873 G4bool crossed = false;
874 G4bool started = true;
875
876 do // Loop checking, 13.08.2015, G.Cosmo
877 {
878 const vector<G4int> &candidates =
879 started ? startingCandidates : fVoxels.GetCandidates(curVoxel);
880 started = false;
881 if (auto candidatesCount = (G4int)candidates.size())
882 {
883 for (G4int i = 0 ; i < candidatesCount; ++i)
884 {
885 G4int candidate = candidates[i];
886 // bits.SetBitNumber(candidate);
887 G4VFacet& facet = *fFacets[candidate];
888
889 crossingO = facet.Intersect(p,v,true,distO,distFromSurfaceO,normalO);
890 crossingI = facet.Intersect(p,v,false,distI,distFromSurfaceI,normalI);
891
892 if (crossingO || crossingI)
893 {
894 crossed = true;
895
896 nearParallel = (crossingO
897 && std::fabs(normalO.dot(v))<dirTolerance)
898 || (crossingI && std::fabs(normalI.dot(v))<dirTolerance);
899 if (!nearParallel)
900 {
901 if (crossingO && distO > 0.0 && distO < distOut)
902 {
903 distOut = distO;
904 }
905 if (crossingI && distI > 0.0 && distI < distIn)
906 {
907 distIn = distI;
908 }
909 }
910 else { break; }
911 }
912 }
913 if (nearParallel) { break; }
914 }
915 else
916 {
917 if (!crossed)
918 {
919 G4int index = fVoxels.GetVoxelsIndex(curVoxel);
920 G4bool inside = fInsides[index];
921 location = inside ? kInside : kOutside;
922 return location;
923 }
924 }
925
926 G4double shift=fVoxels.DistanceToNext(currentPoint, direction, curVoxel);
927 if (shift == kInfinity) { break; }
928
929 currentPoint += direction * (shift + shiftBonus);
930 }
931 while (fVoxels.UpdateCurrentVoxel(currentPoint, direction, curVoxel));
932
933 }
934 while (nearParallel && sm != fMaxTries);
935 //
936 // Here we loop through the facets to find out if there is an intersection
937 // between the ray and that facet. The test if performed separately whether
938 // the ray is entering the facet or exiting.
939 //
940#ifdef G4VERBOSE
941 if (sm == fMaxTries)
942 {
943 //
944 // We've run out of random vector directions. If nTries is set sufficiently
945 // low (nTries <= 0.5*maxTries) then this would indicate that there is
946 // something wrong with geometry.
947 //
948 std::ostringstream message;
949 G4long oldprc = message.precision(16);
950 message << "Cannot determine whether point is inside or outside volume!"
951 << G4endl
952 << "Solid name = " << GetName() << G4endl
953 << "Geometry Type = " << fGeometryType << G4endl
954 << "Number of facets = " << fFacets.size() << G4endl
955 << "Position:" << G4endl << G4endl
956 << "p.x() = " << p.x()/mm << " mm" << G4endl
957 << "p.y() = " << p.y()/mm << " mm" << G4endl
958 << "p.z() = " << p.z()/mm << " mm";
959 message.precision(oldprc);
960 G4Exception("G4TessellatedSolid::Inside()",
961 "GeomSolids1002", JustWarning, message);
962 }
963#endif
964
965 // In the next if-then-elseif G4String the logic is as follows:
966 // (1) You don't hit anything so cannot be inside volume, provided volume
967 // constructed correctly!
968 // (2) Distance to inside (ie. nearest facet such that you enter facet) is
969 // shorter than distance to outside (nearest facet such that you exit
970 // facet) - on condition of safety distance - therefore we're outside.
971 // (3) Distance to outside is shorter than distance to inside therefore
972 // we're inside.
973 //
974 if (distIn == kInfinity && distOut == kInfinity)
975 {
976 location = kOutside;
977 }
978 else if (distIn <= distOut - kCarToleranceHalf)
979 {
980 location = kOutside;
981 }
982 else if (distOut <= distIn - kCarToleranceHalf)
983 {
984 location = kInside;
985 }
986
987 return location;
988}
989
990///////////////////////////////////////////////////////////////////////////////
991//
992EInside G4TessellatedSolid::InsideNoVoxels (const G4ThreeVector &p) const
993{
994 //
995 // First the simple test - check if we're outside of the X-Y-Z extremes
996 // of the tessellated solid.
997 //
998 if (OutsideOfExtent(p, kCarTolerance)) { return kOutside; }
999
1000 const G4double dirTolerance = 1.0E-14;
1001
1002 G4double minDist = kInfinity;
1003 //
1004 // Check if we are close to a surface
1005 //
1006 std::size_t size = fFacets.size();
1007 for (std::size_t i = 0; i < size; ++i)
1008 {
1009 G4VFacet& facet = *fFacets[i];
1010 G4double dist = facet.Distance(p,minDist);
1011 if (dist < minDist) { minDist = dist; }
1012 if (dist <= kCarToleranceHalf)
1013 {
1014 return kSurface;
1015 }
1016 }
1017 //
1018 // The following is something of an adaptation of the method implemented by
1019 // Rickard Holmberg augmented with information from Schneider & Eberly,
1020 // "Geometric Tools for Computer Graphics," pp700-701, 2003. In essence, we're
1021 // trying to determine whether we're inside the volume by projecting a few
1022 // rays and determining if the first surface crossed is has a normal vector
1023 // between 0 to pi/2 (out-going) or pi/2 to pi (in-going). We should also
1024 // avoid rays which are nearly within the plane of the tessellated surface,
1025 // and therefore produce rays randomly. For the moment, this is a bit
1026 // over-engineered (belt-braces-and-ducttape).
1027 //
1028#if G4SPECSDEBUG
1029 G4int nTry = 7;
1030#else
1031 G4int nTry = 3;
1032#endif
1033 G4double distOut = kInfinity;
1034 G4double distIn = kInfinity;
1035 G4double distO = 0.0;
1036 G4double distI = 0.0;
1037 G4double distFromSurfaceO = 0.0;
1038 G4double distFromSurfaceI = 0.0;
1039 G4ThreeVector normalO(0.0,0.0,0.0);
1040 G4ThreeVector normalI(0.0,0.0,0.0);
1041 G4bool crossingO = false;
1042 G4bool crossingI = false;
1043 EInside location = kOutside;
1044 EInside locationprime = kOutside;
1045 G4int sm = 0;
1046
1047 for (G4int i=0; i<nTry; ++i)
1048 {
1049 G4bool nearParallel = false;
1050 do // Loop checking, 13.08.2015, G.Cosmo
1051 {
1052 //
1053 // We loop until we find direction where the vector is not nearly parallel
1054 // to the surface of any facet since this causes ambiguities. The usual
1055 // case is that the angles should be sufficiently different, but there
1056 // are 20 random directions to select from - hopefully sufficient.
1057 //
1058 distOut = distIn = kInfinity;
1059 G4ThreeVector v = fRandir[sm];
1060 sm++;
1061 auto f = fFacets.cbegin();
1062
1063 do // Loop checking, 13.08.2015, G.Cosmo
1064 {
1065 //
1066 // Here we loop through the facets to find out if there is an
1067 // intersection between the ray and that facet. The test if performed
1068 // separately whether the ray is entering the facet or exiting.
1069 //
1070 crossingO = ((*f)->Intersect(p,v,true,distO,distFromSurfaceO,normalO));
1071 crossingI = ((*f)->Intersect(p,v,false,distI,distFromSurfaceI,normalI));
1072 if (crossingO || crossingI)
1073 {
1074 nearParallel = (crossingO && std::fabs(normalO.dot(v))<dirTolerance)
1075 || (crossingI && std::fabs(normalI.dot(v))<dirTolerance);
1076 if (!nearParallel)
1077 {
1078 if (crossingO && distO > 0.0 && distO < distOut) { distOut = distO; }
1079 if (crossingI && distI > 0.0 && distI < distIn) { distIn = distI; }
1080 }
1081 }
1082 } while (!nearParallel && ++f != fFacets.cend());
1083 } while (nearParallel && sm != fMaxTries);
1084
1085#ifdef G4VERBOSE
1086 if (sm == fMaxTries)
1087 {
1088 //
1089 // We've run out of random vector directions. If nTries is set
1090 // sufficiently low (nTries <= 0.5*maxTries) then this would indicate
1091 // that there is something wrong with geometry.
1092 //
1093 std::ostringstream message;
1094 G4long oldprc = message.precision(16);
1095 message << "Cannot determine whether point is inside or outside volume!"
1096 << G4endl
1097 << "Solid name = " << GetName() << G4endl
1098 << "Geometry Type = " << fGeometryType << G4endl
1099 << "Number of facets = " << fFacets.size() << G4endl
1100 << "Position:" << G4endl << G4endl
1101 << "p.x() = " << p.x()/mm << " mm" << G4endl
1102 << "p.y() = " << p.y()/mm << " mm" << G4endl
1103 << "p.z() = " << p.z()/mm << " mm";
1104 message.precision(oldprc);
1105 G4Exception("G4TessellatedSolid::Inside()",
1106 "GeomSolids1002", JustWarning, message);
1107 }
1108#endif
1109 //
1110 // In the next if-then-elseif G4String the logic is as follows:
1111 // (1) You don't hit anything so cannot be inside volume, provided volume
1112 // constructed correctly!
1113 // (2) Distance to inside (ie. nearest facet such that you enter facet) is
1114 // shorter than distance to outside (nearest facet such that you exit
1115 // facet) - on condition of safety distance - therefore we're outside.
1116 // (3) Distance to outside is shorter than distance to inside therefore
1117 // we're inside.
1118 //
1119 if (distIn == kInfinity && distOut == kInfinity)
1120 {
1121 locationprime = kOutside;
1122 }
1123 else if (distIn <= distOut - kCarToleranceHalf)
1124 {
1125 locationprime = kOutside;
1126 }
1127 else if (distOut <= distIn - kCarToleranceHalf)
1128 {
1129 locationprime = kInside;
1130 }
1131
1132 if (i == 0) { location = locationprime; }
1133 }
1134
1135 return location;
1136}
1137
1138///////////////////////////////////////////////////////////////////////////////
1139//
1140// Return index of the facet closest to the point p, normally the point should
1141// be located on the surface. Return -1 if no facet selected.
1142//
1144{
1145 G4int index = -1;
1146
1147 if (fVoxels.GetCountOfVoxels() > 1)
1148 {
1149 vector<G4int> curVoxel(3);
1150 fVoxels.GetVoxel(curVoxel, p);
1151 const vector<G4int> &candidates = fVoxels.GetCandidates(curVoxel);
1152 if (auto limit = (G4int)candidates.size())
1153 {
1154 G4double minDist = kInfinity;
1155 for(G4int i = 0 ; i < limit ; ++i)
1156 {
1157 G4int candidate = candidates[i];
1158 G4VFacet& facet = *fFacets[candidate];
1159 G4double dist = facet.Distance(p, minDist);
1160 if (dist <= kCarToleranceHalf) { return index = candidate; }
1161 if (dist < minDist)
1162 {
1163 minDist = dist;
1164 index = candidate;
1165 }
1166 }
1167 }
1168 }
1169 else
1170 {
1171 G4double minDist = kInfinity;
1172 std::size_t size = fFacets.size();
1173 for (std::size_t i = 0; i < size; ++i)
1174 {
1175 G4VFacet& facet = *fFacets[i];
1176 G4double dist = facet.Distance(p, minDist);
1177 if (dist < minDist)
1178 {
1179 minDist = dist;
1180 index = (G4int)i;
1181 }
1182 }
1183 }
1184 return index;
1185}
1186
1187///////////////////////////////////////////////////////////////////////////////
1188//
1189// Return the outwards pointing unit normal of the shape for the
1190// surface closest to the point at offset p.
1191//
1193 G4ThreeVector& aNormal) const
1194{
1195 G4double minDist;
1196 G4VFacet* facet = nullptr;
1197
1198 if (fVoxels.GetCountOfVoxels() > 1)
1199 {
1200 vector<G4int> curVoxel(3);
1201 fVoxels.GetVoxel(curVoxel, p);
1202 const vector<G4int> &candidates = fVoxels.GetCandidates(curVoxel);
1203 // fVoxels.GetCandidatesVoxelArray(p, candidates, 0);
1204
1205 if (auto limit = (G4int)candidates.size())
1206 {
1207 minDist = kInfinity;
1208 for(G4int i = 0 ; i < limit ; ++i)
1209 {
1210 G4int candidate = candidates[i];
1211 G4VFacet &fct = *fFacets[candidate];
1212 G4double dist = fct.Distance(p,minDist);
1213 if (dist < minDist) { minDist = dist; }
1214 if (dist <= kCarToleranceHalf)
1215 {
1216 aNormal = fct.GetSurfaceNormal();
1217 return true;
1218 }
1219 }
1220 }
1221 minDist = MinDistanceFacet(p, true, facet);
1222 }
1223 else
1224 {
1225 minDist = kInfinity;
1226 std::size_t size = fFacets.size();
1227 for (std::size_t i = 0; i < size; ++i)
1228 {
1229 G4VFacet& f = *fFacets[i];
1230 G4double dist = f.Distance(p, minDist);
1231 if (dist < minDist)
1232 {
1233 minDist = dist;
1234 facet = &f;
1235 }
1236 }
1237 }
1238
1239 if (minDist != kInfinity)
1240 {
1241 if (facet != nullptr) { aNormal = facet->GetSurfaceNormal(); }
1242 return minDist <= kCarToleranceHalf;
1243 }
1244
1245#ifdef G4VERBOSE
1246 std::ostringstream message;
1247 message << "Point p is not on surface !?" << G4endl
1248 << " No facets found for point: " << p << " !" << G4endl
1249 << " Returning approximated value for normal.";
1250
1251 G4Exception("G4TessellatedSolid::SurfaceNormal(p)",
1252 "GeomSolids1002", JustWarning, message );
1253#endif
1254
1255 aNormal = (p.z() > 0 ? G4ThreeVector(0,0,1) : G4ThreeVector(0,0,-1));
1256 return false;
1257}
1258
1259///////////////////////////////////////////////////////////////////////////////
1260//
1261// G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v)
1262//
1263// Return the distance along the normalised vector v to the shape,
1264// from the point at offset p. If there is no intersection, return
1265// kInfinity. The first intersection resulting from 'leaving' a
1266// surface/volume is discarded. Hence, this is tolerant of points on
1267// surface of shape.
1268//
1270G4TessellatedSolid::DistanceToInNoVoxels (const G4ThreeVector& p,
1271 const G4ThreeVector& v,
1272 G4double /*aPstep*/) const
1273{
1274 G4double minDist = kInfinity;
1275 G4double dist = 0.0;
1276 G4double distFromSurface = 0.0;
1277 G4ThreeVector normal;
1278
1279#if G4SPECSDEBUG
1280 if (Inside(p) == kInside )
1281 {
1282 std::ostringstream message;
1283 G4int oldprc = message.precision(16) ;
1284 message << "Point p is already inside!?" << G4endl
1285 << "Position:" << G4endl << G4endl
1286 << " p.x() = " << p.x()/mm << " mm" << G4endl
1287 << " p.y() = " << p.y()/mm << " mm" << G4endl
1288 << " p.z() = " << p.z()/mm << " mm" << G4endl
1289 << "DistanceToOut(p) == " << DistanceToOut(p);
1290 message.precision(oldprc) ;
1291 G4Exception("G4TriangularFacet::DistanceToIn(p,v)",
1292 "GeomSolids1002", JustWarning, message);
1293 }
1294#endif
1295
1296 std::size_t size = fFacets.size();
1297 for (std::size_t i = 0; i < size; ++i)
1298 {
1299 G4VFacet& facet = *fFacets[i];
1300 if (facet.Intersect(p,v,false,dist,distFromSurface,normal))
1301 {
1302 //
1303 // set minDist to the new distance to current facet if distFromSurface
1304 // is in positive direction and point is not at surface. If the point is
1305 // within 0.5*kCarTolerance of the surface, then force distance to be
1306 // zero and leave member function immediately (for efficiency), as
1307 // proposed by & credit to Akira Okumura.
1308 //
1309 if (distFromSurface > kCarToleranceHalf && dist >= 0.0 && dist < minDist)
1310 {
1311 minDist = dist;
1312 }
1313 else
1314 {
1315 if (-kCarToleranceHalf <= dist && dist <= kCarToleranceHalf)
1316 {
1317 return 0.0;
1318 }
1319
1320 if (distFromSurface > -kCarToleranceHalf
1321 && distFromSurface < kCarToleranceHalf)
1322 {
1323 minDist = dist;
1324 }
1325 }
1326 }
1327 }
1328 return minDist;
1329}
1330
1331///////////////////////////////////////////////////////////////////////////////
1332//
1334G4TessellatedSolid::DistanceToOutNoVoxels (const G4ThreeVector& p,
1335 const G4ThreeVector& v,
1336 G4ThreeVector& aNormalVector,
1337 G4bool& aConvex,
1338 G4double /*aPstep*/) const
1339{
1340 G4double minDist = kInfinity;
1341 G4double dist = 0.0;
1342 G4double distFromSurface = 0.0;
1343 G4ThreeVector normal, minNormal;
1344
1345#if G4SPECSDEBUG
1346 if ( Inside(p) == kOutside )
1347 {
1348 std::ostringstream message;
1349 G4int oldprc = message.precision(16) ;
1350 message << "Point p is already outside!?" << G4endl
1351 << "Position:" << G4endl << G4endl
1352 << " p.x() = " << p.x()/mm << " mm" << G4endl
1353 << " p.y() = " << p.y()/mm << " mm" << G4endl
1354 << " p.z() = " << p.z()/mm << " mm" << G4endl
1355 << "DistanceToIn(p) == " << DistanceToIn(p);
1356 message.precision(oldprc) ;
1357 G4Exception("G4TriangularFacet::DistanceToOut(p)",
1358 "GeomSolids1002", JustWarning, message);
1359 }
1360#endif
1361
1362 G4bool isExtreme = false;
1363 std::size_t size = fFacets.size();
1364 for (std::size_t i = 0; i < size; ++i)
1365 {
1366 G4VFacet& facet = *fFacets[i];
1367 if (facet.Intersect(p,v,true,dist,distFromSurface,normal))
1368 {
1369 if (distFromSurface > 0.0 && distFromSurface <= kCarToleranceHalf
1371 {
1372 // We are on a surface. Return zero.
1373 aConvex = (fExtremeFacets.find(&facet) != fExtremeFacets.end());
1374 // Normal(p, aNormalVector);
1375 // aNormalVector = facet.GetSurfaceNormal();
1376 aNormalVector = normal;
1377 return 0.0;
1378 }
1379 if (dist >= 0.0 && dist < minDist)
1380 {
1381 minDist = dist;
1382 minNormal = normal;
1383 isExtreme = (fExtremeFacets.find(&facet) != fExtremeFacets.end());
1384 }
1385 }
1386 }
1387 if (minDist < kInfinity)
1388 {
1389 aNormalVector = minNormal;
1390 aConvex = isExtreme;
1391 return minDist;
1392 }
1393
1394 // No intersection found
1395 aConvex = false;
1396 Normal(p, aNormalVector);
1397 return 0.0;
1398}
1399
1400///////////////////////////////////////////////////////////////////////////////
1401//
1402void G4TessellatedSolid::
1403DistanceToOutCandidates(const std::vector<G4int>& candidates,
1404 const G4ThreeVector& aPoint,
1405 const G4ThreeVector& direction,
1406 G4double& minDist, G4ThreeVector& minNormal,
1407 G4int& minCandidate ) const
1408{
1409 auto candidatesCount = (G4int)candidates.size();
1410 G4double dist = 0.0;
1411 G4double distFromSurface = 0.0;
1412 G4ThreeVector normal;
1413
1414 for (G4int i = 0 ; i < candidatesCount; ++i)
1415 {
1416 G4int candidate = candidates[i];
1417 G4VFacet& facet = *fFacets[candidate];
1418 if (facet.Intersect(aPoint,direction,true,dist,distFromSurface,normal))
1419 {
1420 if (distFromSurface > 0.0 && distFromSurface <= kCarToleranceHalf
1421 && facet.Distance(aPoint,kCarTolerance) <= kCarToleranceHalf)
1422 {
1423 // We are on a surface
1424 //
1425 minDist = 0.0;
1426 minNormal = normal;
1427 minCandidate = candidate;
1428 break;
1429 }
1430 if (dist >= 0.0 && dist < minDist)
1431 {
1432 minDist = dist;
1433 minNormal = normal;
1434 minCandidate = candidate;
1435 }
1436 }
1437 }
1438}
1439
1440///////////////////////////////////////////////////////////////////////////////
1441//
1443G4TessellatedSolid::DistanceToOutCore(const G4ThreeVector& aPoint,
1444 const G4ThreeVector& aDirection,
1445 G4ThreeVector& aNormalVector,
1446 G4bool &aConvex,
1447 G4double aPstep) const
1448{
1449 G4double minDistance;
1450
1451 if (fVoxels.GetCountOfVoxels() > 1)
1452 {
1453 minDistance = kInfinity;
1454
1455 G4ThreeVector currentPoint = aPoint;
1456 G4ThreeVector direction = aDirection.unit();
1457 G4double totalShift = 0.;
1458 vector<G4int> curVoxel(3);
1459 if (!fVoxels.Contains(aPoint)) { return 0.; }
1460
1461 fVoxels.GetVoxel(curVoxel, currentPoint);
1462
1463 G4double shiftBonus = kCarTolerance;
1464
1465 const vector<G4int>* old = nullptr;
1466
1467 G4int minCandidate = -1;
1468 do // Loop checking, 13.08.2015, G.Cosmo
1469 {
1470 const vector<G4int>& candidates = fVoxels.GetCandidates(curVoxel);
1471 if (old == &candidates) { ++old; }
1472 if (old != &candidates && !candidates.empty())
1473 {
1474 DistanceToOutCandidates(candidates, aPoint, direction, minDistance,
1475 aNormalVector, minCandidate);
1476 if (minDistance <= totalShift) { break; }
1477 }
1478
1479 G4double shift=fVoxels.DistanceToNext(currentPoint, direction, curVoxel);
1480 if (shift == kInfinity) { break; }
1481
1482 totalShift += shift;
1483 if (minDistance <= totalShift) { break; }
1484
1485 currentPoint += direction * (shift + shiftBonus);
1486
1487 old = &candidates;
1488 }
1489 while (fVoxels.UpdateCurrentVoxel(currentPoint, direction, curVoxel));
1490
1491 if (minCandidate < 0)
1492 {
1493 // No intersection found
1494 minDistance = 0.;
1495 aConvex = false;
1496 Normal(aPoint, aNormalVector);
1497 }
1498 else
1499 {
1500 aConvex = (fExtremeFacets.find(fFacets[minCandidate])
1501 != fExtremeFacets.end());
1502 }
1503 }
1504 else
1505 {
1506 minDistance = DistanceToOutNoVoxels(aPoint, aDirection, aNormalVector,
1507 aConvex, aPstep);
1508 }
1509 return minDistance;
1510}
1511
1512///////////////////////////////////////////////////////////////////////////////
1513//
1514G4double G4TessellatedSolid::
1515DistanceToInCandidates(const std::vector<G4int>& candidates,
1516 const G4ThreeVector& aPoint,
1517 const G4ThreeVector& direction) const
1518{
1519 auto candidatesCount = (G4int)candidates.size();
1520 G4double dist = 0.0;
1521 G4double distFromSurface = 0.0;
1522 G4ThreeVector normal;
1523
1524 G4double minDistance = kInfinity;
1525 for (G4int i = 0 ; i < candidatesCount; ++i)
1526 {
1527 G4int candidate = candidates[i];
1528 G4VFacet& facet = *fFacets[candidate];
1529 if (facet.Intersect(aPoint,direction,false,dist,distFromSurface,normal))
1530 {
1531 //
1532 // Set minDist to the new distance to current facet if distFromSurface is
1533 // in positive direction and point is not at surface. If the point is
1534 // within 0.5*kCarTolerance of the surface, then force distance to be
1535 // zero and leave member function immediately (for efficiency), as
1536 // proposed by & credit to Akira Okumura.
1537 //
1538 if ( (distFromSurface > kCarToleranceHalf)
1539 && (dist >= 0.0) && (dist < minDistance))
1540 {
1541 minDistance = dist;
1542 }
1543 else
1544 {
1545 if (-kCarToleranceHalf <= dist && dist <= kCarToleranceHalf)
1546 {
1547 return 0.0;
1548 }
1549 if (distFromSurface > -kCarToleranceHalf
1550 && distFromSurface < kCarToleranceHalf)
1551 {
1552 minDistance = dist;
1553 }
1554 }
1555 }
1556 }
1557 return minDistance;
1558}
1559
1560///////////////////////////////////////////////////////////////////////////////
1561//
1563G4TessellatedSolid::DistanceToInCore(const G4ThreeVector& aPoint,
1564 const G4ThreeVector& aDirection,
1565 G4double aPstep) const
1566{
1567 G4double minDistance;
1568
1569 if (fVoxels.GetCountOfVoxels() > 1)
1570 {
1571 minDistance = kInfinity;
1572 G4ThreeVector currentPoint = aPoint;
1573 G4ThreeVector direction = aDirection.unit();
1574 G4double shift = fVoxels.DistanceToFirst(currentPoint, direction);
1575 if (shift == kInfinity) { return shift; }
1576 G4double shiftBonus = kCarTolerance;
1577 if (shift != 0.0) { currentPoint += direction * (shift + shiftBonus); }
1578 // if (!fVoxels.Contains(currentPoint)) return minDistance;
1579 G4double totalShift = shift;
1580
1581 // G4SurfBits exclusion; // (1/*fVoxels.GetBitsPerSlice()*/);
1582 vector<G4int> curVoxel(3);
1583
1584 fVoxels.GetVoxel(curVoxel, currentPoint);
1585 do // Loop checking, 13.08.2015, G.Cosmo
1586 {
1587 const vector<G4int>& candidates = fVoxels.GetCandidates(curVoxel);
1588 if (!candidates.empty())
1589 {
1590 G4double distance=DistanceToInCandidates(candidates, aPoint, direction);
1591 if (minDistance > distance) { minDistance = distance; }
1592 if (distance < totalShift) { break; }
1593 }
1594
1595 shift = fVoxels.DistanceToNext(currentPoint, direction, curVoxel);
1596 if (shift == kInfinity /*|| shift == 0*/) { break; }
1597
1598 totalShift += shift;
1599 if (minDistance < totalShift) { break; }
1600
1601 currentPoint += direction * (shift + shiftBonus);
1602 }
1603 while (fVoxels.UpdateCurrentVoxel(currentPoint, direction, curVoxel));
1604 }
1605 else
1606 {
1607 minDistance = DistanceToInNoVoxels(aPoint, aDirection, aPstep);
1608 }
1609
1610 return minDistance;
1611}
1612
1613///////////////////////////////////////////////////////////////////////////////
1614//
1615G4bool
1616G4TessellatedSolid::CompareSortedVoxel(const std::pair<G4int, G4double>& l,
1617 const std::pair<G4int, G4double>& r)
1618{
1619 return l.second < r.second;
1620}
1621
1622///////////////////////////////////////////////////////////////////////////////
1623//
1625G4TessellatedSolid::MinDistanceFacet(const G4ThreeVector& p,
1626 G4bool simple,
1627 G4VFacet* &minFacet) const
1628{
1629 G4double minDist = kInfinity;
1630
1631 G4int size = fVoxels.GetVoxelBoxesSize();
1632 vector<pair<G4int, G4double> > voxelsSorted(size);
1633
1634 pair<G4int, G4double> info;
1635
1636 for (G4int i = 0; i < size; ++i)
1637 {
1638 const G4VoxelBox& voxelBox = fVoxels.GetVoxelBox(i);
1639
1640 G4ThreeVector pointShifted = p - voxelBox.pos;
1641 G4double safety = fVoxels.MinDistanceToBox(pointShifted, voxelBox.hlen);
1642 info.first = i;
1643 info.second = safety;
1644 voxelsSorted[i] = info;
1645 }
1646
1647 std::sort(voxelsSorted.begin(), voxelsSorted.end(),
1648 &G4TessellatedSolid::CompareSortedVoxel);
1649
1650 for (G4int i = 0; i < size; ++i)
1651 {
1652 const pair<G4int,G4double>& inf = voxelsSorted[i];
1653 G4double dist = inf.second;
1654 if (dist > minDist) { break; }
1655
1656 const vector<G4int>& candidates = fVoxels.GetVoxelBoxCandidates(inf.first);
1657 auto csize = (G4int)candidates.size();
1658 for (G4int j = 0; j < csize; ++j)
1659 {
1660 G4int candidate = candidates[j];
1661 G4VFacet& facet = *fFacets[candidate];
1662 dist = simple ? facet.Distance(p,minDist)
1663 : facet.Distance(p,minDist,false);
1664 if (dist < minDist)
1665 {
1666 minDist = dist;
1667 minFacet = &facet;
1668 }
1669 }
1670 }
1671 return minDist;
1672}
1673
1674///////////////////////////////////////////////////////////////////////////////
1675//
1677 G4bool aAccurate) const
1678{
1679#if G4SPECSDEBUG
1680 if ( Inside(p) == kInside )
1681 {
1682 std::ostringstream message;
1683 G4int oldprc = message.precision(16) ;
1684 message << "Point p is already inside!?" << G4endl
1685 << "Position:" << G4endl << G4endl
1686 << "p.x() = " << p.x()/mm << " mm" << G4endl
1687 << "p.y() = " << p.y()/mm << " mm" << G4endl
1688 << "p.z() = " << p.z()/mm << " mm" << G4endl
1689 << "DistanceToOut(p) == " << DistanceToOut(p);
1690 message.precision(oldprc) ;
1691 G4Exception("G4TriangularFacet::DistanceToIn(p)",
1692 "GeomSolids1002", JustWarning, message);
1693 }
1694#endif
1695
1696 G4double minDist;
1697
1698 if (fVoxels.GetCountOfVoxels() > 1)
1699 {
1700 if (!aAccurate) { return fVoxels.DistanceToBoundingBox(p); }
1701
1702 if (!OutsideOfExtent(p, kCarTolerance))
1703 {
1704 vector<G4int> startingVoxel(3);
1705 fVoxels.GetVoxel(startingVoxel, p);
1706 const vector<G4int> &candidates = fVoxels.GetCandidates(startingVoxel);
1707 if (candidates.empty() && (fInsides.GetNbits() != 0u))
1708 {
1709 G4int index = fVoxels.GetPointIndex(p);
1710 if (fInsides[index]) { return 0.; }
1711 }
1712 }
1713
1714 G4VFacet* facet;
1715 minDist = MinDistanceFacet(p, true, facet);
1716 }
1717 else
1718 {
1719 minDist = kInfinity;
1720 std::size_t size = fFacets.size();
1721 for (std::size_t i = 0; i < size; ++i)
1722 {
1723 G4VFacet& facet = *fFacets[i];
1724 G4double dist = facet.Distance(p,minDist);
1725 if (dist < minDist) { minDist = dist; }
1726 }
1727 }
1728 return minDist;
1729}
1730
1731///////////////////////////////////////////////////////////////////////////////
1732//
1735{
1736#if G4SPECSDEBUG
1737 if ( Inside(p) == kOutside )
1738 {
1739 std::ostringstream message;
1740 G4int oldprc = message.precision(16) ;
1741 message << "Point p is already outside!?" << G4endl
1742 << "Position:" << G4endl << G4endl
1743 << "p.x() = " << p.x()/mm << " mm" << G4endl
1744 << "p.y() = " << p.y()/mm << " mm" << G4endl
1745 << "p.z() = " << p.z()/mm << " mm" << G4endl
1746 << "DistanceToIn(p) == " << DistanceToIn(p);
1747 message.precision(oldprc) ;
1748 G4Exception("G4TriangularFacet::DistanceToOut(p)",
1749 "GeomSolids1002", JustWarning, message);
1750 }
1751#endif
1752
1753 G4double minDist;
1754
1755 if (OutsideOfExtent(p, kCarTolerance)) { return 0.0; }
1756
1757 if (fVoxels.GetCountOfVoxels() > 1)
1758 {
1759 G4VFacet* facet;
1760 minDist = MinDistanceFacet(p, true, facet);
1761 }
1762 else
1763 {
1764 minDist = kInfinity;
1765 G4double dist = 0.0;
1766 std::size_t size = fFacets.size();
1767 for (std::size_t i = 0; i < size; ++i)
1768 {
1769 G4VFacet& facet = *fFacets[i];
1770 dist = facet.Distance(p,minDist);
1771 if (dist < minDist) { minDist = dist; }
1772 }
1773 }
1774 return minDist;
1775}
1776
1777///////////////////////////////////////////////////////////////////////////////
1778//
1779// G4GeometryType GetEntityType() const;
1780//
1781// Provide identification of the class of an object
1782//
1784{
1785 return fGeometryType;
1786}
1787
1788///////////////////////////////////////////////////////////////////////////////
1789//
1790// IsFaceted
1791//
1793{
1794 return true;
1795}
1796
1797///////////////////////////////////////////////////////////////////////////////
1798//
1799std::ostream &G4TessellatedSolid::StreamInfo(std::ostream &os) const
1800{
1801 os << G4endl;
1802 os << "Solid name = " << GetName() << G4endl;
1803 os << "Geometry Type = " << fGeometryType << G4endl;
1804 os << "Number of facets = " << fFacets.size() << G4endl;
1805
1806 std::size_t size = fFacets.size();
1807 for (std::size_t i = 0; i < size; ++i)
1808 {
1809 os << "FACET # = " << i + 1 << G4endl;
1810 G4VFacet &facet = *fFacets[i];
1811 facet.StreamInfo(os);
1812 }
1813 os << G4endl;
1814
1815 return os;
1816}
1817
1818///////////////////////////////////////////////////////////////////////////////
1819//
1820// Make a clone of the object
1821//
1823{
1824 return new G4TessellatedSolid(*this);
1825}
1826
1827///////////////////////////////////////////////////////////////////////////////
1828//
1829// EInside G4TessellatedSolid::Inside (const G4ThreeVector &p) const
1830//
1831// This method must return:
1832// * kOutside if the point at offset p is outside the shape
1833// boundaries plus kCarTolerance/2,
1834// * kSurface if the point is <= kCarTolerance/2 from a surface, or
1835// * kInside otherwise.
1836//
1838{
1839 EInside location;
1840
1841 if (fVoxels.GetCountOfVoxels() > 1)
1842 {
1843 location = InsideVoxels(aPoint);
1844 }
1845 else
1846 {
1847 location = InsideNoVoxels(aPoint);
1848 }
1849 return location;
1850}
1851
1852///////////////////////////////////////////////////////////////////////////////
1853//
1855{
1856 G4ThreeVector n;
1857 Normal(p, n);
1858 return n;
1859}
1860
1861///////////////////////////////////////////////////////////////////////////////
1862//
1863// G4double DistanceToIn(const G4ThreeVector& p)
1864//
1865// Calculate distance to nearest surface of shape from an outside point p. The
1866// distance can be an underestimate.
1867//
1869{
1870 return SafetyFromOutside(p, false);
1871}
1872
1873///////////////////////////////////////////////////////////////////////////////
1874//
1876 const G4ThreeVector& v)const
1877{
1878 G4double dist = DistanceToInCore(p,v,kInfinity);
1879#ifdef G4SPECSDEBUG
1880 if (dist < kInfinity)
1881 {
1882 if (Inside(p + dist*v) != kSurface)
1883 {
1884 std::ostringstream message;
1885 message << "Invalid response from facet in solid '" << GetName() << "',"
1886 << G4endl
1887 << "at point: " << p << "and direction: " << v;
1888 G4Exception("G4TessellatedSolid::DistanceToIn(p,v)",
1889 "GeomSolids1002", JustWarning, message);
1890 }
1891 }
1892#endif
1893 return dist;
1894}
1895
1896///////////////////////////////////////////////////////////////////////////////
1897//
1898// G4double DistanceToOut(const G4ThreeVector& p)
1899//
1900// Calculate distance to nearest surface of shape from an inside
1901// point. The distance can be an underestimate.
1902//
1904{
1905 return SafetyFromInside(p, false);
1906}
1907
1908///////////////////////////////////////////////////////////////////////////////
1909//
1910// G4double DistanceToOut(const G4ThreeVector& p, const G4ThreeVector& v,
1911// const G4bool calcNorm=false,
1912// G4bool *validNorm=0, G4ThreeVector *n=0);
1913//
1914// Return distance along the normalised vector v to the shape, from a
1915// point at an offset p inside or on the surface of the
1916// shape. Intersections with surfaces, when the point is not greater
1917// than kCarTolerance/2 from a surface, must be ignored.
1918// If calcNorm is true, then it must also set validNorm to either
1919// * true, if the solid lies entirely behind or on the exiting
1920// surface. Then it must set n to the outwards normal vector
1921// (the Magnitude of the vector is not defined).
1922// * false, if the solid does not lie entirely behind or on the
1923// exiting surface.
1924// If calcNorm is false, then validNorm and n are unused.
1925//
1927 const G4ThreeVector& v,
1928 const G4bool calcNorm,
1929 G4bool* validNorm,
1930 G4ThreeVector* norm) const
1931{
1932 G4ThreeVector n;
1933 G4bool valid;
1934
1935 G4double dist = DistanceToOutCore(p, v, n, valid);
1936 if (calcNorm)
1937 {
1938 *norm = n;
1939 *validNorm = valid;
1940 }
1941#ifdef G4SPECSDEBUG
1942 if (dist < kInfinity)
1943 {
1944 if (Inside(p + dist*v) != kSurface)
1945 {
1946 std::ostringstream message;
1947 message << "Invalid response from facet in solid '" << GetName() << "',"
1948 << G4endl
1949 << "at point: " << p << "and direction: " << v;
1950 G4Exception("G4TessellatedSolid::DistanceToOut(p,v,..)",
1951 "GeomSolids1002", JustWarning, message);
1952 }
1953 }
1954#endif
1955 return dist;
1956}
1957
1958///////////////////////////////////////////////////////////////////////////////
1959//
1961{
1962 scene.AddSolid (*this);
1963}
1964
1965///////////////////////////////////////////////////////////////////////////////
1966//
1968{
1969 auto nVertices = (G4int)fVertexList.size();
1970 auto nFacets = (G4int)fFacets.size();
1971 auto polyhedron = new G4Polyhedron(nVertices, nFacets);
1972 for (auto i = 0; i < nVertices; ++i)
1973 {
1974 polyhedron->SetVertex(i+1, fVertexList[i]);
1975 }
1976
1977 for (auto i = 0; i < nFacets; ++i)
1978 {
1979 G4VFacet* facet = fFacets[i];
1980 G4int v[4] = {0};
1981 G4int n = facet->GetNumberOfVertices();
1982 if (n > 4) { n = 4; }
1983 for (auto j = 0; j < n; ++j)
1984 {
1985 v[j] = facet->GetVertexIndex(j) + 1;
1986 }
1987 polyhedron->SetFacet(i+1, v[0], v[1], v[2], v[3]);
1988 }
1989 polyhedron->SetReferences();
1990
1991 return polyhedron;
1992}
1993
1994///////////////////////////////////////////////////////////////////////////////
1995//
1996// GetPolyhedron
1997//
1999{
2000 if (fpPolyhedron == nullptr ||
2001 fRebuildPolyhedron ||
2002 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
2003 fpPolyhedron->GetNumberOfRotationSteps())
2004 {
2005 G4AutoLock l(&polyhedronMutex);
2006 delete fpPolyhedron;
2007 fpPolyhedron = CreatePolyhedron();
2008 fRebuildPolyhedron = false;
2009 l.unlock();
2010 }
2011 return fpPolyhedron;
2012}
2013
2014///////////////////////////////////////////////////////////////////////////////
2015//
2016// Get bounding box
2017//
2019 G4ThreeVector& pMax) const
2020{
2021 pMin = fMinExtent;
2022 pMax = fMaxExtent;
2023
2024 // Check correctness of the bounding box
2025 //
2026 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
2027 {
2028 std::ostringstream message;
2029 message << "Bad bounding box (min >= max) for solid: "
2030 << GetName() << " !"
2031 << "\npMin = " << pMin
2032 << "\npMax = " << pMax;
2033 G4Exception("G4TessellatedSolid::BoundingLimits()",
2034 "GeomMgt0001", JustWarning, message);
2035 DumpInfo();
2036 }
2037}
2038
2039///////////////////////////////////////////////////////////////////////////////
2040//
2041// Calculate extent under transform and specified limit
2042//
2043G4bool
2045 const G4VoxelLimits& pVoxelLimit,
2046 const G4AffineTransform& pTransform,
2047 G4double& pMin, G4double& pMax) const
2048{
2049 G4ThreeVector bmin, bmax;
2050
2051 // Check bounding box (bbox)
2052 //
2053 BoundingLimits(bmin,bmax);
2054 G4BoundingEnvelope bbox(bmin,bmax);
2055
2056 // Use simple bounding-box to help in the case of complex meshes
2057 //
2058 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
2059
2060#if 0
2061 // Precise extent computation (disabled by default for this shape)
2062 //
2063 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
2064 {
2065 return (pMin < pMax) ? true : false;
2066 }
2067
2068 // The extent is calculated as cumulative extent of the pyramids
2069 // formed by facets and the center of the bounding box.
2070 //
2071 G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
2072 G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
2073
2074 G4ThreeVectorList base;
2075 G4ThreeVectorList apex(1);
2076 std::vector<const G4ThreeVectorList *> pyramid(2);
2077 pyramid[0] = &base;
2078 pyramid[1] = &apex;
2079 apex[0] = (bmin+bmax)*0.5;
2080
2081 // main loop along facets
2082 pMin = kInfinity;
2083 pMax = -kInfinity;
2084 for (G4int i=0; i<GetNumberOfFacets(); ++i)
2085 {
2086 G4VFacet* facet = GetFacet(i);
2087 if (std::abs((facet->GetSurfaceNormal()).dot(facet->GetVertex(0)-apex[0]))
2088 < kCarToleranceHalf) { continue; }
2089
2090 G4int nv = facet->GetNumberOfVertices();
2091 base.resize(nv);
2092 for (G4int k=0; k<nv; ++k) { base[k] = facet->GetVertex(k); }
2093
2094 G4double emin,emax;
2095 G4BoundingEnvelope benv(pyramid);
2096 if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax))
2097 {
2098 continue;
2099 }
2100 if (emin < pMin) { pMin = emin; }
2101 if (emax > pMax) { pMax = emax; }
2102 if (eminlim > pMin && emaxlim < pMax) { break; } // max possible extent
2103 }
2104 return (pMin < pMax);
2105#endif
2106}
2107
2108///////////////////////////////////////////////////////////////////////////////
2109//
2111{
2112 return fMinExtent.x();
2113}
2114
2115///////////////////////////////////////////////////////////////////////////////
2116//
2118{
2119 return fMaxExtent.x();
2120}
2121
2122///////////////////////////////////////////////////////////////////////////////
2123//
2125{
2126 return fMinExtent.y();
2127}
2128
2129///////////////////////////////////////////////////////////////////////////////
2130//
2132{
2133 return fMaxExtent.y();
2134}
2135
2136///////////////////////////////////////////////////////////////////////////////
2137//
2139{
2140 return fMinExtent.z();
2141}
2142
2143///////////////////////////////////////////////////////////////////////////////
2144//
2146{
2147 return fMaxExtent.z();
2148}
2149
2150///////////////////////////////////////////////////////////////////////////////
2151//
2153{
2154 return { fMinExtent.x(), fMaxExtent.x(),
2155 fMinExtent.y(), fMaxExtent.y(),
2156 fMinExtent.z(), fMaxExtent.z() };
2157}
2158
2159///////////////////////////////////////////////////////////////////////////////
2160//
2162{
2163 if (fCubicVolume != 0.) { return fCubicVolume; }
2164
2165 // For explanation of the following algorithm see:
2166 // https://en.wikipedia.org/wiki/Polyhedron#Volume
2167 // http://wwwf.imperial.ac.uk/~rn/centroid.pdf
2168
2169 G4AutoLock l(&tessMutex);
2170 std::size_t size = fFacets.size();
2171 for (std::size_t i = 0; i < size; ++i)
2172 {
2173 G4VFacet &facet = *fFacets[i];
2174 G4double area = facet.GetArea();
2175 G4ThreeVector unit_normal = facet.GetSurfaceNormal();
2176 fCubicVolume += area * (facet.GetVertex(0).dot(unit_normal));
2177 }
2178 fCubicVolume /= 3.;
2179 l.unlock();
2180
2181 return fCubicVolume;
2182}
2183
2184///////////////////////////////////////////////////////////////////////////////
2185//
2187{
2188 if (fSurfaceArea != 0.) { return fSurfaceArea; }
2189
2190 G4AutoLock l(&tessMutex);
2191 std::size_t size = fFacets.size();
2192 for (std::size_t i = 0; i < size; ++i)
2193 {
2194 G4VFacet &facet = *fFacets[i];
2195 fSurfaceArea += facet.GetArea();
2196 }
2197 l.unlock();
2198
2199 return fSurfaceArea;
2200}
2201
2202///////////////////////////////////////////////////////////////////////////////
2203//
2205{
2206 // Select randomly a facet and return a random point on it
2207
2208 auto i = (G4int)(fFacets.size()*G4QuickRand());
2209 return fFacets[i]->GetPointOnFace();
2210}
2211
2212///////////////////////////////////////////////////////////////////////////////
2213//
2214// SetRandomVectorSet
2215//
2216// This is a set of predefined random vectors (if that isn't a contradition
2217// in terms!) used to generate rays from a user-defined point. The member
2218// function Inside uses these to determine whether the point is inside or
2219// outside of the tessellated solid. All vectors should be unit vectors.
2220//
2221void G4TessellatedSolid::SetRandomVectors ()
2222{
2223 fRandir.resize(20);
2224 fRandir[0] =
2225 G4ThreeVector(-0.9577428892113370, 0.2732676269591740, 0.0897405271949221);
2226 fRandir[1] =
2227 G4ThreeVector(-0.8331264504940770,-0.5162067214954600,-0.1985722492445700);
2228 fRandir[2] =
2229 G4ThreeVector(-0.1516671651108820, 0.9666292616127460, 0.2064580868390110);
2230 fRandir[3] =
2231 G4ThreeVector( 0.6570250350323190,-0.6944539025883300, 0.2933460081893360);
2232 fRandir[4] =
2233 G4ThreeVector(-0.4820456281280320,-0.6331060000098690,-0.6056474264406270);
2234 fRandir[5] =
2235 G4ThreeVector( 0.7629032554236800 , 0.1016854697539910,-0.6384658864065180);
2236 fRandir[6] =
2237 G4ThreeVector( 0.7689540409061150, 0.5034929891988220, 0.3939600142169160);
2238 fRandir[7] =
2239 G4ThreeVector( 0.5765188359255740, 0.5997271636278330,-0.5549354566343150);
2240 fRandir[8] =
2241 G4ThreeVector( 0.6660632777862070,-0.6362809868288380, 0.3892379937580790);
2242 fRandir[9] =
2243 G4ThreeVector( 0.3824415020414780, 0.6541792713761380,-0.6525243125110690);
2244 fRandir[10] =
2245 G4ThreeVector(-0.5107726564526760, 0.6020905056811610, 0.6136760679616570);
2246 fRandir[11] =
2247 G4ThreeVector( 0.7459135439578050, 0.6618796061649330, 0.0743530220183488);
2248 fRandir[12] =
2249 G4ThreeVector( 0.1536405855311580, 0.8117477913978260,-0.5634359711967240);
2250 fRandir[13] =
2251 G4ThreeVector( 0.0744395301705579,-0.8707110101772920,-0.4861286795736560);
2252 fRandir[14] =
2253 G4ThreeVector(-0.1665874645185400, 0.6018553940549240,-0.7810369397872780);
2254 fRandir[15] =
2255 G4ThreeVector( 0.7766902003633100, 0.6014617505959970,-0.1870724331097450);
2256 fRandir[16] =
2257 G4ThreeVector(-0.8710128685847430,-0.1434320216603030,-0.4698551243971010);
2258 fRandir[17] =
2259 G4ThreeVector( 0.8901082092766820,-0.4388411398893870, 0.1229871120030100);
2260 fRandir[18] =
2261 G4ThreeVector(-0.6430417431544370,-0.3295938228697690, 0.6912779675984150);
2262 fRandir[19] =
2263 G4ThreeVector( 0.6331124368380410, 0.6306211461665000, 0.4488714875425340);
2264
2265 fMaxTries = 20;
2266}
2267
2268///////////////////////////////////////////////////////////////////////////////
2269//
2271{
2272 G4int base = sizeof(*this);
2273 base += fVertexList.capacity() * sizeof(G4ThreeVector);
2274 base += fRandir.capacity() * sizeof(G4ThreeVector);
2275
2276 std::size_t limit = fFacets.size();
2277 for (std::size_t i = 0; i < limit; ++i)
2278 {
2279 G4VFacet& facet = *fFacets[i];
2280 base += facet.AllocatedMemory();
2281 }
2282
2283 for (const auto & fExtremeFacet : fExtremeFacets)
2284 {
2285 G4VFacet &facet = *fExtremeFacet;
2286 base += facet.AllocatedMemory();
2287 }
2288 return base;
2289}
2290
2291///////////////////////////////////////////////////////////////////////////////
2292//
2294{
2296 G4int sizeInsides = fInsides.GetNbytes();
2297 G4int sizeVoxels = fVoxels.AllocatedMemory();
2298 size += sizeInsides + sizeVoxels;
2299 return size;
2300}
2301
2302#endif
G4TemplateAutoLock< G4Mutex > G4AutoLock
const G4double kCarTolerance
std::vector< G4ThreeVector > G4ThreeVectorList
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double G4QuickRand(uint32_t seed=0)
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
long G4long
Definition G4Types.hh:87
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4String G4GeometryType
Definition G4VSolid.hh:70
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double z() const
Hep3Vector unit() const
double x() const
double y() const
double dot(const Hep3Vector &) const
void set(double x, double y, double z)
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4BoundingEnvelope is a helper class to facilitate calculation of the extent of a solid within the li...
G4bool BoundingBoxVsVoxelLimits(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4SurfBits provides a simple container of bits, to be used for optimization of tessellated surfaces....
Definition G4SurfBits.hh:57
void SetBitNumber(unsigned int bitnumber, G4bool value=true)
G4TessellatedSolid is a solid defined by a number of facets. It is important that the supplied facets...
virtual G4bool Normal(const G4ThreeVector &p, G4ThreeVector &n) const
G4double GetMinYExtent() const
G4double GetMinZExtent() const
EInside Inside(const G4ThreeVector &p) const override
G4TessellatedSolid & operator=(const G4TessellatedSolid &right)
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
G4Polyhedron * GetPolyhedron() const override
G4TessellatedSolid & operator+=(const G4TessellatedSolid &right)
G4bool AddFacet(G4VFacet *aFacet)
G4int GetNumberOfFacets() const
G4double DistanceToOut(const G4ThreeVector &p) const override
G4ThreeVector GetPointOnSurface() const override
G4double GetMaxYExtent() const
G4double GetMaxZExtent() const
G4double GetMaxXExtent() const
G4bool GetSolidClosed() const
G4VFacet * GetFacet(G4int i) const
virtual G4double SafetyFromInside(const G4ThreeVector &p, G4bool aAccurate=false) const
G4double GetMinXExtent() const
void SetSolidClosed(const G4bool t)
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const override
G4VisExtent GetExtent() const override
G4bool IsFaceted() const override
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
G4VSolid * Clone() const override
G4double GetSurfaceArea() override
G4int GetFacetIndex(const G4ThreeVector &p) const
virtual G4double SafetyFromOutside(const G4ThreeVector &p, G4bool aAccurate=false) const
G4double GetCubicVolume() override
std::ostream & StreamInfo(std::ostream &os) const override
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
G4Polyhedron * CreatePolyhedron() const override
G4GeometryType GetEntityType() const override
G4VFacet is a base class defining the facets which are components of a G4TessellatedSolid shape.
Definition G4VFacet.hh:56
virtual void SetVertexIndex(G4int i, G4int j)=0
virtual G4double Distance(const G4ThreeVector &, G4double minDist)=0
virtual G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal)=0
virtual G4int AllocatedMemory()=0
virtual G4ThreeVector GetCircumcentre() const =0
std::ostream & StreamInfo(std::ostream &os) const
Definition G4VFacet.cc:96
G4bool IsInside(const G4ThreeVector &p) const
Definition G4VFacet.cc:114
virtual G4ThreeVector GetSurfaceNormal() const =0
virtual G4ThreeVector GetVertex(G4int i) const =0
virtual G4double GetArea() const =0
virtual G4int GetNumberOfVertices() const =0
virtual G4int GetVertexIndex(G4int i) const =0
virtual G4VFacet * GetClone()=0
virtual void SetVertices(std::vector< G4ThreeVector > *vertices)=0
virtual G4bool IsDefined() const =0
virtual void AddSolid(const G4Box &)=0
G4String GetName() const
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59
void DumpInfo() const
G4double kCarTolerance
Definition G4VSolid.hh:418
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:108
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
G4double GetMinExtent(const EAxis pAxis) const
G4double GetMaxExtent(const EAxis pAxis) const
G4bool IsEmpty(G4int index) const
G4int GetVoxelsIndex(G4int x, G4int y, G4int z) const
EAxis
Definition geomdefs.hh:54
EInside
Definition geomdefs.hh:67
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68
@ kSurface
Definition geomdefs.hh:69
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4ThreeVector hlen
G4ThreeVector pos