DOLFINx 0.10.0.0
DOLFINx C++ interface
|
Finite element method functionality. More...
Namespaces | |
namespace | sparsitybuild |
Support for building sparsity patterns from degree-of-freedom maps. | |
Classes | |
struct | BasixElementData |
Basix element holder. More... | |
class | Constant |
Constant value which can be attached to a Form. More... | |
class | CoordinateElement |
class | DirichletBC |
class | DofMap |
Degree-of-freedom map. More... | |
class | ElementDofLayout |
class | Expression |
An Expression represents a mathematical expression evaluated at a pre-defined points on a reference cell. More... | |
class | FiniteElement |
Model of a finite element. More... | |
class | Form |
A representation of finite element variational forms. More... | |
class | Function |
class | FunctionSpace |
This class represents a finite element function space defined by a mesh, a finite element, and a local-to-global map of the degrees-of-freedom. More... | |
struct | integral_data |
Represents integral data, containing the kernel, and a list of entities to integrate over and the indicies of the coefficient functions (relative to the Form) active for this integral. More... | |
Concepts | |
concept | MDSpan |
concept | DofTransformKernel |
DOF transform kernel concept. | |
concept | FEkernel |
Finite element cell kernel concept. | |
concept | MDSpan2 |
Concept for mdspan of rank 1 or 2. | |
Enumerations | |
enum class | doftransform { standard = 0 , transpose = 1 , inverse = 2 , inverse_transpose = 3 } |
DOF transformation type. More... | |
enum class | IntegralType : std::int8_t { cell = 0 , exterior_facet = 1 , interior_facet = 2 , vertex = 3 } |
Type of integral. More... | |
Functions | |
template<dolfinx::scalar T, std::floating_point U> | |
void | tabulate_expression (std::span< T > values, const fem::Expression< T, U > &e, md::mdspan< const T, md::dextents< std::size_t, 2 > > coeffs, std::span< const T > constants, const mesh::Mesh< U > &mesh, fem::MDSpan2 auto entities, std::optional< std::pair< std::reference_wrapper< const FiniteElement< U > >, std::size_t > > element) |
Evaluate an Expression on cells or facets. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | tabulate_expression (std::span< T > values, const fem::Expression< T, U > &e, const mesh::Mesh< U > &mesh, fem::MDSpan2 auto entities) |
Evaluate an Expression on cells or facets. | |
template<dolfinx::scalar T> | |
std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > | make_coefficients_span (const std::map< std::pair< IntegralType, int >, std::pair< std::vector< T >, int > > &coeffs) |
Create a map of std::span s from a map of std::vector s. | |
template<dolfinx::scalar T, std::floating_point U> | |
T | assemble_scalar (const Form< T, U > &M, std::span< const T > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > &coefficients) |
Assemble functional into scalar. | |
template<dolfinx::scalar T, std::floating_point U> | |
T | assemble_scalar (const Form< T, U > &M) |
Assemble functional into scalar. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | assemble_vector (std::span< T > b, const Form< T, U > &L, std::span< const T > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > &coefficients) |
Assemble linear form into a vector. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | assemble_vector (std::span< T > b, const Form< T, U > &L) |
Assemble linear form into a vector. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | apply_lifting (std::span< T > b, std::vector< std::optional< std::reference_wrapper< const Form< T, U > > > > a, const std::vector< std::span< const T > > &constants, const std::vector< std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > > &coeffs, const std::vector< std::vector< std::reference_wrapper< const DirichletBC< T, U > > > > &bcs1, const std::vector< std::span< const T > > &x0, T alpha) |
Modify the right-hand side vector to account for constraints (Dirichlet boundary condition constraints). This modification is known as 'lifting'. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | apply_lifting (std::span< T > b, std::vector< std::optional< std::reference_wrapper< const Form< T, U > > > > a, const std::vector< std::vector< std::reference_wrapper< const DirichletBC< T, U > > > > &bcs1, const std::vector< std::span< const T > > &x0, T alpha) |
Modify the right-hand side vector to account for constraints (Dirichlet boundary conditions constraints). This modification is known as 'lifting'. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | assemble_matrix (la::MatSet< T > auto mat_add, const Form< T, U > &a, std::span< const T > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > &coefficients, std::span< const std::int8_t > dof_marker0, std::span< const std::int8_t > dof_marker1) |
Assemble bilinear form into a matrix. Matrix must already be initialised. Does not zero or finalise the matrix. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | assemble_matrix (auto mat_add, const Form< T, U > &a, std::span< const T > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > &coefficients, const std::vector< std::reference_wrapper< const DirichletBC< T, U > > > &bcs) |
Assemble bilinear form into a matrix. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | assemble_matrix (auto mat_add, const Form< T, U > &a, const std::vector< std::reference_wrapper< const DirichletBC< T, U > > > &bcs) |
Assemble bilinear form into a matrix. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | assemble_matrix (auto mat_add, const Form< T, U > &a, std::span< const std::int8_t > dof_marker0, std::span< const std::int8_t > dof_marker1) |
Assemble bilinear form into a matrix. Matrix must already be initialised. Does not zero or finalise the matrix. | |
template<dolfinx::scalar T> | |
void | set_diagonal (auto set_fn, std::span< const std::int32_t > rows, T diagonal=1.0) |
Sets a value to the diagonal of a matrix for specified rows. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | set_diagonal (auto set_fn, const FunctionSpace< U > &V, const std::vector< std::reference_wrapper< const DirichletBC< T, U > > > &bcs, T diagonal=1.0) |
Sets a value to the diagonal of the matrix for rows with a Dirichlet boundary conditions applied. | |
std::vector< std::int32_t > | locate_dofs_topological (const mesh::Topology &topology, const DofMap &dofmap, int dim, std::span< const std::int32_t > entities, bool remote=true) |
Find degrees-of-freedom which belong to the provided mesh entities (topological). | |
std::array< std::vector< std::int32_t >, 2 > | locate_dofs_topological (const mesh::Topology &topology, std::array< std::reference_wrapper< const DofMap >, 2 > dofmaps, int dim, std::span< const std::int32_t > entities, bool remote=true) |
Find degrees-of-freedom which belong to the provided mesh entities (topological). | |
template<std::floating_point T, typename U> | |
std::vector< std::int32_t > | locate_dofs_geometrical (const FunctionSpace< T > &V, U marker_fn) |
Find degrees of freedom whose geometric coordinate is true for the provided marking function. | |
template<std::floating_point T, typename U> | |
std::array< std::vector< std::int32_t >, 2 > | locate_dofs_geometrical (const std::array< std::reference_wrapper< const FunctionSpace< T > >, 2 > &V, U marker_fn) |
template<std::floating_point T, dolfinx::scalar U = T> | |
void | discrete_curl (const FunctionSpace< T > &V0, const FunctionSpace< T > &V1, la::MatSet< U > auto &&mat_set) |
Assemble a discrete curl operator. | |
template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>> | |
void | discrete_gradient (mesh::Topology &topology, std::pair< std::reference_wrapper< const FiniteElement< U > >, std::reference_wrapper< const DofMap > > V0, std::pair< std::reference_wrapper< const FiniteElement< U > >, std::reference_wrapper< const DofMap > > V1, auto &&mat_set) |
Assemble a discrete gradient operator. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | interpolation_matrix (const FunctionSpace< U > &V0, const FunctionSpace< U > &V1, auto &&mat_set) |
Assemble an interpolation operator matrix. | |
graph::AdjacencyList< std::int32_t > | transpose_dofmap (md::mdspan< const std::int32_t, md::dextents< std::size_t, 2 > > dofmap, std::int32_t num_cells) |
Create an adjacency list that maps a global index (process-wise) to the 'unassembled' cell-wise contributions. | |
std::tuple< common::IndexMap, int, std::vector< std::vector< std::int32_t > > > | build_dofmap_data (MPI_Comm comm, const mesh::Topology &topology, const std::vector< ElementDofLayout > &element_dof_layouts, const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &reorder_fn) |
template<typename U, typename V, typename W> | |
BasixElementData (U element, V bs, W symmetry) -> BasixElementData< typename std::remove_cvref< U >::type::scalar_type > | |
Type deduction helper. | |
template<dolfinx::scalar T> | |
std::array< std::vector< std::shared_ptr< const FunctionSpace< T > > >, 2 > | common_function_spaces (const std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace< T > >, 2 > > > &V) |
Extract FunctionSpaces for (0) rows blocks and (1) columns blocks from a rectangular array of (test, trial) space pairs. | |
template<typename U, typename V, typename W> | |
FunctionSpace (U mesh, V element, W dofmap) -> FunctionSpace< typename std::remove_cvref< typename U::element_type >::type::geometry_type::value_type > | |
Type deduction. | |
template<std::floating_point T> | |
std::vector< T > | interpolation_coords (const fem::FiniteElement< T > &element, const mesh::Geometry< T > &geometry, std::span< const std::int32_t > cells) |
Compute the evaluation points in the physical space at which an expression should be computed to interpolate it in a finite element space. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | interpolate (Function< T, U > &u, std::span< const T > f, std::array< std::size_t, 2 > fshape, std::span< const std::int32_t > cells) |
Interpolate an evaluated expression f(x) in a finite element space. | |
template<std::floating_point T> | |
geometry::PointOwnershipData< T > | create_interpolation_data (const mesh::Geometry< T > &geometry0, const FiniteElement< T > &element0, const mesh::Mesh< T > &mesh1, std::span< const std::int32_t > cells, T padding) |
Generate data needed to interpolate finite element fem::Function's across different meshes. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | interpolate (Function< T, U > &u, const Function< T, U > &v, std::span< const std::int32_t > cells, const geometry::PointOwnershipData< U > &interpolation_data) |
Interpolate a finite element Function defined on a mesh to a finite element Function defined on different (non-matching) mesh. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | interpolate (Function< T, U > &u1, std::span< const std::int32_t > cells1, const Function< T, U > &u0, std::span< const std::int32_t > cells0) |
Interpolate from one finite element Function to another Function on the same (sub)mesh. | |
template<dolfinx::scalar T, std::floating_point U> | |
std::pair< std::vector< T >, int > | allocate_coefficient_storage (const Form< T, U > &form, IntegralType integral_type, int id) |
Allocate storage for coefficients of a pair (integral_type, / id) from a Form. | |
template<dolfinx::scalar T, std::floating_point U> | |
std::map< std::pair< IntegralType, int >, std::pair< std::vector< T >, int > > | allocate_coefficient_storage (const Form< T, U > &form) |
Allocate memory for packed coefficients of a Form. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | pack_coefficients (const Form< T, U > &form, std::map< std::pair< IntegralType, int >, std::pair< std::vector< T >, int > > &coeffs) |
Pack coefficients of a Form. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | pack_coefficients (std::vector< std::reference_wrapper< const Function< T, U > > > coeffs, std::span< const int > offsets, fem::MDSpan2 auto entities, std::span< T > c) |
Pack coefficient data over a list of cells or facets. | |
template<typename T> | |
std::vector< T > | pack_constants (std::vector< std::reference_wrapper< const fem::Constant< T > > > c) |
Pack constants of an Expression or Form into a single array ready for assembly. | |
template<typename U> requires std::convertible_to< U, fem::Expression<typename std::decay_t<U>::scalar_type, typename std::decay_t<U>::geometry_type>> or std::convertible_to< U, fem::Form<typename std::decay_t<U>::scalar_type, typename std::decay_t<U>::geometry_type>> | |
std::vector< typename U::scalar_type > | pack_constants (const U &u) |
Pack constants of an Expression or Form into a single array ready for assembly. | |
std::vector< std::int32_t > | compute_integration_domains (IntegralType integral_type, const mesh::Topology &topology, std::span< const std::int32_t > entities) |
Given an integral type and a set of entities, computes and return data for the entities that should be integrated over. | |
template<dolfinx::scalar T, std::floating_point U> | |
std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace< U > >, 2 > > > | extract_function_spaces (const std::vector< std::vector< const Form< T, U > * > > &a) |
Extract test (0) and trial (1) function spaces pairs for each bilinear form for a rectangular array of forms. | |
template<dolfinx::scalar T, std::floating_point U> | |
la::SparsityPattern | create_sparsity_pattern (const Form< T, U > &a) |
Create a sparsity pattern for a given form. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | build_sparsity_pattern (la::SparsityPattern &pattern, const Form< T, U > &a) |
Build a sparsity pattern for a given form. | |
template<std::floating_point T> | |
ElementDofLayout | create_element_dof_layout (const fem::FiniteElement< T > &element, const std::vector< int > &parent_map={}) |
Create an ElementDofLayout from a FiniteElement. | |
DofMap | create_dofmap (MPI_Comm comm, const ElementDofLayout &layout, mesh::Topology &topology, std::function< void(std::span< std::int32_t >, std::uint32_t)> permute_inv, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> reorder_fn) |
Create a dof map on mesh. | |
std::vector< DofMap > | create_dofmaps (MPI_Comm comm, const std::vector< ElementDofLayout > &layouts, mesh::Topology &topology, std::function< void(std::span< std::int32_t >, std::uint32_t)> permute_inv, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> reorder_fn) |
Create a set of dofmaps on a given topology. | |
std::vector< std::string > | get_coefficient_names (const ufcx_form &ufcx_form) |
std::vector< std::string > | get_constant_names (const ufcx_form &ufcx_form) |
Get the name of each constant in a UFC form. | |
template<dolfinx::scalar T, std::floating_point U = scalar_value_t<T>> | |
Form< T, U > | create_form_factory (const std::vector< std::reference_wrapper< const ufcx_form > > &ufcx_forms, const std::vector< std::shared_ptr< const FunctionSpace< U > > > &spaces, const std::vector< std::shared_ptr< const Function< T, U > > > &coefficients, const std::vector< std::shared_ptr< const Constant< T > > > &constants, const std::map< IntegralType, std::vector< std::pair< std::int32_t, std::span< const std::int32_t > > > > &subdomains, const std::map< std::shared_ptr< const mesh::Mesh< U > >, std::span< const std::int32_t > > &entity_maps, std::shared_ptr< const mesh::Mesh< U > > mesh=nullptr) |
Create a Form from UFCx input with coefficients and constants passed in the required order. | |
template<dolfinx::scalar T, std::floating_point U = scalar_value_t<T>> | |
Form< T, U > | create_form (const ufcx_form &ufcx_form, const std::vector< std::shared_ptr< const FunctionSpace< U > > > &spaces, const std::map< std::string, std::shared_ptr< const Function< T, U > > > &coefficients, const std::map< std::string, std::shared_ptr< const Constant< T > > > &constants, const std::map< IntegralType, std::vector< std::pair< std::int32_t, std::span< const std::int32_t > > > > &subdomains, const std::map< std::shared_ptr< const mesh::Mesh< U > >, std::span< const std::int32_t > > &entity_maps, std::shared_ptr< const mesh::Mesh< U > > mesh=nullptr) |
Create a Form from UFC input with coefficients and constants resolved by name. | |
template<dolfinx::scalar T, std::floating_point U = scalar_value_t<T>> | |
Form< T, U > | create_form (ufcx_form *(*fptr)(), const std::vector< std::shared_ptr< const FunctionSpace< U > > > &spaces, const std::map< std::string, std::shared_ptr< const Function< T, U > > > &coefficients, const std::map< std::string, std::shared_ptr< const Constant< T > > > &constants, const std::map< IntegralType, std::vector< std::pair< std::int32_t, std::span< const std::int32_t > > > > &subdomains, const std::map< std::shared_ptr< const mesh::Mesh< U > >, std::span< const std::int32_t > > &entity_maps, std::shared_ptr< const mesh::Mesh< U > > mesh=nullptr) |
Create a Form using a factory function that returns a pointer to a ufcx_form . | |
template<std::floating_point T> | |
FunctionSpace< T > | create_functionspace (std::shared_ptr< mesh::Mesh< T > > mesh, std::shared_ptr< const fem::FiniteElement< T > > e, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> reorder_fn=nullptr) |
NEW Create a function space from a fem::FiniteElement. | |
template<dolfinx::scalar T, std::floating_point U = scalar_value_t<T>> | |
Expression< T, U > | create_expression (const ufcx_expression &e, const std::vector< std::shared_ptr< const Function< T, U > > > &coefficients, const std::vector< std::shared_ptr< const Constant< T > > > &constants, std::shared_ptr< const FunctionSpace< U > > argument_space=nullptr) |
Create Expression from UFC. | |
template<dolfinx::scalar T, std::floating_point U = scalar_value_t<T>> | |
Expression< T, U > | create_expression (const ufcx_expression &e, const std::map< std::string, std::shared_ptr< const Function< T, U > > > &coefficients, const std::map< std::string, std::shared_ptr< const Constant< T > > > &constants, std::shared_ptr< const FunctionSpace< U > > argument_space=nullptr) |
Create Expression from UFC input (with named coefficients and constants). | |
Finite element method functionality.
Classes and algorithms for finite element method spaces and operations.
|
strong |
|
strong |
std::map< std::pair< IntegralType, int >, std::pair< std::vector< T >, int > > allocate_coefficient_storage | ( | const Form< T, U > & | form | ) |
std::pair< std::vector< T >, int > allocate_coefficient_storage | ( | const Form< T, U > & | form, |
IntegralType | integral_type, | ||
int | id ) |
void apply_lifting | ( | std::span< T > | b, |
std::vector< std::optional< std::reference_wrapper< const Form< T, U > > > > | a, | ||
const std::vector< std::span< const T > > & | constants, | ||
const std::vector< std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > > & | coeffs, | ||
const std::vector< std::vector< std::reference_wrapper< const DirichletBC< T, U > > > > & | bcs1, | ||
const std::vector< std::span< const T > > & | x0, | ||
T | alpha ) |
Modify the right-hand side vector to account for constraints (Dirichlet boundary condition constraints). This modification is known as 'lifting'.
Consider the discrete algebraic system
\[\begin{bmatrix} A_{0} & A_{1} \end{bmatrix} \begin{bmatrix} u_{0} \\ u_{1} \end{bmatrix} = b, \]
where \(A_{i}\) is a matrix. Partitioning each vector \(u_{i}\) into 'unknown' ( \(u_{i}^{(0)}\)) and prescribed ( \(u_{i}^{(1)}\)) groups,
\[\begin{bmatrix} A_{0}^{(0)} & A_{0}^{(1)} & A_{1}^{(0)} & A_{1}^{(1)} \end{bmatrix} \begin{bmatrix} u_{0}^{(0)} \\ u_{0}^{(1)} \\ u_{1}^{(0)} \\ u_{1}^{(1)} \end{bmatrix} = b. \]
If \(u_{i}^{(1)} = \alpha(g_{i} - x_{i})\), where \(g_{i}\) is the Dirichlet boundary condition value, \(x_{i}\) is provided and \(\alpha\) is a constant, then
\[\begin{bmatrix} A_{0}^{(0)} & A_{0}^{(1)} & A_{1}^{(0)} & A_{1}^{(1)} \end{bmatrix} \begin{bmatrix} u_{0}^{(0)} \\ \alpha(g_{0} - x_{0}) \\ u_{1}^{(0)} \\ \alpha(g_{1} - x_{1}) \end{bmatrix} = b. \]
Rearranging,
\[\begin{bmatrix} A_{0}^{(0)} & A_{1}^{(0)} \end{bmatrix} \begin{bmatrix} u_{0}^{(0)} \\ u_{1}^{(0)} \end{bmatrix} = b - \alpha A_{0}^{(1)} (g_{0} - x_{0}) - \alpha A_{1}^{(1)} (g_{1} - x_{1}). \]
The modified \(b\) vector is
\[ b \leftarrow b - \alpha A_{0}^{(1)} (g_{0} - x_{0}) - \alpha A_{1}^{(1)} (g_{1} - x_{1}) \]
More generally,
\[ b \leftarrow b - \alpha A_{i}^{(1)} (g_{i} - x_{i}). \]
b
by this function. Use DirichletBC::set to set values in b
.[in,out] | b | The vector to modify inplace. |
[in] | a | List of bilinear forms, where a[i] is the form that generates the matrix \(A_{i}\). All forms in a must share the same test function space. The trial function spaces can differ. |
[in] | constants | Constant data appearing in the forms a . |
[in] | coeffs | Coefficient data appearing in the forms a . |
[in] | x0 | The vector \(x_{i}\) above. If empty it is set to zero. |
[in] | bcs1 | Boundary conditions that provide the \(g_{i}\) values. bcs1[i] is the list of boundary conditions on \(u_{i}\). |
[in] | alpha | Scalar used in the modification of b . |
void apply_lifting | ( | std::span< T > | b, |
std::vector< std::optional< std::reference_wrapper< const Form< T, U > > > > | a, | ||
const std::vector< std::vector< std::reference_wrapper< const DirichletBC< T, U > > > > & | bcs1, | ||
const std::vector< std::span< const T > > & | x0, | ||
T | alpha ) |
Modify the right-hand side vector to account for constraints (Dirichlet boundary conditions constraints). This modification is known as 'lifting'.
See apply_lifting() for a detailed explanation of the lifting. The difference between this function and apply_lifting() is that apply_lifting() requires packed form constant and coefficient data to be passed to the function, whereas this function packs the constant and coefficient form data and then calls apply_lifting().
b
by this function. Use DirichletBC::set to set values in b
.[in,out] | b | The vector to modify inplace. |
[in] | a | List of bilinear forms, where a[i] is the form that generates the matrix \(A_{i}\) (see apply_lifting()). All forms in a must share the same test function space. The trial function spaces can differ. |
[in] | x0 | The vector \(x_{i}\) described in apply_lifting(). If empty it is set to zero. |
[in] | bcs1 | Boundary conditions that provide the \(g_{i}\) values described in apply_lifting(). bcs1[i] is the list of boundary conditions on \(u_{i}\). |
[in] | alpha | Scalar used in the modification of b (see described in apply_lifting()). |
void assemble_matrix | ( | auto | mat_add, |
const Form< T, U > & | a, | ||
const std::vector< std::reference_wrapper< const DirichletBC< T, U > > > & | bcs ) |
Assemble bilinear form into a matrix.
[in] | mat_add | The function for adding values into the matrix. |
[in] | a | The bilinear from to assemble. |
[in] | bcs | Boundary conditions to apply. For boundary condition dofs the row and column are zeroed. The diagonal entry is not set. |
void assemble_matrix | ( | auto | mat_add, |
const Form< T, U > & | a, | ||
std::span< const std::int8_t > | dof_marker0, | ||
std::span< const std::int8_t > | dof_marker1 ) |
Assemble bilinear form into a matrix. Matrix must already be initialised. Does not zero or finalise the matrix.
[in] | mat_add | The function for adding values into the matrix. |
[in] | a | The bilinear form to assemble. |
[in] | dof_marker0 | Boundary condition markers for the rows. If bc[i] is true then rows i in Awill be zeroed. The index i is a local index. |
[in] | dof_marker1 | Boundary condition markers for the columns. If bc[i] is true then rows i in A will be zeroed. The index i is a local index. |
void assemble_matrix | ( | auto | mat_add, |
const Form< T, U > & | a, | ||
std::span< const T > | constants, | ||
const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > & | coefficients, | ||
const std::vector< std::reference_wrapper< const DirichletBC< T, U > > > & | bcs ) |
Assemble bilinear form into a matrix.
[in] | mat_add | The function for adding values into the matrix. |
[in] | a | The bilinear from to assemble. |
[in] | constants | Constants that appear in a . |
[in] | coefficients | Coefficients that appear in a . |
[in] | bcs | Boundary conditions to apply. For boundary condition dofs the row and column are zeroed. The diagonal entry is not set. |
void assemble_matrix | ( | la::MatSet< T > auto | mat_add, |
const Form< T, U > & | a, | ||
std::span< const T > | constants, | ||
const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > & | coefficients, | ||
std::span< const std::int8_t > | dof_marker0, | ||
std::span< const std::int8_t > | dof_marker1 ) |
Assemble bilinear form into a matrix. Matrix must already be initialised. Does not zero or finalise the matrix.
[in] | mat_add | The function for adding values into the matrix. |
[in] | a | The bilinear form to assemble. |
[in] | constants | Constants that appear in a . |
[in] | coefficients | Coefficients that appear in a . |
[in] | dof_marker0 | Boundary condition markers for the rows. If bc[i] is true then rows i in A will be zeroed. The index i is a local index. |
[in] | dof_marker1 | Boundary condition markers for the columns. If bc[i] is true then rows i in A will be zeroed. The index i is a local index. |
T assemble_scalar | ( | const Form< T, U > & | M | ) |
Assemble functional into scalar.
[in] | M | The form (functional) to assemble. |
T assemble_scalar | ( | const Form< T, U > & | M, |
std::span< const T > | constants, | ||
const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > & | coefficients ) |
Assemble functional into scalar.
The caller supplies the form constants and coefficients for this version, which has efficiency benefits if the data can be re-used for multiple calls.
[in] | M | The form (functional) to assemble |
[in] | constants | The constants that appear in M |
[in] | coefficients | The coefficients that appear in M |
void assemble_vector | ( | std::span< T > | b, |
const Form< T, U > & | L ) |
Assemble linear form into a vector.
[in,out] | b | Vector to be assembled. It will not be zeroed before assembly. |
[in] | L | Linear forms to assemble into b. |
void assemble_vector | ( | std::span< T > | b, |
const Form< T, U > & | L, | ||
std::span< const T > | constants, | ||
const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > & | coefficients ) |
Assemble linear form into a vector.
The caller supplies the form constants and coefficients for this version, which has efficiency benefits if the data can be re-used for multiple calls.
[in,out] | b | The vector to be assembled. It will not be zeroed before assembly. |
[in] | L | The linear forms to assemble into b. |
[in] | constants | The constants that appear in L . |
[in] | coefficients | The coefficients that appear in L . |
std::tuple< common::IndexMap, int, std::vector< std::vector< std::int32_t > > > build_dofmap_data | ( | MPI_Comm | comm, |
const mesh::Topology & | topology, | ||
const std::vector< ElementDofLayout > & | element_dof_layouts, | ||
const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> & | reorder_fn ) |
Build dofmap data for elements on a mesh topology
[in] | comm | MPI communicator |
[in] | topology | The mesh topology |
[in] | element_dof_layouts | The element dof layouts for each cell type in topology . |
[in] | reorder_fn | Graph reordering function that is applied to the dofmaps |
void build_sparsity_pattern | ( | la::SparsityPattern & | pattern, |
const Form< T, U > & | a ) |
Build a sparsity pattern for a given form.
[in] | pattern | The sparsity pattern to add to |
[in] | a | A bilinear form |
std::array< std::vector< std::shared_ptr< const FunctionSpace< T > > >, 2 > common_function_spaces | ( | const std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace< T > >, 2 > > > & | V | ) |
Extract FunctionSpaces for (0) rows blocks and (1) columns blocks from a rectangular array of (test, trial) space pairs.
The test space must be the same for each row and the trial spaces must be the same for each column. Raises an exception if there is an inconsistency. e.g. if each form in row i does not have the same test space then an exception is raised.
[in] | V | Vector function spaces for (0) each row block and (1) each column block |
std::vector< std::int32_t > compute_integration_domains | ( | fem::IntegralType | integral_type, |
const mesh::Topology & | topology, | ||
std::span< const std::int32_t > | entities ) |
Given an integral type and a set of entities, computes and return data for the entities that should be integrated over.
This function returns a list data, for each entity in entities
, that is used in assembly. For cell integrals it is simply the cell cell indices. For exterior facet integrals, a list of (cell_index, / local_facet_index)
pairs is returned. For interior facet integrals, a list of (cell_index0, local_facet_index0, cell_index1, / local_facet_index1)
tuples is returned. The data computed by this function is typically used as input to fem::create_form.
[in] | integral_type | Integral type. |
[in] | topology | Mesh topology. |
[in] | entities | List of mesh entities. For integral_type==IntegralType::cell , entities should be cell indices. For other IntegralType , entities should be facet indices. |
fem::DofMap create_dofmap | ( | MPI_Comm | comm, |
const ElementDofLayout & | layout, | ||
mesh::Topology & | topology, | ||
std::function< void(std::span< std::int32_t >, std::uint32_t)> | permute_inv, | ||
std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> | reorder_fn ) |
std::vector< fem::DofMap > create_dofmaps | ( | MPI_Comm | comm, |
const std::vector< ElementDofLayout > & | layouts, | ||
mesh::Topology & | topology, | ||
std::function< void(std::span< std::int32_t >, std::uint32_t)> | permute_inv, | ||
std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> | reorder_fn ) |
Create a set of dofmaps on a given topology.
[in] | comm | MPI communicator |
[in] | layouts | Dof layout on each element type |
[in] | topology | Mesh topology |
[in] | permute_inv | Function to un-permute dofs. nullptr when transformation is not required. |
[in] | reorder_fn | Graph reordering function called on the dofmaps |
Form< T, U > create_form | ( | const ufcx_form & | ufcx_form, |
const std::vector< std::shared_ptr< const FunctionSpace< U > > > & | spaces, | ||
const std::map< std::string, std::shared_ptr< const Function< T, U > > > & | coefficients, | ||
const std::map< std::string, std::shared_ptr< const Constant< T > > > & | constants, | ||
const std::map< IntegralType, std::vector< std::pair< std::int32_t, std::span< const std::int32_t > > > > & | subdomains, | ||
const std::map< std::shared_ptr< const mesh::Mesh< U > >, std::span< const std::int32_t > > & | entity_maps, | ||
std::shared_ptr< const mesh::Mesh< U > > | mesh = nullptr ) |
Create a Form from UFC input with coefficients and constants resolved by name.
[in] | ufcx_form | UFC form |
[in] | spaces | Function spaces for the Form arguments. |
[in] | coefficients | Coefficient fields in the form (by name). |
[in] | constants | Spatial constants in the form (by name). |
[in] | subdomains | Subdomain markers. The data can be computed using fem::compute_integration_domains. |
subdomains
must be sorted by domain id. [in] | entity_maps | The entity maps for the form. Empty for single domain problems. |
[in] | mesh | Mesh of the domain. This is required if the form has no arguments, e.g. a functional. |
Form< T, U > create_form | ( | ufcx_form *(* | fptr )(), |
const std::vector< std::shared_ptr< const FunctionSpace< U > > > & | spaces, | ||
const std::map< std::string, std::shared_ptr< const Function< T, U > > > & | coefficients, | ||
const std::map< std::string, std::shared_ptr< const Constant< T > > > & | constants, | ||
const std::map< IntegralType, std::vector< std::pair< std::int32_t, std::span< const std::int32_t > > > > & | subdomains, | ||
const std::map< std::shared_ptr< const mesh::Mesh< U > >, std::span< const std::int32_t > > & | entity_maps, | ||
std::shared_ptr< const mesh::Mesh< U > > | mesh = nullptr ) |
Create a Form using a factory function that returns a pointer to a ufcx_form
.
Coefficients and constants are resolved by name/string.
[in] | fptr | Pointer to a function returning a pointer to ufcx_form. |
[in] | spaces | Function spaces for the Form arguments. |
[in] | coefficients | Coefficient fields in the form (by name), |
[in] | constants | Spatial constants in the form (by name), |
[in] | subdomains | Subdomain markers. The data can be computed using fem::compute_integration_domains. |
subdomains
must be sorted by domain id. [in] | entity_maps | The entity maps for the form. Empty for single domain problems. |
[in] | mesh | Mesh of the domain. This is required if the form has no arguments, e.g. a functional. |
Form< T, U > create_form_factory | ( | const std::vector< std::reference_wrapper< const ufcx_form > > & | ufcx_forms, |
const std::vector< std::shared_ptr< const FunctionSpace< U > > > & | spaces, | ||
const std::vector< std::shared_ptr< const Function< T, U > > > & | coefficients, | ||
const std::vector< std::shared_ptr< const Constant< T > > > & | constants, | ||
const std::map< IntegralType, std::vector< std::pair< std::int32_t, std::span< const std::int32_t > > > > & | subdomains, | ||
const std::map< std::shared_ptr< const mesh::Mesh< U > >, std::span< const std::int32_t > > & | entity_maps, | ||
std::shared_ptr< const mesh::Mesh< U > > | mesh = nullptr ) |
Create a Form from UFCx input with coefficients and constants passed in the required order.
Use fem::create_form to create a fem::Form with coefficients and constants associated with the name/string.
[in] | ufcx_forms | A list of UFCx forms, one for each cell type. |
[in] | spaces | Vector of function spaces. The number of spaces is equal to the rank of the form. |
[in] | coefficients | Coefficient fields in the form. |
[in] | constants | Spatial constants in the form. |
[in] | subdomains | Subdomain markers. The data can be computed using fem::compute_integration_domains. |
[in] | entity_maps | The entity maps for the form. Empty for single domain problems. |
[in] | mesh | The mesh of the domain. |
subdomains
must be sorted by domain id. geometry::PointOwnershipData< T > create_interpolation_data | ( | const mesh::Geometry< T > & | geometry0, |
const FiniteElement< T > & | element0, | ||
const mesh::Mesh< T > & | mesh1, | ||
std::span< const std::int32_t > | cells, | ||
T | padding ) |
Generate data needed to interpolate finite element fem::Function's across different meshes.
[in] | geometry0 | Mesh geometry of the space to interpolate into. |
[in] | element0 | Element of the space to interpolate into. |
[in] | mesh1 | Mesh of the function to interpolate from. |
[in] | cells | Indices of the cells in the destination mesh on which to interpolate. Should be the same as the list used when calling interpolation_coords . |
[in] | padding | Absolute padding of bounding boxes of all entities on mesh1 . This is used avoid floating point issues when an interpolation point from mesh0 is on the surface of a cell in mesh1 . This parameter can also be used for extrapolation, i.e. if cells in mesh0 is not overlapped by mesh1 . |
padding
to a large value will increase the runtime of this function, as one has to determine what entity is closest if there is no intersection. la::SparsityPattern create_sparsity_pattern | ( | const Form< T, U > & | a | ) |
Create a sparsity pattern for a given form.
[in] | a | A bilinear form |
void discrete_curl | ( | const FunctionSpace< T > & | V0, |
const FunctionSpace< T > & | V1, | ||
la::MatSet< U > auto && | mat_set ) |
Assemble a discrete curl operator.
For vector-valued finite functions \(u \in V_{0} \) and \(v \in V_{1}\), consider the interpolation of the curl of \(u\) in the space \(V_{1}\), i.e. \(\Pi_{V_{1}}: \nabla \times u \rightarrow v\), where \(\Pi_{V_{1}}\) is the interpolation operator associated with \(V_{1}\). This interpolation of \(\nabla \times u\) into \(V_{1}\) is properly posed and exact for specific choices of function spaces. If \(V_{0}\) is a Nédélec ( \(H({\rm curl})\)) space of degree \(k > 1\) and \(V_{1}\) is a Raviart-Thomas ( \(H({\rm div})\)) space of degree of at least \(k - 1\), then the interpolation is exact.
The implementation of this function exploits the result:
\[ \hat{\nabla} \times \psi_{C}(\boldsymbol{u}) = \psi_{D}(\nabla \times \boldsymbol{u}), \]
where \(\psi_{C}\) is the covariant pull-back (to the reference cell) and \(\psi_{D}\) is the contravariant pull-back. See Ern and Guermond (2021), Finite Elements I, Springer Nature, https://doi.org/10.1007/978-3-030-56341-7 [Corollary 9.9 (Commuting properties)]. Consequently, the spaces V0
and V1
must use covariant and contravariant maps, respectively.
This function builds a matrix \(C\) (the 'discrete curl'), which when applied to the degrees-of-freedom of \(u\) gives the degrees-of-freedom of \(v\) such that \(v = \nabla \times u\). If the finite element degree-of-freedom vectors associated with \(u\) and \(v\) are \(a\) and \(b\), respectively, then \(b = C a\), which yields \(v = \Pi_{V} \nabla \times u\). It essentially maps that curl of a function in a degree \(k > 1\) Nédélec space into a degree \(k - 1\) Raviart-Thomas space.
The discerete curl is typically used in constructing algebraic multigrid preconditioners for \(H({\rm div})\) problems, e.g. when using the Hypre Auxiliary-space Divergence Solver (ADS) to solve a mixed Poisson in three-dimension.
V0
and V1
must be vector-valued, in three spatial dimensions, and use covariant and contravariant maps, respectively.T | Scalar type of the mesh and elements. |
U | Scalar type of the matrix being inserted into. This is usually the same as T , but may differ for matrix backends that support only a specific type, e.g. PETSc which supports only one scalar type for a build of PETSc. |
[in] | V0 | Space that \(u\) is from. It must be a covariant Piola mapped element. It is normally an \(H({\rm curl})\)-conforming Nédélec space. |
[in] | V1 | Space that \(v\) is from. It must be a contravariant Piola mapped element. It is normally an \(H({\rm
div})\)-conforming Raviart-Thomas space of one degree lower than V0 . |
[in] | mat_set | A functor that sets (not add) values in a matrix \(C\). |
void discrete_gradient | ( | mesh::Topology & | topology, |
std::pair< std::reference_wrapper< const FiniteElement< U > >, std::reference_wrapper< const DofMap > > | V0, | ||
std::pair< std::reference_wrapper< const FiniteElement< U > >, std::reference_wrapper< const DofMap > > | V1, | ||
auto && | mat_set ) |
Assemble a discrete gradient operator.
The discrete gradient operator \(A\) interpolates the gradient of a Lagrange finite element function in \(V_0 \subset H^1\) into a Nédélec (first kind) space \(V_1 \subset H({\rm curl})\), i.e. \(\nabla V_0 \rightarrow V_1\). If \(u_0\) is the degree-of-freedom vector associated with \(V_0\), then \(u_1=Au_0\) where \(u_1\) is the degrees-of-freedom vector for interpolating function in the \(H({\rm curl})\) space. An example of where discrete gradient operators are used is the creation of algebraic multigrid solvers for \(H({\rm curl})\) and \(H({\rm div})\) problems.
V1
should be used for the rows of the sparsity pattern, V0
for the columns.[in] | topology | Mesh topology |
[in] | V0 | Lagrange element and dofmap for corresponding space to interpolate the gradient from. |
[in] | V1 | Nédélec (first kind) element and dofmap for corresponding space to interpolate into. |
[in] | mat_set | A functor that sets values in a matrix |
std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace< U > >, 2 > > > extract_function_spaces | ( | const std::vector< std::vector< const Form< T, U > * > > & | a | ) |
Extract test (0) and trial (1) function spaces pairs for each bilinear form for a rectangular array of forms.
[in] | a | A rectangular block on bilinear forms. |
a
with a pair of function spaces in each array entry. If a form is null, then the returned function space pair is (null, null). std::vector< std::string > get_coefficient_names | ( | const ufcx_form & | ufcx_form | ) |
Get the name of each coefficient in a UFC form
[in] | ufcx_form | The UFC form |
std::vector< std::string > get_constant_names | ( | const ufcx_form & | ufcx_form | ) |
Get the name of each constant in a UFC form.
[in] | ufcx_form | The UFC form |
void interpolate | ( | Function< T, U > & | u, |
const Function< T, U > & | v, | ||
std::span< const std::int32_t > | cells, | ||
const geometry::PointOwnershipData< U > & | interpolation_data ) |
Interpolate a finite element Function defined on a mesh to a finite element Function defined on different (non-matching) mesh.
T | Function scalar type. |
U | mesh::Mesh geometry scalar type. |
u | Function to interpolate into. |
v | Function to interpolate from. |
cells | Cells indices relative to the mesh associated with u that will be interpolated into. |
interpolation_data | Data required for associating the interpolation points of u with cells in v . This is computed by fem::create_interpolation_data. |
void interpolate | ( | Function< T, U > & | u, |
std::span< const T > | f, | ||
std::array< std::size_t, 2 > | fshape, | ||
std::span< const std::int32_t > | cells ) |
Interpolate an evaluated expression f(x) in a finite element space.
T | Scalar type. |
U | Mesh geometry type. |
[out] | u | Function object to interpolate into. |
[in] | f | Evaluation of the function f(x) at the physical points x given by interpolation_coords . The element used in interpolation_coords should be the same element as associated with u . The shape of f is (value_size, num_points) , with row-major storage. |
[in] | fshape | Shape of f . |
[in] | cells | Indices of the cells in the mesh on which to interpolate. Should be the same as the list of cells used when calling interpolation_coords . |
void interpolate | ( | Function< T, U > & | u1, |
std::span< const std::int32_t > | cells1, | ||
const Function< T, U > & | u0, | ||
std::span< const std::int32_t > | cells0 ) |
Interpolate from one finite element Function to another Function on the same (sub)mesh.
Interpolation can be performed on a subset of mesh cells and Functions may be defined on 'sub-meshes'.
[out] | u1 | Function to interpolate into. |
[in] | cells1 | Cell indices associated with the mesh of u1 that will be interpolated onto. |
[in] | u0 | Function to b interpolated from. |
[in] | cells0 | Cell indices associated with the mesh of u0 that will be interpolated from. If cells1[i] is the index of a cell in the mesh associated with u1 , then cells0[i] is the index of the same cell but in the mesh associated with u0 . cells0 and cells1 must be the same size. |
std::vector< T > interpolation_coords | ( | const fem::FiniteElement< T > & | element, |
const mesh::Geometry< T > & | geometry, | ||
std::span< const std::int32_t > | cells ) |
Compute the evaluation points in the physical space at which an expression should be computed to interpolate it in a finite element space.
[in] | element | Element to be interpolated into. |
[in] | geometry | Mesh geometry. |
[in] | cells | Indices of the cells in the mesh to compute interpolation coordinates for. |
(3, num_points)
and storage is row-major. void interpolation_matrix | ( | const FunctionSpace< U > & | V0, |
const FunctionSpace< U > & | V1, | ||
auto && | mat_set ) |
Assemble an interpolation operator matrix.
The interpolation operator \(A\) interpolates a function in the space \(V_0\) into a space \(V_1\). If \(u_0\) is the degree-of-freedom vector associated with \(V_0\), then the degree-of-freedom vector \(u_1\) for the interpolated function in \(V_1\) is given by \(u_1=Au_0\).
V1
should be used for the rows of the sparsity pattern, V0
for the columns.[in] | V0 | Space to interpolate from. |
[in] | V1 | Space to interpolate to. |
[in] | mat_set | Functor that sets values in a matrix. |
std::vector< std::int32_t > locate_dofs_geometrical | ( | const FunctionSpace< T > & | V, |
U | marker_fn ) |
Find degrees of freedom whose geometric coordinate is true for the provided marking function.
[in] | V | The function (sub)space on which degrees of freedom will be located. |
[in] | marker_fn | Function marking tabulated degrees of freedom |
std::array< std::vector< std::int32_t >, 2 > locate_dofs_geometrical | ( | const std::array< std::reference_wrapper< const FunctionSpace< T > >, 2 > & | V, |
U | marker_fn ) |
Finds degrees of freedom whose geometric coordinate is true for the provided marking function.
[in] | V | The function (sub)space(s) on which degrees of freedom will be located. The spaces must share the same mesh and element type. |
[in] | marker_fn | Function marking tabulated degrees of freedom |
std::vector< std::int32_t > locate_dofs_topological | ( | const mesh::Topology & | topology, |
const DofMap & | dofmap, | ||
int | dim, | ||
std::span< const std::int32_t > | entities, | ||
bool | remote = true ) |
Find degrees-of-freedom which belong to the provided mesh entities (topological).
[in] | topology | Mesh topology. |
[in] | dofmap | Dofmap that associated DOFs with cells. |
[in] | dim | Topological dimension of mesh entities on which degrees-of-freedom will be located |
[in] | entities | Indices of mesh entities. All DOFs associated with the closure of these indices will be returned |
[in] | remote | True to return also "remotely located" degree-of-freedom indices. Remotely located degree-of-freedom indices are local/owned by the current process, but which the current process cannot identify because it does not recognize mesh entity as a marked. For example, a boundary condition dof at a vertex where this process does not have the associated boundary facet. This commonly occurs with partitioned meshes. |
std::array< std::vector< std::int32_t >, 2 > locate_dofs_topological | ( | const mesh::Topology & | topology, |
std::array< std::reference_wrapper< const DofMap >, 2 > | dofmaps, | ||
int | dim, | ||
std::span< const std::int32_t > | entities, | ||
bool | remote = true ) |
Find degrees-of-freedom which belong to the provided mesh entities (topological).
[in] | topology | Mesh topology. |
[in] | dofmaps | The dofmaps. |
[in] | dim | Topological dimension of mesh entities on which degrees-of-freedom will be located |
[in] | entities | Indices of mesh entities. All DOFs associated with the closure of these indices will be returned |
[in] | remote | True to return also "remotely located" degree-of-freedom indices. Remotely located degree-of-freedom indices are local/owned by the current process, but which the current process cannot identify because it does not recognize mesh entity as a marked. For example, a boundary condition dof at a vertex where this process does not have the associated boundary facet. This commonly occurs with partitioned meshes. |
void pack_coefficients | ( | const Form< T, U > & | form, |
std::map< std::pair< IntegralType, int >, std::pair< std::vector< T >, int > > & | coeffs ) |
Pack coefficients of a Form.
[in] | form | Form to pack the coefficients for. |
[in,out] | coeffs | Map from a (integral_type, domain_id) pair to a (coeffs, cstride) pair.
|
void pack_coefficients | ( | std::vector< std::reference_wrapper< const Function< T, U > > > | coeffs, |
std::span< const int > | offsets, | ||
fem::MDSpan2 auto | entities, | ||
std::span< T > | c ) |
Pack coefficient data over a list of cells or facets.
Typically used to prepare coefficient data for an Expression.
T | |
U |
coeffs | Coefficients to pack | |
offsets | Offsets | |
entities | Entities to pack over | |
[out] | c | Packed coefficients. |
std::vector< typename U::scalar_type > pack_constants | ( | const U & | u | ) |
Pack constants of an Expression or Form into a single array ready for assembly.
u | The Expression or Form to pack constant data for. |
std::vector< T > pack_constants | ( | std::vector< std::reference_wrapper< const fem::Constant< T > > > | c | ) |
Pack constants of an Expression or Form into a single array ready for assembly.
c | Constants to pack. |
void set_diagonal | ( | auto | set_fn, |
const FunctionSpace< U > & | V, | ||
const std::vector< std::reference_wrapper< const DirichletBC< T, U > > > & | bcs, | ||
T | diagonal = 1.0 ) |
Sets a value to the diagonal of the matrix for rows with a Dirichlet boundary conditions applied.
This function is typically called after assembly. The assembly function zeroes Dirichlet rows and columns. This function adds the value only to rows that are locally owned, and therefore does not create a need for parallel communication. For block matrices, this function should normally be called only on the diagonal blocks, i.e. blocks for which the test and trial spaces are the same.
[in] | set_fn | The function for setting values to a matrix. |
[in] | V | The function space for the rows and columns of the matrix. It is used to extract only the Dirichlet boundary conditions that are define on V or subspaces of V. |
[in] | bcs | The Dirichlet boundary conditions. |
[in] | diagonal | Value to add to the diagonal for rows with a boundary condition applied. |
void set_diagonal | ( | auto | set_fn, |
std::span< const std::int32_t > | rows, | ||
T | diagonal = 1.0 ) |
Sets a value to the diagonal of a matrix for specified rows.
This function is typically called after assembly. The assembly function zeroes Dirichlet rows and columns. For block matrices, this function should normally be called only on the diagonal blocks, i.e. blocks for which the test and trial spaces are the same.
[in] | set_fn | The function for setting values to a matrix. |
[in] | rows | Row blocks, in local indices, for which to add a value to the diagonal. |
[in] | diagonal | Value to add to the diagonal for the specified rows. |
void tabulate_expression | ( | std::span< T > | values, |
const fem::Expression< T, U > & | e, | ||
const mesh::Mesh< U > & | mesh, | ||
fem::MDSpan2 auto | entities ) |
Evaluate an Expression on cells or facets.
T | Scalar type. |
U | Geometry type |
[in,out] | values | Array to fill with computed values. Row major storage. Sizing should be (num_cells, num_points * value_size * / num_all_argument_dofs columns) . facet index) tuples. Array is flattened per entity. |
[in] | e | Expression to evaluate. |
[in] | mesh | Mesh to compute e on. |
[in] | entities | Mesh entities to evaluate the expression over. For expressions executed on cells, rank is 1 and size is the number of cells. For expressions executed on facets rank is 2, and shape is (num_facets, 2) , where entities[i, 0] is the cell index and entities[i, 1] is the local index of the facet relative to the cell. |
void tabulate_expression | ( | std::span< T > | values, |
const fem::Expression< T, U > & | e, | ||
md::mdspan< const T, md::dextents< std::size_t, 2 > > | coeffs, | ||
std::span< const T > | constants, | ||
const mesh::Mesh< U > & | mesh, | ||
fem::MDSpan2 auto | entities, | ||
std::optional< std::pair< std::reference_wrapper< const FiniteElement< U > >, std::size_t > > | element ) |
Evaluate an Expression on cells or facets.
This function accepts packed coefficient data, which allows it be called without re-packing all coefficient data at each evaluation.
T | Scalar type. |
U | Geometry type |
[in,out] | values | Array to fill with computed values. Shape is (num_entities, num_points, value_size, num_argument_dofs) and storage is row-major. |
[in] | e | Expression to evaluate. |
[in] | coeffs | Packed coefficients for the Expressions. Typically computed using fem::pack_coefficients. |
[in] | constants | Packed constant data. Typically computed using fem::pack_constants. |
[in] | entities | Mesh entities to evaluate the expression over. For cells it is a list of cell indices. For facets is is a list of (cell index, local facet index) index pairs, i.e. entities=[cell0, / facet_local0, cell1, facet_local1, ...] . |
[in] | mesh | Mesh that the Expression is evaluated on. |
[in] | element | Argument element and argument space dimension. |
graph::AdjacencyList< std::int32_t > transpose_dofmap | ( | md::mdspan< const std::int32_t, md::dextents< std::size_t, 2 > > | dofmap, |
std::int32_t | num_cells ) |
Create an adjacency list that maps a global index (process-wise) to the 'unassembled' cell-wise contributions.
It is built from the usual (cell, local index) -> global index dof map. An 'unassembled' vector is the stacked cell contributions, ordered by cell index. If the usual dof map is:
Cell: 0 1 2 3
Global index: [ [0, 3, 5], [3, 2, 4], [4, 3, 2], [2, 1, 0]]
the 'transpose' dof map will be:
Global index: 0 1 2 3 4 5
Unassembled index: [ [0, 11], [10], [4, 8, 9], [1, 3, 7], [5, 6], [2] ]
[in] | dofmap | The standard dof map that for each cell (node) gives the global (process-wise) index of each local (cell-wise) index. |
[in] | num_cells | The number of cells (nodes) in dofmap to consider. The first num_cells are used. This is argument is typically used to exclude ghost cell contributions. |