1
0
mirror of https://git.dev.opencascade.org/repos/occt.git synced 2025-08-09 13:22:24 +03:00

0023981: Wrong section curves

Test case for issue CR23981
Small correction of test case for issue CR23981
This commit is contained in:
ifv
2013-10-03 16:13:10 +04:00
committed by bugmaster
parent b1c5c4e6a6
commit 7c4e9501b4
7 changed files with 220 additions and 247 deletions

View File

@@ -24,7 +24,7 @@
#include <gp_Pnt.hxx> #include <gp_Pnt.hxx>
#include <gp_Trsf.hxx> #include <gp_Trsf.hxx>
#define myInfinite 1.e15 #define myInfinite Precision::Infinite()
static void GetConeApexParam(const gp_Cone& C, Standard_Real& U, Standard_Real& V) static void GetConeApexParam(const gp_Cone& C, Standard_Real& U, Standard_Real& V)
{ {
@@ -282,12 +282,12 @@ TopAbs_State Adaptor3d_TopolTool::Classify(const gp_Pnt2d& P,
if (nbRestr == 4) { if (nbRestr == 4) {
if ((U < Uinf - Tol) || (U > Usup + Tol) || if ((U < Uinf - Tol) || (U > Usup + Tol) ||
(V < Vinf - Tol) || (V > Vsup + Tol)) { (V < Vinf - Tol) || (V > Vsup + Tol)) {
return TopAbs_OUT; return TopAbs_OUT;
} }
if ((Abs(U - Uinf) <= Tol) || (Abs(U - Usup) <= Tol) || if ((Abs(U - Uinf) <= Tol) || (Abs(U - Usup) <= Tol) ||
(Abs(V - Vinf) <= Tol) || (Abs(V - Vsup) <= Tol)) { (Abs(V - Vinf) <= Tol) || (Abs(V - Vsup) <= Tol)) {
return TopAbs_ON; return TopAbs_ON;
} }
return TopAbs_IN; return TopAbs_IN;
} }
@@ -297,100 +297,100 @@ TopAbs_State Adaptor3d_TopolTool::Classify(const gp_Pnt2d& P,
else { else {
Standard_Boolean dansu,dansv,surumin,surumax,survmin,survmax; Standard_Boolean dansu,dansv,surumin,surumax,survmin,survmax;
if (Precision::IsNegativeInfinite(Uinf) && if (Precision::IsNegativeInfinite(Uinf) &&
Precision::IsPositiveInfinite(Usup)) { Precision::IsPositiveInfinite(Usup)) {
dansu = Standard_True; dansu = Standard_True;
surumin = surumax = Standard_False; surumin = surumax = Standard_False;
} }
else if (Precision::IsNegativeInfinite(Uinf)) { else if (Precision::IsNegativeInfinite(Uinf)) {
surumin = Standard_False; surumin = Standard_False;
if (U >= Usup+Tol) { if (U >= Usup+Tol) {
dansu = Standard_False; dansu = Standard_False;
surumax = Standard_False; surumax = Standard_False;
} }
else { else {
dansu = Standard_True; dansu = Standard_True;
surumax = Standard_False; surumax = Standard_False;
if (Abs(U-Usup)<=Tol) { if (Abs(U-Usup)<=Tol) {
surumax = Standard_True; surumax = Standard_True;
} }
} }
} }
else if (Precision::IsPositiveInfinite(Usup)) { else if (Precision::IsPositiveInfinite(Usup)) {
surumax = Standard_False; surumax = Standard_False;
if (U < Uinf-Tol) { if (U < Uinf-Tol) {
dansu = Standard_False; dansu = Standard_False;
surumin = Standard_False; surumin = Standard_False;
} }
else { else {
dansu = Standard_True; dansu = Standard_True;
surumin = Standard_False; surumin = Standard_False;
if (Abs(U-Uinf)<=Tol) { if (Abs(U-Uinf)<=Tol) {
surumin = Standard_True; surumin = Standard_True;
} }
} }
} }
else { else {
if ((U < Uinf - Tol) || (U > Usup + Tol)) { if ((U < Uinf - Tol) || (U > Usup + Tol)) {
surumin = surumax = dansu = Standard_False; surumin = surumax = dansu = Standard_False;
} }
else { else {
dansu = Standard_True; dansu = Standard_True;
surumin = surumax = Standard_False; surumin = surumax = Standard_False;
if (Abs(U-Uinf)<=Tol) { if (Abs(U-Uinf)<=Tol) {
surumin = Standard_True; surumin = Standard_True;
} }
else if (Abs(U-Usup)<=Tol) { else if (Abs(U-Usup)<=Tol) {
surumax = Standard_True; surumax = Standard_True;
} }
} }
} }
if (Precision::IsNegativeInfinite(Vinf) && if (Precision::IsNegativeInfinite(Vinf) &&
Precision::IsPositiveInfinite(Vsup)) { Precision::IsPositiveInfinite(Vsup)) {
dansv = Standard_True; dansv = Standard_True;
survmin = survmax = Standard_False; survmin = survmax = Standard_False;
} }
else if (Precision::IsNegativeInfinite(Vinf)) { else if (Precision::IsNegativeInfinite(Vinf)) {
survmin = Standard_False; survmin = Standard_False;
if (V > Vsup+Tol) { if (V > Vsup+Tol) {
dansv = Standard_False; dansv = Standard_False;
survmax = Standard_False; survmax = Standard_False;
} }
else { else {
dansv = Standard_True; dansv = Standard_True;
survmax = Standard_False; survmax = Standard_False;
if (Abs(V-Vsup)<=Tol) { if (Abs(V-Vsup)<=Tol) {
survmax = Standard_True; survmax = Standard_True;
} }
} }
} }
else if (Precision::IsPositiveInfinite(Vsup)) { else if (Precision::IsPositiveInfinite(Vsup)) {
survmax = Standard_False; survmax = Standard_False;
if (V < Vinf-Tol) { if (V < Vinf-Tol) {
dansv = Standard_False; dansv = Standard_False;
survmin = Standard_False; survmin = Standard_False;
} }
else { else {
dansv = Standard_True; dansv = Standard_True;
survmin = Standard_False; survmin = Standard_False;
if (Abs(V-Vinf)<=Tol) { if (Abs(V-Vinf)<=Tol) {
survmin = Standard_True; survmin = Standard_True;
} }
} }
} }
else { else {
if ((V < Vinf - Tol) || (V > Vsup + Tol)) { if ((V < Vinf - Tol) || (V > Vsup + Tol)) {
survmin = survmax = dansv = Standard_False; survmin = survmax = dansv = Standard_False;
} }
else { else {
dansv = Standard_True; dansv = Standard_True;
survmin = survmax = Standard_False; survmin = survmax = Standard_False;
if (Abs(V-Vinf)<=Tol) { if (Abs(V-Vinf)<=Tol) {
survmin = Standard_True; survmin = Standard_True;
} }
else if (Abs(V-Vsup)<=Tol) { else if (Abs(V-Vsup)<=Tol) {
survmax = Standard_True; survmax = Standard_True;
} }
} }
} }
@@ -1269,12 +1269,12 @@ void Adaptor3d_TopolTool::BSplSamplePnts(const Standard_Real theDefl,
// //
for(j = 0, i = 1; i <= nbsu; ++i) { for(j = 0, i = 1; i <= nbsu; ++i) {
if (bFlag) { if (bFlag) {
myUPars->SetValue(i,anUPars(i)); myUPars->SetValue(i,anUPars(i));
} }
else { else {
if(anUFlg(i)) { if(anUFlg(i)) {
++j; ++j;
myUPars->SetValue(j,anUPars(i)); myUPars->SetValue(j,anUPars(i));
} }
} }
} }
@@ -1293,8 +1293,8 @@ void Adaptor3d_TopolTool::BSplSamplePnts(const Standard_Real theDefl,
} }
else { else {
if(aVFlg(i)) { if(aVFlg(i)) {
++j; ++j;
myVPars->SetValue(j,aVPars(i)); myVPars->SetValue(j,aVPars(i));
} }
} }
} }

View File

@@ -22,6 +22,7 @@
#include <Standard_OutOfRange.hxx> #include <Standard_OutOfRange.hxx>
#include <Standard_NotImplemented.hxx> #include <Standard_NotImplemented.hxx>
#include <ElSLib.hxx> #include <ElSLib.hxx>
static const Standard_Real ExtPElS_MyEps = Epsilon(2. * M_PI);
//============================================================================= //=============================================================================
Extrema_ExtPElS::Extrema_ExtPElS () { myDone = Standard_False; } Extrema_ExtPElS::Extrema_ExtPElS () { myDone = Standard_False; }
@@ -70,6 +71,7 @@ void Extrema_ExtPElS::Perform(const gp_Pnt& P,
if (OPp.Magnitude() < Tol) { return; } if (OPp.Magnitude() < Tol) { return; }
gp_Vec myZ = Pos.XDirection()^Pos.YDirection(); gp_Vec myZ = Pos.XDirection()^Pos.YDirection();
Standard_Real U1 = gp_Vec(Pos.XDirection()).AngleWithRef(OPp,myZ); //-M_PI<U1<M_PI Standard_Real U1 = gp_Vec(Pos.XDirection()).AngleWithRef(OPp,myZ); //-M_PI<U1<M_PI
if (U1 > -ExtPElS_MyEps && U1 < ExtPElS_MyEps) { U1 = 0.; }
Standard_Real U2 = U1 + M_PI; Standard_Real U2 = U1 + M_PI;
if (U1 < 0.) { U1 += 2. * M_PI; } if (U1 < 0.) { U1 += 2. * M_PI; }
@@ -164,6 +166,7 @@ void Extrema_ExtPElS::Perform(const gp_Pnt& P,
Standard_Real B, U1, V1, U2, V2; Standard_Real B, U1, V1, U2, V2;
Standard_Boolean Same = DirZ.Dot(MP) >= 0.0; Standard_Boolean Same = DirZ.Dot(MP) >= 0.0;
U1 = gp_Vec(Pos.XDirection()).AngleWithRef(OPp,myZ); //-M_PI<U1<M_PI U1 = gp_Vec(Pos.XDirection()).AngleWithRef(OPp,myZ); //-M_PI<U1<M_PI
if (U1 > -ExtPElS_MyEps && U1 < ExtPElS_MyEps) { U1 = 0.; }
B = MP.Angle(DirZ); B = MP.Angle(DirZ);
if (!Same) { U1 += M_PI; } if (!Same) { U1 += M_PI; }
U2 = U1 + M_PI; U2 = U1 + M_PI;
@@ -257,6 +260,7 @@ void Extrema_ExtPElS::Perform(const gp_Pnt& P,
else { else {
gp_Vec myZ = Pos.XDirection()^Pos.YDirection(); gp_Vec myZ = Pos.XDirection()^Pos.YDirection();
U1 = gp_Vec(Pos.XDirection()).AngleWithRef(OPp,myZ); U1 = gp_Vec(Pos.XDirection()).AngleWithRef(OPp,myZ);
if (U1 > -ExtPElS_MyEps && U1 < ExtPElS_MyEps) { U1 = 0.; }
U2 = U1 + M_PI; U2 = U1 + M_PI;
if (U1 < 0.) { U1 += 2. * M_PI; } if (U1 < 0.) { U1 += 2. * M_PI; }
V = OP.Angle(OPp); V = OP.Angle(OPp);
@@ -321,6 +325,7 @@ void Extrema_ExtPElS::Perform(const gp_Pnt& P,
gp_Vec myZ = Pos.XDirection()^Pos.YDirection(); gp_Vec myZ = Pos.XDirection()^Pos.YDirection();
Standard_Real U1 = gp_Vec(Pos.XDirection()).AngleWithRef(OPp,myZ); Standard_Real U1 = gp_Vec(Pos.XDirection()).AngleWithRef(OPp,myZ);
if (U1 > -ExtPElS_MyEps && U1 < ExtPElS_MyEps) { U1 = 0.; }
Standard_Real U2 = U1 + M_PI; Standard_Real U2 = U1 + M_PI;
if (U1 < 0.) { U1 += 2. * M_PI; } if (U1 < 0.) { U1 += 2. * M_PI; }
Standard_Real R = sqrt(R2); Standard_Real R = sqrt(R2);
@@ -333,7 +338,10 @@ void Extrema_ExtPElS::Perform(const gp_Pnt& P,
if(O2.SquareDistance(P) < Tol) { return; } if(O2.SquareDistance(P) < Tol) { return; }
Standard_Real V1 = OO1.AngleWithRef(gp_Vec(O1,P),OO1.Crossed(OZ)); Standard_Real V1 = OO1.AngleWithRef(gp_Vec(O1,P),OO1.Crossed(OZ));
if (V1 > -ExtPElS_MyEps && V1 < ExtPElS_MyEps) { V1 = 0.; }
Standard_Real V2 = OO2.AngleWithRef(gp_Vec(P,O2),OO2.Crossed(OZ)); Standard_Real V2 = OO2.AngleWithRef(gp_Vec(P,O2),OO2.Crossed(OZ));
if (V2 > -ExtPElS_MyEps && V2 < ExtPElS_MyEps) { V2 = 0.; }
if (V1 < 0.) { V1 += 2. * M_PI; } if (V1 < 0.) { V1 += 2. * M_PI; }
if (V2 < 0.) { V2 += 2. * M_PI; } if (V2 < 0.) { V2 += 2. * M_PI; }

View File

@@ -336,66 +336,66 @@ Standard_Real ProjectPointOnSurf::LowerDistance() const
} }
// //
if(myApprox2) { if(myApprox2) {
Handle (Geom2d_Curve) C2d; Handle (Geom2d_Curve) C2d;
BuildPCurves(fprm,lprm,Tolpc,myHS2->ChangeSurface().Surface(),newc,C2d); BuildPCurves(fprm,lprm,Tolpc,myHS2->ChangeSurface().Surface(),newc,C2d);
if(Tolpc>myTolReached2d || myTolReached2d==0.) { if(Tolpc>myTolReached2d || myTolReached2d==0.) {
myTolReached2d=Tolpc; myTolReached2d=Tolpc;
} }
// //
slineS2.Append(new Geom2d_TrimmedCurve(C2d,fprm,lprm)); slineS2.Append(new Geom2d_TrimmedCurve(C2d,fprm,lprm));
} }
else { else {
slineS2.Append(H1); slineS2.Append(H1);
} }
} // if (!Precision::IsNegativeInfinite(fprm) && !Precision::IsPositiveInfinite(lprm)) } // if (!Precision::IsNegativeInfinite(fprm) && !Precision::IsPositiveInfinite(lprm))
else { else {
GeomAbs_SurfaceType typS1 = myHS1->Surface().GetType(); GeomAbs_SurfaceType typS1 = myHS1->Surface().GetType();
GeomAbs_SurfaceType typS2 = myHS2->Surface().GetType(); GeomAbs_SurfaceType typS2 = myHS2->Surface().GetType();
if( typS1 == GeomAbs_SurfaceOfExtrusion || if( typS1 == GeomAbs_SurfaceOfExtrusion ||
typS1 == GeomAbs_OffsetSurface || typS1 == GeomAbs_OffsetSurface ||
typS1 == GeomAbs_SurfaceOfRevolution || typS1 == GeomAbs_SurfaceOfRevolution ||
typS2 == GeomAbs_SurfaceOfExtrusion || typS2 == GeomAbs_SurfaceOfExtrusion ||
typS2 == GeomAbs_OffsetSurface || typS2 == GeomAbs_OffsetSurface ||
typS2 == GeomAbs_SurfaceOfRevolution) { typS2 == GeomAbs_SurfaceOfRevolution) {
sline.Append(newc); sline.Append(newc);
slineS1.Append(H1); slineS1.Append(H1);
slineS2.Append(H1); slineS2.Append(H1);
continue; continue;
} }
Standard_Boolean bFNIt, bLPIt; Standard_Boolean bFNIt, bLPIt;
Standard_Real aTestPrm, dT=100.; Standard_Real aTestPrm, dT=100.;
Standard_Real u1, v1, u2, v2, TolX; Standard_Real u1, v1, u2, v2, TolX;
// //
bFNIt=Precision::IsNegativeInfinite(fprm); bFNIt=Precision::IsNegativeInfinite(fprm);
bLPIt=Precision::IsPositiveInfinite(lprm); bLPIt=Precision::IsPositiveInfinite(lprm);
aTestPrm=0.; aTestPrm=0.;
if (bFNIt && !bLPIt) { if (bFNIt && !bLPIt) {
aTestPrm=lprm-dT; aTestPrm=lprm-dT;
} }
else if (!bFNIt && bLPIt) { else if (!bFNIt && bLPIt) {
aTestPrm=fprm+dT; aTestPrm=fprm+dT;
} }
// //
gp_Pnt ptref(newc->Value(aTestPrm)); gp_Pnt ptref(newc->Value(aTestPrm));
// //
TolX = Precision::Confusion(); TolX = Precision::Confusion();
Parameters(myHS1, myHS2, ptref, u1, v1, u2, v2); Parameters(myHS1, myHS2, ptref, u1, v1, u2, v2);
ok = (dom1->Classify(gp_Pnt2d(u1, v1), TolX) != TopAbs_OUT); ok = (dom1->Classify(gp_Pnt2d(u1, v1), TolX) != TopAbs_OUT);
if(ok) { if(ok) {
ok = (dom2->Classify(gp_Pnt2d(u2,v2),TolX) != TopAbs_OUT); ok = (dom2->Classify(gp_Pnt2d(u2,v2),TolX) != TopAbs_OUT);
} }
if (ok) { if (ok) {
sline.Append(newc); sline.Append(newc);
slineS1.Append(H1); slineS1.Append(H1);
slineS2.Append(H1); slineS2.Append(H1);
} }
} }
}// end of for (i=1; i<=myLConstruct.NbParts(); i++) }// end of for (i=1; i<=myLConstruct.NbParts(); i++)
}// case IntPatch_Lin: case IntPatch_Parabola: case IntPatch_Hyperbola: }// case IntPatch_Lin: case IntPatch_Parabola: case IntPatch_Hyperbola:
break; break;
//######################################## //########################################
// Circle and Ellipse // Circle and Ellipse
@@ -952,7 +952,19 @@ Standard_Real ProjectPointOnSurf::LowerDistance() const
mbspc.Degree()); mbspc.Degree());
GeomLib_CheckBSplineCurve Check(BS,myTolCheck,myTolAngCheck); GeomLib_CheckBSplineCurve Check(BS,myTolCheck,myTolAngCheck);
Check.FixTangent(Standard_True,Standard_True); Check.FixTangent(Standard_True,Standard_True);
// //
//Check IsClosed()
Standard_Real aDist = Max(BS->StartPoint().XYZ().SquareModulus(),
BS->EndPoint().XYZ().SquareModulus());
Standard_Real eps = Epsilon(aDist);
if(BS->StartPoint().SquareDistance(BS->EndPoint()) < 2.*eps &&
!BS->IsClosed() && !BS->IsPeriodic())
{
//force Closed()
gp_Pnt aPm((BS->Pole(1).XYZ() + BS->Pole(BS->NbPoles()).XYZ()) / 2.);
BS->SetPole(1, aPm);
BS->SetPole(BS->NbPoles(), aPm);
}
sline.Append(BS); sline.Append(BS);
if(myApprox1) { if(myApprox1) {

View File

@@ -51,6 +51,7 @@
#include <IntPatch_CSFunction.hxx> #include <IntPatch_CSFunction.hxx>
#include <IntPatch_CurvIntSurf.hxx> #include <IntPatch_CurvIntSurf.hxx>
#define myInfinite 1.e15 // the same as was in Adaptor3d_TopolTool
static void Recadre(GeomAbs_SurfaceType typeS1, static void Recadre(GeomAbs_SurfaceType typeS1,
GeomAbs_SurfaceType typeS2, GeomAbs_SurfaceType typeS2,
@@ -558,11 +559,15 @@ void IntPatch_RstInt::PutVertexOnLine (Handle(IntPatch_Line)& L,
Standard_Real VRes = Surf->VResolution(edgeTol); Standard_Real VRes = Surf->VResolution(edgeTol);
IntPatch_HInterTool::Bounds(arc,PFirst,PLast); IntPatch_HInterTool::Bounds(arc,PFirst,PLast);
if (Precision::IsNegativeInfinite(PFirst) || if(Precision::IsNegativeInfinite(PFirst))
Precision::IsPositiveInfinite(PLast)) { PFirst = -myInfinite;
//-- cout<<" IntPatch_RstInt::PutVertexOnLine ---> Restrictions Infinies :"<<endl; if(Precision::IsPositiveInfinite(PLast))
return; PLast = myInfinite;
} //if (Precision::IsNegativeInfinite(PFirst) ||
// Precision::IsPositiveInfinite(PLast)) {
// //-- cout<<" IntPatch_RstInt::PutVertexOnLine ---> Restrictions Infinies :"<<endl;
// return;
//}
Standard_Boolean isVFirst = Standard_False, isVLast = Standard_False; Standard_Boolean isVFirst = Standard_False, isVLast = Standard_False;
gp_Pnt2d p2dFirst,p2dLast; gp_Pnt2d p2dFirst,p2dLast;

View File

@@ -64,17 +64,6 @@ static void BoundedArc (const TheArc& A,
Standard_Boolean& Arcsol, Standard_Boolean& Arcsol,
const Standard_Boolean RecheckOnRegularity); const Standard_Boolean RecheckOnRegularity);
static void InfiniteArc (const TheArc&,
const Handle(TheTopolTool)&,
const Standard_Real,
const Standard_Real,
TheFunction&,
IntStart_SequenceOfPathPoint&,
IntStart_SequenceOfSegment&,
const Standard_Real,
const Standard_Real,
Standard_Boolean&);
static void PointProcess (const gp_Pnt&, static void PointProcess (const gp_Pnt&,
const Standard_Real, const Standard_Real,
const TheArc&, const TheArc&,
@@ -521,7 +510,7 @@ void ComputeBoundsfromInfinite(TheFunction& Func,
// - Inifinies walk. It will take this code // - Inifinies walk. It will take this code
// - With curve surface intersections. // - With curve surface intersections.
NbEchant = 10; NbEchant = 100;
Standard_Real U0 = 0.0; Standard_Real U0 = 0.0;
Standard_Real dU = 0.001; Standard_Real dU = 0.001;
@@ -556,8 +545,8 @@ void ComputeBoundsfromInfinite(TheFunction& Func,
if(Umin>U0) { Umin=U0-10.0; } if(Umin>U0) { Umin=U0-10.0; }
if(Umax<U0) { Umax=U0+10.0; } if(Umax<U0) { Umax=U0+10.0; }
PFin = Umax; PFin = Umax + 10. * (Umax - Umin);
PDeb = Umin; PDeb = Umin - 10. * (Umax - Umin);
} }
else { else {
//-- Possibilite de Arc totalement inclu ds Quad //-- Possibilite de Arc totalement inclu ds Quad
@@ -566,115 +555,6 @@ void ComputeBoundsfromInfinite(TheFunction& Func,
} }
} }
//=======================================================================
//function : InfiniteArc
//purpose :
//=======================================================================
void InfiniteArc (const TheArc& A,
const Handle(TheTopolTool)& Domain,
const Standard_Real Pdeb,
const Standard_Real Pfin,
TheFunction& Func,
IntStart_SequenceOfPathPoint& pnt,
IntStart_SequenceOfSegment& seg,
const Standard_Real TolBoundary,
const Standard_Real TolTangency,
Standard_Boolean& Arcsol)
{
// Find points of solutions and tips bow bow gives a solution.
// The math_FunctionAllRoots function is used. Therefore suitable for
// Beginning of arcs having a point and a closed end point (range
// Parametrage).
Standard_Integer i,Nbi,Nbp;
gp_Pnt ptdeb,ptfin;
Standard_Real pardeb = 0.,parfin = 0.;
Standard_Integer ideb,ifin,range,ranged,rangef;
// Create the Sample Rate (math_FunctionSample or inheriting class)
// Call a math_FunctionAllRoots
Standard_Real EpsX = TheArcTool::Resolution(A,Precision::Confusion());
//@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
//@@@ Tolerance is the asociee al arc (Inconsistency with tracking)
//@@@ (EPSX ~ 1e-5 and ResolutionU and V ~ 1e-9)
//@@@ Vertex is here is not found as a point to stop
//@@@ Wayline
//@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
EpsX = 0.0000000001;
//@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
//@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
//@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
Standard_Integer NbEchant = Func.NbSamples();
//-- Modif 24 Aout 93 -----------------------------
Standard_Real nTolTangency = TolTangency;
if((Pfin - Pdeb) < (TolTangency*10.0)) {
nTolTangency=(Pfin-Pdeb)*0.1;
}
if(EpsX>(nTolTangency+nTolTangency)) {
EpsX = nTolTangency * 0.1;
}
//--------------------------------------------------
// - Plant with a edge with 2 Samples
// - Whose ends are solutions (f = 0)
// - And the derivative is zero or
// - Example: a diameter of a sphere segment
if(NbEchant<3) NbEchant = 3; //-- lbr 19.04.95
//--------------------------------------------------
Standard_Real PDeb = Pdeb;
Standard_Real PFin = Pfin;
ComputeBoundsfromInfinite(Func,PDeb,PFin,NbEchant);
math_FunctionSample Echant(PDeb,PFin,NbEchant);
math_FunctionAllRoots Sol(Func,Echant,EpsX,TolBoundary,nTolTangency);
if (!Sol.IsDone()) {Standard_Failure::Raise();}
Nbp=Sol.NbPoints();
for (i=1; i<=Nbp; i++) {
Standard_Real para = Sol.GetPoint(i);
Standard_Real dist;
if(Func.Value(para,dist)) {
PointProcess(Func.Valpoint(Sol.GetPointState(i)),Sol.GetPoint(i),
A,Domain,pnt,TolBoundary,range);
}
}
// For each interval:
// Process the ends as points
// Add range in the list of segments
Nbi=Sol.NbIntervals();
for (i=1; i<=Nbi; i++) {
IntStart_TheSegment newseg;
newseg.SetValue(A);
// Recover start and end points, and parameter.
Sol.GetInterval(i,pardeb,parfin);
Sol.GetIntervalState(i,ideb,ifin);
ptdeb=Func.Valpoint(ideb);
ptfin=Func.Valpoint(ifin);
PointProcess(ptdeb,pardeb,A,Domain,pnt,TolBoundary,ranged);
newseg.SetLimitPoint(pnt.Value(ranged),Standard_True);
PointProcess(ptfin,parfin,A,Domain,pnt,TolBoundary,rangef);
newseg.SetLimitPoint(pnt.Value(rangef),Standard_False);
seg.Append(newseg);
}
Arcsol=Standard_False;
if (Nbi==1) {
if (pardeb == Pdeb && parfin == Pfin) {
Arcsol=Standard_True;
}
}
}
//======================================================================= //=======================================================================
//function : PointProcess //function : PointProcess
//purpose : //purpose :
@@ -955,11 +835,10 @@ IntStart_SearchOnBoundaries::IntStart_SearchOnBoundaries ()
TheSOBTool::Bounds(A,PDeb,PFin); TheSOBTool::Bounds(A,PDeb,PFin);
if(Precision::IsNegativeInfinite(PDeb) || if(Precision::IsNegativeInfinite(PDeb) ||
Precision::IsPositiveInfinite(PFin)) { Precision::IsPositiveInfinite(PFin)) {
InfiniteArc(A,Domain,PDeb,PFin,Func,spnt,sseg,TolBoundary,TolTangency,Arcsol); Standard_Integer NbEchant;
} ComputeBoundsfromInfinite(Func,PDeb,PFin,NbEchant);
else {
BoundedArc(A,Domain,PDeb,PFin,Func,spnt,sseg,TolBoundary,TolTangency,Arcsol,RecheckOnRegularity);
} }
BoundedArc(A,Domain,PDeb,PFin,Func,spnt,sseg,TolBoundary,TolTangency,Arcsol,RecheckOnRegularity);
all = (all && Arcsol); all = (all && Arcsol);
} }

View File

@@ -41,7 +41,7 @@ prism crg w 0 -$y*2 0 inf
bsection result crg cyl bsection result crg cyl
set length 15.1392 set length 110.093
set 2dviewer 0 set 2dviewer 0
# checksection res # checksection res

69
tests/bugs/moddata_3/bug23981 Executable file
View File

@@ -0,0 +1,69 @@
puts "========="
puts "CR23981"
puts "========="
puts ""
###############################
## Wrong section curves
###############################
restore [locate_data_file bug23981_s1.draw] s1
restore [locate_data_file bug23981_s2.draw] s2
intersect i s1 s2
puts "First test"
dlog reset
dlog on
xdistcs i_1 s1 0 1 100
set Log1 [dlog get]
set List1 [split ${Log1} {TD= \t\n}]
set L1 [llength ${List1}]
set L2 10
set L3 5
set N [expr (${L1} - ${L2})/${L3} + 1]
set Tolerance 1.0e-5
set D_good 0.
for {set i 1} {${i} <= ${N}} {incr i} {
set j1 [expr ${L2} + (${i}-1)*${L3}]
set j2 [expr ${j1} + 2]
set T [lindex ${List1} ${j1}]
set D [lindex ${List1} ${j2}]
#puts "i=${i} j1=${j1} j2=${j2} T=${T} D=${D}"
if { [expr abs(${D} - ${D_good})] > ${Tolerance} } {
puts "Error: i=${i} T=${T} D=${D}"
}
}
puts "Second test"
dlog reset
dlog on
xdistcs i_2 s1 0 1 100
set Log2 [dlog get]
set List2 [split ${Log2} {TD= \t\n}]
set L1 [llength ${List2}]
set L2 10
set L3 5
set N [expr (${L1} - ${L2})/${L3} + 1]
set Tolerance 1.0e-5
set D_good 0.
for {set i 1} {${i} <= ${N}} {incr i} {
set j1 [expr ${L2} + (${i}-1)*${L3}]
set j2 [expr ${j1} + 2]
set T [lindex ${List2} ${j1}]
set D [lindex ${List2} ${j2}]
#puts "i=${i} j1=${j1} j2=${j2} T=${T} D=${D}"
if { [expr abs(${D} - ${D_good})] > ${Tolerance} } {
puts "Error: i=${i} T=${T} D=${D}"
}
}
smallview
fit
set only_screen_axo 1