Salome HOME
22222: [CEA 820] GHS3D in salome 7.2.0 ten times slower than in salome 6.6.0
[modules/smesh.git] / src / StdMeshers / StdMeshers_ProjectionUtils.cxx
index 08136aa73352cfef83e3f6466a0a382761c43574..de59b265423948b4b4eb9698195e9e47fc51640b 100644 (file)
@@ -41,6 +41,7 @@
 #include "SMESH_MesherHelper.hxx"
 #include "SMESH_subMesh.hxx"
 #include "SMESH_subMeshEventListener.hxx"
+#include "SMESH_MeshAlgos.hxx"
 
 #include "utilities.h"
 
@@ -1000,7 +1001,8 @@ bool StdMeshers_ProjectionUtils::FindSubShapeAssociation(const TopoDS_Shape& the
       // get outer edge of theShape1
       TopoDS_Shape wire = OuterShape( face1, TopAbs_WIRE );
       //edge1 = TopoDS::Edge( OuterShape( face1, TopAbs_EDGE ));
-      map<int,TopoDS_Edge> propag_edges; // use map to find the closest propagation edge
+      // use map to find the closest propagation edge
+      map<int, pair< TopoDS_Edge, TopoDS_Edge > > propag_edges;
       for ( TopoDS_Iterator edgeIt( wire ); edgeIt.More(); edgeIt.Next() )
       {
         edge1 = TopoDS::Edge( edgeIt.Value() );
@@ -1009,7 +1011,8 @@ bool StdMeshers_ProjectionUtils::FindSubShapeAssociation(const TopoDS_Shape& the
           edge2 = TopoDS::Edge( exp.Current() );
           pair<int,TopoDS_Edge> step_edge = GetPropagationEdge( theMesh1, edge2, edge1 );
           if ( !step_edge.second.IsNull() ) { // propagation found
-            propag_edges.insert( step_edge );
+            propag_edges.insert( make_pair( step_edge.first,
+                                            ( make_pair( edge1, step_edge.second ))));
             if ( step_edge.first == 1 ) break; // most close found
           }
         }
@@ -1017,7 +1020,8 @@ bool StdMeshers_ProjectionUtils::FindSubShapeAssociation(const TopoDS_Shape& the
       }
       if ( !propag_edges.empty() ) // propagation found
       {
-        edge2 = propag_edges.begin()->second;
+        edge1 = propag_edges.begin()->second.first;
+        edge2 = propag_edges.begin()->second.second;
         TopoDS_Vertex VV1[2], VV2[2];
         TopExp::Vertices( edge1, VV1[0], VV1[1], true );
         TopExp::Vertices( edge2, VV2[0], VV2[1], true );
@@ -1873,12 +1877,12 @@ FindMatchingNodesOnFaces( const TopoDS_Face&     face1,
       TIDSortedElemSet inSet, notInSet;
 
       const SMDS_MeshElement* f1 =
-        SMESH_MeshEditor::FindFaceInSet( vNode, eNode, inSet, notInSet );
+        SMESH_MeshAlgos::FindFaceInSet( vNode, eNode, inSet, notInSet );
       if ( !f1 ) RETURN_BAD_RESULT("The first face on seam not found");
       notInSet.insert( f1 );
 
       const SMDS_MeshElement* f2 =
-        SMESH_MeshEditor::FindFaceInSet( vNode, eNode, inSet, notInSet );
+        SMESH_MeshAlgos::FindFaceInSet( vNode, eNode, inSet, notInSet );
       if ( !f2 ) RETURN_BAD_RESULT("The second face on seam not found");
 
       // select a face with less UV of vNode
@@ -1921,7 +1925,7 @@ FindMatchingNodesOnFaces( const TopoDS_Face&     face1,
       for ( int i = 0; i < nbNodes; ++i ) {
         const SMDS_MeshNode* n1 = faceToKeep->GetNode( i );
         const SMDS_MeshNode* n2 = faceToKeep->GetNode(( i+1 ) % nbNodes );
-        f1 = SMESH_MeshEditor::FindFaceInSet( n1, n2, inSet, notInSet );
+        f1 = SMESH_MeshAlgos::FindFaceInSet( n1, n2, inSet, notInSet );
         if ( f1 )
           elems.insert( f1 );
       }
@@ -2027,8 +2031,8 @@ TopoDS_Shape StdMeshers_ProjectionUtils::OuterShape( const TopoDS_Face& face,
 
 //================================================================================
 /*
- * Check that submesh is computed and try to compute it if is not
- *  \param sm - submesh to compute
+ * Check that sub-mesh is computed and try to compute it if is not
+ *  \param sm - sub-mesh to compute
  *  \param iterationNb - int used to stop infinite recursive call
  *  \retval bool - true if computed
  */
@@ -2043,30 +2047,65 @@ bool StdMeshers_ProjectionUtils::MakeComputed(SMESH_subMesh * sm, const int iter
   if ( sm->IsMeshComputed() )
     return true;
 
-  SMESH_Mesh* mesh = sm->GetFather();
-  SMESH_Gen*   gen = mesh->GetGen();
-  SMESH_Algo* algo = sm->GetAlgo();
+  SMESH_Mesh*   mesh = sm->GetFather();
+  SMESH_Gen*     gen = mesh->GetGen();
+  SMESH_Algo*   algo = sm->GetAlgo();
+  TopoDS_Shape shape = sm->GetSubShape();
   if ( !algo )
   {
-    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;
+    if ( shape.ShapeType() != TopAbs_COMPOUND )
+    {
+      // No algo assigned to a non-compound sub-mesh.
+      // Try to find an all-dimensional algo of an upper dimension
+      int dim = gen->GetShapeDim( shape );
+      for ( ++dim; ( dim <= 3 && !algo ); ++dim )
+      {
+        SMESH_HypoFilter hypoFilter( SMESH_HypoFilter::IsAlgo() );
+        hypoFilter.And( SMESH_HypoFilter::HasDim( dim ));
+        list <const SMESHDS_Hypothesis * > hyps;
+        list< TopoDS_Shape >               assignedTo;
+        int nbAlgos =
+          mesh->GetHypotheses( shape, hypoFilter, hyps, true, &assignedTo );
+        if ( nbAlgos > 1 ) // concurrent algos
+        {
+          list<SMESH_subMesh*> smList; // where an algo is assigned
+          list< TopoDS_Shape >::iterator shapeIt = assignedTo.begin();
+          for ( ; shapeIt != assignedTo.end(); ++shapeIt )
+            smList.push_back( mesh->GetSubMesh( *shapeIt ));
+
+          mesh->SortByMeshOrder( smList );
+          algo  = smList.front()->GetAlgo();
+          shape = smList.front()->GetSubShape();
+        }
+        else if ( nbAlgos == 1 )
+        {
+          algo = (SMESH_Algo*) hyps.front();
+          shape = assignedTo.front();
+        }
+      }
+      if ( !algo )
+        return false;
+    }
+    else
+    {
+      // group
+      bool computed = true;
+      for ( TopoDS_Iterator grMember( shape ); 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_")
-    return gen->Compute( *mesh, sm->GetSubShape() );
+    return gen->Compute( *mesh, shape, /*shapeOnly=*/true );
 
   // try to compute source mesh
 
   const list <const SMESHDS_Hypothesis *> & hyps =
-    algo->GetUsedHypothesis( *mesh, sm->GetSubShape() );
+    algo->GetUsedHypothesis( *mesh, shape );
 
   TopoDS_Shape srcShape;
   SMESH_Mesh* srcMesh = 0;
@@ -2093,16 +2132,16 @@ bool StdMeshers_ProjectionUtils::MakeComputed(SMESH_subMesh * sm, const int iter
     }
   }
   if ( srcShape.IsNull() ) // no projection source defined
-    return gen->Compute( *mesh, sm->GetSubShape() );
+    return gen->Compute( *mesh, shape, /*shapeOnly=*/true );
 
-  if ( srcShape.IsSame( sm->GetSubShape() ))
+  if ( srcShape.IsSame( shape ))
     RETURN_BAD_RESULT("Projection from self");
     
   if ( !srcMesh )
     srcMesh = mesh;
 
   if ( MakeComputed( srcMesh->GetSubMesh( srcShape ), iterationNb + 1 ) &&
-       gen->Compute( *mesh, sm->GetSubShape() ))
+       gen->Compute( *mesh, shape, /*shapeOnly=*/true ))
     return sm->IsMeshComputed();
 
   return false;