G4OCCT 0.1.0
Geant4 interface to Open CASCADE Technology (OCCT) geometry definitions
Loading...
Searching...
No Matches
G4OCCTSolid.cc
Go to the documentation of this file.
1// SPDX-License-Identifier: LGPL-2.1-or-later
2// Copyright (C) 2026 G4OCCT Contributors
3
6
8
9#include <BRepAdaptor_Curve.hxx>
10#include <BRepAdaptor_Surface.hxx>
11#include <BRepBndLib.hxx>
12#include <BRepClass3d_SolidClassifier.hxx>
13#include <BRepExtrema_DistShapeShape.hxx>
14#include <BRepExtrema_SupportType.hxx>
15#include <BRepExtrema_TriangleSet.hxx>
16#include <BRepGProp.hxx>
17#include <BRepLib.hxx>
18#include <BRepLProp_SLProps.hxx>
19#include <BRepMesh_IncrementalMesh.hxx>
20#include <BRep_Builder.hxx>
21#include <BRep_Tool.hxx>
22#include <BRepTools.hxx>
23#include <BRepTools_WireExplorer.hxx>
24#include <BVH_Distance.hxx>
25#include <BVH_Tools.hxx>
26#include <Bnd_Box.hxx>
27#include <ElSLib.hxx>
28#include <GeomAbs_CurveType.hxx>
29#include <GeomAbs_SurfaceType.hxx>
30#include <GProp_GProps.hxx>
31#include <GeomAPI_ProjectPointOnSurf.hxx>
32#include <Geom_Surface.hxx>
33#include <IFSelect_ReturnStatus.hxx>
34#include <IntCurvesFace_Intersector.hxx>
35#include <Poly_Triangulation.hxx>
36#include <Precision.hxx>
37#include <STEPControl_Reader.hxx>
38#include <TopAbs_Orientation.hxx>
39#include <TopAbs_ShapeEnum.hxx>
40#include <TopAbs_State.hxx>
41#include <TopExp_Explorer.hxx>
42#include <TopLoc_Location.hxx>
43#include <TopoDS.hxx>
44#include <TopoDS_Edge.hxx>
45#include <TopoDS_Face.hxx>
46#include <TopoDS_Vertex.hxx>
47#include <TopoDS_Wire.hxx>
48#include <gp_Dir.hxx>
49#include <gp_Lin.hxx>
50#include <gp_Pln.hxx>
51#include <gp_Pnt.hxx>
52#include <gp_Pnt2d.hxx>
53#include <gp_Vec.hxx>
54
55#include <G4BoundingEnvelope.hh>
56#include <G4AffineTransform.hh>
57#include <G4Exception.hh>
58#include <G4GeometryTolerance.hh>
59#include <G4Polyhedron.hh>
60#include <G4TessellatedSolid.hh>
61#include <G4TriangularFacet.hh>
62#include <G4VGraphicsScene.hh>
63#include <G4VisExtent.hh>
64#include <Randomize.hh>
65
66#include <algorithm>
67#include <cmath>
68#include <optional>
69#include <stdexcept>
70#include <string>
71#include <utility>
72
73namespace {
74
79constexpr Standard_Real kRelativeDeflection = 0.01;
80
92class PointToMeshDistance
93 : public BVH_Distance<Standard_Real, 3, BVH_Vec3d, BRepExtrema_TriangleSet> {
94public:
97 Standard_Boolean RejectNode(const BVH_Vec3d& theCornerMin, const BVH_Vec3d& theCornerMax,
98 Standard_Real& theMetric) const override {
99 theMetric =
100 BVH_Tools<Standard_Real, 3>::PointBoxSquareDistance(myObject, theCornerMin, theCornerMax);
101 return RejectMetric(theMetric);
102 }
103
105 Standard_Boolean Accept(const Standard_Integer theIndex, const Standard_Real&) override {
106 BVH_Vec3d v0, v1, v2;
107 myBVHSet->GetVertices(theIndex, v0, v1, v2);
108 const Standard_Real sq =
109 BVH_Tools<Standard_Real, 3>::PointTriangleSquareDistance(myObject, v0, v1, v2);
110 if (sq < myDistance) {
111 myDistance = sq;
112 myBestIndex = theIndex;
113 return Standard_True;
114 }
115 return Standard_False;
116 }
117
121 Standard_Integer BestIndex() const { return myBestIndex; }
122
123private:
124 Standard_Integer myBestIndex{-1};
125};
126
137class TriangleRayCast
138 : public BVH_Traverse<Standard_Real, 3, BRepExtrema_TriangleSet, Standard_Real> {
139public:
142 void SetRay(const BVH_Vec3d& theOrigin, const BVH_Vec3d& theDir, Standard_Real theTolerance) {
143 myOrigin = theOrigin;
144 myDir = theDir;
145 myTolerance = theTolerance;
146
147 // Reset traversal output state so each cast starts from a clean state.
148 myCrossings = 0;
149 myOnSurface = Standard_False;
150 myDegenerate = Standard_False;
151 }
152
154 Standard_Boolean RejectNode(const BVH_Vec3d& theCornerMin, const BVH_Vec3d& theCornerMax,
155 Standard_Real& theMetric) const override {
156 Standard_Real tmin = 0.0;
157 Standard_Real tmax = Precision::Infinite();
158 for (int k = 0; k < 3; ++k) {
159 const Standard_Real dk = (k == 0) ? myDir.x() : (k == 1) ? myDir.y() : myDir.z();
160 const Standard_Real ok = (k == 0) ? myOrigin.x() : (k == 1) ? myOrigin.y() : myOrigin.z();
161 const Standard_Real ck_min = (k == 0) ? theCornerMin.x()
162 : (k == 1) ? theCornerMin.y()
163 : theCornerMin.z();
164 const Standard_Real ck_max = (k == 0) ? theCornerMax.x()
165 : (k == 1) ? theCornerMax.y()
166 : theCornerMax.z();
167 if (std::abs(dk) < Precision::Confusion()) {
168 if (ok < ck_min - myTolerance || ok > ck_max + myTolerance) {
169 return Standard_True;
170 }
171 } else {
172 Standard_Real t1 = (ck_min - ok) / dk;
173 Standard_Real t2 = (ck_max - ok) / dk;
174 if (t1 > t2) {
175 std::swap(t1, t2);
176 }
177 if (t1 > tmin) {
178 tmin = t1;
179 }
180 if (t2 < tmax) {
181 tmax = t2;
182 }
183 if (tmin > tmax + myTolerance) {
184 return Standard_True;
185 }
186 }
187 }
188 theMetric = tmin;
189 return Standard_False;
190 }
191
193 Standard_Boolean Accept(const Standard_Integer theIndex, const Standard_Real&) override {
194 BVH_Vec3d v0, v1, v2;
195 myBVHSet->GetVertices(theIndex, v0, v1, v2);
196 const BVH_Vec3d edge1 = v1 - v0;
197 const BVH_Vec3d edge2 = v2 - v0;
198 const BVH_Vec3d h = BVH_Vec3d::Cross(myDir, edge2);
199 const Standard_Real a = edge1.Dot(h);
200 if (std::abs(a) < 1e-12) {
201 return Standard_False; // ray parallel to triangle plane
202 }
203 const Standard_Real f = 1.0 / a;
204 const BVH_Vec3d s = myOrigin - v0;
205 const Standard_Real u = f * s.Dot(h);
206 if (u < 0.0 || u > 1.0) {
207 return Standard_False;
208 }
209 const BVH_Vec3d q = BVH_Vec3d::Cross(s, edge1);
210 const Standard_Real v = f * myDir.Dot(q);
211 if (v < 0.0 || u + v > 1.0) {
212 return Standard_False;
213 }
214 const Standard_Real t = f * edge2.Dot(q);
215 if (std::abs(t) <= myTolerance) {
216 myOnSurface = Standard_True; // ray origin lies on this triangle
217 return Standard_True;
218 }
219 if (t < -myTolerance) {
220 return Standard_False; // triangle is behind the ray origin
221 }
222 // Forward crossing at t > myTolerance.
223 ++myCrossings;
224 // Proximity to an edge or vertex signals a potentially degenerate parity count.
225 constexpr Standard_Real kEdgeTol = 1e-6;
226 if (u < kEdgeTol || v < kEdgeTol || (1.0 - u - v) < kEdgeTol) {
227 myDegenerate = Standard_True;
228 }
229 return Standard_True;
230 }
231
233 Standard_Boolean RejectMetric(const Standard_Real&) const override { return Standard_False; }
234
236 Standard_Boolean Stop() const override { return Standard_False; }
237
238 int Crossings() const { return myCrossings; }
239 Standard_Boolean OnSurface() const { return myOnSurface; }
240 Standard_Boolean Degenerate() const { return myDegenerate; }
241
242private:
243 BVH_Vec3d myOrigin;
244 BVH_Vec3d myDir;
245 Standard_Real myTolerance{1e-7};
246 int myCrossings{0};
247 Standard_Boolean myOnSurface{Standard_False};
248 Standard_Boolean myDegenerate{Standard_False};
249};
250
252gp_Pnt ToPoint(const G4ThreeVector& point) { return gp_Pnt(point.x(), point.y(), point.z()); }
253
255G4double IntersectionTolerance() {
256 return 0.5 * G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
257}
258
260TopoDS_Vertex MakeVertex(const G4ThreeVector& point) {
261 BRep_Builder builder;
262 TopoDS_Vertex vertex;
263 builder.MakeVertex(vertex, ToPoint(point), IntersectionTolerance());
264 return vertex;
265}
266
268EInside ToG4Inside(const TopAbs_State state) {
269 switch (state) {
270 case TopAbs_IN:
271 return kInside;
272 case TopAbs_ON:
273 return kSurface;
274 case TopAbs_OUT:
275 default:
276 return kOutside;
277 }
278}
279
281G4ThreeVector FallbackNormal() { return G4ThreeVector(0.0, 0.0, 1.0); }
282
288bool PointInPolygon2d(Standard_Real u, Standard_Real v, const std::vector<gp_Pnt2d>& poly) {
289 const std::size_t n = poly.size();
290 int crossings = 0;
291 for (std::size_t i = 0; i < n; ++i) {
292 const gp_Pnt2d& a = poly[i];
293 const gp_Pnt2d& b = poly[(i + 1) % n];
294 const Standard_Real av = a.Y();
295 const Standard_Real bv = b.Y();
296 if ((av <= v && bv > v) || (bv <= v && av > v)) {
297 const Standard_Real uCross = a.X() + (v - av) * (b.X() - a.X()) / (bv - av);
298 if (u < uCross) {
299 ++crossings;
300 }
301 }
302 }
303 return (crossings % 2) == 1;
304}
305
310bool PointOnPolygonBoundary2d(Standard_Real u, Standard_Real v, const std::vector<gp_Pnt2d>& poly,
311 Standard_Real tol) {
312 const Standard_Real tol2 = tol * tol;
313 const std::size_t n = poly.size();
314 for (std::size_t i = 0; i < n; ++i) {
315 const gp_Pnt2d& a = poly[i];
316 const gp_Pnt2d& b = poly[(i + 1) % n];
317 const Standard_Real dx = b.X() - a.X();
318 const Standard_Real dy = b.Y() - a.Y();
319 const Standard_Real len2 = dx * dx + dy * dy;
320 Standard_Real px = 0.0;
321 Standard_Real py = 0.0;
322 if (len2 < 1.0e-20) {
323 px = a.X();
324 py = a.Y();
325 } else {
326 const Standard_Real t_seg =
327 std::max(0.0, std::min(1.0, ((u - a.X()) * dx + (v - a.Y()) * dy) / len2));
328 px = a.X() + t_seg * dx;
329 py = a.Y() + t_seg * dy;
330 }
331 const Standard_Real dist2 = (u - px) * (u - px) + (v - py) * (v - py);
332 if (dist2 <= tol2) {
333 return true;
334 }
335 }
336 return false;
337}
338
357std::optional<Standard_Real>
358RayPlaneFaceHit(const gp_Lin& ray, const gp_Pln& plane, const std::vector<gp_Pnt2d>& uvPoly,
359 Standard_Real tMin, Standard_Real tMax, Standard_Real tolerance,
360 Standard_Real* u_out = nullptr, Standard_Real* v_out = nullptr) {
361 const gp_Dir& lineDir = ray.Direction();
362 const gp_Dir& plnNormal = plane.Axis().Direction();
363 const Standard_Real denom =
364 plnNormal.X() * lineDir.X() + plnNormal.Y() * lineDir.Y() + plnNormal.Z() * lineDir.Z();
365 // Ray parallel to the plane: no intersection (or entire ray lies on plane,
366 // which we treat as no intersection for navigation purposes).
367 if (std::abs(denom) < 1.0e-10) {
368 return std::nullopt;
369 }
370 const gp_Pnt& orig = ray.Location();
371 const gp_Pnt& planePt = plane.Location();
372 const Standard_Real numer = plnNormal.X() * (planePt.X() - orig.X()) +
373 plnNormal.Y() * (planePt.Y() - orig.Y()) +
374 plnNormal.Z() * (planePt.Z() - orig.Z());
375 const Standard_Real t = numer / denom;
376 if (t < tMin || t > tMax) {
377 return std::nullopt;
378 }
379 const gp_Pnt hitPt(orig.X() + t * lineDir.X(), orig.Y() + t * lineDir.Y(),
380 orig.Z() + t * lineDir.Z());
381 Standard_Real u = 0.0;
382 Standard_Real v = 0.0;
383 ElSLib::PlaneParameters(plane.Position(), hitPt, u, v);
384 // Accept the hit when the UV point is strictly inside the polygon, or when it
385 // lies within tolerance of any polygon edge (the TopAbs_ON-equivalent case).
386 // Pure boundary-miss points that are strictly outside are rejected.
387 if (!PointInPolygon2d(u, v, uvPoly) && !PointOnPolygonBoundary2d(u, v, uvPoly, tolerance)) {
388 return std::nullopt;
389 }
390 if (u_out != nullptr) {
391 *u_out = u;
392 }
393 if (v_out != nullptr) {
394 *v_out = v;
395 }
396 return t;
397}
398
402std::optional<G4ThreeVector> TryGetOutwardNormal(const BRepAdaptor_Surface& surface,
403 const TopoDS_Face& face, const Standard_Real u,
404 const Standard_Real v) {
405 Standard_Real adjustedU = u;
406 Standard_Real adjustedV = v;
407 const Standard_Real tolerance = IntersectionTolerance();
408 // Nudge U away from parametric boundaries to avoid seam singularities.
409 // For periodic surfaces this prevents descending into OCCT's 2D seam
410 // classifier; for non-periodic surfaces with degenerate boundary edges
411 // (e.g. the poles of a sphere) this ensures IsNormalDefined() returns true.
412 {
413 const Standard_Real uFirst = surface.FirstUParameter();
414 const Standard_Real uLast = surface.LastUParameter();
415 const Standard_Real uEpsilon = std::min(
416 std::max(surface.UResolution(tolerance), Precision::PConfusion()), 0.5 * (uLast - uFirst));
417 if (std::abs(adjustedU - uFirst) <= uEpsilon) {
418 adjustedU = std::min(uFirst + uEpsilon, uLast);
419 } else if (std::abs(adjustedU - uLast) <= uEpsilon) {
420 adjustedU = std::max(uLast - uEpsilon, uFirst);
421 }
422 }
423 // Nudge V away from parametric boundaries. Non-periodic V boundaries
424 // include degenerate pole edges on spheres (V = ±π/2) and cone apices.
425 {
426 const Standard_Real vFirst = surface.FirstVParameter();
427 const Standard_Real vLast = surface.LastVParameter();
428 const Standard_Real vEpsilon = std::min(
429 std::max(surface.VResolution(tolerance), Precision::PConfusion()), 0.5 * (vLast - vFirst));
430 if (std::abs(adjustedV - vFirst) <= vEpsilon) {
431 adjustedV = std::min(vFirst + vEpsilon, vLast);
432 } else if (std::abs(adjustedV - vLast) <= vEpsilon) {
433 adjustedV = std::max(vLast - vEpsilon, vFirst);
434 }
435 }
436
437 BRepLProp_SLProps props(surface, adjustedU, adjustedV, 1, tolerance);
438 // If IsNormalDefined() fails the initial nudge may still be too small for
439 // OCCT's internal null-derivative threshold (e.g. a sphere pole where
440 // |dP/dU| = R·sin(δ) must exceed IntersectionTolerance()). Retry with
441 // exponentially larger V nudges, walking toward the V midpoint.
442 if (!props.IsNormalDefined()) {
443 const Standard_Real vFirst = surface.FirstVParameter();
444 const Standard_Real vLast = surface.LastVParameter();
445 const Standard_Real vMid = 0.5 * (vFirst + vLast);
446 const bool nearVFirst = (adjustedV < vMid);
447 const Standard_Real vRes = std::max(surface.VResolution(tolerance), Precision::PConfusion());
448 Standard_Real finalRetryV = adjustedV;
449 for (int attempt = 0; attempt < 8 && !props.IsNormalDefined(); ++attempt) {
450 const Standard_Real scale = std::pow(10.0, static_cast<Standard_Real>(attempt));
451 const Standard_Real nudge = scale * vRes;
452 finalRetryV =
453 nearVFirst ? std::min(adjustedV + nudge, vMid) : std::max(adjustedV - nudge, vMid);
454 props = BRepLProp_SLProps(surface, adjustedU, finalRetryV, 1, tolerance);
455 }
456 if (!props.IsNormalDefined()) {
457 return std::nullopt;
458 }
459 // Guard: if the retry had to drift V by more than kMaxRetryVDriftFraction of
460 // the full V range, the normal is sampled from a different surface region
461 // (e.g. the equatorial band instead of the pole on a NURBS ellipsoid), which
462 // would produce a spurious mismatch against the exact analytic normal, so we
463 // treat the normal as unavailable at this point and return std::nullopt for
464 // the caller to handle (e.g. by leaving the normal invalid or using a fallback).
465 constexpr Standard_Real kMaxRetryVDriftFraction = 0.10;
466 if (std::fabs(finalRetryV - adjustedV) > kMaxRetryVDriftFraction * (vLast - vFirst)) {
467 return std::nullopt;
468 }
469 }
470
471 gp_Dir faceNormal = props.Normal();
472 if (face.Orientation() == TopAbs_REVERSED) {
473 faceNormal.Reverse();
474 }
475
476 return G4ThreeVector(faceNormal.X(), faceNormal.Y(), faceNormal.Z());
477}
478
479} // namespace
480
481G4OCCTSolid::G4OCCTSolid(const G4String& name, const TopoDS_Shape& shape)
482 : G4VSolid(name), fShape(shape) {
483 if (fShape.IsNull()) {
484 throw std::invalid_argument("G4OCCTSolid: shape must not be null");
485 }
486 ComputeBounds();
487}
488
490 // Manually clear cached OCCT objects to avoid G4Exception during static destruction.
491 // G4Cache::Destroy() calls G4Exception with FatalException if cache size validation
492 // fails, which can happen when G4SolidStore cleanup occurs during exit() and the
493 // thread-local storage has already been torn down.
494 // Resetting the optional values here releases the OCCT objects before G4Cache
495 // attempts to destroy the thread-local storage.
496 fClassifierCache.Get().classifier.reset();
497 fIntersectorCache.Get().faceIntersectors.clear();
498 fSphereCache.Get().spheres.clear();
499}
500
501// ── G4OCCTSolid static factory ────────────────────────────────────────────────
502
503G4OCCTSolid* G4OCCTSolid::FromSTEP(const G4String& name, const std::string& path) {
504 STEPControl_Reader reader;
505 if (reader.ReadFile(path.c_str()) != IFSelect_RetDone) {
506 throw std::runtime_error("G4OCCTSolid::FromSTEP: failed to read STEP file: " + path);
507 }
508 if (reader.TransferRoots() <= 0) {
509 throw std::runtime_error("G4OCCTSolid::FromSTEP: no roots transferred from: " + path);
510 }
511 TopoDS_Shape shape = reader.OneShape();
512 if (shape.IsNull()) {
513 throw std::runtime_error("G4OCCTSolid::FromSTEP: null shape loaded from: " + path);
514 }
515 return new G4OCCTSolid(name, shape);
516}
517
518// ── G4OCCTSolid private helpers ───────────────────────────────────────────────
519
520void G4OCCTSolid::ComputeBounds() {
521 // Pre-build PCurves for edges on planar faces so BRep_Tool::CurveOnPlane()
522 // returns the stored result on every Inside() query instead of recomputing
523 // via GeomProjLib::ProjectOnPlane each time (~3.3% of total instructions).
524 for (TopExp_Explorer faceEx(fShape, TopAbs_FACE); faceEx.More(); faceEx.Next()) {
525 const TopoDS_Face& face = TopoDS::Face(faceEx.Current());
526 if (BRepAdaptor_Surface(face).GetType() != GeomAbs_Plane) {
527 continue;
528 }
529 for (TopExp_Explorer edgeEx(face, TopAbs_EDGE); edgeEx.More(); edgeEx.Next()) {
530 BRepLib::BuildPCurveForEdgeOnPlane(TopoDS::Edge(edgeEx.Current()), face);
531 }
532 }
533
534 Bnd_Box boundingBox;
535 BRepBndLib::AddOptimal(fShape, boundingBox, /*useTriangulation=*/Standard_False);
536 if (boundingBox.IsVoid()) {
537 throw std::invalid_argument("G4OCCTSolid: shape has no computable bounding box (no geometry)");
538 }
539
540 Standard_Real xMin = 0.0;
541 Standard_Real yMin = 0.0;
542 Standard_Real zMin = 0.0;
543 Standard_Real xMax = 0.0;
544 Standard_Real yMax = 0.0;
545 Standard_Real zMax = 0.0;
546 boundingBox.Get(xMin, yMin, zMin, xMax, yMax, zMax);
547
548 fCachedBounds =
549 AxisAlignedBounds{G4ThreeVector(xMin, yMin, zMin), G4ThreeVector(xMax, yMax, zMax)};
550
551 // Build per-face bounding-box cache for the SurfaceNormal prefilter.
552 fFaceBoundsCache.clear();
553 fFaceAdaptorIndex.clear();
554 G4double maxFaceDiag = 0.0;
555 for (TopExp_Explorer ex(fShape, TopAbs_FACE); ex.More(); ex.Next()) {
556 Bnd_Box faceBox;
557 BRepBndLib::AddOptimal(ex.Current(), faceBox, /*useTriangulation=*/Standard_False);
558 const TopoDS_Face& currentFace = TopoDS::Face(ex.Current());
559 const std::size_t idx = fFaceBoundsCache.size();
560 BRepAdaptor_Surface adaptor(currentFace);
561 std::optional<gp_Pln> maybePlane;
562 std::vector<gp_Pnt2d> uvPolygon;
563 std::optional<G4ThreeVector> outwardNormal;
564 if (adaptor.GetType() == GeomAbs_Plane) {
565 maybePlane = adaptor.Plane();
566 const gp_Ax3& pos = maybePlane->Position();
567 const TopoDS_Wire wire = BRepTools::OuterWire(currentFace);
568 if (!wire.IsNull()) {
569 // Collect vertices of the outer wire in order via BRepTools_WireExplorer.
570 // Only accept faces where every edge is a straight line segment: for faces
571 // with curved edges (arcs, splines) the vertex-only polygon is an
572 // under-approximation and would produce false misses.
573 bool allLinear = true;
574 std::vector<gp_Pnt2d> poly;
575 for (BRepTools_WireExplorer we(wire, currentFace); we.More(); we.Next()) {
576 const BRepAdaptor_Curve ec(we.Current());
577 if (ec.GetType() != GeomAbs_Line) {
578 allLinear = false;
579 break;
580 }
581 const gp_Pnt pt = BRep_Tool::Pnt(we.CurrentVertex());
582 Standard_Real u = 0.0;
583 Standard_Real v = 0.0;
584 ElSLib::PlaneParameters(pos, pt, u, v);
585 poly.emplace_back(u, v);
586 }
587 if (allLinear && poly.size() >= 3) {
588 // Verify the face has no inner wires (holes). The outer-wire-only polygon
589 // cannot represent cutouts: a ray hitting inside the outer boundary but
590 // inside a hole would be incorrectly counted as a face crossing.
591 int wireCount = 0;
592 for (TopExp_Explorer wc(currentFace, TopAbs_WIRE); wc.More(); wc.Next()) {
593 ++wireCount;
594 if (wireCount > 1) {
595 break;
596 }
597 }
598 if (wireCount == 1) {
599 uvPolygon = std::move(poly);
600 // Precompute the constant outward face normal. For planar faces the
601 // normal is uniform; gp_Pln::Axis() points along the plane normal, but
602 // the outward direction depends on the face orientation in the solid.
603 gp_Dir faceNormal = maybePlane->Axis().Direction();
604 if (currentFace.Orientation() == TopAbs_REVERSED) {
605 faceNormal.Reverse();
606 }
607 outwardNormal = G4ThreeVector(faceNormal.X(), faceNormal.Y(), faceNormal.Z());
608 }
609 }
610 }
611 }
612 fFaceBoundsCache.push_back({currentFace, faceBox, std::move(adaptor), std::move(maybePlane),
613 std::move(uvPolygon), std::move(outwardNormal)});
614 fFaceAdaptorIndex[currentFace.TShape().get()].push_back(idx);
615 // Track the largest face bounding-box diagonal to bound the tessellation error.
616 if (!faceBox.IsVoid()) {
617 Standard_Real fx0 = 0.0;
618 Standard_Real fy0 = 0.0;
619 Standard_Real fz0 = 0.0;
620 Standard_Real fx1 = 0.0;
621 Standard_Real fy1 = 0.0;
622 Standard_Real fz1 = 0.0;
623 faceBox.Get(fx0, fy0, fz0, fx1, fy1, fz1);
624 const G4double diag = G4ThreeVector(fx1 - fx0, fy1 - fy0, fz1 - fz0).mag();
625 maxFaceDiag = std::max(maxFaceDiag, diag);
626 }
627 }
628
629 // fBVHDeflection bounds the Hausdorff distance between the analytical surface
630 // and the tessellation built with relative deflection kRelativeDeflection.
631 fBVHDeflection = kRelativeDeflection * maxFaceDiag;
632
633 // Determine if all faces are planar (enables the fast plane-distance path in DistanceToOut).
634 fAllFacesPlanar = std::all_of(fFaceBoundsCache.begin(), fFaceBoundsCache.end(),
635 [](const FaceBounds& fb) { return fb.plane.has_value(); });
636
637 // Ensure a triangulation is present (idempotent if mesh already exists).
638 // Limit lifetime to this scope so mesh resources are released before building
639 // the BVH, preserving the original temporary-object behavior.
640 {
641 [[maybe_unused]] const BRepMesh_IncrementalMesh mesher(fShape, kRelativeDeflection,
642 /*isRelative=*/Standard_True);
643 }
644
645 // Build the BVH-accelerated triangle set over all tessellated faces.
646 BRepExtrema_ShapeList faces;
647 for (TopExp_Explorer ex(fShape, TopAbs_FACE); ex.More(); ex.Next()) {
648 faces.Append(ex.Current());
649 }
650 if (faces.IsEmpty()) {
651 fTriangleSet.Nullify();
652 fFaceDeflections.clear();
653 } else {
654 fTriangleSet = new BRepExtrema_TriangleSet(faces);
655
656 // Build per-face deflection table indexed in the same order as fFaceBoundsCache
657 // (and as the faces list above — both use the same TopExp_Explorer traversal).
658 // At query time, BRepExtrema_TriangleSet::GetFaceID() maps the nearest triangle's
659 // post-BVH-reorder index back to its face index, allowing a per-face lookup here.
660 fFaceDeflections.clear();
661 fFaceDeflections.reserve(fFaceBoundsCache.size());
662 for (const FaceBounds& fb : fFaceBoundsCache) {
663 // Start from the global BVH deflection so that per-face values never
664 // reduce the safety margin below fBVHDeflection, even for faces with
665 // a void bounding box or very small extent.
666 G4double deflection = fBVHDeflection;
667 if (!fb.box.IsVoid()) {
668 Standard_Real fx0 = 0.0, fy0 = 0.0, fz0 = 0.0;
669 Standard_Real fx1 = 0.0, fy1 = 0.0, fz1 = 0.0;
670 fb.box.Get(fx0, fy0, fz0, fx1, fy1, fz1);
671 const G4double faceDiag = G4ThreeVector(fx1 - fx0, fy1 - fy0, fz1 - fz0).mag();
672 deflection = kRelativeDeflection * faceDiag;
673 }
674 fFaceDeflections.push_back(deflection);
675 }
676 }
677
678 ComputeInitialSpheres();
679}
680
681void G4OCCTSolid::ComputeInitialSpheres() {
682 fInitialSpheres.clear();
683
684 const G4double tol = IntersectionTolerance();
685 const G4ThreeVector& bmin = fCachedBounds.min;
686 const G4ThreeVector& bmax = fCachedBounds.max;
687 const G4ThreeVector centre = 0.5 * (bmin + bmax);
688 const G4ThreeVector halfExt = 0.5 * (bmax - bmin);
689
690 // Candidate seed points: AABB centre, 6 axis-offset interior points (each at half
691 // the distance from the AABB centre to the nearest face along that axis), and 8 octant
692 // centres (15 total).
693 std::vector<G4ThreeVector> candidates;
694 candidates.reserve(15);
695 candidates.push_back(centre);
696 for (const G4double s : {-0.5, 0.5}) {
697 candidates.push_back(centre + G4ThreeVector(s * halfExt.x(), 0.0, 0.0));
698 candidates.push_back(centre + G4ThreeVector(0.0, s * halfExt.y(), 0.0));
699 candidates.push_back(centre + G4ThreeVector(0.0, 0.0, s * halfExt.z()));
700 }
701 for (const int sx : {-1, 1}) {
702 for (const int sy : {-1, 1}) {
703 for (const int sz : {-1, 1}) {
704 candidates.push_back(centre + G4ThreeVector(0.75 * sx * halfExt.x(),
705 0.75 * sy * halfExt.y(),
706 0.75 * sz * halfExt.z()));
707 }
708 }
709 }
710
711 // Use a local classifier: construction-time only, no thread-local needed.
712 BRepClass3d_SolidClassifier localClassifier;
713 localClassifier.Load(fShape);
714
715 for (const G4ThreeVector& cand : candidates) {
716 localClassifier.Perform(ToPoint(cand), tol);
717 if (localClassifier.State() != TopAbs_IN) {
718 continue;
719 }
720 G4double d = BVHLowerBoundDistance(cand);
721 if (d >= kInfinity || d <= tol) {
722 // BVH unavailable or too close to surface; try exact distance.
723 const auto match = TryFindClosestFace(fFaceBoundsCache, cand);
724 if (!match.has_value() || match->distance <= tol) {
725 continue;
726 }
727 d = match->distance;
728 }
729 fInitialSpheres.push_back({cand, d});
730 }
731 // Sort descending by radius so the most useful spheres are checked first.
732 std::sort(fInitialSpheres.begin(), fInitialSpheres.end(),
733 [](const InscribedSphere& a, const InscribedSphere& b) { return a.radius > b.radius; });
734}
735
736std::optional<G4OCCTSolid::ClosestFaceMatch>
737G4OCCTSolid::TryFindClosestFace(const std::vector<FaceBounds>& faceBoundsCache,
738 const G4ThreeVector& point, G4double maxDistance) {
739 if (faceBoundsCache.empty()) {
740 return std::nullopt;
741 }
742
743 const gp_Pnt queryPoint = ToPoint(point);
744 const TopoDS_Vertex queryVertex = MakeVertex(point);
745
746 // Build the point box once — it is constant for the entire loop.
747 Bnd_Box queryBox;
748 queryBox.Add(queryPoint);
749
750 std::optional<ClosestFaceMatch> bestMatch;
751 for (std::size_t i = 0; i < faceBoundsCache.size(); ++i) {
752 const FaceBounds& fb = faceBoundsCache[i];
753 // Lower bound: distance from query point to the face's axis-aligned bounding box.
754 // Use maxDistance as the initial pruning threshold (before a bestMatch is found),
755 // tightening to bestMatch->distance once a candidate has been accepted.
756 const G4double threshold = bestMatch.has_value() ? bestMatch->distance : maxDistance;
757 if (threshold < kInfinity && fb.box.Distance(queryBox) > threshold) {
758 continue;
759 }
760
761 // Use BRepExtrema_DistShapeShape for a robust point-to-face distance.
762 // Extrema_ExtPS only finds interior critical points of the distance function
763 // on the (unbounded) surface; for bilinear (degree-1) B-spline faces the
764 // global minimum is often on the face boundary (an edge or vertex), which
765 // Extrema_ExtPS misses. BRepExtrema_DistShapeShape considers the bounded
766 // face including all its edges and vertices, guaranteeing the true minimum.
767 BRepExtrema_DistShapeShape distance(queryVertex, fb.face);
768 if (!distance.IsDone() || distance.NbSolution() == 0) {
769 continue;
770 }
771 const G4double candidateDistance = distance.Value();
772
773 if (bestMatch.has_value() && candidateDistance >= bestMatch->distance) {
774 continue;
775 }
776 ClosestFaceMatch match{.face = fb.face, .distance = candidateDistance, .faceIndex = i};
777 // Extract UV from BRepExtrema when solution is on the face interior.
778 // For edge/vertex solutions ParOnFaceS2 has a different meaning so we
779 // skip caching in those cases and let SurfaceNormal fall back to
780 // GeomAPI_ProjectPointOnSurf.
781 if (distance.NbSolution() > 0 && distance.SupportTypeShape2(1) == BRepExtrema_IsInFace) {
782 Standard_Real u = 0.0;
783 Standard_Real v = 0.0;
784 distance.ParOnFaceS2(1, u, v);
785 match.uv = std::make_pair(u, v);
786 }
787 bestMatch = std::move(match);
788 }
789
790 return bestMatch;
791}
792
793BRepClass3d_SolidClassifier& G4OCCTSolid::GetOrCreateClassifier() const {
794 ClassifierCache& cache = fClassifierCache.Get();
795 const std::uint64_t currentGen = fShapeGeneration.load(std::memory_order_acquire);
796 if (cache.generation != currentGen) {
797 cache.classifier.emplace();
798 cache.classifier->Load(fShape); // O(N_faces) — paid once per thread per shape
799 cache.generation = currentGen;
800 }
801 return *cache.classifier;
802}
803
804G4OCCTSolid::IntersectorCache& G4OCCTSolid::GetOrCreateIntersector() const {
805 IntersectorCache& cache = fIntersectorCache.Get();
806 const std::uint64_t currentGen = fShapeGeneration.load(std::memory_order_acquire);
807 if (cache.generation != currentGen) {
808 const G4double tol = IntersectionTolerance();
809 cache.faceIntersectors.clear();
810 cache.faceIntersectors.reserve(fFaceBoundsCache.size());
811 cache.expandedBoxes.clear();
812 cache.expandedBoxes.reserve(fFaceBoundsCache.size());
813 for (const auto& fb : fFaceBoundsCache) {
814 cache.faceIntersectors.push_back(
815 std::make_unique<IntCurvesFace_Intersector>(fb.face, tol)); // O(1) per face
816 Bnd_Box expanded = fb.box;
817 expanded.Enlarge(tol); // prevent false-positive IsOut for zero-thickness planar face boxes
818 cache.expandedBoxes.push_back(std::move(expanded));
819 }
820 cache.generation = currentGen;
821 }
822 return cache;
823}
824
825G4OCCTSolid::SphereCacheData& G4OCCTSolid::GetOrInitSphereCache() const {
826 SphereCacheData& cache = fSphereCache.Get();
827 const std::uint64_t currentGen = fShapeGeneration.load(std::memory_order_acquire);
828 if (cache.generation != currentGen) {
829 // Re-seed from the construction-time initial spheres (already sorted descending).
830 cache.spheres = fInitialSpheres;
831 cache.generation = currentGen;
832 }
833 return cache;
834}
835
836void G4OCCTSolid::TryInsertSphere(const G4ThreeVector& centre, G4double d) const {
837 if (d <= 0.0) {
838 return;
839 }
840
841 // Enforce a minimum meaningful radius, consistent with ComputeInitialSpheres().
842 const G4double minRadius = IntersectionTolerance();
843 if (d <= minRadius) {
844 return;
845 }
846 SphereCacheData& cache = GetOrInitSphereCache();
847
848 // Quick capacity check: if at capacity and the new radius is no better than
849 // the smallest stored radius, there is nothing to gain.
850 if (cache.spheres.size() >= kMaxInscribedSpheres && d <= cache.spheres.back().radius) {
851 return;
852 }
853
854 // Dominance check: skip if the new sphere is fully contained inside an existing one.
855 // Condition: |centre − cᵢ| + d ≤ rᵢ ↔ |centre − cᵢ|² ≤ (rᵢ − d)² (when rᵢ ≥ d).
856 for (const InscribedSphere& s : cache.spheres) {
857 if (s.radius >= d) {
858 const G4double gap = s.radius - d;
859 if ((centre - s.centre).mag2() <= gap * gap) {
860 return;
861 }
862 }
863 }
864
865 // Insert in sorted position (descending by radius).
866 const InscribedSphere newSphere{centre, d};
867 const auto it = std::lower_bound(
868 cache.spheres.begin(), cache.spheres.end(), newSphere,
869 [](const InscribedSphere& a, const InscribedSphere& b) { return a.radius > b.radius; });
870 cache.spheres.insert(it, newSphere);
871
872 // Evict the smallest sphere if the cache has grown beyond capacity.
873 if (cache.spheres.size() > kMaxInscribedSpheres) {
874 cache.spheres.pop_back();
875 }
876}
877
878// ── G4VSolid pure-virtual implementations ────────────────────────────────────
879
880EInside G4OCCTSolid::Inside(const G4ThreeVector& p) const {
881 const G4double tolerance = IntersectionTolerance();
882 if (p.x() < fCachedBounds.min.x() - tolerance || p.x() > fCachedBounds.max.x() + tolerance ||
883 p.y() < fCachedBounds.min.y() - tolerance || p.y() > fCachedBounds.max.y() + tolerance ||
884 p.z() < fCachedBounds.min.z() - tolerance || p.z() > fCachedBounds.max.z() + tolerance) {
885 return kOutside;
886 }
887
888 // Fast inscribed-sphere check: if p lies strictly inside any cached sphere
889 // (with a tolerance margin), every such sphere is provably interior to the
890 // solid, so we can return kInside immediately without an OCCT classifier call.
891 // The check uses (radius - tolerance) to avoid misclassifying boundary-adjacent
892 // points as kInside — matching the open-ball guarantee from the safety distance.
893 const SphereCacheData& sphereCache = GetOrInitSphereCache();
894 for (const InscribedSphere& s : sphereCache.spheres) {
895 const G4double interiorRadius = s.radius - tolerance;
896 if (interiorRadius > 0.0 && (p - s.centre).mag2() < interiorRadius * interiorRadius) {
897 return kInside;
898 }
899 }
900
901 // Tier-2a: BVH ray-parity test over the pre-built fTriangleSet.
902 // Cast a +Z ray and count crossings via Möller–Trumbore in a single
903 // O(log N_triangles) BVH traversal, replacing the per-face
904 // IntCurvesFace_Intersector calls for curved faces.
905 if (!fTriangleSet.IsNull() && fTriangleSet->Size() > 0) {
906 // Guard: if the point is within the mesh error band (distance-to-surface
907 // not provably greater than tolerance + deflection), the tessellation may
908 // disagree with the exact solid near the boundary. Fall back to the exact
909 // classifier so the classification is always consistent with the analytic shape.
910 if (fBVHDeflection > 0.0) {
911 const G4double bvhLB = BVHLowerBoundDistance(p);
912 if (bvhLB < tolerance) {
913 BRepClass3d_SolidClassifier& classifier = GetOrCreateClassifier();
914 classifier.Perform(ToPoint(p), tolerance);
915 return ToG4Inside(classifier.State());
916 }
917 }
918
919 // Multi-ray majority vote: cast up to three orthogonal rays (+Z, +X, +Y)
920 // through the BVH triangle set and tally inside/outside votes.
921 //
922 // Fast path: if the primary (+Z) ray is non-degenerate and has a non-zero
923 // crossing count, its parity is reliable and we return immediately.
924 //
925 // Slow path (primary degenerate or zero crossings): cast +X and +Y as
926 // tie-breakers. Only non-degenerate rays contribute a vote. A clear
927 // majority (more inside votes than outside, or vice versa) determines the
928 // result. If no majority can be established (all three rays degenerate, or
929 // a genuine 1–1 tie), we fall back to the exact OCCT classifier.
930 //
931 // This eliminates BRepClass3d_SolidClassifier::Perform() calls for all but
932 // the near-surface and all-degenerate cases, collapsing the CSLib_Class2d
933 // and NCollection heap-allocation overhead for complex solids.
934 const BVH_Vec3d bvhOrigin(p.x(), p.y(), p.z());
935 const Standard_Real bvhTol = static_cast<Standard_Real>(tolerance);
936 TriangleRayCast caster;
937 caster.SetBVHSet(fTriangleSet.get());
938
939 // Primary ray: +Z
940 caster.SetRay(bvhOrigin, BVH_Vec3d(0.0, 0.0, 1.0), bvhTol);
941 caster.Select();
942
943 if (caster.OnSurface()) {
944 return kSurface;
945 }
946
947 // Fast path: primary ray is unambiguous.
948 if (!caster.Degenerate() && caster.Crossings() > 0) {
949 return (caster.Crossings() % 2 == 1) ? kInside : kOutside;
950 }
951
952 // Primary is problematic (degenerate or zero crossings). Tally its vote
953 // (excluded if degenerate), then cast two more orthogonal rays.
954 int insideVotes = 0;
955 int outsideVotes = 0;
956 if (!caster.Degenerate()) {
957 if (caster.Crossings() % 2 == 1) {
958 ++insideVotes;
959 } else {
960 ++outsideVotes;
961 }
962 }
963
964 const BVH_Vec3d kExtraRays[2] = {
965 BVH_Vec3d(1.0, 0.0, 0.0),
966 BVH_Vec3d(0.0, 1.0, 0.0),
967 };
968 for (const BVH_Vec3d& dir : kExtraRays) {
969 caster.SetRay(bvhOrigin, dir, bvhTol);
970 caster.Select();
971 if (caster.OnSurface()) {
972 return kSurface;
973 }
974 if (!caster.Degenerate()) {
975 if (caster.Crossings() % 2 == 1) {
976 ++insideVotes;
977 } else {
978 ++outsideVotes;
979 }
980 }
981 }
982
983 if (insideVotes > outsideVotes) {
984 return kInside;
985 }
986 if (outsideVotes > insideVotes) {
987 return kOutside;
988 }
989
990 // No majority (all degenerate, or genuine 1–1 tie): fall back to the exact
991 // OCCT classifier.
992 BRepClass3d_SolidClassifier& classifier = GetOrCreateClassifier();
993 classifier.Perform(ToPoint(p), tolerance);
994 return ToG4Inside(classifier.State());
995 }
996
997 // Tier-2b: fallback when fTriangleSet is unavailable — original per-face
998 // ray-parity test using per-thread IntCurvesFace_Intersector cache.
999 // Cast a ray from p in the canonical +Z direction and count how many faces
1000 // the ray crosses strictly through their interior (TopAbs_IN). Odd = inside.
1001 // TopAbs_ON hits (shared edges/vertices) are excluded from the parity count;
1002 // any such hit signals a degenerate ray and triggers a classifier fallback so
1003 // that altered parity from skipping the crossing does not misclassify the point.
1004 IntersectorCache& cache = GetOrCreateIntersector();
1005 const gp_Lin ray(ToPoint(p), gp_Dir(0.0, 0.0, 1.0));
1006 int crossings = 0;
1007 bool onSurface = false;
1008 bool degenerateRay = false;
1009
1010 for (std::size_t i = 0; i < fFaceBoundsCache.size(); ++i) {
1011 if (cache.expandedBoxes[i].IsOut(ray)) {
1012 continue;
1013 }
1014 const FaceBounds& fb = fFaceBoundsCache[i];
1015 if (!fb.uvPolygon.empty()) {
1016 // Fast path: direct ray–plane intersection without OCCT topology overhead.
1017 // Retrieve the UV hit coordinates so we can detect edge/vertex hits (which
1018 // would be TopAbs_ON in the OCCT intersector) and fall back to the
1019 // classifier just like the degenerate-ray path does for non-planar faces.
1020 Standard_Real u_hit = 0.0;
1021 Standard_Real v_hit = 0.0;
1022 const auto t = RayPlaneFaceHit(ray, *fb.plane, fb.uvPolygon, -tolerance,
1023 Precision::Infinite(), tolerance, &u_hit, &v_hit);
1024 if (t) {
1025 const G4double w = static_cast<G4double>(*t);
1026 if (std::abs(w) <= tolerance) {
1027 onSurface = true;
1028 } else if (w > tolerance) {
1029 // Check whether the hit UV point is on the polygon boundary (within
1030 // tolerance): that corresponds to a TopAbs_ON state in the OCCT
1031 // intersector and requires the same degenerate-ray fallback.
1032 if (PointOnPolygonBoundary2d(u_hit, v_hit, fb.uvPolygon, tolerance)) {
1033 degenerateRay = true;
1034 } else {
1035 ++crossings;
1036 }
1037 }
1038 }
1039 } else {
1040 IntCurvesFace_Intersector& fi = *cache.faceIntersectors[i];
1041 fi.Perform(ray, -tolerance, Precision::Infinite());
1042 if (!fi.IsDone()) {
1043 continue;
1044 }
1045 for (Standard_Integer j = 1; j <= fi.NbPnt(); ++j) {
1046 const G4double w = fi.WParameter(j);
1047 const TopAbs_State state = fi.State(j);
1048 if (std::abs(w) <= tolerance && (state == TopAbs_IN || state == TopAbs_ON)) {
1049 onSurface = true;
1050 } else if (w > tolerance && state == TopAbs_IN) {
1051 ++crossings;
1052 } else if (w > tolerance && state == TopAbs_ON) {
1053 degenerateRay = true;
1054 }
1055 }
1056 }
1057 }
1058
1059 if (onSurface) {
1060 return kSurface;
1061 }
1062 // Safety fallback: use the classifier when the parity count is unreliable —
1063 // either because no TopAbs_IN crossings were found (ray may pass entirely
1064 // through edges/vertices) or because at least one TopAbs_ON edge/vertex hit
1065 // was skipped (altering the effective crossing count).
1066 if (crossings == 0 || degenerateRay) {
1067 BRepClass3d_SolidClassifier& classifier = GetOrCreateClassifier();
1068 classifier.Perform(ToPoint(p), tolerance);
1069 return ToG4Inside(classifier.State());
1070 }
1071 return (crossings % 2 == 1) ? kInside : kOutside;
1072}
1073
1074G4ThreeVector G4OCCTSolid::SurfaceNormal(const G4ThreeVector& p) const {
1075 // Phase 1: all-planar fast path — O(N_faces) dot products, zero OCCT solver calls.
1076 // For planar solids every face's outward normal is constant and precomputed in
1077 // fFaceBoundsCache. The closest face is identified by the minimum plane distance.
1078 if (fAllFacesPlanar) {
1079 const gp_Pnt pt = ToPoint(p);
1080 const FaceBounds* bestFB = nullptr;
1081 G4double bestDist = kInfinity;
1082 for (const FaceBounds& fb : fFaceBoundsCache) {
1083 if (!fb.plane.has_value() || !fb.outwardNormal.has_value()) {
1084 continue;
1085 }
1086 const G4double d = fb.plane->Distance(pt);
1087 if (d < bestDist) {
1088 bestDist = d;
1089 bestFB = &fb;
1090 }
1091 }
1092 if (bestFB) {
1093 return *bestFB->outwardNormal;
1094 }
1095 }
1096
1097 // Phase 2: BVH-seeded face identification for curved (or mixed) solids.
1098 //
1099 // Use face-local extrema + projection instead of a single solid-wide
1100 // BRepExtrema_DistShapeShape(vertex, solid) query. On imported periodic
1101 // surfaces (notably the torus seam hit at (15,0,0) in the geometry-fixture
1102 // tests), the solid-wide path can descend into OCCT's 2D seam classifier and
1103 // either wedge in boundary classification or report an edge-supported
1104 // solution whose direct normal evaluation falls back to +Z. Face-local
1105 // selection still identifies the nearest supporting face, while the explicit
1106 // projection step recovers a usable (u,v) pair that we can nudge off exact
1107 // periodic seams before asking OCCT for derivatives. Revisit this once the
1108 // upstream periodic-classifier fix is available everywhere we build.
1109 //
1110 // BVHLowerBoundDistance gives a lower bound on the true distance to the surface.
1111 // Adding 2×fBVHDeflection (Hausdorff bound between tessellation and analytical
1112 // surface) yields a valid upper bound on the distance to the nearest surface point.
1113 // This seed is passed to TryFindClosestFace to prune faces that are provably
1114 // farther away, reducing the number of BRepExtrema calls for multi-face solids.
1115 // TryFindClosestFace then uses exact BRepExtrema on every candidate, guaranteeing
1116 // the correct nearest face is returned regardless of tessellation approximation
1117 // errors (which could fool a pure-BVH face-identification approach).
1118
1119 // Shared helper: project @p onto @face using a pre-built surface adaptor,
1120 // obtain (u,v), and evaluate the outward normal. Accepting the cached
1121 // BRepAdaptor_Surface avoids reconstructing it on the fly on this hot path.
1122 const auto projectAndGetNormalFallback = [&](const FaceBounds& fb) -> G4ThreeVector {
1123 TopLoc_Location loc;
1124 const Handle(Geom_Surface) surface = BRep_Tool::Surface(fb.face, loc);
1125 if (surface.IsNull()) {
1126 return FallbackNormal();
1127 }
1128 gp_Pnt pLocal = ToPoint(p);
1129 if (!loc.IsIdentity()) {
1130 pLocal.Transform(loc.Transformation().Inverted());
1131 }
1132 GeomAPI_ProjectPointOnSurf projection(pLocal, surface);
1133 if (projection.NbPoints() == 0) {
1134 return FallbackNormal();
1135 }
1136 Standard_Real u = 0.0;
1137 Standard_Real v = 0.0;
1138 projection.LowerDistanceParameters(u, v);
1139 return TryGetOutwardNormal(fb.adaptor, fb.face, u, v).value_or(FallbackNormal());
1140 };
1141
1142 const G4double bvhLB = BVHLowerBoundDistance(p);
1143 const G4double seedDist =
1144 (fBVHDeflection > 0.0 && bvhLB < kInfinity) ? bvhLB + 2.0 * fBVHDeflection : kInfinity;
1145 const auto closestFaceMatch = TryFindClosestFace(fFaceBoundsCache, p, seedDist);
1146 if (!closestFaceMatch.has_value()) {
1147 return FallbackNormal();
1148 }
1149 const FaceBounds& fb = fFaceBoundsCache[closestFaceMatch->faceIndex];
1150
1151 // Fast-path A: planar face with precomputed constant outward normal.
1152 // Applies even in mixed solids (e.g. G4CutTubs) where fAllFacesPlanar is
1153 // false but the closest face happens to be a flat endcap.
1154 if (fb.outwardNormal.has_value()) {
1155 return *fb.outwardNormal;
1156 }
1157
1158 // Fast-path B: reuse UV already computed by BRepExtrema_DistShapeShape
1159 // when the solution is on the face interior, skipping GeomAPI_ProjectPointOnSurf.
1160 if (closestFaceMatch->uv.has_value()) {
1161 const auto [u, v] = *closestFaceMatch->uv;
1162 const auto result = TryGetOutwardNormal(fb.adaptor, fb.face, u, v);
1163 if (result.has_value()) {
1164 return *result;
1165 }
1166 // Fall through to full projection if normal evaluation failed at stored UV
1167 // (e.g. torus seam or degenerate point).
1168 }
1169
1170 // Fallback: full GeomAPI_ProjectPointOnSurf reprojection for edge/vertex
1171 // solutions and cases where TryGetOutwardNormal returned nullopt.
1172 return projectAndGetNormalFallback(fb);
1173}
1174
1175G4double G4OCCTSolid::DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v) const {
1176 const G4double tolerance = IntersectionTolerance();
1177 const gp_Lin ray(ToPoint(p), gp_Dir(v.x(), v.y(), v.z()));
1178
1179 IntersectorCache& cache = GetOrCreateIntersector();
1180
1181 G4double minDistance = kInfinity;
1182 for (std::size_t i = 0; i < fFaceBoundsCache.size(); ++i) {
1183 if (cache.expandedBoxes[i].IsOut(ray)) {
1184 continue;
1185 }
1186 const FaceBounds& fb = fFaceBoundsCache[i];
1187 if (!fb.uvPolygon.empty()) {
1188 // Fast path: bypass IntCurvesFace_Intersector for straight-edge planar faces.
1189 const auto t = RayPlaneFaceHit(ray, *fb.plane, fb.uvPolygon, tolerance, Precision::Infinite(),
1190 tolerance);
1191 if (t && *t < minDistance) {
1192 minDistance = static_cast<G4double>(*t);
1193 }
1194 } else {
1195 IntCurvesFace_Intersector& fi = *cache.faceIntersectors[i];
1196 fi.Perform(ray, tolerance, Precision::Infinite());
1197 if (!fi.IsDone()) {
1198 continue;
1199 }
1200 for (Standard_Integer j = 1; j <= fi.NbPnt(); ++j) {
1201 const G4double w = fi.WParameter(j);
1202 if (w > tolerance && w < minDistance) {
1203 minDistance = w;
1204 }
1205 }
1206 }
1207 }
1208
1209 return minDistance;
1210}
1211
1212G4double G4OCCTSolid::ExactDistanceToIn(const G4ThreeVector& p) const {
1213 BRepClass3d_SolidClassifier& classifier = GetOrCreateClassifier();
1214 classifier.Perform(ToPoint(p), IntersectionTolerance());
1215 if (classifier.State() == TopAbs_IN || classifier.State() == TopAbs_ON) {
1216 return 0.0;
1217 }
1218
1219 const auto match = TryFindClosestFace(fFaceBoundsCache, p);
1220 if (!match.has_value()) {
1221 return kInfinity;
1222 }
1223 return (match->distance <= IntersectionTolerance()) ? 0.0 : match->distance;
1224}
1225
1226G4double G4OCCTSolid::AABBLowerBound(const G4ThreeVector& p) const {
1227 const G4ThreeVector& mn = fCachedBounds.min;
1228 const G4ThreeVector& mx = fCachedBounds.max;
1229 const G4double dx = std::max({0.0, mn.x() - p.x(), p.x() - mx.x()});
1230 const G4double dy = std::max({0.0, mn.y() - p.y(), p.y() - mx.y()});
1231 const G4double dz = std::max({0.0, mn.z() - p.z(), p.z() - mx.z()});
1232 return std::sqrt(dx * dx + dy * dy + dz * dz);
1233}
1234
1235G4double G4OCCTSolid::BVHLowerBoundDistance(const G4ThreeVector& p) const {
1236 if (fTriangleSet.IsNull() || fTriangleSet->Size() == 0) {
1237 return kInfinity;
1238 }
1239 PointToMeshDistance solver;
1240 solver.SetObject(BVH_Vec3d(p.x(), p.y(), p.z()));
1241 solver.SetBVHSet(fTriangleSet.get());
1242 const Standard_Real meshDistSq = solver.ComputeDistance();
1243 if (!solver.IsDone()) {
1244 return kInfinity;
1245 }
1246 // The solver returns a squared distance; take the single sqrt here
1247 // before subtracting the deflection bound (which is in actual-distance space).
1248 const G4double meshDist = std::sqrt(static_cast<G4double>(meshDistSq));
1249
1250 // Use per-face deflection for the face owning the nearest triangle.
1251 // BRepExtrema_TriangleSet::GetFaceID() maps the BVH-reordered triangle index
1252 // back to the face index in fFaceBoundsCache (the same ordering used to build
1253 // the shape list passed to fTriangleSet). This is tighter than fBVHDeflection
1254 // when the nearest triangle belongs to a small face.
1255 G4double deflection = fBVHDeflection;
1256 const Standard_Integer bestIdx = solver.BestIndex();
1257 if (bestIdx >= 0) {
1258 const Standard_Integer faceId = fTriangleSet->GetFaceID(bestIdx);
1259 if (faceId >= 0 && static_cast<std::size_t>(faceId) < fFaceDeflections.size()) {
1260 deflection = fFaceDeflections[static_cast<std::size_t>(faceId)];
1261 }
1262 }
1263
1264 return std::max(0.0, meshDist - deflection);
1265}
1266
1267G4double G4OCCTSolid::PlanarFaceLowerBoundDistance(const G4ThreeVector& p) const {
1268 const gp_Pnt pt = ToPoint(p);
1269 G4double minDist = kInfinity;
1270 for (const FaceBounds& fb : fFaceBoundsCache) {
1271 if (!fb.plane.has_value()) {
1272 continue;
1273 }
1274 const G4double d = static_cast<G4double>(fb.plane->Distance(pt));
1275 if (d < minDist) {
1276 minDist = d;
1277 }
1278 }
1279 return minDist;
1280}
1281
1282G4double G4OCCTSolid::DistanceToIn(const G4ThreeVector& p) const {
1283 // Tier-0: AABB lower bound (O(1)). If the point is clearly outside the shape's
1284 // axis-aligned bounding box, the AABB distance is a guaranteed conservative lower
1285 // bound on the true surface distance and avoids any further OCCT computation.
1286 const G4double aabbDist = AABBLowerBound(p);
1287 if (aabbDist > IntersectionTolerance()) {
1288 return aabbDist;
1289 }
1290
1291 // Tier-1: BVH triangle-mesh lower bound — O(log N_triangles).
1292 // For outside points this is a valid conservative lower bound (as required for
1293 // safety distances). Avoids the expensive per-face BRepExtrema_DistShapeShape
1294 // calls that dominate ExactDistanceToIn for curved surfaces (e.g. ellipsoids).
1295 const G4double bvhDist = BVHLowerBoundDistance(p);
1296 if (bvhDist < kInfinity && bvhDist > IntersectionTolerance()) {
1297 return bvhDist;
1298 }
1299
1300 // Fallback: exact distance (handles points near or inside the surface, including
1301 // interior/surface classification).
1302 return ExactDistanceToIn(p);
1303}
1304
1305G4double G4OCCTSolid::DistanceToOut(const G4ThreeVector& p, const G4ThreeVector& v,
1306 const G4bool calcNorm, G4bool* validNorm,
1307 G4ThreeVector* n) const {
1308 if (validNorm != nullptr) {
1309 *validNorm = false;
1310 }
1311
1312 const G4double tolerance = IntersectionTolerance();
1313 const gp_Lin ray(ToPoint(p), gp_Dir(v.x(), v.y(), v.z()));
1314
1315 IntersectorCache& cache = GetOrCreateIntersector();
1316
1317 G4double minDistance = kInfinity;
1318 std::size_t minFaceIdx = std::numeric_limits<std::size_t>::max();
1319 G4double minU = 0.0;
1320 G4double minV = 0.0;
1321 bool minIsIn = false;
1322 bool minIsFastPath = false;
1323
1324 for (std::size_t i = 0; i < fFaceBoundsCache.size(); ++i) {
1325 if (cache.expandedBoxes[i].IsOut(ray)) {
1326 continue;
1327 }
1328 const FaceBounds& fb = fFaceBoundsCache[i];
1329 if (!fb.uvPolygon.empty()) {
1330 // Fast path: direct ray–plane intersection for straight-edge planar faces.
1331 const auto t = RayPlaneFaceHit(ray, *fb.plane, fb.uvPolygon, tolerance, Precision::Infinite(),
1332 tolerance);
1333 if (t && *t < minDistance) {
1334 minDistance = static_cast<G4double>(*t);
1335 minFaceIdx = i;
1336 minIsIn = true; // mirrors OCCT path: TopAbs_IN and TopAbs_ON both trigger normal
1337 minIsFastPath = true;
1338 }
1339 } else {
1340 IntCurvesFace_Intersector& fi = *cache.faceIntersectors[i];
1341 fi.Perform(ray, tolerance, Precision::Infinite());
1342 if (!fi.IsDone()) {
1343 continue;
1344 }
1345 for (Standard_Integer j = 1; j <= fi.NbPnt(); ++j) {
1346 const G4double w = fi.WParameter(j);
1347 if (w > tolerance && w < minDistance) {
1348 minDistance = w;
1349 minFaceIdx = i;
1350 minU = fi.UParameter(j);
1351 minV = fi.VParameter(j);
1352 minIsIn = (fi.State(j) == TopAbs_IN || fi.State(j) == TopAbs_ON);
1353 minIsFastPath = false;
1354 }
1355 }
1356 }
1357 }
1358
1359 if (minFaceIdx == std::numeric_limits<std::size_t>::max() || minDistance == kInfinity) {
1360 return 0.0;
1361 }
1362
1363 if (calcNorm && validNorm != nullptr && n != nullptr && minIsIn) {
1364 const FaceBounds& fb = fFaceBoundsCache[minFaceIdx];
1365 if (minIsFastPath && fb.outwardNormal.has_value()) {
1366 // Fast path: use the precomputed constant outward normal for planar faces.
1367 *n = *fb.outwardNormal;
1368 *validNorm = true;
1369 } else {
1370 const auto outNorm = TryGetOutwardNormal(fb.adaptor, fb.face, minU, minV);
1371 if (outNorm) {
1372 *n = *outNorm;
1373 *validNorm = true;
1374 }
1375 }
1376 }
1377
1378 return minDistance;
1379}
1380
1381G4double G4OCCTSolid::ExactDistanceToOut(const G4ThreeVector& p) const {
1382 // Use the pre-built per-face AABB cache to prune candidates before calling
1383 // BRepExtrema on individual faces: each query performs O(N_faces) cheap AABB
1384 // checks, but only O(k) faces require expensive extrema evaluations.
1385 const auto match = TryFindClosestFace(fFaceBoundsCache, p);
1386 if (!match.has_value()) {
1387 return 0.0;
1388 }
1389 return (match->distance <= IntersectionTolerance()) ? 0.0 : match->distance;
1390}
1391
1392G4double G4OCCTSolid::DistanceToOut(const G4ThreeVector& p) const {
1393 G4double d;
1394 if (fAllFacesPlanar) {
1395 // All faces are planar: the minimum perpendicular plane distance is exact
1396 // for convex solids and a valid conservative lower bound for all others.
1397 // This avoids the BVH triangle-mesh traversal entirely.
1398 d = PlanarFaceLowerBoundDistance(p);
1399 if (d == kInfinity) {
1400 d = ExactDistanceToOut(p);
1401 }
1402 } else {
1403 // Mixed geometry (some curved faces): the planar lower bound does not cover
1404 // curved faces and is not a valid bound on the true distance. Fall back to
1405 // the BVH lower bound which covers all triangulated faces.
1406 const G4double bvhDist = BVHLowerBoundDistance(p);
1407 d = (bvhDist < kInfinity) ? bvhDist : ExactDistanceToOut(p);
1408 }
1409 // Feed the inscribed-sphere cache: every positive return value proves B(p,d)
1410 // is inside the solid and can accelerate future Inside(p) calls.
1411 TryInsertSphere(p, d);
1412 return d;
1413}
1414
1416 std::unique_lock<std::mutex> lock(fVolumeAreaMutex);
1417 if (!fCachedVolume) {
1418 GProp_GProps props;
1419 BRepGProp::VolumeProperties(fShape, props);
1420 fCachedVolume = props.Mass();
1421 }
1422 return *fCachedVolume;
1423}
1424
1426 std::unique_lock<std::mutex> lock(fVolumeAreaMutex);
1427 if (!fCachedSurfaceArea) {
1428 GProp_GProps props;
1429 BRepGProp::SurfaceProperties(fShape, props);
1430 fCachedSurfaceArea = props.Mass();
1431 }
1432 return *fCachedSurfaceArea;
1433}
1434
1435const G4OCCTSolid::SurfaceSamplingCache& G4OCCTSolid::GetOrBuildSurfaceCache() const {
1436 // Tessellate first, outside the lock: BRepMesh_IncrementalMesh is idempotent
1437 // and calling it from multiple threads simultaneously is safe (extra calls
1438 // after the first are no-ops).
1439 BRepMesh_IncrementalMesh mesher(fShape, kRelativeDeflection, /*isRelative=*/Standard_True);
1440 (void)mesher;
1441
1442 std::unique_lock<std::mutex> lock(fSurfaceCacheMutex);
1443
1444 const std::uint64_t currentGen = fShapeGeneration.load(std::memory_order_acquire);
1445 if (fSurfaceCache.has_value() && fSurfaceCacheGeneration == currentGen) {
1446 return *fSurfaceCache;
1447 }
1448
1449 // Build the cache while holding the lock. Collecting triangle vertices and
1450 // computing areas is fast (just reading the already-computed triangulation)
1451 // so blocking other threads briefly is acceptable.
1452 SurfaceSamplingCache cache;
1453
1454 for (TopExp_Explorer ex(fShape, TopAbs_FACE); ex.More(); ex.Next()) {
1455 const TopoDS_Face& face = TopoDS::Face(ex.Current());
1456 TopLoc_Location loc;
1457 const Handle(Poly_Triangulation) & triangulation = BRep_Tool::Triangulation(face, loc);
1458 if (triangulation.IsNull()) {
1459 continue;
1460 }
1461
1462 // Register this face once in the deduplicated faces array.
1463 const auto faceIndex = static_cast<std::uint32_t>(cache.faces.size());
1464 cache.faces.push_back(face);
1465
1466 const gp_Trsf& transform = loc.Transformation();
1467 const bool reverseWinding = face.Orientation() == TopAbs_REVERSED;
1468
1469 for (Standard_Integer i = 1; i <= triangulation->NbTriangles(); ++i) {
1470 Standard_Integer idx1 = 0;
1471 Standard_Integer idx2 = 0;
1472 Standard_Integer idx3 = 0;
1473 triangulation->Triangle(i).Get(idx1, idx2, idx3);
1474 if (reverseWinding) {
1475 std::swap(idx2, idx3);
1476 }
1477
1478 const gp_Pnt q1 = triangulation->Node(idx1).Transformed(transform);
1479 const gp_Pnt q2 = triangulation->Node(idx2).Transformed(transform);
1480 const gp_Pnt q3 = triangulation->Node(idx3).Transformed(transform);
1481
1482 const G4ThreeVector v1(q1.X(), q1.Y(), q1.Z());
1483 const G4ThreeVector v2(q2.X(), q2.Y(), q2.Z());
1484 const G4ThreeVector v3(q3.X(), q3.Y(), q3.Z());
1485
1486 const G4double area = 0.5 * (v2 - v1).cross(v3 - v1).mag();
1487 if (area > 0.0) {
1488 cache.totalArea += area;
1489 cache.cumulativeAreas.push_back(cache.totalArea);
1490 cache.triangles.push_back({v1, v2, v3, faceIndex});
1491 }
1492 }
1493 }
1494
1495 fSurfaceCache = std::move(cache);
1496 fSurfaceCacheGeneration = currentGen;
1497 return *fSurfaceCache;
1498}
1499
1500G4ThreeVector G4OCCTSolid::GetPointOnSurface() const {
1501 const SurfaceSamplingCache& cache = GetOrBuildSurfaceCache();
1502
1503 if (cache.triangles.empty() || cache.totalArea == 0.0) {
1504 G4ExceptionDescription msg;
1505 msg << "Tessellation of solid \"" << GetName()
1506 << "\" produced no valid triangles; cannot sample a point on the surface.";
1507 G4Exception("G4OCCTSolid::GetPointOnSurface", "GeomMgt1001", FatalException, msg);
1508 return {0.0, 0.0, 0.0}; // unreachable; silences compiler warning
1509 }
1510
1511 // Select a triangle with probability proportional to its area using a
1512 // binary search on the cumulative-area array.
1513 const G4double target = G4UniformRand() * cache.totalArea;
1514 const auto it = std::ranges::lower_bound(cache.cumulativeAreas, target);
1515 const std::size_t idx = std::min(static_cast<std::size_t>(it - cache.cumulativeAreas.begin()),
1516 cache.triangles.size() - 1);
1517 const SurfaceTriangle& chosen = cache.triangles[idx];
1518
1519 // Sample uniformly within the chosen triangle using the standard
1520 // barycentric-coordinate technique: fold the unit square into a triangle
1521 // by reflecting points with r1+r2 > 1.
1522 G4double r1 = G4UniformRand();
1523 G4double r2 = G4UniformRand();
1524 if (r1 + r2 > 1.0) {
1525 r1 = 1.0 - r1;
1526 r2 = 1.0 - r2;
1527 }
1528
1529 const G4ThreeVector tessPoint =
1530 chosen.p1 + r1 * (chosen.p2 - chosen.p1) + r2 * (chosen.p3 - chosen.p1);
1531
1532 // For curved faces the tessellation triangle is a planar chord that lies
1533 // inside the solid. Project the sampled point back to the nearest point on
1534 // the analytical face surface so that the returned point truly lies on the
1535 // boundary and passes the G4PVPlacement::CheckOverlaps() surface test.
1536 const TopoDS_Face& face = cache.faces[chosen.faceIndex];
1537 TopLoc_Location loc;
1538 const Handle(Geom_Surface) geomSurface = BRep_Tool::Surface(face, loc);
1539 if (!geomSurface.IsNull()) {
1540 gp_Pnt tessPointLocal(tessPoint.x(), tessPoint.y(), tessPoint.z());
1541 if (!loc.IsIdentity()) {
1542 tessPointLocal.Transform(loc.Transformation().Inverted());
1543 }
1544 GeomAPI_ProjectPointOnSurf projection(tessPointLocal, geomSurface);
1545 if (projection.NbPoints() > 0) {
1546 gp_Pnt projectedPoint = projection.NearestPoint();
1547 if (!loc.IsIdentity()) {
1548 projectedPoint.Transform(loc.Transformation());
1549 }
1550 return {projectedPoint.X(), projectedPoint.Y(), projectedPoint.Z()};
1551 }
1552 }
1553
1554 return tessPoint;
1555}
1556
1557G4GeometryType G4OCCTSolid::GetEntityType() const { return "G4OCCTSolid"; }
1558
1559G4VisExtent G4OCCTSolid::GetExtent() const {
1560 return {fCachedBounds.min.x(), fCachedBounds.max.x(), fCachedBounds.min.y(),
1561 fCachedBounds.max.y(), fCachedBounds.min.z(), fCachedBounds.max.z()};
1562}
1563
1564void G4OCCTSolid::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const {
1565 pMin = fCachedBounds.min;
1566 pMax = fCachedBounds.max;
1567}
1568
1569G4bool G4OCCTSolid::CalculateExtent(const EAxis pAxis, const G4VoxelLimits& pVoxelLimit,
1570 const G4AffineTransform& pTransform, G4double& pMin,
1571 G4double& pMax) const {
1572 const G4BoundingEnvelope envelope(fCachedBounds.min, fCachedBounds.max);
1573 return envelope.CalculateExtent(pAxis, pVoxelLimit, G4Transform3D(pTransform), pMin, pMax);
1574}
1575
1576void G4OCCTSolid::DescribeYourselfTo(G4VGraphicsScene& scene) const { scene.AddSolid(*this); }
1577
1578G4Polyhedron* G4OCCTSolid::CreatePolyhedron() const {
1579 const auto currentGeneration = fShapeGeneration.load(std::memory_order_acquire);
1580
1581 {
1582 std::unique_lock<std::mutex> lock(fPolyhedronMutex);
1583
1584 // Wait until no other thread is mid-build for a *different* generation, or
1585 // until the cache is already populated for our generation (fast path).
1586 fPolyhedronCV.wait(lock, [this, currentGeneration] {
1587 return !fPolyhedronBuilding || fPolyhedronGeneration == currentGeneration;
1588 });
1589
1590 // Cache hit: return an independent copy to the caller.
1591 if (fCachedPolyhedron && fPolyhedronGeneration == currentGeneration) {
1592 return new G4Polyhedron(*fCachedPolyhedron);
1593 }
1594
1595 // Claim the build slot so other concurrent callers wait for this thread.
1596 fPolyhedronBuilding = true;
1597 }
1598
1599 // ── Tessellation outside the lock ─────────────────────────────────────────
1600 // Use a relative deflection (1 % of each face's bounding-box size) so that
1601 // the mesh density scales with the shape rather than being a fixed world-
1602 // space length that is inappropriate for both very small and very large
1603 // shapes.
1604 BRepMesh_IncrementalMesh mesher(fShape, kRelativeDeflection, /*isRelative=*/Standard_True);
1605 (void)mesher;
1606
1607 G4TessellatedSolid tessellatedSolid(GetName() + "_polyhedron");
1608 G4int facetCount = 0;
1609
1610 for (TopExp_Explorer explorer(fShape, TopAbs_FACE); explorer.More(); explorer.Next()) {
1611 const TopoDS_Face& face = TopoDS::Face(explorer.Current());
1612 TopLoc_Location location;
1613 const Handle(Poly_Triangulation) & triangulation = BRep_Tool::Triangulation(face, location);
1614 if (triangulation.IsNull()) {
1615 continue;
1616 }
1617
1618 const gp_Trsf& transform = location.Transformation();
1619 const bool reverseWinding = face.Orientation() == TopAbs_REVERSED;
1620
1621 for (Standard_Integer triangleIndex = 1; triangleIndex <= triangulation->NbTriangles();
1622 ++triangleIndex) {
1623 Standard_Integer index1 = 0;
1624 Standard_Integer index2 = 0;
1625 Standard_Integer index3 = 0;
1626 triangulation->Triangle(triangleIndex).Get(index1, index2, index3);
1627
1628 if (reverseWinding) {
1629 std::swap(index2, index3);
1630 }
1631
1632 const gp_Pnt point1 = triangulation->Node(index1).Transformed(transform);
1633 const gp_Pnt point2 = triangulation->Node(index2).Transformed(transform);
1634 const gp_Pnt point3 = triangulation->Node(index3).Transformed(transform);
1635
1636 auto* facet =
1637 new G4TriangularFacet(G4ThreeVector(point1.X(), point1.Y(), point1.Z()),
1638 G4ThreeVector(point2.X(), point2.Y(), point2.Z()),
1639 G4ThreeVector(point3.X(), point3.Y(), point3.Z()), ABSOLUTE);
1640 if (!facet->IsDefined() || !tessellatedSolid.AddFacet(facet)) {
1641 delete facet;
1642 continue;
1643 }
1644
1645 ++facetCount;
1646 }
1647 }
1648
1649 // Finalise the tessellated solid and build a polyhedron from it (still
1650 // outside the lock — both objects are local to this thread).
1651 std::unique_ptr<G4Polyhedron> freshPolyhedron;
1652 if (facetCount > 0) {
1653 tessellatedSolid.SetSolidClosed(true);
1654 G4Polyhedron* tmp = tessellatedSolid.GetPolyhedron();
1655 if (tmp != nullptr) {
1656 freshPolyhedron = std::make_unique<G4Polyhedron>(*tmp);
1657 }
1658 }
1659
1660 // ── Write cache and wake waiters ──────────────────────────────────────────
1661 {
1662 std::unique_lock<std::mutex> lock(fPolyhedronMutex);
1663
1664 // Only persist the result if the shape has not been replaced while
1665 // tessellation was running. Re-reading the generation under the mutex
1666 // ensures the check and the write are atomic with respect to SetOCCTShape().
1667 bool cacheWritten = false;
1668 if (freshPolyhedron && fShapeGeneration.load(std::memory_order_acquire) == currentGeneration) {
1669 fCachedPolyhedron = std::make_unique<G4Polyhedron>(*freshPolyhedron);
1670 fPolyhedronGeneration = currentGeneration;
1671 cacheWritten = true;
1672 }
1673
1674 // Clear the build slot *after* the cache write so that threads woken by
1675 // notify_all() see the complete, consistent cache state.
1676 fPolyhedronBuilding = false;
1677 fPolyhedronCV.notify_all();
1678
1679 // Return a copy of the freshly-cached polyhedron when possible; fall back
1680 // to the locally-built one if the shape was replaced mid-tessellation.
1681 if (cacheWritten) {
1682 return new G4Polyhedron(*fCachedPolyhedron);
1683 }
1684 }
1685
1686 return freshPolyhedron ? new G4Polyhedron(*freshPolyhedron) : nullptr;
1687}
1688
1689std::ostream& G4OCCTSolid::StreamInfo(std::ostream& os) const {
1690 os << "-----------------------------------------------------------\n"
1691 << " *** Dump for solid - " << GetName() << " ***\n"
1692 << " ===================================================\n"
1693 << " Solid type: G4OCCTSolid\n"
1694 << " OCCT shape type: " << fShape.ShapeType() << "\n"
1695 << "-----------------------------------------------------------\n";
1696 return os;
1697}
Declaration of G4OCCTSolid.
Geant4 solid wrapping an Open CASCADE Technology (OCCT) TopoDS_Shape.
G4double ExactDistanceToOut(const G4ThreeVector &p) const
EInside Inside(const G4ThreeVector &p) const override
Return kInside, kSurface, or kOutside for point p.
G4double GetCubicVolume() override
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
G4Polyhedron * CreatePolyhedron() const override
Create a polyhedron representation for visualisation.
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
Return the axis-aligned bounding box limits.
G4OCCTSolid(const G4String &name, const TopoDS_Shape &shape)
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
Return the outward unit normal at surface point p.
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const override
Calculate the extent of the solid in the given axis.
~G4OCCTSolid() override
G4double GetSurfaceArea() override
G4GeometryType GetEntityType() const override
Return a string identifying the entity type.
G4ThreeVector GetPointOnSurface() const override
G4VisExtent GetExtent() const override
Return the axis-aligned bounding box extent.
std::ostream & StreamInfo(std::ostream &os) const override
Stream a human-readable description.
static G4OCCTSolid * FromSTEP(const G4String &name, const std::string &path)
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
Describe the solid to the graphics scene.
G4double ExactDistanceToIn(const G4ThreeVector &p) const