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