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
|
"""This module defines the class Argument and a number of related
classes (functions), including TestFunction and TrialFunction."""
# Copyright (C) 2008-2011 Martin Sandve Alnes
#
# This file is part of UFL.
#
# UFL 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.
#
# UFL 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 UFL. If not, see <http://www.gnu.org/licenses/>.
#
# Modified by Anders Logg, 2008-2009.
#
# First added: 2008-03-14
# Last changed: 2011-10-20
from ufl.assertions import ufl_assert
from ufl.common import Counted, product
from ufl.terminal import FormArgument
from ufl.split_functions import split
from ufl.finiteelement import FiniteElementBase
# --- Class representing an argument (basis function) in a form ---
class Argument(FormArgument, Counted):
"""UFL value: Representation of an argument to a form."""
__slots__ = ("_repr", "_element",)
_globalcount = 0
def __init__(self, element, count=None):
FormArgument.__init__(self)
Counted.__init__(self, count, Argument)
ufl_assert(isinstance(element, FiniteElementBase),
"Expecting a FiniteElementBase instance.")
self._element = element
self._repr = "Argument(%r, %r)" % (self._element, self._count)
def reconstruct(self, element=None, count=None):
if element is None or element == self._element:
element = self._element
if count is None or count == self._count:
count = self._count
if count is self._count and element is self._element:
return self
ufl_assert(isinstance(element, FiniteElementBase),
"Expecting an element, not %s" % element)
ufl_assert(isinstance(count, int),
"Expecting an int, not %s" % count)
ufl_assert(element.value_shape() == self._element.value_shape(),
"Cannot reconstruct an Argument with a different value shape.")
return Argument(element, count)
def element(self):
return self._element
def shape(self):
return self._element.value_shape()
def cell(self):
return self._element.cell()
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# TODO: Should in principle do like with Coefficient,
# but that may currently simplify away some arguments
# we want to keep, or?
return False
def __str__(self):
count = str(self._count)
if len(count) == 1:
return "v_%s" % count
else:
return "v_{%s}" % count
def __repr__(self):
return self._repr
def __eq__(self, other):
return isinstance(other, Argument) and self._element == other._element and self._count == other._count
# --- Helper functions for pretty syntax ---
def TestFunction(element):
"""UFL value: Create a test function argument to a form."""
return Argument(element, -2)
def TrialFunction(element):
"""UFL value: Create a trial function argument to a form."""
return Argument(element, -1)
# --- Helper functions for creating subfunctions on mixed elements ---
def Arguments(element):
"""UFL value: Create an Argument in a mixed space, and return a
tuple with the function components corresponding to the subelements."""
return split(Argument(element))
def TestFunctions(element):
"""UFL value: Create a TestFunction in a mixed space, and return a
tuple with the function components corresponding to the subelements."""
return split(TestFunction(element))
def TrialFunctions(element):
"""UFL value: Create a TrialFunction in a mixed space, and return a
tuple with the function components corresponding to the subelements."""
return split(TrialFunction(element))
|