mirror of
https://git.dev.opencascade.org/repos/occt.git
synced 2025-04-03 17:56:21 +03:00
0031995: Modeling Data - Bounding box wrong for face
BndLib/BndLib_AddSurface.cxx : algorithm of building bounding box for BSpline surfaces is improved Test cases are modified according to current state of algorithm Test case bug31995 added
This commit is contained in:
parent
38b336df80
commit
2fda738623
@ -60,6 +60,17 @@ static Standard_Real AdjustExtr(const Adaptor3d_Surface& S,
|
|||||||
const Standard_Real Tol,
|
const Standard_Real Tol,
|
||||||
const Standard_Boolean IsMin);
|
const Standard_Boolean IsMin);
|
||||||
|
|
||||||
|
|
||||||
|
static void ComputePolesIndexes(const TColStd_Array1OfReal &theKnots,
|
||||||
|
const TColStd_Array1OfInteger &theMults,
|
||||||
|
const Standard_Integer theDegree,
|
||||||
|
const Standard_Real theMin,
|
||||||
|
const Standard_Real theMax,
|
||||||
|
const Standard_Integer theMaxPoleIdx,
|
||||||
|
const Standard_Boolean theIsPeriodic,
|
||||||
|
Standard_Integer &theOutMinIdx,
|
||||||
|
Standard_Integer &theOutMaxIdx);
|
||||||
|
|
||||||
//=======================================================================
|
//=======================================================================
|
||||||
//function : Add
|
//function : Add
|
||||||
//purpose :
|
//purpose :
|
||||||
@ -212,35 +223,31 @@ static void TreatInfinitePlane(const gp_Pln &aPlane,
|
|||||||
// theShiftCoeff - shift between flatknots array and poles array.
|
// theShiftCoeff - shift between flatknots array and poles array.
|
||||||
// This vaule should be equal to 1 in case of non periodic BSpline,
|
// This vaule should be equal to 1 in case of non periodic BSpline,
|
||||||
// and (degree + 1) - mults(the lowest index).
|
// and (degree + 1) - mults(the lowest index).
|
||||||
void ComputePolesIndexes(const TColStd_Array1OfReal &theFlatKnots,
|
|
||||||
const Standard_Integer theDegree,
|
void ComputePolesIndexes(const TColStd_Array1OfReal &theKnots,
|
||||||
const Standard_Real theMin,
|
const TColStd_Array1OfInteger &theMults,
|
||||||
const Standard_Real theMax,
|
const Standard_Integer theDegree,
|
||||||
const Standard_Integer theMinIdx,
|
const Standard_Real theMin,
|
||||||
const Standard_Integer theMaxIdx,
|
const Standard_Real theMax,
|
||||||
const Standard_Integer theShiftCoeff,
|
const Standard_Integer theMaxPoleIdx,
|
||||||
Standard_Integer &theOutMinIdx,
|
const Standard_Boolean theIsPeriodic,
|
||||||
Standard_Integer &theOutMaxIdx)
|
Standard_Integer &theOutMinIdx,
|
||||||
|
Standard_Integer &theOutMaxIdx)
|
||||||
{
|
{
|
||||||
// Set initial values for the result indexes to handle situation when requested parameter space
|
BSplCLib::Hunt(theKnots, theMin, theOutMinIdx);
|
||||||
// is slightly greater than B-spline parameter space.
|
theOutMinIdx = Max(theOutMinIdx, theKnots.Lower());
|
||||||
theOutMinIdx = theFlatKnots.Lower();
|
|
||||||
theOutMaxIdx = theFlatKnots.Upper();
|
|
||||||
|
|
||||||
// Compute first and last used flat knots.
|
BSplCLib::Hunt(theKnots, theMax, theOutMaxIdx);
|
||||||
for(Standard_Integer aKnotIdx = theFlatKnots.Lower();
|
theOutMaxIdx++;
|
||||||
aKnotIdx < theFlatKnots.Upper();
|
theOutMaxIdx = Min(theOutMaxIdx, theKnots.Upper());
|
||||||
aKnotIdx++)
|
Standard_Integer mult = theMults(theOutMaxIdx);
|
||||||
{
|
|
||||||
if (theFlatKnots(aKnotIdx) <= theMin)
|
|
||||||
theOutMinIdx = aKnotIdx;
|
|
||||||
|
|
||||||
if (theFlatKnots(theFlatKnots.Upper() - aKnotIdx + theFlatKnots.Lower()) >= theMax)
|
theOutMinIdx = BSplCLib::PoleIndex(theDegree, theOutMinIdx, theIsPeriodic, theMults) + 1;
|
||||||
theOutMaxIdx = theFlatKnots.Upper() - aKnotIdx + theFlatKnots.Lower();
|
theOutMinIdx = Max(theOutMinIdx, 1);
|
||||||
}
|
theOutMaxIdx = BSplCLib::PoleIndex(theDegree, theOutMaxIdx, theIsPeriodic, theMults) + 1;
|
||||||
|
theOutMaxIdx += theDegree - mult;
|
||||||
theOutMinIdx = Max(theOutMinIdx - 2 * theDegree + 2 - theShiftCoeff, theMinIdx);
|
if (!theIsPeriodic)
|
||||||
theOutMaxIdx = Min(theOutMaxIdx - 2 + theDegree + 1 - theShiftCoeff, theMaxIdx);
|
theOutMaxIdx = Min(theOutMaxIdx, theMaxPoleIdx);
|
||||||
}
|
}
|
||||||
|
|
||||||
// Modified by skv - Fri Aug 27 12:29:04 2004 OCC6503 End
|
// Modified by skv - Fri Aug 27 12:29:04 2004 OCC6503 End
|
||||||
@ -336,6 +343,7 @@ void BndLib_AddSurface::Add(const Adaptor3d_Surface& S,
|
|||||||
// Borders of underlying geometry.
|
// Borders of underlying geometry.
|
||||||
Standard_Real anUMinParam = UMin, anUMaxParam = UMax,// BSpline case.
|
Standard_Real anUMinParam = UMin, anUMaxParam = UMax,// BSpline case.
|
||||||
aVMinParam = VMin, aVMaxParam = VMax;
|
aVMinParam = VMin, aVMaxParam = VMax;
|
||||||
|
Handle(Geom_BSplineSurface) aBS;
|
||||||
if (Type == GeomAbs_BezierSurface)
|
if (Type == GeomAbs_BezierSurface)
|
||||||
{
|
{
|
||||||
// Bezier surface:
|
// Bezier surface:
|
||||||
@ -358,8 +366,8 @@ void BndLib_AddSurface::Add(const Adaptor3d_Surface& S,
|
|||||||
// use convex hull algorithm,
|
// use convex hull algorithm,
|
||||||
// if Umin, VMin, Umax, Vmax lies outside then:
|
// if Umin, VMin, Umax, Vmax lies outside then:
|
||||||
// use grid algorithm on analytic continuation (default case).
|
// use grid algorithm on analytic continuation (default case).
|
||||||
S.BSpline()->Bounds(anUMinParam, anUMaxParam, aVMinParam, aVMaxParam);
|
aBS = S.BSpline();
|
||||||
|
aBS->Bounds(anUMinParam, anUMaxParam, aVMinParam, aVMaxParam);
|
||||||
if ( (UMin - anUMinParam) < -PTol ||
|
if ( (UMin - anUMinParam) < -PTol ||
|
||||||
(VMin - aVMinParam) < -PTol ||
|
(VMin - aVMinParam) < -PTol ||
|
||||||
(UMax - anUMaxParam) > PTol ||
|
(UMax - anUMaxParam) > PTol ||
|
||||||
@ -372,94 +380,87 @@ void BndLib_AddSurface::Add(const Adaptor3d_Surface& S,
|
|||||||
|
|
||||||
if (isUseConvexHullAlgorithm)
|
if (isUseConvexHullAlgorithm)
|
||||||
{
|
{
|
||||||
TColgp_Array2OfPnt Tp(1,S.NbUPoles(),1,S.NbVPoles());
|
Standard_Integer aNbUPoles = S.NbUPoles(), aNbVPoles = S.NbVPoles();
|
||||||
Standard_Integer UMinIdx = 0, UMaxIdx = 0;
|
TColgp_Array2OfPnt Tp(1, aNbUPoles, 1, aNbVPoles);
|
||||||
Standard_Integer VMinIdx = 0, VMaxIdx = 0;
|
Standard_Integer UMinIdx = 0, UMaxIdx = 0;
|
||||||
if (Type == GeomAbs_BezierSurface)
|
Standard_Integer VMinIdx = 0, VMaxIdx = 0;
|
||||||
|
Standard_Boolean isUPeriodic = S.IsUPeriodic(), isVPeriodic = S.IsVPeriodic();
|
||||||
|
if (Type == GeomAbs_BezierSurface)
|
||||||
|
{
|
||||||
|
S.Bezier()->Poles(Tp);
|
||||||
|
UMinIdx = 1; UMaxIdx = aNbUPoles;
|
||||||
|
VMinIdx = 1; VMaxIdx = aNbVPoles;
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
aBS->Poles(Tp);
|
||||||
|
|
||||||
|
UMinIdx = 1;
|
||||||
|
UMaxIdx = aNbUPoles;
|
||||||
|
VMinIdx = 1;
|
||||||
|
VMaxIdx = aNbVPoles;
|
||||||
|
|
||||||
|
if (UMin > anUMinParam ||
|
||||||
|
UMax < anUMaxParam)
|
||||||
{
|
{
|
||||||
S.Bezier()->Poles(Tp);
|
TColStd_Array1OfInteger aMults(1, aBS->NbUKnots());
|
||||||
|
TColStd_Array1OfReal aKnots(1, aBS->NbUKnots());
|
||||||
|
aBS->UKnots(aKnots);
|
||||||
|
aBS->UMultiplicities(aMults);
|
||||||
|
|
||||||
UMinIdx = Tp.LowerRow();
|
ComputePolesIndexes(aKnots,
|
||||||
UMaxIdx = Tp.UpperRow();
|
aMults,
|
||||||
VMinIdx = Tp.LowerCol();
|
aBS->UDegree(),
|
||||||
VMaxIdx = Tp.UpperCol();
|
UMin, UMax,
|
||||||
}
|
aNbUPoles,
|
||||||
else
|
isUPeriodic,
|
||||||
{
|
UMinIdx, UMaxIdx); // the Output indexes
|
||||||
S.BSpline()->Poles(Tp);
|
|
||||||
|
|
||||||
UMinIdx = Tp.LowerRow();
|
|
||||||
UMaxIdx = Tp.UpperRow();
|
|
||||||
VMinIdx = Tp.LowerCol();
|
|
||||||
VMaxIdx = Tp.UpperCol();
|
|
||||||
|
|
||||||
if (UMin > anUMinParam ||
|
|
||||||
UMax < anUMaxParam)
|
|
||||||
{
|
|
||||||
Standard_Integer anUFlatKnotsCount = S.BSpline()->NbUPoles() + S.BSpline()->UDegree() + 1;
|
|
||||||
Standard_Integer aShift = 1;
|
|
||||||
|
|
||||||
if (S.BSpline()->IsUPeriodic())
|
|
||||||
{
|
|
||||||
TColStd_Array1OfInteger aMults(1, S.BSpline()->NbUKnots());
|
|
||||||
S.BSpline()->UMultiplicities(aMults);
|
|
||||||
anUFlatKnotsCount = BSplCLib::KnotSequenceLength(aMults, S.BSpline()->UDegree(), Standard_True);
|
|
||||||
|
|
||||||
aShift = S.BSpline()->UDegree() + 1 - S.BSpline()->UMultiplicity(1);
|
|
||||||
}
|
|
||||||
|
|
||||||
TColStd_Array1OfReal anUFlatKnots(1, anUFlatKnotsCount);
|
|
||||||
S.BSpline()->UKnotSequence(anUFlatKnots);
|
|
||||||
|
|
||||||
ComputePolesIndexes(anUFlatKnots,
|
|
||||||
S.BSpline()->UDegree(),
|
|
||||||
UMin, UMax,
|
|
||||||
UMinIdx, UMaxIdx, // Min and Max Indexes
|
|
||||||
aShift,
|
|
||||||
UMinIdx, UMaxIdx); // the Output indexes
|
|
||||||
}
|
|
||||||
|
|
||||||
if (VMin > aVMinParam ||
|
|
||||||
VMax < aVMaxParam)
|
|
||||||
{
|
|
||||||
Standard_Integer anVFlatKnotsCount = S.BSpline()->NbVPoles() + S.BSpline()->VDegree() + 1;
|
|
||||||
Standard_Integer aShift = 1;
|
|
||||||
|
|
||||||
if (S.BSpline()->IsVPeriodic())
|
|
||||||
{
|
|
||||||
TColStd_Array1OfInteger aMults(1, S.BSpline()->NbVKnots());
|
|
||||||
S.BSpline()->VMultiplicities(aMults);
|
|
||||||
anVFlatKnotsCount = BSplCLib::KnotSequenceLength(aMults, S.BSpline()->VDegree(), Standard_True);
|
|
||||||
|
|
||||||
aShift = S.BSpline()->VDegree() + 1 - S.BSpline()->VMultiplicity(1);
|
|
||||||
}
|
|
||||||
|
|
||||||
TColStd_Array1OfReal anVFlatKnots(1, anVFlatKnotsCount);
|
|
||||||
S.BSpline()->VKnotSequence(anVFlatKnots);
|
|
||||||
|
|
||||||
ComputePolesIndexes(anVFlatKnots,
|
|
||||||
S.BSpline()->VDegree(),
|
|
||||||
VMin, VMax,
|
|
||||||
VMinIdx, VMaxIdx, // Min and Max Indexes
|
|
||||||
aShift,
|
|
||||||
VMinIdx, VMaxIdx); // the Output indexes
|
|
||||||
}
|
|
||||||
|
|
||||||
}
|
}
|
||||||
|
|
||||||
// Use poles to build convex hull.
|
if (VMin > aVMinParam ||
|
||||||
for (Standard_Integer i = UMinIdx; i <= UMaxIdx; i++)
|
VMax < aVMaxParam)
|
||||||
{
|
{
|
||||||
for (Standard_Integer j = VMinIdx; j <= VMaxIdx; j++)
|
TColStd_Array1OfInteger aMults(1, aBS->NbVKnots());
|
||||||
{
|
TColStd_Array1OfReal aKnots(1, aBS->NbVKnots());
|
||||||
B.Add(Tp(i,j));
|
aBS->VKnots(aKnots);
|
||||||
}
|
aBS->VMultiplicities(aMults);
|
||||||
|
|
||||||
|
ComputePolesIndexes(aKnots,
|
||||||
|
aMults,
|
||||||
|
aBS->VDegree(),
|
||||||
|
VMin, VMax,
|
||||||
|
aNbVPoles,
|
||||||
|
isVPeriodic,
|
||||||
|
VMinIdx, VMaxIdx); // the Output indexes
|
||||||
}
|
}
|
||||||
|
|
||||||
B.Enlarge(Tol);
|
}
|
||||||
break;
|
|
||||||
|
// Use poles to build convex hull.
|
||||||
|
Standard_Integer ip, jp;
|
||||||
|
for (Standard_Integer i = UMinIdx; i <= UMaxIdx; i++)
|
||||||
|
{
|
||||||
|
ip = i;
|
||||||
|
if (isUPeriodic && ip > aNbUPoles)
|
||||||
|
{
|
||||||
|
ip = ip - aNbUPoles;
|
||||||
|
}
|
||||||
|
for (Standard_Integer j = VMinIdx; j <= VMaxIdx; j++)
|
||||||
|
{
|
||||||
|
jp = j;
|
||||||
|
if (isVPeriodic && jp > aNbVPoles)
|
||||||
|
{
|
||||||
|
jp = jp - aNbVPoles;
|
||||||
|
}
|
||||||
|
B.Add(Tp(ip, jp));
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
B.Enlarge(Tol);
|
||||||
|
break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
Standard_FALLTHROUGH
|
Standard_FALLTHROUGH
|
||||||
default:
|
default:
|
||||||
{
|
{
|
||||||
@ -961,7 +962,7 @@ Standard_Integer NbVSamples(const Adaptor3d_Surface& S,
|
|||||||
{
|
{
|
||||||
const Handle(Geom_BSplineSurface)& BS = S.BSpline();
|
const Handle(Geom_BSplineSurface)& BS = S.BSpline();
|
||||||
N = 2*(BS->VDegree() + 1)*(BS->NbVKnots() - 1) ;
|
N = 2*(BS->VDegree() + 1)*(BS->NbVKnots() - 1) ;
|
||||||
Standard_Real umin, umax, vmin, vmax;
|
Standard_Real umin, umax, vmin, vmax;
|
||||||
BS->Bounds(umin, umax, vmin, vmax);
|
BS->Bounds(umin, umax, vmin, vmax);
|
||||||
Standard_Real dv = (Vmax - Vmin) / (vmax - vmin);
|
Standard_Real dv = (Vmax - Vmin) / (vmax - vmin);
|
||||||
if(dv < .9)
|
if(dv < .9)
|
||||||
|
@ -14,6 +14,6 @@ vdisplay result
|
|||||||
vsetdispmode result 1
|
vsetdispmode result 1
|
||||||
vfit
|
vfit
|
||||||
|
|
||||||
checktrinfo result -tri 2711 -nod 2611
|
checktrinfo result -tri 2722 -nod 2618
|
||||||
|
|
||||||
checkview -display result -2d -path ${imagedir}/${test_image}.png
|
checkview -display result -2d -path ${imagedir}/${test_image}.png
|
||||||
|
@ -12,6 +12,6 @@ vdisplay result
|
|||||||
vviewparams -scale 8.46292 -proj 0.653203 -0.644806 0.396926 -up -0.0109833 0.51609 0.856464 -at 347.559 1026.89 219.262 -eye 2080.75 -684.022 1272.45
|
vviewparams -scale 8.46292 -proj 0.653203 -0.644806 0.396926 -up -0.0109833 0.51609 0.856464 -at 347.559 1026.89 219.262 -eye 2080.75 -684.022 1272.45
|
||||||
|
|
||||||
tricheck result
|
tricheck result
|
||||||
checktrinfo result -tri 11826 -nod 7310 -defl 7.6167024939147652
|
checktrinfo result -tri 11800 -nod 7301 -defl 7.6167024939147652
|
||||||
|
|
||||||
checkview -screenshot -3d -path ${imagedir}/${test_image}.png
|
checkview -screenshot -3d -path ${imagedir}/${test_image}.png
|
||||||
|
@ -13,12 +13,12 @@ bounding result -save Xmin Ymin Zmin Xmax Ymax Zmax
|
|||||||
set tol_abs 1.0e-4
|
set tol_abs 1.0e-4
|
||||||
set tol_rel 1.0e-4
|
set tol_rel 1.0e-4
|
||||||
|
|
||||||
checkreal "Xmin" [dval Xmin] 102.04999989999993 ${tol_abs} ${tol_rel}
|
checkreal "Xmin" [dval Xmin] 107.49999989999992 ${tol_abs} ${tol_rel}
|
||||||
checkreal "Ymin" [dval Ymin] -12.576503364721431 ${tol_abs} ${tol_rel}
|
checkreal "Ymin" [dval Ymin] -12.487307730633276 ${tol_abs} ${tol_rel}
|
||||||
checkreal "Zmin" [dval Zmin] -12.267407382031644 ${tol_abs} ${tol_rel}
|
checkreal "Zmin" [dval Zmin] -12.135364833115547 ${tol_abs} ${tol_rel}
|
||||||
checkreal "Xmax" [dval Xmax] 145.65000009999983 ${tol_abs} ${tol_rel}
|
checkreal "Xmax" [dval Xmax] 140.20000009999987 ${tol_abs} ${tol_rel}
|
||||||
checkreal "Ymax" [dval Ymax] 1.0883692081680807 ${tol_abs} ${tol_rel}
|
checkreal "Ymax" [dval Ymax] -0.8336837432205401 ${tol_abs} ${tol_rel}
|
||||||
checkreal "Zmax" [dval Zmax] 1.4146362604116396 ${tol_abs} ${tol_rel}
|
checkreal "Zmax" [dval Zmax] -0.5631589698265231 ${tol_abs} ${tol_rel}
|
||||||
|
|
||||||
smallview
|
smallview
|
||||||
fit
|
fit
|
||||||
|
23
tests/bugs/moddata_3/bug31995
Normal file
23
tests/bugs/moddata_3/bug31995
Normal file
@ -0,0 +1,23 @@
|
|||||||
|
puts "========="
|
||||||
|
puts "0031995: Modeling Data - Bounding box wrong for face"
|
||||||
|
puts "========="
|
||||||
|
puts ""
|
||||||
|
|
||||||
|
|
||||||
|
restore [locate_data_file bug31995.brep] result
|
||||||
|
|
||||||
|
bounding result -save Xmin Ymin Zmin Xmax Ymax Zmax
|
||||||
|
|
||||||
|
set tol_abs 1.0e-4
|
||||||
|
set tol_rel 1.0e-4
|
||||||
|
|
||||||
|
checkreal "Xmin" [dval Xmin] -0.4992512824357961 ${tol_abs} ${tol_rel}
|
||||||
|
checkreal "Ymin" [dval Ymin] 0.5646627848871493 ${tol_abs} ${tol_rel}
|
||||||
|
checkreal "Zmin" [dval Zmin] 0.5646627848871487 ${tol_abs} ${tol_rel}
|
||||||
|
checkreal "Xmax" [dval Xmax] 0.1535534905932737 ${tol_abs} ${tol_rel}
|
||||||
|
checkreal "Ymax" [dval Ymax] 1.4505295939775213 ${tol_abs} ${tol_rel}
|
||||||
|
checkreal "Zmax" [dval Zmax] 1.4505295939775207 ${tol_abs} ${tol_rel}
|
||||||
|
|
||||||
|
smallview
|
||||||
|
fit
|
||||||
|
checkview -screenshot -2d -path ${imagedir}/${test_image}.png
|
@ -1,4 +1,4 @@
|
|||||||
puts "TODO OCC30286 ALL: Error : The length of result shape is 5496.05, expected 5934.34"
|
puts "TODO OCC30286 ALL: Error : The length of result shape is 5499.57, expected 5934.34"
|
||||||
|
|
||||||
set viewname "vright"
|
set viewname "vright"
|
||||||
set length 5934.34
|
set length 5934.34
|
||||||
|
Loading…
x
Reference in New Issue
Block a user