File indexing completed on 2025-09-17 08:59:37
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040 #ifndef G4UADAPTER_HH
0041 #define G4UADAPTER_HH
0042
0043 #include "G4ThreeVector.hh"
0044 #include "G4VSolid.hh"
0045
0046
0047
0048 #include "G4AffineTransform.hh"
0049 #include "G4VoxelLimits.hh"
0050 #include "G4VGraphicsScene.hh"
0051 #include "G4Polyhedron.hh"
0052 #include "G4VisExtent.hh"
0053 #include "G4BoundingEnvelope.hh"
0054 #include "G4AutoLock.hh"
0055
0056 #include "G4GeomTypes.hh"
0057
0058 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
0059
0060 #include <VecGeom/base/Global.h>
0061 #include <VecGeom/base/Vector3D.h>
0062
0063 class G4VPVParameterisation;
0064
0065 template <class UnplacedVolume_t>
0066 class G4UAdapter : public G4VSolid, protected UnplacedVolume_t
0067 {
0068 public:
0069
0070 using U3Vector = vecgeom::Vector3D<G4double>;
0071
0072 using UnplacedVolume_t::operator delete;
0073 using UnplacedVolume_t::operator new;
0074
0075
0076
0077 G4UAdapter(const G4String& name)
0078 : G4VSolid(name)
0079 { kHalfTolerance = 0.5*kCarTolerance; }
0080
0081 template <typename... T>
0082 G4UAdapter(const G4String& name, const T &... params)
0083 : G4VSolid(name), UnplacedVolume_t(params...)
0084 { kHalfTolerance = 0.5*kCarTolerance; }
0085
0086 virtual ~G4UAdapter();
0087
0088 G4bool operator==(const G4UAdapter& s) const;
0089
0090
0091 virtual G4bool CalculateExtent(const EAxis pAxis,
0092 const G4VoxelLimits& pVoxelLimit,
0093 const G4AffineTransform& pTransform,
0094 G4double& pMin, G4double& pMax) const override;
0095
0096
0097
0098
0099 virtual EInside Inside(const G4ThreeVector& p) const override;
0100
0101
0102
0103
0104 virtual G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const override;
0105
0106
0107
0108 virtual G4double DistanceToIn(const G4ThreeVector& p,
0109 const G4ThreeVector& v) const override;
0110
0111
0112
0113
0114
0115
0116 virtual G4double DistanceToIn(const G4ThreeVector& p) const override;
0117
0118
0119
0120 virtual G4double DistanceToOut(const G4ThreeVector& p,
0121 const G4ThreeVector& v,
0122 const G4bool calcNorm = false,
0123 G4bool* validNorm = 0,
0124 G4ThreeVector* n = 0) const override;
0125
0126
0127
0128
0129
0130
0131
0132
0133
0134
0135
0136
0137
0138
0139
0140
0141 virtual G4double DistanceToOut(const G4ThreeVector& p) const override;
0142
0143
0144
0145 virtual void ComputeDimensions(G4VPVParameterisation* p,
0146 const G4int n,
0147 const G4VPhysicalVolume* pRep) override;
0148
0149
0150
0151 virtual G4double GetCubicVolume() override;
0152
0153
0154
0155
0156
0157
0158 virtual G4double GetSurfaceArea() override;
0159
0160
0161
0162
0163
0164
0165 virtual G4ThreeVector GetPointOnSurface() const override;
0166
0167
0168 virtual G4int GetNumOfConstituents() const override;
0169
0170
0171
0172 virtual G4bool IsFaceted() const override;
0173
0174
0175 virtual G4GeometryType GetEntityType() const override;
0176
0177
0178
0179 virtual G4VSolid* Clone() const override;
0180
0181
0182
0183
0184 virtual std::ostream& StreamInfo(std::ostream& os) const override;
0185
0186
0187 virtual void DescribeYourselfTo(G4VGraphicsScene& scene) const override;
0188
0189
0190
0191 virtual G4VisExtent GetExtent() const override;
0192
0193 virtual G4Polyhedron* CreatePolyhedron() const override;
0194
0195 virtual G4Polyhedron* GetPolyhedron() const override;
0196
0197
0198
0199 G4UAdapter(const G4UAdapter& rhs);
0200 G4UAdapter& operator=(const G4UAdapter& rhs);
0201
0202
0203 public:
0204
0205 vecgeom::Precision
0206 DistanceToOut(U3Vector const& position, U3Vector const& direction,
0207 vecgeom::Precision stepMax = kInfinity) const override
0208 {
0209 return UnplacedVolume_t::DistanceToOut(position, direction, stepMax);
0210 }
0211
0212 vecgeom::EnumInside
0213 Inside(U3Vector const& aPoint) const override
0214 {
0215 return UnplacedVolume_t::Inside(aPoint);
0216 }
0217
0218 vecgeom::Precision
0219 DistanceToIn(U3Vector const& position, U3Vector const& direction,
0220 const vecgeom::Precision step_max = kInfinity) const override
0221 {
0222 return UnplacedVolume_t::DistanceToIn(position, direction, step_max);
0223 }
0224
0225 G4bool Normal(U3Vector const& aPoint, U3Vector& aNormal) const override
0226 {
0227 return UnplacedVolume_t::Normal(aPoint, aNormal);
0228 }
0229
0230 void Extent(U3Vector& aMin, U3Vector& aMax) const override
0231 {
0232 return UnplacedVolume_t::Extent(aMin, aMax);
0233 }
0234
0235 U3Vector SamplePointOnSurface() const override
0236 {
0237 return UnplacedVolume_t::SamplePointOnSurface();
0238 }
0239
0240 protected:
0241
0242 mutable G4bool fRebuildPolyhedron = false;
0243 mutable G4Polyhedron* fPolyhedron = nullptr;
0244
0245 G4double kHalfTolerance;
0246
0247 using UnplacedVolume_t::DistanceToOut;
0248 using UnplacedVolume_t::DistanceToIn;
0249 };
0250
0251
0252
0253 template <class UnplacedVolume_t>
0254 G4UAdapter<UnplacedVolume_t>::~G4UAdapter()
0255 {
0256 delete fPolyhedron; fPolyhedron = nullptr;
0257 }
0258
0259 template <class UnplacedVolume_t>
0260 G4bool G4UAdapter<UnplacedVolume_t>::
0261 operator==(const G4UAdapter& rhs) const
0262 {
0263 return (this == &rhs) ? true : false;
0264 }
0265
0266 template <class UnplacedVolume_t>
0267 G4UAdapter<UnplacedVolume_t>::
0268 G4UAdapter(const G4UAdapter& rhs)
0269 : G4VSolid(rhs), UnplacedVolume_t(rhs)
0270 {
0271 kHalfTolerance = 0.5*kCarTolerance;
0272 }
0273
0274 template <class UnplacedVolume_t>
0275 G4UAdapter<UnplacedVolume_t>& G4UAdapter<UnplacedVolume_t>::
0276 operator=(const G4UAdapter& rhs)
0277 {
0278
0279
0280 if (this == &rhs)
0281 {
0282 return *this;
0283 }
0284
0285
0286
0287 G4VSolid::operator=(rhs);
0288 UnplacedVolume_t::operator=(rhs);
0289
0290
0291
0292 fRebuildPolyhedron = false;
0293 delete fPolyhedron; fPolyhedron = nullptr;
0294 kHalfTolerance = 0.5*kCarTolerance;
0295
0296 return *this;
0297 }
0298
0299 template <class UnplacedVolume_t>
0300 EInside G4UAdapter<UnplacedVolume_t>::
0301 Inside(const G4ThreeVector& p) const
0302 {
0303 U3Vector pt(p.x(), p.y(), p.z());
0304 vecgeom::EnumInside in_temp;
0305 EInside in = kOutside;
0306
0307 in_temp = UnplacedVolume_t::Inside(pt);
0308
0309 if (in_temp == vecgeom::EnumInside::eInside) in = kInside;
0310 else if (in_temp == vecgeom::EnumInside::eSurface) in = kSurface;
0311
0312 return in;
0313 }
0314
0315 template <class UnplacedVolume_t>
0316 G4ThreeVector G4UAdapter<UnplacedVolume_t>::
0317 SurfaceNormal(const G4ThreeVector& pt) const
0318 {
0319 U3Vector p(pt.x(), pt.y(), pt.z());
0320 U3Vector n;
0321 UnplacedVolume_t::Normal(p, n);
0322 return G4ThreeVector(n.x(), n.y(), n.z());
0323 }
0324
0325 template <class UnplacedVolume_t>
0326 G4double G4UAdapter<UnplacedVolume_t>::
0327 DistanceToIn(const G4ThreeVector& pt, const G4ThreeVector& d) const
0328 {
0329 U3Vector p(pt.x(), pt.y(), pt.z());
0330 U3Vector v(d.x(), d.y(), d.z());
0331 G4double dist = UnplacedVolume_t::DistanceToIn(p, v, kInfinity);
0332
0333
0334
0335 if (dist < kHalfTolerance) return 0.0;
0336 return (dist > kInfinity) ? kInfinity : dist;
0337 }
0338
0339 template <class UnplacedVolume_t>
0340 G4double G4UAdapter<UnplacedVolume_t>::
0341 DistanceToIn(const G4ThreeVector& pt) const
0342 {
0343 U3Vector p(pt.x(), pt.y(), pt.z());
0344 G4double dist = UnplacedVolume_t::SafetyToIn(p);
0345
0346
0347
0348 if (dist < kHalfTolerance) return 0.0;
0349 return (dist > kInfinity) ? kInfinity : dist;
0350 }
0351
0352 template <class UnplacedVolume_t>
0353 G4double G4UAdapter<UnplacedVolume_t>::
0354 DistanceToOut(const G4ThreeVector& pt, const G4ThreeVector& d,
0355 const G4bool calcNorm, G4bool* validNorm,
0356 G4ThreeVector* norm) const
0357 {
0358 U3Vector p(pt.x(), pt.y(), pt.z());
0359 U3Vector v(d.x(), d.y(), d.z());
0360
0361 G4double dist = UnplacedVolume_t::DistanceToOut(p, v, kInfinity);
0362 if(calcNorm)
0363 {
0364 *validNorm = UnplacedVolume_t::IsConvex();
0365 U3Vector n, hitpoint = p + dist * v;
0366 UnplacedVolume_t::Normal(hitpoint, n);
0367 norm->set(n.x(), n.y(), n.z());
0368 }
0369
0370
0371
0372 if (dist < kHalfTolerance) return 0.0;
0373 return (dist > kInfinity) ? kInfinity : dist;
0374 }
0375
0376 template <class UnplacedVolume_t>
0377 G4double G4UAdapter<UnplacedVolume_t>::
0378 DistanceToOut(const G4ThreeVector& pt) const
0379 {
0380 U3Vector p(pt.x(), pt.y(), pt.z());
0381 G4double dist = UnplacedVolume_t::SafetyToOut(p);
0382
0383
0384
0385 if (dist < kHalfTolerance) return 0.0;
0386 return (dist > kInfinity) ? kInfinity : dist;
0387 }
0388
0389 template <class UnplacedVolume_t>
0390 G4double G4UAdapter<UnplacedVolume_t>::GetCubicVolume()
0391 {
0392 return UnplacedVolume_t::Capacity();
0393 }
0394
0395 template <class UnplacedVolume_t>
0396 G4double G4UAdapter<UnplacedVolume_t>::GetSurfaceArea()
0397 {
0398 return UnplacedVolume_t::SurfaceArea();
0399 }
0400
0401 template <class UnplacedVolume_t>
0402 G4ThreeVector G4UAdapter<UnplacedVolume_t>::GetPointOnSurface() const
0403 {
0404 U3Vector p = UnplacedVolume_t::SamplePointOnSurface();
0405 return G4ThreeVector(p.x(), p.y(), p.z());
0406 }
0407
0408 template <class UnplacedVolume_t>
0409 G4int G4UAdapter<UnplacedVolume_t>::GetNumOfConstituents() const
0410 {
0411 return 1;
0412 }
0413
0414 template <class UnplacedVolume_t>
0415 G4bool G4UAdapter<UnplacedVolume_t>::IsFaceted() const
0416 {
0417 return false;
0418 }
0419
0420
0421
0422 namespace
0423 {
0424 G4Mutex pMutex = G4MUTEX_INITIALIZER;
0425 }
0426
0427
0428 template <class UnplacedVolume_t>
0429 std::ostream&
0430 operator<<(std::ostream& os, const G4UAdapter<UnplacedVolume_t>& uAdapted)
0431 {
0432 return uAdapted.StreamInfo(os);
0433 }
0434
0435 template <class UnplacedVolume_t>
0436 void G4UAdapter<UnplacedVolume_t>::
0437 ComputeDimensions(G4VPVParameterisation*, const G4int,
0438 const G4VPhysicalVolume*)
0439 {
0440 std::ostringstream message;
0441 message << "Illegal call to G4UAdapter::ComputeDimensions()" << G4endl
0442 << "Method not overloaded by derived class !";
0443 G4Exception("G4UAdapter::ComputeDimensions()", "GeomSolids0003",
0444 FatalException, message);
0445 }
0446
0447 template <class UnplacedVolume_t>
0448 void G4UAdapter<UnplacedVolume_t>::
0449 DescribeYourselfTo(G4VGraphicsScene& scene) const
0450 {
0451 scene.AddSolid(*this);
0452 }
0453
0454 template <class UnplacedVolume_t>
0455 G4GeometryType G4UAdapter<UnplacedVolume_t>::
0456 GetEntityType() const
0457 {
0458
0459 G4String string = "VSolid";
0460 return "G4" + string;
0461 }
0462
0463 template <class UnplacedVolume_t>
0464 std::ostream& G4UAdapter<UnplacedVolume_t>::
0465 StreamInfo(std::ostream& os) const
0466 {
0467 UnplacedVolume_t::Print(os);
0468 return os;
0469 }
0470
0471 template <class UnplacedVolume_t>
0472 G4VSolid* G4UAdapter<UnplacedVolume_t>::Clone() const
0473 {
0474 std::ostringstream message;
0475 message << "Clone() method not implemented for type: "
0476 << GetEntityType() << "!" << G4endl
0477 << "Returning NULL pointer!";
0478 G4Exception("G4UAdapter::Clone()", "GeomSolids1001", JustWarning, message);
0479 return nullptr;
0480 }
0481
0482 template <class UnplacedVolume_t>
0483 G4bool G4UAdapter<UnplacedVolume_t>::CalculateExtent(const EAxis pAxis,
0484 const G4VoxelLimits& pVoxelLimit,
0485 const G4AffineTransform& pTransform,
0486 G4double& pMin, G4double& pMax) const
0487 {
0488 U3Vector vmin, vmax;
0489 UnplacedVolume_t::Extent(vmin,vmax);
0490 G4ThreeVector bmin(vmin.x(),vmin.y(),vmin.z());
0491 G4ThreeVector bmax(vmax.x(),vmax.y(),vmax.z());
0492
0493
0494
0495 if (bmin.x() >= bmax.x() || bmin.y() >= bmax.y() || bmin.z() >= bmax.z())
0496 {
0497 std::ostringstream message;
0498 message << "Bad bounding box (min >= max) for solid: "
0499 << GetName() << " - " << GetEntityType() << " !"
0500 << "\nmin = " << bmin
0501 << "\nmax = " << bmax;
0502 G4Exception("G4UAdapter::CalculateExtent()", "GeomMgt0001",
0503 JustWarning, message);
0504 StreamInfo(G4cout);
0505 }
0506
0507 G4BoundingEnvelope bbox(bmin,bmax);
0508 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
0509 }
0510
0511 template <class UnplacedVolume_t>
0512 G4Polyhedron* G4UAdapter<UnplacedVolume_t>::CreatePolyhedron() const
0513 {
0514
0515
0516 std::ostringstream message;
0517 message << "Visualization not supported for USolid shape "
0518 << GetEntityType() << "... Sorry!" << G4endl;
0519 G4Exception("G4UAdapter::CreatePolyhedron()", "GeomSolids0003",
0520 FatalException, message);
0521 return nullptr;
0522 }
0523
0524 template <class UnplacedVolume_t>
0525 G4Polyhedron* G4UAdapter<UnplacedVolume_t>::GetPolyhedron() const
0526 {
0527 if (!fPolyhedron ||
0528 fRebuildPolyhedron ||
0529 fPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
0530 fPolyhedron->GetNumberOfRotationSteps())
0531 {
0532 G4AutoLock l(&pMutex);
0533 delete fPolyhedron;
0534 fPolyhedron = CreatePolyhedron();
0535 fRebuildPolyhedron = false;
0536 l.unlock();
0537 }
0538 return fPolyhedron;
0539 }
0540
0541 template <class UnplacedVolume_t>
0542 G4VisExtent G4UAdapter<UnplacedVolume_t>::GetExtent() const
0543 {
0544 U3Vector vmin, vmax;
0545 UnplacedVolume_t::Extent(vmin,vmax);
0546 return G4VisExtent(vmin.x(),vmax.x(),
0547 vmin.y(),vmax.y(),
0548 vmin.z(),vmax.z());
0549 }
0550
0551 #endif
0552
0553 #endif