Salome HOME
IPAL54585: Extrusion 3D algo fails with "OCC exception. Standard_NoSuchObject: NColle...
authoreap <eap@opencascade.com>
Mon, 1 Jul 2019 10:44:54 +0000 (13:44 +0300)
committereap <eap@opencascade.com>
Mon, 1 Jul 2019 10:44:54 +0000 (13:44 +0300)
src/SMESH/SMESH_MesherHelper.cxx
src/SMESHUtils/SMESH_Delaunay.cxx
src/SMESHUtils/SMESH_Slot.cxx
src/StdMeshers/StdMeshers_Quadrangle_2D.cxx
src/StdMeshers/StdMeshers_RadialQuadrangle_1D2D.cxx

index d6cd15d3ef8b4402873367a33b9d233304055947..c30234a7e48169024cd38a090184679104d40cd5 100644 (file)
@@ -151,8 +151,6 @@ bool SMESH_MesherHelper::IsQuadraticSubMesh(const TopoDS_Shape& aSh)
   // we can create quadratic elements only if all elements
   // created on sub-shapes of given shape are quadratic
   myCreateQuadratic = true;
-  mySeamShapeIds.clear();
-  myDegenShapeIds.clear();
   TopAbs_ShapeEnum subType( aSh.ShapeType()==TopAbs_FACE ? TopAbs_EDGE : TopAbs_FACE );
   if ( aSh.ShapeType()==TopAbs_COMPOUND )
   {
@@ -741,7 +739,8 @@ gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face&   F,
               if ( !C2d.IsNull() ) {
                 double u = ( V == IthVertex( 0, edge )) ?  f : l;
                 uv = C2d->Value( u );
-                uvOK = true;
+                gp_Pnt p = GetSurface( F )->Value( uv );
+                uvOK = ( p.Distance( BRep_Tool::Pnt( V )) < getFaceMaxTol( F ));
                 break;
               }
             }
index 8ae5f1761a41089737828165354d6c0d61d8311c..8b1f252d6ab05d011d9fa48c569d3e70f87d2d15 100644 (file)
@@ -40,7 +40,6 @@
  *  \param [in] boundaryNodes - vector of nodes of a wire
  *  \param [in] face - the face
  *  \param [in] faceID - the face ID
- *  \param [in] nbNodesToVisit - nb of non-marked nodes on the face
  */
 //================================================================================
 
@@ -56,6 +55,8 @@ SMESH_Delaunay::SMESH_Delaunay(const std::vector< const UVPtStructVec* > & bound
     const int nbDiv = 100;
     const double uRange = surf.LastUParameter() - surf.FirstUParameter();
     const double vRange = surf.LastVParameter() - surf.FirstVParameter();
+    const double uFixed = surf.FirstUParameter() + 0.5 * uRange;
+    const double vFixed = surf.FirstVParameter() + 0.5 * vRange;
     const double dU = uRange / nbDiv;
     const double dV = vRange / nbDiv;
     double u = surf.FirstUParameter(), v = surf.FirstVParameter();
@@ -63,10 +64,10 @@ SMESH_Delaunay::SMESH_Delaunay(const std::vector< const UVPtStructVec* > & bound
     double lenU = 0, lenV = 0;
     for ( ; u < surf.LastUParameter(); u += dU, v += dV )
     {
-      gp_Pnt p1U = surf.Value( u, surf.FirstVParameter() );
+      gp_Pnt p1U = surf.Value( u, vFixed );
       lenU += p1U.Distance( p0U );
       p0U = p1U;
-      gp_Pnt p1V = surf.Value( surf.FirstUParameter(), v );
+      gp_Pnt p1V = surf.Value( uFixed, v );
       lenV += p1V.Distance( p0V );
       p0V = p1V;
     }
@@ -270,6 +271,8 @@ const BRepMesh_Triangle* SMESH_Delaunay::FindTriangle( const gp_XY&
 
 const BRepMesh_Triangle* SMESH_Delaunay::GetTriangleNear( int iBndNode )
 {
+  if ( iBndNode >= _triaDS->NbNodes() )
+    return 0;
   int nodeIDs[3];
   int nbNbNodes = _bndNodes.size();
 #if OCC_VERSION_LARGE <= 0x07030000
index f706f5bceeddf28bea7b575d2cdd3bab0c773f9a..a11555ddd09c26d5b5ee4690c375de0e599cced5 100644 (file)
@@ -127,6 +127,15 @@ namespace
         return 0;
       if ( myFreeEnds.empty() )
       {
+        // remove degenerated cuts
+        // for ( size_t iC1 = 0; iC1 < myCuts.size(); ++iC1 )
+        //   if ( myCuts[ iC1 ][ 0 ].myNode == myCuts[ iC1 ][ 1 ].myNode )
+        //   {
+        //     if ( iC1 < myCuts.size() - 1 )
+        //       myCuts[ iC1 ] = myCuts.back();
+        //     myCuts.pop_back();
+        //   }
+
         int nbShared = 0;
         std::vector< bool > isSharedPnt( myCuts.size() * 2, false );
         for ( size_t iC1 = 0; iC1 < myCuts.size() - 1; ++iC1 )
index 17702dd30d6a0e1f81589de66ac1692a52734594..e416673067eada1a400b8eec05f37ae8744be35b 100644 (file)
@@ -4336,6 +4336,7 @@ void StdMeshers_Quadrangle_2D::smooth (FaceQuadStruct::Ptr quad)
       SMESH_TNodeXYZ a2( quad->UVPt( nbhoriz-1, nbvertic-1 ).node );
       SMESH_TNodeXYZ a3( quad->UVPt( 0,         nbvertic-1 ).node );
 
+      // compute TFI
       for (int i = 1; i < nbhoriz-1; i++)
       {
         SMESH_TNodeXYZ p0( quad->UVPt( i, 0          ).node );
@@ -4347,13 +4348,39 @@ void StdMeshers_Quadrangle_2D::smooth (FaceQuadStruct::Ptr quad)
 
           UVPtStruct& uvp = quad->UVPt( i, j );
 
-          gp_Pnt    p = myHelper->calcTFI(uvp.x,uvp.y, a0,a1,a2,a3, p0,p1,p2,p3);
+          gp_Pnt pnew = myHelper->calcTFI(uvp.x,uvp.y, a0,a1,a2,a3, p0,p1,p2,p3);
+          meshDS->MoveNode( uvp.node, pnew.X(), pnew.Y(), pnew.Z() );
+        }
+      }
+      // project to surface
+      double cellSize;
+      for (int i = 1; i < nbhoriz-1; i++)
+      {
+        for (int j = 1; j < nbvertic-1; j++)
+        {
+          UVPtStruct& uvp = quad->UVPt( i, j );
+          SMESH_NodeXYZ p = uvp.node;
+
+          cellSize = Max( p.SquareDistance( quad->UVPt( i+1, j ).node ),
+                          p.SquareDistance( quad->UVPt( i-1, j ).node ));
+          cellSize = Max( p.SquareDistance( quad->UVPt( i, j+1 ).node ), cellSize );
+          cellSize = Max( p.SquareDistance( quad->UVPt( i, j-1 ).node ), cellSize );
+
           gp_Pnt2d uv = surface->NextValueOfUV( uvp.UV(), p, 10*tol );
           gp_Pnt pnew = surface->Value( uv );
-
-          meshDS->MoveNode( uvp.node, pnew.X(), pnew.Y(), pnew.Z() );
-          uvp.u = uv.X();
-          uvp.v = uv.Y();
+          bool     ok = ( pnew.SquareDistance( p ) < 2 * cellSize );
+          if ( !ok )
+          {
+            uv   = surface->ValueOfUV( p, 10*tol );
+            pnew = surface->Value( uv );
+            ok   = ( pnew.SquareDistance( p ) < 2 * cellSize );
+          }
+          if ( ok )
+          {
+            meshDS->MoveNode( uvp.node, pnew.X(), pnew.Y(), pnew.Z() );
+            uvp.u = uv.X();
+            uvp.v = uv.Y();
+          }
         }
       }
     }
index 8e183e37ea621dd0ff509d988cb62e1203a8f507..740eebaeb0279b0b98428d54417830ce71b94968 100644 (file)
@@ -226,12 +226,19 @@ namespace
     if ( nbWire > 2 || nbEdgesInWire.front() < 1 ) return 0;
 
     // remove degenerated EDGEs
+    TopTools_MapOfShape degenVV;
     list<TopoDS_Edge>::iterator edge = edges.begin();
     while ( edge != edges.end() )
       if ( SMESH_Algo::isDegenerated( *edge ))
+      {
+        degenVV.Add( SMESH_MesherHelper::IthVertex( 0, *edge ));
+        degenVV.Add( SMESH_MesherHelper::IthVertex( 1, *edge ));
         edge = edges.erase( edge );
+      }
       else
+      {
         ++edge;
+      }
     int nbEdges = edges.size();
 
     // find VERTEXes between continues EDGEs
@@ -327,9 +334,25 @@ namespace
       double len1 = (++l2i)->first;
       double len2 = (++l2i)->first;
       if ( len1 - len0 > len2 - len1 )
-        deviation2sideInd.insert( make_pair( 0., len2sideInd.begin()->second ));
+        deviation2sideInd.insert( std::make_pair( 0., len2sideInd.begin()->second ));
       else
-        deviation2sideInd.insert( make_pair( 0., len2sideInd.rbegin()->second ));
+        deviation2sideInd.insert( std::make_pair( 0., len2sideInd.rbegin()->second ));
+    }
+
+    double minDevi = deviation2sideInd.begin()->first;
+    int   iMinCurv = deviation2sideInd.begin()->second;
+    if ( sides.size() == 3 && degenVV.Size() == 1 &&
+         minDevi / sides[ iMinCurv ]->Length() > 1e-3 )
+    {
+      // a triangle with curved sides and a degenerated EDGE (IPAL54585);
+      // use a side opposite to the degenerated EDGE as an elliptic one
+      for ( size_t iS = 0; iS < sides.size(); ++iS )
+        if ( degenVV.Contains( sides[ iS ]->FirstVertex() ))
+        {
+          deviation2sideInd.clear();
+          deviation2sideInd.insert( std::make_pair( 0.,( iS + 1 ) % sides.size() ));
+          break;
+        }
     }
 
     int iCirc = deviation2sideInd.rbegin()->second; 
@@ -994,6 +1017,7 @@ bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh&         aMesh,
   quad->side[1] = linSide1;
   quad->side[2] = StdMeshers_FaceSide::New( circSide.get(), centerNode, &centerUV );
   quad->side[3] = linSide2;
+  quad->face    = F;
 
   myQuadList.push_back( quad );
 
@@ -1004,6 +1028,12 @@ bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh&         aMesh,
   else
     ok = StdMeshers_Quadrangle_2D::computeTriangles( aMesh, F, quad );
 
+  if ( helper.HasDegeneratedEdges() )
+  {
+    StdMeshers_Quadrangle_2D::myNeedSmooth = true;
+    StdMeshers_Quadrangle_2D::smooth( quad );
+  }
+
   StdMeshers_Quadrangle_2D::myHelper = 0;
 
   return ok;