File: test_cell.py

package info (click to toggle)
fenics-basix 0.10.0.post0-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 3,156 kB
  • sloc: cpp: 23,435; python: 10,829; makefile: 43; sh: 26
file content (143 lines) | stat: -rw-r--r-- 4,987 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
# Copyright (c) 2021-2025 Matthew Scroggs, Jørgen S. Dokken
# FEniCS Project
# SPDX-License-Identifier: MIT

import numpy as np
import pytest

import basix

cells = ["triangle", "quadrilateral", "tetrahedron", "hexahedron", "pyramid", "prism"]


@pytest.mark.parametrize("cell", cells)
def test_volume(cell):
    cell_type = getattr(basix.CellType, cell)
    volumes = {
        "triangle": 1 / 2,
        "quadrilateral": 1,
        "tetrahedron": 1 / 6,
        "hexahedron": 1,
        "pyramid": 1 / 3,
        "prism": 1 / 2,
    }
    assert np.isclose(basix.cell.volume(cell_type), volumes[cell])


@pytest.mark.parametrize("cell", cells)
def test_normals(cell):
    cell_type = getattr(basix.CellType, cell)
    normals = basix.cell.facet_normals(cell_type)
    facets = basix.topology(cell_type)[-2]
    geometry = basix.geometry(cell_type)
    for normal, facet in zip(normals, facets):
        assert np.isclose(np.linalg.norm(normal), 1)
        for v in facet[1:]:
            tangent = geometry[v] - geometry[facet[0]]
            assert np.isclose(np.dot(tangent, normal), 0)


@pytest.mark.parametrize("cell", cells)
def test_outward_normals(cell):
    cell_type = getattr(basix.CellType, cell)
    normals = basix.cell.facet_outward_normals(cell_type)
    facets = basix.topology(cell_type)[-2]
    geometry = basix.geometry(cell_type)
    midpoint = sum(geometry) / len(geometry)
    for normal, facet in zip(normals, facets):
        assert np.dot(geometry[facet[0]] - midpoint, normal) > 0


@pytest.mark.parametrize("cell", cells)
def test_facet_orientations(cell):
    cell_type = getattr(basix.CellType, cell)
    normals = basix.cell.facet_normals(cell_type)
    outward_normals = basix.cell.facet_outward_normals(cell_type)
    orientations = basix.cell.facet_orientations(cell_type)
    for n1, n2, orient in zip(normals, outward_normals, orientations):
        if orient:
            assert np.allclose(n1, -n2)
        else:
            assert np.allclose(n1, n2)


@pytest.mark.parametrize("cell", cells)
def test_sub_entity_connectivity(cell):
    cell_type = getattr(basix.CellType, cell)
    connectivity = basix.cell.sub_entity_connectivity(cell_type)
    topology = basix.topology(cell_type)
    assert len(connectivity) == len(topology)
    for dim, entities in enumerate(connectivity):
        assert len(entities) == len(topology[dim])
        for n, entity in enumerate(entities):
            for dim2, connected_entities in enumerate(entity):
                for n2 in connected_entities:
                    if dim > dim2:
                        for i in topology[dim2][n2]:
                            assert i in topology[dim][n]
                    else:
                        for i in topology[dim][n]:
                            assert i in topology[dim2][n2]


def test_sub_entity_type():
    cell_type = basix.CellType.tetrahedron
    for i in range(4):
        assert basix.cell.sub_entity_type(cell_type, 0, i) == basix.CellType.point
    for i in range(6):
        assert basix.cell.sub_entity_type(cell_type, 1, i) == basix.CellType.interval
    for i in range(4):
        assert basix.cell.sub_entity_type(cell_type, 2, i) == basix.CellType.triangle
    assert basix.cell.sub_entity_type(cell_type, 3, 0) == basix.CellType.tetrahedron


def test_facet_jacobians_2D_simplex():
    cell = basix.cell.CellType.triangle
    facet_jacobian = basix.cell.facet_jacobians(cell)

    reference_vertices = np.array([[0.0, 0.0], [1.0, 0.0], [0.0, 1.0]])
    mask = np.zeros(3, dtype=np.bool_)
    for i in range(3):
        mask[:] = True
        mask[i] = False
        facet = reference_vertices[mask]

        reference_facet_jacobian = -facet[0:1, :].T + facet[1:2, :].T
        np.testing.assert_allclose(reference_facet_jacobian, facet_jacobian[i])


def test_facet_jacobians_3D_simplex():
    cell = basix.cell.CellType.tetrahedron
    facet_jacobian = basix.cell.facet_jacobians(cell)

    reference_vertices = np.array(
        [[0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [0.0, 1.0, 0.0], [0.0, 0.0, 1.0]]
    )
    mask = np.zeros(4, dtype=np.bool_)
    for i in range(4):
        mask[:] = True
        mask[i] = False
        facet = reference_vertices[mask]
        reference_facet_jacobian = np.array([-facet[0] + facet[1], -facet[0] + facet[2]]).T
        np.testing.assert_allclose(reference_facet_jacobian, facet_jacobian[i])


@pytest.mark.parametrize(
    "cell",
    [
        basix.cell.CellType.hexahedron,
        basix.cell.CellType.tetrahedron,
        basix.cell.CellType.prism,
        basix.cell.CellType.pyramid,
    ],
)
def test_edge_jacobian_3D_simplex(cell):
    edge_jacobian = basix.cell.edge_jacobians(cell)
    geom = basix.geometry(cell)
    topology = basix.topology(cell)
    edges = topology[1]

    for i, edge in enumerate(edges):
        points = geom[edge]
        reference_edge_jacobian = (points[1:2, :] - points[0:1, :]).T
        np.testing.assert_allclose(reference_edge_jacobian, edge_jacobian[i])