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