DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Form< T, U > Class Template Reference

A representation of finite element variational forms. More...

#include <Form.h>

Public Types

using scalar_type = T
 Scalar type.
using geometry_type = U
 Geometry type.

Public Member Functions

template<typename X>
requires std::is_convertible_v< std::remove_cvref_t<X>, std::map<std::tuple<IntegralType, int, int>, integral_data<scalar_type, geometry_type>>>
 Form (const std::vector< std::shared_ptr< const FunctionSpace< geometry_type > > > &V, X &&integrals, std::shared_ptr< const mesh::Mesh< geometry_type > > mesh, const std::vector< std::shared_ptr< const Function< scalar_type, geometry_type > > > &coefficients, const std::vector< std::shared_ptr< const Constant< scalar_type > > > &constants, bool needs_facet_permutations, const std::vector< std::reference_wrapper< const mesh::EntityMap > > &entity_maps)
 Create a finite element form.
 Form (const Form &form)=delete
 Copy constructor.
 Form (Form &&form)=default
 Move constructor.
virtual ~Form ()=default
 Destructor.
int rank () const
 Rank of the form.
std::shared_ptr< const mesh::Mesh< geometry_type > > mesh () const
 Common mesh for the form (the 'integration domain').
const std::vector< std::shared_ptr< const FunctionSpace< geometry_type > > > & function_spaces () const
 Function spaces for all arguments.
std::function< void(scalar_type *, const scalar_type *, const scalar_type *, const geometry_type *, const int *, const uint8_t *, void *)> kernel (IntegralType type, int id, int kernel_idx) const
 Get the kernel function for an integral.
std::set< IntegralTypeintegral_types () const
 Get types of integrals in the form.
std::vector< int > active_coeffs (IntegralType type, int id) const
 Indices of coefficients that are active for a given integral (kernel).
int num_integrals (IntegralType type, int kernel_idx) const
 Get number of integrals (kernels) for a given integral type and kernel index.
std::span< const std::int32_t > domain (IntegralType type, int idx, int kernel_idx) const
 Mesh entity indices to integrate over for a given integral (kernel).
std::span< const std::int32_t > domain_arg (IntegralType type, int rank, int idx, int kernel_idx) const
 Argument function mesh integration entity indices.
std::span< const std::int32_t > domain_coeff (IntegralType type, int idx, int c) const
 Coefficient function mesh integration entity indices.
const std::vector< std::shared_ptr< const Function< scalar_type, geometry_type > > > & coefficients () const
 Access coefficients.
bool needs_facet_permutations () const
 Get bool indicating whether permutation data needs to be passed into these integrals.
std::vector< int > coefficient_offsets () const
 Offset for each coefficient expansion array on a cell.
const std::vector< std::shared_ptr< const Constant< scalar_type > > > & constants () const
 Access constants.

Detailed Description

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
class dolfinx::fem::Form< T, U >

A representation of finite element variational forms.

A note on the order of trial and test spaces: FEniCS numbers argument spaces starting with the leading dimension of the corresponding tensor (matrix). In other words, the test space is numbered 0 and the trial space is numbered 1. However, in order to have a notation that agrees with most existing finite element literature, in particular

\[ a = a(u, v) \]

the spaces are numbered from right to left

\[ a: V_1 \times V_0 \rightarrow \mathbb{R} \]

This is reflected in the ordering of the spaces that should be supplied to generated subclasses. In particular, when a bilinear form is initialized, it should be initialized as a(V_1, V_0) = / ..., where V_1 is the trial space and V_0 is the test space. However, when a form is initialized by a list of argument spaces (the variable function_spaces in the constructors below), the list of spaces should start with space number 0 (the test space) and then space number 1 (the trial space).

Template Parameters
TScalar type in the form.
UFloat (real) type used for the finite element and geometry.
KernElement kernel.

Constructor & Destructor Documentation

◆ Form()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
template<typename X>
requires std::is_convertible_v< std::remove_cvref_t<X>, std::map<std::tuple<IntegralType, int, int>, integral_data<scalar_type, geometry_type>>>
Form ( const std::vector< std::shared_ptr< const FunctionSpace< geometry_type > > > & V,
X && integrals,
std::shared_ptr< const mesh::Mesh< geometry_type > > mesh,
const std::vector< std::shared_ptr< const Function< scalar_type, geometry_type > > > & coefficients,
const std::vector< std::shared_ptr< const Constant< scalar_type > > > & constants,
bool needs_facet_permutations,
const std::vector< std::reference_wrapper< const mesh::EntityMap > > & entity_maps )
inline

Create a finite element form.

Note
User applications will normally call a factory function rather using this interface directly.
Parameters
[in]VFunction spaces for the form arguments, e.g. test and trial function spaces.
[in]integralsIntegrals in the form, where integrals[IntegralType, i, kernel index] returns the ith integral (integral_data) of type IntegralType with kernel index kernel index. The i-index refers to the position of a kernel when flattened by sorted subdomain ids, sorted by subdomain ids. The subdomain ids can contain duplicate entries referring to different kernels over the same subdomain.
[in]coefficientsCoefficients in the form.
[in]constantsConstants in the form.
[in]meshMesh of the domain to integrate over (the 'integration domain').
[in]needs_facet_permutationsSet to true is any of the integration kernels require cell permutation data.
[in]entity_mapsIf any trial functions, test functions, or coefficients in the form are not defined on mesh (the 'integration domain'),entity_maps must be supplied. For each key (a mesh, which is different to mesh) an array map must be provided which relates the entities in mesh to the entities in the key mesh e.g. for a key/value pair (mesh0, emap) in entity_maps, emap[i] is the entity in mesh0 corresponding to entity i in mesh.
Note
For the single domain case, pass an empty entity_maps.

Member Function Documentation

◆ active_coeffs()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
std::vector< int > active_coeffs ( IntegralType type,
int id ) const
inline

Indices of coefficients that are active for a given integral (kernel).

A form is split into multiple integrals (kernels) and each integral might container only a subset of all coefficients in the form. This function returns an indicator array for a given integral kernel that signifies which coefficients are present.

Parameters
[in]typeIntegral type.
[in]idDomain index (identifier) of the integral.

◆ coefficient_offsets()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
std::vector< int > coefficient_offsets ( ) const
inline

Offset for each coefficient expansion array on a cell.

Used to pack data for multiple coefficients in a flat array. The last entry is the size required to store all coefficients.

◆ domain()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
std::span< const std::int32_t > domain ( IntegralType type,
int idx,
int kernel_idx ) const
inline

Mesh entity indices to integrate over for a given integral (kernel).

These are the entities in the mesh returned by mesh that are integrated over by a given integral (kernel).

  • For IntegralType::cell, returns a list of cell indices.
  • For IntegralType::exterior_facet, returns a list with shape (num_facets, 2), where [cell_index, 0] is the cell index and [cell_index, 1] is the local facet index relative to the cell.
  • For IntegralType::interior_facet the shape is (num_facets, 4), where [cell_index, 0] is one attached cell and [cell_index, 1] is the is the local facet index relative to the cell, and [cell_index, 2] is the other one attached cell and [cell_index, 1] is the is the local facet index relative to this cell. Storage is row-major.
Parameters
[in]typeIntegral type.
[in]idxIntegral index in flattened list of integral kernels. For a form containing two integrals of integral_a and integral_b with subdomain-ids (1, 4) and (3, 4, 5) respectively, the integrals are stored as a flattened list, sorted by sudomain-ids
auto form_integrals = {integral_a, integral_b, integral_a, integral_b,
integral_b}; auto form_integral_ids = {1, 3, 4, 4, 5}.
[in]kernel_idxIndex of the kernel with in the domain (we may have multiple kernels for a given ID in mixed-topology meshes).
Returns
Entity indices in the mesh::Mesh returned by mesh() to integrate over.

◆ domain_arg()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
std::span< const std::int32_t > domain_arg ( IntegralType type,
int rank,
int idx,
int kernel_idx ) const
inline

Argument function mesh integration entity indices.

Integration can be performed over cells/facets involving functions that are defined on different meshes but which share common cells, i.e. meshes can be 'views' into a common mesh. Meshes can share some cells but a common cell will have a different index in each mesh::Mesh. Consider:

auto mesh = this->mesh();
auto entities = this->domain(type, id, kernel_idx);
auto entities0 = this->domain_arg(type, rank, id, kernel_idx);
int rank() const
Rank of the form.
Definition Form.h:357
std::shared_ptr< const mesh::Mesh< geometry_type > > mesh() const
Common mesh for the form (the 'integration domain').
Definition Form.h:361
std::span< const std::int32_t > domain_arg(IntegralType type, int rank, int idx, int kernel_idx) const
Argument function mesh integration entity indices.
Definition Form.h:525
std::span< const std::int32_t > domain(IntegralType type, int idx, int kernel_idx) const
Mesh entity indices to integrate over for a given integral (kernel).
Definition Form.h:482

Assembly is performed over entities, where entities[i] is an entity index (e.g., cell index) in mesh. entities0 holds the corresponding entity indices but in the mesh associated with the argument function (test/trial function) space. entities[i] and entities0[i] point to the same mesh entity, but with respect to different mesh views. In some cases, such as when integrating over the interface between two domains that do not overlap, an entity may exist in one domain but not another. In this case, the entity is marked with -1.

Parameters
[in]typeIntegral type.
[in]rankArgument index, e.g. 0 for the test function space, 1 for the trial function space.
[in]idxIntegral identifier.
[in]kernel_idxKernel index (cell type).
Returns
Entity indices in the argument function space mesh that is integrated over.
  • For cell integrals it has shape (num_cells,).
  • For exterior/interior facet integrals, it has shape (num_facts, 2) (row-major storage), where [i, 0] is the index of a cell and [i, 1] is the local index of the facet relative to the cell.

◆ domain_coeff()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
std::span< const std::int32_t > domain_coeff ( IntegralType type,
int idx,
int c ) const
inline

Coefficient function mesh integration entity indices.

This method is equivalent to domain_arg, but returns mesh entity indices for coefficient Functions.

Parameters
[in]typeIntegral type.
[in]idxIntegral identifier.
[in]cCoefficient index.
Returns
Entity indices in the coefficient function space mesh that is integrated over.
  • For cell integrals it has shape (num_cells,).
  • For exterior/interior facet integrals, it has shape (num_facts, 2) (row-major storage), where [i, 0] is the index of a cell and [i, 1] is the local index of the facet relative to the cell.

◆ function_spaces()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
const std::vector< std::shared_ptr< const FunctionSpace< geometry_type > > > & function_spaces ( ) const
inline

Function spaces for all arguments.

Returns
Function spaces.

◆ integral_types()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
std::set< IntegralType > integral_types ( ) const
inline

Get types of integrals in the form.

Returns
Integrals types.

◆ kernel()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
std::function< void(scalar_type *, const scalar_type *, const scalar_type *, const geometry_type *, const int *, const uint8_t *, void *)> kernel ( IntegralType type,
int id,
int kernel_idx ) const
inline

Get the kernel function for an integral.

Parameters
[in]typeIntegral type.
[in]idIntegral subdomain ID.
[in]kernel_idxIndex of the kernel (we may have multiple kernels for a given ID in mixed-topology meshes).
Returns
Function to call for tabulate_tensor.

◆ mesh()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
std::shared_ptr< const mesh::Mesh< geometry_type > > mesh ( ) const
inline

Common mesh for the form (the 'integration domain').

Returns
The integration domain mesh.

◆ needs_facet_permutations()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
bool needs_facet_permutations ( ) const
inline

Get bool indicating whether permutation data needs to be passed into these integrals.

Returns
True if cell permutation data is required

◆ num_integrals()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
int num_integrals ( IntegralType type,
int kernel_idx ) const
inline

Get number of integrals (kernels) for a given integral type and kernel index.

For a form containing two integrals of integral_a and integral_b with subdomain-ids (1, 4) and (3, 4, 5) respectively, the integrals are stored as a flattened list, sorted by sudomain-ids

auto form_integrals = {integral_a, integral_b, integral_a, integral_b,
integral_b}; auto form_integral_ids = {1, 3, 4, 4, 5}.
Parameters
[in]typeIntegral type.
[in]kernel_idxIndex of the kernel (we may have multiple kernels for a integral type in mixed-topology meshes).

◆ rank()

template<dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
int rank ( ) const
inline

Rank of the form.

bilinear form = 2, linear form = 1, functional = 0, etc.

Returns
The rank of the form.

The documentation for this class was generated from the following files: