File: test_complex.py

package info (click to toggle)
fenics-ufl 2025.2.1-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 2,176 kB
  • sloc: python: 25,267; makefile: 170
file content (220 lines) | stat: -rwxr-xr-x 6,034 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
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
import cmath

import pytest
from utils import LagrangeElement

from ufl import (
    Coefficient,
    FunctionSpace,
    Mesh,
    TestFunction,
    TrialFunction,
    as_tensor,
    as_ufl,
    atan,
    conditional,
    conj,
    cos,
    cosh,
    dot,
    dx,
    exp,
    ge,
    grad,
    gt,
    imag,
    inner,
    le,
    ln,
    lt,
    max_value,
    min_value,
    outer,
    real,
    sin,
    sqrt,
    triangle,
)
from ufl.algebra import Conj, Imag, Real
from ufl.algorithms import estimate_total_polynomial_degree
from ufl.algorithms.apply_algebra_lowering import apply_algebra_lowering
from ufl.algorithms.comparison_checker import ComplexComparisonError, do_comparison_check
from ufl.algorithms.formtransformations import compute_form_adjoint
from ufl.algorithms.remove_complex_nodes import remove_complex_nodes
from ufl.constantvalue import ComplexValue, Zero


def test_conj(self):
    z1 = ComplexValue(1 + 2j)
    z2 = ComplexValue(1 - 2j)

    assert z1 == Conj(z2)
    assert z2 == Conj(z1)


def test_real(self):
    z0 = Zero()
    z1 = as_ufl(1.0)
    z2 = ComplexValue(1j)
    z3 = ComplexValue(1 + 1j)
    assert Real(z1) == z1
    assert Real(z3) == z1
    assert Real(z2) == z0


def test_imag(self):
    z0 = Zero()
    z1 = as_ufl(1.0)
    z2 = as_ufl(1j)
    z3 = ComplexValue(1 + 1j)

    assert Imag(z2) == z1
    assert Imag(z3) == z1
    assert Imag(z1) == z0


def test_compute_form_adjoint(self):
    cell = triangle
    element = LagrangeElement(cell, 1)
    domain = Mesh(LagrangeElement(cell, 1, (2,)))
    space = FunctionSpace(domain, element)

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

    a = inner(grad(u), grad(v)) * dx

    assert compute_form_adjoint(a) == conj(inner(grad(v), grad(u))) * dx


def test_complex_algebra(self):
    z1 = ComplexValue(1j)
    z2 = ComplexValue(1 + 1j)

    # Remember that ufl.algebra functions return ComplexValues, but
    # ufl.mathfunctions return complex Python scalar
    # Any operations with a ComplexValue and a complex Python scalar
    # promote to ComplexValue
    assert z1 * z2 == ComplexValue(-1 + 1j)
    assert z2 / z1 == ComplexValue(1 - 1j)
    assert pow(z2, z1) == ComplexValue((1 + 1j) ** 1j)
    assert sqrt(z2) * as_ufl(1) == ComplexValue(cmath.sqrt(1 + 1j))
    assert (sin(z2) + cosh(z2) - atan(z2)) * z1 == ComplexValue(
        (cmath.sin(1 + 1j) + cmath.cosh(1 + 1j) - cmath.atan(1 + 1j)) * 1j
    )
    assert (abs(z2) - ln(z2)) / exp(z1) == ComplexValue(
        (abs(1 + 1j) - cmath.log(1 + 1j)) / cmath.exp(1j)
    )


def test_automatic_simplification(self):
    cell = triangle
    element = LagrangeElement(cell, 1)
    domain = Mesh(LagrangeElement(cell, 1, (2,)))
    space = FunctionSpace(domain, element)

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

    assert inner(u, v) == u * conj(v)
    assert dot(u, v) == u * v
    assert outer(u, v) == conj(u) * v


def test_apply_algebra_lowering_complex(self):
    cell = triangle
    element = LagrangeElement(cell, 1)
    domain = Mesh(LagrangeElement(cell, 1, (2,)))
    space = FunctionSpace(domain, element)

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

    gv = grad(v)
    gu = grad(u)

    a = dot(gu, gv)
    b = inner(gv, gu)
    c = outer(gu, gv)

    lowered_a = apply_algebra_lowering(a)
    lowered_b = apply_algebra_lowering(b)
    lowered_c = apply_algebra_lowering(c)
    lowered_a_index = lowered_a.index()
    lowered_b_index = lowered_b.index()
    lowered_c_indices = lowered_c.indices()

    assert lowered_a == gu[lowered_a_index] * gv[lowered_a_index]
    assert lowered_b == gv[lowered_b_index] * conj(gu[lowered_b_index])
    assert lowered_c == as_tensor(
        conj(gu[lowered_c_indices[0]]) * gv[lowered_c_indices[1]],
        (lowered_c_indices[0],) + (lowered_c_indices[1],),
    )


def test_remove_complex_nodes(self):
    cell = triangle
    element = LagrangeElement(cell, 1)
    domain = Mesh(LagrangeElement(cell, 1, (2,)))
    space = FunctionSpace(domain, element)

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

    a = conj(v)
    b = real(u)
    c = imag(f)
    d = conj(real(v)) * imag(conj(u))

    assert remove_complex_nodes(a) == v
    assert remove_complex_nodes(b) == u
    with pytest.raises(BaseException):
        remove_complex_nodes(c)
    with pytest.raises(BaseException):
        remove_complex_nodes(d)


def test_comparison_checker(self):
    cell = triangle
    element = LagrangeElement(cell, 1)
    domain = Mesh(LagrangeElement(cell, 1, (2,)))
    space = FunctionSpace(domain, element)

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

    a = conditional(ge(abs(u), imag(v)), u, v)
    b = conditional(le(sqrt(abs(u)), imag(v)), as_ufl(1), as_ufl(1j))
    c = conditional(gt(abs(u), pow(imag(v), 0.5)), sin(u), cos(v))
    d = conditional(lt(as_ufl(-1), as_ufl(1)), u, v)
    e = max_value(as_ufl(0), real(u))
    f = min_value(sin(u), cos(v))
    g = min_value(sin(pow(u, 3)), cos(abs(v)))

    assert do_comparison_check(a) == conditional(ge(real(abs(u)), real(imag(v))), u, v)
    with pytest.raises(ComplexComparisonError):
        b = do_comparison_check(b)
    with pytest.raises(ComplexComparisonError):
        c = do_comparison_check(c)
    assert do_comparison_check(d) == conditional(lt(real(as_ufl(-1)), real(as_ufl(1))), u, v)
    assert do_comparison_check(e) == max_value(real(as_ufl(0)), real(real(u)))
    assert do_comparison_check(f) == min_value(real(sin(u)), real(cos(v)))
    assert do_comparison_check(g) == min_value(real(sin(pow(u, 3))), real(cos(abs(v))))


def test_complex_degree_handling(self):
    cell = triangle
    element = LagrangeElement(cell, 3)
    domain = Mesh(LagrangeElement(cell, 1, (2,)))
    space = FunctionSpace(domain, element)

    v = TestFunction(space)

    a = conj(v)
    b = imag(v)
    c = real(v)

    assert estimate_total_polynomial_degree(a) == 3
    assert estimate_total_polynomial_degree(b) == 3
    assert estimate_total_polynomial_degree(c) == 3