1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212
|
# Copyright (C) 2023 Jorgen Dokken and Matthew Scroggs
#
# This file is part of DOLFINx (https://www.fenicsproject.org)
#
# SPDX-License-Identifier: LGPL-3.0-or-later
from mpi4py import MPI
import numpy as np
import pytest
import basix.ufl
import dolfinx
import ufl
@pytest.mark.parametrize("degree", range(1, 4))
def test_default(degree):
msh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 10, 10)
CG2_vect = dolfinx.fem.functionspace(msh, ("Lagrange", 1))
Qe = basix.ufl.quadrature_element(msh.topology.cell_name(), degree=degree)
Quad = dolfinx.fem.functionspace(msh, Qe)
u = dolfinx.fem.Function(Quad)
v = ufl.TrialFunction(CG2_vect)
dx_m = ufl.Measure(
"dx", domain=msh, metadata={"quadrature_degree": 1, "quadrature_scheme": "default"}
)
ds = ufl.Measure("ds", domain=msh)
residual = u * v * dx_m
vol = dolfinx.fem.form(residual)
residual = v * ds
surf = dolfinx.fem.form(residual)
residual = u * v * dx_m + v * ds
vol_surf = dolfinx.fem.form(residual)
vol_v = dolfinx.fem.assemble_vector(vol)
sur_v = dolfinx.fem.assemble_vector(surf)
vol_surf = dolfinx.fem.assemble_vector(vol_surf)
assert np.allclose(vol_v.array + sur_v.array, vol_surf.array)
def test_points_and_weights():
msh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 10, 10)
CG2_vect = dolfinx.fem.functionspace(msh, ("Lagrange", 1))
Qe = basix.ufl.quadrature_element(
msh.topology.cell_name(),
value_shape=(),
points=np.array([[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1 / 3, 1 / 3]]),
weights=np.array([0.2, 0.2, 0.2, 0.4]),
)
Quad = dolfinx.fem.functionspace(msh, Qe)
u = dolfinx.fem.Function(Quad)
v = ufl.TrialFunction(CG2_vect)
dx_m = ufl.Measure(
"dx", domain=msh, metadata={"quadrature_degree": 1, "quadrature_scheme": "default"}
)
ds = ufl.Measure("ds", domain=msh)
residual = u * v * dx_m
vol = dolfinx.fem.form(residual)
residual = v * ds
surf = dolfinx.fem.form(residual)
residual = u * v * dx_m + v * ds
vol_surf = dolfinx.fem.form(residual)
vol_v = dolfinx.fem.assemble_vector(vol)
sur_v = dolfinx.fem.assemble_vector(surf)
vol_surf = dolfinx.fem.assemble_vector(vol_surf)
assert np.allclose(vol_v.array + sur_v.array, vol_surf.array)
@pytest.mark.parametrize("degree", range(1, 5))
def test_interpolation(degree):
msh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 10, 10)
e = basix.ufl.quadrature_element(msh.topology.cell_name(), degree=degree)
space = dolfinx.fem.functionspace(msh, e)
p4 = dolfinx.fem.functionspace(msh, ("Lagrange", 4))
f_p4 = dolfinx.fem.Function(p4)
f_p4.interpolate(lambda x: x[0] ** 4)
f = dolfinx.fem.Function(space)
f.interpolate(lambda x: x[0] ** 4)
diff = dolfinx.fem.form(ufl.inner(f - f_p4, f - f_p4) * ufl.dx)
error = dolfinx.fem.assemble_scalar(diff)
assert np.isclose(error, 0)
@pytest.mark.parametrize("degree", range(1, 5))
def test_interpolation_blocked(degree):
msh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 10, 10)
e = basix.ufl.quadrature_element(msh.topology.cell_name(), value_shape=(2,), degree=degree)
space = dolfinx.fem.functionspace(msh, e)
p4 = dolfinx.fem.functionspace(msh, ("Lagrange", 4, (2,)))
f_p4 = dolfinx.fem.Function(p4)
f_p4.interpolate(lambda x: ([x[1] ** 4, x[0] ** 3]))
f = dolfinx.fem.Function(space)
f.interpolate(lambda x: ([x[1] ** 4, x[0] ** 3]))
diff = dolfinx.fem.form(ufl.inner(f - f_p4, f - f_p4) * ufl.dx)
error = dolfinx.fem.assemble_scalar(diff)
assert np.isclose(error, 0)
def extract_diagonal(mat):
num_rows = mat._cpp_object.index_map(0).size_local
num_cols = mat._cpp_object.index_map(1).size_local
assert num_rows == num_cols, "Matrix must be square"
bs = mat.block_size[0]
diag = np.empty(num_rows * bs, dtype=mat.data.dtype)
for row, (start, end) in enumerate(zip(mat.indptr[:-1], mat.indptr[1:])):
for i in range(start, end):
if mat.indices[i] == row:
for block in range(bs):
diag[bs * row + block] = mat.data[bs**2 * i + (bs + 1) * block]
return diag
@pytest.mark.parametrize("degree", range(1, 4))
@pytest.mark.parametrize("shape", [(), (1,), (2,), (3,), (4,), (2, 2), (3, 3)])
def test_vector_element(shape, degree):
"""
Compare assembly into a vector with quadrature elements with the diagonal of
an assembled mass matrix with the same quadrature element.
"""
msh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 10, 10)
dx_m = ufl.Measure(
"dx",
domain=msh,
metadata={"quadrature_degree": degree, "quadrature_scheme": "default"},
)
Qe = basix.ufl.quadrature_element(
msh.topology.cell_name(), value_shape=shape, scheme="default", degree=degree
)
Quad = dolfinx.fem.functionspace(msh, Qe)
q_ = ufl.TestFunction(Quad)
dq = ufl.TrialFunction(Quad)
one = dolfinx.fem.Function(Quad)
one.x.array[:] = 1.0
mass_L_form = dolfinx.fem.form(ufl.inner(one, q_) * dx_m)
mass_v = dolfinx.fem.assemble_vector(mass_L_form)
mass_v.scatter_reverse(dolfinx.la.InsertMode.add)
mass_v.scatter_forward()
mass_a_form = dolfinx.fem.form(ufl.inner(dq, q_) * dx_m)
mass_A = dolfinx.fem.assemble_matrix(mass_a_form)
mass_A.scatter_reverse()
num_owned_dofs = Quad.dofmap.index_map.size_local * Quad.dofmap.index_map_bs
assert np.allclose(extract_diagonal(mass_A), mass_v.array[:num_owned_dofs])
@pytest.mark.parametrize("degree", range(1, 4))
def test_quadrature_assembly(degree):
"""
Test quadrature element against assembly with spatial coordinate and a fixed quadrature rule
"""
msh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 5, 7)
dx_m = ufl.Measure(
"dx",
domain=msh,
metadata={"quadrature_degree": degree, "quadrature_scheme": "default"},
)
Qe = basix.ufl.quadrature_element(
msh.topology.cell_name(), value_shape=(), scheme="default", degree=degree
)
Quad = dolfinx.fem.functionspace(msh, Qe)
V = dolfinx.fem.functionspace(msh, ("Lagrange", 1))
v = ufl.TestFunction(V)
def f(x):
return 1 + x[0] + x[1] ** degree
q = dolfinx.fem.Function(Quad)
q.interpolate(f)
L = ufl.inner(q, v) * dx_m
b = dolfinx.fem.assemble_vector(dolfinx.fem.form(L))
b.scatter_reverse(dolfinx.la.InsertMode.add)
b.scatter_forward()
x = ufl.SpatialCoordinate(msh)
L_ref = ufl.inner(f(x), v) * dx_m
b_ref = dolfinx.fem.assemble_vector(dolfinx.fem.form(L_ref))
b_ref.scatter_reverse(dolfinx.la.InsertMode.add)
b_ref.scatter_forward()
np.testing.assert_allclose(b.array, b_ref.array)
|