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
|
"""Unit tests for MeshValueCollection"""
# 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/>.
#
# First added: 2011-03-10
# Last changed: 2011-03-10
import numpy.random
from dolfin import *
def test_assign_2D_cells():
mesh = UnitSquareMesh(3, 3)
ncells = mesh.num_cells()
f = MeshValueCollection("int", mesh, 2)
all_new = True
for cell in cells(mesh):
value = ncells - cell.index()
all_new = all_new and f.set_value(cell.index(), value)
g = MeshValueCollection("int", mesh, 2)
g.assign(f)
assert ncells == f.size()
assert ncells == g.size()
assert all_new
for cell in cells(mesh):
value = ncells - cell.index()
assert value, g.get_value(cell.index() == 0)
old_value = g.get_value(0, 0)
g.set_value(0, 0, old_value + 1)
assert old_value + 1 == g.get_value(0, 0)
def test_assign_2D_facets():
mesh = UnitSquareMesh(3, 3)
mesh.init(2, 1)
ncells = mesh.num_cells()
f = MeshValueCollection("int", mesh, 1)
all_new = True
for cell in cells(mesh):
value = ncells - cell.index()
for i, facet in enumerate(facets(cell)):
all_new = all_new and f.set_value(cell.index(), i, value + i)
g = MeshValueCollection("int", mesh, 1)
g.assign(f)
assert ncells*3 == f.size()
assert ncells*3 == g.size()
assert all_new
for cell in cells(mesh):
value = ncells - cell.index()
for i, facet in enumerate(facets(cell)):
assert value+i == g.get_value(cell.index(), i)
def test_assign_2D_vertices():
mesh = UnitSquareMesh(3, 3)
mesh.init(2, 0)
ncells = mesh.num_cells()
f = MeshValueCollection("int", mesh, 0)
all_new = True
for cell in cells(mesh):
value = ncells - cell.index()
for i, vert in enumerate(vertices(cell)):
all_new = all_new and f.set_value(cell.index(), i, value+i)
g = MeshValueCollection("int", mesh, 0)
g.assign(f)
assert ncells*3 == f.size()
assert ncells*3 == g.size()
assert all_new
for cell in cells(mesh):
value = ncells - cell.index()
for i, vert in enumerate(vertices(cell)):
assert value+i == g.get_value(cell.index(), i)
def test_mesh_function_assign_2D_cells():
mesh = UnitSquareMesh(3, 3)
ncells = mesh.num_cells()
f = MeshFunction("int", mesh, mesh.topology().dim())
for cell in cells(mesh):
f[cell] = ncells - cell.index()
g = MeshValueCollection("int", mesh, 2)
g.assign(f)
assert ncells == f.size()
assert ncells == g.size()
f2 = MeshFunction("int", mesh, g)
for cell in cells(mesh):
value = ncells - cell.index()
assert value == g.get_value(cell.index(), 0)
assert f2[cell] == g.get_value(cell.index(), 0)
h = MeshValueCollection("int", mesh, 2)
global_indices = mesh.topology().global_indices(2)
ncells_global = mesh.num_entities_global(2)
for cell in cells(mesh):
if global_indices[cell.index()] in [5, 8, 10]:
continue
value = ncells_global - global_indices[cell.index()]
h.set_value(cell.index(), int(value))
f3 = MeshFunction("int", mesh, h)
values = f3.array()
values[values > ncells_global] = 0.
info(str(values))
info(str(values.sum()))
assert MPI.sum(mesh.mpi_comm(), values.sum()*1.0) == 140.
def test_mesh_function_assign_2D_facets():
mesh = UnitSquareMesh(3, 3)
mesh.init(1)
f = MeshFunction("int", mesh, mesh.topology().dim()-1, 25)
for cell in cells(mesh):
for i, facet in enumerate(facets(cell)):
assert 25 == f[facet]
g = MeshValueCollection("int", mesh, 1)
g.assign(f)
assert mesh.num_facets() == f.size()
assert mesh.num_cells()*3 == g.size()
for cell in cells(mesh):
for i, facet in enumerate(facets(cell)):
assert 25 == g.get_value(cell.index(), i)
f2 = MeshFunction("int", mesh, g)
for cell in cells(mesh):
for i, facet in enumerate(facets(cell)):
assert f2[facet] == g.get_value(cell.index(), i)
def test_mesh_function_assign_2D_vertices():
mesh = UnitSquareMesh(3, 3)
mesh.init(0)
f = MeshFunction("int", mesh, 0, 25)
g = MeshValueCollection("int", mesh, 0)
g.assign(f)
assert mesh.num_vertices() == f.size()
assert mesh.num_cells()*3 == g.size()
f2 = MeshFunction("int", mesh, g)
for cell in cells(mesh):
for i, vert in enumerate(vertices(cell)):
assert 25 == g.get_value(cell.index(), i)
assert f2[vert] == g.get_value(cell.index(), i)
|