Skip to content

Commit

Permalink
Merge pull request #615 from RemDelaporteMathurin/fix-ufl-cell-vs-bas…
Browse files Browse the repository at this point in the history
…ix-cell

Fix ufl cell vs basix cell
  • Loading branch information
RemDelaporteMathurin authored Oct 20, 2023
2 parents 4d4a65a + 8d621ff commit b26c407
Show file tree
Hide file tree
Showing 3 changed files with 22 additions and 11 deletions.
8 changes: 7 additions & 1 deletion festim/hydrogen_transport_problem.py
Original file line number Diff line number Diff line change
@@ -1,6 +1,7 @@
from dolfinx import fem
from dolfinx.nls.petsc import NewtonSolver
from dolfinx.io import XDMFFile
import basix
import ufl
from mpi4py import MPI
from dolfinx.fem import Function, form, assemble_scalar
Expand Down Expand Up @@ -133,7 +134,12 @@ def defing_export_writers(self):
export.writer.write_mesh(self.mesh.mesh)

def define_function_space(self):
elements = ufl.FiniteElement("CG", self.mesh.mesh.ufl_cell(), 1)
elements = basix.ufl.element(
basix.ElementFamily.P,
self.mesh.mesh.basix_cell(),
1,
basix.LagrangeVariant.equispaced,
)
self.function_space = fem.FunctionSpace(self.mesh.mesh, elements)

def define_markers_and_measures(self):
Expand Down
11 changes: 7 additions & 4 deletions test/benchmark.py
Original file line number Diff line number Diff line change
Expand Up @@ -17,18 +17,16 @@
from dolfinx.nls.petsc import NewtonSolver
from ufl import (
dot,
FiniteElement,
grad,
TestFunction,
exp,
FacetNormal,
dx,
ds,
Cell,
Mesh,
VectorElement,
Measure,
)
import basix
from dolfinx.mesh import create_mesh, meshtags, locate_entities
import numpy as np
import tqdm.autonotebook
Expand All @@ -50,7 +48,12 @@ def fenics_test_permeation_problem(mesh_size=1001):
vdim = my_mesh.topology.dim
n = FacetNormal(my_mesh)

elements = FiniteElement("CG", my_mesh.ufl_cell(), 1)
elements = basix.ufl.element(
basix.ElementFamily.P,
my_mesh.basix_cell(),
1,
basix.LagrangeVariant.equispaced,
)
V = FunctionSpace(my_mesh, elements)
u = Function(V)
u_n = Function(V)
Expand Down
14 changes: 8 additions & 6 deletions test/test_permeation_problem.py
Original file line number Diff line number Diff line change
Expand Up @@ -75,11 +75,12 @@ def test_permeation_problem(mesh_size=1001):
analytical_flux = np.abs(analytical_flux)
flux_values = np.array(np.abs(flux_values))

indices = np.where(analytical_flux > 0.01 * np.max(analytical_flux))
analytical_flux = analytical_flux[indices]
flux_values = flux_values[indices]

relative_error = np.abs((flux_values - analytical_flux) / analytical_flux)

relative_error = relative_error[
np.where(analytical_flux > 0.01 * np.max(analytical_flux))
]
error = relative_error.mean()

assert error < 0.01
Expand Down Expand Up @@ -165,11 +166,12 @@ def test_permeation_problem_multi_volume():
analytical_flux = np.abs(analytical_flux)
flux_values = np.array(np.abs(flux_values))

indices = np.where(analytical_flux > 0.01 * np.max(analytical_flux))
analytical_flux = analytical_flux[indices]
flux_values = flux_values[indices]

relative_error = np.abs((flux_values - analytical_flux) / analytical_flux)

relative_error = relative_error[
np.where(analytical_flux > 0.01 * np.max(analytical_flux))
]
error = relative_error.mean()

assert error < 0.01

0 comments on commit b26c407

Please sign in to comment.