File: test_dofmap.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 (564 lines) | stat: -rw-r--r-- 22,479 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
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
"""Unit tests for the fem interface"""

# Copyright (C) 2009-2014 Garth N. Wells
#
# 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/>.

import pytest
import numpy as np

import sys

from dolfin import *
from dolfin_utils.test import *


xfail = pytest.mark.xfail(strict=True)

@fixture
def mesh():
    return UnitSquareMesh(4, 4)


reorder_dofs = set_parameters_fixture("reorder_dofs_serial", [True, False])


@pytest.mark.parametrize('mesh_factory', [(UnitIntervalMesh, (8,)),
                                          (UnitSquareMesh, (4, 4)),
                                          (UnitCubeMesh, (2, 2, 2)),
                                          # cell.contains(Point) does not work correctly
                                          # for quad/hex cells once it is fixed, this test will pass
                                          pytest.param(((UnitSquareMesh.create, (4, 4, CellType.Type.quadrilateral))), marks=xfail),
                                          pytest.param(((UnitCubeMesh.create, (2, 2, 2, CellType.Type.hexahedron))), marks=xfail)])
def test_tabulate_all_coordinates(mesh_factory):
    func, args = mesh_factory
    mesh = func(*args)
    V = FunctionSpace(mesh, "Lagrange", 1)
    W0 = FiniteElement("Lagrange", mesh.ufl_cell(), 1)
    W1 = VectorElement("Lagrange", mesh.ufl_cell(), 1)
    W = FunctionSpace(mesh, W0*W1)

    D = mesh.geometry().dim()
    V_dofmap = V.dofmap()
    W_dofmap = W.dofmap()

    all_coords_V = V.tabulate_dof_coordinates()
    all_coords_W = W.tabulate_dof_coordinates()
    local_size_V = V_dofmap.ownership_range()[1]-V_dofmap.ownership_range()[0]
    local_size_W = W_dofmap.ownership_range()[1]-W_dofmap.ownership_range()[0]

    all_coords_V = all_coords_V.reshape(local_size_V, D)
    all_coords_W = all_coords_W.reshape(local_size_W, D)

    checked_V = [False]*local_size_V
    checked_W = [False]*local_size_W

    # Check that all coordinates are within the cell it should be
    for cell in cells(mesh):
        dofs_V = V_dofmap.cell_dofs(cell.index())
        for di in dofs_V:
            if di >= local_size_V:
                continue
            assert cell.contains(Point(all_coords_V[di]))
            checked_V[di] = True

        dofs_W = W_dofmap.cell_dofs(cell.index())
        for di in dofs_W:
            if di >= local_size_W:
                continue
            assert cell.contains(Point(all_coords_W[di]))
            checked_W[di] = True

    # Assert that all dofs have been checked by the above
    assert all(checked_V)
    assert all(checked_W)


@pytest.mark.parametrize('mesh_factory', [(UnitSquareMesh, (4, 4)), (UnitSquareMesh.create, (4, 4, CellType.Type.quadrilateral))])
def test_tabulate_dofs(mesh_factory):
    func, args = mesh_factory
    mesh = func(*args)
    W0 = FiniteElement("Lagrange", mesh.ufl_cell(), 1)
    W1 = VectorElement("Lagrange", mesh.ufl_cell(), 1)
    W = FunctionSpace(mesh, W0*W1)

    L0 = W.sub(0)
    L1 = W.sub(1)
    L01 = L1.sub(0)
    L11 = L1.sub(1)

    for i, cell in enumerate(cells(mesh)):
        dofs0 = L0.dofmap().cell_dofs(cell.index())
        dofs1 = L01.dofmap().cell_dofs(cell.index())
        dofs2 = L11.dofmap().cell_dofs(cell.index())
        dofs3 = L1.dofmap().cell_dofs(cell.index())

        assert np.array_equal(dofs0, L0.dofmap().cell_dofs(i))
        assert np.array_equal(dofs1, L01.dofmap().cell_dofs(i))
        assert np.array_equal(dofs2, L11.dofmap().cell_dofs(i))
        assert np.array_equal(dofs3, L1.dofmap().cell_dofs(i))

        assert len(np.intersect1d(dofs0, dofs1)) == 0
        assert len(np.intersect1d(dofs0, dofs2)) == 0
        assert len(np.intersect1d(dofs1, dofs2)) == 0
        assert np.array_equal(np.append(dofs1, dofs2), dofs3)


@pytest.mark.parametrize('mesh_factory', [(UnitSquareMesh, (4, 4)), (UnitSquareMesh.create, (4, 4, CellType.Type.quadrilateral))])
def test_tabulate_coord_periodic(mesh_factory):

    class PeriodicBoundary2(SubDomain):
        def inside(self, x, on_boundary):
            return x[0] < DOLFIN_EPS

        def map(self, x, y):
            y[0] = x[0] - 1.0
            y[1] = x[1]

    # Create periodic boundary condition
    periodic_boundary = PeriodicBoundary2()

    func, args = mesh_factory
    mesh = func(*args)

    V = FiniteElement("Lagrange", mesh.ufl_cell(), 1)
    Q = VectorElement("Lagrange", mesh.ufl_cell(), 1)
    W = V*Q

    V = FunctionSpace(mesh, V, constrained_domain=periodic_boundary)
    W = FunctionSpace(mesh, W, constrained_domain=periodic_boundary)

    L0 = W.sub(0)
    L1 = W.sub(1)
    L01 = L1.sub(0)
    L11 = L1.sub(1)

    sdim = V.element().space_dimension()
    coord0 = np.zeros((sdim, 2), dtype="d")
    coord1 = np.zeros((sdim, 2), dtype="d")
    coord2 = np.zeros((sdim, 2), dtype="d")
    coord3 = np.zeros((sdim, 2), dtype="d")

    for cell in cells(mesh):
        coord0 = V.element().tabulate_dof_coordinates(cell)
        coord1 = L0.element().tabulate_dof_coordinates(cell)
        coord2 = L01.element().tabulate_dof_coordinates(cell)
        coord3 = L11.element().tabulate_dof_coordinates(cell)
        coord4 = L1.element().tabulate_dof_coordinates(cell)

        assert (coord0 == coord1).all()
        assert (coord0 == coord2).all()
        assert (coord0 == coord3).all()
        assert (coord4[:sdim] == coord0).all()
        assert (coord4[sdim:] == coord0).all()


@pytest.mark.parametrize('mesh_factory', [(UnitSquareMesh, (5, 5)), (UnitSquareMesh.create, (5, 5, CellType.Type.quadrilateral))])
def test_tabulate_dofs_periodic(mesh_factory):

    class PeriodicBoundary2(SubDomain):
        def inside(self, x, on_boundary):
            return x[0] < DOLFIN_EPS

        def map(self, x, y):
            y[0] = x[0] - 1.0
            y[1] = x[1]

    func, args = mesh_factory
    mesh = func(*args)

    # Create periodic boundary
    periodic_boundary = PeriodicBoundary2()

    V = FiniteElement("Lagrange", mesh.ufl_cell(), 2)
    Q = VectorElement("Lagrange", mesh.ufl_cell(), 2)
    W = V*Q

    V = FunctionSpace(mesh, V, constrained_domain=periodic_boundary)
    Q = FunctionSpace(mesh, Q, constrained_domain=periodic_boundary)
    W = FunctionSpace(mesh, W, constrained_domain=periodic_boundary)

    L0 = W.sub(0)
    L1 = W.sub(1)
    L01 = L1.sub(0)
    L11 = L1.sub(1)

    # Check dimensions
    assert V.dim() == 110
    assert Q.dim() == 220
    assert L0.dim() == V.dim()
    assert L1.dim() == Q.dim()
    assert L01.dim() == V.dim()
    assert L11.dim() == V.dim()

    for i, cell in enumerate(cells(mesh)):
        dofs0 = L0.dofmap().cell_dofs(cell.index())
        dofs1 = L01.dofmap().cell_dofs(cell.index())
        dofs2 = L11.dofmap().cell_dofs(cell.index())
        dofs3 = L1.dofmap().cell_dofs(cell.index())

        assert np.array_equal(dofs0, L0.dofmap().cell_dofs(i))
        assert np.array_equal(dofs1, L01.dofmap().cell_dofs(i))
        assert np.array_equal(dofs2, L11.dofmap().cell_dofs(i))
        assert np.array_equal(dofs3, L1.dofmap().cell_dofs(i))

        assert len(np.intersect1d(dofs0, dofs1)) == 0
        assert len(np.intersect1d(dofs0, dofs2)) == 0
        assert len(np.intersect1d(dofs1, dofs2)) == 0
        assert np.array_equal(np.append(dofs1, dofs2), dofs3)


@pytest.mark.parametrize('mesh_factory', [(UnitSquareMesh, (3, 3)), (UnitSquareMesh.create, (3, 3, CellType.Type.quadrilateral))])
def test_global_dof_builder(mesh_factory):
    func, args = mesh_factory
    mesh = func(*args)

    V = VectorElement("CG", mesh.ufl_cell(), 1)
    Q = FiniteElement("CG", mesh.ufl_cell(), 1)
    R = FiniteElement("R",  mesh.ufl_cell(), 0)

    W = FunctionSpace(mesh, MixedElement([Q, Q, Q, R]))
    W = FunctionSpace(mesh, MixedElement([Q, Q, R, Q]))
    W = FunctionSpace(mesh, V*R)
    W = FunctionSpace(mesh, R*V)


@pytest.mark.parametrize('mesh_factory', [(UnitSquareMesh, (3, 3)), (UnitSquareMesh.create, (3, 3, CellType.Type.quadrilateral))])
def test_dof_to_vertex_map(mesh_factory, reorder_dofs):
    func, args = mesh_factory
    mesh = func(*args)

    def _test_maps_consistency(space):
        v2d = vertex_to_dof_map(space)
        d2v = dof_to_vertex_map(space)
        assert len(v2d) == len(d2v)
        assert np.all(v2d[d2v] == np.arange(len(v2d)))
        assert np.all(d2v[v2d] == np.arange(len(d2v)))

    # Check for both reordered and UFC ordered dofs
    v = FiniteElement("Lagrange", mesh.ufl_cell(), 1)
    q = VectorElement("Lagrange", mesh.ufl_cell(), 1)
    w = v*q

    V = FunctionSpace(mesh, v)
    Q = FunctionSpace(mesh, q)
    W = FunctionSpace(mesh, w)

    _test_maps_consistency(V)
    _test_maps_consistency(Q)
    _test_maps_consistency(W)

    u = Function(V)
    e = Expression("x[0] + x[1]", degree=1)
    u.interpolate(e)

    vert_values = mesh.coordinates().sum(1)
    func_values = u.vector().get_local(vertex_to_dof_map(V))
    assert round(max(abs(func_values - vert_values)), 7) == 0

    c0 = Constant((1, 2))
    u0 = Function(Q)
    u0.interpolate(c0)

    vert_values = np.zeros(mesh.num_vertices()*2)
    u1 = Function(Q)
    vert_values[::2] = 1
    vert_values[1::2] = 2

    dim = Q.dofmap().index_map().size(IndexMap.MapSize.OWNED)
    u1.vector().set_local(vert_values[dof_to_vertex_map(Q)[:dim]].copy())
    assert round((u0.vector()-u1.vector()).sum() - 0.0, 7) == 0

    W = FunctionSpace(mesh, "DG", 0)
    with pytest.raises(RuntimeError):
        dof_to_vertex_map(W)

    W = FunctionSpace(mesh, q*FiniteElement("R", mesh.ufl_cell(), 0))
    with pytest.raises(RuntimeError):
        dof_to_vertex_map(W)

    W = FunctionSpace(mesh, "CG", 2)
    with pytest.raises(RuntimeError):
        dof_to_vertex_map(W)

    W = VectorFunctionSpace(mesh, "CG", 1)
    with pytest.raises(RuntimeError):
        dof_to_vertex_map(W.sub(0))


def test_entity_dofs(mesh):

    # Test that num entity dofs is correctly wrapped to
    # dolfin::DofMap
    V = FunctionSpace(mesh, "CG", 1)
    assert V.dofmap().num_entity_dofs(0) == 1
    assert V.dofmap().num_entity_dofs(1) == 0
    assert V.dofmap().num_entity_dofs(2) == 0

    V = VectorFunctionSpace(mesh, "CG", 1)
    assert V.dofmap().num_entity_dofs(0) == 2
    assert V.dofmap().num_entity_dofs(1) == 0
    assert V.dofmap().num_entity_dofs(2) == 0

    V = FunctionSpace(mesh, "CG", 2)
    assert V.dofmap().num_entity_dofs(0) == 1
    assert V.dofmap().num_entity_dofs(1) == 1
    assert V.dofmap().num_entity_dofs(2) == 0

    V = FunctionSpace(mesh, "CG", 3)
    assert V.dofmap().num_entity_dofs(0) == 1
    assert V.dofmap().num_entity_dofs(1) == 2
    assert V.dofmap().num_entity_dofs(2) == 1

    V = FunctionSpace(mesh, "DG", 0)
    assert V.dofmap().num_entity_dofs(0) == 0
    assert V.dofmap().num_entity_dofs(1) == 0
    assert V.dofmap().num_entity_dofs(2) == 1

    V = FunctionSpace(mesh, "DG", 1)
    assert V.dofmap().num_entity_dofs(0) == 0
    assert V.dofmap().num_entity_dofs(1) == 0
    assert V.dofmap().num_entity_dofs(2) == 3

    V = VectorFunctionSpace(mesh, "CG", 1)

    # Note this numbering is dependent on FFC and can change This test
    # is here just to check that we get correct numbers mapped from
    # ufc generated code to dolfin
    for i, cdofs in enumerate([[0, 3], [1, 4], [2, 5]]):
        dofs = V.dofmap().tabulate_entity_dofs(0, i)
        assert all(d == cd for d, cd in zip(dofs, cdofs))


@skip_in_parallel
@pytest.mark.parametrize('mesh_factory', [(UnitSquareMesh, (2, 2)), (UnitSquareMesh.create, (2, 2, CellType.Type.quadrilateral))])
def test_entity_closure_dofs(mesh_factory):
    func, args = mesh_factory
    mesh = func(*args)
    tdim = mesh.topology().dim()

    for degree in (1, 2, 3):
        V = FunctionSpace(mesh, "CG", degree)
        for d in range(tdim + 1):
            covered = set()
            covered2 = set()
            all_entities = np.array([entity for entity in range(mesh.num_entities(d))], dtype=np.uintp)
            for entity in all_entities:
                entities = np.array([entity], dtype=np.uintp)
                dofs_on_this_entity = V.dofmap().entity_dofs(mesh, d, entities)
                closure_dofs = V.dofmap().entity_closure_dofs(mesh, d, entities)
                assert len(dofs_on_this_entity) == V.dofmap().num_entity_dofs(d)
                assert len(dofs_on_this_entity) <= len(closure_dofs)
                covered.update(dofs_on_this_entity)
                covered2.update(closure_dofs)
            dofs_on_all_entities = V.dofmap().entity_dofs(mesh, d, all_entities)
            closure_dofs_on_all_entities = V.dofmap().entity_closure_dofs(mesh, d, all_entities)
            assert len(dofs_on_all_entities) == V.dofmap().num_entity_dofs(d) * mesh.num_entities(d)
            assert covered == set(dofs_on_all_entities)
            assert covered2 == set(closure_dofs_on_all_entities)
        d = tdim
        all_cells = np.array([entity for entity in range(mesh.num_entities(d))], dtype=np.uintp)
        assert set(V.dofmap().entity_closure_dofs(mesh, d, all_cells)) == set(range(V.dim()))


def test_clear_sub_map_data_scalar(mesh):
    V = FunctionSpace(mesh, "CG", 2)
    with pytest.raises(ValueError):
        V.sub(1)

    V = VectorFunctionSpace(mesh, "CG", 2)
    V1 = V.sub(1)

    # Clean sub-map data
    V.dofmap().clear_sub_map_data()

    # Can still get previously computed map
    V1 = V.sub(1)

    # New sub-map should throw an error
    with pytest.raises(RuntimeError):
        V.sub(0)


def test_clear_sub_map_data_vector(mesh):
    mesh = UnitSquareMesh(8, 8)
    P1 = FiniteElement("Lagrange", mesh.ufl_cell(), 1)
    W = FunctionSpace(mesh, P1*P1)

    # Check block size
    assert W.dofmap().block_size() == 2

    W.dofmap().clear_sub_map_data()
    with pytest.raises(RuntimeError):
        W0 = W.sub(0)
    with pytest.raises(RuntimeError):
        W1 = W.sub(1)


def test_block_size(mesh):
    meshes = [UnitSquareMesh(8, 8), UnitCubeMesh(4, 4, 4),
              UnitSquareMesh.create(8, 8, CellType.Type.quadrilateral),
              UnitCubeMesh.create(4, 4, 4, CellType.Type.hexahedron)]
    for mesh in meshes:
        P2 = FiniteElement("Lagrange", mesh.ufl_cell(), 2)

        V = FunctionSpace(mesh, P2)
        assert V.dofmap().block_size() == 1

        V = FunctionSpace(mesh, P2*P2)
        assert V.dofmap().block_size() == 2

        for i in range(1, 6):
            W = FunctionSpace(mesh, MixedElement(i*[P2]))
            assert W.dofmap().block_size() == i

        V = VectorFunctionSpace(mesh, "Lagrange", 2)
        assert V.dofmap().block_size() == mesh.geometry().dim()


def test_block_size_real(mesh):
    mesh = UnitIntervalMesh(12)
    V = FiniteElement('DG', mesh.ufl_cell(), 0)
    R = FiniteElement('R',  mesh.ufl_cell(), 0)
    X = FunctionSpace(mesh, V*R)
    assert X.dofmap().block_size() == 1


@skip_in_serial
@pytest.mark.parametrize('mesh_factory', [(UnitIntervalMesh, (8,)),
                                          (UnitSquareMesh, (4, 4)),
                                          (UnitCubeMesh, (2, 2, 2)),
                                          (UnitSquareMesh.create, (4, 4, CellType.Type.quadrilateral)),
                                          (UnitCubeMesh.create, (2, 2, 2, CellType.Type.hexahedron))])
def test_mpi_dofmap_stats(mesh_factory):
    func, args = mesh_factory
    mesh = func(*args)

    V = FunctionSpace(mesh, "CG", 1)
    assert len(V.dofmap().shared_nodes()) > 0
    neighbours = V.dofmap().neighbours()
    for processes in V.dofmap().shared_nodes().values():
        for process in processes:
            assert process in neighbours

    for owner in V.dofmap().off_process_owner():
        assert owner in neighbours

@pytest.mark.parametrize('mesh_factory', [(UnitSquareMesh, (4, 4)), (UnitSquareMesh.create, (4, 4, CellType.Type.quadrilateral))])
def test_local_dimension(mesh_factory):
    func, args = mesh_factory
    mesh = func(*args)

    v = FiniteElement("Lagrange", mesh.ufl_cell(), 1)
    q = VectorElement("Lagrange", mesh.ufl_cell(), 1)
    w = v*q

    V = FunctionSpace(mesh, v)
    Q = FunctionSpace(mesh, q)
    W = FunctionSpace(mesh, w)

    for space in [V, Q, W]:
        dofmap = space.dofmap()
        local_to_global_map = dofmap.tabulate_local_to_global_dofs()
        ownership_range = dofmap.ownership_range()
        dim1 = dofmap.index_map().size(IndexMap.MapSize.OWNED)
        dim2 = dofmap.index_map().size(IndexMap.MapSize.UNOWNED)
        dim3 = dofmap.index_map().size(IndexMap.MapSize.ALL)
        assert dim1 == ownership_range[1] - ownership_range[0]
        assert dim3 == local_to_global_map.size
        assert dim1 + dim2 == dim3
        # with pytest.raises(RuntimeError):
        #    dofmap.index_map().size('foo')


# Failures in FFC on quads/hexes
xfail_ffc = pytest.mark.xfail(raises=Exception, strict=True)

@skip_in_parallel
@pytest.mark.parametrize('space', [
    "FunctionSpace(UnitIntervalMesh.create(10),                                        'P', 1)",
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.triangle),                'P', 1)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.tetrahedron),            'P', 1)",
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.quadrilateral),           'Q', 1)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.hexahedron),             'Q', 1)",
    "FunctionSpace(UnitIntervalMesh.create(10),                                        'P', 2)",
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.triangle),                'P', 2)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.tetrahedron),            'P', 2)",
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.quadrilateral),           'Q', 2)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.hexahedron),             'Q', 2)",
    "FunctionSpace(UnitIntervalMesh.create(10),                                        'P', 3)",
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.triangle),                'P', 3)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.tetrahedron),            'P', 3)",
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.quadrilateral),           'Q', 3)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.hexahedron),             'Q', 3)",
    "FunctionSpace(UnitIntervalMesh.create(10),                                        'DP', 1)",
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.triangle),                'DP', 1)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.tetrahedron),            'DP', 1)",
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.quadrilateral),           'DQ', 1)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.hexahedron),             'DQ', 1)",
    "FunctionSpace(UnitIntervalMesh.create(10),                                        'DP', 2)",
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.triangle),                'DP', 2)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.tetrahedron),            'DP', 2)",
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.quadrilateral),           'DQ', 2)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.hexahedron),             'DQ', 2)",
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.triangle),                'N1curl', 1)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.tetrahedron),            'N1curl', 1)",
    pytest.param(("FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.quadrilateral), 'N1curl', 1)"), marks=xfail_ffc),
    pytest.param(("FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.hexahedron),   'N1curl', 1)"), marks=xfail_ffc),
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.triangle),                'N1curl', 2)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.tetrahedron),            'N1curl', 2)",
    pytest.param(("FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.quadrilateral), 'N1curl', 2)"), marks=xfail_ffc),
    pytest.param(("FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.hexahedron),   'N1curl', 2)"), marks=xfail_ffc),
    "FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.triangle),                'RT', 1)",
    "FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.tetrahedron),            'RT', 1)",
    pytest.param(("FunctionSpace(UnitSquareMesh.create(6, 6, CellType.Type.quadrilateral), 'RT', 1)"), marks=xfail_ffc),
    pytest.param(("FunctionSpace(UnitCubeMesh.create(2, 2, 2, CellType.Type.hexahedron),   'RT', 1)"), marks=xfail_ffc),
])
def test_dofs_dim(space):
    """Test function GenericDofMap::dofs(mesh, dim)"""
    V = eval(space)
    dofmap = V.dofmap()
    mesh = V.mesh()
    for dim in range(0, mesh.topology().dim()):
        edofs = dofmap.dofs(mesh, dim)
        num_mesh_entities = mesh.num_entities(dim)
        dofs_per_entity = dofmap.num_entity_dofs(dim)
        assert len(edofs) == dofs_per_entity*num_mesh_entities


def test_readonly_view_local_to_global_unwoned(mesh):
    """Test that local_to_global_unwoned() returns readonly
    view into the data; in particular test lifetime of data
    owner"""
    V = FunctionSpace(mesh, "P", 1)
    dofmap = V.dofmap()
    index_map = dofmap.index_map()

    rc = sys.getrefcount(dofmap)
    l2gu = dofmap.local_to_global_unowned()
    assert sys.getrefcount(dofmap) == rc + 1 if l2gu.size else rc
    assert not l2gu.flags.writeable
    assert all(l2gu < V.dofmap().global_dimension())
    del l2gu
    assert sys.getrefcount(dofmap) == rc

    rc = sys.getrefcount(index_map)
    l2gu = index_map.local_to_global_unowned()
    assert sys.getrefcount(index_map) == rc + 1 if l2gu.size else rc
    assert not l2gu.flags.writeable
    assert all(l2gu < V.dofmap().global_dimension())
    del l2gu
    assert sys.getrefcount(index_map) == rc