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

0030052: Data Exchange - STEP import missing surfaces

Changes made for #31233 are reverted.

Correction in the ShapeFix_ComposeShell:
 Modification of the method SplitByLine in order to find all points of the intersection for case when initial curve shifts in the positive or negative direction
 to avoid splitting edge when split point lies in the limits of the tolerance of the vertex.
 to correct shift of the pcurves for case when initial shift is more than 1 in the method ShapeFix_ComposeShell::SplitByLine

In the method ShapeFix_Face::FixMissingSeam added removing small edges having length less than working precision and removing wires having area less than working precision

Modified test cases; added test for related issue: bugs step bug31301
This commit is contained in:
gka 2020-02-26 16:35:58 +03:00 committed by bugmaster
parent 22fa1da36e
commit 6a9f983a16
15 changed files with 190 additions and 225 deletions

View File

@ -25,14 +25,12 @@
#include <BRepTools.hxx>
#include <BRepTopAdaptor_FClass2d.hxx>
#include <Extrema_ExtPC2d.hxx>
#include <GCPnts_AbscissaPoint.hxx>
#include <Geom2d_Curve.hxx>
#include <Geom2d_Line.hxx>
#include <Geom2dAdaptor_Curve.hxx>
#include <Geom2dInt_GInter.hxx>
#include <Geom_Curve.hxx>
#include <Geom_ElementarySurface.hxx>
#include <GeomAdaptor_Curve.hxx>
#include <GeomAdaptor_Surface.hxx>
#include <gp_Dir2d.hxx>
#include <gp_Lin2d.hxx>
@ -283,13 +281,6 @@ Standard_Boolean ShapeFix_ComposeShell::Status (const ShapeExtend_Status status)
#define ITP_ENDSEG 32 // stop of tangential segment
#define ITP_TANG 64 // tangential point
enum ShapeFix_MakeCut
{
ShapeFix_NoCut,
ShapeFix_CutHead,
ShapeFix_CutTail
};
//=======================================================================
//function : PointLineDeviation
//purpose : auxilary
@ -818,7 +809,6 @@ static Standard_Real GetGridResolution(const Handle(TColStd_HArray1OfReal) Split
ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wire,
TColStd_SequenceOfInteger& indexes,
const TColStd_SequenceOfReal& values,
const TopTools_MapOfShape& theVerts,
TopTools_SequenceOfShape& vertices,
const TColStd_SequenceOfInteger &SegmentCodes,
const Standard_Boolean isCutByU,
@ -832,7 +822,6 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
ShapeAnalysis_Edge sae;
Standard_Integer start = 1;
TopAbs_Orientation anWireOrient = wire.Orientation();
Standard_Integer aLastSplitted = 0;
gp_Trsf T;
if ( ! myLoc.IsIdentity() ) T = myLoc.Inverted().Transformation();
@ -848,7 +837,6 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
std::cout << "Error: ShapeFix_ComposeShell::SplitWire: edge dismissed" << std::endl;
#endif
i--;
aLastSplitted = 0;
continue;
}
}
@ -867,7 +855,6 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
result.AddEdge ( 0, edge, iumin, iumax, ivmin, ivmax );
if(code!=0 || wire.Orientation()!=TopAbs_EXTERNAL) // pdn 0 code handling for extertnal wires
DefinePatch ( result, code, isCutByU, cutIndex );
aLastSplitted = 0;
continue;
}
//find non-manifold vertices on edge
@ -881,7 +868,6 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
// Collect data on edge
Standard_Real tolEdge = BRep_Tool::Tolerance(edge);
Standard_Real tol = LimitTolerance( tolEdge );
TopoDS_Vertex prevV = sae.FirstVertex(edge);
TopoDS_Vertex lastV = sae.LastVertex(edge);
Standard_Real prevVTol = LimitTolerance( BRep_Tool::Tolerance(prevV) );
@ -961,14 +947,11 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
Standard_Integer NbEdgesStart = result.NbEdges();
Standard_Boolean splitted = Standard_False;
Standard_Real currPar=lastPar; //SK
Standard_Integer aCurSplitted = 0;
for ( Standard_Integer j = start; j <= stop; prevPar = currPar, j++ ) {
if ( ! splitted && j >= stop ) {// no splitting at all
// code = SegmentCodes ( j >1 ? j-1 : SegmentCodes.Length() ); // classification code
break;
}
ShapeFix_MakeCut aCutPart = ShapeFix_NoCut;
Standard_Boolean isCutChecked = Standard_False;
currPar = ( j < stop ? values.Value(j) : lastPar );
//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
@ -1000,83 +983,56 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
currPnt2d = C2d->Value(currPar);
currPnt = myGrid->Value ( currPnt2d );
if ( currPnt.Distance ( lastVPnt ) <= lastVTol &&
lastPnt.Distance ( currPnt ) <= tol &&
// Tolerance is increased to prevent degenerated cuts in cases where all vertex
// 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).
CheckByCurve3d ( lastVPnt, c3d, f3d+(currPar-firstPar)*(l3d-f3d)/span2d,
T, lastVTol + 2 * Precision::Confusion() ) &&
lastPnt.Distance ( myGrid->Value ( C2d->Value(0.5*(currPar+lastPar)) ) ) <= tol )
{
if (theVerts.Contains(lastV))
{
GeomAdaptor_Curve AC(c3d, currPar, (edge.Orientation() == TopAbs_REVERSED) ? firstPar : lastPar);
Standard_Real L = GCPnts_AbscissaPoint::Length(AC, Precision::Confusion());
if (L < BRep_Tool::Tolerance(lastV))
aCutPart = ShapeFix_CutTail;
isCutChecked = Standard_True;
lastPnt.Distance (myGrid->Value (C2d->Value(0.5*(currPar+lastPar)))) <= lastVTol) {
V = lastV;
Standard_Real uRes = myUResolution;
Standard_Real vRes = myVResolution;
if(isCutByU) {
Standard_Real gridRes = GetGridResolution(myGrid->UJointValues(), cutIndex) / lastVTol;
uRes = Min(myUResolution,gridRes);
}
if (aCutPart == ShapeFix_NoCut)
{
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;
}
else {
Standard_Real gridRes = GetGridResolution(myGrid->VJointValues(), cutIndex) / lastVTol;
vRes = Min(myVResolution,gridRes);
}
if (IsCoincided(lastPnt2d, currPnt2d, uRes, vRes, lastVTol) &&
IsCoincided(lastPnt2d, C2d->Value(0.5*(currPar + lastPar)), uRes, vRes, lastVTol))
doCut = Standard_False;
}
else if ( currPnt.Distance ( prevVPnt ) <= prevVTol &&
prevPnt.Distance ( currPnt ) <= tol &&
// Tolerance is increased to prevent degenerated cuts in cases where all vertex
// 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).
CheckByCurve3d ( prevVPnt, c3d, f3d + (currPar - firstPar)*(l3d - f3d) / span2d,
CheckByCurve3d ( prevVPnt, c3d, f3d+(currPar-firstPar)*(l3d-f3d)/span2d,
T, prevVTol + 2 * Precision::Confusion()) &&
prevPnt.Distance ( myGrid->Value ( C2d->Value(0.5*(currPar + prevPar)) ) ) <= tol )
{
if (theVerts.Contains(prevV))
{
GeomAdaptor_Curve AC(c3d, (edge.Orientation() == TopAbs_REVERSED) ? lastPar : firstPar, currPar);
Standard_Real L = GCPnts_AbscissaPoint::Length(AC, Precision::Confusion());
if (L < BRep_Tool::Tolerance(prevV))
aCutPart = ShapeFix_CutHead;
isCutChecked = Standard_True;
prevPnt.Distance (myGrid->Value (C2d->Value(0.5*(currPar+prevPar)))) <= prevVTol) {
V = prevV;
Standard_Real uRes = myUResolution;
Standard_Real vRes = myVResolution;
if(isCutByU) {
Standard_Real gridRes = GetGridResolution(myGrid->UJointValues(), cutIndex) / prevVTol;
uRes = Min(myUResolution,gridRes);
}
if (aCutPart == ShapeFix_NoCut)
{
V = prevV;
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(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
}
else {
Standard_Real gridRes = GetGridResolution(myGrid->VJointValues(), cutIndex) / prevVTol;
vRes = Min(myVResolution,gridRes);
}
if (IsCoincided(prevPnt2d, currPnt2d, uRes, vRes, prevVTol) &&
IsCoincided(prevPnt2d, C2d->Value(0.5*(currPar + prevPar)), uRes, vRes, prevVTol)) {
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
// 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;
}
}
@ -1087,40 +1043,10 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
if ( V.IsNull() ) {
B.MakeVertex ( V, currPnt.Transformed(myLoc.Transformation()), tolEdge );
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 if ( ! doCut ) {
for (; j < stop; j++)
{
vertices.Append(lastV);
aCurSplitted++;
}
for ( ; j < stop; j++ ) vertices.Append ( lastV );
if ( ! splitted ) break; // no splitting at all
currPar = lastPar;
}
@ -1210,21 +1136,8 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
code = ( ( isCutByU == (j == 1) ) ? 1 : 2 );
}
if (aCutPart == ShapeFix_CutHead)
{
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);
}
result.AddEdge ( 0, newEdge, iumin, iumax, ivmin, ivmax );
DefinePatch ( result, code, isCutByU, cutIndex );
// Changing prev parameters
prevV = V;
@ -1232,16 +1145,8 @@ ShapeFix_WireSegment ShapeFix_ComposeShell::SplitWire (ShapeFix_WireSegment &wir
prevVPnt = BRep_Tool::Pnt ( V );
prevPnt = currPnt;
prevPnt2d = currPnt2d;
if (aCutPart == ShapeFix_CutTail)
{
Context()->Replace(lastV, V);
for (; j + 1 < stop; j++) vertices.Append(lastV);
break;
}
}
start = stop;
aLastSplitted = aCurSplitted;
if ( splitted ) {
// record replacement in context
@ -1302,8 +1207,6 @@ Standard_Boolean ShapeFix_ComposeShell::SplitByLine (ShapeFix_WireSegment &wire,
TColStd_SequenceOfReal IntEdgePar; // parameter of intersection point on edge
TColStd_SequenceOfReal IntLinePar; // parameter of intersection point on line
TopTools_MapOfShape exactVertices; // vertices that line passes through
Standard_Boolean isnonmanifold = (wire.Orientation() == TopAbs_INTERNAL);
//gka correction for non-manifold vertices SAMTECH
if(wire.IsVertex()) {
@ -1354,6 +1257,7 @@ Standard_Boolean ShapeFix_ComposeShell::SplitByLine (ShapeFix_WireSegment &wire,
Standard_Real f, l;
Handle(Geom2d_Curve) c2d;
if ( ! sae.PCurve ( E, myFace, c2d, f, l, Standard_False ) ) continue;
Handle(Geom2d_Curve) c2d_sav = c2d;
// get end points
gp_Pnt2d posf = c2d->Value(f), posl = c2d->Value(l);
@ -1384,8 +1288,9 @@ Standard_Boolean ShapeFix_ComposeShell::SplitByLine (ShapeFix_WireSegment &wire,
pppf.SetX ( pppf.X() + shift );
pppl.SetX ( pppl.X() + shift );
}
shiftNext.SetX ( -myUPeriod );
nbIter = (Standard_Integer)( 1 + Abs ( umax + shift - x ) / myUPeriod );
Standard_Real dUmax = umax + shift - x;
shiftNext.SetX (dUmax > 0 ? -myUPeriod : myUPeriod);
nbIter = (Standard_Integer)(1 + Abs (dUmax) / myUPeriod);
shift = ShapeAnalysis::AdjustByPeriod ( posf.X(), x, myUPeriod );
posf.SetX ( posf.X() + shift );
shift = ShapeAnalysis::AdjustByPeriod ( posl.X(), x, myUPeriod );
@ -1401,8 +1306,9 @@ Standard_Boolean ShapeFix_ComposeShell::SplitByLine (ShapeFix_WireSegment &wire,
pppf.SetY ( pppf.Y() + shift );
pppl.SetY ( pppl.Y() + shift );
}
shiftNext.SetY ( -myVPeriod );
nbIter = (Standard_Integer)( 1 + Abs ( umax + shift - y ) / myVPeriod );
Standard_Real dVmax = vmax + shift - y;
shiftNext.SetY (dVmax > 0 ? -myVPeriod : myVPeriod);
nbIter = (Standard_Integer)(1 + Abs (dVmax) / myVPeriod);
shift = ShapeAnalysis::AdjustByPeriod ( posf.Y(), y, myVPeriod );
posf.SetY ( posf.Y() + shift );
shift = ShapeAnalysis::AdjustByPeriod ( posl.Y(), y, myVPeriod );
@ -1420,7 +1326,6 @@ Standard_Boolean ShapeFix_ComposeShell::SplitByLine (ShapeFix_WireSegment &wire,
IntLinePar.Append ( ParamPointsOnLine ( pos, prevPos, line ) ); // !! - maybe compute exactly ?
IntEdgePar.Append ( isreversed ? l : f );
IntEdgeInd.Append ( iedge );
exactVertices.Add (sae.FirstVertex(E));
}
}
@ -1444,8 +1349,11 @@ Standard_Boolean ShapeFix_ComposeShell::SplitByLine (ShapeFix_WireSegment &wire,
Standard_Integer i;
for ( i = 1; i <= Inter.NbPoints(); i++ ) {
IntRes2d_IntersectionPoint IP = Inter.Point (i);
IntLinePar.Append ( IP.ParamOnFirst() );
IntEdgePar.Append ( IP.ParamOnSecond() );
if (IP.TransitionOfSecond().PositionOnCurve() == IntRes2d_Middle || (code != IOR_UNDEF && prevCode != IOR_UNDEF) )
{
IntLinePar.Append (IP.ParamOnFirst());
IntEdgePar.Append (IP.ParamOnSecond());
}
}
for ( i = 1; i <= Inter.NbSegments(); i++ ) {
IntRes2d_IntersectionSegment IS = Inter.Segment (i);
@ -1499,7 +1407,6 @@ Standard_Boolean ShapeFix_ComposeShell::SplitByLine (ShapeFix_WireSegment &wire,
IntLinePar.Append ( ParamPointsOnLine ( pos, firstPos, line ) );
IntEdgePar.Append ( isreversed ? f : l );
IntEdgeInd.Append ( iedge );
exactVertices.Add (sae.LastVertex(E));
}
}
}
@ -1637,7 +1544,7 @@ Standard_Boolean ShapeFix_ComposeShell::SplitByLine (ShapeFix_WireSegment &wire,
//=======================================
// Split edges in the wire by intersection points and fill vertices array
TopTools_SequenceOfShape IntVertices;
wire = SplitWire ( wire, IntEdgeInd, IntEdgePar, exactVertices, IntVertices,
wire = SplitWire ( wire, IntEdgeInd, IntEdgePar, IntVertices,
aNewSegCodes, isCutByU, cutIndex );
// add all data to input arrays
@ -1683,8 +1590,7 @@ void ShapeFix_ComposeShell::SplitByLine (ShapeFix_SequenceOfWireSegment &wires,
// merge null-length tangential segments into one-point tangencies or intersections
for ( i = 1; i < SplitLinePar.Length(); i++ ) {
Standard_Boolean isSameVertex = SplitLineVertex(i).IsSame(SplitLineVertex(i + 1));
if ( Abs ( SplitLinePar(i+1) - SplitLinePar(i) ) > ::Precision::PConfusion() && !isSameVertex) continue;
if ( Abs ( SplitLinePar(i+1) - SplitLinePar(i) ) > ::Precision::PConfusion() && !SplitLineVertex(i).IsSame(SplitLineVertex(i+1)) ) continue;
if ( ( SplitLineCode(i) & ITP_ENDSEG &&
SplitLineCode(i+1) & ITP_BEGSEG ) ||
( SplitLineCode(i) & ITP_BEGSEG &&
@ -1695,17 +1601,6 @@ void ShapeFix_ComposeShell::SplitByLine (ShapeFix_SequenceOfWireSegment &wires,
SplitLineCode.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
@ -1744,8 +1639,19 @@ void ShapeFix_ComposeShell::SplitByLine (ShapeFix_SequenceOfWireSegment &wires,
TopoDS_Shape tmpV2 = Context()->Apply ( SplitLineVertex(i) );
TopoDS_Vertex V1 = TopoDS::Vertex ( tmpV1 );
TopoDS_Vertex V2 = TopoDS::Vertex ( tmpV2 );
// protection against creating null-length edges
if ( SplitLinePar(i) - SplitLinePar(i-1) < ::Precision::PConfusion() ) {
// protection against creating null-length edges or edges lying inside tolerance of vertices
//first and last vertices for split line can not be merged to each other
Standard_Boolean canbeMerged = ( /*myClosedMode &&*/ (i -1 > 1 || i < SplitLinePar.Length()));
Standard_Real aMaxTol = MaxTolerance();
//case when max tolerance is not defined tolerance of vertices will be used as is
if( aMaxTol <= 2. *Precision::Confusion() )
aMaxTol = Precision::Infinite();
Standard_Real aTol1 = Min(BRep_Tool::Tolerance(V1), aMaxTol);
Standard_Real aTol2 = Min(BRep_Tool::Tolerance(V2), aMaxTol);
gp_Pnt aP1 = BRep_Tool::Pnt(V1);
gp_Pnt aP2 = BRep_Tool::Pnt(V2);
Standard_Real aD = aP1.SquareDistance(aP2);
if (SplitLinePar(i) - SplitLinePar(i-1) < ::Precision::PConfusion() || ( canbeMerged && ( aD <= (aTol1 * aTol1) || aD <= (aTol2 * aTol2)))) {// BRepTools::Compare(V1, V2)) ) {
#ifdef OCCT_DEBUG
std::cout << "Info: ShapeFix_ComposeShell::SplitByLine: Short segment ignored" << std::endl;
@ -2257,7 +2163,7 @@ void ShapeFix_ComposeShell::CollectWires (ShapeFix_SequenceOfWireSegment &wires,
// for coincidence (instead of vertex tolerance) in order
// this check to be in agreement with check for position of wire segments
// thus avoiding bad effects on overlapping edges
Standard_Real ctol = Max ( edgeTol, BRep_Tool::Tolerance(lastEdge) );
Standard_Real ctol = Max (edgeTol, BRep_Tool::Tolerance(endV/*lastEdge*/));
Standard_Boolean conn = IsCoincided ( endPnt, lPnt, myUResolution, myVResolution, ctol );
Standard_Real dist = endPnt.SquareDistance ( lPnt );

View File

@ -33,7 +33,6 @@
#include <TColStd_SequenceOfInteger.hxx>
#include <TColStd_SequenceOfReal.hxx>
#include <TopTools_SequenceOfShape.hxx>
#include <TopTools_MapOfShape.hxx>
class ShapeExtend_CompositeSurface;
class ShapeAnalysis_TransferParameters;
class ShapeExtend_WireData;
@ -99,10 +98,7 @@ public:
//! Here face defines both set of wires and way of getting
//! pcurves. Precision is used (together with tolerance of edges)
//! 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'
//! mode which forces ComposeShell to consider
@ -143,8 +139,7 @@ public:
//! and all pcurves on the initial (pseudo)face are reassigned to
//! that surface. If several wires are one inside another, single
//! 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.
Standard_EXPORT void SetTransferParamTool (const Handle(ShapeAnalysis_TransferParameters)& TransferParam);
@ -170,13 +165,7 @@ protected:
//! 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,
//! 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
//! parameters on them. Returns resulting wire and vertices
@ -189,14 +178,7 @@ protected:
//! NOTE: If edge is splitted, it is replaced by wire, and
//! order of edges in the wire corresponds to FORWARD orientation
//! of the edge.
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);
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);
//! Split edges in the wire by cutting line.
//! Wires with FORWARD or REVERSED orientation are considered
@ -209,13 +191,7 @@ protected:
//! Method fills sequences of parameters of intersection points
//! of cutting line with all edges, their types, and corresponding
//! 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.
//! Wires with FORWARD or REVERSED orientation are considered
@ -230,10 +206,7 @@ protected:
//! All modifications (splitting) are recorded in context,
//! except splitting of wires marked as EXTERNAL
//! (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
//! to joints between patches on the composite surface.
@ -258,15 +231,12 @@ protected:
//! taking EXTERNAL as necessary in fork points. Forks are detected
//! by common vertices. In fork point, most left way is seleccted
//! 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
//! (wires) into one or several faces depending on their mutual
//! 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;
TopoDS_Shape myResult;

View File

@ -1863,13 +1863,56 @@ Standard_Boolean ShapeFix_Face::FixMissingSeam()
mySurf = new ShapeAnalysis_Surface ( RTS );
myResult = CompShell.Result();
// if ( myFace.Orientation() == TopAbs_REVERSED ) res.Reverse();
Context()->Replace ( myFace, myResult );
// Remove small wires and / or faces that can be generated by ComposeShell
// (see tests bugs step bug30052_4, de step_3 E6)
Standard_Integer nbFaces = 0;
TopExp_Explorer expF ( myResult, TopAbs_FACE );
for (; expF.More(); expF.Next() )
{
TopoDS_Face aFace = TopoDS::Face(expF.Value());
TopExp_Explorer aExpW(aFace, TopAbs_WIRE);
Standard_Integer nbWires = 0;
for( ;aExpW.More(); aExpW.Next() )
{
ShapeFix_Wire aSfw(TopoDS::Wire(aExpW.Value()), aFace, Precision());
aSfw.SetContext(Context());
if(aSfw.NbEdges())
aSfw.FixSmall (Standard_True, Precision());
if(!aSfw.NbEdges())
{
Context()->Remove(aExpW.Value());
continue;
}
nbWires++;
}
if(!nbWires)
{
Context()->Remove(aFace);
continue;
}
nbFaces++;
}
myResult = Context()->Apply(myResult);
for (TopExp_Explorer exp ( myResult, TopAbs_FACE ); exp.More(); exp.Next() ) {
myFace = TopoDS::Face ( exp.Current() );
myFace = TopoDS::Face ( Context()->Apply(exp.Current() ));
if( myFace.IsNull())
continue;
if(nbFaces > 1)
{
FixSmallAreaWire(Standard_True);
TopoDS_Shape aShape = Context()->Apply(myFace);
if(aShape.IsNull() )
continue;
myFace = TopoDS::Face(aShape);
}
BRepTools::Update(myFace); //:p4
}
myResult = Context()->Apply(myResult);
SendWarning ( Message_Msg ( "FixAdvFace.FixMissingSeam.MSG0" ) );// Missing seam-edge added
return Standard_True;
}

View File

@ -394,9 +394,9 @@ Standard_Boolean ShapeFix_Wire::Perform()
if (myFixTailMode != 0)
{
Fixed |= FixTails();
if (Fixed)
if (FixTails())
{
Fixed =Standard_True;
FixShifted();
}
}

View File

@ -6,4 +6,4 @@ puts ""
pload XDE
stepread [locate_data_file bug30273.stp] res *
checknbshapes res_1 -solid 176 -face 9968 -shape 69497
checknbshapes res_1 -solid 176 -face 9968 -shape 69495

View File

@ -0,0 +1,10 @@
puts "========================"
puts "0030052: Data Exchange - on reading STEP produce invalid shapes"
puts "========================"
stepread [locate_data_file EADE1200A4AQJ.STP] a *
renamevar a_1 result
checkshape result
checkmaxtol result -ref 0.08108
checkview -display result -2d -path ${imagedir}/${test_image}.png

View File

@ -0,0 +1,10 @@
puts "========================"
puts "0030052: Data Exchange -"
puts "========================"
stepread [locate_data_file KM63UTFBHS116.STP] a *
renamevar a_1 result
checkshape result
checkmaxtol result -ref 0.04829
checkview -display result -2d -path ${imagedir}/${test_image}.png

View File

@ -0,0 +1,10 @@
puts "========================"
puts "0030052: Data Exchange "
puts "========================"
stepread [locate_data_file TM711MF160X150R2DHA.STP] a *
renamevar a_1 result
checkshape result
checkmaxtol result -ref 0.051677
checkview -display result -2d -path ${imagedir}/${test_image}.png

View File

@ -0,0 +1,10 @@
puts "========================"
puts "0030052: Data Exchange "
puts "========================"
stepread [locate_data_file C11816W.STP] a *
renamevar a_1 result
checkshape result
checkmaxtol result -ref 0.019014
checkview -display result -2d -path ${imagedir}/${test_image}.png

View File

@ -0,0 +1,11 @@
puts "========================"
puts "0030052: "
puts "========================"
stepread [locate_data_file KNSU040R08X16S100.STP] a *
renamevar a_1 result
checkshape result
checkmaxtol result -ref 0.043776
checkview -display result -2d -path ${imagedir}/${test_image}.png

5
tests/bugs/step/bug31301 Normal file
View File

@ -0,0 +1,5 @@
testreadstep [locate_data_file bug31301.stp] s
testwritestep bug31301_1.stp s
testreadstep bug31301_1.stp s1
checkshape s1 f
checkmaxtol s1 -ref 0.0078

View File

@ -1,16 +1,13 @@
# !!!! This file is generated automatically, do not edit manually! See end script
puts "TODO CR23096 ALL: CHECKSHAPE : Faulty"
set filename r0701_ug.stp
set ref_data {
DATA : Faulties = 0 ( 0 ) Warnings = 0 ( 0 ) Summary = 0 ( 0 )
TPSTAT : Faulties = 0 ( 0 ) Warnings = 4 ( 20 ) Summary = 4 ( 20 )
CHECKSHAPE : Wires = 2 ( 0 ) Faces = 2 ( 0 ) Shells = 0 ( 0 ) Solids = 0 ( 0 )
TPSTAT : Faulties = 0 ( 0 ) Warnings = 3 ( 20 ) Summary = 3 ( 20 )
CHECKSHAPE : Wires = 0 ( 0 ) Faces = 0 ( 0 ) Shells = 0 ( 0 ) Solids = 0 ( 0 )
NBSHAPES : Solid = 1 ( 1 ) Shell = 1 ( 1 ) Face = 26 ( 26 )
STATSHAPE : Solid = 1 ( 1 ) Shell = 1 ( 1 ) Face = 26 ( 26 ) FreeWire = 0 ( 0 )
TOLERANCE : MaxTol = 1.497084727e-005 ( 1.497084651e-005 ) AvgTol = 2.648852507e-006 ( 7.773266296e-006 )
TOLERANCE : MaxTol = 1.497084727e-005 ( 1.497084651e-005 ) AvgTol = 2.422471639e-006 ( 7.787901742e-006 )
LABELS : N0Labels = 1 ( 1 ) N1Labels = 0 ( 0 ) N2Labels = 0 ( 0 ) TotalLabels = 1 ( 1 ) NameLabels = 1 ( 1 ) ColorLabels = 1 ( 1 ) LayerLabels = 0 ( 0 )
PROPS : Centroid = 0 ( 0 ) Volume = 0 ( 0 ) Area = 0 ( 0 )
NCOLORS : NColors = 1 ( 1 )

View File

@ -1,16 +1,13 @@
# !!!! This file is generated automatically, do not edit manually! See end script
puts "TODO CR23096 ALL: CHECKSHAPE : Faulty"
set filename tr10_r0301_ug.stp
set ref_data {
DATA : Faulties = 0 ( 0 ) Warnings = 0 ( 0 ) Summary = 0 ( 0 )
TPSTAT : Faulties = 0 ( 0 ) Warnings = 4 ( 20 ) Summary = 4 ( 20 )
CHECKSHAPE : Wires = 2 ( 0 ) Faces = 2 ( 0 ) Shells = 0 ( 0 ) Solids = 0 ( 0 )
TPSTAT : Faulties = 0 ( 0 ) Warnings = 3 ( 20 ) Summary = 3 ( 20 )
CHECKSHAPE : Wires = 0 ( 0 ) Faces = 0 ( 0 ) Shells = 0 ( 0 ) Solids = 0 ( 0 )
NBSHAPES : Solid = 1 ( 1 ) Shell = 1 ( 1 ) Face = 26 ( 26 )
STATSHAPE : Solid = 1 ( 1 ) Shell = 1 ( 1 ) Face = 26 ( 26 ) FreeWire = 0 ( 0 )
TOLERANCE : MaxTol = 5.560731313e-005 ( 5.560731145e-005 ) AvgTol = 5.452525023e-006 ( 1.059579749e-005 )
TOLERANCE : MaxTol = 5.560731313e-005 ( 5.560731145e-005 ) AvgTol = 5.245002938e-006 ( 1.0648922e-005 )
LABELS : N0Labels = 1 ( 1 ) N1Labels = 0 ( 0 ) N2Labels = 0 ( 0 ) TotalLabels = 1 ( 1 ) NameLabels = 1 ( 1 ) ColorLabels = 1 ( 1 ) LayerLabels = 0 ( 0 )
PROPS : Centroid = 0 ( 0 ) Volume = 0 ( 0 ) Area = 0 ( 0 )
NCOLORS : NColors = 1 ( 1 )

View File

@ -1,17 +1,13 @@
# !!!! This file is generated automatically, do not edit manually! See end script
puts "TODO CR23096 ALL: CHECKSHAPE : Faulty"
puts "TODO CR23096 ALL: TOLERANCE : Faulty"
set filename trj7_pm5-ug-214.stp
set ref_data {
DATA : Faulties = 0 ( 0 ) Warnings = 0 ( 0 ) Summary = 0 ( 0 )
TPSTAT : Faulties = 0 ( 0 ) Warnings = 5 ( 184 ) Summary = 5 ( 184 )
CHECKSHAPE : Wires = 2 ( 0 ) Faces = 2 ( 0 ) Shells = 0 ( 0 ) Solids = 0 ( 0 )
TPSTAT : Faulties = 0 ( 0 ) Warnings = 3 ( 184 ) Summary = 3 ( 184 )
CHECKSHAPE : Wires = 0 ( 0 ) Faces = 0 ( 0 ) Shells = 0 ( 0 ) Solids = 0 ( 0 )
NBSHAPES : Solid = 1 ( 1 ) Shell = 1 ( 1 ) Face = 236 ( 236 )
STATSHAPE : Solid = 1 ( 1 ) Shell = 1 ( 1 ) Face = 236 ( 236 ) FreeWire = 0 ( 0 )
TOLERANCE : MaxTol = 0.9036712735 ( 0.01854047103 ) AvgTol = 0.002850040066 ( 0.0006030086632 )
TOLERANCE : MaxTol = 0.004380869792 ( 0.01854047103 ) AvgTol = 0.0001630183528 ( 0.0006019917644 )
LABELS : N0Labels = 1 ( 1 ) N1Labels = 34 ( 34 ) N2Labels = 0 ( 0 ) TotalLabels = 35 ( 35 ) NameLabels = 1 ( 1 ) ColorLabels = 35 ( 35 ) LayerLabels = 1 ( 1 )
PROPS : Centroid = 0 ( 0 ) Volume = 0 ( 0 ) Area = 0 ( 0 )
NCOLORS : NColors = 11 ( 11 )

View File

@ -9,11 +9,11 @@ set filename Z8M6SAT.stp
set ref_data {
DATA : Faulties = 0 ( 0 ) Warnings = 0 ( 0 ) Summary = 0 ( 0 )
TPSTAT : Faulties = 0 ( 0 ) Warnings = 922 ( 3168 ) Summary = 922 ( 3168 )
CHECKSHAPE : Wires = 45 ( 42 ) Faces = 46 ( 40 ) Shells = 0 ( 2 ) Solids = 0 ( 0 )
NBSHAPES : Solid = 28 ( 28 ) Shell = 768 ( 30 ) Face = 3240 ( 3239 )
STATSHAPE : Solid = 28 ( 28 ) Shell = 768 ( 30 ) Face = 3240 ( 3239 ) FreeWire = 0 ( 0 )
TOLERANCE : MaxTol = 15.00151489 ( 20.46526799 ) AvgTol = 0.02219114253 ( 0.03611458186 )
TPSTAT : Faulties = 0 ( 0 ) Warnings = 911 ( 3168 ) Summary = 911 ( 3168 )
CHECKSHAPE : Wires = 43 ( 40 ) Faces = 44 ( 39 ) Shells = 0 ( 4 ) Solids = 0 ( 0 )
NBSHAPES : Solid = 28 ( 28 ) Shell = 769 ( 32 ) Face = 3241 ( 3241 )
STATSHAPE : Solid = 28 ( 28 ) Shell = 769 ( 32 ) Face = 3241 ( 3241 ) FreeWire = 0 ( 0 )
TOLERANCE : MaxTol = 15.00151489 ( 20.46526799 ) AvgTol = 0.01760149468 ( 0.03612886317 )
LABELS : N0Labels = 3 ( 3 ) N1Labels = 2 ( 2 ) N2Labels = 0 ( 0 ) TotalLabels = 5 ( 5 ) NameLabels = 5 ( 5 ) ColorLabels = 0 ( 0 ) LayerLabels = 0 ( 0 )
PROPS : Centroid = 0 ( 0 ) Volume = 0 ( 0 ) Area = 0 ( 0 )
NCOLORS : NColors = 0 ( 0 )