Skip to content
GitLab
Projects
Groups
Snippets
Help
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
HO_homog
Project overview
Project overview
Details
Activity
Releases
Repository
Repository
Files
Commits
Branches
Tags
Contributors
Graph
Compare
Issues
14
Issues
14
List
Boards
Labels
Service Desk
Milestones
Operations
Operations
Incidents
Analytics
Analytics
Repository
Value Stream
Wiki
Wiki
Snippets
Snippets
Members
Members
Collapse sidebar
Close sidebar
Activity
Graph
Create a new issue
Commits
Issue Boards
Open sidebar
Baptiste Durand
HO_homog
Commits
e3e3a3de
Commit
e3e3a3de
authored
Nov 26, 2019
by
Baptiste Durand
Browse files
Options
Browse Files
Download
Email Patches
Plain Diff
Use of mesh conversion and mesh import tools for the definition of Part instances
parent
27d3edec
Changes
3
Hide whitespace changes
Inline
Side-by-side
Showing
3 changed files
with
87 additions
and
40 deletions
+87
-40
ho_homog/part.py
ho_homog/part.py
+31
-12
ho_homog/toolbox_FEniCS.py
ho_homog/toolbox_FEniCS.py
+3
-0
ho_homog/toolbox_gmsh.py
ho_homog/toolbox_gmsh.py
+53
-28
No files found.
ho_homog/part.py
View file @
e3e3a3de
...
...
@@ -13,6 +13,7 @@ import dolfin as fe
import
ho_homog.materials
as
mat
from
subprocess
import
run
from
pathlib
import
Path
from
ho_homog.toolbox_gmsh
import
msh_conversion
plt
.
ioff
()
...
...
@@ -98,21 +99,39 @@ class FenicsPart(object):
if
not
isinstance
(
mesh_path
,
Path
):
mesh_path
=
Path
(
mesh_path
)
name
=
mesh_path
.
stem
if
mesh_path
.
suffix
==
"xml"
:
mesh
=
fe
.
Mesh
(
str
(
mesh_path
))
elif
mesh_path
.
suffix
==
".xdmf"
:
mesh
=
fe
.
Mesh
()
with
fe
.
XDMFFile
(
str
(
mesh_path
))
as
file_in
:
file_in
.
read
(
mesh
)
else
:
cmd
=
f
"dolfin-convert
{
mesh_path
}
{
mesh_path
.
with_suffix
(
'.xml'
)
}
"
run
(
cmd
,
shell
=
True
,
check
=
True
)
mesh_path
=
mesh_path
.
with_suffix
(
".xml"
)
suffix
=
mesh_path
.
suffix
if
suffix
not
in
fetools
.
SUPPORTED_MESH_SUFFIX
:
# ! TOCOMMIT : UTILISATION DE L'OUTIL DE CONVERSION
mesh_file_paths
=
msh_conversion
(
mesh_path
,
format_
=
".xml"
,
subdomains
=
subdomains_import
)
try
:
mesh_path
=
mesh_file_paths
[
0
]
except
IndexError
as
error
:
mesh_path
=
mesh_file_paths
logger
.
warning
(
error
)
# Each supported mesh format -> one if structure
if
suffix
==
".xml"
:
mesh
=
fe
.
Mesh
(
str
(
mesh_path
))
if
suffix
==
".xdmf"
:
if
subdomains_import
:
mesh
,
subdomains
,
facets
=
fetools
.
xdmf_mesh
(
mesh_path
,
True
)
logger
.
info
(
f
"Import of a mesh from
{
mesh_path
}
file, with subdomains data"
)
else
:
mesh
=
fetools
.
xdmf_mesh
(
mesh_path
,
False
)
logger
.
info
(
f
"Import of a mesh from
{
mesh_path
}
file, without subdomains data"
)
if
subdomains_import
:
# ! Faire adaptation pour xdmf
subdo_path
=
mesh_path
.
with_name
(
name
+
"_physical_region.xml
"
)
facet_path
=
mesh_path
.
with_name
(
name
+
"_facet_region.xml
"
)
subdo_path
=
mesh_path
.
with_name
(
f
"
{
name
}
_physical_region
{
suffix
}
"
)
facet_path
=
mesh_path
.
with_name
(
f
"
{
name
}
_facet_region
{
suffix
}
"
)
if
subdo_path
.
exists
():
subdomains
=
fe
.
MeshFunction
(
"size_t"
,
mesh
,
str
(
subdo_path
))
subdo_val
=
fetools
.
get_MeshFunction_val
(
subdomains
)
...
...
ho_homog/toolbox_FEniCS.py
View file @
e3e3a3de
...
...
@@ -34,6 +34,9 @@ DOLFIN_LU_METHODS = {
}
SUPPORTED_MESH_SUFFIX
=
(
".xml"
,
".xdmf"
)
def
get_MeshFunction_val
(
msh_fctn
):
"""
Get information about the set of values that a given MeshFunction outputs.
...
...
ho_homog/toolbox_gmsh.py
View file @
e3e3a3de
...
...
@@ -43,8 +43,57 @@ def process_gmsh_log(gmsh_log: list, detect_error=True):
raise
AssertionError
(
"Gmsh logging messages signal errors."
)
def
conversion_to_xdmf
(
i_path
,
o_path
,
cell_reg
,
facet_reg
,
dim
,
subdomains
=
False
):
"""Convert a ".msh" mesh generated with Gmsh to a xdmf mesh file.
Parameters
----------
i_path : Path
o_path : Path
Desired path for the converted mesh
cell_reg : Path
Desired path of the extra file for cell regions if subdomains are retained.
facet_reg : Path
Desired path of the extra file for facet regions if subdomains are retained.
dim: int
Geometrical dimension of the mesh (2D or 3D).
subdomains : bool, optional
If True, extra files are created to store information about subdomains.
(default: False)
"""
m
=
meshio
.
read
(
str
(
i_path
))
if
dim
==
2
:
m
.
points
=
m
.
points
[:,
:
2
]
geo_only
=
meshio
.
Mesh
(
points
=
m
.
points
,
cells
=
{
"triangle"
:
m
.
cells
[
"triangle"
]})
cell
=
"triangle"
face
=
"line"
elif
dim
==
3
:
raise
NotImplementedError
(
"3D meshes are not supported yet."
)
else
:
ValueError
meshio
.
write
(
str
(
o_path
),
geo_only
)
if
subdomains
:
cell_funct
=
meshio
.
Mesh
(
points
=
m
.
points
,
cells
=
{
cell
:
m
.
cells
[
cell
]},
cell_data
=
{
cell
:
{
"cell_data"
:
m
.
cell_data
[
cell
][
"gmsh:physical"
]}},
)
meshio
.
write
(
str
(
cell_reg
),
cell_funct
)
facet_funct
=
meshio
.
Mesh
(
points
=
m
.
points
,
cells
=
{
face
:
m
.
cells
[
face
]},
cell_data
=
{
face
:
{
"facet_data"
:
m
.
cell_data
[
face
][
"gmsh:physical"
]}},
)
meshio
.
write
(
str
(
facet_reg
),
facet_funct
)
return
True
def
msh_conversion
(
mesh
,
format_
:
str
=
".xdmf"
,
output_dir
=
None
,
subdomains
:
bool
=
False
,
dim
:
int
=
2
mesh
,
format_
:
str
=
".xdmf"
,
output_dir
=
None
,
subdomains
:
bool
=
False
,
dim
:
int
=
2
,
):
"""
Convert a ".msh" mesh generated with Gmsh to a format suitable for FEniCS.
...
...
@@ -104,33 +153,9 @@ def msh_conversion(
cmd
=
f
"dolfin-convert
{
input_path
}
{
mesh_path
}
"
run
(
cmd
,
shell
=
True
,
check
=
True
)
elif
format_
==
".xdmf"
:
mesh
=
meshio
.
read
(
str
(
input_path
))
if
dim
==
2
:
mesh
.
points
=
mesh
.
points
[:,
:
2
]
geo_only
=
meshio
.
Mesh
(
points
=
mesh
.
points
,
cells
=
{
"triangle"
:
mesh
.
cells
[
"triangle"
]}
)
cell
=
"triangle"
face
=
"line"
elif
dim
==
3
:
raise
NotImplementedError
(
"3D meshes are not supported yet."
)
else
:
ValueError
meshio
.
write
(
str
(
mesh_path
),
geo_only
)
if
subdomains
:
cell_function
=
meshio
.
Mesh
(
points
=
mesh
.
points
,
cells
=
{
cell
:
mesh
.
cells
[
cell
]},
cell_data
=
{
cell
:
{
"cell_data"
:
mesh
.
cell_data
[
cell
][
"gmsh:physical"
]}},
)
meshio
.
write
(
str
(
physical_region
),
cell_function
)
facet_function
=
meshio
.
Mesh
(
points
=
mesh
.
points
,
cells
=
{
face
:
mesh
.
cells
[
face
]},
cell_data
=
{
face
:
{
"facet_data"
:
mesh
.
cell_data
[
face
][
"gmsh:physical"
]}},
)
meshio
.
write
(
str
(
facet_region
),
facet_function
)
conversion_to_xdmf
(
input_path
,
mesh_path
,
physical_region
,
facet_region
,
dim
,
subdomains
)
if
subdomains
:
return
(
mesh_path
,
...
...
Write
Preview
Markdown
is supported
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment