Salome HOME
Copyright update 2022
[modules/smesh.git] / src / StdMeshers / StdMeshers_Import_1D2D.cxx
index e557bc60d6d549114405da755d186a0a3a8061f6..cd2874d4f404213663e9b8516755a20376a6ef53 100644 (file)
@@ -1,4 +1,4 @@
-// Copyright (C) 2007-2014  CEA/DEN, EDF R&D, OPEN CASCADE
+// Copyright (C) 2007-2022  CEA/DEN, EDF R&D, OPEN CASCADE
 //
 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
@@ -20,7 +20,7 @@
 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
 //
 
-//  SMESH SMESH : implementaion of SMESH idl descriptions
+//  SMESH SMESH : implementation of SMESH idl descriptions
 //  File   : StdMeshers_Import_1D2D.cxx
 //  Module : SMESH
 //
 #include "SMESHDS_Group.hxx"
 #include "SMESHDS_Mesh.hxx"
 #include "SMESH_Comment.hxx"
+#include "SMESH_ControlsDef.hxx"
 #include "SMESH_Gen.hxx"
 #include "SMESH_Group.hxx"
 #include "SMESH_Mesh.hxx"
 #include "SMESH_MesherHelper.hxx"
 #include "SMESH_OctreeNode.hxx"
+#include "SMESH_MeshEditor.hxx"
 #include "SMESH_subMesh.hxx"
 
 #include "Utils_SALOME_Exception.hxx"
 #include "utilities.h"
 
+#include <BRepBndLib.hxx>
+#include <BRepClass_FaceClassifier.hxx>
+#include <BRepTools.hxx>
+#include <BRepTopAdaptor_FClass2d.hxx>
 #include <BRep_Builder.hxx>
 #include <BRep_Tool.hxx>
+#include <Bnd_B2d.hxx>
+#include <Bnd_Box.hxx>
+#include <GeomAPI_ProjectPointOnSurf.hxx>
+#include <GeomAdaptor_Surface.hxx>
+#include <Precision.hxx>
 #include <TopExp.hxx>
 #include <TopExp_Explorer.hxx>
 #include <TopoDS.hxx>
 #include <TopoDS_Compound.hxx>
 #include <TopoDS_Edge.hxx>
 #include <TopoDS_Vertex.hxx>
-#include <Precision.hxx>
 
 #include <numeric>
 
@@ -88,15 +98,15 @@ namespace
  */
 //=============================================================================
 
-StdMeshers_Import_1D2D::StdMeshers_Import_1D2D(int hypId, int studyId, SMESH_Gen * gen)
-  :SMESH_2D_Algo(hypId, studyId, gen), _sourceHyp(0)
+StdMeshers_Import_1D2D::StdMeshers_Import_1D2D(int hypId, SMESH_Gen * gen)
+  :SMESH_2D_Algo(hypId, gen), _sourceHyp(0)
 {
-  MESSAGE("StdMeshers_Import_1D2D::StdMeshers_Import_1D2D");
   _name = "Import_1D2D";
   _shapeType = (1 << TopAbs_FACE);
 
   _compatibleHypothesis.push_back("ImportSource2D");
   _requireDiscreteBoundary = false;
+  _supportSubmeshes = true;
 }
 
 //=============================================================================
@@ -180,15 +190,41 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
   SMESHDS_Mesh* tgtMesh = theMesh.GetMeshDS();
 
   const TopoDS_Face& geomFace = TopoDS::Face( theShape );
-  const double faceTol = helper.MaxTolerance( geomFace );
-  const int shapeID = tgtMesh->ShapeToIndex( geomFace );
+  const double  faceTol = helper.MaxTolerance( geomFace );
+  const int     shapeID = tgtMesh->ShapeToIndex( geomFace );
   const bool toCheckOri = (helper.NbAncestors( geomFace, theMesh, TopAbs_SOLID ) == 1 );
 
+
   Handle(Geom_Surface) surface = BRep_Tool::Surface( geomFace );
-  const bool reverse = 
-    ( helper.GetSubShapeOri( tgtMesh->ShapeToMesh(), geomFace) == TopAbs_REVERSED );
+  const bool reverse =
+    ( helper.GetSubShapeOri( tgtMesh->ShapeToMesh(), geomFace ) == TopAbs_REVERSED );
   gp_Pnt p; gp_Vec du, dv;
 
+  // BRepClass_FaceClassifier is most time consuming, so minimize its usage
+  const double clsfTol = 1e2 * BRep_Tool::MaxTolerance( geomFace, TopAbs_VERTEX );
+  BRepTopAdaptor_FClass2d classifier( geomFace, clsfTol ); //Brimless_FaceClassifier classifier;
+  Bnd_B2d bndBox2d;
+  Bnd_Box bndBox3d;
+  {
+    Standard_Real umin,umax,vmin,vmax;
+    BRepTools::UVBounds(geomFace,umin,umax,vmin,vmax);
+    gp_XY pmin( umin,vmin ), pmax( umax,vmax );
+    bndBox2d.Add( pmin );
+    bndBox2d.Add( pmax );
+    if ( helper.HasSeam() )
+    {
+      const int i = helper.GetPeriodicIndex();
+      pmin.SetCoord( i, helper.GetOtherParam( pmin.Coord( i )));
+      pmax.SetCoord( i, helper.GetOtherParam( pmax.Coord( i )));
+      bndBox2d.Add( pmin );
+      bndBox2d.Add( pmax );
+    }
+    bndBox2d.Enlarge( 1e-2 * Sqrt( bndBox2d.SquareExtent() ));
+
+    BRepBndLib::Add( geomFace, bndBox3d );
+    bndBox3d.Enlarge( 1e-2 * sqrt( bndBox3d.SquareExtent() ));
+  }
+
   set<int> subShapeIDs;
   subShapeIDs.insert( shapeID );
 
@@ -212,7 +248,7 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
     existingNodes.insert( n );
   }
 
-  // get EDGESs and their ids and get existing nodes on EDGEs
+  // get EDGEs and their ids and get existing nodes on EDGEs
   vector< TopoDS_Edge > edges;
   for ( exp.Init( theShape, TopAbs_EDGE ); exp.More(); exp.Next() )
   {
@@ -237,13 +273,29 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
   map<TLink, int>::iterator link2Nb;
   double minGroupTol = Precision::Infinite();
 
+  SMESH::Controls::ElementsOnShape onEdgeClassifier;
+  {
+    TopoDS_Compound edgesCompound;
+    BRep_Builder    builder;
+    builder.MakeCompound( edgesCompound );
+    for ( size_t iE = 0; iE < edges.size(); ++iE )
+      builder.Add( edgesCompound, edges[ iE ]);
+    onEdgeClassifier.SetShape( edgesCompound, SMDSAbs_Node );
+  }
+
   // =========================
   // Import faces from groups
   // =========================
 
   StdMeshers_Import_1D::TNodeNodeMap* n2n;
   StdMeshers_Import_1D::TElemElemMap* e2e;
-  vector<const SMDS_MeshNode*> newNodes;
+  StdMeshers_Import_1D::TNodeNodeMap::iterator n2nIt;
+  pair< StdMeshers_Import_1D::TNodeNodeMap::iterator, bool > it_isnew;
+  vector<TopAbs_State>         nodeState;
+  vector<const SMDS_MeshNode*> newNodes; // of a face
+  set   <const SMDS_MeshNode*> bndNodes; // nodes classified ON
+  vector<bool>                 isNodeIn; // nodes classified IN, by node ID
+
   for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
   {
     const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
@@ -256,58 +308,157 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
     const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
     minGroupTol = std::min( groupTol, minGroupTol );
 
+    // clsfTol is 2D tolerance of a probe line
+    //GeomAdaptor_Surface S( surface );
+    // const double clsfTol = Min( S.UResolution( 0.1 * groupTol ), -- issue 0023092
+    //                             S.VResolution( 0.1 * groupTol ));
+    // another idea: try to use max tol of all edges
+    //const double clsfTol = 10 * BRep_Tool::Tolerance( geomFace ); // 0.1 * groupTol;
+
+    onEdgeClassifier.SetMesh( srcMesh->GetMeshDS() );
+    onEdgeClassifier.SetTolerance( groupTol / 10 );
+
+
     SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
-    SMDS_MeshNode *tmpNode = helper.AddNode(0,0,0);
-    gp_XY uv( Precision::Infinite(), Precision::Infinite() );
     while ( srcElems->more() ) // loop on group contents
     {
       const SMDS_MeshElement* face = srcElems->next();
+
+      if ( bndBox3d.IsOut( SMESH_NodeXYZ( face->GetNode(0) )))
+        continue;
+
       // find or create nodes of a new face
-      newNodes.resize( face->NbNodes() );
+      nodeState.resize( face->NbNodes() );
+      newNodes.resize( nodeState.size() );
       newNodes.back() = 0;
       int nbCreatedNodes = 0;
-      SMDS_MeshElement::iterator node = face->begin_nodes();
-      for ( size_t i = 0; i < newNodes.size(); ++i, ++node )
+      bool isOut = false, isIn = false; // if at least one node isIn - do not classify other nodes
+      for ( size_t i = 0; i < newNodes.size(); ++i )
       {
-        StdMeshers_Import_1D::TNodeNodeMap::iterator n2nIt =
-          n2n->insert( make_pair( *node, (SMDS_MeshNode*)0 )).first;
-        if ( n2nIt->second )
+        const SMDS_MeshNode* node = face->GetNode( i );
+        SMESH_NodeXYZ nXYZ = node;
+        nodeState[ i ] = TopAbs_UNKNOWN;
+        newNodes [ i ] = 0;
+
+        it_isnew = n2n->insert( make_pair( node, nullptr ));
+        n2nIt    = it_isnew.first;
+
+        const SMDS_MeshNode* & newNode = n2nIt->second;
+        if ( !it_isnew.second && !newNode )
+          break; // a node is mapped to NULL - it is OUT of the FACE
+
+        if ( newNode )
         {
-          if ( !subShapeIDs.count( n2nIt->second->getshapeId() )) // node already on an EDGE
-            break;
+          if ( !subShapeIDs.count( newNode->getshapeId() ))
+            break; // node is Imported onto other FACE
+          if ( newNode->GetID() < (int) isNodeIn.size() &&
+               isNodeIn[ newNode->GetID() ])
+            isIn = true;
+          if ( !isIn && bndNodes.count( node ))
+            nodeState[ i ] = TopAbs_ON;
         }
         else
         {
-          // find a pre-existing node
+          // find a node pre-existing on EDGE or VERTEX
           dist2foundNodes.clear();
-          existingNodeOcTr.NodesAround( SMESH_TNodeXYZ( *node ), dist2foundNodes, groupTol );
+          existingNodeOcTr.NodesAround( nXYZ, dist2foundNodes, groupTol );
           if ( !dist2foundNodes.empty() )
-            (*n2nIt).second = dist2foundNodes.begin()->second;
+          {
+            newNode = dist2foundNodes.begin()->second;
+            nodeState[ i ] = TopAbs_ON;
+          }
         }
-        if ( !n2nIt->second )
+
+        if ( !newNode )
         {
-          // find out if node lies on theShape
-          tmpNode->setXYZ( (*node)->X(), (*node)->Y(), (*node)->Z());
-          uv.SetCoord( Precision::Infinite(), Precision::Infinite() );
-          if ( helper.CheckNodeUV( geomFace, tmpNode, uv, groupTol, /*force=*/true ))
+          // find out if node lies on the surface of theShape
+          gp_XY uv( Precision::Infinite(), 0 );
+          bool isOutBox = true;
+          isOut = (! helper.CheckNodeUV( geomFace, node, uv, groupTol, /*force=*/true ) ||
+                   ( isOutBox = bndBox2d.IsOut( uv )));
+          //int iCoo;
+          if ( !isOut && !isIn ) // classify
           {
-            SMDS_MeshNode* newNode = tgtMesh->AddNode( (*node)->X(), (*node)->Y(), (*node)->Z());
-            n2nIt->second = newNode;
+            nodeState[i] = classifier.Perform( uv ); //classifier.Perform( geomFace, uv, clsfTol );
+            //nodeState[i] = classifier.State();
+            isOut = ( nodeState[i] == TopAbs_OUT );
+            if (( isOut ) &&
+                ( !isOutBox || helper.IsOnSeam( uv )) &&
+                onEdgeClassifier.IsSatisfy( node ))
+            {
+              // uv.SetCoord( iCoo, helper.GetOtherParam( uv.Coord( iCoo )));
+              // classifier.Perform( geomFace, uv, clsfTol );
+              // nodeState[i] = classifier.State();
+              // isOut = ( nodeState[i] == TopAbs_OUT );
+              nodeState[i] = TopAbs_ON;
+              isOut = false;
+            }
+          }
+          if ( !isOut ) // create a new node
+          {
+            newNode = tgtMesh->AddNode( nXYZ.X(), nXYZ.Y(), nXYZ.Z());
             tgtMesh->SetNodeOnFace( newNode, shapeID, uv.X(), uv.Y() );
             nbCreatedNodes++;
+            if ( newNode->GetID() >= (smIdType) isNodeIn.size() )
+            {
+              isNodeIn.push_back( false ); // allow allocate more than newNode->GetID()
+              isNodeIn.resize( newNode->GetID() + 1, false );
+            }
+            if ( nodeState[i] == TopAbs_ON )
+              bndNodes.insert( node );
+            else if ( nodeState[i] != TopAbs_UNKNOWN )
+              isNodeIn[ newNode->GetID() ] = isIn = true;
           }
         }
-        if ( !(newNodes[i] = n2nIt->second ))
+        if ( !(newNodes[i] = newNode ) || isOut )
           break;
-      }
+
+      } // loop on face nodes
+
       if ( !newNodes.back() )
         continue; // not all nodes of the face lie on theShape
 
+      if ( !isIn ) // if all nodes are on FACE boundary, a mesh face can be OUT
+      {
+        // check state of nodes created for other faces
+        for ( size_t i = 0; i < nodeState.size() && !isIn; ++i )
+        {
+          if ( nodeState[i] != TopAbs_UNKNOWN ) continue;
+          gp_XY uv = helper.GetNodeUV( geomFace, newNodes[i] );
+          nodeState[i] = classifier.Perform( uv ); //geomFace, uv, clsfTol );
+          //nodeState[i] = classifier.State();
+          isIn = ( nodeState[i] == TopAbs_IN );
+        }
+        if ( !isIn ) // classify face center
+        {
+          gp_XYZ gc( 0., 0., 0 );
+          for ( size_t i = 0; i < newNodes.size(); ++i )
+            gc += SMESH_TNodeXYZ( newNodes[i] );
+          gc /= newNodes.size();
+
+          TopLoc_Location loc;
+          GeomAPI_ProjectPointOnSurf& proj = helper.GetProjector( geomFace,
+                                                                  loc,
+                                                                  helper.MaxTolerance( geomFace ));
+          if ( !loc.IsIdentity() ) loc.Transformation().Inverted().Transforms( gc );
+          proj.Perform( gc );
+          if ( !proj.IsDone() || proj.NbPoints() < 1 )
+            continue;
+          Standard_Real U,V;
+          proj.LowerDistanceParameters(U,V);
+          gp_XY uv( U,V );
+          //classifier.Perform( geomFace, uv, clsfTol );
+          TopAbs_State state = classifier.Perform( uv );
+          if ( state != TopAbs_IN )
+            continue;
+        }
+      }
+
       // try to find already created face
       SMDS_MeshElement * newFace = 0;
       if ( nbCreatedNodes == 0 &&
            tgtMesh->FindElement(newNodes, SMDSAbs_Face, /*noMedium=*/false))
-        continue; // repeated face in source groups already created 
+        continue; // repeated face in source groups already created
 
       // check future face orientation
       const int nbCorners = face->NbCornerNodes();
@@ -318,7 +469,7 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
         gp_Vec geomNorm;
         do
         {
-          uv = helper.GetNodeUV( geomFace, newNodes[++iNode] );
+          gp_XY uv = helper.GetNodeUV( geomFace, newNodes[++iNode] );
           surface->D1( uv.X(),uv.Y(), p, du,dv );
           geomNorm = reverse ? dv^du : du^dv;
         }
@@ -336,7 +487,7 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
 
         if ( geomNorm * meshNorm < 0 )
           SMDS_MeshCell::applyInterlace
-            ( SMDS_MeshCell::reverseSmdsOrder( face->GetEntityType() ), newNodes );
+            ( SMDS_MeshCell::reverseSmdsOrder( face->GetEntityType(), newNodes.size() ), newNodes );
       }
 
       // make a new face
@@ -374,16 +525,31 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
           medium = newNodes[i+nbCorners];
         link2Nb = linkCount.insert( make_pair( TLink( n1, n2, medium ), 0)).first;
         ++link2Nb->second;
-        // if ( link2Nb->second == 1 )
-        // {
-        //   // measure link length
-        //   double len2 = SMESH_TNodeXYZ( n1 ).SquareDistance( n2 );
-        //   if ( len2 < minGroupTol )
-        //     minGroupTol = len2;
-        // }
+      }
+    } // loop on group contents
+
+    // Remove OUT nodes from n2n map
+    for ( n2nIt = n2n->begin(); n2nIt != n2n->end(); )
+      if ( !n2nIt->second )
+        n2n->erase( n2nIt++ );
+      else
+        ++n2nIt;
+
+  } // loop on src groups
+
+  // remove free nodes created on EDGEs
+  {
+    set<const SMDS_MeshNode*>::iterator node = bndNodes.begin();
+    for ( ; node != bndNodes.end(); ++node )
+    {
+      n2nIt = n2n->find( *node );
+      const SMDS_MeshNode* newNode = n2nIt->second;
+      if ( newNode && newNode->NbInverseElements() == 0 )
+      {
+        tgtMesh->RemoveFreeNode( newNode, 0, false );
+        n2n->erase( n2nIt );
       }
     }
-    helper.GetMeshDS()->RemoveNode(tmpNode);
   }
 
   // ==========================================================
@@ -402,7 +568,6 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
     // the imported mesh is valid if all external links (encountered once)
     // lie on geom edges
     subShapeIDs.erase( shapeID ); // to contain edges and vertices only
-    double u, f, l;
     for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); ++link2Nb)
     {
       const TLink& link = (*link2Nb).first;
@@ -417,17 +582,14 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
           const SMDS_MeshNode* n = is1stN ? link.node1() : link.node2();
           if ( !subShapeIDs.count( n->getshapeId() )) // n is assigned to FACE
           {
-            for ( size_t iE = 0; iE < edges.size(); ++iE )
-              if ( helper.CheckNodeU( edges[iE], n, u=0, projTol, /*force=*/true ))
-              {
-                BRep_Tool::Range(edges[iE],f,l);
-                if ( Abs(u-f) < 2 * faceTol || Abs(u-l) < 2 * faceTol )
-                  // duplicated node on vertex
-                  return error("Source elements overlap one another");
-                tgtFaceSM->RemoveNode( n, /*isNodeDeleted=*/false );
-                tgtMesh->SetNodeOnEdge( n, edges[iE], u );
-                break;
-              }
+            TopoDS_Shape edge;
+            if ( onEdgeClassifier.IsSatisfy( n, &edge ))
+            {
+              tgtFaceSM->RemoveNode( n );
+              double u, v;
+              onEdgeClassifier.GetParams( u, v );
+              tgtMesh->SetNodeOnEdge( n, TopoDS::Edge(edge), u );
+            }
             nodesOnBoundary = subShapeIDs.count( n->getshapeId());
           }
           if ( nodesOnBoundary )
@@ -475,9 +637,10 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
         {
           edge = tgtMesh->AddEdge( newNodes[0], newNodes[1], link._medium );
 
+          double u;
           TopoDS_Edge geomEdge = TopoDS::Edge(bndShapes.back());
           helper.CheckNodeU( geomEdge, link._medium, u, projTol, /*force=*/true );
-          tgtFaceSM->RemoveNode( link._medium, /*isNodeDeleted=*/false );
+          tgtFaceSM->RemoveNode( link._medium );
           tgtMesh->SetNodeOnEdge( (SMDS_MeshNode*)link._medium, geomEdge, u );
         }
         else
@@ -547,6 +710,7 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
       TUNodeList nodesOnSeam;
       double u = helper.GetNodeU( seamEdge, vertNode );
       nodesOnSeam.push_back( make_pair( u, vertNode ));
+      size_t nbNodesOnSeam = 1;
       TUNodeList::iterator u2nIt = nodesOnSeam.begin();
       for ( ; u2nIt != nodesOnSeam.end(); ++u2nIt )
       {
@@ -560,8 +724,9 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
           {
             const SMDS_MeshNode* n = face->GetNode( i );
             if ( n == startNode || !checkedNodes.insert( n ).second ) continue;
-            if ( helper.CheckNodeU( seamEdge, n, u=0, projTol, /*force=*/true ))
-              nodesOnSeam.push_back( make_pair( u, n ));
+            helper.CheckNodeU( seamEdge, n, u=0, projTol, /*force=*/true );
+            nodesOnSeam.push_back( make_pair( u, n ));
+            ++nbNodesOnSeam;
           }
         }
       }
@@ -570,6 +735,15 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
       map< double, const SMDS_MeshNode* > u2nodeMap;
       for ( u2nIt = nodesOnSeam.begin(); u2nIt != nodesOnSeam.end(); ++u2nIt )
         u2nodeMap.insert( u2nodeMap.end(), *u2nIt );
+      if ( u2nodeMap.size() != nbNodesOnSeam ) // problem with parameters on EDGE
+      {
+        // sort nodes by distance from seamVertex
+        gp_Pnt vertPnt = SMESH_NodeXYZ( vertNode );
+        u2nodeMap.clear();
+        for ( u2nIt = nodesOnSeam.begin(); u2nIt != nodesOnSeam.end(); ++u2nIt )
+          u2nodeMap.insert
+            ({ vertPnt.SquareDistance( SMESH_NodeXYZ( u2nIt->second )), u2nIt->second });
+      }
 
       // create edges
       {
@@ -577,7 +751,7 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
         seamHelper.SetSubShape( edges[ iE ]);
         seamHelper.SetElementsOnShape( true );
 
-        if ( (*checkedFaces.begin())->IsQuadratic() )
+        if ( !checkedFaces.empty() && (*checkedFaces.begin())->IsQuadratic() )
           for ( set< const SMDS_MeshElement* >::iterator fIt = checkedFaces.begin();
                 fIt != checkedFaces.end(); ++fIt )
             seamHelper.AddTLinks( static_cast<const SMDS_MeshFace*>( *fIt ));
@@ -590,7 +764,7 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
           seamHelper.AddEdge( node1, node2 );
           if ( node2->getshapeId() == helper.GetSubShapeID() )
           {
-            tgtFaceSM->RemoveNode( node2, /*isNodeDeleted=*/false );
+            tgtFaceSM->RemoveNode( node2 );
             tgtMesh->SetNodeOnEdge( const_cast<SMDS_MeshNode*>( node2 ), seamEdge, n2->first );
           }
         }
@@ -606,8 +780,7 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
     //   sm->SetIsAlwaysComputed( true );
     sm->ComputeStateEngine(SMESH_subMesh::CHECK_COMPUTE_STATE);
     if ( sm->GetComputeState() != SMESH_subMesh::COMPUTE_OK )
-      return error(SMESH_Comment("Failed to create segments on the edge ")
-                   << tgtMesh->ShapeToIndex( edges[iE ]));
+      return error(SMESH_Comment("Failed to create segments on the edge #") << sm->GetId());
   }
 
   // ============
@@ -661,7 +834,7 @@ bool StdMeshers_Import_1D2D::Evaluate(SMESH_Mesh &         theMesh,
   if ( srcGroups.empty() )
     return error("Invalid source groups");
 
-  vector<int> aVec(SMDSEntity_Last,0);
+  vector<smIdType> aVec(SMDSEntity_Last,0);
 
   bool toCopyMesh, toCopyGroups;
   _sourceHyp->GetCopySourceMesh(toCopyMesh, toCopyGroups);
@@ -701,7 +874,7 @@ bool StdMeshers_Import_1D2D::Evaluate(SMESH_Mesh &         theMesh,
     set<const SMDS_MeshNode* > allNodes;
     gp_XY uv;
     double minGroupTol = 1e100;
-    for ( int iG = 0; iG < srcGroups.size(); ++iG )
+    for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
     {
       const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
       const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
@@ -751,7 +924,7 @@ bool StdMeshers_Import_1D2D::Evaluate(SMESH_Mesh &         theMesh,
     {
       TopoDS_Edge geomEdge = TopoDS::Edge( exp.Current() );
       SMESH_subMesh* sm = theMesh.GetSubMesh( geomEdge );
-      vector<int>& edgeVec = aResMap[sm];
+      vector<smIdType>& edgeVec = aResMap[sm];
       if ( edgeVec.empty() )
       {
         edgeVec.resize(SMDSEntity_Last,0);