Salome HOME
Corrections to handle adapation to shape
authorYoann Audouin <yoann.audouin@edf.fr>
Thu, 6 Oct 2022 06:43:27 +0000 (08:43 +0200)
committerYoann Audouin <yoann.audouin@edf.fr>
Thu, 13 Oct 2022 12:33:36 +0000 (14:33 +0200)
idl/SMESH_Gen.idl
src/SMESHGUI/SMESHGUI_CreateDualMeshOp.cxx
src/SMESH_I/SMESH_Gen_i.cxx
src/SMESH_I/SMESH_Gen_i.hxx
src/SMESH_SWIG/smeshBuilder.py
src/SMESH_SWIG/smesh_tools.py
test/SMESH_create_dual_mesh.py
test/SMESH_create_dual_mesh_adapt.py [new file with mode: 0644]
test/tests.set

index 010525eaf1d8360aeba5abddcc41d55048270973..a5b48aaa6bf75f81c3e09fd8a28488e0f9a5b8ae 100644 (file)
@@ -288,9 +288,11 @@ module SMESH
      * Create a dual mesh of a Tetrahedron mesh
      *  \param mesh - TetraHedron mesh to create dual from
      *  \param meshName - a name of the new mesh
+     *  \param adaptToShape - if True project boundary point on shape
      */
     SMESH_Mesh CreateDualMesh(in SMESH_IDSource mesh,
-                              in string         meshName)
+                              in string         meshName,
+                              in boolean        adaptToShape)
       raises ( SALOME::SALOME_Exception );
 
     /*!
index 449f0a58233eef7bf6b41a815c084acd4c45bd99..4ee81382a7bf2afdc4cf92eb267af4ccf9d916cc 100644 (file)
@@ -45,6 +45,7 @@
 
 // Qt includes
 #include <QLineEdit>
+#include <QCheckBox>
 
 // IDL includes
 #include <SALOMEconfig.h>
@@ -226,10 +227,10 @@ bool SMESHGUI_CreateDualMeshOp::onApply()
   SMESH::SMESH_Gen_var gen = SMESHGUI::GetSMESHGen();
   SMESH::SMESH_Mesh_var newMesh;
   QByteArray newMeshName=myDlg->myMeshName->text().toUtf8();
+  bool adapt_to_shape=myDlg->myProjShape->isChecked();
   try
   {
-    // TODO: change name as previous name + "_dual"
-    newMesh = gen->CreateDualMesh(mesh, newMeshName.constData());
+    newMesh = gen->CreateDualMesh(mesh, newMeshName.constData(), adapt_to_shape);
 
     if ( !newMesh->_is_nil() )
       if ( _PTR(SObject) aSObject = SMESH::ObjectToSObject( newMesh ) )
index 7157f413c0481beecdd7937c58392fa29c799ca3..0fd53cbc7ae57124da6517ee62723282f4d064b9 100644 (file)
@@ -2820,7 +2820,8 @@ SMESH_Gen_i::ConcatenateCommon(const SMESH::ListOfIDSources& theMeshesArray,
 //================================================================================
 
 SMESH::SMESH_Mesh_ptr SMESH_Gen_i::CreateDualMesh(SMESH::SMESH_IDSource_ptr mesh,
-                                                  const char*               meshName)
+                                                  const char*               meshName,
+                                                  CORBA::Boolean            adapt_to_shape)
 {
   Unexpect aCatch(SALOME_SalomeException);
 
@@ -2841,10 +2842,12 @@ SMESH::SMESH_Mesh_ptr SMESH_Gen_i::CreateDualMesh(SMESH::SMESH_IDSource_ptr mesh
   CORBA::String_var mesh_var=GetORB()->object_to_string(mesh);
   std::string mesh_ior = mesh_var.in();
 
-  fs::path tmp_folder = fs::unique_path(fs::path("dual_mesh-%%%%-%%%%"));
+  //temporary folder for the generation of the med file
+  fs::path tmp_folder = fs::temp_directory_path() / fs::unique_path(fs::path("dual_mesh-%%%%"));
   fs::create_directories(tmp_folder);
   fs::path dual_mesh_file = tmp_folder / fs::path("tmp_dual_mesh.med");
-  std::string mesh_name = "MESH";
+  std::string mesh_name = meshName;
+  MESSAGE("Working in folder" + tmp_folder.string());
 
   // Running Python script
   assert(Py_IsInitialized());
@@ -2853,11 +2856,20 @@ SMESH::SMESH_Mesh_ptr SMESH_Gen_i::CreateDualMesh(SMESH::SMESH_IDSource_ptr mesh
 
   std::string cmd="import salome.smesh.smesh_tools as smt";
   PyRun_SimpleString(cmd.c_str());
+  std::string ats;
+  if(adapt_to_shape)
+    ats = "True";
+  else
+    ats = "False";
 
-  cmd = "smt.create_dual_mesh(\"" + mesh_ior + "\", \"" + dual_mesh_file.string() + "\", \"" + mesh_name + "\")";
+  cmd = "smt.smesh_create_dual_mesh(\"" + mesh_ior + "\", \"" +
+        dual_mesh_file.string() + "\", mesh_name=\"" + mesh_name + "\", adapt_to_shape=" + ats + ")";
+  MESSAGE(cmd);
   PyRun_SimpleString(cmd.c_str());
 
   PyGILState_Release(gstate);
+  MESSAGE("Executed python script");
+  MESSAGE("Mesh created in " + dual_mesh_file.string());
 
   // Import created MED
   SMESH::SMESH_Mesh_var newMesh = CreateMesh(GEOM::GEOM_Object::_nil());
@@ -2867,21 +2879,25 @@ SMESH::SMESH_Mesh_ptr SMESH_Gen_i::CreateDualMesh(SMESH::SMESH_IDSource_ptr mesh
   SALOMEDS::SObject_wrap meshSO = ObjectToSObject( newMesh );
   if ( !meshSO->_is_nil() )
   {
-    SetName( meshSO, meshName, "Mesh" );
-    SetPixMap( meshSO, "ICON_SMESH_TREE_MESH");
+    SetName( meshSO, meshName, meshName );
+    SetPixMap( meshSO, "ICON_SMESH_TREE_MESH_IMPORTED");
   }
 
   SMESH_Mesh& newMesh2 = newMesh_i->GetImpl();
 
   newMesh2.MEDToMesh(dual_mesh_file.c_str(), meshName);
 
+  MESSAGE("Imported created MED")
+
   SMESHDS_Mesh* newMeshDS = newMesh_i->GetImpl().GetMeshDS();
 
   newMeshDS->Modified();
 
   *pyDump << newMesh << " = " << this
           << ".CreateDualMesh("
-          << "'" << meshName << "') ";
+          << mesh << ", "
+          << "'" << meshName << "', "
+          << ats << ") ";
 
   return newMesh._retn();
 }
index bd887203878514873887709b6c27d4c7961d9cd1..822e71f2e51c96ca2a0cf6cdd6d8e9f38e2f8b1f 100644 (file)
@@ -255,7 +255,8 @@ public:
 
   // Create dual mesh of a tetrahedron mesh
   SMESH::SMESH_Mesh_ptr CreateDualMesh(SMESH::SMESH_IDSource_ptr meshPart,
-                                       const char*               meshName);
+                                       const char*               meshName,
+                                       CORBA::Boolean            adapt_to_shape);
 
   // Copy a part of mesh
   SMESH::SMESH_Mesh_ptr CopyMesh(SMESH::SMESH_IDSource_ptr meshPart,
index 42bad03cf12a85a0011bb41dd7355fa82892ef57..62db8895183d921e4967b74ec190ee4713171e83 100644 (file)
@@ -773,7 +773,7 @@ class smeshBuilder( SMESH._objref_SMESH_Gen, object ):
         aMesh = Mesh( self, self.geompyD, aSmeshMesh, name=name )
         return aMesh
 
-    def CreateDualMesh( self, mesh, meshName):
+    def CreateDualMesh( self, mesh, meshName, adaptToShape):
         """
         Create a dual of a mesh.
 
@@ -782,13 +782,15 @@ class smeshBuilder( SMESH._objref_SMESH_Gen, object ):
                         :class:`mesh, <SMESH.SMESH_IDSource>`.
 
                 meshName: a name of the new mesh
+                adpatToShape: if true project boundary points on shape
 
         Returns:
                 an instance of class :class:`Mesh`
         """
         if isinstance( mesh, Mesh ):
             mesh = mesh.GetMesh()
-        dualMesh = SMESH._objref_SMESH_Gen.CreateDualMesh(self, mesh, meshName)
+        print("calling createdualmesh from Python")
+        dualMesh = SMESH._objref_SMESH_Gen.CreateDualMesh(self, mesh, meshName, adaptToShape)
         return Mesh(self, self.geompyD, dualMesh)
 
 
index 38cd978e0748d8f879d968544707d91a557d0d71..1a6995eedcfe55d695791988efa21ab142590ac6 100644 (file)
@@ -17,7 +17,7 @@ from salome.smesh import smeshBuilder
 
 smesh = smeshBuilder.New()
 
-def create_dual_mesh(mesh_ior, output_file, adapt_to_shape=True, mesh_name="MESH"):
+def smesh_create_dual_mesh(mesh_ior, output_file, adapt_to_shape=True, mesh_name="MESH"):
     """ Create a dual of the mesh in input_file into output_file
 
     Args:
@@ -25,7 +25,6 @@ def create_dual_mesh(mesh_ior, output_file, adapt_to_shape=True, mesh_name="MESH
         output_file (string): dual mesh file
     """
     # Import mesh from file
-    # mesh = salome.orb.string_to_object(salome.salome_study.myStudy.FindObjectID(mesh_id).GetIOR())
     mesh = salome.orb.string_to_object(mesh_ior)
     if not mesh:
         raise Exception("Could not find mesh using id: ", mesh_ior)
@@ -35,11 +34,14 @@ def create_dual_mesh(mesh_ior, output_file, adapt_to_shape=True, mesh_name="MESH
     # We got a meshProxy so we need to convert pointer to MEDCoupling
     int_ptr = mesh.ExportMEDCoupling(True, True)
     dab = mc.FromPyIntPtrToDataArrayByte(int_ptr)
-    tetras =  mc.MEDFileMesh.New(dab)[0]
+    tetras = mc.MEDFileMesh.New(dab)[0]
     # End of SMESH -> MEDCoupling part for dualmesh
 
     tetras = mc.MEDCoupling1SGTUMesh(tetras)
     polyh = tetras.computeDualMesh()
+    dual_volume_raw = polyh.getMeasureField(True).accumulate()[0]
+
+    # Getting list of new points added on the skin
     skin = tetras.buildUnstructured().computeSkin()
     skin_polyh = polyh.buildUnstructured().computeSkin()
     allNodesOnSkinPolyh = skin_polyh.computeFetchedNodeIds()
@@ -51,24 +53,25 @@ def create_dual_mesh(mesh_ior, output_file, adapt_to_shape=True, mesh_name="MESH
         ptsAddedCoo = ptsAddedMesh.getCoords()
         ptsAddedCooModified = ptsAddedCoo[:]
 
-        # We need the geometry for that
-        # TODO : Loop on faces identify points associated to which face
+        # Matching faces with their ids
         faces = geompy.ExtractShapes(shape, geompy.ShapeType["FACE"], True)
-        #assert( len(faces) == 1 )
-        ## projection des points ajoutés par le dual sur la surface
-        #for i,tup in enumerate(ptsAddedCooModified):
-        #    vertex = geompy.MakeVertex(*tuple(tup))
-        #    prj = geompy.MakeProjection(vertex, faces)
-        #    newCoor = geompy.PointCoordinates( prj )
-        #    ptsAddedCooModified[i] = newCoor
-        ## assign coordinates with projected ones
-        #polyh.getCoords()[ptsAdded] = ptsAddedCooModified
-
-    print("Writing dual mesh in ", output_file)
+        id2face = {}
+        for face in faces:
+            id2face[face.GetSubShapeIndices()[0]] = face
+        print(id2face)
+
+        ## Projecting each points added by the dual mesh on the surface it is
+        # associated with
+        for i, tup in enumerate(ptsAddedCooModified):
+            vertex = geompy.MakeVertex(*tuple(tup))
+            shapes = geompy.GetShapesNearPoint(shape, vertex,
+                                               geompy.ShapeType["FACE"])
+            prj = geompy.MakeProjection(vertex,
+                                        id2face[shapes.GetSubShapeIndices()[0]])
+            new_coor = geompy.PointCoordinates(prj)
+            ptsAddedCooModified[i] = new_coor
+
+        polyh.getCoords()[ptsAdded] = ptsAddedCooModified
+
     polyh.setName(mesh_name)
     polyh.write(output_file)
-
-
-
-
-
index 007dbb8a7337470a78a0fe540c7e56c3e4c60481..08621f30ee061417e0f936f625a29e943b6d1f4f 100644 (file)
@@ -50,7 +50,7 @@ Mesh_1 = smesh.Mesh(Sphere_1,'Mesh_1')
 status = Mesh_1.AddHypothesis( Sphere_1, NETGEN_3D_Parameters_1 )
 NETGEN_1D_2D_3D = Mesh_1.Tetrahedron(algo=smeshBuilder.NETGEN_1D2D3D)
 isDone = Mesh_1.Compute()
-dual_Mesh_1 = smesh.CreateDualMesh( Mesh_1, 'dual_Mesh_1')
+dual_Mesh_1 = smesh.CreateDualMesh( Mesh_1, 'dual_Mesh_1', True)
 
 
 assert(dual_Mesh_1.NbPolyhedrons() > 0)
diff --git a/test/SMESH_create_dual_mesh_adapt.py b/test/SMESH_create_dual_mesh_adapt.py
new file mode 100644 (file)
index 0000000..c95e6a6
--- /dev/null
@@ -0,0 +1,92 @@
+#!/usr/bin/env python
+
+###
+### This file is generated automatically by SALOME v9.9.0 with dump python functionality
+###
+
+import sys
+import salome
+
+salome.salome_init()
+import salome_notebook
+notebook = salome_notebook.NoteBook()
+sys.path.insert(0, r'/home/B61570/work_in_progress/dual_mesh')
+
+###
+### GEOM component
+###
+
+import GEOM
+from salome.geom import geomBuilder
+import math
+import SALOMEDS
+
+
+geompy = geomBuilder.New()
+
+O = geompy.MakeVertex(0, 0, 0)
+OX = geompy.MakeVectorDXDYDZ(1, 0, 0)
+OY = geompy.MakeVectorDXDYDZ(0, 1, 0)
+OZ = geompy.MakeVectorDXDYDZ(0, 0, 1)
+Cylinder_1 = geompy.MakeCylinderRH(100, 400)
+Sphere_1 = geompy.MakeSpherePntR(O, 100)
+Fuse_1 = geompy.MakeFuseList([Cylinder_1, Sphere_1], True, True)
+
+[geomObj_1,geomObj_2,geomObj_3] = geompy.ExtractShapes(Fuse_1, geompy.ShapeType["FACE"], True)
+
+top = geompy.CreateGroup(Fuse_1, geompy.ShapeType["FACE"])
+geompy.UnionIDs(top, geomObj_1.GetSubShapeIndices())
+
+middle = geompy.CreateGroup(Fuse_1, geompy.ShapeType["FACE"])
+geompy.UnionIDs(middle, geomObj_2.GetSubShapeIndices())
+
+bottom = geompy.CreateGroup(Fuse_1, geompy.ShapeType["FACE"])
+geompy.UnionIDs(bottom, geomObj_3.GetSubShapeIndices())
+
+#[top, middle, bottom] = geompy.GetExistingSubObjects(Fuse_1, False)
+
+geompy.addToStudy( O, 'O' )
+geompy.addToStudy( OX, 'OX' )
+geompy.addToStudy( OY, 'OY' )
+geompy.addToStudy( OZ, 'OZ' )
+geompy.addToStudy( Cylinder_1, 'Cylinder_1' )
+geompy.addToStudy( Sphere_1, 'Sphere_1' )
+geompy.addToStudy( Fuse_1, 'Fuse_1' )
+geompy.addToStudyInFather( Fuse_1, top, 'top' )
+geompy.addToStudyInFather( Fuse_1, middle, 'middle' )
+geompy.addToStudyInFather( Fuse_1, bottom, 'bottom' )
+
+###
+### SMESH component
+###
+
+import  SMESH, SALOMEDS
+from salome.smesh import smeshBuilder
+
+smesh = smeshBuilder.New()
+
+NETGEN_3D_Parameters_1 = smesh.CreateHypothesisByAverageLength( 'NETGEN_Parameters', 'NETGENEngine', 50, 0 )
+Mesh_1 = smesh.Mesh(Fuse_1,'Mesh_1')
+status = Mesh_1.AddHypothesis( Fuse_1, NETGEN_3D_Parameters_1 )
+NETGEN_1D_2D_3D = Mesh_1.Tetrahedron(algo=smeshBuilder.NETGEN_1D2D3D)
+top_1 = Mesh_1.GroupOnGeom(top,'top',SMESH.FACE)
+middle_1 = Mesh_1.GroupOnGeom(middle,'middle',SMESH.FACE)
+bottom_1 = Mesh_1.GroupOnGeom(bottom,'bottom',SMESH.FACE)
+isDone = Mesh_1.Compute()
+
+[ top_1, middle_1, bottom_1 ] = Mesh_1.GetGroups()
+
+
+dual_Mesh_1 = smesh.CreateDualMesh(Mesh_1, 'dual_Mesh_1', True)
+dual_Mesh_raw_1 = smesh.CreateDualMesh(Mesh_1, 'dual_Mesh_1', False)
+
+#Comparing volumes
+dual_volume = dual_Mesh_1.GetVolume()
+dual_raw_volume = dual_Mesh_raw_1.GetVolume()
+print("dual_volume: ", dual_volume)
+print("dual_raw_volume: ", dual_raw_volume)
+
+assert (dual_volume >= dual_raw_volume)
+
+if salome.sg.hasDesktop():
+  salome.sg.updateObjBrowser()
index f00bde4c6e4e22001fa81c2761b3426dc8b08f80..977f4a667fda7bc7e95d4ee8499b263bb8139e9b 100644 (file)
@@ -63,6 +63,7 @@ SET(BAD_TESTS
   SMESH_test2.py
   SMESH_test4.py
   SMESH_create_dual_mesh.py
+  SMESH_create_dual_mesh_adapt.py
   )
 IF(NOT WIN32)
   LIST(APPEND BAD_TESTS