]> SALOME platform Git repositories - modules/smesh.git/commitdiff
Salome HOME
Regression of SALOME_TESTS/Grids/smesh/imps_09/K2
authoreap <eap@opencascade.com>
Thu, 24 Apr 2014 16:06:54 +0000 (20:06 +0400)
committereap <eap@opencascade.com>
Thu, 24 Apr 2014 16:06:54 +0000 (20:06 +0400)
src/StdMeshers/StdMeshers_Prism_3D.cxx
src/StdMeshers/StdMeshers_Prism_3D.hxx
src/StdMeshers/StdMeshers_ProjectionUtils.cxx
src/StdMeshers/StdMeshers_ProjectionUtils.hxx

index c43a126c3e2823188a2992bec17ac34dbfbf1e8a..579cbe55011420ac9c6d1725381b2e59bdcaf906 100644 (file)
@@ -157,6 +157,46 @@ namespace {
       return algo;
     }
   };
+  //=======================================================================
+  /*!
+   * \brief Returns already computed EDGEs
+   */
+  void getPrecomputedEdges( SMESH_MesherHelper&    theHelper,
+                            const TopoDS_Shape&    theShape,
+                            vector< TopoDS_Edge >& theEdges)
+  {
+    theEdges.clear();
+
+    SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
+    SMESHDS_SubMesh* sm;
+
+    TopTools_IndexedMapOfShape edges;
+    TopExp::MapShapes( theShape, TopAbs_EDGE, edges );
+    for ( int iE = 1; iE <= edges.Extent(); ++iE )
+    {
+      const TopoDS_Shape edge = edges( iE );
+      if (( ! ( sm = meshDS->MeshElements( edge ))) ||
+          ( sm->NbElements() == 0 ))
+        continue;
+
+      // there must not be FACEs meshed with triangles and sharing a computed EDGE
+      // as the precomputed EDGEs are used for propagation other to 'vertical' EDGEs
+      bool faceFound = false;
+      PShapeIteratorPtr faceIt =
+        theHelper.GetAncestors( edge, *theHelper.GetMesh(), TopAbs_FACE );
+      while ( const TopoDS_Shape* face = faceIt->next() )
+
+        if (( sm = meshDS->MeshElements( *face )) &&
+            ( sm->NbElements() > 0 ) &&
+            ( !theHelper.IsSameElemGeometry( sm, SMDSGeom_QUADRANGLE ) ))
+        {
+          faceFound;
+          break;
+        }
+      if ( !faceFound )
+        theEdges.push_back( TopoDS::Edge( edge ));
+    }
+  }
 
   //================================================================================
   /*!
@@ -623,6 +663,21 @@ bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theSh
              compute( prism ));
   }
 
+  // find propagation chains from already computed EDGEs
+  vector< TopoDS_Edge > computedEdges;
+  getPrecomputedEdges( helper, theShape, computedEdges );
+  myPropagChains = new TopTools_IndexedMapOfShape[ computedEdges.size() + 1 ];
+  SMESHUtils::ArrayDeleter< TopTools_IndexedMapOfShape > pcDel( myPropagChains );
+  for ( size_t i = 0, nb = 0; i < computedEdges.size(); ++i )
+  {
+    StdMeshers_ProjectionUtils::GetPropagationEdge( &theMesh, TopoDS_Edge(),
+                                                    computedEdges[i], myPropagChains + nb );
+    if ( myPropagChains[ nb ].Extent() < 2 ) // an empty map is a termination sign
+      myPropagChains[ nb ].Clear();
+    else
+      nb++;
+  }
+
   TopTools_MapOfShape meshedSolids;
   list< Prism_3D::TPrismTopo > meshedPrism;
   TopTools_ListIteratorOfListOfShape solidIt;
@@ -650,7 +705,11 @@ bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theSh
                !compute( prism ))
             return false;
 
-          meshedFaces.push_front( prism.myTop );
+          SMESHDS_SubMesh* smDS = theMesh.GetMeshDS()->MeshElements( prism.myTop );
+          if ( !myHelper->IsSameElemGeometry( smDS, SMDSGeom_QUADRANGLE ))
+          {
+            meshedFaces.push_front( prism.myTop );
+          }
           meshedPrism.push_back( prism );
         }
       }
@@ -691,6 +750,7 @@ bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theSh
               prism.myBottom  = candidateF;
               mySetErrorToSM = false;
               if ( !myHelper->IsSubShape( candidateF, prismIt->myShape3D ) &&
+                   myHelper->IsSubShape( candidateF, solid ) &&
                    !myHelper->GetMesh()->GetSubMesh( candidateF )->IsMeshComputed() &&
                    initPrism( prism, solid ) &&
                    project2dMesh( prismIt->myBottom, candidateF))
@@ -698,8 +758,12 @@ bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theSh
                 mySetErrorToSM = true;
                 if ( !compute( prism ))
                   return false;
-                meshedFaces.push_front( prism.myTop );
-                meshedFaces.push_front( prism.myBottom );
+                SMESHDS_SubMesh* smDS = theMesh.GetMeshDS()->MeshElements( prism.myTop );
+                if ( !myHelper->IsSameElemGeometry( smDS, SMDSGeom_QUADRANGLE ))
+                {
+                  meshedFaces.push_front( prism.myTop );
+                  meshedFaces.push_front( prism.myBottom );
+                }
                 meshedPrism.push_back( prism );
                 meshedSolids.Add( solid );
               }
@@ -1187,7 +1251,7 @@ 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 );
+  TProjction1dAlgo*      projector1D = TProjction1dAlgo::instance( this );
   StdMeshers_Quadrangle_2D* quadAlgo = TQuadrangleAlgo::instance( this, myHelper );
 
   SMESH_HypoFilter hyp1dFilter( SMESH_HypoFilter::IsAlgo(),/*not=*/true);
@@ -1252,8 +1316,17 @@ bool StdMeshers_Prism_3D::computeWalls(const Prism_3D::TPrismTopo& thePrism)
         SMESH_subMesh*    srcSM = mesh->GetSubMesh( srcE );
         if ( !srcSM->IsMeshComputed() ) {
           DBGOUT( "COMPUTE V edge " << srcSM->GetId() );
-          srcSM->ComputeSubMeshStateEngine( SMESH_subMesh::COMPUTE );
-          srcSM->ComputeStateEngine       ( SMESH_subMesh::COMPUTE );
+          TopoDS_Edge prpgSrcE = findPropagationSource( srcE );
+          if ( !prpgSrcE.IsNull() ) {
+            srcSM->ComputeSubMeshStateEngine( SMESH_subMesh::COMPUTE );
+            projector1D->myHyp.SetSourceEdge( prpgSrcE );
+            projector1D->Compute( *mesh, srcE );
+            srcSM->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
+          }
+          else {
+            srcSM->ComputeSubMeshStateEngine( SMESH_subMesh::COMPUTE );
+            srcSM->ComputeStateEngine       ( SMESH_subMesh::COMPUTE );
+          }
           if ( !srcSM->IsMeshComputed() )
             return toSM( error( "Can't compute 1D mesh" ));
         }
@@ -1464,6 +1537,21 @@ bool StdMeshers_Prism_3D::computeWalls(const Prism_3D::TPrismTopo& thePrism)
   return true;
 }
 
+//=======================================================================
+/*!
+ * \brief Returns a source EDGE of propagation to a given EDGE
+ */
+//=======================================================================
+
+TopoDS_Edge StdMeshers_Prism_3D::findPropagationSource( const TopoDS_Edge& E )
+{
+  for ( size_t i = 0; !myPropagChains[i].IsEmpty(); ++i )
+    if ( myPropagChains[i].Contains( E ))
+      return TopoDS::Edge( myPropagChains[i].FindKey( 1 ));
+
+  return TopoDS_Edge();
+}
+
 //=======================================================================
 //function : Evaluate
 //purpose  : 
index 3a4f61f4dcd49d4bd2e4294be4db7d8f43a841ba..f4e49ce63190253cd9ec38a4abb5defb46ae257d 100644 (file)
@@ -460,6 +460,11 @@ public:
    */
   bool computeWalls(const Prism_3D::TPrismTopo& thePrism);
 
+  /*!
+   * \brief Returns a source EDGE of propagation to a given EDGE
+   */
+  TopoDS_Edge findPropagationSource( const TopoDS_Edge& E );
+
   /*!
    * \brief Find correspondence between bottom and top nodes.
    *  If elements on the bottom and top faces are topologically different,
@@ -514,6 +519,8 @@ private:
   // (the column includes the bottom node)
   TNode2ColumnMap         myBotToColumnMap;
 
+  TopTools_IndexedMapOfShape* myPropagChains;
+
 }; // class StdMeshers_Prism_3D
 
 #endif
index 77d9c1a71b96f07bb0ed84d8b66f33b8324fb3e9..ae6f38419e465580eafdaad0461a7a1eaa107224 100644 (file)
@@ -1613,79 +1613,87 @@ TopoDS_Vertex StdMeshers_ProjectionUtils::GetNextVertex(const TopoDS_Edge&   edg
 /*
  * Return a propagation edge
  *  \param aMesh - mesh
- *  \param theEdge - edge to find by propagation
+ *  \param anEdge - edge to find by propagation
  *  \param fromEdge - start edge for propagation
+ *  \param chain - return, if !NULL, a propagation chain passed till
+ *         anEdge; if anEdge.IsNull() then a full propagation chain is returned;
+ *         fromEdge is the 1st in the chain
  *  \retval pair<int,TopoDS_Edge> - propagation step and found edge
  */
 //================================================================================
 
 pair<int,TopoDS_Edge>
-StdMeshers_ProjectionUtils::GetPropagationEdge( SMESH_Mesh*        aMesh,
-                                                const TopoDS_Edge& theEdge,
-                                                const TopoDS_Edge& fromEdge)
+StdMeshers_ProjectionUtils::GetPropagationEdge( SMESH_Mesh*                 aMesh,
+                                                const TopoDS_Edge&          anEdge,
+                                                const TopoDS_Edge&          fromEdge,
+                                                TopTools_IndexedMapOfShape* chain)
 {
-  TopTools_IndexedMapOfShape aChain;
+  TopTools_IndexedMapOfShape locChain;
+  TopTools_IndexedMapOfShape& aChain = chain ? *chain : locChain;
   int step = 0;
 
+  //TopTools_IndexedMapOfShape checkedWires;
+  BRepTools_WireExplorer aWE;
+  TopoDS_Shape fourEdges[4];
+
   // List of edges, added to chain on the previous cycle pass
   TopTools_ListOfShape listPrevEdges;
-  listPrevEdges.Append(fromEdge);
+  listPrevEdges.Append( fromEdge );
+  aChain.Add( fromEdge );
 
   // Collect all edges pass by pass
-  while (listPrevEdges.Extent() > 0) {
+  while (listPrevEdges.Extent() > 0)
+  {
     step++;
     // List of edges, added to chain on this cycle pass
     TopTools_ListOfShape listCurEdges;
 
     // Find the next portion of edges
     TopTools_ListIteratorOfListOfShape itE (listPrevEdges);
-    for (; itE.More(); itE.Next()) {
-      TopoDS_Shape anE = itE.Value();
+    for (; itE.More(); itE.Next())
+    {
+      const TopoDS_Shape& anE = itE.Value();
 
       // Iterate on faces, having edge <anE>
       TopTools_ListIteratorOfListOfShape itA (aMesh->GetAncestors(anE));
-      for (; itA.More(); itA.Next()) {
-        TopoDS_Shape aW = itA.Value();
+      for (; itA.More(); itA.Next())
+      {
+        const TopoDS_Shape& aW = itA.Value();
 
         // There are objects of different type among the ancestors of edge
-        if (aW.ShapeType() == TopAbs_WIRE) {
-          TopoDS_Shape anOppE;
-
-          BRepTools_WireExplorer aWE (TopoDS::Wire(aW));
-          Standard_Integer nb = 1, found = 0;
-          TopTools_Array1OfShape anEdges (1,4);
-          for (; aWE.More(); aWE.Next(), nb++) {
-            if (nb > 4) {
-              found = 0;
+        if ( aW.ShapeType() == TopAbs_WIRE /*&& checkedWires.Add( aW )*/)
+        {
+          Standard_Integer nb = 0, found = -1;
+          for ( aWE.Init( TopoDS::Wire( aW )); aWE.More(); aWE.Next() ) {
+            if (nb+1 > 4) {
+              found = -1;
               break;
             }
-            anEdges(nb) = aWE.Current();
-            if (anEdges(nb).IsSame(anE)) found = nb;
+            fourEdges[ nb ] = aWE.Current();
+            if ( aWE.Current().IsSame( anE )) found = nb;
+            nb++;
           }
-
-          if (nb == 5 && found > 0) {
+          if (nb == 4 && found >= 0) {
             // Quadrangle face found, get an opposite edge
-            Standard_Integer opp = found + 2;
-            if (opp > 4) opp -= 4;
-            anOppE = anEdges(opp);
+            TopoDS_Shape& anOppE = fourEdges[( found + 2 ) % 4 ];
 
             // add anOppE to aChain if ...
-            if (!aChain.Contains(anOppE)) { // ... anOppE is not in aChain
+            int prevChainSize = aChain.Extent();
+            if ( aChain.Add(anOppE) > prevChainSize ) { // ... anOppE is not in aChain
               // Add found edge to the chain oriented so that to
               // have it co-directed with a forward MainEdge
               TopAbs_Orientation ori = anE.Orientation();
-              if ( anEdges(opp).Orientation() == anEdges(found).Orientation() )
+              if ( anOppE.Orientation() == fourEdges[found].Orientation() )
                 ori = TopAbs::Reverse( ori );
               anOppE.Orientation( ori );
-              if ( anOppE.IsSame( theEdge ))
+              if ( anOppE.IsSame( anEdge ))
                 return make_pair( step, TopoDS::Edge( anOppE ));
-              aChain.Add(anOppE);
               listCurEdges.Append(anOppE);
             }
-          } // if (nb == 5 && found > 0)
+          } // if (nb == 4 && found >= 0)
         } // if (aF.ShapeType() == TopAbs_WIRE)
-      } // for (; itF.More(); itF.Next())
-    } // for (; itE.More(); itE.Next())
+      } // loop on ancestors of anE
+    } // loop on listPrevEdges
 
     listPrevEdges = listCurEdges;
   } // while (listPrevEdges.Extent() > 0)
index 9709cd6c86c53c3214a11345f6f6c4630228cdb2..31afb072f55f7c5ca029bc2f599f83bc329a12b5 100644 (file)
@@ -44,6 +44,7 @@ class SMESH_Hypothesis;
 class SMESH_Mesh;
 class SMESH_subMesh;
 class TopTools_IndexedDataMapOfShapeListOfShape;
+class TopTools_IndexedMapOfShape;
 class TopoDS_Shape;
 
 /*!
@@ -156,11 +157,14 @@ namespace StdMeshers_ProjectionUtils
    * \brief Return an oriented propagation edge
    * \param aMesh - mesh
    * \param fromEdge - start edge for propagation
+   * \param chain - return, if provided, a propagation chain passed till
+   *        anEdge; if anEdge.IsNull() then a full propagation chain is returned
    * \retval pair<int,TopoDS_Edge> - propagation step and found edge
    */
-  std::pair<int,TopoDS_Edge> GetPropagationEdge( SMESH_Mesh*        aMesh,
-                                                 const TopoDS_Edge& anEdge,
-                                                 const TopoDS_Edge& fromEdge);
+  std::pair<int,TopoDS_Edge> GetPropagationEdge( SMESH_Mesh*                 aMesh,
+                                                 const TopoDS_Edge&          anEdge,
+                                                 const TopoDS_Edge&          fromEdge,
+                                                 TopTools_IndexedMapOfShape* chain=0);
 
   /*!
    * \brief Find corresponding nodes on two faces