Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add DoF extraction function #2

Open
wants to merge 1 commit into
base: master
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions fenics_helpers/__init__.py
Original file line number Diff line number Diff line change
@@ -1 +1,2 @@
from . import timestepping
from .dof_value_extraction import extract_dof_values
27 changes: 27 additions & 0 deletions fenics_helpers/dof_value_extraction.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,27 @@
import dolfin
from ufl.core.expr import Expr


def _extract_dof_from_function(f):
# Functions have direct access to their DoF values
try:
values = f.vector().get_local()
# subfunctions need to be interpolated
except RuntimeError:
fspace = f.function_space()
try:
fspace = fspace.collapse()
except RuntimeError:
pass
values = dolfin.interpolate(f, fspace).vector().get_local()

return values


def extract_dof_values(expression):
if isinstance(expression, dolfin.Function):
return _extract_dof_from_function(expression)
if isinstance(expression, Expr):
return dolfin.project(expression).vector().get_local()
else:
raise RuntimeError("Can't extract values from {}".format(type(expression)))
43 changes: 43 additions & 0 deletions tests/test_dof_extraction.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,43 @@
from context import fenics_helpers as fh
import unittest
import dolfin
import numpy as np


class TestDofExtraction(unittest.TestCase):
def setUp(self):
mesh = dolfin.UnitIntervalMesh(100)
P1 = dolfin.FiniteElement("P", dolfin.interval, 1)
fs = dolfin.FunctionSpace(mesh, P1)
self.u = dolfin.Function(fs)

mixed_element = dolfin.MixedElement([P1, P1, P1])
mixed_fs = dolfin.FunctionSpace(mesh, mixed_element)
self.mixed_u = dolfin.Function(mixed_fs)

def test_function(self):
u_vec = fh.extract_dof_values(self.u)
self._is_numpy_array_of_size(u_vec, 101)

def test_subspace_function(self):
u1, u2, u3 = self.mixed_u.split()
u2_values = fh.extract_dof_values(u2)
self._is_numpy_array_of_size(u2_values, 101)

def test_expression(self):
my_expr = self.u ** 2 + self.u + dolfin.sin(self.u) + 73.0
my_expr_values = fh.extract_dof_values(my_expr)
self._is_numpy_array_of_size(my_expr_values, 101)

def test_vector_expression(self):
my_expr = dolfin.as_vector([self.u ** 2, self.u, dolfin.sin(self.u)])
my_expr_values = fh.extract_dof_values(my_expr)
self._is_numpy_array_of_size(my_expr_values, 3 * 101)

def _is_numpy_array_of_size(self, array, size):
self.assertIsInstance(array, np.ndarray)
self.assertEqual(array.size, size)


if __name__ == "__main__":
unittest.main()