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

0021564: Intersection of two planar faces produces curve with too many poles

I ComputePurgedWLine() function:
Excess points in walking line are deleted when:
1) Distance between neighboring points too small.
2) Points lie in one pipe without big jump on chord length.

III
Fixed problem with extremaPC with too close knot distribution to [minParam, maxParam] borders.

IV ApproxInt_Approx.gxx
New division criteria in intersection approximator.

III Test case
Test cases update to the new behavior.
Test case for CR21564

Correction of test cases for issue CR21564
This commit is contained in:
aml 2015-08-13 11:04:03 +03:00 committed by ski
parent 7a324550c8
commit 0cbfb9f151
24 changed files with 481 additions and 99 deletions

View File

@ -332,13 +332,16 @@ void ApproxInt_Approx::Perform(const Handle(TheWLine)& theline,
nbmc++) {
myBezToBSpl.Append(myComputeLineBezier.Value(nbmc));
}
if(imax<indicemax) {
imin = imax;
imax = imin+nbpntbez;
OtherInter = Standard_True;
if((indicemax-imax)<(nbpntbez/2)) {
imax = indicemax;
}
if(imax<indicemax)
{
imin = imax;
imax = imin+nbpntbez;
OtherInter = Standard_True;
if((indicemax-imax)<(nbpntbez/2))
{
imax = indicemax;
}
imax = CorrectFinishIdx(imin, imax, theline);
}
}
}
@ -567,14 +570,17 @@ void ApproxInt_Approx::Perform(const ThePSurface& Surf1,
nbmc++) {
myBezToBSpl.Append(myComputeLineBezier.Value(nbmc));
}
if(imax<indicemax) {
imin = imax;
imax = imin+nbpntbez;
OtherInter = Standard_True;
if((indicemax-imax)<(nbpntbez/2)) {
imax = indicemax;
}
}
if(imax<indicemax)
{
imin = imax;
imax = imin+nbpntbez;
OtherInter = Standard_True;
if((indicemax-imax)<(nbpntbez/2))
{
imax = indicemax;
}
imax = CorrectFinishIdx(imin, imax, theline);
}
}
}
while(OtherInter);
@ -903,13 +909,16 @@ void ApproxInt_Approx::Perform(const ThePSurface& PSurf,
nbmc++) {
myBezToBSpl.Append(myComputeLineBezier.Value(nbmc));
}
if(imax<indicemax) {
imin = imax;
imax = imin+nbpntbez;
OtherInter = Standard_True;
if((indicemax-imax)<(nbpntbez/2)) {
imax = indicemax;
}
if(imax<indicemax)
{
imin = imax;
imax = imin+nbpntbez;
OtherInter = Standard_True;
if((indicemax-imax)<(nbpntbez/2))
{
imax = indicemax;
}
imax = CorrectFinishIdx(imin, imax, theline);
}
}
}
@ -1117,21 +1126,26 @@ void ApproxInt_Approx::Perform(const TheISurface& ISurf,
}
}
OtherInter = Standard_False;
if(myApproxBez) {
for(Standard_Integer nbmc = 1;
nbmc <= myComputeLineBezier.NbMultiCurves() ;
nbmc++) {
myBezToBSpl.Append(myComputeLineBezier.Value(nbmc));
if(myApproxBez)
{
for(Standard_Integer nbmc = 1;
nbmc <= myComputeLineBezier.NbMultiCurves() ;
nbmc++)
{
myBezToBSpl.Append(myComputeLineBezier.Value(nbmc));
}
if(imax<indicemax) {
imin = imax;
imax = imin+nbpntbez;
OtherInter = Standard_True;
if((indicemax-imax)<(nbpntbez/2)) {
imax = indicemax;
}
if(imax<indicemax)
{
imin = imax;
imax = imin+nbpntbez;
OtherInter = Standard_True;
if((indicemax-imax)<(nbpntbez/2))
{
imax = indicemax;
}
imax = CorrectFinishIdx(imin, imax, theline);
}
}
}
}
while(OtherInter);
if(myApproxBez) {
@ -1208,3 +1222,31 @@ const AppParCurves_MultiBSpCurve& ApproxInt_Approx::Value(const Standard_Integer
return(myComputeLine.Value());
}
}
//--------------------------------------------------------------------------------
Standard_Integer ApproxInt_Approx::CorrectFinishIdx(const Standard_Integer theMinIdx,
const Standard_Integer theMaxIdx,
const Handle(TheWLine)& theline)
{
const Standard_Real aNullCoeff = 1.0e-16;
Standard_Real aLimitMaxCoeff = 1.0 / 2500.0;
Standard_Real aDist = theline->Point(theMinIdx).Value().SquareDistance(
theline->Point(theMinIdx + 1).Value());
for(Standard_Integer anIdx = theMinIdx + 1; anIdx < theMaxIdx - 1; anIdx++)
{
Standard_Real aNextDist = theline->Point(anIdx).Value().SquareDistance(
theline->Point(anIdx + 1).Value());
Standard_Real aCoeff = Min (aNextDist, aDist) / Max (aNextDist, aDist);
//
if (aCoeff < aLimitMaxCoeff && // Base criteria.
aNextDist > aDist && // Step increasing.
aNextDist > aNullCoeff && // Avoid separation in case of too small step.
aDist > aNullCoeff) // Usually found when purger not invoked (blend).
{
return anIdx;
}
aDist = aNextDist;
}
return theMaxIdx;
}

View File

@ -90,6 +90,9 @@ protected:
private:
Standard_EXPORT Standard_Integer CorrectFinishIdx(const Standard_Integer theMinIdx,
const Standard_Integer theMaxIdx,
const Handle(BRepApprox_ApproxLine)& theline);
Standard_EXPORT void Perform (const BRepAdaptor_Surface& Surf1, const IntSurf_Quadric& Surf2, const Handle(BRepApprox_ApproxLine)& aLine, const Standard_Boolean ApproxXYZ, const Standard_Boolean ApproxU1V1, const Standard_Boolean ApproxU2V2, const Standard_Integer indicemin, const Standard_Integer indicemax);

View File

@ -105,12 +105,14 @@ void Extrema_GExtPC::Perform(const ThePoint& P)
// Workaround to work with:
// blend, where knots may be moved from param space.
Standard_Real aPeriodJump = 0.0;
// Avoid prolem with too close knots.
const Standard_Real aTolCoeff = (myusup - myuinf) * Precision::PConfusion();
if (TheCurveTool::IsPeriodic(aCurve))
{
Standard_Integer aPeriodShift =
Standard_Integer ((myuinf - aKnots(aFirstIdx)) / TheCurveTool::Period(aCurve));
if (myuinf < aKnots(aFirstIdx))
if (myuinf < aKnots(aFirstIdx) - aTolCoeff)
aPeriodShift--;
aPeriodJump = TheCurveTool::Period(aCurve) * aPeriodShift;
@ -124,7 +126,7 @@ void Extrema_GExtPC::Perform(const ThePoint& P)
for(anIdx = aFirstIdx; anIdx <= aLastIdx; anIdx++)
{
Standard_Real aKnot = aKnots(anIdx) + aPeriodJump;
if (myuinf >= aKnot)
if (myuinf >= aKnot - aTolCoeff)
aFirstUsedKnot = anIdx;
else
break;
@ -133,7 +135,7 @@ void Extrema_GExtPC::Perform(const ThePoint& P)
for(anIdx = aLastIdx; anIdx >= aFirstIdx; anIdx--)
{
Standard_Real aKnot = aKnots(anIdx) + aPeriodJump;
if (myusup <= aKnot)
if (myusup <= aKnot + aTolCoeff)
aLastUsedKnot = anIdx;
else
break;

View File

@ -90,7 +90,10 @@ protected:
private:
Standard_EXPORT Standard_Integer CorrectFinishIdx(const Standard_Integer theMinIdx,
const Standard_Integer theMaxIdx,
const Handle(IntPatch_WLine)& theline);
Standard_EXPORT void Perform (const Handle(Adaptor3d_HSurface)& Surf1, const IntSurf_Quadric& Surf2, const Handle(IntPatch_WLine)& aLine, const Standard_Boolean ApproxXYZ, const Standard_Boolean ApproxU1V1, const Standard_Boolean ApproxU2V2, const Standard_Integer indicemin, const Standard_Integer indicemax);
Standard_EXPORT void Perform (const IntSurf_Quadric& Surf1, const Handle(Adaptor3d_HSurface)& Surf2, const Handle(IntPatch_WLine)& aLine, const Standard_Boolean ApproxXYZ, const Standard_Boolean ApproxU1V1, const Standard_Boolean ApproxU2V2, const Standard_Integer indicemin, const Standard_Integer indicemax);

View File

@ -768,26 +768,38 @@ void IntPatch_RstInt::PutVertexOnLine (const Handle(IntPatch_Line)& L,
}
Standard_Boolean refined = Standard_False;
if (refine[ip]) {
//------------------------------------------------------------------------
//-- On a trouve un point 2d approche Ua,Va intersection de la ligne
//-- de cheminement et de la restriction.
//--
//-- On injecte ce point ds les intersections Courbe-Surface
//--
IntPatch_CSFunction thefunc(OtherSurf,arc,Surf);
// MSV: extend UV bounds to not miss solution near the boundary
Standard_Real margCoef = 0.004;
IntPatch_CurvIntSurf IntCS(U,V,W,thefunc,edgeTol,margCoef);
if (IntCS.IsDone()) {
if (!IntCS.IsEmpty()) {
ptsommet = IntCS.Point();
IntCS.ParameterOnSurface(U2,V2);
paramarc = IntCS.ParameterOnCurve();
refined = Standard_True;
}
}
}
if (refine[ip])
{
//------------------------------------------------------------------------
//-- On a trouve un point 2d approche Ua,Va intersection de la ligne
//-- de cheminement et de la restriction.
//--
//-- On injecte ce point ds les intersections Courbe-Surface
//--
IntPatch_CSFunction thefunc(OtherSurf,arc,Surf);
// MSV: extend UV bounds to not miss solution near the boundary
Standard_Real margCoef = 0.004;
IntPatch_CurvIntSurf IntCS(U,V,W,thefunc,edgeTol,margCoef);
if (IntCS.IsDone())
{
if (!IntCS.IsEmpty())
{
ptsommet = IntCS.Point();
IntCS.ParameterOnSurface(U2,V2);
gp_Pnt anOldPnt, aNewPnt;
OtherSurf->D0(U,V, anOldPnt);
OtherSurf->D0(U2,V2, aNewPnt);
if (anOldPnt.SquareDistance(aNewPnt) < Precision::Confusion()
* Precision::Confusion())
{
U2 = U;
V2 = V;
}
paramarc = IntCS.ParameterOnCurve();
refined = Standard_True;
}
}
}
else {
U2 = U; V2 = V;
paramarc = W;

View File

@ -162,7 +162,9 @@ static
const Standard_Integer ilprm);
static
Handle(IntPatch_WLine) ComputePurgedWLine(const Handle(IntPatch_WLine)& theWLine);
Handle(IntPatch_WLine) ComputePurgedWLine(const Handle(IntPatch_WLine) &theWLine,
const Handle(GeomAdaptor_HSurface) &theS1,
const Handle(GeomAdaptor_HSurface) &theS2);
static
Handle(Geom2d_BSplineCurve) MakeBSpline2d(const Handle(IntPatch_WLine)& theWLine,
@ -680,7 +682,7 @@ void IntTools_FaceFace::Perform(const TopoDS_Face& aF1,
{
const Standard_Real UVMaxStep = 0.001;
const Standard_Real Deflection = (hasCone) ? 0.085 : 0.1;
myIntersector.SetTolerances(TolArc, TolTang, UVMaxStep, Deflection);
myIntersector.SetTolerances(TolArc, TolTang, UVMaxStep, Deflection);
}
if((myHS1->IsUClosed() && !myHS1->IsUPeriodic()) ||
@ -876,7 +878,7 @@ Standard_Real IntTools_FaceFace::ComputeTolerance()
//function :ComputeTolReached3d
//purpose :
//=======================================================================
void IntTools_FaceFace::ComputeTolReached3d()
void IntTools_FaceFace::ComputeTolReached3d()
{
Standard_Integer aNbLin;
GeomAbs_SurfaceType aType1, aType2;
@ -923,9 +925,9 @@ void IntTools_FaceFace::ComputeTolReached3d()
Standard_Real aDMax = ComputeTolerance();
if (aDMax > myTolReached3d)
{
myTolReached3d = aDMax;
myTolReached3d = aDMax;
}
}
}
//=======================================================================
//function : MakeCurve
@ -963,14 +965,12 @@ void IntTools_FaceFace::ComputeTolReached3d()
Handle(IntPatch_Line) anewL;
//
Handle(IntPatch_WLine) aWLine (Handle(IntPatch_WLine)::DownCast(L));
//DumpWLine(aWLine);
anewL = ComputePurgedWLine(aWLine);
anewL = ComputePurgedWLine(aWLine, myHS1, myHS2);
if(anewL.IsNull()) {
return;
}
L = anewL;
//Handle(IntPatch_WLine) aWLineX (Handle(IntPatch_WLine)::DownCast(L));
//DumpWLine(aWLineX);
@ -1884,9 +1884,11 @@ void IntTools_FaceFace::ComputeTolReached3d()
}
//
mySeqOfCurve.Append(aCurve);
}//if(typs1 == GeomAbs_Plane) {
else if(typs2 == GeomAbs_Plane) {
else if(typs2 == GeomAbs_Plane)
{
const AppParCurves_MultiBSpCurve& mbspc = theapp3d.Value(j);
nbpoles = mbspc.NbPoles();
@ -1995,6 +1997,7 @@ void IntTools_FaceFace::ComputeTolReached3d()
} else {
mySeqOfCurve.Append(aCurve);
}
}// else if(typs2 == GeomAbs_Plane)
//
else { //typs2 != GeomAbs_Plane && typs1 != GeomAbs_Plane
@ -2788,15 +2791,53 @@ Standard_Boolean IsDegeneratedZone(const gp_Pnt2d& aP2d,
return !bFlag;
}
// Check if aNextPnt lies inside of tube build on aBasePnt and aBaseVec.
// In 2d space.
static Standard_Boolean IsInsideIn2d(const gp_Pnt2d& aBasePnt,
const gp_Vec2d& aBaseVec,
const gp_Pnt2d& aNextPnt,
const Standard_Real aSquareMaxDist)
{
gp_Vec2d aVec2d(aBasePnt, aNextPnt);
//d*d = (basevec^(nextpnt-basepnt))**2 / basevec**2
Standard_Real aCross = aVec2d.Crossed(aBaseVec);
Standard_Real aSquareDist = aCross * aCross
/ aBaseVec.SquareMagnitude();
return (aSquareDist <= aSquareMaxDist);
}
// Check if aNextPnt lies inside of tube build on aBasePnt and aBaseVec.
// In 3d space.
static Standard_Boolean IsInsideIn3d(const gp_Pnt& aBasePnt,
const gp_Vec& aBaseVec,
const gp_Pnt& aNextPnt,
const Standard_Real aSquareMaxDist)
{
gp_Vec aVec(aBasePnt, aNextPnt);
//d*d = (basevec^(nextpnt-basepnt))**2 / basevec**2
Standard_Real aSquareDist = aVec.CrossSquareMagnitude(aBaseVec)
/ aBaseVec.SquareMagnitude();
return (aSquareDist <= aSquareMaxDist);
}
//=========================================================================
// static function : ComputePurgedWLine
// purpose : Removes equal points (leave one of equal points) from theWLine
// and recompute vertex parameters.
// Removes exceed points using tube criteria:
// delete 7D point if it lies near to expected lines in 2d and 3d.
// Each task (2d, 2d, 3d) have its own tolerance and checked separately.
// Returns new WLine or null WLine if the number
// of the points is less than 2.
//=========================================================================
Handle(IntPatch_WLine) ComputePurgedWLine(const Handle(IntPatch_WLine)& theWLine) {
Handle(IntPatch_WLine) ComputePurgedWLine(const Handle(IntPatch_WLine) &theWLine,
const Handle(GeomAdaptor_HSurface) &theS1,
const Handle(GeomAdaptor_HSurface) &theS2)
{
Standard_Integer i, k, v, nb, nbvtx;
Handle(IntPatch_WLine) aResult;
nbvtx = theWLine->NbVertex();
@ -2820,7 +2861,6 @@ Handle(IntPatch_WLine) ComputePurgedWLine(const Handle(IntPatch_WLine)& theWLine
for(v = 1; v <= nbvtx; v++) {
aLocalWLine->AddVertex(theWLine->Vertex(v));
}
for(i = 1; i <= aLineOn2S->NbPoints(); i++) {
Standard_Integer aStartIndex = i + 1;
Standard_Integer anEndIndex = i + 5;
@ -2838,10 +2878,24 @@ Handle(IntPatch_WLine) ComputePurgedWLine(const Handle(IntPatch_WLine)& theWLine
IntSurf_PntOn2S p1 = aLineOn2S->Value(i);
IntSurf_PntOn2S p2 = aLineOn2S->Value(k);
if(p1.Value().IsEqual(p2.Value(), gp::Resolution())) {
Standard_Real UV[8];
p1.Parameters(UV[0], UV[1], UV[2], UV[3]);
p2.Parameters(UV[4], UV[5], UV[6], UV[7]);
Standard_Real aMax = Abs(UV[0]);
for(Standard_Integer anIdx = 1; anIdx < 8; anIdx++)
{
if (aMax < Abs(UV[anIdx]))
aMax = Abs(UV[anIdx]);
}
if(p1.Value().IsEqual(p2.Value(), gp::Resolution()) ||
Abs(UV[0] - UV[4]) + Abs(UV[1] - UV[5]) < 1.0e-16 * aMax ||
Abs(UV[2] - UV[6]) + Abs(UV[3] - UV[7]) < 1.0e-16 * aMax )
{
aTmpWLine = aLocalWLine;
aLocalWLine = new IntPatch_WLine(aLineOn2S, Standard_False);
for(v = 1; v <= aTmpWLine->NbVertex(); v++) {
IntPatch_Point aVertex = aTmpWLine->Vertex(v);
Standard_Integer avertexindex = (Standard_Integer)aVertex.ParameterOnLine();
@ -2860,7 +2914,211 @@ Handle(IntPatch_WLine) ComputePurgedWLine(const Handle(IntPatch_WLine)& theWLine
}
}
if(aLineOn2S->NbPoints() > 1) {
if (aLineOn2S->NbPoints() <= 2)
{
if (aLineOn2S->NbPoints() == 2)
return aLocalWLine;
else
return aResult;
}
const Standard_Integer aMinNbBadDistr = 15;
const Standard_Integer aNbSingleBezier = 30;
// Avoid purge in case of C0 continuity:
// Intersection approximator may produce invalid curve after purge, example:
// bugs modalg_5 bug24731,
// Do not run purger when base number of points is too small.
if (theS1->UContinuity() == GeomAbs_C0 ||
theS1->VContinuity() == GeomAbs_C0 ||
theS2->UContinuity() == GeomAbs_C0 ||
theS2->VContinuity() == GeomAbs_C0 ||
nb < aNbSingleBezier)
{
return aLocalWLine;
}
// 1 - Delete point.
// 0 - Store point.
// -1 - Vertex point (not delete).
NCollection_Array1<Standard_Integer> aNewPointsHash(1, aLineOn2S->NbPoints());
for(i = 1; i <= aLineOn2S->NbPoints(); i++)
aNewPointsHash.SetValue(i, 0);
for(v = 1; v <= aLocalWLine->NbVertex(); v++)
{
IntPatch_Point aVertex = aLocalWLine->Vertex(v);
Standard_Integer avertexindex = (Standard_Integer)aVertex.ParameterOnLine();
aNewPointsHash.SetValue(avertexindex, -1);
}
// Workaround to handle case of small amount points after purge.
// Test "boolean boptuc_complex B5" and similar.
Standard_Integer aNbPnt = 0;
// Inital computations.
Standard_Real UonS1[3], VonS1[3], UonS2[3], VonS2[3];
aLineOn2S->Value(1).ParametersOnS1(UonS1[0], VonS1[0]);
aLineOn2S->Value(2).ParametersOnS1(UonS1[1], VonS1[1]);
aLineOn2S->Value(1).ParametersOnS2(UonS2[0], VonS2[0]);
aLineOn2S->Value(2).ParametersOnS2(UonS2[1], VonS2[1]);
gp_Pnt2d aBase2dPnt1(UonS1[0], VonS1[0]);
gp_Pnt2d aBase2dPnt2(UonS2[0], VonS2[0]);
gp_Vec2d aBase2dVec1(UonS1[1] - UonS1[0], VonS1[1] - VonS1[0]);
gp_Vec2d aBase2dVec2(UonS2[1] - UonS2[0], VonS2[1] - VonS2[0]);
gp_Pnt aBase3dPnt = aLineOn2S->Value(1).Value();
gp_Vec aBase3dVec(aLineOn2S->Value(1).Value(), aLineOn2S->Value(2).Value());
// Choose base tolerance and scale it to pipe algorithm.
const Standard_Real aBaseTolerance = Precision::Approximation();
Standard_Real aResS1Tol = Min(theS1->UResolution(aBaseTolerance),
theS1->VResolution(aBaseTolerance));
Standard_Real aResS2Tol = Min(theS2->UResolution(aBaseTolerance),
theS2->VResolution(aBaseTolerance));
Standard_Real aTol1 = aResS1Tol * aResS1Tol;
Standard_Real aTol2 = aResS2Tol * aResS2Tol;
Standard_Real aTol3d = aBaseTolerance * aBaseTolerance;
const Standard_Real aLimitCoeff = 0.99 * 0.99;
for(i = 3; i <= aLineOn2S->NbPoints(); i++)
{
Standard_Boolean isDeleteState = Standard_False;
aLineOn2S->Value(i).ParametersOnS1(UonS1[2], VonS1[2]);
aLineOn2S->Value(i).ParametersOnS2(UonS2[2], VonS2[2]);
gp_Pnt2d aPnt2dOnS1(UonS1[2], VonS1[2]);
gp_Pnt2d aPnt2dOnS2(UonS2[2], VonS2[2]);
const gp_Pnt& aPnt3d = aLineOn2S->Value(i).Value();
if (aNewPointsHash(i - 1) != - 1 &&
IsInsideIn2d(aBase2dPnt1, aBase2dVec1, aPnt2dOnS1, aTol1) &&
IsInsideIn2d(aBase2dPnt2, aBase2dVec2, aPnt2dOnS2, aTol2) &&
IsInsideIn3d(aBase3dPnt, aBase3dVec, aPnt3d, aTol3d) )
{
// Handle possible uneven parametrization on one of 2d subspaces.
// Delete point only when expected lengths are close to each other (aLimitCoeff).
// Example:
// c2d1 - line
// c3d - line
// c2d2 - geometrically line, but have uneven parametrization -> c2d2 is bspline.
gp_XY aPntOnS1[2]= { gp_XY(UonS1[1] - UonS1[0], VonS1[1] - VonS1[0])
, gp_XY(UonS1[2] - UonS1[1], VonS1[2] - VonS1[1])};
gp_XY aPntOnS2[2]= { gp_XY(UonS2[1] - UonS2[0], VonS2[1] - VonS2[0])
, gp_XY(UonS2[2] - UonS2[1], VonS2[2] - VonS2[1])};
Standard_Real aStepOnS1 = aPntOnS1[0].SquareModulus() / aPntOnS1[1].SquareModulus();
Standard_Real aStepOnS2 = aPntOnS2[0].SquareModulus() / aPntOnS2[1].SquareModulus();
Standard_Real aStepCoeff = Min(aStepOnS1, aStepOnS2) / Max(aStepOnS1, aStepOnS2);
if (aStepCoeff > aLimitCoeff)
{
// Set hash flag to "Delete" state.
isDeleteState = Standard_True;
aNewPointsHash.SetValue(i - 1, 1);
// Change middle point.
UonS1[1] = UonS1[2];
UonS2[1] = UonS2[2];
VonS1[1] = VonS1[2];
VonS2[1] = VonS2[2];
}
}
if (!isDeleteState)
{
// Compute new pipe parameters.
UonS1[0] = UonS1[1];
VonS1[0] = VonS1[1];
UonS2[0] = UonS2[1];
VonS2[0] = VonS2[1];
UonS1[1] = UonS1[2];
VonS1[1] = VonS1[2];
UonS2[1] = UonS2[2];
VonS2[1] = VonS2[2];
aBase2dPnt1.SetCoord(UonS1[0], VonS1[0]);
aBase2dPnt2.SetCoord(UonS2[0], VonS2[0]);
aBase2dVec1.SetCoord(UonS1[1] - UonS1[0], VonS1[1] - VonS1[0]);
aBase2dVec2.SetCoord(UonS2[1] - UonS2[0], VonS2[1] - VonS2[0]);
aBase3dPnt = aLineOn2S->Value(i - 1).Value();
aBase3dVec = gp_Vec(aLineOn2S->Value(i - 1).Value(), aLineOn2S->Value(i).Value());
aNbPnt++;
}
}
// Workaround to handle case of small amount of points after purge.
// Test "boolean boptuc_complex B5" and similar.
// This is possible since there are at least two points.
if (aNewPointsHash(1) == -1 &&
aNewPointsHash(2) == -1 &&
aNbPnt <= 3)
{
// Delete first.
aNewPointsHash(1) = 1;
}
if (aNewPointsHash(aLineOn2S->NbPoints() - 1) == -1 &&
aNewPointsHash(aLineOn2S->NbPoints() ) == -1 &&
aNbPnt <= 3)
{
// Delete last.
aNewPointsHash(aLineOn2S->NbPoints() ) = 1;
}
// Purgre when too small amount of points left.
if (aNbPnt <= 2)
{
for(i = aNewPointsHash.Lower(); i <= aNewPointsHash.Upper(); i++)
{
if (aNewPointsHash(i) != -1)
{
aNewPointsHash(i) = 1;
}
}
}
// Handle possible bad distribution of points,
// which are will converted into one single bezier curve (less than 30 points).
// Make distribution more even:
// max step will be nearly to 0.1 of param distance.
if (aNbPnt + 2 > aMinNbBadDistr &&
aNbPnt + 2 < aNbSingleBezier )
{
for(Standard_Integer anIdx = 1; anIdx <= 8; anIdx++)
{
Standard_Integer aHashIdx =
Standard_Integer(anIdx * aLineOn2S->NbPoints() / 9);
//Store this point.
aNewPointsHash(aHashIdx) = 0;
}
}
aTmpWLine = aLocalWLine;
Handle(IntSurf_LineOn2S) aPurgedLineOn2S = new IntSurf_LineOn2S();
aLocalWLine = new IntPatch_WLine(aPurgedLineOn2S, Standard_False);
Standard_Integer anOldLineIdx = 1, aVertexIdx = 1;
for(i = 1; i <= aNewPointsHash.Upper(); i++)
{
if (aNewPointsHash(i) == 0)
{
// Store this point.
aPurgedLineOn2S->Add(aLineOn2S->Value(i));
anOldLineIdx++;
}
else if (aNewPointsHash(i) == -1)
{
// Add vertex.
IntPatch_Point aVertex = aTmpWLine->Vertex(aVertexIdx++);
aVertex.SetParameter(anOldLineIdx++);
aLocalWLine->AddVertex(aVertex);
aPurgedLineOn2S->Add(aLineOn2S->Value(i));
}
}
if(aPurgedLineOn2S->NbPoints() > 1) {
aResult = aLocalWLine;
}
return aResult;

View File

@ -7,4 +7,4 @@ restore [locate_data_file buc60462b.brep] b
bsection result a b
set length 261.262
set length 261.53

View File

@ -1,3 +1,4 @@
puts "TODO OCC21564 ALL: Error : The square of result shape is"
puts "=========="
puts "BUC60532"
puts "=========="

View File

@ -33,5 +33,5 @@ if { $MaxFaceTolerance > 1 || $MaxEdgeTolerance > 1 || $MaxVertexTolerance > 1 }
puts "Tolerance of shape is less then 1.0"
}
set square 4.60842e+07
set square 4.10276e+007
set 2dviewer 0

View File

@ -50,5 +50,5 @@ if { $MaxFaceTolerance > 1 || $MaxEdgeTolerance > 1 || $MaxVertexTolerance > 1 }
} else {
puts "Tolerance of shape is less then 1.0"
}
set square 4.52817e+07
set square 3.92639e+007
set 2dviewer 0

View File

@ -61,10 +61,10 @@ regexp { +Face +: +Min +[-0-9.+eE]+ +Max +([-0-9.+eE]+)} $tolerance full MaxFace
regexp { +Edge +: +Min +[-0-9.+eE]+ +Max +([-0-9.+eE]+)} $tolerance full MaxEdgeTolerance
regexp { +Vertex +: +Min +[-0-9.+eE]+ +Max +([-0-9.+eE]+)} $tolerance full MaxVertexTolerance
if { $MaxFaceTolerance > 1 || $MaxEdgeTolerance > 1 || $MaxVertexTolerance > 1 } {
puts "Faulty : Tolerance of shape is more then 1.0"
if { $MaxFaceTolerance > 2 || $MaxEdgeTolerance > 2 || $MaxVertexTolerance > 2 } {
puts "Faulty : Tolerance of shape is more then 2.0"
} else {
puts "Tolerance of shape is less then 1.0"
puts "Tolerance of shape is less then 2.0"
}
set square 2.34799e+007
set square 1.74934e+007
set 2dviewer 0

View File

@ -1,3 +1,6 @@
puts "TODO OCC21564 ALL: The square of result shape is"
puts "TODO OCC21564 ALL: Result shape is WRONG because it must contains"
puts "============"
puts "OCC22557"
puts "============"

View File

@ -1,5 +1,5 @@
puts "TODO OCC25829 ALL: Error : The square of result shape is"
puts "TODO OCC25829 ALL: Faulty shapes in variables faulty_1 to"
puts "TODO OCC25829 Linux: Faulty shapes in variables faulty_1 to"
puts "============"
puts "OCC697"

View File

@ -1,5 +1,5 @@
puts "TODO OCC25829 ALL: Error : The square of result shape is"
puts "TODO OCC25829 ALL: Faulty shapes in variables faulty_1 to"
puts "TODO OCC25829 Linux: Faulty shapes in variables faulty_1 to"
puts "============"
puts "OCC697"

View File

@ -1,5 +1,5 @@
puts "TODO OCC25829 ALL: Error : The square of result shape is"
puts "TODO OCC25829 ALL: Faulty shapes in variables faulty_1 to"
puts "TODO OCC25829 Linux: Faulty shapes in variables faulty_1 to"
puts "============"
puts "OCC697"

View File

@ -1,5 +1,5 @@
puts "TODO OCC25829 ALL: Error : The square of result shape is"
puts "TODO OCC25829 ALL: Faulty shapes in variables faulty_1 to"
puts "TODO OCC25829 Linux: Faulty shapes in variables faulty_1 to"
puts "============"
puts "OCC697"

View File

@ -1,4 +1,3 @@
puts "========"
puts "OCC770"
puts "SAM1636"

View File

@ -0,0 +1,61 @@
puts "================"
puts "OCC21564"
puts "================"
puts ""
#######################################################################
# Intersection of two planar faces produces curve with too many poles
#######################################################################
puts "Construct horizontal plane, convert it to b-spline, and create a face on it"
plane p1 0 0 0 0 0 1
trim p1t p1 -10 10 -10 10
convert bs1 p1t
mkface f1 bs1
puts "Construct vertical plane and face on it"
plane p2 0 0 0 1 0 0
mkface f2 p2 -10 10 -10 10
puts "\nBuild intersection of two faces"
bsection sec f1 f2
puts "Check number of points in resulting curve"
if { ! [regexp {([0-9]+) Poles} [dump sec] str nbp] } {
puts "Error: Could not check number of poles in resulting curve!"
} elseif { $nbp != 2 } {
puts "Error: Intersection curve has too many poles ($nbp while 2 is expected)"
}
puts "\nBuild intersection with approximation of resulting curve"
bsection sec f1 f2 -a
puts "Check number of points in resulting curve"
if { ! [regexp {([0-9]+) Poles} [dump sec] str nbp] } {
puts "Error: Could not check number of poles in resulting curve!"
} elseif { $nbp != 2 } {
puts "Error: Intersection curve has too many poles ($nbp while 2 is expected)"
}
puts "\nNow trying with inclined plane"
plane p2 0 0 0 1 2 0
mkface f2 p2 -15 15 -15 15
puts "\nBuild intersection of two faces"
bsection sec f1 f2
puts "Check number of points in resulting curve"
if { ! [regexp {([0-9]+) Poles} [dump sec] str nbp] } {
puts "Error: Could not check number of poles in resulting curve!"
} elseif { $nbp != 2 } {
puts "Error: Intersection curve has too many poles ($nbp while 2 is expected)"
}
puts "\nBuild intersection with approximation of resulting curve"
bsection sec f1 f2 -a
puts "Check number of points in resulting curve"
if { ! [regexp {([0-9]+) Poles} [dump sec] str nbp] } {
puts "Error: Could not check number of poles in resulting curve!"
} elseif { $nbp != 2 } {
puts "Error: Intersection curve has too many poles ($nbp while 2 is expected)"
}

View File

@ -12,7 +12,7 @@ restore [locate_data_file bug23884_fz124] b2
bop b1 b2
bopfuse result
set square 2415.65
set square 2368.48
set 2dviewer 0

View File

@ -5,7 +5,7 @@ puts ""
###########################################################
# Wrong pcurve of the section curve
###########################################################
set MaxTol 1.9e-5
set MaxTol 1.05e-6
set NbCurv_OK 1
restore [locate_data_file bug24585_b1.brep] b1

View File

@ -6,7 +6,7 @@ puts ""
# Wrong pcurve of the section curve
###########################################################
set MaxTol 5.0e-7
set MaxTol 7.0e-5
set NbCurv_OK 1
restore [locate_data_file bug24612_b1.brep] b1
restore [locate_data_file bug24612_b2.brep] b2

View File

@ -1,4 +1,11 @@
puts "TODO OCC25929 ALL: Error: Tolerance is too big!"
puts "TODO OCC21564 Linux: Error : T=0.464646 D=0.000326627"
puts "TODO OCC21564 Linux: Error : T=0.464646 D=0.00032747"
puts "TODO OCC21564 Windows: Error : T=0.464646 D=0.000326671"
puts "TODO OCC21564 Windows: Error : T=0.464646 D=0.000327516"
puts "========="
puts "CR24915"
puts "========="

View File

@ -1,8 +1,3 @@
puts "TODO OCC25929 ALL: ERROR in bopargcheck res1"
puts "TODO OCC25929 ALL: ERROR. res1 with continuity C0"
puts "TODO OCC25929 ALL: ERROR in bopargcheck res2"
puts "TODO OCC25929 ALL: ERROR. res2 with continuity C0"
puts "========"
puts "OCC26310"
puts "========"

View File

@ -1,7 +1,3 @@
puts "TODO OCC24156 MacOS: \\*\\* Exception \\*\\*.*"
puts "TODO OCC24156 MacOS: An exception was caught"
puts "TODO OCC24156 MacOS: TEST INCOMPLETE"
beziercurve w1 5 0 0 0 20 0 0 20 5 0 25 10 0 10 20 0
mkedge w1 w1
polyline w2 10 20 0 0 10 0 0 0 0