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 14b9f4b0ccf3204a4d3fe1b2b53a61cd63c08c55..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
@@ -40,10 +40,12 @@ extern "C"{
 #include <SMESH_File.hxx>
 #include <SMESH_Gen.hxx>
 #include <SMESH_Group.hxx>
+#include <SMESH_MGLicenseKeyGen.hxx>
 #include <SMESH_Mesh.hxx>
 #include <SMESH_MeshEditor.hxx>
 #include <SMESH_MesherHelper.hxx>
 #include <StdMeshers_FaceSide.hxx>
+#include <StdMeshers_ProjectionUtils.hxx>
 #include <StdMeshers_ViscousLayers2D.hxx>
 
 #include <Basics_Utils.hxx>
@@ -90,13 +92,14 @@ 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>
 #endif
 
+#include <memory>
+#include <functional>
+
 using namespace std;
 
 /* ==================================
@@ -438,9 +441,10 @@ BLSURFPlugin_BLSURF::getProjectionPoint(TopoDS_Face&  theFace,
       {
         // check location on the face
         BRepClass_FaceClassifier FC( face, uv, BRep_Tool::Tolerance( face ));
-        if ( FC.State() == TopAbs_IN )
+        if (( FC.State() == TopAbs_IN ) ||
+            ( FC.State() == TopAbs_ON && theAllowStateON ))
         {
-          if ( !foundFace.IsNull() )
+          if ( !foundFace.IsNull() && !theAllowStateON )
             return myPoint; // thePoint seems to be TopAbs_ON
           foundFace     = face;
           myPoint.uv    = uv.XY();
@@ -511,13 +515,48 @@ TopoDS_Shape BLSURFPlugin_BLSURF::entryToShape(std::string entry)
   return S;
 }
 
-void _createEnforcedVertexOnFace(TopoDS_Face faceShape, gp_Pnt aPnt, BLSURFPlugin_Hypothesis::TEnfVertex *enfVertex)
+// 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,
+                                 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 myPoint =
-    BLSURFPlugin_BLSURF::getProjectionPoint(faceShape,aPnt);
+  BLSURFPlugin_BLSURF::projectionPoint projPnt =
+    BLSURFPlugin_BLSURF::getProjectionPoint( faceShape, aPnt, /*allowStateON=*/true );
   if ( faceShape.IsNull() )
     return;
 
@@ -525,15 +564,16 @@ void _createEnforcedVertexOnFace(TopoDS_Face faceShape, gp_Pnt aPnt, BLSURFPlugi
   enf_coords.push_back(aPnt.Y());
   enf_coords.push_back(aPnt.Z());
 
-  coords.push_back(myPoint.uv.X());
-  coords.push_back(myPoint.uv.Y());
-  coords.push_back(myPoint.xyz.X());
-  coords.push_back(myPoint.xyz.Y());
-  coords.push_back(myPoint.xyz.Z());
+  coords.push_back(projPnt.uv.X());
+  coords.push_back(projPnt.uv.Y());
+  coords.push_back(projPnt.xyz.X());
+  coords.push_back(projPnt.xyz.Y());
+  coords.push_back(projPnt.xyz.Z());
+  coords.push_back(projPnt.state == TopAbs_ON);
 
-  s_coords.push_back(myPoint.xyz.X());
-  s_coords.push_back(myPoint.xyz.Y());
-  s_coords.push_back(myPoint.xyz.Z());
+  s_coords.push_back(projPnt.xyz.X());
+  s_coords.push_back(projPnt.xyz.Y());
+  s_coords.push_back(projPnt.xyz.Z());
 
   // Save pair projected vertex / enf vertex
   EnfVertexCoords2ProjVertex[s_coords] = enf_coords;
@@ -545,16 +585,10 @@ void _createEnforcedVertexOnFace(TopoDS_Face faceShape, gp_Pnt aPnt, BLSURFPlugi
     (*it)->grpName = enfVertex->grpName;
   }
 
-  int key = 0;
-  if (! FacesWithEnforcedVertices.Contains(faceShape)) {
-    key = FacesWithEnforcedVertices.Add(faceShape);
-  }
-  else {
-    key = FacesWithEnforcedVertices.FindIndex(faceShape);
-  }
+  int key = FacesWithEnforcedVertices.Add( faceShape );
 
   // If a node is already created by an attractor, do not create enforced vertex
-  int attractorKey = FacesWithSizeMap.FindIndex(faceShape);
+  int attractorKey = FacesWithSizeMap.FindIndex( faceShape );
   bool sameAttractor = false;
   if (attractorKey >= 0)
     if (FaceId2AttractorCoords.count(attractorKey) > 0)
@@ -563,7 +597,7 @@ void _createEnforcedVertexOnFace(TopoDS_Face faceShape, gp_Pnt aPnt, BLSURFPlugi
 
   if (FaceId2EnforcedVertexCoords.find(key) != FaceId2EnforcedVertexCoords.end()) {
     if (! sameAttractor)
-      FaceId2EnforcedVertexCoords[key].insert(coords); // there should be no redondant coords here (see std::set management)
+      FaceId2EnforcedVertexCoords[key].insert(coords); // there should be no redundant coords here (see std::set management)
   }
   else {
     if (! sameAttractor) {
@@ -577,13 +611,9 @@ void _createEnforcedVertexOnFace(TopoDS_Face faceShape, gp_Pnt aPnt, BLSURFPlugi
 /////////////////////////////////////////////////////////
 void BLSURFPlugin_BLSURF::createEnforcedVertexOnFace(TopoDS_Shape faceShape, BLSURFPlugin_Hypothesis::TEnfVertexList enfVertexList)
 {
-  BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex;
   gp_Pnt aPnt;
-
-  BLSURFPlugin_Hypothesis::TEnfVertexList::const_iterator enfVertexListIt = enfVertexList.begin();
-
-  for( ; enfVertexListIt != enfVertexList.end() ; ++enfVertexListIt ) {
-    enfVertex = *enfVertexListIt;
+  for( BLSURFPlugin_Hypothesis::TEnfVertex* enfVertex : enfVertexList )
+  {
     // Case of manual coords
     if (enfVertex->coords.size() != 0) {
       aPnt.SetCoord(enfVertex->coords[0],enfVertex->coords[1],enfVertex->coords[2]);
@@ -592,7 +622,7 @@ void BLSURFPlugin_BLSURF::createEnforcedVertexOnFace(TopoDS_Shape faceShape, BLS
 
     // Case of geom vertex coords
     if (enfVertex->geomEntry != "") {
-      TopoDS_Shape GeomShape = entryToShape(enfVertex->geomEntry);
+      TopoDS_Shape     GeomShape = entryToShape(enfVertex->geomEntry);
       TopAbs_ShapeEnum GeomType  = GeomShape.ShapeType();
        if (GeomType == TopAbs_VERTEX)
        {
@@ -675,10 +705,10 @@ void createAttractorOnFace(TopoDS_Shape GeomShape, std::string AttractorFunction
   }
 
   // Get the (u,v) values of the attractor on the face
-  BLSURFPlugin_BLSURF::projectionPoint myPoint =
+  BLSURFPlugin_BLSURF::projectionPoint projPnt =
     BLSURFPlugin_BLSURF::getProjectionPoint(TopoDS::Face(GeomShape),gp_Pnt(xa,ya,za));
-  gp_XY    uvPoint = myPoint.uv;
-  gp_XYZ  xyzPoint = myPoint.xyz;
+  gp_XY    uvPoint = projPnt.uv;
+  gp_XYZ  xyzPoint = projPnt.xyz;
   Standard_Real u0 = uvPoint.X();
   Standard_Real v0 = uvPoint.Y();
   Standard_Real x0 = xyzPoint.X();
@@ -1315,7 +1345,7 @@ void BLSURFPlugin_BLSURF::SetParameters(const BLSURFPlugin_Hypothesis* hyp,
             ++theNbAttractors;
         }
         else{
-          MESSAGE("Wrong shape type !!")
+          MESSAGE("Wrong shape type !!");
         }
 
       }
@@ -1375,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;
         }
       }
@@ -1765,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() ));
@@ -2141,14 +2172,19 @@ bool BLSURFPlugin_BLSURF::compute(SMESH_Mesh&         aMesh,
       BLSURFPlugin_Hypothesis::TEnfVertexCoordsList::const_iterator evlIt = evl.begin();
       for (; evlIt != evl.end(); ++evlIt)
       {
+        BLSURFPlugin_Hypothesis::TEnfVertexCoords xyzCoords = { evlIt->at(2),
+                                                                evlIt->at(3),
+                                                                evlIt->at(4)};
+        bool isOnEdge = evlIt->at(5);
+        if ( isOnEdge )
+        {
+          enforcedMesh.AddVertexOnEdge( xyzCoords.data() );
+          continue;
+        }
         double uvCoords[2] = { evlIt->at(0), evlIt->at(1) };
         ienf++;
         cad_point_t* point_p = cad_point_new(fce, ienf, uvCoords);
         int tag = 0;
-        BLSURFPlugin_Hypothesis::TEnfVertexCoords xyzCoords;
-        xyzCoords.push_back(evlIt->at(2));
-        xyzCoords.push_back(evlIt->at(3));
-        xyzCoords.push_back(evlIt->at(4));
         std::map< BLSURFPlugin_Hypothesis::TEnfVertexCoords, BLSURFPlugin_Hypothesis::TEnfVertexList >::const_iterator enfCoordsIt = EnfVertexCoords2EnfVertexList.find(xyzCoords);
         if (enfCoordsIt != EnfVertexCoords2EnfVertexList.end() &&
             !enfCoordsIt->second.empty() )
@@ -2529,6 +2565,9 @@ bool BLSURFPlugin_BLSURF::compute(SMESH_Mesh&         aMesh,
     // set_param(css, "global_physical_size", val_to_string( minFaceSize * 0.5 ).c_str());
     // set_param(css, "max_size",             val_to_string( minFaceSize * 5 ).c_str());
   }
+  std::string errorTxt;
+  if ( !SMESHUtils_MGLicenseKeyGen::SignCAD( c, errorTxt ))
+    return error( "Problem with library SalomeMeshGemsKeyGenerator: " + errorTxt );
 
   // Use the original dcad
   cadsurf_set_dcad(css, dcad);
@@ -2626,18 +2665,22 @@ bool BLSURFPlugin_BLSURF::compute(SMESH_Mesh&         aMesh,
 
   /* retrieve mesh data (see meshgems/mesh.h) */
   integer nv, ne, nt, nq, vtx[4], tag, nb_tag;
-  integer *evedg, *evtri, *evquad, *tags_buff, type;
+  integer type;
   real xyz[3];
 
   mesh_get_vertex_count    (msh, &nv);
   mesh_get_edge_count      (msh, &ne);
   mesh_get_triangle_count  (msh, &nt);
   mesh_get_quadrangle_count(msh, &nq);
+  
+  using deleted_unique_ptr = std::unique_ptr<integer,std::function<void(integer*)>>;
+  
+  deleted_unique_ptr evedg_var((integer *)mesh_calloc_generic_buffer(msh), [](integer *ptr) { free(ptr); });
+  deleted_unique_ptr evtri_var((integer *)mesh_calloc_generic_buffer(msh), [](integer *ptr) { free(ptr); });
+  deleted_unique_ptr evquad_var((integer *)mesh_calloc_generic_buffer(msh), [](integer *ptr) { free(ptr); });
+  deleted_unique_ptr tags_buff_var((integer*)mesh_calloc_generic_buffer(msh), [](integer *ptr) { free(ptr); });
 
-  evedg  = (integer *)mesh_calloc_generic_buffer(msh);
-  evtri  = (integer *)mesh_calloc_generic_buffer(msh);
-  evquad = (integer *)mesh_calloc_generic_buffer(msh);
-  tags_buff = (integer*)mesh_calloc_generic_buffer(msh);
+  integer *evedg(evedg_var.get()),*evtri(evtri_var.get()),*evquad(evquad_var.get()),*tags_buff(tags_buff_var.get());
 
   std::vector<const SMDS_MeshNode*> nodes(nv+1);
   std::vector<bool>                  tags(nv+1);
@@ -2651,6 +2694,7 @@ bool BLSURFPlugin_BLSURF::compute(SMESH_Mesh&         aMesh,
     mesh_get_vertex_tag(msh, iv, &tag);
     // Issue 0020656. Use vertex coordinates
     nodes[iv] = NULL;
+    bool isEnforcedNode = false;
     if ( tag > 0 )
     {
       if ( tag <= pmap.Extent() )
@@ -2664,10 +2708,8 @@ bool BLSURFPlugin_BLSURF::compute(SMESH_Mesh&         aMesh,
           tag = 0; // enforced or attracted vertex
         nodes[iv] = SMESH_Algo::VertexNode( v, meshDS );
       }
-      if ( !nodes[iv] && ( nodes[iv] = enforcedMesh.GetNodeByTag( tag, pmap )))
-      {
-        continue;
-      }
+      if ( !nodes[iv] )
+        isEnforcedNode = ( nodes[iv] = enforcedMesh.GetNodeByTag( tag, pmap ));
     }
     if ( !nodes[iv] )
       nodes[iv] = meshDS->AddNode(xyz[0], xyz[1], xyz[2]);
@@ -2699,6 +2741,8 @@ bool BLSURFPlugin_BLSURF::compute(SMESH_Mesh&         aMesh,
           groupDS->Add( nodes[iv] );
         }
     }
+    if ( isEnforcedNode )
+      continue;
 
     // internal points are tagged to zero
     if ( tag > 0 && tag <= pmap.Extent() )
@@ -2958,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 )
@@ -3061,6 +3113,11 @@ bool BLSURFPlugin_BLSURF::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelpe
       meshgems_mesh_set_quadrangle_vertices( msh, iQ++, nodeIDs );
   }
 
+  std::string errorTxt;
+
+  if ( !SMESHUtils_MGLicenseKeyGen::SignMesh( msh, "cadsurf", errorTxt ))
+    return error( "Problem with library SalomeMeshGemsKeyGenerator: " + errorTxt );
+
   ret = cadsurf_set_mesh(css, msh);
   if ( ret != STATUS_OK ) return error("Pb in cadsurf_set_mesh()");
 
@@ -3481,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;
   }
 }
 
@@ -3663,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 );
+}