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

#include <G4ScaledSolid.hh>

+ Inheritance diagram for G4ScaledSolid:

Public Member Functions

 G4ScaledSolid (const G4String &pName, G4VSolid *pSolid, const G4Scale3D &pScale)
 
virtual ~G4ScaledSolid ()
 
EInside Inside (const G4ThreeVector &p) const
 
void BoundingLimits (G4ThreeVector &pMin, G4ThreeVector &pMax) const
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
G4double DistanceToIn (const G4ThreeVector &p) const
 
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
 
G4double DistanceToOut (const G4ThreeVector &p) const
 
void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
G4double GetCubicVolume ()
 
G4double GetSurfaceArea ()
 
void CleanTransformations ()
 
G4ThreeVector GetPointOnSurface () const
 
G4Scale3D GetScaleTransform () const
 
void SetScaleTransform (const G4Scale3D &scale)
 
G4VSolidGetUnscaledSolid () const
 
G4GeometryType GetEntityType () const
 
G4VSolidClone () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
 G4ScaledSolid (__void__ &)
 
 G4ScaledSolid (const G4ScaledSolid &rhs)
 
G4ScaledSolidoperator= (const G4ScaledSolid &rhs)
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
G4PolyhedronCreatePolyhedron () const
 
G4PolyhedronGetPolyhedron () const
 
- 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
 

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 45 of file G4ScaledSolid.hh.

Constructor & Destructor Documentation

◆ G4ScaledSolid() [1/3]

G4ScaledSolid::G4ScaledSolid ( const G4String pName,
G4VSolid pSolid,
const G4Scale3D pScale 
)

Definition at line 45 of file G4ScaledSolid.cc.

48 : G4VSolid(pName), fPtrSolid(pSolid)
49{
50 fScale = new G4ScaleTransform(pScale);
51}

◆ ~G4ScaledSolid()

G4ScaledSolid::~G4ScaledSolid ( )
virtual

Definition at line 67 of file G4ScaledSolid.cc.

68{
69 delete fpPolyhedron; fpPolyhedron = nullptr;
70 delete fScale; fScale = nullptr;
71}

◆ G4ScaledSolid() [2/3]

G4ScaledSolid::G4ScaledSolid ( __void__ &  a)

Definition at line 58 of file G4ScaledSolid.cc.

59 : G4VSolid(a)
60{
61}

◆ G4ScaledSolid() [3/3]

G4ScaledSolid::G4ScaledSolid ( const G4ScaledSolid rhs)

Definition at line 77 of file G4ScaledSolid.cc.

78 : G4VSolid (rhs), fPtrSolid(rhs.fPtrSolid),
79 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea)
80{
81 fScale = new G4ScaleTransform(*(rhs.fScale));
82}

Member Function Documentation

◆ BoundingLimits()

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

Reimplemented from G4VSolid.

Definition at line 124 of file G4ScaledSolid.cc.

126{
127 G4ThreeVector bmin,bmax;
128 G4ThreeVector scale = fScale->GetScale();
129
130 fPtrSolid->BoundingLimits(bmin,bmax);
131 pMin.set(bmin.x()*scale.x(),bmin.y()*scale.y(),bmin.z()*scale.z());
132 pMax.set(bmax.x()*scale.x(),bmax.y()*scale.y(),bmax.z()*scale.z());
133
134 // Check correctness of the bounding box
135 //
136 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
137 {
138 std::ostringstream message;
139 message << "Bad bounding box (min >= max) for solid: "
140 << GetName() << " !"
141 << "\npMin = " << pMin
142 << "\npMax = " << pMax;
143 G4Exception("G4ScaledSolid::BoundingLimits()", "GeomMgt0001",
144 JustWarning, message);
145 DumpInfo();
146 }
147}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double z() const
double x() const
double y() const
void set(double x, double y, double z)
const G4ThreeVector & GetScale() const
G4String GetName() const
void DumpInfo() const
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4VSolid.cc:653

◆ CalculateExtent()

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

Implements G4VSolid.

Definition at line 154 of file G4ScaledSolid.cc.

159{
160 // Find bounding box of unscaled solid
161 G4ThreeVector bmin,bmax;
162 fPtrSolid->BoundingLimits(bmin,bmax);
163
164 // Set combined transformation
165 G4Transform3D transform3D =
166 G4Transform3D(pTransform.NetRotation().inverse(),
167 pTransform.NetTranslation())*GetScaleTransform();
168
169 // Find extent
170 G4BoundingEnvelope bbox(bmin,bmax);
171 return bbox.CalculateExtent(pAxis,pVoxelLimit,transform3D,pMin,pMax);
172}
HepGeom::Transform3D G4Transform3D
HepRotation inverse() const
G4ThreeVector NetTranslation() const
G4RotationMatrix NetRotation() const
G4Scale3D GetScaleTransform() const

◆ CleanTransformations()

void G4ScaledSolid::CleanTransformations ( )

◆ Clone()

G4VSolid * G4ScaledSolid::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 332 of file G4ScaledSolid.cc.

333{
334 return new G4ScaledSolid(*this);
335}

◆ ComputeDimensions()

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

Reimplemented from G4VSolid.

Definition at line 299 of file G4ScaledSolid.cc.

302{
303 DumpInfo();
304 G4Exception("G4ScaledSolid::ComputeDimensions()",
305 "GeomSolids0001", FatalException,
306 "Method not applicable in this context!");
307}
@ FatalException

◆ CreatePolyhedron()

G4Polyhedron * G4ScaledSolid::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 427 of file G4ScaledSolid.cc.

428{
429 G4Polyhedron* polyhedron = fPtrSolid->CreatePolyhedron();
430 if (polyhedron != nullptr)
431 {
432 polyhedron->Transform(GetScaleTransform());
433 }
434 else
435 {
436 DumpInfo();
437 G4Exception("G4ScaledSolid::CreatePolyhedron()",
438 "GeomSolids2003", JustWarning,
439 "No G4Polyhedron for scaled solid");
440 }
441 return polyhedron;
442}
virtual G4Polyhedron * CreatePolyhedron() const
Definition: G4VSolid.cc:688
HepPolyhedron & Transform(const G4Transform3D &t)

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

void G4ScaledSolid::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 417 of file G4ScaledSolid.cc.

418{
419 scene.AddSolid (*this);
420}
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

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

Implements G4VSolid.

Definition at line 232 of file G4ScaledSolid.cc.

233{
234 // Transform point to unscaled shape frame
235 G4ThreeVector newPoint;
236 fScale->Transform(p, newPoint);
237
238 // Compute unscaled safety, then scale it
239 G4double dist = fPtrSolid->DistanceToIn(newPoint);
240 return fScale->InverseTransformDistance(dist);
241}
double G4double
Definition: G4Types.hh:83
void Transform(const G4ThreeVector &global, G4ThreeVector &local) const
G4double InverseTransformDistance(G4double dist, const G4ThreeVector &dir) const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0

◆ DistanceToIn() [2/2]

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

Implements G4VSolid.

Definition at line 208 of file G4ScaledSolid.cc.

210{
211 // Transform point and direction to unscaled shape frame
212 G4ThreeVector newPoint;
213 fScale->Transform(p, newPoint);
214
215 // Direction is un-normalized after scale transformation
216 G4ThreeVector newDirection;
217 fScale->Transform(v, newDirection);
218 newDirection = newDirection/newDirection.mag();
219
220 // Compute distance in unscaled system
221 G4double dist = fPtrSolid->DistanceToIn(newPoint,newDirection);
222
223 // Return converted distance to global
224 return fScale->InverseTransformDistance(dist, newDirection);
225}
double mag() const

◆ DistanceToOut() [1/2]

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

Implements G4VSolid.

Definition at line 283 of file G4ScaledSolid.cc.

284{
285 // Transform point to unscaled shape frame
286 G4ThreeVector newPoint;
287 fScale->Transform(p, newPoint);
288
289 // Compute unscaled safety, then scale it
290 G4double dist = fPtrSolid->DistanceToOut(newPoint);
291 return fScale->InverseTransformDistance(dist);
292}
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 G4ScaledSolid::DistanceToOut ( const G4ThreeVector p,
const G4ThreeVector v,
const G4bool  calcNorm = false,
G4bool validNorm = nullptr,
G4ThreeVector n = nullptr 
) const
virtual

Implements G4VSolid.

Definition at line 248 of file G4ScaledSolid.cc.

253{
254 // Transform point and direction to unscaled shape frame
255 G4ThreeVector newPoint;
256 fScale->Transform(p, newPoint);
257
258 // Direction is un-normalized after scale transformation
259 G4ThreeVector newDirection;
260 fScale->Transform(v, newDirection);
261 newDirection = newDirection/newDirection.mag();
262
263 // Compute distance in unscaled system
264 G4ThreeVector solNorm;
265 G4double dist = fPtrSolid->DistanceToOut(newPoint,newDirection,
266 calcNorm,validNorm,&solNorm);
267 if(calcNorm)
268 {
269 G4ThreeVector normal;
270 fScale->TransformNormal(solNorm, normal);
271 *n = normal/normal.mag();
272 }
273
274 // Return distance converted to global
275 return fScale->InverseTransformDistance(dist, newDirection);
276}
void TransformNormal(const G4ThreeVector &global, G4ThreeVector &local) const

◆ GetCubicVolume()

G4double G4ScaledSolid::GetCubicVolume ( )
virtual

Reimplemented from G4VSolid.

Definition at line 363 of file G4ScaledSolid.cc.

364{
365 if(fCubicVolume < 0.)
366 {
367 fCubicVolume = fPtrSolid->GetCubicVolume() *
368 fScale->GetScale().x() *
369 fScale->GetScale().y() *
370 fScale->GetScale().z();
371 }
372 return fCubicVolume;
373}
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:176

◆ GetEntityType()

G4GeometryType G4ScaledSolid::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 323 of file G4ScaledSolid.cc.

324{
325 return G4String("G4ScaledSolid");
326}

Referenced by StreamInfo().

◆ GetPointOnSurface()

G4ThreeVector G4ScaledSolid::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 314 of file G4ScaledSolid.cc.

315{
316 return fScale->InverseTransform(fPtrSolid->GetPointOnSurface());
317}
void InverseTransform(const G4ThreeVector &local, G4ThreeVector &global) const
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:140

◆ GetPolyhedron()

G4Polyhedron * G4ScaledSolid::GetPolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 448 of file G4ScaledSolid.cc.

449{
450 if (fpPolyhedron == nullptr ||
451 fRebuildPolyhedron ||
453 fpPolyhedron->GetNumberOfRotationSteps())
454 {
455 fpPolyhedron = CreatePolyhedron();
456 fRebuildPolyhedron = false;
457 }
458 return fpPolyhedron;
459}
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
G4Polyhedron * CreatePolyhedron() const
static G4int GetNumberOfRotationSteps()

◆ GetScaleTransform()

G4Scale3D G4ScaledSolid::GetScaleTransform ( ) const

Definition at line 341 of file G4ScaledSolid.cc.

342{
343 return G4Scale3D(fScale->GetScale().x(),
344 fScale->GetScale().y(),
345 fScale->GetScale().z());
346}
HepGeom::Scale3D G4Scale3D

Referenced by CalculateExtent(), CreatePolyhedron(), and G4GDMLWriteSolids::ScaledWrite().

◆ GetSurfaceArea()

G4double G4ScaledSolid::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 379 of file G4ScaledSolid.cc.

380{
381 if(fSurfaceArea < 0.)
382 {
383 fSurfaceArea = G4VSolid::GetSurfaceArea();
384 }
385 return fSurfaceArea;
386}
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:238

◆ GetUnscaledSolid()

G4VSolid * G4ScaledSolid::GetUnscaledSolid ( ) const

Definition at line 115 of file G4ScaledSolid.cc.

116{
117 return fPtrSolid;
118}

Referenced by G4GDMLWriteSolids::ScaledWrite().

◆ Inside()

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

Implements G4VSolid.

Definition at line 178 of file G4ScaledSolid.cc.

179{
180 return fPtrSolid->Inside(fScale->Transform(p));
181}
virtual EInside Inside(const G4ThreeVector &p) const =0

◆ operator=()

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

Definition at line 88 of file G4ScaledSolid.cc.

89{
90 // Check assignment to self
91 //
92 if (this == &rhs) { return *this; }
93
94 // Copy base class data
95 //
97
98 // Copy data
99 //
100 fPtrSolid = rhs.fPtrSolid;
101 delete fScale;
102 fScale = new G4ScaleTransform(*(rhs.fScale));
103 fCubicVolume = rhs.fCubicVolume;
104 fSurfaceArea = rhs.fSurfaceArea;
105 fRebuildPolyhedron = false;
106 delete fpPolyhedron; fpPolyhedron = nullptr;
107
108 return *this;
109}
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:98

◆ SetScaleTransform()

void G4ScaledSolid::SetScaleTransform ( const G4Scale3D scale)

Definition at line 352 of file G4ScaledSolid.cc.

353{
354 if (fScale != nullptr) { delete fScale; }
355 fScale = new G4ScaleTransform(scale);
356 fRebuildPolyhedron = true;
357}

◆ StreamInfo()

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

Implements G4VSolid.

Definition at line 392 of file G4ScaledSolid.cc.

393{
394 os << "-----------------------------------------------------------\n"
395 << " *** Dump for Scaled solid - " << GetName() << " ***\n"
396 << " ===================================================\n"
397 << " Solid type: " << GetEntityType() << "\n"
398 << " Parameters of constituent solid: \n"
399 << "===========================================================\n";
400 fPtrSolid->StreamInfo(os);
401 os << "===========================================================\n"
402 << " Scaling: \n"
403 << " Scale transformation : \n"
404 << " " << fScale->GetScale().x() << ", "
405 << fScale->GetScale().y() << ", "
406 << fScale->GetScale().z() << "\n"
407 << "===========================================================\n";
408
409 return os;
410}
G4GeometryType GetEntityType() const
virtual std::ostream & StreamInfo(std::ostream &os) const =0

◆ SurfaceNormal()

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

Implements G4VSolid.

Definition at line 188 of file G4ScaledSolid.cc.

189{
190 // Transform point to unscaled shape frame
191 G4ThreeVector newPoint;
192 fScale->Transform(p, newPoint);
193
194 // Compute normal in unscaled frame
195 G4ThreeVector newNormal = fPtrSolid->SurfaceNormal(newPoint);
196 G4ThreeVector normal;
197
198 // Convert normal to scaled frame
199 fScale->InverseTransformNormal(newNormal, normal);
200 return normal/normal.mag();
201}
void InverseTransformNormal(const G4ThreeVector &local, G4ThreeVector &global) const
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0

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