Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ReflectedSolid.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 G4ReflectedSolid class
27//
28// Author: Vladimir Grichine, 23.07.01 ([email protected])
29// --------------------------------------------------------------------
30
31#include "G4ReflectedSolid.hh"
32
33#include <sstream>
34
35#include "G4Point3D.hh"
36#include "G4Vector3D.hh"
37
38#include "G4AffineTransform.hh"
39#include "G4Transform3D.hh"
40#include "G4VoxelLimits.hh"
41
43
44#include "G4VGraphicsScene.hh"
45#include "G4Polyhedron.hh"
46
47/////////////////////////////////////////////////////////////////
48//
49// Constructor using HepTransform3D, in fact HepReflect3D
50
52 G4VSolid* pSolid ,
53 const G4Transform3D& transform )
54 : G4VSolid(pName)
55{
56 fPtrSolid = pSolid;
57 fDirectTransform3D = new G4Transform3D(transform);
58}
59
60///////////////////////////////////////////////////////////////////
61//
62
68
69///////////////////////////////////////////////////////////////////
70//
71
73 : G4VSolid(rhs), fPtrSolid(rhs.fPtrSolid)
74{
76}
77
78///////////////////////////////////////////////////////////////////
79//
80
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}
101
102///////////////////////////////////////////////////////////////////
103//
104
106{
107 return {"G4ReflectedSolid"};
108}
109
111{
112 return this;
113}
114
119
124
125/////////////////////////////////////////////////////////////////////////////
126//
127
132
134{
135 G4Transform3D aTransform = *fDirectTransform3D;
136 return aTransform;
137}
138
140{
141 fDirectTransform3D = &transform;
142 fRebuildPolyhedron = true;
143}
144
145//////////////////////////////////////////////////////////////////////////
146//
147// Get bounding box
148
150 G4ThreeVector& pMax) const
151{
152 fPtrSolid->BoundingLimits(pMin,pMax);
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);
178 G4AffineTransform transform(transform3D.getRotation().inverse(),
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}
206
207//////////////////////////////////////////////////////////////////////////
208//
209// Calculate extent under transform and specified limit
210
211G4bool
213 const G4VoxelLimits& pVoxelLimits,
214 const G4AffineTransform& pTransform,
215 G4double& pMin,
216 G4double& pMax ) const
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);
236 G4AffineTransform transform(transform3D.getRotation().inverse(),
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}
252
253//////////////////////////////////////////////////////////////
254//
255//
256
258{
259 G4ThreeVector newPoint = (*fDirectTransform3D)*G4Point3D(p);
260 return fPtrSolid->Inside(newPoint);
261}
262
263//////////////////////////////////////////////////////////////
264//
265//
266
269{
270 G4ThreeVector newPoint = (*fDirectTransform3D)*G4Point3D(p);
271 G4Vector3D normal = fPtrSolid->SurfaceNormal(newPoint);
272 return (*fDirectTransform3D)*normal;
273}
274
275/////////////////////////////////////////////////////////////
276//
277// The same algorithm as in DistanceToIn(p)
278
281 const G4ThreeVector& v ) const
282{
283 G4ThreeVector newPoint = (*fDirectTransform3D)*G4Point3D(p);
284 G4ThreeVector newDirection = (*fDirectTransform3D)*G4Vector3D(v);
285 return fPtrSolid->DistanceToIn(newPoint,newDirection);
286}
287
288////////////////////////////////////////////////////////
289//
290// Approximate nearest distance from the point p to the intersection of
291// two solids
292
295{
296 G4ThreeVector newPoint = (*fDirectTransform3D)*G4Point3D(p);
297 return fPtrSolid->DistanceToIn(newPoint);
298}
299
300//////////////////////////////////////////////////////////
301//
302// The same algorithm as DistanceToOut(p)
303
306 const G4ThreeVector& v,
307 const G4bool calcNorm,
308 G4bool* validNorm,
309 G4ThreeVector* n ) const
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}
324
325//////////////////////////////////////////////////////////////
326//
327// Inverted algorithm of DistanceToIn(p)
328
331{
332 G4ThreeVector newPoint = (*fDirectTransform3D)*G4Point3D(p);
333 return fPtrSolid->DistanceToOut(newPoint);
334}
335
336//////////////////////////////////////////////////////////////
337//
338//
339
340void
342 const G4int,
343 const G4VPhysicalVolume* )
344{
345 DumpInfo();
346 G4Exception("G4ReflectedSolid::ComputeDimensions()",
347 "GeomMgt0001", FatalException,
348 "Method not applicable in this context!");
349}
350
351//////////////////////////////////////////////////////////////
352//
353// Return volume
354
359
360//////////////////////////////////////////////////////////////
361//
362// Return surface area
363
368
369//////////////////////////////////////////////////////////////
370//
371// Return a point (G4ThreeVector) randomly and uniformly selected
372// on the solid surface
373
379
380//////////////////////////////////////////////////////////////////////////
381//
382// Make a clone of this object
383
385{
386 return new G4ReflectedSolid(*this);
387}
388
389
390//////////////////////////////////////////////////////////////////////////
391//
392// Stream object contents to an output stream
393
394std::ostream& G4ReflectedSolid::StreamInfo(std::ostream& os) const
395{
396 os << "-----------------------------------------------------------\n"
397 << " *** Dump for Reflected solid - " << GetName() << " ***\n"
398 << " ===================================================\n"
399 << " Solid type: " << GetEntityType() << "\n"
400 << " Parameters of constituent solid: \n"
401 << "===========================================================\n";
403 os << "===========================================================\n"
404 << " Transformations: \n"
405 << " Direct transformation - translation : \n"
406 << " " << fDirectTransform3D->getTranslation() << "\n"
407 << " - rotation : \n"
408 << " ";
410 os << "\n"
411 << "===========================================================\n";
412
413 return os;
414}
415
416/////////////////////////////////////////////////
417//
418//
419
420void
422{
423 scene.AddSolid (*this);
424}
425
426////////////////////////////////////////////////////
427//
428//
429
432{
433 G4Polyhedron* polyhedron = fPtrSolid->CreatePolyhedron();
434 if (polyhedron != nullptr)
435 {
436 polyhedron->Transform(*fDirectTransform3D);
437 return polyhedron;
438 }
439 else
440 {
441 std::ostringstream message;
442 message << "Solid - " << GetName()
443 << " - original solid has no" << G4endl
444 << "corresponding polyhedron. Returning NULL!";
445 G4Exception("G4ReflectedSolid::CreatePolyhedron()",
446 "GeomMgt1001", JustWarning, message);
447 return nullptr;
448 }
449}
450
451/////////////////////////////////////////////////////////
452//
453//
454
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
HepGeom::Transform3D G4Transform3D
HepGeom::ReflectZ3D G4ReflectZ3D
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
HepGeom::Vector3D< G4double > G4Vector3D
Definition G4Vector3D.hh:34
#define G4endl
Definition G4ios.hh:67
double z() const
double x() const
double y() const
void set(double x, double y, double z)
std::ostream & print(std::ostream &os) const
Definition RotationIO.cc:17
HepRotation inverse() const
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
G4VSolid * Clone() const override
G4GeometryType GetEntityType() const override
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
G4Transform3D GetDirectTransform3D() const
G4Polyhedron * fpPolyhedron
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const override
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
std::ostream & StreamInfo(std::ostream &os) const override
EInside Inside(const G4ThreeVector &p) const override
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
G4Polyhedron * CreatePolyhedron() const override
G4Polyhedron * GetPolyhedron() const override
void SetDirectTransform3D(G4Transform3D &)
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
G4Transform3D GetTransform3D() const
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
G4VSolid * GetConstituentMovedSolid() const
~G4ReflectedSolid() override
G4double GetCubicVolume() override
G4ReflectedSolid(const G4String &pName, G4VSolid *pSolid, const G4Transform3D &transform)
virtual const G4ReflectedSolid * GetReflectedSolidPtr() const
G4ReflectedSolid & operator=(const G4ReflectedSolid &rhs)
G4ThreeVector GetPointOnSurface() const override
G4Transform3D * fDirectTransform3D
G4double GetSurfaceArea() override
virtual void AddSolid(const G4Box &)=0
G4String GetName() const
virtual std::ostream & StreamInfo(std::ostream &os) const =0
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 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
G4double GetMinZExtent() const
void AddLimit(const EAxis pAxis, const G4double pMin, const G4double pMax)
G4double GetMaxYExtent() const
G4double GetMaxZExtent() const
G4double GetMinYExtent() const
G4double GetMinXExtent() const
G4double GetMaxXExtent() const
double dy() const
double zz() const
double dz() const
CLHEP::HepRotation getRotation() const
double dx() const
CLHEP::Hep3Vector getTranslation() const
Transform3D inverse() const
double xx() const
double yy() const
static G4int GetNumberOfRotationSteps()
HepPolyhedron & Transform(const G4Transform3D &t)
EAxis
Definition geomdefs.hh:54
@ kYAxis
Definition geomdefs.hh:56
@ kXAxis
Definition geomdefs.hh:55
@ kZAxis
Definition geomdefs.hh:57
EInside
Definition geomdefs.hh:67