File: test_forms.py

package info (click to toggle)
fenics-dolfinx 1%3A0.10.0.post5-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 5,952 kB
  • sloc: cpp: 36,535; python: 25,391; makefile: 223; sh: 174; xml: 55
file content (133 lines) | stat: -rw-r--r-- 4,226 bytes parent folder | download
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
"""Tests for DOLFINx integration of various form operations"""

# Copyright (C) 2021 Garth N. Wells
#
# This file is part of DOLFINx (https://www.fenicsproject.org)
#
# SPDX-License-Identifier:    LGPL-3.0-or-later

from mpi4py import MPI

import numpy as np
import pytest

import basix
import basix.ufl
import dolfinx
from dolfinx.fem import IntegralType, extract_function_spaces, form, functionspace
from dolfinx.fem.forms import form_cpp_class
from dolfinx.mesh import create_unit_square
from ufl import Measure, SpatialCoordinate, TestFunction, TrialFunction, dx, inner


def test_extract_forms():
    """Test extraction on unique function spaces for rows and columns of
    a block system"""
    mesh = create_unit_square(MPI.COMM_WORLD, 32, 31)
    V0 = functionspace(mesh, ("Lagrange", 1))
    V1 = functionspace(mesh, ("Lagrange", 2))
    V2 = V0.clone()
    V3 = V1.clone()

    v0, u0 = TestFunction(V0), TrialFunction(V0)
    v1, u1 = TestFunction(V1), TrialFunction(V1)
    v2, u2 = TestFunction(V2), TrialFunction(V2)
    v3, u3 = TestFunction(V3), TrialFunction(V3)

    a = form([[inner(u0, v0) * dx, inner(u1, v1) * dx], [inner(u2, v2) * dx, inner(u3, v3) * dx]])
    with pytest.raises(AssertionError):
        extract_function_spaces(a, 0)
    with pytest.raises(AssertionError):
        extract_function_spaces(a, 1)

    a = form([[inner(u0, v0) * dx, inner(u2, v1) * dx], [inner(u0, v2) * dx, inner(u2, v2) * dx]])
    with pytest.raises(AssertionError):
        extract_function_spaces(a, 0)
    Vc = extract_function_spaces(a, 1)
    assert Vc[0] is V0._cpp_object
    assert Vc[1] is V2._cpp_object

    a = form([[inner(u0, v0) * dx, inner(u1, v0) * dx], [inner(u2, v1) * dx, inner(u3, v1) * dx]])
    Vr = extract_function_spaces(a, 0)
    assert Vr[0] is V0._cpp_object
    assert Vr[1] is V1._cpp_object
    with pytest.raises(AssertionError):
        extract_function_spaces(a, 1)


def test_incorrect_element():
    """Test that an error is raised if an incorrect element is used."""
    mesh = create_unit_square(MPI.COMM_WORLD, 32, 31)
    element = basix.ufl.element(
        "Lagrange",
        "triangle",
        4,
        lagrange_variant=basix.LagrangeVariant.gll_warped,
        dtype=dolfinx.default_real_type,
    )
    incorrect_element = basix.ufl.element(
        "Lagrange",
        "triangle",
        4,
        lagrange_variant=basix.LagrangeVariant.equispaced,
        dtype=dolfinx.default_real_type,
    )

    space = functionspace(mesh, element)
    incorrect_space = functionspace(mesh, incorrect_element)

    u = TrialFunction(space)
    v = TestFunction(space)

    a = inner(u, v) * dx

    dtype = dolfinx.default_scalar_type
    ftype = form_cpp_class(dtype)

    ufcx_form, module, code = dolfinx.jit.ffcx_jit(
        mesh.comm, a, form_compiler_options={"scalar_type": dtype}
    )

    f = ftype(
        [module.ffi.cast("uintptr_t", module.ffi.addressof(ufcx_form))],
        [space._cpp_object, space._cpp_object],
        [],
        [],
        {IntegralType.cell: []},
        [],
        mesh._cpp_object,
    )
    dolfinx.fem.Form(f, ufcx_form, code)

    with pytest.raises(RuntimeError):
        f = ftype(
            [module.ffi.cast("uintptr_t", module.ffi.addressof(ufcx_form))],
            [incorrect_space._cpp_object, incorrect_space._cpp_object],
            [],
            [],
            {IntegralType.cell: []},
            [],
            mesh._cpp_object,
        )
        dolfinx.fem.Form(f, ufcx_form, code)


def test_multiple_measures_one_subdomain_data():
    comm = MPI.COMM_WORLD
    mesh = dolfinx.mesh.create_unit_interval(comm, 10)
    x = SpatialCoordinate(mesh)
    num_cells_local = mesh.topology.index_map(mesh.topology.dim).size_local
    ct = dolfinx.mesh.meshtags(
        mesh,
        mesh.topology.dim,
        np.arange(num_cells_local, dtype=np.int32),
        np.arange(num_cells_local, dtype=np.int32),
    )

    dx = Measure("dx", domain=mesh, subdomain_data=ct)
    dx_stand = Measure("dx", domain=mesh)

    J = dolfinx.fem.form(x[0] ** 2 * dx + x[0] * dx_stand)
    J_local = dolfinx.fem.assemble_scalar(J)
    J_global = comm.allreduce(J_local, op=MPI.SUM)
    assert np.isclose(J_global, 1 / 3 + 1 / 2)