Salome HOME
Fix test MESH_BLSURF_00_A3. Don't add an enforced vertex if it is already on the...
[plugins/blsurfplugin.git] / src / BLSURFPlugin / BLSURFPlugin_BLSURF.cxx
index 8f517c63eb03a501ad1dc99e94610ae60ee7f1fd..94e8065769ef2d2afaf6da97c51990370d6a739c 100644 (file)
@@ -1,4 +1,4 @@
-// Copyright (C) 2007-2021  CEA/DEN, EDF R&D
+// Copyright (C) 2007-2023  CEA, EDF
 //
 // This library is free software; you can redistribute it and/or
 // modify it under the terms of the GNU Lesser General Public
@@ -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>
@@ -516,12 +515,45 @@ TopoDS_Shape BLSURFPlugin_BLSURF::entryToShape(std::string entry)
   return S;
 }
 
+// Check if a point is already on the shape.
+// if faceShape is defined, we use it
+// else we use the sub-shape of the mesh
+bool _doVertexAlreadyExists(const TopoDS_Face&                   faceShape,
+                            const gp_Pnt&                        aPnt)
+{
+  // If the point is already on the shape, no need to add it
+  TopExp_Explorer ex;
+  TopoDS_Shape shapeToCheckVertices(faceShape);
+  bool alreadyExists(false);
+  if (faceShape.IsNull())
+    shapeToCheckVertices = theHelper->GetSubShape();
+  for (ex.Init(shapeToCheckVertices, TopAbs_VERTEX); ex.More(); ex.Next())
+  {
+    TopoDS_Vertex vertex = TopoDS::Vertex(ex.Current());
+    double tol = BRep_Tool::Tolerance( vertex );
+    gp_Pnt p = BRep_Tool::Pnt(vertex);
+    if ( p.IsEqual( aPnt, tol))
+    {
+      MESSAGE("The enforced vertex is already on the shape => No need to add it.");
+      alreadyExists = true;
+      break;
+    }
+  }
+  return alreadyExists;
+}
+
 void _createEnforcedVertexOnFace(TopoDS_Face                          faceShape,
                                  const gp_Pnt&                        aPnt,
-                                 BLSURFPlugin_Hypothesis::TEnfVertex *enfVertex)
+                                 BLSURFPlugin_Hypothesis::TEnfVertex *enfVertex,
+                                 bool checkVertexAlreadyExists=true)
 {
   BLSURFPlugin_Hypothesis::TEnfVertexCoords enf_coords, coords, s_coords;
 
+  // checkVertexAlreadyExists only for non internal vertices.
+  // For them we set checkVertexAlreadyExists=false as we do have to add them
+  if (checkVertexAlreadyExists && _doVertexAlreadyExists(faceShape, aPnt))
+    return;
+
   // Find the face and get the (u,v) values of the enforced vertex on the face
   BLSURFPlugin_BLSURF::projectionPoint projPnt =
     BLSURFPlugin_BLSURF::getProjectionPoint( faceShape, aPnt, /*allowStateON=*/true );
@@ -1313,7 +1345,7 @@ void BLSURFPlugin_BLSURF::SetParameters(const BLSURFPlugin_Hypothesis* hyp,
             ++theNbAttractors;
         }
         else{
-          MESSAGE("Wrong shape type !!")
+          MESSAGE("Wrong shape type !!");
         }
 
       }
@@ -1373,7 +1405,8 @@ void BLSURFPlugin_BLSURF::SetParameters(const BLSURFPlugin_Hypothesis* hyp,
           enfVertex->geomEntry = "";
           enfVertex->grpName = grpName;
           enfVertex->vertex = TopoDS::Vertex( exp_face.Current() );
-          _createEnforcedVertexOnFace( TopoDS::Face(exp.Current()),  aPnt, enfVertex);
+          _createEnforcedVertexOnFace( TopoDS::Face(exp.Current()),  aPnt, enfVertex,
+                                       /*checkVertexAlreadyExists*/false);
           HasSizeMapOnFace = true;
         }
       }
@@ -1763,7 +1796,7 @@ namespace
       TSeg2EdgeMap seg2EdgeMap;
 
       TopoDS_Iterator edgeIt( wire );
-      for ( int iSeg = 1; edgeIt.More(); edgeIt.Next(), ++iSeg )
+      for ( size_t iSeg = 1; edgeIt.More() && iSeg < nodesOfVertices.size(); edgeIt.Next(), ++iSeg )
       {
         SMESH_TLink link( nodesOfVertices[ iSeg-1 ], nodesOfVertices[ iSeg ]);
         TopoDS_Edge edge( TopoDS::Edge( edgeIt.Value() ));
@@ -2969,6 +3002,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 )
@@ -3073,7 +3114,8 @@ bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelpe
   }
 
   std::string errorTxt;
-  if ( !SMESHUtils_MGLicenseKeyGen::SignMesh( msh, errorTxt ))
+
+  if ( !SMESHUtils_MGLicenseKeyGen::SignMesh( msh, "cadsurf", errorTxt ))
     return error( "Problem with library SalomeMeshGemsKeyGenerator: " + errorTxt );
 
   ret = cadsurf_set_mesh(css, msh);
@@ -3496,7 +3538,7 @@ status_t interrupt_cb(integer *interrupt_status, void *user_data)
   else /* you want to stop MG-CADSurf */
   {
     *interrupt_status = INTERRUPT_STOP;
-    return STATUS_ERROR;
+    return STATUS_OK;
   }
 }
 
@@ -3678,3 +3720,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 );
+}