DOLFINx 0.7.3
DOLFINx C++ interface
|
Finite element method functionality. More...
Namespaces | |
namespace | petsc |
Helper functions for assembly into PETSc data structures. | |
namespace | sparsitybuild |
Support for building sparsity patterns from degree-of-freedom maps. | |
Classes | |
class | Constant |
Constant value which can be attached to a Form. More... | |
class | CoordinateElement |
A CoordinateElement manages coordinate mappings for isoparametric cells. More... | |
class | DirichletBC |
Object for setting (strong) Dirichlet boundary conditions. More... | |
class | DofMap |
Degree-of-freedom map. More... | |
class | ElementDofLayout |
The class represents the degree-of-freedom (dofs) for an element. Dofs are associated with a mesh entity. This class also handles sub-space dofs, which are views into the parent dofs. More... | |
class | Expression |
Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell. More... | |
class | FiniteElement |
Finite Element, containing the dof layout on a reference element, and various methods for evaluating and transforming the basis. More... | |
class | Form |
A representation of finite element variational forms. More... | |
class | Function |
This class represents a function \( u_h \) in a finite element function space \( V_h \), given by. More... | |
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... | |
Concepts | |
concept | FEkernel |
Finite element cell kernel concept. | |
Enumerations | |
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::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, const std::vector< std::shared_ptr< 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::shared_ptr< const DirichletBC< T, U > > > > &bcs1, const std::vector< std::span< const T > > &x0, T scale) |
Modify b such that: | |
template<dolfinx::scalar T, std::floating_point U> | |
void | apply_lifting (std::span< T > b, const std::vector< std::shared_ptr< const Form< T, U > > > &a, const std::vector< std::vector< std::shared_ptr< const DirichletBC< T, U > > > > &bcs1, const std::vector< std::span< const T > > &x0, T scale) |
Modify b such that: | |
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::shared_ptr< 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::shared_ptr< 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::shared_ptr< 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. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | set_bc (std::span< T > b, const std::vector< std::shared_ptr< const DirichletBC< T, U > > > &bcs, std::span< const T > x0, T scale=1) |
Set bc values in owned (local) part of the vector, multiplied by 'scale'. The vectors b and x0 must have the same local size. The bcs should be on (sub-)spaces of the form L that b represents. | |
template<dolfinx::scalar T, std::floating_point U> | |
void | set_bc (std::span< T > b, const std::vector< std::shared_ptr< const DirichletBC< T, U > > > &bcs, T scale=1) |
Set bc values in owned (local) part of the vector, multiplied by 'scale'. The bcs should be on (sub-)spaces of the form L that b represents. | |
std::vector< std::int32_t > | locate_dofs_topological (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). Note that degrees-of-freedom for discontinuous elements are associated with the cell even if they may appear to be associated with a facet/edge/vertex. | |
std::array< std::vector< std::int32_t >, 2 > | locate_dofs_topological (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). Note that degrees-of-freedom for discontinuous elements are associated with the cell even if they may appear to be associated with a facet/edge/vertex. | |
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) |
Finds degrees of freedom whose geometric coordinate is true for the provided marking function. | |
template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_type_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 (MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const std::int32_t, MDSPAN_IMPL_STANDARD_NAMESPACE::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::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 an element on a mesh topology. | |
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. 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. | |
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 expression f(x) in a finite element space. | |
template<std::floating_point T> | |
std::tuple< std::vector< std::int32_t >, std::vector< std::int32_t >, std::vector< T >, std::vector< std::int32_t > > | create_nonmatching_meshes_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 discrete functions across different meshes. | |
template<std::floating_point T> | |
std::tuple< std::vector< std::int32_t >, std::vector< std::int32_t >, std::vector< T >, std::vector< std::int32_t > > | create_nonmatching_meshes_interpolation_data (const mesh::Mesh< T > &mesh0, const FiniteElement< T > &element0, const mesh::Mesh< T > &mesh1, T padding) |
Generate data needed to interpolate discrete functions defined on different meshes. Interpolate on all cells in the mesh. | |
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 std::tuple< std::vector< std::int32_t >, std::vector< std::int32_t >, std::vector< U >, std::vector< std::int32_t > > &nmm_interpolation_data={}) |
Interpolate from one finite element Function to another one. | |
std::vector< std::pair< int, std::vector< std::int32_t > > > | compute_integration_domains (IntegralType integral_type, const mesh::Topology &topology, std::span< const std::int32_t > entities, int dim, std::span< const int > values) |
Given an integral type and mesh tag data, compute 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. | |
ElementDofLayout | create_element_dof_layout (const ufcx_dofmap &dofmap, const mesh::CellType cell_type, const std::vector< int > &parent_map={}) |
Create an ElementDofLayout from a ufcx_dofmap. | |
DofMap | create_dofmap (MPI_Comm comm, const ElementDofLayout &layout, mesh::Topology &topology, std::function< void(const std::span< std::int32_t > &, std::uint32_t)> unpermute_dofs, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> reorder_fn) |
Create a dof map on mesh. | |
std::vector< std::string > | get_coefficient_names (const ufcx_form &ufcx_form) |
Get the name of each coefficient in a UFC 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, typename U = dolfinx::scalar_value_type_t<T>> | |
Form< T, U > | create_form (const ufcx_form &ufcx_form, 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::vector< std::int32_t > > > > &subdomains, std::shared_ptr< const mesh::Mesh< U > > mesh=nullptr) |
Create a Form from UFC input. | |
template<dolfinx::scalar T, typename U = dolfinx::scalar_value_type_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::vector< std::int32_t > > > > &subdomains, std::shared_ptr< const mesh::Mesh< U > > mesh=nullptr) |
Create a Form from UFC input. | |
template<dolfinx::scalar T, typename U = dolfinx::scalar_value_type_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::vector< std::int32_t > > > > &subdomains, 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, const basix::FiniteElement< T > &e, const std::vector< std::size_t > &value_shape={}, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> reorder_fn=nullptr) |
Create a function space from a Basix element. | |
template<std::floating_point T> | |
FunctionSpace< T > | create_functionspace (ufcx_function_space *(*fptr)(const char *), const std::string &function_name, std::shared_ptr< mesh::Mesh< T > > mesh, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> reorder_fn=nullptr) |
Create a FunctionSpace from UFC data. | |
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, IntegralType integral_type, int id, std::span< T > c, int cstride) |
Pack coefficients of a Form for a given integral type and domain id. | |
template<dolfinx::scalar T, typename U = dolfinx::scalar_value_type_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_function_space=nullptr) |
Create Expression from UFC. | |
template<dolfinx::scalar T, typename U = dolfinx::scalar_value_type_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_function_space=nullptr) |
Create Expression from UFC input (with named coefficients and constants). | |
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> | |
std::pair< std::vector< T >, int > | pack_coefficients (const Expression< T, U > &e, std::span< const std::int32_t > cells) |
Pack coefficients of a Expression u for a give list of active cells. | |
template<typename U > | |
std::vector< typename U::scalar_type > | pack_constants (const U &u) |
Pack constants of u of generic type U ready for assembly. | |
Finite element method functionality.
Classes and algorithms for finite element method spaces and operations.
|
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, |
const std::vector< std::shared_ptr< 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::shared_ptr< const DirichletBC< T, U > > > > & | bcs1, | ||
const std::vector< std::span< const T > > & | x0, | ||
T | scale | ||
) |
Modify b such that:
b <- b - scale * A_j (g_j - x0_j)
where j is a block (nest) index. For a non-blocked problem j = 0. The boundary conditions bcs1 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.
Ghost contributions are not accumulated (not sent to owner). Caller is responsible for calling VecGhostUpdateBegin/End.
void apply_lifting | ( | std::span< T > | b, |
const std::vector< std::shared_ptr< const Form< T, U > > > & | a, | ||
const std::vector< std::vector< std::shared_ptr< const DirichletBC< T, U > > > > & | bcs1, | ||
const std::vector< std::span< const T > > & | x0, | ||
T | scale | ||
) |
Modify b such that:
b <- b - scale * A_j (g_j - x0_j)
where j is a block (nest) index. For a non-blocked problem j = 0. The boundary conditions bcs1 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.
Ghost contributions are not accumulated (not sent to owner). Caller is responsible for calling VecGhostUpdateBegin/End.
void assemble_matrix | ( | auto | mat_add, |
const Form< T, U > & | a, | ||
const std::vector< std::shared_ptr< 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 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. |
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::shared_ptr< 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 | The vector to be assembled. It will not be zeroed before assembly. |
[in] | L | The 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::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 an element on a mesh topology.
[in] | comm | MPI communicator |
[in] | topology | The mesh topology |
[in] | element_dof_layouts | The element dof layouts for the function space |
[in] | reorder_fn | Graph reordering function that is applied to the dofmap |
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::pair< int, std::vector< std::int32_t > > > compute_integration_domains | ( | fem::IntegralType | integral_type, |
const mesh::Topology & | topology, | ||
std::span< const std::int32_t > | entities, | ||
int | dim, | ||
std::span< const int > | values | ||
) |
Given an integral type and mesh tag data, compute the entities that should be integrated over.
This function returns as list [(id, entities)]
, where entities
are the entities in meshtags
with tag value id
. 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)
.
[in] | integral_type | Integral type |
[in] | topology | Mesh topology |
[in] | entities | List of tagged mesh entities |
[in] | dim | Topological dimension of tagged entities |
[in] | values | Value associated with each entity |
(integral id, entities)
pairs fem::DofMap create_dofmap | ( | MPI_Comm | comm, |
const ElementDofLayout & | layout, | ||
mesh::Topology & | topology, | ||
std::function< void(const std::span< std::int32_t > &, std::uint32_t)> | unpermute_dofs, | ||
std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> | reorder_fn | ||
) |
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::vector< std::int32_t > > > > & | subdomains, | ||
std::shared_ptr< const mesh::Mesh< U > > | mesh = nullptr |
||
) |
Create a Form from UFC input.
[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. |
subdomains
must be sorted by domain id. [in] | mesh | Mesh of the domain. This is required if the form has no arguments, e.g. a functional. |
Form< T, U > create_form | ( | const ufcx_form & | ufcx_form, |
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::vector< std::int32_t > > > > & | subdomains, | ||
std::shared_ptr< const mesh::Mesh< U > > | mesh = nullptr |
||
) |
Create a Form from UFC input.
[in] | ufcx_form | The UFC form |
[in] | spaces | Vector of function spaces |
[in] | coefficients | Coefficient fields in the form |
[in] | constants | Spatial constants in the form |
[in] | subdomains | Subdomain markers |
subdomains
must be sorted by domain id [in] | mesh | The mesh of the domain |
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::vector< std::int32_t > > > > & | subdomains, | ||
std::shared_ptr< const mesh::Mesh< U > > | mesh = nullptr |
||
) |
Create a Form using a factory function that returns a pointer to a ufcx_form.
[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. |
subdomains
must be sorted by domain id. [in] | mesh | Mesh of the domain. This is required if the form has no arguments, e.g. a functional. |
FunctionSpace< T > create_functionspace | ( | std::shared_ptr< mesh::Mesh< T > > | mesh, |
const basix::FiniteElement< T > & | e, | ||
const std::vector< std::size_t > & | value_shape = {} , |
||
std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> | reorder_fn = nullptr |
||
) |
Create a function space from a Basix element.
[in] | mesh | Mesh |
[in] | e | Basix finite element. |
[in] | value_shape | Value shape for 'blocked' elements, e.g. vector-valued Lagrange elements where each component for the vector field is a Lagrange element. For example, a vector-valued element in 3D will have value_shape equal to {3} , and for a second-order tensor element in 2D value_shape equal to {2, 2} . |
[in] | reorder_fn | The graph reordering function to call on the dofmap. If nullptr , the default re-ordering is used. |
FunctionSpace< T > create_functionspace | ( | ufcx_function_space *(*)(const char *) | fptr, |
const std::string & | function_name, | ||
std::shared_ptr< mesh::Mesh< T > > | mesh, | ||
std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> | reorder_fn = nullptr |
||
) |
Create a FunctionSpace from UFC data.
[in] | fptr | Pointer to a ufcx_function_space_create function. |
[in] | function_name | Name of a function whose function space to create. Function name is the name of Python variable for ufl.Coefficient, ufl.TrialFunction or ufl.TestFunction as defined in the UFL file. |
[in] | mesh | Mesh |
[in] | reorder_fn | Graph reordering function to call on the dofmap. If nullptr , the default re-ordering is used. |
std::tuple< std::vector< std::int32_t >, std::vector< std::int32_t >, std::vector< T >, std::vector< std::int32_t > > create_nonmatching_meshes_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 discrete functions 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 fem::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 runtime of this function, as one has to determine what entity is closest if there is no intersection. std::tuple< std::vector< std::int32_t >, std::vector< std::int32_t >, std::vector< T >, std::vector< std::int32_t > > create_nonmatching_meshes_interpolation_data | ( | const mesh::Mesh< T > & | mesh0, |
const FiniteElement< T > & | element0, | ||
const mesh::Mesh< T > & | mesh1, | ||
T | padding | ||
) |
Generate data needed to interpolate discrete functions defined on different meshes. Interpolate on all cells in the mesh.
[in] | mesh0 | Mesh 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] | 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 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_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\), the hen \(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 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 std::tuple< std::vector< std::int32_t >, std::vector< std::int32_t >, std::vector< U >, std::vector< std::int32_t > > & | nmm_interpolation_data = {} |
||
) |
Interpolate from one finite element Function to another one.
[out] | u | The function to interpolate into |
[in] | v | The function to be interpolated |
[in] | cells | List of cell indices to interpolate on |
[in] | nmm_interpolation_data | Auxiliary data to interpolate on nonmatching meshes. This data can be generated with create_nonmatching_meshes_interpolation_data (optional). |
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 expression f(x) in a finite element space.
[out] | u | The Function object to interpolate into |
[in] | f | Evaluation of the function f(x) at the physical points x given by fem::interpolation_coords. The element used in fem::interpolation_coords should be the same element as associated with u . The shape of f should be (value_size, num_points), with row-major storage. |
[in] | fshape | The shape of f . |
[in] | cells | Indices of the cells in the mesh on which to interpolate. Should be the same as the list used when calling fem::interpolation_coords. |
T | Scalar type |
U | Mesh geometry type |
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 | The element to be interpolated into |
[in] | geometry | Mesh geometry |
[in] | cells | Indices of the cells in the mesh to compute interpolation coordinates for |
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 | The space to interpolate from |
[in] | V1 | The space to interpolate to |
[in] | mat_set | A 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 | ( | 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). Note that degrees-of-freedom for discontinuous elements are associated with the cell even if they may appear to be associated with a facet/edge/vertex.
[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 | ( | 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). Note that degrees-of-freedom for discontinuous elements are associated with the cell even if they may appear to be associated with a facet/edge/vertex.
[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. |
std::pair< std::vector< T >, int > pack_coefficients | ( | const Expression< T, U > & | e, |
std::span< const std::int32_t > | cells | ||
) |
Pack coefficients of a Expression u for a give list of active cells.
[in] | e | The Expression |
[in] | cells | A list of active cells |
void pack_coefficients | ( | const Form< T, U > & | form, |
IntegralType | integral_type, | ||
int | id, | ||
std::span< T > | c, | ||
int | cstride | ||
) |
void pack_coefficients | ( | const Form< T, U > & | form, |
std::map< std::pair< IntegralType, int >, std::pair< std::vector< T >, int > > & | coeffs | ||
) |
std::vector< typename U::scalar_type > pack_constants | ( | const U & | u | ) |
Pack constants of u of generic type U ready for assembly.
void set_diagonal | ( | auto | set_fn, |
const FunctionSpace< U > & | V, | ||
const std::vector< std::shared_ptr< 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 | The 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 | The row blocks, in local indices, for which to add a value to the diagonal |
[in] | diagonal | The value to add to the diagonal for the specified rows |
graph::AdjacencyList< std::int32_t > transpose_dofmap | ( | MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const std::int32_t, MDSPAN_IMPL_STANDARD_NAMESPACE::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. |