|
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::map< std::shared_ptr< const mesh::Mesh< geometry_type > >, std::span< const std::int32_t > > &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< IntegralType > | integral_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).
|
|
std::vector< int > | integral_ids (IntegralType type) const |
| Get the IDs for integrals (kernels) for given integral domain type.
|
|
std::span< const std::int32_t > | domain (IntegralType type, int id, 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 id, int kernel_idx) const |
| Argument function mesh integration entity indices.
|
|
std::span< const std::int32_t > | domain_coeff (IntegralType type, int id, 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.
|
|
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
-
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_t<T>>
std::span< const std::int32_t > domain_arg |
( |
IntegralType | type, |
|
|
int | rank, |
|
|
int | id, |
|
|
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 entities = this->
domain(type,
id, kernel_idx);
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.
- Parameters
-
type | Integral type. |
rank | Argument index, e.g. 0 for the test function space, 1 for the trial function space. |
id | Integral domain identifier. |
kernel_idx | Kernel 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.