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