File: test_xdmf_function.py

package info (click to toggle)
fenics-dolfinx 1%3A0.9.0-10
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 5,376 kB
  • sloc: cpp: 33,701; python: 22,338; makefile: 230; sh: 170; xml: 55
file content (409 lines) | stat: -rw-r--r-- 15,611 bytes parent folder | download | duplicates (2)
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
# Copyright (C) 2012-2019 Garth N. Wells
#
# 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

import basix
from dolfinx import default_real_type
from dolfinx.fem import Function, functionspace
from dolfinx.io import XDMFFile
from dolfinx.mesh import CellType, create_unit_cube, create_unit_interval, create_unit_square

# Supported XDMF file encoding
if MPI.COMM_WORLD.size > 1:
    encodings = [XDMFFile.Encoding.HDF5]
else:
    encodings = [XDMFFile.Encoding.HDF5, XDMFFile.Encoding.ASCII]

celltypes_2D = [CellType.triangle, CellType.quadrilateral]
celltypes_3D = [CellType.tetrahedron, CellType.hexahedron]


def mesh_factory(tdim, n):
    if tdim == 1:
        return create_unit_interval(MPI.COMM_WORLD, n)
    elif tdim == 2:
        return create_unit_square(MPI.COMM_WORLD, n, n)
    elif tdim == 3:
        return create_unit_cube(MPI.COMM_WORLD, n, n, n)


# --- Function


@pytest.mark.parametrize("use_pathlib", [True, False])
@pytest.mark.parametrize("encoding", encodings)
@pytest.mark.parametrize("dtype", [np.double, np.complex128])
def test_save_1d_scalar(tempdir, encoding, dtype, use_pathlib):
    xtype = np.real(dtype(0)).dtype
    filename2 = Path(tempdir).joinpath("u1_.xdmf") if use_pathlib else Path(tempdir, "u1_.xdmf")
    mesh = create_unit_interval(MPI.COMM_WORLD, 32, dtype=xtype)
    V = functionspace(mesh, ("Lagrange", 2))
    u = Function(V, dtype=dtype)
    u.x.array[:] = 1.0 + (1j if np.issubdtype(dtype, np.complexfloating) else 0)

    with pytest.raises(RuntimeError):
        with XDMFFile(mesh.comm, filename2, "w", encoding=encoding) as file:
            file.write_mesh(mesh)
            file.write_function(u)

    V1 = functionspace(mesh, ("Lagrange", 1))
    u1 = Function(V1, dtype=dtype)
    u1.interpolate(u)
    with XDMFFile(mesh.comm, filename2, "w", encoding=encoding) as file:
        file.write_mesh(mesh)
        file.write_function(u1)


@pytest.mark.parametrize("cell_type", celltypes_2D)
@pytest.mark.parametrize("encoding", encodings)
@pytest.mark.parametrize("dtype", [np.double, np.complex128])
def test_save_2d_scalar(tempdir, encoding, dtype, cell_type):
    xtype = np.real(dtype(0)).dtype
    filename = Path(tempdir, "u2.xdmf")
    mesh = create_unit_square(MPI.COMM_WORLD, 12, 12, cell_type, dtype=xtype)
    V = functionspace(mesh, ("Lagrange", 2))
    u = Function(V, dtype=dtype)
    u.x.array[:] = 1.0

    V1 = functionspace(mesh, ("Lagrange", 1))
    u1 = Function(V1, dtype=dtype)
    u1.interpolate(u)
    with XDMFFile(mesh.comm, filename, "w", encoding=encoding) as file:
        file.write_mesh(mesh)
        file.write_function(u1)

    # Discontinuous (degree == 0)
    V = functionspace(mesh, ("Discontinuous Lagrange", 0))
    u = Function(V, dtype=dtype)
    with XDMFFile(mesh.comm, filename, "w", encoding=encoding) as file:
        file.write_mesh(mesh)
        file.write_function(u)

    # Discontinuous (degree > 0) should raise exception
    V = functionspace(mesh, ("Discontinuous Lagrange", 1))
    u = Function(V, dtype=dtype)
    with pytest.raises(RuntimeError):
        with XDMFFile(mesh.comm, filename, "w", encoding=encoding) as file:
            file.write_mesh(mesh)
            file.write_function(u)


@pytest.mark.parametrize("cell_type", celltypes_3D)
@pytest.mark.parametrize("encoding", encodings)
@pytest.mark.parametrize("dtype", [np.double, np.complex128])
def test_save_3d_scalar(tempdir, encoding, dtype, cell_type):
    xtype = np.real(dtype(0)).dtype
    filename = Path(tempdir, "u3.xdmf")
    mesh = create_unit_cube(MPI.COMM_WORLD, 4, 3, 4, cell_type, dtype=xtype)
    V = functionspace(mesh, ("Lagrange", 2))
    u = Function(V, dtype=dtype)
    u.x.array[:] = 1.0

    V1 = functionspace(mesh, ("Lagrange", 1))
    u1 = Function(V1, dtype=dtype)
    u1.interpolate(u)
    with XDMFFile(mesh.comm, filename, "w", encoding=encoding) as file:
        file.write_mesh(mesh)
        file.write_function(u1)


@pytest.mark.parametrize("cell_type", celltypes_2D)
@pytest.mark.parametrize("encoding", encodings)
@pytest.mark.parametrize("dtype", [np.double, np.complex128])
def test_save_2d_vector(tempdir, encoding, dtype, cell_type):
    xtype = np.real(dtype(0)).dtype
    filename = Path(tempdir, "u_2dv.xdmf")
    mesh = create_unit_square(MPI.COMM_WORLD, 12, 13, cell_type, dtype=xtype)
    gdim = mesh.geometry.dim

    V = functionspace(mesh, ("Lagrange", 2, (gdim,)))
    u = Function(V, dtype=dtype)
    u.x.array[:] = 1.0 + (1j if np.issubdtype(dtype, np.complexfloating) else 0)

    V1 = functionspace(mesh, ("Lagrange", 1, (gdim,)))
    u1 = Function(V1, dtype=dtype)
    u1.interpolate(u)
    with XDMFFile(mesh.comm, filename, "w", encoding=encoding) as file:
        file.write_mesh(mesh)
        file.write_function(u1)


@pytest.mark.parametrize("cell_type", celltypes_3D)
@pytest.mark.parametrize("encoding", encodings)
@pytest.mark.parametrize("dtype", [np.double, np.complex128])
def test_save_3d_vector(tempdir, encoding, dtype, cell_type):
    xtype = np.real(dtype(0)).dtype
    filename = Path(tempdir, "u_3Dv.xdmf")
    mesh = create_unit_cube(MPI.COMM_WORLD, 2, 2, 2, cell_type, dtype=xtype)
    gdim = mesh.geometry.dim

    u = Function(functionspace(mesh, ("Lagrange", 1, (gdim,))), dtype=dtype)
    u.x.array[:] = 1.0 + (1j if np.issubdtype(dtype, np.complexfloating) else 0)

    V1 = functionspace(mesh, ("Lagrange", 1, (gdim,)))
    u1 = Function(V1, dtype=dtype)
    u1.interpolate(u)
    with XDMFFile(mesh.comm, filename, "w", encoding=encoding) as file:
        file.write_mesh(mesh)
        file.write_function(u1)

    V2 = functionspace(mesh, ("RT", 1))
    u2 = Function(V2, dtype=dtype)
    u2.interpolate(u)
    with pytest.raises(RuntimeError):
        with XDMFFile(mesh.comm, filename, "w", encoding=encoding) as file:
            file.write_mesh(mesh)
            file.write_function(u2)


@pytest.mark.parametrize("cell_type", celltypes_2D)
@pytest.mark.parametrize("encoding", encodings)
@pytest.mark.parametrize("dtype", [np.double, np.complex128])
def test_save_2d_tensor(tempdir, encoding, dtype, cell_type):
    xtype = np.real(dtype(0)).dtype
    filename = Path(tempdir, "tensor.xdmf")
    mesh = create_unit_square(MPI.COMM_WORLD, 16, 16, cell_type, dtype=xtype)
    gdim = mesh.geometry.dim

    u = Function(functionspace(mesh, ("Lagrange", 2, (gdim, gdim))), dtype=dtype)
    u.x.array[:] = 1.0 + (1j if np.issubdtype(dtype, np.complexfloating) else 0)
    u1 = Function(functionspace(mesh, ("Lagrange", 1, (gdim, gdim))), dtype=dtype)
    u1.interpolate(u)
    with XDMFFile(mesh.comm, filename, "w", encoding=encoding) as file:
        file.write_mesh(mesh)
        file.write_function(u1)


@pytest.mark.parametrize("cell_type", celltypes_3D)
@pytest.mark.parametrize("encoding", encodings)
@pytest.mark.parametrize("dtype", [np.float32, np.float64, np.complex64, np.complex128])
def test_save_3d_tensor(tempdir, encoding, dtype, cell_type):
    xtype = np.real(dtype(0)).dtype
    filename = Path(tempdir, "u3t.xdmf")
    mesh = create_unit_cube(MPI.COMM_WORLD, 4, 4, 4, cell_type, dtype=xtype)

    gdim = mesh.geometry.dim
    u = Function(functionspace(mesh, ("Lagrange", 2, (gdim, gdim))), dtype=dtype)
    u.x.array[:] = 1.0 + (1j if np.issubdtype(dtype, np.complexfloating) else 0)

    u1 = Function(functionspace(mesh, ("Lagrange", 1, (gdim, gdim))), dtype=dtype)
    u1.interpolate(u)
    with XDMFFile(mesh.comm, filename, "w", encoding=encoding) as file:
        file.write_mesh(mesh)
        file.write_function(u1)


@pytest.mark.parametrize("cell_type", celltypes_3D)
@pytest.mark.parametrize("encoding", encodings)
@pytest.mark.parametrize("dtype", [np.float32, np.float64, np.complex64, np.complex128])
def test_save_3d_vector_series(tempdir, encoding, dtype, cell_type):
    filename = Path(tempdir, "u_3D.xdmf")
    xtype = np.real(dtype(0)).dtype
    mesh = create_unit_cube(MPI.COMM_WORLD, 2, 2, 2, cell_type, dtype=xtype)
    gdim = mesh.geometry.dim
    u = Function(functionspace(mesh, ("Lagrange", 1, (gdim,))), dtype=dtype)
    with XDMFFile(mesh.comm, filename, "w", encoding=encoding) as file:
        file.write_mesh(mesh)
        u.x.array[:] = 1.0 + (1j if np.issubdtype(dtype, np.complexfloating) else 0)
        file.write_function(u, 0.1)
        u.x.array[:] = 2.0 + (2j if np.issubdtype(dtype, np.complexfloating) else 0)
        file.write_function(u, 0.2)
    with XDMFFile(mesh.comm, filename, "a", encoding=encoding) as file:
        u.x.array[:] = 3.0 + (3j if np.issubdtype(dtype, np.complexfloating) else 0)
        file.write_function(u, 0.3)


def test_higher_order_function(tempdir):
    """Test Function output for higher-order meshes."""
    gmsh = pytest.importorskip("gmsh")
    from dolfinx.io import gmshio

    gmsh.initialize()

    def gmsh_tet_model(order):
        gmsh.option.setNumber("General.Terminal", 0)
        model = gmsh.model()
        comm = MPI.COMM_WORLD
        if comm.rank == 0:
            model.add("Sphere minus box")
            model.setCurrent("Sphere minus box")
            model.occ.addSphere(0, 0, 0, 1)
            model.occ.synchronize()
            volume_entities = [model[1] for model in model.getEntities(3)]
            model.addPhysicalGroup(3, volume_entities, tag=2)
            model.mesh.generate(3)
            gmsh.option.setNumber("General.Terminal", 1)
            model.mesh.setOrder(order)
            gmsh.option.setNumber("General.Terminal", 0)

        msh, _, _ = gmshio.model_to_mesh(model, comm, 0)
        return msh

    def gmsh_hex_model(order):
        model = gmsh.model()
        gmsh.option.setNumber("General.Terminal", 0)
        model.add("Hexahedral mesh")
        model.setCurrent("Hexahedral mesh")
        comm = MPI.COMM_WORLD
        if comm.rank == 0:
            gmsh.option.setNumber("Mesh.RecombinationAlgorithm", 2)
            gmsh.option.setNumber("Mesh.RecombineAll", 2)
            gmsh.option.setNumber("Mesh.CharacteristicLengthFactor", 1)
            circle = model.occ.addDisk(0, 0, 0, 1, 1)
            circle_inner = model.occ.addDisk(0, 0, 0, 0.5, 0.5)
            cut = model.occ.cut([(2, circle)], [(2, circle_inner)])[0]
            extruded_geometry = model.occ.extrude(cut, 0, 0, 0.5, numElements=[5], recombine=True)
            model.occ.synchronize()
            model.addPhysicalGroup(2, [cut[0][1]], tag=1)
            model.setPhysicalName(2, 1, "2D cylinder")
            boundary_entities = model.getEntities(2)
            other_boundary_entities = []
            for entity in boundary_entities:
                if entity != cut[0][1]:
                    other_boundary_entities.append(entity[1])
            model.addPhysicalGroup(2, other_boundary_entities, tag=3)
            model.setPhysicalName(2, 3, "Remaining boundaries")
            model.mesh.generate(3)
            model.mesh.setOrder(order)
            volume_entities = []
            for entity in extruded_geometry:
                if entity[0] == 3:
                    volume_entities.append(entity[1])
            model.addPhysicalGroup(3, volume_entities, tag=1)
            model.setPhysicalName(3, 1, "Mesh volume")

        msh, _, _ = gmshio.model_to_mesh(gmsh.model, comm, 0)
        return msh

    # -- Degree 1 mesh (tet)
    msh = gmsh_tet_model(1)
    gdim = msh.geometry.dim
    assert msh.geometry.cmap.degree == 1

    # Write P1 Function
    u = Function(functionspace(msh, ("Lagrange", 1, (gdim,))))
    filename = Path(tempdir, "u3D_P1.xdmf")
    with XDMFFile(msh.comm, filename, "w") as file:
        file.write_mesh(msh)
        file.write_function(u)

    # Write P2 Function (exception expected)
    u = Function(functionspace(msh, ("Lagrange", 2, (gdim,))))
    filename = Path(tempdir, "u3D_P2.xdmf")
    with pytest.raises(RuntimeError):
        with XDMFFile(msh.comm, filename, "w") as file:
            file.write_mesh(msh)
            file.write_function(u)

    # -- Degree 2 mesh (tet)
    msh = gmsh_tet_model(2)
    gdim = msh.geometry.dim
    assert msh.geometry.cmap.degree == 2

    # Write P1 Function (exception expected)
    u = Function(functionspace(msh, ("Lagrange", 1, (gdim,))))
    with pytest.raises(RuntimeError):
        filename = Path(tempdir, "u3D_P1.xdmf")
        with XDMFFile(msh.comm, filename, "w") as file:
            file.write_mesh(msh)
            file.write_function(u)

    # Write P2 Function
    u = Function(functionspace(msh, ("Lagrange", 2, (gdim,))))
    filename = Path(tempdir, "u3D_P2.xdmf")
    with XDMFFile(msh.comm, filename, "w") as file:
        file.write_mesh(msh)
        file.write_function(u)

    # -- Degree 3 mesh (tet)
    # NOTE: XDMF/ParaView does not support TETRAHEDRON_20
    msh = gmsh_tet_model(3)
    gdim = msh.geometry.dim
    assert msh.geometry.cmap.degree == 3

    # Write P2 Function (exception expected)
    u = Function(functionspace(msh, ("Lagrange", 2, (gdim,))))
    with pytest.raises(RuntimeError):
        filename = Path(tempdir, "u3D_P3.xdmf")
        with XDMFFile(msh.comm, filename, "w") as file:
            file.write_mesh(msh)
            file.write_function(u)

    # Write P3 GLL Function (exception expected)
    ufl_e = basix.ufl.element(
        basix.ElementFamily.P,
        basix.CellType.tetrahedron,
        3,
        basix.LagrangeVariant.gll_warped,
        dtype=default_real_type,
    )
    u = Function(functionspace(msh, ufl_e))
    with pytest.raises(RuntimeError):
        filename = Path(tempdir, "u3D_P3.xdmf")
        with XDMFFile(msh.comm, filename, "w") as file:
            file.write_mesh(msh)
            file.write_function(u)

    # Write P3 equispaced Function
    ufl_e = basix.ufl.element(
        basix.ElementFamily.P,
        basix.CellType.tetrahedron,
        3,
        basix.LagrangeVariant.equispaced,
        dtype=default_real_type,
    )
    u1 = Function(functionspace(msh, ufl_e))
    u1.interpolate(u)
    filename = Path(tempdir, "u3D_P3.xdmf")
    with XDMFFile(msh.comm, filename, "w") as file:
        file.write_mesh(msh)
        file.write_function(u1)

    # --  Degree 2 mesh (hex)
    msh = gmsh_hex_model(2)
    gdim = msh.geometry.dim
    assert msh.geometry.cmap.degree == 2

    # Write Q1 Function (exception expected)
    u = Function(functionspace(msh, ("Lagrange", 1, (gdim,))))
    with pytest.raises(RuntimeError):
        filename = Path(tempdir, "u3D_Q1.xdmf")
        with XDMFFile(msh.comm, filename, "w") as file:
            file.write_mesh(msh)
            file.write_function(u)

    # Write Q2 Function
    u = Function(functionspace(msh, ("Lagrange", 2, (gdim,))))
    filename = Path(tempdir, "u3D_Q2.xdmf")
    with XDMFFile(msh.comm, filename, "w") as file:
        file.write_mesh(msh)
        file.write_function(u)

    # TODO: Higher-order gmsh hex meshes not yet supported by DOLFINx
    #
    # # Degree 3 mesh (hex)
    # msh = gmsh_hex_model(3)
    # assert msh.geometry.cmap.degree == 3
    # gdim = msh.geometry.dim
    # u = Function(functionspace(msh, ("Lagrange", 1, (gdim,))))
    # with pytest.raises(RuntimeError):
    #     filename = Path(tempdir, "u3D_Q1.xdmf")
    #     with XDMFFile(msh.comm, filename, "w") as file:
    #         file.write_mesh(msh)
    #         file.write_function(u)
    # u = Function(functionspace(msh, ("Lagrange", 2, (gdim,))))
    # filename = Path(tempdir, "u3D_Q2.xdmf")
    # with XDMFFile(msh.comm, filename, "w") as file:
    #     file.write_mesh(msh)
    #     file.write_function(u)

    gmsh.finalize()