dolfinx.fem
Tools for assembling and manipulating finite element forms.
Functions
Create a sparsity pattern from a bilinear form. |
- class dolfinx.fem.Constant(domain, c: Union[ndarray, Sequence, float])[source]
Bases:
Constant
A constant with respect to a domain.
- Parameters
domain – DOLFINx or UFL mesh
c – Value of the constant.
- property dtype: dtype
- property value
The value of the constant
- class dolfinx.fem.DirichletBCMetaClass(value: Union[Function, Constant, numpy.ndarray], dofs: numpy.typing.ArrayLike, V: dolfinx.fem.FunctionSpace = None)[source]
Bases:
object
Representation of Dirichlet boundary condition which is imposed on a linear system.
Notes
Dirichlet boundary conditions should normally be constructed using
fem.dirichletbc()
and not using this class initialiser. This class is combined with different base classes that depend on the scalar type of the boundary condition.- Parameters
value – Lifted boundary values function.
dofs – Local indices of degrees of freedom in function space to which boundary condition applies. Expects array of size (number of dofs, 2) if function space of the problem,
V
is passed. Otherwise assumes function space of the problem is the same of function space of boundary values function.V – Function space of a problem to which boundary conditions are applied.
- property function_space: FunctionSpace
The function space on which the boundary condition is defined
- property g
The boundary condition value(s)
- class dolfinx.fem.DofMap(dofmap: DofMap)[source]
Bases:
object
Degree-of-freedom map
This class handles the mapping of degrees of freedom. It builds a dof map based on a ufcx_dofmap on a specific mesh.
- property bs
Returns the block size of the dofmap
- cell_dofs(cell_index: int)[source]
Cell local-global dof map
- Parameters
cell – The cell index
- Returns
Local-global dof map for the cell (using process-local indices)
- property dof_layout
Layout of dofs on an element
- property index_map
Index map that described the parallel distribution of the dofmap
- property index_map_bs
Block size of the index map
- property list
Adjacency list with dof indices for each cell
- class dolfinx.fem.Expression(ufl_expression: ~ufl.core.expr.Expr, X: ~numpy.ndarray, form_compiler_params: dict = {}, jit_params: dict = {}, dtype=<class 'numpy.float64'>)[source]
Bases:
object
Create DOLFINx Expression.
Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell. This class closely follows the concept of a UFC Expression.
This functionality can be used to evaluate a gradient of a Function at the quadrature points in all cells. This evaluated gradient can then be used as input to a non-FEniCS function that calculates a material constitutive model.
- Parameters
ufl_expression – Pure UFL expression
X – Array of points of shape (num_points, tdim) on the reference element.
form_compiler_params – Parameters used in FFCx compilation of this Expression. Run
ffcx --help
in the commandline to see all available options.jit_params – Parameters controlling JIT compilation of C code.
Notes
This wrapper is responsible for the FFCx compilation of the UFL Expr and attaching the correct data to the underlying C++ Expression.
- property argument_function_space: Optional[FunctionSpace]
The argument function space if expression has argument
- property code: str
C code strings
- property dtype: dtype
- eval(cells: ndarray, values: Optional[ndarray] = None) ndarray [source]
Evaluate Expression in cells. Values should have shape (cells.shape[0], num_points * value_size * num_all_argument_dofs). If values is not passed then a new array will be allocated.
- property ufcx_expression
The compiled ufcx_expression object
- property ufl_expression
Original UFL Expression
- property value_size: int
Value size of the expression
- class dolfinx.fem.FormMetaClass(form, V: list[dolfinx.cpp.fem.FunctionSpace], coeffs, constants, subdomains: dict[dolfinx.cpp.mesh.MeshTags_int32, Optional[Any]], mesh: Mesh, ffi, code)[source]
Bases:
object
A finite element form
Notes
Forms should normally be constructed using
forms.form()
and not using this class initialiser. This class is combined with different base classes that depend on the scalar type used in the Form.- Parameters
form – Compiled UFC form
V – The argument function spaces
coeffs – Finite element coefficients that appear in the form
constants – Constants appearing in the form
subdomains – Subdomains for integrals
mesh – The mesh that the form is deined on
- property code: str
C code strings
- property dtype: dtype
dtype of this form
- property function_spaces: List[FunctionSpace]
Function spaces on which this form is defined
- property integral_types
Integral types in the form
- property ufcx_form
The compiled ufcx_form object
- class dolfinx.fem.Function(V: ~dolfinx.fem.function.FunctionSpace, x: ~typing.Optional[~dolfinx.la.VectorMetaClass] = None, name: ~typing.Optional[str] = None, dtype: ~numpy.dtype = <class 'numpy.float64'>)[source]
Bases:
Coefficient
A finite element function that is represented by a function space (domain, element and dofmap) and a vector holding the degrees-of-freedom
Initialize a finite element Function.
- Parameters
V – The function space that the Function is defined on.
x – Function degree-of-freedom vector. Typically required only when reading a saved Function from file.
name – Function name.
dtype – Scalar type.
- copy() Function [source]
Return a copy of the Function. The FunctionSpace is shared and the degree-of-freedom vector is copied.
- property dtype: dtype
- eval(x: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], cells: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], u=None) ndarray [source]
Evaluate Function at points x, where x has shape (num_points, 3), and cells has shape (num_points,) and cell[i] is the index of the cell containing point x[i]. If the cell index is negative the point is ignored.
- property function_space: FunctionSpace
The FunctionSpace that the Function is defined on
- interpolate(u: Union[Callable, Expression, Function], cells: Optional[ndarray] = None) None [source]
Interpolate an expression
- Parameters
u – The function, Expression or Function to interpolate.
cells – The cells to interpolate over. If None then all cells are interpolated over.
- property name: str
Name of the Function.
- split() tuple[dolfinx.fem.function.Function, ...] [source]
Extract any sub functions.
A sub function can be extracted from a discrete function that is in a mixed, vector, or tensor FunctionSpace. The sub function resides in the subspace of the mixed space.
- Parameters
subspaces. (Function space) –
- sub(i: int) Function [source]
Return a sub function.
- Parameters
i – The index of the sub-function to extract.
Note
The sub functions are numbered i = 0..N-1, where N is the total number of sub spaces.
- property vector
PETSc vector holding the degrees-of-freedom.
- property x
Vector holding the degrees-of-freedom.
- class dolfinx.fem.FunctionSpace(mesh: Union[None, Mesh], element: Union[ufl.FiniteElementBase, ElementMetaData, Tuple[str, int]], cppV: Optional[_cpp.fem.FunctionSpace] = None, form_compiler_params: dict[str, Any] = {}, jit_params: dict[str, Any] = {})[source]
Bases:
FunctionSpace
A space on which Functions (fields) can be defined.
Create a finite element function space.
- clone() FunctionSpace [source]
Return a new FunctionSpace \(W\) which shares data with this FunctionSpace \(V\), but with a different unique integer ID.
This function is helpful for defining mixed problems and using blocked linear algebra. For example, a matrix block defined on the spaces \(V \times W\) where, \(V\) and \(W\) are defined on the same finite element and mesh can be identified as an off-diagonal block whereas the \(V \times V\) and \(V \times V\) matrices can be identified as diagonal blocks. This is relevant for the handling of boundary conditions.
- collapse() tuple[dolfinx.fem.function.FunctionSpace, numpy.ndarray] [source]
Collapse a subspace and return a new function space and a map from new to old dofs.
- Returns
The new function space and the map from new to old degrees-of-freedom.
- contains(V) bool [source]
Check if a space is contained in, or is the same as (identity), this space.
- Parameters
V – The space to check to for inclusion.
- Returns
True is
V
is contained in, or is the same as, this space
- property element
- property num_sub_spaces: int
Number of sub spaces
- sub(i: int) FunctionSpace [source]
Return the i-th sub space.
- Parameters
i – The subspace index
- Returns
A subspace
- class dolfinx.fem.IntegralType(self: dolfinx.cpp.fem.IntegralType, value: int)
Bases:
pybind11_object
Members:
cell
exterior_facet
interior_facet
vertex
- cell = <IntegralType.cell: 0>
- exterior_facet = <IntegralType.exterior_facet: 1>
- interior_facet = <IntegralType.interior_facet: 2>
- property name
- property value
- vertex = <IntegralType.vertex: 3>
- dolfinx.fem.TensorFunctionSpace(mesh: Mesh, element: ElementMetaData, shape=None, symmetry: Optional[bool] = None, restriction=None) FunctionSpace [source]
Create tensor finite element (composition of scalar elements) function space.
- dolfinx.fem.VectorFunctionSpace(mesh: Mesh, element: Union[ElementMetaData, Tuple[str, int]], dim=None, restriction=None) FunctionSpace [source]
Create vector finite element (composition of scalar elements) function space.
- dolfinx.fem.apply_lifting(b: ndarray, a: List[FormMetaClass], bcs: List[List[DirichletBCMetaClass]], x0: Optional[List[ndarray]] = None, scale: float = 1.0, constants=None, coeffs=None) None [source]
Modify RHS vector b for lifting of Dirichlet boundary conditions.
It modifies b such that:
\[b \leftarrow b - \text{scale} * A_j (g_j - x0_j)\]where j is a block (nest) index. For a non-blocked problem j = 0. The boundary conditions bcs are on the trial spaces V_j. The forms in [a] must have the same test space as L (from which b was built), but the trial space may differ. If x0 is not supplied, then it is treated as zero.
Note
Ghost contributions are not accumulated (not sent to owner). Caller is responsible for calling VecGhostUpdateBegin/End.
- dolfinx.fem.assemble_matrix(a: Any, bcs: List[DirichletBCMetaClass] = None, diagonal: float = 1.0, constants=None, coeffs=None)[source]
- dolfinx.fem.assemble_matrix(A: MatrixCSRMetaClass, a: Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128], bcs: Optional[List[DirichletBCMetaClass]] = None, diagonal: float = 1.0, constants=None, coeffs=None) MatrixCSRMetaClass
- dolfinx.fem.assemble_matrix(a: Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128], bcs: Optional[List[DirichletBCMetaClass]] = None, diagonal: float = 1.0, constants=None, coeffs=None) MatrixCSRMetaClass
- dolfinx.fem.assemble_scalar(M: FormMetaClass, constants=None, coeffs=None)[source]
Assemble functional. The returned value is local and not accumulated across processes.
- Parameters
M – The functional to compute.
constants – Constants that appear in the form. If not provided, any required constants will be computed.
coeffs – Coefficients that appear in the form. If not provided, any required coefficients will be computed.
- Returns
The computed scalar on the calling rank.
Note
Passing constants and coefficients is a performance optimisation for when a form is assembled multiple times and when (some) constants and coefficients are unchanged.
To compute the functional value on the whole domain, the output of this function is typically summed across all MPI ranks.
- dolfinx.fem.assemble_vector(L: Any, constants=None, coeffs=None)[source]
- dolfinx.fem.assemble_vector(L: Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128], constants=None, coeffs=None) VectorMetaClass
- dolfinx.fem.assemble_vector(b: ndarray, L: FormMetaClass, constants=None, coeffs=None)
- dolfinx.fem.bcs_by_block(spaces: Iterable[Optional[FunctionSpace]], bcs: Iterable[DirichletBCMetaClass]) List[List[DirichletBCMetaClass]] [source]
Arrange Dirichlet boundary conditions by the function space that they constrain.
Given a sequence of function spaces spaces and a sequence of DirichletBC objects bcs, return a list where the ith entry is the list of DirichletBC objects whose space is contained in space[i].
- dolfinx.fem.create_sparsity_pattern(a: FormMetaClass)[source]
Create a sparsity pattern from a bilinear form.
- Parameters
a – The bilinear form to build a sparsity pattern for.
- Returns
Sparsity pattern for the form
a
.
- dolfinx.fem.dirichletbc(value: Union[Function, Constant, np.ndarray], dofs: numpy.typing.NDArray[np.int32], V: dolfinx.fem.FunctionSpace = None) DirichletBCMetaClass [source]
Create a representation of Dirichlet boundary condition which is imposed on a linear system.
- Parameters
value – Lifted boundary values function. It must have a
dtype
property. –
dofs – Local indices of degrees of freedom in function space to which boundary condition applies. Expects array of size (number of dofs, 2) if function space of the problem,
V
, is passed. Otherwise assumes function space of the problem is the same of function space of boundary values function.V – Function space of a problem to which boundary conditions are applied.
- Returns
A representation of the boundary condition for modifying linear systems.
- dolfinx.fem.extract_function_spaces(forms: Union[Iterable[FormMetaClass], Iterable[Iterable[FormMetaClass]]], index: int = 0) Iterable[Union[None, function.FunctionSpace]] [source]
Extract common function spaces from an array of forms. If forms is a list of linear form, this function returns of list of the corresponding test functions. If forms is a 2D array of bilinear forms, for index=0 the list common test function space for each row is returned, and if index=1 the common trial function spaces for each column are returned.
- dolfinx.fem.form(form: ~typing.Union[~ufl.form.Form, ~typing.Iterable[~ufl.form.Form]], dtype: ~numpy.dtype = <class 'numpy.float64'>, form_compiler_params: dict = {}, jit_params: dict = {})[source]
Create a DOLFINx Form or an array of Forms
- Parameters
- Returns
Compiled finite element Form
Notes
This function is responsible for the compilation of a UFL form (using FFCx) and attaching coefficients and domains specific data to the underlying C++ form. It dynamically create a
Form
instance with an appropriate base class for the scalar type, e.g. _cpp.fem.Form_float64.
- dolfinx.fem.locate_dofs_geometrical(V: Union[FunctionSpace, Iterable[FunctionSpace]], marker: Callable) ndarray [source]
Locate degrees-of-freedom geometrically using a marker function.
- Parameters
V – Function space(s) in which to search for degree-of-freedom indices.
marker – A function that takes an array of points
x
with shape(gdim, num_points)
and returns an array of booleans of lengthnum_points
, evaluating toTrue
for entities whose degree-of-freedom should be returned.
- Returns
An array of degree-of-freedom indices (local to the process) for degrees-of-freedom whose coordinate evaluates to True for the marker function.
If
V
is a list of two function spaces, then a 2-D array of shape (number of dofs, 2) is returned.Returned degree-of-freedom indices are unique and ordered by the first column.
- dolfinx.fem.locate_dofs_topological(V: Union[FunctionSpace, Iterable[FunctionSpace]], entity_dim: int, entities: ndarray[Any, dtype[int32]], remote: bool = True) ndarray [source]
Locate degrees-of-freedom belonging to mesh entities topologically.
- Parameters
V – Function space(s) in which to search for degree-of-freedom indices.
entity_dim – Topological dimension of entities where degrees-of-freedom are located.
entities – Indices of mesh entities of dimension
entity_dim
where degrees-of-freedom are located.remote – True to return also “remotely located” degree-of-freedom indices.
- Returns
An array of degree-of-freedom indices (local to the process) for degrees-of-freedom topologically belonging to mesh entities.
If
V
is a list of two function spaces, then a 2-D array of shape (number of dofs, 2) is returned.Returned degree-of-freedom indices are unique and ordered by the first column.
- dolfinx.fem.set_bc(b: ndarray, bcs: List[DirichletBCMetaClass], x0: Optional[ndarray] = None, scale: float = 1.0) None [source]
Insert boundary condition values into vector. Only local (owned) entries are set, hence communication after calling this function is not required unless ghost entries need to be updated to the boundary condition value.