FunctionSpace¶
-
class
dolfin.cpp.function.FunctionSpace(*args)¶ Bases:
dolfin.cpp.common.Variable,dolfin.cpp.function.HierarchicalFunctionSpaceThis class represents a finite element function space defined by a mesh, a finite element, and a local-to-global mapping of the degrees of freedom (dofmap).
Overloaded versions
FunctionSpace(mesh, element, dofmap)
Create function space for given mesh, element and dofmap (shared data)
- Arguments
- mesh (
Mesh) The mesh.
- element (
FiniteElement) The element.
- dofmap (
GenericDofMap) The dofmap.
- mesh (
FunctionSpace(mesh)
Create empty function space for later initialization. This constructor is intended for use by any sub-classes which need to construct objects before the initialisation of the base class. Data can be attached to the base class using FunctionSpace::attach(...).
- Arguments
- mesh (
Mesh) The mesh.
- mesh (
FunctionSpace(V)
Copy constructor
- Arguments
- V (
FunctionSpace) The object to be copied.
- V (
-
assign()¶ Assignment operator
- Arguments
- V (
FunctionSpace) - Another function space.
- V (
-
child()¶ Return the child FunctionSpace in the hierarchy
-
collapse()¶ Overloaded versions
collapse()
Collapse a subspace and return a new function space
- Returns
FunctionSpaceThe new function space.
collapse(collapsed_dofs)
Collapse a subspace and return a new function space and a map from new to old dofs
- Arguments
- collapsed_dofs (std::unordered_map<std::size_t, std::size_t>)
The map from new to old dofs.
- Returns
FunctionSpaceThe new function space.
-
component()¶ - Return component w.r.t. to root superspace, i.e.
- W.sub(1).sub(0) == [1, 0].
- Returns
- numpy.array(int)
- The component (w.r.t to root superspace).
-
contains()¶ Check whether V is subspace of this, or this itself
- Arguments
- V (
FunctionSpace) - The space to be tested for inclusion.
- V (
- Returns
- bool
- True if V is contained or equal to this.
-
dim()¶ Return dimension of function space
- Returns
- int
- The dimension of the function space.
-
dofmap()¶ Return dofmap
- Returns
GenericDofMap- The dofmap.
-
element()¶ Return finite element
- Returns
FiniteElement- The finite element.
-
extract_sub_space()¶ Extract subspace for component
- Arguments
- component (numpy.array(int))
- The component.
- Returns
FunctionSpace- The subspace.
-
has_cell()¶ Check if function space has given cell
- Arguments
- cell (
Cell) - The cell.
- cell (
- Returns
- bool
- True if the function space has the given cell.
-
has_element()¶ Check if function space has given element
- Arguments
- element (
FiniteElement) - The finite element.
- element (
- Returns
- bool
- True if the function space has the given element.
-
interpolate()¶ Interpolate function v into function space, returning the vector of expansion coefficients
- Arguments
- expansion_coefficients (
GenericVector) - The expansion coefficients.
- v (
GenericFunction) - The function to be interpolated.
- expansion_coefficients (
-
leaf_node()¶ Return the finest FunctionSpace in hierarchy
-
mesh()¶ Return mesh
- Returns
Mesh- The mesh.
-
parent()¶ Return the parent FunctionSpace in the hierarchy
-
print_dofmap()¶ Print dofmap (useful for debugging)
-
root_node()¶ Return the coarsest FunctionSpace in hierarchy
-
set_x()¶ Set dof entries in vector to value*x[i], where [x][i] is the coordinate of the dof spatial coordinate. Parallel layout of vector must be consistent with dof map range This function is typically used to construct the null space of a matrix operator, e.g. rigid body rotations.
- Arguments
- vector (
GenericVector) - The vector to set.
- value (float)
- The value to multiply to coordinate by.
- component (int)
- The coordinate index.
- mesh (
Mesh) - The mesh.
- vector (
-
sub()¶ Overloaded versions
sub(component)
Extract subspace for component
- Arguments
- component (int)
Index of the subspace.
- Returns
FunctionSpaceThe subspace.
sub(component)
Extract subspace for component
- Arguments
- component (numpy.array(int))
The component.
- Returns
FunctionSpaceThe subspace.
-
tabulate_dof_coordinates()¶ Tabulate the coordinates of all dofs on this process. This function is typically used by preconditioners that require the spatial coordinates of dofs, for example for re-partitioning or nullspace computations.
- Arguments
- mesh (
Mesh) - The mesh.
- mesh (
- Returns
- numpy.array(float)
- The dof coordinates (x0, y0, x1, y1, . . .)
-
thisown¶ The membership flag