[FEM] proper support for Gmsh HighOrder optimization

There are 4 optimizers available while we only supported 1. This PR adds support for all.
This commit is contained in:
donovaly
2021-04-03 03:55:15 +02:00
committed by Bernd Hahnebach
parent 3a36602ae8
commit cc0af9b4dd
4 changed files with 49 additions and 19 deletions

View File

@@ -125,6 +125,21 @@ class GmshTools():
else:
self.algorithm3D = "1"
# HighOrderOptimize
optimizers = self.mesh_obj.HighOrderOptimize
if optimizers == "None":
self.HighOrderOptimize = "0"
elif optimizers == "Optimization":
self.HighOrderOptimize = "1"
elif optimizers == "Elastic+Optimization":
self.HighOrderOptimize = "2"
elif optimizers == "Elastic":
self.HighOrderOptimize = "3"
elif optimizers == "Fast Curving":
self.HighOrderOptimize = "4"
else:
self.HighOrderOptimize = "0"
# mesh groups
if self.mesh_obj.GroupsOfNodes is True:
self.group_nodes_export = True
@@ -767,16 +782,11 @@ class GmshTools():
else:
geo.write("Mesh.OptimizeNetgen = 0;\n")
# higher order mesh optimizing
if hasattr(self.mesh_obj, "HighOrderOptimize") and self.mesh_obj.HighOrderOptimize is True:
geo.write(
"Mesh.HighOrderOptimize = 1; // for more HighOrderOptimize "
"parameter check http://gmsh.info/doc/texinfo/gmsh.html\n"
)
else:
geo.write(
"Mesh.HighOrderOptimize = 0; // for more HighOrderOptimize "
"parameter check http://gmsh.info/doc/texinfo/gmsh.html\n"
)
geo.write(
"// High-order meshes optimization (0=none, 1=optimization, 2=elastic+optimization, "
"3=elastic, 4=fast curving)\n"
)
geo.write("Mesh.HighOrderOptimize = " + self.HighOrderOptimize + ";\n")
geo.write("\n")
geo.write("// mesh order\n")
@@ -850,7 +860,7 @@ class GmshTools():
# some useful information
geo.write("// " + "*" * 70 + "\n")
geo.write("// Gmsh documentation:\n")
geo.write("// http://gmsh.info/doc/texinfo/gmsh.html#Mesh\n")
geo.write("// https://gmsh.info/doc/texinfo/gmsh.html#Mesh\n")
geo.write("//\n")
geo.write(
"// We do not check if something went wrong, like negative "

View File

@@ -60,16 +60,34 @@ class MeshGmsh(base_fempythonobject.BaseFemPythonObject):
"R-tree",
"HXT"
]
known_mesh_HighOrderOptimizers = [
"None",
"Optimization",
"Elastic+Optimization",
"Elastic",
"Fast curving"
]
def __init__(self, obj):
super(MeshGmsh, self).__init__(obj)
self.add_properties(obj)
def onDocumentRestored(self, obj):
# HighOrderOptimize was once App::PropertyBool, so check this
HighOrderOptimizer = ""
if obj.HighOrderOptimize is True:
HighOrderOptimizer = "Optimization"
obj.removeProperty("HighOrderOptimize")
elif obj.HighOrderOptimize is False:
HighOrderOptimizer = "None"
obj.removeProperty("HighOrderOptimize")
self.add_properties(obj)
# refresh the list of known 3D algorithms for existing meshes
# since some algos are meanwhile deprecated and new algos are available
obj.Algorithm3D = MeshGmsh.known_mesh_algorithm_3D
# write the stored HighOrderOptimizer
if HighOrderOptimizer:
obj.HighOrderOptimize = HighOrderOptimizer
def add_properties(self, obj):
if not hasattr(obj, "MeshBoundaryLayerList"):
@@ -166,12 +184,13 @@ class MeshGmsh(base_fempythonobject.BaseFemPythonObject):
if not hasattr(obj, "HighOrderOptimize"):
obj.addProperty(
"App::PropertyBool",
"App::PropertyEnumeration",
"HighOrderOptimize",
"FEM Gmsh Mesh Params",
"Optimize high order meshes"
"Optimization of high order meshes"
)
obj.HighOrderOptimize = False
obj.HighOrderOptimize = MeshGmsh.known_mesh_HighOrderOptimizers
obj.HighOrderOptimize = "None"
if not hasattr(obj, "RecombineAll"):
obj.addProperty(
@@ -232,7 +251,7 @@ class MeshGmsh(base_fempythonobject.BaseFemPythonObject):
"mesh algorithm 2D"
)
obj.Algorithm2D = MeshGmsh.known_mesh_algorithm_2D
obj.Algorithm2D = "Automatic" # ?
obj.Algorithm2D = "Automatic"
if not hasattr(obj, "Algorithm3D"):
obj.addProperty(
@@ -242,7 +261,7 @@ class MeshGmsh(base_fempythonobject.BaseFemPythonObject):
"mesh algorithm 3D"
)
obj.Algorithm3D = MeshGmsh.known_mesh_algorithm_3D
obj.Algorithm3D = "Automatic" # ?
obj.Algorithm3D = "Automatic"
if not hasattr(obj, "GroupsOfNodes"):
obj.addProperty(

View File

@@ -19,7 +19,8 @@ Mesh.MeshSizeFromCurvature = 12; // number of elements per 2*pi radians, 0 to de
// optimize the mesh
Mesh.Optimize = 1;
Mesh.OptimizeNetgen = 0;
Mesh.HighOrderOptimize = 0; // for more HighOrderOptimize parameter check http://gmsh.info/doc/texinfo/gmsh.html
// High-order meshes optimization (0=none, 1=optimization, 2=elastic+optimization, 3=elastic, 4=fast curving)
Mesh.HighOrderOptimize = 0;
// mesh order
Mesh.ElementOrder = 2;
@@ -45,7 +46,7 @@ Save "/tmp/tmpjVhNNb.unv";
// **********************************************************************
// Gmsh documentation:
// http://gmsh.info/doc/texinfo/gmsh.html#Mesh
// https://gmsh.info/doc/texinfo/gmsh.html#Mesh
//
// We do not check if something went wrong, like negative jacobians etc. You can run Gmsh manually yourself:
//

View File

@@ -514,7 +514,7 @@ bool Mesh2ShapeGmsh::writeProject(QString& inpFile, QString& outFile)
<< "// optimize the mesh\n"
<< "Mesh.Optimize = 1;\n"
<< "Mesh.OptimizeNetgen = 0;\n"
<< "// for more HighOrderOptimize parameter check http://gmsh.info/doc/texinfo/gmsh.html\n"
<< "// High-order meshes optimization (0=none, 1=optimization, 2=elastic+optimization, 3=elastic, 4=fast curving)\n"
<< "Mesh.HighOrderOptimize = 0;\n\n"
<< "// mesh order\n"
<< "Mesh.ElementOrder = 2;\n"