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
|
# Copyright (C) 2007-2024 Anders Logg Garth N. Wells, Marie E. Rognes, Lizao Li, Matthew Scroggs
# This test was originally part of FFCx
# Copyright (C) 2007-2017 Anders Logg and Garth N. Wells
#
# This file is part of FFCx.
#
# FFCx is free software: you can redistribute it and/or modify
# it under the terms of the GNU Lesser General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# FFCx is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public License
# along with FFCx. If not, see <http://www.gnu.org/licenses/>.
#
# Modified by Marie E. Rognes, 2010
# Modified by Lizao Li, 2016
import random
import numpy as np
import pytest
import basix
import basix.ufl
def random_point(cell):
vertices = basix.geometry(basix.CellType[cell])
w = [random.random() for _ in vertices]
return sum(v * i for v, i in zip(vertices, w)) / sum(w)
@pytest.mark.parametrize(
"family, cell, degree, dim",
[
("Lagrange", "triangle", 1, 3),
("Lagrange", "triangle", 2, 6),
("Lagrange", "triangle", 3, 10),
("DG", "triangle", 1, 3),
("DG", "triangle", 2, 6),
("DG", "triangle", 3, 10),
("Lagrange", "quadrilateral", 1, 4),
("Lagrange", "quadrilateral", 2, 9),
("Lagrange", "quadrilateral", 3, 16),
("Regge", "triangle", 0, 3),
("Regge", "triangle", 1, 9),
("Regge", "triangle", 2, 18),
("Regge", "triangle", 3, 30),
("HHJ", "triangle", 0, 3),
("HHJ", "triangle", 1, 9),
("HHJ", "triangle", 2, 18),
("HHJ", "triangle", 3, 30),
],
)
def test_dimension(family, cell, degree, dim):
e = basix.ufl.element(family, cell, degree)
assert e.dim == dim
@pytest.mark.parametrize(
"family, cell, degree, functions",
[
("Lagrange", "interval", 1, [lambda x: 1 - x[0], lambda x: x[0]]),
(
"Lagrange",
"triangle",
1,
[lambda x: 1 - x[0] - x[1], lambda x: x[0], lambda x: x[1]],
),
(
"Lagrange",
"tetrahedron",
1,
[
lambda x: 1 - x[0] - x[1] - x[2],
lambda x: x[0],
lambda x: x[1],
lambda x: x[2],
],
),
(
"Lagrange",
"quadrilateral",
1,
[
lambda x: (1 - x[0]) * (1 - x[1]),
lambda x: x[0] * (1 - x[1]),
lambda x: (1 - x[0]) * x[1],
lambda x: x[0] * x[1],
],
),
(
"Lagrange",
"hexahedron",
1,
[
lambda x: (1 - x[0]) * (1 - x[1]) * (1 - x[2]),
lambda x: x[0] * (1 - x[1]) * (1 - x[2]),
lambda x: (1 - x[0]) * x[1] * (1 - x[2]),
lambda x: x[0] * x[1] * (1 - x[2]),
lambda x: (1 - x[0]) * (1 - x[1]) * x[2],
lambda x: x[0] * (1 - x[1]) * x[2],
lambda x: (1 - x[0]) * x[1] * x[2],
lambda x: x[0] * x[1] * x[2],
],
),
(
"Brezzi-Douglas-Marini",
"triangle",
1,
[
lambda x: (-x[0], -x[1]),
lambda x: (3**0.5 * x[0], -(3**0.5) * x[1]),
lambda x: (x[0] - 1, x[1]),
lambda x: (3**0.5 * (1 - x[0] - 2 * x[1]), 3**0.5 * x[1]),
lambda x: (-x[0], 1 - x[1]),
lambda x: (-(3**0.5) * x[0], 3**0.5 * (2 * x[0] + x[1] - 1)),
],
),
(
"Raviart-Thomas",
"triangle",
1,
[
lambda x: (-x[0], -x[1]),
lambda x: (x[0] - 1, x[1]),
lambda x: (-x[0], 1 - x[1]),
],
),
(
"Raviart-Thomas",
"tetrahedron",
1,
[
lambda x: (2**0.5 * x[0], 2**0.5 * x[1], 2**0.5 * x[2]),
lambda x: (2**0.5 - 2**0.5 * x[0], -(2**0.5) * x[1], -(2**0.5) * x[2]),
lambda x: (2**0.5 * x[0], 2**0.5 * x[1] - 2**0.5, 2**0.5 * x[2]),
lambda x: (-(2**0.5) * x[0], -(2**0.5) * x[1], 2**0.5 - 2**0.5 * x[2]),
],
),
(
"N1curl",
"triangle",
1,
[
lambda x: (-x[1], x[0]),
lambda x: (x[1], 1 - x[0]),
lambda x: (1.0 - x[1], x[0]),
],
),
(
"N1curl",
"tetrahedron",
1,
[
lambda x: (0.0, -x[2], x[1]),
lambda x: (-x[2], 0.0, x[0]),
lambda x: (-x[1], x[0], 0.0),
lambda x: (x[2], x[2], 1.0 - x[0] - x[1]),
lambda x: (x[1], 1.0 - x[0] - x[2], x[1]),
lambda x: (1.0 - x[1] - x[2], x[0], x[0]),
],
),
],
)
def test_values(family, cell, degree, functions):
# Create element
e = basix.ufl.element(family, cell, degree)
# Get some points and check basis function values at points
points = [random_point(cell) for i in range(5)]
tables = e.tabulate(0, np.array(points, dtype=np.float64))[0]
for x, t in zip(points, tables):
for i, f in enumerate(functions):
assert np.allclose(t[i :: len(functions)], f(x))
def test_hash():
e0 = basix.create_element(basix.ElementFamily.P, basix.CellType.interval, 1)
e1 = basix.create_element(basix.ElementFamily.P, basix.CellType.interval, 1)
e2 = basix.create_element(
basix.ElementFamily.P, basix.CellType.interval, 1, dof_ordering=[0, 1]
)
e3 = basix.create_element(basix.ElementFamily.P, basix.CellType.interval, 2)
e4 = basix.create_element(basix.ElementFamily.P, basix.CellType.triangle, 1)
wcoeffs = np.eye(3)
z = np.zeros((0, 2))
x = [
[np.array([[0.0, 0.0]]), np.array([[1.0, 0.0]]), np.array([[0.0, 1.0]])],
[z, z, z],
[z],
[],
]
z = np.zeros((0, 1, 0, 1))
M = [
[np.array([[[[1.0]]]]), np.array([[[[1.0]]]]), np.array([[[[1.0]]]])],
[z, z, z],
[z],
[],
]
e5 = basix.create_custom_element(
basix.CellType.triangle,
[],
wcoeffs,
x,
M,
0,
basix.MapType.L2Piola,
basix.SobolevSpace.L2,
False,
1,
1,
basix.PolysetType.standard,
)
e6 = basix.create_element(basix.ElementFamily.P, basix.CellType.quadrilateral, 2)
e7 = basix.create_tp_element(basix.ElementFamily.P, basix.CellType.quadrilateral, 2)
assert hash(e0) == hash(e1) == hash(e2)
different_elements = [e2, e3, e4, e5, e6, e7]
for i, d0 in enumerate(different_elements):
for d1 in different_elements[:i]:
assert hash(d0) != hash(d1)
|