Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MultiUnion Class Reference

An instance of G4MultiUnion constitutes a grouping of several solids. The constituent solids are stored with their respective location in a node instance. An instance of G4MultiUnion is subsequently composed of one or several nodes. More...

#include <G4MultiUnion.hh>

Inheritance diagram for G4MultiUnion:

Public Member Functions

 G4MultiUnion ()
 G4MultiUnion (const G4String &name)
 ~G4MultiUnion () override=default
 G4MultiUnion (__void__ &)
void AddNode (G4VSolid &solid, const G4Transform3D &trans)
void AddNode (G4VSolid *solid, const G4Transform3D &trans)
 G4MultiUnion (const G4MultiUnion &rhs)
G4MultiUnionoperator= (const G4MultiUnion &rhs)
const G4Transform3DGetTransformation (G4int index) const
G4VSolidGetSolid (G4int index) const
G4int GetNumberOfSolids () const
EInside Inside (const G4ThreeVector &aPoint) const override
G4double DistanceToIn (const G4ThreeVector &aPoint) const override
G4double DistanceToOut (const G4ThreeVector &aPoint) const override
void SetAccurateSafety (G4bool flag)
G4double DistanceToIn (const G4ThreeVector &aPoint, const G4ThreeVector &aDirection) const override
G4double DistanceToOut (const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *aNormalVector=nullptr) const override
G4double DistanceToInNoVoxels (const G4ThreeVector &aPoint, const G4ThreeVector &aDirection) const
G4double DistanceToOutVoxels (const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const
G4double DistanceToOutNoVoxels (const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const
G4ThreeVector SurfaceNormal (const G4ThreeVector &aPoint) const override
void Extent (EAxis aAxis, G4double &aMin, G4double &aMax) const
void BoundingLimits (G4ThreeVector &aMin, G4ThreeVector &aMax) const override
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const override
G4double GetCubicVolume () override
G4double GetSurfaceArea () override
G4int GetNumOfConstituents () const override
G4bool IsFaceted () const override
G4VSolidClone () const override
G4GeometryType GetEntityType () const override
void Voxelize ()
G4VoxelizerGetVoxels () const
std::ostream & StreamInfo (std::ostream &os) const override
G4ThreeVector GetPointOnSurface () const override
void DescribeYourselfTo (G4VGraphicsScene &scene) const override
G4PolyhedronCreatePolyhedron () const override
G4PolyhedronGetPolyhedron () const override
Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
virtual ~G4VSolid ()
 G4VSolid (const G4VSolid &rhs)
G4VSolidoperator= (const G4VSolid &rhs)
G4bool operator== (const G4VSolid &s) const
G4String GetName () const
void SetName (const G4String &name)
G4double GetTolerance () const
virtual void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
void DumpInfo () const
virtual G4VisExtent GetExtent () const
virtual const G4VSolidGetConstituentSolid (G4int no) const
virtual G4VSolidGetConstituentSolid (G4int no)
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 G4VSolid (__void__ &)
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
G4double EstimateSurfaceArea (G4int nStat, G4double epsilon) const

Friends

class G4Voxelizer

Additional Inherited Members

Protected Member Functions inherited from G4VSolid
void CalculateClippedPolygonExtent (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
void ClipCrossSection (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
void ClipBetweenSections (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
void ClipPolygon (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
Protected Attributes inherited from G4VSolid
G4double kCarTolerance

Detailed Description

An instance of G4MultiUnion constitutes a grouping of several solids. The constituent solids are stored with their respective location in a node instance. An instance of G4MultiUnion is subsequently composed of one or several nodes.

Definition at line 61 of file G4MultiUnion.hh.

Constructor & Destructor Documentation

◆ G4MultiUnion() [1/4]

G4MultiUnion::G4MultiUnion ( )

Empty default constructor.

Definition at line 57 of file G4MultiUnion.cc.

58 : G4VSolid("")
59{
60}
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59

Referenced by Clone(), G4MultiUnion(), operator=(), and ~G4MultiUnion().

◆ G4MultiUnion() [2/4]

G4MultiUnion::G4MultiUnion ( const G4String & name)

Constructor assigning a name and initialising components.

Definition at line 63 of file G4MultiUnion.cc.

64 : G4VSolid(name)
65{
66 SetName(name);
67 fSolids.clear();
68 fTransformObjs.clear();
70}
G4double GetRadialTolerance() const
static G4GeometryTolerance * GetInstance()
void SetName(const G4String &name)
Definition G4VSolid.cc:126

◆ ~G4MultiUnion()

G4MultiUnion::~G4MultiUnion ( )
overridedefault

Default destructor.

◆ G4MultiUnion() [3/4]

G4MultiUnion::G4MultiUnion ( __void__ & a)

Fake default constructor for usage restricted to direct object persistency for clients requiring preallocation of memory for persistifiable objects.

Definition at line 103 of file G4MultiUnion.cc.

104 : G4VSolid(a)
105{
106}

◆ G4MultiUnion() [4/4]

G4MultiUnion::G4MultiUnion ( const G4MultiUnion & rhs)

Copy constructor and assignment operator.

Definition at line 94 of file G4MultiUnion.cc.

95 : G4VSolid(rhs), fCubicVolume(rhs.fCubicVolume),
96 fSurfaceArea(rhs.fSurfaceArea),
97 kRadTolerance(rhs.kRadTolerance), fAccurate(rhs.fAccurate)
98{
99}

Member Function Documentation

◆ AddNode() [1/2]

void G4MultiUnion::AddNode ( G4VSolid & solid,
const G4Transform3D & trans )

Methods to build the multiple union by adding nodes (by pointer or ref).

Parameters
[in]solidThe solid to be added to the structure.
[in]transThe 3D transformation relative to the structure.

Definition at line 73 of file G4MultiUnion.cc.

74{
75 fSolids.push_back(&solid);
76 fTransformObjs.push_back(trans); // Store a local copy of transformations
77}

Referenced by G4tgbVolume::FindOrConstructG4Solid(), and G4GDMLReadSolids::MultiUnionNodeRead().

◆ AddNode() [2/2]

void G4MultiUnion::AddNode ( G4VSolid * solid,
const G4Transform3D & trans )

Definition at line 80 of file G4MultiUnion.cc.

81{
82 fSolids.push_back(solid);
83 fTransformObjs.push_back(trans); // Store a local copy of transformations
84}

◆ BoundingLimits()

void G4MultiUnion::BoundingLimits ( G4ThreeVector & aMin,
G4ThreeVector & aMax ) const
overridevirtual

Computes the bounding limits of the solid.

Parameters
[out]aMinThe minimum bounding limit point.
[out]aMaxThe maximum bounding limit point.

Reimplemented from G4VSolid.

Definition at line 583 of file G4MultiUnion.cc.

585{
586 Extent(kXAxis, aMin[0], aMax[0]);
587 Extent(kYAxis, aMin[1], aMax[1]);
588 Extent(kZAxis, aMin[2], aMax[2]);
589}
void Extent(EAxis aAxis, G4double &aMin, G4double &aMax) const
@ kYAxis
Definition geomdefs.hh:56
@ kXAxis
Definition geomdefs.hh:55
@ kZAxis
Definition geomdefs.hh:57

Referenced by CalculateExtent().

◆ CalculateExtent()

G4bool G4MultiUnion::CalculateExtent ( const EAxis pAxis,
const G4VoxelLimits & pVoxelLimit,
const G4AffineTransform & pTransform,
G4double & pMin,
G4double & pMax ) const
overridevirtual

Calculates the minimum and maximum extent of a solid, when under the specified transform, and within the specified limits.

Parameters
[in]pAxisThe axis along which compute the extent.
[in]pVoxelLimitThe limiting space dictated by voxels.
[in]pTransformThe internal transformation applied to the solid.
[out]pMinThe minimum extent value.
[out]pMaxThe maximum extent value.
Returns
True if the solid is intersected by the extent region.

Implements G4VSolid.

Definition at line 593 of file G4MultiUnion.cc.

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}
CLHEP::Hep3Vector G4ThreeVector
void BoundingLimits(G4ThreeVector &aMin, G4ThreeVector &aMax) const override

◆ Clone()

G4VSolid * G4MultiUnion::Clone ( ) const
overridevirtual

Returns a new allocated clone of the multi-union structure.

Reimplemented from G4VSolid.

Definition at line 87 of file G4MultiUnion.cc.

88{
89 return new G4MultiUnion(*this);
90}

◆ CreatePolyhedron()

G4Polyhedron * G4MultiUnion::CreatePolyhedron ( ) const
overridevirtual

Creates a Polyhedron used for Visualisation. It is the caller's responsibility to delete it. A null pointer means "not created".

Reimplemented from G4VSolid.

Definition at line 969 of file G4MultiUnion.cc.

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}
HepGeom::Transform3D G4Transform3D
int G4int
Definition G4Types.hh:85
static G4VBooleanProcessor * GetExternalBooleanProcessor()
const G4Transform3D & GetTransformation(G4int index) const
G4int GetNumberOfSolids() const
G4VSolid * GetSolid(G4int index) const
virtual G4PolyhedronArbitrary * Process(const G4VSolid *)
bool execute(HepPolyhedron &)
void push_back(Operation, const HepPolyhedron &)

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

void G4MultiUnion::DescribeYourselfTo ( G4VGraphicsScene & scene) const
overridevirtual

Methods for creating graphical representations (i.e. for visualisation).

Implements G4VSolid.

Definition at line 963 of file G4MultiUnion.cc.

964{
965 scene.AddSolid (*this);
966}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

G4double G4MultiUnion::DistanceToIn ( const G4ThreeVector & p) const
overridevirtual

Calculates the distance to the nearest surface of a shape from an outside point. The distance can be an underestimate.

Parameters
[in]pThe point at offset p.
Returns
The safety distance to enter the shape.

Implements G4VSolid.

Definition at line 723 of file G4MultiUnion.cc.

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}
double G4double
Definition G4Types.hh:83
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0

◆ DistanceToIn() [2/2]

G4double G4MultiUnion::DistanceToIn ( const G4ThreeVector & aPoint,
const G4ThreeVector & aDirection ) const
overridevirtual

Returns the distance along the normalised vector "aDirection" to the shape, from the point at offset "aPoint". If there is no intersection, return kInfinity. The first intersection resulting from leaving a surface/volume is discarded. Hence, it is tolerant of points on the surface of the shape.

Implements G4VSolid.

Definition at line 187 of file G4MultiUnion.cc.

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}
Hep3Vector unit() const

◆ DistanceToInNoVoxels()

G4double G4MultiUnion::DistanceToInNoVoxels ( const G4ThreeVector & aPoint,
const G4ThreeVector & aDirection ) const

Methods to compute the distance to enter/exit a volume, given point and direction, in presence of voxels-based optimisation structure or not.

Definition at line 128 of file G4MultiUnion.cc.

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}

◆ DistanceToOut() [1/2]

G4double G4MultiUnion::DistanceToOut ( const G4ThreeVector & p) const
overridevirtual

Calculates the distance to the nearest surface of a shape from an inside point 'p'. The distance can be an underestimate.

Parameters
[in]pThe point at offset p.
Returns
The safety distance to exit the shape.

Implements G4VSolid.

Definition at line 686 of file G4MultiUnion.cc.

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}
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
@ kInside
Definition geomdefs.hh:70

◆ DistanceToOut() [2/2]

G4double G4MultiUnion::DistanceToOut ( const G4ThreeVector & aPoint,
const G4ThreeVector & aDirection,
const G4bool calcNorm = false,
G4bool * validNorm = nullptr,
G4ThreeVector * aNormalVector = nullptr ) const
overridevirtual

Computes distance from a point presumably inside the solid to the solid surface. Ignores first surface along each axis systematically (for points inside or outside. Early returns zero in case the second surface is behind the starting point. The normal vector to the crossed surface is always filled. In the case the considered point is located inside the G4MultiUnion structure, it acts as follows:

  • investigation of the candidates for the passed point
  • progressive moving of the point towards the surface, along the provided direction
  • processing of the normal.
    Parameters
    [in]aPointThe reference point in space.
    [in]aDirectionThe normalised direction.
    [in]calcNormFlag unused.
    [out]validNormUnused.
    [out]aNormalVectorThe exiting outwards normal vector (undefined Magnitude).
    Returns
    The distance value to exit a volume.

Implements G4VSolid.

Definition at line 272 of file G4MultiUnion.cc.

277{
278 return DistanceToOutVoxels(aPoint, aDirection, aNormal);
279}
G4double DistanceToOutVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const

◆ DistanceToOutNoVoxels()

G4double G4MultiUnion::DistanceToOutNoVoxels ( const G4ThreeVector & aPoint,
const G4ThreeVector & aDirection,
G4ThreeVector * aNormalVector ) const

Definition at line 221 of file G4MultiUnion.cc.

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}
EInside
Definition geomdefs.hh:67
@ kOutside
Definition geomdefs.hh:68

◆ DistanceToOutVoxels()

G4double G4MultiUnion::DistanceToOutVoxels ( const G4ThreeVector & aPoint,
const G4ThreeVector & aDirection,
G4ThreeVector * aNormalVector ) const

Definition at line 282 of file G4MultiUnion.cc.

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}
bool G4bool
Definition G4Types.hh:86

Referenced by DistanceToOut().

◆ Extent()

void G4MultiUnion::Extent ( EAxis aAxis,
G4double & aMin,
G4double & aMax ) const

Determines the bounding box for the considered instance of G4MultiUnion.

Parameters
[in]aAxisThe axis along which computing the extent.
[out]aMinThe minimum bounding limit point.
[out]aMaxThe maximum bounding limit point.

Definition at line 523 of file G4MultiUnion.cc.

524{
525 // Determines the bounding box for the considered instance of G4MultiUnion
526
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}
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition G4VSolid.cc:691
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

Referenced by BoundingLimits().

◆ GetCubicVolume()

G4double G4MultiUnion::GetCubicVolume ( )
overridevirtual

Returns an estimate of the structure capacity or surface area.

Reimplemented from G4VSolid.

Definition at line 938 of file G4MultiUnion.cc.

939{
940 if (fCubicVolume == 0)
941 {
942 G4AutoLock l(&munionMutex);
943 fCubicVolume = EstimateCubicVolume(1000000, 0.001);
944 l.unlock();
945 }
946 return fCubicVolume;
947}
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4double EstimateCubicVolume(G4int nStat, G4double epsilon) const
Definition G4VSolid.cc:229

◆ GetEntityType()

G4GeometryType G4MultiUnion::GetEntityType ( ) const
inlineoverridevirtual

Returns the type ID, "G4MultiUnion" of the solid.

Implements G4VSolid.

Definition at line 231 of file G4MultiUnion.hh.

231{ return "G4MultiUnion"; }

◆ GetNumberOfSolids()

G4int G4MultiUnion::GetNumberOfSolids ( ) const
inline

◆ GetNumOfConstituents()

G4int G4MultiUnion::GetNumOfConstituents ( ) const
overridevirtual

Returns the number of solids part of the structure.

Reimplemented from G4VSolid.

Definition at line 780 of file G4MultiUnion.cc.

781{
782 G4int num = 0;
783 for (const auto solid : fSolids)
784 {
785 num += solid->GetNumOfConstituents();
786 }
787 return num;
788}

◆ GetPointOnSurface()

G4ThreeVector G4MultiUnion::GetPointOnSurface ( ) const
overridevirtual

Returns a point (G4ThreeVector) randomly and uniformly generated on the surface of a solid.

Reimplemented from G4VSolid.

Definition at line 921 of file G4MultiUnion.cc.

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}
G4double G4QuickRand(uint32_t seed=0)
long G4long
Definition G4Types.hh:87
EInside Inside(const G4ThreeVector &aPoint) const override
virtual G4ThreeVector GetPointOnSurface() const
Definition G4VSolid.cc:151
@ kSurface
Definition geomdefs.hh:69

◆ GetPolyhedron()

G4Polyhedron * G4MultiUnion::GetPolyhedron ( ) const
overridevirtual

Smart access function - creates on request and stores for future access. A null pointer means "not available".

Reimplemented from G4VSolid.

Definition at line 1004 of file G4MultiUnion.cc.

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}
G4Polyhedron * CreatePolyhedron() const override

◆ GetSolid()

G4VSolid * G4MultiUnion::GetSolid ( G4int index) const
inline

◆ GetSurfaceArea()

G4double G4MultiUnion::GetSurfaceArea ( )
overridevirtual

Returns an estimation of the solid surface area in internal units. This method may be overloaded by derived classes to compute the exact geometrical quantity for solids where this is possible, or anyway to cache the computed value. Note: the computed value is NOT cached.

Reimplemented from G4VSolid.

Definition at line 950 of file G4MultiUnion.cc.

951{
952 if (fSurfaceArea == 0)
953 {
954 G4AutoLock l(&munionMutex);
955 fSurfaceArea = EstimateSurfaceArea(1000000, -1.);
956 l.unlock();
957 }
958 return fSurfaceArea;
959}
G4double EstimateSurfaceArea(G4int nStat, G4double epsilon) const
Definition G4VSolid.cc:290

◆ GetTransformation()

const G4Transform3D & G4MultiUnion::GetTransformation ( G4int index) const
inline

Accessors to retrieve a transformation or a solid, given an index and the total number of solids in the structure.

Referenced by CreatePolyhedron(), G4tgbGeometryDumper::DumpMultiUnionVolume(), Extent(), and G4GDMLWriteSolids::MultiUnionWrite().

◆ GetVoxels()

G4Voxelizer & G4MultiUnion::GetVoxels ( ) const
inline

Returns the xoxelised optimisation structure.

◆ Inside()

EInside G4MultiUnion::Inside ( const G4ThreeVector & aPoint) const
overridevirtual

Returns if the given point "aPoint" is inside or not the solid.

Implements G4VSolid.

Definition at line 490 of file G4MultiUnion.cc.

491{
492 return InsideWithExclusion(aPoint);
493}

Referenced by GetPointOnSurface().

◆ IsFaceted()

G4bool G4MultiUnion::IsFaceted ( ) const
overridevirtual

Returns false if any of the solids part of the structure is not faceted.

Reimplemented from G4VSolid.

Definition at line 791 of file G4MultiUnion.cc.

792{
793 for (const auto solid : fSolids)
794 {
795 if (!solid->IsFaceted()) { return false; }
796 }
797 return true;
798}

◆ operator=()

G4MultiUnion & G4MultiUnion::operator= ( const G4MultiUnion & rhs)

Definition at line 110 of file G4MultiUnion.cc.

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}
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:108

◆ SetAccurateSafety()

void G4MultiUnion::SetAccurateSafety ( G4bool flag)
inline

◆ StreamInfo()

std::ostream & G4MultiUnion::StreamInfo ( std::ostream & os) const
overridevirtual

Streams the object contents to an output stream.

Implements G4VSolid.

Definition at line 895 of file G4MultiUnion.cc.

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}
G4String GetName() const
virtual std::ostream & StreamInfo(std::ostream &os) const =0
CLHEP::HepRotation getRotation() const
CLHEP::Hep3Vector getTranslation() const

◆ SurfaceNormal()

G4ThreeVector G4MultiUnion::SurfaceNormal ( const G4ThreeVector & aPoint) const
overridevirtual

Returns the outwards pointing unit normal of the shape for the surface closest to the point at offset "aPoint".

Implements G4VSolid.

Definition at line 609 of file G4MultiUnion.cc.

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}
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0

◆ Voxelize()

void G4MultiUnion::Voxelize ( )

Finalises and prepares for use, creating the optimisation structure for all solids in the structure. It must be called once before navigation use.

Definition at line 801 of file G4MultiUnion.cc.

802{
803 fVoxels.Voxelize(fSolids, fTransformObjs);
804}

Referenced by G4tgbVolume::FindOrConstructG4Solid(), and G4GDMLReadSolids::MultiUnionRead().

◆ G4Voxelizer

friend class G4Voxelizer
friend

Definition at line 63 of file G4MultiUnion.hh.

Referenced by G4Voxelizer, and GetVoxels().


The documentation for this class was generated from the following files: