60 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
62 #include <base/Global.h>
63 #include <base/Vector3D.h>
67 template <
class UnplacedVolume_t>
68 class G4UAdapter :
public G4VSolid,
protected UnplacedVolume_t
72 typedef vecgeom::Vector3D<G4double> U3Vector;
74 using UnplacedVolume_t::operator
delete;
75 using UnplacedVolume_t::operator
new;
80 :
G4VSolid(name), fRebuildPolyhedron(false), fPolyhedron(0)
83 template <
typename... T>
85 :
G4VSolid(name), UnplacedVolume_t(params...),
86 fRebuildPolyhedron(false), fPolyhedron(0)
89 virtual ~G4UAdapter();
125 const G4bool calcNorm =
false,
180 virtual std::ostream&
StreamInfo(std::ostream& os)
const override;
197 G4UAdapter(__void__&);
202 G4UAdapter(
const G4UAdapter& rhs);
203 G4UAdapter&
operator=(
const G4UAdapter& rhs);
210 vecgeom::Precision stepMax =
kInfinity)
const override
212 return UnplacedVolume_t::DistanceToOut(position, direction, stepMax);
216 Inside(U3Vector
const &aPoint)
const override
218 return UnplacedVolume_t::Inside(aPoint);
223 const vecgeom::Precision step_max =
kInfinity)
const override
225 return UnplacedVolume_t::DistanceToIn(position, direction, step_max);
228 G4bool Normal(U3Vector
const &aPoint, U3Vector &aNormal)
const override
230 return UnplacedVolume_t::Normal(aPoint, aNormal);
233 void Extent(U3Vector &aMin, U3Vector &aMax)
const override
235 return UnplacedVolume_t::Extent(aMin, aMax);
238 U3Vector SamplePointOnSurface()
const override
240 return UnplacedVolume_t::SamplePointOnSurface();
245 mutable G4bool fRebuildPolyhedron;
250 using UnplacedVolume_t::DistanceToOut;
251 using UnplacedVolume_t::DistanceToIn;
256 template <
class UnplacedVolume_t>
257 G4UAdapter<UnplacedVolume_t>::G4UAdapter(__void__&
a)
258 :
G4VSolid(a), UnplacedVolume_t(*this),
259 fRebuildPolyhedron(false), fPolyhedron(0),
264 template <
class UnplacedVolume_t>
265 G4UAdapter<UnplacedVolume_t>::~G4UAdapter()
267 delete fPolyhedron; fPolyhedron = 0;
270 template <
class UnplacedVolume_t>
274 return (
this == &rhs) ?
true :
false;
277 template <
class UnplacedVolume_t>
278 G4UAdapter<UnplacedVolume_t>::
279 G4UAdapter(
const G4UAdapter& rhs)
280 :
G4VSolid(rhs), UnplacedVolume_t(rhs),
281 fRebuildPolyhedron(false), fPolyhedron(0)
286 template <
class UnplacedVolume_t>
287 G4UAdapter<UnplacedVolume_t>& G4UAdapter<UnplacedVolume_t>::
288 operator=(
const G4UAdapter& rhs)
300 UnplacedVolume_t::operator=(rhs);
304 fRebuildPolyhedron =
false;
305 delete fPolyhedron; fPolyhedron = 0;
311 template <
class UnplacedVolume_t>
312 EInside G4UAdapter<UnplacedVolume_t>::
315 U3Vector
pt(p.
x(), p.
y(), p.
z());
316 vecgeom::EnumInside in_temp;
319 in_temp = UnplacedVolume_t::Inside(
pt);
321 if (in_temp == vecgeom::EnumInside::eInside) in =
kInside;
322 else if (in_temp == vecgeom::EnumInside::eSurface) in =
kSurface;
327 template <
class UnplacedVolume_t>
331 U3Vector
p(pt.
x(), pt.
y(), pt.
z());
333 UnplacedVolume_t::Normal(p, n);
337 template <
class UnplacedVolume_t>
338 G4double G4UAdapter<UnplacedVolume_t>::
341 U3Vector
p(pt.
x(), pt.
y(), pt.
z());
342 U3Vector v(d.
x(), d.
y(), d.
z());
347 if (dist < kHalfTolerance)
return 0.0;
351 template <
class UnplacedVolume_t>
352 G4double G4UAdapter<UnplacedVolume_t>::
355 U3Vector
p(pt.
x(), pt.
y(), pt.
z());
356 G4double dist = UnplacedVolume_t::SafetyToIn(p);
360 if (dist < kHalfTolerance)
return 0.0;
364 template <
class UnplacedVolume_t>
365 G4double G4UAdapter<UnplacedVolume_t>::
370 U3Vector
p(pt.
x(), pt.
y(), pt.
z());
371 U3Vector v(d.
x(), d.
y(), d.
z());
376 if (UnplacedVolume_t::IsConvex())
378 U3Vector
n, hitpoint = p + dist * v;
379 UnplacedVolume_t::Normal(hitpoint, n);
381 norm->
set(n.x(), n.y(), n.z());
391 if (dist < kHalfTolerance)
return 0.0;
395 template <
class UnplacedVolume_t>
396 G4double G4UAdapter<UnplacedVolume_t>::
399 U3Vector
p(pt.
x(), pt.
y(), pt.
z());
400 G4double dist = UnplacedVolume_t::SafetyToOut(p);
404 if (dist < kHalfTolerance)
return 0.0;
408 template <
class UnplacedVolume_t>
409 G4double G4UAdapter<UnplacedVolume_t>::GetCubicVolume()
411 return UnplacedVolume_t::Capacity();
414 template <
class UnplacedVolume_t>
415 G4double G4UAdapter<UnplacedVolume_t>::GetSurfaceArea()
417 return UnplacedVolume_t::SurfaceArea();
420 template <
class UnplacedVolume_t>
421 G4ThreeVector G4UAdapter<UnplacedVolume_t>::GetPointOnSurface()
const
423 U3Vector p = UnplacedVolume_t::SamplePointOnSurface();;
434 template <
class UnplacedVolume_t>
435 void G4UAdapter<UnplacedVolume_t>::
440 message <<
"Illegal call to G4UAdapter::ComputeDimensions()" <<
G4endl
441 <<
"Method not overloaded by derived class !";
442 G4Exception(
"G4UAdapter::ComputeDimensions()",
"GeomSolids0003",
446 template <
class UnplacedVolume_t>
447 void G4UAdapter<UnplacedVolume_t>::
453 template <
class UnplacedVolume_t>
455 GetEntityType()
const
459 return "G4" + string;
462 template <
class UnplacedVolume_t>
463 std::ostream& G4UAdapter<UnplacedVolume_t>::
464 StreamInfo(std::ostream& os)
const
470 template <
class UnplacedVolume_t>
471 G4VSolid* G4UAdapter<UnplacedVolume_t>::Clone()
const
474 message <<
"Clone() method not implemented for type: "
475 << GetEntityType() <<
"!" <<
G4endl
476 <<
"Returning NULL pointer!";
481 template <
class UnplacedVolume_t>
482 G4bool G4UAdapter<UnplacedVolume_t>::CalculateExtent(
const EAxis pAxis,
488 UnplacedVolume_t::Extent(vmin,vmax);
494 if (bmin.x() >= bmax.x() || bmin.y() >= bmax.y() || bmin.z() >= bmax.z())
497 message <<
"Bad bounding box (min >= max) for solid: "
498 << GetName() <<
" - " << GetEntityType() <<
" !"
499 <<
"\nmin = " << bmin
500 <<
"\nmax = " << bmax;
501 G4Exception(
"G4UAdapter::CalculateExtent()",
"GeomMgt0001",
507 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
510 template <
class UnplacedVolume_t>
511 G4Polyhedron* G4UAdapter<UnplacedVolume_t>::CreatePolyhedron()
const
516 message <<
"Visualization not supported for USolid shape "
517 << GetEntityType() <<
"... Sorry!" <<
G4endl;
518 G4Exception(
"G4UAdapter::CreatePolyhedron()",
"GeomSolids0003",
523 template <
class UnplacedVolume_t>
524 G4Polyhedron* G4UAdapter<UnplacedVolume_t>::GetPolyhedron()
const
527 fRebuildPolyhedron ||
528 fPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
529 fPolyhedron->GetNumberOfRotationSteps())
533 fPolyhedron = CreatePolyhedron();
534 fRebuildPolyhedron =
false;
540 template <
class UnplacedVolume_t>
541 G4VisExtent G4UAdapter<UnplacedVolume_t>::GetExtent()
const
544 UnplacedVolume_t::Extent(vmin,vmax);
550 #endif // G4GEOM_USE_USOLIDS
552 #endif // G4UADAPTER_HH
void set(double x, double y, double z)
CLHEP::Hep3Vector G4ThreeVector
std::vector< ExP01TrackerHit * > a
virtual void AddSolid(const G4Box &)=0
static const G4double kInfinity
virtual G4VisExtent GetExtent() const
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
virtual G4Polyhedron * GetPolyhedron() const
void message(RunManager *runmanager)
virtual G4double GetSurfaceArea()
virtual G4ThreeVector GetPointOnSurface() const
virtual G4VSolid * Clone() const
G4bool operator==(const G4VSolid &s) const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const =0
const G4double kCarTolerance
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
virtual EInside Inside(const G4ThreeVector &p) const =0
G4VSolid & operator=(const G4VSolid &rhs)
#define G4MUTEX_INITIALIZER
virtual std::ostream & StreamInfo(std::ostream &os) const =0
void Print(G4Element &ele)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
virtual G4double GetCubicVolume()
G4GLOB_DLL std::ostream G4cout
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
virtual G4Polyhedron * CreatePolyhedron() const
virtual G4GeometryType GetEntityType() const =0
bool operator==(const HepRotation &r, const HepLorentzRotation <)