File indexing completed on 2025-01-18 09:59:15
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 G4GeometryType GetEntityType() const override;
0169
0170
0171
0172 virtual G4VSolid* Clone() const override;
0173
0174
0175
0176
0177 virtual std::ostream& StreamInfo(std::ostream& os) const override;
0178
0179
0180 virtual void DescribeYourselfTo(G4VGraphicsScene& scene) const override;
0181
0182
0183
0184 virtual G4VisExtent GetExtent() const override;
0185
0186 virtual G4Polyhedron* CreatePolyhedron() const override;
0187
0188 virtual G4Polyhedron* GetPolyhedron() const override;
0189
0190
0191
0192 public:
0193
0194 G4UAdapter(__void__&);
0195
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(__void__& a)
0255 : G4VSolid(a), UnplacedVolume_t(*this),
0256 kHalfTolerance(0.5*kCarTolerance)
0257 {
0258 }
0259
0260 template <class UnplacedVolume_t>
0261 G4UAdapter<UnplacedVolume_t>::~G4UAdapter()
0262 {
0263 delete fPolyhedron; fPolyhedron = nullptr;
0264 }
0265
0266 template <class UnplacedVolume_t>
0267 G4bool G4UAdapter<UnplacedVolume_t>::
0268 operator==(const G4UAdapter& rhs) const
0269 {
0270 return (this == &rhs) ? true : false;
0271 }
0272
0273 template <class UnplacedVolume_t>
0274 G4UAdapter<UnplacedVolume_t>::
0275 G4UAdapter(const G4UAdapter& rhs)
0276 : G4VSolid(rhs), UnplacedVolume_t(rhs)
0277 {
0278 kHalfTolerance = 0.5*kCarTolerance;
0279 }
0280
0281 template <class UnplacedVolume_t>
0282 G4UAdapter<UnplacedVolume_t>& G4UAdapter<UnplacedVolume_t>::
0283 operator=(const G4UAdapter& rhs)
0284 {
0285
0286
0287 if (this == &rhs)
0288 {
0289 return *this;
0290 }
0291
0292
0293
0294 G4VSolid::operator=(rhs);
0295 UnplacedVolume_t::operator=(rhs);
0296
0297
0298
0299 fRebuildPolyhedron = false;
0300 delete fPolyhedron; fPolyhedron = nullptr;
0301 kHalfTolerance = 0.5*kCarTolerance;
0302
0303 return *this;
0304 }
0305
0306 template <class UnplacedVolume_t>
0307 EInside G4UAdapter<UnplacedVolume_t>::
0308 Inside(const G4ThreeVector& p) const
0309 {
0310 U3Vector pt(p.x(), p.y(), p.z());
0311 vecgeom::EnumInside in_temp;
0312 EInside in = kOutside;
0313
0314 in_temp = UnplacedVolume_t::Inside(pt);
0315
0316 if (in_temp == vecgeom::EnumInside::eInside) in = kInside;
0317 else if (in_temp == vecgeom::EnumInside::eSurface) in = kSurface;
0318
0319 return in;
0320 }
0321
0322 template <class UnplacedVolume_t>
0323 G4ThreeVector G4UAdapter<UnplacedVolume_t>::
0324 SurfaceNormal(const G4ThreeVector& pt) const
0325 {
0326 U3Vector p(pt.x(), pt.y(), pt.z());
0327 U3Vector n;
0328 UnplacedVolume_t::Normal(p, n);
0329 return G4ThreeVector(n.x(), n.y(), n.z());
0330 }
0331
0332 template <class UnplacedVolume_t>
0333 G4double G4UAdapter<UnplacedVolume_t>::
0334 DistanceToIn(const G4ThreeVector& pt, const G4ThreeVector& d) const
0335 {
0336 U3Vector p(pt.x(), pt.y(), pt.z());
0337 U3Vector v(d.x(), d.y(), d.z());
0338 G4double dist = UnplacedVolume_t::DistanceToIn(p, v, kInfinity);
0339
0340
0341
0342 if (dist < kHalfTolerance) return 0.0;
0343 return (dist > kInfinity) ? kInfinity : dist;
0344 }
0345
0346 template <class UnplacedVolume_t>
0347 G4double G4UAdapter<UnplacedVolume_t>::
0348 DistanceToIn(const G4ThreeVector& pt) const
0349 {
0350 U3Vector p(pt.x(), pt.y(), pt.z());
0351 G4double dist = UnplacedVolume_t::SafetyToIn(p);
0352
0353
0354
0355 if (dist < kHalfTolerance) return 0.0;
0356 return (dist > kInfinity) ? kInfinity : dist;
0357 }
0358
0359 template <class UnplacedVolume_t>
0360 G4double G4UAdapter<UnplacedVolume_t>::
0361 DistanceToOut(const G4ThreeVector& pt, const G4ThreeVector& d,
0362 const G4bool calcNorm, G4bool* validNorm,
0363 G4ThreeVector* norm) const
0364 {
0365 U3Vector p(pt.x(), pt.y(), pt.z());
0366 U3Vector v(d.x(), d.y(), d.z());
0367
0368 G4double dist = UnplacedVolume_t::DistanceToOut(p, v, kInfinity);
0369 if(calcNorm)
0370 {
0371 *validNorm = UnplacedVolume_t::IsConvex();
0372 U3Vector n, hitpoint = p + dist * v;
0373 UnplacedVolume_t::Normal(hitpoint, n);
0374 norm->set(n.x(), n.y(), n.z());
0375 }
0376
0377
0378
0379 if (dist < kHalfTolerance) return 0.0;
0380 return (dist > kInfinity) ? kInfinity : dist;
0381 }
0382
0383 template <class UnplacedVolume_t>
0384 G4double G4UAdapter<UnplacedVolume_t>::
0385 DistanceToOut(const G4ThreeVector& pt) const
0386 {
0387 U3Vector p(pt.x(), pt.y(), pt.z());
0388 G4double dist = UnplacedVolume_t::SafetyToOut(p);
0389
0390
0391
0392 if (dist < kHalfTolerance) return 0.0;
0393 return (dist > kInfinity) ? kInfinity : dist;
0394 }
0395
0396 template <class UnplacedVolume_t>
0397 G4double G4UAdapter<UnplacedVolume_t>::GetCubicVolume()
0398 {
0399 return UnplacedVolume_t::Capacity();
0400 }
0401
0402 template <class UnplacedVolume_t>
0403 G4double G4UAdapter<UnplacedVolume_t>::GetSurfaceArea()
0404 {
0405 return UnplacedVolume_t::SurfaceArea();
0406 }
0407
0408 template <class UnplacedVolume_t>
0409 G4ThreeVector G4UAdapter<UnplacedVolume_t>::GetPointOnSurface() const
0410 {
0411 U3Vector p = UnplacedVolume_t::SamplePointOnSurface();
0412 return G4ThreeVector(p.x(), p.y(), p.z());
0413 }
0414
0415
0416
0417 namespace
0418 {
0419 G4Mutex pMutex = G4MUTEX_INITIALIZER;
0420 }
0421
0422
0423 template <class UnplacedVolume_t>
0424 std::ostream&
0425 operator<<(std::ostream& os, const G4UAdapter<UnplacedVolume_t>& uAdapted)
0426 {
0427 return uAdapted.StreamInfo(os);
0428 }
0429
0430 template <class UnplacedVolume_t>
0431 void G4UAdapter<UnplacedVolume_t>::
0432 ComputeDimensions(G4VPVParameterisation*, const G4int,
0433 const G4VPhysicalVolume*)
0434 {
0435 std::ostringstream message;
0436 message << "Illegal call to G4UAdapter::ComputeDimensions()" << G4endl
0437 << "Method not overloaded by derived class !";
0438 G4Exception("G4UAdapter::ComputeDimensions()", "GeomSolids0003",
0439 FatalException, message);
0440 }
0441
0442 template <class UnplacedVolume_t>
0443 void G4UAdapter<UnplacedVolume_t>::
0444 DescribeYourselfTo(G4VGraphicsScene& scene) const
0445 {
0446 scene.AddSolid(*this);
0447 }
0448
0449 template <class UnplacedVolume_t>
0450 G4GeometryType G4UAdapter<UnplacedVolume_t>::
0451 GetEntityType() const
0452 {
0453
0454 G4String string = "VSolid";
0455 return "G4" + string;
0456 }
0457
0458 template <class UnplacedVolume_t>
0459 std::ostream& G4UAdapter<UnplacedVolume_t>::
0460 StreamInfo(std::ostream& os) const
0461 {
0462 UnplacedVolume_t::Print(os);
0463 return os;
0464 }
0465
0466 template <class UnplacedVolume_t>
0467 G4VSolid* G4UAdapter<UnplacedVolume_t>::Clone() const
0468 {
0469 std::ostringstream message;
0470 message << "Clone() method not implemented for type: "
0471 << GetEntityType() << "!" << G4endl
0472 << "Returning NULL pointer!";
0473 G4Exception("G4UAdapter::Clone()", "GeomSolids1001", JustWarning, message);
0474 return nullptr;
0475 }
0476
0477 template <class UnplacedVolume_t>
0478 G4bool G4UAdapter<UnplacedVolume_t>::CalculateExtent(const EAxis pAxis,
0479 const G4VoxelLimits& pVoxelLimit,
0480 const G4AffineTransform& pTransform,
0481 G4double& pMin, G4double& pMax) const
0482 {
0483 U3Vector vmin, vmax;
0484 UnplacedVolume_t::Extent(vmin,vmax);
0485 G4ThreeVector bmin(vmin.x(),vmin.y(),vmin.z());
0486 G4ThreeVector bmax(vmax.x(),vmax.y(),vmax.z());
0487
0488
0489
0490 if (bmin.x() >= bmax.x() || bmin.y() >= bmax.y() || bmin.z() >= bmax.z())
0491 {
0492 std::ostringstream message;
0493 message << "Bad bounding box (min >= max) for solid: "
0494 << GetName() << " - " << GetEntityType() << " !"
0495 << "\nmin = " << bmin
0496 << "\nmax = " << bmax;
0497 G4Exception("G4UAdapter::CalculateExtent()", "GeomMgt0001",
0498 JustWarning, message);
0499 StreamInfo(G4cout);
0500 }
0501
0502 G4BoundingEnvelope bbox(bmin,bmax);
0503 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
0504 }
0505
0506 template <class UnplacedVolume_t>
0507 G4Polyhedron* G4UAdapter<UnplacedVolume_t>::CreatePolyhedron() const
0508 {
0509
0510
0511 std::ostringstream message;
0512 message << "Visualization not supported for USolid shape "
0513 << GetEntityType() << "... Sorry!" << G4endl;
0514 G4Exception("G4UAdapter::CreatePolyhedron()", "GeomSolids0003",
0515 FatalException, message);
0516 return nullptr;
0517 }
0518
0519 template <class UnplacedVolume_t>
0520 G4Polyhedron* G4UAdapter<UnplacedVolume_t>::GetPolyhedron() const
0521 {
0522 if (!fPolyhedron ||
0523 fRebuildPolyhedron ||
0524 fPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
0525 fPolyhedron->GetNumberOfRotationSteps())
0526 {
0527 G4AutoLock l(&pMutex);
0528 delete fPolyhedron;
0529 fPolyhedron = CreatePolyhedron();
0530 fRebuildPolyhedron = false;
0531 l.unlock();
0532 }
0533 return fPolyhedron;
0534 }
0535
0536 template <class UnplacedVolume_t>
0537 G4VisExtent G4UAdapter<UnplacedVolume_t>::GetExtent() const
0538 {
0539 U3Vector vmin, vmax;
0540 UnplacedVolume_t::Extent(vmin,vmax);
0541 return G4VisExtent(vmin.x(),vmax.x(),
0542 vmin.y(),vmax.y(),
0543 vmin.z(),vmax.z());
0544 }
0545
0546 #endif
0547
0548 #endif