dolfinx.fem
Tools for assembling and manipulating finite element forms.
Functions
|
Create a sparse matrix that is compatible with a given bilinear form. |
|
Given an integral type and a set of entities compute integration entities. |
|
Create a finite element function space. |
Create a sparsity pattern from a bilinear form. |
|
|
Assemble a discrete gradient operator. |
|
Assemble functional. |
Assemble bilinear form into a matrix. |
|
|
|
|
Modify RHS vector b for lifting of Dirichlet boundary conditions. |
|
Create a representation of Dirichlet boundary condition which is imposed on a linear system. |
|
Arrange Dirichlet boundary conditions by the function space that they constrain. |
|
Create a Form or an array of Forms. |
Create a Vector that is compatible with a given linear form |
|
|
Locate degrees-of-freedom geometrically using a marker function. |
|
Locate degrees-of-freedom belonging to mesh entities topologically. |
|
Extract common function spaces from an array of forms. |
|
Generate data needed to interpolate discrete functions across different meshes. |
Create a Lagrange CoordinateElement from element metadata. |
|
|
Return the wrapped C++ class of a variational form of a specific scalar type. |
|
Create a Form object from a data-independent compiled form. |
|
Compile UFL form without associated DOLFINx data. |
|
Insert boundary condition values into vector. |
Classes
|
A constant with respect to a domain. |
|
Create a DOLFINx Expression. |
|
A finite element function that is represented by a function space (domain, element and dofmap) and a vector holding the degrees-of-freedom. |
|
Data for representing a finite element |
|
A space on which Functions (fields) can be defined. |
|
Representation of Dirichlet boundary condition which is imposed on a linear system. |
|
Degree-of-freedom map |
|
A finite element form. |
|
|
|
Coordinate element describing the geometry map for mesh cells. |
- class dolfinx.fem.Constant(domain, c: ndarray | Sequence | floating | complexfloating)[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.CoordinateElement(cmap: CoordinateElement_float32 | CoordinateElement_float64)[source]
Bases:
object
Coordinate element describing the geometry map for mesh cells.
Create a coordinate map element.
Note
This initialiser is for internal use and not normally called in user code. Use
coordinate_element()
to create a CoordinateElement.- Parameters:
cmap – A C++ CoordinateElement.
- create_dof_layout() ElementDofLayout [source]
Compute and return the dof layout
- property degree: int
Polynomial degree of the coordinate element.
- property dim: int
Dimension of the coordinate element space.
This is number of basis functions that span the coordinate space, e.g., for a linear triangle cell the dimension will be 3.
- property dtype: dtype
Scalar type for the coordinate element.
- pull_back(x: ndarray[Any, dtype[float32]] | ndarray[Any, dtype[float64]], cell_geometry: ndarray[Any, dtype[float32]] | ndarray[Any, dtype[float64]]) ndarray[Any, dtype[float32]] | ndarray[Any, dtype[float64]] [source]
Pull points on the physical cell back to the reference cell.
For non-affine cells, the pull-back is a nonlinear operation.
- Parameters:
x – Physical coordinates to pull back to the reference cells,
shape=(num_points, geometrical_dimension)
.cell_geometry – Physical coordinates describing the cell, shape
(num_of_geometry_basis_functions, geometrical_dimension)
They can be created by accessing geometry.x[geometry.dofmap.cell_dofs(i)],
- Returns:
Reference coordinates of the physical points
x
.
- push_forward(X: ndarray[Any, dtype[float32]] | ndarray[Any, dtype[float64]], cell_geometry: ndarray[Any, dtype[float32]] | ndarray[Any, dtype[float64]]) ndarray[Any, dtype[float32]] | ndarray[Any, dtype[float64]] [source]
Push points on the reference cell forward to the physical cell.
- Parameters:
X – Coordinates of points on the reference cell,
shape=(num_points, topological_dimension)
.cell_geometry – Coordinate ‘degrees-of-freedom’ (nodes) of the cell,
shape=(num_geometry_basis_functions, geometrical_dimension)
. Can be created by accessinggeometry.x[geometry.dofmap.cell_dofs(i)]
,
- Returns:
Physical coordinates of the points reference points
X
.
- property variant: int
Lagrange variant of the coordinate element.
Note
The return type is an integer. A Basix enum can be created using
basix.LagrangeVariant(value)
.
- class dolfinx.fem.DirichletBC(bc)[source]
Bases:
object
Representation of Dirichlet boundary condition which is imposed on a linear system.
Note
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.
- dof_indices() tuple[ndarray[Any, dtype[int32]], int] [source]
Access dof indices (local indices, unrolled), including ghosts, to which a Dirichlet condition is applied, and the index to the first non-owned (ghost) index. The array of indices is sorted.
Note
The returned array is read-only.
- Returns:
Sorted array of dof indices (unrolled) and index to the first entry in the dof index array that is not owned. Entries dofs[:pos] are owned and entries dofs[pos:] are ghosts.
- property function_space: FunctionSpace
The function space on which the boundary condition is defined
- set(x: ndarray[Any, dtype[_ScalarType_co]], x0: ndarray[Any, dtype[int32]] | None = None, alpha: float = 1) None [source]
Set entries in an array that are constrained by Dirichlet boundary conditions.
Entries in
x
that are constrained by a Dirichlet boundary conditions are set toalpha * (x_bc - x0)
, wherex_bc
is the (interpolated) boundary condition value.For elements with point-wise evaluated degrees-of-freedom, e.g. Lagrange elements,
x_bc
is the value of the boundary condition at the degree-of-freedom. For elements with moment degrees-of-freedom,x_bc
is the value of the boundary condition interpolated into the finite element space.If x includes ghosted entries (entries available on the calling rank but owned by another rank), ghosted entries constrained by a Dirichlet condition will also be set.
- Parameters:
x – Array to modify for Dirichlet boundary conditions.
x0 – Optional array used in computing the value to set. If not provided it is treated as zero.
alpha – Scaling factor.
- 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 FiniteElement 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.ElementMetaData(family: str, degree: int, shape: Optional[tuple[int, ...]] = None, symmetry: Optional[bool] = None)[source]
Bases:
NamedTuple
Data for representing a finite element
- Parameters:
family – Element type.
degree – Polynomial degree of the element.
shape – Shape for vector/tensor valued elements that are constructed from blocked scalar elements (e.g., Lagrange).
symmetry – Symmetry option for blocked tensor elements.
Create new instance of ElementMetaData(family, degree, shape, symmetry)
- degree: int
Alias for field number 1
- family: str
Alias for field number 0
- shape: tuple[int, ...] | None
Alias for field number 2
- symmetry: bool | None
Alias for field number 3
- class dolfinx.fem.Expression(e: ufl.core.expr.Expr, X: np.ndarray, comm: Optional[_MPI.Comm] = None, form_compiler_options: Optional[dict] = None, jit_options: Optional[dict] = None, dtype: Optional[npt.DTypeLike] = None)[source]
Bases:
object
Create a 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:
e – UFL expression.
X – Array of points of shape
(num_points, tdim)
on the reference element.comm – Communicator that the Expression is defined on.
form_compiler_options – Options used in FFCx compilation of this Expression. Run
ffcx --help
in the commandline to see all available options.jit_options – Options controlling JIT compilation of C code.
Note
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: FunctionSpace | None
The argument function space if expression has argument
- property code: str
C code strings
- property dtype: dtype
- eval(mesh: Mesh, entities: np.ndarray, values: Optional[np.ndarray] = None) np.ndarray [source]
Evaluate Expression on entities.
- Parameters:
mesh – Mesh to evaluate Expression on.
entities – Either an array of cells (index local to process) or an array of integral tuples (cell index, local facet index). The array is flattened.
values – Array to fill with evaluated values. If
None
, storage will be allocated. Otherwise must have shape(num_entities, num_points * value_size * num_all_argument_dofs)
- Returns:
Expression evaluated at points for entities.
- 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.Form(form: Form_complex64 | Form_complex128 | Form_float32 | Form_float64, ufcx_form=None, code: str | None = None, module: ModuleType | None = None)[source]
Bases:
object
A finite element form.
Note
Forms should normally be constructed using
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 form object.
ufcx_form – UFCx form.
code – Form C++ code.
module – CFFI module.
- property code: str | None
C code strings.
- property dtype: dtype
Scalar type 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 mesh: Mesh_float32 | Mesh_float64
Mesh on which this form is defined.
- property module: ModuleType | None
The CFFI module
- property rank: int
- property ufcx_form
The compiled ufcx_form object.
- class dolfinx.fem.Function(*args, **kw)[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. Is not set, the DOLFINx default scalar type is used.
- copy() Function [source]
Create a copy of the Function.
The function space is shared and the degree-of-freedom vector is copied.
- Returns:
A new Function with a copy of the degree-of-freedom vector.
- property dtype: dtype
- eval(x: Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], cells: Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], u=None) ndarray [source]
Evaluate Function at points x.
Points 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(u0: Callable | Expression | Function, cells0: ndarray | None = None, cells1: ndarray | None = None) None [source]
Interpolate an expression.
- Parameters:
u0 – Callable function, Expression or Function to interpolate.
cells0 – Cells in mesh associated with
u0
to interpolate over. IfNone
then all cells are interpolated over.cells1 – Cells in the mesh associated with
self
to interpolate over. IfNone
, then taken to be the same cells ascells0
. Ifcells1
is notNone
, then it must have the same length ascells0
.
- interpolate_nonmatching(u0: Function, cells: ndarray[Any, dtype[int32]], interpolation_data: PointOwnershipData) None [source]
Interpolate a Function defined on one mesh to a function defined on a different mesh.
- Parameters:
u0 – The Function to interpolate.
cells – The cells to interpolate over. If None then all cells are interpolated over.
interpolation_data – Data needed to interpolate functions defined on other meshes. Created by
dolfinx.fem.create_interpolation_data()
.
- property name: str
Name of the Function.
- split() tuple[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.
- Returns:
First level of subspaces of the function space.
- sub(i: int) Function [source]
Return a sub-function (a view into the Function).
Sub-functions are indexed i = 0, …, N-1, where N is the number of sub-spaces.
- Parameters:
i – Index of the sub-function to extract.
- Returns:
A view into the parent Function.
Note
If the sub-Function is re-used, for performance reasons the returned Function should be stored by the caller to avoid repeated re-computation of the subspace.
- class dolfinx.fem.FunctionSpace(mesh: Mesh, element: ufl.FiniteElementBase, cppV: Union[_cpp.fem.FunctionSpace_float32, _cpp.fem.FunctionSpace_float64])[source]
Bases:
FunctionSpace
A space on which Functions (fields) can be defined.
Create a finite element function space.
Note
This initialiser is for internal use and not normally called in user code. Use
functionspace()
to create a function space.- Parameters:
mesh – Mesh that space is defined on.
element – UFL finite element.
cppV – Compiled C++ function space.
- clone() FunctionSpace [source]
Create 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.
- Returns:
A new function space that shares data
- collapse() tuple[FunctionSpace, ndarray] [source]
Collapse a subspace and return a new function space and a map from new to old dofs.
- Returns:
A 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`` if
V
is contained in, or is the same as, this space.
- property element: FiniteElement_float32 | FiniteElement_float64
Function space finite element.
- property num_sub_spaces: int
Number of sub spaces.
- sub(i: int) FunctionSpace [source]
Return the i-th sub space.
- Parameters:
i – Index of the subspace to extract.
- Returns:
A subspace.
Note
If the subspace is re-used, for performance reasons the returned subspace should be stored by the caller to avoid repeated re-computation of the subspace.
- class dolfinx.fem.IntegralType(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)
Bases:
Enum
- cell = 0
- exterior_facet = 1
- interior_facet = 2
- vertex = 3
- dolfinx.fem.apply_lifting(b: ndarray, a: list[Form], bcs: list[list[DirichletBC]], x0: list[ndarray] | None = None, alpha: float = 1, 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[DirichletBC] | None = None, diagonal: float = 1.0, constants=None, coeffs=None, block_mode: BlockMode | None = None)[source]
- dolfinx.fem.assemble_matrix(A: MatrixCSR, a: Form, bcs: list[DirichletBC] | None = None, diagonal: float = 1.0, constants=None, coeffs=None) MatrixCSR
Assemble bilinear form into a matrix.
- Parameters:
a – The bilinear form assemble.
bcs – Boundary conditions that affect the assembled matrix. Degrees-of-freedom constrained by a boundary condition will have their rows/columns zeroed and the value
diagonal
set on on the matrix diagonal.diagonal – Value to set on the matrix diagonal for Dirichlet boundary condition constrained degrees-of-freedom belonging to the same trial and test space.
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.
- block_mode: Block size mode for the returned space matrix. If
None
, default is used.
- Returns:
Matrix representation of the bilinear form
a
.
Note
The returned matrix is not finalised, i.e. ghost values are not accumulated.
- dolfinx.fem.assemble_scalar(M: Form, 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: Form, constants=None, coeffs=None) Vector
- dolfinx.fem.assemble_vector(b: ndarray, L: Form, constants=None, coeffs=None)
- dolfinx.fem.bcs_by_block(spaces: Iterable[FunctionSpace | None], bcs: Iterable[DirichletBC]) list[list[DirichletBC]] [source]
Arrange Dirichlet boundary conditions by the function space that they constrain.
Given a sequence of function spaces
spaces
and a sequence of DirichletBC objectsbcs
, return a list where the ith entry is the list of DirichletBC objects whose space is contained inspace[i]
.
- dolfinx.fem.compile_form(comm: MPI.Intracomm, form: ufl.Form, form_compiler_options: typing.Optional[dict] = {'scalar_type': <class 'numpy.float64'>}, jit_options: typing.Optional[dict] = None) CompiledForm [source]
Compile UFL form without associated DOLFINx data.
- dolfinx.fem.compute_integration_domains(integral_type: IntegralType, topology: Topology, entities: ndarray)[source]
Given an integral type and a set of entities compute integration entities.
This function returns a list [(id, entities)]. For cell integrals entities are the cell indices. For exterior facet integrals, entities is a list of (cell_index, local_facet_index) pairs. For interior facet integrals, entities is a list of (cell_index0, local_facet_index0, cell_index1, local_facet_index1). id refers to the subdomain id used in the definition of the integration measures of the variational form.
Note
Owned mesh entities only are returned. Ghost entities are not included.
Note
For facet integrals, the topology facet-to-cell and cell-to-facet connectivity must be computed before calling this function.
- Parameters:
integral_type – Integral type.
topology – Mesh topology.
entities – List of mesh entities. For
integral_type==IntegralType.cell
,entities
should be cell indices. For otherIntegralType``s, ``entities
should be facet indices.
- Returns:
List of integration entities.
- dolfinx.fem.coordinate_element(celltype: ~dolfinx.cpp.mesh.CellType, degree: int, variant=0, dtype: ~numpy.dtype[~typing.Any] | None | type[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | str | tuple[~typing.Any, int] | tuple[~typing.Any, ~typing.SupportsIndex | ~collections.abc.Sequence[~typing.SupportsIndex]] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | tuple[~typing.Any, ~typing.Any] = <class 'numpy.float64'>)[source]
- dolfinx.fem.coordinate_element(e: FiniteElement)
Create a Lagrange CoordinateElement from element metadata.
Coordinate elements are typically used to create meshes.
- Parameters:
celltype – Cell shape
degree – Polynomial degree of the coordinate element map.
variant – Basix Lagrange variant (affects node placement).
dtype – Scalar type for the coordinate element.
- Returns:
A coordinate element.
- dolfinx.fem.create_form(form: CompiledForm, function_spaces: list[function.FunctionSpace], mesh: Mesh, subdomains: dict[IntegralType, list[tuple[int, np.ndarray]]], coefficient_map: dict[ufl.Function, function.Function], constant_map: dict[ufl.Constant, function.Constant], entity_maps: dict[Mesh, np.typing.NDArray[np.int32]] | None = None) Form [source]
Create a Form object from a data-independent compiled form.
- Parameters:
form – Compiled ufl form
function_spaces – List of function spaces associated with the form. Should match the number of arguments in the form.
mesh – Mesh to associate form with
subdomains – A map from integral type to a list of pairs, where each pair corresponds to a subdomain id and the set of of integration entities to integrate over. Can be computed with {py:func}`dolfinx.fem.compute_integration_domains`.
coefficient_map – Map from UFL coefficient to function with data.
constant_map – Map from UFL constant to constant with data.
entity_map – A map where each key corresponds to a mesh different to the integration domain mesh. The value of the map is an array of integers, where the i-th entry is the entity in the key mesh.
- dolfinx.fem.create_interpolation_data(V_to: FunctionSpace, V_from: FunctionSpace, cells: ndarray[Any, dtype[int32]], padding: float = 1e-14) PointOwnershipData [source]
Generate data needed to interpolate discrete functions across different meshes.
- Parameters:
V_to – Function space to interpolate into
V_from – Function space to interpolate from
cells – Indices of the cells associated with V_to on which to interpolate into.
padding – Absolute padding of bounding boxes of all entities on mesh_to
- Returns:
Data needed to interpolation functions defined on function spaces on the meshes.
- dolfinx.fem.create_matrix(a: Form, block_mode: BlockMode | None = None) MatrixCSR [source]
Create a sparse matrix that is compatible with a given bilinear form.
- Parameters:
a – Bilinear form to assemble.
block_mode – Block mode of the CSR matrix. If
None
, default is used.
- Returns:
Assembled sparse matrix.
- dolfinx.fem.create_sparsity_pattern(a: Form)[source]
Create a sparsity pattern from a bilinear form.
- Parameters:
a – Bilinear form to build a sparsity pattern for.
- Returns:
Sparsity pattern for the form
a
.
Note
The pattern is not finalised, i.e. the caller is responsible for calling
assemble
on the sparsity pattern.
- dolfinx.fem.create_vector(L: Form) Vector [source]
Create a Vector that is compatible with a given linear form
- dolfinx.fem.dirichletbc(value: Union[Function, Constant, np.ndarray], dofs: npt.NDArray[np.int32], V: Optional[dolfinx.fem.FunctionSpace] = None) DirichletBC [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.discrete_gradient(space0: FunctionSpace, space1: FunctionSpace) MatrixCSR [source]
Assemble a discrete gradient operator.
The discrete gradient operator interpolates the gradient of a H1 finite element function into a H(curl) space. It is assumed that the H1 space uses an identity map and the H(curl) space uses a covariant Piola map.
- Parameters:
space0 – H1 space to interpolate the gradient from
space1 – H(curl) space to interpolate into
- Returns:
Discrete gradient operator
- dolfinx.fem.extract_function_spaces(forms: Union[Iterable[Form], Iterable[Iterable[Form]]], 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, typing.Iterable[ufl.Form]], dtype: npt.DTypeLike = <class 'numpy.float64'>, form_compiler_options: typing.Optional[dict] = None, jit_options: typing.Optional[dict] = None, entity_maps: typing.Optional[dict[Mesh, np.typing.NDArray[np.int32]]] = None)[source]
Create a Form or an array of Forms.
- Parameters:
form – A UFL form or list(s) of UFL forms.
dtype – Scalar type to use for the compiled form.
form_compiler_options – See
ffcx_jit
jit_options – See
ffcx_jit
.entity_maps – If any trial functions, test functions, or coefficients in the form are not defined over the same mesh as the integration domain, entity_maps must be supplied. For each key (a mesh, different to the integration domain mesh) a map should be provided relating the entities in the integration domain mesh to the entities in the key mesh e.g. for a key-value pair (msh, emap) in entity_maps, emap[i] is the entity in msh corresponding to entity i in the integration domain mesh.
- Returns:
Compiled finite element Form.
Note
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.form_cpp_class(dtype: dtype[Any] | None | type[Any] | _SupportsDType[dtype[Any]] | str | tuple[Any, int] | tuple[Any, SupportsIndex | Sequence[SupportsIndex]] | list[Any] | _DTypeDict | tuple[Any, Any]) Form_float32 | Form_float64 | Form_complex64 | Form_complex128 [source]
Return the wrapped C++ class of a variational form of a specific scalar type.
- Parameters:
dtype – Scalar type of the required form class.
- Returns:
Wrapped C++ form class of the requested type.
Note
This function is for advanced usage, typically when writing custom kernels using Numba or C.
- dolfinx.fem.functionspace(mesh: Mesh, element: Union[ufl.FiniteElementBase, ElementMetaData, tuple[str, int, tuple, bool]], form_compiler_options: Optional[dict[str, Any]] = None, jit_options: Optional[dict[str, Any]] = None) FunctionSpace [source]
Create a finite element function space.
- Parameters:
mesh – Mesh that space is defined on.
element – Finite element description.
form_compiler_options – Options passed to the form compiler.
jit_options – Options controlling just-in-time compilation.
- Returns:
A function space.
- dolfinx.fem.locate_dofs_geometrical(V: 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: 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[DirichletBC], x0: ndarray | None = None, scale: float = 1) None [source]
Insert boundary condition values into vector.
Note
This function is deprecated.
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.