Salome HOME
IPAL53561: Extrusion 3D algo is not applicable to a prism with a seam on top/bottom
authoreap <eap@opencascade.com>
Wed, 7 Sep 2016 11:43:02 +0000 (14:43 +0300)
committereap <eap@opencascade.com>
Wed, 7 Sep 2016 11:43:02 +0000 (14:43 +0300)
src/SMESH_SWIG/StdMeshersBuilder.py
src/StdMeshers/CMakeLists.txt
src/StdMeshers/StdMeshers_Prism_3D.cxx
src/StdMeshers/StdMeshers_Prism_3D.hxx
src/StdMeshers/StdMeshers_ProjectionUtils.cxx

index 35e43e551ba6d5b8ce79acfb63de41e21cfb4fec..203a9ff480e4b8483270e56d4f1f5fc3da6dd428 100644 (file)
@@ -21,6 +21,8 @@
 # @package StdMeshersBuilder
 # Python API for the standard meshing plug-in module.
 
 # @package StdMeshersBuilder
 # Python API for the standard meshing plug-in module.
 
+LIBRARY = "libStdMeshersEngine.so"
+
 from salome.smesh.smesh_algorithm import Mesh_Algorithm
 import StdMeshers
 
 from salome.smesh.smesh_algorithm import Mesh_Algorithm
 import StdMeshers
 
@@ -961,10 +963,8 @@ class StdMeshersBuilder_Prism3D(Mesh_Algorithm):
         shape = geom
         if not shape:
             shape = mesh.geom
         shape = geom
         if not shape:
             shape = mesh.geom
-        from salome.geom import geomBuilder
-        nbSolids = len( geomBuilder.geom.SubShapeAll( shape, geomBuilder.geomBuilder.ShapeType["SOLID"] ))
-        nbShells = len( geomBuilder.geom.SubShapeAll( shape, geomBuilder.geomBuilder.ShapeType["SHELL"] ))
-        if nbSolids == 0 or nbSolids == nbShells:
+        isRadial = mesh.smeshpyD.IsApplicable("RadialPrism_3D", LIBRARY, shape, False )
+        if not isRadial:
             self.Create(mesh, geom, "Prism_3D")
             pass
         else:
             self.Create(mesh, geom, "Prism_3D")
             pass
         else:
index 3a1348297710287166b4abdab4fde5ddc5b750a2..c0097c64a27f2b02d3039ab977164d89605cbf9e 100644 (file)
@@ -29,7 +29,6 @@ INCLUDE_DIRECTORIES(
   ${Boost_INCLUDE_DIRS}
   ${VTK_INCLUDE_DIRS}
   ${KERNEL_INCLUDE_DIRS}
   ${Boost_INCLUDE_DIRS}
   ${VTK_INCLUDE_DIRS}
   ${KERNEL_INCLUDE_DIRS}
-  ${GUI_INCLUDE_DIRS}
   ${GEOM_INCLUDE_DIRS}
   ${PROJECT_SOURCE_DIR}/src/SMESHUtils
   ${PROJECT_SOURCE_DIR}/src/SMESH
   ${GEOM_INCLUDE_DIRS}
   ${PROJECT_SOURCE_DIR}/src/SMESHUtils
   ${PROJECT_SOURCE_DIR}/src/SMESH
index 5734c02e99354a81c2242fb804dc4edc558be3b0..9aed55e3e88f7ea1e16b69652f2b98737fc31c31 100644 (file)
@@ -944,11 +944,10 @@ bool StdMeshers_Prism_3D::getWallFaces( Prism_3D::TPrismTopo & thePrism,
   std::list< int >::iterator     nbE = thePrism.myNbEdgesInWires.begin();
   std::list< int > nbQuadsPerWire;
   int iE = 0;
   std::list< int >::iterator     nbE = thePrism.myNbEdgesInWires.begin();
   std::list< int > nbQuadsPerWire;
   int iE = 0;
-  double f,l;
   while ( edge != thePrism.myBottomEdges.end() )
   {
     ++iE;
   while ( edge != thePrism.myBottomEdges.end() )
   {
     ++iE;
-    if ( BRep_Tool::Curve( *edge, f,l ).IsNull() )
+    if ( SMESH_Algo::isDegenerated( *edge ))
     {
       edge = thePrism.myBottomEdges.erase( edge );
       --iE;
     {
       edge = thePrism.myBottomEdges.erase( edge );
       --iE;
@@ -956,12 +955,14 @@ bool StdMeshers_Prism_3D::getWallFaces( Prism_3D::TPrismTopo & thePrism,
     }
     else
     {
     }
     else
     {
+      bool hasWallFace = false;
       TopTools_ListIteratorOfListOfShape faceIt( edgeToFaces.FindFromKey( *edge ));
       for ( ; faceIt.More(); faceIt.Next() )
       {
         const TopoDS_Face& face = TopoDS::Face( faceIt.Value() );
         if ( !thePrism.myBottom.IsSame( face ))
         {
       TopTools_ListIteratorOfListOfShape faceIt( edgeToFaces.FindFromKey( *edge ));
       for ( ; faceIt.More(); faceIt.Next() )
       {
         const TopoDS_Face& face = TopoDS::Face( faceIt.Value() );
         if ( !thePrism.myBottom.IsSame( face ))
         {
+          hasWallFace = true;
           Prism_3D::TQuadList quadList( 1, quadAlgo->CheckNbEdges( *mesh, face ));
           if ( !quadList.back() )
             return toSM( error(TCom("Side face #") << shapeID( face )
           Prism_3D::TQuadList quadList( 1, quadAlgo->CheckNbEdges( *mesh, face ));
           if ( !quadList.back() )
             return toSM( error(TCom("Side face #") << shapeID( face )
@@ -980,7 +981,16 @@ bool StdMeshers_Prism_3D::getWallFaces( Prism_3D::TPrismTopo & thePrism,
           break;
         }
       }
           break;
         }
       }
-      ++edge;
+      if ( hasWallFace )
+      {
+        ++edge;
+      }
+      else // seam edge (IPAL53561)
+      {
+        edge = thePrism.myBottomEdges.erase( edge );
+        --iE;
+        --(*nbE);
+      }
     }
     if ( iE == *nbE )
     {
     }
     if ( iE == *nbE )
     {
@@ -1210,7 +1220,8 @@ bool StdMeshers_Prism_3D::compute(const Prism_3D::TPrismTopo& thePrism)
     for ( ; bot_column != myBotToColumnMap.end(); ++bot_column )
     {
       const Prism_3D::TNode& tBotNode = bot_column->first; // bottom TNode
     for ( ; bot_column != myBotToColumnMap.end(); ++bot_column )
     {
       const Prism_3D::TNode& tBotNode = bot_column->first; // bottom TNode
-      if ( tBotNode.GetPositionType() != SMDS_TOP_FACE )
+      if ( tBotNode.GetPositionType() != SMDS_TOP_FACE &&
+           myBlock.HasNodeColumn( tBotNode.myNode ))
         continue; // node is not inside the FACE
 
       // column nodes; middle part of the column are zero pointers
         continue; // node is not inside the FACE
 
       // column nodes; middle part of the column are zero pointers
@@ -1314,16 +1325,17 @@ bool StdMeshers_Prism_3D::compute(const Prism_3D::TPrismTopo& thePrism)
     for ( int i = 0; i < nbNodes; ++i )
     {
       const SMDS_MeshNode* n = face->GetNode( i );
     for ( int i = 0; i < nbNodes; ++i )
     {
       const SMDS_MeshNode* n = face->GetNode( i );
-      if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ) {
+      columns[ i ] = NULL;
+
+      if ( n->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE )
+        columns[ i ] = myBlock.GetNodeColumn( n );
+
+      if ( !columns[ i ] )
+      {
         TNode2ColumnMap::iterator bot_column = myBotToColumnMap.find( n );
         if ( bot_column == myBotToColumnMap.end() )
         TNode2ColumnMap::iterator bot_column = myBotToColumnMap.find( n );
         if ( bot_column == myBotToColumnMap.end() )
-          return toSM( error(TCom("No nodes found above node ") << n->GetID() ));
-        columns[ i ] = & bot_column->second;
-      }
-      else {
-        columns[ i ] = myBlock.GetNodeColumn( n );
-        if ( !columns[ i ] )
           return toSM( error(TCom("No side nodes found above node ") << n->GetID() ));
           return toSM( error(TCom("No side nodes found above node ") << n->GetID() ));
+        columns[ i ] = & bot_column->second;
       }
     }
     // create prisms
       }
     }
     // create prisms
@@ -2074,7 +2086,8 @@ bool StdMeshers_Prism_3D::assocOrProjBottom2Top( const gp_Trsf & bottomToTopTrsf
   {
     const SMDS_MeshNode* botNode = bN_tN->first;
     const SMDS_MeshNode* topNode = bN_tN->second;
   {
     const SMDS_MeshNode* botNode = bN_tN->first;
     const SMDS_MeshNode* topNode = bN_tN->second;
-    if ( botNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE )
+    if ( botNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE &&
+         myBlock.HasNodeColumn( botNode ))
       continue; // wall columns are contained in myBlock
     // create node column
     Prism_3D::TNode bN( botNode );
       continue; // wall columns are contained in myBlock
     // create node column
     Prism_3D::TNode bN( botNode );
@@ -2523,13 +2536,15 @@ namespace // utils used by StdMeshers_Prism_3D::IsApplicable()
   struct EdgeWithNeighbors
   {
     TopoDS_Edge _edge;
   struct EdgeWithNeighbors
   {
     TopoDS_Edge _edge;
-    int         _iL, _iR;
-    EdgeWithNeighbors(const TopoDS_Edge& E, int iE, int nbE, int shift = 0 ):
-      _edge( E ),
+    int         _iBase;   /* index in a WIRE with non-base EDGEs excluded */
+    int         _iL, _iR; /* used to connect edges in a base FACE */
+    bool        _isBase;  /* is used in a base FACE */
+    EdgeWithNeighbors(const TopoDS_Edge& E, int iE, int nbE, int shift, bool isBase ):
+      _edge( E ), _iBase( iE + shift ),
       _iL( SMESH_MesherHelper::WrapIndex( iE-1, nbE ) + shift ),
       _iL( SMESH_MesherHelper::WrapIndex( iE-1, nbE ) + shift ),
-      _iR( SMESH_MesherHelper::WrapIndex( iE+1, nbE ) + shift )
+      _iR( SMESH_MesherHelper::WrapIndex( iE+1, nbE ) + shift ),
+      _isBase( isBase )
     {
     {
-      //_edge.Orientation( TopAbs_FORWARD ); // for operator==() to work
     }
     EdgeWithNeighbors() {}
     bool IsInternal() const { return !_edge.IsNull() && _edge.Orientation() == TopAbs_INTERNAL; }
     }
     EdgeWithNeighbors() {}
     bool IsInternal() const { return !_edge.IsNull() && _edge.Orientation() == TopAbs_INTERNAL; }
@@ -2571,13 +2586,29 @@ namespace // utils used by StdMeshers_Prism_3D::IsApplicable()
       return false;
     }
   };
       return false;
     }
   };
+  //--------------------------------------------------------------------------------
+  /*!
+   * \brief Return another faces sharing an edge
+   */
+  const TopoDS_Face & getAnotherFace( const TopoDS_Face& face,
+                                      const TopoDS_Edge& edge,
+                                      TopTools_IndexedDataMapOfShapeListOfShape& facesOfEdge)
+  {
+    TopTools_ListIteratorOfListOfShape faceIt( facesOfEdge.FindFromKey( edge ));
+    for ( ; faceIt.More(); faceIt.Next() )
+      if ( !face.IsSame( faceIt.Value() ))
+        return TopoDS::Face( faceIt.Value() );
+    return face;
+  }
+
   //--------------------------------------------------------------------------------
   /*!
    * \brief Return ordered edges of a face
    */
   //--------------------------------------------------------------------------------
   /*!
    * \brief Return ordered edges of a face
    */
-  bool getEdges( const TopoDS_Face&            face,
-                 vector< EdgeWithNeighbors > & edges,
-                 const bool                    noHolesAllowed)
+  bool getEdges( const TopoDS_Face&                         face,
+                 vector< EdgeWithNeighbors > &              edges,
+                 TopTools_IndexedDataMapOfShapeListOfShape& facesOfEdge,
+                 const bool                                 noHolesAllowed)
   {
     TopoDS_Face f = face;
     if ( f.Orientation() != TopAbs_FORWARD &&
   {
     TopoDS_Face f = face;
     if ( f.Orientation() != TopAbs_FORWARD &&
@@ -2589,26 +2620,40 @@ namespace // utils used by StdMeshers_Prism_3D::IsApplicable()
     if ( nbW > 1 && noHolesAllowed )
       return false;
 
     if ( nbW > 1 && noHolesAllowed )
       return false;
 
-    int iE, nbTot = 0;
+    int iE, nbTot = 0, nbBase, iBase;
     list< TopoDS_Edge >::iterator   e = ee.begin();
     list< int         >::iterator nbE = nbEdgesInWires.begin();
     for ( ; nbE != nbEdgesInWires.end(); ++nbE )
       for ( iE = 0; iE < *nbE; ++e, ++iE )
     list< TopoDS_Edge >::iterator   e = ee.begin();
     list< int         >::iterator nbE = nbEdgesInWires.begin();
     for ( ; nbE != nbEdgesInWires.end(); ++nbE )
       for ( iE = 0; iE < *nbE; ++e, ++iE )
-        if ( SMESH_Algo::isDegenerated( *e ))
+        if ( SMESH_Algo::isDegenerated( *e )) // degenerated EDGE is never used
         {
           e = --ee.erase( e );
           --(*nbE);
           --iE;
         }
 
         {
           e = --ee.erase( e );
           --(*nbE);
           --iE;
         }
 
+    vector<int> isBase;
     edges.clear();
     e = ee.begin();
     for ( nbE = nbEdgesInWires.begin(); nbE != nbEdgesInWires.end(); ++nbE )
     {
     edges.clear();
     e = ee.begin();
     for ( nbE = nbEdgesInWires.begin(); nbE != nbEdgesInWires.end(); ++nbE )
     {
-      for ( iE = 0; iE < *nbE; ++e, ++iE )
-        edges.push_back( EdgeWithNeighbors( *e, iE, *nbE, nbTot ));
-      nbTot += *nbE;
+      nbBase = 0;
+      isBase.resize( *nbE );
+      list< TopoDS_Edge >::iterator eIt = e;
+      for ( iE = 0; iE < *nbE; ++eIt, ++iE )
+      {
+        isBase[ iE ] = ( getAnotherFace( face, *eIt, facesOfEdge ) != face );
+        nbBase += isBase[ iE ];
+      }
+      for ( iBase = 0, iE = 0; iE < *nbE; ++e, ++iE )
+      {
+        edges.push_back( EdgeWithNeighbors( *e, iBase, nbBase, nbTot, isBase[ iE ] ));
+        iBase += isBase[ iE ];
+      }
+      nbTot += nbBase;
     }
     }
+    if ( nbTot == 0 )
+      return false;
 
     // IPAL53099. Set correct neighbors to INTERNAL EDGEs, which can be connected to
     // EDGEs of the outer WIRE but this fact can't be detected by their order.
 
     // IPAL53099. Set correct neighbors to INTERNAL EDGEs, which can be connected to
     // EDGEs of the outer WIRE but this fact can't be detected by their order.
@@ -2623,6 +2668,9 @@ namespace // utils used by StdMeshers_Prism_3D::IsApplicable()
         bool isConnectOk = ( vv[0].IsSame( vv[1] ));
         if ( !isConnectOk )
         {
         bool isConnectOk = ( vv[0].IsSame( vv[1] ));
         if ( !isConnectOk )
         {
+          edges[ iFirst ]._iL = edges[ iFirst ]._iBase; // connect to self
+          edges[ iLast  ]._iR = edges[ iLast ]._iBase;
+
           // look for an EDGE of the outer WIRE connected to vv
           TopoDS_Vertex v0, v1;
           for ( iE = 0; iE < nbEdgesInWires.front(); ++iE )
           // look for an EDGE of the outer WIRE connected to vv
           TopoDS_Vertex v0, v1;
           for ( iE = 0; iE < nbEdgesInWires.front(); ++iE )
@@ -2630,9 +2678,9 @@ namespace // utils used by StdMeshers_Prism_3D::IsApplicable()
             v0 = SMESH_MesherHelper::IthVertex( 0, edges[ iE ]._edge );
             v1 = SMESH_MesherHelper::IthVertex( 1, edges[ iE ]._edge );
             if ( vv[0].IsSame( v0 ) || vv[0].IsSame( v1 ))
             v0 = SMESH_MesherHelper::IthVertex( 0, edges[ iE ]._edge );
             v1 = SMESH_MesherHelper::IthVertex( 1, edges[ iE ]._edge );
             if ( vv[0].IsSame( v0 ) || vv[0].IsSame( v1 ))
-              edges[ iFirst ]._iL = iE;
+              edges[ iFirst ]._iL = edges[ iE ]._iBase;
             if ( vv[1].IsSame( v0 ) || vv[1].IsSame( v1 ))
             if ( vv[1].IsSame( v0 ) || vv[1].IsSame( v1 ))
-              edges[ iLast ]._iR = iE;
+              edges[ iLast ]._iR = edges[ iE ]._iBase;
           }
         }
         iFirst += *nbE;
           }
         }
         iFirst += *nbE;
@@ -2640,21 +2688,7 @@ namespace // utils used by StdMeshers_Prism_3D::IsApplicable()
     }
     return edges.size();
   }
     }
     return edges.size();
   }
-  //--------------------------------------------------------------------------------
-  /*!
-   * \brief Return another faces sharing an edge
-   */
-  const TopoDS_Face & getAnotherFace( const TopoDS_Face& face,
-                                      const TopoDS_Edge& edge,
-                                      TopTools_IndexedDataMapOfShapeListOfShape& facesOfEdge)
-  {
-    TopTools_ListIteratorOfListOfShape faceIt( facesOfEdge.FindFromKey( edge ));
-    for ( ; faceIt.More(); faceIt.Next() )
-      if ( !face.IsSame( faceIt.Value() ))
-        return TopoDS::Face( faceIt.Value() );
-    return face;
-  }
-
+  
   //--------------------------------------------------------------------------------
   /*!
    * \brief Return number of faces sharing given edges
   //--------------------------------------------------------------------------------
   /*!
    * \brief Return number of faces sharing given edges
@@ -2694,10 +2728,10 @@ bool StdMeshers_Prism_3D::IsApplicable(const TopoDS_Shape & shape, bool toCheckA
     // check nb shells
     TopoDS_Shape shell;
     TopExp_Explorer shExp( sExp.Current(), TopAbs_SHELL );
     // check nb shells
     TopoDS_Shape shell;
     TopExp_Explorer shExp( sExp.Current(), TopAbs_SHELL );
-    if ( shExp.More() ) {
+    while ( shExp.More() ) {
       shell = shExp.Current();
       shExp.Next();
       shell = shExp.Current();
       shExp.Next();
-      if ( shExp.More() )
+      if ( shExp.More() && BRep_Tool::IsClosed( shExp.Current() ))
         shell.Nullify();
     }
     if ( shell.IsNull() ) {
         shell.Nullify();
     }
     if ( shell.IsNull() ) {
@@ -2706,7 +2740,7 @@ bool StdMeshers_Prism_3D::IsApplicable(const TopoDS_Shape & shape, bool toCheckA
     }
     // get all faces
     TopTools_IndexedMapOfShape allFaces;
     }
     // get all faces
     TopTools_IndexedMapOfShape allFaces;
-    TopExp::MapShapes( shell, TopAbs_FACE, allFaces );
+    TopExp::MapShapes( sExp.Current(), TopAbs_FACE, allFaces );
     if ( allFaces.Extent() < 3 ) {
       if ( toCheckAll ) return false;
       continue;
     if ( allFaces.Extent() < 3 ) {
       if ( toCheckAll ) return false;
       continue;
@@ -2723,7 +2757,7 @@ bool StdMeshers_Prism_3D::IsApplicable(const TopoDS_Shape & shape, bool toCheckA
       }
     }
 #ifdef _DEBUG_
       }
     }
 #ifdef _DEBUG_
-    TopTools_IndexedMapOfShape allShapes;
+    TopTools_IndexedMapOfShape allShapes; // usage: allShapes.FindIndex( s )
     TopExp::MapShapes( shape, allShapes );
 #endif
 
     TopExp::MapShapes( shape, allShapes );
 #endif
 
@@ -2750,23 +2784,32 @@ bool StdMeshers_Prism_3D::IsApplicable(const TopoDS_Shape & shape, bool toCheckA
 
       TEdgeWithNeighborsVec& botEdges = faceEdgesVec[ iF ];
       if ( botEdges.empty() )
 
       TEdgeWithNeighborsVec& botEdges = faceEdgesVec[ iF ];
       if ( botEdges.empty() )
-        if ( !getEdges( botF, botEdges, /*noHoles=*/false ))
+        if ( !getEdges( botF, botEdges, facesOfEdge, /*noHoles=*/false ))
           break;
           break;
-      if ( allFaces.Extent()-1 <= (int) botEdges.size() )
+
+      int nbBase = 0;
+      for ( size_t iS = 0; iS < botEdges.size(); ++iS )
+        nbBase += botEdges[ iS ]._isBase;
+
+      if ( allFaces.Extent()-1 <= nbBase )
         continue; // all faces are adjacent to botF - no top FACE
 
       // init data of side FACEs
       sides.clear();
         continue; // all faces are adjacent to botF - no top FACE
 
       // init data of side FACEs
       sides.clear();
-      sides.resize( botEdges.size() );
-      for ( size_t iS = 0; iS < botEdges.size(); ++iS )
+      sides.resize( nbBase );
+      size_t iS = 0;
+      for ( size_t iE = 0; iE < botEdges.size(); ++iE )
       {
       {
-        sides[ iS ]._topEdge    = botEdges[ iS ]._edge;
+        if ( !botEdges[ iE ]._isBase )
+          continue;
+        sides[ iS ]._topEdge    = botEdges[ iE ]._edge;
         sides[ iS ]._face       = botF;
         sides[ iS ]._face       = botF;
-        sides[ iS ]._leftSide   = & sides[ botEdges[ iS ]._iR ];
-        sides[ iS ]._rightSide  = & sides[ botEdges[ iS ]._iL ];
-        sides[ iS ]._isInternal = botEdges[ iS ].IsInternal();
+        sides[ iS ]._leftSide   = & sides[ botEdges[ iE ]._iR ];
+        sides[ iS ]._rightSide  = & sides[ botEdges[ iE ]._iL ];
+        sides[ iS ]._isInternal = botEdges[ iE ].IsInternal();
         sides[ iS ]._faces      = & facesOfSide[ iS ];
         sides[ iS ]._faces->Clear();
         sides[ iS ]._faces      = & facesOfSide[ iS ];
         sides[ iS ]._faces->Clear();
+        ++iS;
       }
 
       bool isOK = true; // ok for a current botF
       }
 
       bool isOK = true; // ok for a current botF
@@ -2861,7 +2904,7 @@ bool StdMeshers_Prism_3D::IsApplicable(const TopoDS_Shape & shape, bool toCheckA
             int faceID  = allFaces.FindIndex( side._face );
             side._edges = & faceEdgesVec[ faceID ];
             if ( side._edges->empty() )
             int faceID  = allFaces.FindIndex( side._face );
             side._edges = & faceEdgesVec[ faceID ];
             if ( side._edges->empty() )
-              if ( !getEdges( side._face, * side._edges, /*noHoles=*/true ))
+              if ( !getEdges( side._face, * side._edges, facesOfEdge, /*noHoles=*/true ))
                 break;
             const int nbE = side._edges->size();
             if ( nbE >= 4 )
                 break;
             const int nbE = side._edges->size();
             if ( nbE >= 4 )
index c113fd1e867a519e52bb987062ac0b578628bb09..ddfbff96a13a5c917b54248cfedad4d506334b9b 100644 (file)
@@ -185,6 +185,16 @@ class STDMESHERS_EXPORT StdMeshers_PrismAsBlock: public SMESH_Block
     return col_frw.first;
   }
 
     return col_frw.first;
   }
 
+  /*!
+   * \brief Return pointer to column of nodes
+    * \param node - bottom node from which the returned column goes up
+    * \retval const TNodeColumn* - the found column
+   */
+  bool HasNodeColumn(const SMDS_MeshNode* node) const
+  {
+    return myShapeIndex2ColumnMap.count( node->getshapeId() );
+  }
+
   /*!
    * \brief Return transformations to get coordinates of nodes of each internal layer
    *        by nodes of the bottom. Layer is a set of nodes at a certain step
   /*!
    * \brief Return transformations to get coordinates of nodes of each internal layer
    *        by nodes of the bottom. Layer is a set of nodes at a certain step
index d81ade633a2008dcad61c8935e6dddbdc839b790..c4b082adcc93c72d16d06cb72c166f5698a06ac7 100644 (file)
@@ -414,9 +414,8 @@ namespace {
                      const gp_Pnt2d&    uv,
                      const double&      tol2d )
   {
                      const gp_Pnt2d&    uv,
                      const double&      tol2d )
   {
-    TopoDS_Vertex VV[2];
-    TopExp::Vertices( edge, VV[0], VV[1], true);
-    gp_Pnt2d v1UV = BRep_Tool::Parameters( VV[vIndex], face);
+    TopoDS_Vertex V = SMESH_MesherHelper::IthVertex( vIndex, edge, /*CumOri=*/true );
+    gp_Pnt2d v1UV = BRep_Tool::Parameters( V, face);
     double dist2d = v1UV.Distance( uv );
     return dist2d < tol2d;
   }
     double dist2d = v1UV.Distance( uv );
     return dist2d < tol2d;
   }
@@ -1111,8 +1110,8 @@ bool StdMeshers_ProjectionUtils::FindSubShapeAssociation(const TopoDS_Shape& the
         for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
         {
           InsertAssociation( *eIt1, *eIt2, theMap );
         for ( ; eIt1 != edges1.end(); ++eIt1, ++eIt2 )
         {
           InsertAssociation( *eIt1, *eIt2, theMap );
-          VV1[0] = TopExp::FirstVertex( *eIt1, true );
-          VV2[0] = TopExp::FirstVertex( *eIt2, true );
+          VV1[0] = SMESH_MesherHelper::IthVertex( 0, *eIt1, true );
+          VV2[0] = SMESH_MesherHelper::IthVertex( 0, *eIt2, true );
           InsertAssociation( VV1[0], VV2[0], theMap );
         }
         InsertAssociation( theShape1, theShape2, theMap );
           InsertAssociation( VV1[0], VV2[0], theMap );
         }
         InsertAssociation( theShape1, theShape2, theMap );