Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4ReflectedSolid Class Reference

#include <G4ReflectedSolid.hh>

Inheritance diagram for G4ReflectedSolid:
G4VSolid

Public Member Functions

void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
 
G4VSolidClone () const
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
G4PolyhedronCreatePolyhedron () const
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
G4double DistanceToIn (const G4ThreeVector &p) const
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
G4double DistanceToOut (const G4ThreeVector &p) const
 
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
 
void DumpInfo () const
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 
 G4ReflectedSolid (const G4ReflectedSolid &rhs)
 
 G4ReflectedSolid (const G4String &pName, G4VSolid *pSolid, const G4Transform3D &transform)
 
G4VSolidGetConstituentMovedSolid () const
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
virtual G4double GetCubicVolume ()
 
G4Transform3D GetDirectTransform3D () const
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
virtual G4GeometryType GetEntityType () const
 
virtual G4VisExtent GetExtent () const
 
G4String GetName () const
 
G4ThreeVector GetPointOnSurface () const
 
G4PolyhedronGetPolyhedron () const
 
virtual G4ReflectedSolidGetReflectedSolidPtr ()
 
virtual const G4ReflectedSolidGetReflectedSolidPtr () const
 
virtual G4double GetSurfaceArea ()
 
G4double GetTolerance () const
 
G4Transform3D GetTransform3D () const
 
EInside Inside (const G4ThreeVector &p) const
 
G4ReflectedSolidoperator= (const G4ReflectedSolid &rhs)
 
G4bool operator== (const G4VSolid &s) const
 
void SetDirectTransform3D (G4Transform3D &)
 
void SetName (const G4String &name)
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
virtual ~G4ReflectedSolid ()
 

Protected Member Functions

void CalculateClippedPolygonExtent (G4ThreeVectorList &pPolygon, 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 ClipCrossSection (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

G4Transform3DfDirectTransform3D = nullptr
 
G4PolyhedronfpPolyhedron = nullptr
 
G4VSolidfPtrSolid = nullptr
 
G4bool fRebuildPolyhedron = false
 
G4double kCarTolerance
 

Private Member Functions

void ClipPolygonToSimpleLimits (G4ThreeVectorList &pPolygon, G4ThreeVectorList &outputPolygon, const G4VoxelLimits &pVoxelLimit) const
 

Private Attributes

G4String fshapeName
 

Detailed Description

Definition at line 42 of file G4ReflectedSolid.hh.

Constructor & Destructor Documentation

◆ G4ReflectedSolid() [1/2]

G4ReflectedSolid::G4ReflectedSolid ( const G4String pName,
G4VSolid pSolid,
const G4Transform3D transform 
)

Definition at line 51 of file G4ReflectedSolid.cc.

54 : G4VSolid(pName)
55{
56 fPtrSolid = pSolid;
58}
HepGeom::Transform3D G4Transform3D
G4Transform3D * fDirectTransform3D
G4VSolid(const G4String &name)
Definition: G4VSolid.cc:57
G4bool transform(G4String &input, const G4String &type)

References fDirectTransform3D, fPtrSolid, and G4coutFormatters::anonymous_namespace{G4coutFormatters.cc}::transform().

Referenced by Clone().

◆ ~G4ReflectedSolid()

G4ReflectedSolid::~G4ReflectedSolid ( )
virtual

Definition at line 63 of file G4ReflectedSolid.cc.

64{
65 delete fDirectTransform3D; fDirectTransform3D = nullptr;
66 delete fpPolyhedron; fpPolyhedron = nullptr;
67}
G4Polyhedron * fpPolyhedron

References fDirectTransform3D, and fpPolyhedron.

◆ G4ReflectedSolid() [2/2]

G4ReflectedSolid::G4ReflectedSolid ( const G4ReflectedSolid rhs)

Definition at line 72 of file G4ReflectedSolid.cc.

References fDirectTransform3D.

Member Function Documentation

◆ BoundingLimits()

void G4ReflectedSolid::BoundingLimits ( G4ThreeVector pMin,
G4ThreeVector pMax 
) const
virtual

Reimplemented from G4VSolid.

Definition at line 149 of file G4ReflectedSolid.cc.

151{
153 G4double xmin = pMin.x(), ymin = pMin.y(), zmin = pMin.z();
154 G4double xmax = pMax.x(), ymax = pMax.y(), zmax = pMax.z();
158
159 if (std::abs(xx) == 1 && std::abs(yy) == 1 && std::abs(zz) == 1)
160 {
161 // Special case of reflection in axis and pure translation
162 //
163 if (xx == -1) { G4double tmp = -xmin; xmin = -xmax; xmax = tmp; }
164 if (yy == -1) { G4double tmp = -ymin; ymin = -ymax; ymax = tmp; }
165 if (zz == -1) { G4double tmp = -zmin; zmin = -zmax; zmax = tmp; }
166 xmin += fDirectTransform3D->dx();
167 xmax += fDirectTransform3D->dx();
168 ymin += fDirectTransform3D->dy();
169 ymax += fDirectTransform3D->dy();
170 zmin += fDirectTransform3D->dz();
171 zmax += fDirectTransform3D->dz();
172 }
173 else
174 {
175 // Use additional reflection in Z to set up affine transformation
176 //
177 G4Transform3D transform3D = G4ReflectZ3D()*(*fDirectTransform3D);
179 transform3D.getTranslation());
180
181 // Find bounding box
182 //
183 G4VoxelLimits unLimit;
184 fPtrSolid->CalculateExtent(kXAxis,unLimit,transform,xmin,xmax);
185 fPtrSolid->CalculateExtent(kYAxis,unLimit,transform,ymin,ymax);
186 fPtrSolid->CalculateExtent(kZAxis,unLimit,transform,zmin,zmax);
187 }
188
189 pMin.set(xmin,ymin,-zmax);
190 pMax.set(xmax,ymax,-zmin);
191
192 // Check correctness of the bounding box
193 //
194 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
195 {
196 std::ostringstream message;
197 message << "Bad bounding box (min >= max) for solid: "
198 << GetName() << " !"
199 << "\npMin = " << pMin
200 << "\npMax = " << pMax;
201 G4Exception("G4ReflectedSolid::BoundingLimits()", "GeomMgt0001",
202 JustWarning, message);
203 DumpInfo();
204 }
205}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static const G4double pMax
static const G4double pMin
HepGeom::ReflectZ3D G4ReflectZ3D
double G4double
Definition: G4Types.hh:83
HepRotation inverse() const
G4String GetName() const
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
void DumpInfo() const
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4VSolid.cc:665
double dy() const
Definition: Transform3D.h:287
double zz() const
Definition: Transform3D.h:281
double dz() const
Definition: Transform3D.h:290
CLHEP::HepRotation getRotation() const
double dx() const
Definition: Transform3D.h:284
CLHEP::Hep3Vector getTranslation() const
double xx() const
Definition: Transform3D.h:257
double yy() const
Definition: Transform3D.h:269
@ kYAxis
Definition: geomdefs.hh:56
@ kXAxis
Definition: geomdefs.hh:55
@ kZAxis
Definition: geomdefs.hh:57

References G4VSolid::BoundingLimits(), G4VSolid::CalculateExtent(), G4VSolid::DumpInfo(), HepGeom::Transform3D::dx(), HepGeom::Transform3D::dy(), HepGeom::Transform3D::dz(), fDirectTransform3D, fPtrSolid, G4Exception(), G4VSolid::GetName(), HepGeom::Transform3D::getRotation(), HepGeom::Transform3D::getTranslation(), CLHEP::HepRotation::inverse(), JustWarning, kXAxis, kYAxis, kZAxis, pMax, pMin, G4coutFormatters::anonymous_namespace{G4coutFormatters.cc}::transform(), HepGeom::Transform3D::xx(), HepGeom::Transform3D::yy(), and HepGeom::Transform3D::zz().

◆ CalculateClippedPolygonExtent()

void G4VSolid::CalculateClippedPolygonExtent ( G4ThreeVectorList pPolygon,
const G4VoxelLimits pVoxelLimit,
const EAxis  pAxis,
G4double pMin,
G4double pMax 
) const
protectedinherited

Definition at line 489 of file G4VSolid.cc.

494{
495 G4int noLeft,i;
496 G4double component;
497
498 ClipPolygon(pPolygon,pVoxelLimit,pAxis);
499 noLeft = pPolygon.size();
500
501 if ( noLeft )
502 {
503 for (i=0; i<noLeft; ++i)
504 {
505 component = pPolygon[i].operator()(pAxis);
506
507 if (component < pMin)
508 {
509 pMin = component;
510 }
511 if (component > pMax)
512 {
513 pMax = component;
514 }
515 }
516 }
517}
int G4int
Definition: G4Types.hh:85
void ClipPolygon(G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
Definition: G4VSolid.cc:539

References G4VSolid::ClipPolygon(), pMax, and pMin.

Referenced by G4VSolid::ClipBetweenSections(), and G4VSolid::ClipCrossSection().

◆ CalculateExtent()

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

Implements G4VSolid.

Definition at line 212 of file G4ReflectedSolid.cc.

217{
218 // Separation of transformations. Calculation of the extent is done
219 // in a reflection of the global space. In such way, the voxel is
220 // reflected, but the solid is transformed just by G4AffineTransform.
221 // It allows one to use CalculateExtent() of the solid.
222
223 // Reflect voxel limits in Z
224 //
225 G4VoxelLimits limits;
226 limits.AddLimit(kXAxis, pVoxelLimits.GetMinXExtent(),
227 pVoxelLimits.GetMaxXExtent());
228 limits.AddLimit(kYAxis, pVoxelLimits.GetMinYExtent(),
229 pVoxelLimits.GetMaxYExtent());
230 limits.AddLimit(kZAxis,-pVoxelLimits.GetMaxZExtent(),
231 -pVoxelLimits.GetMinZExtent());
232
233 // Set affine transformation
234 //
235 G4Transform3D transform3D = G4ReflectZ3D()*pTransform*(*fDirectTransform3D);
237 transform3D.getTranslation());
238
239 // Find extent
240 //
241 if (!fPtrSolid->CalculateExtent(pAxis, limits, transform, pMin, pMax))
242 {
243 return false;
244 }
245 if (pAxis == kZAxis)
246 {
247 G4double tmp= -pMin; pMin= -pMax; pMax= tmp;
248 }
249
250 return true;
251}
void AddLimit(const EAxis pAxis, const G4double pMin, const G4double pMax)

References G4VoxelLimits::AddLimit(), G4VSolid::CalculateExtent(), fPtrSolid, G4VoxelLimits::GetMaxXExtent(), G4VoxelLimits::GetMaxYExtent(), G4VoxelLimits::GetMaxZExtent(), G4VoxelLimits::GetMinXExtent(), G4VoxelLimits::GetMinYExtent(), G4VoxelLimits::GetMinZExtent(), HepGeom::Transform3D::getRotation(), HepGeom::Transform3D::getTranslation(), CLHEP::HepRotation::inverse(), kXAxis, kYAxis, kZAxis, pMax, pMin, and G4coutFormatters::anonymous_namespace{G4coutFormatters.cc}::transform().

◆ ClipBetweenSections()

void G4VSolid::ClipBetweenSections ( G4ThreeVectorList pVertices,
const G4int  pSectionIndex,
const G4VoxelLimits pVoxelLimit,
const EAxis  pAxis,
G4double pMin,
G4double pMax 
) const
protectedinherited

Definition at line 444 of file G4VSolid.cc.

449{
450 G4ThreeVectorList polygon;
451 polygon.reserve(4);
452 polygon.push_back((*pVertices)[pSectionIndex]);
453 polygon.push_back((*pVertices)[pSectionIndex+4]);
454 polygon.push_back((*pVertices)[pSectionIndex+5]);
455 polygon.push_back((*pVertices)[pSectionIndex+1]);
456 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
457 polygon.clear();
458
459 polygon.push_back((*pVertices)[pSectionIndex+1]);
460 polygon.push_back((*pVertices)[pSectionIndex+5]);
461 polygon.push_back((*pVertices)[pSectionIndex+6]);
462 polygon.push_back((*pVertices)[pSectionIndex+2]);
463 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
464 polygon.clear();
465
466 polygon.push_back((*pVertices)[pSectionIndex+2]);
467 polygon.push_back((*pVertices)[pSectionIndex+6]);
468 polygon.push_back((*pVertices)[pSectionIndex+7]);
469 polygon.push_back((*pVertices)[pSectionIndex+3]);
470 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
471 polygon.clear();
472
473 polygon.push_back((*pVertices)[pSectionIndex+3]);
474 polygon.push_back((*pVertices)[pSectionIndex+7]);
475 polygon.push_back((*pVertices)[pSectionIndex+4]);
476 polygon.push_back((*pVertices)[pSectionIndex]);
477 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
478 return;
479}
std::vector< G4ThreeVector > G4ThreeVectorList
void CalculateClippedPolygonExtent(G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:489

References G4VSolid::CalculateClippedPolygonExtent(), pMax, and pMin.

◆ ClipCrossSection()

void G4VSolid::ClipCrossSection ( G4ThreeVectorList pVertices,
const G4int  pSectionIndex,
const G4VoxelLimits pVoxelLimit,
const EAxis  pAxis,
G4double pMin,
G4double pMax 
) const
protectedinherited

Definition at line 414 of file G4VSolid.cc.

419{
420
421 G4ThreeVectorList polygon;
422 polygon.reserve(4);
423 polygon.push_back((*pVertices)[pSectionIndex]);
424 polygon.push_back((*pVertices)[pSectionIndex+1]);
425 polygon.push_back((*pVertices)[pSectionIndex+2]);
426 polygon.push_back((*pVertices)[pSectionIndex+3]);
427 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
428 return;
429}

References G4VSolid::CalculateClippedPolygonExtent(), pMax, and pMin.

◆ ClipPolygon()

void G4VSolid::ClipPolygon ( G4ThreeVectorList pPolygon,
const G4VoxelLimits pVoxelLimit,
const EAxis  pAxis 
) const
protectedinherited

Definition at line 539 of file G4VSolid.cc.

542{
543 G4ThreeVectorList outputPolygon;
544
545 if ( pVoxelLimit.IsLimited() )
546 {
547 if (pVoxelLimit.IsXLimited() ) // && pAxis != kXAxis)
548 {
549 G4VoxelLimits simpleLimit1;
550 simpleLimit1.AddLimit(kXAxis,pVoxelLimit.GetMinXExtent(),kInfinity);
551 ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
552
553 pPolygon.clear();
554
555 if ( !outputPolygon.size() ) return;
556
557 G4VoxelLimits simpleLimit2;
558 simpleLimit2.AddLimit(kXAxis,-kInfinity,pVoxelLimit.GetMaxXExtent());
559 ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
560
561 if ( !pPolygon.size() ) return;
562 else outputPolygon.clear();
563 }
564 if ( pVoxelLimit.IsYLimited() ) // && pAxis != kYAxis)
565 {
566 G4VoxelLimits simpleLimit1;
567 simpleLimit1.AddLimit(kYAxis,pVoxelLimit.GetMinYExtent(),kInfinity);
568 ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
569
570 // Must always clear pPolygon - for clip to simpleLimit2 and in case of
571 // early exit
572
573 pPolygon.clear();
574
575 if ( !outputPolygon.size() ) return;
576
577 G4VoxelLimits simpleLimit2;
578 simpleLimit2.AddLimit(kYAxis,-kInfinity,pVoxelLimit.GetMaxYExtent());
579 ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
580
581 if ( !pPolygon.size() ) return;
582 else outputPolygon.clear();
583 }
584 if ( pVoxelLimit.IsZLimited() ) // && pAxis != kZAxis)
585 {
586 G4VoxelLimits simpleLimit1;
587 simpleLimit1.AddLimit(kZAxis,pVoxelLimit.GetMinZExtent(),kInfinity);
588 ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
589
590 // Must always clear pPolygon - for clip to simpleLimit2 and in case of
591 // early exit
592
593 pPolygon.clear();
594
595 if ( !outputPolygon.size() ) return;
596
597 G4VoxelLimits simpleLimit2;
598 simpleLimit2.AddLimit(kZAxis,-kInfinity,pVoxelLimit.GetMaxZExtent());
599 ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
600
601 // Return after final clip - no cleanup
602 }
603 }
604}
void ClipPolygonToSimpleLimits(G4ThreeVectorList &pPolygon, G4ThreeVectorList &outputPolygon, const G4VoxelLimits &pVoxelLimit) const
Definition: G4VSolid.cc:612
G4bool IsYLimited() const
G4double GetMinZExtent() const
G4bool IsXLimited() const
G4double GetMaxYExtent() const
G4double GetMaxZExtent() const
G4double GetMinYExtent() const
G4double GetMinXExtent() const
G4bool IsZLimited() const
G4bool IsLimited() const
G4double GetMaxXExtent() const
static const G4double kInfinity
Definition: geomdefs.hh:41

References G4VoxelLimits::AddLimit(), G4VSolid::ClipPolygonToSimpleLimits(), G4VoxelLimits::GetMaxXExtent(), G4VoxelLimits::GetMaxYExtent(), G4VoxelLimits::GetMaxZExtent(), G4VoxelLimits::GetMinXExtent(), G4VoxelLimits::GetMinYExtent(), G4VoxelLimits::GetMinZExtent(), G4VoxelLimits::IsLimited(), G4VoxelLimits::IsXLimited(), G4VoxelLimits::IsYLimited(), G4VoxelLimits::IsZLimited(), kInfinity, kXAxis, kYAxis, and kZAxis.

Referenced by G4VSolid::CalculateClippedPolygonExtent().

◆ ClipPolygonToSimpleLimits()

void G4VSolid::ClipPolygonToSimpleLimits ( G4ThreeVectorList pPolygon,
G4ThreeVectorList outputPolygon,
const G4VoxelLimits pVoxelLimit 
) const
privateinherited

Definition at line 612 of file G4VSolid.cc.

615{
616 G4int i;
617 G4int noVertices=pPolygon.size();
618 G4ThreeVector vEnd,vStart;
619
620 for (i = 0 ; i < noVertices ; ++i )
621 {
622 vStart = pPolygon[i];
623 if ( i == noVertices-1 ) vEnd = pPolygon[0];
624 else vEnd = pPolygon[i+1];
625
626 if ( pVoxelLimit.Inside(vStart) )
627 {
628 if (pVoxelLimit.Inside(vEnd))
629 {
630 // vStart and vEnd inside -> output end point
631 //
632 outputPolygon.push_back(vEnd);
633 }
634 else
635 {
636 // vStart inside, vEnd outside -> output crossing point
637 //
638 pVoxelLimit.ClipToLimits(vStart,vEnd);
639 outputPolygon.push_back(vEnd);
640 }
641 }
642 else
643 {
644 if (pVoxelLimit.Inside(vEnd))
645 {
646 // vStart outside, vEnd inside -> output inside section
647 //
648 pVoxelLimit.ClipToLimits(vStart,vEnd);
649 outputPolygon.push_back(vStart);
650 outputPolygon.push_back(vEnd);
651 }
652 else // Both point outside -> no output
653 {
654 // outputPolygon.push_back(vStart);
655 // outputPolygon.push_back(vEnd);
656 }
657 }
658 }
659}
G4bool ClipToLimits(G4ThreeVector &pStart, G4ThreeVector &pEnd) const
G4bool Inside(const G4ThreeVector &pVec) const

References G4VoxelLimits::ClipToLimits(), and G4VoxelLimits::Inside().

Referenced by G4VSolid::ClipPolygon().

◆ Clone()

G4VSolid * G4ReflectedSolid::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 366 of file G4ReflectedSolid.cc.

367{
368 return new G4ReflectedSolid(*this);
369}
G4ReflectedSolid(const G4String &pName, G4VSolid *pSolid, const G4Transform3D &transform)

References G4ReflectedSolid().

◆ ComputeDimensions()

void G4ReflectedSolid::ComputeDimensions ( G4VPVParameterisation p,
const G4int  n,
const G4VPhysicalVolume pRep 
)
virtual

Reimplemented from G4VSolid.

Definition at line 341 of file G4ReflectedSolid.cc.

344{
345 DumpInfo();
346 G4Exception("G4ReflectedSolid::ComputeDimensions()",
347 "GeomMgt0001", FatalException,
348 "Method not applicable in this context!");
349}
@ FatalException

References G4VSolid::DumpInfo(), FatalException, and G4Exception().

◆ CreatePolyhedron()

G4Polyhedron * G4ReflectedSolid::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 413 of file G4ReflectedSolid.cc.

414{
415 G4Polyhedron* polyhedron = fPtrSolid->CreatePolyhedron();
416 if (polyhedron != nullptr)
417 {
418 polyhedron->Transform(*fDirectTransform3D);
419 return polyhedron;
420 }
421 else
422 {
423 std::ostringstream message;
424 message << "Solid - " << GetName()
425 << " - original solid has no" << G4endl
426 << "corresponding polyhedron. Returning NULL!";
427 G4Exception("G4ReflectedSolid::CreatePolyhedron()",
428 "GeomMgt1001", JustWarning, message);
429 return nullptr;
430 }
431}
#define G4endl
Definition: G4ios.hh:57
virtual G4Polyhedron * CreatePolyhedron() const
Definition: G4VSolid.cc:700
HepPolyhedron & Transform(const G4Transform3D &t)

References G4VSolid::CreatePolyhedron(), fDirectTransform3D, fPtrSolid, G4endl, G4Exception(), G4VSolid::GetName(), JustWarning, and HepPolyhedron::Transform().

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

void G4ReflectedSolid::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 403 of file G4ReflectedSolid.cc.

404{
405 scene.AddSolid (*this);
406}
virtual void AddSolid(const G4Box &)=0

References G4VGraphicsScene::AddSolid().

◆ DistanceToIn() [1/2]

G4double G4ReflectedSolid::DistanceToIn ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 294 of file G4ReflectedSolid.cc.

295{
296 G4ThreeVector newPoint = (*fDirectTransform3D)*G4Point3D(p);
297 return fPtrSolid->DistanceToIn(newPoint);
298}
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:34
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0

References G4VSolid::DistanceToIn(), and fPtrSolid.

◆ DistanceToIn() [2/2]

G4double G4ReflectedSolid::DistanceToIn ( const G4ThreeVector p,
const G4ThreeVector v 
) const
virtual

Implements G4VSolid.

Definition at line 280 of file G4ReflectedSolid.cc.

282{
283 G4ThreeVector newPoint = (*fDirectTransform3D)*G4Point3D(p);
284 G4ThreeVector newDirection = (*fDirectTransform3D)*G4Vector3D(v);
285 return fPtrSolid->DistanceToIn(newPoint,newDirection);
286}
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:34

References G4VSolid::DistanceToIn(), and fPtrSolid.

◆ DistanceToOut() [1/2]

G4double G4ReflectedSolid::DistanceToOut ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 330 of file G4ReflectedSolid.cc.

331{
332 G4ThreeVector newPoint = (*fDirectTransform3D)*G4Point3D(p);
333 return fPtrSolid->DistanceToOut(newPoint);
334}
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0

References G4VSolid::DistanceToOut(), and fPtrSolid.

◆ DistanceToOut() [2/2]

G4double G4ReflectedSolid::DistanceToOut ( const G4ThreeVector p,
const G4ThreeVector v,
const G4bool  calcNorm = false,
G4bool validNorm = nullptr,
G4ThreeVector n = nullptr 
) const
virtual

Implements G4VSolid.

Definition at line 305 of file G4ReflectedSolid.cc.

310{
311 G4ThreeVector solNorm;
312
313 G4ThreeVector newPoint = (*fDirectTransform3D)*G4Point3D(p);
314 G4ThreeVector newDirection = (*fDirectTransform3D)*G4Vector3D(v);
315
316 G4double dist = fPtrSolid->DistanceToOut(newPoint, newDirection,
317 calcNorm, validNorm, &solNorm);
318 if(calcNorm)
319 {
320 *n = (*fDirectTransform3D)*G4Vector3D(solNorm);
321 }
322 return dist;
323}

References G4VSolid::DistanceToOut(), fPtrSolid, and CLHEP::detail::n.

◆ DumpInfo()

void G4VSolid::DumpInfo ( ) const
inlineinherited

Referenced by G4Cons::ApproxSurfaceNormal(), G4CutTubs::ApproxSurfaceNormal(), G4Sphere::ApproxSurfaceNormal(), G4Torus::ApproxSurfaceNormal(), G4Tubs::ApproxSurfaceNormal(), BoundingLimits(), G4DisplacedSolid::BoundingLimits(), G4IntersectionSolid::BoundingLimits(), G4ScaledSolid::BoundingLimits(), G4SubtractionSolid::BoundingLimits(), G4UnionSolid::BoundingLimits(), G4Box::BoundingLimits(), G4Cons::BoundingLimits(), G4CutTubs::BoundingLimits(), G4Orb::BoundingLimits(), G4Para::BoundingLimits(), G4Sphere::BoundingLimits(), G4Torus::BoundingLimits(), G4Trap::BoundingLimits(), G4Trd::BoundingLimits(), G4Tubs::BoundingLimits(), G4EllipticalCone::BoundingLimits(), G4ExtrudedSolid::BoundingLimits(), G4GenericPolycone::BoundingLimits(), G4GenericTrap::BoundingLimits(), G4Hype::BoundingLimits(), G4Paraboloid::BoundingLimits(), G4Polycone::BoundingLimits(), G4Polyhedra::BoundingLimits(), G4TessellatedSolid::BoundingLimits(), G4TwistedTubs::BoundingLimits(), G4ParameterisationBoxX::ComputeDimensions(), G4ParameterisationBoxY::ComputeDimensions(), G4ParameterisationBoxZ::ComputeDimensions(), G4ParameterisationConsRho::ComputeDimensions(), G4ParameterisationConsPhi::ComputeDimensions(), G4ParameterisationConsZ::ComputeDimensions(), G4ParameterisationParaX::ComputeDimensions(), G4ParameterisationParaY::ComputeDimensions(), G4ParameterisationParaZ::ComputeDimensions(), G4ParameterisationPolyconeRho::ComputeDimensions(), G4ParameterisationPolyconePhi::ComputeDimensions(), G4ParameterisationPolyconeZ::ComputeDimensions(), G4ParameterisationPolyhedraRho::ComputeDimensions(), G4ParameterisationPolyhedraPhi::ComputeDimensions(), G4ParameterisationPolyhedraZ::ComputeDimensions(), G4ParameterisationTrdX::ComputeDimensions(), G4ParameterisationTrdY::ComputeDimensions(), G4ParameterisationTrdZ::ComputeDimensions(), G4ParameterisationTubsRho::ComputeDimensions(), G4ParameterisationTubsPhi::ComputeDimensions(), G4ParameterisationTubsZ::ComputeDimensions(), ComputeDimensions(), G4DisplacedSolid::ComputeDimensions(), G4ScaledSolid::ComputeDimensions(), G4ParameterisedNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4DisplacedSolid::CreatePolyhedron(), G4ScaledSolid::CreatePolyhedron(), G4SubtractionSolid::DistanceToIn(), G4Box::DistanceToOut(), G4Orb::DistanceToOut(), G4Para::DistanceToOut(), G4Trap::DistanceToOut(), G4Trd::DistanceToOut(), G4Paraboloid::DistanceToOut(), G4VTwistedFaceted::DistanceToOut(), G4Cons::DistanceToOut(), G4CutTubs::DistanceToOut(), G4Sphere::DistanceToOut(), G4Torus::DistanceToOut(), G4Tubs::DistanceToOut(), G4Ellipsoid::DistanceToOut(), G4EllipticalCone::DistanceToOut(), G4EllipticalTube::DistanceToOut(), G4GenericTrap::DistanceToOut(), export_G4VSolid(), G4Polycone::G4Polycone(), G4Polyhedra::G4Polyhedra(), G4BooleanSolid::GetConstituentSolid(), G4NavigationLogger::PostComputeStepLog(), G4Box::SurfaceNormal(), G4Para::SurfaceNormal(), G4Trap::SurfaceNormal(), G4Trd::SurfaceNormal(), G4Ellipsoid::SurfaceNormal(), G4EllipticalCone::SurfaceNormal(), G4EllipticalTube::SurfaceNormal(), G4ExtrudedSolid::SurfaceNormal(), and G4Tet::SurfaceNormal().

◆ EstimateCubicVolume()

G4double G4VSolid::EstimateCubicVolume ( G4int  nStat,
G4double  epsilon 
) const
inherited

Definition at line 203 of file G4VSolid.cc.

204{
205 G4int iInside=0;
206 G4double px,py,pz,minX,maxX,minY,maxY,minZ,maxZ,volume,halfepsilon;
208 EInside in;
209
210 // values needed for CalculateExtent signature
211
212 G4VoxelLimits limit; // Unlimited
213 G4AffineTransform origin;
214
215 // min max extents of pSolid along X,Y,Z
216
217 CalculateExtent(kXAxis,limit,origin,minX,maxX);
218 CalculateExtent(kYAxis,limit,origin,minY,maxY);
219 CalculateExtent(kZAxis,limit,origin,minZ,maxZ);
220
221 // limits
222
223 if(nStat < 100) nStat = 100;
224 if(epsilon > 0.01) epsilon = 0.01;
225 halfepsilon = 0.5*epsilon;
226
227 for(auto i = 0; i < nStat; ++i )
228 {
229 px = minX-halfepsilon+(maxX-minX+epsilon)*G4QuickRand();
230 py = minY-halfepsilon+(maxY-minY+epsilon)*G4QuickRand();
231 pz = minZ-halfepsilon+(maxZ-minZ+epsilon)*G4QuickRand();
232 p = G4ThreeVector(px,py,pz);
233 in = Inside(p);
234 if(in != kOutside) ++iInside;
235 }
236 volume = (maxX-minX+epsilon)*(maxY-minY+epsilon)
237 * (maxZ-minZ+epsilon)*iInside/nStat;
238 return volume;
239}
G4double epsilon(G4double density, G4double temperature)
static const G4int maxZ
G4double G4QuickRand()
Definition: G4QuickRand.hh:34
CLHEP::Hep3Vector G4ThreeVector
virtual EInside Inside(const G4ThreeVector &p) const =0
EInside
Definition: geomdefs.hh:67
@ kOutside
Definition: geomdefs.hh:68

References G4VSolid::CalculateExtent(), epsilon(), G4QuickRand(), G4VSolid::Inside(), kOutside, kXAxis, kYAxis, kZAxis, and maxZ.

Referenced by G4VSolid::GetCubicVolume(), G4BooleanSolid::GetCubicVolume(), and G4VCSGfaceted::GetCubicVolume().

◆ EstimateSurfaceArea()

G4double G4VSolid::EstimateSurfaceArea ( G4int  nStat,
G4double  ell 
) const
inherited

Definition at line 265 of file G4VSolid.cc.

266{
267 static const G4double s2 = 1./std::sqrt(2.);
268 static const G4double s3 = 1./std::sqrt(3.);
269 static const G4ThreeVector directions[64] =
270 {
271 G4ThreeVector( 0, 0, 0), G4ThreeVector( -1, 0, 0), // ( , , ) ( -, , )
272 G4ThreeVector( 1, 0, 0), G4ThreeVector( -1, 0, 0), // ( +, , ) (-+, , )
273 G4ThreeVector( 0, -1, 0), G4ThreeVector(-s2,-s2, 0), // ( , -, ) ( -, -, )
274 G4ThreeVector( s2, -s2, 0), G4ThreeVector( 0, -1, 0), // ( +, -, ) (-+, -, )
275
276 G4ThreeVector( 0, 1, 0), G4ThreeVector( -s2, s2, 0), // ( , +, ) ( -, +, )
277 G4ThreeVector( s2, s2, 0), G4ThreeVector( 0, 1, 0), // ( +, +, ) (-+, +, )
278 G4ThreeVector( 0, -1, 0), G4ThreeVector( -1, 0, 0), // ( ,-+, ) ( -,-+, )
279 G4ThreeVector( 1, 0, 0), G4ThreeVector( -1, 0, 0), // ( +,-+, ) (-+,-+, )
280
281 G4ThreeVector( 0, 0, -1), G4ThreeVector(-s2, 0,-s2), // ( , , -) ( -, , -)
282 G4ThreeVector( s2, 0,-s2), G4ThreeVector( 0, 0, -1), // ( +, , -) (-+, , -)
283 G4ThreeVector( 0,-s2,-s2), G4ThreeVector(-s3,-s3,-s3), // ( , -, -) ( -, -, -)
284 G4ThreeVector( s3,-s3,-s3), G4ThreeVector( 0,-s2,-s2), // ( +, -, -) (-+, -, -)
285
286 G4ThreeVector( 0, s2,-s2), G4ThreeVector(-s3, s3,-s3), // ( , +, -) ( -, +, -)
287 G4ThreeVector( s3, s3,-s3), G4ThreeVector( 0, s2,-s2), // ( +, +, -) (-+, +, -)
288 G4ThreeVector( 0, 0, -1), G4ThreeVector(-s2, 0,-s2), // ( ,-+, -) ( -,-+, -)
289 G4ThreeVector( s2, 0,-s2), G4ThreeVector( 0, 0, -1), // ( +,-+, -) (-+,-+, -)
290
291 G4ThreeVector( 0, 0, 1), G4ThreeVector(-s2, 0, s2), // ( , , +) ( -, , +)
292 G4ThreeVector( s2, 0, s2), G4ThreeVector( 0, 0, 1), // ( +, , +) (-+, , +)
293 G4ThreeVector( 0,-s2, s2), G4ThreeVector(-s3,-s3, s3), // ( , -, +) ( -, -, +)
294 G4ThreeVector( s3,-s3, s3), G4ThreeVector( 0,-s2, s2), // ( +, -, +) (-+, -, +)
295
296 G4ThreeVector( 0, s2, s2), G4ThreeVector(-s3, s3, s3), // ( , +, +) ( -, +, +)
297 G4ThreeVector( s3, s3, s3), G4ThreeVector( 0, s2, s2), // ( +, +, +) (-+, +, +)
298 G4ThreeVector( 0, 0, 1), G4ThreeVector(-s2, 0, s2), // ( ,-+, +) ( -,-+, +)
299 G4ThreeVector( s2, 0, s2), G4ThreeVector( 0, 0, 1), // ( +,-+, +) (-+,-+, +)
300
301 G4ThreeVector( 0, 0, -1), G4ThreeVector( -1, 0, 0), // ( , ,-+) ( -, ,-+)
302 G4ThreeVector( 1, 0, 0), G4ThreeVector( -1, 0, 0), // ( +, ,-+) (-+, ,-+)
303 G4ThreeVector( 0, -1, 0), G4ThreeVector(-s2,-s2, 0), // ( , -,-+) ( -, -,-+)
304 G4ThreeVector( s2, -s2, 0), G4ThreeVector( 0, -1, 0), // ( +, -,-+) (-+, -,-+)
305
306 G4ThreeVector( 0, 1, 0), G4ThreeVector( -s2, s2, 0), // ( , +,-+) ( -, +,-+)
307 G4ThreeVector( s2, s2, 0), G4ThreeVector( 0, 1, 0), // ( +, +,-+) (-+, +,-+)
308 G4ThreeVector( 0, -1, 0), G4ThreeVector( -1, 0, 0), // ( ,-+,-+) ( -,-+,-+)
309 G4ThreeVector( 1, 0, 0), G4ThreeVector( -1, 0, 0), // ( +,-+,-+) (-+,-+,-+)
310 };
311
312 G4ThreeVector bmin, bmax;
313 BoundingLimits(bmin, bmax);
314
315 G4double dX = bmax.x() - bmin.x();
316 G4double dY = bmax.y() - bmin.y();
317 G4double dZ = bmax.z() - bmin.z();
318
319 // Define statistics and shell thickness
320 //
321 G4int npoints = (nstat < 1000) ? 1000 : nstat;
322 G4double coeff = 0.5 / std::cbrt(G4double(npoints));
323 G4double eps = (ell > 0) ? ell : coeff * std::min(std::min(dX, dY), dZ);
324 G4double del = 1.8 * eps; // shold be more than sqrt(3.)
325
326 G4double minX = bmin.x() - eps;
327 G4double minY = bmin.y() - eps;
328 G4double minZ = bmin.z() - eps;
329
330 G4double dd = 2. * eps;
331 dX += dd;
332 dY += dd;
333 dZ += dd;
334
335 // Calculate surface area
336 //
337 G4int icount = 0;
338 for(auto i = 0; i < npoints; ++i)
339 {
340 G4double px = minX + dX*G4QuickRand();
341 G4double py = minY + dY*G4QuickRand();
342 G4double pz = minZ + dZ*G4QuickRand();
343 G4ThreeVector p = G4ThreeVector(px, py, pz);
344 EInside in = Inside(p);
345 G4double dist = 0;
346 if (in == kInside)
347 {
348 if (DistanceToOut(p) >= eps) continue;
349 G4int icase = 0;
350 if (Inside(G4ThreeVector(px-del, py, pz)) != kInside) icase += 1;
351 if (Inside(G4ThreeVector(px+del, py, pz)) != kInside) icase += 2;
352 if (Inside(G4ThreeVector(px, py-del, pz)) != kInside) icase += 4;
353 if (Inside(G4ThreeVector(px, py+del, pz)) != kInside) icase += 8;
354 if (Inside(G4ThreeVector(px, py, pz-del)) != kInside) icase += 16;
355 if (Inside(G4ThreeVector(px, py, pz+del)) != kInside) icase += 32;
356 if (icase == 0) continue;
357 G4ThreeVector v = directions[icase];
358 dist = DistanceToOut(p, v);
359 G4ThreeVector n = SurfaceNormal(p + v*dist);
360 dist *= v.dot(n);
361 }
362 else if (in == kOutside)
363 {
364 if (DistanceToIn(p) >= eps) continue;
365 G4int icase = 0;
366 if (Inside(G4ThreeVector(px-del, py, pz)) != kOutside) icase += 1;
367 if (Inside(G4ThreeVector(px+del, py, pz)) != kOutside) icase += 2;
368 if (Inside(G4ThreeVector(px, py-del, pz)) != kOutside) icase += 4;
369 if (Inside(G4ThreeVector(px, py+del, pz)) != kOutside) icase += 8;
370 if (Inside(G4ThreeVector(px, py, pz-del)) != kOutside) icase += 16;
371 if (Inside(G4ThreeVector(px, py, pz+del)) != kOutside) icase += 32;
372 if (icase == 0) continue;
373 G4ThreeVector v = directions[icase];
374 dist = DistanceToIn(p, v);
375 if (dist == kInfinity) continue;
376 G4ThreeVector n = SurfaceNormal(p + v*dist);
377 dist *= -(v.dot(n));
378 }
379 if (dist < eps) ++icount;
380 }
381 return dX*dY*dZ*icount/npoints/dd;
382}
static const G4double eps
double z() const
double x() const
double y() const
double dot(const Hep3Vector &) const
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
@ kInside
Definition: geomdefs.hh:70
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References G4VSolid::BoundingLimits(), G4VSolid::DistanceToIn(), G4VSolid::DistanceToOut(), CLHEP::Hep3Vector::dot(), eps, G4QuickRand(), G4VSolid::Inside(), kInfinity, kInside, kOutside, G4INCL::Math::min(), CLHEP::detail::n, G4VSolid::SurfaceNormal(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by G4VSolid::GetSurfaceArea(), G4MultiUnion::GetSurfaceArea(), and G4VCSGfaceted::GetSurfaceArea().

◆ GetConstituentMovedSolid()

G4VSolid * G4ReflectedSolid::GetConstituentMovedSolid ( ) const

Definition at line 120 of file G4ReflectedSolid.cc.

121{
122 return fPtrSolid;
123}

References fPtrSolid.

Referenced by G4tgbGeometryDumper::DumpSolid().

◆ GetConstituentSolid() [1/2]

G4VSolid * G4VSolid::GetConstituentSolid ( G4int  no)
virtualinherited

Reimplemented in G4BooleanSolid.

Definition at line 170 of file G4VSolid.cc.

171{ return nullptr; }

◆ GetConstituentSolid() [2/2]

const G4VSolid * G4VSolid::GetConstituentSolid ( G4int  no) const
virtualinherited

Reimplemented in G4BooleanSolid.

Definition at line 167 of file G4VSolid.cc.

168{ return nullptr; }

Referenced by G4BooleanSolid::StackPolyhedron().

◆ GetCubicVolume()

G4double G4VSolid::GetCubicVolume ( )
virtualinherited

◆ GetDirectTransform3D()

G4Transform3D G4ReflectedSolid::GetDirectTransform3D ( ) const

Definition at line 133 of file G4ReflectedSolid.cc.

134{
135 G4Transform3D aTransform = *fDirectTransform3D;
136 return aTransform;
137}

References fDirectTransform3D.

◆ GetDisplacedSolidPtr() [1/2]

G4DisplacedSolid * G4VSolid::GetDisplacedSolidPtr ( )
virtualinherited

Reimplemented in G4DisplacedSolid.

Definition at line 176 of file G4VSolid.cc.

177{ return nullptr; }

◆ GetDisplacedSolidPtr() [2/2]

const G4DisplacedSolid * G4VSolid::GetDisplacedSolidPtr ( ) const
virtualinherited

Reimplemented in G4DisplacedSolid.

Definition at line 173 of file G4VSolid.cc.

174{ return nullptr; }

◆ GetEntityType()

G4GeometryType G4ReflectedSolid::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 105 of file G4ReflectedSolid.cc.

106{
107 return G4String("G4ReflectedSolid");
108}

Referenced by StreamInfo().

◆ GetExtent()

G4VisExtent G4VSolid::GetExtent ( ) const
virtualinherited

Reimplemented in G4Box, G4Orb, G4Sphere, G4Ellipsoid, G4EllipticalCone, G4EllipticalTube, G4GenericTrap, G4Hype, G4TessellatedSolid, G4Tet, G4TwistedTubs, G4VCSGfaceted, and G4VTwistedFaceted.

Definition at line 682 of file G4VSolid.cc.

683{
684 G4VisExtent extent;
685 G4VoxelLimits voxelLimits; // Defaults to "infinite" limits.
686 G4AffineTransform affineTransform;
687 G4double vmin, vmax;
688 CalculateExtent(kXAxis,voxelLimits,affineTransform,vmin,vmax);
689 extent.SetXmin (vmin);
690 extent.SetXmax (vmax);
691 CalculateExtent(kYAxis,voxelLimits,affineTransform,vmin,vmax);
692 extent.SetYmin (vmin);
693 extent.SetYmax (vmax);
694 CalculateExtent(kZAxis,voxelLimits,affineTransform,vmin,vmax);
695 extent.SetZmin (vmin);
696 extent.SetZmax (vmax);
697 return extent;
698}
void SetYmin(G4double ymin)
Definition: G4VisExtent.hh:114
void SetYmax(G4double ymax)
Definition: G4VisExtent.hh:116
void SetXmax(G4double xmax)
Definition: G4VisExtent.hh:112
void SetXmin(G4double xmin)
Definition: G4VisExtent.hh:110
void SetZmax(G4double zmax)
Definition: G4VisExtent.hh:120
void SetZmin(G4double zmin)
Definition: G4VisExtent.hh:118

References G4VSolid::CalculateExtent(), kXAxis, kYAxis, kZAxis, G4VisExtent::SetXmax(), G4VisExtent::SetXmin(), G4VisExtent::SetYmax(), G4VisExtent::SetYmin(), G4VisExtent::SetZmax(), and G4VisExtent::SetZmin().

Referenced by G4tgbVolume::BuildSolidForDivision(), G4BoundingExtentScene::ProcessVolume(), G4BoundingSphereScene::ProcessVolume(), and G4VisCommandsTouchable::SetNewValue().

◆ GetName()

G4String G4VSolid::GetName ( ) const
inlineinherited

Referenced by G4GMocrenFileSceneHandler::AddDetector(), G4HepRepFileSceneHandler::AddHepRepInstance(), G4GMocrenFileSceneHandler::AddPrimitive(), G4HepRepFileSceneHandler::AddSolid(), G4GMocrenFileSceneHandler::AddSolid(), G4VtkSceneHandler::AddSolid(), G4GDMLWriteSolids::AddSolid(), G4NavigationLogger::AlongComputeStepLog(), G4GDMLWriteSolids::BooleanWrite(), BoundingLimits(), G4DisplacedSolid::BoundingLimits(), G4IntersectionSolid::BoundingLimits(), G4ScaledSolid::BoundingLimits(), G4SubtractionSolid::BoundingLimits(), G4UnionSolid::BoundingLimits(), G4Box::BoundingLimits(), G4Cons::BoundingLimits(), G4CutTubs::BoundingLimits(), G4Orb::BoundingLimits(), G4Para::BoundingLimits(), G4Sphere::BoundingLimits(), G4Torus::BoundingLimits(), G4Trap::BoundingLimits(), G4Trd::BoundingLimits(), G4Tubs::BoundingLimits(), G4EllipticalCone::BoundingLimits(), G4ExtrudedSolid::BoundingLimits(), G4GenericPolycone::BoundingLimits(), G4GenericTrap::BoundingLimits(), G4Hype::BoundingLimits(), G4Paraboloid::BoundingLimits(), G4Polycone::BoundingLimits(), G4Polyhedra::BoundingLimits(), G4TessellatedSolid::BoundingLimits(), G4TwistedTubs::BoundingLimits(), G4GDMLWriteSolids::BoxWrite(), G4ExtrudedSolid::CalculateExtent(), G4GenericPolycone::CalculateExtent(), G4Polycone::CalculateExtent(), G4Polyhedra::CalculateExtent(), G4NavigationLogger::CheckDaughterEntryPoint(), G4VDivisionParameterisation::CheckNDivAndWidth(), G4VDivisionParameterisation::CheckOffset(), G4GenericTrap::CheckOrder(), G4Para::CheckParameters(), G4Trap::CheckParameters(), G4Trd::CheckParameters(), G4Ellipsoid::CheckParameters(), G4EllipticalTube::CheckParameters(), G4ParameterisationPolyconeRho::CheckParametersValidity(), G4ParameterisationPolyconeZ::CheckParametersValidity(), G4ParameterisationPolyhedraRho::CheckParametersValidity(), G4ParameterisationPolyhedraPhi::CheckParametersValidity(), G4ParameterisationPolyhedraZ::CheckParametersValidity(), G4PhantomParameterisation::CheckVoxelsFillContainer(), G4GenericTrap::ComputeIsTwisted(), G4VoxelNavigation::ComputeSafety(), G4VoxelSafety::ComputeSafety(), G4NavigationLogger::ComputeSafetyLog(), G4ParameterisedNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4GDMLWriteSolids::ConeWrite(), G4Polyhedra::Create(), G4GenericPolycone::Create(), G4Polycone::Create(), G4PhysicalVolumeModel::CreateCurrentAttValues(), CreatePolyhedron(), G4ReflectionFactory::CreateReflectedLV(), G4GenericTrap::CreateTessellatedSolid(), G4GDMLWriteSolids::CutTubeWrite(), G4SolidStore::DeRegister(), G4PhysicalVolumeModel::DescribeSolid(), G4SubtractionSolid::DistanceToIn(), G4Paraboloid::DistanceToIn(), G4TessellatedSolid::DistanceToIn(), G4Box::DistanceToOut(), G4Orb::DistanceToOut(), G4Para::DistanceToOut(), G4Trap::DistanceToOut(), G4Trd::DistanceToOut(), G4EllipticalCone::DistanceToOut(), G4TessellatedSolid::DistanceToOut(), G4Ellipsoid::DistanceToOut(), G4EllipticalTube::DistanceToOut(), G4tgbGeometryDumper::DumpMultiUnionVolume(), G4tgbGeometryDumper::DumpScaledVolume(), G4tgbGeometryDumper::DumpSolid(), G4GDMLWriteSolids::ElconeWrite(), G4GDMLWriteSolids::EllipsoidWrite(), G4GDMLWriteSolids::EltubeWrite(), G4PVDivision::ErrorInAxis(), G4ReplicatedSlice::ErrorInAxis(), export_G4VSolid(), G4Box::G4Box(), G4Cons::G4Cons(), G4CutTubs::G4CutTubs(), G4EllipticalCone::G4EllipticalCone(), G4Hype::G4Hype(), G4Para::G4Para(), G4Paraboloid::G4Paraboloid(), G4Polycone::G4Polycone(), G4Polyhedra::G4Polyhedra(), G4Sphere::G4Sphere(), G4Tet::G4Tet(), G4Trap::G4Trap(), G4Tubs::G4Tubs(), G4VParameterisationCons::G4VParameterisationCons(), G4VParameterisationPara::G4VParameterisationPara(), G4VParameterisationPolycone::G4VParameterisationPolycone(), G4VParameterisationPolyhedra::G4VParameterisationPolyhedra(), G4VParameterisationTrd::G4VParameterisationTrd(), G4VTwistedFaceted::G4VTwistedFaceted(), G4GDMLWriteSolids::GenericPolyconeWrite(), G4GDMLWriteSolids::GenTrapWrite(), G4Navigator::GetGlobalExitNormal(), G4Navigator::GetLocalExitNormal(), G4ITNavigator1::GetLocalExitNormal(), G4ITNavigator2::GetLocalExitNormal(), G4BooleanSolid::GetPointOnSurface(), G4PhantomParameterisation::GetReplicaNo(), G4GDMLWriteSolids::HypeWrite(), G4TessellatedSolid::InsideNoVoxels(), G4TessellatedSolid::InsideVoxels(), G4ITNavigator1::LocateGlobalPointAndSetup(), G4ITNavigator2::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointAndSetup(), G4GenericTrap::MakeDownFacet(), G4Trap::MakePlanes(), G4GenericTrap::MakeUpFacet(), G4GDMLWriteSolids::MultiUnionWrite(), G4GDMLWriteSolids::OrbWrite(), G4GDMLWriteSolids::ParaboloidWrite(), G4GDMLWriteParamvol::ParametersWrite(), G4GDMLWriteSolids::ParaWrite(), G4GDMLWriteSolids::PolyconeWrite(), G4GDMLWriteSolids::PolyhedraWrite(), G4NavigationLogger::PostComputeStepLog(), G4NavigationLogger::PreComputeStepLog(), G4NavigationLogger::PrintDaughterLog(), G4PseudoScene::ProcessVolume(), G4SolidStore::Register(), G4tgbVolumeMgr::RegisterMe(), G4NavigationLogger::ReportOutsideMother(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4VSceneHandler::RequestPrimitives(), G4GenericPolycone::Reset(), G4Polyhedra::Reset(), G4VoxelSafety::SafetyForVoxelNode(), G4GDMLWriteSolids::ScaledWrite(), G4Torus::SetAllParameters(), G4Tet::SetBoundingLimits(), G4Polycone::SetOriginalParameters(), G4Polyhedra::SetOriginalParameters(), G4TessellatedSolid::SetSolidClosed(), G4Tet::SetVertices(), G4Box::SetXHalfLength(), G4Box::SetYHalfLength(), G4Box::SetZHalfLength(), G4GDMLWriteSolids::SphereWrite(), G4BooleanSolid::StackPolyhedron(), StreamInfo(), G4BooleanSolid::StreamInfo(), G4DisplacedSolid::StreamInfo(), G4MultiUnion::StreamInfo(), G4ScaledSolid::StreamInfo(), G4Box::StreamInfo(), G4Cons::StreamInfo(), G4CSGSolid::StreamInfo(), G4CutTubs::StreamInfo(), G4Orb::StreamInfo(), G4Para::StreamInfo(), G4Sphere::StreamInfo(), G4Torus::StreamInfo(), G4Trap::StreamInfo(), G4Trd::StreamInfo(), G4Tubs::StreamInfo(), G4Ellipsoid::StreamInfo(), G4EllipticalCone::StreamInfo(), G4EllipticalTube::StreamInfo(), G4ExtrudedSolid::StreamInfo(), G4GenericPolycone::StreamInfo(), G4GenericTrap::StreamInfo(), G4Hype::StreamInfo(), G4Paraboloid::StreamInfo(), G4Polycone::StreamInfo(), G4Polyhedra::StreamInfo(), G4TessellatedSolid::StreamInfo(), G4Tet::StreamInfo(), G4TwistedBox::StreamInfo(), G4TwistedTrap::StreamInfo(), G4TwistedTrd::StreamInfo(), G4TwistedTubs::StreamInfo(), G4VCSGfaceted::StreamInfo(), G4VTwistedFaceted::StreamInfo(), G4GDMLRead::StripNames(), SubstractSolids(), G4UnionSolid::SurfaceNormal(), G4Box::SurfaceNormal(), G4Para::SurfaceNormal(), G4Trap::SurfaceNormal(), G4Trd::SurfaceNormal(), G4Ellipsoid::SurfaceNormal(), G4EllipticalCone::SurfaceNormal(), G4EllipticalTube::SurfaceNormal(), G4ExtrudedSolid::SurfaceNormal(), G4Tet::SurfaceNormal(), G4GDMLWriteSolids::TessellatedWrite(), G4GDMLWriteSolids::TetWrite(), G4GDMLWriteSolids::TorusWrite(), G4GDMLWriteSolids::TrapWrite(), G4GDMLWriteStructure::TraverseVolumeTree(), G4GDMLWriteSolids::TrdWrite(), G4GDMLWriteSolids::TubeWrite(), G4GDMLWriteSolids::TwistedboxWrite(), G4GDMLWriteSolids::TwistedtrapWrite(), G4GDMLWriteSolids::TwistedtrdWrite(), G4GDMLWriteSolids::TwistedtubsWrite(), G4PhysicalVolumeModel::VisitGeometryAndGetVisReps(), and G4GDMLWriteSolids::XtruWrite().

◆ GetPointOnSurface()

G4ThreeVector G4ReflectedSolid::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 356 of file G4ReflectedSolid.cc.

357{
359 return (*fDirectTransform3D)*G4Point3D(p);
360}
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152

References fDirectTransform3D, fPtrSolid, and G4VSolid::GetPointOnSurface().

◆ GetPolyhedron()

G4Polyhedron * G4ReflectedSolid::GetPolyhedron ( ) const
virtual

◆ GetReflectedSolidPtr() [1/2]

G4ReflectedSolid * G4ReflectedSolid::GetReflectedSolidPtr ( )
virtual

Definition at line 115 of file G4ReflectedSolid.cc.

116{
117 return this;
118}

◆ GetReflectedSolidPtr() [2/2]

const G4ReflectedSolid * G4ReflectedSolid::GetReflectedSolidPtr ( ) const
virtual

Definition at line 110 of file G4ReflectedSolid.cc.

111{
112 return this;
113}

◆ GetSurfaceArea()

G4double G4VSolid::GetSurfaceArea ( )
virtualinherited

◆ GetTolerance()

G4double G4VSolid::GetTolerance ( ) const
inlineinherited

◆ GetTransform3D()

G4Transform3D G4ReflectedSolid::GetTransform3D ( ) const

Definition at line 128 of file G4ReflectedSolid.cc.

129{
130 return fDirectTransform3D->inverse();
131}
Transform3D inverse() const
Definition: Transform3D.cc:141

References fDirectTransform3D, and HepGeom::Transform3D::inverse().

◆ Inside()

EInside G4ReflectedSolid::Inside ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 257 of file G4ReflectedSolid.cc.

258{
259 G4ThreeVector newPoint = (*fDirectTransform3D)*G4Point3D(p);
260 return fPtrSolid->Inside(newPoint);
261}

References fPtrSolid, and G4VSolid::Inside().

◆ operator=()

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

Definition at line 81 of file G4ReflectedSolid.cc.

82{
83 // Check assignment to self
84 //
85 if (this == &rhs) { return *this; }
86
87 // Copy base class data
88 //
90
91 // Copy data
92 //
93 fPtrSolid = rhs.fPtrSolid;
94 delete fDirectTransform3D;
96 fRebuildPolyhedron = false;
97 delete fpPolyhedron; fpPolyhedron = nullptr;
98
99 return *this;
100}
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:107

References fDirectTransform3D, fpPolyhedron, fPtrSolid, fRebuildPolyhedron, and G4VSolid::operator=().

◆ operator==()

G4bool G4VSolid::operator== ( const G4VSolid s) const
inlineinherited

◆ SetDirectTransform3D()

void G4ReflectedSolid::SetDirectTransform3D ( G4Transform3D transform)

◆ SetName()

void G4VSolid::SetName ( const G4String name)
inherited

Definition at line 127 of file G4VSolid.cc.

128{
131}
void SetMapValid(G4bool val)
Definition: G4SolidStore.hh:76
static G4SolidStore * GetInstance()
G4String fshapeName
Definition: G4VSolid.hh:312
const char * name(G4int ptype)

References G4VSolid::fshapeName, G4SolidStore::GetInstance(), G4InuclParticleNames::name(), and G4SolidStore::SetMapValid().

Referenced by export_G4VSolid(), G4MultiUnion::G4MultiUnion(), and G4GDMLRead::StripNames().

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 376 of file G4ReflectedSolid.cc.

377{
378 os << "-----------------------------------------------------------\n"
379 << " *** Dump for Reflected solid - " << GetName() << " ***\n"
380 << " ===================================================\n"
381 << " Solid type: " << GetEntityType() << "\n"
382 << " Parameters of constituent solid: \n"
383 << "===========================================================\n";
385 os << "===========================================================\n"
386 << " Transformations: \n"
387 << " Direct transformation - translation : \n"
388 << " " << fDirectTransform3D->getTranslation() << "\n"
389 << " - rotation : \n"
390 << " ";
392 os << "\n"
393 << "===========================================================\n";
394
395 return os;
396}
std::ostream & print(std::ostream &os) const
Definition: RotationIO.cc:17
virtual G4GeometryType GetEntityType() const
virtual std::ostream & StreamInfo(std::ostream &os) const =0

References fDirectTransform3D, fPtrSolid, GetEntityType(), G4VSolid::GetName(), HepGeom::Transform3D::getRotation(), HepGeom::Transform3D::getTranslation(), CLHEP::HepRotation::print(), and G4VSolid::StreamInfo().

◆ SurfaceNormal()

G4ThreeVector G4ReflectedSolid::SurfaceNormal ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 268 of file G4ReflectedSolid.cc.

269{
270 G4ThreeVector newPoint = (*fDirectTransform3D)*G4Point3D(p);
272 return (*fDirectTransform3D)*normal;
273}
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:79

References fDirectTransform3D, fPtrSolid, CLHEP::normal(), and G4VSolid::SurfaceNormal().

Field Documentation

◆ fDirectTransform3D

G4Transform3D* G4ReflectedSolid::fDirectTransform3D = nullptr
protected

◆ fpPolyhedron

G4Polyhedron* G4ReflectedSolid::fpPolyhedron = nullptr
mutableprotected

Definition at line 125 of file G4ReflectedSolid.hh.

Referenced by GetPolyhedron(), operator=(), and ~G4ReflectedSolid().

◆ fPtrSolid

G4VSolid* G4ReflectedSolid::fPtrSolid = nullptr
protected

◆ fRebuildPolyhedron

G4bool G4ReflectedSolid::fRebuildPolyhedron = false
mutableprotected

Definition at line 124 of file G4ReflectedSolid.hh.

Referenced by GetPolyhedron(), operator=(), and SetDirectTransform3D().

◆ fshapeName

G4String G4VSolid::fshapeName
privateinherited

Definition at line 312 of file G4VSolid.hh.

Referenced by G4VSolid::operator=(), and G4VSolid::SetName().

◆ kCarTolerance

G4double G4VSolid::kCarTolerance
protectedinherited

Definition at line 299 of file G4VSolid.hh.

Referenced by G4TessellatedSolid::AddFacet(), G4Polycone::CalculateExtent(), G4Polyhedra::CalculateExtent(), G4Tet::CheckDegeneracy(), G4Para::CheckParameters(), G4Trd::CheckParameters(), G4Ellipsoid::CheckParameters(), G4EllipticalTube::CheckParameters(), G4GenericTrap::ComputeIsTwisted(), G4Polyhedra::Create(), G4GenericPolycone::Create(), G4Polycone::Create(), G4CutTubs::CreatePolyhedron(), G4TessellatedSolid::CreateVertexList(), G4VCSGfaceted::DistanceTo(), G4Sphere::DistanceToIn(), G4Ellipsoid::DistanceToIn(), G4Hype::DistanceToIn(), G4Paraboloid::DistanceToIn(), G4VCSGfaceted::DistanceToIn(), G4TessellatedSolid::DistanceToInCore(), G4Cons::DistanceToOut(), G4CutTubs::DistanceToOut(), G4Sphere::DistanceToOut(), G4Torus::DistanceToOut(), G4Tubs::DistanceToOut(), G4GenericTrap::DistanceToOut(), G4Hype::DistanceToOut(), G4Paraboloid::DistanceToOut(), G4VCSGfaceted::DistanceToOut(), G4TessellatedSolid::DistanceToOutCandidates(), G4TessellatedSolid::DistanceToOutCore(), G4TessellatedSolid::DistanceToOutNoVoxels(), G4GenericTrap::DistToPlane(), G4GenericTrap::DistToTriangle(), G4Box::G4Box(), G4Cons::G4Cons(), G4CutTubs::G4CutTubs(), G4EllipticalCone::G4EllipticalCone(), G4ExtrudedSolid::G4ExtrudedSolid(), G4GenericTrap::G4GenericTrap(), G4Hype::G4Hype(), G4Para::G4Para(), G4Sphere::G4Sphere(), G4Tet::G4Tet(), G4Trap::G4Trap(), G4Tubs::G4Tubs(), G4UnionSolid::G4UnionSolid(), G4VSolid::G4VSolid(), G4VTwistedFaceted::G4VTwistedFaceted(), G4GenericPolycone::GetPointOnSurface(), G4Polycone::GetPointOnSurface(), G4UnionSolid::Init(), G4Orb::Initialize(), G4TessellatedSolid::Initialize(), G4SubtractionSolid::Inside(), G4Hype::Inside(), G4Paraboloid::Inside(), G4VCSGfaceted::Inside(), G4VTwistedFaceted::Inside(), G4TessellatedSolid::InsideNoVoxels(), G4GenericTrap::InsidePolygone(), G4TessellatedSolid::InsideVoxels(), G4CutTubs::IsCrossingCutPlanes(), G4GenericTrap::IsSegCrossingZ(), G4Trap::MakePlane(), G4GenericTrap::NormalToPlane(), G4VSolid::operator=(), G4TessellatedSolid::SafetyFromInside(), G4TessellatedSolid::SafetyFromOutside(), G4Torus::SetAllParameters(), G4Polycone::SetOriginalParameters(), G4Polyhedra::SetOriginalParameters(), G4Box::SetXHalfLength(), G4Box::SetYHalfLength(), G4Box::SetZHalfLength(), G4Torus::SurfaceNormal(), G4GenericTrap::SurfaceNormal(), and G4Paraboloid::SurfaceNormal().


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