Salome HOME
0020832: EDF 1359 SMESH : Automatic meshing of boundary layers
authoreap <eap@opencascade.com>
Tue, 18 Jan 2011 11:01:28 +0000 (11:01 +0000)
committereap <eap@opencascade.com>
Tue, 18 Jan 2011 11:01:28 +0000 (11:01 +0000)
+ SMESH_ProxyMesh.hxx

src/SMESH/Makefile.am
src/SMESH/SMESH_ProxyMesh.cxx [new file with mode: 0644]
src/SMESH/SMESH_ProxyMesh.hxx [new file with mode: 0644]

index 47543642883ce6afd8481d75c44aa445beda1390..ac1138b4985b90fc36c827287a1c059177787bd1 100644 (file)
@@ -22,7 +22,6 @@
 #  Author : Paul RASCLE, EDF
 #  Modified by : Alexander BORODIN (OCN) - autotools usage
 #  Module : SMESH
-#  $Header$
 #
 include $(top_srcdir)/adm_local/unix/make_common_starter.am
 
@@ -53,7 +52,8 @@ salomeinclude_HEADERS = \
        SMESH_Comment.hxx \
        SMESH_ComputeError.hxx \
        SMESH_File.hxx \
-       SMESH_SMESH.hxx
+       SMESH_SMESH.hxx \
+       SMESH_ProxyMesh.hxx
 
 # Libraries targets
 
@@ -78,7 +78,8 @@ dist_libSMESHimpl_la_SOURCES = \
        SMESH_MesherHelper.cxx \
        SMESH_Octree.cxx \
        SMESH_OctreeNode.cxx \
-       SMESH_File.cxx
+       SMESH_File.cxx \
+       SMESH_ProxyMesh.cxx
 
 # additionnal information to compile and link file
 libSMESHimpl_la_CPPFLAGS = \
diff --git a/src/SMESH/SMESH_ProxyMesh.cxx b/src/SMESH/SMESH_ProxyMesh.cxx
new file mode 100644 (file)
index 0000000..aba3c8e
--- /dev/null
@@ -0,0 +1,545 @@
+//  Copyright (C) 2007-2010  CEA/DEN, EDF R&D, OPEN CASCADE
+//
+//  This library is free software; you can redistribute it and/or
+//  modify it under the terms of the GNU Lesser General Public
+//  License as published by the Free Software Foundation; either
+//  version 2.1 of the License.
+//
+//  This library is distributed in the hope that it will be useful,
+//  but WITHOUT ANY WARRANTY; without even the implied warranty of
+//  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+//  Lesser General Public License for more details.
+//
+//  You should have received a copy of the GNU Lesser General Public
+//  License along with this library; if not, write to the Free Software
+//  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
+//
+//  See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
+//
+// File      : SMESH_ProxyMesh.cxx
+// Created   : Thu Dec  2 12:32:53 2010
+// Author    : Edward AGAPOV (eap)
+
+#include "SMESH_ProxyMesh.hxx"
+
+#include "SMDS_IteratorOnIterators.hxx"
+#include "SMDS_SetIterator.hxx"
+#include "SMESH_MesherHelper.hxx"
+
+#include <TopTools_ListIteratorOfListOfShape.hxx>
+#include <TopExp.hxx>
+#include <TopTools_IndexedMapOfShape.hxx>
+
+//================================================================================
+/*!
+ * \brief Constructor; mesh must be set by a descendant class
+ */
+//================================================================================
+
+SMESH_ProxyMesh::SMESH_ProxyMesh():_mesh(0)
+{
+}
+//================================================================================
+/*!
+ * \brief Make a proxy mesh from components. Components become empty
+ */
+//================================================================================
+
+SMESH_ProxyMesh::SMESH_ProxyMesh(vector<SMESH_ProxyMesh::Ptr>& components):
+  _mesh(0)
+{
+  if ( components.empty() ) return;
+
+  for ( unsigned i = 0; i < components.size(); ++i )
+  {
+    SMESH_ProxyMesh* m = components[i].get();
+    if ( !m ) continue;
+
+    takeTmpElemsInMesh( m );
+
+    if ( !_mesh ) _mesh = m->_mesh;
+    if ( _allowedTypes.empty() ) _allowedTypes = m->_allowedTypes;
+
+    if ( _subMeshes.size() < m->_subMeshes.size() )
+      _subMeshes.resize( m->_subMeshes.size(), 0 );
+    for ( unsigned j = 0; j < m->_subMeshes.size(); ++j )
+    {
+      if ( !m->_subMeshes[j] ) continue;
+      if ( _subMeshes[j] )
+      {
+        // unite 2 sub-meshes
+        set< const SMDS_MeshElement * > elems( _subMeshes[j]->_elements.begin(),
+                                               _subMeshes[j]->_elements.end());
+        elems.insert( m->_subMeshes[j]->_elements.begin(),
+                      m->_subMeshes[j]->_elements.end());
+        _subMeshes[j]->_elements.assign( elems.begin(), elems.end() );
+        m->_subMeshes[j]->_elements.clear();
+
+        if ( !_subMeshes[j]->_n2n )
+          _subMeshes[j]->_n2n = m->_subMeshes[j]->_n2n, m->_subMeshes[j]->_n2n = 0;
+
+        else if ( _subMeshes[j]->_n2n && m->_subMeshes[j]->_n2n )
+          _subMeshes[j]->_n2n->insert( m->_subMeshes[j]->_n2n->begin(),
+                                       m->_subMeshes[j]->_n2n->end());
+      }
+      else
+      {
+        _subMeshes[j] = m->_subMeshes[j];
+        m->_subMeshes[j] = 0;
+      }
+    }
+  }
+}
+
+//================================================================================
+/*!
+ * \brief Destructor deletes proxy submeshes and tmp elemens
+ */
+//================================================================================
+
+SMESH_ProxyMesh::~SMESH_ProxyMesh()
+{
+  for ( unsigned i = 0; i < _subMeshes.size(); ++i )
+    delete _subMeshes[i];
+  _subMeshes.clear();
+
+  set< const SMDS_MeshElement* >::iterator i = _elemsInMesh.begin();
+  for ( ; i != _elemsInMesh.end(); ++i )
+    GetMeshDS()->RemoveFreeElement( *i, 0 );
+  _elemsInMesh.clear();
+}
+
+//================================================================================
+/*!
+ * \brief Returns index of a shape
+ */
+//================================================================================
+
+int SMESH_ProxyMesh::shapeIndex(const TopoDS_Shape& shape) const
+{
+  return ( shape.IsNull() || !_mesh->HasShapeToMesh() ? 0 : GetMeshDS()->ShapeToIndex(shape));
+}
+
+//================================================================================
+/*!
+ * \brief Returns the submesh of a shape; it can be a proxy sub-mesh
+ */
+//================================================================================
+
+const SMESHDS_SubMesh* SMESH_ProxyMesh::GetSubMesh(const TopoDS_Shape& shape) const
+{
+  const SMESHDS_SubMesh* sm = 0;
+
+  int i = shapeIndex(shape);
+  if ( i < _subMeshes.size() )
+    sm = _subMeshes[i];
+  if ( !sm )
+    sm = GetMeshDS()->MeshElements( i );
+
+  return sm;
+}
+
+//================================================================================
+/*!
+ * \brief Returns the proxy sub-mesh of a shape; it can be NULL
+ */
+//================================================================================
+
+const SMESH_ProxyMesh::SubMesh*
+SMESH_ProxyMesh::GetProxySubMesh(const TopoDS_Shape& shape) const
+{
+  int i = shapeIndex(shape);
+  return i < _subMeshes.size() ? _subMeshes[i] : 0;
+}
+
+//================================================================================
+/*!
+ * \brief Returns the proxy node of a node; the input node is returned if no proxy exists
+ */
+//================================================================================
+
+const SMDS_MeshNode* SMESH_ProxyMesh::GetProxyNode( const SMDS_MeshNode* node ) const
+{
+  const SMDS_MeshNode* proxy = node;
+  if ( node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE )
+  {
+    if ( const SubMesh* proxySM = findProxySubMesh( node->getshapeId() ))
+      proxy = proxySM->GetProxyNode( node );
+  }
+  else
+  {
+    TopoDS_Shape shape = SMESH_MesherHelper::GetSubShapeByNode( node, GetMeshDS());
+    TopTools_ListIteratorOfListOfShape ancIt;
+    if ( !shape.IsNull() ) ancIt.Initialize( _mesh->GetAncestors( shape ));
+    for ( ; ancIt.More() && proxy == node; ancIt.Next() )
+      if ( const SubMesh* proxySM = findProxySubMesh( shapeIndex(ancIt.Value())))
+        proxy = proxySM->GetProxyNode( node );
+  }
+  return proxy;
+}
+
+namespace
+{
+  //================================================================================
+  /*!
+   * \brief Iterator filtering elements by type
+   */
+  //================================================================================
+
+  class TFilteringIterator : public SMDS_ElemIterator
+  {
+    SMDS_ElemIteratorPtr        _iter;
+    const SMDS_MeshElement *    _curElem;
+    vector< SMDSAbs_EntityType> _okTypes;
+  public:
+    TFilteringIterator( const vector< SMDSAbs_EntityType>& okTypes,
+                        const SMDS_ElemIteratorPtr&        elemIterator)
+      :_iter(elemIterator), _curElem(0), _okTypes(okTypes)
+    {
+      next();
+    }
+    virtual bool more()
+    {
+      return _curElem;
+    }
+    virtual const SMDS_MeshElement* next()
+    {
+      const SMDS_MeshElement* res = _curElem;
+      _curElem = 0;
+      while ( _iter->more() && !_curElem )
+      {
+        _curElem = _iter->next();
+        if ( find( _okTypes.begin(), _okTypes.end(), _curElem->GetEntityType()) == _okTypes.end())
+          _curElem = 0;
+      }
+      return res;
+    }
+  };
+}
+
+//================================================================================
+/*!
+ * \brief Returns iterator on all faces on the shape taking into account substitutions
+ */
+//================================================================================
+
+SMDS_ElemIteratorPtr SMESH_ProxyMesh::GetFaces(const TopoDS_Shape& shape) const
+{
+  if ( !_mesh->HasShapeToMesh() )
+    return SMDS_ElemIteratorPtr();
+
+  _subContainer.RemoveAllSubmeshes();
+
+  TopTools_IndexedMapOfShape FF;
+  TopExp::MapShapes( shape, TopAbs_FACE, FF );
+  for ( int i = 1; i <= FF.Extent(); ++i )
+    if ( const SMESHDS_SubMesh* sm = GetSubMesh( FF(i)))
+      _subContainer.AddSubMesh( sm );
+
+  return _subContainer.SMESHDS_SubMesh::GetElements();
+}
+
+//================================================================================
+/*!
+ * \brief Returns iterator on all faces of the mesh taking into account substitutions
+ * To be used in case of mesh without shape
+ */
+//================================================================================
+
+SMDS_ElemIteratorPtr SMESH_ProxyMesh::GetFaces() const
+{
+  if ( _mesh->HasShapeToMesh() )
+    return SMDS_ElemIteratorPtr();
+
+  _subContainer.RemoveAllSubmeshes();
+  for ( unsigned i = 0; i < _subMeshes.size(); ++i )
+    if ( _subMeshes[i] )
+      _subContainer.AddSubMesh( _subMeshes[i] );
+
+  if ( _subContainer.NbSubMeshes() == 0 ) // no elements substituted
+    return GetMeshDS()->elementsIterator(SMDSAbs_Face);
+
+  // if _allowedTypes is empty, only elements from _subMeshes are returned,...
+  SMDS_ElemIteratorPtr proxyIter = _subContainer.SMESHDS_SubMesh::GetElements();
+  if ( _allowedTypes.empty() || NbFaces() == _mesh->NbFaces() )
+    return proxyIter;
+
+  // ... else elements filtered using allowedTypes are additionally returned
+  SMDS_ElemIteratorPtr facesIter = GetMeshDS()->elementsIterator(SMDSAbs_Face);
+  SMDS_ElemIteratorPtr filterIter( new TFilteringIterator( _allowedTypes, facesIter ));
+  vector< SMDS_ElemIteratorPtr > iters(2);
+  iters[0] = proxyIter;
+  iters[1] = filterIter;
+    
+  typedef vector< SMDS_ElemIteratorPtr > TElemIterVector;
+  typedef SMDS_IteratorOnIterators<const SMDS_MeshElement *, TElemIterVector> TItersIter;
+  return SMDS_ElemIteratorPtr( new TItersIter( iters ));
+}
+
+//================================================================================
+/*!
+ * \brief Return total nb of faces taking into account substitutions
+ */
+//================================================================================
+
+int SMESH_ProxyMesh::NbFaces() const
+{
+  int nb = 0;
+  if ( _mesh->HasShapeToMesh() )
+  {
+    TopTools_IndexedMapOfShape FF;
+    TopExp::MapShapes( _mesh->GetShapeToMesh(), TopAbs_FACE, FF );
+    for ( int i = 1; i <= FF.Extent(); ++i )
+      if ( const SMESHDS_SubMesh* sm = GetSubMesh( FF(i)))
+        nb += sm->NbElements();
+  }
+  else
+  {
+    if ( _subMeshes.empty() )
+      return GetMeshDS()->NbFaces();
+
+    for ( unsigned i = 0; i < _subMeshes.size(); ++i )
+      if ( _subMeshes[i] )
+        nb += _subMeshes[i]->NbElements();
+
+    // if _allowedTypes is empty, only elements from _subMeshes are returned,
+    // else elements filtered using allowedTypes are additionally returned
+    if ( !_allowedTypes.empty() )
+    {
+      for ( int t = SMDSEntity_Triangle; t <= SMDSEntity_Quad_Quadrangle; ++t )
+      {
+        bool allowed =
+          ( find( _allowedTypes.begin(), _allowedTypes.end(), t ) != _allowedTypes.end() );
+        if ( allowed )
+          nb += GetMeshDS()->GetMeshInfo().NbEntities( SMDSAbs_EntityType( t ));
+      }
+    }
+  }
+  return nb;
+}
+
+//================================================================================
+/*!
+ * \brief Returns a proxy sub-mesh; it is created if not yet exists
+ */
+//================================================================================
+
+SMESH_ProxyMesh::SubMesh* SMESH_ProxyMesh::getProxySubMesh(int index)
+{
+  if ( int(_subMeshes.size()) <= index )
+    _subMeshes.resize( index+1, 0 );
+  if ( !_subMeshes[index] )
+    _subMeshes[index] = new SubMesh( index );
+  return _subMeshes[index];
+}
+
+//================================================================================
+/*!
+ * \brief Returns a proxy sub-mesh; it is created if not yet exists
+ */
+//================================================================================
+
+SMESH_ProxyMesh::SubMesh* SMESH_ProxyMesh::getProxySubMesh(const TopoDS_Shape& shape)
+{
+  return getProxySubMesh( shapeIndex( shape ));
+}
+
+//================================================================================
+/*!
+ * \brief Returns a proxy sub-mesh
+ */
+//================================================================================
+
+SMESH_ProxyMesh::SubMesh* SMESH_ProxyMesh::findProxySubMesh(int shapeIndex) const
+{
+  return shapeIndex < int(_subMeshes.size()) ? _subMeshes[shapeIndex] : 0;
+}
+
+//================================================================================
+/*!
+ * \brief Returns mesh DS
+ */
+//================================================================================
+
+SMESHDS_Mesh* SMESH_ProxyMesh::GetMeshDS() const
+{
+  return (SMESHDS_Mesh*)( _mesh ? _mesh->GetMeshDS() : 0 );
+}
+
+//================================================================================
+/*!
+ * \brief Move proxy sub-mesh from other proxy mesh to this, returns true if sub-mesh found
+ */
+//================================================================================
+
+bool SMESH_ProxyMesh::takeProxySubMesh( const TopoDS_Shape&   shape,
+                                             SMESH_ProxyMesh* proxyMesh )
+{
+  if ( proxyMesh && proxyMesh->_mesh == _mesh )
+  {
+    int iS = shapeIndex( shape );
+    if ( SubMesh* sm = proxyMesh->findProxySubMesh( iS ))
+    {
+      if ( iS >= int(_subMeshes.size()) )
+        _subMeshes.resize( iS + 1, 0 );
+      _subMeshes[iS] = sm;
+      proxyMesh->_subMeshes[iS] = 0;
+      return true;
+    }
+  }
+  return false;
+}
+
+//================================================================================
+/*!
+ * \brief Move tmp elements residing the _mesh from other proxy mesh to this
+ */
+//================================================================================
+
+void SMESH_ProxyMesh::takeTmpElemsInMesh( SMESH_ProxyMesh* proxyMesh )
+{
+  if ( proxyMesh )
+  {
+    _elemsInMesh.insert( proxyMesh->_elemsInMesh.begin(),
+                         proxyMesh->_elemsInMesh.end());
+    proxyMesh->_elemsInMesh.clear();
+  }
+}
+
+//================================================================================
+/*!
+ * \brief Removes tmp faces from the _mesh
+ */
+//================================================================================
+
+void SMESH_ProxyMesh::removeTmpElement( const SMDS_MeshElement* face )
+{
+  if ( face && face->GetID() > 0 )
+  {
+    set< const SMDS_MeshElement* >::iterator i =  _elemsInMesh.find( face );
+    if ( i != _elemsInMesh.end() )
+    {
+      GetMeshDS()->RemoveFreeElement( face, 0 );
+      _elemsInMesh.erase( i );
+    }
+  }
+  else
+  {
+    delete face;
+  }
+}
+
+//================================================================================
+/*!
+ * \brief Stores tmp element residing the _mesh
+ */
+//================================================================================
+
+void SMESH_ProxyMesh::storeTmpElement( const SMDS_MeshElement* face )
+{
+  _elemsInMesh.insert( face );
+}
+
+//================================================================================
+/*!
+ * \brief Set node-node correspondence
+ */
+//================================================================================
+
+void SMESH_ProxyMesh::setNode2Node(const SMDS_MeshNode* srcNode,
+                                        const SMDS_MeshNode* proxyNode,
+                                        const SubMesh*       subMesh)
+{
+  SubMesh* sm = const_cast<SubMesh*>( subMesh );
+  if ( !subMesh->_n2n )
+    sm->_n2n = new TN2NMap;
+  sm->_n2n->insert( make_pair( srcNode, proxyNode ));
+}
+
+//================================================================================
+/*!
+ * \brief Return true if the element is a temporary one
+ */
+//================================================================================
+
+bool SMESH_ProxyMesh::IsTemporary(const SMDS_MeshElement* elem ) const
+{
+  return ( elem->GetID() < 1 ) || _elemsInMesh.count( elem );
+}
+
+//================================================================================
+/*!
+ * \brief Return a proxy node or an input node
+ */
+//================================================================================
+
+const SMDS_MeshNode* SMESH_ProxyMesh::SubMesh::GetProxyNode( const SMDS_MeshNode* n ) const
+{
+  TN2NMap::iterator n2n;
+  if ( _n2n && ( n2n = _n2n->find( n )) != _n2n->end())
+    return n2n->second;
+  return n;
+}
+
+//================================================================================
+/*!
+ * \brief Deletes temporary elements
+ */
+//================================================================================
+
+void SMESH_ProxyMesh::SubMesh::Clear()
+{
+  for ( unsigned i = 0; i < _elements.size(); ++i )
+    if ( _elements[i]->GetID() < 0 )
+      delete _elements[i];
+  _elements.clear();
+  if ( _n2n )
+    delete _n2n, _n2n = 0;
+}
+
+//================================================================================
+/*!
+ * \brief Return number of elements in a proxy submesh
+ */
+//================================================================================
+
+int SMESH_ProxyMesh::SubMesh::NbElements() const
+{
+  return _elements.size();
+}
+
+//================================================================================
+/*!
+ * \brief Return elements of a proxy submesh
+ */
+//================================================================================
+
+SMDS_ElemIteratorPtr SMESH_ProxyMesh::SubMesh::GetElements() const
+{
+  return SMDS_ElemIteratorPtr
+    ( new SMDS_ElementVectorIterator( _elements.begin(), _elements.end() ));
+}
+
+//================================================================================
+/*!
+ * \brief Store an element
+ */
+//================================================================================
+
+void SMESH_ProxyMesh::SubMesh::AddElement(const SMDS_MeshElement * e)
+{
+  _elements.push_back( e );
+}
+
+//================================================================================
+/*!
+ * \brief Check presence of element inside it-self
+ */
+//================================================================================
+
+bool SMESH_ProxyMesh::SubMesh::Contains(const SMDS_MeshElement * ME) const
+{
+  if ( ME->GetType() != SMDSAbs_Node )
+    return find( _elements.begin(), _elements.end(), ME ) != _elements.end();
+  return false;
+}
diff --git a/src/SMESH/SMESH_ProxyMesh.hxx b/src/SMESH/SMESH_ProxyMesh.hxx
new file mode 100644 (file)
index 0000000..3bcd18d
--- /dev/null
@@ -0,0 +1,173 @@
+//  Copyright (C) 2007-2010  CEA/DEN, EDF R&D, OPEN CASCADE
+//
+//  This library is free software; you can redistribute it and/or
+//  modify it under the terms of the GNU Lesser General Public
+//  License as published by the Free Software Foundation; either
+//  version 2.1 of the License.
+//
+//  This library is distributed in the hope that it will be useful,
+//  but WITHOUT ANY WARRANTY; without even the implied warranty of
+//  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+//  Lesser General Public License for more details.
+//
+//  You should have received a copy of the GNU Lesser General Public
+//  License along with this library; if not, write to the Free Software
+//  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
+//
+//  See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
+//
+
+// File      : SMESH_ProxyMesh.hxx
+// Created   : Thu Dec  2 10:05:35 2010
+// Author    : Edward AGAPOV (eap)
+
+#ifndef __SMESH_ProxyMesh_HXX__
+#define __SMESH_ProxyMesh_HXX__
+
+#include "SMESH_SMESH.hxx"
+
+#include "SMDS_MeshElement.hxx"
+#include "SMESHDS_SubMesh.hxx"
+
+#include <TopoDS_Shape.hxx>
+
+#include <map>
+#include <vector>
+#include <boost/shared_ptr.hpp>
+
+class SMDS_MeshNode;
+class SMESHDS_Mesh;
+class SMESH_Mesh;
+
+/*!
+ * \brief Container of mesh faces substituting other faces in the input mesh of 3D algorithm
+ */
+class SMESH_EXPORT SMESH_ProxyMesh
+{
+public:
+
+  typedef boost::shared_ptr<SMESH_ProxyMesh> Ptr;
+
+  typedef std::map<const SMDS_MeshNode*, const SMDS_MeshNode*, TIDCompare > TN2NMap;
+
+  //--------------------------------------------------------------------------------
+  /*!
+   * \brief Proxy sub-mesh
+   */
+  class SubMesh : public SMESHDS_SubMesh
+  {
+  public:
+
+    const TN2NMap* GetNodeNodeMap() const { return _n2n; }
+    const SMDS_MeshNode* GetProxyNode( const SMDS_MeshNode* n ) const;
+    virtual void AddElement(const SMDS_MeshElement * e);
+    virtual int NbElements() const;
+    virtual SMDS_ElemIteratorPtr GetElements() const;
+    virtual void Clear();
+    virtual bool Contains(const SMDS_MeshElement * ME) const;
+
+    template< class ITERATOR >
+    void ChangeElements( ITERATOR it, ITERATOR end )
+    {
+      // change SubMesh contents without deleting tmp faces
+      // for which the caller is responsible
+      _elements.clear();
+      while ( it != end ) _elements.push_back( *it++ );
+    }
+    SubMesh(int index=0):SMESHDS_SubMesh(0,index),_n2n(0) {}
+    ~SubMesh() { Clear(); }
+
+  private:
+    std::vector<const SMDS_MeshElement *> _elements;
+    TN2NMap*                              _n2n;
+    friend class SMESH_ProxyMesh;
+  };
+  //--------------------------------------------------------------------------------
+  // Public interface
+
+  SMESH_ProxyMesh();
+  SMESH_ProxyMesh(std::vector<SMESH_ProxyMesh::Ptr>& components);
+  SMESH_ProxyMesh(const SMESH_Mesh& mesh) { _mesh = &mesh; }
+  virtual ~SMESH_ProxyMesh();
+
+  // Returns the submesh of a face; it can be a proxy sub-mesh
+  const SMESHDS_SubMesh* GetSubMesh(const TopoDS_Shape& face) const;
+
+  // Returns the proxy sub-mesh of a face; it can be NULL
+  const SubMesh* GetProxySubMesh(const TopoDS_Shape& face) const;
+
+  // Returns the proxy node of a node; the input node is returned if no proxy exists
+  const SMDS_MeshNode* GetProxyNode( const SMDS_MeshNode* node ) const;
+
+  // Returns iterator on all faces of the mesh taking into account substitutions
+  // To be used in case of mesh without shape
+  SMDS_ElemIteratorPtr GetFaces() const;
+
+  // Returns iterator on all faces on the face taking into account substitutions
+  SMDS_ElemIteratorPtr GetFaces(const TopoDS_Shape& face) const;
+
+  // Return total nb of faces taking into account substitutions
+  int NbFaces() const;
+
+  bool IsTemporary(const SMDS_MeshElement* elem ) const;
+
+
+
+  const SMESH_Mesh* GetMesh() const { return _mesh; }
+
+  SMESHDS_Mesh* GetMeshDS() const;
+
+  //--------------------------------------------------------------------------------
+  // Interface for descendants
+ protected:
+
+  void setMesh(const SMESH_Mesh& mesh) { _mesh = &mesh; }
+
+  int shapeIndex(const TopoDS_Shape& shape) const;
+
+  // returns a proxy sub-mesh; zero index is for the case of mesh w/o shape
+  SubMesh* findProxySubMesh(int shapeIndex=0) const;
+
+  // returns a proxy sub-mesh; it is created if not yet exists
+  SubMesh* getProxySubMesh(int shapeIndex);
+
+  // returns a proxy sub-mesh; it is created if not yet exists
+  SubMesh* getProxySubMesh(const TopoDS_Shape& shape=TopoDS_Shape());
+
+  // move proxy sub-mesh from other proxy mesh to this, returns true if sub-mesh found
+  bool takeProxySubMesh( const TopoDS_Shape& shape, SMESH_ProxyMesh* proxyMesh );
+
+  // move tmp elements residing the _mesh from other proxy mesh to this
+  void takeTmpElemsInMesh( SMESH_ProxyMesh* proxyMesh );
+
+  // removes tmp faces from the _mesh
+  void removeTmpElement( const SMDS_MeshElement* face );
+
+  // stores tmp element residing the _mesh
+  void storeTmpElement( const SMDS_MeshElement* face );
+
+  // store node-node correspondence
+  void setNode2Node(const SMDS_MeshNode* srcNode,
+                    const SMDS_MeshNode* proxyNode,
+                    const SubMesh*       subMesh);
+
+  // types of elements needed to implement NbFaces() and GetFaces();
+  // if _allowedTypes is empty, only elements from _subMeshes are returned,
+  // else elements of _mesh filtered using allowedTypes are additionally returned
+  std::vector< SMDSAbs_EntityType> _allowedTypes;
+
+ private:
+
+  const SMESH_Mesh*       _mesh;
+
+  // proxy sub-meshes; index in vector == shapeIndex(shape)
+  std::vector< SubMesh* > _subMeshes;
+
+  // tmp elements residing the _mesh, to be deleted at destruction
+  std::set< const SMDS_MeshElement* > _elemsInMesh;
+
+  // Complex submesh used to iterate over elements in other sub-meshes
+  mutable SubMesh _subContainer;
+};
+
+#endif