ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4TessellatedSolid.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4TessellatedSolid.cc
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration 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 
60 #include "geomdefs.hh"
61 #include "Randomize.hh"
62 #include "G4SystemOfUnits.hh"
63 #include "G4PhysicalConstants.hh"
64 #include "G4GeometryTolerance.hh"
65 #include "G4VoxelLimits.hh"
66 #include "G4AffineTransform.hh"
67 #include "G4BoundingEnvelope.hh"
68 
69 #include "G4PolyhedronArbitrary.hh"
70 #include "G4VGraphicsScene.hh"
71 #include "G4VisExtent.hh"
72 
73 #include "G4AutoLock.hh"
74 
75 namespace
76 {
77  G4Mutex polyhedronMutex = G4MUTEX_INITIALIZER;
78 }
79 
80 using namespace std;
81 
83 //
84 // Standard contructor has blank name and defines no fFacets.
85 //
87 {
88  Initialize();
89 }
90 
92 //
93 // Alternative constructor. Simple define name and geometry type - no fFacets
94 // to detine.
95 //
97  : G4VSolid(name)
98 {
99  Initialize();
100 }
101 
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 
117 {
118  DeleteObjects();
119 }
120 
122 //
123 // Copy constructor.
124 //
126  : G4VSolid(ts)
127 {
128  Initialize();
129 
130  CopyObjects(ts);
131 }
132 
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 
155 //
157 {
159 
160  fRebuildPolyhedron = false; fpPolyhedron = nullptr;
161  fCubicVolume = 0.; fSurfaceArea = 0.;
162 
163  fGeometryType = "G4TessellatedSolid";
164  fSolidClosed = false;
165 
168 
170 }
171 
173 //
175 {
176  G4int size = fFacets.size();
177  for (G4int i = 0; i < size; ++i) { delete fFacets[i]; }
178  fFacets.clear();
179  delete fpPolyhedron; fpPolyhedron = nullptr;
180 }
181 
183 //
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 
193  G4int n = ts.GetNumberOfFacets();
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 
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();
223  value.id = fFacetList.size();
224  value.mag2 = p.x() + p.y() + p.z();
225 
226  G4bool found = false;
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 
278 //
279 G4int 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  G4int cc = 0, nz = 0;
288 
289  vector<G4int> candidates;
290 
291  while (!pos.empty()) // Loop checking, 13.08.2015, G.Cosmo
292  {
293  xyz = pos.top();
294  pos.pop();
295  G4int index = fVoxels.GetVoxelsIndex(xyz);
296  if (!checked[index])
297  {
298  checked.SetBitNumber(index, true);
299  ++cc;
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  else
325  {
326  ++nz;
327  }
328  }
329  }
330  return filled;
331 }
332 
334 //
336 {
337  vector<G4int> voxel(3), maxVoxels(3);
338  for (auto i = 0; i <= 2; ++i) maxVoxels[i] = fVoxels.GetBoundary(i).size();
339  G4int size = maxVoxels[0] * maxVoxels[1] * maxVoxels[2];
340 
341  G4SurfBits checked(size-1);
342  fInsides.Clear();
343  fInsides.ResetBitNumber(size-1);
344 
346  for (voxel[2] = 0; voxel[2] < maxVoxels[2] - 1; ++voxel[2])
347  {
348  for (voxel[1] = 0; voxel[1] < maxVoxels[1] - 1; ++voxel[1])
349  {
350  for (voxel[0] = 0; voxel[0] < maxVoxels[0] - 1; ++voxel[0])
351  {
352  G4int index = fVoxels.GetVoxelsIndex(voxel);
353  if (!checked[index] && fVoxels.IsEmpty(index))
354  {
355  for (auto i = 0; i <= 2; ++i)
356  {
357  point[i] = fVoxels.GetBoundary(i)[voxel[i]];
358  }
359  G4bool inside = (G4bool) (InsideNoVoxels(point) == kInside);
360  SetAllUsingStack(voxel, maxVoxels, inside, checked);
361  }
362  else checked.SetBitNumber(index);
363  }
364  }
365  }
366 }
367 
369 //
371 {
372 #ifdef G4SPECSDEBUG
373  G4cout << "Voxelizing..." << G4endl;
374 #endif
376 
377  if (fVoxels.Empty().GetNbits())
378  {
379 #ifdef G4SPECSDEBUG
380  G4cout << "Precalculating Insides..." << G4endl;
381 #endif
383  }
384 }
385 
387 //
388 // Compute extremeFacets, i.e. find those facets that have surface
389 // planes that bound the volume.
390 // Note that this is going to reject concaved surfaces as being extreme. Also
391 // note that if the vertex is on the facet, displacement is zero, so IsInside
392 // returns true. So will this work?? Need non-equality
393 // "G4bool inside = displacement < 0.0;"
394 // or
395 // "G4bool inside = displacement <= -0.5*kCarTolerance"
396 // (Notes from PT 13/08/2007).
397 //
399 {
400  G4int size = fFacets.size();
401  for (G4int j = 0; j < size; ++j)
402  {
403  G4VFacet &facet = *fFacets[j];
404 
405  G4bool isExtreme = true;
406  G4int vsize = fVertexList.size();
407  for (G4int i=0; i < vsize; ++i)
408  {
409  if (!facet.IsInside(fVertexList[i]))
410  {
411  isExtreme = false;
412  break;
413  }
414  }
415  if (isExtreme) fExtremeFacets.insert(&facet);
416  }
417 }
418 
420 //
422 {
423  // The algorithm:
424  // we will have additional vertexListSorted, where all the items will be
425  // sorted by magnitude of vertice vector.
426  // New candidate for fVertexList - we will determine the position fo first
427  // item which would be within its magnitude - 0.5*kCarTolerance.
428  // We will go trough until we will reach > +0.5 kCarTolerance.
429  // Comparison (q-p).mag() < 0.5*kCarTolerance will be made.
430  // They can be just stored in std::vector, with custom insertion based
431  // on binary search.
432 
433  set<G4VertexInfo,G4VertexComparator> vertexListSorted;
434  set<G4VertexInfo,G4VertexComparator>::iterator begin
435  = vertexListSorted.begin(), end = vertexListSorted.end(), pos, it;
438 
439  fVertexList.clear();
440  G4int size = fFacets.size();
441 
442  G4double kCarTolerance24 = kCarTolerance * kCarTolerance / 4.0;
443  G4double kCarTolerance3 = 3 * kCarTolerance;
444  vector<G4int> newIndex(100);
445 
446  for (G4int k = 0; k < size; ++k)
447  {
448  G4VFacet &facet = *fFacets[k];
449  G4int max = facet.GetNumberOfVertices();
450 
451  for (G4int i = 0; i < max; ++i)
452  {
453  p = facet.GetVertex(i);
454  value.id = fVertexList.size();
455  value.mag2 = p.x() + p.y() + p.z();
456 
457  G4bool found = false;
458  G4int id = 0;
460  {
461  pos = vertexListSorted.lower_bound(value);
462  it = pos;
463  while (it != end) // Loop checking, 13.08.2015, G.Cosmo
464  {
465  id = (*it).id;
466  G4ThreeVector q = fVertexList[id];
467  G4double dif = (q-p).mag2();
468  found = (dif < kCarTolerance24);
469  if (found) break;
470  dif = q.x() + q.y() + q.z() - value.mag2;
471  if (dif > kCarTolerance3) break;
472  it++;
473  }
474 
475  if (!found && (fVertexList.size() > 1))
476  {
477  it = pos;
478  while (it != begin) // Loop checking, 13.08.2015, G.Cosmo
479  {
480  --it;
481  id = (*it).id;
482  G4ThreeVector q = fVertexList[id];
483  G4double dif = (q-p).mag2();
484  found = (dif < kCarTolerance24);
485  if (found) break;
486  dif = value.mag2 - (q.x() + q.y() + q.z());
487  if (dif > kCarTolerance3) break;
488  }
489  }
490  }
491 
492  if (!found)
493  {
494 #ifdef G4SPECSDEBUG
495  G4cout << p.x() << ":" << p.y() << ":" << p.z() << G4endl;
496  G4cout << "Adding new vertex #" << i << " of facet " << k
497  << " id " << value.id << G4endl;
498  G4cout << "===" << G4endl;
499 #endif
500  fVertexList.push_back(p);
501  vertexListSorted.insert(value);
502  begin = vertexListSorted.begin();
503  end = vertexListSorted.end();
504  newIndex[i] = value.id;
505  //
506  // Now update the maximum x, y and z limits of the volume.
507  //
508  if (value.id == 0) fMinExtent = fMaxExtent = p;
509  else
510  {
511  if (p.x() > fMaxExtent.x()) fMaxExtent.setX(p.x());
512  else if (p.x() < fMinExtent.x()) fMinExtent.setX(p.x());
513  if (p.y() > fMaxExtent.y()) fMaxExtent.setY(p.y());
514  else if (p.y() < fMinExtent.y()) fMinExtent.setY(p.y());
515  if (p.z() > fMaxExtent.z()) fMaxExtent.setZ(p.z());
516  else if (p.z() < fMinExtent.z()) fMinExtent.setZ(p.z());
517  }
518  }
519  else
520  {
521 #ifdef G4SPECSDEBUG
522  G4cout << p.x() << ":" << p.y() << ":" << p.z() << G4endl;
523  G4cout << "Vertex #" << i << " of facet " << k
524  << " found, redirecting to " << id << G4endl;
525  G4cout << "===" << G4endl;
526 #endif
527  newIndex[i] = id;
528  }
529  }
530  // only now it is possible to change vertices pointer
531  //
532  facet.SetVertices(&fVertexList);
533  for (G4int i = 0; i < max; ++i)
534  facet.SetVertexIndex(i,newIndex[i]);
535  }
536  vector<G4ThreeVector>(fVertexList).swap(fVertexList);
537 
538 #ifdef G4SPECSDEBUG
539  G4double previousValue = 0.;
540  for (auto res=vertexListSorted.cbegin(); res!=vertexListSorted.cend(); ++res)
541  {
542  G4int id = (*res).id;
543  G4ThreeVector vec = fVertexList[id];
544  G4double mvalue = vec.x() + vec.y() + vec.z();
545  if (previousValue && (previousValue - 1e-9 > mvalue))
546  G4cout << "Error in CreateVertexList: previousValue " << previousValue
547  << " is smaller than mvalue " << mvalue << G4endl;
548  previousValue = mvalue;
549  }
550 #endif
551 }
552 
554 //
556 {
558  G4int with = AllocatedMemory();
559  G4double ratio = (G4double) with / without;
560  G4cout << "G4TessellatedSolid - Allocated memory without voxel overhead "
561  << without << "; with " << with << "; ratio: " << ratio << G4endl;
562 }
563 
565 //
567 {
568  if (t)
569  {
570 #ifdef G4SPECSDEBUG
571  G4cout << "Creating vertex list..." << G4endl;
572 #endif
574 
575 #ifdef G4SPECSDEBUG
576  G4cout << "Setting extreme facets..." << G4endl;
577 #endif
579 
580 #ifdef G4SPECSDEBUG
581  G4cout << "Voxelizing..." << G4endl;
582 #endif
583  Voxelize();
584 
585 #ifdef G4SPECSDEBUG
587 #endif
588 
589  }
590  fSolidClosed = t;
591 }
592 
594 //
595 // GetSolidClosed
596 //
597 // Used to determine whether the solid is closed to adding further fFacets.
598 //
600 {
601  return fSolidClosed;
602 }
603 
605 //
606 // operator+=
607 //
608 // This operator allows the user to add two tessellated solids together, so
609 // that the solid on the left then includes all of the facets in the solid
610 // on the right. Note that copies of the facets are generated, rather than
611 // using the original facet set of the solid on the right.
612 //
615 {
616  G4int size = right.GetNumberOfFacets();
617  for (G4int i = 0; i < size; ++i)
618  AddFacet(right.GetFacet(i)->GetClone());
619 
620  return *this;
621 }
622 
624 //
625 // GetNumberOfFacets
626 //
628 {
629  return fFacets.size();
630 }
631 
633 //
635 {
636  //
637  // First the simple test - check if we're outside of the X-Y-Z extremes
638  // of the tessellated solid.
639  //
641  return kOutside;
642 
643  vector<G4int> startingVoxel(3);
644  fVoxels.GetVoxel(startingVoxel, p);
645 
646  const G4double dirTolerance = 1.0E-14;
647 
648  const vector<G4int> &startingCandidates =
649  fVoxels.GetCandidates(startingVoxel);
650  G4int limit = startingCandidates.size();
651  if (limit == 0 && fInsides.GetNbits())
652  {
653  G4int index = fVoxels.GetPointIndex(p);
654  EInside location = fInsides[index] ? kInside : kOutside;
655  return location;
656  }
657 
658  G4double minDist = kInfinity;
659 
660  for(G4int i = 0; i < limit; ++i)
661  {
662  G4int candidate = startingCandidates[i];
663  G4VFacet &facet = *fFacets[candidate];
664  G4double dist = facet.Distance(p,minDist);
665  if (dist < minDist) minDist = dist;
666  if (dist <= kCarToleranceHalf)
667  return kSurface;
668  }
669 
670  // The following is something of an adaptation of the method implemented by
671  // Rickard Holmberg augmented with information from Schneider & Eberly,
672  // "Geometric Tools for Computer Graphics," pp700-701, 2003. In essence,
673  // we're trying to determine whether we're inside the volume by projecting
674  // a few rays and determining if the first surface crossed is has a normal
675  // vector between 0 to pi/2 (out-going) or pi/2 to pi (in-going).
676  // We should also avoid rays which are nearly within the plane of the
677  // tessellated surface, and therefore produce rays randomly.
678  // For the moment, this is a bit over-engineered (belt-braces-and-ducttape).
679  //
680  G4double distOut = kInfinity;
681  G4double distIn = kInfinity;
682  G4double distO = 0.0;
683  G4double distI = 0.0;
684  G4double distFromSurfaceO = 0.0;
685  G4double distFromSurfaceI = 0.0;
686  G4ThreeVector normalO, normalI;
687  G4bool crossingO = false;
688  G4bool crossingI = false;
689  EInside location = kOutside;
690  G4int sm = 0;
691 
692  G4bool nearParallel = false;
693  do // Loop checking, 13.08.2015, G.Cosmo
694  {
695  // We loop until we find direction where the vector is not nearly parallel
696  // to the surface of any facet since this causes ambiguities. The usual
697  // case is that the angles should be sufficiently different, but there
698  // are 20 random directions to select from - hopefully sufficient.
699  //
700  distOut = distIn = kInfinity;
701  const G4ThreeVector& v = fRandir[sm];
702  ++sm;
703  //
704  // This code could be voxelized by the same algorithm, which is used for
705  // DistanceToOut(). We will traverse through fVoxels. we will call
706  // intersect only for those, which would be candidates and was not
707  // checked before.
708  //
709  G4ThreeVector currentPoint = p;
710  G4ThreeVector direction = v.unit();
711  // G4SurfBits exclusion(fVoxels.GetBitsPerSlice());
712  vector<G4int> curVoxel(3);
713  curVoxel = startingVoxel;
714  G4double shiftBonus = kCarTolerance;
715 
716  G4bool crossed = false;
717  G4bool started = true;
718 
719  do // Loop checking, 13.08.2015, G.Cosmo
720  {
721  const vector<G4int> &candidates =
722  started ? startingCandidates : fVoxels.GetCandidates(curVoxel);
723  started = false;
724  if (G4int candidatesCount = candidates.size())
725  {
726  for (G4int i = 0 ; i < candidatesCount; ++i)
727  {
728  G4int candidate = candidates[i];
729  // bits.SetBitNumber(candidate);
730  G4VFacet& facet = *fFacets[candidate];
731 
732  crossingO = facet.Intersect(p,v,true,distO,distFromSurfaceO,normalO);
733  crossingI = facet.Intersect(p,v,false,distI,distFromSurfaceI,normalI);
734 
735  if (crossingO || crossingI)
736  {
737  crossed = true;
738 
739  nearParallel = (crossingO
740  && std::fabs(normalO.dot(v))<dirTolerance)
741  || (crossingI && std::fabs(normalI.dot(v))<dirTolerance);
742  if (!nearParallel)
743  {
744  if (crossingO && distO > 0.0 && distO < distOut)
745  distOut = distO;
746  if (crossingI && distI > 0.0 && distI < distIn)
747  distIn = distI;
748  }
749  else break;
750  }
751  }
752  if (nearParallel) break;
753  }
754  else
755  {
756  if (!crossed)
757  {
758  G4int index = fVoxels.GetVoxelsIndex(curVoxel);
759  G4bool inside = fInsides[index];
760  location = inside ? kInside : kOutside;
761  return location;
762  }
763  }
764 
765  G4double shift=fVoxels.DistanceToNext(currentPoint, direction, curVoxel);
766  if (shift == kInfinity) break;
767 
768  currentPoint += direction * (shift + shiftBonus);
769  }
770  while (fVoxels.UpdateCurrentVoxel(currentPoint, direction, curVoxel));
771 
772  }
773  while (nearParallel && sm != fMaxTries);
774  //
775  // Here we loop through the facets to find out if there is an intersection
776  // between the ray and that facet. The test if performed separately whether
777  // the ray is entering the facet or exiting.
778  //
779 #ifdef G4VERBOSE
780  if (sm == fMaxTries)
781  {
782  //
783  // We've run out of random vector directions. If nTries is set sufficiently
784  // low (nTries <= 0.5*maxTries) then this would indicate that there is
785  // something wrong with geometry.
786  //
787  std::ostringstream message;
788  G4int oldprc = message.precision(16);
789  message << "Cannot determine whether point is inside or outside volume!"
790  << G4endl
791  << "Solid name = " << GetName() << G4endl
792  << "Geometry Type = " << fGeometryType << G4endl
793  << "Number of facets = " << fFacets.size() << G4endl
794  << "Position:" << G4endl << G4endl
795  << "p.x() = " << p.x()/mm << " mm" << G4endl
796  << "p.y() = " << p.y()/mm << " mm" << G4endl
797  << "p.z() = " << p.z()/mm << " mm";
798  message.precision(oldprc);
799  G4Exception("G4TessellatedSolid::Inside()",
800  "GeomSolids1002", JustWarning, message);
801  }
802 #endif
803 
804  // In the next if-then-elseif G4String the logic is as follows:
805  // (1) You don't hit anything so cannot be inside volume, provided volume
806  // constructed correctly!
807  // (2) Distance to inside (ie. nearest facet such that you enter facet) is
808  // shorter than distance to outside (nearest facet such that you exit
809  // facet) - on condition of safety distance - therefore we're outside.
810  // (3) Distance to outside is shorter than distance to inside therefore
811  // we're inside.
812  //
813  if (distIn == kInfinity && distOut == kInfinity)
814  location = kOutside;
815  else if (distIn <= distOut - kCarToleranceHalf)
816  location = kOutside;
817  else if (distOut <= distIn - kCarToleranceHalf)
818  location = kInside;
819 
820  return location;
821 }
822 
824 //
826 {
827  //
828  // First the simple test - check if we're outside of the X-Y-Z extremes
829  // of the tessellated solid.
830  //
832  return kOutside;
833 
834  const G4double dirTolerance = 1.0E-14;
835 
836  G4double minDist = kInfinity;
837  //
838  // Check if we are close to a surface
839  //
840  G4int size = fFacets.size();
841  for (G4int i = 0; i < size; ++i)
842  {
843  G4VFacet& facet = *fFacets[i];
844  G4double dist = facet.Distance(p,minDist);
845  if (dist < minDist) minDist = dist;
846  if (dist <= kCarToleranceHalf)
847  {
848  return kSurface;
849  }
850  }
851  //
852  // The following is something of an adaptation of the method implemented by
853  // Rickard Holmberg augmented with information from Schneider & Eberly,
854  // "Geometric Tools for Computer Graphics," pp700-701, 2003. In essence, we're
855  // trying to determine whether we're inside the volume by projecting a few
856  // rays and determining if the first surface crossed is has a normal vector
857  // between 0 to pi/2 (out-going) or pi/2 to pi (in-going). We should also
858  // avoid rays which are nearly within the plane of the tessellated surface,
859  // and therefore produce rays randomly. For the moment, this is a bit
860  // over-engineered (belt-braces-and-ducttape).
861  //
862 #if G4SPECSDEBUG
863  G4int nTry = 7;
864 #else
865  G4int nTry = 3;
866 #endif
867  G4double distOut = kInfinity;
868  G4double distIn = kInfinity;
869  G4double distO = 0.0;
870  G4double distI = 0.0;
871  G4double distFromSurfaceO = 0.0;
872  G4double distFromSurfaceI = 0.0;
873  G4ThreeVector normalO(0.0,0.0,0.0);
874  G4ThreeVector normalI(0.0,0.0,0.0);
875  G4bool crossingO = false;
876  G4bool crossingI = false;
877  EInside location = kOutside;
878  EInside locationprime = kOutside;
879  G4int sm = 0;
880 
881  for (G4int i=0; i<nTry; ++i)
882  {
883  G4bool nearParallel = false;
884  do // Loop checking, 13.08.2015, G.Cosmo
885  {
886  //
887  // We loop until we find direction where the vector is not nearly parallel
888  // to the surface of any facet since this causes ambiguities. The usual
889  // case is that the angles should be sufficiently different, but there
890  // are 20 random directions to select from - hopefully sufficient.
891  //
892  distOut = distIn = kInfinity;
894  sm++;
895  vector<G4VFacet*>::const_iterator f = fFacets.begin();
896 
897  do // Loop checking, 13.08.2015, G.Cosmo
898  {
899  //
900  // Here we loop through the facets to find out if there is an
901  // intersection between the ray and that facet. The test if performed
902  // separately whether the ray is entering the facet or exiting.
903  //
904  crossingO = ((*f)->Intersect(p,v,true,distO,distFromSurfaceO,normalO));
905  crossingI = ((*f)->Intersect(p,v,false,distI,distFromSurfaceI,normalI));
906  if (crossingO || crossingI)
907  {
908  nearParallel = (crossingO && std::fabs(normalO.dot(v))<dirTolerance)
909  || (crossingI && std::fabs(normalI.dot(v))<dirTolerance);
910  if (!nearParallel)
911  {
912  if (crossingO && distO > 0.0 && distO < distOut) distOut = distO;
913  if (crossingI && distI > 0.0 && distI < distIn) distIn = distI;
914  }
915  }
916  } while (!nearParallel && ++f != fFacets.end());
917  } while (nearParallel && sm != fMaxTries);
918 
919 #ifdef G4VERBOSE
920  if (sm == fMaxTries)
921  {
922  //
923  // We've run out of random vector directions. If nTries is set
924  // sufficiently low (nTries <= 0.5*maxTries) then this would indicate
925  // that there is something wrong with geometry.
926  //
927  std::ostringstream message;
928  G4int oldprc = message.precision(16);
929  message << "Cannot determine whether point is inside or outside volume!"
930  << G4endl
931  << "Solid name = " << GetName() << G4endl
932  << "Geometry Type = " << fGeometryType << G4endl
933  << "Number of facets = " << fFacets.size() << G4endl
934  << "Position:" << G4endl << G4endl
935  << "p.x() = " << p.x()/mm << " mm" << G4endl
936  << "p.y() = " << p.y()/mm << " mm" << G4endl
937  << "p.z() = " << p.z()/mm << " mm";
938  message.precision(oldprc);
939  G4Exception("G4TessellatedSolid::Inside()",
940  "GeomSolids1002", JustWarning, message);
941  }
942 #endif
943  //
944  // In the next if-then-elseif G4String the logic is as follows:
945  // (1) You don't hit anything so cannot be inside volume, provided volume
946  // constructed correctly!
947  // (2) Distance to inside (ie. nearest facet such that you enter facet) is
948  // shorter than distance to outside (nearest facet such that you exit
949  // facet) - on condition of safety distance - therefore we're outside.
950  // (3) Distance to outside is shorter than distance to inside therefore
951  // we're inside.
952  //
953  if (distIn == kInfinity && distOut == kInfinity)
954  locationprime = kOutside;
955  else if (distIn <= distOut - kCarToleranceHalf)
956  locationprime = kOutside;
957  else if (distOut <= distIn - kCarToleranceHalf)
958  locationprime = kInside;
959 
960  if (i == 0) location = locationprime;
961  }
962 
963  return location;
964 }
965 
967 //
968 // Return the outwards pointing unit normal of the shape for the
969 // surface closest to the point at offset p.
970 //
972  G4ThreeVector& aNormal) const
973 {
974  G4double minDist;
975  G4VFacet* facet = nullptr;
976 
977  if (fVoxels.GetCountOfVoxels() > 1)
978  {
979  vector<G4int> curVoxel(3);
980  fVoxels.GetVoxel(curVoxel, p);
981  const vector<G4int> &candidates = fVoxels.GetCandidates(curVoxel);
982  // fVoxels.GetCandidatesVoxelArray(p, candidates, 0);
983 
984  if (G4int limit = candidates.size())
985  {
986  minDist = kInfinity;
987  for(G4int i = 0 ; i < limit ; ++i)
988  {
989  G4int candidate = candidates[i];
990  G4VFacet &fct = *fFacets[candidate];
991  G4double dist = fct.Distance(p,minDist);
992  if (dist < minDist) minDist = dist;
993  if (dist <= kCarToleranceHalf)
994  {
995  aNormal = fct.GetSurfaceNormal();
996  return true;
997  }
998  }
999  }
1000  minDist = MinDistanceFacet(p, true, facet);
1001  }
1002  else
1003  {
1004  minDist = kInfinity;
1005  G4int size = fFacets.size();
1006  for (G4int i = 0; i < size; ++i)
1007  {
1008  G4VFacet& f = *fFacets[i];
1009  G4double dist = f.Distance(p, minDist);
1010  if (dist < minDist)
1011  {
1012  minDist = dist;
1013  facet = &f;
1014  }
1015  }
1016  }
1017 
1018  if (minDist != kInfinity)
1019  {
1020  if (facet) { aNormal = facet->GetSurfaceNormal(); }
1021  return minDist <= kCarToleranceHalf;
1022  }
1023  else
1024  {
1025 #ifdef G4VERBOSE
1026  std::ostringstream message;
1027  message << "Point p is not on surface !?" << G4endl
1028  << " No facets found for point: " << p << " !" << G4endl
1029  << " Returning approximated value for normal.";
1030 
1031  G4Exception("G4TessellatedSolid::SurfaceNormal(p)",
1032  "GeomSolids1002", JustWarning, message );
1033 #endif
1034  aNormal = (p.z() > 0 ? G4ThreeVector(0,0,1) : G4ThreeVector(0,0,-1));
1035  return false;
1036  }
1037 }
1038 
1040 //
1041 // G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v)
1042 //
1043 // Return the distance along the normalised vector v to the shape,
1044 // from the point at offset p. If there is no intersection, return
1045 // kInfinity. The first intersection resulting from 'leaving' a
1046 // surface/volume is discarded. Hence, this is tolerant of points on
1047 // surface of shape.
1048 //
1049 G4double
1051  const G4ThreeVector& v,
1052  G4double /*aPstep*/) const
1053 {
1054  G4double minDist = kInfinity;
1055  G4double dist = 0.0;
1056  G4double distFromSurface = 0.0;
1058 
1059 #if G4SPECSDEBUG
1060  if (Inside(p) == kInside )
1061  {
1062  std::ostringstream message;
1063  G4int oldprc = message.precision(16) ;
1064  message << "Point p is already inside!?" << G4endl
1065  << "Position:" << G4endl << G4endl
1066  << " p.x() = " << p.x()/mm << " mm" << G4endl
1067  << " p.y() = " << p.y()/mm << " mm" << G4endl
1068  << " p.z() = " << p.z()/mm << " mm" << G4endl
1069  << "DistanceToOut(p) == " << DistanceToOut(p);
1070  message.precision(oldprc) ;
1071  G4Exception("G4TriangularFacet::DistanceToIn(p,v)",
1072  "GeomSolids1002", JustWarning, message);
1073  }
1074 #endif
1075 
1076  G4int size = fFacets.size();
1077  for (G4int i = 0; i < size; ++i)
1078  {
1079  G4VFacet& facet = *fFacets[i];
1080  if (facet.Intersect(p,v,false,dist,distFromSurface,normal))
1081  {
1082  //
1083  // set minDist to the new distance to current facet if distFromSurface
1084  // is in positive direction and point is not at surface. If the point is
1085  // within 0.5*kCarTolerance of the surface, then force distance to be
1086  // zero and leave member function immediately (for efficiency), as
1087  // proposed by & credit to Akira Okumura.
1088  //
1089  if (distFromSurface > kCarToleranceHalf && dist >= 0.0 && dist < minDist)
1090  {
1091  minDist = dist;
1092  }
1093  else
1094  {
1095  if (-kCarToleranceHalf <= dist && dist <= kCarToleranceHalf)
1096  {
1097  return 0.0;
1098  }
1099  else
1100  {
1101  if (distFromSurface > -kCarToleranceHalf
1102  && distFromSurface < kCarToleranceHalf)
1103  {
1104  minDist = dist;
1105  }
1106  }
1107  }
1108  }
1109  }
1110  return minDist;
1111 }
1112 
1114 //
1115 G4double
1117  const G4ThreeVector& v,
1118  G4ThreeVector& aNormalVector,
1119  G4bool& aConvex,
1120  G4double /*aPstep*/) const
1121 {
1122  G4double minDist = kInfinity;
1123  G4double dist = 0.0;
1124  G4double distFromSurface = 0.0;
1125  G4ThreeVector normal, minNormal;
1126 
1127 #if G4SPECSDEBUG
1128  if ( Inside(p) == kOutside )
1129  {
1130  std::ostringstream message;
1131  G4int oldprc = message.precision(16) ;
1132  message << "Point p is already outside!?" << G4endl
1133  << "Position:" << G4endl << G4endl
1134  << " p.x() = " << p.x()/mm << " mm" << G4endl
1135  << " p.y() = " << p.y()/mm << " mm" << G4endl
1136  << " p.z() = " << p.z()/mm << " mm" << G4endl
1137  << "DistanceToIn(p) == " << DistanceToIn(p);
1138  message.precision(oldprc) ;
1139  G4Exception("G4TriangularFacet::DistanceToOut(p)",
1140  "GeomSolids1002", JustWarning, message);
1141  }
1142 #endif
1143 
1144  G4bool isExtreme = false;
1145  G4int size = fFacets.size();
1146  for (G4int i = 0; i < size; ++i)
1147  {
1148  G4VFacet& facet = *fFacets[i];
1149  if (facet.Intersect(p,v,true,dist,distFromSurface,normal))
1150  {
1151  if (distFromSurface > 0.0 && distFromSurface <= kCarToleranceHalf
1152  && facet.Distance(p,kCarTolerance) <= kCarToleranceHalf)
1153  {
1154  // We are on a surface. Return zero.
1155  aConvex = (fExtremeFacets.find(&facet) != fExtremeFacets.end());
1156  // Normal(p, aNormalVector);
1157  // aNormalVector = facet.GetSurfaceNormal();
1158  aNormalVector = normal;
1159  return 0.0;
1160  }
1161  if (dist >= 0.0 && dist < minDist)
1162  {
1163  minDist = dist;
1164  minNormal = normal;
1165  isExtreme = (fExtremeFacets.find(&facet) != fExtremeFacets.end());
1166  }
1167  }
1168  }
1169  if (minDist < kInfinity)
1170  {
1171  aNormalVector = minNormal;
1172  aConvex = isExtreme;
1173  return minDist;
1174  }
1175  else
1176  {
1177  // No intersection found
1178  aConvex = false;
1179  Normal(p, aNormalVector);
1180  return 0.0;
1181  }
1182 }
1183 
1185 //
1187 DistanceToOutCandidates(const std::vector<G4int>& candidates,
1188  const G4ThreeVector& aPoint,
1189  const G4ThreeVector& direction,
1190  G4double& minDist, G4ThreeVector& minNormal,
1191  G4int& minCandidate ) const
1192 {
1193  G4int candidatesCount = candidates.size();
1194  G4double dist = 0.0;
1195  G4double distFromSurface = 0.0;
1197 
1198  for (G4int i = 0 ; i < candidatesCount; ++i)
1199  {
1200  G4int candidate = candidates[i];
1201  G4VFacet& facet = *fFacets[candidate];
1202  if (facet.Intersect(aPoint,direction,true,dist,distFromSurface,normal))
1203  {
1204  if (distFromSurface > 0.0 && distFromSurface <= kCarToleranceHalf
1205  && facet.Distance(aPoint,kCarTolerance) <= kCarToleranceHalf)
1206  {
1207  // We are on a surface
1208  //
1209  minDist = 0.0;
1210  minNormal = normal;
1211  minCandidate = candidate;
1212  break;
1213  }
1214  if (dist >= 0.0 && dist < minDist)
1215  {
1216  minDist = dist;
1217  minNormal = normal;
1218  minCandidate = candidate;
1219  }
1220  }
1221  }
1222 }
1223 
1225 //
1226 G4double
1228  const G4ThreeVector& aDirection,
1229  G4ThreeVector& aNormalVector,
1230  G4bool &aConvex,
1231  G4double aPstep) const
1232 {
1233  G4double minDistance;
1234 
1235  if (fVoxels.GetCountOfVoxels() > 1)
1236  {
1237  minDistance = kInfinity;
1238 
1239  G4ThreeVector currentPoint = aPoint;
1240  G4ThreeVector direction = aDirection.unit();
1241  G4double totalShift = 0.;
1242  vector<G4int> curVoxel(3);
1243  if (!fVoxels.Contains(aPoint)) return 0.;
1244 
1245  fVoxels.GetVoxel(curVoxel, currentPoint);
1246 
1247  G4double shiftBonus = kCarTolerance;
1248 
1249  const vector<G4int>* old = nullptr;
1250 
1251  G4int minCandidate = -1;
1252  do // Loop checking, 13.08.2015, G.Cosmo
1253  {
1254  const vector<G4int>& candidates = fVoxels.GetCandidates(curVoxel);
1255  if (old == &candidates)
1256  ++old;
1257  if (old != &candidates && candidates.size())
1258  {
1259  DistanceToOutCandidates(candidates, aPoint, direction, minDistance,
1260  aNormalVector, minCandidate);
1261  if (minDistance <= totalShift) break;
1262  }
1263 
1264  G4double shift=fVoxels.DistanceToNext(currentPoint, direction, curVoxel);
1265  if (shift == kInfinity) break;
1266 
1267  totalShift += shift;
1268  if (minDistance <= totalShift) break;
1269 
1270  currentPoint += direction * (shift + shiftBonus);
1271 
1272  old = &candidates;
1273  }
1274  while (fVoxels.UpdateCurrentVoxel(currentPoint, direction, curVoxel));
1275 
1276  if (minCandidate < 0)
1277  {
1278  // No intersection found
1279  minDistance = 0.;
1280  aConvex = false;
1281  Normal(aPoint, aNormalVector);
1282  }
1283  else
1284  {
1285  aConvex = (fExtremeFacets.find(fFacets[minCandidate])
1286  != fExtremeFacets.end());
1287  }
1288  }
1289  else
1290  {
1291  minDistance = DistanceToOutNoVoxels(aPoint, aDirection, aNormalVector,
1292  aConvex, aPstep);
1293  }
1294  return minDistance;
1295 }
1296 
1298 //
1300 DistanceToInCandidates(const std::vector<G4int>& candidates,
1301  const G4ThreeVector& aPoint,
1302  const G4ThreeVector& direction) const
1303 {
1304  G4int candidatesCount = candidates.size();
1305  G4double dist = 0.0;
1306  G4double distFromSurface = 0.0;
1308 
1309  G4double minDistance = kInfinity;
1310  for (G4int i = 0 ; i < candidatesCount; ++i)
1311  {
1312  G4int candidate = candidates[i];
1313  G4VFacet& facet = *fFacets[candidate];
1314  if (facet.Intersect(aPoint,direction,false,dist,distFromSurface,normal))
1315  {
1316  //
1317  // Set minDist to the new distance to current facet if distFromSurface is
1318  // in positive direction and point is not at surface. If the point is
1319  // within 0.5*kCarTolerance of the surface, then force distance to be
1320  // zero and leave member function immediately (for efficiency), as
1321  // proposed by & credit to Akira Okumura.
1322  //
1323  if ( (distFromSurface > kCarToleranceHalf)
1324  && (dist >= 0.0) && (dist < minDistance))
1325  {
1326  minDistance = dist;
1327  }
1328  else
1329  {
1330  if (-kCarToleranceHalf <= dist && dist <= kCarToleranceHalf)
1331  {
1332  return 0.0;
1333  }
1334  else if (distFromSurface > -kCarToleranceHalf
1335  && distFromSurface < kCarToleranceHalf)
1336  {
1337  minDistance = dist;
1338  }
1339  }
1340  }
1341  }
1342  return minDistance;
1343 }
1344 
1346 //
1347 G4double
1349  const G4ThreeVector& aDirection,
1350  G4double aPstep) const
1351 {
1352  G4double minDistance;
1353 
1354  if (fVoxels.GetCountOfVoxels() > 1)
1355  {
1356  minDistance = kInfinity;
1357  G4ThreeVector currentPoint = aPoint;
1358  G4ThreeVector direction = aDirection.unit();
1359  G4double shift = fVoxels.DistanceToFirst(currentPoint, direction);
1360  if (shift == kInfinity) return shift;
1361  G4double shiftBonus = kCarTolerance;
1362  if (shift)
1363  currentPoint += direction * (shift + shiftBonus);
1364  // if (!fVoxels.Contains(currentPoint)) return minDistance;
1365  G4double totalShift = shift;
1366 
1367  // G4SurfBits exclusion; // (1/*fVoxels.GetBitsPerSlice()*/);
1368  vector<G4int> curVoxel(3);
1369 
1370  fVoxels.GetVoxel(curVoxel, currentPoint);
1371  do // Loop checking, 13.08.2015, G.Cosmo
1372  {
1373  const vector<G4int>& candidates = fVoxels.GetCandidates(curVoxel);
1374  if (candidates.size())
1375  {
1376  G4double distance=DistanceToInCandidates(candidates, aPoint, direction);
1377  if (minDistance > distance) minDistance = distance;
1378  if (distance < totalShift) break;
1379  }
1380 
1381  shift = fVoxels.DistanceToNext(currentPoint, direction, curVoxel);
1382  if (shift == kInfinity /*|| shift == 0*/) break;
1383 
1384  totalShift += shift;
1385  if (minDistance < totalShift) break;
1386 
1387  currentPoint += direction * (shift + shiftBonus);
1388  }
1389  while (fVoxels.UpdateCurrentVoxel(currentPoint, direction, curVoxel));
1390  }
1391  else
1392  {
1393  minDistance = DistanceToInNoVoxels(aPoint, aDirection, aPstep);
1394  }
1395 
1396  return minDistance;
1397 }
1398 
1400 //
1401 G4bool
1402 G4TessellatedSolid::CompareSortedVoxel(const std::pair<G4int, G4double>& l,
1403  const std::pair<G4int, G4double>& r)
1404 {
1405  return l.second < r.second;
1406 }
1407 
1409 //
1410 G4double
1412  G4bool simple,
1413  G4VFacet* &minFacet) const
1414 {
1415  G4double minDist = kInfinity;
1416 
1417  G4int size = fVoxels.GetVoxelBoxesSize();
1418  vector<pair<G4int, G4double> > voxelsSorted(size);
1419 
1420  pair<G4int, G4double> info;
1421 
1422  for (G4int i = 0; i < size; ++i)
1423  {
1424  const G4VoxelBox& voxelBox = fVoxels.GetVoxelBox(i);
1425 
1426  G4ThreeVector pointShifted = p - voxelBox.pos;
1427  G4double safety = fVoxels.MinDistanceToBox(pointShifted, voxelBox.hlen);
1428  info.first = i;
1429  info.second = safety;
1430  voxelsSorted[i] = info;
1431  }
1432 
1433  std::sort(voxelsSorted.begin(), voxelsSorted.end(),
1435 
1436  for (G4int i = 0; i < size; ++i)
1437  {
1438  const pair<G4int,G4double>& inf = voxelsSorted[i];
1439  G4double dist = inf.second;
1440  if (dist > minDist) break;
1441 
1442  const vector<G4int>& candidates = fVoxels.GetVoxelBoxCandidates(inf.first);
1443  G4int csize = candidates.size();
1444  for (G4int j = 0; j < csize; ++j)
1445  {
1446  G4int candidate = candidates[j];
1447  G4VFacet& facet = *fFacets[candidate];
1448  dist = simple ? facet.Distance(p,minDist)
1449  : facet.Distance(p,minDist,false);
1450  if (dist < minDist)
1451  {
1452  minDist = dist;
1453  minFacet = &facet;
1454  }
1455  }
1456  }
1457  return minDist;
1458 }
1459 
1461 //
1463  G4bool aAccurate) const
1464 {
1465 #if G4SPECSDEBUG
1466  if ( Inside(p) == kInside )
1467  {
1468  std::ostringstream message;
1469  G4int oldprc = message.precision(16) ;
1470  message << "Point p is already inside!?" << G4endl
1471  << "Position:" << G4endl << G4endl
1472  << "p.x() = " << p.x()/mm << " mm" << G4endl
1473  << "p.y() = " << p.y()/mm << " mm" << G4endl
1474  << "p.z() = " << p.z()/mm << " mm" << G4endl
1475  << "DistanceToOut(p) == " << DistanceToOut(p);
1476  message.precision(oldprc) ;
1477  G4Exception("G4TriangularFacet::DistanceToIn(p)",
1478  "GeomSolids1002", JustWarning, message);
1479  }
1480 #endif
1481 
1482  G4double minDist;
1483 
1484  if (fVoxels.GetCountOfVoxels() > 1)
1485  {
1486  if (!aAccurate)
1487  return fVoxels.DistanceToBoundingBox(p);
1488 
1489  if (!OutsideOfExtent(p, kCarTolerance))
1490  {
1491  vector<G4int> startingVoxel(3);
1492  fVoxels.GetVoxel(startingVoxel, p);
1493  const vector<G4int> &candidates = fVoxels.GetCandidates(startingVoxel);
1494  if (candidates.size() == 0 && fInsides.GetNbits())
1495  {
1496  G4int index = fVoxels.GetPointIndex(p);
1497  if (fInsides[index]) return 0.;
1498  }
1499  }
1500 
1501  G4VFacet* facet;
1502  minDist = MinDistanceFacet(p, true, facet);
1503  }
1504  else
1505  {
1506  minDist = kInfinity;
1507  G4int size = fFacets.size();
1508  for (G4int i = 0; i < size; ++i)
1509  {
1510  G4VFacet& facet = *fFacets[i];
1511  G4double dist = facet.Distance(p,minDist);
1512  if (dist < minDist) minDist = dist;
1513  }
1514  }
1515  return minDist;
1516 }
1517 
1519 //
1520 G4double
1522 {
1523 #if G4SPECSDEBUG
1524  if ( Inside(p) == kOutside )
1525  {
1526  std::ostringstream message;
1527  G4int oldprc = message.precision(16) ;
1528  message << "Point p is already outside!?" << G4endl
1529  << "Position:" << G4endl << G4endl
1530  << "p.x() = " << p.x()/mm << " mm" << G4endl
1531  << "p.y() = " << p.y()/mm << " mm" << G4endl
1532  << "p.z() = " << p.z()/mm << " mm" << G4endl
1533  << "DistanceToIn(p) == " << DistanceToIn(p);
1534  message.precision(oldprc) ;
1535  G4Exception("G4TriangularFacet::DistanceToOut(p)",
1536  "GeomSolids1002", JustWarning, message);
1537  }
1538 #endif
1539 
1540  G4double minDist;
1541 
1542  if (OutsideOfExtent(p, kCarTolerance)) return 0.0;
1543 
1544  if (fVoxels.GetCountOfVoxels() > 1)
1545  {
1546  G4VFacet* facet;
1547  minDist = MinDistanceFacet(p, true, facet);
1548  }
1549  else
1550  {
1551  minDist = kInfinity;
1552  G4double dist = 0.0;
1553  G4int size = fFacets.size();
1554  for (G4int i = 0; i < size; ++i)
1555  {
1556  G4VFacet& facet = *fFacets[i];
1557  dist = facet.Distance(p,minDist);
1558  if (dist < minDist) minDist = dist;
1559  }
1560  }
1561  return minDist;
1562 }
1563 
1565 //
1566 // G4GeometryType GetEntityType() const;
1567 //
1568 // Provide identification of the class of an object
1569 //
1571 {
1572  return fGeometryType;
1573 }
1574 
1576 //
1577 std::ostream &G4TessellatedSolid::StreamInfo(std::ostream &os) const
1578 {
1579  os << G4endl;
1580  os << "Solid name = " << GetName() << G4endl;
1581  os << "Geometry Type = " << fGeometryType << G4endl;
1582  os << "Number of facets = " << fFacets.size() << G4endl;
1583 
1584  G4int size = fFacets.size();
1585  for (G4int i = 0; i < size; ++i)
1586  {
1587  os << "FACET # = " << i + 1 << G4endl;
1588  G4VFacet &facet = *fFacets[i];
1589  facet.StreamInfo(os);
1590  }
1591  os << G4endl;
1592 
1593  return os;
1594 }
1595 
1597 //
1598 // Make a clone of the object
1599 //
1601 {
1602  return new G4TessellatedSolid(*this);
1603 }
1604 
1606 //
1607 // EInside G4TessellatedSolid::Inside (const G4ThreeVector &p) const
1608 //
1609 // This method must return:
1610 // * kOutside if the point at offset p is outside the shape
1611 // boundaries plus kCarTolerance/2,
1612 // * kSurface if the point is <= kCarTolerance/2 from a surface, or
1613 // * kInside otherwise.
1614 //
1616 {
1617  EInside location;
1618 
1619  if (fVoxels.GetCountOfVoxels() > 1)
1620  {
1621  location = InsideVoxels(aPoint);
1622  }
1623  else
1624  {
1625  location = InsideNoVoxels(aPoint);
1626  }
1627  return location;
1628 }
1629 
1631 //
1633 {
1634  G4ThreeVector n;
1635  Normal(p, n);
1636  return n;
1637 }
1638 
1640 //
1641 // G4double DistanceToIn(const G4ThreeVector& p)
1642 //
1643 // Calculate distance to nearest surface of shape from an outside point p. The
1644 // distance can be an underestimate.
1645 //
1647 {
1648  return SafetyFromOutside(p, false);
1649 }
1650 
1652 //
1654  const G4ThreeVector& v)const
1655 {
1656  G4double dist = DistanceToInCore(p,v,kInfinity);
1657 #ifdef G4SPECSDEBUG
1658  if (dist < kInfinity)
1659  {
1660  if (Inside(p + dist*v) != kSurface)
1661  {
1662  std::ostringstream message;
1663  message << "Invalid response from facet in solid '" << GetName() << "',"
1664  << G4endl
1665  << "at point: " << p << "and direction: " << v;
1666  G4Exception("G4TessellatedSolid::DistanceToIn(p,v)",
1667  "GeomSolids1002", JustWarning, message);
1668  }
1669  }
1670 #endif
1671  return dist;
1672 }
1673 
1675 //
1676 // G4double DistanceToOut(const G4ThreeVector& p)
1677 //
1678 // Calculate distance to nearest surface of shape from an inside
1679 // point. The distance can be an underestimate.
1680 //
1682 {
1683  return SafetyFromInside(p, false);
1684 }
1685 
1687 //
1688 // G4double DistanceToOut(const G4ThreeVector& p, const G4ThreeVector& v,
1689 // const G4bool calcNorm=false,
1690 // G4bool *validNorm=0, G4ThreeVector *n=0);
1691 //
1692 // Return distance along the normalised vector v to the shape, from a
1693 // point at an offset p inside or on the surface of the
1694 // shape. Intersections with surfaces, when the point is not greater
1695 // than kCarTolerance/2 from a surface, must be ignored.
1696 // If calcNorm is true, then it must also set validNorm to either
1697 // * true, if the solid lies entirely behind or on the exiting
1698 // surface. Then it must set n to the outwards normal vector
1699 // (the Magnitude of the vector is not defined).
1700 // * false, if the solid does not lie entirely behind or on the
1701 // exiting surface.
1702 // If calcNorm is false, then validNorm and n are unused.
1703 //
1705  const G4ThreeVector& v,
1706  const G4bool calcNorm,
1707  G4bool* validNorm,
1708  G4ThreeVector* norm) const
1709 {
1710  G4ThreeVector n;
1711  G4bool valid;
1712 
1713  G4double dist = DistanceToOutCore(p, v, n, valid);
1714  if (calcNorm)
1715  {
1716  *norm = n;
1717  *validNorm = valid;
1718  }
1719 #ifdef G4SPECSDEBUG
1720  if (dist < kInfinity)
1721  {
1722  if (Inside(p + dist*v) != kSurface)
1723  {
1724  std::ostringstream message;
1725  message << "Invalid response from facet in solid '" << GetName() << "',"
1726  << G4endl
1727  << "at point: " << p << "and direction: " << v;
1728  G4Exception("G4TessellatedSolid::DistanceToOut(p,v,..)",
1729  "GeomSolids1002", JustWarning, message);
1730  }
1731  }
1732 #endif
1733  return dist;
1734 }
1735 
1737 //
1739 {
1740  scene.AddSolid (*this);
1741 }
1742 
1744 //
1746 {
1747  G4int nVertices = fVertexList.size();
1748  G4int nFacets = fFacets.size();
1749  G4PolyhedronArbitrary* polyhedron =
1750  new G4PolyhedronArbitrary (nVertices, nFacets);
1751  for (auto v= fVertexList.cbegin(); v!=fVertexList.cend(); ++v)
1752  {
1753  polyhedron->AddVertex(*v);
1754  }
1755 
1756  G4int size = fFacets.size();
1757  for (G4int i = 0; i < size; ++i)
1758  {
1759  G4VFacet* facet = fFacets[i];
1760  G4int v[4] = {0};
1761  G4int n = facet->GetNumberOfVertices();
1762  if (n > 4) n = 4;
1763  for (G4int j=0; j<n; ++j)
1764  {
1765  G4int k = facet->GetVertexIndex(j);
1766  v[j] = k+1;
1767  }
1768  polyhedron->AddFacet(v[0],v[1],v[2],v[3]);
1769  }
1770  polyhedron->SetReferences();
1771 
1772  return (G4Polyhedron*) polyhedron;
1773 }
1774 
1776 //
1777 // GetPolyhedron
1778 //
1780 {
1781  if (fpPolyhedron == nullptr ||
1785  {
1786  G4AutoLock l(&polyhedronMutex);
1787  delete fpPolyhedron;
1789  fRebuildPolyhedron = false;
1790  l.unlock();
1791  }
1792  return fpPolyhedron;
1793 }
1794 
1796 //
1797 // Get bounding box
1798 //
1800  G4ThreeVector& pMax) const
1801 {
1802  pMin = fMinExtent;
1803  pMax = fMaxExtent;
1804 
1805  // Check correctness of the bounding box
1806  //
1807  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
1808  {
1809  std::ostringstream message;
1810  message << "Bad bounding box (min >= max) for solid: "
1811  << GetName() << " !"
1812  << "\npMin = " << pMin
1813  << "\npMax = " << pMax;
1814  G4Exception("G4TessellatedSolid::BoundingLimits()",
1815  "GeomMgt0001", JustWarning, message);
1816  DumpInfo();
1817  }
1818 }
1819 
1821 //
1822 // Calculate extent under transform and specified limit
1823 //
1824 G4bool
1826  const G4VoxelLimits& pVoxelLimit,
1827  const G4AffineTransform& pTransform,
1828  G4double& pMin, G4double& pMax) const
1829 {
1830  G4ThreeVector bmin, bmax;
1831 
1832  // Check bounding box (bbox)
1833  //
1834  BoundingLimits(bmin,bmax);
1835  G4BoundingEnvelope bbox(bmin,bmax);
1836 
1837  // Use simple bounding-box to help in the case of complex meshes
1838  //
1839  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
1840 
1841 #if 0
1842  // Precise extent computation (disabled by default for this shape)
1843  //
1844  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
1845  {
1846  return (pMin < pMax) ? true : false;
1847  }
1848 
1849  // The extent is calculated as cumulative extent of the pyramids
1850  // formed by facets and the center of the bounding box.
1851  //
1852  G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
1853  G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
1854 
1856  G4ThreeVectorList apex(1);
1857  std::vector<const G4ThreeVectorList *> pyramid(2);
1858  pyramid[0] = &base;
1859  pyramid[1] = &apex;
1860  apex[0] = (bmin+bmax)*0.5;
1861 
1862  // main loop along facets
1863  pMin = kInfinity;
1864  pMax = -kInfinity;
1865  for (G4int i=0; i<GetNumberOfFacets(); ++i)
1866  {
1867  G4VFacet* facet = GetFacet(i);
1868  if (std::abs((facet->GetSurfaceNormal()).dot(facet->GetVertex(0)-apex[0]))
1869  < kCarToleranceHalf) continue;
1870 
1871  G4int nv = facet->GetNumberOfVertices();
1872  base.resize(nv);
1873  for (G4int k=0; k<nv; ++k) { base[k] = facet->GetVertex(k); }
1874 
1875  G4double emin,emax;
1876  G4BoundingEnvelope benv(pyramid);
1877  if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
1878  if (emin < pMin) pMin = emin;
1879  if (emax > pMax) pMax = emax;
1880  if (eminlim > pMin && emaxlim < pMax) break; // max possible extent
1881  }
1882  return (pMin < pMax);
1883 #endif
1884 }
1885 
1887 //
1889 {
1890  return fMinExtent.x();
1891 }
1892 
1894 //
1896 {
1897  return fMaxExtent.x();
1898 }
1899 
1901 //
1903 {
1904  return fMinExtent.y();
1905 }
1906 
1908 //
1910 {
1911  return fMaxExtent.y();
1912 }
1913 
1915 //
1917 {
1918  return fMinExtent.z();
1919 }
1920 
1922 //
1924 {
1925  return fMaxExtent.z();
1926 }
1927 
1929 //
1931 {
1932  return G4VisExtent (fMinExtent.x(), fMaxExtent.x(),
1933  fMinExtent.y(), fMaxExtent.y(),
1934  fMinExtent.z(), fMaxExtent.z());
1935 }
1936 
1938 //
1940 {
1941  if (fCubicVolume != 0.) return fCubicVolume;
1942 
1943  // For explanation of the following algorithm see:
1944  // https://en.wikipedia.org/wiki/Polyhedron#Volume
1945  // http://wwwf.imperial.ac.uk/~rn/centroid.pdf
1946 
1947  G4int size = fFacets.size();
1948  for (G4int i = 0; i < size; ++i)
1949  {
1950  G4VFacet &facet = *fFacets[i];
1951  G4double area = facet.GetArea();
1952  G4ThreeVector unit_normal = facet.GetSurfaceNormal();
1953  fCubicVolume += area * (facet.GetVertex(0).dot(unit_normal));
1954  }
1955  fCubicVolume /= 3.;
1956  return fCubicVolume;
1957 }
1958 
1960 //
1962 {
1963  if (fSurfaceArea != 0.) return fSurfaceArea;
1964 
1965  G4int size = fFacets.size();
1966  for (G4int i = 0; i < size; ++i)
1967  {
1968  G4VFacet &facet = *fFacets[i];
1969  fSurfaceArea += facet.GetArea();
1970  }
1971  return fSurfaceArea;
1972 }
1973 
1975 //
1977 {
1978  // Select randomly a facet and return a random point on it
1979 
1980  G4int i = (G4int) G4RandFlat::shoot(0., fFacets.size());
1981  return fFacets[i]->GetPointOnFace();
1982 }
1983 
1985 //
1986 // SetRandomVectorSet
1987 //
1988 // This is a set of predefined random vectors (if that isn't a contradition
1989 // in terms!) used to generate rays from a user-defined point. The member
1990 // function Inside uses these to determine whether the point is inside or
1991 // outside of the tessellated solid. All vectors should be unit vectors.
1992 //
1994 {
1995  fRandir.resize(20);
1996  fRandir[0] =
1997  G4ThreeVector(-0.9577428892113370, 0.2732676269591740, 0.0897405271949221);
1998  fRandir[1] =
1999  G4ThreeVector(-0.8331264504940770,-0.5162067214954600,-0.1985722492445700);
2000  fRandir[2] =
2001  G4ThreeVector(-0.1516671651108820, 0.9666292616127460, 0.2064580868390110);
2002  fRandir[3] =
2003  G4ThreeVector( 0.6570250350323190,-0.6944539025883300, 0.2933460081893360);
2004  fRandir[4] =
2005  G4ThreeVector(-0.4820456281280320,-0.6331060000098690,-0.6056474264406270);
2006  fRandir[5] =
2007  G4ThreeVector( 0.7629032554236800 , 0.1016854697539910,-0.6384658864065180);
2008  fRandir[6] =
2009  G4ThreeVector( 0.7689540409061150, 0.5034929891988220, 0.3939600142169160);
2010  fRandir[7] =
2011  G4ThreeVector( 0.5765188359255740, 0.5997271636278330,-0.5549354566343150);
2012  fRandir[8] =
2013  G4ThreeVector( 0.6660632777862070,-0.6362809868288380, 0.3892379937580790);
2014  fRandir[9] =
2015  G4ThreeVector( 0.3824415020414780, 0.6541792713761380,-0.6525243125110690);
2016  fRandir[10] =
2017  G4ThreeVector(-0.5107726564526760, 0.6020905056811610, 0.6136760679616570);
2018  fRandir[11] =
2019  G4ThreeVector( 0.7459135439578050, 0.6618796061649330, 0.0743530220183488);
2020  fRandir[12] =
2021  G4ThreeVector( 0.1536405855311580, 0.8117477913978260,-0.5634359711967240);
2022  fRandir[13] =
2023  G4ThreeVector( 0.0744395301705579,-0.8707110101772920,-0.4861286795736560);
2024  fRandir[14] =
2025  G4ThreeVector(-0.1665874645185400, 0.6018553940549240,-0.7810369397872780);
2026  fRandir[15] =
2027  G4ThreeVector( 0.7766902003633100, 0.6014617505959970,-0.1870724331097450);
2028  fRandir[16] =
2029  G4ThreeVector(-0.8710128685847430,-0.1434320216603030,-0.4698551243971010);
2030  fRandir[17] =
2031  G4ThreeVector( 0.8901082092766820,-0.4388411398893870, 0.1229871120030100);
2032  fRandir[18] =
2033  G4ThreeVector(-0.6430417431544370,-0.3295938228697690, 0.6912779675984150);
2034  fRandir[19] =
2035  G4ThreeVector( 0.6331124368380410, 0.6306211461665000, 0.4488714875425340);
2036 
2037  fMaxTries = 20;
2038 }
2039 
2041 //
2043 {
2044  G4int base = sizeof(*this);
2045  base += fVertexList.capacity() * sizeof(G4ThreeVector);
2046  base += fRandir.capacity() * sizeof(G4ThreeVector);
2047 
2048  G4int limit = fFacets.size();
2049  for (G4int i = 0; i < limit; ++i)
2050  {
2051  G4VFacet& facet = *fFacets[i];
2052  base += facet.AllocatedMemory();
2053  }
2054 
2055  for (auto it = fExtremeFacets.cbegin(); it != fExtremeFacets.cend(); ++it)
2056  {
2057  G4VFacet &facet = *(*it);
2058  base += facet.AllocatedMemory();
2059  }
2060  return base;
2061 }
2062 
2064 //
2066 {
2068  G4int sizeInsides = fInsides.GetNbytes();
2069  G4int sizeVoxels = fVoxels.AllocatedMemory();
2070  size += sizeInsides + sizeVoxels;
2071  return size;
2072 }
2073 
2074 #endif