Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ScaledSolid.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// Implementation for G4ScaledSolid class
27//
28// 27.10.15 G.Cosmo: created, based on implementation also provided in Root
29// --------------------------------------------------------------------
30
31#include "G4ScaledSolid.hh"
32#include "G4BoundingEnvelope.hh"
33
35
36#include "G4ScaleTransform.hh"
37
38#include "G4VGraphicsScene.hh"
39#include "G4Polyhedron.hh"
40
42//
43// Constructor
44//
46 G4VSolid* pSolid,
47 const G4Scale3D& pScale )
48 : G4VSolid(pName), fPtrSolid(pSolid)
49{
50 fScale = new G4ScaleTransform(pScale);
51}
52
54//
55// Fake default constructor - sets only member data and allocates memory
56// for usage restricted to object persistency
57//
59 : G4VSolid(a)
60{
61}
62
64//
65// Destructor
66//
68{
69 delete fpPolyhedron; fpPolyhedron = nullptr;
70 delete fScale; fScale = nullptr;
71}
72
74//
75// Copy constructor
76//
78 : G4VSolid (rhs), fPtrSolid(rhs.fPtrSolid),
79 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea)
80{
81 fScale = new G4ScaleTransform(*(rhs.fScale));
82}
83
85//
86// Assignment operator
87//
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));
105 fRebuildPolyhedron = false;
106 delete fpPolyhedron; fpPolyhedron = nullptr;
107
108 return *this;
109}
110
112//
113// Return original solid not scaled
114//
116{
117 return fPtrSolid;
118}
119
121//
122// Get bounding box
123
125 G4ThreeVector& pMax) const
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}
148
150//
151// Calculate extent under transform and specified limit
152//
153G4bool
155 const G4VoxelLimits& pVoxelLimit,
156 const G4AffineTransform& pTransform,
157 G4double& pMin,
158 G4double& pMax ) const
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}
173
175//
176// Inside
177//
179{
180 return fPtrSolid->Inside(fScale->Transform(p));
181}
182
184//
185// SurfaceNormal
186//
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);
197
198 // Convert normal to scaled frame
200 return normal/normal.mag();
201}
202
204//
205// The same algorithm as in DistanceToIn(p)
206//
209 const G4ThreeVector& v ) const
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}
226
228//
229// Approximate nearest distance from the point p to the solid from outside
230//
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}
242
244//
245// The same algorithm as DistanceToOut(p)
246//
249 const G4ThreeVector& v,
250 const G4bool calcNorm,
251 G4bool *validNorm,
252 G4ThreeVector *n ) const
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 {
270 fScale->TransformNormal(solNorm, normal);
271 *n = normal.unit();
272 }
273
274 // Return distance converted to global
275 return fScale->InverseTransformDistance(dist, newDirection);
276}
277
279//
280// Approximate nearest distance from the point p to the solid from inside
281//
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}
293
295//
296// ComputeDimensions
297//
298void
300 const G4int,
301 const G4VPhysicalVolume* )
302{
303 DumpInfo();
304 G4Exception("G4ScaledSolid::ComputeDimensions()",
305 "GeomSolids0001", FatalException,
306 "Method not applicable in this context!");
307}
308
310//
311// Returns a point (G4ThreeVector) randomly and uniformly selected
312// on the solid surface
313//
315{
317}
318
320//
321// Return object type name
322//
324{
325 return G4String("G4ScaledSolid");
326}
327
329//
330// Make a clone of the object
331//
333{
334 return new G4ScaledSolid(*this);
335}
336
338//
339// Returning the scaling transformation
340//
342{
343 return G4Scale3D(fScale->GetScale().x(),
344 fScale->GetScale().y(),
345 fScale->GetScale().z());
346}
347
349//
350// Setting the scaling transformation
351//
353{
354 if (fScale != nullptr) { delete fScale; }
355 fScale = new G4ScaleTransform(scale);
356 fRebuildPolyhedron = true;
357}
358
360//
361// Get volume of the scaled solid
362//
364{
365 if(fCubicVolume < 0.)
366 {
368 fScale->GetScale().x() *
369 fScale->GetScale().y() *
370 fScale->GetScale().z();
371 }
372 return fCubicVolume;
373}
374
376//
377// Get estimated surface area of the scaled solid
378//
380{
381 if(fSurfaceArea < 0.)
382 {
384 }
385 return fSurfaceArea;
386}
387
389//
390// Stream object contents to an output stream
391//
392std::ostream& G4ScaledSolid::StreamInfo(std::ostream& os) const
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";
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}
411
413//
414// DescribeYourselfTo
415//
416void
418{
419 scene.AddSolid (*this);
420}
421
423//
424// CreatePolyhedron
425//
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}
443
445//
446// GetPolyhedron
447//
449{
450 if (fpPolyhedron == nullptr ||
454 {
456 fRebuildPolyhedron = false;
457 }
458 return fpPolyhedron;
459}
@ JustWarning
@ FatalException
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::Transform3D G4Transform3D
HepGeom::Scale3D G4Scale3D
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
double z() const
double x() const
double y() const
double mag() const
HepRotation inverse() const
G4ThreeVector NetTranslation() const
G4RotationMatrix NetRotation() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
void TransformNormal(const G4ThreeVector &global, G4ThreeVector &local) const
void InverseTransformNormal(const G4ThreeVector &local, G4ThreeVector &global) const
void InverseTransform(const G4ThreeVector &local, G4ThreeVector &global) const
const G4ThreeVector & GetScale() const
void Transform(const G4ThreeVector &global, G4ThreeVector &local) const
G4double InverseTransformDistance(G4double dist, const G4ThreeVector &dir) const
G4VSolid * GetUnscaledSolid() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
G4double GetCubicVolume()
G4ScaledSolid & operator=(const G4ScaledSolid &rhs)
G4double GetSurfaceArea()
G4VSolid * Clone() const
G4Scale3D GetScaleTransform() const
G4Polyhedron * fpPolyhedron
std::ostream & StreamInfo(std::ostream &os) const
G4ThreeVector GetPointOnSurface() const
G4double fSurfaceArea
G4VSolid * fPtrSolid
EInside Inside(const G4ThreeVector &p) const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4double fCubicVolume
void SetScaleTransform(const G4Scale3D &scale)
G4Polyhedron * CreatePolyhedron() const
G4ScaledSolid(const G4String &pName, G4VSolid *pSolid, const G4Scale3D &pScale)
G4GeometryType GetEntityType() const
G4bool fRebuildPolyhedron
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
G4ScaleTransform * fScale
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
virtual ~G4ScaledSolid()
G4Polyhedron * GetPolyhedron() const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
virtual void AddSolid(const G4Box &)=0
G4String GetName() const
virtual std::ostream & StreamInfo(std::ostream &os) const =0
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
void DumpInfo() const
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4VSolid.cc:665
virtual G4Polyhedron * CreatePolyhedron() const
Definition: G4VSolid.cc:700
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:107
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:250
static G4int GetNumberOfRotationSteps()
HepPolyhedron & Transform(const G4Transform3D &t)
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:79