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
|
import math
import pytest
from numpy import ndindex, reshape
from utils import LagrangeElement
from ufl import (
Coefficient,
FunctionSpace,
Mesh,
TestFunction,
TrialFunction,
VectorConstant,
acos,
as_tensor,
as_ufl,
asin,
atan,
cos,
cosh,
dx,
exp,
i,
j,
ln,
max_value,
min_value,
outer,
sin,
sinh,
tan,
tanh,
triangle,
)
from ufl.algorithms import compute_form_data
from ufl.constantvalue import Zero
from ufl.core.multiindex import FixedIndex, Index, MultiIndex, indices
from ufl.indexed import Indexed
from ufl.tensors import ComponentTensor, ListTensor
def xtest_zero_times_argument(self):
# FIXME: Allow zero forms
element = LagrangeElement(triangle, 1)
domain = Mesh(LagrangeElement(triangle, 1, (2,)))
space = FunctionSpace(domain, element)
v = TestFunction(space)
u = TrialFunction(space)
L = 0 * v * dx
a = 0 * (u * v) * dx
b = (0 * u) * v * dx
assert len(compute_form_data(L).arguments) == 1
assert len(compute_form_data(a).arguments) == 2
assert len(compute_form_data(b).arguments) == 2
def test_divisions(self):
element = LagrangeElement(triangle, 1)
domain = Mesh(LagrangeElement(triangle, 1, (2,)))
space = FunctionSpace(domain, element)
f = Coefficient(space)
# Test simplification of division by 1
a = f
b = f / 1
assert a == b
# Test simplification of division by 1.0
a = f
b = f / 1.0
assert a == b
# Test simplification of division by of zero by something
a = 0 / f
b = 0 * f
assert a == b
# Test simplification of division by self (this simplification has been disabled)
# a = f/f
# b = 1
# assert a == b
def test_products(self):
element = LagrangeElement(triangle, 1)
domain = Mesh(LagrangeElement(triangle, 1, (2,)))
space = FunctionSpace(domain, element)
f = Coefficient(space)
g = Coefficient(space)
# Test simplification of literal multiplication
assert f * 0 == as_ufl(0)
assert 0 * f == as_ufl(0)
assert 1 * f == f
assert f * 1 == f
assert as_ufl(2) * as_ufl(3) == as_ufl(6)
assert as_ufl(2.0) * as_ufl(3.0) == as_ufl(6.0)
# Test reordering of operands
assert f * g == g * f
# Test simplification of self-multiplication (this simplification has been disabled)
# assert f*f == f**2
def test_sums(self):
element = LagrangeElement(triangle, 1)
domain = Mesh(LagrangeElement(triangle, 1, (2,)))
space = FunctionSpace(domain, element)
f = Coefficient(space)
g = Coefficient(space)
# Test reordering of operands
assert f + g == g + f
# Test adding zero
assert f + 0 == f
assert 0 + f == f
# Test collapsing of basic sum (this simplification has been disabled)
# assert f + f == 2 * f
# Test reordering of operands and collapsing sum
a = f + g + f # not collapsed, but ordered
b = g + f + f # not collapsed, but ordered
c = (g + f) + f # not collapsed, but ordered
d = f + (f + g) # not collapsed, but ordered
assert a == b
assert a == c
assert a == d
# Test reordering of operands and collapsing sum
a = f + f + g # collapsed
b = g + (f + f) # collapsed
assert a == b
def test_mathfunctions(self):
for a in (0.1, 0.3, 0.9):
assert math.sin(a) == sin(a)
assert math.cos(a) == cos(a)
assert math.tan(a) == tan(a)
assert math.sinh(a) == sinh(a)
assert math.cosh(a) == cosh(a)
assert math.tanh(a) == tanh(a)
assert math.asin(a) == asin(a)
assert math.acos(a) == acos(a)
assert math.atan(a) == atan(a)
assert math.exp(a) == exp(a)
assert math.log(a) == ln(a)
# TODO: Implement automatic simplification of conditionals?
assert a == float(max_value(a, a - 1))
# TODO: Implement automatic simplification of conditionals?
assert a == float(min_value(a, a + 1))
def test_indexing(self):
domain = Mesh(LagrangeElement(triangle, 1, (2,)))
u = VectorConstant(domain)
v = VectorConstant(domain)
A = outer(u, v)
A2 = as_tensor(A[i, j], (i, j))
assert A2 == A
Bij = u[i] * v[j]
Bij2 = as_tensor(Bij, (i, j))[i, j]
as_tensor(Bij, (i, j))
assert Bij2 == Bij
@pytest.mark.parametrize("shape", [(3,), (3, 2)], ids=("vector", "matrix"))
def test_tensor_from_indexed(self, shape):
element = LagrangeElement(triangle, 1, shape)
domain = Mesh(LagrangeElement(triangle, 1, (2,)))
space = FunctionSpace(domain, element)
f = Coefficient(space)
assert as_tensor(reshape([f[i] for i in ndindex(f.ufl_shape)], f.ufl_shape).tolist()) is f
def test_nested_indexed(self):
# Test that a nested Indexed expression simplifies to the existing Indexed object
shape = (2,)
element = LagrangeElement(triangle, 1, shape)
domain = Mesh(LagrangeElement(triangle, 1, (2,)))
space = FunctionSpace(domain, element)
f = Coefficient(space)
comps = tuple(f[i] for i in range(2))
assert all(isinstance(c, Indexed) for c in comps)
expr = as_tensor(list(reversed(comps)))
multiindex = MultiIndex((FixedIndex(0),))
assert Indexed(expr, multiindex) is expr[0]
assert Indexed(expr, multiindex) is comps[1]
def test_repeated_indexing(self):
# Test that an Indexed with repeated indices does not contract indices
shape = (2, 2)
element = LagrangeElement(triangle, 1, shape)
domain = Mesh(LagrangeElement(triangle, 1, (2,)))
space = FunctionSpace(domain, element)
x = Coefficient(space)
C = as_tensor([x, x])
fi = FixedIndex(0)
i = Index()
ii = MultiIndex((fi, i, i))
expr = Indexed(C, ii)
assert i.count() in expr.ufl_free_indices
assert isinstance(expr, Indexed)
B, jj = expr.ufl_operands
assert B is x
assert tuple(jj) == tuple(ii[1:])
def test_untangle_indexed_component_tensor(self):
shape = (2, 2, 2, 2)
element = LagrangeElement(triangle, 1, shape)
domain = Mesh(LagrangeElement(triangle, 1, (2,)))
space = FunctionSpace(domain, element)
C = Coefficient(space)
r = len(shape)
kk = indices(r)
# Untangle as_tensor(C[kk], kk) -> C
B = as_tensor(Indexed(C, MultiIndex(kk)), kk)
assert B is C
# Untangle as_tensor(C[kk], jj)[ii] -> C[ll]
jj = kk[2:]
A = as_tensor(Indexed(C, MultiIndex(kk)), jj)
assert A is not C
ii = kk
expr = Indexed(A, MultiIndex(ii))
assert isinstance(expr, Indexed)
B, ll = expr.ufl_operands
assert B is C
rep = dict(zip(jj, ii))
expected = tuple(rep.get(k, k) for k in kk)
assert tuple(ll) == expected
def test_simplify_indexed(self):
element = LagrangeElement(triangle, 1, (3,))
domain = Mesh(LagrangeElement(triangle, 1, (2,)))
space = FunctionSpace(domain, element)
u = Coefficient(space)
z = Zero(())
i = Index()
j = Index()
# ListTensor
lt = ListTensor(z, z, u[1])
assert Indexed(lt, MultiIndex((FixedIndex(2),))) == u[1]
# ListTensor -- nested
l0 = ListTensor(z, u[1], z)
l1 = ListTensor(z, z, u[2])
l2 = ListTensor(u[0], z, z)
ll = ListTensor(l0, l1, l2)
assert Indexed(ll, MultiIndex((FixedIndex(1), FixedIndex(2)))) == u[2]
assert Indexed(ll, MultiIndex((FixedIndex(2), i))) == l2[i]
# ComponentTensor + ListTensor
c = ComponentTensor(Indexed(ll, MultiIndex((i, j))), MultiIndex((j, i)))
assert Indexed(c, MultiIndex((FixedIndex(1), FixedIndex(2)))) == l2[1]
|