FunctionSpace¶
-
class
dolfin.functions.functionspace.
FunctionSpace
(*args, **kwargs)¶ Bases:
ufl.functionspace.FunctionSpace
,dolfin.cpp.function.FunctionSpace
Base class for all function spaces.
Create finite element function space.
Overloaded versions
FunctionSpace(mesh, element, constrained_domain=None)
Create function space on given mesh for given finite element
- FunctionSpace(mesh, family, degree, form_degree=None,
constrained_domain=None, restriction=None)
Convenience function for specifying element by string, degree, etc.
FunctionSpace(cppV)
Internal library function for wrapping :py:class`cpp.FunctionSpace <dolfin.cpp.FunctionSpace>` into
FunctionSpace
. No copying is done.
- Arguments
- mesh (
Mesh
) - the mesh
- element (
ufl.FiniteElement
) - the element
- family (string)
- specification of the element family, see below for alternatives.
- degree (int)
- the degree of the element.
- form_degree (int)
- the form degree (FEEC notation, used when field is viewed as k-form)
- constrained_domain
- constrained subdomain with map function.
- restriction
- restriction of the element (e.g. to cell facets).
- mesh (
Which families and degrees that are supported is determined by the form compiler used to generate the element, but typical families include
Name Usage Argyris* “ARG” Arnold-Winther* “AW” Brezzi-Douglas-Fortin-Marini* “BDFM” Brezzi-Douglas-Marini “BDM” Bubble “B” Crouzeix-Raviart “CR” Discontinuous Lagrange “DG” Discontinuous Raviart-Thomas “DRT” Hermite* “HER” Lagrange “CG” Mardal-Tai-Winther* “MTW” Morley* “MOR” Nedelec 1st kind H(curl) “N1curl” Nedelec 2nd kind H(curl) “N2curl” Quadrature “Quadrature” Raviart-Thomas “RT” Real “R” *only partly supported.
- Examples of usage
To define a discrete function space over e.g. the unit square:
mesh = UnitSquare(32,32) V = FunctionSpace(mesh, "CG", 1)
Here,
"CG"
stands for Continuous Galerkin, implying the standard Lagrange family of elements. Instead of"CG"
, we could have written"Lagrange"
. With degree 1, we get the linear Lagrange element. Other examples include:# Discontinuous element of degree 0 V = FunctionSpace(mesh, "DG", 0) # Brezzi-Douglas-Marini element of degree 2 W = FunctionSpace(mesh, "BDM", 2) # Real element with one global degree of freedom R = FunctionSpace(mesh, "R", 0)
-
cell
()¶ Return the UFL cell.
-
collapse
(collapsed_dofs=False)¶ Collapse a subspace and return a new function space and a map from new to old dofs
- Arguments
- collapsed_dofs (bool)
- Return the map from new to old dofs
- Returns
- _FunctionSpace_
- The new function space.
- dict
- The map from new to old dofs (optional)
-
dolfin_element
()¶ Return the DOLFIN element.
-
extract_sub_space
(component)¶ Extract subspace for component
- Arguments
- component (numpy.array(uint))
- The component.
- Returns
- _FunctionSpace_
- The subspace.
-
num_sub_spaces
()¶ Return the number of sub spaces
-
split
()¶ Split a mixed functionspace into its sub spaces
-
sub
(i)¶ Return the i-th sub space
-
ufl_cell
()¶ Return the UFL cell.
-
ufl_function_space
()¶