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

#include <G4MultiUnion.hh>

+ Inheritance diagram for G4MultiUnion:

Public Member Functions

 G4MultiUnion ()
 
 G4MultiUnion (const G4String &name)
 
 ~G4MultiUnion ()
 
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
 
EInside InsideIterator (const G4ThreeVector &aPoint) const
 
G4double DistanceToIn (const G4ThreeVector &aPoint) const
 
G4double DistanceToOut (const G4ThreeVector &aPoint) const
 
void SetAccurateSafety (G4bool flag)
 
G4double DistanceToIn (const G4ThreeVector &aPoint, const G4ThreeVector &aDirection) const
 
G4double DistanceToOut (const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *aNormalVector=nullptr) const
 
G4double DistanceToInNoVoxels (const G4ThreeVector &aPoint, const G4ThreeVector &aDirection) const
 
G4double DistanceToOutVoxels (const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const
 
G4double DistanceToOutVoxelsCore (const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector, G4bool &aConvex, std::vector< G4int > &candidates) const
 
G4double DistanceToOutNoVoxels (const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &aPoint) const
 
void Extent (EAxis aAxis, G4double &aMin, G4double &aMax) const
 
void BoundingLimits (G4ThreeVector &aMin, G4ThreeVector &aMax) const
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
 
G4double GetCubicVolume ()
 
G4double GetSurfaceArea ()
 
G4VSolidClone () const
 
G4GeometryType GetEntityType () const
 
void Voxelize ()
 
EInside InsideNoVoxels (const G4ThreeVector &aPoint) const
 
G4VoxelizerGetVoxels () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4ThreeVector GetPointOnSurface () const
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
G4PolyhedronCreatePolyhedron () const
 
G4PolyhedronGetPolyhedron () const
 
 G4MultiUnion (__void__ &)
 
- Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
 
virtual ~G4VSolid ()
 
G4bool operator== (const G4VSolid &s) const
 
G4String GetName () const
 
void SetName (const G4String &name)
 
G4double GetTolerance () const
 
virtual void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const
 
virtual G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
 
virtual EInside Inside (const G4ThreeVector &p) const =0
 
virtual G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const =0
 
virtual G4double DistanceToIn (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
 
virtual G4double DistanceToOut (const G4ThreeVector &p) const =0
 
virtual void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
virtual G4double GetCubicVolume ()
 
virtual G4double GetSurfaceArea ()
 
virtual G4GeometryType GetEntityType () const =0
 
virtual G4ThreeVector GetPointOnSurface () const
 
virtual G4VSolidClone () const
 
virtual std::ostream & StreamInfo (std::ostream &os) const =0
 
void DumpInfo () const
 
virtual void DescribeYourselfTo (G4VGraphicsScene &scene) const =0
 
virtual G4VisExtent GetExtent () const
 
virtual G4PolyhedronCreatePolyhedron () const
 
virtual G4PolyhedronGetPolyhedron () const
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
 G4VSolid (__void__ &)
 
 G4VSolid (const G4VSolid &rhs)
 
G4VSolidoperator= (const G4VSolid &rhs)
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) 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

Definition at line 53 of file G4MultiUnion.hh.

Constructor & Destructor Documentation

◆ G4MultiUnion() [1/4]

G4MultiUnion::G4MultiUnion ( )
inline

Definition at line 59 of file G4MultiUnion.hh.

59: G4VSolid("") {}

Referenced by Clone().

◆ G4MultiUnion() [2/4]

G4MultiUnion::G4MultiUnion ( const G4String name)

Definition at line 54 of file G4MultiUnion.cc.

55 : G4VSolid(name)
56{
57 SetName(name);
58 fSolids.clear();
59 fTransformObjs.clear();
61}
G4double GetRadialTolerance() const
static G4GeometryTolerance * GetInstance()
void SetName(const G4String &name)
Definition: G4VSolid.cc:127

◆ ~G4MultiUnion()

G4MultiUnion::~G4MultiUnion ( )

Definition at line 64 of file G4MultiUnion.cc.

65{
66}

◆ G4MultiUnion() [3/4]

G4MultiUnion::G4MultiUnion ( const G4MultiUnion rhs)

Definition at line 90 of file G4MultiUnion.cc.

91 : G4VSolid(rhs), fCubicVolume(rhs.fCubicVolume),
92 fSurfaceArea(rhs.fSurfaceArea),
93 kRadTolerance(rhs.kRadTolerance), fAccurate(rhs.fAccurate)
94{
95}

◆ G4MultiUnion() [4/4]

G4MultiUnion::G4MultiUnion ( __void__ &  a)

Definition at line 99 of file G4MultiUnion.cc.

100 : G4VSolid(a)
101{
102}

Member Function Documentation

◆ AddNode() [1/2]

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

Definition at line 69 of file G4MultiUnion.cc.

70{
71 fSolids.push_back(&solid);
72 fTransformObjs.push_back(trans); // Store a local copy of transformations
73}

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

◆ AddNode() [2/2]

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

Definition at line 76 of file G4MultiUnion.cc.

77{
78 fSolids.push_back(solid);
79 fTransformObjs.push_back(trans); // Store a local copy of transformations
80}

◆ BoundingLimits()

void G4MultiUnion::BoundingLimits ( G4ThreeVector aMin,
G4ThreeVector aMax 
) const
virtual

Reimplemented from G4VSolid.

Definition at line 620 of file G4MultiUnion.cc.

622{
623 Extent(kXAxis, aMin[0], aMax[0]);
624 Extent(kYAxis, aMin[1], aMax[1]);
625 Extent(kZAxis, aMin[2], aMax[2]);
626}
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(), and GetCubicVolume().

◆ CalculateExtent()

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

Implements G4VSolid.

Definition at line 630 of file G4MultiUnion.cc.

634{
635 G4ThreeVector bmin, bmax;
636
637 // Get bounding box
638 BoundingLimits(bmin,bmax);
639
640 // Find extent
641 G4BoundingEnvelope bbox(bmin,bmax);
642 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
643}
void BoundingLimits(G4ThreeVector &aMin, G4ThreeVector &aMax) const

◆ Clone()

G4VSolid * G4MultiUnion::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 83 of file G4MultiUnion.cc.

84{
85 return new G4MultiUnion(*this);
86}

◆ CreatePolyhedron()

G4Polyhedron * G4MultiUnion::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 972 of file G4MultiUnion.cc.

973{
974 HepPolyhedronProcessor processor;
976
977 G4VSolid* solidA = GetSolid(0);
978 const G4Transform3D transform0=GetTransformation(0);
979 G4DisplacedSolid dispSolidA("placedA",solidA,transform0);
980
981 G4Polyhedron* top = new G4Polyhedron(*dispSolidA.GetPolyhedron());
982
983 for(G4int i=1; i<GetNumberOfSolids(); ++i)
984 {
985 G4VSolid* solidB = GetSolid(i);
986 const G4Transform3D transform=GetTransformation(i);
987 G4DisplacedSolid dispSolidB("placedB",solidB,transform);
988 G4Polyhedron* operand = dispSolidB.GetPolyhedron();
989 processor.push_back (operation, *operand);
990 }
991
992 if (processor.execute(*top)) { return top; }
993 else { return 0; }
994}
int G4int
Definition: G4Types.hh:85
const G4Transform3D & GetTransformation(G4int index) const
G4int GetNumberOfSolids() const
G4VSolid * GetSolid(G4int index) const
bool execute(HepPolyhedron &)
void push_back(Operation, const HepPolyhedron &)

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

void G4MultiUnion::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 966 of file G4MultiUnion.cc.

967{
968 scene.AddSolid (*this);
969}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

G4double G4MultiUnion::DistanceToIn ( const G4ThreeVector aPoint) const
virtual

Implements G4VSolid.

Definition at line 764 of file G4MultiUnion.cc.

765{
766 // Estimates the isotropic safety from a point outside the current solid to
767 // any of its surfaces. The algorithm may be accurate or should provide a fast
768 // underestimate.
769
770 if (!fAccurate) { return fVoxels.DistanceToBoundingBox(point); }
771
772 const std::vector<G4VoxelBox>& boxes = fVoxels.GetBoxes();
773 G4double safetyMin = kInfinity;
774 G4ThreeVector localPoint;
775
776 std::size_t numNodes = fSolids.size();
777 for (std::size_t j = 0; j < numNodes; ++j)
778 {
779 G4ThreeVector dxyz;
780 if (j > 0)
781 {
782 const G4ThreeVector& pos = boxes[j].pos;
783 const G4ThreeVector& hlen = boxes[j].hlen;
784 for (auto i = 0; i <= 2; ++i)
785 // distance to middle point - hlength => distance from point to border
786 // of x,y,z
787 if ((dxyz[i] = std::abs(point[i] - pos[i]) - hlen[i]) > safetyMin)
788 continue;
789
790 G4double d2xyz = 0.;
791 for (auto i = 0; i <= 2; ++i)
792 if (dxyz[i] > 0) d2xyz += dxyz[i] * dxyz[i];
793
794 // minimal distance is at least this, but could be even higher. therefore,
795 // we can stop if previous was already lower, let us check if it does any
796 // chance to be better tha previous values...
797 if (d2xyz >= safetyMin * safetyMin)
798 {
799 continue;
800 }
801 }
802 const G4Transform3D& transform = fTransformObjs[j];
803 localPoint = GetLocalPoint(transform, point);
804 G4VSolid& solid = *fSolids[j];
805
806 G4double safety = solid.DistanceToIn(localPoint);
807 if (safety <= 0) return safety;
808 // it was detected, that the point is not located outside
809 if (safetyMin > safety) safetyMin = safety;
810 }
811 return safetyMin;
812}
double G4double
Definition: G4Types.hh:83
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
G4double DistanceToBoundingBox(const G4ThreeVector &point) const
const std::vector< G4VoxelBox > & GetBoxes() const

◆ DistanceToIn() [2/2]

G4double G4MultiUnion::DistanceToIn ( const G4ThreeVector aPoint,
const G4ThreeVector aDirection 
) const
virtual

Implements G4VSolid.

Definition at line 212 of file G4MultiUnion.cc.

214{
215 G4double minDistance = kInfinity;
216 G4ThreeVector direction = aDirection.unit();
217 G4double shift = fVoxels.DistanceToFirst(aPoint, direction);
218 if (shift == kInfinity) return shift;
219
220 G4ThreeVector currentPoint = aPoint;
221 if (shift) currentPoint += direction * shift;
222
223 G4SurfBits exclusion(fVoxels.GetBitsPerSlice());
224 std::vector<G4int> candidates, curVoxel(3);
225 fVoxels.GetVoxel(curVoxel, currentPoint);
226
227 do
228 {
229 {
230 if (fVoxels.GetCandidatesVoxelArray(curVoxel, candidates, &exclusion))
231 {
232 G4double distance = DistanceToInCandidates(aPoint, direction,
233 candidates, exclusion);
234 if (minDistance > distance) minDistance = distance;
235 if (distance < shift) break;
236 }
237 }
238 shift = fVoxels.DistanceToNext(aPoint, direction, curVoxel);
239 }
240 while (minDistance > shift);
241
242 return minDistance;
243}
Hep3Vector unit() const
void GetVoxel(std::vector< G4int > &curVoxel, const G4ThreeVector &point) const
G4int GetBitsPerSlice() const
G4double DistanceToFirst(const G4ThreeVector &point, const G4ThreeVector &direction) const
G4int GetCandidatesVoxelArray(const G4ThreeVector &point, std::vector< G4int > &list, G4SurfBits *crossed=nullptr) const
Definition: G4Voxelizer.cc:966
G4double DistanceToNext(const G4ThreeVector &point, const G4ThreeVector &direction, std::vector< G4int > &curVoxel) const

◆ DistanceToInNoVoxels()

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

Definition at line 153 of file G4MultiUnion.cc.

155{
156 G4ThreeVector direction = aDirection.unit();
157 G4ThreeVector localPoint, localDirection;
158 G4double minDistance = kInfinity;
159
160 std::size_t numNodes = fSolids.size();
161 for (std::size_t i = 0 ; i < numNodes ; ++i)
162 {
163 G4VSolid& solid = *fSolids[i];
164 const G4Transform3D& transform = fTransformObjs[i];
165
166 localPoint = GetLocalPoint(transform, aPoint);
167 localDirection = GetLocalVector(transform, direction);
168
169 G4double distance = solid.DistanceToIn(localPoint, localDirection);
170 if (minDistance > distance) minDistance = distance;
171 }
172 return minDistance;
173}

◆ DistanceToOut() [1/2]

G4double G4MultiUnion::DistanceToOut ( const G4ThreeVector aPoint) const
virtual

Implements G4VSolid.

Definition at line 727 of file G4MultiUnion.cc.

728{
729 // Estimates isotropic distance to the surface of the solid. This must
730 // be either accurate or an underestimate.
731 // Two modes: - default/fast mode, sacrificing accuracy for speed
732 // - "precise" mode, requests accurate value if available.
733
734 std::vector<G4int> candidates;
735 G4ThreeVector localPoint;
736 G4double safetyMin = kInfinity;
737
738 // In general, the value return by DistanceToIn(p) will not be the exact
739 // but only an undervalue (cf. overlaps)
740 fVoxels.GetCandidatesVoxelArray(point, candidates);
741
742 std::size_t limit = candidates.size();
743 for (std::size_t i = 0; i < limit; ++i)
744 {
745 G4int candidate = candidates[i];
746
747 // The coordinates of the point are modified so as to fit the intrinsic
748 // solid local frame:
749 const G4Transform3D& transform = fTransformObjs[candidate];
750 localPoint = GetLocalPoint(transform, point);
751 G4VSolid& solid = *fSolids[candidate];
752 if (solid.Inside(localPoint) == EInside::kInside)
753 {
754 G4double safety = solid.DistanceToOut(localPoint);
755 if (safetyMin > safety) safetyMin = safety;
756 }
757 }
758 if (safetyMin == kInfinity) safetyMin = 0; // we are not inside
759
760 return safetyMin;
761}
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

◆ DistanceToOut() [2/2]

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

Implements G4VSolid.

Definition at line 297 of file G4MultiUnion.cc.

302{
303 return DistanceToOutVoxels(aPoint, aDirection, aNormal);
304}
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 246 of file G4MultiUnion.cc.

249{
250 // Computes distance from a point presumably outside the solid to the solid
251 // surface. Ignores first surface if the point is actually inside.
252 // Early return infinity in case the safety to any surface is found greater
253 // than the proposed step aPstep.
254 // The normal vector to the crossed surface is filled only in case the box
255 // is crossed, otherwise aNormal->IsNull() is true.
256
257 // algorithm:
258 G4ThreeVector direction = aDirection.unit();
259 G4ThreeVector localPoint, localDirection;
260 G4int ignoredSolid = -1;
261 G4double resultDistToOut = 0;
262 G4ThreeVector currentPoint = aPoint;
263
264 G4int numNodes = (G4int)fSolids.size();
265 for (auto i = 0; i < numNodes; ++i)
266 {
267 if (i != ignoredSolid)
268 {
269 G4VSolid& solid = *fSolids[i];
270 const G4Transform3D& transform = fTransformObjs[i];
271 localPoint = GetLocalPoint(transform, currentPoint);
272 localDirection = GetLocalVector(transform, direction);
273 EInside location = solid.Inside(localPoint);
274 if (location != EInside::kOutside)
275 {
276 G4double distance = solid.DistanceToOut(localPoint, localDirection,
277 aNormal);
278 if (distance < kInfinity)
279 {
280 if (resultDistToOut == kInfinity) resultDistToOut = 0;
281 if (distance > 0)
282 {
283 currentPoint = GetGlobalPoint(transform, localPoint
284 + distance*localDirection);
285 resultDistToOut += distance;
286 ignoredSolid = i; // skip the solid which we have just left
287 i = -1; // force the loop to continue from 0
288 }
289 }
290 }
291 }
292 }
293 return resultDistToOut;
294}
EInside
Definition: geomdefs.hh:67

◆ DistanceToOutVoxels()

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

Definition at line 307 of file G4MultiUnion.cc.

310{
311 // Computes distance from a point presumably inside the solid to the solid
312 // surface. Ignores first surface along each axis systematically (for points
313 // inside or outside. Early returns zero in case the second surface is behind
314 // the starting point.
315 // o The proposed step is ignored.
316 // o The normal vector to the crossed surface is always filled.
317
318 // In the case the considered point is located inside the G4MultiUnion
319 // structure, the treatments are as follows:
320 // - investigation of the candidates for the passed point
321 // - progressive moving of the point towards the surface, along the
322 // passed direction
323 // - processing of the normal
324
325 G4ThreeVector direction = aDirection.unit();
326 std::vector<G4int> candidates;
327 G4double distance = 0;
328 std::size_t numNodes = 2*fSolids.size();
329 std::size_t count=0;
330
331 if (fVoxels.GetCandidatesVoxelArray(aPoint, candidates))
332 {
333 // For normal case for which we presume the point is inside
334 G4ThreeVector localPoint, localDirection, localNormal;
335 G4ThreeVector currentPoint = aPoint;
336 G4SurfBits exclusion(fVoxels.GetBitsPerSlice());
337 G4bool notOutside;
338 G4ThreeVector maxNormal;
339
340 do
341 {
342 notOutside = false;
343
344 G4double maxDistance = -kInfinity;
345 G4int maxCandidate = 0;
346 G4ThreeVector maxLocalPoint;
347
348 std::size_t limit = candidates.size();
349 for (std::size_t i = 0 ; i < limit ; ++i)
350 {
351 G4int candidate = candidates[i];
352 // ignore the current component (that you just got out of) since
353 // numerically the propagated point will be on its surface
354
355 G4VSolid& solid = *fSolids[candidate];
356 const G4Transform3D& transform = fTransformObjs[candidate];
357
358 // The coordinates of the point are modified so as to fit the
359 // intrinsic solid local frame:
360 localPoint = GetLocalPoint(transform, currentPoint);
361
362 // DistanceToOut at least for Trd sometimes return non-zero value
363 // even from points that are outside. Therefore, this condition
364 // must currently be here, otherwise it would not work.
365 // But it means it would be slower.
366
367 if (solid.Inside(localPoint) != EInside::kOutside)
368 {
369 notOutside = true;
370
371 localDirection = GetLocalVector(transform, direction);
372
373 // propagate with solid.DistanceToOut
374 G4double shift = solid.DistanceToOut(localPoint, localDirection,
375 false, 0, &localNormal);
376 if (maxDistance < shift)
377 {
378 maxDistance = shift;
379 maxCandidate = candidate;
380 maxNormal = localNormal;
381 }
382 }
383 }
384
385 if (notOutside)
386 {
387 const G4Transform3D& transform = fTransformObjs[maxCandidate];
388
389 // convert from local normal
390 if (aNormal) *aNormal = GetGlobalVector(transform, maxNormal);
391
392 distance += maxDistance;
393 currentPoint += maxDistance * direction;
394 if(maxDistance == 0.) ++count;
395
396 // the current component will be ignored
397 exclusion.SetBitNumber(maxCandidate);
398 EInside location = InsideWithExclusion(currentPoint, &exclusion);
399
400 // perform a Inside
401 // it should be excluded current solid from checking
402 // we have to collect the maximum distance from all given candidates.
403 // such "maximum" candidate should be then used for finding next
404 // candidates
405 if (location == EInside::kOutside)
406 {
407 // else return cumulated distances to outside of the traversed
408 // components
409 break;
410 }
411 // if inside another component, redo 1 to 3 but add the next
412 // DistanceToOut on top of the previous.
413
414 // and fill the candidates for the corresponding voxel (just
415 // exiting current component along direction)
416 candidates.clear();
417
418 fVoxels.GetCandidatesVoxelArray(currentPoint, candidates, &exclusion);
419 exclusion.ResetBitNumber(maxCandidate);
420 }
421 }
422 while ((notOutside) && (count < numNodes));
423 }
424
425 return distance;
426}
bool G4bool
Definition: G4Types.hh:86

Referenced by DistanceToOut().

◆ DistanceToOutVoxelsCore()

G4double G4MultiUnion::DistanceToOutVoxelsCore ( const G4ThreeVector aPoint,
const G4ThreeVector aDirection,
G4ThreeVector aNormalVector,
G4bool aConvex,
std::vector< G4int > &  candidates 
) const

◆ Extent()

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

Definition at line 555 of file G4MultiUnion.cc.

556{
557 // Determines the bounding box for the considered instance of "UMultipleUnion"
559
560 G4int numNodes = (G4int)fSolids.size();
561 for (auto i = 0 ; i < numNodes ; ++i)
562 {
563 G4VSolid& solid = *fSolids[i];
564 G4Transform3D transform = GetTransformation(i);
565 solid.BoundingLimits(min, max);
566
567 TransformLimits(min, max, transform);
568
569 if (i == 0)
570 {
571 switch (aAxis)
572 {
573 case kXAxis:
574 aMin = min.x();
575 aMax = max.x();
576 break;
577 case kYAxis:
578 aMin = min.y();
579 aMax = max.y();
580 break;
581 case kZAxis:
582 aMin = min.z();
583 aMax = max.z();
584 break;
585 default:
586 break;
587 }
588 }
589 else
590 {
591 // Determine the min/max on the considered axis:
592 switch (aAxis)
593 {
594 case kXAxis:
595 if (min.x() < aMin)
596 aMin = min.x();
597 if (max.x() > aMax)
598 aMax = max.x();
599 break;
600 case kYAxis:
601 if (min.y() < aMin)
602 aMin = min.y();
603 if (max.y() > aMax)
604 aMax = max.y();
605 break;
606 case kZAxis:
607 if (min.z() < aMin)
608 aMin = min.z();
609 if (max.z() > aMax)
610 aMax = max.z();
611 break;
612 default:
613 break;
614 }
615 }
616 }
617}
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4VSolid.cc:665
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 ( )
virtual

Reimplemented from G4VSolid.

Definition at line 123 of file G4MultiUnion.cc.

124{
125 // Computes the cubic volume of the "G4MultiUnion" structure using
126 // random points
127
128 if (!fCubicVolume)
129 {
130 G4ThreeVector extentMin, extentMax, d, p, point;
131 G4int inside = 0, generated;
132 BoundingLimits(extentMin, extentMax);
133 d = (extentMax - extentMin) / 2.;
134 p = (extentMax + extentMin) / 2.;
135 G4ThreeVector left = p - d;
136 G4ThreeVector length = d * 2;
137 for (generated = 0; generated < 10000; ++generated)
138 {
140 point = left + G4ThreeVector(length.x()*rvec.x(),
141 length.y()*rvec.y(),
142 length.z()*rvec.z());
143 if (Inside(point) != EInside::kOutside) ++inside;
144 }
145 G4double vbox = (2 * d.x()) * (2 * d.y()) * (2 * d.z());
146 fCubicVolume = inside * vbox / generated;
147 }
148 return fCubicVolume;
149}
CLHEP::Hep3Vector G4ThreeVector
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
double x() const
double y() const
EInside Inside(const G4ThreeVector &aPoint) const

◆ GetEntityType()

G4GeometryType G4MultiUnion::GetEntityType ( ) const
inlinevirtual

Implements G4VSolid.

Definition at line 121 of file G4MultiUnion.hh.

121{ return "G4MultiUnion"; }

◆ GetNumberOfSolids()

G4int G4MultiUnion::GetNumberOfSolids ( ) const
inline

Definition at line 204 of file G4MultiUnion.hh.

205{
206 return G4int(fSolids.size());
207}

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

◆ GetPointOnSurface()

G4ThreeVector G4MultiUnion::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 945 of file G4MultiUnion.cc.

946{
947 G4ThreeVector point;
948
949 G4long size = fSolids.size();
950
951 do
952 {
953 G4long rnd = G4RandFlat::shootInt(G4long(0), size);
954 G4VSolid& solid = *fSolids[rnd];
955 point = solid.GetPointOnSurface();
956 const G4Transform3D& transform = fTransformObjs[rnd];
957 point = GetGlobalPoint(transform, point);
958 }
959 while (Inside(point) != EInside::kSurface);
960
961 return point;
962}
long G4long
Definition: G4Types.hh:87
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152

◆ GetPolyhedron()

G4Polyhedron * G4MultiUnion::GetPolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 997 of file G4MultiUnion.cc.

998{
999 if (fpPolyhedron == nullptr ||
1000 fRebuildPolyhedron ||
1002 fpPolyhedron->GetNumberOfRotationSteps())
1003 {
1004 G4AutoLock l(&polyhedronMutex);
1005 delete fpPolyhedron;
1006 fpPolyhedron = CreatePolyhedron();
1007 fRebuildPolyhedron = false;
1008 l.unlock();
1009 }
1010 return fpPolyhedron;
1011}
G4Polyhedron * CreatePolyhedron() const
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
static G4int GetNumberOfRotationSteps()

◆ GetSolid()

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

Definition at line 198 of file G4MultiUnion.hh.

199{
200 return fSolids[index];
201}

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

◆ GetSurfaceArea()

G4double G4MultiUnion::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 815 of file G4MultiUnion.cc.

816{
817 if (!fSurfaceArea)
818 {
819 fSurfaceArea = EstimateSurfaceArea(1000000, 0.001);
820 }
821 return fSurfaceArea;
822}
G4double EstimateSurfaceArea(G4int nStat, G4double ell) const
Definition: G4VSolid.cc:265

◆ GetTransformation()

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

Definition at line 192 of file G4MultiUnion.hh.

193{
194 return fTransformObjs[index];
195}

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

◆ GetVoxels()

G4Voxelizer & G4MultiUnion::GetVoxels ( ) const
inline

Definition at line 186 of file G4MultiUnion.hh.

187{
188 return (G4Voxelizer&)fVoxels;
189}

◆ Inside()

EInside G4MultiUnion::Inside ( const G4ThreeVector aPoint) const
virtual

Implements G4VSolid.

Definition at line 507 of file G4MultiUnion.cc.

508{
509 // Classify point location with respect to solid:
510 // o eInside - inside the solid
511 // o eSurface - close to surface within tolerance
512 // o eOutside - outside the solid
513
514 // Hitherto, it is considered that only parallelepipedic nodes can be
515 // added to the container
516
517 // Implementation using voxelisation techniques:
518 // ---------------------------------------------
519
520 // return InsideIterator(aPoint);
521
522 EInside location = InsideWithExclusion(aPoint);
523 return location;
524}

Referenced by GetCubicVolume(), and GetPointOnSurface().

◆ InsideIterator()

EInside G4MultiUnion::InsideIterator ( const G4ThreeVector aPoint) const

◆ InsideNoVoxels()

EInside G4MultiUnion::InsideNoVoxels ( const G4ThreeVector aPoint) const

Definition at line 527 of file G4MultiUnion.cc.

528{
529 G4ThreeVector localPoint;
530 EInside location = EInside::kOutside;
531 G4int countSurface = 0;
532
533 G4int numNodes = (G4int)fSolids.size();
534 for (auto i = 0 ; i < numNodes ; ++i)
535 {
536 G4VSolid& solid = *fSolids[i];
537 G4Transform3D transform = GetTransformation(i);
538
539 // The coordinates of the point are modified so as to fit the
540 // intrinsic solid local frame:
541 localPoint = GetLocalPoint(transform, aPoint);
542
543 location = solid.Inside(localPoint);
544
545 if (location == EInside::kSurface)
546 ++countSurface;
547
548 if (location == EInside::kInside) return EInside::kInside;
549 }
550 if (countSurface != 0) return EInside::kSurface;
551 return EInside::kOutside;
552}

◆ operator=()

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

Definition at line 106 of file G4MultiUnion.cc.

107{
108 // Check assignment to self
109 //
110 if (this == &rhs)
111 {
112 return *this;
113 }
114
115 // Copy base class data
116 //
118
119 return *this;
120}
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:107

◆ SetAccurateSafety()

void G4MultiUnion::SetAccurateSafety ( G4bool  flag)
inline

Definition at line 210 of file G4MultiUnion.hh.

211{
212 fAccurate = flag;
213}

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 919 of file G4MultiUnion.cc.

920{
921 G4long oldprc = os.precision(16);
922 os << "-----------------------------------------------------------\n"
923 << " *** Dump for solid - " << GetName() << " ***\n"
924 << " ===================================================\n"
925 << " Solid type: G4MultiUnion\n"
926 << " Parameters: \n";
927 std::size_t numNodes = fSolids.size();
928 for (std::size_t i = 0 ; i < numNodes ; ++i)
929 {
930 G4VSolid& solid = *fSolids[i];
931 solid.StreamInfo(os);
932 const G4Transform3D& transform = fTransformObjs[i];
933 os << " Translation is " << transform.getTranslation() << " \n";
934 os << " Rotation is :" << " \n";
935 os << " " << transform.getRotation() << "\n";
936 }
937 os << " \n"
938 << "-----------------------------------------------------------\n";
939 os.precision(oldprc);
940
941 return os;
942}
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
virtual

Implements G4VSolid.

Definition at line 646 of file G4MultiUnion.cc.

647{
648 // Computes the localNormal on a surface and returns it as a unit vector.
649 // Must return a valid vector. (even if the point is not on the surface).
650 //
651 // On an edge or corner, provide an average localNormal of all facets within
652 // tolerance
653 // NOTE: the tolerance value used in here is not yet the global surface
654 // tolerance - we will have to revise this value - TODO
655
656 std::vector<G4int> candidates;
657 G4ThreeVector localPoint, normal, localNormal;
658 G4double safety = kInfinity;
659 G4int node = 0;
660
661 ///////////////////////////////////////////////////////////////////////////
662 // Important comment: Cases for which the point is located on an edge or
663 // on a vertice remain to be treated
664
665 // determine weather we are in voxel area
666 if (fVoxels.GetCandidatesVoxelArray(aPoint, candidates))
667 {
668 std::size_t limit = candidates.size();
669 for (std::size_t i = 0 ; i < limit ; ++i)
670 {
671 G4int candidate = candidates[i];
672 const G4Transform3D& transform = fTransformObjs[candidate];
673
674 // The coordinates of the point are modified so as to fit the intrinsic
675 // solid local frame:
676 localPoint = GetLocalPoint(transform, aPoint);
677 G4VSolid& solid = *fSolids[candidate];
678 EInside location = solid.Inside(localPoint);
679
680 if (location == EInside::kSurface)
681 {
682 // normal case when point is on surface, we pick first solid
683 normal = GetGlobalVector(transform, solid.SurfaceNormal(localPoint));
684 return normal.unit();
685 }
686 else
687 {
688 // collect the smallest safety and remember solid node
689 G4double s = (location == EInside::kInside)
690 ? solid.DistanceToOut(localPoint)
691 : solid.DistanceToIn(localPoint);
692 if (s < safety)
693 {
694 safety = s;
695 node = candidate;
696 }
697 }
698 }
699 // on none of the solids, the point was not on the surface
700 G4VSolid& solid = *fSolids[node];
701 const G4Transform3D& transform = fTransformObjs[node];
702 localPoint = GetLocalPoint(transform, aPoint);
703
704 normal = GetGlobalVector(transform, solid.SurfaceNormal(localPoint));
705 return normal.unit();
706 }
707 else
708 {
709 // for the case when point is certainly outside:
710
711 // find a solid in union with the smallest safety
712 node = SafetyFromOutsideNumberNode(aPoint, safety);
713 G4VSolid& solid = *fSolids[node];
714
715 const G4Transform3D& transform = fTransformObjs[node];
716 localPoint = GetLocalPoint(transform, aPoint);
717
718 // evaluate normal for point at this found solid
719 // and transform multi-union coordinates
720 normal = GetGlobalVector(transform, solid.SurfaceNormal(localPoint));
721
722 return normal.unit();
723 }
724}
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0

◆ Voxelize()

void G4MultiUnion::Voxelize ( )

Definition at line 825 of file G4MultiUnion.cc.

826{
827 fVoxels.Voxelize(fSolids, fTransformObjs);
828}
void Voxelize(std::vector< G4VSolid * > &solids, std::vector< G4Transform3D > &transforms)
Definition: G4Voxelizer.cc:706

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

Friends And Related Function Documentation

◆ G4Voxelizer

friend class G4Voxelizer
friend

Definition at line 55 of file G4MultiUnion.hh.


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