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
|
import itertools
from mpi4py import MPI
import basix
import basix.ufl
import dolfinx
import numpy as np
import pytest
dtypes = [np.float64, np.float32] # Mesh geometry dtypes
write_comm = [MPI.COMM_SELF, MPI.COMM_WORLD] # Communicators for creating mesh
simplex_two_dim = itertools.product(dtypes, [dolfinx.mesh.CellType.triangle], write_comm)
simplex_three_dim = itertools.product(dtypes, [dolfinx.mesh.CellType.tetrahedron], write_comm)
non_simplex_two_dim = itertools.product(dtypes, [dolfinx.mesh.CellType.quadrilateral], write_comm)
non_simplex_three_dim = itertools.product(dtypes, [dolfinx.mesh.CellType.hexahedron], write_comm)
@pytest.fixture(params=simplex_two_dim, scope="module")
def simplex_mesh_2D(request):
dtype, cell_type, write_comm = request.param
mesh = dolfinx.mesh.create_unit_square(write_comm, 10, 10, cell_type=cell_type, dtype=dtype)
return mesh
@pytest.fixture(params=simplex_three_dim, scope="module")
def simplex_mesh_3D(request):
dtype, cell_type, write_comm = request.param
mesh = dolfinx.mesh.create_unit_cube(write_comm, 5, 5, 5, cell_type=cell_type, dtype=dtype)
return mesh
@pytest.fixture(params=non_simplex_two_dim, scope="module")
def non_simplex_mesh_2D(request):
dtype, cell_type, write_comm = request.param
mesh = dolfinx.mesh.create_unit_square(write_comm, 10, 10, cell_type=cell_type, dtype=dtype)
return mesh
@pytest.fixture(params=non_simplex_three_dim, scope="module")
def non_simplex_mesh_3D(request):
dtype, cell_type, write_comm = request.param
mesh = dolfinx.mesh.create_unit_cube(write_comm, 5, 5, 5, cell_type=cell_type, dtype=dtype)
return mesh
@pytest.mark.parametrize("is_complex", [True, False])
@pytest.mark.parametrize("family", ["N1curl", "RT"])
@pytest.mark.parametrize("degree", [1, 4])
@pytest.mark.parametrize("read_comm", [MPI.COMM_SELF, MPI.COMM_WORLD])
def test_read_write_2D(
read_comm, family, degree, is_complex, simplex_mesh_2D, get_dtype, write_function, read_function
):
mesh = simplex_mesh_2D
f_dtype = get_dtype(mesh.geometry.x.dtype, is_complex)
el = basix.ufl.element(family, mesh.ufl_cell().cellname(), degree, dtype=mesh.geometry.x.dtype)
def f(x):
values = np.empty((2, x.shape[1]), dtype=f_dtype)
values[0] = np.full(x.shape[1], np.pi) + x[0]
values[1] = x[1]
if is_complex:
values[0] += 2j * x[1]
values[1] += 2j * x[0]
return values
fname = write_function(mesh, el, f, f_dtype)
MPI.COMM_WORLD.Barrier()
read_function(read_comm, el, f, fname, f_dtype)
@pytest.mark.parametrize("is_complex", [True, False])
@pytest.mark.parametrize("family", ["N1curl", "RT"])
@pytest.mark.parametrize("degree", [1, 4])
@pytest.mark.parametrize("read_comm", [MPI.COMM_SELF, MPI.COMM_WORLD])
def test_read_write_3D(
read_comm, family, degree, is_complex, simplex_mesh_3D, get_dtype, write_function, read_function
):
mesh = simplex_mesh_3D
f_dtype = get_dtype(mesh.geometry.x.dtype, is_complex)
el = basix.ufl.element(family, mesh.ufl_cell().cellname(), degree, dtype=mesh.geometry.x.dtype)
def f(x):
values = np.empty((3, x.shape[1]), dtype=f_dtype)
values[0] = np.full(x.shape[1], np.pi)
values[1] = x[1] + 2 * x[0]
values[2] = np.cos(x[2])
if is_complex:
values[0] += 2j * x[2]
values[1] += 2j * np.cos(x[2])
return values
fname = write_function(mesh, el, f, dtype=f_dtype)
MPI.COMM_WORLD.Barrier()
read_function(read_comm, el, f, fname, dtype=f_dtype)
@pytest.mark.parametrize("is_complex", [True, False])
@pytest.mark.parametrize("family", ["RTCF"])
@pytest.mark.parametrize("degree", [1, 2, 3])
@pytest.mark.parametrize("read_comm", [MPI.COMM_SELF, MPI.COMM_WORLD])
def test_read_write_2D_quad(
read_comm,
family,
degree,
is_complex,
non_simplex_mesh_2D,
get_dtype,
write_function,
read_function,
):
mesh = non_simplex_mesh_2D
f_dtype = get_dtype(mesh.geometry.x.dtype, is_complex)
el = basix.ufl.element(family, mesh.ufl_cell().cellname(), degree, dtype=mesh.geometry.x.dtype)
def f(x):
values = np.empty((2, x.shape[1]), dtype=f_dtype)
values[0] = np.full(x.shape[1], np.pi)
values[1] = x[1] + 2 * x[0]
if is_complex:
values[0] += 2j * x[2]
values[1] += 2j * np.cos(x[2])
return values
hash = write_function(mesh, el, f, f_dtype)
MPI.COMM_WORLD.Barrier()
read_function(read_comm, el, f, hash, f_dtype)
@pytest.mark.parametrize("is_complex", [True, False])
@pytest.mark.parametrize("family", ["NCF"])
@pytest.mark.parametrize("degree", [1, 4])
@pytest.mark.parametrize("read_comm", [MPI.COMM_SELF, MPI.COMM_WORLD])
def test_read_write_hex(
read_comm,
family,
degree,
is_complex,
non_simplex_mesh_3D,
get_dtype,
write_function,
read_function,
):
mesh = non_simplex_mesh_3D
f_dtype = get_dtype(mesh.geometry.x.dtype, is_complex)
el = basix.ufl.element(family, mesh.ufl_cell().cellname(), degree, dtype=mesh.geometry.x.dtype)
def f(x):
values = np.empty((3, x.shape[1]), dtype=f_dtype)
values[0] = np.full(x.shape[1], np.pi) + x[0]
values[1] = np.cos(x[2])
values[2] = x[0]
if is_complex:
values[0] += 2j * x[2]
values[2] -= 1j * x[1]
return values
hash = write_function(mesh, el, f, dtype=f_dtype)
MPI.COMM_WORLD.Barrier()
read_function(read_comm, el, f, hash, dtype=f_dtype)
@pytest.mark.parametrize("is_complex", [True, False])
@pytest.mark.parametrize("family", ["RTCF"])
@pytest.mark.parametrize("degree", [1, 2, 3])
@pytest.mark.parametrize("read_comm", [MPI.COMM_SELF, MPI.COMM_WORLD])
def test_read_write_multiple(
read_comm,
family,
degree,
is_complex,
non_simplex_mesh_2D,
get_dtype,
write_function,
read_function,
):
mesh = non_simplex_mesh_2D
f_dtype = get_dtype(mesh.geometry.x.dtype, is_complex)
el = basix.ufl.element(family, mesh.ufl_cell().cellname(), degree, dtype=mesh.geometry.x.dtype)
def f(x):
values = np.empty((2, x.shape[1]), dtype=f_dtype)
values[0] = np.full(x.shape[1], np.pi)
values[1] = x[1] + 2 * x[0]
if is_complex:
values[0] -= 2j * x[2]
values[1] += 2j * np.cos(x[2])
return values
def g(x):
values = np.empty((2, x.shape[1]), dtype=f_dtype)
values[0] = 2 * x[1]
values[1] = 3 * x[0]
if is_complex:
values[0] += 3j * x[0]
values[1] += 2j * x[0] * x[1]
return values
hash_f = write_function(mesh, el, f, dtype=f_dtype, name="f", append=False)
hash_g = write_function(mesh, el, g, dtype=f_dtype, name="g", append=True)
assert hash_f == hash_g
MPI.COMM_WORLD.Barrier()
read_function(read_comm, el, f, hash_f, dtype=f_dtype, name="f")
read_function(read_comm, el, g, hash_g, dtype=f_dtype, name="g")
|