Salome HOME
Fix MinDistance for node-group (SALOME_TESTS/Grids/smesh/imps_09/K0)
[modules/smesh.git] / src / SMESHUtils / SMESH_MeshAlgos.cxx
index 8464e245fe197e6765f100a08498d730211522bb..a2c9c10368fe357f5ee0c6011b6ac255774f8ff5 100644 (file)
@@ -1,4 +1,4 @@
-// Copyright (C) 2007-2016  CEA/DEN, EDF R&D, OPEN CASCADE
+// Copyright (C) 2007-2019  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
 // Created   : Tue Apr 30 18:00:36 2013
 // Author    : Edward AGAPOV (eap)
 
-// This file holds some low level algorithms extracted from SMESH_MeshEditor
+// Initially this file held some low level algorithms extracted from SMESH_MeshEditor
 // to make them accessible from Controls package
 
 #include "SMESH_MeshAlgos.hxx"
 
+#include "ObjectPool.hxx"
 #include "SMDS_FaceOfNodes.hxx"
 #include "SMDS_LinearEdge.hxx"
 #include "SMDS_Mesh.hxx"
@@ -35,6 +36,8 @@
 #include "SMDS_VolumeTool.hxx"
 #include "SMESH_OctreeNode.hxx"
 
+#include <Utils_SALOME_Exception.hxx>
+
 #include <GC_MakeSegment.hxx>
 #include <GeomAPI_ExtremaCurveCurve.hxx>
 #include <Geom_Line.hxx>
 #include <IntAna_Quadric.hxx>
 #include <gp_Lin.hxx>
 #include <gp_Pln.hxx>
+#include <NCollection_DataMap.hxx>
 
 #include <limits>
 #include <numeric>
 
-using namespace std;
+#include <boost/container/flat_set.hpp>
 
 //=======================================================================
 /*!
@@ -67,7 +71,7 @@ struct SMESH_NodeSearcherImpl: public SMESH_NodeSearcher
 
     TIDSortedNodeSet nodes;
     if ( theMesh ) {
-      SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator(/*idInceasingOrder=*/true);
+      SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator();
       while ( nIt->more() )
         nodes.insert( nodes.end(), nIt->next() );
     }
@@ -108,21 +112,22 @@ struct SMESH_NodeSearcherImpl: public SMESH_NodeSearcher
    */
   const SMDS_MeshNode* FindClosestTo( const gp_Pnt& thePnt )
   {
-    map<double, const SMDS_MeshNode*> dist2Nodes;
+    std::map<double, const SMDS_MeshNode*> dist2Nodes;
     myOctreeNode->NodesAround( thePnt.Coord(), dist2Nodes, myHalfLeafSize );
     if ( !dist2Nodes.empty() )
       return dist2Nodes.begin()->second;
-    list<const SMDS_MeshNode*> nodes;
+
+    std::vector<const SMDS_MeshNode*> nodes;
     //myOctreeNode->NodesAround( &tgtNode, &nodes, myHalfLeafSize );
 
     double minSqDist = DBL_MAX;
     if ( nodes.empty() )  // get all nodes of OctreeNode's closest to thePnt
     {
       // sort leafs by their distance from thePnt
-      typedef map< double, SMESH_OctreeNode* > TDistTreeMap;
+      typedef std::multimap< double, SMESH_OctreeNode* > TDistTreeMap;
       TDistTreeMap treeMap;
-      list< SMESH_OctreeNode* > treeList;
-      list< SMESH_OctreeNode* >::iterator trIt;
+      std::list< SMESH_OctreeNode* > treeList;
+      std::list< SMESH_OctreeNode* >::iterator trIt;
       treeList.push_back( myOctreeNode );
 
       gp_XYZ pointNode( thePnt.X(), thePnt.Y(), thePnt.Z() );
@@ -141,9 +146,7 @@ struct SMESH_NodeSearcherImpl: public SMESH_NodeSearcher
         {
           const Bnd_B3d& box = *tree->getBox();
           double sqDist = thePnt.SquareDistance( 0.5 * ( box.CornerMin() + box.CornerMax() ));
-          pair<TDistTreeMap::iterator,bool> it_in = treeMap.insert( make_pair( sqDist, tree ));
-          if ( !it_in.second ) // not unique distance to box center
-            treeMap.insert( it_in.first, make_pair( sqDist + 1e-13*treeMap.size(), tree ));
+          treeMap.insert( std::make_pair( sqDist, tree ));
         }
       }
       // find distance after which there is no sense to check tree's
@@ -160,17 +163,17 @@ struct SMESH_NodeSearcherImpl: public SMESH_NodeSearcher
         if ( sqDist_tree->first > sqLimit )
           break;
         SMESH_OctreeNode* tree = sqDist_tree->second;
-        tree->NodesAround( tree->GetNodeIterator()->next(), &nodes );
+        tree->AllNodesAround( tree->GetNodeIterator()->next(), &nodes );
       }
     }
     // find closest among nodes
     minSqDist = DBL_MAX;
     const SMDS_MeshNode* closestNode = 0;
-    list<const SMDS_MeshNode*>::iterator nIt = nodes.begin();
-    for ( ; nIt != nodes.end(); ++nIt ) {
-      double sqDist = thePnt.SquareDistance( SMESH_TNodeXYZ( *nIt ) );
+    for ( size_t i = 0; i < nodes.size(); ++i )
+    {
+      double sqDist = thePnt.SquareDistance( SMESH_NodeXYZ( nodes[ i ]));
       if ( minSqDist > sqDist ) {
-        closestNode = *nIt;
+        closestNode = nodes[ i ];
         minSqDist = sqDist;
       }
     }
@@ -179,7 +182,7 @@ struct SMESH_NodeSearcherImpl: public SMESH_NodeSearcher
 
   //---------------------------------------------------------------------
   /*!
-   * \brief Finds nodes located within a tolerance near a point 
+   * \brief Finds nodes located within a tolerance near a point
    */
   int FindNearPoint(const gp_Pnt&                        point,
                     const double                         tolerance,
@@ -224,16 +227,18 @@ namespace // Utils used in SMESH_ElementSearcherImpl::FindElementsByPoint()
   {
   public:
 
+    typedef boost::container::flat_set< const SMDS_MeshElement*, TIDCompare > TElemSeq;
+
     ElementBndBoxTree(const SMDS_Mesh&     mesh,
                       SMDSAbs_ElementType  elemType,
                       SMDS_ElemIteratorPtr theElemIt = SMDS_ElemIteratorPtr(),
                       double               tolerance = NodeRadius );
-    void prepare(); // !!!call it before calling the following methods!!!
-    void getElementsNearPoint( const gp_Pnt& point, vector<const SMDS_MeshElement*>& foundElems );
-    void getElementsNearLine ( const gp_Ax1& line,  vector<const SMDS_MeshElement*>& foundElems);
-    void getElementsInSphere ( const gp_XYZ&                    center,
-                               const double                     radius,
-                               vector<const SMDS_MeshElement*>& foundElems);
+    void getElementsNearPoint( const gp_Pnt& point, TElemSeq& foundElems );
+    void getElementsNearLine ( const gp_Ax1& line,  TElemSeq& foundElems );
+    void getElementsInBox    ( const Bnd_B3d& box,  TElemSeq& foundElems );
+    void getElementsInSphere ( const gp_XYZ& center, const double radius, TElemSeq& foundElems );
+    ElementBndBoxTree* getLeafAtPoint( const gp_XYZ& point );
+    int  getNbElements();
 
   protected:
     ElementBndBoxTree() {}
@@ -245,18 +250,16 @@ namespace // Utils used in SMESH_ElementSearcherImpl::FindElementsByPoint()
     struct ElementBox : public Bnd_B3d
     {
       const SMDS_MeshElement* _element;
-      bool                    _isMarked;
       void init(const SMDS_MeshElement* elem, double tolerance);
     };
-    vector< ElementBox* > _elements;
+    std::vector< ElementBox* > _elements;
 
     typedef ObjectPool< ElementBox > TElementBoxPool;
 
     //!< allocator of ElementBox's and SMESH_TreeLimit
     struct LimitAndPool : public SMESH_TreeLimit
     {
-      TElementBoxPool            _elBoPool;
-      std::vector< ElementBox* > _markedElems;
+      TElementBoxPool _elBoPool;
       LimitAndPool():SMESH_TreeLimit( MaxLevel, /*minSize=*/0. ) {}
     };
     LimitAndPool* getLimitAndPool() const
@@ -283,6 +286,11 @@ namespace // Utils used in SMESH_ElementSearcherImpl::FindElementsByPoint()
 
     TElementBoxPool& elBoPool = getLimitAndPool()->_elBoPool;
 
+#ifdef _DEBUG_
+    if ( theElemIt && !theElemIt->more() )
+      std::cout << "WARNING: ElementBndBoxTree constructed on empty iterator!" << std::endl;
+#endif
+
     SMDS_ElemIteratorPtr elemIt = theElemIt ? theElemIt : mesh.elementsIterator( elemType );
     while ( elemIt->more() )
     {
@@ -337,56 +345,27 @@ namespace // Utils used in SMESH_ElementSearcherImpl::FindElementsByPoint()
     }
   }
 
-  //================================================================================
-  /*!
-   * \brief Un-mark all elements
-   */
-  //================================================================================
-
-  void ElementBndBoxTree::prepare()
-  {
-    // TElementBoxPool& elBoPool = getElementBoxPool();
-    // for ( size_t i = 0; i < elBoPool.nbElements(); ++i )
-    //   const_cast< ElementBox* >( elBoPool[ i ])->_isMarked = false;
-  }
-
   //================================================================================
   /*!
    * \brief Return elements which can include the point
    */
   //================================================================================
 
-  void ElementBndBoxTree::getElementsNearPoint( const gp_Pnt&                    point,
-                                                vector<const SMDS_MeshElement*>& foundElems)
+  void ElementBndBoxTree::getElementsNearPoint( const gp_Pnt& point, TElemSeq& foundElems)
   {
     if ( getBox()->IsOut( point.XYZ() ))
       return;
 
     if ( isLeaf() )
     {
-      LimitAndPool* pool = getLimitAndPool();
-
       for ( size_t i = 0; i < _elements.size(); ++i )
-        if ( !_elements[i]->IsOut( point.XYZ() ) &&
-             !_elements[i]->_isMarked )
-        {
-          foundElems.push_back( _elements[i]->_element );
-          _elements[i]->_isMarked = true;
-          pool->_markedElems.push_back( _elements[i] );
-        }
+        if ( !_elements[i]->IsOut( point.XYZ() ))
+          foundElems.insert( _elements[i]->_element );
     }
     else
     {
       for (int i = 0; i < 8; i++)
         ((ElementBndBoxTree*) myChildren[i])->getElementsNearPoint( point, foundElems );
-
-      if ( level() == 0 )
-      {
-        LimitAndPool* pool = getLimitAndPool();
-        for ( size_t i = 0; i < pool->_markedElems.size(); ++i )
-          pool->_markedElems[i]->_isMarked = false;
-        pool->_markedElems.clear();
-      }
     }
   }
 
@@ -396,37 +375,21 @@ namespace // Utils used in SMESH_ElementSearcherImpl::FindElementsByPoint()
    */
   //================================================================================
 
-  void ElementBndBoxTree::getElementsNearLine( const gp_Ax1&                    line,
-                                               vector<const SMDS_MeshElement*>& foundElems)
+  void ElementBndBoxTree::getElementsNearLine( const gp_Ax1& line, TElemSeq& foundElems )
   {
     if ( getBox()->IsOut( line ))
       return;
 
     if ( isLeaf() )
     {
-      LimitAndPool* pool = getLimitAndPool();
-
       for ( size_t i = 0; i < _elements.size(); ++i )
-        if ( !_elements[i]->IsOut( line ) &&
-             !_elements[i]->_isMarked )
-        {
-          foundElems.push_back( _elements[i]->_element );
-          _elements[i]->_isMarked = true;
-          pool->_markedElems.push_back( _elements[i] );
-        }
+        if ( !_elements[i]->IsOut( line ) )
+          foundElems.insert( _elements[i]->_element );
     }
     else
     {
       for (int i = 0; i < 8; i++)
         ((ElementBndBoxTree*) myChildren[i])->getElementsNearLine( line, foundElems );
-
-      if ( level() == 0 )
-      {
-        LimitAndPool* pool = getLimitAndPool();
-        for ( size_t i = 0; i < pool->_markedElems.size(); ++i )
-          pool->_markedElems[i]->_isMarked = false;
-        pool->_markedElems.clear();
-      }
     }
   }
 
@@ -436,41 +399,95 @@ namespace // Utils used in SMESH_ElementSearcherImpl::FindElementsByPoint()
    */
   //================================================================================
 
-  void ElementBndBoxTree::getElementsInSphere ( const gp_XYZ&                    center,
-                                                const double                     radius,
-                                                vector<const SMDS_MeshElement*>& foundElems)
+  void ElementBndBoxTree::getElementsInSphere ( const gp_XYZ& center,
+                                                const double  radius,
+                                                TElemSeq&     foundElems)
   {
     if ( getBox()->IsOut( center, radius ))
       return;
 
     if ( isLeaf() )
     {
-      LimitAndPool* pool = getLimitAndPool();
-
       for ( size_t i = 0; i < _elements.size(); ++i )
-        if ( !_elements[i]->IsOut( center, radius ) &&
-             !_elements[i]->_isMarked )
-        {
-          foundElems.push_back( _elements[i]->_element );
-          _elements[i]->_isMarked = true;
-          pool->_markedElems.push_back( _elements[i] );
-        }
+        if ( !_elements[i]->IsOut( center, radius ))
+          foundElems.insert( _elements[i]->_element );
     }
     else
     {
       for (int i = 0; i < 8; i++)
         ((ElementBndBoxTree*) myChildren[i])->getElementsInSphere( center, radius, foundElems );
+    }
+  }
 
-      if ( level() == 0 )
-      {
-        LimitAndPool* pool = getLimitAndPool();
-        for ( size_t i = 0; i < pool->_markedElems.size(); ++i )
-          pool->_markedElems[i]->_isMarked = false;
-        pool->_markedElems.clear();
-      }
+  //================================================================================
+  /*!
+   * \brief Return elements from leaves intersecting the box
+   */
+  //================================================================================
+
+  void ElementBndBoxTree::getElementsInBox( const Bnd_B3d& box,  TElemSeq& foundElems )
+  {
+    if ( getBox()->IsOut( box ))
+      return;
+
+    if ( isLeaf() )
+    {
+      for ( size_t i = 0; i < _elements.size(); ++i )
+        if ( !_elements[i]->IsOut( box ))
+          foundElems.insert( _elements[i]->_element );
+    }
+    else
+    {
+      for (int i = 0; i < 8; i++)
+        ((ElementBndBoxTree*) myChildren[i])->getElementsInBox( box, foundElems );
     }
   }
 
+  //================================================================================
+  /*!
+   * \brief Return a leaf including a point
+   */
+  //================================================================================
+
+  ElementBndBoxTree* ElementBndBoxTree::getLeafAtPoint( const gp_XYZ& point )
+  {
+    if ( getBox()->IsOut( point ))
+      return 0;
+
+    if ( isLeaf() )
+    {
+      return this;
+    }
+    else
+    {
+      for (int i = 0; i < 8; i++)
+        if ( ElementBndBoxTree* l = ((ElementBndBoxTree*) myChildren[i])->getLeafAtPoint( point ))
+          return l;
+    }
+    return 0;
+  }
+
+  //================================================================================
+  /*!
+   * \brief Return number of elements
+   */
+  //================================================================================
+
+  int ElementBndBoxTree::getNbElements()
+  {
+    int nb = 0;
+    if ( isLeaf() )
+    {
+      nb = _elements.size();
+    }
+    else
+    {
+      for (int i = 0; i < 8; i++)
+        nb += ((ElementBndBoxTree*) myChildren[i])->getNbElements();
+    }
+    return nb;
+  }
+
   //================================================================================
   /*!
    * \brief Construct the element box
@@ -480,7 +497,6 @@ namespace // Utils used in SMESH_ElementSearcherImpl::FindElementsByPoint()
   void ElementBndBoxTree::ElementBox::init(const SMDS_MeshElement* elem, double tolerance)
   {
     _element  = elem;
-    _isMarked = false;
     SMDS_ElemIteratorPtr nIt = elem->nodesIterator();
     while ( nIt->more() )
       Add( SMESH_NodeXYZ( nIt->next() ));
@@ -502,15 +518,15 @@ SMESH_ElementSearcher::~SMESH_ElementSearcher()
 
 struct SMESH_ElementSearcherImpl: public SMESH_ElementSearcher
 {
-  SMDS_Mesh*                   _mesh;
-  SMDS_ElemIteratorPtr         _meshPartIt;
-  ElementBndBoxTree*           _ebbTree      [SMDSAbs_NbElementTypes];
-  int                          _ebbTreeHeight[SMDSAbs_NbElementTypes];
-  SMESH_NodeSearcherImpl*      _nodeSearcher;
-  SMDSAbs_ElementType          _elementType;
-  double                       _tolerance;
-  bool                         _outerFacesFound;
-  set<const SMDS_MeshElement*> _outerFaces; // empty means "no internal faces at all"
+  SMDS_Mesh*                        _mesh;
+  SMDS_ElemIteratorPtr              _meshPartIt;
+  ElementBndBoxTree*                _ebbTree      [SMDSAbs_NbElementTypes];
+  int                               _ebbTreeHeight[SMDSAbs_NbElementTypes];
+  SMESH_NodeSearcherImpl*           _nodeSearcher;
+  SMDSAbs_ElementType               _elementType;
+  double                            _tolerance;
+  bool                              _outerFacesFound;
+  std::set<const SMDS_MeshElement*> _outerFaces; // empty means "no internal faces at all"
 
   SMESH_ElementSearcherImpl( SMDS_Mesh&           mesh,
                              double               tol=-1,
@@ -532,20 +548,26 @@ struct SMESH_ElementSearcherImpl: public SMESH_ElementSearcher
     }
     if ( _nodeSearcher ) delete _nodeSearcher; _nodeSearcher = 0;
   }
-  virtual int FindElementsByPoint(const gp_Pnt&                      point,
-                                  SMDSAbs_ElementType                type,
-                                  vector< const SMDS_MeshElement* >& foundElements);
+  virtual int FindElementsByPoint(const gp_Pnt&                           point,
+                                  SMDSAbs_ElementType                     type,
+                                  std::vector< const SMDS_MeshElement* >& foundElements);
   virtual TopAbs_State GetPointState(const gp_Pnt& point);
   virtual const SMDS_MeshElement* FindClosestTo( const gp_Pnt&       point,
                                                  SMDSAbs_ElementType type );
 
-  void GetElementsNearLine( const gp_Ax1&                      line,
-                            SMDSAbs_ElementType                type,
-                            vector< const SMDS_MeshElement* >& foundElems);
-  void GetElementsInSphere( const gp_XYZ&                      center,
-                            const double                       radius,
-                            SMDSAbs_ElementType                type,
-                            vector< const SMDS_MeshElement* >& foundElems);
+  virtual void GetElementsNearLine( const gp_Ax1&                           line,
+                                    SMDSAbs_ElementType                     type,
+                                    std::vector< const SMDS_MeshElement* >& foundElems);
+  virtual void GetElementsInSphere( const gp_XYZ&                           center,
+                                    const double                            radius,
+                                    SMDSAbs_ElementType                     type,
+                                    std::vector< const SMDS_MeshElement* >& foundElems);
+  virtual void GetElementsInBox( const Bnd_B3d&                          box,
+                                 SMDSAbs_ElementType                     type,
+                                 std::vector< const SMDS_MeshElement* >& foundElems);
+  virtual gp_XYZ Project(const gp_Pnt&            point,
+                         SMDSAbs_ElementType      type,
+                         const SMDS_MeshElement** closestElem);
   double getTolerance();
   bool getIntersParamOnLine(const gp_Lin& line, const SMDS_MeshElement* face,
                             const double tolerance, double & param);
@@ -633,7 +655,7 @@ double SMESH_ElementSearcherImpl::getTolerance()
         while ( nodeIt->more() )
         {
           double dist = n1.Distance( static_cast<const SMDS_MeshNode*>( nodeIt->next() ));
-          elemSize = max( dist, elemSize );
+          elemSize = std::max( dist, elemSize );
         }
       }
       _tolerance = 1e-4 * elemSize;
@@ -691,10 +713,10 @@ void SMESH_ElementSearcherImpl::findOuterBoundary(const SMDS_MeshElement* outerF
   // and BTW find out if there are internal faces at all.
 
   // checked links and links where outer boundary meets internal one
-  set< SMESH_TLink > visitedLinks, seamLinks;
+  std::set< SMESH_TLink > visitedLinks, seamLinks;
 
   // links to treat with already visited faces sharing them
-  list < TFaceLink > startLinks;
+  std::list < TFaceLink > startLinks;
 
   // load startLinks with the first outerFace
   startLinks.push_back( TFaceLink( outerFace->GetNode(0), outerFace->GetNode(1), outerFace));
@@ -736,8 +758,8 @@ void SMESH_ElementSearcherImpl::findOuterBoundary(const SMDS_MeshElement* outerF
         // direction from the link inside outerFace
         gp_Vec dirInOF = gp_Vec( ofNorm ) ^ n1n2;
         // sort all other faces by angle with the dirInOF
-        map< double, const SMDS_MeshElement* > angle2Face;
-        set< const SMDS_MeshElement*, TIDCompare >::const_iterator face = faces.begin();
+        std::map< double, const SMDS_MeshElement* > angle2Face;
+        std::set< const SMDS_MeshElement*, TIDCompare >::const_iterator face = faces.begin();
         for ( ; face != faces.end(); ++face )
         {
           if ( *face == outerFace ) continue;
@@ -746,7 +768,7 @@ void SMESH_ElementSearcherImpl::findOuterBoundary(const SMDS_MeshElement* outerF
           gp_Vec dirInF = gp_Vec( fNorm ) ^ n1n2;
           double angle = dirInOF.AngleWithRef( dirInF, n1n2 );
           if ( angle < 0 ) angle += 2. * M_PI;
-          angle2Face.insert( make_pair( angle, *face ));
+          angle2Face.insert( std::make_pair( angle, *face ));
         }
         if ( !angle2Face.empty() )
           outerFace2 = angle2Face.begin()->second;
@@ -791,9 +813,9 @@ void SMESH_ElementSearcherImpl::findOuterBoundary(const SMDS_MeshElement* outerF
 //=======================================================================
 
 int SMESH_ElementSearcherImpl::
-FindElementsByPoint(const gp_Pnt&                      point,
-                    SMDSAbs_ElementType                type,
-                    vector< const SMDS_MeshElement* >& foundElements)
+FindElementsByPoint(const gp_Pnt&                           point,
+                    SMDSAbs_ElementType                     type,
+                    std::vector< const SMDS_MeshElement* >& foundElements)
 {
   foundElements.clear();
   _elementType = type;
@@ -834,13 +856,9 @@ FindElementsByPoint(const gp_Pnt&                      point,
     {
       _ebbTree[_elementType] = new ElementBndBoxTree( *_mesh, type, _meshPartIt, tolerance );
     }
-    else
-    {
-      _ebbTree[ type ]->prepare();
-    }
-    vector< const SMDS_MeshElement* > suspectElems;
+    ElementBndBoxTree::TElemSeq suspectElems;
     _ebbTree[ type ]->getElementsNearPoint( point, suspectElems );
-    vector< const SMDS_MeshElement* >::iterator elem = suspectElems.begin();
+    ElementBndBoxTree::TElemSeq::iterator elem = suspectElems.begin();
     for ( ; elem != suspectElems.end(); ++elem )
       if ( !SMESH_MeshAlgos::IsOut( *elem, point, tolerance ))
         foundElements.push_back( *elem );
@@ -863,15 +881,15 @@ SMESH_ElementSearcherImpl::FindClosestTo( const gp_Pnt&       point,
   const SMDS_MeshElement* closestElem = 0;
   _elementType = type;
 
-  if ( type == SMDSAbs_Face || type == SMDSAbs_Volume )
+  if ( type == SMDSAbs_Face ||
+       type == SMDSAbs_Volume ||
+       type == SMDSAbs_Edge )
   {
     ElementBndBoxTree*& ebbTree = _ebbTree[ type ];
     if ( !ebbTree )
       ebbTree = new ElementBndBoxTree( *_mesh, type, _meshPartIt );
-    else
-      ebbTree->prepare();
 
-    vector<const SMDS_MeshElement*> suspectElems;
+    ElementBndBoxTree::TElemSeq suspectElems;
     ebbTree->getElementsNearPoint( point, suspectElems );
 
     if ( suspectElems.empty() && ebbTree->maxSize() > 0 )
@@ -883,28 +901,27 @@ SMESH_ElementSearcherImpl::FindClosestTo( const gp_Pnt&       point,
         radius = point.Distance( boxCenter ) - 0.5 * ebbTree->maxSize();
       if ( radius < 0 )
         radius = ebbTree->maxSize() / pow( 2., getTreeHeight()) / 2;
-      while ( suspectElems.empty() )
+      while ( suspectElems.empty() && radius < 1e100 )
       {
-        ebbTree->prepare();
         ebbTree->getElementsInSphere( point.XYZ(), radius, suspectElems );
         radius *= 1.1;
       }
     }
     double minDist = std::numeric_limits<double>::max();
-    multimap< double, const SMDS_MeshElement* > dist2face;
-    vector<const SMDS_MeshElement*>::iterator elem = suspectElems.begin();
+    std::multimap< double, const SMDS_MeshElement* > dist2face;
+    ElementBndBoxTree::TElemSeq::iterator elem = suspectElems.begin();
     for ( ; elem != suspectElems.end(); ++elem )
     {
       double dist = SMESH_MeshAlgos::GetDistance( *elem, point );
       if ( dist < minDist + 1e-10)
       {
         minDist = dist;
-        dist2face.insert( dist2face.begin(), make_pair( dist, *elem ));
+        dist2face.insert( dist2face.begin(), std::make_pair( dist, *elem ));
       }
     }
     if ( !dist2face.empty() )
     {
-      multimap< double, const SMDS_MeshElement* >::iterator d2f = dist2face.begin();
+      std::multimap< double, const SMDS_MeshElement* >::iterator d2f = dist2face.begin();
       closestElem = d2f->second;
       // if there are several elements at the same distance, select one
       // with GC closest to the point
@@ -956,8 +973,6 @@ TopAbs_State SMESH_ElementSearcherImpl::GetPointState(const gp_Pnt& point)
   ElementBndBoxTree*& ebbTree = _ebbTree[ SMDSAbs_Face ];
   if ( !ebbTree )
     ebbTree = new ElementBndBoxTree( *_mesh, _elementType, _meshPartIt );
-  else
-    ebbTree->prepare();
 
   // Algo: analyse transition of a line starting at the point through mesh boundary;
   // try three lines parallel to axis of the coordinate system and perform rough
@@ -965,22 +980,21 @@ TopAbs_State SMESH_ElementSearcherImpl::GetPointState(const gp_Pnt& point)
 
   const int nbAxes = 3;
   gp_Dir axisDir[ nbAxes ] = { gp::DX(), gp::DY(), gp::DZ() };
-  map< double, TInters >   paramOnLine2TInters[ nbAxes ];
-  list< TInters > tangentInters[ nbAxes ]; // of faces whose plane includes the line
-  multimap< int, int > nbInt2Axis; // to find the simplest case
+  std::map< double, TInters >   paramOnLine2TInters[ nbAxes ];
+  std::list< TInters > tangentInters[ nbAxes ]; // of faces whose plane includes the line
+  std::multimap< int, int > nbInt2Axis; // to find the simplest case
   for ( int axis = 0; axis < nbAxes; ++axis )
   {
     gp_Ax1 lineAxis( point, axisDir[axis]);
     gp_Lin line    ( lineAxis );
 
-    vector<const SMDS_MeshElement*> suspectFaces; // faces possibly intersecting the line
-    if ( axis > 0 ) ebbTree->prepare();
+    ElementBndBoxTree::TElemSeq suspectFaces; // faces possibly intersecting the line
     ebbTree->getElementsNearLine( lineAxis, suspectFaces );
 
     // Intersect faces with the line
 
-    map< double, TInters > & u2inters = paramOnLine2TInters[ axis ];
-    vector<const SMDS_MeshElement*>::iterator face = suspectFaces.begin();
+    std::map< double, TInters > & u2inters = paramOnLine2TInters[ axis ];
+    ElementBndBoxTree::TElemSeq::iterator face = suspectFaces.begin();
     for ( ; face != suspectFaces.end(); ++face )
     {
       // get face plane
@@ -1001,7 +1015,7 @@ TopAbs_State SMESH_ElementSearcherImpl::GetPointState(const gp_Pnt& point)
         double tol = 1e-4 * Sqrt( fNorm.Modulus() );
         gp_Pnt intersectionPoint = intersection.Point(1);
         if ( !SMESH_MeshAlgos::IsOut( *face, intersectionPoint, tol ))
-          u2inters.insert(make_pair( intersection.ParamOnConic(1), TInters( *face, fNorm )));
+          u2inters.insert( std::make_pair( intersection.ParamOnConic(1), TInters( *face, fNorm )));
       }
     }
     // Analyse intersections roughly
@@ -1025,7 +1039,7 @@ TopAbs_State SMESH_ElementSearcherImpl::GetPointState(const gp_Pnt& point)
     if ( nbIntBeforePoint == 1 || nbIntAfterPoint == 1 )
       return TopAbs_IN;
 
-    nbInt2Axis.insert( make_pair( min( nbIntBeforePoint, nbIntAfterPoint ), axis ));
+    nbInt2Axis.insert( std::make_pair( std::min( nbIntBeforePoint, nbIntAfterPoint ), axis ));
 
     if ( _outerFacesFound ) break; // pass to thorough analysis
 
@@ -1039,29 +1053,29 @@ TopAbs_State SMESH_ElementSearcherImpl::GetPointState(const gp_Pnt& point)
 
   for ( int hasPositionInfo = _outerFacesFound; hasPositionInfo < 2; ++hasPositionInfo )
   {
-    multimap< int, int >::const_iterator nb_axis = nbInt2Axis.begin();
+    std::multimap< int, int >::const_iterator nb_axis = nbInt2Axis.begin();
     for ( ; nb_axis != nbInt2Axis.end(); ++nb_axis )
     {
       int axis = nb_axis->second;
-      map< double, TInters > & u2inters = paramOnLine2TInters[ axis ];
+      std::map< double, TInters > & u2inters = paramOnLine2TInters[ axis ];
 
       gp_Ax1 lineAxis( point, axisDir[axis]);
       gp_Lin line    ( lineAxis );
 
       // add tangent intersections to u2inters
       double param;
-      list< TInters >::const_iterator tgtInt = tangentInters[ axis ].begin();
+      std::list< TInters >::const_iterator tgtInt = tangentInters[ axis ].begin();
       for ( ; tgtInt != tangentInters[ axis ].end(); ++tgtInt )
         if ( getIntersParamOnLine( line, tgtInt->_face, tolerance, param ))
-          u2inters.insert(make_pair( param, *tgtInt ));
+          u2inters.insert( std::make_pair( param, *tgtInt ));
       tangentInters[ axis ].clear();
 
       // Count intersections before and after the point excluding touching ones.
       // If hasPositionInfo we count intersections of outer boundary only
 
       int nbIntBeforePoint = 0, nbIntAfterPoint = 0;
-      double f = numeric_limits<double>::max(), l = -numeric_limits<double>::max();
-      map< double, TInters >::iterator u_int1 = u2inters.begin(), u_int2 = u_int1;
+      double f = std::numeric_limits<double>::max(), l = -std::numeric_limits<double>::max();
+      std::map< double, TInters >::iterator u_int1 = u2inters.begin(), u_int2 = u_int1;
       bool ok = ! u_int1->second._coincides;
       while ( ok && u_int1 != u2inters.end() )
       {
@@ -1110,8 +1124,8 @@ TopAbs_State SMESH_ElementSearcherImpl::GetPointState(const gp_Pnt& point)
           // decide if we skipped a touching intersection
           if ( nbSamePnt + nbTgt > 0 )
           {
-            double minDot = numeric_limits<double>::max(), maxDot = -numeric_limits<double>::max();
-            map< double, TInters >::iterator u_int = u_int1;
+            double minDot = std::numeric_limits<double>::max(), maxDot = -minDot;
+            std::map< double, TInters >::iterator u_int = u_int1;
             for ( ; u_int != u_int2; ++u_int )
             {
               if ( u_int->second._coincides ) continue;
@@ -1164,7 +1178,7 @@ TopAbs_State SMESH_ElementSearcherImpl::GetPointState(const gp_Pnt& point)
     if ( !hasPositionInfo )
     {
       // gather info on faces position - is face in the outer boundary or not
-      map< double, TInters > & u2inters = paramOnLine2TInters[ 0 ];
+      std::map< double, TInters > & u2inters = paramOnLine2TInters[ 0 ];
       findOuterBoundary( u2inters.begin()->second._face );
     }
 
@@ -1179,18 +1193,20 @@ TopAbs_State SMESH_ElementSearcherImpl::GetPointState(const gp_Pnt& point)
  */
 //=======================================================================
 
-void SMESH_ElementSearcherImpl::GetElementsNearLine( const gp_Ax1&                      line,
-                                                     SMDSAbs_ElementType                type,
-                                                     vector< const SMDS_MeshElement* >& foundElems)
+void SMESH_ElementSearcherImpl::
+GetElementsNearLine( const gp_Ax1&                           line,
+                     SMDSAbs_ElementType                     type,
+                     std::vector< const SMDS_MeshElement* >& foundElems)
 {
   _elementType = type;
   ElementBndBoxTree*& ebbTree = _ebbTree[ type ];
   if ( !ebbTree )
     ebbTree = new ElementBndBoxTree( *_mesh, _elementType, _meshPartIt );
-  else
-    ebbTree->prepare();
 
-  ebbTree->getElementsNearLine( line, foundElems );
+  ElementBndBoxTree::TElemSeq elems;
+  ebbTree->getElementsNearLine( line, elems );
+
+  foundElems.insert( foundElems.end(), elems.begin(), elems.end() );
 }
 
 //=======================================================================
@@ -1199,19 +1215,135 @@ void SMESH_ElementSearcherImpl::GetElementsNearLine( const gp_Ax1&
  */
 //=======================================================================
 
-void SMESH_ElementSearcherImpl::GetElementsInSphere( const gp_XYZ&                      center,
-                                                     const double                       radius,
-                                                     SMDSAbs_ElementType                type,
-                                                     vector< const SMDS_MeshElement* >& foundElems)
+void SMESH_ElementSearcherImpl::
+GetElementsInSphere( const gp_XYZ&                           center,
+                     const double                            radius,
+                     SMDSAbs_ElementType                     type,
+                     std::vector< const SMDS_MeshElement* >& foundElems)
 {
   _elementType = type;
   ElementBndBoxTree*& ebbTree = _ebbTree[ type ];
   if ( !ebbTree )
     ebbTree = new ElementBndBoxTree( *_mesh, _elementType, _meshPartIt );
-  else
-    ebbTree->prepare();
 
-  ebbTree->getElementsInSphere( center, radius, foundElems );
+  ElementBndBoxTree::TElemSeq elems;
+  ebbTree->getElementsInSphere( center, radius, elems );
+
+  foundElems.insert( foundElems.end(), elems.begin(), elems.end() );
+}
+
+//=======================================================================
+/*
+ * Return elements whose bounding box intersects a given bounding box
+ */
+//=======================================================================
+
+void SMESH_ElementSearcherImpl::
+GetElementsInBox( const Bnd_B3d&                          box,
+                  SMDSAbs_ElementType                     type,
+                  std::vector< const SMDS_MeshElement* >& foundElems)
+{
+  _elementType = type;
+  ElementBndBoxTree*& ebbTree = _ebbTree[ type ];
+  if ( !ebbTree )
+    ebbTree = new ElementBndBoxTree( *_mesh, _elementType, _meshPartIt, getTolerance() );
+
+  ElementBndBoxTree::TElemSeq elems;
+  ebbTree->getElementsInBox( box, elems );
+
+  foundElems.insert( foundElems.end(), elems.begin(), elems.end() );
+}
+
+//=======================================================================
+/*
+ * \brief Return a projection of a given point to a mesh.
+ *        Optionally return the closest element
+ */
+//=======================================================================
+
+gp_XYZ SMESH_ElementSearcherImpl::Project(const gp_Pnt&            point,
+                                          SMDSAbs_ElementType      type,
+                                          const SMDS_MeshElement** closestElem)
+{
+  _elementType = type;
+  if ( _mesh->GetMeshInfo().NbElements( _elementType ) == 0 )
+    throw SALOME_Exception( LOCALIZED( "No elements of given type in the mesh" ));
+
+  ElementBndBoxTree*& ebbTree = _ebbTree[ _elementType ];
+  if ( !ebbTree )
+    ebbTree = new ElementBndBoxTree( *_mesh, _elementType, _meshPartIt );
+
+  gp_XYZ p = point.XYZ();
+  ElementBndBoxTree* ebbLeaf = ebbTree->getLeafAtPoint( p );
+  const Bnd_B3d* box = ebbLeaf ? ebbLeaf->getBox() : ebbTree->getBox();
+  gp_XYZ pMin = box->CornerMin(), pMax = box->CornerMax();
+  double radius = Precision::Infinite();
+  if ( ebbLeaf || !box->IsOut( p ))
+  {
+    for ( int i = 1; i <= 3; ++i )
+    {
+      double d = 0.5 * ( pMax.Coord(i) - pMin.Coord(i) );
+      if ( d > Precision::Confusion() )
+        radius = Min( d, radius );
+    }
+    if ( !ebbLeaf )
+      radius /= ebbTree->getHeight( /*full=*/true );
+  }
+  else // p outside of box
+  {
+    for ( int i = 1; i <= 3; ++i )
+    {
+      double d = 0;
+      if ( point.Coord(i) < pMin.Coord(i) )
+        d = pMin.Coord(i) - point.Coord(i);
+      else if ( point.Coord(i) > pMax.Coord(i) )
+        d = point.Coord(i) - pMax.Coord(i);
+      if ( d > Precision::Confusion() )
+        radius = Min( d, radius );
+    }
+  }
+
+  ElementBndBoxTree::TElemSeq elems;
+  ebbTree->getElementsInSphere( p, radius, elems );
+  while ( elems.empty() && radius < 1e100 )
+  {
+    radius *= 1.1;
+    ebbTree->getElementsInSphere( p, radius, elems );
+  }
+  gp_XYZ proj, bestProj;
+  const SMDS_MeshElement* elem = 0;
+  double minDist = Precision::Infinite();
+  ElementBndBoxTree::TElemSeq::iterator e = elems.begin();
+  for ( ; e != elems.end(); ++e )
+  {
+    double d = SMESH_MeshAlgos::GetDistance( *e, point, &proj );
+    if ( d < minDist )
+    {
+      bestProj = proj;
+      elem = *e;
+      minDist = d;
+    }
+  }
+  if ( minDist > radius )
+  {
+    ElementBndBoxTree::TElemSeq elems2;
+    ebbTree->getElementsInSphere( p, minDist, elems2 );
+    for ( e = elems2.begin(); e != elems2.end(); ++e )
+    {
+      if ( elems.count( *e ))
+        continue;
+      double d = SMESH_MeshAlgos::GetDistance( *e, point, &proj );
+      if ( d < minDist )
+      {
+        bestProj = proj;
+        elem = *e;
+        minDist = d;
+      }
+    }
+  }
+  if ( closestElem ) *closestElem = elem;
+
+  return bestProj;
 }
 
 //=======================================================================
@@ -1229,9 +1361,9 @@ bool SMESH_MeshAlgos::IsOut( const SMDS_MeshElement* element, const gp_Pnt& poin
 
   // get ordered nodes
 
-  vector< SMESH_TNodeXYZ > xyz; xyz.reserve( element->NbNodes()+1 );
+  std::vector< SMESH_TNodeXYZ > xyz; xyz.reserve( element->NbNodes()+1 );
 
-  SMDS_ElemIteratorPtr nodeIt = element->interlacedNodesElemIterator();
+  SMDS_NodeIteratorPtr nodeIt = element->interlacedNodesIterator();
   for ( int i = 0; nodeIt->more(); ++i )
     xyz.push_back( SMESH_TNodeXYZ( nodeIt->next() ));
 
@@ -1294,7 +1426,7 @@ bool SMESH_MeshAlgos::IsOut( const SMDS_MeshElement* element, const gp_Pnt& poin
     if ( n2pSize * n2pSize > fNormSize * 100 ) return true; // point is very far
     plnNorm /= n2pSize;
     // for each node of the face, compute its signed distance to the cutting plane
-    vector<double> dist( nbNodes + 1);
+    std::vector<double> dist( nbNodes + 1);
     for ( i = 0; i < nbNodes; ++i )
     {
       gp_Vec n2p( xyz[i], point );
@@ -1397,7 +1529,8 @@ namespace
   //            .  RIGHT    .
   //            .           .
   enum PositionName { POS_LEFT = 1, POS_VERTEX = 2, POS_RIGHT = 4, //POS_ON = 8,
-                      POS_ALL = POS_LEFT | POS_RIGHT | POS_VERTEX };
+                      POS_ALL = POS_LEFT | POS_RIGHT | POS_VERTEX,
+                      POS_MAX = POS_RIGHT };
   struct PointPos
   {
     PositionName _name;
@@ -1414,9 +1547,9 @@ namespace
 
   //================================================================================
   /*!
-   * \brief Return of a point relative to a segment
+   * \brief Return position of a point relative to a segment
    *  \param point2D      - the point to analyze position of
-   *  \param xyVec        - end points of segments
+   *  \param segEnds      - end points of segments
    *  \param index0       - 0-based index of the first point of segment
    *  \param posToFindOut - flags of positions to detect
    *  \retval PointPos - point position
@@ -1461,17 +1594,19 @@ namespace
 //=======================================================================
 
 double SMESH_MeshAlgos::GetDistance( const SMDS_MeshElement* elem,
-                                     const gp_Pnt&           point )
+                                     const gp_Pnt&           point,
+                                     gp_XYZ*                 closestPnt )
 {
   switch ( elem->GetType() )
   {
   case SMDSAbs_Volume:
-    return GetDistance( dynamic_cast<const SMDS_MeshVolume*>( elem ), point);
+    return GetDistance( static_cast<const SMDS_MeshVolume*>( elem ), point, closestPnt );
   case SMDSAbs_Face:
-    return GetDistance( dynamic_cast<const SMDS_MeshFace*>( elem ), point);
+    return GetDistance( static_cast<const SMDS_MeshFace*>( elem ), point, closestPnt );
   case SMDSAbs_Edge:
-    return GetDistance( dynamic_cast<const SMDS_MeshEdge*>( elem ), point);
+    return GetDistance( static_cast<const SMDS_MeshEdge*>( elem ), point, closestPnt );
   case SMDSAbs_Node:
+    if ( closestPnt ) *closestPnt = SMESH_TNodeXYZ( elem );
     return point.Distance( SMESH_TNodeXYZ( elem ));
   default:;
   }
@@ -1487,15 +1622,41 @@ double SMESH_MeshAlgos::GetDistance( const SMDS_MeshElement* elem,
 //=======================================================================
 
 double SMESH_MeshAlgos::GetDistance( const SMDS_MeshFace* face,
-                                     const gp_Pnt&        point )
+                                     const gp_Pnt&        point,
+                                     gp_XYZ*              closestPnt )
 {
-  double badDistance = -1;
+  const double badDistance = -1;
   if ( !face ) return badDistance;
 
+  int nbCorners = face->NbCornerNodes();
+  if ( nbCorners > 3 )
+  {
+    std::vector< const SMDS_MeshNode* > nodes;
+    int nbTria = SMESH_MeshAlgos::Triangulate().GetTriangles( face, nodes );
+
+    double minDist = Precision::Infinite();
+    gp_XYZ cp;
+    for ( int i = 0; i < 3 * nbTria; i += 3 )
+    {
+      SMDS_FaceOfNodes triangle( nodes[i], nodes[i+1], nodes[i+2] );
+      double dist = GetDistance( &triangle, point, closestPnt );
+      if ( dist < minDist )
+      {
+        minDist = dist;
+        if ( closestPnt )
+          cp = *closestPnt;
+      }
+    }
+
+    if ( closestPnt )
+      *closestPnt = cp;
+    return minDist;
+  }
+
   // coordinates of nodes (medium nodes, if any, ignored)
   typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
-  vector<gp_XYZ> xyz( TXyzIterator( face->nodesIterator()), TXyzIterator() );
-  xyz.resize( face->NbCornerNodes()+1 );
+  std::vector<gp_XYZ> xyz( TXyzIterator( face->nodesIterator()), TXyzIterator() );
+  xyz.resize( 4 );
 
   // transformation to get xyz[0] lies on the origin, xyz[1] lies on the Z axis,
   // and xyz[2] lies in the XZ plane. This is to pass to 2D space on XZ plane.
@@ -1518,8 +1679,8 @@ double SMESH_MeshAlgos::GetDistance( const SMDS_MeshFace* face,
   trsf.SetTransformation( tgtCS );
 
   // move all the nodes to 2D
-  vector<gp_XY> xy( xyz.size() );
-  for ( size_t i = 0;i < xyz.size()-1; ++i )
+  std::vector<gp_XY> xy( xyz.size() );
+  for ( size_t i = 0; i < 3; ++i )
   {
     gp_XYZ p3d = xyz[i];
     trsf.Transforms( p3d );
@@ -1533,46 +1694,64 @@ double SMESH_MeshAlgos::GetDistance( const SMDS_MeshFace* face,
   trsf.Transforms( tmpPnt );
   gp_XY point2D( tmpPnt.X(), tmpPnt.Z() );
 
-  // loop on segments of the face to analyze point position ralative to the face
-  set< PointPos > pntPosSet;
+  // loop on edges of the face to analyze point position ralative to the face
+  std::vector< PointPos > pntPosByType[ POS_MAX + 1 ];
   for ( size_t i = 1; i < xy.size(); ++i )
   {
     PointPos pos = getPointPosition( point2D, &xy[0], i-1 );
-    pntPosSet.insert( pos );
+    pntPosByType[ pos._name ].push_back( pos );
   }
 
   // compute distance
-  PointPos pos = *pntPosSet.begin();
-  // cout << "Face " << face->GetID() << " DIST: ";
-  switch ( pos._name )
-  {
-  case POS_LEFT: {
-    // point is most close to a segment
-    gp_Vec p0p1( point, xyz[ pos._index ] );
-    gp_Vec p1p2( xyz[ pos._index ], xyz[ pos._index+1 ]); // segment vector
-    p1p2.Normalize();
-    double projDist = p0p1 * p1p2; // distance projected to the segment
-    gp_Vec projVec = p1p2 * projDist;
-    gp_Vec distVec = p0p1 - projVec;
-    // cout << distVec.Magnitude()  << ", SEG " << face->GetNode(pos._index)->GetID()
-    //      << " - " << face->GetNodeWrap(pos._index+1)->GetID() << endl;
-    return distVec.Magnitude();
-  }
-  case POS_RIGHT: {
-    // point is inside the face
-    double distToFacePlane = tmpPnt.Y();
-    // cout << distToFacePlane << ", INSIDE " << endl;
-    return Abs( distToFacePlane );
+
+  double dist = badDistance;
+
+  if ( pntPosByType[ POS_LEFT ].size() > 0 ) // point is most close to an edge
+  {
+    PointPos& pos = pntPosByType[ POS_LEFT ][0];
+
+    gp_Vec edge( xyz[ pos._index ], xyz[ pos._index+1 ]);
+    gp_Vec n1p ( xyz[ pos._index ], point  );
+    double u = ( edge * n1p ) / edge.SquareMagnitude(); // param [0,1] on the edge
+    gp_XYZ proj = xyz[ pos._index ] + u * edge.XYZ(); // projection on the edge
+    dist = point.Distance( proj );
+    if ( closestPnt ) *closestPnt = proj;
   }
-  case POS_VERTEX: {
-    // point is most close to a node
-    gp_Vec distVec( point, xyz[ pos._index ]);
-    // cout << distVec.Magnitude()  << " VERTEX " << face->GetNode(pos._index)->GetID() << endl;
-    return distVec.Magnitude();
+
+  else if ( pntPosByType[ POS_RIGHT ].size() >= 2 ) // point is inside the face
+  {
+    dist = Abs( tmpPnt.Y() );
+    if ( closestPnt )
+    {
+      if ( dist < std::numeric_limits<double>::min() ) {
+        *closestPnt = point.XYZ();
+      }
+      else {
+        tmpPnt.SetY( 0 );
+        trsf.Inverted().Transforms( tmpPnt );
+        *closestPnt = tmpPnt;
+      }
+    }
   }
-  default:;
+
+  else if ( pntPosByType[ POS_VERTEX ].size() > 0 ) // point is most close to a node
+  {
+    double minDist2 = Precision::Infinite();
+    for ( size_t i = 0; i < pntPosByType[ POS_VERTEX ].size(); ++i )
+    {
+      PointPos& pos = pntPosByType[ POS_VERTEX ][i];
+
+      double d2 = point.SquareDistance( xyz[ pos._index ]);
+      if ( minDist2 > d2 )
+      {
+        if ( closestPnt ) *closestPnt = xyz[ pos._index ];
+        minDist2 = d2;
+      }
+    }
+    dist = Sqrt( minDist2 );
   }
-  return badDistance;
+
+  return dist;
 }
 
 //=======================================================================
@@ -1581,32 +1760,42 @@ double SMESH_MeshAlgos::GetDistance( const SMDS_MeshFace* face,
  */
 //=======================================================================
 
-double SMESH_MeshAlgos::GetDistance( const SMDS_MeshEdge* seg, const gp_Pnt& point )
+double SMESH_MeshAlgos::GetDistance( const SMDS_MeshEdge* seg,
+                                     const gp_Pnt&        point,
+                                     gp_XYZ*              closestPnt )
 {
   double dist = Precision::Infinite();
   if ( !seg ) return dist;
 
   int i = 0, nbNodes = seg->NbNodes();
 
-  vector< SMESH_TNodeXYZ > xyz( nbNodes );
-  SMDS_ElemIteratorPtr nodeIt = seg->interlacedNodesElemIterator();
-  while ( nodeIt->more() )
-    xyz[ i++ ].Set( nodeIt->next() );
+  std::vector< SMESH_TNodeXYZ > xyz( nbNodes );
+  for ( SMDS_NodeIteratorPtr nodeIt = seg->interlacedNodesIterator(); nodeIt->more(); i++ )
+    xyz[ i ].Set( nodeIt->next() );
 
   for ( i = 1; i < nbNodes; ++i )
   {
     gp_Vec edge( xyz[i-1], xyz[i] );
     gp_Vec n1p ( xyz[i-1], point  );
-    double u = ( edge * n1p ) / edge.SquareMagnitude(); // param [0,1] on the edge
+    double d, u = ( edge * n1p ) / edge.SquareMagnitude(); // param [0,1] on the edge
     if ( u <= 0. ) {
-      dist = Min( dist, n1p.SquareMagnitude() );
+      if (( d = n1p.SquareMagnitude() ) < dist ) {
+        dist = d;
+        if ( closestPnt ) *closestPnt = xyz[i-1];
+      }
     }
     else if ( u >= 1. ) {
-      dist = Min( dist, point.SquareDistance( xyz[i] ));
+      if (( d = point.SquareDistance( xyz[i] )) < dist ) {
+        dist = d;
+        if ( closestPnt ) *closestPnt = xyz[i];
+      }
     }
     else {
-      gp_XYZ proj = ( 1. - u ) * xyz[i-1] + u * xyz[i]; // projection of the point on the edge
-      dist = Min( dist, point.SquareDistance( proj ));
+      gp_XYZ proj = xyz[i-1] + u * edge.XYZ(); // projection of the point on the edge
+      if (( d = point.SquareDistance( proj )) < dist ) {
+        dist = d;
+        if ( closestPnt ) *closestPnt = proj;
+      }
     }
   }
   return Sqrt( dist );
@@ -1620,7 +1809,9 @@ double SMESH_MeshAlgos::GetDistance( const SMDS_MeshEdge* seg, const gp_Pnt& poi
  */
 //=======================================================================
 
-double SMESH_MeshAlgos::GetDistance( const SMDS_MeshVolume* volume, const gp_Pnt& point )
+double SMESH_MeshAlgos::GetDistance( const SMDS_MeshVolume* volume,
+                                     const gp_Pnt&          point,
+                                     gp_XYZ*                closestPnt )
 {
   SMDS_VolumeTool vTool( volume );
   vTool.SetExternalNormal();
@@ -1628,6 +1819,8 @@ double SMESH_MeshAlgos::GetDistance( const SMDS_MeshVolume* volume, const gp_Pnt
 
   double n[3], bc[3];
   double minDist = 1e100, dist;
+  gp_XYZ closeP = point.XYZ();
+  bool isOut = false;
   for ( int iF = 0; iF < vTool.NbFaces(); ++iF )
   {
     // skip a facet with normal not "looking at" the point
@@ -1635,7 +1828,7 @@ double SMESH_MeshAlgos::GetDistance( const SMDS_MeshVolume* volume, const gp_Pnt
          !vTool.GetFaceBaryCenter( iF, bc[0], bc[1], bc[2] ))
       continue;
     gp_XYZ bcp = point.XYZ() - gp_XYZ( bc[0], bc[1], bc[2] );
-    if ( gp_XYZ( n[0], n[1], n[2] ) * bcp < 1e-6 )
+    if ( gp_XYZ( n[0], n[1], n[2] ) * bcp < -1e-12 )
       continue;
 
     // find distance to a facet
@@ -1644,23 +1837,34 @@ double SMESH_MeshAlgos::GetDistance( const SMDS_MeshVolume* volume, const gp_Pnt
     case 3:
     {
       SMDS_FaceOfNodes tmpFace( nodes[0], nodes[ 1*iQ ], nodes[ 2*iQ ] );
-      dist = GetDistance( &tmpFace, point );
+      dist = GetDistance( &tmpFace, point, closestPnt );
       break;
     }
     case 4:
     {
       SMDS_FaceOfNodes tmpFace( nodes[0], nodes[ 1*iQ ], nodes[ 2*iQ ], nodes[ 3*iQ ]);
-      dist = GetDistance( &tmpFace, point );
+      dist = GetDistance( &tmpFace, point, closestPnt );
       break;
     }
     default:
-      vector<const SMDS_MeshNode *> nvec( nodes, nodes + vTool.NbFaceNodes( iF ));
+      std::vector<const SMDS_MeshNode *> nvec( nodes, nodes + vTool.NbFaceNodes( iF ));
       SMDS_PolygonalFaceOfNodes tmpFace( nvec );
-      dist = GetDistance( &tmpFace, point );
+      dist = GetDistance( &tmpFace, point, closestPnt );
+    }
+    if ( dist < minDist )
+    {
+      minDist = dist;
+      isOut = true;
+      if ( closestPnt ) closeP = *closestPnt;
     }
-    minDist = Min( minDist, dist );
   }
-  return minDist;
+  if ( isOut )
+  {
+    if ( closestPnt ) *closestPnt = closeP;
+    return minDist;
+  }
+
+  return 0; // point is inside the volume
 }
 
 //================================================================================
@@ -1691,7 +1895,7 @@ void SMESH_MeshAlgos::GetBarycentricCoords( const gp_XY& p,
   const double t11 = T22, t12 = -T12, t21 = -T21, t22 = T11;
   // vector
   const double r11 = p.X()-t2.X(), r12 = p.Y()-t2.Y();
-  // barycentric coordinates: mutiply matrix by vector
+  // barycentric coordinates: multiply matrix by vector
   bc0 = (t11 * r11 + t12 * r12)/Tdet;
   bc1 = (t21 * r11 + t22 * r12)/Tdet;
 }
@@ -1737,7 +1941,7 @@ SMESH_MeshAlgos::FindFaceInSet(const SMDS_MeshNode*    n1,
     if ( !face && elem->IsQuadratic())
     {
       // analysis for quadratic elements using all nodes
-      SMDS_ElemIteratorPtr anIter = elem->interlacedNodesElemIterator();
+      SMDS_NodeIteratorPtr anIter = elem->interlacedNodesIterator();
       const SMDS_MeshNode* prevN = static_cast<const SMDS_MeshNode*>( anIter->next() );
       for ( i1 = -1, i2 = 0; anIter->more() && !face; i1++, i2++ )
       {
@@ -1748,7 +1952,7 @@ SMESH_MeshAlgos::FindFaceInSet(const SMDS_MeshNode*    n1,
         }
         else if ( n2 == prevN && n1 == n )
         {
-          face = elem; swap( i1, i2 );
+          face = elem; std::swap( i1, i2 );
         }
         prevN = n;
       }
@@ -1759,6 +1963,222 @@ SMESH_MeshAlgos::FindFaceInSet(const SMDS_MeshNode*    n1,
   return face;
 }
 
+//================================================================================
+/*!
+ * Return sharp edges of faces and non-manifold ones. Optionally adds existing edges.
+ */
+//================================================================================
+
+std::vector< SMESH_MeshAlgos::Edge >
+SMESH_MeshAlgos::FindSharpEdges( SMDS_Mesh* theMesh,
+                                 double     theAngle,
+                                 bool       theAddExisting )
+{
+  std::vector< Edge > resultEdges;
+  if ( !theMesh ) return resultEdges;
+
+  typedef std::pair< bool, const SMDS_MeshNode* >                            TIsSharpAndMedium;
+  typedef NCollection_DataMap< SMESH_TLink, TIsSharpAndMedium, SMESH_TLink > TLinkSharpMap;
+
+  TLinkSharpMap linkIsSharp( theMesh->NbFaces() );
+  TIsSharpAndMedium sharpMedium( true, 0 );
+  bool                 & isSharp = sharpMedium.first;
+  const SMDS_MeshNode* & nMedium = sharpMedium.second;
+
+  if ( theAddExisting )
+  {
+    for ( SMDS_EdgeIteratorPtr edgeIt = theMesh->edgesIterator(); edgeIt->more(); )
+    {
+      const SMDS_MeshElement* edge = edgeIt->next();
+      nMedium = ( edge->IsQuadratic() ) ? edge->GetNode(2) : 0;
+      linkIsSharp.Bind( SMESH_TLink( edge->GetNode(0), edge->GetNode(1)), sharpMedium );
+    }
+  }
+
+  // check angles between face normals
+
+  const double angleCos = Cos( theAngle * M_PI / 180. ), angleCos2 = angleCos * angleCos;
+  gp_XYZ norm1, norm2;
+  std::vector< const SMDS_MeshNode* > faceNodes, linkNodes(2);
+  std::vector<const SMDS_MeshElement *> linkFaces;
+
+  int nbSharp = linkIsSharp.Extent();
+  for ( SMDS_FaceIteratorPtr faceIt = theMesh->facesIterator(); faceIt->more(); )
+  {
+    const SMDS_MeshElement* face = faceIt->next();
+    size_t             nbCorners = face->NbCornerNodes();
+
+    faceNodes.assign( face->begin_nodes(), face->end_nodes() );
+    if ( faceNodes.size() == nbCorners )
+      faceNodes.resize( nbCorners * 2, 0 );
+
+    const SMDS_MeshNode* nPrev = faceNodes[ nbCorners-1 ];
+    for ( size_t i = 0; i < nbCorners; ++i )
+    {
+      SMESH_TLink link( nPrev, faceNodes[i] );
+      if ( !linkIsSharp.IsBound( link ))
+      {
+        linkNodes[0] = link.node1();
+        linkNodes[1] = link.node2();
+        linkFaces.clear();
+        theMesh->GetElementsByNodes( linkNodes, linkFaces, SMDSAbs_Face );
+
+        isSharp = false;
+        if ( linkFaces.size() > 2 )
+        {
+          isSharp = true;
+        }
+        else if ( linkFaces.size() == 2 &&
+                  FaceNormal( linkFaces[0], norm1, /*normalize=*/false ) &&
+                  FaceNormal( linkFaces[1], norm2, /*normalize=*/false ))
+        {
+          double dot = norm1 * norm2; // == cos * |norm1| * |norm2|
+          if (( dot < 0 ) == ( angleCos < 0 ))
+          {
+            double cos2 = dot * dot / norm1.SquareModulus() / norm2.SquareModulus();
+            isSharp = ( angleCos < 0 ) ? ( cos2 > angleCos2 ) : ( cos2 < angleCos2 );
+          }
+          else
+          {
+            isSharp = ( angleCos > 0 );
+          }
+        }
+        nMedium = faceNodes[( i-1+nbCorners ) % nbCorners + nbCorners ];
+
+        linkIsSharp.Bind( link, sharpMedium );
+        nbSharp += isSharp;
+      }
+
+      nPrev = faceNodes[i];
+    }
+  }
+
+  resultEdges.resize( nbSharp );
+  TLinkSharpMap::Iterator linkIsSharpIter( linkIsSharp );
+  for ( int i = 0; linkIsSharpIter.More() && i < nbSharp; linkIsSharpIter.Next() )
+  {
+    const SMESH_TLink&                link = linkIsSharpIter.Key();
+    const TIsSharpAndMedium& isSharpMedium = linkIsSharpIter.Value();
+    if ( isSharpMedium.first )
+    {
+      Edge & edge  = resultEdges[ i++ ];
+      edge._node1  = link.node1();
+      edge._node2  = link.node2();
+      edge._medium = isSharpMedium.second;
+    }
+  }
+
+  return resultEdges;
+}
+
+//================================================================================
+/*!
+ * Distribute all faces of the mesh between groups using given edges as group boundaries
+ */
+//================================================================================
+
+std::vector< std::vector< const SMDS_MeshElement* > >
+SMESH_MeshAlgos::SeparateFacesByEdges( SMDS_Mesh* theMesh, const std::vector< Edge >& theEdges )
+{
+  std::vector< std::vector< const SMDS_MeshElement* > > groups;
+  if ( !theMesh ) return groups;
+
+  // build map of face edges (SMESH_TLink) and their faces
+
+  typedef std::vector< const SMDS_MeshElement* >                    TFaceVec;
+  typedef NCollection_DataMap< SMESH_TLink, TFaceVec, SMESH_TLink > TFacesByLinks;
+  TFacesByLinks facesByLink( theMesh->NbFaces() );
+
+  std::vector< const SMDS_MeshNode* > faceNodes;
+  for ( SMDS_FaceIteratorPtr faceIt = theMesh->facesIterator(); faceIt->more(); )
+  {
+    const SMDS_MeshElement* face = faceIt->next();
+    size_t             nbCorners = face->NbCornerNodes();
+
+    faceNodes.assign( face->begin_nodes(), face->end_nodes() );
+    faceNodes.resize( nbCorners + 1 );
+    faceNodes[ nbCorners ] = faceNodes[0];
+
+    face->setIsMarked( false );
+
+    for ( size_t i = 0; i < nbCorners; ++i )
+    {
+      SMESH_TLink link( faceNodes[i], faceNodes[i+1] );
+      TFaceVec* linkFaces = facesByLink.ChangeSeek( link );
+      if ( !linkFaces )
+      {
+        linkFaces = facesByLink.Bound( link, TFaceVec() );
+        linkFaces->reserve(2);
+      }
+      linkFaces->push_back( face );
+    }
+  }
+
+  // remove the given edges from facesByLink map
+
+  for ( size_t i = 0; i < theEdges.size(); ++i )
+  {
+    SMESH_TLink link( theEdges[i]._node1, theEdges[i]._node2 );
+    facesByLink.UnBind( link );
+  }
+
+  // faces connected via links of facesByLink map form a group
+
+  while ( !facesByLink.IsEmpty() )
+  {
+    groups.push_back( TFaceVec() );
+    TFaceVec & group = groups.back();
+
+    group.push_back( TFacesByLinks::Iterator( facesByLink ).Value()[0] );
+    group.back()->setIsMarked( true );
+
+    for ( size_t iF = 0; iF < group.size(); ++iF )
+    {
+      const SMDS_MeshElement* face = group[iF];
+      size_t             nbCorners = face->NbCornerNodes();
+      faceNodes.assign( face->begin_nodes(), face->end_nodes() );
+      faceNodes.resize( nbCorners + 1 );
+      faceNodes[ nbCorners ] = faceNodes[0];
+
+      for ( size_t iN = 0; iN < nbCorners; ++iN )
+      {
+        SMESH_TLink link( faceNodes[iN], faceNodes[iN+1] );
+        if ( const TFaceVec* faces = facesByLink.Seek( link ))
+        {
+          const TFaceVec& faceNeighbors = *faces;
+          for ( size_t i = 0; i < faceNeighbors.size(); ++i )
+            if ( !faceNeighbors[i]->isMarked() )
+            {
+              group.push_back( faceNeighbors[i] );
+              faceNeighbors[i]->setIsMarked( true );
+            }
+          facesByLink.UnBind( link );
+        }
+      }
+    }
+  }
+
+  // find faces that are alone in its group; they were not in facesByLink
+
+  int nbInGroups = 0;
+  for ( size_t i = 0; i < groups.size(); ++i )
+    nbInGroups += groups[i].size();
+  if ( nbInGroups < theMesh->NbFaces() )
+  {
+    for ( SMDS_FaceIteratorPtr faceIt = theMesh->facesIterator(); faceIt->more(); )
+    {
+      const SMDS_MeshElement* face = faceIt->next();
+      if ( !face->isMarked() )
+      {
+        groups.push_back( TFaceVec() );
+        groups.back().push_back( face );
+      }
+    }
+  }
+
+  return groups;
+}
+
 //================================================================================
 /*!
  * \brief Calculate normal of a mesh face
@@ -1783,7 +2203,7 @@ bool SMESH_MeshAlgos::FaceNormal(const SMDS_MeshElement* F, gp_XYZ& normal, bool
     normal += ( p[2] - p[1] ) ^ ( p[0] - p[1] );
   }
   double size2 = normal.SquareModulus();
-  bool ok = ( size2 > numeric_limits<double>::min() * numeric_limits<double>::min());
+  bool ok = ( size2 > std::numeric_limits<double>::min() * std::numeric_limits<double>::min());
   if ( normalized && ok )
     normal /= sqrt( size2 );
 
@@ -1795,15 +2215,232 @@ bool SMESH_MeshAlgos::FaceNormal(const SMDS_MeshElement* F, gp_XYZ& normal, bool
 //purpose  : Return nodes common to two elements
 //=======================================================================
 
-vector< const SMDS_MeshNode*> SMESH_MeshAlgos::GetCommonNodes(const SMDS_MeshElement* e1,
-                                                              const SMDS_MeshElement* e2)
+std::vector< const SMDS_MeshNode*> SMESH_MeshAlgos::GetCommonNodes(const SMDS_MeshElement* e1,
+                                                                   const SMDS_MeshElement* e2)
 {
-  vector< const SMDS_MeshNode*> common;
+  std::vector< const SMDS_MeshNode*> common;
   for ( int i = 0 ; i < e1->NbNodes(); ++i )
     if ( e2->GetNodeIndex( e1->GetNode( i )) >= 0 )
       common.push_back( e1->GetNode( i ));
   return common;
 }
+//================================================================================
+/*!
+ * \brief Return true if node1 encounters first in the face and node2, after
+ */
+//================================================================================
+
+bool SMESH_MeshAlgos::IsRightOrder( const SMDS_MeshElement* face,
+                                    const SMDS_MeshNode*    node0,
+                                    const SMDS_MeshNode*    node1 )
+{
+  int i0 = face->GetNodeIndex( node0 );
+  int i1 = face->GetNodeIndex( node1 );
+  if ( face->IsQuadratic() )
+  {
+    if ( face->IsMediumNode( node0 ))
+    {
+      i0 -= ( face->NbNodes()/2 - 1 );
+      i1 *= 2;
+    }
+    else
+    {
+      i1 -= ( face->NbNodes()/2 - 1 );
+      i0 *= 2;
+    }
+  }
+  int diff = i1 - i0;
+  return ( diff == 1 ) || ( diff == -face->NbNodes()+1 );
+}
+
+//=======================================================================
+/*!
+ * \brief Partition given 1D elements into groups of contiguous edges.
+ *        A node where number of meeting edges != 2 is a group end.
+ *        An optional startNode is used to orient groups it belongs to.
+ * \return a list of edge groups and a list of corresponding node groups.
+ *         If a group is closed, the first and last nodes of the group are same.
+ */
+//=======================================================================
+
+void SMESH_MeshAlgos::Get1DBranches( SMDS_ElemIteratorPtr theEdgeIt,
+                                     TElemGroupVector&    theEdgeGroups,
+                                     TNodeGroupVector&    theNodeGroups,
+                                     const SMDS_MeshNode* theStartNode )
+{
+  if ( !theEdgeIt )
+    return;
+
+  // build map of nodes and their adjacent edges
+
+  typedef std::vector< const SMDS_MeshNode* >                                 TNodeVec;
+  typedef std::vector< const SMDS_MeshElement* >                              TEdgeVec;
+  typedef NCollection_DataMap< const SMDS_MeshNode*, TEdgeVec, SMESH_Hasher > TEdgesByNodeMap;
+  TEdgesByNodeMap edgesByNode;
+
+  while ( theEdgeIt->more() )
+  {
+    const SMDS_MeshElement* edge = theEdgeIt->next();
+    if ( edge->GetType() != SMDSAbs_Edge )
+      continue;
+
+    const SMDS_MeshNode* nodes[2] = { edge->GetNode(0), edge->GetNode(1) };
+    for ( int i = 0; i < 2; ++i )
+    {
+      TEdgeVec* nodeEdges = edgesByNode.ChangeSeek( nodes[i] );
+      if ( !nodeEdges )
+      {
+        nodeEdges = edgesByNode.Bound( nodes[i], TEdgeVec() );
+        nodeEdges->reserve(2);
+      }
+      nodeEdges->push_back( edge );
+    }
+  }
+
+  if ( edgesByNode.IsEmpty() )
+    return;
+
+
+  // build edge branches
+
+  TElemGroupVector branches(2);
+  TNodeGroupVector nodeBranches(2);
+
+  while ( !edgesByNode.IsEmpty() )
+  {
+    if ( !theStartNode || !edgesByNode.IsBound( theStartNode ))
+    {
+      theStartNode = TEdgesByNodeMap::Iterator( edgesByNode ).Key();
+    }
+
+    size_t nbBranches = 0;
+    bool startIsBranchEnd = false;
+
+    while ( edgesByNode.IsBound( theStartNode ))
+    {
+      // initialize a new branch
+
+      ++nbBranches;
+      if ( branches.size() < nbBranches )
+      {
+        branches.push_back   ( TEdgeVec() );
+        nodeBranches.push_back( TNodeVec() );
+      }
+      TEdgeVec & branch     = branches    [ nbBranches - 1 ];
+      TNodeVec & nodeBranch = nodeBranches[ nbBranches - 1 ];
+      branch.clear();
+      nodeBranch.clear();
+      {
+        TEdgeVec& edges = edgesByNode( theStartNode );
+        startIsBranchEnd = ( edges.size() != 2 );
+
+        int nbEdges = 0;
+        const SMDS_MeshElement* startEdge = 0;
+        for ( size_t i = 0; i < edges.size(); ++i )
+        {
+          if ( !startEdge && edges[i] )
+          {
+            startEdge = edges[i];
+            edges[i] = 0;
+          }
+          nbEdges += bool( edges[i] );
+        }
+        if ( nbEdges == 0 )
+          edgesByNode.UnBind( theStartNode );
+        if ( !startEdge )
+          continue;
+
+        branch.push_back( startEdge );
+
+        nodeBranch.push_back( theStartNode );
+        nodeBranch.push_back( branch.back()->GetNode(0) );
+        if ( nodeBranch.back() == theStartNode )
+          nodeBranch.back() = branch.back()->GetNode(1);
+      }
+
+      // fill the branch
+
+      bool isBranchEnd = false;
+      TEdgeVec* edgesPtr;
+
+      while (( !isBranchEnd ) && ( edgesPtr = edgesByNode.ChangeSeek( nodeBranch.back() )))
+      {
+        TEdgeVec& edges = *edgesPtr;
+
+        isBranchEnd = ( edges.size() != 2 );
+
+        const SMDS_MeshNode* lastNode = nodeBranch.back();
+
+        switch ( edges.size() )
+        {
+        case 1:
+          edgesByNode.UnBind( lastNode );
+          break;
+
+        case 2:
+        {
+          if ( const SMDS_MeshElement* nextEdge = edges[ edges[0] == branch.back() ])
+          {
+            branch.push_back( nextEdge );
+
+            const SMDS_MeshNode* nextNode = nextEdge->GetNode(0);
+            if ( nodeBranch.back() == nextNode )
+              nextNode = nextEdge->GetNode(1);
+            nodeBranch.push_back( nextNode );
+          }
+          edgesByNode.UnBind( lastNode );
+          break;
+        }
+
+        default:
+          int nbEdges = 0;
+          for ( size_t i = 0; i < edges.size(); ++i )
+          {
+            if ( edges[i] == branch.back() )
+              edges[i] = 0;
+            nbEdges += bool( edges[i] );
+          }
+          if ( nbEdges == 0 )
+            edgesByNode.UnBind( lastNode );
+        }
+      }
+    } // while ( edgesByNode.IsBound( theStartNode ))
+
+
+    // put the found branches to the result
+
+    if ( nbBranches == 2 && !startIsBranchEnd ) // join two branches starting at the same node
+    {
+      std::reverse( nodeBranches[0].begin(), nodeBranches[0].end() );
+      nodeBranches[0].pop_back();
+      nodeBranches[0].reserve( nodeBranches[0].size() + nodeBranches[1].size() );
+      nodeBranches[0].insert( nodeBranches[0].end(),
+                              nodeBranches[1].begin(), nodeBranches[1].end() );
+
+      std::reverse( branches[0].begin(), branches[0].end() );
+      branches[0].reserve( branches[0].size() + branches[1].size() );
+      branches[0].insert( branches[0].end(), branches[1].begin(), branches[1].end() );
+
+      nodeBranches[1].clear();
+      branches[1].clear();
+    }
+
+    for ( size_t i = 0; i < nbBranches; ++i )
+    {
+      if ( branches[i].empty() )
+        continue;
+
+      theEdgeGroups.push_back( TEdgeVec() );
+      theEdgeGroups.back().swap( branches[i] );
+
+      theNodeGroups.push_back( TNodeVec() );
+      theNodeGroups.back().swap( nodeBranches[i] );
+    }
+
+  } // while ( !edgesByNode.IsEmpty() )
+
+  return;
+}
 
 //=======================================================================
 /*!