Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/df/d02/classdolfinx_1_1fem_1_1Form.html
DOLFINx  0.5.1
DOLFINx C++ interface
Classes | Public Types | Public Member Functions | List of all members
Form< T > Class Template Reference

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

#include <Form.h>

Public Types

using scalar_type = T
 Scalar type (T)
 

Public Member Functions

 Form (const std::vector< std::shared_ptr< const FunctionSpace >> &function_spaces, const std::map< IntegralType, std::pair< std::vector< std::pair< int, std::function< void(T *, const T *, const T *, const scalar_value_type_t *, const int *, const std::uint8_t *)>>>, const mesh::MeshTags< int > * >> &integrals, const std::vector< std::shared_ptr< const Function< T >>> &coefficients, const std::vector< std::shared_ptr< const Constant< T >>> &constants, bool needs_facet_permutations, const std::shared_ptr< const mesh::Mesh > &mesh=nullptr)
 Create a finite element form. More...
 
 Form (const Form &form)=delete
 Copy constructor.
 
 Form (Form &&form)=default
 Move constructor.
 
virtual ~Form ()=default
 Destructor.
 
int rank () const
 Rank of the form (bilinear form = 2, linear form = 1, functional = 0, etc) More...
 
std::shared_ptr< const mesh::Meshmesh () const
 Extract common mesh for the form. More...
 
const std::vector< std::shared_ptr< const FunctionSpace > > & function_spaces () const
 Return function spaces for all arguments. More...
 
const std::function< void(T *, const T *, const T *, const scalar_value_type_t *, const int *, const std::uint8_t *)> & kernel (IntegralType type, int i) const
 Get the function for 'kernel' for integral i of given type. More...
 
std::set< IntegralTypeintegral_types () const
 Get types of integrals in the form. More...
 
int num_integrals (IntegralType type) const
 Number of integrals of given type. More...
 
std::vector< int > integral_ids (IntegralType type) const
 Get the IDs for integrals (kernels) for given integral type. The IDs correspond to the domain IDs which the integrals are defined for in the form. ID=-1 is the default integral over the whole domain. More...
 
const std::vector< std::int32_t > & cell_domains (int i) const
 Get the list of cell indices for the ith integral (kernel) for the cell domain type. More...
 
const std::vector< std::int32_t > & exterior_facet_domains (int i) const
 Get the list of (cell_index, local_facet_index) pairs for the ith integral (kernel) for the exterior facet domain type. More...
 
const std::vector< std::int32_t > & interior_facet_domains (int i) const
 Get the list of (cell_index_0, local_facet_index_0, cell_index_1, local_facet_index_1) quadruplets for the ith integral (kernel) for the interior facet domain type. More...
 
const std::vector< std::shared_ptr< const Function< T > > > & coefficients () const
 Access coefficients.
 
bool needs_facet_permutations () const
 Get bool indicating whether permutation data needs to be passed into these integrals. More...
 
std::vector< int > coefficient_offsets () const
 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.
 
const std::vector< std::shared_ptr< const Constant< T > > > & constants () const
 Access constants.
 

Detailed Description

template<typename T>
class dolfinx::fem::Form< T >

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).

Constructor & Destructor Documentation

◆ Form()

Form ( const std::vector< std::shared_ptr< const FunctionSpace >> &  function_spaces,
const std::map< IntegralType, std::pair< std::vector< std::pair< int, std::function< void(T *, const T *, const T *, const scalar_value_type_t *, const int *, const std::uint8_t *)>>>, const mesh::MeshTags< int > * >> &  integrals,
const std::vector< std::shared_ptr< const Function< T >>> &  coefficients,
const std::vector< std::shared_ptr< const Constant< T >>> &  constants,
bool  needs_facet_permutations,
const std::shared_ptr< const mesh::Mesh > &  mesh = nullptr 
)
inline

Create a finite element form.

Note
User applications will normally call a fem::Form builder function rather using this interfcae directly.
Parameters
[in]function_spacesFunction spaces for the form arguments
[in]integralsThe integrals in the form. The first key is the domain type. For each key there is a pair (list[domain id, integration kernel], domain markers).
[in]coefficients
[in]constantsConstants in the Form
[in]needs_facet_permutationsSet to true is any of the integration kernels require cell permutation data
[in]meshThe mesh of the domain. This is required when there are not argument functions from which the mesh can be extracted, e.g. for functionals

Member Function Documentation

◆ cell_domains()

const std::vector<std::int32_t>& cell_domains ( int  i) const
inline

Get the list of cell indices for the ith integral (kernel) for the cell domain type.

Parameters
[in]iIntegral ID, i.e. (sub)domain index
Returns
List of active cell entities for the given integral (kernel)

◆ exterior_facet_domains()

const std::vector<std::int32_t>& exterior_facet_domains ( int  i) const
inline

Get the list of (cell_index, local_facet_index) pairs for the ith integral (kernel) for the exterior facet domain type.

Parameters
[in]iIntegral ID, i.e. (sub)domain index
Returns
List of (cell_index, local_facet_index) pairs. This data is flattened with row-major layout, shape=(num_facets, 2)

◆ function_spaces()

const std::vector<std::shared_ptr<const FunctionSpace> >& function_spaces ( ) const
inline

Return function spaces for all arguments.

Returns
Function spaces

◆ integral_ids()

std::vector<int> integral_ids ( IntegralType  type) const
inline

Get the IDs for integrals (kernels) for given integral type. The IDs correspond to the domain IDs which the integrals are defined for in the form. ID=-1 is the default integral over the whole domain.

Parameters
[in]typeIntegral type
Returns
List of IDs for given integral type

◆ integral_types()

std::set<IntegralType> integral_types ( ) const
inline

Get types of integrals in the form.

Returns
Integrals types

◆ interior_facet_domains()

const std::vector<std::int32_t>& interior_facet_domains ( int  i) const
inline

Get the list of (cell_index_0, local_facet_index_0, cell_index_1, local_facet_index_1) quadruplets for the ith integral (kernel) for the interior facet domain type.

Parameters
[in]iIntegral ID, i.e. (sub)domain index
Returns
List of tuples of the form (cell_index_0, local_facet_index_0, cell_index_1, local_facet_index_1). This data is flattened with row-major layout, shape=(num_facets, 4)

◆ kernel()

const std::function<void(T*, const T*, const T*, const scalar_value_type_t*, const int*, const std::uint8_t*)>& kernel ( IntegralType  type,
int  i 
) const
inline

Get the function for 'kernel' for integral i of given type.

Parameters
[in]typeIntegral type
[in]iDomain index
Returns
Function to call for tabulate_tensor

◆ mesh()

std::shared_ptr<const mesh::Mesh> mesh ( ) const
inline

Extract common mesh for the form.

Returns
The mesh

◆ needs_facet_permutations()

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()

int num_integrals ( IntegralType  type) const
inline

Number of integrals of given type.

Parameters
[in]typeIntegral type
Returns
Number of integrals

◆ rank()

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: