58 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
60 #include <base/Global.h>
61 #include <base/Vector3D.h>
65 template <
class UnplacedVolume_t>
66 class G4UAdapter :
public G4VSolid,
protected UnplacedVolume_t
70 using U3Vector = vecgeom::Vector3D<G4double>;
72 using UnplacedVolume_t::operator
delete;
73 using UnplacedVolume_t::operator
new;
81 template <
typename...
T>
83 :
G4VSolid(name), UnplacedVolume_t(params...)
86 virtual ~G4UAdapter();
122 const G4bool calcNorm =
false,
177 virtual std::ostream&
StreamInfo(std::ostream& os)
const override;
194 G4UAdapter(__void__&);
199 G4UAdapter(
const G4UAdapter& rhs);
200 G4UAdapter&
operator=(
const G4UAdapter& rhs);
207 vecgeom::Precision stepMax =
kInfinity)
const override
209 return UnplacedVolume_t::DistanceToOut(position, direction, stepMax);
213 Inside(U3Vector
const& aPoint)
const override
215 return UnplacedVolume_t::Inside(aPoint);
220 const vecgeom::Precision step_max =
kInfinity)
const override
222 return UnplacedVolume_t::DistanceToIn(position, direction, step_max);
225 G4bool Normal(U3Vector
const& aPoint, U3Vector& aNormal)
const override
227 return UnplacedVolume_t::Normal(aPoint, aNormal);
230 void Extent(U3Vector& aMin, U3Vector& aMax)
const override
232 return UnplacedVolume_t::Extent(aMin, aMax);
235 U3Vector SamplePointOnSurface()
const override
237 return UnplacedVolume_t::SamplePointOnSurface();
242 mutable G4bool fRebuildPolyhedron =
false;
247 using UnplacedVolume_t::DistanceToOut;
248 using UnplacedVolume_t::DistanceToIn;
253 template <
class UnplacedVolume_t>
254 G4UAdapter<UnplacedVolume_t>::G4UAdapter(__void__&
a)
255 :
G4VSolid(a), UnplacedVolume_t(*this),
260 template <
class UnplacedVolume_t>
261 G4UAdapter<UnplacedVolume_t>::~G4UAdapter()
263 delete fPolyhedron; fPolyhedron =
nullptr;
266 template <
class UnplacedVolume_t>
270 return (
this == &rhs) ?
true :
false;
273 template <
class UnplacedVolume_t>
274 G4UAdapter<UnplacedVolume_t>::
275 G4UAdapter(
const G4UAdapter& rhs)
276 :
G4VSolid(rhs), UnplacedVolume_t(rhs)
281 template <
class UnplacedVolume_t>
282 G4UAdapter<UnplacedVolume_t>& G4UAdapter<UnplacedVolume_t>::
283 operator=(
const G4UAdapter& rhs)
295 UnplacedVolume_t::operator=(rhs);
299 fRebuildPolyhedron =
false;
300 delete fPolyhedron; fPolyhedron =
nullptr;
306 template <
class UnplacedVolume_t>
307 EInside G4UAdapter<UnplacedVolume_t>::
310 U3Vector
pt(p.
x(), p.
y(), p.
z());
311 vecgeom::EnumInside in_temp;
314 in_temp = UnplacedVolume_t::Inside(
pt);
316 if (in_temp == vecgeom::EnumInside::eInside) in =
kInside;
317 else if (in_temp == vecgeom::EnumInside::eSurface) in =
kSurface;
322 template <
class UnplacedVolume_t>
326 U3Vector
p(pt.
x(), pt.
y(), pt.
z());
328 UnplacedVolume_t::Normal(p, n);
332 template <
class UnplacedVolume_t>
333 G4double G4UAdapter<UnplacedVolume_t>::
336 U3Vector
p(pt.
x(), pt.
y(), pt.
z());
337 U3Vector
v(d.
x(), d.
y(), d.
z());
342 if (dist < kHalfTolerance)
return 0.0;
346 template <
class UnplacedVolume_t>
347 G4double G4UAdapter<UnplacedVolume_t>::
350 U3Vector
p(pt.
x(), pt.
y(), pt.
z());
351 G4double dist = UnplacedVolume_t::SafetyToIn(p);
355 if (dist < kHalfTolerance)
return 0.0;
359 template <
class UnplacedVolume_t>
360 G4double G4UAdapter<UnplacedVolume_t>::
365 U3Vector
p(pt.
x(), pt.
y(), pt.
z());
366 U3Vector
v(d.
x(), d.
y(), d.
z());
371 if (UnplacedVolume_t::IsConvex())
373 U3Vector
n, hitpoint = p + dist *
v;
374 UnplacedVolume_t::Normal(hitpoint, n);
376 norm->
set(n.x(), n.y(), n.z());
386 if (dist < kHalfTolerance)
return 0.0;
390 template <
class UnplacedVolume_t>
391 G4double G4UAdapter<UnplacedVolume_t>::
394 U3Vector
p(pt.
x(), pt.
y(), pt.
z());
395 G4double dist = UnplacedVolume_t::SafetyToOut(p);
399 if (dist < kHalfTolerance)
return 0.0;
403 template <
class UnplacedVolume_t>
404 G4double G4UAdapter<UnplacedVolume_t>::GetCubicVolume()
406 return UnplacedVolume_t::Capacity();
409 template <
class UnplacedVolume_t>
410 G4double G4UAdapter<UnplacedVolume_t>::GetSurfaceArea()
412 return UnplacedVolume_t::SurfaceArea();
415 template <
class UnplacedVolume_t>
416 G4ThreeVector G4UAdapter<UnplacedVolume_t>::GetPointOnSurface()
const
418 U3Vector p = UnplacedVolume_t::SamplePointOnSurface();
429 template <
class UnplacedVolume_t>
430 void G4UAdapter<UnplacedVolume_t>::
435 message <<
"Illegal call to G4UAdapter::ComputeDimensions()" <<
G4endl
436 <<
"Method not overloaded by derived class !";
437 G4Exception(
"G4UAdapter::ComputeDimensions()",
"GeomSolids0003",
441 template <
class UnplacedVolume_t>
442 void G4UAdapter<UnplacedVolume_t>::
448 template <
class UnplacedVolume_t>
450 GetEntityType()
const
454 return "G4" + string;
457 template <
class UnplacedVolume_t>
458 std::ostream& G4UAdapter<UnplacedVolume_t>::
459 StreamInfo(std::ostream& os)
const
465 template <
class UnplacedVolume_t>
466 G4VSolid* G4UAdapter<UnplacedVolume_t>::Clone()
const
469 message <<
"Clone() method not implemented for type: "
470 << GetEntityType() <<
"!" <<
G4endl
471 <<
"Returning NULL pointer!";
476 template <
class UnplacedVolume_t>
477 G4bool G4UAdapter<UnplacedVolume_t>::CalculateExtent(
const EAxis pAxis,
483 UnplacedVolume_t::Extent(vmin,vmax);
489 if (bmin.x() >= bmax.x() || bmin.y() >= bmax.y() || bmin.z() >= bmax.z())
492 message <<
"Bad bounding box (min >= max) for solid: "
493 << GetName() <<
" - " << GetEntityType() <<
" !"
494 <<
"\nmin = " << bmin
495 <<
"\nmax = " << bmax;
496 G4Exception(
"G4UAdapter::CalculateExtent()",
"GeomMgt0001",
502 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
505 template <
class UnplacedVolume_t>
506 G4Polyhedron* G4UAdapter<UnplacedVolume_t>::CreatePolyhedron()
const
511 message <<
"Visualization not supported for USolid shape "
512 << GetEntityType() <<
"... Sorry!" <<
G4endl;
513 G4Exception(
"G4UAdapter::CreatePolyhedron()",
"GeomSolids0003",
518 template <
class UnplacedVolume_t>
519 G4Polyhedron* G4UAdapter<UnplacedVolume_t>::GetPolyhedron()
const
522 fRebuildPolyhedron ||
523 fPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
524 fPolyhedron->GetNumberOfRotationSteps())
528 fPolyhedron = CreatePolyhedron();
529 fRebuildPolyhedron =
false;
535 template <
class UnplacedVolume_t>
536 G4VisExtent G4UAdapter<UnplacedVolume_t>::GetExtent()
const
539 UnplacedVolume_t::Extent(vmin,vmax);
545 #endif // G4GEOM_USE_USOLIDS
547 #endif // G4UADAPTER_HH