Salome HOME
54122: Bad quality prismatic mesh
[modules/smesh.git] / src / StdMeshers / StdMeshers_Prism_3D.cxx
index 5e253829ab1af52a83a2b8a25489a086df85af4a..0ee4041ce5ef74b936d21ad7b4fc5273cea105ff 100644 (file)
@@ -51,7 +51,7 @@
 #include <Geom2d_Line.hxx>
 #include <GeomLib_IsPlanarSurface.hxx>
 #include <Geom_Curve.hxx>
-#include <TColStd_DataMapOfIntegerInteger.hxx>
+#include <Standard_ErrorHandler.hxx>
 #include <TopExp.hxx>
 #include <TopExp_Explorer.hxx>
 #include <TopTools_ListIteratorOfListOfShape.hxx>
@@ -529,6 +529,27 @@ namespace {
     return nbSides;
   }
 
+  //================================================================================
+  /*!
+   * \brief Set/get wire index to FaceQuadStruct
+   */
+  //================================================================================
+
+  void setWireIndex( TFaceQuadStructPtr& quad, int iWire )
+  {
+    quad->iSize = iWire;
+  }
+  int getWireIndex( const TFaceQuadStructPtr& quad )
+  {
+    return quad->iSize;
+  }
+
+  //================================================================================
+  /*!
+   * \brief Print Python commands adding given points to a mesh
+   */
+  //================================================================================
+
   void pointsToPython(const std::vector<gp_XYZ>& p)
   {
 #ifdef _DEBUG_
@@ -559,6 +580,7 @@ StdMeshers_Prism_3D::StdMeshers_Prism_3D(int hypId, int studyId, SMESH_Gen* gen)
 
   //myProjectTriangles       = false;
   mySetErrorToSM           = true;  // to pass an error to a sub-mesh of a current solid or not
+  myPrevBottomSM           = 0;     // last treated bottom sub-mesh with a suitable algorithm
 }
 
 //================================================================================
@@ -581,39 +603,6 @@ bool StdMeshers_Prism_3D::CheckHypothesis(SMESH_Mesh&                          a
                                           const TopoDS_Shape&                  aShape,
                                           SMESH_Hypothesis::Hypothesis_Status& aStatus)
 {
-  // Check shape geometry
-/*  PAL16229
-  aStatus = SMESH_Hypothesis::HYP_BAD_GEOMETRY;
-
-  // find not quadrangle faces
-  list< TopoDS_Shape > notQuadFaces;
-  int nbEdge, nbWire, nbFace = 0;
-  TopExp_Explorer exp( aShape, TopAbs_FACE );
-  for ( ; exp.More(); exp.Next() ) {
-    ++nbFace;
-    const TopoDS_Shape& face = exp.Current();
-    nbEdge = NSProjUtils::Count( face, TopAbs_EDGE, 0 );
-    nbWire = NSProjUtils::Count( face, TopAbs_WIRE, 0 );
-    if (  nbEdge!= 4 || nbWire!= 1 ) {
-      if ( !notQuadFaces.empty() ) {
-        if ( NSProjUtils::Count( notQuadFaces.back(), TopAbs_EDGE, 0 ) != nbEdge ||
-             NSProjUtils::Count( notQuadFaces.back(), TopAbs_WIRE, 0 ) != nbWire )
-          RETURN_BAD_RESULT("Different not quad faces");
-      }
-      notQuadFaces.push_back( face );
-    }
-  }
-  if ( !notQuadFaces.empty() )
-  {
-    if ( notQuadFaces.size() != 2 )
-      RETURN_BAD_RESULT("Bad nb not quad faces: " << notQuadFaces.size());
-
-    // check total nb faces
-    nbEdge = NSProjUtils::Count( notQuadFaces.back(), TopAbs_EDGE, 0 );
-    if ( nbFace != nbEdge + 2 )
-      RETURN_BAD_RESULT("Bad nb of faces: " << nbFace << " but must be " << nbEdge + 2);
-  }
-*/
   // no hypothesis
   aStatus = SMESH_Hypothesis::HYP_OK;
   return true;
@@ -628,6 +617,7 @@ bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theSh
 {
   SMESH_MesherHelper helper( theMesh );
   myHelper = &helper;
+  myPrevBottomSM = 0;
 
   int nbSolids = helper.Count( theShape, TopAbs_SOLID, /*skipSame=*/false );
   if ( nbSolids < 1 )
@@ -944,7 +934,7 @@ bool StdMeshers_Prism_3D::getWallFaces( Prism_3D::TPrismTopo & thePrism,
   list< TopoDS_Edge >::iterator edge = thePrism.myBottomEdges.begin();
   std::list< int >::iterator     nbE = thePrism.myNbEdgesInWires.begin();
   std::list< int > nbQuadsPerWire;
-  int iE = 0;
+  int iE = 0, iWire = 0;
   while ( edge != thePrism.myBottomEdges.end() )
   {
     ++iE;
@@ -978,7 +968,10 @@ bool StdMeshers_Prism_3D::getWallFaces( Prism_3D::TPrismTopo & thePrism,
                 return toSM( error(TCom("Composite 'horizontal' edges are not supported")));
           }
           if ( faceMap.Add( face ))
+          {
+            setWireIndex( quadList.back(), iWire ); // for use in makeQuadsForOutInProjection()
             thePrism.myWallQuads.push_back( quadList );
+          }
           break;
         }
       }
@@ -996,6 +989,7 @@ bool StdMeshers_Prism_3D::getWallFaces( Prism_3D::TPrismTopo & thePrism,
     if ( iE == *nbE )
     {
       iE = 0;
+      ++iWire;
       ++nbE;
       int nbQuadPrev = std::accumulate( nbQuadsPerWire.begin(), nbQuadsPerWire.end(), 0 );
       nbQuadsPerWire.push_back( thePrism.myWallQuads.size() - nbQuadPrev );
@@ -1130,13 +1124,8 @@ bool StdMeshers_Prism_3D::compute(const Prism_3D::TPrismTopo& thePrism)
     return toSM( error( SMESH_ComputeError::New(COMPERR_CANCELED)));
 
   // Assure the bottom is meshed
-  SMESH_subMesh * botSM = myHelper->GetMesh()->GetSubMesh( thePrism.myBottom );
-  if (( botSM->IsEmpty() ) &&
-      ( ! botSM->GetAlgo() ||
-        ! _gen->Compute( *botSM->GetFather(), botSM->GetSubShape(), /*shapeOnly=*/true )))
-    return error( COMPERR_BAD_INPUT_MESH,
-                  TCom( "No mesher defined to compute the base face #")
-                  << shapeID( thePrism.myBottom ));
+  if ( !computeBase( thePrism ))
+    return false;
 
   // Make all side FACEs of thePrism meshed with quads
   if ( !computeWalls( thePrism ))
@@ -1161,7 +1150,8 @@ bool StdMeshers_Prism_3D::compute(const Prism_3D::TPrismTopo& thePrism)
   // else if ( !trsf.empty() )
   //   bottomToTopTrsf = trsf.back();
 
-  // To compute coordinates of a node inside a block, it is necessary to know
+  // To compute coordinates of a node inside a block using "block approach",
+  // it is necessary to know
   // 1. normalized parameters of the node by which
   // 2. coordinates of node projections on all block sub-shapes are computed
 
@@ -1186,6 +1176,9 @@ bool StdMeshers_Prism_3D::compute(const Prism_3D::TPrismTopo& thePrism)
   {
     // use transformation (issue 0020680, IPAL0052499) or a "straight line" approach
     StdMeshers_Sweeper sweeper;
+    sweeper.myHelper  = myHelper;
+    sweeper.myBotFace = thePrism.myBottom;
+    sweeper.myTopFace = thePrism.myTop;
 
     // load boundary nodes into sweeper
     bool dummy;
@@ -1222,15 +1215,15 @@ bool StdMeshers_Prism_3D::compute(const Prism_3D::TPrismTopo& thePrism)
     {
       double tol = getSweepTolerance( thePrism );
       bool allowHighBndError = !isSimpleBottom( thePrism );
-      myUseBlock = !sweeper.ComputeNodes( *myHelper, tol, allowHighBndError );
+      myUseBlock = !sweeper.ComputeNodesByTrsf( tol, allowHighBndError );
     }
     else if ( sweeper.CheckSameZ() )
     {
-      myUseBlock = !sweeper.ComputeNodesOnStraightSameZ( *myHelper );
+      myUseBlock = !sweeper.ComputeNodesOnStraightSameZ();
     }
     else
     {
-      myUseBlock = !sweeper.ComputeNodesOnStraight( *myHelper, thePrism.myBottom, thePrism.myTop );
+      myUseBlock = !sweeper.ComputeNodesOnStraight();
     }
     myHelper->SetElementsOnShape( false );
   }
@@ -1384,6 +1377,75 @@ bool StdMeshers_Prism_3D::compute(const Prism_3D::TPrismTopo& thePrism)
   return true;
 }
 
+//=======================================================================
+//function : computeBase
+//purpose  : Compute the base face of a prism
+//=======================================================================
+
+bool StdMeshers_Prism_3D::computeBase(const Prism_3D::TPrismTopo& thePrism)
+{
+  SMESH_Mesh*     mesh = myHelper->GetMesh();
+  SMESH_subMesh* botSM = mesh->GetSubMesh( thePrism.myBottom );
+  if (( botSM->IsEmpty() ) &&
+      ( ! botSM->GetAlgo() ||
+        ! _gen->Compute( *botSM->GetFather(), botSM->GetSubShape(), /*shapeOnly=*/true )))
+  {
+    // find any applicable algorithm assigned to any FACE of the main shape
+    std::vector< TopoDS_Shape > faces;
+    if ( myPrevBottomSM &&
+         myPrevBottomSM->GetAlgo()->IsApplicableToShape( thePrism.myBottom, /*all=*/false ))
+      faces.push_back( myPrevBottomSM->GetSubShape() );
+
+    TopExp_Explorer faceIt( mesh->GetShapeToMesh(), TopAbs_FACE );
+    for ( ; faceIt.More(); faceIt.Next() )
+      faces.push_back( faceIt.Current() );
+
+    faces.push_back( TopoDS_Shape() ); // to try quadrangle algorithm
+
+    SMESH_Algo* algo = 0;
+    for ( size_t i = 0; i < faces.size() &&  botSM->IsEmpty(); ++i )
+    {
+      if ( faces[i].IsNull() ) algo = TQuadrangleAlgo::instance( this, myHelper );
+      else                     algo = mesh->GetSubMesh( faces[i] )->GetAlgo();
+      if ( algo && algo->IsApplicableToShape( thePrism.myBottom, /*all=*/false ))
+      {
+        // try to compute the bottom FACE
+        if ( algo->NeedDiscreteBoundary() )
+        {
+          // compute sub-shapes
+          SMESH_subMeshIteratorPtr smIt = botSM->getDependsOnIterator(false,false);
+          bool subOK = true;
+          while ( smIt->more() && subOK )
+          {
+            SMESH_subMesh* sub = smIt->next();
+            sub->ComputeStateEngine( SMESH_subMesh::COMPUTE );
+            subOK = sub->IsMeshComputed();
+          }
+          if ( !subOK )
+            continue;
+        }
+        try {
+          OCC_CATCH_SIGNALS;
+          algo->InitComputeError();
+          algo->Compute( *mesh, botSM->GetSubShape() );
+        }
+        catch (...) {
+        }
+      }
+    }
+  }
+
+  if ( botSM->IsEmpty() )
+    return error( COMPERR_BAD_INPUT_MESH,
+                  TCom( "No mesher defined to compute the base face #")
+                  << shapeID( thePrism.myBottom ));
+
+  if ( botSM->GetAlgo() )
+    myPrevBottomSM = botSM;
+
+  return true;
+}
+
 //=======================================================================
 //function : computeWalls
 //purpose  : Compute 2D mesh on walls FACEs of a prism
@@ -1414,6 +1476,7 @@ bool StdMeshers_Prism_3D::computeWalls(const Prism_3D::TPrismTopo& thePrism)
     for ( ; quad != thePrism.myWallQuads[iW].end(); ++quad )
     {
       StdMeshers_FaceSidePtr lftSide = (*quad)->side[ QUAD_LEFT_SIDE ];
+      lftSide->Reverse(); // to go up
       for ( int i = 0; i < lftSide->NbEdges(); ++i )
       {
         ++wgt[ iW ];
@@ -1441,6 +1504,11 @@ bool StdMeshers_Prism_3D::computeWalls(const Prism_3D::TPrismTopo& thePrism)
   for ( size_t iW = 0; iW != nbWalls; ++iW )
     wgt2quad.insert( make_pair( wgt[ iW ], iW ));
 
+  // artificial quads to do outer <-> inner wall projection
+  std::map< int, FaceQuadStruct > iW2oiQuads;
+  std::map< int, FaceQuadStruct >::iterator w2oiq;
+  makeQuadsForOutInProjection( thePrism, wgt2quad, iW2oiQuads );
+
   // Project 'vertical' EDGEs, from left to right
   multimap< int, int >::reverse_iterator w2q = wgt2quad.rbegin();
   for ( ; w2q != wgt2quad.rend(); ++w2q )
@@ -1457,10 +1525,25 @@ bool StdMeshers_Prism_3D::computeWalls(const Prism_3D::TPrismTopo& thePrism)
       if ( swapLeftRight )
         std::swap( lftSide, rgtSide );
 
+      bool isArtificialQuad = (( w2oiq = iW2oiQuads.find( iW )) != iW2oiQuads.end() );
+      if ( isArtificialQuad )
+      {
+        // reset sides to perform the outer <-> inner projection
+        FaceQuadStruct& oiQuad = w2oiq->second;
+        rgtSide = oiQuad.side[ QUAD_RIGHT_SIDE ];
+        lftSide = oiQuad.side[ QUAD_LEFT_SIDE ];
+        iW2oiQuads.erase( w2oiq );
+      }
+
       // assure that all the source (left) EDGEs are meshed
       int nbSrcSegments = 0;
       for ( int i = 0; i < lftSide->NbEdges(); ++i )
       {
+        if ( isArtificialQuad )
+        {
+          nbSrcSegments = lftSide->NbPoints()-1;
+          continue;
+        }
         const TopoDS_Edge& srcE = lftSide->Edge(i);
         SMESH_subMesh*    srcSM = mesh->GetSubMesh( srcE );
         if ( !srcSM->IsMeshComputed() ) {
@@ -1536,7 +1619,7 @@ bool StdMeshers_Prism_3D::computeWalls(const Prism_3D::TPrismTopo& thePrism)
         const UVPtStructVec&  srcNodeStr = lftSide->GetUVPtStruct();
         if ( srcNodeStr.size() == 0 )
           return toSM( error( TCom("Invalid node positions on edge #") <<
-                              shapeID( lftSide->Edge(0) )));
+                              lftSide->EdgeID(0) ));
         vector< SMDS_MeshNode* > newNodes( srcNodeStr.size() );
         for ( int is2ndV = 0; is2ndV < 2; ++is2ndV )
         {
@@ -1549,7 +1632,7 @@ bool StdMeshers_Prism_3D::computeWalls(const Prism_3D::TPrismTopo& thePrism)
 
         // compute nodes on target EDGEs
         DBGOUT( "COMPUTE V edge (proj) " << shapeID( lftSide->Edge(0)));
-        rgtSide->Reverse(); // direct it same as the lftSide
+        //rgtSide->Reverse(); // direct it same as the lftSide
         myHelper->SetElementsOnShape( false ); // myHelper holds the prism shape
         TopoDS_Edge tgtEdge;
         for ( size_t iN = 1; iN < srcNodeStr.size()-1; ++iN ) // add nodes
@@ -1718,9 +1801,8 @@ bool StdMeshers_Prism_3D::computeWalls(const Prism_3D::TPrismTopo& thePrism)
 }
 
 //=======================================================================
-/*!
- * \brief Returns a source EDGE of propagation to a given EDGE
- */
+//function : findPropagationSource
+//purpose  : Returns a source EDGE of propagation to a given EDGE
 //=======================================================================
 
 TopoDS_Edge StdMeshers_Prism_3D::findPropagationSource( const TopoDS_Edge& E )
@@ -1733,9 +1815,93 @@ TopoDS_Edge StdMeshers_Prism_3D::findPropagationSource( const TopoDS_Edge& E )
   return TopoDS_Edge();
 }
 
+//=======================================================================
+//function : makeQuadsForOutInProjection
+//purpose  : Create artificial wall quads for vertical projection between
+//           the outer and inner walls
+//=======================================================================
+
+void StdMeshers_Prism_3D::makeQuadsForOutInProjection( const Prism_3D::TPrismTopo& thePrism,
+                                                       multimap< int, int >&       wgt2quad,
+                                                       map< int, FaceQuadStruct >& iQ2oiQuads)
+{
+  if ( thePrism.NbWires() <= 1 )
+    return;
+
+  std::set< int > doneWires; // processed wires
+
+  SMESH_Mesh*      mesh = myHelper->GetMesh();
+  const bool  isForward = true;
+  const bool skipMedium = myHelper->GetIsQuadratic();
+
+  // make a source side for all projections
+
+  multimap< int, int >::reverse_iterator w2q = wgt2quad.rbegin();
+  const int iQuad = w2q->second;
+  const int iWire = getWireIndex( thePrism.myWallQuads[ iQuad ].front() );
+  doneWires.insert( iWire );
+
+  UVPtStructVec srcNodes;
+
+  Prism_3D::TQuadList::const_iterator quad = thePrism.myWallQuads[ iQuad ].begin();
+  for ( ; quad != thePrism.myWallQuads[ iQuad ].end(); ++quad )
+  {
+    StdMeshers_FaceSidePtr lftSide = (*quad)->side[ QUAD_LEFT_SIDE ];
+
+    // assure that all the source (left) EDGEs are meshed
+    for ( int i = 0; i < lftSide->NbEdges(); ++i )
+    {
+      const TopoDS_Edge& srcE = lftSide->Edge(i);
+      SMESH_subMesh*    srcSM = mesh->GetSubMesh( srcE );
+      if ( !srcSM->IsMeshComputed() ) {
+        srcSM->ComputeSubMeshStateEngine( SMESH_subMesh::COMPUTE );
+        srcSM->ComputeStateEngine       ( SMESH_subMesh::COMPUTE );
+      }
+      if ( !srcSM->IsMeshComputed() )
+        return;
+    }
+    const UVPtStructVec& subNodes = lftSide->GetUVPtStruct();
+    UVPtStructVec::const_iterator subBeg = subNodes.begin(), subEnd = subNodes.end();
+    if ( !srcNodes.empty() ) ++subBeg;
+    srcNodes.insert( srcNodes.end(), subBeg, subEnd );
+  }
+  StdMeshers_FaceSidePtr srcSide = StdMeshers_FaceSide::New( srcNodes );
+
+  // make the quads
+
+  list< TopoDS_Edge > sideEdges;
+  TopoDS_Face face;
+  for ( ++w2q; w2q != wgt2quad.rend(); ++w2q )
+  {
+    const int                  iQuad = w2q->second;
+    const Prism_3D::TQuadList& quads = thePrism.myWallQuads[ iQuad ];
+    const int                  iWire = getWireIndex( quads.front() );
+    if ( !doneWires.insert( iWire ).second )
+      continue;
+
+    sideEdges.clear();
+    for ( quad = quads.begin(); quad != quads.end(); ++quad )
+    {
+      StdMeshers_FaceSidePtr lftSide = (*quad)->side[ QUAD_LEFT_SIDE ];
+      for ( int i = 0; i < lftSide->NbEdges(); ++i )
+        sideEdges.push_back( lftSide->Edge( i ));
+      face = lftSide->Face();
+    }
+    StdMeshers_FaceSidePtr tgtSide =
+      StdMeshers_FaceSide::New( face, sideEdges, mesh, isForward, skipMedium, myHelper );
+
+    FaceQuadStruct& newQuad = iQ2oiQuads[ iQuad ];
+    newQuad.side.resize( 4 );
+    newQuad.side[ QUAD_LEFT_SIDE  ] = srcSide;
+    newQuad.side[ QUAD_RIGHT_SIDE ] = tgtSide;
+
+    wgt2quad.insert( *w2q ); // to process this quad after processing the newQuad
+  }
+}
+
 //=======================================================================
 //function : Evaluate
-//purpose  : 
+//purpose  :
 //=======================================================================
 
 bool StdMeshers_Prism_3D::Evaluate(SMESH_Mesh&         theMesh,
@@ -2297,7 +2463,7 @@ bool StdMeshers_Prism_3D::projectBottomToTop( const gp_Trsf &             bottom
 
   // Check the projected mesh
 
-  if ( thePrism.myNbEdgesInWires.size() > 1 && // there are holes
+  if ( thePrism.NbWires() > 1 && // there are holes
        topHelper.IsDistorted2D( topSM, /*checkUV=*/false ))
   {
     SMESH_MeshEditor editor( topHelper.GetMesh() );
@@ -4792,13 +4958,13 @@ void StdMeshers_Sweeper::applyBoundaryError(const vector< gp_XYZ >& bndPoints,
 
 //================================================================================
 /*!
- * \brief Create internal nodes of the prism
+ * \brief Create internal nodes of the prism by computing an affine transformation
+ *        from layer to layer
  */
 //================================================================================
 
-bool StdMeshers_Sweeper::ComputeNodes( SMESH_MesherHelper& helper,
-                                       const double        tol,
-                                       const bool          allowHighBndError)
+bool StdMeshers_Sweeper::ComputeNodesByTrsf( const double tol,
+                                             const bool   allowHighBndError)
 {
   const size_t zSize = myBndColumns[0]->size();
   const size_t zSrc = 0, zTgt = zSize-1;
@@ -5065,7 +5231,7 @@ bool StdMeshers_Sweeper::ComputeNodes( SMESH_MesherHelper& helper,
     for ( size_t z = zSrc + 1; z < zTgt; ++z ) // vertical loop on layers
     {
       const gp_XYZ & xyz = intPntsOfLayer[ z ][ iP ];
-      if ( !( nodeCol[ z ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() )))
+      if ( !( nodeCol[ z ] = myHelper->AddNode( xyz.X(), xyz.Y(), xyz.Z() )))
         return false;
     }
   }
@@ -5137,7 +5303,7 @@ bool StdMeshers_Sweeper::CheckSameZ()
  */
 //================================================================================
 
-bool StdMeshers_Sweeper::ComputeNodesOnStraightSameZ( SMESH_MesherHelper& helper )
+bool StdMeshers_Sweeper::ComputeNodesOnStraightSameZ()
 {
   TZColumn& z = myZColumns[0];
 
@@ -5149,7 +5315,7 @@ bool StdMeshers_Sweeper::ComputeNodesOnStraightSameZ( SMESH_MesherHelper& helper
     for ( size_t iZ = 0; iZ < z.size(); ++iZ )
     {
       gp_XYZ p = n0 * ( 1 - z[iZ] ) + n1 * z[iZ];
-      nodes[ iZ+1 ] = helper.AddNode( p.X(), p.Y(), p.Z() );
+      nodes[ iZ+1 ] = myHelper->AddNode( p.X(), p.Y(), p.Z() );
     }
   }
 
@@ -5163,144 +5329,51 @@ bool StdMeshers_Sweeper::ComputeNodesOnStraightSameZ( SMESH_MesherHelper& helper
  */
 //================================================================================
 
-bool StdMeshers_Sweeper::ComputeNodesOnStraight( SMESH_MesherHelper& helper,
-                                                 const TopoDS_Face&  botFace,
-                                                 const TopoDS_Face&  topFace )
+bool StdMeshers_Sweeper::ComputeNodesOnStraight()
 {
-  // get data to create a Morph
-  UVPtStructVec botUV( myBndColumns.size() + 1 );
-  UVPtStructVec topUV( myBndColumns.size() + 1 );
-  for ( size_t i = 0; i < myBndColumns.size(); ++i )
-  {
-    TNodeColumn& nodes = *myBndColumns[i];
-    botUV[i].node = nodes[0];
-    botUV[i].SetUV( helper.GetNodeUV( botFace, nodes[0] ));
-    topUV[i].node = nodes.back();
-    topUV[i].SetUV( helper.GetNodeUV( topFace, nodes.back() ));
-    botUV[i].node->setIsMarked( true );
-  }
-  botUV.back() = botUV[0];
-  topUV.back() = topUV[0];
-
-  TopoDS_Edge dummyE;
-  TSideVector botWires( 1, StdMeshers_FaceSide::New( botUV, botFace, dummyE, helper.GetMesh() ));
-  TSideVector topWires( 1, StdMeshers_FaceSide::New( topUV, topFace, dummyE, helper.GetMesh() ));
-
-  // use Morph to make delauney mesh on the FACEs. Locating of a node within a
-  // delauney triangle will be used to get a weighted Z.
-  NSProjUtils::Morph botDelauney( botWires );
-  NSProjUtils::Morph topDelauney( topWires );
-
-  if ( helper.GetIsQuadratic() )
-  {
-    // mark all medium nodes of faces on botFace to avoid their treating
-    SMESHDS_SubMesh* smDS = helper.GetMeshDS()->MeshElements( botFace );
-    SMDS_ElemIteratorPtr eIt = smDS->GetElements();
-    while ( eIt->more() )
-    {
-      const SMDS_MeshElement* e = eIt->next();
-      for ( int i = e->NbCornerNodes(), nb = e->NbNodes(); i < nb; ++i )
-        e->GetNode( i )->setIsMarked( true );
-    }
-  }
-
-  // map to get a node column by a bottom node
-  TColStd_DataMapOfIntegerInteger iNode2iCol( myIntColumns.size() );
-
-  // un-mark nodes to treat (internal bottom nodes); later we will mark treated nodes
-  for ( size_t i = 0; i < myIntColumns.size(); ++i )
-  {
-    const SMDS_MeshNode* botNode = myIntColumns[i]->front();
-    botNode->setIsMarked( false );
-    iNode2iCol.Bind( botNode->GetID(), i );
-  }
+  prepareTopBotDelaunay();
 
-  const int botFaceID = helper.GetMesh()->GetSubMesh( botFace )->GetId();
   const SMDS_MeshNode     *botNode, *topNode;
-  const BRepMesh_Triangle *botTria, *topTria;
+  const BRepMesh_Triangle *topTria;
   double botBC[3], topBC[3]; // barycentric coordinates
   int    botTriaNodes[3], topTriaNodes[3];
   bool   checkUV = true;
 
-  // a queue of bottom nodes with starting delauney triangles
-  NSProjUtils::Morph::TNodeTriaList botNoTriQueue;
+  int nbInternalNodes = myIntColumns.size();
+  myBotDelaunay->InitTraversal( nbInternalNodes );
 
-  size_t iBndN = 1; // index of a bottom boundary node
-  int nbNodesToProcess = myIntColumns.size();
-  while ( nbNodesToProcess > 0 )
+  while (( botNode = myBotDelaunay->NextNode( botBC, botTriaNodes )))
   {
-    while ( !botNoTriQueue.empty() ) // treat all nodes in the queue
-    {
-      botNode = botNoTriQueue.front().first;
-      botTria = botNoTriQueue.front().second;
-      botNoTriQueue.pop_front();
-      if ( botNode->isMarked() )
-        continue;
-      --nbNodesToProcess;
-      botNode->setIsMarked( true );
-
-      TNodeColumn* column = myIntColumns[ iNode2iCol( botNode->GetID() )];
-
-      // find a delauney triangle containing the botNode
-      gp_XY botUV = helper.GetNodeUV( botFace, botNode, NULL, &checkUV );
-      botUV  *= botDelauney.GetScale();
-      botTria = botDelauney.FindTriangle( botUV, botTria, botBC, botTriaNodes );
-      if ( !botTria )
-        return false;
-
-      // find a delauney triangle containing the topNode
-      topNode = column->back();
-      gp_XY topUV = helper.GetNodeUV( topFace, topNode, NULL, &checkUV );
-      topUV *= topDelauney.GetScale();
-      // get a starting triangle basing on that top and bot boundary nodes have same index
-      topTria = topDelauney.GetTriangleNear( botTriaNodes[0] );
-      topTria = topDelauney.FindTriangle( topUV, topTria, topBC, topTriaNodes );
-      if ( !topTria )
-        return false;
-
-      // create nodes along a line
-      SMESH_NodeXYZ botP( botNode ), topP( topNode);
-      for ( size_t iZ = 0; iZ < myZColumns[0].size(); ++iZ )
-      {
-        // use barycentric coordinates as weight of Z of boundary columns
-        double botZ = 0, topZ = 0;
-        for ( int i = 0; i < 3; ++i )
-        {
-          botZ += botBC[i] * myZColumns[ botTriaNodes[i]-1 ][ iZ ];
-          topZ += topBC[i] * myZColumns[ topTriaNodes[i]-1 ][ iZ ];
-        }
-        double rZ = double( iZ + 1 ) / ( myZColumns[0].size() + 1 );
-        double z = botZ * ( 1 - rZ ) + topZ * rZ;
-        gp_XYZ p = botP * ( 1 - z  ) + topP * z;
-        (*column)[ iZ+1 ] = helper.AddNode( p.X(), p.Y(), p.Z() );
-      }
-
-      // add neighbor nodes to the queue
-      botDelauney.AddCloseNodes( botNode, botTria, botFaceID, botNoTriQueue );
-    }
+    TNodeColumn* column = myIntColumns[ myNodeID2ColID( botNode->GetID() )];
+
+    // find a Delaunay triangle containing the topNode
+    topNode = column->back();
+    gp_XY topUV = myHelper->GetNodeUV( myTopFace, topNode, NULL, &checkUV );
+    // get a starting triangle basing on that top and bot boundary nodes have same index
+    topTria = myTopDelaunay->GetTriangleNear( botTriaNodes[0] );
+    topTria = myTopDelaunay->FindTriangle( topUV, topTria, topBC, topTriaNodes );
+    if ( !topTria )
+      return false;
 
-    if ( nbNodesToProcess > 0 ) // fill the queue
+    // create nodes along a line
+    SMESH_NodeXYZ botP( botNode ), topP( topNode);
+    for ( size_t iZ = 0; iZ < myZColumns[0].size(); ++iZ )
     {
-      // assure that all bot nodes are visited
-      for ( ; iBndN-1 < myBndColumns.size() &&  botNoTriQueue.empty();  ++iBndN )
+      // use barycentric coordinates as weight of Z of boundary columns
+      double botZ = 0, topZ = 0;
+      for ( int i = 0; i < 3; ++i )
       {
-        botTria = botDelauney.GetTriangleNear( iBndN );
-        const SMDS_MeshNode*  bndNode = botDelauney.GetBndNodes()[ iBndN ];
-        botDelauney.AddCloseNodes( bndNode, botTria, botFaceID, botNoTriQueue );
-      }
-      if ( botNoTriQueue.empty() )
-      {
-        for ( size_t i = 0; i < myIntColumns.size(); ++i )
-        {
-          botNode = myIntColumns[i]->front();
-          if ( !botNode->isMarked() )
-            botNoTriQueue.push_back( make_pair( botNode, botTria ));
-        }
+        botZ += botBC[i] * myZColumns[ botTriaNodes[i] ][ iZ ];
+        topZ += topBC[i] * myZColumns[ topTriaNodes[i] ][ iZ ];
       }
+      double rZ = double( iZ + 1 ) / ( myZColumns[0].size() + 1 );
+      double z = botZ * ( 1 - rZ ) + topZ * rZ;
+      gp_XYZ p = botP * ( 1 - z  ) + topP * z;
+      (*column)[ iZ+1 ] = myHelper->AddNode( p.X(), p.Y(), p.Z() );
     }
   }
 
-  return true;
+  return myBotDelaunay->NbVisitedNodes() == nbInternalNodes;
 }
 
 //================================================================================
@@ -5326,3 +5399,58 @@ void StdMeshers_Sweeper::fillZColumn( TZColumn&    zColumn,
     zColumn[i] = ( line * vec ) / len2; // param [0,1] on the line
   }
 }
+
+//================================================================================
+/*!
+ * \brief Initialize *Delaunay members
+ */
+//================================================================================
+
+void StdMeshers_Sweeper::prepareTopBotDelaunay()
+{
+  UVPtStructVec botUV( myBndColumns.size() );
+  UVPtStructVec topUV( myBndColumns.size() );
+  for ( size_t i = 0; i < myBndColumns.size(); ++i )
+  {
+    TNodeColumn& nodes = *myBndColumns[i];
+    botUV[i].node = nodes[0];
+    botUV[i].SetUV( myHelper->GetNodeUV( myBotFace, nodes[0] ));
+    topUV[i].node = nodes.back();
+    topUV[i].SetUV( myHelper->GetNodeUV( myTopFace, nodes.back() ));
+    botUV[i].node->setIsMarked( true );
+  }
+  TopoDS_Edge dummyE;
+  SMESH_Mesh* mesh = myHelper->GetMesh();
+  TSideVector botWires( 1, StdMeshers_FaceSide::New( botUV, myBotFace, dummyE, mesh ));
+  TSideVector topWires( 1, StdMeshers_FaceSide::New( topUV, myTopFace, dummyE, mesh ));
+
+  // Delaunay mesh on the FACEs.
+  bool checkUV = false;
+  myBotDelaunay.reset( new NSProjUtils::Delaunay( botWires, checkUV ));
+  myTopDelaunay.reset( new NSProjUtils::Delaunay( topWires, checkUV ));
+
+  if ( myHelper->GetIsQuadratic() )
+  {
+    // mark all medium nodes of faces on botFace to avoid their treating
+    SMESHDS_SubMesh* smDS = myHelper->GetMeshDS()->MeshElements( myBotFace );
+    SMDS_ElemIteratorPtr eIt = smDS->GetElements();
+    while ( eIt->more() )
+    {
+      const SMDS_MeshElement* e = eIt->next();
+      for ( int i = e->NbCornerNodes(), nb = e->NbNodes(); i < nb; ++i )
+        e->GetNode( i )->setIsMarked( true );
+    }
+  }
+
+  // map to get a node column by a bottom node
+  myNodeID2ColID.Clear(/*doReleaseMemory=*/false);
+  myNodeID2ColID.ReSize( myIntColumns.size() );
+
+  // un-mark nodes to treat (internal bottom nodes) to be returned by myBotDelaunay
+  for ( size_t i = 0; i < myIntColumns.size(); ++i )
+  {
+    const SMDS_MeshNode* botNode = myIntColumns[i]->front();
+    botNode->setIsMarked( false );
+    myNodeID2ColID.Bind( botNode->GetID(), i );
+  }
+}