File: test_function_space.py

package info (click to toggle)
dolfin 2018.1.0.post1-16
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 28,764 kB
  • sloc: xml: 104,040; cpp: 98,856; python: 22,511; makefile: 204; sh: 182
file content (246 lines) | stat: -rwxr-xr-x 6,739 bytes parent folder | download | duplicates (2)
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
b"""Unit tests for the FunctionSpace class"""

# Copyright (C) 2011 Johan Hake
#
# This file is part of DOLFIN.
#
# DOLFIN 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.
#
# DOLFIN 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 DOLFIN. If not, see <http://www.gnu.org/licenses/>.
#
# Modified by Oeyvind Evju 2013
#
# First added:  2011-09-21
# Last changed: 2013-10-11

import pytest
from dolfin import *
from ufl.log import UFLException
from dolfin_utils.test import fixture


@fixture
def mesh():
    return UnitCubeMesh(8, 8, 8)


@fixture
def V(mesh):
    return FunctionSpace(mesh, 'CG', 1)


@fixture
def W(mesh):
    return VectorFunctionSpace(mesh, 'CG', 1)


@fixture
def Q(mesh):
    W = VectorElement('CG', mesh.ufl_cell(), 1)
    V = FiniteElement('CG', mesh.ufl_cell(), 1)
    return FunctionSpace(mesh, W*V)


@fixture
def f(V):
    return Function(V)


@fixture
def V2(f):
    return f.function_space()


@fixture
def g(W):
    return Function(W)


@fixture
def W2(g):
    return g.function_space()


def test_python_interface(V, V2, W, W2, Q):
    # Test Python interface of cpp generated FunctionSpace
    assert isinstance(V, FunctionSpace)
    assert isinstance(W, FunctionSpace)
    assert isinstance(V2, FunctionSpace)
    assert isinstance(W2, FunctionSpace)

    assert V.ufl_cell() == V2.ufl_cell()
    assert W.ufl_cell() == W2.ufl_cell()
    assert V.dolfin_element().signature() == V2.dolfin_element().signature()
    assert W.dolfin_element().signature() == W2.dolfin_element().signature()
    assert V.ufl_element() == V2.ufl_element()
    assert W.ufl_element() == W2.ufl_element()
    assert W.id() == W2.id()
    assert V.id() == V2.id()


def test_component(V, W, Q):
    assert not W.component()
    assert not V.component()
    assert W.sub(0).component()[0] == 0
    assert W.sub(1).component()[0] == 1
    assert Q.sub(0).component()[0] == 0
    assert Q.sub(1).component()[0] == 1


def test_equality(V, V2, W, W2):
    assert V == V
    assert V == V2
    assert W == W
    assert W == W2


def test_inclusion(V, Q):
    assert V.contains(V)
    assert not Q.contains(V)

    assert Q.contains(Q)
    assert Q.contains(Q.sub(0))
    assert Q.contains(Q.sub(1))
    assert Q.contains(Q.sub(0).sub(0))
    assert Q.contains(Q.sub(0).sub(1))
    assert Q.contains(Q.extract_sub_space((0, 0)))
    assert Q.contains(Q.extract_sub_space((0, 1)))
    assert Q.contains(Q.extract_sub_space((1,)))

    assert not Q.sub(0).contains(Q)
    assert Q.sub(0).contains(Q.sub(0))
    assert not Q.sub(0).contains(Q.sub(1))
    assert Q.sub(0).contains(Q.sub(0).sub(0))
    assert Q.sub(0).contains(Q.sub(0).sub(1))
    assert Q.sub(0).contains(Q.extract_sub_space((0, 0)))
    assert Q.sub(0).contains(Q.extract_sub_space((0, 1)))
    assert not Q.sub(0).contains(Q.extract_sub_space((1,)))

    assert not Q.sub(1).contains(Q)
    assert not Q.sub(1).contains(Q.sub(0))
    assert Q.sub(1).contains(Q.sub(1))
    assert not Q.sub(1).contains(Q.sub(0).sub(0))
    assert not Q.sub(1).contains(Q.sub(0).sub(1))
    assert not Q.sub(1).contains(Q.extract_sub_space((0, 0)))
    assert not Q.sub(1).contains(Q.extract_sub_space((0, 1)))
    assert Q.sub(1).contains(Q.extract_sub_space((1,)))

    assert not Q.sub(0).sub(0).contains(Q)
    assert not Q.sub(0).sub(0).contains(Q.sub(0))
    assert not Q.sub(0).sub(0).contains(Q.sub(1))
    assert Q.sub(0).sub(0).contains(Q.sub(0).sub(0))
    assert not Q.sub(0).sub(0).contains(Q.sub(0).sub(1))
    assert Q.sub(0).sub(0).contains(Q.extract_sub_space((0, 0)))
    assert not Q.sub(0).sub(0).contains(Q.extract_sub_space((0, 1)))
    assert not Q.sub(0).sub(0).contains(Q.extract_sub_space((1,)))


def test_boundary(mesh):
    bmesh = BoundaryMesh(mesh, "exterior")
    Vb = FunctionSpace(bmesh, "DG", 0)
    Wb = VectorFunctionSpace(bmesh, "CG", 1)
    assert Vb.dim() == 768
    assert Wb.dim() == 1158


def test_not_equal(W, V, W2, V2):
    assert W != V
    assert W2 != V2


def test_sub_equality(W, Q):
    assert W.sub(0) == W.sub(0)
    assert W.sub(0) != W.sub(1)
    assert W.sub(0) == W.extract_sub_space([0])
    assert W.sub(1) == W.extract_sub_space([1])
    assert Q.sub(0) == Q.extract_sub_space([0])


def test_in_operator(f, g, V, V2, W, W2, mesh):
    assert f in V
    assert f in V2
    assert g in W
    assert g in W2
    with pytest.raises(RuntimeError):
        mesh in V


def test_collapse(W, V):
    Vs = W.sub(2)
    with pytest.raises(RuntimeError):
        Function(Vs)
    assert Vs.dofmap().cell_dofs(0)[0] != V.dofmap().cell_dofs(0)[0]

    # Collapse the space it should now be the same as V
    Vc, dofmap_new_old = Vs.collapse(True)
    assert Vc.dofmap().cell_dofs(0)[0] == V.dofmap().cell_dofs(0)[0]
    f0 = Function(V)
    f1 = Function(Vc)
    assert len(f0.vector()) == len(f1.vector())


def test_argument_equality(mesh, V, V2, W, W2):
    """Placed this test here because it's mainly about detecting differing
    function spaces.

    """
    mesh2 = UnitCubeMesh(8, 8, 8)
    V3 = FunctionSpace(mesh2, 'CG', 1)
    W3 = VectorFunctionSpace(mesh2, 'CG', 1)

    for TF in (TestFunction, TrialFunction):
        v = TF(V)
        v2 = TF(V2)
        v3 = TF(V3)
        assert v == v2
        assert v2 == v
        assert V != V3
        assert V2 != V3
        assert not v == v3
        assert not v2 == v3
        assert v != v3
        assert v2 != v3
        assert v != v3
        assert v2 != v3

        w = TF(W)
        w2 = TF(W2)
        w3 = TF(W3)
        assert w == w2
        assert w2 == w
        assert w != w3
        assert w2 != w3

        assert v != w
        assert w != v

        s1 = set((v, w))
        s2 = set((v2, w2))
        s3 = set((v, v2, w, w2))
        assert len(s1) == 2
        assert len(s2) == 2
        assert len(s3) == 2
        assert s1 == s2
        assert s1 == s3
        assert s2 == s3

        # Test that the dolfin implementation of Argument.__eq__
        # is triggered when comparing ufl expressions
        assert grad(v) == grad(v2)
        assert grad(v) != grad(v3)


def test_cell_mismatch(mesh):
    """Test that cell mismatch raises early enough from UFL"""
    element = FiniteElement("P", triangle, 1)
    with pytest.raises(UFLException):
        FunctionSpace(mesh, element)