1
0
mirror of https://git.dev.opencascade.org/repos/occt.git synced 2025-04-16 10:08:36 +03:00

0028903: BRepOffset_MakeOffset produces invalid shape (thickshell) in Intersection mode

1. Method BRepOffset_Tool::Inter3D is modified: now selection of proper edges is performed here, they are not concatenated into one edge if they go through a vertex on a boundary.

2. Method BRepOffset_Inter3d::ConnexIntByInt is modified: selection of edges is eliminated.

3. Method BRepOffset_Inter2d::ConnexIntByInt is corrected to be able to process seam edges correct.
This commit is contained in:
jgv 2018-01-23 17:54:06 +03:00 committed by bugmaster
parent f84d6446a7
commit bad76cfc7a
17 changed files with 636 additions and 684 deletions

View File

@ -88,6 +88,15 @@ void BRepAlgo_Loop::Init(const TopoDS_Face& F)
static void Bubble(const TopoDS_Edge& E, static void Bubble(const TopoDS_Edge& E,
TopTools_SequenceOfShape& Seq) TopTools_SequenceOfShape& Seq)
{ {
//Remove duplicates
for (Standard_Integer i = 1; i < Seq.Length(); i++)
for (Standard_Integer j = i+1; j <= Seq.Length(); j++)
if (Seq(i) == Seq(j))
{
Seq.Remove(j);
j--;
}
Standard_Boolean Invert = Standard_True; Standard_Boolean Invert = Standard_True;
Standard_Integer NbPoints = Seq.Length(); Standard_Integer NbPoints = Seq.Length();
Standard_Real U1,U2; Standard_Real U1,U2;

View File

@ -112,6 +112,28 @@ static TopoDS_Vertex CommonVertex(TopoDS_Edge& E1,
return V; return V;
} }
static Standard_Boolean IsOrientationChanged(TopTools_IndexedMapOfShape& theMap,
const TopoDS_Edge& theEdge)
{
Standard_Boolean IsOrChanged = Standard_False;
if (!theMap.Contains(theEdge))
theMap.Add(theEdge);
else
{
Standard_Integer anInd = theMap.FindIndex(theEdge);
const TopoDS_Shape& anEdge = theMap(anInd);
if (theEdge.Orientation() != anEdge.Orientation())
{
theMap.Substitute( anInd, theEdge );
IsOrChanged = Standard_True;
}
}
return IsOrChanged;
}
//======================================================================= //=======================================================================
//function : Store //function : Store
//purpose : Store the vertices <theLV> into AsDes for the edge <theEdge>. //purpose : Store the vertices <theLV> into AsDes for the edge <theEdge>.
@ -1480,6 +1502,9 @@ void BRepOffset_Inter2d::ConnexIntByInt
if (!wexp.More()) if (!wexp.More())
continue; // Protection from case when explorer does not contain edges. continue; // Protection from case when explorer does not contain edges.
CurE = FirstE = wexp.Current(); CurE = FirstE = wexp.Current();
TopTools_IndexedMapOfShape Edges;
Standard_Boolean ToReverse1, ToReverse2;
ToReverse1 = IsOrientationChanged(Edges, CurE);
while (!end) { while (!end) {
wexp.Next(); wexp.Next();
if (wexp.More()) { if (wexp.More()) {
@ -1490,6 +1515,8 @@ void BRepOffset_Inter2d::ConnexIntByInt
} }
if (CurE.IsSame(NextE)) continue; if (CurE.IsSame(NextE)) continue;
ToReverse2 = IsOrientationChanged(Edges, NextE);
TopoDS_Vertex Vref = CommonVertex(CurE, NextE); TopoDS_Vertex Vref = CommonVertex(CurE, NextE);
gp_Pnt Pref = BRep_Tool::Pnt(Vref); gp_Pnt Pref = BRep_Tool::Pnt(Vref);
@ -1515,6 +1542,9 @@ void BRepOffset_Inter2d::ConnexIntByInt
else if (Build.IsBound(NextE) && MES.IsBound(CEO)) { else if (Build.IsBound(NextE) && MES.IsBound(CEO)) {
NE1 = Build(NextE); NE1 = Build(NextE);
NE2 = MES(CEO); NE2 = MES(CEO);
Standard_Boolean Tmp = ToReverse1;
ToReverse1 = ToReverse2;
ToReverse2 = Tmp;
} }
else { else {
DoInter = 0; DoInter = 0;
@ -1526,9 +1556,17 @@ void BRepOffset_Inter2d::ConnexIntByInt
Standard_Boolean bCoincide; Standard_Boolean bCoincide;
TopExp_Explorer Exp1, Exp2; TopExp_Explorer Exp1, Exp2;
for (Exp1.Init(NE1, TopAbs_EDGE); Exp1.More(); Exp1.Next()) { for (Exp1.Init(NE1, TopAbs_EDGE); Exp1.More(); Exp1.Next()) {
const TopoDS_Edge& aE1 = TopoDS::Edge(Exp1.Current()); TopoDS_Edge aE1 = TopoDS::Edge(Exp1.Current());
for (Exp2.Init(NE2, TopAbs_EDGE); Exp2.More(); Exp2.Next()) { for (Exp2.Init(NE2, TopAbs_EDGE); Exp2.More(); Exp2.Next()) {
const TopoDS_Edge& aE2 = TopoDS::Edge(Exp2.Current()); TopoDS_Edge aE2 = TopoDS::Edge(Exp2.Current());
//Correct orientation of edges
if (ToReverse1)
aE1.Reverse();
if (ToReverse2)
aE2.Reverse();
//////////////////////////////
RefEdgeInter(FIO, BAsurf, aE1, aE2, AsDes2d, RefEdgeInter(FIO, BAsurf, aE1, aE2, AsDes2d,
Tol, Standard_True, Pref, theDMVV, bCoincide); Tol, Standard_True, Pref, theDMVV, bCoincide);
} }
@ -1553,6 +1591,7 @@ void BRepOffset_Inter2d::ConnexIntByInt
} }
} }
CurE = NextE; CurE = NextE;
ToReverse1 = ToReverse2;
} }
} }
} }

View File

@ -100,57 +100,6 @@ static void ExtentEdge(const TopoDS_Face& /*F*/,
} }
//=======================================================================
//function : SelectEdge
//purpose :
//=======================================================================
static void SelectEdge (const TopoDS_Shape& theS,
TopTools_ListOfShape& theLE)
{
Standard_Real aT1, aT2, aDist, aDistMin;
TopExp_Explorer aExp;
TopTools_ListIteratorOfListOfShape aIt;
GeomAPI_ProjectPointOnCurve aProjPC;
gp_Pnt aPE1, aPE2;
TopoDS_Edge aRE;
//
aDistMin = RealLast();
//
aIt.Initialize(theLE);
for (; aIt.More(); aIt.Next()) {
const TopoDS_Edge& aE = *(TopoDS_Edge*)&aIt.Value();
//
const Handle(Geom_Curve)& aC = BRep_Tool::Curve(aE, aT1, aT2);
//
aProjPC.Init(aC, aT1, aT2);
aPE1 = aC->Value(aT1);
aPE2 = aC->Value(aT2);
//
aDist = 0.;
aExp.Init(theS, TopAbs_VERTEX);
for (; aExp.More(); aExp.Next()) {
const TopoDS_Vertex& aV = *(TopoDS_Vertex*)&aExp.Current();
const gp_Pnt aP = BRep_Tool::Pnt(aV);
//
aProjPC.Perform(aP);
if (aProjPC.NbPoints()) {
aDist += aProjPC.LowerDistance();
}
else {
aDist += Min(aP.Distance(aPE1), aP.Distance(aPE2));
}
}
//
if (aDist < aDistMin) {
aDistMin = aDist;
aRE = aE;
}
}
//
theLE.Clear();
theLE.Append(aRE);
}
//======================================================================= //=======================================================================
//function : CompletInt //function : CompletInt
//purpose : //purpose :
@ -617,11 +566,6 @@ void BRepOffset_Inter3d::ConnexIntByInt
if (!IsDone(NF1,NF2)) { if (!IsDone(NF1,NF2)) {
TopTools_ListOfShape LInt1,LInt2; TopTools_ListOfShape LInt1,LInt2;
BRepOffset_Tool::Inter3D (NF1,NF2,LInt1,LInt2,CurSide,E,bEdge); BRepOffset_Tool::Inter3D (NF1,NF2,LInt1,LInt2,CurSide,E,bEdge);
if (LInt1.Extent() > 1) {
// intersection is in seceral edges (free sewing)
SelectEdge(aS, LInt1);
SelectEdge(aS, LInt2);
}
SetDone(NF1,NF2); SetDone(NF1,NF2);
if (!LInt1.IsEmpty()) { if (!LInt1.IsEmpty()) {
Store (NF1,NF2,LInt1,LInt2); Store (NF1,NF2,LInt1,LInt2);

View File

@ -320,6 +320,26 @@ static Standard_Real ComputeMaxDist(const gp_Pln& thePlane,
static void CorrectSolid(TopoDS_Solid& theSol, TopTools_ListOfShape& theSolList); static void CorrectSolid(TopoDS_Solid& theSol, TopTools_ListOfShape& theSolList);
//--------------------------------------------------------------------- //---------------------------------------------------------------------
static TopAbs_Orientation OrientationOfEdgeInFace(const TopoDS_Edge& theEdge,
const TopoDS_Face& theFace)
{
TopAbs_Orientation anOr = TopAbs_EXTERNAL;
TopExp_Explorer Explo(theFace, TopAbs_EDGE);
for (; Explo.More(); Explo.Next())
{
const TopoDS_Shape& anEdge = Explo.Current();
if (anEdge.IsSame(theEdge))
{
anOr = anEdge.Orientation();
break;
}
}
return anOr;
}
// //
static Standard_Boolean FindParameter(const TopoDS_Vertex& V, static Standard_Boolean FindParameter(const TopoDS_Vertex& V,
const TopoDS_Edge& E, const TopoDS_Edge& E,
@ -480,7 +500,7 @@ static void GetEdgePoints(const TopoDS_Edge& anEdge,
//======================================================================= //=======================================================================
static void FillContours(const TopoDS_Shape& aShape, static void FillContours(const TopoDS_Shape& aShape,
const BRepOffset_Analyse& Analyser, const BRepOffset_Analyse& Analyser,
TopTools_DataMapOfShapeListOfShape& Contours, TopTools_IndexedDataMapOfShapeListOfShape& Contours,
TopTools_DataMapOfShapeShape& MapEF) TopTools_DataMapOfShapeShape& MapEF)
{ {
TopTools_ListOfShape Edges; TopTools_ListOfShape Edges;
@ -533,7 +553,7 @@ static void FillContours(const TopoDS_Shape& aShape,
break; break;
} }
} }
Contours.Bind(StartVertex, aContour); Contours.Add(StartVertex, aContour);
} }
} }
@ -2497,18 +2517,17 @@ static void UpdateInitOffset (BRepAlgo_Image& myInitOffset,
//======================================================================= //=======================================================================
void BRepOffset_MakeOffset::MakeMissingWalls () void BRepOffset_MakeOffset::MakeMissingWalls ()
{ {
TopTools_DataMapOfShapeListOfShape Contours; //Start vertex + list of connected edges (free boundary) TopTools_IndexedDataMapOfShapeListOfShape Contours; //Start vertex + list of connected edges (free boundary)
TopTools_DataMapOfShapeShape MapEF; //Edges of contours: edge + face TopTools_DataMapOfShapeShape MapEF; //Edges of contours: edge + face
Standard_Real OffsetVal = Abs(myOffset); Standard_Real OffsetVal = Abs(myOffset);
FillContours(myShape, myAnalyse, Contours, MapEF); FillContours(myShape, myAnalyse, Contours, MapEF);
TopTools_DataMapIteratorOfDataMapOfShapeListOfShape iter(Contours); for (Standard_Integer ic = 1; ic <= Contours.Extent(); ic++)
for (; iter.More(); iter.Next())
{ {
TopoDS_Vertex StartVertex = TopoDS::Vertex(iter.Key()); TopoDS_Vertex StartVertex = TopoDS::Vertex(Contours.FindKey(ic));
TopoDS_Edge StartEdge; TopoDS_Edge StartEdge;
const TopTools_ListOfShape& aContour = iter.Value(); const TopTools_ListOfShape& aContour = Contours(ic);
TopTools_ListIteratorOfListOfShape itl(aContour); TopTools_ListIteratorOfListOfShape itl(aContour);
Standard_Boolean FirstStep = Standard_True; Standard_Boolean FirstStep = Standard_True;
TopoDS_Edge PrevEdge; TopoDS_Edge PrevEdge;
@ -2517,6 +2536,7 @@ void BRepOffset_MakeOffset::MakeMissingWalls ()
for (; itl.More(); itl.Next()) for (; itl.More(); itl.Next())
{ {
TopoDS_Edge anEdge = TopoDS::Edge(itl.Value()); TopoDS_Edge anEdge = TopoDS::Edge(itl.Value());
TopoDS_Face aFaceOfEdge = TopoDS::Face(MapEF(anEdge));
// Check for offset existence. // Check for offset existence.
if (!myInitOffsetEdge.HasImage(anEdge)) if (!myInitOffsetEdge.HasImage(anEdge))
@ -2856,6 +2876,12 @@ void BRepOffset_MakeOffset::MakeMissingWalls ()
} }
BRepLib::SameParameter(NewFace); BRepLib::SameParameter(NewFace);
BRepTools::Update(NewFace); BRepTools::Update(NewFace);
//Check orientation
TopAbs_Orientation anOr = OrientationOfEdgeInFace(anEdge, aFaceOfEdge);
TopAbs_Orientation OrInNewFace = OrientationOfEdgeInFace(anEdge, NewFace);
if (OrInNewFace != TopAbs::Reverse(anOr))
NewFace.Reverse();
///////////////////
myWalls.Append(NewFace); myWalls.Append(NewFace);
if (ArcOnV2) if (ArcOnV2)
{ {
@ -2938,7 +2964,9 @@ void BRepOffset_MakeOffset::MakeShells ()
TopTools_ListIteratorOfListOfShape it(R); TopTools_ListIteratorOfListOfShape it(R);
// //
for (; it.More(); it.Next()) { for (; it.More(); it.Next()) {
const TopoDS_Shape& aF = it.Value(); TopoDS_Shape aF = it.Value();
if (myThickening) //offsetted faces must change their orientations
aF.Reverse();
// //
TopTools_ListOfShape Image; TopTools_ListOfShape Image;
myImageOffset.LastImage(aF,Image); myImageOffset.LastImage(aF,Image);

View File

@ -19,6 +19,7 @@
#include <BndLib_Add3dCurve.hxx> #include <BndLib_Add3dCurve.hxx>
#include <BOPAlgo_PaveFiller.hxx> #include <BOPAlgo_PaveFiller.hxx>
#include <BOPDS_DS.hxx> #include <BOPDS_DS.hxx>
#include <BOPTools_AlgoTools.hxx>
#include <BOPTools_AlgoTools2D.hxx> #include <BOPTools_AlgoTools2D.hxx>
#include <BRep_CurveRepresentation.hxx> #include <BRep_CurveRepresentation.hxx>
#include <BRep_ListIteratorOfListOfCurveRepresentation.hxx> #include <BRep_ListIteratorOfListOfCurveRepresentation.hxx>
@ -39,7 +40,6 @@
#include <BRepLib_MakeFace.hxx> #include <BRepLib_MakeFace.hxx>
#include <BRepLib_MakePolygon.hxx> #include <BRepLib_MakePolygon.hxx>
#include <BRepLib_MakeVertex.hxx> #include <BRepLib_MakeVertex.hxx>
#include <BRepLib_MakeWire.hxx>
#include <BRepOffset_Analyse.hxx> #include <BRepOffset_Analyse.hxx>
#include <BRepOffset_Interval.hxx> #include <BRepOffset_Interval.hxx>
#include <BRepOffset_ListOfInterval.hxx> #include <BRepOffset_ListOfInterval.hxx>
@ -53,6 +53,7 @@
#include <ElSLib.hxx> #include <ElSLib.hxx>
#include <Extrema_ExtPC.hxx> #include <Extrema_ExtPC.hxx>
#include <Extrema_ExtPC2d.hxx> #include <Extrema_ExtPC2d.hxx>
#include <BRepExtrema_DistShapeShape.hxx>
#include <GCPnts_AbscissaPoint.hxx> #include <GCPnts_AbscissaPoint.hxx>
#include <GCPnts_QuasiUniformDeflection.hxx> #include <GCPnts_QuasiUniformDeflection.hxx>
#include <GCPnts_UniformAbscissa.hxx> #include <GCPnts_UniformAbscissa.hxx>
@ -137,7 +138,6 @@ const Standard_Real TheInfini = 1.e+7;
#ifdef DRAW #ifdef DRAW
#include <DBRep.hxx> #include <DBRep.hxx>
#include <Geom2d_Conic.hxx> #include <Geom2d_Conic.hxx>
#include <Geom_ElementarySurface.hxx>
#include <Geom_BoundedCurve.hxx> #include <Geom_BoundedCurve.hxx>
Standard_Boolean AffichInter = Standard_False; Standard_Boolean AffichInter = Standard_False;
static Standard_Integer NbNewEdges = 1; static Standard_Integer NbNewEdges = 1;
@ -158,6 +158,8 @@ static
TopTools_ListOfShape& theL1, TopTools_ListOfShape& theL1,
TopTools_ListOfShape& theL2); TopTools_ListOfShape& theL2);
static void UpdateVertexTolerances(const TopoDS_Face& theFace);
inline inline
Standard_Boolean IsInf(const Standard_Real theVal); Standard_Boolean IsInf(const Standard_Real theVal);
@ -198,7 +200,6 @@ TopAbs_Orientation BRepOffset_Tool::OriEdgeInFace (const TopoDS_Edge& E,
throw Standard_ConstructionError("BRepOffset_Tool::OriEdgeInFace"); throw Standard_ConstructionError("BRepOffset_Tool::OriEdgeInFace");
} }
//======================================================================= //=======================================================================
//function : FindPeriod //function : FindPeriod
//purpose : //purpose :
@ -799,6 +800,31 @@ void BRepOffset_Tool::PipeInter(const TopoDS_Face& F1,
} }
} }
//=======================================================================
//function : IsAutonomVertex
//purpose : Checks wether a vertex is "autonom" or not
//=======================================================================
static Standard_Boolean IsAutonomVertex(const TopoDS_Shape& theVertex,
const BOPDS_PDS& thePDS,
const TopoDS_Face& theFace1,
const TopoDS_Face& theFace2)
{
Standard_Integer nV = thePDS->Index(theVertex);
Standard_Integer nF [2];
nF[0] = thePDS->Index(theFace1);
nF[1] = thePDS->Index(theFace2);
for (Standard_Integer i = 0; i < 2; i++)
{
const BOPDS_FaceInfo& aFaceInfo = thePDS->FaceInfo(nF[i]);
const TColStd_MapOfInteger& IndMap = aFaceInfo.VerticesOn();
if (IndMap.Contains(nV))
return Standard_False;
}
return Standard_True;
}
//======================================================================= //=======================================================================
//function : IsAutonomVertex //function : IsAutonomVertex
@ -874,8 +900,7 @@ static Standard_Boolean IsAutonomVertex(const TopoDS_Shape& aVertex,
//======================================================================= //=======================================================================
static Standard_Boolean AreConnex(const TopoDS_Wire& W1, static Standard_Boolean AreConnex(const TopoDS_Wire& W1,
const TopoDS_Wire& W2, const TopoDS_Wire& W2)
const BOPDS_PDS& pDS)
{ {
TopoDS_Vertex V11, V12, V21, V22; TopoDS_Vertex V11, V12, V21, V22;
TopExp::Vertices( W1, V11, V12 ); TopExp::Vertices( W1, V11, V12 );
@ -883,55 +908,7 @@ static Standard_Boolean AreConnex(const TopoDS_Wire& W1,
if (V11.IsSame(V21) || V11.IsSame(V22) || if (V11.IsSame(V21) || V11.IsSame(V22) ||
V12.IsSame(V21) || V12.IsSame(V22)) V12.IsSame(V21) || V12.IsSame(V22))
{
Standard_Boolean isCV1 = V11.IsSame(V21) || V11.IsSame(V22);
Standard_Boolean isCV2 = V12.IsSame(V21) || V12.IsSame(V22);
if (isCV1 && !IsAutonomVertex(V11, pDS))
{
if (!isCV2)
return Standard_False;
if (!IsAutonomVertex(V12, pDS))
return Standard_False;
}
if (!isCV1 && !IsAutonomVertex(V12, pDS))
return Standard_False;
TopoDS_Vertex CV = (V11.IsSame(V21) || V11.IsSame(V22))? V11 : V12;
TopoDS_Edge E1, E2;
TopoDS_Iterator itw( W1 );
for (; itw.More(); itw.Next())
{
E1 = TopoDS::Edge(itw.Value());
TopoDS_Vertex V1, V2;
TopExp::Vertices( E1, V1, V2 );
if (V1.IsSame(CV) || V2.IsSame(CV))
break;
}
itw.Initialize( W2 );
E2 = TopoDS::Edge(itw.Value());
Standard_Real f, l;
Handle(Geom_Curve) C1 = BRep_Tool::Curve( E1, f, l );
if (C1->IsInstance(STANDARD_TYPE(Geom_TrimmedCurve)))
C1 = Handle(Geom_TrimmedCurve)::DownCast (C1)->BasisCurve();
Handle(Geom_Curve) C2 = BRep_Tool::Curve( E2, f, l );
if (C2->IsInstance(STANDARD_TYPE(Geom_TrimmedCurve)))
C2 = Handle(Geom_TrimmedCurve)::DownCast (C2)->BasisCurve();
if (C1->IsInstance(STANDARD_TYPE(Geom_Line)) &&
C2->IsInstance(STANDARD_TYPE(Geom_Line)))
{
Handle(Geom_Line) L1 = Handle(Geom_Line)::DownCast (C1);
gp_Ax1 Axis1 = L1->Position();
Handle(Geom_Line) L2 = Handle(Geom_Line)::DownCast (C2);
gp_Ax1 Axis2 = L2->Position();
if (! Axis1.IsParallel( Axis2, Precision::Angular() ))
return Standard_False;
}
return Standard_True; return Standard_True;
}
return Standard_False; return Standard_False;
} }
@ -1190,6 +1167,8 @@ static TopoDS_Edge Glue(const TopoDS_Edge& E1,
const Standard_Boolean addPCurve2, const Standard_Boolean addPCurve2,
const Standard_Real theGlueTol) const Standard_Real theGlueTol)
{ {
TopoDS_Edge newEdge;
Standard_Real Tol = 1.e-7; Standard_Real Tol = 1.e-7;
GeomAbs_Shape Continuity = GeomAbs_C1; GeomAbs_Shape Continuity = GeomAbs_C1;
Standard_Integer MaxDeg = 14; Standard_Integer MaxDeg = 14;
@ -1226,7 +1205,8 @@ static TopoDS_Edge Glue(const TopoDS_Edge& E1,
Handle(Geom_TrimmedCurve) TC1 = new Geom_TrimmedCurve( C1, first1, last1 ); Handle(Geom_TrimmedCurve) TC1 = new Geom_TrimmedCurve( C1, first1, last1 );
Handle(Geom_TrimmedCurve) TC2 = new Geom_TrimmedCurve( C2, first2, last2 ); Handle(Geom_TrimmedCurve) TC2 = new Geom_TrimmedCurve( C2, first2, last2 );
GeomConvert_CompCurveToBSplineCurve Concat( TC1 ); GeomConvert_CompCurveToBSplineCurve Concat( TC1 );
Concat.Add( TC2, theGlueTol, After ); if (!Concat.Add( TC2, theGlueTol, After ))
return newEdge;
newCurve = Concat.BSplineCurve(); newCurve = Concat.BSplineCurve();
if (newCurve->Continuity() < GeomAbs_C1) if (newCurve->Continuity() < GeomAbs_C1)
{ {
@ -1238,7 +1218,6 @@ static TopoDS_Edge Glue(const TopoDS_Edge& E1,
lparam = newCurve->LastParameter(); lparam = newCurve->LastParameter();
} }
TopoDS_Edge newEdge;
BRep_Builder BB; BRep_Builder BB;
if (IsCanonic) if (IsCanonic)
@ -1263,85 +1242,22 @@ static TopoDS_Edge Glue(const TopoDS_Edge& E1,
return newEdge; return newEdge;
} }
//=======================================================================
//function : FindNewVerticesOnBoundsOfFace
//purpose :
//=======================================================================
static void FindNewVerticesOnBoundsOfFace(const BOPDS_PDS& pDS,
const TopoDS_Face& aFace,
TopTools_DataMapOfShapeShape& VEmap)
{
TopTools_IndexedMapOfShape OldVertices;
TopExp::MapShapes( aFace, TopAbs_VERTEX, OldVertices );
BOPDS_ListIteratorOfListOfPaveBlock aItLPB;
TopoDS_Vertex V1, V2;
TopExp_Explorer Explo( aFace, TopAbs_EDGE );
for (; Explo.More(); Explo.Next()) {
const TopoDS_Shape& aE = Explo.Current();
Standard_Integer nE = pDS->Index(aE);
//
const BOPDS_ListOfPaveBlock& aLPB = pDS->PaveBlocks(nE);
aItLPB.Initialize(aLPB);
for (; aItLPB.More(); aItLPB.Next()) {
const Handle(BOPDS_PaveBlock)& aPB = aItLPB.Value();
const TopoDS_Edge& aESp = *(TopoDS_Edge*)&pDS->Shape(aPB->Edge());
//
TopExp::Vertices( aESp, V1, V2 );
if (!OldVertices.Contains( V1 )) {
VEmap.Bind( V1, aE );
}
//
if (!OldVertices.Contains( V2 )) {
VEmap.Bind( V2, aE );
}
}
}
}
//======================================================================= //=======================================================================
//function : CheckIntersFF //function : CheckIntersFF
//purpose : //purpose :
//======================================================================= //=======================================================================
static Standard_Boolean CheckIntersFF(const BOPDS_PDS& pDS, static void CheckIntersFF(const BOPDS_PDS& pDS,
const TopoDS_Edge& RefEdge, const TopoDS_Edge& RefEdge,
const TopoDS_Face& F1,
const TopoDS_Face& F2,
TopTools_IndexedMapOfShape& TrueEdges) TopTools_IndexedMapOfShape& TrueEdges)
{ {
Standard_Boolean isEl1 = Standard_False, isEl2 = Standard_False;
Standard_Boolean isPlane1 = Standard_False, isPlane2 = Standard_False;
Handle(Geom_Surface) aSurf = BRep_Tool::Surface(F1);
if (aSurf->IsInstance(STANDARD_TYPE(Geom_RectangularTrimmedSurface)))
aSurf = Handle(Geom_RectangularTrimmedSurface)::DownCast (aSurf)->BasisSurface();
if (aSurf->IsInstance(STANDARD_TYPE(Geom_Plane)))
isPlane1 = Standard_True;
else if (aSurf->IsKind(STANDARD_TYPE(Geom_ElementarySurface)))
isEl1 = Standard_True;
aSurf = BRep_Tool::Surface(F2);
if (aSurf->IsInstance(STANDARD_TYPE(Geom_RectangularTrimmedSurface)))
aSurf = Handle(Geom_RectangularTrimmedSurface)::DownCast (aSurf)->BasisSurface();
if (aSurf->IsInstance(STANDARD_TYPE(Geom_Plane)))
isPlane2 = Standard_True;
else if (aSurf->IsKind(STANDARD_TYPE(Geom_ElementarySurface)))
isEl2 = Standard_True;
if (isPlane1 || isPlane2)
return Standard_True;
if (isEl1 && isEl2)
return Standard_True;
BOPDS_VectorOfInterfFF& aFFs = pDS->InterfFF(); BOPDS_VectorOfInterfFF& aFFs = pDS->InterfFF();
Standard_Integer aNb = aFFs.Length(); Standard_Integer aNb = aFFs.Length();
Standard_Integer i, j, nbe = 0; Standard_Integer i, j, nbe = 0;
TopTools_SequenceOfShape Edges; TopoDS_Compound Edges;
BRep_Builder BB;
BB.MakeCompound(Edges);
for (i = 0; i < aNb; ++i) for (i = 0; i < aNb; ++i)
{ {
@ -1362,143 +1278,46 @@ static Standard_Boolean CheckIntersFF(const BOPDS_PDS& pDS,
const Handle(BOPDS_PaveBlock)& aPB = aPBIt.Value(); const Handle(BOPDS_PaveBlock)& aPB = aPBIt.Value();
Standard_Integer nSect = aPB->Edge(); Standard_Integer nSect = aPB->Edge();
const TopoDS_Edge& anEdge = *(TopoDS_Edge*)&pDS->Shape(nSect); const TopoDS_Edge& anEdge = *(TopoDS_Edge*)&pDS->Shape(nSect);
Edges.Append( anEdge ); BB.Add(Edges, anEdge);
nbe++; nbe++;
} }
} }
} }
if (nbe <= 1) if (nbe == 0)
return Standard_True; return;
//define tangents of RefEdge on start and on end TopTools_ListOfShape CompList;
BRepAdaptor_Curve cref(RefEdge); BOPTools_AlgoTools::MakeConnexityBlocks(Edges, TopAbs_VERTEX, TopAbs_EDGE, CompList);
gp_Vec RefTangFirst = cref.DN(cref.FirstParameter(), 1);
gp_Vec RefTangLast = cref.DN(cref.LastParameter(), 1);
//find the start edge and take it from Edges TopoDS_Shape NearestCompound;
TopoDS_Edge StartEdge; //, StartEonF1, StartEonF2, EndEonF1, EndEonF2; if (CompList.Extent() == 1)
NearestCompound = CompList.First();
TopTools_DataMapOfShapeShape VEmapF1, VEmapF2;
FindNewVerticesOnBoundsOfFace( pDS, F1, VEmapF1 );
FindNewVerticesOnBoundsOfFace( pDS, F2, VEmapF2 );
Standard_Real AngTol = 0.1;
Standard_Boolean V1onBound = Standard_False;
Standard_Boolean V2onBound = Standard_False;
TopoDS_Vertex V1, V2, Vcur;
gp_Vec TangFirst, TangLast, TangCur;
for (i = 1; i <= Edges.Length(); i++)
{
StartEdge = TopoDS::Edge(Edges(i));
TopExp::Vertices( StartEdge, V1, V2 );
V1onBound = VEmapF1.IsBound(V1) || VEmapF2.IsBound(V1); // && ?
V2onBound = VEmapF1.IsBound(V2) || VEmapF2.IsBound(V2); // && ?
if (V1onBound || V2onBound)
{
BRepAdaptor_Curve CE(StartEdge);
TangFirst = CE.DN( CE.FirstParameter(), 1 );
TangLast = CE.DN( CE.LastParameter(), 1 );
if (V1onBound)
{
if (TangFirst.IsParallel( RefTangFirst, AngTol ) ||
TangFirst.IsParallel( RefTangLast, AngTol ))
break; //start edged found
}
else if (V2onBound)
{
if (TangLast.IsParallel( RefTangLast, AngTol ) ||
TangLast.IsParallel( RefTangFirst, AngTol ))
break; //start edged found
}
}
}
if (i > Edges.Length()) //start edged not found
return Standard_False;
if (V1onBound && V2onBound)
{
if ((TangFirst.IsParallel(RefTangFirst,AngTol) && TangLast.IsParallel(RefTangLast,AngTol)) ||
(TangFirst.IsParallel(RefTangLast,AngTol) && TangLast.IsParallel(RefTangFirst,AngTol)))
{
TrueEdges.Add( Edges(i) );
return Standard_True;
}
else else
return Standard_False; {
} BRepAdaptor_Curve BAcurve(RefEdge);
gp_Pnt Pref = BAcurve.Value((BAcurve.FirstParameter()+BAcurve.LastParameter())/2);
TopoDS_Vertex Vref = BRepLib_MakeVertex(Pref);
Standard_Real MinDist = RealLast();
TopTools_ListIteratorOfListOfShape itl(CompList);
for (; itl.More(); itl.Next())
{
const TopoDS_Shape& aCompound = itl.Value();
//StartEonF1 = (V1onBound)? VEmapF1( V1 ) : VEmapF1( V2 ); BRepExtrema_DistShapeShape Projector(Vref, aCompound);
//StartEonF2 = (V1onBound)? VEmapF2( V1 ) : VEmapF2( V2 ); if (!Projector.IsDone() || Projector.NbSolution() == 0)
continue;
TrueEdges.Add( Edges(i) ); Standard_Real aDist = Projector.Value();
Edges.Remove(i); if (aDist < MinDist)
Vcur = (V1onBound)? V2 : V1;
TangCur = (V1onBound)? TangLast : TangFirst;
if (V2onBound)
TangCur.Reverse();
//build the chain from StartEdge till the opposite bound of face
for (;;)
{ {
TColStd_SequenceOfInteger Candidates; MinDist = aDist;
for (i = 1; i <= Edges.Length(); i++) NearestCompound = aCompound;
{
TopoDS_Edge anEdge = TopoDS::Edge(Edges(i));
TopExp::Vertices( anEdge, V1, V2 );
if (V1.IsSame(Vcur) || V2.IsSame(Vcur))
Candidates.Append(i);
}
if (Candidates.IsEmpty())
{
TrueEdges.Clear();
return Standard_False;
}
Standard_Integer minind = 1;
if (Candidates.Length() > 1)
{
Standard_Real MinAngle = RealLast();
for (i = 1; i <= Candidates.Length(); i++)
{
TopoDS_Edge anEdge = TopoDS::Edge(Edges(Candidates(i)));
TopExp::Vertices( anEdge, V1, V2 );
Standard_Boolean ConnectByFirst = (Vcur.IsSame(V1))? Standard_True : Standard_False;
BRepAdaptor_Curve CE(anEdge);
gp_Vec aTang = (ConnectByFirst)?
CE.DN( CE.FirstParameter(), 1 ) : CE.DN( CE.LastParameter(), 1 );
if (!ConnectByFirst)
aTang.Reverse();
Standard_Real anAngle = TangCur.Angle(aTang);
if (anAngle < MinAngle)
{
MinAngle = anAngle;
minind = i;
} }
} }
} }
TopoDS_Edge CurEdge = TopoDS::Edge(Edges(Candidates(minind)));
TrueEdges.Add( CurEdge );
Edges.Remove(Candidates(minind));
TopExp::Vertices( CurEdge, V1, V2 );
Standard_Boolean ConnectByFirst = (Vcur.IsSame(V1))? Standard_True : Standard_False;
Vcur = (ConnectByFirst)? V2 : V1;
BRepAdaptor_Curve CE(CurEdge);
TangCur = (ConnectByFirst)? CE.DN( CE.LastParameter(), 1 ) : CE.DN( CE.FirstParameter(), 1 );
if (!ConnectByFirst)
TangCur.Reverse();
//check if Vcur is on bounds of faces
if (VEmapF1.IsBound(Vcur) || VEmapF2.IsBound(Vcur))
break;
} //end of for (;;)
if (TangCur.IsParallel( RefTangFirst, AngTol ) || TopExp::MapShapes(NearestCompound, TopAbs_EDGE, TrueEdges);
TangCur.IsParallel( RefTangLast, AngTol ))
return Standard_True;
TrueEdges.Clear();
return Standard_False;
} }
//======================================================================= //=======================================================================
@ -1513,32 +1332,42 @@ static TopoDS_Edge AssembleEdge(const BOPDS_PDS& pDS,
const Standard_Boolean addPCurve2, const Standard_Boolean addPCurve2,
const TopTools_SequenceOfShape& EdgesForConcat) const TopTools_SequenceOfShape& EdgesForConcat)
{ {
TopoDS_Edge NullEdge;
TopoDS_Edge CurEdge = TopoDS::Edge( EdgesForConcat(1) ); TopoDS_Edge CurEdge = TopoDS::Edge( EdgesForConcat(1) );
Standard_Real aGlueTol = Precision::Confusion(); Standard_Real aGlueTol = Precision::Confusion();
for (Standard_Integer j = 2; j <= EdgesForConcat.Length(); j++) for (Standard_Integer j = 2; j <= EdgesForConcat.Length(); j++)
{ {
TopoDS_Edge anEdge = TopoDS::Edge( EdgesForConcat(j) ); TopoDS_Edge anEdge = TopoDS::Edge( EdgesForConcat(j) );
Standard_Boolean After = Standard_False; Standard_Boolean After = Standard_False;
TopoDS_Vertex Vfirst, Vlast; TopoDS_Vertex Vfirst, Vlast;
if (AreClosed( CurEdge, anEdge )) Standard_Boolean AreClosedWire = AreClosed( CurEdge, anEdge );
if (AreClosedWire)
{ {
TopoDS_Vertex V1, V2; TopoDS_Vertex V1, V2;
TopExp::Vertices( CurEdge, V1, V2 ); TopExp::Vertices( CurEdge, V1, V2 );
if (IsAutonomVertex( V1, pDS )) Standard_Boolean IsAutonomV1 = IsAutonomVertex( V1, pDS, F1, F2 );
Standard_Boolean IsAutonomV2 = IsAutonomVertex( V2, pDS, F1, F2 );
if (IsAutonomV1)
{ {
After = Standard_False; After = Standard_False;
Vfirst = Vlast = V2; Vfirst = Vlast = V2;
} }
else else if (IsAutonomV2)
{ {
After = Standard_True; After = Standard_True;
Vfirst = Vlast = V1; Vfirst = Vlast = V1;
} }
else
return NullEdge;
} }
else else
{ {
TopoDS_Vertex CV, V11, V12, V21, V22; TopoDS_Vertex CV, V11, V12, V21, V22;
TopExp::CommonVertex( CurEdge, anEdge, CV ); TopExp::CommonVertex( CurEdge, anEdge, CV );
Standard_Boolean IsAutonomCV = IsAutonomVertex( CV, pDS, F1, F2 );
if (IsAutonomCV)
{
aGlueTol = BRep_Tool::Tolerance(CV); aGlueTol = BRep_Tool::Tolerance(CV);
TopExp::Vertices( CurEdge, V11, V12 ); TopExp::Vertices( CurEdge, V11, V12 );
TopExp::Vertices( anEdge, V21, V22 ); TopExp::Vertices( anEdge, V21, V22 );
@ -1562,10 +1391,16 @@ static TopoDS_Edge AssembleEdge(const BOPDS_PDS& pDS,
Vfirst = V11; Vfirst = V11;
Vlast = V21; Vlast = V21;
} }
}
else
return NullEdge;
} //end of else (open wire) } //end of else (open wire)
TopoDS_Edge NewEdge = Glue(CurEdge, anEdge, Vfirst, Vlast, After, TopoDS_Edge NewEdge = Glue(CurEdge, anEdge, Vfirst, Vlast, After,
F1, addPCurve1, F2, addPCurve2, aGlueTol); F1, addPCurve1, F2, addPCurve2, aGlueTol);
if (NewEdge.IsNull())
return NullEdge;
else
CurEdge = NewEdge; CurEdge = NewEdge;
} //end of for (Standard_Integer j = 2; j <= EdgesForConcat.Length(); j++) } //end of for (Standard_Integer j = 2; j <= EdgesForConcat.Length(); j++)
@ -1610,49 +1445,38 @@ void BRepOffset_Tool::Inter3D(const TopoDS_Face& F1,
} }
} }
} }
//
TopoDS_Face cpF1=F1;
TopoDS_Face cpF2=F2;
// create 3D curves on faces // create 3D curves on faces
BRepLib::BuildCurves3d(cpF1); BRepLib::BuildCurves3d(F1);
BRepLib::BuildCurves3d(cpF2); BRepLib::BuildCurves3d(F2);
UpdateVertexTolerances(F1);
UpdateVertexTolerances(F2);
BOPAlgo_PaveFiller aPF1, aPF2; BOPAlgo_PaveFiller aPF;
TopTools_ListOfShape aLS; TopTools_ListOfShape aLS;
aLS.Append(cpF1); aLS.Append(F1);
aLS.Append(cpF2); aLS.Append(F2);
aPF1.SetArguments(aLS); aPF.SetArguments(aLS);
// //
aPF1.Perform(); aPF.Perform();
BOPAlgo_PaveFiller *pPF = &aPF1;
aLS.Clear();
TopTools_IndexedMapOfShape TrueEdges; TopTools_IndexedMapOfShape TrueEdges;
if (IsRefEdgeDefined && !CheckIntersFF( pPF->PDS(), RefEdge, cpF1, cpF2, TrueEdges )) if (IsRefEdgeDefined)
{ CheckIntersFF( aPF.PDS(), RefEdge, TrueEdges );
cpF1 = F2;
cpF2 = F1;
aLS.Append(cpF1);
aLS.Append(cpF2);
aPF2.SetArguments(aLS);
aPF2.Perform();
pPF = &aPF2;
CheckIntersFF( pPF->PDS(), RefEdge, cpF1, cpF2, TrueEdges );
}
Standard_Boolean addPCurve1 = 1; Standard_Boolean addPCurve1 = 1;
Standard_Boolean addPCurve2 = 1; Standard_Boolean addPCurve2 = 1;
const BOPDS_PDS& pDS = pPF->PDS(); const BOPDS_PDS& pDS = aPF.PDS();
BOPDS_VectorOfInterfFF& aFFs=pDS->InterfFF(); BOPDS_VectorOfInterfFF& aFFs=pDS->InterfFF();
Standard_Integer aNb = aFFs.Length(); Standard_Integer aNb = aFFs.Length();
Standard_Integer i = 0, j = 0, k; Standard_Integer i = 0, j = 0, k;
// Store Result // Store Result
L1.Clear(); L2.Clear(); L1.Clear(); L2.Clear();
TopAbs_Orientation O1,O2; TopAbs_Orientation O1,O2;
BRep_Builder BB;
// //
const Handle(IntTools_Context)& aContext = pPF->Context(); const Handle(IntTools_Context)& aContext = aPF.Context();
// //
for (i = 0; i < aNb; i++) { for (i = 0; i < aNb; i++) {
BOPDS_InterfFF& aFFi=aFFs(i); BOPDS_InterfFF& aFFi=aFFs(i);
@ -1681,39 +1505,38 @@ void BRepOffset_Tool::Inter3D(const TopoDS_Face& F1,
if(!aC3DE.IsNull()) if(!aC3DE.IsNull())
aC3DETrim = new Geom_TrimmedCurve(aC3DE, f, l); aC3DETrim = new Geom_TrimmedCurve(aC3DE, f, l);
BRep_Builder aBB;
Standard_Real aTolEdge = BRep_Tool::Tolerance(anEdge); Standard_Real aTolEdge = BRep_Tool::Tolerance(anEdge);
if (!BOPTools_AlgoTools2D::HasCurveOnSurface(anEdge, cpF1)) { if (!BOPTools_AlgoTools2D::HasCurveOnSurface(anEdge, F1)) {
Handle(Geom2d_Curve) aC2d = aBC.Curve().FirstCurve2d(); Handle(Geom2d_Curve) aC2d = aBC.Curve().FirstCurve2d();
if(!aC3DETrim.IsNull()) { if(!aC3DETrim.IsNull()) {
Handle(Geom2d_Curve) aC2dNew; Handle(Geom2d_Curve) aC2dNew;
if(aC3DE->IsPeriodic()) { if(aC3DE->IsPeriodic()) {
BOPTools_AlgoTools2D::AdjustPCurveOnFace(cpF1, f, l, aC2d, aC2dNew, aContext); BOPTools_AlgoTools2D::AdjustPCurveOnFace(F1, f, l, aC2d, aC2dNew, aContext);
} }
else { else {
BOPTools_AlgoTools2D::AdjustPCurveOnFace(cpF1, aC3DETrim, aC2d, aC2dNew, aContext); BOPTools_AlgoTools2D::AdjustPCurveOnFace(F1, aC3DETrim, aC2d, aC2dNew, aContext);
} }
aC2d = aC2dNew; aC2d = aC2dNew;
} }
aBB.UpdateEdge(anEdge, aC2d, cpF1, aTolEdge); BB.UpdateEdge(anEdge, aC2d, F1, aTolEdge);
} }
if (!BOPTools_AlgoTools2D::HasCurveOnSurface(anEdge, cpF2)) { if (!BOPTools_AlgoTools2D::HasCurveOnSurface(anEdge, F2)) {
Handle(Geom2d_Curve) aC2d = aBC.Curve().SecondCurve2d(); Handle(Geom2d_Curve) aC2d = aBC.Curve().SecondCurve2d();
if(!aC3DETrim.IsNull()) { if(!aC3DETrim.IsNull()) {
Handle(Geom2d_Curve) aC2dNew; Handle(Geom2d_Curve) aC2dNew;
if(aC3DE->IsPeriodic()) { if(aC3DE->IsPeriodic()) {
BOPTools_AlgoTools2D::AdjustPCurveOnFace(cpF2, f, l, aC2d, aC2dNew, aContext); BOPTools_AlgoTools2D::AdjustPCurveOnFace(F2, f, l, aC2d, aC2dNew, aContext);
} }
else { else {
BOPTools_AlgoTools2D::AdjustPCurveOnFace(cpF2, aC3DETrim, aC2d, aC2dNew, aContext); BOPTools_AlgoTools2D::AdjustPCurveOnFace(F2, aC3DETrim, aC2d, aC2dNew, aContext);
} }
aC2d = aC2dNew; aC2d = aC2dNew;
} }
aBB.UpdateEdge(anEdge, aC2d, cpF2, aTolEdge); BB.UpdateEdge(anEdge, aC2d, F2, aTolEdge);
} }
OrientSection (anEdge, F1, F2, O1, O2); OrientSection (anEdge, F1, F2, O1, O2);
@ -1740,7 +1563,7 @@ void BRepOffset_Tool::Inter3D(const TopoDS_Face& F1,
Standard_Real aSameParTol = Precision::Confusion(); Standard_Real aSameParTol = Precision::Confusion();
Standard_Boolean isEl1 = Standard_False, isEl2 = Standard_False; Standard_Boolean isEl1 = Standard_False, isEl2 = Standard_False;
Handle(Geom_Surface) aSurf = BRep_Tool::Surface(cpF1); Handle(Geom_Surface) aSurf = BRep_Tool::Surface(F1);
if (aSurf->IsInstance(STANDARD_TYPE(Geom_RectangularTrimmedSurface))) if (aSurf->IsInstance(STANDARD_TYPE(Geom_RectangularTrimmedSurface)))
aSurf = Handle(Geom_RectangularTrimmedSurface)::DownCast (aSurf)->BasisSurface(); aSurf = Handle(Geom_RectangularTrimmedSurface)::DownCast (aSurf)->BasisSurface();
if (aSurf->IsInstance(STANDARD_TYPE(Geom_Plane))) if (aSurf->IsInstance(STANDARD_TYPE(Geom_Plane)))
@ -1748,7 +1571,7 @@ void BRepOffset_Tool::Inter3D(const TopoDS_Face& F1,
else if (aSurf->IsKind(STANDARD_TYPE(Geom_ElementarySurface))) else if (aSurf->IsKind(STANDARD_TYPE(Geom_ElementarySurface)))
isEl1 = Standard_True; isEl1 = Standard_True;
aSurf = BRep_Tool::Surface(cpF2); aSurf = BRep_Tool::Surface(F2);
if (aSurf->IsInstance(STANDARD_TYPE(Geom_RectangularTrimmedSurface))) if (aSurf->IsInstance(STANDARD_TYPE(Geom_RectangularTrimmedSurface)))
aSurf = Handle(Geom_RectangularTrimmedSurface)::DownCast (aSurf)->BasisSurface(); aSurf = Handle(Geom_RectangularTrimmedSurface)::DownCast (aSurf)->BasisSurface();
if (aSurf->IsInstance(STANDARD_TYPE(Geom_Plane))) if (aSurf->IsInstance(STANDARD_TYPE(Geom_Plane)))
@ -1764,13 +1587,16 @@ void BRepOffset_Tool::Inter3D(const TopoDS_Face& F1,
{ {
for (i = TrueEdges.Extent(); i >= 1; i--) for (i = TrueEdges.Extent(); i >= 1; i--)
EdgesForConcat.Append( TrueEdges(i) ); EdgesForConcat.Append( TrueEdges(i) );
TopoDS_Edge theEdge = TopoDS_Edge AssembledEdge =
AssembleEdge( pDS, cpF1, cpF2, addPCurve1, addPCurve2, EdgesForConcat ); AssembleEdge( pDS, F1, F2, addPCurve1, addPCurve2, EdgesForConcat );
eseq.Append(theEdge); if (AssembledEdge.IsNull())
for (i = TrueEdges.Extent(); i >= 1; i--)
eseq.Append( TrueEdges(i) );
else
eseq.Append(AssembledEdge);
} }
else else
{ {
TopTools_SequenceOfShape wseq; TopTools_SequenceOfShape wseq;
TopTools_SequenceOfShape edges; TopTools_SequenceOfShape edges;
TopTools_ListIteratorOfListOfShape itl(L1); TopTools_ListIteratorOfListOfShape itl(L1);
@ -1779,12 +1605,14 @@ void BRepOffset_Tool::Inter3D(const TopoDS_Face& F1,
while (!edges.IsEmpty()) while (!edges.IsEmpty())
{ {
TopoDS_Edge anEdge = TopoDS::Edge( edges.First() ); TopoDS_Edge anEdge = TopoDS::Edge( edges.First() );
TopoDS_Wire aWire = BRepLib_MakeWire( anEdge ), resWire; TopoDS_Wire aWire, resWire;
BB.MakeWire(aWire);
BB.Add( aWire, anEdge );
TColStd_SequenceOfInteger Candidates; TColStd_SequenceOfInteger Candidates;
for (k = 1; k <= wseq.Length(); k++) for (k = 1; k <= wseq.Length(); k++)
{ {
resWire = TopoDS::Wire(wseq(k)); resWire = TopoDS::Wire(wseq(k));
if (AreConnex( resWire, aWire, pDS )) if (AreConnex( resWire, aWire ))
{ {
Candidates.Append( 1 ); Candidates.Append( 1 );
break; break;
@ -1800,10 +1628,10 @@ void BRepOffset_Tool::Inter3D(const TopoDS_Face& F1,
for (j = 2; j <= edges.Length(); j++) for (j = 2; j <= edges.Length(); j++)
{ {
anEdge = TopoDS::Edge( edges(j) ); anEdge = TopoDS::Edge( edges(j) );
//if (anEdge.IsSame(edges(Candidates.First()))) aWire.Nullify();
//continue; BB.MakeWire(aWire);
aWire = BRepLib_MakeWire( anEdge ); BB.Add( aWire, anEdge );
if (AreConnex( resWire, aWire, pDS )) if (AreConnex( resWire, aWire ))
Candidates.Append( j ); Candidates.Append( j );
} }
Standard_Integer minind = 1; Standard_Integer minind = 1;
@ -1821,8 +1649,8 @@ void BRepOffset_Tool::Inter3D(const TopoDS_Face& F1,
} }
} }
} }
TopoDS_Wire NewWire = BRepLib_MakeWire( resWire, TopoDS::Edge(edges(Candidates(minind))) ); BB.Add( resWire, TopoDS::Edge(edges(Candidates(minind))) );
wseq(k) = NewWire; wseq(k) = resWire;
edges.Remove(Candidates(minind)); edges.Remove(Candidates(minind));
} }
} //end of while (!edges.IsEmpty()) } //end of while (!edges.IsEmpty())
@ -1905,9 +1733,13 @@ void BRepOffset_Tool::Inter3D(const TopoDS_Face& F1,
aLocalEdgesForConcat.Append( Wexp.Current() ); aLocalEdgesForConcat.Append( Wexp.Current() );
} }
TopoDS_Edge theEdge = TopoDS_Edge AssembledEdge =
AssembleEdge( pDS, cpF1, cpF2, addPCurve1, addPCurve2, aLocalEdgesForConcat ); AssembleEdge( pDS, F1, F2, addPCurve1, addPCurve2, aLocalEdgesForConcat );
eseq.Append( theEdge ); if (AssembledEdge.IsNull())
for (j = aLocalEdgesForConcat.Length(); j >= 1; j--)
eseq.Append( aLocalEdgesForConcat(j) );
else
eseq.Append( AssembledEdge );
} }
} //end of else (when TrueEdges is empty) } //end of else (when TrueEdges is empty)
@ -1917,6 +1749,7 @@ void BRepOffset_Tool::Inter3D(const TopoDS_Face& F1,
L2.Clear(); L2.Clear();
for (i = 1; i <= eseq.Length(); i++) for (i = 1; i <= eseq.Length(); i++)
{ {
TopoDS_Shape aShape = eseq(i);
TopoDS_Edge anEdge = TopoDS::Edge(eseq(i)); TopoDS_Edge anEdge = TopoDS::Edge(eseq(i));
BRepLib::SameParameter(anEdge, aSameParTol, Standard_True); BRepLib::SameParameter(anEdge, aSameParTol, Standard_True);
Standard_Real EdgeTol = BRep_Tool::Tolerance(anEdge); Standard_Real EdgeTol = BRep_Tool::Tolerance(anEdge);
@ -4155,3 +3988,50 @@ Standard_Boolean IsInf(const Standard_Real theVal)
{ {
return (theVal > TheInfini*0.9); return (theVal > TheInfini*0.9);
} }
static void UpdateVertexTolerances(const TopoDS_Face& theFace)
{
BRep_Builder BB;
TopTools_IndexedDataMapOfShapeListOfShape VEmap;
TopExp::MapShapesAndAncestors(theFace, TopAbs_VERTEX, TopAbs_EDGE, VEmap);
for (Standard_Integer i = 1; i <= VEmap.Extent(); i++)
{
const TopoDS_Vertex& aVertex = TopoDS::Vertex(VEmap.FindKey(i));
const TopTools_ListOfShape& Elist = VEmap(i);
gp_Pnt PntVtx = BRep_Tool::Pnt(aVertex);
TopTools_ListIteratorOfListOfShape itl(Elist);
for (; itl.More(); itl.Next())
{
const TopoDS_Edge& anEdge = TopoDS::Edge(itl.Value());
TopoDS_Vertex V1, V2;
TopExp::Vertices(anEdge, V1, V2);
Standard_Real fpar, lpar;
BRep_Tool::Range(anEdge, fpar, lpar);
Standard_Real aParam = (V1.IsSame(aVertex))? fpar : lpar;
if (!BRep_Tool::Degenerated(anEdge))
{
BRepAdaptor_Curve BAcurve(anEdge);
gp_Pnt aPnt = BAcurve.Value(aParam);
Standard_Real aDist = PntVtx.Distance(aPnt);
BB.UpdateVertex(aVertex, aDist);
if (V1.IsSame(V2))
{
aPnt = BAcurve.Value(lpar);
aDist = PntVtx.Distance(aPnt);
BB.UpdateVertex(aVertex, aDist);
}
}
BRepAdaptor_Curve BAcurveonsurf(anEdge, theFace);
gp_Pnt aPnt = BAcurveonsurf.Value(aParam);
Standard_Real aDist = PntVtx.Distance(aPnt);
BB.UpdateVertex(aVertex, aDist);
if (V1.IsSame(V2))
{
aPnt = BAcurveonsurf.Value(lpar);
aDist = PntVtx.Distance(aPnt);
BB.UpdateVertex(aVertex, aDist);
}
}
}
}

View File

@ -60,7 +60,11 @@ public:
//! <E> is a section between <F1> and <F2>. Computes //! <E> is a section between <F1> and <F2>. Computes
//! <O1> the orientation of <E> in <F1> influenced by <F2>. //! <O1> the orientation of <E> in <F1> influenced by <F2>.
//! idem for <O2>. //! idem for <O2>.
Standard_EXPORT static void OrientSection (const TopoDS_Edge& E, const TopoDS_Face& F1, const TopoDS_Face& F2, TopAbs_Orientation& O1, TopAbs_Orientation& O2); Standard_EXPORT static void OrientSection (const TopoDS_Edge& E,
const TopoDS_Face& F1,
const TopoDS_Face& F2,
TopAbs_Orientation& O1,
TopAbs_Orientation& O2);
//! Looks for the common Vertices and Edges between faces <theF1> and <theF2>.<br> //! Looks for the common Vertices and Edges between faces <theF1> and <theF2>.<br>
//! Returns TRUE if common shapes have been found.<br> //! Returns TRUE if common shapes have been found.<br>
@ -83,42 +87,85 @@ public:
//! edges solution are stored in <LInt1> with the //! edges solution are stored in <LInt1> with the
//! orientation on <F1>, the sames edges are stored in //! orientation on <F1>, the sames edges are stored in
//! <Lint2> with the orientation on <F2>. //! <Lint2> with the orientation on <F2>.
Standard_EXPORT static void Inter3D (const TopoDS_Face& F1, const TopoDS_Face& F2, TopTools_ListOfShape& LInt1, TopTools_ListOfShape& LInt2, const TopAbs_State Side, const TopoDS_Edge& RefEdge, const Standard_Boolean IsRefEdgeDefined = Standard_False); Standard_EXPORT static void Inter3D (const TopoDS_Face& F1,
const TopoDS_Face& F2,
TopTools_ListOfShape& LInt1,
TopTools_ListOfShape& LInt2,
const TopAbs_State Side,
const TopoDS_Edge& RefEdge,
const Standard_Boolean IsRefEdgeDefined = Standard_False);
//! Find if the edges <Edges> of the face <F2> are on //! Find if the edges <Edges> of the face <F2> are on
//! the face <F1>. //! the face <F1>.
//! Set in <LInt1> <LInt2> the updated edges. //! Set in <LInt1> <LInt2> the updated edges.
//! If all the edges are computed, returns true. //! If all the edges are computed, returns true.
Standard_EXPORT static Standard_Boolean TryProject (const TopoDS_Face& F1, const TopoDS_Face& F2, const TopTools_ListOfShape& Edges, TopTools_ListOfShape& LInt1, TopTools_ListOfShape& LInt2, const TopAbs_State Side, const Standard_Real TolConf); Standard_EXPORT static Standard_Boolean TryProject (const TopoDS_Face& F1,
const TopoDS_Face& F2,
const TopTools_ListOfShape& Edges,
TopTools_ListOfShape& LInt1,
TopTools_ListOfShape& LInt2,
const TopAbs_State Side,
const Standard_Real TolConf);
Standard_EXPORT static void PipeInter (const TopoDS_Face& F1, const TopoDS_Face& F2, TopTools_ListOfShape& LInt1, TopTools_ListOfShape& LInt2, const TopAbs_State Side); Standard_EXPORT static void PipeInter (const TopoDS_Face& F1,
const TopoDS_Face& F2,
TopTools_ListOfShape& LInt1,
TopTools_ListOfShape& LInt2,
const TopAbs_State Side);
Standard_EXPORT static void Inter2d (const TopoDS_Face& F, const TopoDS_Edge& E1, const TopoDS_Edge& E2, TopTools_ListOfShape& LV, const Standard_Real Tol); Standard_EXPORT static void Inter2d (const TopoDS_Face& F,
const TopoDS_Edge& E1,
const TopoDS_Edge& E2,
TopTools_ListOfShape& LV,
const Standard_Real Tol);
Standard_EXPORT static void InterOrExtent (const TopoDS_Face& F1, const TopoDS_Face& F2, TopTools_ListOfShape& LInt1, TopTools_ListOfShape& LInt2, const TopAbs_State Side); Standard_EXPORT static void InterOrExtent (const TopoDS_Face& F1,
const TopoDS_Face& F2,
TopTools_ListOfShape& LInt1,
TopTools_ListOfShape& LInt2,
const TopAbs_State Side);
Standard_EXPORT static void CheckBounds (const TopoDS_Face& F, const BRepOffset_Analyse& Analyse, Standard_Boolean& enlargeU, Standard_Boolean& enlargeVfirst, Standard_Boolean& enlargeVlast); Standard_EXPORT static void CheckBounds (const TopoDS_Face& F,
const BRepOffset_Analyse& Analyse,
Standard_Boolean& enlargeU,
Standard_Boolean& enlargeVfirst,
Standard_Boolean& enlargeVlast);
//! if <ChangeGeom> is TRUE , the surface can be //! if <ChangeGeom> is TRUE , the surface can be
//! changed . //! changed .
//! if <UpdatePCurve> is TRUE, update the pcurves of the //! if <UpdatePCurve> is TRUE, update the pcurves of the
//! edges of <F> on the new surface.if the surface has been changed, //! edges of <F> on the new surface.if the surface has been changed,
//! Returns True if The Surface of <NF> has changed. //! Returns True if The Surface of <NF> has changed.
Standard_EXPORT static Standard_Boolean EnLargeFace (const TopoDS_Face& F, TopoDS_Face& NF, const Standard_Boolean ChangeGeom, const Standard_Boolean UpDatePCurve = Standard_False, const Standard_Boolean enlargeU = Standard_True, const Standard_Boolean enlargeVfirst = Standard_True, const Standard_Boolean enlargeVlast = Standard_True); Standard_EXPORT static Standard_Boolean EnLargeFace (const TopoDS_Face& F,
TopoDS_Face& NF,
const Standard_Boolean ChangeGeom,
const Standard_Boolean UpDatePCurve = Standard_False,
const Standard_Boolean enlargeU = Standard_True,
const Standard_Boolean enlargeVfirst = Standard_True,
const Standard_Boolean enlargeVlast = Standard_True);
Standard_EXPORT static void ExtentFace (const TopoDS_Face& F, TopTools_DataMapOfShapeShape& ConstShapes, TopTools_DataMapOfShapeShape& ToBuild, const TopAbs_State Side, const Standard_Real TolConf, TopoDS_Face& NF); Standard_EXPORT static void ExtentFace (const TopoDS_Face& F,
TopTools_DataMapOfShapeShape& ConstShapes,
TopTools_DataMapOfShapeShape& ToBuild,
const TopAbs_State Side,
const Standard_Real TolConf,
TopoDS_Face& NF);
//! Via the wire explorer store in <NOnV1> for //! Via the wire explorer store in <NOnV1> for
//! an Edge <E> of <W> his Edge neighbour on the first //! an Edge <E> of <W> his Edge neighbour on the first
//! vertex <V1> of <E>. //! vertex <V1> of <E>.
//! Store in NOnV2 the Neighbour of <E>on the last //! Store in NOnV2 the Neighbour of <E>on the last
//! vertex <V2> of <E>. //! vertex <V2> of <E>.
Standard_EXPORT static void BuildNeighbour (const TopoDS_Wire& W, const TopoDS_Face& F, TopTools_DataMapOfShapeShape& NOnV1, TopTools_DataMapOfShapeShape& NOnV2); Standard_EXPORT static void BuildNeighbour (const TopoDS_Wire& W,
const TopoDS_Face& F,
TopTools_DataMapOfShapeShape& NOnV1,
TopTools_DataMapOfShapeShape& NOnV2);
//! Store in MVE for a vertex <V> in <S> the incident //! Store in MVE for a vertex <V> in <S> the incident
//! edges <E> in <S>. //! edges <E> in <S>.
//! An Edge is Store only one Time for a vertex. //! An Edge is Store only one Time for a vertex.
Standard_EXPORT static void MapVertexEdges (const TopoDS_Shape& S, TopTools_DataMapOfShapeListOfShape& MVE); Standard_EXPORT static void MapVertexEdges (const TopoDS_Shape& S,
TopTools_DataMapOfShapeListOfShape& MVE);
//! Remove the non valid part of an offsetshape //! Remove the non valid part of an offsetshape
//! 1 - Remove all the free boundary and the faces //! 1 - Remove all the free boundary and the faces
@ -126,9 +173,14 @@ public:
//! 2 - Remove all the shapes not valid in the result //! 2 - Remove all the shapes not valid in the result
//! (according to the side of offseting) //! (according to the side of offseting)
//! in this verion only the first point is implemented. //! in this verion only the first point is implemented.
Standard_EXPORT static TopoDS_Shape Deboucle3D (const TopoDS_Shape& S, const TopTools_MapOfShape& Boundary); Standard_EXPORT static TopoDS_Shape Deboucle3D (const TopoDS_Shape& S,
const TopTools_MapOfShape& Boundary);
Standard_EXPORT static void CorrectOrientation (const TopoDS_Shape& SI, const TopTools_IndexedMapOfShape& NewEdges, Handle(BRepAlgo_AsDes)& AsDes, BRepAlgo_Image& InitOffset, const Standard_Real Offset); Standard_EXPORT static void CorrectOrientation (const TopoDS_Shape& SI,
const TopTools_IndexedMapOfShape& NewEdges,
Handle(BRepAlgo_AsDes)& AsDes,
BRepAlgo_Image& InitOffset,
const Standard_Real Offset);
Standard_EXPORT static Standard_Real Gabarit (const Handle(Geom_Curve)& aCurve); Standard_EXPORT static Standard_Real Gabarit (const Handle(Geom_Curve)& aCurve);

View File

@ -1,8 +1,4 @@
puts "TODO OCC25925 ALL: ERROR. offsetperform operation not done." puts "TODO OCC25925 ALL: Faulty shapes in variables faulty_1 to faulty_"
puts "TODO OCC25925 ALL: Error: The command cannot be built"
puts "TODO OCC25925 ALL: Faulty OCC5805 : result is not Closed shape"
puts "TODO OCC25925 ALL: Tcl Exception: result is not a topological shape"
puts "TODO OCC25925 ALL: TEST INCOMPLETE"
puts "============" puts "============"
puts "OCC5805" puts "OCC5805"

View File

@ -1,6 +1,6 @@
puts "TODO OCC24862 ALL: Error : is WRONG because number of" puts "TODO OCC23068 ALL: Error : is WRONG because number of"
puts "TODO OCC24862 ALL: Error : The area of result shape is" puts "TODO OCC23068 ALL: Error : The area of result shape is"
puts "TODO OCC24682 ALL: Faulty shapes in variables faulty_1 to faulty_" puts "TODO OCC23068 ALL: Faulty shapes in variables faulty_1 to faulty_"
puts "============" puts "============"
puts "OCC5805" puts "OCC5805"

View File

@ -1,7 +1,4 @@
puts "TODO OCC25925 ALL: ERROR. offsetperform operation not done." puts "TODO OCC25925 ALL: Faulty shapes in variables faulty_1 to faulty_"
puts "TODO OCC25925 ALL: Error: The command cannot be built"
puts "TODO OCC25925 ALL: Tcl Exception: result is not a topological shape"
puts "TODO OCC25925 ALL: TEST INCOMPLETE"
puts "============" puts "============"
puts "OCC5805" puts "OCC5805"
@ -41,8 +38,6 @@ catch { OFFSETSHAPE $distance {s_2} $calcul $type }
checkprops result -s 495.635 checkprops result -s 495.635
checkshape result checkshape result
checknbshapes result -vertex 2 -edge 3 -wire 3 -face 3 -shell 1 -solid 1 -compsolid 0 -compound 0 -shape 13
set index [lsearch [whatis s] Closed] set index [lsearch [whatis s] Closed]
if {$index == -1} { if {$index == -1} {
puts "Faulty ${BugNumber} : s is not Closed shape" puts "Faulty ${BugNumber} : s is not Closed shape"

View File

@ -1,6 +1,7 @@
puts "TODO OCC24862 ALL: Error : is WRONG because number of" puts "TODO OCC23068 ALL: ERROR. offsetperform operation not done."
puts "TODO OCC24862 ALL: Error : The area of result shape is" puts "TODO OCC25925 ALL: Error: The command cannot be built"
puts "TODO OCC24682 ALL: Faulty shapes in variables faulty_1 to faulty_" puts "TODO OCC25925 ALL: Tcl Exception: result is not a topological shape"
puts "TODO OCC25925 ALL: TEST INCOMPLETE"
puts "============" puts "============"
puts "OCC5805" puts "OCC5805"

View File

@ -8,20 +8,17 @@ puts ""
psphere a 100 psphere a 100
explode a f explode a f
thickshell r a_1 10 i 1.e-7 thickshell result a_1 10 i 1.e-7
explode r donly result
set bug_info [whatis r] checkshape result
if {$bug_info != "r is a shape SOLID FORWARD Free Modified Closed\n"} {
puts "ERROR: OCC26233 is reproduced. Shape r has incorrect characteristics." checknbshapes result -solid 1 -shell 2 -face 2 -wire 2 -edge 6 -vertex 4 -shape 17
set tolres [checkmaxtol result]
if { ${tolres} > 2.e-7} {
puts "Error: bad tolerance of result"
} }
set bug_info [whatis r_1] checkprops result -v 1.38649e+006
if {$bug_info != "r_1 is a shape SHELL FORWARD Modified Orientable Closed\n"} {
puts "ERROR: OCC26233 is reproduced. Shape r_1 has incorrect characteristics."
}
set bug_info [whatis r_2]
if {$bug_info != "r_2 is a shape SHELL REVERSED Modified Orientable Closed\n"} {
puts "ERROR: OCC26233 is reproduced. Shape r_2 has incorrect characteristics."
}

View File

@ -0,0 +1,20 @@
puts "============"
puts "OCC28903"
puts "============"
puts ""
##################################################################################
# BRepOffset_MakeOffset produces invalid shape (thickshell) in Intersection mode
##################################################################################
restore [locate_data_file bug28903_Fuse_3.brep] a
thickshell result a 10 i
donly result
checkshape result
checknbshapes result -solid 1 -shell 1 -face 6 -wire 7 -edge 12 -vertex 7 -shape 34
checkmaxtol result -min_tol 0.0015
checkprops result -v 1.1845e+006

View File

@ -1,8 +1,6 @@
puts "TODO OCC23748 ALL: ERROR. offsetperform operation not done." puts "TODO OCC23748 ALL: ERROR. offsetperform operation not done."
puts "TODO OCC23748 ALL: Error: The command cannot be built"
puts "TODO OCC23748 ALL: Error : The offset cannot be built." puts "TODO OCC23748 ALL: Error : The offset cannot be built."
psphere s 15 270 psphere s 15 270
OFFSETSHAPE 1 {s_2} $calcul $type OFFSETSHAPE 1 {s_2} $calcul $type
checkprops result -v 0

View File

@ -2,4 +2,4 @@ ptorus s 20 5 270
OFFSETSHAPE 1 {s_2} $calcul $type OFFSETSHAPE 1 {s_2} $calcul $type
checkprops result -v -296.088 checkprops result -v 3370.13

View File

@ -1,4 +1,3 @@
puts "TODO OCC23068 ALL: Faulty shapes in variables faulty_1 to faulty_2 "
ptorus s 20 5 270 ptorus s 20 5 270
OFFSETSHAPE 1 {s_2 s_3} $calcul $type OFFSETSHAPE 1 {s_2 s_3} $calcul $type

View File

@ -1,8 +1,8 @@
puts "TODO OCC24156 MacOS: \\*\\* Exception \\*\\*.*" puts "TODO OCC24156 MacOS: \\*\\* Exception \\*\\*.*"
puts "TODO OCC24156 MacOS: An exception was caught" puts "TODO OCC24156 MacOS: An exception was caught"
puts "TODO OCC24156 MacOS: TEST INCOMPLETE" puts "TODO OCC24156 MacOS: TEST INCOMPLETE"
puts "TODO OCC23068 ALL: Error : The volume of result shape is" puts "TODO OCC23068 ALL: Faulty shapes in variables faulty_1 to faulty_"
puts "TODO OCC25406 ALL: Error: bsection of the result and s is not equal to zero" puts "TODO OCC23068 ALL: Error : The area of face result_3 of the resulting shape is negative."
ellipse w1 0 0 0 15 10 ellipse w1 0 0 0 15 10
mkedge w1 w1 0 pi/2 mkedge w1 w1 0 pi/2
@ -13,5 +13,3 @@ mkplane w w
revol s w 0 0 0 0 0 1 90 revol s w 0 0 0 0 0 1 90
OFFSETSHAPE 1 {} $calcul $type OFFSETSHAPE 1 {} $calcul $type
checkprops result -v 0

View File

@ -1,7 +1,3 @@
puts "TODO OCC26578 All:An exception was caught"
puts "TODO OCC26578 All:\\*\\* Exception \\*\\*"
puts "TODO OCC26578 All:TEST INCOMPLETE"
restore [locate_data_file bug26663_test_offset_J9.brep] s restore [locate_data_file bug26663_test_offset_J9.brep] s
OFFSETSHAPE ${off_param} {} ${calcul} ${type} OFFSETSHAPE ${off_param} {} ${calcul} ${type}
checknbshapes result -ref [lrange [nbshapes s] 8 19]