]> SALOME platform Git repositories - plugins/blsurfplugin.git/commitdiff
Salome HOME
bos #16292 [CEA 6719] MGCADSurf: option SetEnforced1D mesh
authoreap <eap@opencascade.com>
Thu, 24 Feb 2022 10:05:01 +0000 (13:05 +0300)
committereap <eap@opencascade.com>
Thu, 24 Feb 2022 10:05:01 +0000 (13:05 +0300)
 Add dependency of an enforced mesh

src/BLSURFPlugin/BLSURFPlugin_BLSURF.cxx
src/BLSURFPlugin/BLSURFPlugin_BLSURF.hxx
src/BLSURFPlugin/BLSURFPlugin_EnforcedMesh1D.cxx

index abb3dfb8d4897f689c9e160faee62a5229b8dab4..cd4eef4631f4467d9fea118f63552ea404589e11 100644 (file)
@@ -45,6 +45,7 @@ extern "C"{
 #include <SMESH_MeshEditor.hxx>
 #include <SMESH_MesherHelper.hxx>
 #include <StdMeshers_FaceSide.hxx>
+#include <StdMeshers_ProjectionUtils.hxx>
 #include <StdMeshers_ViscousLayers2D.hxx>
 
 #include <Basics_Utils.hxx>
@@ -91,8 +92,6 @@ extern "C"{
 #include <TopoDS_Wire.hxx>
 #include <gp_Pnt.hxx>
 #include <gp_Pnt2d.hxx>
-#include <gp_XY.hxx>
-#include <gp_XYZ.hxx>
 
 #ifndef WIN32
 #include <fenv.h>
@@ -2969,6 +2968,14 @@ bool BLSURFPlugin_BLSURF::compute(SMESH_Mesh&         aMesh,
     error(_comment);
   }
 
+  // Set event listeners
+  if ( _hypothesis )
+    for ( int iF = 1; iF <= fmap.Size(); ++iF )
+    {
+      const TopoDS_Shape& face = fmap( iF );
+      SetEventListener( aMesh.GetSubMesh( face ));
+    }
+
   // Issue 0019864. On DebianSarge, FE signals do not obey to OSD::SetSignal(false)
 #ifndef WIN32
   if ( oldFEFlags > 0 )
@@ -3678,3 +3685,52 @@ void BLSURFPlugin_BLSURF::FillEntryToShape( const BLSURFPlugin_Hypothesis*
         entryToShape.insert({ entry, shape });
     }
 }
+
+//================================================================================
+/*!
+ * \brief Sets event listener to submeshes if enforced mesh is defined
+ * \param subMesh - submesh where algo is set
+ *
+ * This method is called when a submesh gets HYP_OK algo_state.
+ * After being set, event listener is notified on each event of a submesh.
+ * By default none listener is set
+ */
+//================================================================================
+
+void BLSURFPlugin_BLSURF::SetEventListener(SMESH_subMesh* faceSubMesh)
+{
+  if ( !_hypothesis )
+    return;
+
+  for ( const BLSURFPlugin_Hypothesis::EnforcedMesh& enfMesh : _hypothesis->GetEnforcedMeshes() )
+  {
+    SMESH_Mesh* mesh1D;
+    _hypothesis->GetEnforcedSegments( enfMesh, mesh1D );
+    if ( !mesh1D )
+      continue;
+
+    TopExp_Explorer edgeExp( mesh1D->GetShapeToMesh(), TopAbs_EDGE );
+    if ( edgeExp.More() )
+      StdMeshers_ProjectionUtils::SetEventListener( faceSubMesh,
+                                                    edgeExp.Current(),
+                                                    mesh1D );
+    // StdMeshers_ProjectionUtils::SetEventListener( faceSubMesh,
+    //                                               mesh1D->GetShapeToMesh(),
+    //                                               mesh1D );
+  }
+  return;
+}
+
+//================================================================================
+/*!
+ * \brief Allow algo to do something after persistent restoration
+ * \param subMesh - restored submesh
+ *
+ * This method is called only if a submesh has HYP_OK algo_state.
+ */
+//================================================================================
+
+void BLSURFPlugin_BLSURF::SubmeshRestored(SMESH_subMesh* subMesh)
+{
+  SetEventListener( subMesh );
+}
index 936bdd3dab7fa1e739fc63abd2ed86a8e802695a..4c045a4acd58f06cb673188aebfaa43b1affeeef 100644 (file)
@@ -98,6 +98,11 @@ public:
   static void FillEntryToShape( const BLSURFPlugin_Hypothesis*          hyp,
                                 std::map< std::string, TopoDS_Shape > & s2eMap );
 
+  virtual void SetEventListener(SMESH_subMesh* subMesh) override;
+  virtual void SubmeshRestored(SMESH_subMesh* subMesh) override;
+
+
+
   // List of ids
   typedef std::vector<int> TListOfIDs;
 
index 45a302b4a1161eec26d1f10cc138777727eccdf6..1268ccb2568e17550c5df7f0e711b3eb4b0b3c29 100644 (file)
@@ -101,8 +101,6 @@ namespace
    *  \param [in] braNodes - nodes of the branch
    *  \param [in] nodeIndex - index of a node of the branch
    *  \param [inout] mesh - mesh holding the nodes and segments
-   * 
-   * 
    */
   //================================================================================
 
@@ -336,7 +334,7 @@ copyEnforcedMesh( const BLSURFPlugin_Hypothesis::EnforcedMesh& theEnfMesh,
   SMESH_Mesh* mesh1D;
   SMDS_ElemIteratorPtr segIt = theHyp->GetEnforcedSegments( theEnfMesh, mesh1D );
   if ( !segIt->more() )
-    return;
+    throw SALOME_Exception("No edges in an enforced mesh");
 
   // setup predicates to detect nodes on FACE boundary
   setupPredicates( theShape );