Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Mesh properties fdim and vdim #599

Merged
merged 6 commits into from
Oct 13, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
16 changes: 8 additions & 8 deletions festim/hydrogen_transport_problem.py
Original file line number Diff line number Diff line change
Expand Up @@ -120,12 +120,8 @@ def define_markers_and_measures(self):

dofs_facets, tags_facets = [], []

# TODO this should be a property of mesh
fdim = self.mesh.mesh.topology.dim - 1
vdim = self.mesh.mesh.topology.dim

# find all cells in domain and mark them as 0
num_cells = self.mesh.mesh.topology.index_map(vdim).size_local
num_cells = self.mesh.mesh.topology.index_map(self.mesh.vdim).size_local
mesh_cell_indices = np.arange(num_cells, dtype=np.int32)
tags_volumes = np.full(num_cells, 0, dtype=np.int32)

Expand All @@ -137,17 +133,21 @@ def define_markers_and_measures(self):
if isinstance(sub_dom, F.VolumeSubdomain1D):
# find all cells in subdomain and mark them as sub_dom.id
self.volume_subdomains.append(sub_dom)
entities = sub_dom.locate_subdomain_entities(self.mesh.mesh, vdim)
entities = sub_dom.locate_subdomain_entities(
self.mesh.mesh, self.mesh.vdim
)
tags_volumes[entities] = sub_dom.id

# dofs and tags need to be in np.in32 format for meshtags
dofs_facets = np.array(dofs_facets, dtype=np.int32)
tags_facets = np.array(tags_facets, dtype=np.int32)

# define mesh tags
self.facet_meshtags = meshtags(self.mesh.mesh, fdim, dofs_facets, tags_facets)
self.facet_meshtags = meshtags(
self.mesh.mesh, self.mesh.fdim, dofs_facets, tags_facets
)
self.volume_meshtags = meshtags(
self.mesh.mesh, vdim, mesh_cell_indices, tags_volumes
self.mesh.mesh, self.mesh.vdim, mesh_cell_indices, tags_volumes
)

# define measures
Expand Down
15 changes: 11 additions & 4 deletions festim/mesh/mesh.py
Original file line number Diff line number Diff line change
@@ -1,15 +1,14 @@
import ufl


class Mesh:
"""
Mesh class

Args:
mesh (dolfinx.mesh.Mesh, optional): the mesh. Defaults to None.
mesh (dolfinx.mesh.Mesh, optional): the mesh. Defaults to None.

Attributes:
mesh (dolfinx.mesh.Mesh): the mesh
vdim (int): the dimension of the mesh cells
fdim (int): the dimension of the mesh facets
"""

def __init__(self, mesh=None):
Expand All @@ -25,3 +24,11 @@ def __init__(self, mesh=None):
self.mesh.topology.create_connectivity(
self.mesh.topology.dim - 1, self.mesh.topology.dim
)

@property
def vdim(self):
return self.mesh.topology.dim

@property
def fdim(self):
return self.mesh.topology.dim - 1
4 changes: 2 additions & 2 deletions festim/mesh/mesh_1d.py
Original file line number Diff line number Diff line change
@@ -1,4 +1,4 @@
from dolfinx import fem, mesh
from dolfinx import mesh
from mpi4py import MPI
import ufl
import numpy as np
Expand All @@ -10,7 +10,7 @@ class Mesh1D(Mesh):
1D Mesh

Args:
vertices (list): the mesh x-coordinates (m)
vertices (list): the mesh x-coordinates (m)

Attributes:
vertices (list): the mesh x-coordinates (m)
Expand Down
38 changes: 38 additions & 0 deletions test/test_mesh.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,38 @@
import festim as F
from dolfinx import mesh as fenics_mesh
from mpi4py import MPI
import pytest

mesh_1D = fenics_mesh.create_unit_interval(MPI.COMM_WORLD, 10)
mesh_2D = fenics_mesh.create_unit_square(MPI.COMM_WORLD, 10, 10)
mesh_3D = fenics_mesh.create_unit_cube(MPI.COMM_WORLD, 10, 10, 10)


@pytest.mark.parametrize("mesh", [mesh_1D, mesh_2D, mesh_3D])
def test_get_fdim(mesh):
my_mesh = F.Mesh(mesh)

assert my_mesh.fdim == mesh.topology.dim - 1


def test_fdim_changes_when_mesh_changes():
my_mesh = F.Mesh()

for mesh in [mesh_1D, mesh_2D, mesh_3D]:
my_mesh.mesh = mesh
assert my_mesh.fdim == mesh.topology.dim - 1


@pytest.mark.parametrize("mesh", [mesh_1D, mesh_2D, mesh_3D])
def test_get_vdim(mesh):
my_mesh = F.Mesh(mesh)

assert my_mesh.vdim == mesh.topology.dim


def test_vdim_changes_when_mesh_changes():
my_mesh = F.Mesh()

for mesh in [mesh_1D, mesh_2D, mesh_3D]:
my_mesh.mesh = mesh
assert my_mesh.vdim == mesh.topology.dim
Loading