Geant4 10.7.0
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, 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)

◆ ~G4MultiUnion()

G4MultiUnion::~G4MultiUnion ( )

Definition at line 64 of file G4MultiUnion.cc.

65{
66}

◆ G4MultiUnion() [3/4]

G4MultiUnion::G4MultiUnion ( const G4MultiUnion rhs)

Definition at line 83 of file G4MultiUnion.cc.

84 : G4VSolid(rhs), fCubicVolume(rhs.fCubicVolume),
85 fSurfaceArea(rhs.fSurfaceArea),
86 kRadTolerance(rhs.kRadTolerance), fAccurate(rhs.fAccurate)
87{
88}

◆ G4MultiUnion() [4/4]

G4MultiUnion::G4MultiUnion ( __void__ &  a)

Definition at line 92 of file G4MultiUnion.cc.

93 : G4VSolid(a)
94{
95}

Member Function Documentation

◆ AddNode()

void G4MultiUnion::AddNode ( G4VSolid solid,
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 G4GDMLReadSolids::MultiUnionNodeRead().

◆ BoundingLimits()

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

Reimplemented from G4VSolid.

Definition at line 609 of file G4MultiUnion.cc.

611{
612 Extent(kXAxis, aMin[0], aMax[0]);
613 Extent(kYAxis, aMin[1], aMax[1]);
614 Extent(kZAxis, aMin[2], aMax[2]);
615}
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 619 of file G4MultiUnion.cc.

623{
624 G4ThreeVector bmin, bmax;
625
626 // Get bounding box
627 BoundingLimits(bmin,bmax);
628
629 // Find extent
630 G4BoundingEnvelope bbox(bmin,bmax);
631 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
632}
void BoundingLimits(G4ThreeVector &aMin, G4ThreeVector &aMax) const

◆ Clone()

G4VSolid * G4MultiUnion::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 76 of file G4MultiUnion.cc.

77{
78 return new G4MultiUnion(*this);
79}

◆ CreatePolyhedron()

G4Polyhedron * G4MultiUnion::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 961 of file G4MultiUnion.cc.

962{
965
966 G4VSolid* solidA = GetSolid(0);
967 const G4Transform3D transform0=GetTransformation(0);
968 G4DisplacedSolid dispSolidA("placedA",solidA,transform0);
969
970 G4Polyhedron* top = new G4Polyhedron(*dispSolidA.GetPolyhedron());
971
972 for(G4int i=1; i<GetNumberOfSolids(); ++i)
973 {
974 G4VSolid* solidB = GetSolid(i);
975 const G4Transform3D transform=GetTransformation(i);
976 G4DisplacedSolid dispSolidB("placedB",solidB,transform);
977 G4Polyhedron* operand = dispSolidB.GetPolyhedron();
978 processor.push_back (operation, *operand);
979 }
980
981 if (processor.execute(*top)) { return top; }
982 else { return 0; }
983}
int G4int
Definition: G4Types.hh:85
const G4Transform3D & GetTransformation(G4int index) const
G4int GetNumberOfSolids() const
G4VSolid * GetSolid(G4int index) const
#define processor
Definition: xmlparse.cc:617

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

void G4MultiUnion::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 955 of file G4MultiUnion.cc.

956{
957 scene.AddSolid (*this);
958}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 753 of file G4MultiUnion.cc.

754{
755 // Estimates the isotropic safety from a point outside the current solid to
756 // any of its surfaces. The algorithm may be accurate or should provide a fast
757 // underestimate.
758
759 if (!fAccurate) { return fVoxels.DistanceToBoundingBox(point); }
760
761 const std::vector<G4VoxelBox>& boxes = fVoxels.GetBoxes();
762 G4double safetyMin = kInfinity;
763 G4ThreeVector localPoint;
764
765 G4int numNodes = fSolids.size();
766 for (G4int j = 0; j < numNodes; ++j)
767 {
768 G4ThreeVector dxyz;
769 if (j > 0)
770 {
771 const G4ThreeVector& pos = boxes[j].pos;
772 const G4ThreeVector& hlen = boxes[j].hlen;
773 for (auto i = 0; i <= 2; ++i)
774 // distance to middle point - hlength => distance from point to border
775 // of x,y,z
776 if ((dxyz[i] = std::abs(point[i] - pos[i]) - hlen[i]) > safetyMin)
777 continue;
778
779 G4double d2xyz = 0.;
780 for (auto i = 0; i <= 2; ++i)
781 if (dxyz[i] > 0) d2xyz += dxyz[i] * dxyz[i];
782
783 // minimal distance is at least this, but could be even higher. therefore,
784 // we can stop if previous was already lower, let us check if it does any
785 // chance to be better tha previous values...
786 if (d2xyz >= safetyMin * safetyMin)
787 {
788 continue;
789 }
790 }
791 const G4Transform3D& transform = fTransformObjs[j];
792 localPoint = GetLocalPoint(transform, point);
793 G4VSolid& solid = *fSolids[j];
794
795 G4double safety = solid.DistanceToIn(localPoint);
796 if (safety <= 0) return safety;
797 // it was detected, that the point is not located outside
798 if (safetyMin > safety) safetyMin = safety;
799 }
800 return safetyMin;
801}
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 205 of file G4MultiUnion.cc.

207{
208 G4double minDistance = kInfinity;
209 G4ThreeVector direction = aDirection.unit();
210 G4double shift = fVoxels.DistanceToFirst(aPoint, direction);
211 if (shift == kInfinity) return shift;
212
213 G4ThreeVector currentPoint = aPoint;
214 if (shift) currentPoint += direction * shift;
215
216 G4SurfBits exclusion(fVoxels.GetBitsPerSlice());
217 std::vector<G4int> candidates, curVoxel(3);
218 fVoxels.GetVoxel(curVoxel, currentPoint);
219
220 do
221 {
222 {
223 if (fVoxels.GetCandidatesVoxelArray(curVoxel, candidates, &exclusion))
224 {
225 G4double distance = DistanceToInCandidates(aPoint, direction,
226 candidates, exclusion);
227 if (minDistance > distance) minDistance = distance;
228 if (distance < shift) break;
229 }
230 }
231 shift = fVoxels.DistanceToNext(aPoint, direction, curVoxel);
232 }
233 while (minDistance > shift);
234
235 return minDistance;
236}
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:978
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 146 of file G4MultiUnion.cc.

148{
149 G4ThreeVector direction = aDirection.unit();
150 G4ThreeVector localPoint, localDirection;
151 G4double minDistance = kInfinity;
152
153 G4int numNodes = fSolids.size();
154 for (G4int i = 0 ; i < numNodes ; ++i)
155 {
156 G4VSolid& solid = *fSolids[i];
157 const G4Transform3D& transform = fTransformObjs[i];
158
159 localPoint = GetLocalPoint(transform, aPoint);
160 localDirection = GetLocalVector(transform, direction);
161
162 G4double distance = solid.DistanceToIn(localPoint, localDirection);
163 if (minDistance > distance) minDistance = distance;
164 }
165 return minDistance;
166}

◆ DistanceToOut() [1/2]

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

Implements G4VSolid.

Definition at line 716 of file G4MultiUnion.cc.

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

295{
296 return DistanceToOutVoxels(aPoint, aDirection, aNormal);
297}
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 239 of file G4MultiUnion.cc.

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

◆ DistanceToOutVoxels()

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

Definition at line 300 of file G4MultiUnion.cc.

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

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

117{
118 // Computes the cubic volume of the "G4MultiUnion" structure using
119 // random points
120
121 if (!fCubicVolume)
122 {
123 G4ThreeVector extentMin, extentMax, d, p, point;
124 G4int inside = 0, generated;
125 BoundingLimits(extentMin, extentMax);
126 d = (extentMax - extentMin) / 2.;
127 p = (extentMax + extentMin) / 2.;
128 G4ThreeVector left = p - d;
129 G4ThreeVector length = d * 2;
130 for (generated = 0; generated < 10000; ++generated)
131 {
133 point = left + G4ThreeVector(length.x()*rvec.x(),
134 length.y()*rvec.y(),
135 length.z()*rvec.z());
136 if (Inside(point) != EInside::kOutside) ++inside;
137 }
138 G4double vbox = (2 * d.x()) * (2 * d.y()) * (2 * d.z());
139 fCubicVolume = inside * vbox / generated;
140 }
141 return fCubicVolume;
142}
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 120 of file G4MultiUnion.hh.

120{ return "G4MultiUnion"; }

◆ GetNumberOfSolids()

G4int G4MultiUnion::GetNumberOfSolids ( ) const
inline

Definition at line 203 of file G4MultiUnion.hh.

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

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

◆ GetPointOnSurface()

G4ThreeVector G4MultiUnion::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 934 of file G4MultiUnion.cc.

935{
936 G4ThreeVector point;
937
938 G4long size = fSolids.size();
939
940 do
941 {
942 G4long rnd = G4RandFlat::shootInt(G4long(0), size);
943 G4VSolid& solid = *fSolids[rnd];
944 point = solid.GetPointOnSurface();
945 const G4Transform3D& transform = fTransformObjs[rnd];
946 point = GetGlobalPoint(transform, point);
947 }
948 while (Inside(point) != EInside::kSurface);
949
950 return point;
951}
long G4long
Definition: G4Types.hh:87
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:140

◆ GetPolyhedron()

G4Polyhedron * G4MultiUnion::GetPolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 986 of file G4MultiUnion.cc.

987{
988 if (fpPolyhedron == nullptr ||
989 fRebuildPolyhedron ||
991 fpPolyhedron->GetNumberOfRotationSteps())
992 {
993 G4AutoLock l(&polyhedronMutex);
994 delete fpPolyhedron;
995 fpPolyhedron = CreatePolyhedron();
996 fRebuildPolyhedron = false;
997 l.unlock();
998 }
999 return fpPolyhedron;
1000}
G4Polyhedron * CreatePolyhedron() const
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
static G4int GetNumberOfRotationSteps()

◆ GetSolid()

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

Definition at line 197 of file G4MultiUnion.hh.

198{
199 return fSolids[index];
200}

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

◆ GetSurfaceArea()

G4double G4MultiUnion::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 804 of file G4MultiUnion.cc.

805{
806 if (!fSurfaceArea)
807 {
808 fSurfaceArea = EstimateSurfaceArea(1000000, 0.001);
809 }
810 return fSurfaceArea;
811}
G4double EstimateSurfaceArea(G4int nStat, G4double ell) const
Definition: G4VSolid.cc:253

◆ GetTransformation()

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

Definition at line 191 of file G4MultiUnion.hh.

192{
193 return fTransformObjs[index];
194}

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

◆ GetVoxels()

G4Voxelizer & G4MultiUnion::GetVoxels ( ) const
inline

Definition at line 185 of file G4MultiUnion.hh.

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

◆ Inside()

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

Implements G4VSolid.

Definition at line 496 of file G4MultiUnion.cc.

497{
498 // Classify point location with respect to solid:
499 // o eInside - inside the solid
500 // o eSurface - close to surface within tolerance
501 // o eOutside - outside the solid
502
503 // Hitherto, it is considered that only parallelepipedic nodes can be
504 // added to the container
505
506 // Implementation using voxelisation techniques:
507 // ---------------------------------------------
508
509 // return InsideIterator(aPoint);
510
511 EInside location = InsideWithExclusion(aPoint);
512 return location;
513}

Referenced by GetCubicVolume(), and GetPointOnSurface().

◆ InsideIterator()

EInside G4MultiUnion::InsideIterator ( const G4ThreeVector aPoint) const

◆ InsideNoVoxels()

EInside G4MultiUnion::InsideNoVoxels ( const G4ThreeVector aPoint) const

Definition at line 516 of file G4MultiUnion.cc.

517{
518 G4ThreeVector localPoint;
519 EInside location = EInside::kOutside;
520 G4int countSurface = 0;
521
522 G4int numNodes = fSolids.size();
523 for (G4int i = 0 ; i < numNodes ; ++i)
524 {
525 G4VSolid& solid = *fSolids[i];
526 G4Transform3D transform = GetTransformation(i);
527
528 // The coordinates of the point are modified so as to fit the
529 // intrinsic solid local frame:
530 localPoint = GetLocalPoint(transform, aPoint);
531
532 location = solid.Inside(localPoint);
533
534 if (location == EInside::kSurface)
535 ++countSurface;
536
537 if (location == EInside::kInside) return EInside::kInside;
538 }
539 if (countSurface != 0) return EInside::kSurface;
540 return EInside::kOutside;
541}

◆ operator=()

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

Definition at line 99 of file G4MultiUnion.cc.

100{
101 // Check assignment to self
102 //
103 if (this == &rhs)
104 {
105 return *this;
106 }
107
108 // Copy base class data
109 //
111
112 return *this;
113}
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:98

◆ SetAccurateSafety()

void G4MultiUnion::SetAccurateSafety ( G4bool  flag)
inline

Definition at line 209 of file G4MultiUnion.hh.

210{
211 fAccurate = flag;
212}

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 908 of file G4MultiUnion.cc.

909{
910 G4int oldprc = os.precision(16);
911 os << "-----------------------------------------------------------\n"
912 << " *** Dump for solid - " << GetName() << " ***\n"
913 << " ===================================================\n"
914 << " Solid type: G4MultiUnion\n"
915 << " Parameters: \n";
916 G4int numNodes = fSolids.size();
917 for (G4int i = 0 ; i < numNodes ; ++i)
918 {
919 G4VSolid& solid = *fSolids[i];
920 solid.StreamInfo(os);
921 const G4Transform3D& transform = fTransformObjs[i];
922 os << " Translation is " << transform.getTranslation() << " \n";
923 os << " Rotation is :" << " \n";
924 os << " " << transform.getRotation() << "\n";
925 }
926 os << " \n"
927 << "-----------------------------------------------------------\n";
928 os.precision(oldprc);
929
930 return os;
931}
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 635 of file G4MultiUnion.cc.

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

◆ Voxelize()

void G4MultiUnion::Voxelize ( )

Definition at line 814 of file G4MultiUnion.cc.

815{
816 fVoxels.Voxelize(fSolids, fTransformObjs);
817}
void Voxelize(std::vector< G4VSolid * > &solids, std::vector< G4Transform3D > &transforms)
Definition: G4Voxelizer.cc:718

Referenced by 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: