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

0027299: Incorrect result of the normal projection algorithm

Curve splitting is added to handle seam passing by initial curve.
test cases are added.

Minor corrections.
This commit is contained in:
aml 2016-03-25 13:10:12 +03:00 committed by bugmaster
parent 58e14d59e7
commit 5333268def
8 changed files with 408 additions and 55 deletions

View File

@ -194,12 +194,21 @@ void Extrema_GenExtCC::Perform()
Curve2 &C2 = *(Curve2*)myC[1]; Curve2 &C2 = *(Curve2*)myC[1];
Standard_Integer aNbInter[2]; Standard_Integer aNbInter[2];
aNbInter[0] = C1.NbIntervals(GeomAbs_C2); GeomAbs_Shape aContinuity = GeomAbs_C2;
aNbInter[1] = C2.NbIntervals(GeomAbs_C2); aNbInter[0] = C1.NbIntervals(aContinuity);
aNbInter[1] = C2.NbIntervals(aContinuity);
if (aNbInter[0] * aNbInter[1] > 100)
{
aContinuity = GeomAbs_C1;
aNbInter[0] = C1.NbIntervals(aContinuity);
aNbInter[1] = C2.NbIntervals(aContinuity);
}
TColStd_Array1OfReal anIntervals1(1, aNbInter[0] + 1); TColStd_Array1OfReal anIntervals1(1, aNbInter[0] + 1);
TColStd_Array1OfReal anIntervals2(1, aNbInter[1] + 1); TColStd_Array1OfReal anIntervals2(1, aNbInter[1] + 1);
C1.Intervals(anIntervals1, GeomAbs_C2); C1.Intervals(anIntervals1, aContinuity);
C2.Intervals(anIntervals2, GeomAbs_C2); C2.Intervals(anIntervals2, aContinuity);
// Lipchitz constant approximation. // Lipchitz constant approximation.
Standard_Real aLC = 9.0; // Default value. Standard_Real aLC = 9.0; // Default value.
@ -226,6 +235,7 @@ void Extrema_GenExtCC::Perform()
Extrema_GlobOptFuncCCC2 aFunc (C1, C2); Extrema_GlobOptFuncCCC2 aFunc (C1, C2);
math_GlobOptMin aFinder(&aFunc, myLowBorder, myUppBorder, aLC); math_GlobOptMin aFinder(&aFunc, myLowBorder, myUppBorder, aLC);
aFinder.SetLipConstState(isConstLockedFlag); aFinder.SetLipConstState(isConstLockedFlag);
aFinder.SetContinuity(aContinuity == GeomAbs_C2 ? 2 : 1);
Standard_Real aDiscTol = 1.0e-2; Standard_Real aDiscTol = 1.0e-2;
Standard_Real aValueTol = 1.0e-2; Standard_Real aValueTol = 1.0e-2;
Standard_Real aSameTol = myCurveMinTol / (aDiscTol); Standard_Real aSameTol = myCurveMinTol / (aDiscTol);

View File

@ -15,6 +15,8 @@
// commercial license or contractual agreement. // commercial license or contractual agreement.
#include <algorithm>
#include <Adaptor2d_HCurve2d.hxx> #include <Adaptor2d_HCurve2d.hxx>
#include <Adaptor3d_HCurve.hxx> #include <Adaptor3d_HCurve.hxx>
#include <Adaptor3d_HSurface.hxx> #include <Adaptor3d_HSurface.hxx>
@ -38,6 +40,11 @@
#include <Standard_NotImplemented.hxx> #include <Standard_NotImplemented.hxx>
#include <Standard_OutOfRange.hxx> #include <Standard_OutOfRange.hxx>
#include <TColgp_HSequenceOfPnt.hxx> #include <TColgp_HSequenceOfPnt.hxx>
#include <Adaptor3d_CurveOnSurface.hxx>
#include <Geom2d_Line.hxx>
#include <Geom2dAdaptor_HCurve.hxx>
#include <Extrema_ExtCC.hxx>
#include <NCollection_Vector.hxx>
#define FuncTol 1.e-10 #define FuncTol 1.e-10
@ -64,6 +71,59 @@ static void ResultChron( OSD_Chronometer & ch, Standard_Real & time)
} }
#endif #endif
// Structure to perform splits computation.
// This structure is not thread-safe since operations under mySplits should be performed in a critical section.
// myPeriodicDir - 0 for U periodicity and 1 for V periodicity.
struct SplitDS
{
SplitDS(const Handle(Adaptor3d_HCurve) &theCurve,
const Handle(Adaptor3d_HSurface) &theSurface,
NCollection_Vector<Standard_Real> &theSplits)
: myCurve(theCurve),
mySurface(theSurface),
mySplits(theSplits)
{ }
// Assignment operator is forbidden.
void operator=(const SplitDS &theSplitDS);
const Handle(Adaptor3d_HCurve) myCurve;
const Handle(Adaptor3d_HSurface) mySurface;
NCollection_Vector<Standard_Real> &mySplits;
Standard_Real myPerMinParam;
Standard_Real myPerMaxParam;
Standard_Integer myPeriodicDir;
Extrema_ExtCC *myExtCC;
Extrema_ExtPS *myExtPS;
};
//! Compute split points in the parameter space of the curve.
static void BuildCurveSplits(const Handle(Adaptor3d_HCurve) &theCurve,
const Handle(Adaptor3d_HSurface) &theSurface,
const Standard_Real theTolU,
const Standard_Real theTolV,
NCollection_Vector<Standard_Real> &theSplits);
//! Perform splitting on a specified direction. Sub-method in BuildCurveSplits.
static void SplitOnDirection(SplitDS & theSplitDS);
//! Perform recursive search of the split points.
static void FindSplitPoint(SplitDS & theSplitDS,
const Standard_Real theMinParam,
const Standard_Real theMaxParam);
//=======================================================================
//function : Comparator
//purpose : used in sort algorithm
//=======================================================================
inline Standard_Boolean Comparator(const Standard_Real theA,
const Standard_Real theB)
{
return theA < theB;
}
//======================================================================= //=======================================================================
//function : d1 //function : d1
@ -564,34 +624,37 @@ ProjLib_CompProjectedCurve::ProjLib_CompProjectedCurve
void ProjLib_CompProjectedCurve::Init() void ProjLib_CompProjectedCurve::Init()
{ {
myTabInt.Nullify(); myTabInt.Nullify();
NCollection_Vector<Standard_Real> aSplits;
aSplits.Clear();
Standard_Real Tol;// Tolerance for ExactBound Standard_Real Tol;// Tolerance for ExactBound
Standard_Integer i, Nend = 0; Standard_Integer i, Nend = 0, aSplitIdx = 0;
Standard_Boolean FromLastU=Standard_False; Standard_Boolean FromLastU = Standard_False,
isSplitsComputed = Standard_False;
//new part (to discard far solutions) const Standard_Real aTol3D = Precision::Confusion();
Standard_Real TolC = Precision::Confusion(), TolS = Precision::Confusion(); Extrema_ExtCS CExt(myCurve->Curve(), mySurface->Surface(), aTol3D, aTol3D);
Extrema_ExtCS CExt(myCurve->Curve(),
mySurface->Surface(),
TolC,
TolS);
if (CExt.IsDone() && CExt.NbExt()) if (CExt.IsDone() && CExt.NbExt())
{ {
// Search for the minimum solution // Search for the minimum solution.
Nend = CExt.NbExt(); // Avoid usage of extrema result that can be wrong for extrusion.
if(myMaxDist > 0 && if(myMaxDist > 0 &&
// Avoid usage of extrema result that can be wrong for extrusion
mySurface->GetType() != GeomAbs_SurfaceOfExtrusion) mySurface->GetType() != GeomAbs_SurfaceOfExtrusion)
{ {
Standard_Real min_val2; Standard_Real min_val2;
min_val2 = CExt.SquareDistance(1); min_val2 = CExt.SquareDistance(1);
Nend = CExt.NbExt();
for(i = 2; i <= Nend; i++) for(i = 2; i <= Nend; i++)
if (CExt.SquareDistance(i) < min_val2) min_val2 = CExt.SquareDistance(i); {
if (CExt.SquareDistance(i) < min_val2)
min_val2 = CExt.SquareDistance(i);
}
if (min_val2 > myMaxDist * myMaxDist) if (min_val2 > myMaxDist * myMaxDist)
return; return; // No near solution -> exit.
} }
} }
// end of new part
Standard_Real FirstU, LastU, Step, SearchStep, WalkStep, t; Standard_Real FirstU, LastU, Step, SearchStep, WalkStep, t;
@ -605,12 +668,9 @@ void ProjLib_CompProjectedCurve::Init()
SearchStep = 10*MinStep; SearchStep = 10*MinStep;
Step = SearchStep; Step = SearchStep;
//Initialization of aPrjPS gp_Pnt2d aLowBorder(mySurface->FirstUParameter(),mySurface->FirstVParameter());
Standard_Real Uinf = mySurface->FirstUParameter(); gp_Pnt2d aUppBorder(mySurface->LastUParameter(), mySurface->LastVParameter());
Standard_Real Usup = mySurface->LastUParameter(); gp_Pnt2d aTol(myTolU, myTolV);
Standard_Real Vinf = mySurface->FirstVParameter();
Standard_Real Vsup = mySurface->LastVParameter();
ProjLib_PrjResolve aPrjPS(myCurve->Curve(), mySurface->Surface(), 1); ProjLib_PrjResolve aPrjPS(myCurve->Curve(), mySurface->Surface(), 1);
t = FirstU; t = FirstU;
@ -648,10 +708,8 @@ void ProjLib_CompProjectedCurve::Init()
ParT=P1.Parameter(); ParT=P1.Parameter();
P2.Parameter(ParU, ParV); P2.Parameter(ParU, ParV);
aPrjPS.Perform(ParT, ParU, ParV, gp_Pnt2d(myTolU, myTolV), aPrjPS.Perform(ParT, ParU, ParV, aTol, aLowBorder, aUppBorder, FuncTol, Standard_True);
gp_Pnt2d(mySurface->FirstUParameter(),mySurface->FirstVParameter()),
gp_Pnt2d(mySurface->LastUParameter(), mySurface->LastVParameter()),
FuncTol, Standard_True);
if ( aPrjPS.IsDone() && P1.Parameter() > Max(FirstU,t-Step+Precision::PConfusion()) if ( aPrjPS.IsDone() && P1.Parameter() > Max(FirstU,t-Step+Precision::PConfusion())
&& P1.Parameter() <= t) && P1.Parameter() <= t)
{ {
@ -684,38 +742,37 @@ void ProjLib_CompProjectedCurve::Init()
gp_Vec2d D; gp_Vec2d D;
if ((mySurface->IsUPeriodic() && if ((mySurface->IsUPeriodic() &&
Abs(Usup - Uinf - mySurface->UPeriod()) < Precision::Confusion()) || Abs(aUppBorder.X() - aLowBorder.X() - mySurface->UPeriod()) < Precision::Confusion()) ||
(mySurface->IsVPeriodic() && (mySurface->IsVPeriodic() &&
Abs(Vsup - Vinf - mySurface->VPeriod()) < Precision::Confusion())) Abs(aUppBorder.Y() - aLowBorder.Y() - mySurface->VPeriod()) < Precision::Confusion()))
{ {
if((Abs(U - Uinf) < mySurface->UResolution(Precision::PConfusion())) && if((Abs(U - aLowBorder.X()) < mySurface->UResolution(Precision::PConfusion())) &&
mySurface->IsUPeriodic()) mySurface->IsUPeriodic())
{ {
d1(t, U, V, D, myCurve, mySurface); d1(t, U, V, D, myCurve, mySurface);
if (D.X() < 0 ) U = Usup; if (D.X() < 0 ) U = aUppBorder.X();
} }
else if((Abs(U - Usup) < mySurface->UResolution(Precision::PConfusion())) && else if((Abs(U - aUppBorder.X()) < mySurface->UResolution(Precision::PConfusion())) &&
mySurface->IsUPeriodic()) mySurface->IsUPeriodic())
{ {
d1(t, U, V, D, myCurve, mySurface); d1(t, U, V, D, myCurve, mySurface);
if (D.X() > 0) U = Uinf; if (D.X() > 0) U = aLowBorder.X();
} }
if((Abs(V - Vinf) < mySurface->VResolution(Precision::PConfusion())) && if((Abs(V - aLowBorder.Y()) < mySurface->VResolution(Precision::PConfusion())) &&
mySurface->IsVPeriodic()) mySurface->IsVPeriodic())
{ {
d1(t, U, V, D, myCurve, mySurface); d1(t, U, V, D, myCurve, mySurface);
if (D.Y() < 0) V = Vsup; if (D.Y() < 0) V = aUppBorder.Y();
} }
else if((Abs(V - Vsup) <= mySurface->VResolution(Precision::PConfusion())) && else if((Abs(V - aUppBorder.Y()) <= mySurface->VResolution(Precision::PConfusion())) &&
mySurface->IsVPeriodic()) mySurface->IsVPeriodic())
{ {
d1(t, U, V, D, myCurve, mySurface); d1(t, U, V, D, myCurve, mySurface);
if (D.Y() > 0) V = Vinf; if (D.Y() > 0) V = aLowBorder.Y();
} }
} }
if (myMaxDist > 0) if (myMaxDist > 0)
{ {
// Here we are going to stop if the distance between projection and // Here we are going to stop if the distance between projection and
@ -765,7 +822,6 @@ void ProjLib_CompProjectedCurve::Init()
} }
if (!new_part) break; if (!new_part) break;
//We have found a new continuous part //We have found a new continuous part
Handle(TColgp_HSequenceOfPnt) hSeq = new TColgp_HSequenceOfPnt(); Handle(TColgp_HSequenceOfPnt) hSeq = new TColgp_HSequenceOfPnt();
mySequence->Append(hSeq); mySequence->Append(hSeq);
@ -790,9 +846,7 @@ void ProjLib_CompProjectedCurve::Init()
if (t > LastU) t = LastU; if (t > LastU) t = LastU;
Standard_Real prevStep = Step; Standard_Real prevStep = Step;
Standard_Real U0, V0; Standard_Real U0, V0;
gp_Pnt2d aLowBorder(mySurface->FirstUParameter(),mySurface->FirstVParameter());
gp_Pnt2d aUppBorder(mySurface->LastUParameter(), mySurface->LastVParameter());
gp_Pnt2d aTol(myTolU, myTolV);
//Here we are trying to prolong continuous part //Here we are trying to prolong continuous part
while (t <= LastU && new_part) while (t <= LastU && new_part)
{ {
@ -865,6 +919,33 @@ void ProjLib_CompProjectedCurve::Init()
// Check for possible local traps. // Check for possible local traps.
UpdateTripleByTrapCriteria(Triple); UpdateTripleByTrapCriteria(Triple);
// Protection from case when the whole curve lies on a seam.
if (!isSplitsComputed)
{
Standard_Boolean isUPossible = Standard_False;
if (mySurface->IsUPeriodic() &&
(Abs(Triple.Y() - mySurface->FirstUParameter() ) > Precision::PConfusion() &&
Abs(Triple.Y() - mySurface->LastUParameter() ) > Precision::PConfusion()))
{
isUPossible = Standard_True;
}
Standard_Boolean isVPossible = Standard_False;
if (mySurface->IsVPeriodic() &&
(Abs(Triple.Z() - mySurface->FirstVParameter() ) > Precision::PConfusion() &&
Abs(Triple.Z() - mySurface->LastVParameter() ) > Precision::PConfusion()))
{
isVPossible = Standard_True;
}
if (isUPossible || isVPossible)
{
// When point is good conditioned.
BuildCurveSplits(myCurve, mySurface, myTolU, myTolV, aSplits);
isSplitsComputed = Standard_True;
}
}
if((Triple.X() - mySequence->Value(myNbCurves)->Value(mySequence->Value(myNbCurves)->Length()).X()) > 1.e-10) if((Triple.X() - mySequence->Value(myNbCurves)->Value(mySequence->Value(myNbCurves)->Length()).X()) > 1.e-10)
mySequence->Value(myNbCurves)->Append(Triple); mySequence->Value(myNbCurves)->Append(Triple);
if (t == LastU) {t = LastU + 1; break;}//return; if (t == LastU) {t = LastU + 1; break;}//return;
@ -882,9 +963,32 @@ void ProjLib_CompProjectedCurve::Init()
Step = Step + LastU - t; Step = Step + LastU - t;
t = LastU; t = LastU;
} }
// We assume at least one point of cache inside of a split.
const Standard_Integer aSize = aSplits.Size();
for(Standard_Integer anIdx = aSplitIdx; anIdx < aSize; ++anIdx)
{
const Standard_Real aParam = aSplits(anIdx);
if (Abs(aParam - Triple.X() ) < Precision::PConfusion())
{
// The current point is equal to a split point.
new_part = Standard_False;
// Move split index to avoid check of the whole list.
++aSplitIdx;
break;
}
else if (aParam < t + Precision::PConfusion() )
{
// The next point crosses the split point.
t = aParam;
Step = t - prevTriple.X();
}
} // for(Standard_Integer anIdx = aSplitIdx; anIdx < aSize; ++anIdx)
} }
} }
} }
// Sequence post-proceeding. // Sequence post-proceeding.
Standard_Integer j; Standard_Integer j;
@ -1667,3 +1771,160 @@ void ProjLib_CompProjectedCurve::UpdateTripleByTrapCriteria(gp_Pnt &thePoint) co
thePoint.SetY(U); thePoint.SetY(U);
thePoint.SetZ(V); thePoint.SetZ(V);
} }
//=======================================================================
//function : BuildCurveSplits
//purpose :
//=======================================================================
void BuildCurveSplits(const Handle(Adaptor3d_HCurve) &theCurve,
const Handle(Adaptor3d_HSurface) &theSurface,
const Standard_Real theTolU,
const Standard_Real theTolV,
NCollection_Vector<Standard_Real> &theSplits)
{
SplitDS aDS(theCurve, theSurface, theSplits);
Extrema_ExtPS anExtPS;
anExtPS.Initialize(theSurface->Surface(),
theSurface->FirstUParameter(), theSurface->LastUParameter(),
theSurface->FirstVParameter(), theSurface->LastVParameter(),
theTolU, theTolV);
aDS.myExtPS = &anExtPS;
if (theSurface->IsUPeriodic())
{
aDS.myPeriodicDir = 0;
SplitOnDirection(aDS);
}
if (theSurface->IsVPeriodic())
{
aDS.myPeriodicDir = 1;
SplitOnDirection(aDS);
}
std::sort(aDS.mySplits.begin(), aDS.mySplits.end(), Comparator);
}
//=======================================================================
//function : SplitOnDirection
//purpose : This method compute points in the parameter space of the curve
// on which curve should be split since period jump is happen.
//=======================================================================
void SplitOnDirection(SplitDS & theSplitDS)
{
// Algorithm:
// Create 3D curve which is correspond to the periodic bound in 2d space.
// Run curve / curve extrema and run extrema point / surface to check that
// the point will be projected to the periodic bound.
// In this method assumed that the points cannot be closer to each other that 1% of the parameter space.
gp_Pnt2d aStartPnt(theSplitDS.mySurface->FirstUParameter(), theSplitDS.mySurface->FirstVParameter());
gp_Dir2d aDir(theSplitDS.myPeriodicDir, (Standard_Integer)!theSplitDS.myPeriodicDir);
theSplitDS.myPerMinParam = !theSplitDS.myPeriodicDir ? theSplitDS.mySurface->FirstUParameter():
theSplitDS.mySurface->FirstVParameter();
theSplitDS.myPerMaxParam = !theSplitDS.myPeriodicDir ? theSplitDS.mySurface->LastUParameter():
theSplitDS.mySurface->LastVParameter();
Standard_Real aLast2DParam = theSplitDS.myPeriodicDir ?
theSplitDS.mySurface->LastUParameter() - theSplitDS.mySurface->FirstUParameter():
theSplitDS.mySurface->LastVParameter() - theSplitDS.mySurface->FirstVParameter();
// Create line which is represent periodic border.
Handle(Geom2d_Curve) aC2GC = new Geom2d_Line(aStartPnt, aDir);
Handle(Geom2dAdaptor_HCurve) aC = new Geom2dAdaptor_HCurve(aC2GC, 0, aLast2DParam);
Adaptor3d_CurveOnSurface aCOnS(aC, theSplitDS.mySurface);
Extrema_ExtCC anExtCC;
anExtCC.SetCurve(1, aCOnS);
anExtCC.SetCurve(2, theSplitDS.myCurve->Curve());
anExtCC.SetSingleSolutionFlag(Standard_True); // Search only one solution since multiple invocations are needed.
anExtCC.SetRange(1, 0, aLast2DParam);
theSplitDS.myExtCC = &anExtCC;
FindSplitPoint(theSplitDS,
theSplitDS.myCurve->FirstParameter(), // Initial curve range.
theSplitDS.myCurve->LastParameter());
}
//=======================================================================
//function : FindSplitPoint
//purpose :
//=======================================================================
void FindSplitPoint(SplitDS &theSplitDS,
const Standard_Real theMinParam,
const Standard_Real theMaxParam)
{
// Make extrema copy to avoid dependencies between different levels of the recursion.
Extrema_ExtCC anExtCC(*theSplitDS.myExtCC);
anExtCC.SetRange(2, theMinParam, theMaxParam);
anExtCC.Perform();
if (anExtCC.IsDone())
{
const Standard_Integer aNbExt = anExtCC.NbExt();
for (Standard_Integer anIdx = 1; anIdx <= aNbExt; ++anIdx)
{
Extrema_POnCurv aPOnC1, aPOnC2;
anExtCC.Points(anIdx, aPOnC1, aPOnC2);
theSplitDS.myExtPS->Perform(aPOnC2.Value());
if (!theSplitDS.myExtPS->IsDone())
return;
// Find point with the minimal Euclidean distance to avoid
// false positive points detection.
Standard_Integer aMinIdx = -1;
Standard_Real aMinSqDist = RealLast();
const Standard_Integer aNbPext = theSplitDS.myExtPS->NbExt();
for(Standard_Integer aPIdx = 1; aPIdx <= aNbPext; ++aPIdx)
{
const Standard_Real aCurrSqDist = theSplitDS.myExtPS->SquareDistance(aPIdx);
if (aCurrSqDist < aMinSqDist)
{
aMinSqDist = aCurrSqDist;
aMinIdx = aPIdx;
}
}
// Check that is point will be projected to the periodic border.
const Extrema_POnSurf &aPOnS = theSplitDS.myExtPS->Point(aMinIdx);
Standard_Real U, V, aProjParam;
aPOnS.Parameter(U, V);
aProjParam = theSplitDS.myPeriodicDir ? V : U;
if (Abs(aProjParam - theSplitDS.myPerMinParam) < Precision::PConfusion() ||
Abs(aProjParam - theSplitDS.myPerMaxParam) < Precision::PConfusion() )
{
const Standard_Real aParam = aPOnC2.Parameter();
const Standard_Real aCFParam = theSplitDS.myCurve->FirstParameter();
const Standard_Real aCLParam = theSplitDS.myCurve->LastParameter();
if (aParam > aCFParam + Precision::PConfusion() &&
aParam < aCLParam - Precision::PConfusion() )
{
// Add only inner points.
theSplitDS.mySplits.Append(aParam);
}
const Standard_Real aDeltaCoeff = 0.01;
const Standard_Real aDelta = (theMaxParam - theMinParam +
aCLParam - aCFParam) * aDeltaCoeff;
if (aParam - aDelta > theMinParam + Precision::PConfusion())
{
FindSplitPoint(theSplitDS,
theMinParam, aParam - aDelta); // Curve parameters.
}
if (aParam + aDelta < theMaxParam - Precision::PConfusion())
{
FindSplitPoint(theSplitDS,
aParam + aDelta, theMaxParam); // Curve parameters.
}
}
} // for (Standard_Integer anIdx = 1; anIdx <= aNbExt; ++anIdx)
}
}

View File

@ -141,7 +141,7 @@ public:
//! Returns the parameters corresponding to //! Returns the parameters corresponding to
//! S discontinuities. //! S discontinuities.
//! //!
//! The array must provide enough room to accomodate //! The array must provide enough room to accommodate
//! for the parameters. i.e. T.Length() > NbIntervals() //! for the parameters. i.e. T.Length() > NbIntervals()
Standard_EXPORT void Intervals (TColStd_Array1OfReal& T, const GeomAbs_Shape S) const Standard_OVERRIDE; Standard_EXPORT void Intervals (TColStd_Array1OfReal& T, const GeomAbs_Shape S) const Standard_OVERRIDE;
@ -167,7 +167,7 @@ protected:
private: private:
// This method performs check possibility of optimization traps and tries to go out from them. //! This method performs check possibility of optimization traps and tries to go out from them.
//@return thePoint - input / corrected point. //@return thePoint - input / corrected point.
Standard_EXPORT void UpdateTripleByTrapCriteria(gp_Pnt &thePoint) const; Standard_EXPORT void UpdateTripleByTrapCriteria(gp_Pnt &thePoint) const;
@ -186,8 +186,6 @@ private:
Handle(TColStd_HArray1OfBoolean) mySnglPnts; Handle(TColStd_HArray1OfBoolean) mySnglPnts;
Handle(TColStd_HArray1OfReal) myMaxDistance; Handle(TColStd_HArray1OfReal) myMaxDistance;
Handle(TColStd_HArray1OfReal) myTabInt; Handle(TColStd_HArray1OfReal) myTabInt;
}; };

View File

@ -48,7 +48,8 @@ math_GlobOptMin::math_GlobOptMin(math_MultipleVarFunction* theFunc,
myMaxV(1, myN), myMaxV(1, myN),
myExpandCoeff(1, myN), myExpandCoeff(1, myN),
myCellSize(0, myN - 1), myCellSize(0, myN - 1),
myFilter(theFunc->NbVariables()) myFilter(theFunc->NbVariables()),
myCont(2)
{ {
Standard_Integer i; Standard_Integer i;
@ -278,7 +279,8 @@ Standard_Boolean math_GlobOptMin::computeLocalExtremum(const math_Vector& thePnt
Standard_Integer i; Standard_Integer i;
//Newton method //Newton method
if (dynamic_cast<math_MultipleVarFunctionWithHessian*>(myFunc)) if (myCont >= 2 &&
dynamic_cast<math_MultipleVarFunctionWithHessian*>(myFunc))
{ {
math_MultipleVarFunctionWithHessian* aTmp = math_MultipleVarFunctionWithHessian* aTmp =
dynamic_cast<math_MultipleVarFunctionWithHessian*> (myFunc); dynamic_cast<math_MultipleVarFunctionWithHessian*> (myFunc);
@ -295,7 +297,8 @@ Standard_Boolean math_GlobOptMin::computeLocalExtremum(const math_Vector& thePnt
} else } else
// BFGS method used. // BFGS method used.
if (dynamic_cast<math_MultipleVarFunctionWithGradient*>(myFunc)) if (myCont >= 1 &&
dynamic_cast<math_MultipleVarFunctionWithGradient*>(myFunc))
{ {
math_MultipleVarFunctionWithGradient* aTmp = math_MultipleVarFunctionWithGradient* aTmp =
dynamic_cast<math_MultipleVarFunctionWithGradient*> (myFunc); dynamic_cast<math_MultipleVarFunctionWithGradient*> (myFunc);

View File

@ -152,6 +152,12 @@ public:
//! Get functional minimal value. //! Get functional minimal value.
Standard_EXPORT Standard_Real GetFunctionalMinimalValue(); Standard_EXPORT Standard_Real GetFunctionalMinimalValue();
//! Get continuity of local borders splits.
inline Standard_Integer GetContinuity() const {return myCont; }
//! Set continuity of local borders splits.
inline void SetContinuity(const Standard_Integer theCont) { myCont = theCont; }
private: private:
// Compute cell size. // Compute cell size.
@ -223,6 +229,9 @@ private:
Standard_Boolean isFirstCellFilterInvoke; Standard_Boolean isFirstCellFilterInvoke;
NCollection_CellFilter<NCollection_CellFilter_Inspector> myFilter; NCollection_CellFilter<NCollection_CellFilter_Inspector> myFilter;
// Continuity of local borders.
Standard_Integer myCont;
Standard_Real myF; // Current value of Global optimum. Standard_Real myF; // Current value of Global optimum.
}; };

View File

@ -0,0 +1,23 @@
puts "================"
puts "0027299"
puts "================"
puts ""
##############################################################
# Incorrect result of the normal projection algorithm
# Exception during the exectuion
##############################################################
restore [locate_data_file bug27299_1.brep] aShape
explode aShape
nproject result aShape_1 aShape_2
# Result length check.
checkprops result -l 800
# Visual check.
donly result
smallview
fit
display aShape_2
checkview -screenshot -2d -path ${imagedir}/${test_image}.png

View File

@ -0,0 +1,26 @@
puts "================"
puts "0027299"
puts "================"
puts ""
##############################################################
# Incorrect result of the normal projection algorithm
# Exception during the exectuion
##############################################################
restore [locate_data_file bug27299_1.brep] aShape
explode aShape
# To make task non-symmetry.
ttranslate aShape_1 0 0 135.123
nproject result aShape_1 aShape_2
# Result length check.
checkprops result -l 800
# Visual check.
donly result
smallview
fit
display aShape_2
checkview -screenshot -2d -path ${imagedir}/${test_image}.png

View File

@ -0,0 +1,23 @@
puts "================"
puts "0027299"
puts "================"
puts ""
##############################################################
# Incorrect result of the normal projection algorithm
# Exception during the exectuion
##############################################################
restore [locate_data_file bug27299_2.brep] aShape
explode aShape
nproject result aShape_1 aShape_2
# Result length check.
checkprops result -l 10086
# Visual check.
donly result
smallview
fit
display aShape_2
checkview -screenshot -2d -path ${imagedir}/${test_image}.png