]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
Merge branch 'V9_6_BR'
authorvsr <vsr@opencascade.com>
Wed, 11 Nov 2020 15:57:11 +0000 (18:57 +0300)
committervsr <vsr@opencascade.com>
Wed, 11 Nov 2020 15:57:46 +0000 (18:57 +0300)
CMakeLists.txt
doc/salome/gui/SMESH/input/define_mesh_by_script.rst
src/SMESH/SMESH_Algo.cxx
src/SMESH/SMESH_MeshEditor.cxx
src/SMESH_I/SMESH_MeshEditor_i.cxx
src/SMESH_I/SMESH_Mesh_i.cxx
src/SMESH_SWIG/smeshBuilder.py
src/StdMeshers/StdMeshers_Import_1D.cxx
src/StdMeshers/StdMeshers_Prism_3D.cxx

index 6e6839abd3aeef2b325a8e7fba5c51fef88879ea..adfdf638c5a2854f0191ecadd1cc440668758614 100644 (file)
@@ -27,7 +27,7 @@ INCLUDE(CMakeDependentOption)
 STRING(TOUPPER ${PROJECT_NAME} PROJECT_NAME_UC)
 
 SET(${PROJECT_NAME_UC}_MAJOR_VERSION 9)
-SET(${PROJECT_NAME_UC}_MINOR_VERSION 5)
+SET(${PROJECT_NAME_UC}_MINOR_VERSION 6)
 SET(${PROJECT_NAME_UC}_PATCH_VERSION 0)
 SET(${PROJECT_NAME_UC}_VERSION
   ${${PROJECT_NAME_UC}_MAJOR_VERSION}.${${PROJECT_NAME_UC}_MINOR_VERSION}.${${PROJECT_NAME_UC}_PATCH_VERSION})
index 92387778c8580b7397244ee6b596519b744c3d6b..43d7f89afbc0f0a92215b1dd3f0146c7c26b7bb4 100644 (file)
@@ -15,7 +15,9 @@ meshes and to create 2D mesh by your python code. For this, you
 #. run your python code, which creates a 2D mesh,
 #. invoke **Compute** command, which computes a 3D mesh.
 
-.. warning:: **Use Edges to be Created Manually** and **Use Faces to be Created Manually** algorithms should be assigned *before* mesh generation by the Python code.
+  .. warning::
+     * **Use Edges to be Created Manually** and **Use Faces to be Created Manually** algorithms should be assigned *before* mesh generation by the Python code.
+     * Nodes and elements created in your script must be assigned to geometry entities by calling *SetMeshElementOnShape*, *SetNodeOnVertex*, *SetNodeOnEdge* etc. in order to be used by an algorithm of upper dimension.
 
 Consider trying a sample script demonstrating the usage of :ref:`Use Faces to be Created Manually <tui_use_existing_faces>` algorithm for construction of a 2D mesh using Python commands.
 
index d4930818e22508fbb5daed1a4ff2242a91702e11..bdc93b1a12322f164fcac8eec8c1864ff88dbd70 100644 (file)
@@ -950,10 +950,9 @@ void SMESH_Algo::InitComputeError()
 {
   _error = COMPERR_OK;
   _comment.clear();
-  list<const SMDS_MeshElement*>::iterator elem = _badInputElements.begin();
-  for ( ; elem != _badInputElements.end(); ++elem )
-    if ( (*elem)->GetID() < 1 )
-      delete *elem;
+  for ( const SMDS_MeshElement* & elem : _badInputElements )
+    if ( !elem->IsNull() && elem->GetID() < 1 )
+      delete elem;
   _badInputElements.clear();
   _mesh = 0;
 
@@ -1015,6 +1014,7 @@ void SMESH_Algo::addBadInputElements(const SMESHDS_SubMesh* sm,
       SMDS_ElemIteratorPtr eIt = sm->GetElements();
       while ( eIt->more() ) addBadInputElement( eIt->next() );
     }
+    _mesh = sm->GetParent();
   }
 }
 
index ae3b70321c91196a98ba745dfa86b1466f83824d..1d679c2bc0dd3cf03287b93241a0e2edb9280205 100644 (file)
@@ -1504,9 +1504,12 @@ void SMESH_MeshEditor::QuadTo4Tri (TIDSortedElemSet & theElems)
   SMESH_MesherHelper helper( *GetMesh() );
   helper.SetElementsOnShape( true );
 
-  SMDS_ElemIteratorPtr faceIt;
-  if ( theElems.empty() ) faceIt = GetMeshDS()->elementsIterator(SMDSAbs_Face);
-  else                    faceIt = SMESHUtils::elemSetIterator( theElems );
+  // get standalone groups of faces
+  vector< SMDS_MeshGroup* > allFaceGroups, faceGroups;
+  for ( SMESHDS_GroupBase* grBase : GetMeshDS()->GetGroups() )
+    if ( SMESHDS_Group* group = dynamic_cast<SMESHDS_Group*>( grBase ))
+      if ( group->GetType() == SMDSAbs_Face && !group->IsEmpty() )
+        allFaceGroups.push_back( & group->SMDSGroup() );
 
   bool   checkUV;
   gp_XY  uv [9]; uv[8] = gp_XY(0,0);
@@ -1517,6 +1520,10 @@ void SMESH_MeshEditor::QuadTo4Tri (TIDSortedElemSet & theElems)
   Handle(Geom_Surface)           surface;
   TopLoc_Location                loc;
 
+  SMDS_ElemIteratorPtr faceIt;
+  if ( theElems.empty() ) faceIt = GetMeshDS()->elementsIterator(SMDSAbs_Face);
+  else                    faceIt = SMESHUtils::elemSetIterator( theElems );
+
   while ( faceIt->more() )
   {
     const SMDS_MeshElement* quad = faceIt->next();
@@ -1593,13 +1600,19 @@ void SMESH_MeshEditor::QuadTo4Tri (TIDSortedElemSet & theElems)
       myLastCreatedNodes.push_back( nCentral );
     }
 
-    // create 4 triangles
-
     helper.SetIsQuadratic  ( nodes.size() > 4 );
     helper.SetIsBiQuadratic( nodes.size() == 9 );
     if ( helper.GetIsQuadratic() )
       helper.AddTLinks( static_cast< const SMDS_MeshFace*>( quad ));
 
+    // select groups to update
+    faceGroups.clear();
+    for ( SMDS_MeshGroup* group : allFaceGroups )
+      if ( group->Remove( quad ))
+        faceGroups.push_back( group );
+
+    // create 4 triangles
+
     GetMeshDS()->RemoveFreeElement( quad, subMeshDS, /*fromGroups=*/false );
 
     for ( int i = 0; i < 4; ++i )
@@ -1607,8 +1620,9 @@ void SMESH_MeshEditor::QuadTo4Tri (TIDSortedElemSet & theElems)
       SMDS_MeshElement* tria = helper.AddFace( nodes[ i ],
                                                nodes[(i+1)%4],
                                                nCentral );
-      ReplaceElemInGroups( tria, quad, GetMeshDS() );
       myLastCreatedElems.push_back( tria );
+      for ( SMDS_MeshGroup* group : faceGroups )
+        group->Add( tria );
     }
   }
 }
index 16cc20053d24e05aa276992c8d0193d0c6a23c46..3330950a2acf244f0ebda8d92d64dfc560b91a93 100644 (file)
@@ -1966,6 +1966,8 @@ void SMESH_MeshEditor_i::QuadTo4Tri (SMESH::SMESH_IDSource_ptr theObject)
   getEditor().QuadTo4Tri( faces );
   TPythonDump() << this << ".QuadTo4Tri( " << theObject << " )";
 
+  declareMeshModified( /*isReComputeSafe=*/false );
+
   SMESH_CATCH( SMESH::throwCorbaException );
 }
 
index 8a4313de59b324bccd92a775f058d82c58e3fcad..8455f6b4681a3886ab21329db943221ea44f882f 100644 (file)
@@ -649,11 +649,12 @@ SMESH_Mesh_i::AddHypothesis(GEOM::GEOM_Object_ptr       aSubShape,
   throw(SALOME::SALOME_Exception)
 {
   Unexpect aCatch(SALOME_SalomeException);
+
+  const int prevNbMeshEnt = NbNodes() + NbElements();
+
   if ( _preMeshInfo )
     _preMeshInfo->ForgetOrLoad();
 
-  const int prevNbMeshEnt = _impl->NbNodes() + _impl->GetMeshDS()->NbElements();
-
   std::string error;
   SMESH_Hypothesis::Hypothesis_Status status = addHypothesis( aSubShape, anHyp, &error );
   anErrorText = error.c_str();
index b40bc6d1eea57f812cdb8ec5e10a3d48313cd31e..c99e958214e6c0d16b649cb24eeb24c9e4cf9d11 100644 (file)
@@ -4450,7 +4450,7 @@ class Mesh(metaclass = MeshMeta):
 
         Returns:
              A list of edge groups and a list of corresponding node groups,
-             where the group is a list of IDs of edges or elements, like follows
+             where the group is a list of IDs of edges or nodes, like follows
              [[[branch_edges_1],[branch_edges_2]], [[branch_nodes_1],[branch_nodes_2]]].
              If a group is closed, the first and last nodes of the group are same.
         """
index e3d59a79c8e1fc3195d5319fe7b3da575481db24..af0ea24bd113dc76c3358d4bc0c2b354850f60a5 100644 (file)
@@ -200,22 +200,37 @@ namespace // INTERNAL STUFF
     if ( getBox()->IsOut( point ))
       return false;
 
+    bool ok = false;
+    double dist2, param;
+    distance2 = Precision::Infinite();
+
     if ( isLeaf() )
     {
       for ( size_t i = 0; i < _segments.size(); ++i )
         if ( !_segments[i].IsOut( point ) &&
-             _segments[i].IsOn( point, distance2, u ))
-          return true;
+             _segments[i].IsOn( point, dist2, param ) &&
+             dist2 < distance2 )
+        {
+          distance2 = dist2;
+          u         = param;
+          ok        = true;
+        }
+      return ok;
     }
     else
     {
       for (int i = 0; i < 8; i++)
-        if (((CurveProjector*) myChildren[i])->IsOnCurve( point, distance2, u ))
-          return true;
+        if (((CurveProjector*) myChildren[i])->IsOnCurve( point, dist2, param ) &&
+            dist2 < distance2 )
+        {
+          distance2 = dist2;
+          u         = param;
+          ok        = true;
+        }
     }
-    return false;
+    return ok;
   }
-  
+
   //================================================================================
   /*!
    * \brief Initialize
@@ -233,13 +248,13 @@ namespace // INTERNAL STUFF
     _pLast   = pl;
     _curve   = curve;
     _length2 = pf.SquareDistance( pl );
+    _line.SetLocation( pf );
+    _line.SetDirection( gp_Vec( pf, pl ));
     _chord2  = Max( _line.     SquareDistance( curve->Value( uf + 0.25 * ( ul - uf ))),
                     Max( _line.SquareDistance( curve->Value( uf + 0.5  * ( ul - uf ))),
                          _line.SquareDistance( curve->Value( uf + 0.75 * ( ul - uf )))));
     _chord2  = Max( tol, _chord2 );
     _chord   = Sqrt( _chord2 );
-    _line.SetLocation( pf );
-    _line.SetDirection( gp_Vec( pf, pl ));
 
     Bnd_Box bb;
     BndLib_Add3dCurve::Add( GeomAdaptor_Curve( curve, uf, ul ), tol, bb );
@@ -268,12 +283,12 @@ namespace // INTERNAL STUFF
       gp_Vec edge( _pFirst, _pLast );
       gp_Vec n1p ( _pFirst, point  );
       u = ( edge * n1p ) / _length2; // param [0,1] on the edge
-      if ( u < 0 )
+      if ( u < 0. )
       {
         if ( _pFirst.SquareDistance( point ) > _chord2 )
           return false;
       }
-      else if ( u > _chord )
+      else if ( u > 1. )
       {
         if ( _pLast.SquareDistance( point ) > _chord2 )
           return false;
index 4f4cf213b8effd493797feba0d5f712d85f4411c..e6db9e5014f8fe75ce40c6642fbbe689733d409b 100644 (file)
@@ -89,20 +89,48 @@ enum { ID_BOT_FACE = SMESH_Block::ID_Fxy0,
        BOTTOM_EDGE = 0, TOP_EDGE, V0_EDGE, V1_EDGE, // edge IDs in face
        NB_WALL_FACES = 4 }; //
 
-namespace {
-
+namespace
+{
+  //=======================================================================
+  /*!
+   * \brief Auxiliary mesh
+   */
+  struct TmpMesh: public SMESH_Mesh
+  {
+    TmpMesh() {
+      _isShapeToMesh = (_id = 0);
+      _myMeshDS  = new SMESHDS_Mesh( _id, true );
+    }
+  };
   //=======================================================================
   /*!
    * \brief Quadrangle algorithm
    */
-  struct TQuadrangleAlgo : public StdMeshers_Quadrangle_2D
+  class TQuadrangleAlgo : public StdMeshers_Quadrangle_2D
   {
+    typedef NCollection_DataMap< TopoDS_Face, FaceQuadStruct::Ptr > TFace2QuadMap;
+    TFace2QuadMap myFace2QuadMap;
+
     TQuadrangleAlgo(SMESH_Gen* gen)
       : StdMeshers_Quadrangle_2D( gen->GetANewId(), gen)
     {
     }
-    static StdMeshers_Quadrangle_2D* instance( SMESH_Algo*         fatherAlgo,
-                                               SMESH_MesherHelper* helper=0)
+  public:
+
+    //================================================================================
+    // Clear data of TQuadrangleAlgo at destruction
+    struct Cleaner
+    {
+      TQuadrangleAlgo* myAlgo;
+
+      Cleaner(TQuadrangleAlgo* algo): myAlgo( algo ){}
+      ~Cleaner() { myAlgo->reset(); }
+    };
+
+    //================================================================================
+    // Return TQuadrangleAlgo singleton
+    static TQuadrangleAlgo* instance( SMESH_Algo*         fatherAlgo,
+                                      SMESH_MesherHelper* helper=0)
     {
       static TQuadrangleAlgo* algo = new TQuadrangleAlgo( fatherAlgo->GetGen() );
       if ( helper &&
@@ -118,7 +146,123 @@ namespace {
 
       return algo;
     }
+
+    //================================================================================
+    // Clear collected data
+    void reset()
+    {
+      myFace2QuadMap.Clear();
+      StdMeshers_Quadrangle_2D::myQuadList.clear();
+      StdMeshers_Quadrangle_2D::myHelper = nullptr;
+    }
+
+    //================================================================================
+    /*!
+     * \brief Return FaceQuadStruct if a given FACE can be meshed by StdMeshers_Quadrangle_2D
+     */
+    FaceQuadStruct::Ptr CheckNbEdges(SMESH_Mesh&         theMesh,
+                                     const TopoDS_Shape& theShape )
+    {
+      const TopoDS_Face& face = TopoDS::Face( theShape );
+      if ( myFace2QuadMap.IsBound( face ))
+        return myFace2QuadMap.Find( face );
+
+      FaceQuadStruct::Ptr &  resultQuad = * myFace2QuadMap.Bound( face, FaceQuadStruct::Ptr() );
+
+      FaceQuadStruct::Ptr quad =
+        StdMeshers_Quadrangle_2D::CheckNbEdges( theMesh, face, /*considerMesh=*/false, myHelper );
+      if ( quad )
+      {
+        // check if the quadrangle mesh would be valid
+
+        // check existing 1D mesh
+        // int nbSegments[4], i = 0;
+        // for ( FaceQuadStruct::Side & side : quad->side )
+        //   nbSegments[ i++ ] = side.grid->NbSegments();
+        // if ( nbSegments[0] > 0 && nbSegments[2] > 0 && nbSegments[0] != nbSegments[2] ||
+        //      nbSegments[1] > 0 && nbSegments[3] > 0 && nbSegments[1] != nbSegments[3] )
+        //   return resultQuad;
+
+        int nbEdges = 0;
+        for ( FaceQuadStruct::Side & side : quad->side )
+          nbEdges += side.grid->NbEdges();
+        if ( nbEdges == 4 )
+          return resultQuad = quad;
+
+        TmpMesh mesh;
+        mesh.ShapeToMesh( face );
+        SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
+        SMESH_MesherHelper helper( mesh );
+        helper.SetSubShape( face );
+        helper.SetElementsOnShape( true );
+
+        // create nodes on all VERTEX'es
+        for ( TopExp_Explorer vert( face, TopAbs_VERTEX ); vert.More(); vert.Next() )
+          mesh.GetSubMesh( vert.Current() )->ComputeStateEngine( SMESH_subMesh::COMPUTE );
+
+        FaceQuadStruct::Ptr tmpQuad( new FaceQuadStruct() );
+        tmpQuad->side.resize( 4 );
+
+        // divide quad sides into halves at least
+        const SMDS_MeshNode* node;
+        for ( int iDir = 0; iDir < 2; ++iDir )
+        {
+          StdMeshers_FaceSidePtr sides[2] = { quad->side[iDir], quad->side[iDir+2] };
+          std::map< double, const SMDS_MeshNode* > nodes[2];
+          for ( int iS : { 0, 1 } )
+          {
+            node = SMESH_Algo::VertexNode( sides[iS]->FirstVertex(), meshDS );
+            nodes[iS].insert( std::make_pair( 0, node ));
+            double curLen = 0;
+            for ( int iE = 1; iE < sides[iS]->NbEdges(); ++iE )
+            {
+              curLen += sides[iS]->EdgeLength( iE - 1 );
+              double u = curLen / sides[iS]->Length();
+              node = SMESH_Algo::VertexNode( sides[iS]->FirstVertex( iE ), meshDS );
+              nodes[iS  ].insert( std::make_pair( u, node ));
+              nodes[1-iS].insert( std::make_pair( u, nullptr ));
+            }
+            nodes[iS].insert( std::make_pair( 0.5, nullptr ));
+            node = SMESH_Algo::VertexNode( sides[iS]->LastVertex(), meshDS );
+            nodes[iS].insert( std::make_pair( 1, node ));
+          }
+
+          for ( int iS : { 0, 1 } )
+          {
+            UVPtStructVec sideNodes;
+            sideNodes.reserve( nodes[ iS ].size() );
+            for ( auto & u_node : nodes[ iS ])
+            {
+              if ( !u_node.second )
+              {
+                gp_Pnt p = sides[iS]->Value3d( u_node.first );
+                u_node.second = meshDS->AddNode( p.X(), p.Y(), p.Z() );
+                TopoDS_Edge edge;
+                double param = sides[iS]->Parameter( u_node.first, edge );
+                meshDS->SetNodeOnEdge( u_node.second, edge, param );
+              }
+              sideNodes.push_back( u_node.second );
+              sideNodes.back().SetUV( helper.GetNodeUV( face, u_node.second ));
+            }
+            tmpQuad->side[ iS ? iDir+2 : iDir ] = StdMeshers_FaceSide::New( sideNodes, face );
+          }
+        }
+        StdMeshers_Quadrangle_2D::myCheckOri = true;
+        StdMeshers_Quadrangle_2D::myQuadList.clear();
+        StdMeshers_Quadrangle_2D::myQuadList.push_back( tmpQuad );
+        StdMeshers_Quadrangle_2D::myHelper = &helper;
+        if ( StdMeshers_Quadrangle_2D::computeQuadDominant( mesh, face, tmpQuad ) &&
+             StdMeshers_Quadrangle_2D::check())
+        {
+          resultQuad = quad;
+        }
+        StdMeshers_Quadrangle_2D::myQuadList.clear();
+        StdMeshers_Quadrangle_2D::myHelper = nullptr;
+      }
+      return resultQuad;
+    }
   };
+
   //=======================================================================
   /*!
    * \brief Algorithm projecting 1D mesh
@@ -400,7 +544,7 @@ namespace {
 
   int removeQuasiQuads(list< SMESH_subMesh* >&   notQuadSubMesh,
                        SMESH_MesherHelper*       helper,
-                       StdMeshers_Quadrangle_2D* quadAlgo)
+                       TQuadrangleAlgo*          quadAlgo)
   {
     int nbRemoved = 0;
     //SMESHDS_Mesh* mesh = notQuadSubMesh.front()->GetFather()->GetMeshDS();
@@ -528,6 +672,22 @@ namespace {
     return nbSides;
   }
 
+  //================================================================================
+  /*!
+   * \brief Count EDGEs ignoring degenerated ones
+   */
+  //================================================================================
+
+  int CountEdges( const TopoDS_Face& face )
+  {
+    int nbE = 0;
+    for ( TopExp_Explorer edgeExp( face, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
+      if ( !SMESH_Algo::isDegenerated( TopoDS::Edge( edgeExp.Current() )))
+        ++nbE;
+
+    return nbE;
+  }
+
   //================================================================================
   /*!
    * \brief Set/get wire index to FaceQuadStruct
@@ -559,6 +719,7 @@ namespace {
     }
 #endif
   }
+
 } // namespace
 
 //=======================================================================
@@ -617,6 +778,7 @@ bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theSh
   SMESH_MesherHelper helper( theMesh );
   myHelper = &helper;
   myPrevBottomSM = 0;
+  TQuadrangleAlgo::Cleaner( TQuadrangleAlgo::instance( this ));
 
   int nbSolids = helper.Count( theShape, TopAbs_SOLID, /*skipSame=*/false );
   if ( nbSolids < 1 )
@@ -628,7 +790,6 @@ bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theSh
   // look for meshed FACEs ("source" FACEs) that must be prism bottoms
   list< TopoDS_Face > meshedFaces, notQuadMeshedFaces, notQuadFaces;
   const bool meshHasQuads = ( theMesh.NbQuadrangles() > 0 );
-  //StdMeshers_Quadrangle_2D* quadAlgo = TQuadrangleAlgo::instance( this );
   for ( int iF = 1; iF <= faceToSolids.Extent(); ++iF )
   {
     const TopoDS_Face& face = TopoDS::Face( faceToSolids.FindKey( iF ));
@@ -715,7 +876,8 @@ bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theSh
             return false;
 
           SMESHDS_SubMesh* smDS = theMesh.GetMeshDS()->MeshElements( prism.myTop );
-          if ( !myHelper->IsSameElemGeometry( smDS, SMDSGeom_QUADRANGLE ))
+          if ( !myHelper->IsSameElemGeometry( smDS, SMDSGeom_QUADRANGLE ) ||
+               !myHelper->IsStructured( theMesh.GetSubMesh( prism.myTop )))
           {
             meshedFaces.push_front( prism.myTop );
           }
@@ -919,7 +1081,7 @@ bool StdMeshers_Prism_3D::getWallFaces( Prism_3D::TPrismTopo & thePrism,
 
   SMESH_Mesh* mesh = myHelper->GetMesh();
 
-  StdMeshers_Quadrangle_2D* quadAlgo = TQuadrangleAlgo::instance( this, myHelper );
+  TQuadrangleAlgo* quadAlgo = TQuadrangleAlgo::instance( this, myHelper );
 
   TopTools_MapOfShape faceMap;
   TopTools_IndexedDataMapOfShapeListOfShape edgeToFaces;
@@ -1465,8 +1627,8 @@ bool StdMeshers_Prism_3D::computeWalls(const Prism_3D::TPrismTopo& thePrism)
   SMESHDS_Mesh* meshDS = myHelper->GetMeshDS();
   DBGOUT( endl << "COMPUTE Prism " << meshDS->ShapeToIndex( thePrism.myShape3D ));
 
-  TProjction1dAlgo*      projector1D = TProjction1dAlgo::instance( this );
-  StdMeshers_Quadrangle_2D* quadAlgo = TQuadrangleAlgo::instance( this, myHelper );
+  TProjction1dAlgo* projector1D = TProjction1dAlgo::instance( this );
+  TQuadrangleAlgo*     quadAlgo = TQuadrangleAlgo::instance( this, myHelper );
 
   // SMESH_HypoFilter hyp1dFilter( SMESH_HypoFilter::IsAlgo(),/*not=*/true);
   // hyp1dFilter.And( SMESH_HypoFilter::HasDim( 1 ));
@@ -2674,6 +2836,9 @@ bool StdMeshers_Prism_3D::allVerticalEdgesStraight( const Prism_3D::TPrismTopo&
 bool StdMeshers_Prism_3D::project2dMesh(const TopoDS_Face& theSrcFace,
                                         const TopoDS_Face& theTgtFace)
 {
+  if ( CountEdges( theSrcFace ) != CountEdges( theTgtFace ))
+    return false;
+
   TProjction2dAlgo* projector2D = TProjction2dAlgo::instance( this );
   projector2D->myHyp.SetSourceFace( theSrcFace );
   bool ok = projector2D->Compute( *myHelper->GetMesh(), theTgtFace );