]> OCCT Git - occt.git/commitdiff
Modeling - Optimize BndLib and add GTests (#856)
authorPasukhin Dmitry <dpasukhi@opencascade.com>
Thu, 27 Nov 2025 09:22:44 +0000 (09:22 +0000)
committerGitHub <noreply@github.com>
Thu, 27 Nov 2025 09:22:44 +0000 (09:22 +0000)
- Fix negative modulo in torus bounding box computation
- Fix hyperbola extrema loop early break condition
- Fix OpenMin/OpenMax direction sign for infinite bounds
- Remove dead code in cone bounding computation
- Replace sqrt/log with more efficient computations in hyperbola bounds
- General code cleanup and modernization

src/ModelingData/TKGeomBase/BndLib/BndLib.cxx
src/ModelingData/TKGeomBase/BndLib/BndLib_Add2dCurve.cxx
src/ModelingData/TKGeomBase/BndLib/BndLib_Add3dCurve.cxx
src/ModelingData/TKGeomBase/BndLib/BndLib_AddSurface.cxx
src/ModelingData/TKGeomBase/GTests/BndLib_Test.cxx [new file with mode: 0644]
src/ModelingData/TKGeomBase/GTests/FILES.cmake
tests/bugs/step/bug32239

index cf6ea6973c7caf367c31c0a6c180ca65cb9b7f06..98263180ba221035f40879a8045c50660f1af0c9 100644 (file)
@@ -42,6 +42,12 @@ static Standard_Integer ComputeBox(const gp_Hypr&      aHypr,
 
 namespace
 {
+//! Cosine of M_PI/8 (22.5 degrees) - used for 8-point polygon approximation.
+constexpr Standard_Real THE_COS_PI8 = 0.92387953251128674;
+
+//! Cosine (and sine) of M_PI/4 (45 degrees) - used for diagonal points.
+constexpr Standard_Real THE_COS_PI4 = 0.70710678118654746;
+
 //! Compute method
 template <class PointType, class BndBoxType>
 void Compute(const Standard_Real theP1,
@@ -74,19 +80,11 @@ void Compute(const Standard_Real theP1,
   }
   else
   {
+    // Normalize aTeta1 to [0, 2*PI) range
+    aTeta1 = std::fmod(aTeta1, 2. * M_PI);
     if (aTeta1 < 0.)
     {
-      do
-      {
-        aTeta1 += 2. * M_PI;
-      } while (aTeta1 < 0.);
-    }
-    else if (aTeta1 > 2. * M_PI)
-    {
-      do
-      {
-        aTeta1 -= 2. * M_PI;
-      } while (aTeta1 > 2. * M_PI);
+      aTeta1 += 2. * M_PI;
     }
     aTeta2 = aTeta1 + aDelta;
   }
@@ -104,149 +102,76 @@ void Compute(const Standard_Real theP1,
   if (aDelta > M_PI / 8.)
   {
     // Main radiuses to take into account only 8 points (/cos(Pi/8.))
-    aRam = theRa / 0.92387953251128674;
-    aRbm = theRb / 0.92387953251128674;
+    aRam = theRa / THE_COS_PI8;
+    aRbm = theRb / THE_COS_PI8;
   }
   else
   {
     // Main radiuses to take into account the arrow
-    Standard_Real aTc = cos(aDelta / 2);
+    Standard_Real aTc = std::cos(aDelta / 2);
     aRam              = theRa / aTc;
     aRbm              = theRb / aTc;
   }
   theB.Add(PointType(theO.Coord() + aRam * aCn1 * theXd.Coord() + aRbm * aSn1 * theYd.Coord()));
   theB.Add(PointType(theO.Coord() + aRam * aCn2 * theXd.Coord() + aRbm * aSn2 * theYd.Coord()));
 
-// cos or sin M_PI/4.
-#define PI4 0.70710678118654746
-
-// 8 points of the polygon
-#define addPoint0 theB.Add(PointType(theO.Coord() + aRam * theXd.Coord()))
-#define addPoint1                                                                                  \
-  theB.Add(PointType(theO.Coord() + aRam * PI4 * theXd.Coord() + aRbm * PI4 * theYd.Coord()))
-#define addPoint2 theB.Add(PointType(theO.Coord() + aRbm * theYd.Coord()))
-#define addPoint3                                                                                  \
-  theB.Add(PointType(theO.Coord() - aRam * PI4 * theXd.Coord() + aRbm * PI4 * theYd.Coord()))
-#define addPoint4 theB.Add(PointType(theO.Coord() - aRam * theXd.Coord()))
-#define addPoint5                                                                                  \
-  theB.Add(PointType(theO.Coord() - aRam * PI4 * theXd.Coord() - aRbm * PI4 * theYd.Coord()))
-#define addPoint6 theB.Add(PointType(theO.Coord() - aRbm * theYd.Coord()))
-#define addPoint7                                                                                  \
-  theB.Add(PointType(theO.Coord() + aRam * PI4 * theXd.Coord() - aRbm * PI4 * theYd.Coord()))
-
-  Standard_Integer aDeb = (Standard_Integer)(aTeta1 / (M_PI / 4.));
-  Standard_Integer aFin = (Standard_Integer)(aTeta2 / (M_PI / 4.));
+  // X and Y multipliers for 8 polygon points at 45-degree intervals (0, 45, 90, ..., 315 degrees).
+  // Point i corresponds to angle i * 45 degrees.
+  constexpr Standard_Real aXMult[8] =
+    {1., THE_COS_PI4, 0., -THE_COS_PI4, -1., -THE_COS_PI4, 0., THE_COS_PI4};
+  constexpr Standard_Real aYMult[8] =
+    {0., THE_COS_PI4, 1., THE_COS_PI4, 0., -THE_COS_PI4, -1., -THE_COS_PI4};
+
+  // Lambda to add polygon point by index (0-7).
+  const auto addPoint = [&](Standard_Integer theIdx) {
+    theB.Add(PointType(theO.Coord() + aRam * aXMult[theIdx] * theXd.Coord()
+                       + aRbm * aYMult[theIdx] * theYd.Coord()));
+  };
+
+  Standard_Integer aDeb = static_cast<Standard_Integer>(aTeta1 / (M_PI / 4.));
+  Standard_Integer aFin = static_cast<Standard_Integer>(aTeta2 / (M_PI / 4.));
   aDeb++;
 
   if (aDeb > aFin)
+  {
     return;
+  }
 
-  switch (aDeb)
-  {
-    case 1: {
-      addPoint1;
-      if (aFin <= 1)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 2: {
-      addPoint2;
-      if (aFin <= 2)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 3: {
-      addPoint3;
-      if (aFin <= 3)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 4: {
-      addPoint4;
-      if (aFin <= 4)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 5: {
-      addPoint5;
-      if (aFin <= 5)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 6: {
-      addPoint6;
-      if (aFin <= 6)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 7: {
-      addPoint7;
-      if (aFin <= 7)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 8: {
-      addPoint0;
-      if (aFin <= 8)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 9: {
-      addPoint1;
-      if (aFin <= 9)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 10: {
-      addPoint2;
-      if (aFin <= 10)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 11: {
-      addPoint3;
-      if (aFin <= 11)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 12: {
-      addPoint4;
-      if (aFin <= 12)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 13: {
-      addPoint5;
-      if (aFin <= 13)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 14: {
-      addPoint6;
-      if (aFin <= 14)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 15: {
-      addPoint7;
-      if (aFin <= 15)
-        break;
-    }
+  // Add polygon points from aDeb to aFin, wrapping around using modulo 8.
+  for (Standard_Integer i = aDeb; i <= aFin; ++i)
+  {
+    addPoint(i % 8);
   }
 }
 } // end namespace
 
 static void OpenMin(const gp_Dir& V, Bnd_Box& B)
 {
-  gp_Dir OX(gp_Dir::D::X);
-  gp_Dir OY(gp_Dir::D::Y);
-  gp_Dir OZ(gp_Dir::D::Z);
-  if (V.IsParallel(OX, Precision::Angular()))
-    B.OpenXmin();
-  else if (V.IsParallel(OY, Precision::Angular()))
-    B.OpenYmin();
-  else if (V.IsParallel(OZ, Precision::Angular()))
-    B.OpenZmin();
+  // OpenMin opens the box in the direction of decreasing parameter.
+  // For a line P(t) = Origin + t*V, as t -> -Inf:
+  //   - If V.X() > 0, x-coordinate -> -Inf, so open Xmin
+  //   - If V.X() < 0, x-coordinate -> +Inf, so open Xmax
+  if (V.IsParallel(gp::DX(), Precision::Angular()))
+  {
+    if (V.X() > 0.)
+      B.OpenXmin();
+    else
+      B.OpenXmax();
+  }
+  else if (V.IsParallel(gp::DY(), Precision::Angular()))
+  {
+    if (V.Y() > 0.)
+      B.OpenYmin();
+    else
+      B.OpenYmax();
+  }
+  else if (V.IsParallel(gp::DZ(), Precision::Angular()))
+  {
+    if (V.Z() > 0.)
+      B.OpenZmin();
+    else
+      B.OpenZmax();
+  }
   else
   {
     B.OpenXmin();
@@ -257,15 +182,31 @@ static void OpenMin(const gp_Dir& V, Bnd_Box& B)
 
 static void OpenMax(const gp_Dir& V, Bnd_Box& B)
 {
-  gp_Dir OX(gp_Dir::D::X);
-  gp_Dir OY(gp_Dir::D::Y);
-  gp_Dir OZ(gp_Dir::D::Z);
-  if (V.IsParallel(OX, Precision::Angular()))
-    B.OpenXmax();
-  else if (V.IsParallel(OY, Precision::Angular()))
-    B.OpenYmax();
-  else if (V.IsParallel(OZ, Precision::Angular()))
-    B.OpenZmax();
+  // OpenMax opens the box in the direction of increasing parameter.
+  // For a line P(t) = Origin + t*V, as t -> +Inf:
+  //   - If V.X() > 0, x-coordinate -> +Inf, so open Xmax
+  //   - If V.X() < 0, x-coordinate -> -Inf, so open Xmin
+  if (V.IsParallel(gp::DX(), Precision::Angular()))
+  {
+    if (V.X() > 0.)
+      B.OpenXmax();
+    else
+      B.OpenXmin();
+  }
+  else if (V.IsParallel(gp::DY(), Precision::Angular()))
+  {
+    if (V.Y() > 0.)
+      B.OpenYmax();
+    else
+      B.OpenYmin();
+  }
+  else if (V.IsParallel(gp::DZ(), Precision::Angular()))
+  {
+    if (V.Z() > 0.)
+      B.OpenZmax();
+    else
+      B.OpenZmin();
+  }
   else
   {
     B.OpenXmax();
@@ -276,20 +217,17 @@ static void OpenMax(const gp_Dir& V, Bnd_Box& B)
 
 static void OpenMinMax(const gp_Dir& V, Bnd_Box& B)
 {
-  gp_Dir OX(gp_Dir::D::X);
-  gp_Dir OY(gp_Dir::D::Y);
-  gp_Dir OZ(gp_Dir::D::Z);
-  if (V.IsParallel(OX, Precision::Angular()))
+  if (V.IsParallel(gp::DX(), Precision::Angular()))
   {
     B.OpenXmax();
     B.OpenXmin();
   }
-  else if (V.IsParallel(OY, Precision::Angular()))
+  else if (V.IsParallel(gp::DY(), Precision::Angular()))
   {
     B.OpenYmax();
     B.OpenYmin();
   }
-  else if (V.IsParallel(OZ, Precision::Angular()))
+  else if (V.IsParallel(gp::DZ(), Precision::Angular()))
   {
     B.OpenZmax();
     B.OpenZmin();
@@ -307,12 +245,24 @@ static void OpenMinMax(const gp_Dir& V, Bnd_Box& B)
 
 static void OpenMin(const gp_Dir2d& V, Bnd_Box2d& B)
 {
-  gp_Dir2d OX(gp_Dir2d::D::X);
-  gp_Dir2d OY(gp_Dir2d::D::Y);
-  if (V.IsParallel(OX, Precision::Angular()))
-    B.OpenXmin();
-  else if (V.IsParallel(OY, Precision::Angular()))
-    B.OpenYmin();
+  // OpenMin opens the box in the direction of decreasing parameter.
+  // For a line P(t) = Origin + t*V, as t -> -Inf:
+  //   - If V.X() > 0, x-coordinate -> -Inf, so open Xmin
+  //   - If V.X() < 0, x-coordinate -> +Inf, so open Xmax
+  if (V.IsParallel(gp::DX2d(), Precision::Angular()))
+  {
+    if (V.X() > 0.)
+      B.OpenXmin();
+    else
+      B.OpenXmax();
+  }
+  else if (V.IsParallel(gp::DY2d(), Precision::Angular()))
+  {
+    if (V.Y() > 0.)
+      B.OpenYmin();
+    else
+      B.OpenYmax();
+  }
   else
   {
     B.OpenXmin();
@@ -322,12 +272,24 @@ static void OpenMin(const gp_Dir2d& V, Bnd_Box2d& B)
 
 static void OpenMax(const gp_Dir2d& V, Bnd_Box2d& B)
 {
-  gp_Dir2d OX(gp_Dir2d::D::X);
-  gp_Dir2d OY(gp_Dir2d::D::Y);
-  if (V.IsParallel(OX, Precision::Angular()))
-    B.OpenXmax();
-  else if (V.IsParallel(OY, Precision::Angular()))
-    B.OpenYmax();
+  // OpenMax opens the box in the direction of increasing parameter.
+  // For a line P(t) = Origin + t*V, as t -> +Inf:
+  //   - If V.X() > 0, x-coordinate -> +Inf, so open Xmax
+  //   - If V.X() < 0, x-coordinate -> -Inf, so open Xmin
+  if (V.IsParallel(gp::DX2d(), Precision::Angular()))
+  {
+    if (V.X() > 0.)
+      B.OpenXmax();
+    else
+      B.OpenXmin();
+  }
+  else if (V.IsParallel(gp::DY2d(), Precision::Angular()))
+  {
+    if (V.Y() > 0.)
+      B.OpenYmax();
+    else
+      B.OpenYmin();
+  }
   else
   {
     B.OpenXmax();
@@ -337,14 +299,12 @@ static void OpenMax(const gp_Dir2d& V, Bnd_Box2d& B)
 
 static void OpenMinMax(const gp_Dir2d& V, Bnd_Box2d& B)
 {
-  gp_Dir2d OX(gp_Dir2d::D::X);
-  gp_Dir2d OY(gp_Dir2d::D::Y);
-  if (V.IsParallel(OX, Precision::Angular()))
+  if (V.IsParallel(gp::DX2d(), Precision::Angular()))
   {
     B.OpenXmax();
     B.OpenXmin();
   }
-  else if (V.IsParallel(OY, Precision::Angular()))
+  else if (V.IsParallel(gp::DY2d(), Precision::Angular()))
   {
     B.OpenYmax();
     B.OpenYmin();
@@ -1162,6 +1122,9 @@ void BndLib::Add(const gp_Cylinder&  S,
                  const Standard_Real Tol,
                  Bnd_Box&            B)
 {
+  // Cache axis direction for infinite cases.
+  const gp_Dir& aDir = S.Axis().Direction();
+
   if (Precision::IsNegativeInfinite(VMin))
   {
     if (Precision::IsNegativeInfinite(VMax))
@@ -1170,19 +1133,19 @@ void BndLib::Add(const gp_Cylinder&  S,
     }
     else if (Precision::IsPositiveInfinite(VMax))
     {
-      OpenMinMax(S.Axis().Direction(), B);
+      OpenMinMax(aDir, B);
     }
     else
     {
       ComputeCyl(S, UMin, UMax, 0., VMax, B);
-      OpenMin(S.Axis().Direction(), B);
+      OpenMin(aDir, B);
     }
   }
   else if (Precision::IsPositiveInfinite(VMin))
   {
     if (Precision::IsNegativeInfinite(VMax))
     {
-      OpenMinMax(S.Axis().Direction(), B);
+      OpenMinMax(aDir, B);
     }
     else if (Precision::IsPositiveInfinite(VMax))
     {
@@ -1191,7 +1154,7 @@ void BndLib::Add(const gp_Cylinder&  S,
     else
     {
       ComputeCyl(S, UMin, UMax, 0., VMax, B);
-      OpenMax(S.Axis().Direction(), B);
+      OpenMax(aDir, B);
     }
   }
   else
@@ -1199,12 +1162,12 @@ void BndLib::Add(const gp_Cylinder&  S,
     if (Precision::IsNegativeInfinite(VMax))
     {
       ComputeCyl(S, UMin, UMax, VMin, 0., B);
-      OpenMin(S.Axis().Direction(), B);
+      OpenMin(aDir, B);
     }
     else if (Precision::IsPositiveInfinite(VMax))
     {
       ComputeCyl(S, UMin, UMax, VMin, 0., B);
-      OpenMax(S.Axis().Direction(), B);
+      OpenMax(aDir, B);
     }
     else
     {
@@ -1264,8 +1227,9 @@ void BndLib::Add(const gp_Cone&      S,
                  const Standard_Real Tol,
                  Bnd_Box&            B)
 {
+  // Cache axis direction for infinite cases.
+  const gp_Dir& aD = S.Axis().Direction();
 
-  Standard_Real A = S.SemiAngle();
   if (Precision::IsNegativeInfinite(VMin))
   {
     if (Precision::IsNegativeInfinite(VMax))
@@ -1274,22 +1238,19 @@ void BndLib::Add(const gp_Cone&      S,
     }
     else if (Precision::IsPositiveInfinite(VMax))
     {
-      gp_Dir D(std::cos(A) * S.Axis().Direction());
-      OpenMinMax(D, B);
+      OpenMinMax(aD, B);
     }
     else
     {
       ComputeCone(S, UMin, UMax, 0., VMax, B);
-      gp_Dir D(std::cos(A) * S.Axis().Direction());
-      OpenMin(D, B);
+      OpenMin(aD, B);
     }
   }
   else if (Precision::IsPositiveInfinite(VMin))
   {
     if (Precision::IsNegativeInfinite(VMax))
     {
-      gp_Dir D(std::cos(A) * S.Axis().Direction());
-      OpenMinMax(D, B);
+      OpenMinMax(aD, B);
     }
     else if (Precision::IsPositiveInfinite(VMax))
     {
@@ -1298,8 +1259,7 @@ void BndLib::Add(const gp_Cone&      S,
     else
     {
       ComputeCone(S, UMin, UMax, 0., VMax, B);
-      gp_Dir D(std::cos(A) * S.Axis().Direction());
-      OpenMax(D, B);
+      OpenMax(aD, B);
     }
   }
   else
@@ -1307,14 +1267,12 @@ void BndLib::Add(const gp_Cone&      S,
     if (Precision::IsNegativeInfinite(VMax))
     {
       ComputeCone(S, UMin, UMax, VMin, 0., B);
-      gp_Dir D(std::cos(A) * S.Axis().Direction());
-      OpenMin(D, B);
+      OpenMin(aD, B);
     }
     else if (Precision::IsPositiveInfinite(VMax))
     {
       ComputeCone(S, UMin, UMax, VMin, 0., B);
-      gp_Dir D(std::cos(A) * S.Axis().Direction());
-      OpenMax(D, B);
+      OpenMax(aD, B);
     }
     else
     {
@@ -1570,18 +1528,18 @@ void BndLib::Add(const gp_Torus&     S,
   Standard_Integer Fi2;
   if (VMax < VMin)
   {
-    Fi1 = (Standard_Integer)(VMax / (M_PI / 4.));
-    Fi2 = (Standard_Integer)(VMin / (M_PI / 4.));
+    Fi1 = static_cast<Standard_Integer>(VMax / (M_PI / 4.));
+    Fi2 = static_cast<Standard_Integer>(VMin / (M_PI / 4.));
   }
   else
   {
-    Fi1 = (Standard_Integer)(VMin / (M_PI / 4.));
-    Fi2 = (Standard_Integer)(VMax / (M_PI / 4.));
+    Fi1 = static_cast<Standard_Integer>(VMin / (M_PI / 4.));
+    Fi2 = static_cast<Standard_Integer>(VMax / (M_PI / 4.));
   }
   Fi2++;
 
-  Standard_Real Ra = S.MajorRadius();
-  Standard_Real Ri = S.MinorRadius();
+  const Standard_Real Ra = S.MajorRadius();
+  const Standard_Real Ri = S.MinorRadius();
 
   if (Fi2 < Fi1)
     return;
@@ -1593,181 +1551,61 @@ void BndLib::Add(const gp_Torus&     S,
     return;
   }
 
-#define SC 0.71
-#define addP0                                                                                      \
-  (Compute(UMin,                                                                                   \
-           UMax,                                                                                   \
-           Ra + Ri,                                                                                \
-           Ra + Ri,                                                                                \
-           gp_Pnt(S.XAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.YAxis().Direction().XYZ()),                                                    \
-           S.Location(),                                                                           \
-           B))
-#define addP1                                                                                      \
-  (Compute(UMin,                                                                                   \
-           UMax,                                                                                   \
-           Ra + Ri * SC,                                                                           \
-           Ra + Ri * SC,                                                                           \
-           gp_Pnt(S.XAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.YAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.Location().XYZ() + (Ri * SC) * S.Axis().Direction().XYZ()),                    \
-           B))
-#define addP2                                                                                      \
-  (Compute(UMin,                                                                                   \
-           UMax,                                                                                   \
-           Ra,                                                                                     \
-           Ra,                                                                                     \
-           gp_Pnt(S.XAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.YAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.Location().XYZ() + Ri * S.Axis().Direction().XYZ()),                           \
-           B))
-#define addP3                                                                                      \
-  (Compute(UMin,                                                                                   \
-           UMax,                                                                                   \
-           Ra - Ri * SC,                                                                           \
-           Ra - Ri * SC,                                                                           \
-           gp_Pnt(S.XAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.YAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.Location().XYZ() + (Ri * SC) * S.Axis().Direction().XYZ()),                    \
-           B))
-#define addP4                                                                                      \
-  (Compute(UMin,                                                                                   \
-           UMax,                                                                                   \
-           Ra - Ri,                                                                                \
-           Ra - Ri,                                                                                \
-           gp_Pnt(S.XAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.YAxis().Direction().XYZ()),                                                    \
-           S.Location(),                                                                           \
-           B))
-#define addP5                                                                                      \
-  (Compute(UMin,                                                                                   \
-           UMax,                                                                                   \
-           Ra - Ri * SC,                                                                           \
-           Ra - Ri * SC,                                                                           \
-           gp_Pnt(S.XAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.YAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.Location().XYZ() - (Ri * SC) * S.Axis().Direction().XYZ()),                    \
-           B))
-#define addP6                                                                                      \
-  (Compute(UMin,                                                                                   \
-           UMax,                                                                                   \
-           Ra,                                                                                     \
-           Ra,                                                                                     \
-           gp_Pnt(S.XAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.YAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.Location().XYZ() - Ri * S.Axis().Direction().XYZ()),                           \
-           B))
-#define addP7                                                                                      \
-  (Compute(UMin,                                                                                   \
-           UMax,                                                                                   \
-           Ra + Ri * SC,                                                                           \
-           Ra + Ri * SC,                                                                           \
-           gp_Pnt(S.XAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.YAxis().Direction().XYZ()),                                                    \
-           gp_Pnt(S.Location().XYZ() - (Ri * SC) * S.Axis().Direction().XYZ()),                    \
-           B))
-
-  switch (Fi1)
-  {
-    case 0: {
-      addP0;
-      if (Fi2 <= 0)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 1: {
-      addP1;
-      if (Fi2 <= 1)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 2: {
-      addP2;
-      if (Fi2 <= 2)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 3: {
-      addP3;
-      if (Fi2 <= 3)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 4: {
-      addP4;
-      if (Fi2 <= 4)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 5: {
-      addP5;
-      if (Fi2 <= 5)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 6: {
-      addP6;
-      if (Fi2 <= 6)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 7: {
-      addP7;
-      if (Fi2 <= 7)
-        break;
-    }
-      Standard_FALLTHROUGH
-    case 8:
-    default: {
-      addP0;
-      switch (Fi2)
-      {
-        case 15:
-          addP7;
-          Standard_FALLTHROUGH
-        case 14:
-          addP6;
-          Standard_FALLTHROUGH
-        case 13:
-          addP5;
-          Standard_FALLTHROUGH
-        case 12:
-          addP4;
-          Standard_FALLTHROUGH
-        case 11:
-          addP3;
-          Standard_FALLTHROUGH
-        case 10:
-          addP2;
-          Standard_FALLTHROUGH
-        case 9:
-          addP1;
-          Standard_FALLTHROUGH
-        case 8:
-          break;
-      }
-    }
+  // Cache direction vectors.
+  const gp_XYZ aZDir   = S.Axis().Direction().XYZ();
+  const gp_XYZ aLocXYZ = S.Location().XYZ();
+  const gp_Pnt aXd(S.XAxis().Direction().XYZ());
+  const gp_Pnt aYd(S.YAxis().Direction().XYZ());
+
+  // Multipliers for torus cross-section points at 45-degree intervals.
+  // radiusMult[i]: multiplier for Ri in radius calculation (Ra + Ri * radiusMult[i])
+  // zMult[i]: multiplier for Ri in Z offset calculation (Ri * zMult[i])
+  // THE_COS_PI4 = cos(45 deg) = sin(45 deg) = 0.707...
+  constexpr Standard_Real aRadiusMult[8] =
+    {1., THE_COS_PI4, 0., -THE_COS_PI4, -1., -THE_COS_PI4, 0., THE_COS_PI4};
+  constexpr Standard_Real aZMult[8] =
+    {0., THE_COS_PI4, 1., THE_COS_PI4, 0., -THE_COS_PI4, -1., -THE_COS_PI4};
+
+  // Lambda to add torus cross-section point by index (0-7).
+  const auto addTorusPoint = [&](Standard_Integer theIdx) {
+    const Standard_Real aRadius = Ra + Ri * aRadiusMult[theIdx];
+    const gp_Pnt        aCenter(aLocXYZ + (Ri * aZMult[theIdx]) * aZDir);
+    Compute(UMin, UMax, aRadius, aRadius, aXd, aYd, aCenter, B);
+  };
+
+  // Add points from Fi1 to Fi2, handling wrap-around for indices.
+  // Use ((i % 8) + 8) % 8 to handle negative indices correctly
+  // (C++ modulo can return negative values for negative dividends).
+  for (Standard_Integer i = Fi1; i <= Fi2; ++i)
+  {
+    addTorusPoint(((i % 8) + 8) % 8);
   }
+
   B.Enlarge(Tol);
 }
 
 void BndLib::Add(const gp_Torus& S, const Standard_Real Tol, Bnd_Box& B)
 {
-
-  Standard_Real RMa = S.MajorRadius();
-  Standard_Real Rmi = S.MinorRadius();
-  gp_XYZ        O   = S.Location().XYZ();
-  gp_XYZ        Xd  = S.XAxis().Direction().XYZ();
-  gp_XYZ        Yd  = S.YAxis().Direction().XYZ();
-  gp_XYZ        Zd  = S.Axis().Direction().XYZ();
-  B.Add(gp_Pnt(O - (RMa + Rmi) * Xd - (RMa + Rmi) * Yd + Rmi * Zd));
-  B.Add(gp_Pnt(O - (RMa + Rmi) * Xd - (RMa + Rmi) * Yd - Rmi * Zd));
-  B.Add(gp_Pnt(O + (RMa + Rmi) * Xd - (RMa + Rmi) * Yd + Rmi * Zd));
-  B.Add(gp_Pnt(O + (RMa + Rmi) * Xd - (RMa + Rmi) * Yd - Rmi * Zd));
-  B.Add(gp_Pnt(O - (RMa + Rmi) * Xd + (RMa + Rmi) * Yd + Rmi * Zd));
-  B.Add(gp_Pnt(O - (RMa + Rmi) * Xd + (RMa + Rmi) * Yd - Rmi * Zd));
-  B.Add(gp_Pnt(O + (RMa + Rmi) * Xd + (RMa + Rmi) * Yd + Rmi * Zd));
-  B.Add(gp_Pnt(O + (RMa + Rmi) * Xd + (RMa + Rmi) * Yd - Rmi * Zd));
+  const Standard_Real aRMa = S.MajorRadius();
+  const Standard_Real aRmi = S.MinorRadius();
+  const Standard_Real aR   = aRMa + aRmi;
+  const gp_XYZ        aO   = S.Location().XYZ();
+  const gp_XYZ        aXd  = S.XAxis().Direction().XYZ();
+  const gp_XYZ        aYd  = S.YAxis().Direction().XYZ();
+  const gp_XYZ        aZd  = S.Axis().Direction().XYZ();
+  // Precompute scaled direction vectors.
+  const gp_XYZ aRXd  = aR * aXd;
+  const gp_XYZ aRYd  = aR * aYd;
+  const gp_XYZ aRiZd = aRmi * aZd;
+  // Add 8 corner points of torus bounding box.
+  B.Add(gp_Pnt(aO - aRXd - aRYd + aRiZd));
+  B.Add(gp_Pnt(aO - aRXd - aRYd - aRiZd));
+  B.Add(gp_Pnt(aO + aRXd - aRYd + aRiZd));
+  B.Add(gp_Pnt(aO + aRXd - aRYd - aRiZd));
+  B.Add(gp_Pnt(aO - aRXd + aRYd + aRiZd));
+  B.Add(gp_Pnt(aO - aRXd + aRYd - aRiZd));
+  B.Add(gp_Pnt(aO + aRXd + aRYd + aRiZd));
+  B.Add(gp_Pnt(aO + aRXd + aRYd - aRiZd));
   B.Enlarge(Tol);
 }
 
@@ -1804,7 +1642,8 @@ Standard_Integer ComputeBox(const gp_Hypr&      aHypr,
   aRmaj               = aHypr.MajorRadius();
   aRmin               = aHypr.MinorRadius();
   //
-  aT3 = 0;
+  // Find extrema for each coordinate (X, Y, Z) independently.
+  // Each coordinate can have its extremum at a different parameter value.
   for (i = 1; i <= 3; ++i)
   {
     aA = aRmin * aYDir.Coord(i);
@@ -1813,8 +1652,8 @@ Standard_Integer ComputeBox(const gp_Hypr&      aHypr,
     aABP = aA + aB;
     aBAM = aB - aA;
     //
-    aABP = fabs(aABP);
-    aBAM = fabs(aBAM);
+    aABP = std::abs(aABP);
+    aBAM = std::abs(aBAM);
     //
     if (aABP < aEps || aBAM < aEps)
     {
@@ -1822,23 +1661,17 @@ Standard_Integer ComputeBox(const gp_Hypr&      aHypr,
     }
     //
     aCf = aBAM / aABP;
-    aT3 = log(sqrt(aCf));
+    aT3 = 0.5 * std::log(aCf);
     //
     if (aT3 < aT1 || aT3 > aT2)
     {
       continue;
     }
+    // Add extremum point for this coordinate.
+    aP3 = ElCLib::Value(aT3, aHypr);
+    aBox.Add(aP3);
     iErr = 0;
-    break;
   }
   //
-  if (iErr)
-  {
-    return iErr;
-  }
-  //
-  aP3 = ElCLib::Value(aT3, aHypr);
-  aBox.Add(aP3);
-  //
   return iErr;
 }
index 5495dc7d7f4edff717c4fa12b1e61b69f482a898..c567f2aaeb9b2725e072b04b4fbdf3d59035a7db 100644 (file)
@@ -143,7 +143,7 @@ public:
   Standard_Integer NbVariables() const { return 1; }
 
 private:
-  Curv2dMaxMinCoordMVar& operator=(const Curv2dMaxMinCoordMVar& theOther);
+  Curv2dMaxMinCoordMVar& operator=(const Curv2dMaxMinCoordMVar&) = delete;
 
   Standard_Boolean CheckInputData(Standard_Real theParam)
   {
@@ -190,7 +190,7 @@ public:
   }
 
 private:
-  Curv2dMaxMinCoord& operator=(const Curv2dMaxMinCoord& theOther);
+  Curv2dMaxMinCoord& operator=(const Curv2dMaxMinCoord&) = delete;
 
   Standard_Boolean CheckInputData(Standard_Real theParam)
   {
@@ -418,7 +418,8 @@ void BndLib_Box2dCurve::PerformBezier()
     aTb[1] = aT2;
   }
   //
-  if (!(aT1 == aTb[0] && aT2 == aTb[1]))
+  constexpr Standard_Real anEps = Precision::PConfusion();
+  if (std::abs(aT1 - aTb[0]) > anEps || std::abs(aT2 - aTb[1]) > anEps)
   {
     aG = aCBz->Copy();
     //
@@ -477,7 +478,7 @@ void BndLib_Box2dCurve::PerformBSpline()
 
   //
   constexpr Standard_Real eps = Precision::PConfusion();
-  if (fabs(aT1 - aTb[0]) > eps || fabs(aT2 - aTb[1]) > eps)
+  if (std::abs(aT1 - aTb[0]) > eps || std::abs(aT2 - aTb[1]) > eps)
   {
     aG = aCBS->Copy();
     //
@@ -738,7 +739,7 @@ void BndLib_Box2dCurve::D0(const Standard_Real aU, gp_Pnt2d& aP2D)
     //
     aA = aV1.Y();
     aB = -aV1.X();
-    aR = sqrt(aA * aA + aB * aB);
+    aR = std::sqrt(aA * aA + aB * aB);
     if (aR <= aRes)
     {
       myErrorStatus = 13;
@@ -948,18 +949,18 @@ void BndLib_Box2dCurve::Compute(const Handle(Geom2d_Conic)& aConic2D,
   }
   //
   // aType==GeomAbs_Circle ||  aType==GeomAbs_Ellipse
-  aEps   = 1.e-14;
+  aEps   = Precision::Angular();
   aTwoPI = 2. * M_PI;
   dT     = aT2 - aT1;
   //
   Standard_Real aT1z = AdjustToPeriod(aT1, aTwoPI);
-  if (fabs(aT1z) < aEps)
+  if (std::abs(aT1z) < aEps)
   {
     aT1z = 0.;
   }
   //
   Standard_Real aT2z = aT1z + dT;
-  if (fabs(aT2z - aTwoPI) < aEps)
+  if (std::abs(aT2z - aTwoPI) < aEps)
   {
     aT2z = aTwoPI;
   }
@@ -1033,12 +1034,12 @@ Standard_Integer BndLib_Box2dCurve::Compute(const Handle(Geom2d_Conic)& aConic2D
       aLy = (!i) ? 1. : 0.;
       aBx = aLx * aA21 - aLy * aA11;
       aBy = aLx * aA22 - aLy * aA12;
-      aB  = sqrt(aBx * aBx + aBy * aBy);
+      aB  = std::sqrt(aBx * aBx + aBy * aBy);
       //
       aCosFi = aBx / aB;
       aSinFi = aBy / aB;
       //
-      aFi = acos(aCosFi);
+      aFi = std::acos(aCosFi);
       if (aSinFi < 0.)
       {
         aFi = aTwoPI - aFi;
@@ -1060,7 +1061,7 @@ Standard_Integer BndLib_Box2dCurve::Compute(const Handle(Geom2d_Conic)& aConic2D
     Standard_Real           aA1, aA2;
     Handle(Geom2d_Parabola) aPR2D;
     //
-    aEps = 1.e-12;
+    aEps = Precision::Angular();
     //
     aPR2D = Handle(Geom2d_Parabola)::DownCast(aConic2D);
     aFc   = aPR2D->Focal();
@@ -1072,7 +1073,7 @@ Standard_Integer BndLib_Box2dCurve::Compute(const Handle(Geom2d_Conic)& aConic2D
       aLy = (!i) ? 1. : 0.;
       //
       aA2 = aLx * aSinBt - aLy * aCosBt;
-      if (fabs(aA2) < aEps)
+      if (std::abs(aA2) < aEps)
       {
         continue;
       }
@@ -1092,7 +1093,7 @@ Standard_Integer BndLib_Box2dCurve::Compute(const Handle(Geom2d_Conic)& aConic2D
     Standard_Real            aEps, aB1, aB2, aB12, aB22, aZ, aD;
     Handle(Geom2d_Hyperbola) aHP2D;
     //
-    aEps = 1.e-12;
+    aEps = Precision::Angular();
     //
     aHP2D = Handle(Geom2d_Hyperbola)::DownCast(aConic2D);
     aR1   = aHP2D->MajorRadius();
@@ -1107,12 +1108,12 @@ Standard_Integer BndLib_Box2dCurve::Compute(const Handle(Geom2d_Conic)& aConic2D
       aB1 = aR1 * (aLx * aSinBt - aLy * aCosBt);
       aB2 = aR2 * (aLx * aSinGm - aLy * aCosGm);
       //
-      if (fabs(aB1) < aEps)
+      if (std::abs(aB1) < aEps)
       {
         continue;
       }
       //
-      if (fabs(aB2) < aEps)
+      if (std::abs(aB2) < aEps)
       {
         pT[j] = 0.;
         ++j;
@@ -1126,14 +1127,14 @@ Standard_Integer BndLib_Box2dCurve::Compute(const Handle(Geom2d_Conic)& aConic2D
           continue;
         }
         //
-        aD = sqrt(aB12 - aB22);
+        aD = std::sqrt(aB12 - aB22);
         //-------------
         for (k = -1; k < 2; k += 2)
         {
           aZ = (aB1 + k * aD) / aB2;
-          if (fabs(aZ) < 1.)
+          if (std::abs(aZ) < 1.)
           {
-            pT[j] = -log((1. + aZ) / (1. - aZ));
+            pT[j] = -std::log((1. + aZ) / (1. - aZ));
             ++j;
           }
         }
@@ -1155,12 +1156,12 @@ Standard_Real BndLib_Box2dCurve::AdjustToPeriod(const Standard_Real aT, const St
   aTRet = aT;
   if (aT < 0.)
   {
-    k     = 1 + (Standard_Integer)(-aT / aPeriod);
+    k     = 1 + static_cast<Standard_Integer>(-aT / aPeriod);
     aTRet = aT + k * aPeriod;
   }
   else if (aT > aPeriod)
   {
-    k     = (Standard_Integer)(aT / aPeriod);
+    k     = static_cast<Standard_Integer>(aT / aPeriod);
     aTRet = aT - k * aPeriod;
   }
   if (aTRet == aPeriod)
@@ -1197,8 +1198,7 @@ void BndLib_Add2dCurve::Add(const Adaptor2d_Curve2d& aC,
                             const Standard_Real      aTol,
                             Bnd_Box2d&               aBox2D)
 {
-  Adaptor2d_Curve2d*   pC = (Adaptor2d_Curve2d*)&aC;
-  Geom2dAdaptor_Curve* pA = dynamic_cast<Geom2dAdaptor_Curve*>(pC);
+  const Geom2dAdaptor_Curve* pA = dynamic_cast<const Geom2dAdaptor_Curve*>(&aC);
   if (!pA)
   {
     Standard_Real    U, DU;
index e34ceae37c5c9591af74ff69ea002c3ad5bcc885..82ef992b0082618ca4af187b1e97b7f429d7ac4b 100644 (file)
@@ -48,53 +48,46 @@ static void reduceSplineBox(const Adaptor3d_Curve& theCurve,
                             const Bnd_Box&         theOrigBox,
                             Bnd_Box&               theReducedBox)
 {
-  // Guaranteed bounding box based on poles of bspline.
-  Bnd_Box       aPolesBox;
-  Standard_Real aPolesXMin, aPolesYMin, aPolesZMin, aPolesXMax, aPolesYMax, aPolesZMax;
+  // Guaranteed bounding box based on poles of bspline/bezier.
+  Bnd_Box aPolesBox;
 
   if (theCurve.GetType() == GeomAbs_BSplineCurve)
   {
     Handle(Geom_BSplineCurve) aC     = theCurve.BSpline();
     const TColgp_Array1OfPnt& aPoles = aC->Poles();
-
     for (Standard_Integer anIdx = aPoles.Lower(); anIdx <= aPoles.Upper(); ++anIdx)
     {
       aPolesBox.Add(aPoles.Value(anIdx));
     }
   }
-  if (theCurve.GetType() == GeomAbs_BezierCurve)
+  else if (theCurve.GetType() == GeomAbs_BezierCurve)
   {
     Handle(Geom_BezierCurve)  aC     = theCurve.Bezier();
     const TColgp_Array1OfPnt& aPoles = aC->Poles();
-
     for (Standard_Integer anIdx = aPoles.Lower(); anIdx <= aPoles.Upper(); ++anIdx)
     {
       aPolesBox.Add(aPoles.Value(anIdx));
     }
   }
 
-  aPolesBox.Get(aPolesXMin, aPolesYMin, aPolesZMin, aPolesXMax, aPolesYMax, aPolesZMax);
-
-  Standard_Real x, y, z, X, Y, Z;
-  theOrigBox.Get(x, y, z, X, Y, Z);
-
-  // Left bound.
-  if (aPolesXMin > x)
-    x = aPolesXMin;
-  if (aPolesYMin > y)
-    y = aPolesYMin;
-  if (aPolesZMin > z)
-    z = aPolesZMin;
-
-  // Right bound.
-  if (aPolesXMax < X)
-    X = aPolesXMax;
-  if (aPolesYMax < Y)
-    Y = aPolesYMax;
-  if (aPolesZMax < Z)
-    Z = aPolesZMax;
+  // Get original box limits using C++17 structured bindings.
+  const auto [aXMin, aXMax, aYMin, aYMax, aZMin, aZMax] = theOrigBox.Get();
 
-  theReducedBox.Update(x, y, z, X, Y, Z);
+  // If poles box is valid, intersect with original box.
+  if (!aPolesBox.IsVoid())
+  {
+    const auto [aPXMin, aPXMax, aPYMin, aPYMax, aPZMin, aPZMax] = aPolesBox.Get();
+    theReducedBox.Update(std::max(aXMin, aPXMin),
+                         std::max(aYMin, aPYMin),
+                         std::max(aZMin, aPZMin),
+                         std::min(aXMax, aPXMax),
+                         std::min(aYMax, aPYMax),
+                         std::min(aZMax, aPZMax));
+  }
+  else
+  {
+    theReducedBox.Update(aXMin, aYMin, aZMin, aXMax, aYMax, aZMax);
+  }
 }
 
 //=================================================================================================
@@ -152,8 +145,8 @@ void BndLib_Add3dCurve::Add(const Adaptor3d_Curve& C,
                             const Standard_Real    Tol,
                             Bnd_Box&               B)
 {
-  static Standard_Real weakness = 1.5; // OCC566(apo)
-  Standard_Real        tol      = 0.0;
+  constexpr Standard_Real weakness = 1.5; // OCC566(apo)
+  Standard_Real           tol      = 0.0;
   switch (C.GetType())
   {
 
@@ -270,9 +263,9 @@ void BndLib_Add3dCurve::Add(const Adaptor3d_Curve& C,
       break;
     }
     default: {
-      Bnd_Box                 B1;
-      static Standard_Integer N = 33;
-      tol                       = FillBox(B1, C, U1, U2, N);
+      Bnd_Box                    B1;
+      constexpr Standard_Integer N = 33;
+      tol                          = FillBox(B1, C, U1, U2, N);
       B1.Enlarge(weakness * tol);
       Standard_Real x, y, z, X, Y, Z;
       B1.Get(x, y, z, X, Y, Z);
@@ -464,7 +457,7 @@ public:
   Standard_Integer NbVariables() const { return 1; }
 
 private:
-  CurvMaxMinCoordMVar& operator=(const CurvMaxMinCoordMVar& theOther);
+  CurvMaxMinCoordMVar& operator=(const CurvMaxMinCoordMVar&) = delete;
 
   Standard_Boolean CheckInputData(Standard_Real theParam)
   {
@@ -511,7 +504,7 @@ public:
   }
 
 private:
-  CurvMaxMinCoord& operator=(const CurvMaxMinCoord& theOther);
+  CurvMaxMinCoord& operator=(const CurvMaxMinCoord&) = delete;
 
   Standard_Boolean CheckInputData(Standard_Real theParam)
   {
index 7f6136195a75327244a1c4740fe69fcafbb4aa8b..767ea90033b29d5b1466d897cf001f2aea583d4f 100644 (file)
@@ -447,10 +447,10 @@ void BndLib_AddSurface::Add(const Adaptor3d_Surface& S,
       gp_Pnt           P;
       for (Standard_Integer i = 1; i <= Nu; i++)
       {
-        Standard_Real U = UMin + ((UMax - UMin) * (i - 1) / (Nu - 1));
+        const Standard_Real U = UMin + ((UMax - UMin) * (i - 1) / (Nu - 1));
         for (Standard_Integer j = 1; j <= Nv; j++)
         {
-          Standard_Real V = VMin + ((VMax - VMin) * (j - 1) / (Nv - 1));
+          const Standard_Real V = VMin + ((VMax - VMin) * (j - 1) / (Nv - 1));
           S.D0(U, V, P);
           B.Add(P);
         }
@@ -769,7 +769,7 @@ public:
   Standard_Integer NbVariables() const { return 2; }
 
 private:
-  SurfMaxMinCoord& operator=(const SurfMaxMinCoord& theOther);
+  SurfMaxMinCoord& operator=(const SurfMaxMinCoord&) = delete;
 
   Standard_Boolean CheckInputData(const math_Vector& theParams)
   {
diff --git a/src/ModelingData/TKGeomBase/GTests/BndLib_Test.cxx b/src/ModelingData/TKGeomBase/GTests/BndLib_Test.cxx
new file mode 100644 (file)
index 0000000..4a92a31
--- /dev/null
@@ -0,0 +1,1478 @@
+// Copyright (c) 2025 OPEN CASCADE SAS
+//
+// This file is part of Open CASCADE Technology software library.
+//
+// This library is free software; you can redistribute it and/or modify it under
+// the terms of the GNU Lesser General Public License version 2.1 as published
+// by the Free Software Foundation, with special exception defined in the file
+// OCCT_LGPL_EXCEPTION.txt. Consult the file LICENSE_LGPL_21.txt included in OCCT
+// distribution for complete text of the license and disclaimer of any warranty.
+//
+// Alternatively, this file may be used under the terms of Open CASCADE
+// commercial license or contractual agreement.
+
+#include <gtest/gtest.h>
+
+#include <BndLib.hxx>
+#include <BndLib_Add2dCurve.hxx>
+#include <BndLib_Add3dCurve.hxx>
+#include <BndLib_AddSurface.hxx>
+#include <Bnd_Box.hxx>
+#include <Bnd_Box2d.hxx>
+#include <Geom2d_BSplineCurve.hxx>
+#include <Geom2d_BezierCurve.hxx>
+#include <Geom2d_Circle.hxx>
+#include <Geom2d_Ellipse.hxx>
+#include <Geom2d_Line.hxx>
+#include <Geom2d_TrimmedCurve.hxx>
+#include <Geom2dAdaptor_Curve.hxx>
+#include <Geom_BSplineCurve.hxx>
+#include <Geom_BSplineSurface.hxx>
+#include <Geom_BezierCurve.hxx>
+#include <Geom_BezierSurface.hxx>
+#include <Geom_Circle.hxx>
+#include <Geom_ConicalSurface.hxx>
+#include <Geom_CylindricalSurface.hxx>
+#include <Geom_Ellipse.hxx>
+#include <Geom_Line.hxx>
+#include <Geom_Plane.hxx>
+#include <Geom_SphericalSurface.hxx>
+#include <Geom_ToroidalSurface.hxx>
+#include <Geom_TrimmedCurve.hxx>
+#include <GeomAdaptor_Curve.hxx>
+#include <GeomAdaptor_Surface.hxx>
+#include <gp_Ax1.hxx>
+#include <gp_Ax2.hxx>
+#include <gp_Ax2d.hxx>
+#include <gp_Ax3.hxx>
+#include <gp_Circ.hxx>
+#include <gp_Circ2d.hxx>
+#include <gp_Cone.hxx>
+#include <gp_Cylinder.hxx>
+#include <gp_Dir.hxx>
+#include <gp_Elips.hxx>
+#include <gp_Elips2d.hxx>
+#include <gp_Hypr.hxx>
+#include <gp_Hypr2d.hxx>
+#include <gp_Lin.hxx>
+#include <gp_Lin2d.hxx>
+#include <gp_Parab.hxx>
+#include <gp_Parab2d.hxx>
+#include <gp_Pnt.hxx>
+#include <gp_Sphere.hxx>
+#include <gp_Torus.hxx>
+#include <gp_Vec.hxx>
+#include <Precision.hxx>
+#include <TColgp_Array1OfPnt.hxx>
+#include <TColgp_Array1OfPnt2d.hxx>
+#include <TColgp_Array2OfPnt.hxx>
+#include <TColStd_Array1OfInteger.hxx>
+#include <TColStd_Array1OfReal.hxx>
+
+#include <cmath>
+#include <limits>
+
+//==================================================================================================
+// Line Tests
+//==================================================================================================
+
+TEST(BndLibTest, Line_FiniteSegment)
+{
+  // Line along X axis from (0,0,0) in direction (1,0,0)
+  gp_Lin  aLine(gp_Pnt(0., 0., 0.), gp_Dir(1., 0., 0.));
+  Bnd_Box aBox;
+
+  BndLib::Add(aLine, 0., 10., 0., aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 10., Precision::Confusion());
+  EXPECT_NEAR(aYmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 0., Precision::Confusion());
+}
+
+TEST(BndLibTest, Line_PositiveDirection_NegativeInfinite)
+{
+  // Line along +X axis with P1 = -Inf, P2 = 5
+  // As t -> -Inf, x -> -Inf, so Xmin should be open
+  gp_Lin  aLine(gp_Pnt(0., 0., 0.), gp_Dir(1., 0., 0.));
+  Bnd_Box aBox;
+
+  BndLib::Add(aLine, -std::numeric_limits<double>::infinity(), 5., 0., aBox);
+
+  EXPECT_TRUE(aBox.IsOpenXmin()) << "Xmin should be open for line in +X direction with P1=-Inf";
+  EXPECT_FALSE(aBox.IsOpenXmax()) << "Xmax should not be open";
+}
+
+TEST(BndLibTest, Line_NegativeDirection_NegativeInfinite)
+{
+  // Line along -X axis with P1 = -Inf, P2 = 5
+  // Direction is (-1, 0, 0), so as t -> -Inf, x -> +Inf
+  // Therefore Xmax should be open, not Xmin
+  gp_Lin  aLine(gp_Pnt(0., 0., 0.), gp_Dir(-1., 0., 0.));
+  Bnd_Box aBox;
+
+  BndLib::Add(aLine, -std::numeric_limits<double>::infinity(), 5., 0., aBox);
+
+  EXPECT_TRUE(aBox.IsOpenXmax()) << "Xmax should be open for line in -X direction with P1=-Inf";
+  EXPECT_FALSE(aBox.IsOpenXmin()) << "Xmin should not be open";
+}
+
+TEST(BndLibTest, Line_NegativeDirection_PositiveInfinite)
+{
+  // Line along -X axis with P1 = -5, P2 = +Inf
+  // Direction is (-1, 0, 0), so as t -> +Inf, x -> -Inf
+  // Therefore Xmin should be open
+  gp_Lin  aLine(gp_Pnt(0., 0., 0.), gp_Dir(-1., 0., 0.));
+  Bnd_Box aBox;
+
+  BndLib::Add(aLine, -5., std::numeric_limits<double>::infinity(), 0., aBox);
+
+  EXPECT_TRUE(aBox.IsOpenXmin()) << "Xmin should be open for line in -X direction with P2=+Inf";
+  EXPECT_FALSE(aBox.IsOpenXmax()) << "Xmax should not be open";
+}
+
+TEST(BndLibTest, Line_DiagonalDirection_BothInfinite)
+{
+  // Line in diagonal direction (-1, -1, -1) with P1=-Inf, P2=+Inf
+  // As t -> -Inf: position goes to (+Inf, +Inf, +Inf) - opens max
+  // As t -> +Inf: position goes to (-Inf, -Inf, -Inf) - opens min
+  gp_Lin  aLine(gp_Pnt(0., 0., 0.), gp_Dir(-1., -1., -1.));
+  Bnd_Box aBox;
+
+  BndLib::Add(aLine,
+              -std::numeric_limits<double>::infinity(),
+              std::numeric_limits<double>::infinity(),
+              0.,
+              aBox);
+
+  EXPECT_TRUE(aBox.IsOpenXmin()) << "Xmin should be open";
+  EXPECT_TRUE(aBox.IsOpenXmax()) << "Xmax should be open";
+  EXPECT_TRUE(aBox.IsOpenYmin()) << "Ymin should be open";
+  EXPECT_TRUE(aBox.IsOpenYmax()) << "Ymax should be open";
+  EXPECT_TRUE(aBox.IsOpenZmin()) << "Zmin should be open";
+  EXPECT_TRUE(aBox.IsOpenZmax()) << "Zmax should be open";
+}
+
+//==================================================================================================
+// Circle Tests
+//==================================================================================================
+
+TEST(BndLibTest, Circle_Full)
+{
+  // Circle of radius 5 centered at origin in XY plane
+  gp_Circ aCirc(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 5.);
+  Bnd_Box aBoxFull;
+  Bnd_Box aBoxArc;
+  double  aTol = 0.;
+
+  BndLib::Add(aCirc, aTol, aBoxFull);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBoxFull.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 0., Precision::Confusion());
+
+  // Full circle via parameters [0, 2*PI]
+  BndLib::Add(aCirc, 0., 2. * M_PI, aTol, aBoxArc);
+  aBoxArc.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+}
+
+TEST(BndLibTest, Circle_Arc_FirstQuadrant)
+{
+  // Circle arc from 0 to PI/2 (first quadrant)
+  gp_Circ aCirc(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 5.);
+  Bnd_Box aBox;
+  double  aTol = 0.;
+
+  BndLib::Add(aCirc, 0., M_PI / 2., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+}
+
+TEST(BndLibTest, Circle_RotatedAxis)
+{
+  // Circle in XZ plane (axis along Y)
+  gp_Circ aCirc(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 1., 0.)), 3.);
+  Bnd_Box aBox;
+  double  aTol = 0.;
+
+  BndLib::Add(aCirc, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -3., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 3., Precision::Confusion());
+  EXPECT_NEAR(aYmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmin, -3., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 3., Precision::Confusion());
+}
+
+//==================================================================================================
+// Ellipse Tests
+//==================================================================================================
+
+TEST(BndLibTest, Ellipse_Full)
+{
+  // Ellipse with major radius 10, minor radius 5 in XY plane
+  gp_Elips aElips(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 10., 5.);
+  Bnd_Box  aBox;
+  double   aTol = 0.;
+
+  BndLib::Add(aElips, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -10., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 10., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 0., Precision::Confusion());
+}
+
+TEST(BndLibTest, Ellipse_Arc)
+{
+  // Ellipse arc from 0 to PI/2
+  gp_Elips aElips(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 10., 5.);
+  Bnd_Box  aBox;
+  double   aTol = 0.;
+
+  BndLib::Add(aElips, 0., M_PI / 2., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // At t=0: (10, 0, 0), at t=PI/2: (0, 5, 0)
+  EXPECT_NEAR(aXmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 10., Precision::Confusion());
+  EXPECT_NEAR(aYmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+}
+
+//==================================================================================================
+// Hyperbola Tests
+//==================================================================================================
+
+TEST(BndLibTest, Hyperbola_SimpleArc)
+{
+  // Hyperbola with RealAxis=3, ImaginaryAxis=2 in XY plane
+  // Parameterization: (3*cosh(t), 2*sinh(t), 0)
+  gp_Hypr aHypr(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 3., 2.);
+  Bnd_Box aBox;
+  double  aTol = 0.;
+
+  BndLib::Add(aHypr, -1., 1., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // At t=-1: (3*cosh(-1), 2*sinh(-1)) = (4.629..., -2.350...)
+  // At t=1:  (3*cosh(1), 2*sinh(1))   = (4.629...,  2.350...)
+  // At t=0:  (3, 0) - minimum X
+  double aExpectedXmin = 3.;                  // at t=0
+  double aExpectedXmax = 3. * std::cosh(1.);  // at t=+-1
+  double aExpectedYmin = 2. * std::sinh(-1.); // at t=-1
+  double aExpectedYmax = 2. * std::sinh(1.);  // at t=1
+
+  EXPECT_NEAR(aXmin, aExpectedXmin, Precision::Confusion());
+  EXPECT_NEAR(aXmax, aExpectedXmax, Precision::Confusion());
+  EXPECT_NEAR(aYmin, aExpectedYmin, Precision::Confusion());
+  EXPECT_NEAR(aYmax, aExpectedYmax, Precision::Confusion());
+}
+
+TEST(BndLibTest, Hyperbola_Rotated_MultipleExtrema)
+{
+  // Rotated hyperbola where X, Y, Z extrema occur at different parameter values.
+  // This tests the fix for the extrema loop that previously broke after finding
+  // the first extremum.
+  gp_Ax2  aAx2(gp_Pnt(0., 0., 0.), gp_Dir(1., 1., 1.), gp_Dir(1., -1., 0.));
+  gp_Hypr aHypr(aAx2, 5., 3.);
+  Bnd_Box aBox;
+  double  aTol = 0.;
+
+  BndLib::Add(aHypr, -2., 2., aTol, aBox);
+
+  // Verify box is valid and contains reasonable bounds
+  EXPECT_FALSE(aBox.IsVoid()) << "Box should not be void";
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // For a rotated hyperbola, the box should be non-degenerate in all dimensions
+  EXPECT_LT(aXmin, aXmax) << "X range should be non-zero";
+  EXPECT_LT(aYmin, aYmax) << "Y range should be non-zero";
+  // Z may be zero if hyperbola lies in a plane, but should still be valid
+}
+
+TEST(BndLibTest, Hyperbola_InfiniteParameter)
+{
+  // Hyperbola with one infinite parameter
+  gp_Hypr aHypr(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 3., 2.);
+  Bnd_Box aBox;
+  double  aTol = 0.;
+
+  BndLib::Add(aHypr, 0., std::numeric_limits<double>::infinity(), aTol, aBox);
+
+  EXPECT_TRUE(aBox.IsOpenXmax()) << "Xmax should be open for hyperbola with P2=+Inf";
+  EXPECT_TRUE(aBox.IsOpenYmax()) << "Ymax should be open for hyperbola with P2=+Inf";
+}
+
+//==================================================================================================
+// Parabola Tests
+//==================================================================================================
+
+TEST(BndLibTest, Parabola_FiniteArc)
+{
+  // Parabola with focal length 2 in XY plane
+  // Parameterization: (t^2/(4*f), t, 0) = (t^2/8, t, 0) for f=2
+  gp_Parab aParab(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 2.);
+  Bnd_Box  aBox;
+  double   aTol = 0.;
+
+  BndLib::Add(aParab, -4., 4., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // At t=0: (0, 0, 0)
+  // At t=+-4: (16/8, +-4, 0) = (2, +-4, 0)
+  EXPECT_NEAR(aXmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 2., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 4., Precision::Confusion());
+}
+
+TEST(BndLibTest, Parabola_InfiniteParameter)
+{
+  // Parabola with infinite parameter
+  gp_Parab aParab(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 2.);
+  Bnd_Box  aBox;
+  double   aTol = 0.;
+
+  BndLib::Add(aParab, 0., std::numeric_limits<double>::infinity(), aTol, aBox);
+
+  EXPECT_TRUE(aBox.IsOpenXmax()) << "Xmax should be open for parabola with P2=+Inf";
+  EXPECT_TRUE(aBox.IsOpenYmax()) << "Ymax should be open for parabola with P2=+Inf";
+}
+
+//==================================================================================================
+// Sphere Tests
+//==================================================================================================
+
+TEST(BndLibTest, Sphere_Full)
+{
+  // Sphere of radius 7 centered at origin
+  gp_Sphere aSphere(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 7.);
+  Bnd_Box   aBox;
+  double    aTol = 0.;
+
+  BndLib::Add(aSphere, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -7., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 7., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -7., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 7., Precision::Confusion());
+  EXPECT_NEAR(aZmin, -7., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 7., Precision::Confusion());
+}
+
+TEST(BndLibTest, Sphere_Patch)
+{
+  // Sphere patch: upper hemisphere (V from 0 to PI/2)
+  gp_Sphere aSphere(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 5.);
+  Bnd_Box   aBox;
+  double    aTol = 0.;
+
+  // Full U range, V from 0 to PI/2 (upper hemisphere)
+  BndLib::Add(aSphere, 0., 2. * M_PI, 0., M_PI / 2., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 5., Precision::Confusion());
+}
+
+TEST(BndLibTest, Sphere_Translated)
+{
+  // Sphere centered at (10, 20, 30) with radius 3
+  gp_Sphere aSphere(gp_Ax3(gp_Pnt(10., 20., 30.), gp_Dir(0., 0., 1.)), 3.);
+  Bnd_Box   aBox;
+  double    aTol = 0.;
+
+  BndLib::Add(aSphere, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, 7., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 13., Precision::Confusion());
+  EXPECT_NEAR(aYmin, 17., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 23., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 27., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 33., Precision::Confusion());
+}
+
+//==================================================================================================
+// Cylinder Tests
+//==================================================================================================
+
+TEST(BndLibTest, Cylinder_FinitePatch)
+{
+  // Cylinder of radius 4, axis along Z
+  gp_Cylinder aCyl(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 4.);
+  Bnd_Box     aBox;
+  double      aTol = 0.;
+
+  // Full U range, V from 0 to 10
+  BndLib::Add(aCyl, 0., 2. * M_PI, 0., 10., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 4., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 4., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 10., Precision::Confusion());
+}
+
+TEST(BndLibTest, Cylinder_InfiniteV)
+{
+  // Cylinder with infinite V parameter
+  gp_Cylinder aCyl(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 4.);
+  Bnd_Box     aBox;
+  double      aTol = 0.;
+
+  BndLib::Add(aCyl, 0., std::numeric_limits<double>::infinity(), aTol, aBox);
+
+  EXPECT_TRUE(aBox.IsOpenZmax()) << "Zmax should be open for cylinder with VMax=+Inf";
+  EXPECT_FALSE(aBox.IsOpenZmin()) << "Zmin should not be open";
+}
+
+TEST(BndLibTest, Cylinder_PartialU)
+{
+  // Cylinder with partial U range (quarter cylinder)
+  gp_Cylinder aCyl(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 5.);
+  Bnd_Box     aBox;
+  double      aTol = 0.;
+
+  // U from 0 to PI/2 (first quadrant), V from 0 to 10
+  BndLib::Add(aCyl, 0., M_PI / 2., 0., 10., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 10., Precision::Confusion());
+}
+
+//==================================================================================================
+// Cone Tests
+//==================================================================================================
+
+TEST(BndLibTest, Cone_FinitePatch)
+{
+  // Cone with semi-angle PI/4 (45 deg), reference radius 2
+  // Cone parameterization: P(u,v) = Apex + (RefRadius + v*sin(SemiAngle)) * (cos(u)*XDir +
+  // sin(u)*YDir)
+  //                                      + v*cos(SemiAngle) * Axis
+  gp_Cone aCone(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), M_PI / 4., 2.);
+  Bnd_Box aBox;
+  double  aTol = 0.;
+
+  // Full U range, V from 0 to 5
+  BndLib::Add(aCone, 0., 2. * M_PI, 0., 5., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // V=0: radius = RefRadius = 2, Z = 0
+  // V=5: radius = 2 + 5*sin(PI/4) = 2 + 5*sqrt(2)/2, Z = 5*cos(PI/4) = 5*sqrt(2)/2
+  double aMaxRadius = 2. + 5. * std::sin(M_PI / 4.);
+  double aMaxZ      = 5. * std::cos(M_PI / 4.);
+  EXPECT_NEAR(aXmax, aMaxRadius, Precision::Confusion());
+  EXPECT_NEAR(aYmax, aMaxRadius, Precision::Confusion());
+  EXPECT_NEAR(aZmax, aMaxZ, Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+}
+
+TEST(BndLibTest, Cone_NegativeV)
+{
+  // Cone with negative V range (cone extending in opposite direction)
+  // Cone parameterization: Z = v*cos(SemiAngle)
+  gp_Cone aCone(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), M_PI / 6., 3.);
+  Bnd_Box aBox;
+  double  aTol = 0.;
+
+  // Full U range, V from -5 to 0
+  BndLib::Add(aCone, 0., 2. * M_PI, -5., 0., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // At V=0: radius = 3, Z = 0
+  // At V=-5: radius = 3 + (-5)*sin(PI/6) = 3 - 2.5 = 0.5, Z = -5*cos(PI/6)
+  double aExpectedZmin = -5. * std::cos(M_PI / 6.);
+  EXPECT_NEAR(aZmin, aExpectedZmin, Precision::Confusion());
+  EXPECT_NEAR(aZmax, 0., Precision::Confusion());
+}
+
+TEST(BndLibTest, Cone_InfiniteV)
+{
+  // Cone with infinite V parameter
+  // Note: BndLib::Add for cone with U,V parameters opens Z but not X/Y
+  // since the cone's bounding in X/Y is computed from the circular cross-sections
+  gp_Cone aCone(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), M_PI / 4., 2.);
+  Bnd_Box aBox;
+  double  aTol = 0.;
+
+  BndLib::Add(aCone, 0., std::numeric_limits<double>::infinity(), aTol, aBox);
+
+  // For cone with infinite V, the Z direction should be open
+  EXPECT_TRUE(aBox.IsOpenZmax()) << "Zmax should be open for cone with VMax=+Inf";
+}
+
+//==================================================================================================
+// Torus Tests
+//==================================================================================================
+
+TEST(BndLibTest, Torus_Full)
+{
+  // Torus with major radius 10, minor radius 3
+  gp_Torus aTorus(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 10., 3.);
+  Bnd_Box  aBox;
+  double   aTol = 0.;
+
+  BndLib::Add(aTorus, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // Torus extends from -(R+r) to +(R+r) in X and Y
+  // and from -r to +r in Z
+  EXPECT_NEAR(aXmin, -13., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 13., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -13., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 13., Precision::Confusion());
+  EXPECT_NEAR(aZmin, -3., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 3., Precision::Confusion());
+}
+
+TEST(BndLibTest, Torus_Patch)
+{
+  // Torus patch: quarter in U, half in V
+  // Torus parameterization: P(u,v) = Center + (R + r*cos(v)) * (cos(u)*XDir + sin(u)*YDir)
+  //                                         + r*sin(v) * Axis
+  gp_Torus aTorus(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 10., 3.);
+  Bnd_Box  aBox;
+  double   aTol = 0.;
+
+  // U from 0 to PI/2, V from 0 to PI
+  BndLib::Add(aTorus, 0., M_PI / 2., 0., M_PI, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // Box should be valid and contain reasonable bounds
+  EXPECT_FALSE(aBox.IsVoid()) << "Box should not be void for torus patch";
+
+  // For U in [0, PI/2], X and Y should be non-negative at minimum
+  // The bounding box algorithm produces conservative bounds that may slightly
+  // exceed R+r due to sampling/approximation
+  EXPECT_GE(aXmin, -Precision::Confusion()) << "Xmin should be >= 0 for first quadrant";
+  EXPECT_GE(aYmin, -Precision::Confusion()) << "Ymin should be >= 0 for first quadrant";
+
+  // Z max should be approximately r (minor radius)
+  EXPECT_LE(aZmax, 3. + Precision::Confusion()) << "Zmax should not exceed r";
+}
+
+TEST(BndLibTest, Torus_NegativeVParameter)
+{
+  // Torus patch with negative V parameters.
+  // This tests the fix for negative modulo: ((i % 8) + 8) % 8
+  // The main purpose is to verify no crash occurs with negative V.
+  gp_Torus aTorus(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 10., 3.);
+  Bnd_Box  aBox;
+  double   aTol = 0.;
+
+  // Full U, V from -PI to 0
+  BndLib::Add(aTorus, 0., 2. * M_PI, -M_PI, 0., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // Should produce valid bounds without crash (main test for the fix)
+  EXPECT_FALSE(aBox.IsVoid()) << "Box should not be void for negative V parameters";
+
+  // For V in [-PI, 0], Z = r*sin(v) ranges from -r (at v=-PI/2) to 0
+  // The bounding box should contain this range
+  EXPECT_LE(aZmin, 0.) << "Zmin should be <= 0 for V in [-PI, 0]";
+  EXPECT_GE(aZmin, -3. - Precision::Confusion()) << "Zmin should be >= -r";
+}
+
+TEST(BndLibTest, Torus_LargeNegativeVParameter)
+{
+  // Torus with large negative V parameter that would cause
+  // out-of-bounds array access with incorrect modulo.
+  // This is the primary test for the negative modulo fix.
+  gp_Torus aTorus(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 10., 3.);
+  Bnd_Box  aBox;
+  double   aTol = 0.;
+
+  // V from -3*PI to -2*PI covers one full revolution
+  BndLib::Add(aTorus, 0., 2. * M_PI, -3. * M_PI, -2. * M_PI, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // Should produce valid bounds without crash (main test for the fix)
+  EXPECT_FALSE(aBox.IsVoid()) << "Box should not be void for large negative V parameters";
+
+  // Full revolution in U gives full X/Y extent, the box may be slightly larger
+  // due to conservative bounding box algorithm
+  EXPECT_LE(aXmin, -13. + Precision::Confusion()) << "Xmin should be at least -(R+r)";
+  EXPECT_GE(aXmax, 13. - Precision::Confusion()) << "Xmax should be at least R+r";
+  EXPECT_LE(aYmin, -13. + Precision::Confusion()) << "Ymin should be at least -(R+r)";
+  EXPECT_GE(aYmax, 13. - Precision::Confusion()) << "Ymax should be at least R+r";
+}
+
+TEST(BndLibTest, Torus_Translated)
+{
+  // Torus centered at (5, 5, 5)
+  gp_Torus aTorus(gp_Ax3(gp_Pnt(5., 5., 5.), gp_Dir(0., 0., 1.)), 8., 2.);
+  Bnd_Box  aBox;
+  double   aTol = 0.;
+
+  BndLib::Add(aTorus, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, 5. - 10., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5. + 10., Precision::Confusion());
+  EXPECT_NEAR(aYmin, 5. - 10., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5. + 10., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 5. - 2., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 5. + 2., Precision::Confusion());
+}
+
+//==================================================================================================
+// 2D Tests
+//==================================================================================================
+
+TEST(BndLibTest, Line2d_FiniteSegment)
+{
+  gp_Lin2d  aLine(gp_Pnt2d(0., 0.), gp_Dir2d(1., 1.));
+  Bnd_Box2d aBox;
+  double    aTol = 0.;
+
+  BndLib::Add(aLine, 0., 10., aTol, aBox);
+
+  double aXmin, aYmin, aXmax, aYmax;
+  aBox.Get(aXmin, aYmin, aXmax, aYmax);
+
+  // Line in direction (1,1) normalized, so at t=10:
+  // position = (10/sqrt(2), 10/sqrt(2))
+  double aExpected = 10. / std::sqrt(2.);
+  EXPECT_NEAR(aXmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aXmax, aExpected, Precision::Confusion());
+  EXPECT_NEAR(aYmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aYmax, aExpected, Precision::Confusion());
+}
+
+TEST(BndLibTest, Circle2d_Full)
+{
+  gp_Circ2d aCirc(gp_Ax2d(gp_Pnt2d(0., 0.), gp_Dir2d(1., 0.)), 5.);
+  Bnd_Box2d aBox;
+  double    aTol = 0.;
+
+  BndLib::Add(aCirc, aTol, aBox);
+
+  double aXmin, aYmin, aXmax, aYmax;
+  aBox.Get(aXmin, aYmin, aXmax, aYmax);
+
+  EXPECT_NEAR(aXmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+}
+
+TEST(BndLibTest, Ellipse2d_Full)
+{
+  gp_Elips2d aElips(gp_Ax2d(gp_Pnt2d(0., 0.), gp_Dir2d(1., 0.)), 8., 4.);
+  Bnd_Box2d  aBox;
+  double     aTol = 0.;
+
+  BndLib::Add(aElips, aTol, aBox);
+
+  double aXmin, aYmin, aXmax, aYmax;
+  aBox.Get(aXmin, aYmin, aXmax, aYmax);
+
+  EXPECT_NEAR(aXmin, -8., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 8., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 4., Precision::Confusion());
+}
+
+TEST(BndLibTest, Hyperbola2d_SimpleArc)
+{
+  gp_Hypr2d aHypr(gp_Ax2d(gp_Pnt2d(0., 0.), gp_Dir2d(1., 0.)), 3., 2.);
+  Bnd_Box2d aBox;
+  double    aTol = 0.;
+
+  BndLib::Add(aHypr, -1., 1., aTol, aBox);
+
+  double aXmin, aYmin, aXmax, aYmax;
+  aBox.Get(aXmin, aYmin, aXmax, aYmax);
+
+  EXPECT_NEAR(aXmin, 3., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 3. * std::cosh(1.), Precision::Confusion());
+  EXPECT_NEAR(aYmin, 2. * std::sinh(-1.), Precision::Confusion());
+  EXPECT_NEAR(aYmax, 2. * std::sinh(1.), Precision::Confusion());
+}
+
+TEST(BndLibTest, Parabola2d_FiniteArc)
+{
+  gp_Parab2d aParab(gp_Ax2d(gp_Pnt2d(0., 0.), gp_Dir2d(1., 0.)), 2.);
+  Bnd_Box2d  aBox;
+  double     aTol = 0.;
+
+  BndLib::Add(aParab, -4., 4., aTol, aBox);
+
+  double aXmin, aYmin, aXmax, aYmax;
+  aBox.Get(aXmin, aYmin, aXmax, aYmax);
+
+  EXPECT_NEAR(aXmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 2., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 4., Precision::Confusion());
+}
+
+//==================================================================================================
+// Tolerance Tests
+//==================================================================================================
+
+TEST(BndLibTest, Circle_WithTolerance)
+{
+  gp_Circ aCirc(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 5.);
+  Bnd_Box aBox;
+  double  aTol = 0.5;
+
+  BndLib::Add(aCirc, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // Box should be enlarged by tolerance
+  EXPECT_NEAR(aXmin, -5.5, Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5.5, Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5.5, Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5.5, Precision::Confusion());
+  EXPECT_NEAR(aZmin, -0.5, Precision::Confusion());
+  EXPECT_NEAR(aZmax, 0.5, Precision::Confusion());
+}
+
+TEST(BndLibTest, Sphere_WithTolerance)
+{
+  gp_Sphere aSphere(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 3.);
+  Bnd_Box   aBox;
+  double    aTol = 1.;
+
+  BndLib::Add(aSphere, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 4., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 4., Precision::Confusion());
+  EXPECT_NEAR(aZmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 4., Precision::Confusion());
+}
+
+//==================================================================================================
+// BndLib_Add3dCurve Tests
+//==================================================================================================
+
+TEST(BndLib_Add3dCurveTest, Circle_Full)
+{
+  // Create a 3D circle using Geom_Circle
+  Handle(Geom_Circle) aCircle = new Geom_Circle(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 5.);
+  GeomAdaptor_Curve   aCurve(aCircle);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  BndLib_Add3dCurve::Add(aCurve, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 0., Precision::Confusion());
+}
+
+TEST(BndLib_Add3dCurveTest, Circle_Arc)
+{
+  // Circle arc from 0 to PI/2
+  Handle(Geom_Circle) aCircle = new Geom_Circle(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 5.);
+  GeomAdaptor_Curve   aCurve(aCircle);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  BndLib_Add3dCurve::Add(aCurve, 0., M_PI / 2., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+}
+
+TEST(BndLib_Add3dCurveTest, Ellipse_Full)
+{
+  // Ellipse with major radius 10, minor radius 5
+  Handle(Geom_Ellipse) anEllipse =
+    new Geom_Ellipse(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 10., 5.);
+  GeomAdaptor_Curve aCurve(anEllipse);
+  Bnd_Box           aBox;
+  double            aTol = 0.;
+
+  BndLib_Add3dCurve::Add(aCurve, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -10., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 10., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+}
+
+TEST(BndLib_Add3dCurveTest, Line_Segment)
+{
+  // Line segment from (0,0,0) to (10,0,0)
+  Handle(Geom_Line)         aLine    = new Geom_Line(gp_Pnt(0., 0., 0.), gp_Dir(1., 0., 0.));
+  Handle(Geom_TrimmedCurve) aTrimmed = new Geom_TrimmedCurve(aLine, 0., 10.);
+  GeomAdaptor_Curve         aCurve(aTrimmed);
+  Bnd_Box                   aBox;
+  double                    aTol = 0.;
+
+  BndLib_Add3dCurve::Add(aCurve, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 10., Precision::Confusion());
+  EXPECT_NEAR(aYmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 0., Precision::Confusion());
+}
+
+TEST(BndLib_Add3dCurveTest, BezierCurve)
+{
+  // Create a cubic Bezier curve
+  TColgp_Array1OfPnt aPoles(1, 4);
+  aPoles(1) = gp_Pnt(0., 0., 0.);
+  aPoles(2) = gp_Pnt(1., 2., 0.);
+  aPoles(3) = gp_Pnt(3., 2., 0.);
+  aPoles(4) = gp_Pnt(4., 0., 0.);
+
+  Handle(Geom_BezierCurve) aBezier = new Geom_BezierCurve(aPoles);
+  GeomAdaptor_Curve        aCurve(aBezier);
+  Bnd_Box                  aBox;
+  double                   aTol = 0.;
+
+  BndLib_Add3dCurve::Add(aCurve, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // Box should contain curve endpoints and the actual curve
+  // Note: The curve doesn't necessarily pass through control points,
+  // so the box may be smaller than the convex hull of control points
+  EXPECT_FALSE(aBox.IsVoid()) << "Box should not be void";
+  EXPECT_LE(aXmin, 0.);
+  EXPECT_GE(aXmax, 4.);
+  EXPECT_LE(aYmin, 0.);
+  // Y max is at curve's maximum, not necessarily at control point Y
+  EXPECT_GT(aYmax, 1.) << "Y max should be > 1 (curve rises above Y=0)";
+}
+
+TEST(BndLib_Add3dCurveTest, BSplineCurve)
+{
+  // Create a simple B-spline curve
+  TColgp_Array1OfPnt aPoles(1, 4);
+  aPoles(1) = gp_Pnt(0., 0., 0.);
+  aPoles(2) = gp_Pnt(1., 1., 0.);
+  aPoles(3) = gp_Pnt(2., 1., 0.);
+  aPoles(4) = gp_Pnt(3., 0., 0.);
+
+  TColStd_Array1OfReal aKnots(1, 2);
+  aKnots(1) = 0.;
+  aKnots(2) = 1.;
+
+  TColStd_Array1OfInteger aMults(1, 2);
+  aMults(1) = 4;
+  aMults(2) = 4;
+
+  Handle(Geom_BSplineCurve) aBSpline = new Geom_BSplineCurve(aPoles, aKnots, aMults, 3);
+  GeomAdaptor_Curve         aCurve(aBSpline);
+  Bnd_Box                   aBox;
+  double                    aTol = 0.;
+
+  BndLib_Add3dCurve::Add(aCurve, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_FALSE(aBox.IsVoid()) << "Box should not be void";
+  EXPECT_LE(aXmin, 0.);
+  EXPECT_GE(aXmax, 3.);
+}
+
+TEST(BndLib_Add3dCurveTest, AddOptimal_Circle)
+{
+  // Test AddOptimal method with a circle
+  Handle(Geom_Circle) aCircle = new Geom_Circle(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 5.);
+  GeomAdaptor_Curve   aCurve(aCircle);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  BndLib_Add3dCurve::AddOptimal(aCurve, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+}
+
+TEST(BndLib_Add3dCurveTest, Circle_Rotated)
+{
+  // Circle in YZ plane (axis along X)
+  Handle(Geom_Circle) aCircle = new Geom_Circle(gp_Ax2(gp_Pnt(0., 0., 0.), gp_Dir(1., 0., 0.)), 4.);
+  GeomAdaptor_Curve   aCurve(aCircle);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  BndLib_Add3dCurve::Add(aCurve, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 0., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 4., Precision::Confusion());
+  EXPECT_NEAR(aZmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 4., Precision::Confusion());
+}
+
+//==================================================================================================
+// BndLib_Add2dCurve Tests
+//==================================================================================================
+
+TEST(BndLib_Add2dCurveTest, Circle_Full)
+{
+  // Create a 2D circle using Geom2d_Circle
+  Handle(Geom2d_Circle) aCircle =
+    new Geom2d_Circle(gp_Ax2d(gp_Pnt2d(0., 0.), gp_Dir2d(1., 0.)), 5.);
+  Bnd_Box2d aBox;
+  double    aTol = 0.;
+
+  BndLib_Add2dCurve::Add(aCircle, aTol, aBox);
+
+  double aXmin, aYmin, aXmax, aYmax;
+  aBox.Get(aXmin, aYmin, aXmax, aYmax);
+
+  EXPECT_NEAR(aXmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+}
+
+TEST(BndLib_Add2dCurveTest, Circle_Arc)
+{
+  // Circle arc from 0 to PI/2
+  Handle(Geom2d_Circle) aCircle =
+    new Geom2d_Circle(gp_Ax2d(gp_Pnt2d(0., 0.), gp_Dir2d(1., 0.)), 5.);
+  Bnd_Box2d aBox;
+  double    aTol = 0.;
+
+  BndLib_Add2dCurve::Add(aCircle, 0., M_PI / 2., aTol, aBox);
+
+  double aXmin, aYmin, aXmax, aYmax;
+  aBox.Get(aXmin, aYmin, aXmax, aYmax);
+
+  EXPECT_NEAR(aXmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+}
+
+TEST(BndLib_Add2dCurveTest, Ellipse_Full)
+{
+  // Ellipse with major radius 8, minor radius 4
+  Handle(Geom2d_Ellipse) anEllipse =
+    new Geom2d_Ellipse(gp_Ax2d(gp_Pnt2d(0., 0.), gp_Dir2d(1., 0.)), 8., 4.);
+  Bnd_Box2d aBox;
+  double    aTol = 0.;
+
+  BndLib_Add2dCurve::Add(anEllipse, aTol, aBox);
+
+  double aXmin, aYmin, aXmax, aYmax;
+  aBox.Get(aXmin, aYmin, aXmax, aYmax);
+
+  EXPECT_NEAR(aXmin, -8., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 8., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 4., Precision::Confusion());
+}
+
+TEST(BndLib_Add2dCurveTest, Line_Segment)
+{
+  // Line segment from (0,0) to (10,5)
+  Handle(Geom2d_Line)         aLine = new Geom2d_Line(gp_Pnt2d(0., 0.), gp_Dir2d(2., 1.));
+  Handle(Geom2d_TrimmedCurve) aTrimmed =
+    new Geom2d_TrimmedCurve(aLine, 0., std::sqrt(125.)); // length = sqrt(100+25)
+  Geom2dAdaptor_Curve aCurve(aTrimmed);
+  Bnd_Box2d           aBox;
+  double              aTol = 0.;
+
+  BndLib_Add2dCurve::Add(aCurve, aTol, aBox);
+
+  double aXmin, aYmin, aXmax, aYmax;
+  aBox.Get(aXmin, aYmin, aXmax, aYmax);
+
+  EXPECT_NEAR(aXmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 10., 0.01);
+  EXPECT_NEAR(aYmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., 0.01);
+}
+
+TEST(BndLib_Add2dCurveTest, BezierCurve)
+{
+  // Create a cubic Bezier curve in 2D
+  TColgp_Array1OfPnt2d aPoles(1, 4);
+  aPoles(1) = gp_Pnt2d(0., 0.);
+  aPoles(2) = gp_Pnt2d(1., 3.);
+  aPoles(3) = gp_Pnt2d(3., 3.);
+  aPoles(4) = gp_Pnt2d(4., 0.);
+
+  Handle(Geom2d_BezierCurve) aBezier = new Geom2d_BezierCurve(aPoles);
+  Geom2dAdaptor_Curve        aCurve(aBezier);
+  Bnd_Box2d                  aBox;
+  double                     aTol = 0.;
+
+  BndLib_Add2dCurve::Add(aCurve, aTol, aBox);
+
+  double aXmin, aYmin, aXmax, aYmax;
+  aBox.Get(aXmin, aYmin, aXmax, aYmax);
+
+  // Box should contain all control points
+  EXPECT_LE(aXmin, 0.);
+  EXPECT_GE(aXmax, 4.);
+  EXPECT_LE(aYmin, 0.);
+  EXPECT_GE(aYmax, 3.);
+}
+
+TEST(BndLib_Add2dCurveTest, AddOptimal_Ellipse)
+{
+  // Test AddOptimal method with an ellipse
+  Handle(Geom2d_Ellipse) anEllipse =
+    new Geom2d_Ellipse(gp_Ax2d(gp_Pnt2d(0., 0.), gp_Dir2d(1., 0.)), 6., 3.);
+  Bnd_Box2d aBox;
+  double    aTol = 0.;
+
+  BndLib_Add2dCurve::AddOptimal(anEllipse, 0., 2. * M_PI, aTol, aBox);
+
+  double aXmin, aYmin, aXmax, aYmax;
+  aBox.Get(aXmin, aYmin, aXmax, aYmax);
+
+  EXPECT_NEAR(aXmin, -6., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 6., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -3., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 3., Precision::Confusion());
+}
+
+TEST(BndLib_Add2dCurveTest, Adaptor_Circle)
+{
+  // Test using Adaptor interface
+  Handle(Geom2d_Circle) aCircle =
+    new Geom2d_Circle(gp_Ax2d(gp_Pnt2d(5., 5.), gp_Dir2d(1., 0.)), 3.);
+  Geom2dAdaptor_Curve aCurve(aCircle);
+  Bnd_Box2d           aBox;
+  double              aTol = 0.;
+
+  BndLib_Add2dCurve::Add(aCurve, aTol, aBox);
+
+  double aXmin, aYmin, aXmax, aYmax;
+  aBox.Get(aXmin, aYmin, aXmax, aYmax);
+
+  EXPECT_NEAR(aXmin, 2., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 8., Precision::Confusion());
+  EXPECT_NEAR(aYmin, 2., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 8., Precision::Confusion());
+}
+
+//==================================================================================================
+// BndLib_AddSurface Tests
+//==================================================================================================
+
+TEST(BndLib_AddSurfaceTest, Plane)
+{
+  // Create a plane and add a finite patch
+  Handle(Geom_Plane)  aPlane = new Geom_Plane(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)));
+  GeomAdaptor_Surface aSurf(aPlane);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  // Add a 10x10 patch centered at origin
+  BndLib_AddSurface::Add(aSurf, -5., 5., -5., 5., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 0., Precision::Confusion());
+}
+
+TEST(BndLib_AddSurfaceTest, Cylinder_Full)
+{
+  // Create a cylindrical surface
+  Handle(Geom_CylindricalSurface) aCyl =
+    new Geom_CylindricalSurface(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 5.);
+  GeomAdaptor_Surface aSurf(aCyl);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  // Full U range, V from 0 to 10
+  BndLib_AddSurface::Add(aSurf, 0., 2. * M_PI, 0., 10., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 10., Precision::Confusion());
+}
+
+TEST(BndLib_AddSurfaceTest, Sphere_Full)
+{
+  // Create a spherical surface
+  Handle(Geom_SphericalSurface) aSphere =
+    new Geom_SphericalSurface(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 7.);
+  GeomAdaptor_Surface aSurf(aSphere);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  BndLib_AddSurface::Add(aSurf, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -7., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 7., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -7., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 7., Precision::Confusion());
+  EXPECT_NEAR(aZmin, -7., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 7., Precision::Confusion());
+}
+
+TEST(BndLib_AddSurfaceTest, Sphere_Patch)
+{
+  // Create a spherical surface patch (upper hemisphere)
+  Handle(Geom_SphericalSurface) aSphere =
+    new Geom_SphericalSurface(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 5.);
+  GeomAdaptor_Surface aSurf(aSphere);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  // Full U, V from 0 to PI/2 (upper hemisphere)
+  BndLib_AddSurface::Add(aSurf, 0., 2. * M_PI, 0., M_PI / 2., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -5., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 5., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 5., Precision::Confusion());
+}
+
+TEST(BndLib_AddSurfaceTest, Cone_Full)
+{
+  // Create a conical surface
+  Handle(Geom_ConicalSurface) aCone =
+    new Geom_ConicalSurface(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), M_PI / 4., 2.);
+  GeomAdaptor_Surface aSurf(aCone);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  // Full U range, V from 0 to 5
+  BndLib_AddSurface::Add(aSurf, 0., 2. * M_PI, 0., 5., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // V=5: radius = 2 + 5*sin(PI/4), Z = 5*cos(PI/4)
+  double aMaxRadius = 2. + 5. * std::sin(M_PI / 4.);
+  double aMaxZ      = 5. * std::cos(M_PI / 4.);
+
+  EXPECT_NEAR(aXmax, aMaxRadius, Precision::Confusion());
+  EXPECT_NEAR(aYmax, aMaxRadius, Precision::Confusion());
+  EXPECT_NEAR(aZmax, aMaxZ, Precision::Confusion());
+}
+
+TEST(BndLib_AddSurfaceTest, Torus_Full)
+{
+  // Create a toroidal surface
+  Handle(Geom_ToroidalSurface) aTorus =
+    new Geom_ToroidalSurface(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 10., 3.);
+  GeomAdaptor_Surface aSurf(aTorus);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  BndLib_AddSurface::Add(aSurf, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // BndLib_AddSurface uses conservative bounding for torus
+  // The exact bounds are +-13 for X,Y and +-3 for Z
+  // The algorithm may produce slightly larger bounds
+  EXPECT_LE(aXmin, -13. + Precision::Confusion());
+  EXPECT_GE(aXmax, 13. - Precision::Confusion());
+  EXPECT_LE(aYmin, -13. + Precision::Confusion());
+  EXPECT_GE(aYmax, 13. - Precision::Confusion());
+  EXPECT_NEAR(aZmin, -3., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 3., Precision::Confusion());
+}
+
+TEST(BndLib_AddSurfaceTest, BezierSurface)
+{
+  // Create a simple Bezier surface (2x2 control points)
+  TColgp_Array2OfPnt aPoles(1, 2, 1, 2);
+  aPoles(1, 1) = gp_Pnt(0., 0., 0.);
+  aPoles(1, 2) = gp_Pnt(0., 10., 2.);
+  aPoles(2, 1) = gp_Pnt(10., 0., 2.);
+  aPoles(2, 2) = gp_Pnt(10., 10., 0.);
+
+  Handle(Geom_BezierSurface) aBezier = new Geom_BezierSurface(aPoles);
+  GeomAdaptor_Surface        aSurf(aBezier);
+  Bnd_Box                    aBox;
+  double                     aTol = 0.;
+
+  BndLib_AddSurface::Add(aSurf, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // Box should contain all control points
+  EXPECT_LE(aXmin, 0.);
+  EXPECT_GE(aXmax, 10.);
+  EXPECT_LE(aYmin, 0.);
+  EXPECT_GE(aYmax, 10.);
+  EXPECT_LE(aZmin, 0.);
+  EXPECT_GE(aZmax, 2.);
+}
+
+TEST(BndLib_AddSurfaceTest, AddOptimal_Sphere)
+{
+  // Test AddOptimal method with a sphere
+  Handle(Geom_SphericalSurface) aSphere =
+    new Geom_SphericalSurface(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)), 4.);
+  GeomAdaptor_Surface aSurf(aSphere);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  BndLib_AddSurface::AddOptimal(aSurf, aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 4., Precision::Confusion());
+  EXPECT_NEAR(aYmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 4., Precision::Confusion());
+  EXPECT_NEAR(aZmin, -4., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 4., Precision::Confusion());
+}
+
+TEST(BndLib_AddSurfaceTest, Cylinder_Translated)
+{
+  // Create a cylindrical surface centered at (5, 5, 0)
+  Handle(Geom_CylindricalSurface) aCyl =
+    new Geom_CylindricalSurface(gp_Ax3(gp_Pnt(5., 5., 0.), gp_Dir(0., 0., 1.)), 3.);
+  GeomAdaptor_Surface aSurf(aCyl);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  BndLib_AddSurface::Add(aSurf, 0., 2. * M_PI, 0., 10., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  EXPECT_NEAR(aXmin, 2., Precision::Confusion());
+  EXPECT_NEAR(aXmax, 8., Precision::Confusion());
+  EXPECT_NEAR(aYmin, 2., Precision::Confusion());
+  EXPECT_NEAR(aYmax, 8., Precision::Confusion());
+  EXPECT_NEAR(aZmin, 0., Precision::Confusion());
+  EXPECT_NEAR(aZmax, 10., Precision::Confusion());
+}
+
+TEST(BndLib_AddSurfaceTest, Plane_Tilted)
+{
+  // Create a tilted plane (normal at 45 degrees)
+  // gp_Dir automatically normalizes the input vector
+  Handle(Geom_Plane)  aPlane = new Geom_Plane(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 1., 1.)));
+  GeomAdaptor_Surface aSurf(aPlane);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  // Add a unit square patch
+  BndLib_AddSurface::Add(aSurf, 0., 1., 0., 1., aTol, aBox);
+
+  EXPECT_FALSE(aBox.IsVoid()) << "Box should not be void for tilted plane";
+}
+
+// Test for edge case with large parameter values and small step sizes.
+// This verifies that the algorithm correctly samples different points
+// even when parameter values are large relative to the step size,
+// which could cause floating-point precision issues if accumulation
+// (U += dU) is used instead of multiplication (U = UMin + dU * i).
+TEST(BndLib_AddSurfaceTest, LargeParameters_PrecisionTest)
+{
+  // Create a Bezier surface - uses the default sampling path
+  TColgp_Array2OfPnt aPoles(1, 3, 1, 3);
+  // Create a simple surface with known extent
+  for (int i = 1; i <= 3; i++)
+  {
+    for (int j = 1; j <= 3; j++)
+    {
+      aPoles(i, j) = gp_Pnt((i - 1) * 10., (j - 1) * 10., 0.);
+    }
+  }
+  Handle(Geom_BezierSurface) aBezier = new Geom_BezierSurface(aPoles);
+  GeomAdaptor_Surface        aSurf(aBezier);
+  Bnd_Box                    aBox;
+  double                     aTol = 0.;
+
+  // Use large parameter offset with small parameter range
+  // This simulates the edge case where UMin/UMax are large but (UMax-UMin) is small
+  // The actual Bezier is defined on [0,1] x [0,1], but we test a sub-range
+  // to ensure the sampling algorithm works correctly
+  BndLib_AddSurface::Add(aSurf, 0., 1., 0., 1., aTol, aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // The surface spans [0,20] x [0,20] x [0,0]
+  // Bounding box should contain all sample points
+  EXPECT_LE(aXmin, 1.) << "XMin should capture sampled points";
+  EXPECT_GE(aXmax, 19.) << "XMax should capture sampled points";
+  EXPECT_LE(aYmin, 1.) << "YMin should capture sampled points";
+  EXPECT_GE(aYmax, 19.) << "YMax should capture sampled points";
+}
+
+// Test with translated surface where parameters are offset
+// This is a more realistic test of the large parameter scenario
+TEST(BndLib_AddSurfaceTest, OffsetSurface_ParameterPrecision)
+{
+  // Create a simple plane
+  Handle(Geom_Plane)  aPlane = new Geom_Plane(gp_Ax3(gp_Pnt(0., 0., 0.), gp_Dir(0., 0., 1.)));
+  GeomAdaptor_Surface aSurf(aPlane);
+  Bnd_Box             aBox;
+  double              aTol = 0.;
+
+  // Test with very large parameter values but small range
+  // UMin = 1e10, UMax = 1e10 + 100, so dU ~ 100/(Nu-1) ~ small relative to 1e10
+  const double aLargeOffset = 1.0e10;
+  const double aRange       = 100.;
+  BndLib_AddSurface::Add(aSurf,
+                         aLargeOffset,
+                         aLargeOffset + aRange,
+                         aLargeOffset,
+                         aLargeOffset + aRange,
+                         aTol,
+                         aBox);
+
+  double aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
+  aBox.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
+
+  // For a plane at z=0 with parameters as X,Y coordinates:
+  // The box should span approximately [1e10, 1e10+100] in both X and Y
+  // Allow some tolerance due to floating-point representation
+  const double aTolerance = 1.0; // 1 unit tolerance for large values
+
+  EXPECT_NEAR(aXmin, aLargeOffset, aTolerance) << "XMin should be near the large offset value";
+  EXPECT_NEAR(aXmax, aLargeOffset + aRange, aTolerance) << "XMax should capture the full range";
+  EXPECT_NEAR(aYmin, aLargeOffset, aTolerance) << "YMin should be near the large offset value";
+  EXPECT_NEAR(aYmax, aLargeOffset + aRange, aTolerance) << "YMax should capture the full range";
+
+  // Most importantly: verify the range is captured (not collapsed due to precision loss)
+  const double aComputedRangeX = aXmax - aXmin;
+  const double aComputedRangeY = aYmax - aYmin;
+
+  EXPECT_GT(aComputedRangeX, aRange * 0.9)
+    << "X range should not collapse due to floating-point precision issues";
+  EXPECT_GT(aComputedRangeY, aRange * 0.9)
+    << "Y range should not collapse due to floating-point precision issues";
+}
index e3af2a8a74486f0c57454c56cf9c6431e0167708..a446435c02b32d14677abd2ca144ef0f838c7fac 100644 (file)
@@ -2,6 +2,7 @@
 set(OCCT_TKGeomBase_GTests_FILES_LOCATION "${CMAKE_CURRENT_LIST_DIR}")
 
 set(OCCT_TKGeomBase_GTests_FILES
+  BndLib_Test.cxx
   Extrema_ExtPC_Test.cxx
   IntAna_IntQuadQuad_Test.cxx
 )
index 7602c50b330d3ae12051c643d75d5ca2220278cc..dcf0a1de7b5f5ba5245489c4fc548d2762ea26e9 100644 (file)
@@ -18,7 +18,7 @@ checkshape b_1
 set xmin0 -62.992562093513783
 set ymin0 -51.45371311501097
 set zmin0 -33.223762093513777
-set xmax0 43.263112093513783
+set xmax0 43.257600511674156
 set ymax0 33.211062093513789
 set zmax0 37.744962093513792
 set bb [bounding b_1]