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