Loading [MathJax]/extensions/tex2jax.js
DOLFINx 0.9.0
DOLFINx C++ interface
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages Concepts
utils.h File Reference

Functions supporting finite element method operations. More...

Go to the source code of this file.

Namespaces

namespace  dolfinx
 Top-level namespace.
 
namespace  dolfinx::common
 Miscellaneous classes, functions and types.
 
namespace  dolfinx::fem
 Finite element method functionality.
 

Concepts

concept  dolfinx::fem::impl::FetchCells
 Concepts for function that returns cell index.
 

Functions

template<int num_cells>
std::array< std::int32_t, 2 *num_cells > get_cell_facet_pairs (std::int32_t f, std::span< const std::int32_t > cells, const graph::AdjacencyList< std::int32_t > &c_to_f)
 
std::vector< std::int32_t > compute_integration_domains (IntegralType integral_type, const mesh::Topology &topology, std::span< const std::int32_t > entities, int dim)
 Given an integral type and a set of entities, 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.
 
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< DofMapcreate_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_type_t<T>>
Form< T, U > create_form_factory (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::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_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::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_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::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, 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<int _bs, dolfinx::scalar T>
void pack (std::span< T > coeffs, std::int32_t cell, int bs, std::span< const T > v, std::span< const std::uint32_t > cell_info, const DofMap &dofmap, auto transform)
 Pack a single coefficient for a single cell.
 
template<dolfinx::scalar T, std::floating_point U>
void pack_coefficient_entity (std::span< T > c, int cstride, const Function< T, U > &u, std::span< const std::uint32_t > cell_info, std::span< const std::int32_t > entities, std::size_t estride, FetchCells auto &&fetch_cells, std::int32_t offset)
 Pack a single coefficient for a set of active entities.
 
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, std::floating_point U = 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, std::floating_point U = 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 > entities, std::size_t estride)
 Pack coefficients of a Expression u for a give list of active entities.
 
template<typename U >
std::vector< typename U::scalar_type > pack_constants (const U &u)
 Pack constants of u into a single array ready for assembly.
 

Detailed Description

Functions supporting finite element method operations.

Function Documentation

◆ get_cell_facet_pairs()

template<int num_cells>
std::array< std::int32_t, 2 *num_cells > get_cell_facet_pairs ( std::int32_t f,
std::span< const std::int32_t > cells,
const graph::AdjacencyList< std::int32_t > & c_to_f )

Helper function to get an array of of (cell, local_facet) pairs corresponding to a given facet index.

Parameters
[in]fFacet index
[in]cellsList of cells incident to the facet
[in]c_to_fCell to facet connectivity
Returns
Vector of (cell, local_facet) pairs

◆ pack_coefficient_entity()

template<dolfinx::scalar T, std::floating_point U>
void pack_coefficient_entity ( std::span< T > c,
int cstride,
const Function< T, U > & u,
std::span< const std::uint32_t > cell_info,
std::span< const std::int32_t > entities,
std::size_t estride,
FetchCells auto && fetch_cells,
std::int32_t offset )

Pack a single coefficient for a set of active entities.

Parameters
[out]cCoefficient to be packed.
[in]cstrideTotal number of coefficient values to pack for each entity.
[in]uFunction to extract coefficient data from.
[in]cell_infoArray of bytes describing which transformation has to be applied on the cell to map it to the reference element.
[in]entitiesSet of active entities.
[in]estrideStride for each entity in active entities.
[in]fetch_cellsFunction that fetches the cell index for an entity in active_entities.
[in]offsetThe offset for c.