Salome HOME
IPAL52499: Prismatic mesh is not computed on a prismatic shape
[modules/smesh.git] / src / StdMeshers / StdMeshers_FaceSide.cxx
index d1d20245d4858158eaccc04221b357aefba97717..1424e4226c40d4b92ffce8f8355cbbc73e6f79c1 100644 (file)
@@ -1,4 +1,4 @@
-// Copyright (C) 2007-2015  CEA/DEN, EDF R&D, OPEN CASCADE
+// Copyright (C) 2007-2016  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
 
 #include "utilities.h"
 
+using namespace std;
+
 //================================================================================
 /*!
  * \brief Constructor of a side of one edge
   * \param theFace - the face
   * \param theEdge - the edge
- */
 */
 //================================================================================
 
 StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face&   theFace,
@@ -69,11 +71,15 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face&   theFace,
                                          SMESH_Mesh*          theMesh,
                                          const bool           theIsForward,
                                          const bool           theIgnoreMediumNodes,
+                                         SMESH_MesherHelper*  theFaceHelper,
                                          SMESH_ProxyMesh::Ptr theProxyMesh)
 {
-  list<TopoDS_Edge> edges(1,theEdge);
-  *this = StdMeshers_FaceSide( theFace, edges, theMesh, theIsForward,
-                               theIgnoreMediumNodes, theProxyMesh );
+  std::list<TopoDS_Edge> edges(1,theEdge);
+  StdMeshers_FaceSide tmp( theFace, edges, theMesh, theIsForward,
+                           theIgnoreMediumNodes, theFaceHelper, theProxyMesh );
+  *this = tmp;
+
+  tmp.myHelper = NULL;
 }
 
 //================================================================================
@@ -82,12 +88,13 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face&   theFace,
  */
 //================================================================================
 
-StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face&   theFace,
-                                         list<TopoDS_Edge>&   theEdges,
-                                         SMESH_Mesh*          theMesh,
-                                         const bool           theIsForward,
-                                         const bool           theIgnoreMediumNodes,
-                                         SMESH_ProxyMesh::Ptr theProxyMesh)
+StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face&            theFace,
+                                         const std::list<TopoDS_Edge>& theEdges,
+                                         SMESH_Mesh*                   theMesh,
+                                         const bool                    theIsForward,
+                                         const bool                    theIgnoreMediumNodes,
+                                         SMESH_MesherHelper*           theFaceHelper,
+                                         SMESH_ProxyMesh::Ptr          theProxyMesh)
 {
   int nbEdges = theEdges.size();
   myEdge.resize      ( nbEdges );
@@ -106,13 +113,19 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face&   theFace,
   myMissingVertexNodes = false;
   myIgnoreMediumNodes  = theIgnoreMediumNodes;
   myDefaultPnt2d       = gp_Pnt2d( 1e+100, 1e+100 );
+  myHelper             = NULL;
   if ( !myProxyMesh ) myProxyMesh.reset( new SMESH_ProxyMesh( *theMesh ));
+  if ( theFaceHelper && theFaceHelper->GetSubShape() == myFace )
+  {
+    myHelper           = new SMESH_MesherHelper( * myProxyMesh->GetMesh() );
+    myHelper->CopySubShapeInfo( *theFaceHelper );
+  }
   if ( nbEdges == 0 ) return;
 
   SMESHDS_Mesh* meshDS = myProxyMesh->GetMeshDS();
 
   int nbDegen = 0;
-  list<TopoDS_Edge>::iterator edge = theEdges.begin();
+  std::list<TopoDS_Edge>::const_iterator edge = theEdges.begin();
   for ( int index = 0; edge != theEdges.end(); ++index, ++edge )
   {
     int i = theIsForward ? index : nbEdges-index-1;
@@ -142,7 +155,6 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face&   theFace,
         double p4 = myFirst[i]+(myLast[i]-myFirst[i])/4.;
         double d2 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p2 );
         double d4 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p4 );
-        //cout<<"len = "<<len<<"  d2 = "<<d2<<"  fabs(2*d2/len-1.0) = "<<fabs(2*d2/len-1.0)<<endl;
         myIsUniform[i] = !( fabs(2*d2/myEdgeLength[i]-1.0) > 0.01 || fabs(2*d4/d2-1.0) > 0.01 );
         Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],d2,d4);
         myC3dAdaptor[i].Load( C3d, d2,d4 );
@@ -194,12 +206,15 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const StdMeshers_FaceSide*  theSide,
                                          const double                theUFirst,
                                          const double                theULast)
 {
+  myEdge.resize        ( 1 );
+  myEdgeID.resize      ( 1, 0 );
   myC2d.push_back      ( theC2d );
+  myC3dAdaptor.resize  ( 1 );
   myFirst.push_back    ( theUFirst );
   myLast.push_back     ( theULast );
   myNormPar.push_back  ( 1. );
   myIsUniform.push_back( true );
-  myEdgeID.push_back   ( 0 );
+  myHelper       = NULL;
   myLength       = 0;
   myProxyMesh    = theSide->myProxyMesh;
   myDefaultPnt2d = *thePnt2d1;
@@ -230,9 +245,11 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(const StdMeshers_FaceSide*  theSide,
 //================================================================================
 
 StdMeshers_FaceSide::StdMeshers_FaceSide(UVPtStructVec&     theSideNodes,
-                                         const TopoDS_Face& theFace)
+                                         const TopoDS_Face& theFace,
+                                         const TopoDS_Edge& theEdge,
+                                         SMESH_Mesh*        theMesh)
 {
-  myEdge.resize( 1 );
+  myEdge.resize( 1, theEdge );
   myEdgeID.resize( 1, -1 );
   myC2d.resize( 1 );
   myC3dAdaptor.resize( 1 );
@@ -242,8 +259,20 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(UVPtStructVec&     theSideNodes,
   myIsUniform.resize( 1, 1 );
   myMissingVertexNodes = myIgnoreMediumNodes = false;
   myDefaultPnt2d.SetCoord( 1e100, 1e100 );
+  if ( theMesh ) myProxyMesh.reset( new SMESH_ProxyMesh( *theMesh ));
+  if ( !theEdge.IsNull() )
+  {
+    if ( theMesh ) myEdgeID[0] = theMesh->GetMeshDS()->ShapeToIndex( theEdge );
+    if ( theFace.IsNull() )
+      BRep_Tool::Range( theEdge, myFirst[0], myLast[0] );
+    else
+      myC2d[0] = BRep_Tool::CurveOnSurface( theEdge, theFace, myFirst[0], myLast[0] );
+    if ( theEdge.Orientation() == TopAbs_REVERSED )
+      std::swap( myFirst[0], myLast[0] );
+  }
 
   myFace       = theFace;
+  myHelper     = NULL;
   myPoints     = theSideNodes;
   myNbPonits   = myPoints.size();
   myNbSegments = myNbPonits + 1;
@@ -296,28 +325,38 @@ StdMeshers_FaceSide::StdMeshers_FaceSide(UVPtStructVec&     theSideNodes,
   myEdgeLength.resize( 1, myLength );
 }
 
+//================================================================================
+/*!
+ * \brief Destructor
+ */
+//================================================================================
+
+StdMeshers_FaceSide::~StdMeshers_FaceSide()
+{
+  delete myHelper; myHelper = NULL;
+}
+
 //================================================================================
 /*
  * Return info on nodes on the side
  */
 //================================================================================
 
-const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool   isXConst,
-                                                             double constValue) const
+const std::vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool   isXConst,
+                                                                  double constValue) const
 {
   if ( myPoints.empty() )
   {
     if ( NbEdges() == 0 ) return myPoints;
 
     StdMeshers_FaceSide* me = const_cast< StdMeshers_FaceSide* >( this );
-    SMESHDS_Mesh*    meshDS = myProxyMesh->GetMeshDS();
-    SMESH_MesherHelper eHelper( *myProxyMesh->GetMesh() );
-    SMESH_MesherHelper fHelper( *myProxyMesh->GetMesh() );
-    fHelper.SetSubShape( myFace );
-    bool paramOK;
+    bool paramOK = true;
     double eps = 1e-100;
 
-    // sort nodes of all edges putting them into a map
+    SMESH_MesherHelper  eHelper( *myProxyMesh->GetMesh() );
+    SMESH_MesherHelper& fHelper = *FaceHelper();
+
+    // sort nodes of all edges by putting them into a map
 
     map< double, const SMDS_MeshNode*>            u2node;
     vector< pair< double, const SMDS_MeshNode*> > u2nodeVec;
@@ -349,19 +388,19 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool   isXConst,
       if ( node ) // nodes on internal vertices may be missing
       {
         if ( vertexNodes.insert( node ).second ||
-             fHelper.IsRealSeam( node->getshapeId() ))
+             fHelper.IsRealSeam  ( node->getshapeId() ) ||
+             fHelper.IsDegenShape( node->getshapeId() ))
           u2node.insert( u2node.end(), make_pair( prevNormPar, node ));
       }
       else if ( iE == 0 )
       {
-        if ( nodes.empty() ) {
-          for ( ++iE; iE < myEdge.size(); ++iE )
-            if (( node = VertexNode( iE ))) {
-              u2node.insert( make_pair( prevNormPar, node ));
-              break;
-            }
-          --iE;
-        }
+        for ( ++iE; iE < myEdge.size(); ++iE )
+          if (( node = VertexNode( iE ))) {
+            u2node.insert( make_pair( prevNormPar, node ));
+            break;
+          }
+        --iE;
+
         if ( !node )
           return myPoints;
         vertexNodes.insert( node );
@@ -410,26 +449,34 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool   isXConst,
       if ( u2node.empty() ) return myPoints;
 
       const SMDS_MeshNode* node;
-      if ( IsClosed() )
+      if ( IsClosed() && !proxySubMesh[0] )
         node = u2node.begin()->second;
       else
       {
         node = VertexNode( iE );
-        while ( !node && iE > 0 )
-          node = VertexNode( --iE );
-        if ( !node )
-          return myPoints;
+        if ( myProxyMesh->GetMesh()->HasModificationsToDiscard() )
+          while ( !node && iE > 1 ) // check intermediate VERTEXes
+            node = VertexNode( --iE );
+      }
+      if ( node )
+      {
+        if ( u2node.rbegin()->second == node &&
+             !fHelper.IsRealSeam  ( node->getshapeId() ) &&
+             !fHelper.IsDegenShape( node->getshapeId() ))
+          u2node.erase( --u2node.end() );
+
+        u2node.insert( u2node.end(), make_pair( 1., node ));
       }
-      if ( u2node.rbegin()->second == node )
-        u2node.erase( --u2node.end() );
-      u2node.insert( u2node.end(), make_pair( 1., node ));
     }
 
-    if ( u2node.size() + nbProxyNodes != myNbPonits &&
-         u2node.size() + nbProxyNodes != NbPoints( /*update=*/true ))
+    if ((int) u2node.size() + nbProxyNodes != myNbPonits &&
+        (int) u2node.size() + nbProxyNodes != NbPoints( /*update=*/true ))
+    {
+      return myPoints;
+    }
+    if (( u2node.size() > 0 ) &&
+        ( u2node.begin()->first < 0 || u2node.rbegin()->first > 1 ))
     {
-      MESSAGE("Wrong node parameters on edges, u2node.size():"
-              <<u2node.size()<<" !=  myNbPonits:"<<myNbPonits);
       return myPoints;
     }
 
@@ -446,25 +493,39 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool   isXConst,
       if ( proxySubMesh[ iE ] ) // copy data from a proxy sub-mesh
       {
         const UVPtStructVec& edgeUVPtStruct = proxySubMesh[iE]->GetUVPtStructVec();
-        std::copy( edgeUVPtStruct.begin(), edgeUVPtStruct.end(), & points[iPt] );
+        UVPtStruct* pointsPtr = & points[iPt];
+        std::copy( edgeUVPtStruct.begin(), edgeUVPtStruct.end(), pointsPtr );
         // check orientation
         double du1 = edgeUVPtStruct.back().param - edgeUVPtStruct[0].param;
         double du2 = myLast[iE] - myFirst[iE];
         if ( du1 * du2 < 0 )
         {
-          std::reverse( & points[iPt], & points[iPt + edgeUVPtStruct.size()]);
+          std::reverse( pointsPtr, pointsPtr + edgeUVPtStruct.size());
           for ( size_t i = 0; i < edgeUVPtStruct.size(); ++i )
-            points[iPt+i].normParam = 1. - points[iPt+i].normParam;
+            pointsPtr[i].normParam = 1. - pointsPtr[i].normParam;
         }
         // update normalized params
         if ( myEdge.size() > 1 ) {
-          for ( size_t i = 0; i < edgeUVPtStruct.size(); ++i, ++iPt )
+          for ( size_t i = 0; i < edgeUVPtStruct.size(); ++i )
           {
-            UVPtStruct & uvPt = points[iPt];
+            UVPtStruct & uvPt = pointsPtr[i];
             uvPt.normParam    = prevNormPar + uvPt.normParam * paramSize;
             uvPt.x = uvPt.y   = uvPt.normParam;
           }
-          --iPt; // to point to the 1st VERTEX of the next EDGE
+          iPt += edgeUVPtStruct.size() - 1; // to point to the 1st VERTEX of the next EDGE
+        }
+        // update UV on a seam EDGE
+        if ( fHelper.IsRealSeam( myEdgeID[ iE ]))
+        {
+          // check if points lye on the EDGE
+          const UVPtStruct& pm = edgeUVPtStruct[ edgeUVPtStruct.size()/2 ];
+          gp_Pnt         pNode = SMESH_TNodeXYZ( pm.node );
+          gp_Pnt         pCurv = myC3dAdaptor[ iE ].Value( pm.param );
+          double           tol = BRep_Tool::Tolerance( myEdge[ iE ]) * 10;
+          bool   isPointOnEdge = ( pNode.SquareDistance( pCurv ) < tol * tol );
+          if ( isPointOnEdge )
+            for ( size_t i = 0; i < edgeUVPtStruct.size(); ++i )
+              pointsPtr[i].SetUV( myC2d[ iE ]->Value( pointsPtr[i].param ).XY() );
         }
       }
       else
@@ -472,7 +533,7 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool   isXConst,
         for ( ; u_node != u2node.end(); ++u_node, ++iPt )
         {
           if ( myNormPar[ iE ]-eps < u_node->first )
-            break; // u_node is at VERTEX of the next EDGE 
+            break; // u_node is at VERTEX of the next EDGE
 
           UVPtStruct & uvPt = points[iPt];
           uvPt.node       = u_node->second;
@@ -510,9 +571,9 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool   isXConst,
 
     // set <constValue>
     if ( isXConst )
-      for ( iPt = 0; iPt < points.size(); ++iPt ) points[ iPt ].x = constValue;
+      for ( iPt = 0; iPt < (int)points.size(); ++iPt ) points[ iPt ].x = constValue;
     else
-      for ( iPt = 0; iPt < points.size(); ++iPt ) points[ iPt ].y = constValue;
+      for ( iPt = 0; iPt < (int)points.size(); ++iPt ) points[ iPt ].y = constValue;
 
   } // if ( myPoints.empty())
 
@@ -531,8 +592,8 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::SimulateUVPtStruct(int    nbSeg,
                                                                   bool   isXConst,
                                                                   double constValue) const
 {
-  if ( myFalsePoints.empty() ) {
-
+  if ( myFalsePoints.empty() )
+  {
     if ( NbEdges() == 0 ) return myFalsePoints;
 
     vector<uvPtStruct>* points = const_cast<vector<uvPtStruct>*>( &myFalsePoints );
@@ -540,28 +601,29 @@ const vector<UVPtStruct>& StdMeshers_FaceSide::SimulateUVPtStruct(int    nbSeg,
 
     int EdgeIndex = 0;
     double prevNormPar = 0, paramSize = myNormPar[ EdgeIndex ];
-    for (int i = 0 ; i < myFalsePoints.size(); ++i ) {
+    gp_Pnt2d p;
+    for ( size_t i = 0 ; i < myFalsePoints.size(); ++i )
+    {
       double normPar = double(i) / double(nbSeg);
       UVPtStruct & uvPt = (*points)[i];
       uvPt.node = 0;
       uvPt.x = uvPt.y = uvPt.param = uvPt.normParam = normPar;
       if ( isXConst ) uvPt.x = constValue;
       else            uvPt.y = constValue;
-      if ( myNormPar[ EdgeIndex ] < normPar ) {
+      if ( myNormPar[ EdgeIndex ] < normPar )
+      {
         prevNormPar = myNormPar[ EdgeIndex ];
         ++EdgeIndex;
         paramSize = myNormPar[ EdgeIndex ] - prevNormPar;
       }
       double r = ( normPar - prevNormPar )/ paramSize;
       uvPt.param = myFirst[EdgeIndex] * ( 1 - r ) + myLast[EdgeIndex] * r;
-      if ( !myC2d[ EdgeIndex ].IsNull() ) {
-        gp_Pnt2d p = myC2d[ EdgeIndex ]->Value( uvPt.param );
-        uvPt.u = p.X();
-        uvPt.v = p.Y();
-      }
-      else {
-        uvPt.u = uvPt.v = 1e+100;
-      }
+      if ( !myC2d[ EdgeIndex ].IsNull() )
+        p = myC2d[ EdgeIndex ]->Value( uvPt.param );
+      else
+        p = Value2d( normPar );
+      uvPt.u = p.X();
+      uvPt.v = p.Y();
     }
   }
   return myFalsePoints;
@@ -579,17 +641,16 @@ std::vector<const SMDS_MeshNode*> StdMeshers_FaceSide::GetOrderedNodes(int theEd
   {
     if ( NbEdges() == 0 ) return resultNodes;
 
-    SMESHDS_Mesh* meshDS = myProxyMesh->GetMeshDS();
-    SMESH_MesherHelper eHelper( *myProxyMesh->GetMesh() );
-    SMESH_MesherHelper fHelper( *myProxyMesh->GetMesh() );
-    fHelper.SetSubShape( myFace );
+    //SMESHDS_Mesh* meshDS = myProxyMesh->GetMeshDS();
+    SMESH_MesherHelper  eHelper( *myProxyMesh->GetMesh() );
+    SMESH_MesherHelper& fHelper = * FaceHelper();
     bool paramOK = true;
 
     // Sort nodes of all edges putting them into a map
 
-    map< double, const SMDS_MeshNode*>        u2node;
-    vector<const SMDS_MeshNode*>              nodes;
-    set<const SMDS_MeshNode*>                 vertexNodes;
+    map< double, const SMDS_MeshNode*> u2node;
+    vector<const SMDS_MeshNode*>       nodes;
+    set<const SMDS_MeshNode*>          vertexNodes;
     int iE = 0, iEnd = myEdge.size();
     if ( theEdgeInd >= 0 )
     {
@@ -613,7 +674,8 @@ std::vector<const SMDS_MeshNode*> StdMeshers_FaceSide::GetOrderedNodes(int theEd
       const SMDS_MeshNode* node = VertexNode( iE );
       if ( node ) { // nodes on internal vertices may be missing
         if ( vertexNodes.insert( node ).second ||
-             fHelper.IsRealSeam( node->getshapeId() ))
+             fHelper.IsRealSeam  ( node->getshapeId() ) ||
+             fHelper.IsDegenShape( node->getshapeId() ))
           u2node.insert( make_pair( prevNormPar, node ));
       }
       else if ( iE == 0 )
@@ -667,20 +729,26 @@ std::vector<const SMDS_MeshNode*> StdMeshers_FaceSide::GetOrderedNodes(int theEd
         if ( !node )
           return resultNodes;
       }
-      if ( u2node.rbegin()->second == node )
+      if ( u2node.rbegin()->second == node &&
+           !fHelper.IsRealSeam  ( node->getshapeId() ) &&
+           !fHelper.IsDegenShape( node->getshapeId() ))
         u2node.erase( --u2node.end() );
+
       u2node.insert( u2node.end(), make_pair( 1., node ));
     }
 
     // Fill the result vector
 
-    if ( u2node.size() == myNbPonits )
+    if ( theEdgeInd < 0 &&
+         (int) u2node.size() != myNbPonits &&
+         (int) u2node.size() != NbPoints( /*update=*/true ))
     {
-      resultNodes.reserve( u2node.size() );
-      map< double, const SMDS_MeshNode*>::iterator u2n = u2node.begin();
-      for ( ; u2n != u2node.end(); ++u2n )
-        resultNodes.push_back( u2n->second );
+      u2node.clear();
     }
+    resultNodes.reserve( u2node.size() );
+    map< double, const SMDS_MeshNode*>::iterator u2n = u2node.begin();
+    for ( ; u2n != u2node.end(); ++u2n )
+      resultNodes.push_back( u2n->second );
   }
   else
   {
@@ -800,7 +868,7 @@ const SMDS_MeshNode* StdMeshers_FaceSide::VertexNode(std::size_t i, bool* isMove
     n = SMESH_Algo::VertexNode( V, sm, myProxyMesh->GetMesh(), /*checkV=*/false );
 
     if (( !n ) &&
-        (( i > 0 && i < NbEdges() ) || IsClosed() ))
+        (( i > 0 && (int) i < NbEdges() ) || IsClosed() ))
     {
       iE = SMESH_MesherHelper::WrapIndex( int(i)-1, NbEdges() );
       sm = myProxyMesh->GetMeshDS()->MeshElements( myEdgeID[ iE ]);
@@ -953,7 +1021,7 @@ int StdMeshers_FaceSide::NbPoints(const bool update) const
     {
       if ( const SMESHDS_SubMesh* sm = myProxyMesh->GetSubMesh( Edge(i) ))
       {
-        if ( sm->NbNodes() == sm->NbElements() - 1 )
+        if ( sm->NbNodes() == sm->NbElements()-1 || sm->NbElements() == 0 )
         {
           me->myNbPonits += sm->NbNodes();
           if ( myIgnoreMediumNodes && sm->IsQuadratic() )
@@ -969,16 +1037,16 @@ int StdMeshers_FaceSide::NbPoints(const bool update) const
       }
     }
 
-    SMESH_MesherHelper helper( *myProxyMesh->GetMesh() );
-    helper.SetSubShape( myFace );
+    SMESH_MesherHelper* helper = FaceHelper();
 
     std::set< const SMDS_MeshNode* > vNodes;
-    for ( int i = 0; i <= NbEdges(); ++i ) // nb VERTEXes is more than NbEdges() if !IsClosed()
+    const int nbV = NbEdges() + !IsClosed();
+    for ( int i = 0; i < nbV; ++i )
       if ( const SMDS_MeshNode* n = VertexNode( i ))
       {
         if ( !vNodes.insert( n ).second &&
-             helper.IsRealSeam( n->getshapeId() ) &&
-             i < NbEdges())
+             ( helper->IsRealSeam  ( n->getshapeId() ) ||
+               helper->IsDegenShape( n->getshapeId() )))
           me->myNbPonits++;
       }
       else
@@ -1040,7 +1108,7 @@ void StdMeshers_FaceSide::dump(const char* msg) const
   if (msg) MESSAGE ( std::endl << msg );
   MESSAGE_BEGIN ("NB EDGES: "<< myEdge.size() );
   MESSAGE_ADD ( "nbPoints: "<<myNbPonits<<" vecSize: " << myPoints.size()<<" "<<myFalsePoints.size() );
-  for ( int i=0; i<myEdge.size(); ++i)
+  for ( size_t i = 0; i < myEdge.size(); ++i )
   {
     MESSAGE_ADD ( "\t"<<i+1 );
     MESSAGE_ADD ( "\tEDGE: " );
@@ -1050,17 +1118,17 @@ void StdMeshers_FaceSide::dump(const char* msg) const
     else {
       TopAbs::Print(myEdge[i].Orientation(),cout)<<" "<<myEdge[i].TShape().operator->()<<endl;
       MESSAGE_ADD ( "\tV1: " << TopExp::FirstVertex( myEdge[i], 1).TShape().operator->()
-                 << "  V2: " << TopExp::LastVertex( myEdge[i], 1).TShape().operator->() );
+                    << "  V2: " << TopExp::LastVertex( myEdge[i], 1).TShape().operator->() );
     }
     MESSAGE_ADD ( "\tC2d: ");
-    
+
     if (myC2d[i].IsNull()) {
       MESSAGE_ADD ( "NULL" );
     }
     else {
       MESSAGE_ADD ( myC2d[i].operator->() );
     }
-      
+
     MESSAGE_ADD ( "\tF: "<<myFirst[i]<< " L: "<< myLast[i] );
     MESSAGE_END ( "\tnormPar: "<<myNormPar[i]<<endl );
   }
@@ -1102,7 +1170,7 @@ BRepAdaptor_CompCurve* StdMeshers_FaceSide::GetCurve3d() const
   TopoDS_Wire aWire;
   BRep_Builder aBuilder;
   aBuilder.MakeWire(aWire);
-  for ( int i=0; i<myEdge.size(); ++i )
+  for ( size_t i = 0; i < myEdge.size(); ++i )
     aBuilder.Add( aWire, myEdge[i] );
 
   if ( myEdge.size() == 2 && IsClosed() )
@@ -1146,7 +1214,7 @@ gp_Pnt2d StdMeshers_FaceSide::Value2d(double U) const
     int i = U * double( myPoints.size()-1 );
     while ( i > 0 && myPoints[ i ].normParam > U )
       --i;
-    while ( i+1 < myPoints.size() && myPoints[ i+1 ].normParam < U )
+    while ( i+1 < (int)myPoints.size() && myPoints[ i+1 ].normParam < U )
       ++i;
     double r = (( U - myPoints[ i ].normParam ) /
                 ( myPoints[ i+1 ].normParam - myPoints[ i ].normParam ));
@@ -1195,9 +1263,14 @@ TSideVector StdMeshers_FaceSide::GetFaceWires(const TopoDS_Face&   theFace,
                                               SMESH_Mesh &         theMesh,
                                               const bool           theIgnoreMediumNodes,
                                               TError &             theError,
+                                              SMESH_MesherHelper*  theFaceHelper,
                                               SMESH_ProxyMesh::Ptr theProxyMesh,
                                               const bool           theCheckVertexNodes)
 {
+  SMESH_MesherHelper helper( theMesh );
+  if ( theFaceHelper && theFaceHelper->GetSubShape() == theFace )
+    helper.CopySubShapeInfo( *theFaceHelper );
+
   list< TopoDS_Edge > edges, internalEdges;
   list< int > nbEdgesInWires;
   int nbWires = SMESH_Block::GetOrderedEdges (theFace, edges, nbEdgesInWires);
@@ -1241,7 +1314,7 @@ TSideVector StdMeshers_FaceSide::GetFaceWires(const TopoDS_Face&   theFace,
 
     StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( theFace, wireEdges, &theMesh,
                                                          /*isForward=*/true, theIgnoreMediumNodes,
-                                                         theProxyMesh );
+                                                         &helper, theProxyMesh );
     wires[ iW ] = StdMeshers_FaceSidePtr( wire );
     from = to;
   }
@@ -1249,7 +1322,7 @@ TSideVector StdMeshers_FaceSide::GetFaceWires(const TopoDS_Face&   theFace,
   {
     StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( theFace, internalEdges.back(), &theMesh,
                                                          /*isForward=*/true, theIgnoreMediumNodes,
-                                                         theProxyMesh );
+                                                         &helper, theProxyMesh );
     wires.push_back( StdMeshers_FaceSidePtr( wire ));
     internalEdges.pop_back();
   }
@@ -1305,3 +1378,20 @@ bool StdMeshers_FaceSide::IsClosed() const
 {
   return myEdge.empty() ? false : FirstVertex().IsSame( LastVertex() );
 }
+
+//================================================================================
+/*!
+ * \brief Return a helper initialized with the FACE
+ */
+//================================================================================
+
+SMESH_MesherHelper* StdMeshers_FaceSide::FaceHelper() const
+{
+  StdMeshers_FaceSide* me = const_cast< StdMeshers_FaceSide* >( this );
+  if ( !myHelper && myProxyMesh )
+  {
+    me->myHelper = new SMESH_MesherHelper( *myProxyMesh->GetMesh() );
+    me->myHelper->SetSubShape( myFace );
+  }
+  return me->myHelper;
+}