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
|
# Copyright (C) 20222 Garth N. Wells
#
# This file is part of DOLFINx (https://www.fenicsproject.org)
#
# SPDX-License-Identifier: LGPL-3.0-or-later
"""Unit tests for MatrixCSR"""
from mpi4py import MPI
import numpy as np
import pytest
from dolfinx import cpp as _cpp
from dolfinx import la
from dolfinx.fem import functionspace
from dolfinx.mesh import create_unit_square
def test_create_matrix_csr():
"""Test creation of CSR matrix with specified types"""
mesh = create_unit_square(MPI.COMM_WORLD, 10, 11)
V = functionspace(mesh, ("Lagrange", 1))
map = V.dofmap.index_map
bs = V.dofmap.index_map_bs
pattern = _cpp.la.SparsityPattern(mesh.comm, [map, map], [bs, bs])
rows = np.arange(0, bs * map.size_local)
cols = np.arange(0, bs * map.size_local)
pattern.insert(rows, cols)
pattern.finalize()
A = la.matrix_csr(pattern)
assert A.data.dtype == np.float64
A = la.matrix_csr(pattern, dtype=np.float64)
assert A.data.dtype == np.float64
A = la.matrix_csr(pattern, dtype=np.float32)
assert A.data.dtype == np.float32
A = la.matrix_csr(pattern, dtype=np.complex128)
assert A.data.dtype == np.complex128
dense = A.to_dense()
assert np.allclose(dense, np.zeros(dense.shape, dtype=np.complex128))
@pytest.mark.parametrize(
"dtype",
[
np.float32,
np.float64,
np.complex64,
np.complex128,
],
)
def test_matvec(dtype):
mesh = create_unit_square(MPI.COMM_WORLD, 3, 3)
imap = mesh.topology.index_map(0)
sp = _cpp.la.SparsityPattern(mesh.comm, [imap, imap], [1, 1])
rows = np.arange(0, imap.size_local)
cols = np.arange(0, imap.size_local)
sp.insert(rows, cols)
sp.finalize()
# Identity
A = la.matrix_csr(sp, dtype=dtype)
for i in range(imap.size_local):
A.add([2.0], [i], [i])
A.scatter_reverse()
b = la.vector(imap, dtype=dtype)
u = la.vector(imap, dtype=dtype)
b.array[:] = 1.0
A.mult(b, u)
u.scatter_forward()
assert np.allclose(u.array[: imap.size_local], 2.0)
@pytest.mark.parametrize(
"dtype",
[
np.float32,
np.float64,
np.complex64,
np.complex128,
np.int8,
np.int32,
np.int64,
],
)
def test_create_vector(dtype):
"""Test creation of a distributed vector"""
mesh = create_unit_square(MPI.COMM_WORLD, 5, 5)
im = mesh.topology.index_map(0)
for bs in range(1, 4):
x = la.vector(im, bs=bs, dtype=dtype)
assert x.array.dtype == dtype
assert x.array.size == bs * (im.size_local + im.num_ghosts)
def xfail_norm_of_integral_type_vector(dtype):
return pytest.param(
dtype,
marks=pytest.mark.xfail(
reason="Norm of vector of integers not implemented", strict=True, raises=TypeError
),
)
@pytest.mark.parametrize(
"dtype",
[
np.float32,
np.float64,
np.complex64,
np.complex128,
xfail_norm_of_integral_type_vector(np.int8),
xfail_norm_of_integral_type_vector(np.int32),
xfail_norm_of_integral_type_vector(np.int64),
],
)
@pytest.mark.parametrize(
"norm_type",
[
la.Norm.l1,
la.Norm.l2,
la.Norm.linf,
pytest.param(
la.Norm.frobenius,
marks=pytest.mark.xfail(reason="Norm type not supported for vector", strict=True),
),
],
)
def test_vector_norm(dtype, norm_type):
mesh = create_unit_square(MPI.COMM_WORLD, 5, 5)
im = mesh.topology.index_map(0)
x = la.vector(im, dtype=dtype)
x.array[:] = 0.0
normed_value = la.norm(x, norm_type)
assert np.isclose(normed_value, 0.0)
|