mirror of
https://git.dev.opencascade.org/repos/occt.git
synced 2025-04-04 18:06:22 +03:00
0026184: GeomAPI_ExtremaCurveCurve hangs on parallel b-spline curves
Class CellFilterNDim added. Class CellFulterNDim used in GlobOptMin to improve performance in case of many solutions. Memory leak eliminated. Test cases for issue CR26184 Small correction of test cases for issue CR26184
This commit is contained in:
parent
f0e3a4bac7
commit
4b65fc7750
@ -70,6 +70,7 @@ NCollection_SparseArray.hxx
|
||||
NCollection_SparseArrayBase.hxx
|
||||
NCollection_SparseArrayBase.cxx
|
||||
NCollection_CellFilter.hxx
|
||||
NCollection_CellFilterNDim.hxx
|
||||
NCollection_Handle.hxx
|
||||
NCollection_Handle.cxx
|
||||
|
||||
|
404
src/NCollection/NCollection_CellFilterNDim.hxx
Normal file
404
src/NCollection/NCollection_CellFilterNDim.hxx
Normal file
@ -0,0 +1,404 @@
|
||||
// Created on: 2015-06-17
|
||||
// Created by: Alexander Malyshev
|
||||
// Copyright (c) 2007-2015 OPEN CASCADE SAS
|
||||
//
|
||||
// This file is part of Open CASCADE Technology software library.
|
||||
//
|
||||
// This library is free software; you can redistribute it and/or modify it under
|
||||
// the terms of the GNU Lesser General Public License version 2.1 as published
|
||||
// by the Free Software Foundation, with special exception defined in the file
|
||||
// OCCT_LGPL_EXCEPTION.txt. Consult the file LICENSE_LGPL_21.txt included in OCCT
|
||||
// distribution for complete text of the license and disclaimer of any warranty.
|
||||
//
|
||||
// Alternatively, this file may be used under the terms of Open CASCADE
|
||||
// commercial license or contractual agreement.
|
||||
|
||||
#ifndef NCollection_CellFilterNDim_HeaderFile
|
||||
#define NCollection_CellFilterNDim_HeaderFile
|
||||
|
||||
#include <NCollection_CellFilter.hxx>
|
||||
#include <Standard_Real.hxx>
|
||||
#include <NCollection_List.hxx>
|
||||
#include <NCollection_Map.hxx>
|
||||
#include <NCollection_DataMap.hxx>
|
||||
#include <NCollection_IncAllocator.hxx>
|
||||
#include <NCollection_TypeDef.hxx>
|
||||
#include <NCollection_Array1.hxx>
|
||||
|
||||
/**
|
||||
* A data structure for sorting geometric objects (called targets) in
|
||||
* n-dimensional space into cells, with associated algorithm for fast checking
|
||||
* of coincidence (overlapping, intersection, etc.) with other objects
|
||||
* (called here bullets).
|
||||
*
|
||||
* Purpose of this class is to add possibility to work with CellFilter with unknown
|
||||
dimension count at compilation time.
|
||||
*
|
||||
* For more details look at base class NCollection_CellFilter.
|
||||
*
|
||||
*/
|
||||
|
||||
template <class Inspector> class NCollection_CellFilterNDim
|
||||
{
|
||||
public:
|
||||
typedef TYPENAME Inspector::Target Target;
|
||||
typedef TYPENAME Inspector::Point Point;
|
||||
|
||||
public:
|
||||
|
||||
//! Constructor; initialized by dimension count and cell size.
|
||||
//!
|
||||
//! Note: the cell size must be ensured to be greater than
|
||||
//! maximal co-ordinate of the involved points divided by INT_MAX,
|
||||
//! in order to avoid integer overflow of cell index.
|
||||
//!
|
||||
//! By default cell size is 0, which is invalid; thus if default
|
||||
//! constructor is used, the tool must be initialized later with
|
||||
//! appropriate cell size by call to Reset()
|
||||
NCollection_CellFilterNDim (const Standard_Integer theDim,
|
||||
const Standard_Real theCellSize = 0,
|
||||
const Handle(NCollection_IncAllocator)& theAlloc = 0)
|
||||
: myCellSize(0, theDim - 1)
|
||||
{
|
||||
myDim = theDim;
|
||||
Reset (theCellSize, theAlloc);
|
||||
}
|
||||
|
||||
//! Constructor; initialized by dimension count and cell sizes along each dimension.
|
||||
//! Note: the cell size in each dimension must be ensured to be greater than
|
||||
//! maximal co-ordinate of the involved points by this dimension divided by INT_MAX,
|
||||
//! in order to avoid integer overflow of cell index.
|
||||
NCollection_CellFilterNDim (const Standard_Integer theDim,
|
||||
const NCollection_Array1<Standard_Real> theCellSize,
|
||||
const Handle(NCollection_IncAllocator)& theAlloc = 0)
|
||||
: myCellSize(0, theDim - 1)
|
||||
{
|
||||
myDim = theDim;
|
||||
Reset (theCellSize, theAlloc);
|
||||
}
|
||||
|
||||
//! Clear the data structures, set new cell size and allocator
|
||||
void Reset (Standard_Real theCellSize,
|
||||
const Handle(NCollection_IncAllocator)& theAlloc=0)
|
||||
{
|
||||
for (int i=0; i < myDim; i++)
|
||||
myCellSize(i) = theCellSize;
|
||||
resetAllocator ( theAlloc );
|
||||
}
|
||||
|
||||
//! Clear the data structures and set new cell sizes and allocator
|
||||
void Reset (NCollection_Array1<Standard_Real> theCellSize,
|
||||
const Handle(NCollection_IncAllocator)& theAlloc=0)
|
||||
{
|
||||
myCellSize = theCellSize;
|
||||
resetAllocator ( theAlloc );
|
||||
}
|
||||
|
||||
//! Adds a target object for further search at a point (into only one cell)
|
||||
void Add (const Target& theTarget, const Point &thePnt)
|
||||
{
|
||||
Cell aCell (thePnt, myCellSize);
|
||||
add (aCell, theTarget);
|
||||
}
|
||||
|
||||
//! Adds a target object for further search in the range of cells
|
||||
//! defined by two points (the first point must have all co-ordinates equal or
|
||||
//! less than the same co-ordinate of the second point)
|
||||
void Add (const Target& theTarget,
|
||||
const Point &thePntMin, const Point &thePntMax)
|
||||
{
|
||||
// get cells range by minimal and maximal co-ordinates
|
||||
Cell aCellMin (thePntMin, myCellSize);
|
||||
Cell aCellMax (thePntMax, myCellSize);
|
||||
Cell aCell = aCellMin;
|
||||
// add object recursively into all cells in range
|
||||
iterateAdd (myDim-1, aCell, aCellMin, aCellMax, theTarget);
|
||||
}
|
||||
|
||||
//! Find a target object at a point and remove it from the structures.
|
||||
//! For usage of this method "operator ==" should be defined for Target.
|
||||
void Remove (const Target& theTarget, const Point &thePnt)
|
||||
{
|
||||
Cell aCell (thePnt, myCellSize);
|
||||
remove (aCell, theTarget);
|
||||
}
|
||||
|
||||
//! Find a target object in the range of cells defined by two points and
|
||||
//! remove it from the structures
|
||||
//! (the first point must have all co-ordinates equal or
|
||||
//! less than the same co-ordinate of the second point).
|
||||
//! For usage of this method "operator ==" should be defined for Target.
|
||||
void Remove (const Target& theTarget,
|
||||
const Point &thePntMin, const Point &thePntMax)
|
||||
{
|
||||
// get cells range by minimal and maximal co-ordinates
|
||||
Cell aCellMin (thePntMin, myCellSize);
|
||||
Cell aCellMax (thePntMax, myCellSize);
|
||||
Cell aCell = aCellMin;
|
||||
// remove object recursively from all cells in range
|
||||
iterateRemove (myDim-1, aCell, aCellMin, aCellMax, theTarget);
|
||||
}
|
||||
|
||||
//! Inspect all targets in the cell corresponding to the given point
|
||||
void Inspect (const Point& thePnt, Inspector &theInspector)
|
||||
{
|
||||
Cell aCell (thePnt, myCellSize);
|
||||
inspect (aCell, theInspector);
|
||||
}
|
||||
|
||||
//! Inspect all targets in the cells range limited by two given points
|
||||
//! (the first point must have all co-ordinates equal or
|
||||
//! less than the same co-ordinate of the second point)
|
||||
void Inspect (const Point& thePntMin, const Point& thePntMax,
|
||||
Inspector &theInspector)
|
||||
{
|
||||
// get cells range by minimal and maximal co-ordinates
|
||||
Cell aCellMin (thePntMin, myCellSize);
|
||||
Cell aCellMax (thePntMax, myCellSize);
|
||||
Cell aCell = aCellMin;
|
||||
// inspect object recursively into all cells in range
|
||||
iterateInspect (myDim-1, aCell,
|
||||
aCellMin, aCellMax, theInspector);
|
||||
}
|
||||
|
||||
#if defined(__SUNPRO_CC) && (__SUNPRO_CC <= 0x530)
|
||||
public: // work-around against obsolete SUN WorkShop 5.3 compiler
|
||||
#else
|
||||
protected:
|
||||
#endif
|
||||
|
||||
/**
|
||||
* Auxiliary class for storing points belonging to the cell as the list
|
||||
*/
|
||||
struct ListNode
|
||||
{
|
||||
ListNode()
|
||||
{
|
||||
// Empty constructor is forbidden.
|
||||
Standard_NoSuchObject::Raise("NCollection_CellFilterNDim::ListNode()");
|
||||
}
|
||||
|
||||
Target Object;
|
||||
ListNode *Next;
|
||||
};
|
||||
|
||||
/**
|
||||
* Auxilary structure representing a cell in the space.
|
||||
* Cells are stored in the map, each cell contains list of objects
|
||||
* that belong to that cell.
|
||||
*/
|
||||
struct Cell
|
||||
{
|
||||
public:
|
||||
|
||||
//! Constructor; computes cell indices
|
||||
Cell (const Point& thePnt,
|
||||
const NCollection_Array1<Standard_Real> theCellSize)
|
||||
: index(0, theCellSize.Size() - 1),
|
||||
Objects(0)
|
||||
{
|
||||
for (int i=0; i < theCellSize.Size(); i++)
|
||||
{
|
||||
Standard_Real val = (Standard_Real)(Inspector::Coord(i, thePnt) / theCellSize(i));
|
||||
//If the value of index is greater than
|
||||
//INT_MAX it is decreased correspondingly for the value of INT_MAX. If the value
|
||||
//of index is less than INT_MIN it is increased correspondingly for the absolute
|
||||
//value of INT_MIN.
|
||||
index(i) = long((val > INT_MAX - 1) ? fmod(val, (Standard_Real) INT_MAX)
|
||||
: (val < INT_MIN + 1) ? fmod(val, (Standard_Real) INT_MIN)
|
||||
: val);
|
||||
}
|
||||
}
|
||||
|
||||
//! Copy constructor: ensure that list is not deleted twice
|
||||
Cell (const Cell& theOther)
|
||||
: index(0, theOther.index.Size() - 1)
|
||||
{
|
||||
(*this) = theOther;
|
||||
}
|
||||
|
||||
//! Assignment operator: ensure that list is not deleted twice
|
||||
void operator = (const Cell& theOther)
|
||||
{
|
||||
index = theOther.index;
|
||||
Objects = theOther.Objects;
|
||||
((Cell&)theOther).Objects = 0;
|
||||
}
|
||||
|
||||
//! Destructor; calls destructors for targets contained in the list
|
||||
~Cell ()
|
||||
{
|
||||
for ( ListNode* aNode = Objects; aNode; aNode = aNode->Next )
|
||||
aNode->Object.~Target();
|
||||
// note that list nodes need not to be freed, since IncAllocator is used
|
||||
Objects = 0;
|
||||
}
|
||||
|
||||
//! Compare cell with other one
|
||||
Standard_Boolean IsEqual (const Cell& theOther) const
|
||||
{
|
||||
Standard_Integer myDim = theOther.index.Size();
|
||||
for (int i=0; i < myDim; i++)
|
||||
if ( index(i) != theOther.index(i) ) return Standard_False;
|
||||
return Standard_True;
|
||||
}
|
||||
|
||||
//! Compute hash code
|
||||
Standard_Integer HashCode (const Standard_Integer theUpper) const
|
||||
{
|
||||
// number of bits per each dimension in the hash code
|
||||
Standard_Integer myDim = index.Size();
|
||||
const Standard_Size aShiftBits = (BITS(long)-1) / myDim;
|
||||
long aCode=0;
|
||||
for (int i=0; i < myDim; i++)
|
||||
aCode = ( aCode << aShiftBits ) ^ index(i);
|
||||
return (unsigned)aCode % theUpper;
|
||||
}
|
||||
|
||||
public:
|
||||
NCollection_Array1<long> index;
|
||||
ListNode *Objects;
|
||||
};
|
||||
|
||||
// definition of global functions is needed for map
|
||||
friend Standard_Integer HashCode (const Cell &aCell, const Standard_Integer theUpper)
|
||||
{ return aCell.HashCode(theUpper); }
|
||||
friend Standard_Boolean IsEqual (const Cell &aCell1, const Cell &aCell2)
|
||||
{ return aCell1.IsEqual(aCell2); }
|
||||
|
||||
protected:
|
||||
|
||||
//! Reset allocator to the new one
|
||||
void resetAllocator (const Handle(NCollection_IncAllocator)& theAlloc)
|
||||
{
|
||||
if ( theAlloc.IsNull() )
|
||||
myAllocator = new NCollection_IncAllocator;
|
||||
else
|
||||
myAllocator = theAlloc;
|
||||
myCells.Clear ( myAllocator );
|
||||
}
|
||||
|
||||
//! Add a new target object into the specified cell
|
||||
void add (const Cell& theCell, const Target& theTarget)
|
||||
{
|
||||
// add a new cell or get reference to existing one
|
||||
Cell& aMapCell = (Cell&)myCells.Added (theCell);
|
||||
|
||||
// create a new list node and add it to the beginning of the list
|
||||
ListNode* aNode = (ListNode*)myAllocator->Allocate(sizeof(ListNode));
|
||||
new (&aNode->Object) Target (theTarget);
|
||||
aNode->Next = aMapCell.Objects;
|
||||
aMapCell.Objects = aNode;
|
||||
}
|
||||
|
||||
//! Internal addition function, performing iteration for adjacent cells
|
||||
//! by one dimension; called recursively to cover all dimensions
|
||||
void iterateAdd (int idim, Cell &theCell,
|
||||
const Cell& theCellMin, const Cell& theCellMax,
|
||||
const Target& theTarget)
|
||||
{
|
||||
int start = theCellMin.index(idim);
|
||||
int end = theCellMax.index(idim);
|
||||
for (int i=start; i <= end; i++) {
|
||||
theCell.index(idim) = i;
|
||||
if ( idim ) // recurse
|
||||
iterateAdd (idim-1, theCell, theCellMin, theCellMax, theTarget);
|
||||
else // add to this cell
|
||||
add (theCell, theTarget);
|
||||
}
|
||||
}
|
||||
|
||||
//! Remove the target object from the specified cell
|
||||
void remove (const Cell& theCell, const Target& theTarget)
|
||||
{
|
||||
// check if any objects are recorded in that cell
|
||||
if ( ! myCells.Contains (theCell) )
|
||||
return;
|
||||
|
||||
// iterate by objects in the cell and check each
|
||||
Cell& aMapCell = (Cell&)myCells.Added (theCell);
|
||||
ListNode* aNode = aMapCell.Objects;
|
||||
ListNode* aPrev = NULL;
|
||||
while (aNode)
|
||||
{
|
||||
ListNode* aNext = aNode->Next;
|
||||
if (Inspector::IsEqual (aNode->Object, theTarget))
|
||||
{
|
||||
aNode->Object.~Target();
|
||||
(aPrev ? aPrev->Next : aMapCell.Objects) = aNext;
|
||||
// note that aNode itself need not to be freed, since IncAllocator is used
|
||||
}
|
||||
else
|
||||
aPrev = aNode;
|
||||
aNode = aNext;
|
||||
}
|
||||
}
|
||||
|
||||
//! Internal removal function, performing iteration for adjacent cells
|
||||
//! by one dimension; called recursively to cover all dimensions
|
||||
void iterateRemove (int idim, Cell &theCell,
|
||||
const Cell& theCellMin, const Cell& theCellMax,
|
||||
const Target& theTarget)
|
||||
{
|
||||
int start = theCellMin.index(idim);
|
||||
int end = theCellMax.index(idim);
|
||||
for (int i=start; i <= end; i++) {
|
||||
theCell.index(idim) = i;
|
||||
if ( idim ) // recurse
|
||||
iterateRemove (idim-1, theCell, theCellMin, theCellMax, theTarget);
|
||||
else // remove from this cell
|
||||
remove (theCell, theTarget);
|
||||
}
|
||||
}
|
||||
|
||||
//! Inspect the target objects in the specified cell.
|
||||
void inspect (const Cell& theCell, Inspector& theInspector)
|
||||
{
|
||||
// check if any objects are recorded in that cell
|
||||
if ( ! myCells.Contains (theCell) )
|
||||
return;
|
||||
|
||||
// iterate by objects in the cell and check each
|
||||
Cell& aMapCell = (Cell&)myCells.Added (theCell);
|
||||
ListNode* aNode = aMapCell.Objects;
|
||||
ListNode* aPrev = NULL;
|
||||
while(aNode) {
|
||||
ListNode* aNext = aNode->Next;
|
||||
NCollection_CellFilter_Action anAction =
|
||||
theInspector.Inspect (aNode->Object);
|
||||
// delete items requested to be purged
|
||||
if ( anAction == CellFilter_Purge ) {
|
||||
aNode->Object.~Target();
|
||||
(aPrev ? aPrev->Next : aMapCell.Objects) = aNext;
|
||||
// note that aNode itself need not to be freed, since IncAllocator is used
|
||||
}
|
||||
else
|
||||
aPrev = aNode;
|
||||
aNode = aNext;
|
||||
}
|
||||
}
|
||||
|
||||
//! Inspect the target objects in the specified range of the cells
|
||||
void iterateInspect (int idim, Cell &theCell,
|
||||
const Cell& theCellMin, const Cell& theCellMax,
|
||||
Inspector& theInspector)
|
||||
{
|
||||
int start = theCellMin.index(idim);
|
||||
int end = theCellMax.index(idim);
|
||||
for (int i=start; i <= end; i++) {
|
||||
theCell.index(idim) = i;
|
||||
if ( idim ) // recurse
|
||||
iterateInspect (idim-1, theCell, theCellMin, theCellMax, theInspector);
|
||||
else // inspect this cell
|
||||
inspect (theCell, theInspector);
|
||||
}
|
||||
}
|
||||
|
||||
protected:
|
||||
Standard_Integer myDim;
|
||||
Handle(NCollection_BaseAllocator) myAllocator;
|
||||
NCollection_Map<Cell> myCells;
|
||||
NCollection_Array1<Standard_Real> myCellSize;
|
||||
};
|
||||
|
||||
#endif
|
@ -2273,7 +2273,7 @@ static Standard_Integer OCC25004 (Draw_Interpretor& theDI,
|
||||
Standard_Integer /*theNArg*/,
|
||||
const char** /*theArgs*/)
|
||||
{
|
||||
math_MultipleVarFunction* aFunc = new BraninFunction();
|
||||
BraninFunction aFunc;
|
||||
|
||||
math_Vector aLower(1,2), aUpper(1,2);
|
||||
aLower(1) = -5;
|
||||
@ -2296,7 +2296,7 @@ static Standard_Integer OCC25004 (Draw_Interpretor& theDI,
|
||||
aCurrPnt1(1) = aLower(1) + (aUpper(1) - aLower(1)) * (i - 1) / (aGridOrder - 1.0);
|
||||
aCurrPnt1(2) = aLower(2) + (aUpper(2) - aLower(2)) * (j - 1) / (aGridOrder - 1.0);
|
||||
|
||||
aFunc->Value(aCurrPnt1, aFuncValues(idx));
|
||||
aFunc.Value(aCurrPnt1, aFuncValues(idx));
|
||||
idx++;
|
||||
}
|
||||
}
|
||||
@ -2327,7 +2327,7 @@ static Standard_Integer OCC25004 (Draw_Interpretor& theDI,
|
||||
aLipConst = C;
|
||||
}
|
||||
|
||||
math_GlobOptMin aFinder(aFunc, aLower, aUpper, aLipConst);
|
||||
math_GlobOptMin aFinder(&aFunc, aLower, aUpper, aLipConst);
|
||||
aFinder.Perform();
|
||||
//(-pi , 12.275), (pi , 2.275), (9.42478, 2.475)
|
||||
|
||||
|
@ -1,6 +1,6 @@
|
||||
// Created on: 2014-01-20
|
||||
// Created by: Alexaner Malyshev
|
||||
// Copyright (c) 2014-2014 OPEN CASCADE SAS
|
||||
// Copyright (c) 2014-2015 OPEN CASCADE SAS
|
||||
//
|
||||
// This file is part of Open CASCADE Technology software library.
|
||||
//
|
||||
@ -45,7 +45,9 @@ math_GlobOptMin::math_GlobOptMin(math_MultipleVarFunction* theFunc,
|
||||
myTmp(1, myN),
|
||||
myV(1, myN),
|
||||
myMaxV(1, myN),
|
||||
myExpandCoeff(1, myN)
|
||||
myExpandCoeff(1, myN),
|
||||
myCellSize(0, myN - 1),
|
||||
myFilter(theFunc->NbVariables())
|
||||
{
|
||||
Standard_Integer i;
|
||||
|
||||
@ -78,6 +80,11 @@ math_GlobOptMin::math_GlobOptMin(math_MultipleVarFunction* theFunc,
|
||||
myTol = theDiscretizationTol;
|
||||
mySameTol = theSameTol;
|
||||
|
||||
const Standard_Integer aMaxSquareSearchSol = 200;
|
||||
Standard_Integer aSolNb = Standard_Integer(Pow(3.0, Standard_Real(myN)));
|
||||
myMinCellFilterSol = Max(2 * aSolNb, aMaxSquareSearchSol);
|
||||
initCellSize();
|
||||
|
||||
myDone = Standard_False;
|
||||
}
|
||||
|
||||
@ -122,6 +129,8 @@ void math_GlobOptMin::SetGlobalParams(math_MultipleVarFunction* theFunc,
|
||||
myTol = theDiscretizationTol;
|
||||
mySameTol = theSameTol;
|
||||
|
||||
initCellSize();
|
||||
|
||||
myDone = Standard_False;
|
||||
}
|
||||
|
||||
@ -238,6 +247,7 @@ void math_GlobOptMin::Perform(const Standard_Boolean isFindSingleSolution)
|
||||
myE3 = - maxLength * myTol * myC / 4.0;
|
||||
}
|
||||
|
||||
isFirstCellFilterInvoke = Standard_True;
|
||||
computeGlobalExtremum(myN);
|
||||
|
||||
myDone = Standard_True;
|
||||
@ -398,7 +408,6 @@ void math_GlobOptMin::computeGlobalExtremum(Standard_Integer j)
|
||||
Standard_Boolean isInside = Standard_False;
|
||||
Standard_Real r;
|
||||
|
||||
|
||||
for(myX(j) = myA(j) + myE1; myX(j) < myB(j) + myE1; myX(j) += myV(j))
|
||||
{
|
||||
if (myX(j) > myB(j))
|
||||
@ -444,6 +453,8 @@ void math_GlobOptMin::computeGlobalExtremum(Standard_Integer j)
|
||||
for(i = 1; i <= myN; i++)
|
||||
myY.Append(aStepBestPoint(i));
|
||||
mySolCount++;
|
||||
|
||||
isFirstCellFilterInvoke = Standard_True;
|
||||
}
|
||||
|
||||
aRealStep = myE2 + Abs(myF - d) / myC;
|
||||
@ -500,20 +511,55 @@ Standard_Boolean math_GlobOptMin::isStored(const math_Vector& thePnt)
|
||||
math_Vector aTol(1, myN);
|
||||
aTol = (myB - myA) * mySameTol;
|
||||
|
||||
for(i = 0; i < mySolCount; i++)
|
||||
// C1 * n^2 = C2 * 3^dim * n
|
||||
if (mySolCount < myMinCellFilterSol)
|
||||
{
|
||||
isSame = Standard_True;
|
||||
for(j = 1; j <= myN; j++)
|
||||
for(i = 0; i < mySolCount; i++)
|
||||
{
|
||||
if ((Abs(thePnt(j) - myY(i * myN + j))) > aTol(j))
|
||||
isSame = Standard_True;
|
||||
for(j = 1; j <= myN; j++)
|
||||
{
|
||||
isSame = Standard_False;
|
||||
break;
|
||||
if ((Abs(thePnt(j) - myY(i * myN + j))) > aTol(j))
|
||||
{
|
||||
isSame = Standard_False;
|
||||
break;
|
||||
}
|
||||
}
|
||||
if (isSame == Standard_True)
|
||||
return Standard_True;
|
||||
}
|
||||
}
|
||||
else
|
||||
{
|
||||
NCollection_CellFilter_NDimInspector anInspector(myN, Precision::PConfusion());
|
||||
if (isFirstCellFilterInvoke)
|
||||
{
|
||||
myFilter.Reset(myCellSize);
|
||||
|
||||
// Copy initial data into cell filter.
|
||||
for(Standard_Integer aSolIdx = 0; aSolIdx < mySolCount; aSolIdx++)
|
||||
{
|
||||
math_Vector aVec(1, myN);
|
||||
for(Standard_Integer aSolDim = 1; aSolDim <= myN; aSolDim++)
|
||||
aVec(aSolDim) = myY(aSolIdx * myN + aSolDim);
|
||||
|
||||
myFilter.Add(aVec, aVec);
|
||||
}
|
||||
}
|
||||
if (isSame == Standard_True)
|
||||
return Standard_True;
|
||||
|
||||
isFirstCellFilterInvoke = Standard_False;
|
||||
|
||||
math_Vector aLow(1, myN), anUp(1, myN);
|
||||
anInspector.Shift(thePnt, myCellSize, aLow, anUp);
|
||||
|
||||
anInspector.ClearFind();
|
||||
anInspector.SetCurrent(thePnt);
|
||||
myFilter.Inspect(aLow, anUp, anInspector);
|
||||
if (!anInspector.isFind())
|
||||
{
|
||||
// Point is out of close cells, add new one.
|
||||
myFilter.Add(thePnt, thePnt);
|
||||
}
|
||||
}
|
||||
return Standard_False;
|
||||
}
|
||||
@ -556,3 +602,16 @@ void math_GlobOptMin::Points(const Standard_Integer theIndex, math_Vector& theSo
|
||||
for(j = 1; j <= myN; j++)
|
||||
theSol(j) = myY((theIndex - 1) * myN + j);
|
||||
}
|
||||
|
||||
//=======================================================================
|
||||
//function : initCellSize
|
||||
//purpose :
|
||||
//=======================================================================
|
||||
void math_GlobOptMin::initCellSize()
|
||||
{
|
||||
for(Standard_Integer anIdx = 1; anIdx <= myN; anIdx++)
|
||||
{
|
||||
myCellSize(anIdx - 1) = (myGlobB(anIdx) - myGlobA(anIdx))
|
||||
* Precision::PConfusion() / (2.0 * Sqrt(2.0));
|
||||
}
|
||||
}
|
||||
|
@ -1,6 +1,6 @@
|
||||
// Created on: 2014-01-20
|
||||
// Created by: Alexaner Malyshev
|
||||
// Copyright (c) 2014-2014 OPEN CASCADE SAS
|
||||
// Copyright (c) 2014-2015 OPEN CASCADE SAS
|
||||
//
|
||||
// This file is part of Open CASCADE Technology software library.
|
||||
//
|
||||
@ -16,10 +16,86 @@
|
||||
#ifndef _math_GlobOptMin_HeaderFile
|
||||
#define _math_GlobOptMin_HeaderFile
|
||||
|
||||
#include <gp_Pnt.hxx>
|
||||
#include <gp_Pnt2d.hxx>
|
||||
#include <NCollection_CellFilterNDim.hxx>
|
||||
#include <math_MultipleVarFunction.hxx>
|
||||
#include <NCollection_Sequence.hxx>
|
||||
#include <Standard_Type.hxx>
|
||||
|
||||
class NCollection_CellFilter_NDimInspector
|
||||
{
|
||||
public:
|
||||
|
||||
//! Points and target type
|
||||
typedef math_Vector Point;
|
||||
typedef math_Vector Target;
|
||||
|
||||
NCollection_CellFilter_NDimInspector(const Standard_Integer theDim,
|
||||
const Standard_Real theTol)
|
||||
: myCurrent(1, theDim)
|
||||
{
|
||||
myTol = theTol * theTol;
|
||||
myIsFind = Standard_False;
|
||||
Dimension = theDim;
|
||||
}
|
||||
|
||||
//! Access to co-ordinate
|
||||
static Standard_Real Coord (int i, const Point &thePnt)
|
||||
{
|
||||
return thePnt(i + 1);
|
||||
}
|
||||
|
||||
//! Auxiliary method to shift point by each coordinate on given value;
|
||||
//! useful for preparing a points range for Inspect with tolerance
|
||||
void Shift (const Point& thePnt,
|
||||
const NCollection_Array1<Standard_Real> &theTol,
|
||||
Point& theLowPnt,
|
||||
Point& theUppPnt) const
|
||||
{
|
||||
for(Standard_Integer anIdx = 1; anIdx <= Dimension; anIdx++)
|
||||
{
|
||||
theLowPnt(anIdx) = thePnt(anIdx) - theTol(anIdx - 1);
|
||||
theUppPnt(anIdx) = thePnt(anIdx) + theTol(anIdx - 1);
|
||||
}
|
||||
}
|
||||
|
||||
void ClearFind()
|
||||
{
|
||||
myIsFind = Standard_False;
|
||||
}
|
||||
|
||||
Standard_Boolean isFind()
|
||||
{
|
||||
return myIsFind;
|
||||
}
|
||||
|
||||
//! Set current point to search for coincidence
|
||||
void SetCurrent (const math_Vector& theCurPnt)
|
||||
{
|
||||
myCurrent = theCurPnt;
|
||||
}
|
||||
|
||||
//! Implementation of inspection method
|
||||
NCollection_CellFilter_Action Inspect (const Target& theObject)
|
||||
{
|
||||
Standard_Real aSqDist = (myCurrent - theObject).Norm2();
|
||||
|
||||
if(aSqDist < myTol)
|
||||
{
|
||||
myIsFind = Standard_True;
|
||||
}
|
||||
|
||||
return CellFilter_Keep;
|
||||
}
|
||||
|
||||
private:
|
||||
Standard_Real myTol;
|
||||
math_Vector myCurrent;
|
||||
Standard_Boolean myIsFind;
|
||||
Standard_Integer Dimension;
|
||||
};
|
||||
|
||||
//! This class represents Evtushenko's algorithm of global optimization based on nonuniform mesh.<br>
|
||||
//! Article: Yu. Evtushenko. Numerical methods for finding global extreme (case of a non-uniform mesh). <br>
|
||||
//! U.S.S.R. Comput. Maths. Math. Phys., Vol. 11, N 6, pp. 38-54.
|
||||
@ -69,6 +145,9 @@ public:
|
||||
|
||||
private:
|
||||
|
||||
// Compute cell size.
|
||||
void initCellSize();
|
||||
|
||||
math_GlobOptMin & operator = (const math_GlobOptMin & theOther);
|
||||
|
||||
Standard_Boolean computeLocalExtremum(const math_Vector& thePnt, Standard_Real& theVal, math_Vector& theOutPnt);
|
||||
@ -119,6 +198,11 @@ private:
|
||||
math_Vector myMaxV; // Max Steps array.
|
||||
math_Vector myExpandCoeff; // Define expand coefficient between neighboring indiced dimensions.
|
||||
|
||||
NCollection_Array1<Standard_Real> myCellSize;
|
||||
Standard_Integer myMinCellFilterSol;
|
||||
Standard_Boolean isFirstCellFilterInvoke;
|
||||
NCollection_CellFilterNDim<NCollection_CellFilter_NDimInspector> myFilter;
|
||||
|
||||
Standard_Real myF; // Current value of Global optimum.
|
||||
};
|
||||
|
||||
|
26
tests/bugs/fclasses/bug26184_1
Normal file
26
tests/bugs/fclasses/bug26184_1
Normal file
@ -0,0 +1,26 @@
|
||||
puts "============"
|
||||
puts "OCC26184"
|
||||
puts "============"
|
||||
puts ""
|
||||
#######################################################################
|
||||
# GeomAPI_ExtremaCurveCurve hangs on parallel b-spline curves
|
||||
#######################################################################
|
||||
|
||||
restore [locate_data_file bug26184_Curve_Extrema_1_12971.brep] a1
|
||||
restore [locate_data_file bug26184_Curve_Extrema_2_12971.brep] a2
|
||||
|
||||
mkcurve c1 a1
|
||||
mkcurve c2 a2
|
||||
|
||||
cpulimit 20
|
||||
|
||||
dchrono h reset; dchrono h start
|
||||
extrema c1 c2
|
||||
dchrono h stop; dchrono h show
|
||||
|
||||
regexp {CPU user time: (\d*)} [dchrono h show] dummy sec
|
||||
if {$sec > 10} {
|
||||
puts "Error: too long computation time $sec seconds"
|
||||
} else {
|
||||
puts "Computation time is OK"
|
||||
}
|
26
tests/bugs/fclasses/bug26184_2
Normal file
26
tests/bugs/fclasses/bug26184_2
Normal file
@ -0,0 +1,26 @@
|
||||
puts "============"
|
||||
puts "OCC26184"
|
||||
puts "============"
|
||||
puts ""
|
||||
#######################################################################
|
||||
# GeomAPI_ExtremaCurveCurve hangs on parallel b-spline curves
|
||||
#######################################################################
|
||||
|
||||
restore [locate_data_file bug26184_Curve_Extrema_1_13767.brep] a1
|
||||
restore [locate_data_file bug26184_Curve_Extrema_2_13767.brep] a2
|
||||
|
||||
mkcurve c1 a1
|
||||
mkcurve c2 a2
|
||||
|
||||
cpulimit 20
|
||||
|
||||
dchrono h reset; dchrono h start
|
||||
extrema c1 c2
|
||||
dchrono h stop; dchrono h show
|
||||
|
||||
regexp {CPU user time: (\d*)} [dchrono h show] dummy sec
|
||||
if {$sec > 10} {
|
||||
puts "Error: too long computation time $sec seconds"
|
||||
} else {
|
||||
puts "Computation time is OK"
|
||||
}
|
@ -24,7 +24,7 @@ if { [regexp {Debug mode} [dversion]] } {
|
||||
set max_t_1 20
|
||||
} else {
|
||||
if { [regexp {Windows} [dversion]] } {
|
||||
set max_t_1 15
|
||||
set max_t_1 9
|
||||
} else {
|
||||
set max_t_1 8
|
||||
}
|
||||
|
@ -1,5 +1,5 @@
|
||||
puts "TODO OCC12345 ALL: Faulty OCC565: function intersection works wrongly with trimmed Surfaces"
|
||||
puts "TODO OCC12345 Windows Linux MacOS: Faulty OCC565: function intersection works wrongly with infinite Surfaces"
|
||||
puts "TODO ?OCC12345 ALL: Faulty OCC565: function intersection works wrongly with infinite Surfaces"
|
||||
|
||||
puts "========"
|
||||
puts "OCC567"
|
||||
|
Loading…
x
Reference in New Issue
Block a user