DOLFIN-X
DOLFIN-X C++ interface
Public Member Functions | List of all members
dolfinx::fem::Form Class Reference

Base class for variational forms. More...

#include <Form.h>

Public Member Functions

 Form (const std::vector< std::shared_ptr< const function::FunctionSpace >> &function_spaces, const FormIntegrals &integrals, const FormCoefficients &coefficients, const std::vector< std::pair< std::string, std::shared_ptr< const function::Constant >>> constants)
 Create form. More...
 
 Form (const std::vector< std::shared_ptr< const function::FunctionSpace >> &function_spaces)
 Create form (no UFC integrals). Integrals can be attached later using FormIntegrals::set_cell_tabulate_tensor. More...
 
 Form (Form &&form)=default
 Move constructor.
 
virtual ~Form ()=default
 Destructor.
 
int rank () const
 Rank of form (bilinear form = 2, linear form = 1, functional = 0, etc) More...
 
void set_coefficients (std::map< std::size_t, std::shared_ptr< const function::Function >> coefficients)
 Set coefficient with given number (shared pointer version) More...
 
void set_coefficients (std::map< std::string, std::shared_ptr< const function::Function >> coefficients)
 Set coefficient with given name (shared pointer version) More...
 
int original_coefficient_position (int i) const
 Return original coefficient position for each coefficient (0 <= i < n) More...
 
void set_constants (std::map< std::string, std::shared_ptr< const function::Constant >> constants)
 Set constants based on their names. More...
 
void set_constants (std::vector< std::shared_ptr< const function::Constant >> constants)
 Set constants based on their order (without names) More...
 
bool all_constants_set () const
 Check if all constants associated with the form have been set. More...
 
std::set< std::string > get_unset_constants () const
 Return names of any constants that have not been set. More...
 
void set_mesh (std::shared_ptr< const mesh::Mesh > mesh)
 Set mesh, necessary for functionals when there are no function spaces. More...
 
std::shared_ptr< const mesh::Meshmesh () const
 Extract common mesh from form. More...
 
std::shared_ptr< const function::FunctionSpacefunction_space (int i) const
 Return function space for given argument. More...
 
void set_tabulate_tensor (FormIntegrals::Type type, int i, std::function< void(PetscScalar *, const PetscScalar *, const PetscScalar *, const double *, const int *, const std::uint8_t *, const std::uint32_t)> fn)
 Register the function for 'tabulate_tensor' for cell integral i.
 
void set_cell_domains (const mesh::MeshTags< int > &cell_domains)
 Set cell domains. More...
 
void set_exterior_facet_domains (const mesh::MeshTags< int > &exterior_facet_domains)
 Set exterior facet domains. More...
 
void set_interior_facet_domains (const mesh::MeshTags< int > &interior_facet_domains)
 Set interior facet domains. More...
 
void set_vertex_domains (const mesh::MeshTags< int > &vertex_domains)
 Set vertex domains. More...
 
FormCoefficientscoefficients ()
 Access coefficients.
 
const FormCoefficientscoefficients () const
 Access coefficients.
 
const FormIntegralsintegrals () const
 Access form integrals.
 
const std::vector< std::pair< std::string, std::shared_ptr< const function::Constant > > > & constants () const
 Access constants. More...
 

Detailed Description

Base class for 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() [1/2]

Form::Form ( const std::vector< std::shared_ptr< const function::FunctionSpace >> &  function_spaces,
const FormIntegrals integrals,
const FormCoefficients coefficients,
const std::vector< std::pair< std::string, std::shared_ptr< const function::Constant >>>  constants 
)

Create form.

Parameters
[in]function_spacesFunction Spaces
[in]integrals
[in]coefficients
[in]constantsVector of pairs (name, constant). The index in the vector is the position of the constant in the original (nonsimplified) form.

◆ Form() [2/2]

Form::Form ( const std::vector< std::shared_ptr< const function::FunctionSpace >> &  function_spaces)
explicit

Create form (no UFC integrals). Integrals can be attached later using FormIntegrals::set_cell_tabulate_tensor.

Warning
Experimental
Parameters
[in]function_spacesVector of function spaces

Member Function Documentation

◆ all_constants_set()

bool Form::all_constants_set ( ) const

Check if all constants associated with the form have been set.

Returns
True if all Form constants have been set

◆ constants()

const std::vector< std::pair< std::string, std::shared_ptr< const function::Constant > > > & Form::constants ( ) const

Access constants.

Returns
Vector of attached constants with their names. Names are used to set constants in user's c++ code. Index in the vector is the position of the constant in the original (nonsimplified) form.

◆ function_space()

std::shared_ptr< const function::FunctionSpace > Form::function_space ( int  i) const

Return function space for given argument.

Parameters
[in]iIndex of the argument
Returns
Function space

◆ get_unset_constants()

std::set< std::string > Form::get_unset_constants ( ) const

Return names of any constants that have not been set.

Returns
Names of unset constants

◆ mesh()

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

Extract common mesh from form.

Returns
The mesh

◆ original_coefficient_position()

int Form::original_coefficient_position ( int  i) const

Return original coefficient position for each coefficient (0 <= i < n)

Returns
The position of coefficient i in original ufl form coefficients.

◆ rank()

int Form::rank ( ) const

Rank of form (bilinear form = 2, linear form = 1, functional = 0, etc)

Returns
The rank of the form

◆ set_cell_domains()

void Form::set_cell_domains ( const mesh::MeshTags< int > &  cell_domains)

Set cell domains.

Parameters
[in]cell_domainsThe cell domains

◆ set_coefficients() [1/2]

void Form::set_coefficients ( std::map< std::size_t, std::shared_ptr< const function::Function >>  coefficients)

Set coefficient with given number (shared pointer version)

Parameters
[in]coefficientsMap from coefficient index to the coefficient

◆ set_coefficients() [2/2]

void Form::set_coefficients ( std::map< std::string, std::shared_ptr< const function::Function >>  coefficients)

Set coefficient with given name (shared pointer version)

Parameters
[in]coefficientsMap from coefficient name to the coefficient

◆ set_constants() [1/2]

void Form::set_constants ( std::map< std::string, std::shared_ptr< const function::Constant >>  constants)

Set constants based on their names.

This method is used in command-line workflow, when users set constants to the form in cpp file.

Names of the constants must agree with their names in UFL file.

◆ set_constants() [2/2]

void Form::set_constants ( std::vector< std::shared_ptr< const function::Constant >>  constants)

Set constants based on their order (without names)

This method is used in Python workflow, when constants are automatically attached to the form based on their order in the original form.

The order of constants must match their order in original ufl Form.

◆ set_exterior_facet_domains()

void Form::set_exterior_facet_domains ( const mesh::MeshTags< int > &  exterior_facet_domains)

Set exterior facet domains.

Parameters
[in]exterior_facet_domainsThe exterior facet domains

◆ set_interior_facet_domains()

void Form::set_interior_facet_domains ( const mesh::MeshTags< int > &  interior_facet_domains)

Set interior facet domains.

Parameters
[in]interior_facet_domainsThe interior facet domains

◆ set_mesh()

void Form::set_mesh ( std::shared_ptr< const mesh::Mesh mesh)

Set mesh, necessary for functionals when there are no function spaces.

Parameters
[in]meshThe mesh

◆ set_vertex_domains()

void Form::set_vertex_domains ( const mesh::MeshTags< int > &  vertex_domains)

Set vertex domains.

Parameters
[in]vertex_domainsThe vertex domains.

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