|
template<typename X >
requires std::is_convertible_v< std::remove_cvref_t<X>, std::map<IntegralType, std::vector<integral_data< scalar_type, geometry_type>>>> |
| Form (const std::vector< std::shared_ptr< const FunctionSpace< geometry_type > > > &V, X &&integrals, 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::map< std::shared_ptr< const mesh::Mesh< geometry_type > >, std::span< const std::int32_t > > &entity_maps, std::shared_ptr< const mesh::Mesh< geometry_type > > mesh=nullptr) |
| 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 |
| Extract common mesh for the form.
|
|
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 *)> | kernel (IntegralType type, int i) const |
| Get the kernel function for integral i on given domain type.
|
|
std::set< IntegralType > | integral_types () const |
| Get types of integrals in the form.
|
|
int | num_integrals (IntegralType type) const |
| Number of integrals on given domain type.
|
|
std::vector< int > | active_coeffs (IntegralType type, std::size_t i) const |
| Indices of coefficients that are active for a given integral (kernel).
|
|
std::vector< int > | integral_ids (IntegralType type) const |
| Get the IDs for integrals (kernels) for given integral type.
|
|
std::span< const std::int32_t > | domain (IntegralType type, int i) const |
| Get the list of mesh entity indices for the ith integral (kernel) of a given type.
|
|
std::vector< std::int32_t > | domain (IntegralType type, int i, const mesh::Mesh< geometry_type > &mesh) const |
| Compute the list of entity indices in mesh for the ith integral (kernel) of a given type (i.e. cell, exterior facet, or interior facet).
|
|
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.
|
|
template<
dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_type_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
-
T | Scalar type in the form. |
U | Float (real) type used for the finite element and geometry. |
Kern | Element kernel. |
template<
dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_type_t<T>>
std::span< const std::int32_t > domain |
( |
IntegralType | type, |
|
|
int | i ) const |
|
inline |
Get the list of mesh entity indices for the ith integral (kernel) of a given type.
For IntegralType::cell, returns a list of cell indices.
For IntegralType::exterior_facet, returns a list of (cell_index, local_facet_index) pairs. Data is flattened with row-major layout, shape=(num_facets, 2)
.
For IntegralType::interior_facet, returns list of tuples of the form (cell_index_0, local_facet_index_0, cell_index_1, / local_facet_index_1)
. Data is flattened with row-major layout, shape=(num_facets, 4)
.
- Parameters
-
[in] | type | Integral domain type. |
[in] | i | Integral ID, i.e. (sub)domain index. |
- Returns
- List of active entities for the given integral (kernel).