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

0030489: Modeling Algorithms - BRepBuilderAPI_GTransform hangs

Approx_ComputeCLine.gxx - criterium to stop interval cutting is increased.
ProjLib_ProjectedCurve.cxx - fix regression for bugs modalg_5 bug25886
Other tests are modified according to current state of algorithm
This commit is contained in:
ifv 2019-03-12 16:47:21 +03:00 committed by apn
parent 2328cae25d
commit afb3647b34
8 changed files with 224 additions and 122 deletions

View File

@ -34,14 +34,14 @@
//======================================================================= //=======================================================================
Approx_ComputeCLine::Approx_ComputeCLine Approx_ComputeCLine::Approx_ComputeCLine
(const MultiLine& Line, (const MultiLine& Line,
const Standard_Integer degreemin, const Standard_Integer degreemin,
const Standard_Integer degreemax, const Standard_Integer degreemax,
const Standard_Real Tolerance3d, const Standard_Real Tolerance3d,
const Standard_Real Tolerance2d, const Standard_Real Tolerance2d,
const Standard_Boolean cutting, const Standard_Boolean cutting,
const AppParCurves_Constraint FirstC, const AppParCurves_Constraint FirstC,
const AppParCurves_Constraint LastC) const AppParCurves_Constraint LastC)
{ {
mydegremin = degreemin; mydegremin = degreemin;
mydegremax = degreemax; mydegremax = degreemax;
@ -61,13 +61,13 @@ Approx_ComputeCLine::Approx_ComputeCLine
//======================================================================= //=======================================================================
Approx_ComputeCLine::Approx_ComputeCLine Approx_ComputeCLine::Approx_ComputeCLine
(const Standard_Integer degreemin, (const Standard_Integer degreemin,
const Standard_Integer degreemax, const Standard_Integer degreemax,
const Standard_Real Tolerance3d, const Standard_Real Tolerance3d,
const Standard_Real Tolerance2d, const Standard_Real Tolerance2d,
const Standard_Boolean cutting, const Standard_Boolean cutting,
const AppParCurves_Constraint FirstC, const AppParCurves_Constraint FirstC,
const AppParCurves_Constraint LastC) const AppParCurves_Constraint LastC)
{ {
alldone = Standard_False; alldone = Standard_False;
mydegremin = degreemin; mydegremin = degreemin;
@ -88,21 +88,22 @@ Approx_ComputeCLine::Approx_ComputeCLine
void Approx_ComputeCLine::Perform(const MultiLine& Line) void Approx_ComputeCLine::Perform(const MultiLine& Line)
{ {
Standard_Real UFirst, ULast; Standard_Real UFirst, ULast;
Standard_Boolean Finish = Standard_False, Standard_Boolean Finish = Standard_False,
begin = Standard_True, Ok = Standard_False; begin = Standard_True, Ok = Standard_False;
Standard_Real thetol3d = Precision::Confusion(), thetol2d = Precision::Confusion(); Standard_Real thetol3d = Precision::Confusion(), thetol2d = Precision::Confusion();
UFirst = Line.FirstParameter(); UFirst = Line.FirstParameter();
ULast = Line.LastParameter(); ULast = Line.LastParameter();
Standard_Real TolU = Max((ULast-UFirst)*1.e-05, Precision::PApproximation()); Standard_Real TolU = Max((ULast - UFirst)*1.e-03, Precision::Confusion());
Standard_Real myfirstU = UFirst; Standard_Real myfirstU = UFirst;
Standard_Real mylastU = ULast; Standard_Real mylastU = ULast;
Standard_Integer aMaxSegments = 0; Standard_Integer aMaxSegments = 0;
Standard_Integer aMaxSegments1 = myMaxSegments - 1; Standard_Integer aMaxSegments1 = myMaxSegments - 1;
Standard_Integer aNbCut = 0, aNbImp = 0, aNbComp = 5;
if (!mycut) if (!mycut)
{ {
alldone = Compute(Line, UFirst, ULast, thetol3d, thetol2d); alldone = Compute(Line, UFirst, ULast, thetol3d, thetol2d);
if (!alldone) if (!alldone)
{ {
tolreached = Standard_False; tolreached = Standard_False;
myfirstparam.Append(UFirst); myfirstparam.Append(UFirst);
@ -112,25 +113,27 @@ void Approx_ComputeCLine::Perform(const MultiLine& Line)
Tolers2d.Append(currenttol2d); Tolers2d.Append(currenttol2d);
} }
} }
else else
{ {
// previous decision to be taken if we get worse with next cut (eap) // previous decision to be taken if we get worse with next cut (eap)
AppParCurves_MultiCurve KeptMultiCurve; AppParCurves_MultiCurve KeptMultiCurve;
Standard_Real KeptUfirst = 0., KeptUlast = 0., KeptT3d = RealLast(), KeptT2d = 0.; Standard_Real KeptUfirst = 0., KeptUlast = 0., KeptT3d = RealLast(), KeptT2d = 0.;
while (!Finish) while (!Finish)
{ {
// Gestion du decoupage de la multiline pour approximer: // Gestion du decoupage de la multiline pour approximer:
if (!begin) if (!begin)
{ {
if (Ok) if (Ok)
{ {
// Calcul de la partie a approximer. // Calcul de la partie a approximer.
myfirstU = mylastU; myfirstU = mylastU;
mylastU = ULast; mylastU = ULast;
if (Abs(ULast-myfirstU) <= RealEpsilon() aNbCut = 0;
aNbImp = 0;
if (Abs(ULast - myfirstU) <= RealEpsilon()
|| aMaxSegments >= myMaxSegments) || aMaxSegments >= myMaxSegments)
{ {
Finish = Standard_True; Finish = Standard_True;
@ -147,50 +150,59 @@ void Approx_ComputeCLine::Perform(const MultiLine& Line)
if ((thetol3d + thetol2d) < (KeptT3d + KeptT2d)) if ((thetol3d + thetol2d) < (KeptT3d + KeptT2d))
{ {
KeptMultiCurve = TheMultiCurve; KeptMultiCurve = TheMultiCurve;
KeptUfirst = myfirstU; KeptUfirst = myfirstU;
KeptUlast = mylastU; KeptUlast = mylastU;
KeptT3d = thetol3d; KeptT3d = thetol3d;
KeptT2d = thetol2d; KeptT2d = thetol2d;
aNbImp++;
} }
// cut an interval // cut an interval
mylastU = (myfirstU + mylastU)/2; mylastU = (myfirstU + mylastU) / 2;
aNbCut++;
} }
} }
// Calcul des parametres sur ce nouvel intervalle. // Calcul des parametres sur ce nouvel intervalle.
Ok = Compute(Line, myfirstU, mylastU, thetol3d, thetol2d); Ok = Compute(Line, myfirstU, mylastU, thetol3d, thetol2d);
if(Ok) if (Ok)
{ {
aMaxSegments++; aMaxSegments++;
} }
//cout << myfirstU << " - " << mylastU << " tol : " << thetol3d << " " << thetol2d << endl; //cout << myfirstU << " - " << mylastU << " tol : " << thetol3d << " " << thetol2d << endl;
Standard_Boolean aStopCutting = Standard_False;
// is new decision better? if (aNbCut >= aNbComp)
if (!Ok && (Abs(myfirstU-mylastU) <= TolU || aMaxSegments >= aMaxSegments1))
{ {
Ok = Standard_True; // stop interval cutting, approx the rest part if (aNbCut > aNbImp)
{
if ((thetol3d + thetol2d) < (KeptT3d + KeptT2d)) aStopCutting = Standard_True;
{
KeptMultiCurve = TheMultiCurve;
KeptUfirst = myfirstU;
KeptUlast = mylastU;
KeptT3d = thetol3d;
KeptT2d = thetol2d;
}
mylastU = KeptUlast;
tolreached = Standard_False; // helas
myMultiCurves.Append(KeptMultiCurve);
aMaxSegments++;
Tolers3d.Append (KeptT3d);
Tolers2d.Append (KeptT2d);
myfirstparam.Append (KeptUfirst);
mylastparam.Append (KeptUlast);
} }
}
// is new decision better?
if (!Ok && (Abs(myfirstU - mylastU) <= TolU || aMaxSegments >= aMaxSegments1 || aStopCutting ))
{
Ok = Standard_True; // stop interval cutting, approx the rest part
if ((thetol3d + thetol2d) < (KeptT3d + KeptT2d))
{
KeptMultiCurve = TheMultiCurve;
KeptUfirst = myfirstU;
KeptUlast = mylastU;
KeptT3d = thetol3d;
KeptT2d = thetol2d;
}
mylastU = KeptUlast;
tolreached = Standard_False; // helas
myMultiCurves.Append(KeptMultiCurve);
aMaxSegments++;
Tolers3d.Append(KeptT3d);
Tolers2d.Append(KeptT2d);
myfirstparam.Append(KeptUfirst);
mylastparam.Append(KeptUlast);
}
begin = Standard_False; begin = Standard_False;
} // while (!Finish) } // while (!Finish)
@ -225,10 +237,10 @@ const
//======================================================================= //=======================================================================
Standard_Boolean Approx_ComputeCLine::Compute(const MultiLine& Line, Standard_Boolean Approx_ComputeCLine::Compute(const MultiLine& Line,
const Standard_Real Ufirst, const Standard_Real Ufirst,
const Standard_Real Ulast, const Standard_Real Ulast,
Standard_Real& TheTol3d, Standard_Real& TheTol3d,
Standard_Real& TheTol2d) Standard_Real& TheTol2d)
{ {
@ -243,14 +255,14 @@ Standard_Boolean Approx_ComputeCLine::Compute(const MultiLine& Line,
if (mydone) { if (mydone) {
LSquare.Error(Fv, TheTol3d, TheTol2d); LSquare.Error(Fv, TheTol3d, TheTol2d);
if (TheTol3d <= mytol3d && TheTol2d <= mytol2d) { if (TheTol3d <= mytol3d && TheTol2d <= mytol2d) {
// Stockage de la multicurve approximee. // Stockage de la multicurve approximee.
tolreached = Standard_True; tolreached = Standard_True;
myMultiCurves.Append(LSquare.Value()); myMultiCurves.Append(LSquare.Value());
myfirstparam.Append(Ufirst); myfirstparam.Append(Ufirst);
mylastparam.Append(Ulast); mylastparam.Append(Ulast);
Tolers3d.Append(TheTol3d); Tolers3d.Append(TheTol3d);
Tolers2d.Append(TheTol2d); Tolers2d.Append(TheTol2d);
return Standard_True; return Standard_True;
} }
} }
if (deg == mydegremax) { if (deg == mydegremax) {
@ -258,7 +270,7 @@ Standard_Boolean Approx_ComputeCLine::Compute(const MultiLine& Line,
currenttol3d = TheTol3d; currenttol3d = TheTol3d;
currenttol2d = TheTol2d; currenttol2d = TheTol2d;
} }
} }
return Standard_False; return Standard_False;
} }
@ -270,11 +282,11 @@ Standard_Boolean Approx_ComputeCLine::Compute(const MultiLine& Line,
//======================================================================= //=======================================================================
void Approx_ComputeCLine::Parameters(const Standard_Integer Index, void Approx_ComputeCLine::Parameters(const Standard_Integer Index,
Standard_Real& firstpar, Standard_Real& firstpar,
Standard_Real& lastpar) const Standard_Real& lastpar) const
{ {
firstpar = myfirstparam.Value(Index); firstpar = myfirstparam.Value(Index);
lastpar = mylastparam.Value(Index); lastpar = mylastparam.Value(Index);
} }
//======================================================================= //=======================================================================
@ -283,7 +295,7 @@ void Approx_ComputeCLine::Parameters(const Standard_Integer Index,
//======================================================================= //=======================================================================
void Approx_ComputeCLine::SetDegrees(const Standard_Integer degreemin, void Approx_ComputeCLine::SetDegrees(const Standard_Integer degreemin,
const Standard_Integer degreemax) const Standard_Integer degreemax)
{ {
mydegremin = degreemin; mydegremin = degreemin;
mydegremax = degreemax; mydegremax = degreemax;
@ -295,7 +307,7 @@ void Approx_ComputeCLine::SetDegrees(const Standard_Integer degreemin,
//======================================================================= //=======================================================================
void Approx_ComputeCLine::SetTolerances(const Standard_Real Tolerance3d, void Approx_ComputeCLine::SetTolerances(const Standard_Real Tolerance3d,
const Standard_Real Tolerance2d) const Standard_Real Tolerance2d)
{ {
mytol3d = Tolerance3d; mytol3d = Tolerance3d;
mytol2d = Tolerance2d; mytol2d = Tolerance2d;
@ -307,10 +319,10 @@ void Approx_ComputeCLine::SetTolerances(const Standard_Real Tolerance3d,
//======================================================================= //=======================================================================
void Approx_ComputeCLine::SetConstraints(const AppParCurves_Constraint FirstC, void Approx_ComputeCLine::SetConstraints(const AppParCurves_Constraint FirstC,
const AppParCurves_Constraint LastC) const AppParCurves_Constraint LastC)
{ {
myfirstC = FirstC; myfirstC = FirstC;
mylastC = LastC; mylastC = LastC;
} }
//======================================================================= //=======================================================================
@ -318,7 +330,7 @@ void Approx_ComputeCLine::SetConstraints(const AppParCurves_Constraint FirstC,
//purpose : Changes the max number of segments, which is allowed for cutting. //purpose : Changes the max number of segments, which is allowed for cutting.
//======================================================================= //=======================================================================
void Approx_ComputeCLine:: SetMaxSegments(const Standard_Integer theMaxSegments) void Approx_ComputeCLine::SetMaxSegments(const Standard_Integer theMaxSegments)
{ {
myMaxSegments = theMaxSegments; myMaxSegments = theMaxSegments;
} }
@ -351,8 +363,8 @@ const {
//======================================================================= //=======================================================================
void Approx_ComputeCLine::Error(const Standard_Integer Index, void Approx_ComputeCLine::Error(const Standard_Integer Index,
Standard_Real& tol3d, Standard_Real& tol3d,
Standard_Real& tol2d) const Standard_Real& tol2d) const
{ {
tol3d = Tolers3d.Value(Index); tol3d = Tolers3d.Value(Index);
tol2d = Tolers2d.Value(Index); tol2d = Tolers2d.Value(Index);

View File

@ -42,6 +42,7 @@
#include <TColgp_HArray1OfPnt2d.hxx> #include <TColgp_HArray1OfPnt2d.hxx>
#include <TColStd_HArray1OfReal.hxx> #include <TColStd_HArray1OfReal.hxx>
#include <Geom2dConvert_CompCurveToBSplineCurve.hxx> #include <Geom2dConvert_CompCurveToBSplineCurve.hxx>
#include <Geom2dConvert.hxx>
#include <TColStd_Array1OfReal.hxx> #include <TColStd_Array1OfReal.hxx>
#include <TColStd_Array1OfInteger.hxx> #include <TColStd_Array1OfInteger.hxx>
#include <TColgp_Array1OfPnt2d.hxx> #include <TColgp_Array1OfPnt2d.hxx>
@ -55,6 +56,8 @@
#include <GeomLib.hxx> #include <GeomLib.hxx>
#include <Extrema_ExtPC.hxx> #include <Extrema_ExtPC.hxx>
#include <NCollection_DataMap.hxx> #include <NCollection_DataMap.hxx>
#include <ElSLib.hxx>
#include <ElCLib.hxx>
//======================================================================= //=======================================================================
//function : ComputeTolU //function : ComputeTolU
//purpose : //purpose :
@ -146,18 +149,19 @@ static Standard_Boolean IsoIsDeg (const Adaptor3d_Surface& S,
//======================================================================= //=======================================================================
static void TrimC3d(Handle(Adaptor3d_HCurve)& myCurve, static void TrimC3d(Handle(Adaptor3d_HCurve)& myCurve,
Standard_Boolean* IsTrimmed, Standard_Boolean* IsTrimmed,
const Standard_Real dt, const Standard_Real dt,
const gp_Pnt& Pole, const gp_Pnt& Pole,
Standard_Integer* SingularCase, Standard_Integer* SingularCase,
const Standard_Integer NumberOfSingularCase) const Standard_Integer NumberOfSingularCase,
const Standard_Real TolConf)
{ {
Standard_Real f = myCurve->FirstParameter(); Standard_Real f = myCurve->FirstParameter();
Standard_Real l = myCurve->LastParameter(); Standard_Real l = myCurve->LastParameter();
gp_Pnt P = myCurve->Value(f); gp_Pnt P = myCurve->Value(f);
if(P.Distance(Pole) < Precision::Confusion()) { if(P.Distance(Pole) <= TolConf) {
IsTrimmed[0] = Standard_True; IsTrimmed[0] = Standard_True;
f = f+dt; f = f+dt;
myCurve = myCurve->Trim(f, l, Precision::Confusion()); myCurve = myCurve->Trim(f, l, Precision::Confusion());
@ -165,7 +169,7 @@ static void TrimC3d(Handle(Adaptor3d_HCurve)& myCurve,
} }
P = myCurve->Value(l); P = myCurve->Value(l);
if(P.Distance(Pole) < Precision::Confusion()) { if(P.Distance(Pole) <= TolConf) {
IsTrimmed[1] = Standard_True; IsTrimmed[1] = Standard_True;
l = l-dt; l = l-dt;
myCurve = myCurve->Trim(f, l, Precision::Confusion()); myCurve = myCurve->Trim(f, l, Precision::Confusion());
@ -232,17 +236,24 @@ static void ExtendC2d (Handle(Geom2d_BSplineCurve)& aRes,
} }
} }
gp_Lin2d BoundLin(thePole, theBoundDir); //one of the bounds of rectangle gp_Lin2d BoundLin(thePole, theBoundDir); //one of the bounds of rectangle
Standard_Real ParOnLin = 0.;
if (theBoundDir.IsParallel(aDBnd, 100.*Precision::Angular()))
{
ParOnLin = ElCLib::Parameter(aLin, thePole);
}
else
{
Standard_Real U1x = BoundLin.Direction().X();
Standard_Real U1y = BoundLin.Direction().Y();
Standard_Real U2x = aLin.Direction().X();
Standard_Real U2y = aLin.Direction().Y();
Standard_Real Uo21x = aLin.Location().X() - BoundLin.Location().X();
Standard_Real Uo21y = aLin.Location().Y() - BoundLin.Location().Y();
Standard_Real U1x = BoundLin.Direction().X(); Standard_Real D = U1y*U2x - U1x*U2y;
Standard_Real U1y = BoundLin.Direction().Y();
Standard_Real U2x = aLin.Direction().X(); ParOnLin = (Uo21y * U1x - Uo21x * U1y) / D; //parameter of intersection point
Standard_Real U2y = aLin.Direction().Y(); }
Standard_Real Uo21x = aLin.Location().X() - BoundLin.Location().X();
Standard_Real Uo21y = aLin.Location().Y() - BoundLin.Location().Y();
Standard_Real D = U1y*U2x-U1x*U2y;
Standard_Real ParOnLin = (Uo21y * U1x - Uo21x * U1y)/D; //parameter of intersection point
Handle(Geom2d_Line) aSegLine = new Geom2d_Line(aLin); Handle(Geom2d_Line) aSegLine = new Geom2d_Line(aLin);
aSegment = (FirstOrLast == 0)? aSegment = (FirstOrLast == 0)?
@ -392,6 +403,16 @@ void ProjLib_ProjectedCurve::Perform(const Handle(Adaptor3d_HCurve)& C)
GeomAbs_SurfaceType SType = mySurface->GetType(); GeomAbs_SurfaceType SType = mySurface->GetType();
GeomAbs_CurveType CType = myCurve->GetType(); GeomAbs_CurveType CType = myCurve->GetType();
Standard_Boolean isAnalyticalSurf = Standard_True; Standard_Boolean isAnalyticalSurf = Standard_True;
Standard_Boolean IsTrimmed[2] = { Standard_False, Standard_False };
Standard_Integer SingularCase[2];
const Standard_Real eps = 0.01;
Standard_Real TolConf = Precision::Confusion();
Standard_Real dt = (LastPar - FirstPar) * eps;
Standard_Real U1 = 0.0, U2 = 0.0, V1 = 0.0, V2 = 0.0;
U1 = mySurface->FirstUParameter();
U2 = mySurface->LastUParameter();
V1 = mySurface->FirstVParameter();
V2 = mySurface->LastVParameter();
switch (SType) switch (SType)
{ {
@ -429,6 +450,28 @@ void ProjLib_ProjectedCurve::Perform(const Handle(Adaptor3d_HCurve)& C)
// periodique en V !) // periodique en V !)
P.SetInBounds(myCurve->FirstParameter()); P.SetInBounds(myCurve->FirstParameter());
} }
else
{
const Standard_Real Vmax = M_PI / 2.;
const Standard_Real Vmin = -Vmax;
const Standard_Real minang = 1.e-5 * M_PI;
gp_Sphere aSph = mySurface->Sphere();
Standard_Real anR = aSph.Radius();
Standard_Real f = myCurve->FirstParameter();
Standard_Real l = myCurve->LastParameter();
gp_Pnt Pf = myCurve->Value(f);
gp_Pnt Pl = myCurve->Value(l);
gp_Pnt aLoc = aSph.Position().Location();
Standard_Real maxdist = Max(Pf.Distance(aLoc), Pl.Distance(aLoc));
TolConf = Max(anR * minang, Abs(anR - maxdist));
//Surface has pole at V = Vmin and Vmax
gp_Pnt Pole = mySurface->Value(U1, Vmin);
TrimC3d(myCurve, IsTrimmed, dt, Pole, SingularCase, 3, TolConf);
Pole = mySurface->Value(U1, Vmax);
TrimC3d(myCurve, IsTrimmed, dt, Pole, SingularCase, 4, TolConf);
}
myResult = P; myResult = P;
} }
break; break;
@ -445,15 +488,11 @@ void ProjLib_ProjectedCurve::Perform(const Handle(Adaptor3d_HCurve)& C)
case GeomAbs_BSplineSurface: case GeomAbs_BSplineSurface:
{ {
isAnalyticalSurf = Standard_False; isAnalyticalSurf = Standard_False;
Standard_Boolean IsTrimmed[2] = {Standard_False, Standard_False}; Standard_Real f, l;
Standard_Integer SingularCase[2];
Standard_Real f, l, dt;
const Standard_Real eps = 0.01;
f = myCurve->FirstParameter(); f = myCurve->FirstParameter();
l = myCurve->LastParameter(); l = myCurve->LastParameter();
dt = (l - f) * eps; dt = (l - f) * eps;
Standard_Real U1 = 0.0, U2=0.0, V1=0.0, V2=0.0;
const Adaptor3d_Surface& S = mySurface->Surface(); const Adaptor3d_Surface& S = mySurface->Surface();
U1 = S.FirstUParameter(); U1 = S.FirstUParameter();
U2 = S.LastUParameter(); U2 = S.LastUParameter();
@ -464,28 +503,28 @@ void ProjLib_ProjectedCurve::Perform(const Handle(Adaptor3d_HCurve)& C)
{ {
//Surface has pole at U = Umin //Surface has pole at U = Umin
gp_Pnt Pole = mySurface->Value(U1, V1); gp_Pnt Pole = mySurface->Value(U1, V1);
TrimC3d(myCurve, IsTrimmed, dt, Pole, SingularCase, 1); TrimC3d(myCurve, IsTrimmed, dt, Pole, SingularCase, 1, TolConf);
} }
if(IsoIsDeg(S, U2, GeomAbs_IsoU, 0., myTolerance)) if(IsoIsDeg(S, U2, GeomAbs_IsoU, 0., myTolerance))
{ {
//Surface has pole at U = Umax //Surface has pole at U = Umax
gp_Pnt Pole = mySurface->Value(U2, V1); gp_Pnt Pole = mySurface->Value(U2, V1);
TrimC3d(myCurve, IsTrimmed, dt, Pole, SingularCase, 2); TrimC3d(myCurve, IsTrimmed, dt, Pole, SingularCase, 2, TolConf);
} }
if(IsoIsDeg(S, V1, GeomAbs_IsoV, 0., myTolerance)) if(IsoIsDeg(S, V1, GeomAbs_IsoV, 0., myTolerance))
{ {
//Surface has pole at V = Vmin //Surface has pole at V = Vmin
gp_Pnt Pole = mySurface->Value(U1, V1); gp_Pnt Pole = mySurface->Value(U1, V1);
TrimC3d(myCurve, IsTrimmed, dt, Pole, SingularCase, 3); TrimC3d(myCurve, IsTrimmed, dt, Pole, SingularCase, 3, TolConf);
} }
if(IsoIsDeg(S, V2, GeomAbs_IsoV, 0., myTolerance)) if(IsoIsDeg(S, V2, GeomAbs_IsoV, 0., myTolerance))
{ {
//Surface has pole at V = Vmax //Surface has pole at V = Vmax
gp_Pnt Pole = mySurface->Value(U1, V2); gp_Pnt Pole = mySurface->Value(U1, V2);
TrimC3d(myCurve, IsTrimmed, dt, Pole, SingularCase, 4); TrimC3d(myCurve, IsTrimmed, dt, Pole, SingularCase, 4, TolConf);
} }
ProjLib_ComputeApproxOnPolarSurface polar; ProjLib_ComputeApproxOnPolarSurface polar;
@ -531,10 +570,9 @@ void ProjLib_ProjectedCurve::Perform(const Handle(Adaptor3d_HCurve)& C)
default: default:
{ {
isAnalyticalSurf = Standard_False; isAnalyticalSurf = Standard_False;
Standard_Boolean IsTrimmed[2] = {Standard_False, Standard_False};
Standard_Real Vsingular[2] = {0.0 , 0.0}; //for surfaces of revolution Standard_Real Vsingular[2] = {0.0 , 0.0}; //for surfaces of revolution
Standard_Real f = 0.0, l = 0.0, dt = 0.0; Standard_Real f = 0.0, l = 0.0;
const Standard_Real eps = 0.01; dt = 0.0;
if(mySurface->GetType() == GeomAbs_SurfaceOfRevolution) if(mySurface->GetType() == GeomAbs_SurfaceOfRevolution)
{ {
@ -710,26 +748,59 @@ void ProjLib_ProjectedCurve::Perform(const Handle(Adaptor3d_HCurve)& C)
Comp.SetDegree(myDegMin, myDegMax); Comp.SetDegree(myDegMin, myDegMax);
Comp.SetMaxSegments(myMaxSegments); Comp.SetMaxSegments(myMaxSegments);
Comp.SetBndPnt(myBndPnt); Comp.SetBndPnt(myBndPnt);
Comp.Perform( myCurve, mySurface); Comp.Perform(myCurve, mySurface);
if (Comp.Bezier().IsNull() && Comp.BSpline().IsNull()) if (Comp.Bezier().IsNull() && Comp.BSpline().IsNull())
return; // advanced projector has been failed too return; // advanced projector has been failed too
myResult.Done(); myResult.Done();
Handle(Geom2d_BSplineCurve) aRes;
// set the type if (Comp.BSpline().IsNull())
if ( SType == GeomAbs_Plane && CType == GeomAbs_BezierCurve)
{ {
myResult.SetType(GeomAbs_BezierCurve); aRes = Geom2dConvert::CurveToBSplineCurve(Comp.Bezier());
myResult.SetBezier(Comp.Bezier()) ;
} }
else else
{ {
aRes = Comp.BSpline();
}
if ((IsTrimmed[0] || IsTrimmed[1]))
{
if (IsTrimmed[0])
{
//Add segment before start of curve
Standard_Real f = myCurve->FirstParameter();
ExtendC2d(aRes, f, -dt, U1, U2, V1, V2, 0, SingularCase[0]);
}
if (IsTrimmed[1])
{
//Add segment after end of curve
Standard_Real l = myCurve->LastParameter();
ExtendC2d(aRes, l, dt, U1, U2, V1, V2, 1, SingularCase[1]);
}
Handle(Geom2d_Curve) NewCurve2d;
GeomLib::SameRange(Precision::PConfusion(), aRes,
aRes->FirstParameter(), aRes->LastParameter(),
FirstPar, LastPar, NewCurve2d);
aRes = Handle(Geom2d_BSplineCurve)::DownCast(NewCurve2d);
myResult.SetBSpline(aRes);
myResult.SetType(GeomAbs_BSplineCurve); myResult.SetType(GeomAbs_BSplineCurve);
myResult.SetBSpline(Comp.BSpline()) ; }
else
{
// set the type
if (SType == GeomAbs_Plane && CType == GeomAbs_BezierCurve)
{
myResult.SetType(GeomAbs_BezierCurve);
myResult.SetBezier(Comp.Bezier());
}
else
{
myResult.SetType(GeomAbs_BSplineCurve);
myResult.SetBSpline(Comp.BSpline());
}
} }
// set the periodicity flag // set the periodicity flag
if (SType == GeomAbs_Plane && if (SType == GeomAbs_Plane &&
CType == GeomAbs_BSplineCurve && CType == GeomAbs_BSplineCurve &&
myCurve->IsPeriodic() ) myCurve->IsPeriodic())
{ {
myResult.SetPeriodic(); myResult.SetPeriodic();
} }

View File

@ -6,8 +6,8 @@ puts ""
# Wrong result obtained by projection algorithm # Wrong result obtained by projection algorithm
################################################# #################################################
set ok_len_c3x "3.28347" set ok_len_c3x "1.57079"
set ok_len_c5x "3.28346" set ok_len_c5x "1.57079"
smallview -2D- smallview -2D-

View File

@ -26,7 +26,7 @@ regexp {Tolerance Reached=+([-0-9.+eE]+)\n+([-0-9.+eE]+)} ${log} full Tolerance_
set GoodNbCurv 2 set GoodNbCurv 2
set expected_Tolerance_Reached 1.2482990218170969e-007 set expected_Tolerance_Reached 1.042095984078466e-05
set tol_abs_Tolerance_Reached 1.0e-7 set tol_abs_Tolerance_Reached 1.0e-7
set tol_rel_Tolerance_Reached 0.0 set tol_rel_Tolerance_Reached 0.0
checkreal "Tolerance Reached" ${Tolerance_Reached} ${expected_Tolerance_Reached} ${tol_abs_Tolerance_Reached} ${tol_rel_Tolerance_Reached} checkreal "Tolerance Reached" ${Tolerance_Reached} ${expected_Tolerance_Reached} ${tol_abs_Tolerance_Reached} ${tol_rel_Tolerance_Reached}

View File

@ -31,7 +31,7 @@ regexp {Tolerance Reached=+([-0-9.+eE]+)\n+([-0-9.+eE]+)} ${log} full Tolerance_
set GoodNbCurv 2 set GoodNbCurv 2
set expected_Tolerance_Reached 1.2482990218170969e-007 set expected_Tolerance_Reached 1.0420959841458885e-05
set tol_abs_Tolerance_Reached 1.0e-7 set tol_abs_Tolerance_Reached 1.0e-7
set tol_rel_Tolerance_Reached 0.0 set tol_rel_Tolerance_Reached 0.0
checkreal "Tolerance Reached" ${Tolerance_Reached} ${expected_Tolerance_Reached} ${tol_abs_Tolerance_Reached} ${tol_rel_Tolerance_Reached} checkreal "Tolerance Reached" ${Tolerance_Reached} ${expected_Tolerance_Reached} ${tol_abs_Tolerance_Reached} ${tol_rel_Tolerance_Reached}

View File

@ -26,6 +26,6 @@ checkreal "Reached tolerance" ${Tolerance} 5.8654166482879483e-009 1.e-7 0
set bop_info_2d [bopcurves f1 f2 -2d] set bop_info_2d [bopcurves f1 f2 -2d]
regexp {Tolerance Reached=([-0-9.+eE]+)} $bop_info_2d full Tolerance_2d regexp {Tolerance Reached=([-0-9.+eE]+)} $bop_info_2d full Tolerance_2d
checkreal "Reached tolerance" ${Tolerance_2d} 1.4569392656749484e-008 1.e-7 0 checkreal "Reached tolerance" ${Tolerance_2d} 1.4915699300398263e-07 1.e-7 0
checkview -screenshot -2d -path ${imagedir}/${test_image}.png checkview -screenshot -2d -path ${imagedir}/${test_image}.png

View File

@ -26,6 +26,6 @@ checkreal "Reached tolerance" ${Tolerance} 1.2530391548405894e-008 1.e-7 0
set bop_info_2d [bopcurves f1 f2 -2d] set bop_info_2d [bopcurves f1 f2 -2d]
regexp {Tolerance Reached=([-0-9.+eE]+)} $bop_info_2d full Tolerance_2d regexp {Tolerance Reached=([-0-9.+eE]+)} $bop_info_2d full Tolerance_2d
checkreal "Reached tolerance" ${Tolerance_2d} 1.4134494834137484e-005 0 1.e-2 checkreal "Reached tolerance" ${Tolerance_2d} 0.00040497924613267194 0 1.e-2
checkview -screenshot -2d -path ${imagedir}/${test_image}.png checkview -screenshot -2d -path ${imagedir}/${test_image}.png

View File

@ -0,0 +1,19 @@
puts "============"
puts "BUC30489"
puts "============"
puts ""
###############################
## 0030489: Modeling Algorithms - BRepBuilderAPI_GTransform hangs
###############################
restore [locate_data_file bug30489.brep] a
dchrono h restart
scalexyz r a 1. 1. 1.
dchrono h stop counter scalexyz
checkprops r -s 2.86706e+007
checkshape r
checkview -display r -2d -path ${imagedir}/${test_image}.png