Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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#include "G4AutoLock.hh"
41
42namespace
43{
45}
46
47///////////////////////////////////////////////////////////////////
48//
49// Constructor
50//
52 G4VSolid* pSolid,
53 const G4Scale3D& pScale )
54 : G4VSolid(pName), fPtrSolid(pSolid)
55{
56 fScale = new G4ScaleTransform(pScale);
57}
58
59///////////////////////////////////////////////////////////////////
60//
61// Fake default constructor - sets only member data and allocates memory
62// for usage restricted to object persistency
63//
65 : G4VSolid(a)
66{
67}
68
69///////////////////////////////////////////////////////////////////
70//
71// Destructor
72//
74{
75 delete fpPolyhedron; fpPolyhedron = nullptr;
76 delete fScale; fScale = nullptr;
77}
78
79///////////////////////////////////////////////////////////////
80//
81// Copy constructor
82//
84 : G4VSolid (rhs), fPtrSolid(rhs.fPtrSolid),
85 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea)
86{
87 fScale = new G4ScaleTransform(*(rhs.fScale));
88}
89
90///////////////////////////////////////////////////////////////
91//
92// Assignment operator
93//
95{
96 // Check assignment to self
97 //
98 if (this == &rhs) { return *this; }
99
100 // Copy base class data
101 //
103
104 // Copy data
105 //
106 fPtrSolid = rhs.fPtrSolid;
107 delete fScale;
108 fScale = new G4ScaleTransform(*(rhs.fScale));
109 fCubicVolume = rhs.fCubicVolume;
110 fSurfaceArea = rhs.fSurfaceArea;
111 fRebuildPolyhedron = false;
112 delete fpPolyhedron; fpPolyhedron = nullptr;
113
114 return *this;
115}
116
117//////////////////////////////////////////////////////////////////////////
118//
119// Return original solid not scaled
120//
122{
123 return fPtrSolid;
124}
125
126//////////////////////////////////////////////////////////////////////////
127//
128// Get bounding box
129//
131 G4ThreeVector& pMax) const
132{
133 G4ThreeVector bmin,bmax;
134 G4ThreeVector scale = fScale->GetScale();
135
136 fPtrSolid->BoundingLimits(bmin,bmax);
137 pMin.set(bmin.x()*scale.x(),bmin.y()*scale.y(),bmin.z()*scale.z());
138 pMax.set(bmax.x()*scale.x(),bmax.y()*scale.y(),bmax.z()*scale.z());
139
140 // Check correctness of the bounding box
141 //
142 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
143 {
144 std::ostringstream message;
145 message << "Bad bounding box (min >= max) for solid: "
146 << GetName() << " !"
147 << "\npMin = " << pMin
148 << "\npMax = " << pMax;
149 G4Exception("G4ScaledSolid::BoundingLimits()", "GeomMgt0001",
150 JustWarning, message);
151 DumpInfo();
152 }
153}
154
155//////////////////////////////////////////////////////////////////////////
156//
157// Calculate extent under transform and specified limit
158//
159G4bool
161 const G4VoxelLimits& pVoxelLimit,
162 const G4AffineTransform& pTransform,
163 G4double& pMin,
164 G4double& pMax ) const
165{
166 // Find bounding box of unscaled solid
167 G4ThreeVector bmin,bmax;
168 fPtrSolid->BoundingLimits(bmin,bmax);
169
170 // Set combined transformation
171 G4Transform3D transform3D =
172 G4Transform3D(pTransform.NetRotation().inverse(),
173 pTransform.NetTranslation())*GetScaleTransform();
174
175 // Find extent
176 G4BoundingEnvelope bbox(bmin,bmax);
177 return bbox.CalculateExtent(pAxis,pVoxelLimit,transform3D,pMin,pMax);
178}
179
180/////////////////////////////////////////////////////
181//
182// Inside
183//
185{
186 return fPtrSolid->Inside(fScale->Transform(p));
187}
188
189//////////////////////////////////////////////////////////////
190//
191// SurfaceNormal
192//
195{
196 // Transform point to unscaled shape frame
197 G4ThreeVector newPoint;
198 fScale->Transform(p, newPoint);
199
200 // Compute normal in unscaled frame
201 G4ThreeVector newNormal = fPtrSolid->SurfaceNormal(newPoint);
202 G4ThreeVector normal;
203
204 // Convert normal to scaled frame
205 fScale->InverseTransformNormal(newNormal, normal);
206 return normal/normal.mag();
207}
208
209/////////////////////////////////////////////////////////////
210//
211// The same algorithm as in DistanceToIn(p)
212//
215 const G4ThreeVector& v ) const
216{
217 // Transform point and direction to unscaled shape frame
218 G4ThreeVector newPoint;
219 fScale->Transform(p, newPoint);
220
221 // Direction is un-normalized after scale transformation
222 G4ThreeVector newDirection;
223 fScale->Transform(v, newDirection);
224 newDirection = newDirection/newDirection.mag();
225
226 // Compute distance in unscaled system
227 G4double dist = fPtrSolid->DistanceToIn(newPoint,newDirection);
228
229 // Return converted distance to global
230 return fScale->InverseTransformDistance(dist, newDirection);
231}
232
233////////////////////////////////////////////////////////
234//
235// Approximate nearest distance from the point p to the solid from outside
236//
239{
240 // Transform point to unscaled shape frame
241 G4ThreeVector newPoint;
242 fScale->Transform(p, newPoint);
243
244 // Compute unscaled safety, then scale it
245 G4double dist = fPtrSolid->DistanceToIn(newPoint);
246 return fScale->InverseTransformDistance(dist);
247}
248
249//////////////////////////////////////////////////////////
250//
251// The same algorithm as DistanceToOut(p)
252//
255 const G4ThreeVector& v,
256 const G4bool calcNorm,
257 G4bool *validNorm,
258 G4ThreeVector *n ) const
259{
260 // Transform point and direction to unscaled shape frame
261 G4ThreeVector newPoint;
262 fScale->Transform(p, newPoint);
263
264 // Direction is un-normalized after scale transformation
265 G4ThreeVector newDirection;
266 fScale->Transform(v, newDirection);
267 newDirection = newDirection/newDirection.mag();
268
269 // Compute distance in unscaled system
270 G4ThreeVector solNorm;
271 G4double dist = fPtrSolid->DistanceToOut(newPoint,newDirection,
272 calcNorm,validNorm,&solNorm);
273 if(calcNorm)
274 {
275 G4ThreeVector normal;
276 fScale->TransformNormal(solNorm, normal);
277 *n = normal.unit();
278 }
279
280 // Return distance converted to global
281 return fScale->InverseTransformDistance(dist, newDirection);
282}
283
284//////////////////////////////////////////////////////////////
285//
286// Approximate nearest distance from the point p to the solid from inside
287//
290{
291 // Transform point to unscaled shape frame
292 G4ThreeVector newPoint;
293 fScale->Transform(p, newPoint);
294
295 // Compute unscaled safety, then scale it
296 G4double dist = fPtrSolid->DistanceToOut(newPoint);
297 return fScale->InverseTransformDistance(dist);
298}
299
300//////////////////////////////////////////////////////////////
301//
302// ComputeDimensions
303//
304void
306 const G4int,
307 const G4VPhysicalVolume* )
308{
309 DumpInfo();
310 G4Exception("G4ScaledSolid::ComputeDimensions()",
311 "GeomSolids0001", FatalException,
312 "Method not applicable in this context!");
313}
314
315//////////////////////////////////////////////////////////////////////////
316//
317// Returns a point (G4ThreeVector) randomly and uniformly selected
318// on the solid surface
319//
321{
322 return fScale->InverseTransform(fPtrSolid->GetPointOnSurface());
323}
324//////////////////////////////////////////////////////////////////////////
325//
326// Return the number of constituents used for construction of the solid
327//
329{
330 return fPtrSolid->GetNumOfConstituents();
331}
332
333//////////////////////////////////////////////////////////////////////////
334//
335// Return true if the solid has only planar faces
336//
338{
339 return fPtrSolid->IsFaceted();
340}
341
342//////////////////////////////////////////////////////////////////////////
343//
344// Return object type name
345//
347{
348 return {"G4ScaledSolid"};
349}
350
351//////////////////////////////////////////////////////////////////////////
352//
353// Make a clone of the object
354//
356{
357 return new G4ScaledSolid(*this);
358}
359
360//////////////////////////////////////////////////////////////////////////
361//
362// Returning the scaling transformation
363//
365{
366 return { fScale->GetScale().x(),
367 fScale->GetScale().y(),
368 fScale->GetScale().z() };
369}
370
371//////////////////////////////////////////////////////////////////////////
372//
373// Setting the scaling transformation
374//
376{
377 delete fScale;
378 fScale = new G4ScaleTransform(scale);
379 fRebuildPolyhedron = true;
380}
381
382//////////////////////////////////////////////////////////////////////////
383//
384// Get volume of the scaled solid
385//
387{
388 if(fCubicVolume < 0.)
389 {
390 G4RecursiveAutoLock l(&scaledMutex);
391 fCubicVolume = fPtrSolid->GetCubicVolume() *
392 fScale->GetScale().x() *
393 fScale->GetScale().y() *
394 fScale->GetScale().z();
395 l.unlock();
396 }
397 return fCubicVolume;
398}
399
400//////////////////////////////////////////////////////////////////////////
401//
402// Get estimated surface area of the scaled solid
403//
405{
406 if(fSurfaceArea < 0.)
407 {
408 G4RecursiveAutoLock l(&scaledMutex);
409 fSurfaceArea = G4VSolid::GetSurfaceArea();
410 l.unlock();
411 }
412 return fSurfaceArea;
413}
414
415//////////////////////////////////////////////////////////////////////////
416//
417// Stream object contents to an output stream
418//
419std::ostream& G4ScaledSolid::StreamInfo(std::ostream& os) const
420{
421 os << "-----------------------------------------------------------\n"
422 << " *** Dump for Scaled solid - " << GetName() << " ***\n"
423 << " ===================================================\n"
424 << " Solid type: " << GetEntityType() << "\n"
425 << " Parameters of constituent solid: \n"
426 << "===========================================================\n";
427 fPtrSolid->StreamInfo(os);
428 os << "===========================================================\n"
429 << " Scaling: \n"
430 << " Scale transformation : \n"
431 << " " << fScale->GetScale().x() << ", "
432 << fScale->GetScale().y() << ", "
433 << fScale->GetScale().z() << "\n"
434 << "===========================================================\n";
435
436 return os;
437}
438
439//////////////////////////////////////////////////////////////////////////
440//
441// DescribeYourselfTo
442//
443void
445{
446 scene.AddSolid (*this);
447}
448
449//////////////////////////////////////////////////////////////////////////
450//
451// CreatePolyhedron
452//
455{
456 G4Polyhedron* polyhedron = fPtrSolid->CreatePolyhedron();
457 if (polyhedron != nullptr)
458 {
459 polyhedron->Transform(GetScaleTransform());
460 }
461 else
462 {
463 DumpInfo();
464 G4Exception("G4ScaledSolid::CreatePolyhedron()",
465 "GeomSolids2003", JustWarning,
466 "No G4Polyhedron for scaled solid");
467 }
468 return polyhedron;
469}
470
471//////////////////////////////////////////////////////////////////////////
472//
473// GetPolyhedron
474//
476{
477 if (fpPolyhedron == nullptr ||
478 fRebuildPolyhedron ||
479 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
480 fpPolyhedron->GetNumberOfRotationSteps())
481 {
482 fpPolyhedron = CreatePolyhedron();
483 fRebuildPolyhedron = false;
484 }
485 return fpPolyhedron;
486}
G4TemplateAutoLock< G4RecursiveMutex > G4RecursiveAutoLock
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
#define G4MUTEX_INITIALIZER
std::recursive_mutex G4RecursiveMutex
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Transform3D G4Transform3D
HepGeom::Scale3D G4Scale3D
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4String G4GeometryType
Definition G4VSolid.hh:70
double z() const
Hep3Vector unit() const
double x() const
double y() const
double mag() const
void set(double x, double y, double z)
HepRotation inverse() const
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4ThreeVector NetTranslation() const
G4RotationMatrix NetRotation() const
G4BoundingEnvelope is a helper class to facilitate calculation of the extent of a solid within the li...
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4ScaleTransform is a class for geometric scaling transformations. It supports efficient arbitrary tr...
G4Polyhedron * CreatePolyhedron() const override
G4VSolid * GetUnscaledSolid() const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
G4ThreeVector GetPointOnSurface() const override
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
G4ScaledSolid & operator=(const G4ScaledSolid &rhs)
G4bool IsFaceted() const override
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
G4Scale3D GetScaleTransform() const
G4double GetCubicVolume() override
G4int GetNumOfConstituents() const override
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
void SetScaleTransform(const G4Scale3D &scale)
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const override
G4ScaledSolid(const G4String &pName, G4VSolid *pSolid, const G4Scale3D &pScale)
EInside Inside(const G4ThreeVector &p) const override
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
G4Polyhedron * GetPolyhedron() const override
G4VSolid * Clone() const override
G4double GetSurfaceArea() override
std::ostream & StreamInfo(std::ostream &os) const override
~G4ScaledSolid() override
G4GeometryType GetEntityType() const override
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
virtual void AddSolid(const G4Box &)=0
G4VPVParameterisation ia an abstract base class for Parameterisation, able to compute the transformat...
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4String GetName() const
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59
void DumpInfo() const
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:108
virtual G4double GetSurfaceArea()
Definition G4VSolid.cc:275
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
HepPolyhedron & Transform(const G4Transform3D &t)
EAxis
Definition geomdefs.hh:54
EInside
Definition geomdefs.hh:67