Skip to content

Commit 2f7a203

Browse files
committed
Mesh: extend P1 cell_sizes to high-order meshes
1 parent abd9e35 commit 2f7a203

File tree

2 files changed

+33
-3
lines changed

2 files changed

+33
-3
lines changed

firedrake/mesh.py

+21-2
Original file line numberDiff line numberDiff line change
@@ -2399,10 +2399,29 @@ def cell_sizes(self):
23992399
24002400
This is computed by the :math:`L^2` projection of the local mesh element size."""
24012401
from firedrake.ufl_expr import CellSize
2402+
from firedrake.function import Function
24022403
from firedrake.functionspace import FunctionSpace
24032404
from firedrake.projection import project
2404-
P1 = FunctionSpace(self, "Lagrange", 1)
2405-
return project(CellSize(self), P1)
2405+
2406+
if self.topological_dimension() == 0:
2407+
# On vertex-only meshes we define the cell sizes as 1
2408+
P0 = FunctionSpace(self, "DG", 0)
2409+
return Function(P0).assign(1)
2410+
2411+
if self.ufl_coordinate_element().degree() > 1:
2412+
# We need a P1 mesh, as CellSize is not defined on high-order meshes
2413+
VectorP1 = self.coordinates.function_space().reconstruct(degree=1)
2414+
mesh = Mesh(Function(VectorP1).interpolate(self.coordinates))
2415+
else:
2416+
mesh = self
2417+
2418+
P1 = FunctionSpace(mesh, "Lagrange", 1)
2419+
h = project(CellSize(mesh), P1)
2420+
2421+
if P1.mesh() is not self:
2422+
# Transfer the Function on the P1 mesh into the original mesh
2423+
h = Function(P1.reconstruct(mesh=self), val=h.dat)
2424+
return h
24062425

24072426
def clear_cell_sizes(self):
24082427
"""Reset the :attr:`cell_sizes` field on this mesh geometry.

tests/firedrake/regression/test_interpolate_zany.py

+12-1
Original file line numberDiff line numberDiff line change
@@ -51,7 +51,7 @@ def expect(V, which):
5151
return a + b
5252

5353

54-
def test_interpolate(V, mesh, which, expect, tolerance):
54+
def test_interpolate_zany_into_cg(V, mesh, which, expect, tolerance):
5555
degree = V.ufl_element().degree()
5656
Vcg = FunctionSpace(mesh, "P", degree)
5757

@@ -71,3 +71,14 @@ def test_interpolate(V, mesh, which, expect, tolerance):
7171
g.interpolate(a + b)
7272

7373
assert numpy.allclose(norm(g - expect), 0, atol=tolerance)
74+
75+
76+
def test_high_order_mesh_cell_sizes():
77+
msh1 = UnitSquareMesh(2, 2)
78+
h1 = msh1.cell_sizes
79+
80+
P2 = msh1.coordinates.function_space().reconstruct(degree=2)
81+
msh2 = Mesh(Function(P2).interpolate(msh1.coordinates))
82+
h2 = msh2.cell_sizes
83+
84+
assert numpy.allclose(h1.dat.data, h2.dat.data)

0 commit comments

Comments
 (0)