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

0026441: Modeling Algorithms - BRepOffset_MakeOffset affects original shape

BRepOffset_MakeOffset.cxx - tolerance control for building planar faces is implemented,
                            updating tolerance for initial entities is avoided

BRepAlgo_Loop.cxx - "total" setting tolerance 0.001 is removing

BRepLib.cxx - checking of "locked" for vertex is removing in static function UpdShTol

QABugs_20.cxx - add new command OCC26441 for checking tolerance differenses between two "identical" shapes

tests/bugs/modalg_7/bug30054 - case now is "BAD", because really result shape is invalid: many faces has not closed wires with huge gaps between ends of edges. Result was "OK" only because tolerances of vertices were increased by algorithm to cover all gaps.

tests/bugs/modalg_8/bug26441 - new test case added

Other test: B3, C8, A7, C8: they were "BAD" and now are "BAD", only some problems are changed.
This commit is contained in:
ifv 2022-09-30 15:45:10 +03:00 committed by Vadim Glukhikh
parent 8175a70c4e
commit 91a2f58f8f
11 changed files with 333 additions and 98 deletions

View File

@ -42,12 +42,14 @@
#include <TopTools_SequenceOfShape.hxx>
#include <stdio.h>
//#define OCCT_DEBUG_ALGO
//#define DRAW
#ifdef DRAW
#include <DBRep.hxx>
#pragma comment(lib,"TKDraw")
#endif
#ifdef OCCT_DEBUG_ALGO
Standard_Boolean AffichLoop = Standard_False;
Standard_Boolean AffichLoop = Standard_True;
Standard_Integer NbLoops = 0;
Standard_Integer NbWires = 1;
static char* name = new char[100];
@ -58,7 +60,8 @@ static char* name = new char[100];
//purpose :
//=======================================================================
BRepAlgo_Loop::BRepAlgo_Loop()
BRepAlgo_Loop::BRepAlgo_Loop():
myTolConf (0.001)
{
}
@ -185,7 +188,6 @@ static TopoDS_Vertex UpdateClosedEdge(const TopoDS_Edge& E,
Standard_Boolean OnStart = 0, OnEnd = 0;
//// modified by jgv, 13.04.04 for OCC5634 ////
TopExp::Vertices (E,V1,V2);
//Standard_Real Tol = Precision::Confusion();
Standard_Real Tol = BRep_Tool::Tolerance( V1 );
///////////////////////////////////////////////
@ -427,13 +429,12 @@ static void StoreInMVE (const TopoDS_Face& F,
TopoDS_Edge& E,
TopTools_IndexedDataMapOfShapeListOfShape& MVE,
Standard_Boolean& YaCouture,
TopTools_DataMapOfShapeShape& VerticesForSubstitute )
TopTools_DataMapOfShapeShape& VerticesForSubstitute,
const Standard_Real theTolConf)
{
TopoDS_Vertex V1, V2, V;
TopTools_ListOfShape Empty;
Standard_Real Tol = 0.001; //5.e-05; //5.e-07;
// gp_Pnt P1, P2, P;
gp_Pnt P1, P;
BRep_Builder BB;
for (Standard_Integer iV = 1; iV <= MVE.Extent(); iV++)
@ -449,7 +450,7 @@ static void StoreInMVE (const TopoDS_Face& F,
{
V1 = TopoDS::Vertex( itl.Value() );
P1 = BRep_Tool::Pnt( V1 );
if (P.IsEqual( P1, Tol ) && !V.IsSame(V1))
if (P.IsEqual( P1, theTolConf ) && !V.IsSame(V1))
{
V.Orientation( V1.Orientation() );
if (VerticesForSubstitute.IsBound( V1 ))
@ -574,7 +575,7 @@ void BRepAlgo_Loop::Perform()
TopoDS_Edge& E = TopoDS::Edge(itl1.Value());
if (!Emap.Add(E))
continue;
StoreInMVE(myFace,E,MVE,YaCouture,myVerticesForSubstitute);
StoreInMVE(myFace,E,MVE,YaCouture,myVerticesForSubstitute, myTolConf);
}
}
}
@ -586,7 +587,7 @@ void BRepAlgo_Loop::Perform()
for (itl.Initialize(myConstEdges); itl.More(); itl.Next()) {
TopoDS_Edge& E = TopoDS::Edge(itl.Value());
if (DejaVu.Add(E))
StoreInMVE(myFace,E,MVE,YaCouture,myVerticesForSubstitute);
StoreInMVE(myFace,E,MVE,YaCouture,myVerticesForSubstitute, myTolConf);
}
#ifdef DRAW
@ -626,42 +627,42 @@ void BRepAlgo_Loop::Perform()
//--------------------------------
RemovePendingEdges(MVE);
if (MVE.Extent() == 0) break;
if (MVE.Extent() == 0) break;
//--------------------------------
// Start edge.
//--------------------------------
EF = CE = TopoDS::Edge(MVE(1).First());
TopExp::Vertices(CE,V1,V2);
TopExp::Vertices(CE, V1, V2);
//--------------------------------
// VF vertex start of new wire
//--------------------------------
if (CE.Orientation() == TopAbs_FORWARD) { CV = VF = V1;}
else { CV = VF = V2;}
if (CE.Orientation() == TopAbs_FORWARD) { CV = VF = V1; }
else { CV = VF = V2; }
if (!MVE.Contains(CV)) continue;
TopTools_ListOfShape& aListEdges = MVE.ChangeFromKey(CV);
for ( itl.Initialize(aListEdges); itl.More(); itl.Next()) {
for (itl.Initialize(aListEdges); itl.More(); itl.Next()) {
if (itl.Value().IsEqual(CE)) {
aListEdges.Remove(itl);
break;
aListEdges.Remove(itl);
break;
}
}
End = Standard_False;
End = Standard_False;
while (!End) {
//-------------------------------
// Construction of a wire.
//-------------------------------
TopExp::Vertices(CE,V1,V2);
TopExp::Vertices(CE, V1, V2);
if (!CV.IsSame(V1)) CV = V1; else CV = V2;
B.Add (NW,CE);
B.Add(NW, CE);
UsedEdges.Add(CE);
if (!MVE.Contains(CV) || MVE.FindFromKey(CV).IsEmpty()) {
End = Standard_True;
}
else {
End = !SelectEdge(myFace,CE,CV,NE,MVE.ChangeFromKey(CV));
End = !SelectEdge(myFace, CE, CV, NE, MVE.ChangeFromKey(CV));
if (!End) {
CE = NE;
if (MVE.FindFromKey(CV).IsEmpty())
@ -672,35 +673,41 @@ void BRepAlgo_Loop::Perform()
//--------------------------------------------------
// Add new wire to the set of wires
//------------------------------------------------
Standard_Real Tol = 0.001; //5.e-05; //5.e-07;
TopExp_Explorer explo( NW, TopAbs_VERTEX );
for (; explo.More(); explo.Next())
{
const TopoDS_Vertex& aV = TopoDS::Vertex( explo.Current() );
Handle(BRep_TVertex)& TV = *((Handle(BRep_TVertex)*) &(aV).TShape());
TV->Tolerance( Tol );
TV->Modified( Standard_True );
}
for (explo.Init( NW, TopAbs_EDGE ); explo.More(); explo.Next())
{
const TopoDS_Edge& aE = TopoDS::Edge( explo.Current() );
Handle(BRep_TEdge)& TE = *((Handle(BRep_TEdge)*) &(aE).TShape());
TE->Tolerance( Tol );
TE->Modified( Standard_True );
}
if (VF.IsSame(CV) && SamePnt2d(VF,EF,CE,myFace))
if (VF.IsSame(CV))
{
NW.Closed (Standard_True);
myNewWires.Append (NW);
if (SamePnt2d(VF, EF, CE, myFace))
{
NW.Closed(Standard_True);
myNewWires.Append(NW);
}
else if(BRep_Tool::Tolerance(VF) < myTolConf)
{
BRep_Builder aBB;
aBB.UpdateVertex(VF, myTolConf);
if (SamePnt2d(VF, EF, CE, myFace))
{
NW.Closed(Standard_True);
myNewWires.Append(NW);
}
#ifdef OCCT_DEBUG_ALGO
else
{
std::cout << "BRepAlgo_Loop: Open Wire" << std::endl;
if (AffichLoop)
std::cout << "OpenWire is : NW_" << NbLoops << "_" << NbWires << std::endl;
}
#endif
}
}
#ifdef OCCT_DEBUG_ALGO
else {
std::cout <<"BRepAlgo_Loop: Open Wire"<<std::endl;
std::cout << "BRepAlgo_Loop: Open Wire" << std::endl;
if (AffichLoop)
std::cout << "OpenWire is : NW_"<<NbLoops<<"_"<<NbWires<<std::endl;
}
std::cout << "OpenWire is : NW_" << NbLoops << "_" << NbWires << std::endl;
}
#endif
#ifdef DRAW
if (AffichLoop) {
sprintf(name,"NW_%d_%d",NbLoops,NbWires++);
@ -777,8 +784,6 @@ void BRepAlgo_Loop::CutEdge (const TopoDS_Edge& E,
VF = TopoDS::Vertex(aLocalV);
aLocalV = VCEI.Oriented(TopAbs_REVERSED);
VL = TopoDS::Vertex(aLocalV);
// VF = TopoDS::Vertex(VCEI.Oriented(TopAbs_FORWARD));
// VL = TopoDS::Vertex(VCEI.Oriented(TopAbs_REVERSED));
}
SV.Prepend(VF);
SV.Append(VL);
@ -813,13 +818,9 @@ void BRepAlgo_Loop::CutEdge (const TopoDS_Edge& E,
B.Add (NewEdge,aLocalEdge);
aLocalEdge = V2.Oriented(TopAbs_REVERSED);
B.Add (TopoDS::Edge(NewEdge),aLocalEdge);
// B.Add (NewEdge,V1.Oriented(TopAbs_FORWARD));
// B.Add (NewEdge,V2.Oriented(TopAbs_REVERSED));
if (V1.IsSame(VF))
U1 = f;
else
// U1=BRep_Tool::Parameter
// (TopoDS::Vertex(V1.Oriented(TopAbs_INTERNAL)),WE);
{
TopoDS_Shape aLocalV = V1.Oriented(TopAbs_INTERNAL);
U1=BRep_Tool::Parameter(TopoDS::Vertex(aLocalV),WE);
@ -830,8 +831,6 @@ void BRepAlgo_Loop::CutEdge (const TopoDS_Edge& E,
{
TopoDS_Shape aLocalV = V2.Oriented(TopAbs_INTERNAL);
U2=BRep_Tool::Parameter(TopoDS::Vertex(aLocalV),WE);
// U2=BRep_Tool::Parameter
// (TopoDS::Vertex(V2.Oriented(TopAbs_INTERNAL)),WE);
}
B.Range (TopoDS::Edge(NewEdge),U1,U2);
#ifdef DRAW

View File

@ -86,8 +86,17 @@ public:
Standard_EXPORT void VerticesForSubstitute (TopTools_DataMapOfShapeShape& VerVerMap);
//! Set maximal tolerance used for comparing distaces between vertices.
void SetTolConf(const Standard_Real theTolConf)
{
myTolConf = theTolConf;
}
//! Get maximal tolerance used for comparing distaces between vertices.
Standard_Real GetTolConf() const
{
return myTolConf;
}
protected:
@ -108,6 +117,7 @@ private:
TopTools_DataMapOfShapeListOfShape myCutEdges;
TopTools_DataMapOfShapeShape myVerticesForSubstitute;
BRepAlgo_Image myImageVV;
Standard_Real myTolConf;
};

View File

@ -802,7 +802,7 @@ static void GetEdgeTol(const TopoDS_Edge& theEdge,
}
if(temp > d2) d2 = temp;
}
d2 = 1.5*sqrt(d2);
d2 = 1.05*sqrt(d2);
theEdTol = d2;
}
@ -884,10 +884,6 @@ static void UpdShTol(const TopTools_DataMapOfShapeReal& theShToTol,
case TopAbs_VERTEX:
{
const Handle(BRep_TVertex)& aTV = *((Handle(BRep_TVertex)*)&aNsh.TShape());
//
if(aTV->Locked())
throw TopoDS_LockedShape("BRep_Builder::UpdateVertex");
//
if (theVForceUpdate)
aTV->Tolerance(aTol);
else
@ -1709,8 +1705,8 @@ static void InternalUpdateTolerances(const TopoDS_Shape& theOldShape,
for (iCur=1; iCur<=nbV; iCur++) {
tol=0;
const TopoDS_Vertex& V = TopoDS::Vertex(parents.FindKey(iCur));
Bnd_Box box;
box.Add(BRep_Tool::Pnt(V));
gp_Pnt aPV = BRep_Tool::Pnt(V);
Standard_Real aMaxDist = 0.;
gp_Pnt p3d;
for (lConx.Initialize(parents(iCur)); lConx.More(); lConx.Next()) {
const TopoDS_Edge& E = TopoDS::Edge(lConx.Value());
@ -1732,8 +1728,10 @@ static void InternalUpdateTolerances(const TopoDS_Shape& theOldShape,
if (!C.IsNull()) { // edge non degenerated
p3d = C->Value(par);
p3d.Transform(L.Transformation());
box.Add(p3d);
}
Standard_Real aDist = p3d.SquareDistance(aPV);
if (aDist > aMaxDist)
aMaxDist = aDist;
}
}
else if (cr->IsCurveOnSurface()) {
const Handle(Geom_Surface)& Su = cr->Surface();
@ -1745,21 +1743,22 @@ static void InternalUpdateTolerances(const TopoDS_Shape& theOldShape,
gp_Pnt2d p2d = PC->Value(par);
p3d = Su->Value(p2d.X(),p2d.Y());
p3d.Transform(L.Transformation());
box.Add(p3d);
Standard_Real aDist = p3d.SquareDistance(aPV);
if (aDist > aMaxDist)
aMaxDist = aDist;
if (!PC2.IsNull()) {
p2d = PC2->Value(par);
p3d = Su->Value(p2d.X(),p2d.Y());
p3d.Transform(L.Transformation());
box.Add(p3d);
aDist = p3d.SquareDistance(aPV);
if (aDist > aMaxDist)
aMaxDist = aDist;
}
}
itcr.Next();
}
}
Standard_Real aXmin, aYmin, aZmin, aXmax, aYmax, aZmax;
box.Get(aXmin, aYmin, aZmin, aXmax, aYmax, aZmax);
aXmax -= aXmin; aYmax -= aYmin; aZmax -= aZmin;
tol = Max(tol,sqrt(aXmax*aXmax+aYmax*aYmax+aZmax*aZmax));
tol = Max(tol, sqrt(aMaxDist));
tol += 2.*Epsilon(tol);
//
Standard_Real aVTol = BRep_Tool::Tolerance(V);

View File

@ -336,8 +336,10 @@ static BRepOffset_Error checkSinglePoint(const Standard_Real theUParam,
const NCollection_Vector<gp_Pnt>& theBadPoints);
//---------------------------------------------------------------------
static void UpdateTolerance ( TopoDS_Shape& myShape,
const TopTools_IndexedMapOfShape& myFaces);
static void UpdateTolerance ( TopoDS_Shape& theShape,
const TopTools_IndexedMapOfShape& theFaces,
const TopoDS_Shape& theInitShape);
static Standard_Real ComputeMaxDist(const gp_Pln& thePlane,
const Handle(Geom_Curve)& theCrv,
const Standard_Real theFirst,
@ -1036,8 +1038,16 @@ void BRepOffset_MakeOffset::MakeOffsetShape(const Message_ProgressRange& theRang
// MAJ Tolerance edge and Vertex
// ----------------------------
if (!myOffsetShape.IsNull()) {
UpdateTolerance (myOffsetShape,myFaces);
BRepLib::UpdateTolerances( myOffsetShape );
if (myThickening)
{
UpdateTolerance(myOffsetShape, myFaces, myShape);
}
else
{
TopoDS_Shape aDummy;
UpdateTolerance(myOffsetShape, myFaces, aDummy);
}
BRepLib::UpdateTolerances(myOffsetShape);
}
CorrectConicalFaces();
@ -3165,13 +3175,38 @@ void BRepOffset_MakeOffset::MakeMissingWalls (const Message_ProgressRange& theRa
} //if both edges are arcs of circles
if (NewFace.IsNull())
{
BRepLib_MakeFace MF(theWire, Standard_True); //Only plane
if (MF.Error() == BRepLib_FaceDone)
Standard_Real anEdgeTol = BRep_Tool::Tolerance(anEdge);
//Tolerances of input shape should not be increased by BRepLib_MakeFace
BRepLib_FindSurface aFindPlane(theWire, anEdgeTol, Standard_True); //only plane
IsPlanar = Standard_False;
if(aFindPlane.Found() && aFindPlane.ToleranceReached() <= anEdgeTol)
{
NewFace = MF.Face();
IsPlanar = Standard_True;
Standard_Real f, l;
Handle(Geom_Curve) aGC = BRep_Tool::Curve(anEdge, f, l);
Handle(Geom_Plane) aPln = Handle(Geom_Plane)::DownCast(aFindPlane.Surface());
Standard_Real aMaxDist = ComputeMaxDist(aPln->Pln(), aGC, f, l);
if (aMaxDist <= anEdgeTol)
{
BRepLib_MakeFace MF(aPln->Pln(), theWire);
if (MF.IsDone())
{
NewFace = MF.Face();
TopoDS_Iterator anItE(theWire);
for (; anItE.More(); anItE.Next())
{
const TopoDS_Edge& anE = TopoDS::Edge(anItE.Value());
if (anE.IsSame(anEdge))
continue;
aGC = BRep_Tool::Curve(anE, f, l);
aMaxDist = ComputeMaxDist(aPln->Pln(), aGC, f, l);
BB.UpdateEdge(anE, aMaxDist);
}
IsPlanar = Standard_True;
}
}
}
else //Extrusion (by thrusections)
//
if(!IsPlanar) //Extrusion (by thrusections)
{
Handle(Geom_Curve) EdgeCurve = BRep_Tool::Curve(anEdge, fpar, lpar);
Handle(Geom_TrimmedCurve) TrEdgeCurve =
@ -3185,7 +3220,6 @@ void BRepOffset_MakeOffset::MakeMissingWalls (const Message_ProgressRange& theRa
ThrusecGenerator.AddCurve( TrOffsetCurve );
ThrusecGenerator.Perform( Precision::PConfusion() );
theSurf = ThrusecGenerator.Surface();
//theSurf = new Geom_SurfaceOfLinearExtrusion( TrOffsetCurve, OffsetDir );
Standard_Real Uf, Ul, Vf, Vl;
theSurf->Bounds(Uf, Ul, Vf, Vl);
TopLoc_Location Loc;
@ -3272,8 +3306,14 @@ void BRepOffset_MakeOffset::MakeMissingWalls (const Message_ProgressRange& theRa
BB.Range( anE3, FirstPar, LastPar );
}
}
BRepLib::SameParameter(NewFace);
BRepTools::Update(NewFace);
if (!IsPlanar)
{
// For planar faces these operations are useless,
// because there are no curves on surface
BRepLib::SameParameter(NewFace);
BRepTools::Update(NewFace);
}
//Check orientation
TopAbs_Orientation anOr = OrientationOfEdgeInFace(anEdge, aFaceOfEdge);
TopAbs_Orientation OrInNewFace = OrientationOfEdgeInFace(anEdge, NewFace);
@ -3781,6 +3821,7 @@ void BRepOffset_MakeOffset::EncodeRegularity ()
#endif
}
//=======================================================================
//function : ComputeMaxDist
//purpose :
@ -3807,13 +3848,15 @@ Standard_Real ComputeMaxDist(const gp_Pln& thePlane,
}
return sqrt(aMaxDist)*1.05;
}
//=======================================================================
//function : UpDateTolerance
//function : UpdateTolerance
//purpose :
//=======================================================================
void UpdateTolerance (TopoDS_Shape& S,
const TopTools_IndexedMapOfShape& Faces)
const TopTools_IndexedMapOfShape& Faces,
const TopoDS_Shape& theInitShape)
{
BRep_Builder B;
TopTools_MapOfShape View;
@ -3829,12 +3872,31 @@ void UpdateTolerance (TopoDS_Shape& S,
}
}
Standard_Real Tol;
TopExp_Explorer ExpF;
for (ExpF.Init(S, TopAbs_FACE); ExpF.More(); ExpF.Next())
// The edges of initial shape are not modified
TopTools_MapOfShape aMapInitF;
if (!theInitShape.IsNull())
{
const TopoDS_Shape& F = ExpF.Current();
if (Faces.Contains(F))
TopExp_Explorer anExpF(theInitShape, TopAbs_FACE);
for (; anExpF.More(); anExpF.Next()) {
aMapInitF.Add(anExpF.Current());
TopExp_Explorer anExpE;
for (anExpE.Init(anExpF.Current(), TopAbs_EDGE); anExpE.More(); anExpE.Next()) {
View.Add(anExpE.Current());
TopoDS_Iterator anItV(anExpE.Current());
for (; anItV.More(); anItV.Next())
{
View.Add(anItV.Value());
}
}
}
}
Standard_Real Tol;
TopExp_Explorer anExpF(S, TopAbs_FACE);
for (; anExpF.More(); anExpF.Next())
{
const TopoDS_Shape& F = anExpF.Current();
if (Faces.Contains(F) || aMapInitF.Contains(F))
{
continue;
}
@ -3843,6 +3905,7 @@ void UpdateTolerance (TopoDS_Shape& S,
for (Exp.Init(F, TopAbs_EDGE); Exp.More(); Exp.Next()) {
TopoDS_Edge E = TopoDS::Edge(Exp.Current());
Standard_Boolean isUpdated = Standard_False;
Standard_Real aCurrTol = BRep_Tool::Tolerance(E);
if (aBAS.GetType() == GeomAbs_Plane)
{
//Edge does not seem to have pcurve on plane,
@ -3850,17 +3913,22 @@ void UpdateTolerance (TopoDS_Shape& S,
Standard_Real aFirst, aLast;
Handle(Geom_Curve) aCrv = BRep_Tool::Curve(E, aFirst, aLast);
Standard_Real aMaxDist = ComputeMaxDist(aBAS.Plane(), aCrv, aFirst, aLast);
E.Locked (Standard_False);
B.UpdateEdge(E, aMaxDist);
isUpdated = Standard_True;
if (aMaxDist > aCurrTol)
{
B.UpdateEdge(E, aMaxDist);
isUpdated = Standard_True;
}
}
if (View.Add(E))
{
E.Locked(Standard_False);
BRepCheck_Edge EdgeCorrector(E);
Tol = EdgeCorrector.Tolerance();
B.UpdateEdge(E, Tol);
isUpdated = Standard_True;
if (Tol > aCurrTol)
{
B.UpdateEdge(E, Tol);
isUpdated = Standard_True;
}
}
if (isUpdated)
{
@ -3869,11 +3937,11 @@ void UpdateTolerance (TopoDS_Shape& S,
TopExp::Vertices(E, V[0], V[1]);
for (Standard_Integer i = 0; i <= 1; i++) {
V[i].Locked(Standard_False);
if (View.Add(V[i])) {
Handle(BRep_TVertex) TV = Handle(BRep_TVertex)::DownCast(V[i].TShape());
TV->Tolerance(0.);
BRepCheck_Vertex VertexCorrector(V[i]);
V[i].Locked (Standard_False);
B.UpdateVertex(V[i], VertexCorrector.Tolerance());
// use the occasion to clean the vertices.
(TV->ChangePoints()).Clear();

View File

@ -4359,6 +4359,120 @@ static Standard_Integer QACheckBends(Draw_Interpretor& theDI,
return 0;
}
static Standard_Integer OCC26441(Draw_Interpretor& theDi, Standard_Integer theNbArgs, const char** theArgVec)
{
if (theNbArgs < 3)
{
theDi << "Syntax error: wrong number of arguments!\n";
return 1;
}
const TopoDS_Shape& aShape = DBRep::Get(theArgVec[1]);
if (aShape.IsNull())
{
theDi << " Null Shape is not allowed here\n";
return 1;
}
const TopoDS_Shape& aRefShape = DBRep::Get(theArgVec[2]);
if (aRefShape.IsNull())
{
theDi << " Null Shape is not allowed here\n";
return 1;
}
if(aShape.ShapeType() != aRefShape.ShapeType())
{
theDi << " Shape types are not the same\n";
return 1;
}
Standard_Real anEps = Precision::Confusion();
if (theNbArgs > 3)
{
anEps = Draw::Atof(theArgVec[3]);
}
Standard_Boolean isAllDiff = Standard_False;
if (theNbArgs > 4)
{
Standard_Integer Inc = Draw::Atoi(theArgVec[4]);
if (Inc > 0)
isAllDiff = Standard_True;
}
BRep_Builder aBB;
TopExp_Explorer anExp, anExpRef;
Standard_Real aMaxE = 0., aMaxV = 0.;
TopTools_MapOfShape aChecked;
TopoDS_Vertex aV[2], aRefV[2];
//Checking edge and vertex tolerances
TopoDS_Compound aBadEdges;
aBB.MakeCompound(aBadEdges);
TopoDS_Compound aBadVerts;
aBB.MakeCompound(aBadVerts);
anExp.Init(aShape, TopAbs_EDGE);
anExpRef.Init(aRefShape, TopAbs_EDGE);
for (; anExpRef.More(); anExpRef.Next())
{
const TopoDS_Edge& aRefE = TopoDS::Edge(anExpRef.Current());
if (!anExp.More())
{
theDi << " Different number of edges\n";
return 1;
}
const TopoDS_Edge& anE = TopoDS::Edge(anExp.Current());
if (!aChecked.Add(anE))
continue;
Standard_Real aTolRef = BRep_Tool::Tolerance(aRefE);
Standard_Real aTol = BRep_Tool::Tolerance(anE);
Standard_Real aDiff = aTol - aTolRef;
if (isAllDiff && aDiff < 0)
aDiff = -aDiff;
if (aDiff > anEps)
{
if (aDiff > aMaxE)
aMaxE = aDiff;
aBB.Add(aBadEdges, anE);
}
TopExp::Vertices(aRefE, aRefV[0], aRefV[1]);
TopExp::Vertices(anE, aV[0], aV[1]);
for (int i = 0; i < 2; ++i)
{
if (aRefV[i].IsNull())
continue;
if (!aChecked.Add(aV[i]))
continue;
aTolRef = BRep_Tool::Tolerance(aRefV[i]);
aTol = BRep_Tool::Tolerance(aV[i]);
aDiff = aTol - aTolRef;
if (aDiff > anEps)
{
if (aDiff > aMaxV)
aMaxV = aDiff;
aBB.Add(aBadVerts, aV[i]);
}
}
anExp.Next();
}
if (aMaxE > anEps)
{
theDi << " Maximal difference for edges : " << aMaxE << "\n";
DBRep::Set("BadEdges", aBadEdges);
}
if (aMaxV > anEps)
{
theDi << " Maximal difference for vertices : " << aMaxV << "\n";
DBRep::Set("BadVerts", aBadVerts);
}
return 0;
}
void QABugs::Commands_20(Draw_Interpretor& theCommands) {
@ -4466,6 +4580,10 @@ void QABugs::Commands_20(Draw_Interpretor& theCommands) {
"QACheckBends curve [CosMaxAngle [theNbPoints]]",
__FILE__,
QACheckBends, group);
theCommands.Add("OCC26441",
"OCC26441 shape ref_shape [tol [all_diff 0/1]] \nif all_diff = 0, only icreasing tolerances is considered" ,
__FILE__,
OCC26441, group);
return;
}

View File

@ -2,7 +2,7 @@ puts "================================================="
puts "0030054: BRepOffset_MakeOffset fails to build joints in intersection mode"
puts "================================================="
puts ""
puts "TODO OCC33166 ALL: Faulty shapes"
restore [locate_data_file bug30054.brep] a
thickshell result a 1 i

View File

@ -0,0 +1,40 @@
puts "============"
puts "0026441: Modeling Algorithms - BRepOffset_MakeOffset affects original shape"
puts "============"
puts ""
pload QAcommands
restore [locate_data_file bug26440_plate.brep] sh
tcopy sh sh_ref
thickshell result sh 160
checkprops result -s 2.40831e+07
checkshape result
set nbshapes_expected "
Number of shapes in shape
VERTEX : 196
EDGE : 308
WIRE : 110
FACE : 110
SHELL : 1
SOLID : 1
COMPSOLID : 0
COMPOUND : 0
SHAPE : 726
"
checknbshapes result -ref ${nbshapes_expected} -t -m "solid construction"
if { [info exist BadEdges] } {
unset BadEdges
}
if { [info exist BadVerts] } {
unset BadVerts
}
OCC26441 sh sh_ref 1.e-7 1
if { [isdraw BadEdges] || [isdraw BadVerts]} {
puts "Error: tolerances of some subshapes of initial shape are changed"
}

View File

@ -1,6 +1,6 @@
puts "TODO OCC23748 ALL: ERROR. offsetperform operation not done."
puts "TODO OCC23748 ALL: Error: The command cannot be built"
puts "TODO OCC26556 ALL: Error : The offset cannot be built."
puts "TODO OCC26556 ALL: Error : The volume"
puts "TODO OCC26556 ALL: Faulty shapes"
puts "TODO OCC26556 ALL: Error : The area"
pcone s 5 0 12 270

View File

@ -1,5 +1,6 @@
puts "TODO OCC25406 ALL: Error : The volume of result shape is"
puts "TODO OCC25406 ALL: Error: bsection of the result and s is not equal to zero"
puts "TODO OCC25406 ALL: Faulty shapes"
ellipse w1 0 0 0 15 10
mkedge w1 w1 0 pi/2

View File

@ -1,5 +1,4 @@
puts "TODO OCC23068 ALL: Error : The volume of result shape "
puts "TODO OCC25406 ALL: Error: bsection of the result and s is not equal to zero"
pcone s 5 0 12 270

View File

@ -1,5 +1,6 @@
puts "TODO OCC23068 ALL: Error : The volume of result shape is"
puts "TODO OCC23068 ALL: Error: bsection of the result and s is not equal to zero"
puts "TODO OCC23068 ALL: Faulty shapes"
ellipse w1 0 0 0 15 10
mkedge w1 w1 0 pi/2