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

0031233: Reading SAT files produces invalid shapes

Added removal of overlapping "tails" while splitting wires with a seam edge in ShapeFix_ComposeShell::SplitWire.
This commit is contained in:
anv 2020-01-21 14:47:38 +03:00 committed by bugmaster
parent e2550e48f1
commit 8b3fbdef34
3 changed files with 183 additions and 46 deletions

View File

@ -25,12 +25,14 @@
#include <BRepTools.hxx> #include <BRepTools.hxx>
#include <BRepTopAdaptor_FClass2d.hxx> #include <BRepTopAdaptor_FClass2d.hxx>
#include <Extrema_ExtPC2d.hxx> #include <Extrema_ExtPC2d.hxx>
#include <GCPnts_AbscissaPoint.hxx>
#include <Geom2d_Curve.hxx> #include <Geom2d_Curve.hxx>
#include <Geom2d_Line.hxx> #include <Geom2d_Line.hxx>
#include <Geom2dAdaptor_Curve.hxx> #include <Geom2dAdaptor_Curve.hxx>
#include <Geom2dInt_GInter.hxx> #include <Geom2dInt_GInter.hxx>
#include <Geom_Curve.hxx> #include <Geom_Curve.hxx>
#include <Geom_ElementarySurface.hxx> #include <Geom_ElementarySurface.hxx>
#include <GeomAdaptor_Curve.hxx>
#include <GeomAdaptor_Surface.hxx> #include <GeomAdaptor_Surface.hxx>
#include <gp_Dir2d.hxx> #include <gp_Dir2d.hxx>
#include <gp_Lin2d.hxx> #include <gp_Lin2d.hxx>
@ -281,6 +283,13 @@ Standard_Boolean ShapeFix_ComposeShell::Status (const ShapeExtend_Status status)
#define ITP_ENDSEG 32 // stop of tangential segment #define ITP_ENDSEG 32 // stop of tangential segment
#define ITP_TANG 64 // tangential point #define ITP_TANG 64 // tangential point
enum ShapeFix_MakeCut
{
ShapeFix_NoCut,
ShapeFix_CutHead,
ShapeFix_CutTail
};
//======================================================================= //=======================================================================
//function : PointLineDeviation //function : PointLineDeviation
//purpose : auxilary //purpose : auxilary
@ -809,6 +818,7 @@ static Standard_Real GetGridResolution(const Handle(TColStd_HArray1OfReal) Split
ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wire, ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wire,
TColStd_SequenceOfInteger& indexes, TColStd_SequenceOfInteger& indexes,
const TColStd_SequenceOfReal& values, const TColStd_SequenceOfReal& values,
const TopTools_MapOfShape& theVerts,
TopTools_SequenceOfShape& vertices, TopTools_SequenceOfShape& vertices,
const TColStd_SequenceOfInteger &SegmentCodes, const TColStd_SequenceOfInteger &SegmentCodes,
const Standard_Boolean isCutByU, const Standard_Boolean isCutByU,
@ -822,6 +832,7 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
ShapeAnalysis_Edge sae; ShapeAnalysis_Edge sae;
Standard_Integer start = 1; Standard_Integer start = 1;
TopAbs_Orientation anWireOrient = wire.Orientation(); TopAbs_Orientation anWireOrient = wire.Orientation();
Standard_Integer aLastSplitted = 0;
gp_Trsf T; gp_Trsf T;
if ( ! myLoc.IsIdentity() ) T = myLoc.Inverted().Transformation(); if ( ! myLoc.IsIdentity() ) T = myLoc.Inverted().Transformation();
@ -837,6 +848,7 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
std::cout << "Error: ShapeFix_ComposeShell::SplitWire: edge dismissed" << std::endl; std::cout << "Error: ShapeFix_ComposeShell::SplitWire: edge dismissed" << std::endl;
#endif #endif
i--; i--;
aLastSplitted = 0;
continue; continue;
} }
} }
@ -855,6 +867,7 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
result.AddEdge ( 0, edge, iumin, iumax, ivmin, ivmax ); result.AddEdge ( 0, edge, iumin, iumax, ivmin, ivmax );
if(code!=0 || wire.Orientation()!=TopAbs_EXTERNAL) // pdn 0 code handling for extertnal wires if(code!=0 || wire.Orientation()!=TopAbs_EXTERNAL) // pdn 0 code handling for extertnal wires
DefinePatch ( result, code, isCutByU, cutIndex ); DefinePatch ( result, code, isCutByU, cutIndex );
aLastSplitted = 0;
continue; continue;
} }
//find non-manifold vertices on edge //find non-manifold vertices on edge
@ -948,11 +961,14 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
Standard_Integer NbEdgesStart = result.NbEdges(); Standard_Integer NbEdgesStart = result.NbEdges();
Standard_Boolean splitted = Standard_False; Standard_Boolean splitted = Standard_False;
Standard_Real currPar=lastPar; //SK Standard_Real currPar=lastPar; //SK
Standard_Integer aCurSplitted = 0;
for ( Standard_Integer j = start; j <= stop; prevPar = currPar, j++ ) { for ( Standard_Integer j = start; j <= stop; prevPar = currPar, j++ ) {
if ( ! splitted && j >= stop ) {// no splitting at all if ( ! splitted && j >= stop ) {// no splitting at all
// code = SegmentCodes ( j >1 ? j-1 : SegmentCodes.Length() ); // classification code // code = SegmentCodes ( j >1 ? j-1 : SegmentCodes.Length() ); // classification code
break; break;
} }
ShapeFix_MakeCut aCutPart = ShapeFix_NoCut;
Standard_Boolean isCutChecked = Standard_False;
currPar = ( j < stop ? values.Value(j) : lastPar ); currPar = ( j < stop ? values.Value(j) : lastPar );
//fix for case when pcurve is periodic and first parameter of edge is more than 2P //fix for case when pcurve is periodic and first parameter of edge is more than 2P
//method ShapeBuild_Edge::CopyRanges shift pcurve to range 0-2P and parameters of cutting //method ShapeBuild_Edge::CopyRanges shift pcurve to range 0-2P and parameters of cutting
@ -990,51 +1006,77 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
// Doubled to prevent edge being fully covered by its vertices tolerance (invalid edge). // Doubled to prevent edge being fully covered by its vertices tolerance (invalid edge).
CheckByCurve3d ( lastVPnt, c3d, f3d+(currPar-firstPar)*(l3d-f3d)/span2d, CheckByCurve3d ( lastVPnt, c3d, f3d+(currPar-firstPar)*(l3d-f3d)/span2d,
T, lastVTol + 2 * Precision::Confusion() ) && T, lastVTol + 2 * Precision::Confusion() ) &&
lastPnt.Distance ( myGrid->Value ( C2d->Value(0.5*(currPar+lastPar)) ) ) <= tol ) { lastPnt.Distance ( myGrid->Value ( C2d->Value(0.5*(currPar+lastPar)) ) ) <= tol )
V = lastV; {
Standard_Real uRes = myUResolution; if (theVerts.Contains(lastV))
Standard_Real vRes = myVResolution; {
if(isCutByU) { GeomAdaptor_Curve AC(c3d, currPar, (edge.Orientation() == TopAbs_REVERSED) ? firstPar : lastPar);
Standard_Real gridRes = GetGridResolution(myGrid->UJointValues(),cutIndex)/tol; Standard_Real L = GCPnts_AbscissaPoint::Length(AC, Precision::Confusion());
uRes = Min(myUResolution,gridRes); if (L < BRep_Tool::Tolerance(lastV))
aCutPart = ShapeFix_CutTail;
isCutChecked = Standard_True;
} }
else { if (aCutPart == ShapeFix_NoCut)
Standard_Real gridRes = GetGridResolution(myGrid->VJointValues(),cutIndex)/tol; {
vRes = Min(myVResolution,gridRes); V = lastV;
Standard_Real uRes = myUResolution;
Standard_Real vRes = myVResolution;
if (isCutByU) {
Standard_Real gridRes = GetGridResolution(myGrid->UJointValues(), cutIndex) / tol;
uRes = Min(myUResolution, gridRes);
}
else {
Standard_Real gridRes = GetGridResolution(myGrid->VJointValues(), cutIndex) / tol;
vRes = Min(myVResolution, gridRes);
}
if (IsCoincided(lastPnt2d, currPnt2d, uRes, vRes, tol) &&
IsCoincided(lastPnt2d, C2d->Value(0.5*(currPar + lastPar)), uRes, vRes, tol))
{
doCut = Standard_False;
}
} }
if ( IsCoincided ( lastPnt2d, currPnt2d, uRes, vRes, tol ) &&
IsCoincided ( lastPnt2d, C2d->Value(0.5*(currPar+lastPar)), uRes, vRes, tol ) )
doCut = Standard_False;
} }
else if ( currPnt.Distance ( prevVPnt ) <= prevVTol && else if ( currPnt.Distance ( prevVPnt ) <= prevVTol &&
prevPnt.Distance ( currPnt ) <= tol && prevPnt.Distance ( currPnt ) <= tol &&
// Tolerance is increased to prevent degenerated cuts in cases where all vertex // Tolerance is increased to prevent degenerated cuts in cases where all vertex
// tolerance is covered by distance of the edge curve from vertex point. // tolerance is covered by distance of the edge curve from vertex point.
// Doubled to prevent edge being fully covered by its vertices tolerance (invalid edge). // Doubled to prevent edge being fully covered by its vertices tolerance (invalid edge).
CheckByCurve3d ( prevVPnt, c3d, f3d+(currPar-firstPar)*(l3d-f3d)/span2d, CheckByCurve3d ( prevVPnt, c3d, f3d + (currPar - firstPar)*(l3d - f3d) / span2d,
T, prevVTol + 2 * Precision::Confusion()) && T, prevVTol + 2 * Precision::Confusion()) &&
prevPnt.Distance ( myGrid->Value ( C2d->Value(0.5*(currPar+prevPar)) ) ) <= tol ) { prevPnt.Distance ( myGrid->Value ( C2d->Value(0.5*(currPar + prevPar)) ) ) <= tol )
V = prevV; {
Standard_Real uRes = myUResolution; if (theVerts.Contains(prevV))
Standard_Real vRes = myVResolution; {
if(isCutByU) { GeomAdaptor_Curve AC(c3d, (edge.Orientation() == TopAbs_REVERSED) ? lastPar : firstPar, currPar);
Standard_Real gridRes = GetGridResolution(myGrid->UJointValues(),cutIndex)/tol; Standard_Real L = GCPnts_AbscissaPoint::Length(AC, Precision::Confusion());
uRes = Min(myUResolution,gridRes); if (L < BRep_Tool::Tolerance(prevV))
aCutPart = ShapeFix_CutHead;
isCutChecked = Standard_True;
} }
else { if (aCutPart == ShapeFix_NoCut)
Standard_Real gridRes = GetGridResolution(myGrid->VJointValues(),cutIndex)/tol; {
vRes = Min(myVResolution,gridRes); V = prevV;
} Standard_Real uRes = myUResolution;
if ( IsCoincided ( prevPnt2d, currPnt2d, uRes, vRes, tol ) && Standard_Real vRes = myVResolution;
IsCoincided ( prevPnt2d, C2d->Value(0.5*(currPar+prevPar)), uRes, vRes, tol ) ) { if (isCutByU) {
vertices.Append ( prevV ); Standard_Real gridRes = GetGridResolution(myGrid->UJointValues(), cutIndex) / tol;
code = SegmentCodes ( j ); // classification code - update for next segment uRes = Min(myUResolution, gridRes);
continue; // no splitting at this point, go to next one }
else {
Standard_Real gridRes = GetGridResolution(myGrid->VJointValues(), cutIndex) / tol;
vRes = Min(myVResolution, gridRes);
}
if (IsCoincided(prevPnt2d, currPnt2d, uRes, vRes, tol) &&
IsCoincided(prevPnt2d, C2d->Value(0.5*(currPar + prevPar)), uRes, vRes, tol)) {
vertices.Append(prevV);
code = SegmentCodes(j); // classification code - update for next segment
continue; // no splitting at this point, go to next one
}
} }
} }
//:abv 28.05.02: OCC320 Sample_2: if maxtol = 1e-7, the vertex tolerance //:abv 28.05.02: OCC320 Sample_2: if maxtol = 1e-7, the vertex tolerance
// is actually ignored - protect against new vertex on degenerated edge // is actually ignored - protect against new vertex on degenerated edge
else if ( BRep_Tool::Degenerated(edge) && prevV.IsSame(lastV) ) { else if (BRep_Tool::Degenerated(edge) && prevV.IsSame(lastV)) {
V = prevV; V = prevV;
} }
} }
@ -1045,10 +1087,40 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
if ( V.IsNull() ) { if ( V.IsNull() ) {
B.MakeVertex ( V, currPnt.Transformed(myLoc.Transformation()), tolEdge ); B.MakeVertex ( V, currPnt.Transformed(myLoc.Transformation()), tolEdge );
vertices.Append ( V ); vertices.Append ( V );
if (!isCutChecked)
{
Standard_Real aVTol;
if (theVerts.Contains(prevV))
{
aVTol = BRep_Tool::Tolerance(prevV);
if (currPnt.SquareDistance(prevVPnt) < aVTol * aVTol)
{
GeomAdaptor_Curve AC(c3d, (edge.Orientation() == TopAbs_REVERSED) ? lastPar : firstPar, currPar);
Standard_Real L = GCPnts_AbscissaPoint::Length(AC, Precision::Confusion());
if (L < aVTol)
aCutPart = ShapeFix_CutHead;
}
}
else if (theVerts.Contains(lastV))
{
aVTol = BRep_Tool::Tolerance(lastV);
if (currPnt.SquareDistance(lastVPnt) < aVTol * aVTol)
{
GeomAdaptor_Curve AC(c3d, currPar, (edge.Orientation() == TopAbs_REVERSED) ? firstPar : lastPar);
Standard_Real L = GCPnts_AbscissaPoint::Length(AC, Precision::Confusion());
if (L < aVTol)
aCutPart = ShapeFix_CutTail;
}
}
}
} }
// else adjusted to end, fill all resting vertices // else adjusted to end, fill all resting vertices
else if ( ! doCut ) { else if ( ! doCut ) {
for ( ; j < stop; j++ ) vertices.Append ( lastV ); for (; j < stop; j++)
{
vertices.Append(lastV);
aCurSplitted++;
}
if ( ! splitted ) break; // no splitting at all if ( ! splitted ) break; // no splitting at all
currPar = lastPar; currPar = lastPar;
} }
@ -1138,8 +1210,21 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
code = ( ( isCutByU == (j == 1) ) ? 1 : 2 ); code = ( ( isCutByU == (j == 1) ) ? 1 : 2 );
} }
result.AddEdge ( 0, newEdge, iumin, iumax, ivmin, ivmax ); if (aCutPart == ShapeFix_CutHead)
DefinePatch ( result, code, isCutByU, cutIndex ); {
V.Orientation(TopAbs_FORWARD);
Context()->Replace(prevV, V);
// also replce this vertice in the sequence
for (Standard_Integer iV = 1; aLastSplitted > 0; aLastSplitted--, iV++)
{
vertices.ChangeValue(vertices.Length() - iV) = V;
}
}
else
{
result.AddEdge(0, newEdge, iumin, iumax, ivmin, ivmax);
DefinePatch(result, code, isCutByU, cutIndex);
}
// Changing prev parameters // Changing prev parameters
prevV = V; prevV = V;
@ -1147,8 +1232,16 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
prevVPnt = BRep_Tool::Pnt ( V ); prevVPnt = BRep_Tool::Pnt ( V );
prevPnt = currPnt; prevPnt = currPnt;
prevPnt2d = currPnt2d; prevPnt2d = currPnt2d;
if (aCutPart == ShapeFix_CutTail)
{
Context()->Replace(lastV, V);
for (; j + 1 < stop; j++) vertices.Append(lastV);
break;
}
} }
start = stop; start = stop;
aLastSplitted = aCurSplitted;
if ( splitted ) { if ( splitted ) {
// record replacement in context // record replacement in context
@ -1209,6 +1302,8 @@ Standard_Boolean ShapeFix_ComposeShell::SplitByLine (ShapeFix_WireSegment &wire,
TColStd_SequenceOfReal IntEdgePar; // parameter of intersection point on edge TColStd_SequenceOfReal IntEdgePar; // parameter of intersection point on edge
TColStd_SequenceOfReal IntLinePar; // parameter of intersection point on line TColStd_SequenceOfReal IntLinePar; // parameter of intersection point on line
TopTools_MapOfShape exactVertices; // vertices that line passes through
Standard_Boolean isnonmanifold = (wire.Orientation() == TopAbs_INTERNAL); Standard_Boolean isnonmanifold = (wire.Orientation() == TopAbs_INTERNAL);
//gka correction for non-manifold vertices SAMTECH //gka correction for non-manifold vertices SAMTECH
if(wire.IsVertex()) { if(wire.IsVertex()) {
@ -1325,6 +1420,7 @@ Standard_Boolean ShapeFix_ComposeShell::SplitByLine (ShapeFix_WireSegment &wire,
IntLinePar.Append ( ParamPointsOnLine ( pos, prevPos, line ) ); // !! - maybe compute exactly ? IntLinePar.Append ( ParamPointsOnLine ( pos, prevPos, line ) ); // !! - maybe compute exactly ?
IntEdgePar.Append ( isreversed ? l : f ); IntEdgePar.Append ( isreversed ? l : f );
IntEdgeInd.Append ( iedge ); IntEdgeInd.Append ( iedge );
exactVertices.Add (sae.FirstVertex(E));
} }
} }
@ -1403,6 +1499,7 @@ Standard_Boolean ShapeFix_ComposeShell::SplitByLine (ShapeFix_WireSegment &wire,
IntLinePar.Append ( ParamPointsOnLine ( pos, firstPos, line ) ); IntLinePar.Append ( ParamPointsOnLine ( pos, firstPos, line ) );
IntEdgePar.Append ( isreversed ? f : l ); IntEdgePar.Append ( isreversed ? f : l );
IntEdgeInd.Append ( iedge ); IntEdgeInd.Append ( iedge );
exactVertices.Add (sae.LastVertex(E));
} }
} }
} }
@ -1540,7 +1637,7 @@ Standard_Boolean ShapeFix_ComposeShell::SplitByLine (ShapeFix_WireSegment &wire,
//======================================= //=======================================
// Split edges in the wire by intersection points and fill vertices array // Split edges in the wire by intersection points and fill vertices array
TopTools_SequenceOfShape IntVertices; TopTools_SequenceOfShape IntVertices;
wire = SplitWire ( wire, IntEdgeInd, IntEdgePar, IntVertices, wire = SplitWire ( wire, IntEdgeInd, IntEdgePar, exactVertices, IntVertices,
aNewSegCodes, isCutByU, cutIndex ); aNewSegCodes, isCutByU, cutIndex );
// add all data to input arrays // add all data to input arrays
@ -1586,7 +1683,8 @@ void ShapeFix_ComposeShell::SplitByLine (ShapeFix_SequenceOfWireSegment &wires,
// merge null-length tangential segments into one-point tangencies or intersections // merge null-length tangential segments into one-point tangencies or intersections
for ( i = 1; i < SplitLinePar.Length(); i++ ) { for ( i = 1; i < SplitLinePar.Length(); i++ ) {
if ( Abs ( SplitLinePar(i+1) - SplitLinePar(i) ) > ::Precision::PConfusion() && !SplitLineVertex(i).IsSame(SplitLineVertex(i+1)) ) continue; Standard_Boolean isSameVertex = SplitLineVertex(i).IsSame(SplitLineVertex(i + 1));
if ( Abs ( SplitLinePar(i+1) - SplitLinePar(i) ) > ::Precision::PConfusion() && !isSameVertex) continue;
if ( ( SplitLineCode(i) & ITP_ENDSEG && if ( ( SplitLineCode(i) & ITP_ENDSEG &&
SplitLineCode(i+1) & ITP_BEGSEG ) || SplitLineCode(i+1) & ITP_BEGSEG ) ||
( SplitLineCode(i) & ITP_BEGSEG && ( SplitLineCode(i) & ITP_BEGSEG &&
@ -1597,6 +1695,17 @@ void ShapeFix_ComposeShell::SplitByLine (ShapeFix_SequenceOfWireSegment &wires,
SplitLineCode.Remove(i+1); SplitLineCode.Remove(i+1);
SplitLineVertex.Remove(i+1); SplitLineVertex.Remove(i+1);
} }
else if (isSameVertex &&
((SplitLineCode (i) & ITP_TANG &&
SplitLineCode (i + 1) & ITP_INTER) ||
(SplitLineCode (i) & ITP_INTER &&
SplitLineCode (i + 1) & ITP_TANG)))
{
SplitLineCode.SetValue(i, IOR_BOTH | ITP_INTER );
SplitLinePar.Remove(i + 1);
SplitLineCode.Remove(i + 1);
SplitLineVertex.Remove(i + 1);
}
} }
// go along line, split it by intersection points and create edges // go along line, split it by intersection points and create edges

View File

@ -33,6 +33,7 @@
#include <TColStd_SequenceOfInteger.hxx> #include <TColStd_SequenceOfInteger.hxx>
#include <TColStd_SequenceOfReal.hxx> #include <TColStd_SequenceOfReal.hxx>
#include <TopTools_SequenceOfShape.hxx> #include <TopTools_SequenceOfShape.hxx>
#include <TopTools_MapOfShape.hxx>
class ShapeExtend_CompositeSurface; class ShapeExtend_CompositeSurface;
class ShapeAnalysis_TransferParameters; class ShapeAnalysis_TransferParameters;
class ShapeExtend_WireData; class ShapeExtend_WireData;
@ -98,7 +99,10 @@ public:
//! Here face defines both set of wires and way of getting //! Here face defines both set of wires and way of getting
//! pcurves. Precision is used (together with tolerance of edges) //! pcurves. Precision is used (together with tolerance of edges)
//! for handling subtle cases, such as tangential intersections. //! for handling subtle cases, such as tangential intersections.
Standard_EXPORT void Init (const Handle(ShapeExtend_CompositeSurface)& Grid, const TopLoc_Location& L, const TopoDS_Face& Face, const Standard_Real Prec); Standard_EXPORT void Init (const Handle(ShapeExtend_CompositeSurface)& Grid,
const TopLoc_Location& L,
const TopoDS_Face& Face,
const Standard_Real Prec);
//! Returns (modifiable) flag for special 'closed' //! Returns (modifiable) flag for special 'closed'
//! mode which forces ComposeShell to consider //! mode which forces ComposeShell to consider
@ -139,7 +143,8 @@ public:
//! and all pcurves on the initial (pseudo)face are reassigned to //! and all pcurves on the initial (pseudo)face are reassigned to
//! that surface. If several wires are one inside another, single //! that surface. If several wires are one inside another, single
//! face is created. //! face is created.
Standard_EXPORT void DispatchWires (TopTools_SequenceOfShape& faces, ShapeFix_SequenceOfWireSegment& wires) const; Standard_EXPORT void DispatchWires (TopTools_SequenceOfShape& faces,
ShapeFix_SequenceOfWireSegment& wires) const;
//! Sets tool for transfer parameters from 3d to 2d and vice versa. //! Sets tool for transfer parameters from 3d to 2d and vice versa.
Standard_EXPORT void SetTransferParamTool (const Handle(ShapeAnalysis_TransferParameters)& TransferParam); Standard_EXPORT void SetTransferParamTool (const Handle(ShapeAnalysis_TransferParameters)& TransferParam);
@ -165,7 +170,13 @@ protected:
//! between two intersections: tells if segment is on left or right side //! between two intersections: tells if segment is on left or right side
//! of cutting line, or tangent to it (by several points recomputed to 3d, //! of cutting line, or tangent to it (by several points recomputed to 3d,
//! distance is compared with tolerance of corresponding edge). //! distance is compared with tolerance of corresponding edge).
Standard_EXPORT Standard_Integer ComputeCode (const Handle(ShapeExtend_WireData)& wire, const gp_Lin2d& line, const Standard_Integer begInd, const Standard_Integer endInd, const Standard_Real begPar, const Standard_Real endPar, const Standard_Boolean IsInternal = Standard_False); Standard_EXPORT Standard_Integer ComputeCode (const Handle(ShapeExtend_WireData)& wire,
const gp_Lin2d& line,
const Standard_Integer begInd,
const Standard_Integer endInd,
const Standard_Real begPar,
const Standard_Real endPar,
const Standard_Boolean IsInternal = Standard_False);
//! Splits edges in the wire by given indices of edges and //! Splits edges in the wire by given indices of edges and
//! parameters on them. Returns resulting wire and vertices //! parameters on them. Returns resulting wire and vertices
@ -178,7 +189,14 @@ protected:
//! NOTE: If edge is splitted, it is replaced by wire, and //! NOTE: If edge is splitted, it is replaced by wire, and
//! order of edges in the wire corresponds to FORWARD orientation //! order of edges in the wire corresponds to FORWARD orientation
//! of the edge. //! of the edge.
Standard_EXPORT ShapeFix_WireSegment SplitWire (ShapeFix_WireSegment& wire, TColStd_SequenceOfInteger& indexes, const TColStd_SequenceOfReal& values, TopTools_SequenceOfShape& vertices, const TColStd_SequenceOfInteger& segcodes, const Standard_Boolean cutbyu, const Standard_Integer cutindex); Standard_EXPORT ShapeFix_WireSegment SplitWire (ShapeFix_WireSegment& wire,
TColStd_SequenceOfInteger& indexes,
const TColStd_SequenceOfReal& values,
const TopTools_MapOfShape& theVerts,
TopTools_SequenceOfShape& vertices,
const TColStd_SequenceOfInteger& segcodes,
const Standard_Boolean cutbyu,
const Standard_Integer cutindex);
//! Split edges in the wire by cutting line. //! Split edges in the wire by cutting line.
//! Wires with FORWARD or REVERSED orientation are considered //! Wires with FORWARD or REVERSED orientation are considered
@ -191,7 +209,13 @@ protected:
//! Method fills sequences of parameters of intersection points //! Method fills sequences of parameters of intersection points
//! of cutting line with all edges, their types, and corresponding //! of cutting line with all edges, their types, and corresponding
//! vertices (including ones created during splitting edges). //! vertices (including ones created during splitting edges).
Standard_EXPORT Standard_Boolean SplitByLine (ShapeFix_WireSegment& wire, const gp_Lin2d& line, const Standard_Boolean cutbyu, const Standard_Integer cutindex, TColStd_SequenceOfReal& SplitLinePar, TColStd_SequenceOfInteger& SplitLineCode, TopTools_SequenceOfShape& SplitLineVertex); Standard_EXPORT Standard_Boolean SplitByLine (ShapeFix_WireSegment& wire,
const gp_Lin2d& line,
const Standard_Boolean cutbyu,
const Standard_Integer cutindex,
TColStd_SequenceOfReal& SplitLinePar,
TColStd_SequenceOfInteger& SplitLineCode,
TopTools_SequenceOfShape& SplitLineVertex);
//! Split edges in the sequence of wires by cutting line. //! Split edges in the sequence of wires by cutting line.
//! Wires with FORWARD or REVERSED orientation are considered //! Wires with FORWARD or REVERSED orientation are considered
@ -206,7 +230,10 @@ protected:
//! All modifications (splitting) are recorded in context, //! All modifications (splitting) are recorded in context,
//! except splitting of wires marked as EXTERNAL //! except splitting of wires marked as EXTERNAL
//! (they are supposed to be former cutting lines). //! (they are supposed to be former cutting lines).
Standard_EXPORT void SplitByLine (ShapeFix_SequenceOfWireSegment& seqw, const gp_Lin2d& line, const Standard_Boolean cutbyu, const Standard_Integer cutindex); Standard_EXPORT void SplitByLine (ShapeFix_SequenceOfWireSegment& seqw,
const gp_Lin2d& line,
const Standard_Boolean cutbyu,
const Standard_Integer cutindex);
//! Split initial set of (closed) wires by grid of lines corresponding //! Split initial set of (closed) wires by grid of lines corresponding
//! to joints between patches on the composite surface. //! to joints between patches on the composite surface.
@ -231,12 +258,15 @@ protected:
//! taking EXTERNAL as necessary in fork points. Forks are detected //! taking EXTERNAL as necessary in fork points. Forks are detected
//! by common vertices. In fork point, most left way is seleccted //! by common vertices. In fork point, most left way is seleccted
//! among all possible ways. //! among all possible ways.
Standard_EXPORT void CollectWires (ShapeFix_SequenceOfWireSegment& wires, ShapeFix_SequenceOfWireSegment& seqw); Standard_EXPORT void CollectWires (ShapeFix_SequenceOfWireSegment& wires,
ShapeFix_SequenceOfWireSegment& seqw);
//! Creates new faces on one path of grid. It dispatches given loops //! Creates new faces on one path of grid. It dispatches given loops
//! (wires) into one or several faces depending on their mutual //! (wires) into one or several faces depending on their mutual
//! position. //! position.
Standard_EXPORT void MakeFacesOnPatch (TopTools_SequenceOfShape& faces, const Handle(Geom_Surface)& surf, TopTools_SequenceOfShape& loops) const; Standard_EXPORT void MakeFacesOnPatch (TopTools_SequenceOfShape& faces,
const Handle(Geom_Surface)& surf,
TopTools_SequenceOfShape& loops) const;
TopAbs_Orientation myOrient; TopAbs_Orientation myOrient;
TopoDS_Shape myResult; TopoDS_Shape myResult;

View File

@ -1,5 +1,3 @@
puts "TODO All: Faulty shapes in variables faulty_1 to faulty_"
puts "========" puts "========"
puts "OCC22919" puts "OCC22919"
puts "========" puts "========"