File: test_ghost_mesh.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 (135 lines) | stat: -rw-r--r-- 4,834 bytes parent folder | download | duplicates (5)
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
"Unit tests for ghosted meshes"

# Copyright (C) 2016 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
from dolfin import *
import os

from dolfin_utils.test import (fixture, skip_in_parallel, xfail_in_parallel, cd_tempdir,
                               pushpop_parameters)


# See https://bitbucket.org/fenics-project/dolfin/issues/579
def xtest_ghost_vertex_1d(pushpop_parameters):
    parameters["ghost_mode"] = "shared_vertex"
    mesh = UnitIntervalMesh(20)
    #print("Test: {}".format(MPI.sum(mesh.mpi_comm(), mesh.num_cells())))


def xtest_ghost_facet_1d(pushpop_parameters):
    parameters["ghost_mode"] = "shared_facet"
    mesh = UnitIntervalMesh(20)


def test_ghost_2d(pushpop_parameters):
    modes = ["shared_vertex", "shared_facet"]
    for mode in modes:
        parameters["ghost_mode"] = mode
        N = 8
        num_cells = 128

        mesh = UnitSquareMesh(N, N)
        if MPI.size(mesh.mpi_comm()) > 1:
            assert MPI.sum(mesh.mpi_comm(), mesh.num_cells()) > num_cells

        parameters["reorder_cells_gps"] = True
        parameters["reorder_vertices_gps"] = False
        mesh = UnitSquareMesh(N, N)
        if MPI.size(mesh.mpi_comm()) > 1:
            assert MPI.sum(mesh.mpi_comm(), mesh.num_cells()) > num_cells

        parameters["reorder_cells_gps"] = True
        parameters["reorder_vertices_gps"] = True
        mesh = UnitSquareMesh(N, N)
        if MPI.size(mesh.mpi_comm()) > 1:
            assert MPI.sum(mesh.mpi_comm(), mesh.num_cells()) > num_cells

        parameters["reorder_cells_gps"] = False
        parameters["reorder_vertices_gps"] = True
        mesh = UnitSquareMesh(N, N)
        if MPI.size(mesh.mpi_comm()) > 1:
            assert MPI.sum(mesh.mpi_comm(), mesh.num_cells()) > num_cells


def test_ghost_3d(pushpop_parameters):
    modes = ["shared_vertex", "shared_facet"]
    for mode in modes:
        parameters["ghost_mode"] = mode
        N = 2
        num_cells = 48

        mesh = UnitCubeMesh(N, N, N)
        if MPI.size(mesh.mpi_comm()) > 1:
            assert MPI.sum(mesh.mpi_comm(), mesh.num_cells()) > num_cells

        parameters["reorder_cells_gps"] = True
        parameters["reorder_vertices_gps"] = False
        mesh = UnitCubeMesh(N, N, N)
        if MPI.size(mesh.mpi_comm()) > 1:
            assert MPI.sum(mesh.mpi_comm(), mesh.num_cells()) > num_cells

        parameters["reorder_cells_gps"] = True
        parameters["reorder_vertices_gps"] = True
        mesh = UnitCubeMesh(N, N, N)
        if MPI.size(mesh.mpi_comm()) > 1:
            assert MPI.sum(mesh.mpi_comm(), mesh.num_cells()) > num_cells

        parameters["reorder_cells_gps"] = False
        parameters["reorder_vertices_gps"] = True
        mesh = UnitCubeMesh(N, N, N)
        if MPI.size(mesh.mpi_comm()) > 1:
            assert MPI.sum(mesh.mpi_comm(), mesh.num_cells()) > num_cells


@pytest.mark.parametrize('gmode', ['shared_vertex', 'shared_facet', 'none'])
def test_ghost_connectivities(gmode, pushpop_parameters):
    parameters['ghost_mode'] = gmode

    # Ghosted mesh
    meshG = UnitSquareMesh(MPI.comm_world, 4, 4)
    meshG.init(1, 2)

    # Reference mesh, not ghosted, not parallel
    meshR = UnitSquareMesh(MPI.comm_self, 4, 4)
    meshR.init(1, 2)

    # Create reference mapping from facet midpoint to cell midpoint
    reference = {}
    for facet in facets(meshR):
        fidx = facet.index()
        facet_mp = tuple(facet.midpoint()[:])
        reference[facet_mp] = []
        for cidx in meshR.topology()(1, 2)(fidx):
            cell = Cell(meshR, cidx)
            cell_mp = tuple(cell.midpoint()[:])
            reference[facet_mp].append(cell_mp)

    # Loop through ghosted mesh and check connectivities
    allowable_cell_indices = [cell.index() for cell in cells(meshG, 'all')]
    for facet in facets(meshG, 'regular'):
        fidx = facet.index()
        facet_mp = tuple(facet.midpoint()[:])
        assert facet_mp in reference

        for cidx in meshG.topology()(1, 2)(fidx):
            assert cidx in allowable_cell_indices
            cell = Cell(meshG, cidx)
            cell_mp = tuple(cell.midpoint()[:])
            assert cell_mp in reference[facet_mp]