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

0027552: Wire creation fails depending on the order of edges

1) BRepBuilderAPI_MakeWire::Add (const TopTools_ListOfShape &L) method have been completely rewritten. The order of edges is not significant now.
2) The geometric proximity of free vertices from already existing wire and from input list of edges are also have been taken into account. If such vertices are coincident with each other then they are fused into the one. The original wire remains untouched topologically (yet the tolerances and points can be modified).
3) UBTreeFiller is used to speed up the process of picking of coincident vertices.
4) BRepLib now contains the 'new' method - BoundingVertex(..). The implemenation of this method are taken from BOPTools_AlgoTools::MakeVertex(..).
5) The '-unsorted' argument have been added to 'wire' command.

Conflicts:
	src/QABugs/QABugs_20.cxx

Add missing include.

Eliminate warning.
This commit is contained in:
isn 2016-06-20 15:58:05 +03:00 committed by bugmaster
parent 01b5b3df55
commit 07ef8bdfa2
12 changed files with 750 additions and 342 deletions

View File

@ -107,24 +107,7 @@ static
const BOPTools_ListOfCoupleOfShape& theLCS,
const gp_Pnt& aP);
//=======================================================================
// function: BOPTools_AlgoTools_ComparePoints
// purpose: implementation of IsLess() function for two points
//=======================================================================
struct BOPTools_AlgoTools_ComparePoints {
bool operator()(const gp_Pnt& theP1, const gp_Pnt& theP2)
{
for (Standard_Integer i = 1; i <= 3; ++i) {
if (theP1.Coord(i) < theP2.Coord(i)) {
return Standard_True;
}
else if (theP1.Coord(i) > theP2.Coord(i)) {
return Standard_False;
}
}
return Standard_False;
}
};
//=======================================================================
// function: MakeConnexityBlocks
@ -1514,105 +1497,16 @@ Standard_Integer BOPTools_AlgoTools::ComputeVV(const TopoDS_Vertex& aV1,
void BOPTools_AlgoTools::MakeVertex(const BOPCol_ListOfShape& aLV,
TopoDS_Vertex& aVnew)
{
Standard_Integer aNb;
//
aNb=aLV.Extent();
if (!aNb) {
return;
}
//
else if (aNb==1) {
Standard_Integer aNb = aLV.Extent();
if (aNb == 1)
aVnew=*((TopoDS_Vertex*)(&aLV.First()));
return;
}
//
else if (aNb==2) {
Standard_Integer m, n;
Standard_Real aR[2], dR, aD, aEps;
TopoDS_Vertex aV[2];
gp_Pnt aP[2];
else if (aNb > 1)
{
Standard_Real aNTol;
gp_Pnt aNC;
BRepLib::BoundingVertex(aLV, aNC, aNTol);
BRep_Builder aBB;
//
aEps=RealEpsilon();
for (m=0; m<aNb; ++m) {
aV[m]=(!m)?
*((TopoDS_Vertex*)(&aLV.First())):
*((TopoDS_Vertex*)(&aLV.Last()));
aP[m]=BRep_Tool::Pnt(aV[m]);
aR[m]=BRep_Tool::Tolerance(aV[m]);
}
//
m=0; // max R
n=1; // min R
if (aR[0]<aR[1]) {
m=1;
n=0;
}
//
dR=aR[m]-aR[n]; // dR >= 0.
gp_Vec aVD(aP[m], aP[n]);
aD=aVD.Magnitude();
//
if (aD<=dR || aD<aEps) {
aBB.MakeVertex (aVnew, aP[m], aR[m]);
}
else {
Standard_Real aRr;
gp_XYZ aXYZr;
gp_Pnt aPr;
//
aRr=0.5*(aR[m]+aR[n]+aD);
aXYZr=0.5*(aP[m].XYZ()+aP[n].XYZ()-aVD.XYZ()*(dR/aD));
aPr.SetXYZ(aXYZr);
//
aBB.MakeVertex (aVnew, aPr, aRr);
}
return;
}// else if (aNb==2) {
//
else { // if (aNb>2)
// compute the point
//
// issue 0027540 - sum of doubles may depend on the order
// of addition, thus sort the coordinates for stable result
Standard_Integer i;
NCollection_Array1<gp_Pnt> aPoints(0, aNb-1);
BOPCol_ListIteratorOfListOfShape aIt(aLV);
for (i = 0; aIt.More(); aIt.Next(), ++i) {
const TopoDS_Vertex& aVi = *((TopoDS_Vertex*)(&aIt.Value()));
gp_Pnt aPi = BRep_Tool::Pnt(aVi);
aPoints(i) = aPi;
}
//
std::sort(aPoints.begin(), aPoints.end(), BOPTools_AlgoTools_ComparePoints());
//
gp_XYZ aXYZ(0., 0., 0.);
for (i = 0; i < aNb; ++i) {
aXYZ += aPoints(i).XYZ();
}
aXYZ.Divide((Standard_Real)aNb);
//
gp_Pnt aP(aXYZ);
//
// compute the tolerance for the new vertex
Standard_Real aTi, aDi, aDmax;
//
aDmax=-1.;
aIt.Initialize(aLV);
for (; aIt.More(); aIt.Next()) {
TopoDS_Vertex& aVi=*((TopoDS_Vertex*)(&aIt.Value()));
gp_Pnt aPi=BRep_Tool::Pnt(aVi);
aTi=BRep_Tool::Tolerance(aVi);
aDi=aP.SquareDistance(aPi);
aDi=sqrt(aDi);
aDi=aDi+aTi;
if (aDi > aDmax) {
aDmax=aDi;
}
}
//
BRep_Builder aBB;
aBB.MakeVertex (aVnew, aP, aDmax);
aBB.MakeVertex(aVnew, aNC, aNTol);
}
}
//=======================================================================

View File

@ -90,10 +90,32 @@
#include <TShort_HArray1OfShortReal.hxx>
#include <TColgp_Array1OfXY.hxx>
#include <algorithm>
// TODO - not thread-safe static variables
static Standard_Real thePrecision = Precision::Confusion();
static Handle(Geom_Plane) thePlane;
//=======================================================================
// function: BRepLib_ComparePoints
// purpose: implementation of IsLess() function for two points
//=======================================================================
struct BRepLib_ComparePoints {
bool operator()(const gp_Pnt& theP1, const gp_Pnt& theP2)
{
for (Standard_Integer i = 1; i <= 3; ++i) {
if (theP1.Coord(i) < theP2.Coord(i)) {
return Standard_True;
}
else if (theP1.Coord(i) > theP2.Coord(i)) {
return Standard_False;
}
}
return Standard_False;
}
};
//=======================================================================
//function : Precision
//purpose :
@ -2059,3 +2081,108 @@ void BRepLib::ReverseSortFaces (const TopoDS_Shape& Sh,
}
//=======================================================================
// function: BoundingVertex
// purpose :
//=======================================================================
void BRepLib::BoundingVertex(const NCollection_List<TopoDS_Shape>& theLV,
gp_Pnt& theNewCenter, Standard_Real& theNewTol)
{
Standard_Integer aNb;
//
aNb=theLV.Extent();
if (aNb < 2) {
return;
}
//
else if (aNb==2) {
Standard_Integer m, n;
Standard_Real aR[2], dR, aD, aEps;
TopoDS_Vertex aV[2];
gp_Pnt aP[2];
//
aEps=RealEpsilon();
for (m=0; m<aNb; ++m) {
aV[m]=(!m)?
*((TopoDS_Vertex*)(&theLV.First())):
*((TopoDS_Vertex*)(&theLV.Last()));
aP[m]=BRep_Tool::Pnt(aV[m]);
aR[m]=BRep_Tool::Tolerance(aV[m]);
}
//
m=0; // max R
n=1; // min R
if (aR[0]<aR[1]) {
m=1;
n=0;
}
//
dR=aR[m]-aR[n]; // dR >= 0.
gp_Vec aVD(aP[m], aP[n]);
aD=aVD.Magnitude();
//
if (aD<=dR || aD<aEps) {
theNewCenter = aP[m];
theNewTol = aR[m];
}
else {
Standard_Real aRr;
gp_XYZ aXYZr;
gp_Pnt aPr;
//
aRr=0.5*(aR[m]+aR[n]+aD);
aXYZr=0.5*(aP[m].XYZ()+aP[n].XYZ()-aVD.XYZ()*(dR/aD));
aPr.SetXYZ(aXYZr);
//
theNewCenter = aPr;
theNewTol = aRr;
//aBB.MakeVertex (aVnew, aPr, aRr);
}
return;
}// else if (aNb==2) {
//
else { // if (aNb>2)
// compute the point
//
// issue 0027540 - sum of doubles may depend on the order
// of addition, thus sort the coordinates for stable result
Standard_Integer i;
NCollection_Array1<gp_Pnt> aPoints(0, aNb-1);
NCollection_List<TopoDS_Edge>::Iterator aIt(theLV);
for (i = 0; aIt.More(); aIt.Next(), ++i) {
const TopoDS_Vertex& aVi = *((TopoDS_Vertex*)(&aIt.Value()));
gp_Pnt aPi = BRep_Tool::Pnt(aVi);
aPoints(i) = aPi;
}
//
std::sort(aPoints.begin(), aPoints.end(), BRepLib_ComparePoints());
//
gp_XYZ aXYZ(0., 0., 0.);
for (i = 0; i < aNb; ++i) {
aXYZ += aPoints(i).XYZ();
}
aXYZ.Divide((Standard_Real)aNb);
//
gp_Pnt aP(aXYZ);
//
// compute the tolerance for the new vertex
Standard_Real aTi, aDi, aDmax;
//
aDmax=-1.;
aIt.Initialize(theLV);
for (; aIt.More(); aIt.Next()) {
TopoDS_Vertex& aVi=*((TopoDS_Vertex*)(&aIt.Value()));
gp_Pnt aPi=BRep_Tool::Pnt(aVi);
aTi=BRep_Tool::Tolerance(aVi);
aDi=aP.SquareDistance(aPi);
aDi=sqrt(aDi);
aDi=aDi+aTi;
if (aDi > aDmax) {
aDmax=aDi;
}
}
//
theNewCenter = aP;
theNewTol = aDmax;
}
}

View File

@ -26,6 +26,7 @@
#include <GeomAbs_Shape.hxx>
#include <Standard_Integer.hxx>
#include <TopTools_ListOfShape.hxx>
#include <NCollection_List.hxx>
class Geom_Plane;
class TopoDS_Edge;
class TopoDS_Shape;
@ -184,8 +185,12 @@ public:
//! Returns TRUE if any correction is done.
Standard_EXPORT static Standard_Boolean EnsureNormalConsistency (const TopoDS_Shape& S, const Standard_Real theAngTol = 0.001, const Standard_Boolean ForceComputeNormals = Standard_False);
//! Calculates the bounding sphere around the set of vertexes from the theLV list.
//! Returns the center (theNewCenter) and the radius (theNewTol) of this sphere.
//! This can be used to construct the new vertex which covers the given set of
//! other vertices.
Standard_EXPORT static void BoundingVertex(const NCollection_List<TopoDS_Shape>& theLV,
gp_Pnt& theNewCenter, Standard_Real& theNewTol);
protected:

View File

@ -138,15 +138,24 @@ void BRepLib_MakeWire::Add(const TopoDS_Wire& W)
}
}
//=======================================================================
//function : Add
//purpose :
//=======================================================================
void BRepLib_MakeWire::Add(const TopoDS_Edge& E)
{
Add(E, Standard_True);
}
//=======================================================================
//function : Add
//purpose :
// PMN 19/03/1998 For the Problem of performance TopExp::Vertices are not used on wire
// PMN 10/09/1998 In case if the wire is previously closed (or degenerated)
// TopExp::Vertices is used to reduce the ambiguity.
// IsCheckGeometryProximity flag : If true => check for the geometry proximity of vertices
//=======================================================================
void BRepLib_MakeWire::Add(const TopoDS_Edge& E)
void BRepLib_MakeWire::Add(const TopoDS_Edge& E, Standard_Boolean IsCheckGeometryProximity)
{
Standard_Boolean forward = Standard_False;
@ -219,15 +228,13 @@ void BRepLib_MakeWire::Add(const TopoDS_Edge& E)
}
}
}
else {
else if (IsCheckGeometryProximity)
{
// search if there is a similar vertex in the edge
gp_Pnt PE = BRep_Tool::Pnt(VE);
// Standard_Boolean newvertex = Standard_False;
TopTools_MapIteratorOfMapOfShape itm(myVertices);
while (itm.More()) {
const TopoDS_Vertex& VW = TopoDS::Vertex(itm.Key());
for (Standard_Integer i = 1; i <= myVertices.Extent(); i++) {
const TopoDS_Vertex& VW = TopoDS::Vertex(myVertices.FindKey(i));
gp_Pnt PW = BRep_Tool::Pnt(VW);
Standard_Real l = PE.Distance(PW);
@ -259,7 +266,6 @@ void BRepLib_MakeWire::Add(const TopoDS_Edge& E)
}
break;
}
itm.Next();
}
if (copyedge) {
connected = Standard_True;
@ -289,10 +295,8 @@ void BRepLib_MakeWire::Add(const TopoDS_Edge& E)
gp_Pnt PE = BRep_Tool::Pnt(VE);
Standard_Boolean newvertex = Standard_False;
TopTools_MapIteratorOfMapOfShape itm(myVertices);
while (itm.More()) {
const TopoDS_Vertex& VW = TopoDS::Vertex(itm.Key());
for (Standard_Integer i = 1; i <= myVertices.Extent(); i++) {
const TopoDS_Vertex& VW = TopoDS::Vertex(myVertices.FindKey(i));
gp_Pnt PW = BRep_Tool::Pnt(VW);
Standard_Real l = PE.Distance(PW), tolE, tolW;
tolW = BRep_Tool::Tolerance(VW);
@ -320,7 +324,6 @@ void BRepLib_MakeWire::Add(const TopoDS_Edge& E)
break;
}
itm.Next();
}
if (!newvertex) {
myVertices.Add(VE);

View File

@ -24,15 +24,18 @@
#include <BRepLib_WireError.hxx>
#include <TopoDS_Edge.hxx>
#include <TopoDS_Vertex.hxx>
#include <TopTools_MapOfShape.hxx>
#include <TopTools_IndexedMapOfShape.hxx>
#include <BRepLib_MakeShape.hxx>
#include <TopTools_ListOfShape.hxx>
#include <NCollection_Map.hxx>
#include <Bnd_Box.hxx>
#include <NCollection_UBTree.hxx>
class StdFail_NotDone;
class TopoDS_Edge;
class TopoDS_Wire;
class TopoDS_Vertex;
//! Provides methods to build wires.
//!
//! A wire may be built :
@ -75,6 +78,7 @@ class TopoDS_Vertex;
//! MW.Add(anEdge);
//!
//! TopoDS_Wire W = MW;
class BRepLib_MakeWire : public BRepLib_MakeShape
{
public:
@ -126,6 +130,58 @@ public:
//! Returns the last connecting vertex.
Standard_EXPORT const TopoDS_Vertex& Vertex() const;
private:
class BRepLib_BndBoxVertexSelector : public NCollection_UBTree <Standard_Integer,Bnd_Box>::Selector
{
public:
BRepLib_BndBoxVertexSelector(const TopTools_IndexedMapOfShape& theMapOfShape)
: BRepLib_BndBoxVertexSelector::Selector(), myMapOfShape (theMapOfShape)
{}
Standard_Boolean Reject (const Bnd_Box& theBox) const
{
return theBox.IsOut(myVBox);
}
Standard_Boolean Accept (const Standard_Integer& theObj);
void SetCurrentVertex (const gp_Pnt& theP, Standard_Real theTol,
Standard_Integer theVInd);
const NCollection_List<Standard_Integer>& GetResultInds () const
{
return myResultInd;
}
void ClearResInds()
{
myResultInd.Clear();
}
private:
BRepLib_BndBoxVertexSelector(const BRepLib_BndBoxVertexSelector& );
BRepLib_BndBoxVertexSelector& operator=(const BRepLib_BndBoxVertexSelector& );
const TopTools_IndexedMapOfShape& myMapOfShape; //vertices
gp_Pnt myP;
Standard_Real mySTol;
Standard_Integer myVInd;
Bnd_Box myVBox;
NCollection_List<Standard_Integer> myResultInd;
};
void CollectCoincidentVertices(const TopTools_ListOfShape& theL,
NCollection_List<NCollection_List<TopoDS_Vertex>>& theGrVL);
void CreateNewVertices(const NCollection_List<NCollection_List<TopoDS_Vertex>>& theGrVL,
NCollection_DataMap<TopoDS_Vertex, TopoDS_Vertex>& theO2NV);
void CreateNewListOfEdges(const TopTools_ListOfShape& theL,
const NCollection_DataMap<TopoDS_Vertex, TopoDS_Vertex>& theO2NV,
TopTools_ListOfShape& theNewEList);
void Add(const TopoDS_Edge& E, Standard_Boolean IsCheckGeometryProximity);
@ -133,21 +189,16 @@ protected:
private:
BRepLib_WireError myError;
TopoDS_Edge myEdge;
TopoDS_Vertex myVertex;
TopTools_MapOfShape myVertices;
TopTools_IndexedMapOfShape myVertices;
TopoDS_Vertex FirstVertex;
TopoDS_Vertex VF;
TopoDS_Vertex VL;
};

View File

@ -36,6 +36,9 @@
#include <TopTools_MapIteratorOfMapOfShape.hxx>
#include <TopTools_MapOfOrientedShape.hxx>
#include <TopTools_MapOfShape.hxx>
#include <NCollection_UBTreeFiller.hxx>
#include <BRepBndLib.hxx>
#include <BRepLib_MakeVertex.hxx>
//=======================================================================
//function : Add
@ -44,199 +47,352 @@
void BRepLib_MakeWire::Add(const TopTools_ListOfShape& L)
{
myError = BRepLib_WireDone;
if (!myShape.IsNull()) myShape.Closed(Standard_False);
Standard_Integer aLSize = 0;
Standard_Integer aRefSize = L.Size();
if (!L.IsEmpty())
{
///
NCollection_List<NCollection_List<TopoDS_Vertex>> aGrVL;
TopTools_IndexedDataMapOfShapeListOfShape aMapVE;
if (!L.IsEmpty()) {
CollectCoincidentVertices(L, aGrVL);
NCollection_DataMap<TopoDS_Vertex, TopoDS_Vertex> anO2NV;
CreateNewVertices(aGrVL, anO2NV);
TopTools_ListOfShape aNewEList;
CreateNewListOfEdges(L, anO2NV, aNewEList);
///
TopExp::MapShapesAndAncestors(myShape, TopAbs_VERTEX, TopAbs_EDGE, aMapVE);
TopTools_MapOfShape aProcessedEdges;
TopExp_Explorer anExp;
TopTools_ListOfShape anActEdges, aNeighEdges;
if (myEdge.IsNull())
{
//take the first edge from the list and add it
const TopoDS_Edge& aFE = TopoDS::Edge(aNewEList.First());
Add(aFE);
aProcessedEdges.Add(aFE);
anActEdges.Append(aFE);
aLSize++;
}
else
{
//existing edges are already connected
for (anExp.Init(myShape, TopAbs_EDGE); anExp.More(); anExp.Next())
{
const TopoDS_Shape& aCSh = anExp.Current();
aProcessedEdges.Add(aCSh);
anActEdges.Append(aCSh);
}
}
TopTools_ListIteratorOfListOfShape anItL1, anItL2;
for (anItL1.Initialize(aNewEList); anItL1.More(); anItL1.Next())
TopExp::MapShapesAndAncestors(anItL1.Value(), TopAbs_VERTEX, TopAbs_EDGE, aMapVE);
while (!anActEdges.IsEmpty())
{
anItL2.Initialize(anActEdges);
for (;anItL2.More(); anItL2.Next())
{
const TopoDS_Shape& aCE = anItL2.Value();
anExp.Init(aCE, TopAbs_VERTEX);
for (;anExp.More(); anExp.Next())
{
const TopoDS_Shape& aCV = anExp.Current();
for (anItL1.Initialize(aMapVE.FindFromKey(aCV)); anItL1.More(); anItL1.Next())
{
const TopoDS_Shape& aNE = anItL1.Value(); //neighbor edge
if (!aProcessedEdges.Contains(aNE))
{
Add(TopoDS::Edge(aNE), Standard_False);
aNeighEdges.Append(aNE);
aProcessedEdges.Add(aNE);
aLSize++;
}
}
}
}
anActEdges.Clear();
anActEdges.Append(aNeighEdges);
}
}
if (aLSize == aRefSize)
Done();
else
{
NotDone();
TopTools_MapOfShape mapLocale;
mapLocale.Assign(myVertices);
TopTools_DataMapOfShapeShape toCopy;
TopTools_ListOfShape toAdd, nlist, rlist;
BRep_Builder BB;
myError = BRepLib_DisconnectedWire;
}
}
TopExp_Explorer exv;
TopTools_MapIteratorOfMapOfShape itMS;
TopTools_ListIteratorOfListOfShape itList(L);
for (;itList.More(); itList.Next()) {
const TopoDS_Edge& curEd=TopoDS::Edge(itList.Value());
if (!curEd.IsNull()) {
rlist.Clear();
nlist.Clear();
Standard_Boolean copEd=Standard_False;
if (myEdge.IsNull()) {
Add(curEd);
if (!VF.IsNull()) mapLocale.Add(VF);
if (!VL.IsNull()) mapLocale.Add(VL);
NotDone();
continue;
}
for (exv.Init(curEd, TopAbs_VERTEX); exv.More(); exv.Next()) {
const TopoDS_Vertex& edVer=TopoDS::Vertex(exv.Current());
rlist.Prepend(edVer);
nlist.Prepend(edVer);
if (!mapLocale.Contains(edVer)) {
Standard_Boolean notYetFound = Standard_True;
Standard_Real gap=BRep_Tool::Tolerance(edVer);
gp_Pnt pVer=BRep_Tool::Pnt(edVer);
for (itMS.Initialize(mapLocale); itMS.More(); itMS.Next()) {
notYetFound=Standard_True;
const TopoDS_Vertex& refVer=TopoDS::Vertex(itMS.Key());
gap +=BRep_Tool::Tolerance(refVer);
if (pVer.Distance(BRep_Tool::Pnt(TopoDS::Vertex(refVer))) <= gap) {
nlist.RemoveFirst();
nlist.Prepend(refVer.Oriented(edVer.Orientation()));
copEd=Standard_True;
notYetFound=Standard_False;
break;
}
}
if (notYetFound) mapLocale.Add(edVer);
}
}
if (copEd) {
TopoDS_Shape aLocalShape = curEd.EmptyCopied();
TopoDS_Edge newEd=TopoDS::Edge(aLocalShape);
// TopoDS_Edge newEd=TopoDS::Edge(curEd.EmptyCopied());
BB.Transfert(curEd, newEd);
TopTools_ListIteratorOfListOfShape itV(nlist);
for (; itV.More(); itV.Next()) {
BB.Add(newEd, itV.Value());
BB.Transfert(curEd, newEd, TopoDS::Vertex(rlist.First()), TopoDS::Vertex(itV.Value()));
rlist.RemoveFirst();
}
toAdd.Append(newEd);
}
else {
toAdd.Append(curEd);
}
}
}
if (!toAdd.IsEmpty()) {
TopoDS_Compound comp;
BB.MakeCompound(comp);
TopTools_MapIteratorOfMapOfOrientedShape itMOS;
TopTools_MapOfOrientedShape theEdges;
for (itList.Initialize(toAdd); itList.More(); itList.Next()) {
BB.Add(comp, itList.Value());
theEdges.Add(itList.Value());
}
TopTools_IndexedDataMapOfShapeListOfShape lesMeres;
TopExp::MapShapesAndAncestors(comp, TopAbs_VERTEX, TopAbs_EDGE, lesMeres);
TopoDS_Vertex vf, vl;
TopoDS_Shape theKey;
Standard_Boolean usedVertex;
Standard_Boolean closedEdge = Standard_False;
Standard_Integer vvInd, lastInd;
do {
if (!VL.IsNull() && lesMeres.Contains(VL)) {
if (!VF.IsNull()) closedEdge=VF.IsSame(VL);
usedVertex=Standard_True;
for (itList.Initialize(lesMeres.FindFromKey(VL)); itList.More(); itList.Next()) {
if (theEdges.Contains(itList.Value())) {
usedVertex=Standard_False;
theEdges.Remove(itList.Value());
TopExp::Vertices(TopoDS::Edge(itList.Value()), vf,vl);
if (vf.IsSame(VL)) {
BB.Add(myShape, itList.Value());
myVertices.Add(vl);
VL=vl;
}
else {
if (closedEdge) {
BB.Add(myShape, itList.Value());
VF=vf;
}
else {
BB.Add(myShape, itList.Value().Reversed());
vf.Reverse();
VL=vf;
}
myVertices.Add(vf);
}
}
}
if (usedVertex) {
lastInd=lesMeres.Extent();
vvInd=lesMeres.FindIndex(VL);
if (vvInd != lastInd) {
theKey=lesMeres.FindKey(lastInd);
nlist=lesMeres.FindFromIndex(lastInd);
}
lesMeres.RemoveLast();
if (vvInd != lastInd) {
lesMeres.Substitute(vvInd, theKey, nlist);
}
}
}
else if (!VF.IsNull() && lesMeres.Contains(VF)) {
usedVertex=Standard_True;
for (itList.Initialize(lesMeres.FindFromKey(VF)); itList.More(); itList.Next()) {
if (theEdges.Contains(itList.Value())) {
usedVertex=Standard_False;
theEdges.Remove(itList.Value());
TopExp::Vertices(TopoDS::Edge(itList.Value()), vf,vl);
if (vl.IsSame(VF)) {
BB.Add(myShape, itList.Value());
myVertices.Add(vf);
VF=vf;
}
else {
BB.Add(myShape, itList.Value().Reversed());
vl.Reverse();
myVertices.Add(vl);
VF=vl;
}
}
}
if (usedVertex) {
lastInd=lesMeres.Extent();
vvInd=lesMeres.FindIndex(VF);
if (vvInd != lastInd) {
theKey=lesMeres.FindKey(lastInd);
nlist=lesMeres.FindFromIndex(lastInd);
}
lesMeres.RemoveLast();
if (vvInd != lastInd) {
lesMeres.Substitute(vvInd, theKey, nlist);
}
}
}
else {
if (theEdges.Extent()>0) {
Standard_Boolean noCandidat=Standard_True;
for (itMOS.Initialize(theEdges); itMOS.More(); itMOS.Next()) {
TopExp::Vertices(TopoDS::Edge(itMOS.Key()), vf,vl);
if (myVertices.Contains(vl)) {
if (myError==BRepLib_WireDone) myError = BRepLib_NonManifoldWire;
BB.Add(myShape, itMOS.Key());
myVertices.Add(vf);
VF=vf;
noCandidat=Standard_False;
break;
}
else if (myVertices.Contains(vf)) {
if (myError==BRepLib_WireDone) myError = BRepLib_NonManifoldWire;
BB.Add(myShape, itMOS.Key());
myVertices.Add(vl);
VL=vl;
noCandidat=Standard_False;
break;
}
}
if (noCandidat) {
theEdges.Clear();
// Some Edges are not connected to first edge and the diagnosis is as follows
// but the "Maker" is Done() because otherwise it is not possible to return the constructed connected part...
myError=BRepLib_DisconnectedWire;
}
else theEdges.Remove(itMOS.Key());
}
}
} while (theEdges.Extent()>0);
}
//=======================================================================
//function : Accept
//purpose :
//=======================================================================
Standard_Boolean BRepLib_MakeWire::BRepLib_BndBoxVertexSelector::
Accept (const Standard_Integer& theObj)
{
if (theObj > myMapOfShape.Extent())
return Standard_False;
const TopoDS_Vertex& aV = TopoDS::Vertex(myMapOfShape(theObj));
if (theObj == myVInd)
return Standard_False;
gp_Pnt aVPnt = BRep_Tool::Pnt(aV);
Standard_Real aTolV = BRep_Tool::Tolerance(aV);
Standard_Real aL = myP.SquareDistance(aVPnt);
if (aL < Max(aTolV*aTolV, mySTol))
{
myResultInd.Append(theObj);
return Standard_True;
}
if (!VF.IsNull() && !VL.IsNull() && VF.IsSame(VL))
myShape.Closed(Standard_True);
Done();
return Standard_False;
}
//=======================================================================
//function : SetCurrentVertex
//purpose :
//=======================================================================
void BRepLib_MakeWire::BRepLib_BndBoxVertexSelector::
SetCurrentVertex (const gp_Pnt& theP, Standard_Real theTol,
Standard_Integer theVInd)
{
myP = theP;
myVBox.Add(myP);
myVBox.Enlarge(theTol);
mySTol = theTol*theTol;
myVInd = theVInd;
}
//=======================================================================
//function : CollectCoincidentVertices
//purpose :
//=======================================================================
void BRepLib_MakeWire::CollectCoincidentVertices(const TopTools_ListOfShape& theL,
NCollection_List<NCollection_List<TopoDS_Vertex>>& theGrVL)
{
TopTools_IndexedMapOfShape anAllV;
TopTools_ListIteratorOfListOfShape anItL;
TopTools_IndexedDataMapOfShapeListOfShape aMV2EL;
TopExp::MapShapesAndAncestors(myShape, TopAbs_VERTEX, TopAbs_EDGE, aMV2EL);
TopExp_Explorer exp;
for (anItL.Initialize(theL); anItL.More(); anItL.Next())
TopExp::MapShapesAndAncestors(anItL.Value(), TopAbs_VERTEX, TopAbs_EDGE, aMV2EL);
for (int i = 1; i <= aMV2EL.Extent(); i++)
if (aMV2EL(i).Extent() == 1)
anAllV.Add(aMV2EL.FindKey(i));
//aV2CV : vertex <-> its coincident vertices
NCollection_DataMap<TopoDS_Vertex, NCollection_Map<TopoDS_Vertex>> aV2CV;
NCollection_UBTree <Standard_Integer, Bnd_Box> aTree;
NCollection_UBTreeFiller <Standard_Integer, Bnd_Box> aTreeFiller (aTree);
NCollection_Map<TopoDS_Vertex> aNonGroupedV;
/// add vertices from anAllV to treefiller
for (Standard_Integer i = 1; i <= anAllV.Extent(); i++)
{
const TopoDS_Shape& aSh = anAllV(i);
Bnd_Box aBB;
BRepBndLib::Add(aSh, aBB);
aTreeFiller.Add(i, aBB);
}
aTreeFiller.Fill();
BRepLib_BndBoxVertexSelector aSelector(anAllV);
Standard_Integer aNbColl = 0;
NCollection_List<Standard_Integer>::Iterator itI;
for (Standard_Integer i = 1; i <= anAllV.Extent(); i++ )
{
const TopoDS_Vertex& aV = TopoDS::Vertex(anAllV(i));
if (myVertices.Contains(aV))
continue;
aSelector.SetCurrentVertex(BRep_Tool::Pnt(aV), BRep_Tool::Tolerance(aV), i );
aNbColl = aTree.Select(aSelector);
if (aNbColl > 0)
{
const NCollection_List<Standard_Integer>& aResInds = aSelector.GetResultInds();
NCollection_Map<TopoDS_Vertex>* aVM =
aV2CV.Bound(aV, NCollection_Map<TopoDS_Vertex>());
for (itI.Initialize(aResInds); itI.More(); itI.Next() )
{
const TopoDS_Vertex& aCV = TopoDS::Vertex(anAllV(itI.Value()));
aVM->Add(aCV);
if (myVertices.Contains(aCV))
{
if (aV2CV.IsBound(aCV))
aV2CV(aCV).Add(aV);
else
{
aV2CV.Bound(aCV, NCollection_Map<TopoDS_Vertex>())->Add(aV);
aNonGroupedV.Add(aCV);
}
}
}
if (!aVM->IsEmpty())
aNonGroupedV.Add(aV); //vertexes to be grouped; store only coincident vertices
}
aSelector.ClearResInds();
}
/// group the coincident vertices
NCollection_Map<TopoDS_Vertex>::Iterator itMV;
NCollection_List<TopoDS_Vertex> aStartV, aCurrentV, anOneGrV;
NCollection_List<TopoDS_Vertex>::Iterator itLV;
Standard_Boolean IsStartNewGroup = Standard_True;
while(!aNonGroupedV.IsEmpty() || !IsStartNewGroup)
//exit only if there are no nongrouped vertices
//and the last group are fully are constructed
{
if (IsStartNewGroup)
{
//start list of vertices is empty => append one from aNonGroupedV
// and remove it from it (i.e. mark as grouped)
itMV.Initialize(aNonGroupedV);
const TopoDS_Vertex& aCurV = itMV.Value();
aStartV.Append(aCurV);
aNonGroupedV.Remove(aCurV);
}
itLV.Init(aStartV);
for (;itLV.More();itLV.Next())
{
const TopoDS_Vertex& aSV = itLV.Value();
anOneGrV.Append(aSV);
itMV.Initialize(aV2CV(aSV));
for (;itMV.More();itMV.Next())
{
const TopoDS_Vertex& aCV = itMV.Value();
if (aNonGroupedV.Contains(aCV))
{
aCurrentV.Append(aCV);
aNonGroupedV.Remove(aCV);
}
}
}
aStartV.Clear();
aStartV.Append(aCurrentV);
IsStartNewGroup = aStartV.IsEmpty();
if (IsStartNewGroup && !anOneGrV.IsEmpty())
{
theGrVL.Append(anOneGrV);
anOneGrV.Clear();
}
}
}
//=======================================================================
//function : CreateNewVertices
//purpose :
//=======================================================================
void BRepLib_MakeWire::CreateNewVertices(const NCollection_List<NCollection_List<TopoDS_Vertex>>& theGrVL,
NCollection_DataMap<TopoDS_Vertex, TopoDS_Vertex>& theO2NV)
{
//map [old vertex => new vertex]
//note that already existing shape (i.e. the original ones)
//shouldnt be modified on the topological level
NCollection_List<NCollection_List<TopoDS_Vertex>>::Iterator itLLV;
NCollection_List<TopoDS_Vertex>::Iterator itLV;
BRep_Builder aBB;
itLLV.Init(theGrVL);
for (;itLLV.More();itLLV.Next())
{
TopoDS_Vertex aNewV;
NCollection_List<TopoDS_Shape> aValList;
const NCollection_List<TopoDS_Vertex>& aVal = itLLV.Value();
itLV.Initialize(aVal);
Standard_Real aNewTol = 0;
gp_Pnt aNewC;
for (;itLV.More();itLV.Next())
{
const TopoDS_Vertex& aVV = itLV.Value();
aValList.Append(aVV);
if (myVertices.Contains(aVV))
aNewV = aVV;
}
BRepLib::BoundingVertex(aValList, aNewC, aNewTol);
if (aNewV.IsNull())
{
//vertices from the original shape isnt found in this group
//create the new vertex
aNewV = BRepLib_MakeVertex(aNewC);
aBB.UpdateVertex(aNewV, aNewTol);
}
else
//update already existing vertex
aBB.UpdateVertex(aNewV, gp_Pnt(aNewC), aNewTol);
//fill the map of old->new vertices
itLV.Initialize(aVal);
for (;itLV.More();itLV.Next())
{
const TopoDS_Vertex& aVV = itLV.Value();
theO2NV.Bind(aVV, aNewV);
}
}
}
//=======================================================================
//function : CreateNewListOfEdges
//purpose :
//=======================================================================
void BRepLib_MakeWire::CreateNewListOfEdges(const TopTools_ListOfShape& theL,
const NCollection_DataMap<TopoDS_Vertex, TopoDS_Vertex>& theO2NV,
TopTools_ListOfShape& theNewEList)
{
///create the new list (theNewEList) from the input list L
Standard_Boolean IsNewEdge;
NCollection_List<TopoDS_Vertex> aVList;
TopExp_Explorer exp;
BRep_Builder aBB;
TopTools_ListIteratorOfListOfShape anItL;
for (anItL.Initialize(theL); anItL.More(); anItL.Next())
{
IsNewEdge = Standard_False;
aVList.Clear();
const TopoDS_Edge& aCE = TopoDS::Edge(anItL.Value());
exp.Init(aCE, TopAbs_VERTEX);
for (;exp.More(); exp.Next())
{
const TopoDS_Vertex& aVal = TopoDS::Vertex(exp.Current());
if (theO2NV.IsBound(aVal))
{
IsNewEdge = Standard_True;
//append the new vertex
aVList.Append(TopoDS::Vertex(theO2NV(aVal).Oriented(aVal.Orientation())));
}
else
aVList.Append(aVal);
}
if (IsNewEdge)
{
TopoDS_Shape NewE = aCE.EmptyCopied();
NCollection_List<TopoDS_Edge>::Iterator it(aVList);
for (; it.More(); it.Next())
aBB.Add(NewE, it.Value());
theNewEList.Append(TopoDS::Edge(NewE));
}
else
theNewEList.Append(aCE);
}
}

View File

@ -224,22 +224,56 @@ static Standard_Integer wire(Draw_Interpretor& di, Standard_Integer n, const cha
if (n < 3) return 1;
Standard_Integer i;
BRepBuilderAPI_MakeWire MW;
for (i = 2; i < n; i ++) {
TopoDS_Shape S = DBRep::Get(a[i]);
if (S.IsNull()) continue;
if (S.ShapeType() == TopAbs_EDGE)
MW.Add(TopoDS::Edge(S));
else if (S.ShapeType() == TopAbs_WIRE)
MW.Add(TopoDS::Wire(S));
else
continue;
Standard_Boolean IsUnsorted = !strcmp(a[2], "-unsorted");
if (!IsUnsorted)
for (i = 2; i < n; i ++) {
TopoDS_Shape S = DBRep::Get(a[i]);
if (S.IsNull()) continue;
if (S.ShapeType() == TopAbs_EDGE)
MW.Add(TopoDS::Edge(S));
else if (S.ShapeType() == TopAbs_WIRE)
MW.Add(TopoDS::Wire(S));
else
continue;
}
else
{
TopTools_ListOfShape aLE;
for (i = 3; i < n; i ++)
{
TopoDS_Shape S = DBRep::Get(a[i]);
TopExp_Explorer Exp(S, TopAbs_EDGE);
for (;Exp.More();Exp.Next())
{
const TopoDS_Edge& anE = TopoDS::Edge(Exp.Current());
if (!anE.IsNull())
aLE.Append(anE);
}
}
MW.Add(aLE);
}
if (!MW.IsDone()) {
//cout << "Wire not done" << endl;
di << "Wire not done\n";
return 0;
di << "Wire not done with an error:\n";
switch (MW.Error())
{
case BRepBuilderAPI_EmptyWire:
di << "BRepBuilderAPI_EmptyWire\n";
break;
case BRepBuilderAPI_DisconnectedWire:
di << "BRepBuilderAPI_DisconnectedWire\n";
break;
case BRepBuilderAPI_NonManifoldWire:
di << "BRepBuilderAPI_NonManifoldWire\n";
break;
default:
break;
}
}
DBRep::Set(a[1],MW);
else
DBRep::Set(a[1],MW);
return 0;
}
@ -1864,7 +1898,7 @@ void BRepTest::CurveCommands(Draw_Interpretor& theCommands)
polyvertex,g);
theCommands.Add("wire",
"wire wirename e1/w1 [e2/w2 ...]",__FILE__,
"wire wirename [-unsorted] e1/w1 [e2/w2 ...]",__FILE__,
wire,g);
theCommands.Add("profile",

View File

@ -33,6 +33,7 @@
#include <TopExp_Explorer.hxx>
#include <TopExp.hxx>
#include <TopoDS.hxx>
#include <TopTools_MapOfShape.hxx>
#include <ChFi2d_FilletAPI.hxx>
#include <ChFi2d_ChamferAPI.hxx>

View File

@ -2079,6 +2079,40 @@ static Standard_Integer OCC26270(Draw_Interpretor& theDI,
return 0;
}
#include "BRepBuilderAPI_MakeWire.hxx"
#include "BRepBuilderAPI_MakeEdge.hxx"
static Standard_Integer OCC27552(Draw_Interpretor&,
Standard_Integer,
const char ** )
{
BRep_Builder BB;
TopoDS_Vertex V1, V2, V3;
TopoDS_Edge E1, E2;
BB.MakeVertex(V1, gp_Pnt(0,0,0), 0.1);
BB.MakeVertex(V2, gp_Pnt(5,0,0), 0.1);
BB.MakeVertex(V3, gp_Pnt(10,0,0), 0.1);
E1 = BRepBuilderAPI_MakeEdge(V1, V2).Edge();
E2 = BRepBuilderAPI_MakeEdge(V2, V3).Edge();
BRepBuilderAPI_MakeWire MW;
MW.Add(E1);
MW.Add(E2);
TopoDS_Vertex V4, V5, V6, V7;
TopoDS_Edge E3, E4;
BB.MakeVertex(V4, gp_Pnt(10,0+0.05,0), 0.07);
BB.MakeVertex(V5, gp_Pnt(10,0-0.05,0), 0.07);
BB.MakeVertex(V6, gp_Pnt(10,0+2,0), 0.07);
BB.MakeVertex(V7, gp_Pnt(10,0-2,0), 0.07);
E3 = BRepBuilderAPI_MakeEdge(V4, V6).Edge();
E4 = BRepBuilderAPI_MakeEdge(V5, V7).Edge();
TopTools_ListOfShape LLE;
LLE.Append(E3);
LLE.Append(E4);
MW.Add(LLE);
TopoDS_Shape W = MW.Wire();
DBRep::Set("outw", W);
return 0;
}
void QABugs::Commands_20(Draw_Interpretor& theCommands) {
const char *group = "QABugs";
@ -2097,5 +2131,6 @@ void QABugs::Commands_20(Draw_Interpretor& theCommands) {
theCommands.Add ("OCC26747_3", "OCC26747_3 result", __FILE__, OCC26747_3, group);
theCommands.Add ("OCC27357", "OCC27357", __FILE__, OCC27357, group);
theCommands.Add("OCC26270", "OCC26270 shape result", __FILE__, OCC26270, group);
theCommands.Add ("OCC27552", "OCC27552", __FILE__, OCC27552, group);
return;
}

View File

@ -0,0 +1,43 @@
puts "=========="
puts "OCC27552"
puts "=========="
puts ""
#######################################
# Wire creation fails depending on the order of edges
#######################################
vertex v1 0 0 0
vertex v2 -100 0 100
vertex v3 100 0 100
vertex v4 0 0 -100
vertex v5 -100 0 -200
edge e1 v2 v1
edge e2 v3 v1
edge e3 v1 v4
edge e4 v4 v5
wire w1 -unsorted e1 e2 e3 e4
checkshape w1
checknbshapes w1 -vertex 5 -edge 4
wire w2 -unsorted e1 e4 e3 e2
checkshape w2
checknbshapes w2 -vertex 5 -edge 4
wire w3 -unsorted e2 e1 e3 e4
checkshape w3
checknbshapes w3 -vertex 5 -edge 4
wire w4 -unsorted e4 e3 e2 e1
checkshape w4
checknbshapes w4 -vertex 5 -edge 4
wire w5 -unsorted e4 e3 e1 e2
checkshape w5
checknbshapes w5 -vertex 5 -edge 4
wire w6 -unsorted e3 e1 e2 e4
checkshape w6
checknbshapes w6 -vertex 5 -edge 4

View File

@ -0,0 +1,43 @@
puts "=========="
puts "OCC27552"
puts "=========="
puts ""
#######################################
# Wire creation fails depending on the order of edges
#######################################
#4 vertices in row
point p1 0 0 0
point p2 5 0 0
point p4 10 0 0
point p3 12.5 0 0
point p5 20 0 0
vertex v1 p1
vertex v2 p2
vertex v3 p3
vertex v4 p4
vertex v5 p5
settolerance v1 3
settolerance v2 6
settolerance v3 3
settolerance v4 1
point d 0 10 0
vertex dv d
edge e1 v1 dv
edge e2 v2 dv
edge e3 v3 dv
edge e4 v4 dv
edge e5 dv v5
wire w2 -unsorted e1 e2 e3 e4 e5
checkshape w2
checknbshapes w2 -vertex 3 -edge 5
nexplode w2 v
checkmaxtol w2_2 -ref 9.875

View File

@ -0,0 +1,16 @@
puts "=========="
puts "OCC27552"
puts "=========="
puts ""
#######################################
# Wire creation fails depending on the order of edges
#######################################
pload QAcommands
OCC27552
#outw is output wire
checkshape outw
checknbshapes outw -vertex 5 -edge 4
nexplode outw v
checkmaxtol outw_4 -ref 0.12