1
0
mirror of https://git.dev.opencascade.org/repos/occt.git synced 2025-04-10 18:51:21 +03:00

0026699: Wrong section curves

1. Algorithm of Restriction line processing has been improved in IntTools_FaceFace.cxx file.
2. Algorithm of checking, if Restriction line and Walking line are coincided has been improved in IntPatch_ImpPrmIntersection.cxx file.
3. Algorithm of extending check if starting point of Walking line is a tangent point has been added.

Small correction of some test cases.
Creation of test case for issue #0026699.

Small correction of test case for issue CR26699
This commit is contained in:
nbv 2015-10-28 10:10:37 +03:00 committed by bugmaster
parent 5344638378
commit 71958f7db2
7 changed files with 157 additions and 86 deletions

View File

@ -1351,7 +1351,8 @@ void IntPatch_ImpPrmIntersection::Perform (const Handle(Adaptor3d_HSurface)& Sur
}// if (NbSegm) }// if (NbSegm)
// //
// on traite les restrictions de la surface implicite // on traite les restrictions de la surface implicite
for (Standard_Integer i=1; i<=slin.Length(); i++)
for (Standard_Integer i=1, aNbLin = slin.Length(); i<=aNbLin; i++)
{ {
Handle(IntPatch_Line)& aL = slin(i); Handle(IntPatch_Line)& aL = slin(i);
@ -1359,7 +1360,19 @@ void IntPatch_ImpPrmIntersection::Perform (const Handle(Adaptor3d_HSurface)& Sur
IntPatch_RstInt::PutVertexOnLine(aL,Surf1,D1,Surf2,Standard_True,TolTang); IntPatch_RstInt::PutVertexOnLine(aL,Surf1,D1,Surf2,Standard_True,TolTang);
else else
IntPatch_RstInt::PutVertexOnLine(aL,Surf2,D2,Surf1,Standard_False,TolTang); IntPatch_RstInt::PutVertexOnLine(aL,Surf2,D2,Surf1,Standard_False,TolTang);
if(aL->ArcType() == IntPatch_Walking)
{
const Handle(IntPatch_WLine) aWL = Handle(IntPatch_WLine)::DownCast(aL);
slin.Append(aWL);
slin.Remove(i);
i--;
aNbLin--;
} }
}
// Now slin is filled as follows: lower indices correspond to Restriction line,
// after (higher indices) - only Walking-line.
const Standard_Real aUPeriodOfSurf1 = Surf1->IsUPeriodic() ? Surf1->UPeriod() : 0.0, const Standard_Real aUPeriodOfSurf1 = Surf1->IsUPeriodic() ? Surf1->UPeriod() : 0.0,
aUPeriodOfSurf2 = Surf2->IsUPeriodic() ? Surf2->UPeriod() : 0.0, aUPeriodOfSurf2 = Surf2->IsUPeriodic() ? Surf2->UPeriod() : 0.0,
@ -1378,21 +1391,13 @@ void IntPatch_ImpPrmIntersection::Perform (const Handle(Adaptor3d_HSurface)& Sur
Handle(IntPatch_RLine) aRL1 = Handle(IntPatch_RLine)::DownCast(aL1); Handle(IntPatch_RLine) aRL1 = Handle(IntPatch_RLine)::DownCast(aL1);
Handle(IntPatch_RLine) aRL2 = Handle(IntPatch_RLine)::DownCast(aL2); Handle(IntPatch_RLine) aRL2 = Handle(IntPatch_RLine)::DownCast(aL2);
if(aRL1.IsNull() && aRL2.IsNull()) if(aRL1.IsNull())
{//If Walking-Walking {//If Walking-Walking
continue; continue;
} }
else if(aRL1.IsNull())
{// i-th line is not restriction,
// but j-th is restriction
slin.Append(aL1);
slin.SetValue(i, aL2);
slin.Remove(j);
j--;
continue;
}
//Here aL1 (i-th line) is Restriction-line and aL2 (j-th line) is not Restriction //Here aL1 (i-th line) is Restriction-line and aL2 (j-th line) is
//Restriction or Walking
if(aBRL.IsVoid()) if(aBRL.IsVoid())
{//Fill aBRL {//Fill aBRL

View File

@ -657,7 +657,7 @@ void IntTools_FaceFace::Perform(const TopoDS_Face& aF1,
// //
const Standard_Integer aNbLinIntersector = myIntersector.NbLines(); const Standard_Integer aNbLinIntersector = myIntersector.NbLines();
for (Standard_Integer i=1; i <= aNbLinIntersector; ++i) { for (Standard_Integer i=1; i <= aNbLinIntersector; ++i) {
MakeCurve(i, dom1, dom2); MakeCurve(i, dom1, dom2, TolArc);
} }
// //
ComputeTolReached3d(); ComputeTolReached3d();
@ -843,7 +843,8 @@ Standard_Real IntTools_FaceFace::ComputeTolerance()
//======================================================================= //=======================================================================
void IntTools_FaceFace::MakeCurve(const Standard_Integer Index, void IntTools_FaceFace::MakeCurve(const Standard_Integer Index,
const Handle(Adaptor3d_TopolTool)& dom1, const Handle(Adaptor3d_TopolTool)& dom1,
const Handle(Adaptor3d_TopolTool)& dom2) const Handle(Adaptor3d_TopolTool)& dom2,
const Standard_Real theToler)
{ {
Standard_Boolean bDone, rejectSurface, reApprox, bAvoidLineConstructor; Standard_Boolean bDone, rejectSurface, reApprox, bAvoidLineConstructor;
Standard_Boolean ok, bPCurvesOk; Standard_Boolean ok, bPCurvesOk;
@ -2050,16 +2051,27 @@ Standard_Real IntTools_FaceFace::ComputeTolerance()
GeomInt_IntSS:: GeomInt_IntSS::
TrimILineOnSurfBoundaries(aC2d1, aC2d2, aBox1, aBox2, anArrayOfParameters); TrimILineOnSurfBoundaries(aC2d1, aC2d2, aBox1, aBox2, anArrayOfParameters);
//Intersect with true boundaries. After that, enlarge bounding-boxes in order to
//correct definition, if point on curve is inscribed in the box.
aBox1.Enlarge(theToler);
aBox2.Enlarge(theToler);
const Standard_Integer aNbIntersSolutionsm1 = anArrayOfParameters.Length() - 1; const Standard_Integer aNbIntersSolutionsm1 = anArrayOfParameters.Length() - 1;
//Trim RLine found. //Trim RLine found.
for(Standard_Integer anInd = 0; anInd < aNbIntersSolutionsm1; anInd++) for(Standard_Integer anInd = 0; anInd < aNbIntersSolutionsm1; anInd++)
{ {
const Standard_Real aParF = anArrayOfParameters(anInd), Standard_Real &aParF = anArrayOfParameters(anInd),
aParL = anArrayOfParameters(anInd+1); &aParL = anArrayOfParameters(anInd+1);
if((aParL - aParF) <= Precision::PConfusion()) if((aParL - aParF) <= Precision::PConfusion())
{
//In order to more precise extending to the boundaries of source curves.
if(anInd < aNbIntersSolutionsm1-1)
aParL = aParF;
continue; continue;
}
const Standard_Real aPar = 0.5*(aParF + aParL); const Standard_Real aPar = 0.5*(aParF + aParL);
gp_Pnt2d aPt; gp_Pnt2d aPt;

View File

@ -117,7 +117,10 @@ public:
protected: protected:
Standard_EXPORT void MakeCurve (const Standard_Integer Index, const Handle(Adaptor3d_TopolTool)& D1, const Handle(Adaptor3d_TopolTool)& D2); Standard_EXPORT void MakeCurve (const Standard_Integer Index,
const Handle(Adaptor3d_TopolTool)& D1,
const Handle(Adaptor3d_TopolTool)& D2,
const Standard_Real theToler);
Standard_EXPORT void ComputeTolReached3d(); Standard_EXPORT void ComputeTolReached3d();

View File

@ -21,6 +21,54 @@ OSD_Chronometer Chronrsnld;
#include <NCollection_IncAllocator.hxx> #include <NCollection_IncAllocator.hxx>
#include <Precision.hxx> #include <Precision.hxx>
//==================================================================================
// function : IsTangentExtCheck
// purpose : Additional check if the point (theU, theV) in parametric surface
// is a tangent point.
// If that is TRUE then we can go along any (!) direction in order to
// obtain next intersection point. At present, this case is difficult
// for processing. Therefore, we will return an empty intersection line
// in this case.
//==================================================================================
static Standard_Boolean IsTangentExtCheck(TheIWFunction& theFunc,
const Standard_Real theU,
const Standard_Real theV,
const Standard_Real theStepU,
const Standard_Real theStepV,
const Standard_Real theUinf,
const Standard_Real theUsup,
const Standard_Real theVinf,
const Standard_Real theVsup)
{
//Factor 2.0 is chosen because we compare distance(s) between TWO faces
const Standard_Real aTol = 2.0*Precision::Confusion();
const Standard_Integer aNbItems = 4;
const Standard_Real aParU[aNbItems] = { Min(theU + theStepU, theUsup),
Max(theU - theStepU, theUinf),
theU, theU};
const Standard_Real aParV[aNbItems] = { theV, theV,
Min(theV + theStepV, theVsup),
Max(theV - theStepV, theVinf)};
math_Vector aX(1, 2), aVal(1, 1);
for(Standard_Integer i = 0; i < aNbItems; i++)
{
aX.Value(1) = aParU[i];
aX.Value(2) = aParV[i];
if(!theFunc.Value(aX, aVal))
continue;
if(aVal(1) > aTol)
return Standard_False;
}
return Standard_True;
}
IntWalk_IWalking::IntWalk_IWalking (const Standard_Real Epsilon, IntWalk_IWalking::IntWalk_IWalking (const Standard_Real Epsilon,
const Standard_Real Deflection, const Standard_Real Deflection,
const Standard_Real Increment ) : const Standard_Real Increment ) :
@ -85,7 +133,6 @@ void IntWalk_IWalking::Perform(const ThePOPIterator& Pnts1,
const Standard_Boolean Reversed) const Standard_Boolean Reversed)
{ {
Standard_Integer I; Standard_Integer I;
Standard_Boolean Rajout = Standard_False; Standard_Boolean Rajout = Standard_False;
Standard_Integer nbPnts1 = Pnts1.Length(); Standard_Integer nbPnts1 = Pnts1.Length();
@ -95,6 +142,23 @@ void IntWalk_IWalking::Perform(const ThePOPIterator& Pnts1,
Clear(); Clear();
reversed = Reversed; reversed = Reversed;
Um = ThePSurfaceTool::FirstUParameter(Caro);
Vm = ThePSurfaceTool::FirstVParameter(Caro);
UM = ThePSurfaceTool::LastUParameter(Caro);
VM = ThePSurfaceTool::LastVParameter(Caro);
if (UM < Um) {
Standard_Real utemp = UM;
UM = Um;
Um = utemp;
}
if (VM < Vm) {
Standard_Real vtemp = VM;
VM = Vm;
Vm = vtemp;
}
const Standard_Real aStepU = pas*(UM-Um), aStepV = pas*(VM-Vm);
// Loading of etat1 and etat2 as well as ustart and vstart. // Loading of etat1 and etat2 as well as ustart and vstart.
@ -132,30 +196,18 @@ void IntWalk_IWalking::Perform(const ThePOPIterator& Pnts1,
wd2.reserve (nbPnts2); wd2.reserve (nbPnts2);
for (I = 1; I <= nbPnts2; I++) { for (I = 1; I <= nbPnts2; I++) {
IntWalk_WalkingData aWD2; IntWalk_WalkingData aWD2;
aWD2.etat = 1;
const IntSurf_InteriorPoint& anIP = Pnts2.Value(I);
ThePointOfLoopTool::Value2d(anIP, aWD2.ustart, aWD2.vstart);
if (!IsTangentExtCheck(Func, aWD2.ustart, aWD2.vstart, aStepU, aStepV, Um, UM, Vm, VM))
aWD2.etat = 13; aWD2.etat = 13;
ThePointOfLoopTool::Value2d(Pnts2.Value(I), aWD2.ustart, aWD2.vstart);
wd2.push_back (aWD2); wd2.push_back (aWD2);
} }
tolerance(1) = ThePSurfaceTool::UResolution(Caro,Precision::Confusion()); tolerance(1) = ThePSurfaceTool::UResolution(Caro,Precision::Confusion());
tolerance(2) = ThePSurfaceTool::VResolution(Caro,Precision::Confusion()); tolerance(2) = ThePSurfaceTool::VResolution(Caro,Precision::Confusion());
Um = ThePSurfaceTool::FirstUParameter(Caro);
Vm = ThePSurfaceTool::FirstVParameter(Caro);
UM = ThePSurfaceTool::LastUParameter(Caro);
VM = ThePSurfaceTool::LastVParameter(Caro);
if (UM < Um) {
Standard_Real utemp = UM;
UM = Um;
Um = utemp;
}
if (VM < Vm) {
Standard_Real vtemp = VM;
VM = Vm;
Vm = vtemp;
}
Func.Set(Caro); Func.Set(Caro);
// calculation of all open lines // calculation of all open lines

View File

@ -23,10 +23,8 @@ set log [bopcurves f1 f2 -2d]
regexp {Tolerance Reached=+([-0-9.+eE]+)\n+([-0-9.+eE]+)} ${log} full Toler NbCurv regexp {Tolerance Reached=+([-0-9.+eE]+)\n+([-0-9.+eE]+)} ${log} full Toler NbCurv
#This value must be equal to the analogical value in bug25292_31 and bug25292_32 of "bugs modalg_5" testgrid.
set MaxTol 1.e-7 set MaxTol 1.e-7
#This value must be equal to the analogical value in bug25292_31 and bug25292_32 of "bugs modalg_5" testgrid.
set GoodNbCurv 1 set GoodNbCurv 1
if {${Toler} > ${MaxTol}} { if {${Toler} > ${MaxTol}} {
@ -44,30 +42,11 @@ mksurface s2 f2
erase s1 s2 erase s1 s2
for {set i 1} {$i <= ${NbCurv}} {incr i} { for {set i 1} {$i <= ${NbCurv}} {incr i} {
set log [dump c_$i] bounds c_$i U1 U2
set dumptrimres [regexp {Trimmed curve\nParameters : +([-0-9.+eE]+) +([-0-9.+eE]+)} ${log} full U1 U2]
if {${dumptrimres} == 0} { dump U1 U2
regexp {Degree +([-0-9.+eE]+), +([-0-9.+eE]+) Poles, +([-0-9.+eE]+)} ${log} full Degree Poles KnotsPoles
puts "Degree=${Degree}" if {[dval $U2-$U1] < 1.0e-20} {
puts "Poles=${Poles}"
puts "KnotsPoles=${KnotsPoles}"
puts ""
set Knot 1
set exp_string "Knots :\n\n +${Knot} : +(\[-0-9.+eE\]+) +(\[-0-9.+eE\]+)"
regexp ${exp_string} ${log} full U1 Mult1
set Knot ${KnotsPoles}
set exp_string " +${Knot} : +(\[-0-9.+eE\]+) +(\[-0-9.+eE\]+)"
regexp ${exp_string} ${log} full U2 Mult2
}
puts "U1=${U1}"
puts "U2=${U2}"
if {[expr {$U2 - $U1}] < 1.0e-20} {
puts "Error: Wrong curve's range!" puts "Error: Wrong curve's range!"
} }

View File

@ -27,7 +27,7 @@ mkface ff1 s1
donly ff1 f2 donly ff1 f2
############################# #############################
set log [bopcurves ff1 f2] set log [bopcurves ff1 f2 -2d]
############################# #############################
regexp {Tolerance Reached=+([-0-9.+eE]+)\n+([-0-9.+eE]+)} ${log} full Toler NbCurv regexp {Tolerance Reached=+([-0-9.+eE]+)\n+([-0-9.+eE]+)} ${log} full Toler NbCurv
@ -48,30 +48,11 @@ if {${NbCurv} != ${GoodNbCurv}} {
#------------- #-------------
for {set i 1} {$i <= ${NbCurv}} {incr i} { for {set i 1} {$i <= ${NbCurv}} {incr i} {
set log [dump c_$i] bounds c_$i U1 U2
set dumptrimres [regexp {Trimmed curve\nParameters : +([-0-9.+eE]+) +([-0-9.+eE]+)} ${log} full U1 U2]
if {${dumptrimres} == 0} { dump U1 U2
regexp {Degree +([-0-9.+eE]+), +([-0-9.+eE]+) Poles, +([-0-9.+eE]+)} ${log} full Degree Poles KnotsPoles
puts "Degree=${Degree}" if {[dval $U2-$U1] < 1.0e-20} {
puts "Poles=${Poles}"
puts "KnotsPoles=${KnotsPoles}"
puts ""
set Knot 1
set exp_string "Knots :\n\n +${Knot} : +(\[-0-9.+eE\]+) +(\[-0-9.+eE\]+)"
regexp ${exp_string} ${log} full U1 Mult1
set Knot ${KnotsPoles}
set exp_string " +${Knot} : +(\[-0-9.+eE\]+) +(\[-0-9.+eE\]+)"
regexp ${exp_string} ${log} full U2 Mult2
}
puts "U1=${U1}"
puts "U2=${U2}"
if {[expr {$U2 - $U1}] < 1.0e-20} {
puts "Error: Wrong curve's range!" puts "Error: Wrong curve's range!"
} }

View File

@ -0,0 +1,39 @@
puts "================"
puts "OCC26699"
puts "================"
puts ""
#######################################################################
# Wrong section curves
#######################################################################
set MaxTol 1.e-7
set GoodNbCurv 2
restore [locate_data_file bug26699_f1.brep] f1
restore [locate_data_file bug26699_f2.brep] f2
set log [bopcurves f1 f2 -2d]
regexp {Tolerance Reached=+([-0-9.+eE]+)\n+([-0-9.+eE]+)} ${log} full Toler NbCurv
if {${Toler} > ${MaxTol}} {
puts "Error: Tolerance is too big!"
}
if {${NbCurv} != ${GoodNbCurv}} {
puts "Error: Curve Number is bad!"
}
set expL1 3.0
set expL2 3.0
regexp {The length c_1 is ([-0-9.+eE]+)} [length c_1] full ll1
regexp {The length c_2 is ([-0-9.+eE]+)} [length c_2] full ll2
checkreal "length c_1 " ${ll1} $expL1 0.0 1.0e-6
checkreal "length c_2 " ${ll2} $expL1 0.0 1.0e-6
axo
donly f* c_*
fit
set only_screen_axo 1