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
|
# Copyright (C) 2021 Matthew W. Scroggs and Jack Hale
#
# This file is part of DOLFINx (https://www.fenicsproject.org)
#
# SPDX-License-Identifier: LGPL-3.0-or-later
"""Test that matrices are symmetric."""
from mpi4py import MPI
import numpy as np
import pytest
import basix
import dolfinx
import ufl
from basix.ufl import element, mixed_element
from dolfinx.fem import form, functionspace
from dolfinx.mesh import CellType, create_unit_cube, create_unit_square
from ufl import grad, inner
def check_symmetry(A, tol):
Ad = A.to_dense()
assert np.allclose(Ad, Ad.T, atol=tol)
def run_symmetry_test(cell_type, e, form_f):
dtype = np.float64
if cell_type == CellType.triangle or cell_type == CellType.quadrilateral:
mesh = create_unit_square(MPI.COMM_WORLD, 2, 2, cell_type, dtype=dtype)
else:
mesh = create_unit_cube(MPI.COMM_WORLD, 2, 2, 2, cell_type, dtype=dtype)
space = functionspace(mesh, e)
u = ufl.TrialFunction(space)
v = ufl.TestFunction(space)
f = form(form_f(u, v), dtype=dtype)
A = dolfinx.fem.assemble_matrix(f)
A.scatter_reverse()
tol = np.sqrt(np.finfo(dtype).eps)
check_symmetry(A, tol)
parametrize_elements = pytest.mark.parametrize(
"cell_type, family",
[
(CellType.triangle, "Lagrange"),
(CellType.triangle, "N1curl"),
(CellType.triangle, "RT"),
(CellType.triangle, "Regge"),
(CellType.quadrilateral, "Lagrange"),
(CellType.quadrilateral, "RTCE"),
(CellType.quadrilateral, "RTCF"),
(CellType.tetrahedron, "Lagrange"),
(CellType.tetrahedron, "N1curl"),
(CellType.tetrahedron, "RT"),
(CellType.tetrahedron, "Regge"),
(CellType.hexahedron, "Lagrange"),
(CellType.hexahedron, "NCE"),
(CellType.hexahedron, "NCF"),
],
)
parametrize_lagrange_elements = pytest.mark.parametrize(
"cell_type, family",
[
(CellType.triangle, "Lagrange"),
(CellType.quadrilateral, "Lagrange"),
(CellType.tetrahedron, "Lagrange"),
(CellType.hexahedron, "Lagrange"),
],
)
@pytest.mark.skip_in_parallel
@parametrize_elements
@pytest.mark.parametrize("order", range(1, 2))
def test_mass_matrix_dx(cell_type, family, order):
run_symmetry_test(cell_type, (family, order), lambda u, v: inner(u, v) * ufl.dx)
@pytest.mark.skip_in_parallel
@parametrize_lagrange_elements
@pytest.mark.parametrize("order", range(1, 2))
def test_stiffness_matrix_dx(cell_type, family, order):
run_symmetry_test(cell_type, (family, order), lambda u, v: inner(grad(u), grad(v)) * ufl.dx)
@pytest.mark.skip_in_parallel
@parametrize_elements
@pytest.mark.parametrize("order", range(1, 2))
def test_mass_matrix_ds(cell_type, family, order):
run_symmetry_test(cell_type, (family, order), lambda u, v: inner(u, v) * ufl.ds)
@pytest.mark.skip_in_parallel
@parametrize_lagrange_elements
@pytest.mark.parametrize("order", range(1, 2))
def test_stiffness_matrix_ds(cell_type, family, order):
run_symmetry_test(cell_type, (family, order), lambda u, v: inner(grad(u), grad(v)) * ufl.ds)
@pytest.mark.skip_in_parallel
@parametrize_elements
@pytest.mark.parametrize("order", range(1, 2))
@pytest.mark.parametrize("sign", ["+", "-"])
def test_mass_matrix_dS(cell_type, family, order, sign):
run_symmetry_test(cell_type, (family, order), lambda u, v: inner(u, v)(sign) * ufl.dS)
@pytest.mark.skip_in_parallel
@parametrize_lagrange_elements
@pytest.mark.parametrize("order", range(1, 2))
@pytest.mark.parametrize("sign", ["+", "-"])
def test_stiffness_matrix_dS(cell_type, family, order, sign):
run_symmetry_test(
cell_type, (family, order), lambda u, v: inner(grad(u), grad(v))(sign) * ufl.dS
)
@pytest.mark.skip_in_parallel
@pytest.mark.parametrize(
"cell_type",
[CellType.triangle, CellType.quadrilateral, CellType.tetrahedron, CellType.hexahedron],
)
@pytest.mark.parametrize("sign", ["+", "-"])
@pytest.mark.parametrize("order", range(1, 2))
@pytest.mark.parametrize("dtype", [np.float32, np.float64])
def test_mixed_element_form(cell_type, sign, order, dtype):
if cell_type == CellType.triangle or cell_type == CellType.quadrilateral:
mesh = create_unit_square(MPI.COMM_WORLD, 2, 2, cell_type, dtype=dtype)
else:
mesh = create_unit_cube(MPI.COMM_WORLD, 2, 2, 2, cell_type, dtype=dtype)
U_el = mixed_element(
[
element(basix.ElementFamily.P, cell_type.name, order, dtype=dtype),
element(basix.ElementFamily.N1E, cell_type.name, order, dtype=dtype),
]
)
U = functionspace(mesh, U_el)
u, p = ufl.TrialFunctions(U)
v, q = ufl.TestFunctions(U)
f = form(inner(u, v) * ufl.dx + inner(p, q)(sign) * ufl.dS, dtype=dtype)
A = dolfinx.fem.assemble_matrix(f)
A.scatter_reverse()
tol = np.sqrt(np.finfo(dtype).eps)
check_symmetry(A, tol)
@pytest.mark.skip_in_parallel
@pytest.mark.parametrize("cell_type", [CellType.triangle, CellType.quadrilateral])
@pytest.mark.parametrize("sign", ["+", "-"])
@pytest.mark.parametrize("order", range(1, 2))
@pytest.mark.parametrize("dtype", [np.float32, np.float64])
def test_mixed_element_vector_element_form(cell_type, sign, order, dtype):
if cell_type == CellType.triangle or cell_type == CellType.quadrilateral:
mesh = create_unit_square(MPI.COMM_WORLD, 2, 2, cell_type, dtype=dtype)
else:
mesh = create_unit_cube(MPI.COMM_WORLD, 2, 2, 2, cell_type, dtype=dtype)
U_el = mixed_element(
[
element(
basix.ElementFamily.P,
cell_type.name,
order,
shape=(mesh.geometry.dim,),
dtype=dtype,
),
element(basix.ElementFamily.N1E, cell_type.name, order, dtype=dtype),
]
)
U = functionspace(mesh, U_el)
u, p = ufl.TrialFunctions(U)
v, q = ufl.TestFunctions(U)
f = form(inner(u, v) * ufl.dx + inner(p, q)(sign) * ufl.dS, dtype=dtype)
A = dolfinx.fem.assemble_matrix(f)
A.scatter_reverse()
tol = np.sqrt(np.finfo(dtype).eps)
check_symmetry(A, tol)
|