Salome HOME
0020222: Quandrangle_2D meshing fail
[modules/smesh.git] / src / StdMeshers / StdMeshers_ProjectionUtils.cxx
index 3b69a201f4612ce648f92389c883be8da7489828..1fa7c45249d03f142d0be502c17a07df444ecba8 100644 (file)
@@ -1,32 +1,29 @@
-//  SMESH SMESH : idl implementation based on 'SMESH' unit's calsses
+//  Copyright (C) 2007-2008  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
 //
-//  Copyright (C) 2003  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
-//  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS 
-// 
-//  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
+//  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
+//
+//  SMESH SMESH : idl implementation based on 'SMESH' unit's calsses
 // File      : StdMeshers_ProjectionUtils.cxx
 // Created   : Fri Oct 27 10:24:28 2006
 // Author    : Edward AGAPOV (eap)
-
-using namespace std;
-
+//
 #include "StdMeshers_ProjectionUtils.hxx"
 
 #include "StdMeshers_ProjectionSource1D.hxx"
@@ -66,6 +63,8 @@ using namespace std;
 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
 #include <TopTools_DataMapIteratorOfDataMapOfShapeListOfShape.hxx>
 
+using namespace std;
+
 
 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
 #define SHOW_VERTEX(v,msg) // { \
@@ -86,7 +85,26 @@ using namespace std;
 //     cout << endl;\
 //   }
 
+#define HERE StdMeshers_ProjectionUtils
+
 namespace {
+
+  //================================================================================
+  /*!
+   * \brief Write shape for debug purposes
+   */
+  //================================================================================
+
+  bool _StoreBadShape(const TopoDS_Shape& shape)
+  {
+#ifdef _DEBUG_
+    const char* type[] ={"COMPOUND","COMPSOLID","SOLID","SHELL","FACE","WIRE","EDGE","VERTEX"};
+    BRepTools::Write( shape, SMESH_Comment("/tmp/") << type[shape.ShapeType()] << "_"
+                      << shape.TShape().operator->() << ".brep");
+#endif
+    return false;
+  }
+  
   //================================================================================
   /*!
    * \brief Reverse order of edges in a list and their orientation
@@ -163,8 +181,7 @@ namespace {
     if ( nbEdges == 2 && IsPropagationPossible( theMesh1, theMesh2 ) )
     {
       list< TopoDS_Edge >::iterator eIt2 = ++edges2.begin(); // 2nd edge of the 2nd face
-      TopoDS_Edge edge2 =
-        StdMeshers_ProjectionUtils::GetPropagationEdge( theMesh1, *eIt2, edges1.front() ).second;
+      TopoDS_Edge edge2 = HERE::GetPropagationEdge( theMesh1, *eIt2, edges1.front() ).second;
       if ( !edge2.IsNull() ) { // propagation found for the second edge
         Reverse( edges2, nbEdges );
         return true;
@@ -203,7 +220,135 @@ namespace {
     }
     return TopoDS_Shape();
   }
-}
+
+  //================================================================================
+  /*!
+   * \brief Find association of groups at top and bottom of prism
+   */
+  //================================================================================
+
+  bool AssocGroupsByPropagation(const TopoDS_Shape&   theGroup1,
+                                const TopoDS_Shape&   theGroup2,
+                                SMESH_Mesh&           theMesh,
+                                HERE::TShapeShapeMap& theMap)
+  {
+    // If groups are on top and bottom of prism then we can associate
+    // them using "vertical" (or "side") edges and faces of prism since
+    // they connect corresponding vertices and edges of groups.
+
+    TopTools_IndexedMapOfShape subshapes1, subshapes2;
+    TopExp::MapShapes( theGroup1, subshapes1 );
+    TopExp::MapShapes( theGroup2, subshapes2 );
+    TopTools_ListIteratorOfListOfShape ancestIt;
+
+    // Iterate on vertices of group1 to find corresponding vertices in group2
+    // and associate adjacent edges and faces
+
+    TopTools_MapOfShape verticShapes;
+    TopExp_Explorer vExp1( theGroup1, TopAbs_VERTEX );
+    for ( ; vExp1.More(); vExp1.Next() )
+    {
+      const TopoDS_Vertex& v1 = TopoDS::Vertex( vExp1.Current() );
+      if ( theMap.IsBound( v1 )) continue; // already processed
+
+      // Find "vertical" edge ending in v1 and whose other vertex belongs to group2
+      TopoDS_Shape verticEdge, v2;
+      ancestIt.Initialize( theMesh.GetAncestors( v1 ));
+      for ( ; verticEdge.IsNull() && ancestIt.More(); ancestIt.Next() )
+      {
+        if ( ancestIt.Value().ShapeType() != TopAbs_EDGE ) continue;
+        v2 = HERE::GetNextVertex( TopoDS::Edge( ancestIt.Value() ), v1 );
+        if ( subshapes2.Contains( v2 ))
+          verticEdge = ancestIt.Value();
+      }
+      if ( verticEdge.IsNull() )
+        return false;
+
+      HERE::InsertAssociation( v1, v2, theMap);
+
+      // Associate edges by vertical faces sharing the found vertical edge
+      ancestIt.Initialize( theMesh.GetAncestors( verticEdge ) );
+      for ( ; ancestIt.More(); ancestIt.Next() )
+      {
+        if ( ancestIt.Value().ShapeType() != TopAbs_FACE ) continue;
+        if ( !verticShapes.Add( ancestIt.Value() )) continue;
+        const TopoDS_Face& face = TopoDS::Face( ancestIt.Value() );
+
+        // get edges of the face
+        TopoDS_Edge edgeGr1, edgeGr2, verticEdge2;
+        list< TopoDS_Edge > edges;    list< int > nbEdgesInWire;
+        SMESH_Block::GetOrderedEdges( face, v1, edges, nbEdgesInWire);
+        if ( nbEdgesInWire.front() != 4 )
+          return _StoreBadShape( face );
+        list< TopoDS_Edge >::iterator edge = edges.begin();
+        if ( verticEdge.IsSame( *edge )) {
+          edgeGr2     = *(++edge);
+          verticEdge2 = *(++edge);
+          edgeGr1     = *(++edge);
+        } else {
+          edgeGr1     = *(edge++);
+          verticEdge2 = *(edge++);
+          edgeGr2     = *(edge++);
+        }
+
+        HERE::InsertAssociation( edgeGr1, edgeGr2.Reversed(), theMap);
+      }
+    }
+
+    // Associate faces
+    TopoDS_Iterator gr1It( theGroup1 );
+    if ( gr1It.Value().ShapeType() == TopAbs_FACE )
+    {
+      // find a boundary edge of group1 to start from
+      TopoDS_Shape bndEdge;
+      TopExp_Explorer edgeExp1( theGroup1, TopAbs_EDGE );
+      for ( ; bndEdge.IsNull() && edgeExp1.More(); edgeExp1.Next())
+        if ( HERE::IsBoundaryEdge( TopoDS::Edge( edgeExp1.Current()), theGroup1, theMesh ))
+          bndEdge = edgeExp1.Current();
+      if ( bndEdge.IsNull() )
+        return false;
+
+      list< TopoDS_Shape > edges(1, bndEdge);
+      list< TopoDS_Shape >::iterator edge1 = edges.begin();
+      for ( ; edge1 != edges.end(); ++edge1 )
+      {
+        // there must be one or zero not associated faces between ancestors of edge
+        // belonging to theGroup1
+        TopoDS_Shape face1;
+        ancestIt.Initialize( theMesh.GetAncestors( *edge1 ) );
+        for ( ; ancestIt.More() && face1.IsNull(); ancestIt.Next() ) {
+          if ( ancestIt.Value().ShapeType() == TopAbs_FACE &&
+               !theMap.IsBound( ancestIt.Value() ) &&
+               subshapes1.Contains( ancestIt.Value() ))
+            face1 = ancestIt.Value();
+
+          // add edges of face1 to start searching for adjacent faces from
+          for ( TopExp_Explorer e(face1, TopAbs_EDGE); e.More(); e.Next())
+            if ( !edge1->IsSame( e.Current() ))
+              edges.push_back( e.Current() );
+        }
+        if ( !face1.IsNull() ) {
+          // find the corresponding face of theGroup2
+          TopoDS_Shape edge2 = theMap( *edge1 );
+          TopoDS_Shape face2;
+          ancestIt.Initialize( theMesh.GetAncestors( edge2 ) );
+          for ( ; ancestIt.More() && face2.IsNull(); ancestIt.Next() ) {
+            if ( ancestIt.Value().ShapeType() == TopAbs_FACE &&
+                 !theMap.IsBound( ancestIt.Value() ) &&
+                 subshapes2.Contains( ancestIt.Value() ))
+              face2 = ancestIt.Value();
+          }
+          if ( face2.IsNull() )
+            return false;
+
+          HERE::InsertAssociation( face1, face2, theMap);
+        }
+      }
+    }
+    return true;
+  }
+
+} // namespace
 
 //=======================================================================
 /*!
@@ -312,7 +457,7 @@ bool StdMeshers_ProjectionUtils::FindSubShapeAssociation(const TopoDS_Shape& the
     case TopAbs_SOLID: {
       // ----------------------------------------------------------------------
       TopoDS_Vertex VV1[2], VV2[2];
-      // find a not closed edge of shape1 both vertices of which are associated
+      // try to find a not closed edge of shape1 both vertices of which are associated
       TopoDS_Edge edge1;
       TopExp_Explorer exp ( theShape1, TopAbs_EDGE );
       for ( ; VV2[ 1 ].IsNull() && exp.More(); exp.Next() ) {
@@ -326,6 +471,7 @@ bool StdMeshers_ProjectionUtils::FindSubShapeAssociation(const TopoDS_Shape& the
       }
       if ( VV2[ 1 ].IsNull() ) // 2 bound vertices not found
         RETURN_BAD_RESULT("2 bound vertices not found" );
+      // get an edge2 of theShape2 corresponding to edge1
       TopoDS_Edge edge2 = GetEdgeByVertices( theMesh2, VV2[ 0 ], VV2[ 1 ]);
       if ( edge2.IsNull() )
         RETURN_BAD_RESULT("GetEdgeByVertices() failed");
@@ -773,6 +919,77 @@ bool StdMeshers_ProjectionUtils::FindSubShapeAssociation(const TopoDS_Shape& the
     }
     break; // try by vertex closeness
   }
+  case TopAbs_COMPOUND: {
+    // ----------------------------------------------------------------------
+    if ( IsPropagationPossible( theMesh1, theMesh2 )) {
+
+      // try to accosiate all using propagation
+      if ( AssocGroupsByPropagation( theShape1, theShape2, *theMesh1, theMap ))
+        return true;
+
+      // find a boundary edge for theShape1
+      TopoDS_Edge E;
+      for(TopExp_Explorer exp(theShape1, TopAbs_EDGE); exp.More(); exp.Next() ) {
+        E = TopoDS::Edge( exp.Current() );
+        if ( IsBoundaryEdge( E, theShape1, *theMesh1 ))
+          break;
+        else
+          E.Nullify();
+      }
+      if ( E.IsNull() )
+        break; // try by vertex closeness
+
+      // find association for vertices of edge E
+      TopoDS_Vertex VV1[2], VV2[2];
+      for(TopExp_Explorer eexp(E, TopAbs_VERTEX); eexp.More(); eexp.Next()) {
+        TopoDS_Vertex V1 = TopoDS::Vertex( eexp.Current() );
+        // look for an edge ending in E whose one vertex is in theShape1
+        // and the other, in theShape2
+        const TopTools_ListOfShape& Ancestors = theMesh1->GetAncestors(V1);
+        TopTools_ListIteratorOfListOfShape ita(Ancestors);
+        for(; ita.More(); ita.Next()) {
+          if( ita.Value().ShapeType() != TopAbs_EDGE ) continue;
+          TopoDS_Edge edge = TopoDS::Edge(ita.Value());
+          bool FromShape1 = false;
+          for(TopExp_Explorer expe(theShape1, TopAbs_EDGE); expe.More(); expe.Next() ) {
+            if(edge.IsSame(expe.Current())) {
+              FromShape1 = true;
+              break;
+            }
+          }
+          if(!FromShape1) {
+            // is it an edge between theShape1 and theShape2?
+            TopExp_Explorer expv(edge, TopAbs_VERTEX);
+            TopoDS_Vertex V2 = TopoDS::Vertex( expv.Current() );
+            if(V2.IsSame(V1)) {
+              expv.Next();
+              V2 = TopoDS::Vertex( expv.Current() );
+            }
+            bool FromShape2 = false;
+            for ( expv.Init( theShape2, TopAbs_VERTEX ); expv.More(); expv.Next()) {
+              if ( V2.IsSame( expv.Current() )) {
+                FromShape2 = true;
+                break;
+              }
+            }
+            if ( FromShape2 ) {
+              if ( VV1[0].IsNull() )
+                VV1[0] = V1, VV2[0] = V2;
+              else
+                VV1[1] = V1, VV2[1] = V2;
+              break; // from loop on ancestors of V1
+            }
+          }
+        }
+      }
+      if ( !VV1[1].IsNull() ) {
+        InsertAssociation( VV1[0], VV2[0], theMap, bidirect);
+        InsertAssociation( VV1[1], VV2[1], theMap, bidirect);
+        return FindSubShapeAssociation( theShape1, theMesh1, theShape2, theMesh2, theMap);
+      }
+    }
+    break; // try by vertex closeness
+  }
   default:;
   }
 
@@ -881,7 +1098,7 @@ bool StdMeshers_ProjectionUtils::FindSubShapeAssociation(const TopoDS_Shape& the
  * \param face1 - face 1
  * \param VV1 - vertices of face 1
  * \param face2 - face 2
- * \param VV2 - vertices of face 2 associated with oned of face 1
+ * \param VV2 - vertices of face 2 associated with ones of face 1
  * \param edges1 - out list of edges of face 1
  * \param edges2 - out list of edges of face 2
  * \retval int - nb of edges in an outer wire in a success case, else zero
@@ -1285,7 +1502,7 @@ FindMatchingNodesOnFaces( const TopoDS_Face&     face1,
     }
     if ( nodesOfFaces )
     {
-      if ( BRep_Tool::IsClosed( e2, face2 )) {
+      if ( helper2.IsRealSeam( e2 )) {
         seam1 = e1; seam2 = e2;
       }
       else {
@@ -1357,7 +1574,7 @@ FindMatchingNodesOnFaces( const TopoDS_Face&     face1,
     const TopoDS_Face &             face = is2 ? face2 : face1;
     SMDS_ElemIteratorPtr eIt = sm->GetElements();
 
-    if ( !helper->IsSeamShape( is2 ? edge2 : edge1 ))
+    if ( !helper->IsRealSeam( is2 ? edge2 : edge1 ))
     {
       while ( eIt->more() ) elems.insert( eIt->next() );
     }
@@ -1442,7 +1659,7 @@ FindMatchingNodesOnFaces( const TopoDS_Face&     face1,
 
   // On a sphere, add matching nodes on the edge
 
-  if ( helper1.IsSeamShape( edge1 ))
+  if ( helper1.IsRealSeam( edge1 ))
   {
     // sort nodes on edges by param on edge
     map< double, const SMDS_MeshNode* > u2nodesMaps[2];
@@ -1562,7 +1779,17 @@ bool StdMeshers_ProjectionUtils::MakeComputed(SMESH_subMesh * sm, const int iter
   SMESH_Gen* gen   = mesh->GetGen();
   SMESH_Algo* algo = gen->GetAlgo( *mesh, sm->GetSubShape() );
   if ( !algo )
-    RETURN_BAD_RESULT("No algo assigned to submesh " << sm->GetId());
+  {
+    if ( sm->GetSubShape().ShapeType() != TopAbs_COMPOUND )
+      RETURN_BAD_RESULT("No algo assigned to submesh " << sm->GetId());
+    // group
+    bool computed = true;
+    for ( TopoDS_Iterator grMember( sm->GetSubShape() ); grMember.More(); grMember.Next())
+      if ( SMESH_subMesh* grSub = mesh->GetSubMesh( grMember.Value() ))
+        if ( !MakeComputed( grSub, iterationNb + 1 ))
+          computed = false;
+    return computed;
+  }
 
   string algoType = algo->GetName();
   if ( algoType.substr(0, 11) != "Projection_")
@@ -1606,16 +1833,19 @@ bool StdMeshers_ProjectionUtils::MakeComputed(SMESH_subMesh * sm, const int iter
   if ( !srcMesh )
     srcMesh = mesh;
 
-  return MakeComputed( srcMesh->GetSubMesh( srcShape ), iterationNb + 1 );
+  if ( MakeComputed( srcMesh->GetSubMesh( srcShape ), iterationNb + 1 ))
+    return gen->Compute( *mesh, sm->GetSubShape() );
+
+  return false;
 }
 
 //================================================================================
-  /*!
  * \brief Count nb of subshapes
   * \param shape - the shape
   * \param type - the type of subshapes to count
   * \retval int - the calculated number
  */
+/*!
+ * \brief Count nb of subshapes
+ * \param shape - the shape
+ * \param type - the type of subshapes to count
+ * \retval int - the calculated number
+ */
 //================================================================================
 
 int StdMeshers_ProjectionUtils::Count(const TopoDS_Shape&    shape,
@@ -1635,6 +1865,34 @@ int StdMeshers_ProjectionUtils::Count(const TopoDS_Shape&    shape,
   }
 }
 
+//================================================================================
+/*!
+ * \brief Return true if edge is a boundary of edgeContainer
+ */
+//================================================================================
+
+bool StdMeshers_ProjectionUtils::IsBoundaryEdge(const TopoDS_Edge&  edge,
+                                                const TopoDS_Shape& edgeContainer,
+                                                SMESH_Mesh&         mesh)
+{
+  TopTools_IndexedMapOfShape facesOfEdgeContainer, facesNearEdge;
+  TopExp::MapShapes( edgeContainer, TopAbs_FACE, facesOfEdgeContainer );
+
+  const TopTools_ListOfShape& EAncestors = mesh.GetAncestors(edge);
+  TopTools_ListIteratorOfListOfShape itea(EAncestors);
+  for(; itea.More(); itea.Next()) {
+    if( itea.Value().ShapeType() == TopAbs_FACE &&
+        facesOfEdgeContainer.Contains( itea.Value() ))
+    {
+      facesNearEdge.Add( itea.Value() );
+      if ( facesNearEdge.Extent() > 1 )
+        return false;
+    }
+  }
+  return ( facesNearEdge.Extent() == 1 );
+}
+
+
 namespace {
 
   SMESH_subMeshEventListener* GetSrcSubMeshListener();