DOLFINx  0.5.1
DOLFINx C++ interface
Namespaces | Classes | Enumerations | Functions
dolfinx::fem Namespace Reference

Finite element method functionality. More...

Namespaces

 petsc
 Helper functions for assembly into PETSc data structures.
 
 sparsitybuild
 Support for building sparsity patterns from degree-of-freedom maps.
 

Classes

class  DirichletBC
 Object for setting (strong) Dirichlet boundary conditions. More...
 
class  Form
 A representation of finite element variational forms. More...
 
class  Constant
 Constant value which can be attached to a Form. Constants may be scalar (rank 0), vector (rank 1), or tensor-valued. More...
 
class  CoordinateElement
 A CoordinateElement manages coordinate mappings for isoparametric cells. 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. This class closely follows the concept of a UFC Expression. More...
 
class  FiniteElement
 Finite Element, containing the dof layout on a reference element, and various methods for evaluating and transforming the basis. 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 (dofmap). More...
 

Enumerations

enum class  IntegralType : std::int8_t { cell = 0 , exterior_facet = 1 , interior_facet = 2 , vertex = 3 }
 Type of integral. More...
 

Functions

template<typename T >
std::map< std::pair< dolfinx::fem::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 std::span from a map of std::vector.
 
template<typename T >
assemble_scalar (const Form< T > &M, const 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. More...
 
template<typename T >
assemble_scalar (const Form< T > &M)
 Assemble functional into scalar. More...
 
template<typename T >
void assemble_vector (std::span< T > b, const Form< T > &L, const 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. More...
 
template<typename T >
void assemble_vector (std::span< T > b, const Form< T > &L)
 Assemble linear form into a vector. More...
 
template<typename T >
void apply_lifting (std::span< T > b, const std::vector< std::shared_ptr< const Form< T >>> &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 >>>> &bcs1, const std::vector< std::span< const T >> &x0, double scale)
 Modify b such that: More...
 
template<typename T >
void apply_lifting (std::span< T > b, const std::vector< std::shared_ptr< const Form< T >>> &a, const std::vector< std::vector< std::shared_ptr< const DirichletBC< T >>>> &bcs1, const std::vector< std::span< const T >> &x0, double scale)
 Modify b such that: More...
 
template<typename T , typename U >
void assemble_matrix (U mat_add, const Form< T > &a, const 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 >>> &bcs)
 Assemble bilinear form into a matrix. More...
 
template<typename T , typename U >
void assemble_matrix (U mat_add, const Form< T > &a, const std::vector< std::shared_ptr< const DirichletBC< T >>> &bcs)
 Assemble bilinear form into a matrix. More...
 
template<typename T , typename U >
void assemble_matrix (U mat_add, const Form< T > &a, const std::span< const T > &constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int >> &coefficients, const std::span< const std::int8_t > &dof_marker0, const 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. More...
 
template<typename T , typename U >
void assemble_matrix (U mat_add, const Form< T > &a, const std::span< const std::int8_t > &dof_marker0, const 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. More...
 
template<typename T , typename U >
void set_diagonal (U set_fn, const std::span< const std::int32_t > &rows, T diagonal=1.0)
 Sets a value to the diagonal of a matrix for specified rows. It 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. More...
 
template<typename T , typename U >
void set_diagonal (U set_fn, const fem::FunctionSpace &V, const std::vector< std::shared_ptr< const DirichletBC< T >>> &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. More...
 
template<typename T >
void set_bc (std::span< T > b, const std::vector< std::shared_ptr< const DirichletBC< T >>> &bcs, const std::span< const T > &x0, double scale=1.0)
 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<typename T >
void set_bc (std::span< T > b, const std::vector< std::shared_ptr< const DirichletBC< T >>> &bcs, double scale=1.0)
 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 (const FunctionSpace &V, int dim, const 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. More...
 
std::array< std::vector< std::int32_t >, 2 > locate_dofs_topological (const std::array< std::reference_wrapper< const FunctionSpace >, 2 > &V, int dim, const 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. More...
 
std::vector< std::int32_t > locate_dofs_geometrical (const FunctionSpace &V, const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> &marker_fn)
 Finds degrees of freedom whose geometric coordinate is true for the provided marking function. More...
 
std::array< std::vector< std::int32_t >, 2 > locate_dofs_geometrical (const std::array< std::reference_wrapper< const FunctionSpace >, 2 > &V, const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> &marker_fn)
 Finds degrees of freedom whose geometric coordinate is true for the provided marking function. More...
 
template<typename T , typename U >
void discrete_gradient (const FunctionSpace &V0, const FunctionSpace &V1, U &&mat_set)
 Assemble a discrete gradient operator. More...
 
template<typename T , typename U >
void interpolation_matrix (const FunctionSpace &V0, const FunctionSpace &V1, U &&mat_set)
 Assemble an interpolation operator matrix. More...
 
graph::AdjacencyList< std::int32_t > transpose_dofmap (const graph::AdjacencyList< std::int32_t > &dofmap, std::int32_t num_cells)
 Create an adjacency list that maps a global index (process-wise) to the 'unassembled' cell-wise contributions. More...
 
std::tuple< common::IndexMap, int, graph::AdjacencyList< std::int32_t > > build_dofmap_data (MPI_Comm comm, const mesh::Topology &topology, const ElementDofLayout &element_dof_layout, const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &reorder_fn)
 Build dofmap data for an element on a mesh topology. More...
 
std::array< std::vector< std::shared_ptr< const FunctionSpace > >, 2 > common_function_spaces (const std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace >, 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. More...
 
std::vector< double > interpolation_coords (const FiniteElement &element, const mesh::Mesh &mesh, 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. More...
 
template<typename T >
void interpolate (Function< T > &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. More...
 
template<typename T >
void interpolate (Function< T > &u, const Function< T > &v, const std::span< const std::int32_t > &cells)
 Interpolate from one finite element Function to another on the same mesh. More...
 
template<typename T >
std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace >, 2 > > > extract_function_spaces (const std::vector< std::vector< const Form< T > * >> &a)
 Extract test (0) and trial (1) function spaces pairs for each bilinear form for a rectangular array of forms. More...
 
la::SparsityPattern create_sparsity_pattern (const mesh::Topology &topology, const std::array< std::reference_wrapper< const DofMap >, 2 > &dofmaps, const std::set< IntegralType > &integrals)
 Create a sparsity pattern for a given form. More...
 
template<typename T >
la::SparsityPattern create_sparsity_pattern (const Form< T > &a)
 Create a sparsity pattern for a given form. More...
 
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, const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &reorder_fn, const FiniteElement &element)
 Create a dof map on mesh. More...
 
std::vector< std::string > get_coefficient_names (const ufcx_form &ufcx_form)
 Get the name of each coefficient in a UFC form. More...
 
std::vector< std::string > get_constant_names (const ufcx_form &ufcx_form)
 Get the name of each constant in a UFC form. More...
 
template<typename T >
Form< T > create_form (const ufcx_form &ufcx_form, const std::vector< std::shared_ptr< const FunctionSpace >> &spaces, const std::vector< std::shared_ptr< const Function< T >>> &coefficients, const std::vector< std::shared_ptr< const Constant< T >>> &constants, const std::map< IntegralType, const mesh::MeshTags< int > * > &subdomains, const std::shared_ptr< const mesh::Mesh > &mesh=nullptr)
 Create a Form from UFC input. More...
 
template<typename T >
Form< T > create_form (const ufcx_form &ufcx_form, const std::vector< std::shared_ptr< const FunctionSpace >> &spaces, const std::map< std::string, std::shared_ptr< const Function< T >>> &coefficients, const std::map< std::string, std::shared_ptr< const Constant< T >>> &constants, const std::map< IntegralType, const mesh::MeshTags< int > * > &subdomains, const std::shared_ptr< const mesh::Mesh > &mesh=nullptr)
 Create a Form from UFC input. More...
 
template<typename T >
Form< T > create_form (ufcx_form *(*fptr)(), const std::vector< std::shared_ptr< const FunctionSpace >> &spaces, const std::map< std::string, std::shared_ptr< const Function< T >>> &coefficients, const std::map< std::string, std::shared_ptr< const Constant< T >>> &constants, const std::map< IntegralType, const mesh::MeshTags< int > * > &subdomains, const std::shared_ptr< const mesh::Mesh > &mesh=nullptr)
 Create a Form using a factory function that returns a pointer to a ufcx_form. More...
 
FunctionSpace create_functionspace (const std::shared_ptr< mesh::Mesh > &mesh, const basix::FiniteElement &e, int bs, const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &reorder_fn=nullptr)
 Create a FunctionSpace from a Basix element. More...
 
FunctionSpace create_functionspace (ufcx_function_space *(*fptr)(const char *), const std::string &function_name, const std::shared_ptr< mesh::Mesh > &mesh, const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &reorder_fn=nullptr)
 Create a FunctionSpace from UFC data. More...
 
template<typename T >
std::pair< std::vector< T >, int > allocate_coefficient_storage (const Form< T > &form, IntegralType integral_type, int id)
 Allocate storage for coefficients of a pair (integral_type, id) from a fem::Form form. More...
 
template<typename T >
std::map< std::pair< IntegralType, int >, std::pair< std::vector< T >, int > > allocate_coefficient_storage (const Form< T > &form)
 Allocate memory for packed coefficients of a Form. More...
 
template<typename T >
void pack_coefficients (const Form< T > &form, IntegralType integral_type, int id, const std::span< T > &c, int cstride)
 Pack coefficients of a Form for a given integral type and domain id. More...
 
template<typename T >
Expression< T > create_expression (const ufcx_expression &expression, const std::vector< std::shared_ptr< const Function< T >>> &coefficients, const std::vector< std::shared_ptr< const Constant< T >>> &constants, const std::shared_ptr< const mesh::Mesh > &mesh=nullptr, const std::shared_ptr< const FunctionSpace > &argument_function_space=nullptr)
 Create Expression from UFC.
 
template<typename T >
Expression< T > create_expression (const ufcx_expression &expression, const std::map< std::string, std::shared_ptr< const Function< T >>> &coefficients, const std::map< std::string, std::shared_ptr< const Constant< T >>> &constants, const std::shared_ptr< const mesh::Mesh > &mesh=nullptr, const std::shared_ptr< const FunctionSpace > &argument_function_space=nullptr)
 Create Expression from UFC input (with named coefficients and constants)
 
template<typename T >
void pack_coefficients (const Form< T > &form, std::map< std::pair< IntegralType, int >, std::pair< std::vector< T >, int >> &coeffs)
 Pack coefficients of a Form. More...
 
template<typename T >
std::pair< std::vector< T >, int > pack_coefficients (const Expression< T > &u, const std::span< const std::int32_t > &cells)
 Pack coefficients of a Expression u for a give list of active cells. More...
 
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. More...
 

Detailed Description

Finite element method functionality.

Classes and algorithms for finite element method spaces and operations.

Enumeration Type Documentation

◆ IntegralType

enum IntegralType : std::int8_t
strong

Type of integral.

Enumerator
cell 

Cell.

exterior_facet 

Exterior facet.

interior_facet 

Interior facet.

vertex 

Vertex.

Function Documentation

◆ allocate_coefficient_storage() [1/2]

std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, int> > dolfinx::fem::allocate_coefficient_storage ( const Form< T > &  form)

Allocate memory for packed coefficients of a Form.

Parameters
[in]formThe Form
Returns
A map from a form (integral_type, domain_id) pair to a (coeffs, cstride) pair

◆ allocate_coefficient_storage() [2/2]

std::pair<std::vector<T>, int> dolfinx::fem::allocate_coefficient_storage ( const Form< T > &  form,
IntegralType  integral_type,
int  id 
)

Allocate storage for coefficients of a pair (integral_type, id) from a fem::Form form.

Parameters
[in]formThe Form
[in]integral_typeType of integral
[in]idThe id of the integration domain
Returns
A storage container and the column stride

◆ apply_lifting() [1/2]

void dolfinx::fem::apply_lifting ( std::span< T >  b,
const std::vector< std::shared_ptr< const Form< T >>> &  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 >>>> &  bcs1,
const std::vector< std::span< const T >> &  x0,
double  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.

◆ apply_lifting() [2/2]

void dolfinx::fem::apply_lifting ( std::span< T >  b,
const std::vector< std::shared_ptr< const Form< T >>> &  a,
const std::vector< std::vector< std::shared_ptr< const DirichletBC< T >>>> &  bcs1,
const std::vector< std::span< const T >> &  x0,
double  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.

◆ assemble_matrix() [1/4]

void dolfinx::fem::assemble_matrix ( mat_add,
const Form< T > &  a,
const std::span< const std::int8_t > &  dof_marker0,
const 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.

Parameters
[in]mat_addThe function for adding values into the matrix
[in]aThe bilinear form to assemble
[in]dof_marker0Boundary 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_marker1Boundary 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.

◆ assemble_matrix() [2/4]

void dolfinx::fem::assemble_matrix ( mat_add,
const Form< T > &  a,
const std::span< const T > &  constants,
const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int >> &  coefficients,
const std::span< const std::int8_t > &  dof_marker0,
const 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.

Parameters
[in]mat_addThe function for adding values into the matrix
[in]aThe bilinear form to assemble
[in]constantsConstants that appear in a
[in]coefficientsCoefficients that appear in a
[in]dof_marker0Boundary 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_marker1Boundary 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.

◆ assemble_matrix() [3/4]

void dolfinx::fem::assemble_matrix ( mat_add,
const Form< T > &  a,
const 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 >>> &  bcs 
)

Assemble bilinear form into a matrix.

Parameters
[in]mat_addThe function for adding values into the matrix
[in]aThe bilinear from to assemble
[in]constantsConstants that appear in a
[in]coefficientsCoefficients that appear in a
[in]bcsBoundary conditions to apply. For boundary condition dofs the row and column are zeroed. The diagonal entry is not set.

◆ assemble_matrix() [4/4]

void dolfinx::fem::assemble_matrix ( mat_add,
const Form< T > &  a,
const std::vector< std::shared_ptr< const DirichletBC< T >>> &  bcs 
)

Assemble bilinear form into a matrix.

Parameters
[in]mat_addThe function for adding values into the matrix
[in]aThe bilinear from to assemble
[in]bcsBoundary conditions to apply. For boundary condition dofs the row and column are zeroed. The diagonal entry is not set.

◆ assemble_scalar() [1/2]

T dolfinx::fem::assemble_scalar ( const Form< T > &  M)

Assemble functional into scalar.

Note
Caller is responsible for accumulation across processes.
Parameters
[in]MThe form (functional) to assemble
Returns
The contribution to the form (functional) from the local process

◆ assemble_scalar() [2/2]

T dolfinx::fem::assemble_scalar ( const Form< T > &  M,
const 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.

Note
Caller is responsible for accumulation across processes.
Parameters
[in]MThe form (functional) to assemble
[in]constantsThe constants that appear in M
[in]coefficientsThe coefficients that appear in M
Returns
The contribution to the form (functional) from the local process

◆ assemble_vector() [1/2]

void dolfinx::fem::assemble_vector ( std::span< T >  b,
const Form< T > &  L 
)

Assemble linear form into a vector.

Parameters
[in,out]bThe vector to be assembled. It will not be zeroed before assembly.
[in]LThe linear forms to assemble into b

◆ assemble_vector() [2/2]

void dolfinx::fem::assemble_vector ( std::span< T >  b,
const Form< T > &  L,
const 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.

Parameters
[in,out]bThe vector to be assembled. It will not be zeroed before assembly.
[in]LThe linear forms to assemble into b
[in]constantsThe constants that appear in L
[in]coefficientsThe coefficients that appear in L

◆ build_dofmap_data()

std::tuple< common::IndexMap, int, graph::AdjacencyList< std::int32_t > > build_dofmap_data ( MPI_Comm  comm,
const mesh::Topology topology,
const ElementDofLayout element_dof_layout,
const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &  reorder_fn 
)

Build dofmap data for an element on a mesh topology.

Parameters
[in]commMPI communicator
[in]topologyThe mesh topology
[in]element_dof_layoutThe element dof layout for the function space
[in]reorder_fnGraph reordering function that is applied to the dofmap
Returns
The index map and local to global DOF data for the DOF map

◆ common_function_spaces()

std::array< std::vector< std::shared_ptr< const FunctionSpace > >, 2 > common_function_spaces ( const std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace >, 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.

Parameters
[in]VVector function spaces for (0) each row block and (1) each column block

◆ create_dofmap()

fem::DofMap create_dofmap ( MPI_Comm  comm,
const ElementDofLayout layout,
mesh::Topology topology,
const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &  reorder_fn,
const FiniteElement element 
)

Create a dof map on mesh.

Parameters
[in]commMPI communicator
[in]layoutThe dof layout on an element
[in]topologyThe mesh topology
[in]elementThe finite element
[in]reorder_fnThe graph reordering function called on the dofmap
Returns
A new dof map

◆ create_form() [1/3]

Form<T> dolfinx::fem::create_form ( const ufcx_form &  ufcx_form,
const std::vector< std::shared_ptr< const FunctionSpace >> &  spaces,
const std::map< std::string, std::shared_ptr< const Function< T >>> &  coefficients,
const std::map< std::string, std::shared_ptr< const Constant< T >>> &  constants,
const std::map< IntegralType, const mesh::MeshTags< int > * > &  subdomains,
const std::shared_ptr< const mesh::Mesh > &  mesh = nullptr 
)

Create a Form from UFC input.

Parameters
[in]ufcx_formThe UFC form
[in]spacesThe function spaces for the Form arguments
[in]coefficientsCoefficient fields in the form (by name)
[in]constantsSpatial constants in the form (by name)
[in]subdomainsSubdomain makers
[in]meshThe mesh of the domain. This is required if the form has no arguments, e.g. a functional
Returns
A Form

◆ create_form() [2/3]

Form<T> dolfinx::fem::create_form ( const ufcx_form &  ufcx_form,
const std::vector< std::shared_ptr< const FunctionSpace >> &  spaces,
const std::vector< std::shared_ptr< const Function< T >>> &  coefficients,
const std::vector< std::shared_ptr< const Constant< T >>> &  constants,
const std::map< IntegralType, const mesh::MeshTags< int > * > &  subdomains,
const std::shared_ptr< const mesh::Mesh > &  mesh = nullptr 
)

Create a Form from UFC input.

Parameters
[in]ufcx_formThe UFC form
[in]spacesVector of function spaces
[in]coefficientsCoefficient fields in the form
[in]constantsSpatial constants in the form
[in]subdomainsSubdomain markers
[in]meshThe mesh of the domain

◆ create_form() [3/3]

Form<T> dolfinx::fem::create_form ( ufcx_form *(*)()  fptr,
const std::vector< std::shared_ptr< const FunctionSpace >> &  spaces,
const std::map< std::string, std::shared_ptr< const Function< T >>> &  coefficients,
const std::map< std::string, std::shared_ptr< const Constant< T >>> &  constants,
const std::map< IntegralType, const mesh::MeshTags< int > * > &  subdomains,
const std::shared_ptr< const mesh::Mesh > &  mesh = nullptr 
)

Create a Form using a factory function that returns a pointer to a ufcx_form.

Parameters
[in]fptrpointer to a function returning a pointer to ufcx_form
[in]spacesThe function spaces for the Form arguments
[in]coefficientsCoefficient fields in the form (by name)
[in]constantsSpatial constants in the form (by name)
[in]subdomainsSubdomain markers
[in]meshThe mesh of the domain. This is required if the form has no arguments, e.g. a functional.
Returns
A Form

◆ create_functionspace() [1/2]

fem::FunctionSpace create_functionspace ( const std::shared_ptr< mesh::Mesh > &  mesh,
const basix::FiniteElement &  e,
int  bs,
const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &  reorder_fn = nullptr 
)

Create a FunctionSpace from a Basix element.

Parameters
[in]meshMesh
[in]eBasix finite element
[in]bsThe block size, e.g. 3 for a 'vector' Lagrange element in 3D
[in]reorder_fnThe graph reordering function to call on the dofmap. If nullptr, the default re-ordering is used.
Returns
The created function space

◆ create_functionspace() [2/2]

fem::FunctionSpace create_functionspace ( ufcx_function_space *(*)(const char *)  fptr,
const std::string &  function_name,
const std::shared_ptr< mesh::Mesh > &  mesh,
const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &  reorder_fn = nullptr 
)

Create a FunctionSpace from UFC data.

Parameters
[in]fptrFunction Pointer to a ufcx_function_space_create function
[in]function_nameName 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]meshMesh
[in]reorder_fnThe graph reordering function to call on the dofmap. If nullptr, the default re-ordering is used.
Returns
The created function space

◆ create_sparsity_pattern() [1/2]

la::SparsityPattern dolfinx::fem::create_sparsity_pattern ( const Form< T > &  a)

Create a sparsity pattern for a given form.

Note
The pattern is not finalised, i.e. the caller is responsible for calling SparsityPattern::assemble.
Parameters
[in]aA bilinear form
Returns
The corresponding sparsity pattern

◆ create_sparsity_pattern() [2/2]

la::SparsityPattern create_sparsity_pattern ( const mesh::Topology topology,
const std::array< std::reference_wrapper< const DofMap >, 2 > &  dofmaps,
const std::set< IntegralType > &  integrals 
)

Create a sparsity pattern for a given form.

Note
The pattern is not finalised, i.e. the caller is responsible for calling SparsityPattern::assemble.

◆ discrete_gradient()

void dolfinx::fem::discrete_gradient ( const FunctionSpace V0,
const FunctionSpace V1,
U &&  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.

Note
The sparsity pattern for a discrete operator can be initialised using sparsitybuild::cells. The space V1 should be used for the rows of the sparsity pattern, V0 for the columns.
Warning
This function relies on the user supplying appropriate input and output spaces. See parameter descriptions.
Parameters
[in]V0A Lagrange space to interpolate the gradient from
[in]V1A Nédélec (first kind) space to interpolate into
[in]mat_setA functor that sets values in a matrix

◆ extract_function_spaces()

std::vector<std::vector<std::array<std::shared_ptr<const FunctionSpace>, 2> > > dolfinx::fem::extract_function_spaces ( const std::vector< std::vector< const Form< T > * >> &  a)

Extract test (0) and trial (1) function spaces pairs for each bilinear form for a rectangular array of forms.

Parameters
[in]aA rectangular block on bilinear forms
Returns
Rectangular array of the same shape as 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).

◆ get_coefficient_names()

std::vector< std::string > get_coefficient_names ( const ufcx_form &  ufcx_form)

Get the name of each coefficient in a UFC form.

Parameters
[in]ufcx_formThe UFC form return The name of each coefficient

◆ get_constant_names()

std::vector< std::string > get_constant_names ( const ufcx_form &  ufcx_form)

Get the name of each constant in a UFC form.

Parameters
[in]ufcx_formThe UFC form
Returns
The name of each constant

◆ interpolate() [1/2]

void dolfinx::fem::interpolate ( Function< T > &  u,
const Function< T > &  v,
const std::span< const std::int32_t > &  cells 
)

Interpolate from one finite element Function to another on the same mesh.

Parameters
[out]uThe function to interpolate into
[in]vThe function to be interpolated
[in]cellsList of cell indices to interpolate on

◆ interpolate() [2/2]

void dolfinx::fem::interpolate ( Function< T > &  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.

Parameters
[out]uThe function to interpolate into
[in]fEvaluation 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]fshapeThe shape of f.
[in]cellsIndices of the cells in the mesh on which to interpolate. Should be the same as the list used when calling fem::interpolation_coords.

◆ interpolation_coords()

std::vector< double > interpolation_coords ( const FiniteElement element,
const mesh::Mesh mesh,
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.

Parameters
[in]elementThe element to be interpolated into
[in]meshThe domain
[in]cellsIndices of the cells in the mesh to compute interpolation coordinates for
Returns
The coordinates in the physical space at which to evaluate an expression. The shape is (3, num_points) and storage is row-major.

◆ interpolation_matrix()

void dolfinx::fem::interpolation_matrix ( const FunctionSpace V0,
const FunctionSpace V1,
U &&  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\).

Note
The sparsity pattern for a discrete operator can be initialised using sparsitybuild::cells. The space V1 should be used for the rows of the sparsity pattern, V0 for the columns.
Parameters
[in]V0The space to interpolate from
[in]V1The space to interpolate to
[in]mat_setA functor that sets values in a matrix

◆ locate_dofs_geometrical() [1/2]

std::vector< std::int32_t > locate_dofs_geometrical ( const FunctionSpace V,
const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> &  marker_fn 
)

Finds degrees of freedom whose geometric coordinate is true for the provided marking function.

Attention
This function is slower than the topological version
Parameters
[in]VThe function (sub)space on which degrees of freedom will be located.
[in]marker_fnFunction marking tabulated degrees of freedom
Returns
Array of DOF index blocks (local to the MPI rank) in the space V. The array uses the block size of the dofmap associated with V.

◆ locate_dofs_geometrical() [2/2]

std::array< std::vector< std::int32_t >, 2 > locate_dofs_geometrical ( const std::array< std::reference_wrapper< const FunctionSpace >, 2 > &  V,
const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> &  marker_fn 
)

Finds degrees of freedom whose geometric coordinate is true for the provided marking function.

Attention
This function is slower than the topological version
Parameters
[in]VThe function (sub)space(s) on which degrees of freedom will be located. The spaces must share the same mesh and element type.
[in]marker_fnFunction marking tabulated degrees of freedom
Returns
Array of DOF indices (local to the MPI rank) in the spaces V[0] and V[1]. The array[0](i) entry is the DOF index in the space V[0] and array[1](i) is the correspinding DOF entry in the space V[1]. The returned dofs are 'unrolled', i.e. block size = 1.

◆ locate_dofs_topological() [1/2]

std::vector< std::int32_t > locate_dofs_topological ( const FunctionSpace V,
int  dim,
const 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.

Parameters
[in]VThe function (sub)space on which degrees-of-freedom (DOFs) will be located.
[in]dimTopological dimension of mesh entities on which degrees-of-freedom will be located
[in]entitiesIndices of mesh entities. All DOFs associated with the closure of these indices will be returned
[in]remoteTrue 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.
Returns
Array of DOF index blocks (local to the MPI rank) in the space V. The array uses the block size of the dofmap associated with V.

◆ locate_dofs_topological() [2/2]

std::array< std::vector< std::int32_t >, 2 > locate_dofs_topological ( const std::array< std::reference_wrapper< const FunctionSpace >, 2 > &  V,
int  dim,
const 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.

Parameters
[in]VThe function (sub)spaces on which degrees-of-freedom (DOFs) will be located. The spaces must share the same mesh and element type.
[in]dimTopological dimension of mesh entities on which degrees-of-freedom will be located
[in]entitiesIndices of mesh entities. All DOFs associated with the closure of these indices will be returned
[in]remoteTrue 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.
Returns
Array of DOF indices (local to the MPI rank) in the spaces V[0] and V[1]. The array[0](i) entry is the DOF index in the space V[0] and array[1](i) is the corresponding DOF entry in the space V[1]. The returned dofs are 'unrolled', i.e. block size = 1.

◆ pack_coefficients() [1/3]

std::pair<std::vector<T>, int> dolfinx::fem::pack_coefficients ( const Expression< T > &  u,
const std::span< const std::int32_t > &  cells 
)

Pack coefficients of a Expression u for a give list of active cells.

Parameters
[in]uThe Expression
[in]cellsA list of active cells
Returns
A pair of the form (coeffs, cstride)

◆ pack_coefficients() [2/3]

void dolfinx::fem::pack_coefficients ( const Form< T > &  form,
IntegralType  integral_type,
int  id,
const std::span< T > &  c,
int  cstride 
)

Pack coefficients of a Form for a given integral type and domain id.

Parameters
[in]formThe Form
[in]integral_typeType of integral
[in]idThe id of the integration domain
[in]cThe coefficient array
[in]cstrideThe coefficient stride

◆ pack_coefficients() [3/3]

void dolfinx::fem::pack_coefficients ( const Form< T > &  form,
std::map< std::pair< IntegralType, int >, std::pair< std::vector< T >, int >> &  coeffs 
)

Pack coefficients of a Form.

Warning
This is subject to change
Parameters
[in]formThe Form
[in]coeffsA map from a (integral_type, domain_id) pair to a (coeffs, cstride) pair

◆ pack_constants()

std::vector<typename U::scalar_type> dolfinx::fem::pack_constants ( const U &  u)

Pack constants of u of generic type U ready for assembly.

Warning
This function is subject to change

◆ set_diagonal() [1/2]

void dolfinx::fem::set_diagonal ( set_fn,
const fem::FunctionSpace V,
const std::vector< std::shared_ptr< const DirichletBC< T >>> &  bcs,
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.

Parameters
[in]set_fnThe function for setting values to a matrix
[in]VThe 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]bcsThe Dirichlet boundary condtions
[in]diagonalThe value to add to the diagonal for rows with a boundary condition applied

◆ set_diagonal() [2/2]

void dolfinx::fem::set_diagonal ( set_fn,
const std::span< const std::int32_t > &  rows,
diagonal = 1.0 
)

Sets a value to the diagonal of a matrix for specified rows. It 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.

Parameters
[in]set_fnThe function for setting values to a matrix
[in]rowsThe row blocks, in local indices, for which to add a value to the diagonal
[in]diagonalThe value to add to the diagonal for the specified rows

◆ transpose_dofmap()

graph::AdjacencyList< std::int32_t > transpose_dofmap ( const graph::AdjacencyList< std::int32_t > &  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] ]

Parameters
[in]dofmapThe standard dof map that for each cell (node) gives the global (process-wise) index of each local (cell-wise) index.
[in]num_cellsThe 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.
Returns
Map from global (process-wise) index to positions in an unaassembled array. The links for each node are sorted.