Dual spaces¶
Mathematical background¶
If \(V\) is a vector space, then define the (anti-)dual space \(V^*\) to be the space of bounded conjugate linear functionals \(V \to K\). Therefore, the dual space is the space containing functions from \(V\) to \(K\), where \(K\) can either be \(\mathbb{R}\) or \(\mathbb{C}\).
If \(\{\phi_i\}\) is a basis for a finite vector space \(V\) then there exists \(\{\phi_i^*\}\) a basis for \(V^*\) such that:
The basis \(\{\phi_i^*\}\) is termed the dual basis. Where it is necessary to make the distinction, we will refer to the space to which a dual space is dual as the primal space and its basis as the primal basis.
Since UFL function spaces are finite-dimensional Hilbert spaces which result from the discretisation of infinite-dimensional Hilbert spaces, all of the function spaces with which we are concerned are reflexive, ie \((V^*)^*\) is isomorphic to V under the canonical map. That is, we can identify \((V^*)^*\) and \(V\):
A form defined over an unknown \(a\) in the primal space \(V\) is a known object in the dual space. For example:
with basis coefficients \(h_i = \displaystyle\int_\Omega \phi_i \text{ d}x\).
Dual objects in UFL¶
For an arbitrary FunctionSpace
, V
, the corresponding dual space \(V^*\) can be obtained by calling the dual()
method:
from firedrake import *
mesh = UnitSquareMesh(10, 10)
V = FunctionSpace(mesh, "CG", 1)
V_star = V.dual()
A Coefficient
defines a known function c
in V
. A Function
is a subclass of Coefficient
.
Consequently,
c = Function(V)
f_0 = c * dx
is a symbolic expression for the integral of c
over the domain and represents a scalar value. f_0
is a Python object of type Form
, once assembled, it is a scalar object.
Conversely, Argument
defines a placeholder symbol a
for an unknown function in V
. TestFunction
and TrialFunction
are syntactic sugar for Argument(V, 0)
and Argument(V, 1)
respectively.
a = TrialFunction(V)
f_1 = a * dx
represents the integration of the unknown function a
over the domain. It’s therefore a linear 1-form, or a function in the dual space \(V^* = V \rightarrow K\). f_1
is also a Python object of type Form
. When assembled, it is an object of type Cofunction
:
cf = assemble(f_1) # type Cofunction
cf
is a known object in the dual space, and the dual equivalent of Coefficient
. The more consistent name Cocoefficient
was rejected as confusing and risible. Cofunction
objects can be combined with symbolic Form
objects:
v = TestFunction(V)
a = v * dx
b = assemble(a)
res = a + b
c = assemble(res)
Furthermore, we will want to express unknown objects in the dual space. For example, in order to represent interpolation from a space \(U\) to a space \(V\), it is convienent to reframe this as a problem involving the dual space:
Using the reflexivity of the function space \(U\). This form therefore has two arguments, one in the primal space \(V\) and one in the dual space \(U^*\). Therefore, we need to represent arguments in the dual space - we will call these coarguments. The details of interpolation will be discussed in its own section.
A Coargument
can be constructed by either calling Argument
on a dual space object or calling Coargument
on a dual space.
v = Argument(V, 1) # type Argument
u = Argument(V.dual(), 2) # type Coargument
w = Coargument(V.dual(), 3) # type Coargument
There is a further dual-related type avalilable in UFL. In Cofunction
, we have represented an assembled 1-form. However, commonly we also assemble 2-forms. Matrix
allows an analogous use, and assembled 2-forms can be naturally combined with 2-forms that have not yet been assembled:
mesh = UnitSquareMesh(10,10)
V = FunctionSpace(mesh, "Lagrange", 1)
u = TrialFunction(V)
v = TestFunction(V)
a = (u*v - inner(grad(u),grad(v)) ) * dx
M = assemble(a) # type Matrix
res = assemble(M + a)
Operations supported symbolically, such as the adjoint and action, are also supported on the dual space equivalent.
mesh = UnitSquareMesh(10,10)
V = FunctionSpace(mesh, "Lagrange", 1)
U = FunctionSpace(mesh, "Lagrange", 1)
u = TrialFunction(U)
v = TestFunction(V)
a = u * v * dx
a = assemble(a) # type Matrix
adj = adjoint(a)
b = Matrix(V, U.dual())
u = Coefficient(U)
u_a = Argument(U, 0)
u_form = u_a * dx
primal_action = action(a, u)
dual_action = action(b, u_form)
In summary, this table describes the dual types corresponding to primal finite element spaces, and to known and unknown functions in those spaces:
Primal quantity |
Dual quantity |
---|---|