Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4MultiUnion.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. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// Implementation of G4MultiUnion class
27//
28// 19.10.12 M.Gayer - Original implementation from USolids module
29// 06.04.17 G.Cosmo - Adapted implementation in Geant4 for VecGeom migration
30// --------------------------------------------------------------------
31
32#include <iostream>
33#include <sstream>
34
35#include "G4MultiUnion.hh"
36#include "Randomize.hh"
38#include "G4BoundingEnvelope.hh"
39#include "G4AffineTransform.hh"
40#include "G4DisplacedSolid.hh"
41
42#include "G4VGraphicsScene.hh"
43#include "G4Polyhedron.hh"
45
46#include "G4AutoLock.hh"
47
48namespace
49{
50 G4Mutex polyhedronMutex = G4MUTEX_INITIALIZER;
51}
52
53//______________________________________________________________________________
55 : G4VSolid(name)
56{
57 SetName(name);
58 fSolids.clear();
59 fTransformObjs.clear();
61}
62
63//______________________________________________________________________________
65{
66}
67
68//______________________________________________________________________________
70{
71 fSolids.push_back(&solid);
72 fTransformObjs.push_back(trans); // Store a local copy of transformations
73}
74
75//______________________________________________________________________________
77{
78 fSolids.push_back(solid);
79 fTransformObjs.push_back(trans); // Store a local copy of transformations
80}
81
82//______________________________________________________________________________
84{
85 return new G4MultiUnion(*this);
86}
87
88// Copy constructor
89//______________________________________________________________________________
91 : G4VSolid(rhs), fCubicVolume(rhs.fCubicVolume),
92 fSurfaceArea(rhs.fSurfaceArea),
93 kRadTolerance(rhs.kRadTolerance), fAccurate(rhs.fAccurate)
94{
95}
96
97// Fake default constructor for persistency
98//______________________________________________________________________________
100 : G4VSolid(a)
101{
102}
103
104// Assignment operator
105//______________________________________________________________________________
107{
108 // Check assignment to self
109 //
110 if (this == &rhs)
111 {
112 return *this;
113 }
114
115 // Copy base class data
116 //
118
119 return *this;
120}
121
122//______________________________________________________________________________
124{
125 // Computes the cubic volume of the "G4MultiUnion" structure using
126 // random points
127
128 if (!fCubicVolume)
129 {
130 G4ThreeVector extentMin, extentMax, d, p, point;
131 G4int inside = 0, generated;
132 BoundingLimits(extentMin, extentMax);
133 d = (extentMax - extentMin) / 2.;
134 p = (extentMax + extentMin) / 2.;
135 G4ThreeVector left = p - d;
136 G4ThreeVector length = d * 2;
137 for (generated = 0; generated < 10000; ++generated)
138 {
140 point = left + G4ThreeVector(length.x()*rvec.x(),
141 length.y()*rvec.y(),
142 length.z()*rvec.z());
143 if (Inside(point) != EInside::kOutside) ++inside;
144 }
145 G4double vbox = (2 * d.x()) * (2 * d.y()) * (2 * d.z());
146 fCubicVolume = inside * vbox / generated;
147 }
148 return fCubicVolume;
149}
150
151//______________________________________________________________________________
154 const G4ThreeVector& aDirection) const
155{
156 G4ThreeVector direction = aDirection.unit();
157 G4ThreeVector localPoint, localDirection;
158 G4double minDistance = kInfinity;
159
160 std::size_t numNodes = fSolids.size();
161 for (std::size_t i = 0 ; i < numNodes ; ++i)
162 {
163 G4VSolid& solid = *fSolids[i];
164 const G4Transform3D& transform = fTransformObjs[i];
165
166 localPoint = GetLocalPoint(transform, aPoint);
167 localDirection = GetLocalVector(transform, direction);
168
169 G4double distance = solid.DistanceToIn(localPoint, localDirection);
170 if (minDistance > distance) minDistance = distance;
171 }
172 return minDistance;
173}
174
175//______________________________________________________________________________
176G4double G4MultiUnion::DistanceToInCandidates(const G4ThreeVector& aPoint,
177 const G4ThreeVector& direction,
178 std::vector<G4int>& candidates,
179 G4SurfBits& bits) const
180{
181 std::size_t candidatesCount = candidates.size();
182 G4ThreeVector localPoint, localDirection;
183
184 G4double minDistance = kInfinity;
185 for (std::size_t i = 0 ; i < candidatesCount; ++i)
186 {
187 G4int candidate = candidates[i];
188 G4VSolid& solid = *fSolids[candidate];
189 const G4Transform3D& transform = fTransformObjs[candidate];
190
191 localPoint = GetLocalPoint(transform, aPoint);
192 localDirection = GetLocalVector(transform, direction);
193 G4double distance = solid.DistanceToIn(localPoint, localDirection);
194 if (minDistance > distance) minDistance = distance;
195 bits.SetBitNumber(candidate);
196 if (minDistance == 0) break;
197 }
198 return minDistance;
199}
200
201// Algorithm note: we have to look also for all other objects in next voxels,
202// if the distance is not shorter ... we have to do it because,
203// for example for objects which starts in first voxel in which they
204// do not collide with direction line, but in second it collides...
205// The idea of crossing voxels would be still applicable,
206// because this way we could exclude from the testing such solids,
207// which were found that obviously are not good candidates, because
208// they would return infinity
209// But if distance is smaller than the shift to next voxel, we can return
210// it immediately
211//______________________________________________________________________________
213 const G4ThreeVector& aDirection) const
214{
215 G4double minDistance = kInfinity;
216 G4ThreeVector direction = aDirection.unit();
217 G4double shift = fVoxels.DistanceToFirst(aPoint, direction);
218 if (shift == kInfinity) return shift;
219
220 G4ThreeVector currentPoint = aPoint;
221 if (shift) currentPoint += direction * shift;
222
223 G4SurfBits exclusion(fVoxels.GetBitsPerSlice());
224 std::vector<G4int> candidates, curVoxel(3);
225 fVoxels.GetVoxel(curVoxel, currentPoint);
226
227 do
228 {
229 {
230 if (fVoxels.GetCandidatesVoxelArray(curVoxel, candidates, &exclusion))
231 {
232 G4double distance = DistanceToInCandidates(aPoint, direction,
233 candidates, exclusion);
234 if (minDistance > distance) minDistance = distance;
235 if (distance < shift) break;
236 }
237 }
238 shift = fVoxels.DistanceToNext(aPoint, direction, curVoxel);
239 }
240 while (minDistance > shift);
241
242 return minDistance;
243}
244
245//______________________________________________________________________________
247 const G4ThreeVector& aDirection,
248 G4ThreeVector* aNormal) const
249{
250 // Computes distance from a point presumably outside the solid to the solid
251 // surface. Ignores first surface if the point is actually inside.
252 // Early return infinity in case the safety to any surface is found greater
253 // than the proposed step aPstep.
254 // The normal vector to the crossed surface is filled only in case the box
255 // is crossed, otherwise aNormal->IsNull() is true.
256
257 // algorithm:
258 G4ThreeVector direction = aDirection.unit();
259 G4ThreeVector localPoint, localDirection;
260 G4int ignoredSolid = -1;
261 G4double resultDistToOut = 0;
262 G4ThreeVector currentPoint = aPoint;
263
264 G4int numNodes = (G4int)fSolids.size();
265 for (auto i = 0; i < numNodes; ++i)
266 {
267 if (i != ignoredSolid)
268 {
269 G4VSolid& solid = *fSolids[i];
270 const G4Transform3D& transform = fTransformObjs[i];
271 localPoint = GetLocalPoint(transform, currentPoint);
272 localDirection = GetLocalVector(transform, direction);
273 EInside location = solid.Inside(localPoint);
274 if (location != EInside::kOutside)
275 {
276 G4double distance = solid.DistanceToOut(localPoint, localDirection,
277 aNormal);
278 if (distance < kInfinity)
279 {
280 if (resultDistToOut == kInfinity) resultDistToOut = 0;
281 if (distance > 0)
282 {
283 currentPoint = GetGlobalPoint(transform, localPoint
284 + distance*localDirection);
285 resultDistToOut += distance;
286 ignoredSolid = i; // skip the solid which we have just left
287 i = -1; // force the loop to continue from 0
288 }
289 }
290 }
291 }
292 }
293 return resultDistToOut;
294}
295
296//______________________________________________________________________________
298 const G4ThreeVector& aDirection,
299 const G4bool /* calcNorm */,
300 G4bool* /* validNorm */,
301 G4ThreeVector* aNormal) const
302{
303 return DistanceToOutVoxels(aPoint, aDirection, aNormal);
304}
305
306//______________________________________________________________________________
308 const G4ThreeVector& aDirection,
309 G4ThreeVector* aNormal) const
310{
311 // Computes distance from a point presumably inside the solid to the solid
312 // surface. Ignores first surface along each axis systematically (for points
313 // inside or outside. Early returns zero in case the second surface is behind
314 // the starting point.
315 // o The proposed step is ignored.
316 // o The normal vector to the crossed surface is always filled.
317
318 // In the case the considered point is located inside the G4MultiUnion
319 // structure, the treatments are as follows:
320 // - investigation of the candidates for the passed point
321 // - progressive moving of the point towards the surface, along the
322 // passed direction
323 // - processing of the normal
324
325 G4ThreeVector direction = aDirection.unit();
326 std::vector<G4int> candidates;
327 G4double distance = 0;
328 std::size_t numNodes = 2*fSolids.size();
329 std::size_t count=0;
330
331 if (fVoxels.GetCandidatesVoxelArray(aPoint, candidates))
332 {
333 // For normal case for which we presume the point is inside
334 G4ThreeVector localPoint, localDirection, localNormal;
335 G4ThreeVector currentPoint = aPoint;
336 G4SurfBits exclusion(fVoxels.GetBitsPerSlice());
337 G4bool notOutside;
338 G4ThreeVector maxNormal;
339
340 do
341 {
342 notOutside = false;
343
344 G4double maxDistance = -kInfinity;
345 G4int maxCandidate = 0;
346 G4ThreeVector maxLocalPoint;
347
348 std::size_t limit = candidates.size();
349 for (std::size_t i = 0 ; i < limit ; ++i)
350 {
351 G4int candidate = candidates[i];
352 // ignore the current component (that you just got out of) since
353 // numerically the propagated point will be on its surface
354
355 G4VSolid& solid = *fSolids[candidate];
356 const G4Transform3D& transform = fTransformObjs[candidate];
357
358 // The coordinates of the point are modified so as to fit the
359 // intrinsic solid local frame:
360 localPoint = GetLocalPoint(transform, currentPoint);
361
362 // DistanceToOut at least for Trd sometimes return non-zero value
363 // even from points that are outside. Therefore, this condition
364 // must currently be here, otherwise it would not work.
365 // But it means it would be slower.
366
367 if (solid.Inside(localPoint) != EInside::kOutside)
368 {
369 notOutside = true;
370
371 localDirection = GetLocalVector(transform, direction);
372
373 // propagate with solid.DistanceToOut
374 G4double shift = solid.DistanceToOut(localPoint, localDirection,
375 false, 0, &localNormal);
376 if (maxDistance < shift)
377 {
378 maxDistance = shift;
379 maxCandidate = candidate;
380 maxNormal = localNormal;
381 }
382 }
383 }
384
385 if (notOutside)
386 {
387 const G4Transform3D& transform = fTransformObjs[maxCandidate];
388
389 // convert from local normal
390 if (aNormal) *aNormal = GetGlobalVector(transform, maxNormal);
391
392 distance += maxDistance;
393 currentPoint += maxDistance * direction;
394 if(maxDistance == 0.) ++count;
395
396 // the current component will be ignored
397 exclusion.SetBitNumber(maxCandidate);
398 EInside location = InsideWithExclusion(currentPoint, &exclusion);
399
400 // perform a Inside
401 // it should be excluded current solid from checking
402 // we have to collect the maximum distance from all given candidates.
403 // such "maximum" candidate should be then used for finding next
404 // candidates
405 if (location == EInside::kOutside)
406 {
407 // else return cumulated distances to outside of the traversed
408 // components
409 break;
410 }
411 // if inside another component, redo 1 to 3 but add the next
412 // DistanceToOut on top of the previous.
413
414 // and fill the candidates for the corresponding voxel (just
415 // exiting current component along direction)
416 candidates.clear();
417
418 fVoxels.GetCandidatesVoxelArray(currentPoint, candidates, &exclusion);
419 exclusion.ResetBitNumber(maxCandidate);
420 }
421 }
422 while ((notOutside) && (count < numNodes));
423 }
424
425 return distance;
426}
427
428//______________________________________________________________________________
429EInside G4MultiUnion::InsideWithExclusion(const G4ThreeVector& aPoint,
430 G4SurfBits* exclusion) const
431{
432 // Classify point location with respect to solid:
433 // o eInside - inside the solid
434 // o eSurface - close to surface within tolerance
435 // o eOutside - outside the solid
436
437 // Hitherto, it is considered that only parallelepipedic nodes
438 // can be added to the container
439
440 // Implementation using voxelisation techniques:
441 // ---------------------------------------------
442
443 G4ThreeVector localPoint;
444 EInside location = EInside::kOutside;
445
446 std::vector<G4int> candidates;
447 std::vector<G4MultiUnionSurface> surfaces;
448
449 // TODO: test if it works well and if so measure performance
450 // TODO: getPointIndex should not be used, instead GetVoxel + GetVoxelsIndex
451 // should be used
452 // TODO: than pass result to GetVoxel further to GetCandidatesVoxelArray
453 // TODO: eventually GetVoxel should be inlined here, early exit if any
454 // binary search is -1
455
456 G4int limit = fVoxels.GetCandidatesVoxelArray(aPoint, candidates, exclusion);
457 for (G4int i = 0 ; i < limit ; ++i)
458 {
459 G4int candidate = candidates[i];
460 G4VSolid& solid = *fSolids[candidate];
461 const G4Transform3D& transform = fTransformObjs[candidate];
462
463 // The coordinates of the point are modified so as to fit the intrinsic
464 // solid local frame:
465 localPoint = GetLocalPoint(transform, aPoint);
466 location = solid.Inside(localPoint);
467 if (location == EInside::kInside) return EInside::kInside;
468 else if (location == EInside::kSurface)
469 {
470 G4MultiUnionSurface surface;
471 surface.point = localPoint;
472 surface.solid = &solid;
473 surfaces.push_back(surface);
474 }
475 }
476
477 ///////////////////////////////////////////////////////////////////////////
478 // Important comment: When two solids touch each other along a flat
479 // surface, the surface points will be considered as kSurface, while points
480 // located around will correspond to kInside (cf. G4UnionSolid)
481
482 std::size_t size = surfaces.size();
483
484 if (size == 0)
485 {
486 return EInside::kOutside;
487 }
488
489 for (std::size_t i = 0; i < size - 1; ++i)
490 {
491 G4MultiUnionSurface& left = surfaces[i];
492 for (std::size_t j = i + 1; j < size; ++j)
493 {
494 G4MultiUnionSurface& right = surfaces[j];
495 G4ThreeVector n, n2;
496 n = left.solid->SurfaceNormal(left.point);
497 n2 = right.solid->SurfaceNormal(right.point);
498 if ((n + n2).mag2() < 1000 * kRadTolerance)
499 return EInside::kInside;
500 }
501 }
502
503 return EInside::kSurface;
504}
505
506//______________________________________________________________________________
508{
509 // Classify point location with respect to solid:
510 // o eInside - inside the solid
511 // o eSurface - close to surface within tolerance
512 // o eOutside - outside the solid
513
514 // Hitherto, it is considered that only parallelepipedic nodes can be
515 // added to the container
516
517 // Implementation using voxelisation techniques:
518 // ---------------------------------------------
519
520 // return InsideIterator(aPoint);
521
522 EInside location = InsideWithExclusion(aPoint);
523 return location;
524}
525
526//______________________________________________________________________________
528{
529 G4ThreeVector localPoint;
530 EInside location = EInside::kOutside;
531 G4int countSurface = 0;
532
533 G4int numNodes = (G4int)fSolids.size();
534 for (auto i = 0 ; i < numNodes ; ++i)
535 {
536 G4VSolid& solid = *fSolids[i];
537 G4Transform3D transform = GetTransformation(i);
538
539 // The coordinates of the point are modified so as to fit the
540 // intrinsic solid local frame:
541 localPoint = GetLocalPoint(transform, aPoint);
542
543 location = solid.Inside(localPoint);
544
545 if (location == EInside::kSurface)
546 ++countSurface;
547
548 if (location == EInside::kInside) return EInside::kInside;
549 }
550 if (countSurface != 0) return EInside::kSurface;
551 return EInside::kOutside;
552}
553
554//______________________________________________________________________________
555void G4MultiUnion::Extent(EAxis aAxis, G4double& aMin, G4double& aMax) const
556{
557 // Determines the bounding box for the considered instance of "UMultipleUnion"
558 G4ThreeVector min, max;
559
560 G4int numNodes = (G4int)fSolids.size();
561 for (auto i = 0 ; i < numNodes ; ++i)
562 {
563 G4VSolid& solid = *fSolids[i];
564 G4Transform3D transform = GetTransformation(i);
565 solid.BoundingLimits(min, max);
566
567 TransformLimits(min, max, transform);
568
569 if (i == 0)
570 {
571 switch (aAxis)
572 {
573 case kXAxis:
574 aMin = min.x();
575 aMax = max.x();
576 break;
577 case kYAxis:
578 aMin = min.y();
579 aMax = max.y();
580 break;
581 case kZAxis:
582 aMin = min.z();
583 aMax = max.z();
584 break;
585 default:
586 break;
587 }
588 }
589 else
590 {
591 // Determine the min/max on the considered axis:
592 switch (aAxis)
593 {
594 case kXAxis:
595 if (min.x() < aMin)
596 aMin = min.x();
597 if (max.x() > aMax)
598 aMax = max.x();
599 break;
600 case kYAxis:
601 if (min.y() < aMin)
602 aMin = min.y();
603 if (max.y() > aMax)
604 aMax = max.y();
605 break;
606 case kZAxis:
607 if (min.z() < aMin)
608 aMin = min.z();
609 if (max.z() > aMax)
610 aMax = max.z();
611 break;
612 default:
613 break;
614 }
615 }
616 }
617}
618
619//______________________________________________________________________________
621 G4ThreeVector& aMax) const
622{
623 Extent(kXAxis, aMin[0], aMax[0]);
624 Extent(kYAxis, aMin[1], aMax[1]);
625 Extent(kZAxis, aMin[2], aMax[2]);
626}
627
628//______________________________________________________________________________
629G4bool
631 const G4VoxelLimits& pVoxelLimit,
632 const G4AffineTransform& pTransform,
633 G4double& pMin, G4double& pMax) const
634{
635 G4ThreeVector bmin, bmax;
636
637 // Get bounding box
638 BoundingLimits(bmin,bmax);
639
640 // Find extent
641 G4BoundingEnvelope bbox(bmin,bmax);
642 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
643}
644
645//______________________________________________________________________________
647{
648 // Computes the localNormal on a surface and returns it as a unit vector.
649 // Must return a valid vector. (even if the point is not on the surface).
650 //
651 // On an edge or corner, provide an average localNormal of all facets within
652 // tolerance
653 // NOTE: the tolerance value used in here is not yet the global surface
654 // tolerance - we will have to revise this value - TODO
655
656 std::vector<G4int> candidates;
657 G4ThreeVector localPoint, normal, localNormal;
658 G4double safety = kInfinity;
659 G4int node = 0;
660
661 ///////////////////////////////////////////////////////////////////////////
662 // Important comment: Cases for which the point is located on an edge or
663 // on a vertice remain to be treated
664
665 // determine weather we are in voxel area
666 if (fVoxels.GetCandidatesVoxelArray(aPoint, candidates))
667 {
668 std::size_t limit = candidates.size();
669 for (std::size_t i = 0 ; i < limit ; ++i)
670 {
671 G4int candidate = candidates[i];
672 const G4Transform3D& transform = fTransformObjs[candidate];
673
674 // The coordinates of the point are modified so as to fit the intrinsic
675 // solid local frame:
676 localPoint = GetLocalPoint(transform, aPoint);
677 G4VSolid& solid = *fSolids[candidate];
678 EInside location = solid.Inside(localPoint);
679
680 if (location == EInside::kSurface)
681 {
682 // normal case when point is on surface, we pick first solid
683 normal = GetGlobalVector(transform, solid.SurfaceNormal(localPoint));
684 return normal.unit();
685 }
686 else
687 {
688 // collect the smallest safety and remember solid node
689 G4double s = (location == EInside::kInside)
690 ? solid.DistanceToOut(localPoint)
691 : solid.DistanceToIn(localPoint);
692 if (s < safety)
693 {
694 safety = s;
695 node = candidate;
696 }
697 }
698 }
699 // on none of the solids, the point was not on the surface
700 G4VSolid& solid = *fSolids[node];
701 const G4Transform3D& transform = fTransformObjs[node];
702 localPoint = GetLocalPoint(transform, aPoint);
703
704 normal = GetGlobalVector(transform, solid.SurfaceNormal(localPoint));
705 return normal.unit();
706 }
707 else
708 {
709 // for the case when point is certainly outside:
710
711 // find a solid in union with the smallest safety
712 node = SafetyFromOutsideNumberNode(aPoint, safety);
713 G4VSolid& solid = *fSolids[node];
714
715 const G4Transform3D& transform = fTransformObjs[node];
716 localPoint = GetLocalPoint(transform, aPoint);
717
718 // evaluate normal for point at this found solid
719 // and transform multi-union coordinates
720 normal = GetGlobalVector(transform, solid.SurfaceNormal(localPoint));
721
722 return normal.unit();
723 }
724}
725
726//______________________________________________________________________________
728{
729 // Estimates isotropic distance to the surface of the solid. This must
730 // be either accurate or an underestimate.
731 // Two modes: - default/fast mode, sacrificing accuracy for speed
732 // - "precise" mode, requests accurate value if available.
733
734 std::vector<G4int> candidates;
735 G4ThreeVector localPoint;
736 G4double safetyMin = kInfinity;
737
738 // In general, the value return by DistanceToIn(p) will not be the exact
739 // but only an undervalue (cf. overlaps)
740 fVoxels.GetCandidatesVoxelArray(point, candidates);
741
742 std::size_t limit = candidates.size();
743 for (std::size_t i = 0; i < limit; ++i)
744 {
745 G4int candidate = candidates[i];
746
747 // The coordinates of the point are modified so as to fit the intrinsic
748 // solid local frame:
749 const G4Transform3D& transform = fTransformObjs[candidate];
750 localPoint = GetLocalPoint(transform, point);
751 G4VSolid& solid = *fSolids[candidate];
752 if (solid.Inside(localPoint) == EInside::kInside)
753 {
754 G4double safety = solid.DistanceToOut(localPoint);
755 if (safetyMin > safety) safetyMin = safety;
756 }
757 }
758 if (safetyMin == kInfinity) safetyMin = 0; // we are not inside
759
760 return safetyMin;
761}
762
763//______________________________________________________________________________
765{
766 // Estimates the isotropic safety from a point outside the current solid to
767 // any of its surfaces. The algorithm may be accurate or should provide a fast
768 // underestimate.
769
770 if (!fAccurate) { return fVoxels.DistanceToBoundingBox(point); }
771
772 const std::vector<G4VoxelBox>& boxes = fVoxels.GetBoxes();
773 G4double safetyMin = kInfinity;
774 G4ThreeVector localPoint;
775
776 std::size_t numNodes = fSolids.size();
777 for (std::size_t j = 0; j < numNodes; ++j)
778 {
779 G4ThreeVector dxyz;
780 if (j > 0)
781 {
782 const G4ThreeVector& pos = boxes[j].pos;
783 const G4ThreeVector& hlen = boxes[j].hlen;
784 for (auto i = 0; i <= 2; ++i)
785 // distance to middle point - hlength => distance from point to border
786 // of x,y,z
787 if ((dxyz[i] = std::abs(point[i] - pos[i]) - hlen[i]) > safetyMin)
788 continue;
789
790 G4double d2xyz = 0.;
791 for (auto i = 0; i <= 2; ++i)
792 if (dxyz[i] > 0) d2xyz += dxyz[i] * dxyz[i];
793
794 // minimal distance is at least this, but could be even higher. therefore,
795 // we can stop if previous was already lower, let us check if it does any
796 // chance to be better tha previous values...
797 if (d2xyz >= safetyMin * safetyMin)
798 {
799 continue;
800 }
801 }
802 const G4Transform3D& transform = fTransformObjs[j];
803 localPoint = GetLocalPoint(transform, point);
804 G4VSolid& solid = *fSolids[j];
805
806 G4double safety = solid.DistanceToIn(localPoint);
807 if (safety <= 0) return safety;
808 // it was detected, that the point is not located outside
809 if (safetyMin > safety) safetyMin = safety;
810 }
811 return safetyMin;
812}
813
814//______________________________________________________________________________
816{
817 if (!fSurfaceArea)
818 {
819 fSurfaceArea = EstimateSurfaceArea(1000000, 0.001);
820 }
821 return fSurfaceArea;
822}
823
824//______________________________________________________________________________
826{
827 fVoxels.Voxelize(fSolids, fTransformObjs);
828}
829
830//______________________________________________________________________________
831G4int G4MultiUnion::SafetyFromOutsideNumberNode(const G4ThreeVector& aPoint,
832 G4double& safetyMin) const
833{
834 // Method returning the closest node from a point located outside a
835 // G4MultiUnion.
836 // This is used to compute the normal in the case no candidate has been found.
837
838 const std::vector<G4VoxelBox>& boxes = fVoxels.GetBoxes();
839 safetyMin = kInfinity;
840 std::size_t safetyNode = 0;
841 G4ThreeVector localPoint;
842
843 std::size_t numNodes = fSolids.size();
844 for (std::size_t i = 0; i < numNodes; ++i)
845 {
846 G4double d2xyz = 0.;
847 G4double dxyz0 = std::abs(aPoint.x() - boxes[i].pos.x()) - boxes[i].hlen.x();
848 if (dxyz0 > safetyMin) continue;
849 G4double dxyz1 = std::abs(aPoint.y() - boxes[i].pos.y()) - boxes[i].hlen.y();
850 if (dxyz1 > safetyMin) continue;
851 G4double dxyz2 = std::abs(aPoint.z() - boxes[i].pos.z()) - boxes[i].hlen.z();
852 if (dxyz2 > safetyMin) continue;
853
854 if (dxyz0 > 0) d2xyz += dxyz0 * dxyz0;
855 if (dxyz1 > 0) d2xyz += dxyz1 * dxyz1;
856 if (dxyz2 > 0) d2xyz += dxyz2 * dxyz2;
857 if (d2xyz >= safetyMin * safetyMin) continue;
858
859 G4VSolid& solid = *fSolids[i];
860 const G4Transform3D& transform = fTransformObjs[i];
861 localPoint = GetLocalPoint(transform, aPoint);
862 fAccurate = true;
863 G4double safety = solid.DistanceToIn(localPoint);
864 fAccurate = false;
865 if (safetyMin > safety)
866 {
867 safetyMin = safety;
868 safetyNode = i;
869 }
870 }
871 return (G4int)safetyNode;
872}
873
874//______________________________________________________________________________
875void G4MultiUnion::TransformLimits(G4ThreeVector& min, G4ThreeVector& max,
876 const G4Transform3D& transformation) const
877{
878 // The goal of this method is to convert the quantities min and max
879 // (representing the bounding box of a given solid in its local frame)
880 // to the main frame, using "transformation"
881
882 G4ThreeVector vertices[8] = // Detemination of the vertices thanks to
883 { // the extension of each solid:
884 G4ThreeVector(min.x(), min.y(), min.z()), // 1st vertice:
885 G4ThreeVector(min.x(), max.y(), min.z()), // 2nd vertice:
886 G4ThreeVector(max.x(), max.y(), min.z()),
887 G4ThreeVector(max.x(), min.y(), min.z()),
888 G4ThreeVector(min.x(), min.y(), max.z()),
889 G4ThreeVector(min.x(), max.y(), max.z()),
890 G4ThreeVector(max.x(), max.y(), max.z()),
891 G4ThreeVector(max.x(), min.y(), max.z())
892 };
893
894 min.set(kInfinity,kInfinity,kInfinity);
895 max.set(-kInfinity,-kInfinity,-kInfinity);
896
897 // Loop on th vertices
898 G4int limit = sizeof(vertices) / sizeof(G4ThreeVector);
899 for (G4int i = 0 ; i < limit; ++i)
900 {
901 // From local frame to the global one:
902 // Current positions on the three axis:
903 G4ThreeVector current = GetGlobalPoint(transformation, vertices[i]);
904
905 // If need be, replacement of the min & max values:
906 if (current.x() > max.x()) max.setX(current.x());
907 if (current.x() < min.x()) min.setX(current.x());
908
909 if (current.y() > max.y()) max.setY(current.y());
910 if (current.y() < min.y()) min.setY(current.y());
911
912 if (current.z() > max.z()) max.setZ(current.z());
913 if (current.z() < min.z()) min.setZ(current.z());
914 }
915}
916
917// Stream object contents to an output stream
918//______________________________________________________________________________
919std::ostream& G4MultiUnion::StreamInfo(std::ostream& os) const
920{
921 G4long oldprc = os.precision(16);
922 os << "-----------------------------------------------------------\n"
923 << " *** Dump for solid - " << GetName() << " ***\n"
924 << " ===================================================\n"
925 << " Solid type: G4MultiUnion\n"
926 << " Parameters: \n";
927 std::size_t numNodes = fSolids.size();
928 for (std::size_t i = 0 ; i < numNodes ; ++i)
929 {
930 G4VSolid& solid = *fSolids[i];
931 solid.StreamInfo(os);
932 const G4Transform3D& transform = fTransformObjs[i];
933 os << " Translation is " << transform.getTranslation() << " \n";
934 os << " Rotation is :" << " \n";
935 os << " " << transform.getRotation() << "\n";
936 }
937 os << " \n"
938 << "-----------------------------------------------------------\n";
939 os.precision(oldprc);
940
941 return os;
942}
943
944//______________________________________________________________________________
946{
947 G4ThreeVector point;
948
949 G4long size = fSolids.size();
950
951 do
952 {
953 G4long rnd = G4RandFlat::shootInt(G4long(0), size);
954 G4VSolid& solid = *fSolids[rnd];
955 point = solid.GetPointOnSurface();
956 const G4Transform3D& transform = fTransformObjs[rnd];
957 point = GetGlobalPoint(transform, point);
958 }
959 while (Inside(point) != EInside::kSurface);
960
961 return point;
962}
963
964//______________________________________________________________________________
965void
967{
968 scene.AddSolid (*this);
969}
970
971//______________________________________________________________________________
973{
974 HepPolyhedronProcessor processor;
976
977 G4VSolid* solidA = GetSolid(0);
978 const G4Transform3D transform0=GetTransformation(0);
979 G4DisplacedSolid dispSolidA("placedA",solidA,transform0);
980
981 G4Polyhedron* top = new G4Polyhedron(*dispSolidA.GetPolyhedron());
982
983 for(G4int i=1; i<GetNumberOfSolids(); ++i)
984 {
985 G4VSolid* solidB = GetSolid(i);
986 const G4Transform3D transform=GetTransformation(i);
987 G4DisplacedSolid dispSolidB("placedB",solidB,transform);
988 G4Polyhedron* operand = dispSolidB.GetPolyhedron();
989 processor.push_back (operation, *operand);
990 }
991
992 if (processor.execute(*top)) { return top; }
993 else { return 0; }
994}
995
996//______________________________________________________________________________
998{
999 if (fpPolyhedron == nullptr ||
1000 fRebuildPolyhedron ||
1002 fpPolyhedron->GetNumberOfRotationSteps())
1003 {
1004 G4AutoLock l(&polyhedronMutex);
1005 delete fpPolyhedron;
1006 fpPolyhedron = CreatePolyhedron();
1007 fRebuildPolyhedron = false;
1008 l.unlock();
1009 }
1010 return fpPolyhedron;
1011}
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
std::mutex G4Mutex
Definition: G4Threading.hh:81
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
long G4long
Definition: G4Types.hh:87
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
Hep3Vector unit() const
double x() const
double y() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4Polyhedron * GetPolyhedron() const
G4double GetRadialTolerance() const
static G4GeometryTolerance * GetInstance()
const G4Transform3D & GetTransformation(G4int index) const
G4int GetNumberOfSolids() const
G4double DistanceToIn(const G4ThreeVector &aPoint) const
G4double GetSurfaceArea()
G4double DistanceToOutNoVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const
G4double GetCubicVolume()
EInside InsideNoVoxels(const G4ThreeVector &aPoint) const
G4double DistanceToInNoVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection) const
EInside Inside(const G4ThreeVector &aPoint) const
G4double DistanceToOut(const G4ThreeVector &aPoint) const
G4Polyhedron * CreatePolyhedron() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
void BoundingLimits(G4ThreeVector &aMin, G4ThreeVector &aMax) const
G4Polyhedron * GetPolyhedron() const
G4double DistanceToOutVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const
std::ostream & StreamInfo(std::ostream &os) const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
void AddNode(G4VSolid &solid, const G4Transform3D &trans)
Definition: G4MultiUnion.cc:69
void Extent(EAxis aAxis, G4double &aMin, G4double &aMax) const
G4VSolid * Clone() const
Definition: G4MultiUnion.cc:83
G4ThreeVector GetPointOnSurface() const
G4VSolid * GetSolid(G4int index) const
G4ThreeVector SurfaceNormal(const G4ThreeVector &aPoint) const
G4MultiUnion & operator=(const G4MultiUnion &rhs)
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
void ResetBitNumber(unsigned int bitnumber)
Definition: G4SurfBits.hh:154
void SetBitNumber(unsigned int bitnumber, G4bool value=true)
Definition: G4SurfBits.hh:114
virtual void AddSolid(const G4Box &)=0
G4double EstimateSurfaceArea(G4int nStat, G4double ell) const
Definition: G4VSolid.cc:265
G4String GetName() const
virtual std::ostream & StreamInfo(std::ostream &os) const =0
virtual EInside Inside(const G4ThreeVector &p) const =0
void SetName(const G4String &name)
Definition: G4VSolid.cc:127
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4VSolid.cc:665
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:107
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
G4double DistanceToBoundingBox(const G4ThreeVector &point) const
const std::vector< G4VoxelBox > & GetBoxes() const
void GetVoxel(std::vector< G4int > &curVoxel, const G4ThreeVector &point) const
G4int GetBitsPerSlice() const
G4double DistanceToFirst(const G4ThreeVector &point, const G4ThreeVector &direction) const
G4int GetCandidatesVoxelArray(const G4ThreeVector &point, std::vector< G4int > &list, G4SurfBits *crossed=nullptr) const
Definition: G4Voxelizer.cc:966
G4double DistanceToNext(const G4ThreeVector &point, const G4ThreeVector &direction, std::vector< G4int > &curVoxel) const
void Voxelize(std::vector< G4VSolid * > &solids, std::vector< G4Transform3D > &transforms)
Definition: G4Voxelizer.cc:706
CLHEP::HepRotation getRotation() const
CLHEP::Hep3Vector getTranslation() const
bool execute(HepPolyhedron &)
void push_back(Operation, const HepPolyhedron &)
static G4int GetNumberOfRotationSteps()
EAxis
Definition: geomdefs.hh:54
@ kYAxis
Definition: geomdefs.hh:56
@ kXAxis
Definition: geomdefs.hh:55
@ kZAxis
Definition: geomdefs.hh:57
EInside
Definition: geomdefs.hh:67
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments