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 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452
|
# Copyright (C) 2011-2021 Garth N. Wells and Jørgen S. Dokken
#
# This file is part of DOLFINx (https://www.fenicsproject.org)
#
# SPDX-License-Identifier: LGPL-3.0-or-later
from pathlib import Path
from mpi4py import MPI
import numpy as np
import pytest
from numpy.testing import assert_array_equal
import ufl
from basix.ufl import element, mixed_element
from dolfinx import default_real_type
from dolfinx.fem import Function, functionspace
from dolfinx.io import VTKFile
from dolfinx.io.utils import cell_perm_vtk # F401
from dolfinx.mesh import (
CellType,
create_mesh,
create_unit_cube,
create_unit_interval,
create_unit_square,
)
from dolfinx.plot import vtk_mesh
cell_types_2D = [CellType.triangle, CellType.quadrilateral]
cell_types_3D = [CellType.tetrahedron, CellType.hexahedron]
def test_save_1d_mesh_subdir(tempdir):
filename = Path(tempdir, "mesh.pvd")
mesh = create_unit_interval(MPI.COMM_WORLD, 32)
with VTKFile(MPI.COMM_WORLD, filename, "w") as vtk:
vtk.write_mesh(mesh)
vtk.write_mesh(mesh, 1)
@pytest.mark.parametrize("cell_type", cell_types_2D)
def test_save_2d_mesh(tempdir, cell_type):
mesh = create_unit_square(MPI.COMM_WORLD, 32, 32, cell_type=cell_type)
filename = Path(tempdir, f"mesh_{cell_type.name}.pvd")
with VTKFile(MPI.COMM_WORLD, filename, "w") as vtk:
vtk.write_mesh(mesh, 0.0)
vtk.write_mesh(mesh, 2.0)
@pytest.mark.parametrize("cell_type", cell_types_3D)
def test_save_3d_mesh(tempdir, cell_type):
mesh = create_unit_cube(MPI.COMM_WORLD, 8, 8, 8, cell_type=cell_type)
filename = Path(tempdir, f"mesh_{cell_type.name}.pvd")
with VTKFile(MPI.COMM_WORLD, filename, "w") as vtk:
vtk.write_mesh(mesh, 0.0)
vtk.write_mesh(mesh, 2.0)
def test_save_1d_scalar(tempdir):
mesh = create_unit_interval(MPI.COMM_WORLD, 32)
u = Function(functionspace(mesh, ("Lagrange", 2)))
u.interpolate(lambda x: x[0])
filename = Path(tempdir, "u.pvd")
with VTKFile(MPI.COMM_WORLD, filename, "w") as vtk:
vtk.write_function(u, 0.0)
@pytest.mark.parametrize("cell_type", cell_types_2D)
def test_save_2d_scalar(tempdir, cell_type):
mesh = create_unit_square(MPI.COMM_WORLD, 16, 16, cell_type=cell_type)
u = Function(functionspace(mesh, ("Lagrange", 2)))
u.x.array[:] = 1.0
filename = Path(tempdir, "u.pvd")
with VTKFile(MPI.COMM_WORLD, filename, "w") as vtk:
vtk.write_function(u, 0.0)
vtk.write_function(u, 1.0)
@pytest.mark.parametrize("cell_type", cell_types_3D)
def test_save_3d_scalar(tempdir, cell_type):
mesh = create_unit_cube(MPI.COMM_WORLD, 8, 8, 8, cell_type=cell_type)
u = Function(functionspace(mesh, ("Lagrange", 2)))
u.x.array[:] = 1.0
filename = Path(tempdir, "u.pvd")
with VTKFile(MPI.COMM_WORLD, filename, "w") as vtk:
vtk.write_function(u, 0.0)
vtk.write_function(u, 1.0)
def test_save_1d_vector(tempdir):
mesh = create_unit_interval(MPI.COMM_WORLD, 32)
def f(x):
vals = np.zeros((2, x.shape[1]))
vals[0] = x[0]
vals[1] = 2 * x[0] * x[0]
return vals
e = element("Lagrange", mesh.basix_cell(), 2, shape=(2,), dtype=default_real_type)
u = Function(functionspace(mesh, e))
u.interpolate(f)
filename = Path(tempdir, "u.pvd")
with VTKFile(MPI.COMM_WORLD, filename, "w") as vtk:
vtk.write_function(u, 0.0)
@pytest.mark.parametrize("cell_type", cell_types_2D)
def test_save_2d_vector(tempdir, cell_type):
mesh = create_unit_square(MPI.COMM_WORLD, 16, 16, cell_type=cell_type)
gdim = mesh.geometry.dim
u = Function(functionspace(mesh, ("Lagrange", 1, (gdim,))))
def f(x):
vals = np.zeros((2, x.shape[1]))
vals[0] = x[0]
vals[1] = 2 * x[0] * x[1]
return vals
u.interpolate(f)
filename = Path(tempdir, "u.pvd")
with VTKFile(MPI.COMM_WORLD, filename, "w") as vtk:
vtk.write_function(u, 0.0)
vtk.write_function(u, 1.0)
@pytest.mark.skip_in_parallel
def test_save_2d_vector_CG2(tempdir):
points = np.array(
[[0, 0], [1, 0], [1, 2], [0, 2], [1 / 2, 0], [1, 1], [1 / 2, 2], [0, 1], [1 / 2, 1]],
dtype=default_real_type,
)
points = np.array(
[[0, 0], [1, 0], [0, 2], [0.5, 1], [0, 1], [0.5, 0], [1, 2], [0.5, 2], [1, 1]],
dtype=default_real_type,
)
cells = np.array([[0, 1, 2, 3, 4, 5], [1, 6, 2, 7, 3, 8]])
domain = ufl.Mesh(element("Lagrange", "triangle", 2, shape=(2,), dtype=default_real_type))
mesh = create_mesh(MPI.COMM_WORLD, cells, points, domain)
gdim = mesh.geometry.dim
u = Function(functionspace(mesh, ("Lagrange", 2, (gdim,))))
u.interpolate(lambda x: np.vstack((x[0], x[1])))
filename = Path(tempdir, "u.pvd")
with VTKFile(mesh.comm, filename, "w") as vtk:
vtk.write_function(u, 0.0)
def test_save_vtk_mixed(tempdir):
mesh = create_unit_cube(MPI.COMM_WORLD, 3, 3, 3)
P2 = element(
"Lagrange", mesh.basix_cell(), 1, shape=(mesh.geometry.dim,), dtype=default_real_type
)
P1 = element("Lagrange", mesh.basix_cell(), 1, dtype=default_real_type)
W = functionspace(mesh, mixed_element([P2, P1]))
V1 = functionspace(mesh, P1)
V2 = functionspace(mesh, P2)
U = Function(W)
U.sub(0).interpolate(lambda x: np.vstack((x[0], 0.2 * x[1], np.zeros_like(x[0]))))
U.sub(1).interpolate(lambda x: 0.5 * x[0])
U1, U2 = Function(V1), Function(V2)
U1.interpolate(U.sub(1))
U2.interpolate(U.sub(0))
U2.name = "u"
U1.name = "p"
filename = Path(tempdir, "u.pvd")
with VTKFile(mesh.comm, filename, "w") as vtk:
vtk.write_function([U2, U1], 0.0)
with VTKFile(mesh.comm, filename, "w") as vtk:
vtk.write_function([U1, U2], 0.0)
Up = U.sub(1)
Up.name = "psub"
with pytest.raises(RuntimeError):
with VTKFile(mesh.comm, filename, "w") as vtk:
vtk.write_function([U2, Up, U1], 0)
with pytest.raises(RuntimeError):
with VTKFile(mesh.comm, filename, "w") as vtk:
vtk.write_function([U.sub(i) for i in range(W.num_sub_spaces)], 0)
@pytest.mark.parametrize("cell_type", cell_types_2D)
def test_save_vector_element(tempdir, cell_type):
mesh = create_unit_square(MPI.COMM_WORLD, 16, 16, cell_type=cell_type)
u = Function(functionspace(mesh, ("RT", 1)))
def f(x):
vals = np.zeros((2, x.shape[1]))
vals[0] = x[0]
vals[1] = 2 * x[0] * x[1]
return vals
u.interpolate(f)
filename = Path(tempdir, "u.pvd")
with pytest.raises(RuntimeError):
with VTKFile(MPI.COMM_WORLD, filename, "w") as vtk:
vtk.write_function(u, 0.0)
vtk.write_function(u, 1.0)
def test_save_vtk_cell_point(tempdir):
"""Test writing cell-wise and point-wise data"""
mesh = create_unit_cube(MPI.COMM_WORLD, 3, 3, 3)
P2 = element("Lagrange", mesh.basix_cell(), 1, shape=(3,), dtype=default_real_type)
P1 = element("Discontinuous Lagrange", mesh.basix_cell(), 0, dtype=default_real_type)
V2, V1 = functionspace(mesh, P2), functionspace(mesh, P1)
U2, U1 = Function(V2), Function(V1)
U2.interpolate(lambda x: np.vstack((x[0], 0.2 * x[1], np.zeros_like(x[0]))))
U1.interpolate(lambda x: 0.5 * x[0])
U2.name = "A"
U1.name = "B"
filename = Path(tempdir, "u.pvd")
with VTKFile(mesh.comm, filename, "w") as vtk:
vtk.write_function([U2, U1], 0.0)
with VTKFile(mesh.comm, filename, "w") as vtk:
vtk.write_function((U1, U2), 0.0)
def test_save_1d_tensor(tempdir):
mesh = create_unit_interval(MPI.COMM_WORLD, 32)
e = element("Lagrange", mesh.basix_cell(), 2, shape=(2, 2), dtype=default_real_type)
u = Function(functionspace(mesh, e))
u.x.array[:] = 1.0
filename = Path(tempdir, "u.pvd")
with VTKFile(mesh.comm, filename, "w") as vtk:
vtk.write_function(u, 0.0)
def test_save_2d_tensor(tempdir):
mesh = create_unit_square(MPI.COMM_WORLD, 16, 16)
gdim = mesh.geometry.dim
u = Function(functionspace(mesh, ("Lagrange", 2, (gdim, gdim))))
u.x.array[:] = 1.0
filename = Path(tempdir, "u.pvd")
with VTKFile(mesh.comm, filename, "w") as vtk:
vtk.write_function(u, 0.0)
u.x.array[:] = 2.0
vtk.write_function(u, 1.0)
def test_save_3d_tensor(tempdir):
mesh = create_unit_cube(MPI.COMM_WORLD, 8, 8, 8)
gdim = mesh.geometry.dim
u = Function(functionspace(mesh, ("Lagrange", 2, (gdim, gdim))))
u.x.array[:] = 1.0
filename = Path(tempdir, "u.pvd")
with VTKFile(mesh.comm, filename, "w") as vtk:
vtk.write_function(u, 0.0)
def test_triangle_perm_vtk():
higher_order_triangle_perm = {
10: np.array([0, 1, 2, 5, 6, 8, 7, 3, 4, 9]),
15: np.array([0, 1, 2, 6, 7, 8, 11, 10, 9, 3, 4, 5, 12, 13, 14]),
21: np.array([0, 1, 2, 7, 8, 9, 10, 14, 13, 12, 11, 3, 4, 5, 6, 15, 18, 16, 20, 19, 17]),
28: np.array(
[
0,
1,
2,
8,
9,
10,
11,
12,
17,
16,
15,
14,
13,
3,
4,
5,
6,
7,
18,
21,
22,
19,
26,
27,
23,
25,
24,
20,
]
),
36: np.array(
[
0,
1,
2,
9,
10,
11,
12,
13,
14,
20,
19,
18,
17,
16,
15,
3,
4,
5,
6,
7,
8,
21,
24,
25,
26,
22,
32,
33,
34,
27,
31,
35,
28,
30,
29,
23,
]
),
45: np.array(
[
0,
1,
2,
10,
11,
12,
13,
14,
15,
16,
23,
22,
21,
20,
19,
18,
17,
3,
4,
5,
6,
7,
8,
9,
24,
27,
28,
29,
30,
25,
38,
39,
42,
40,
31,
37,
44,
43,
32,
36,
41,
33,
35,
34,
26,
]
),
55: np.array(
[
0,
1,
2,
11,
12,
13,
14,
15,
16,
17,
18,
26,
25,
24,
23,
22,
21,
20,
19,
3,
4,
5,
6,
7,
8,
9,
10,
27,
30,
31,
32,
33,
34,
28,
44,
45,
48,
49,
46,
35,
43,
53,
54,
50,
36,
42,
52,
51,
37,
41,
47,
38,
40,
39,
29,
]
),
}
for p_test, v_test in higher_order_triangle_perm.items():
v = cell_perm_vtk(CellType.triangle, p_test)
assert_array_equal(v, v_test)
def test_vtk_mesh():
comm = MPI.COMM_WORLD
mesh = create_unit_square(comm, 2 * comm.size, 2 * comm.size)
V = functionspace(mesh, ("Lagrange", 1))
vtk_mesh(V)
|