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