DOLFIN-X
DOLFIN-X C++ interface
|
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::Mesh > | mesh () const |
Extract common mesh from form. More... | |
std::shared_ptr< const function::FunctionSpace > | function_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... | |
FormCoefficients & | coefficients () |
Access coefficients. | |
const FormCoefficients & | coefficients () const |
Access coefficients. | |
const FormIntegrals & | integrals () const |
Access form integrals. | |
const std::vector< std::pair< std::string, std::shared_ptr< const function::Constant > > > & | constants () const |
Access constants. More... | |
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).
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.
[in] | function_spaces | Function Spaces |
[in] | integrals | |
[in] | coefficients | |
[in] | constants | Vector of pairs (name, constant). The index in the vector is the position of the constant in the original (nonsimplified) form. |
|
explicit |
Create form (no UFC integrals). Integrals can be attached later using FormIntegrals::set_cell_tabulate_tensor.
[in] | function_spaces | Vector of function spaces |
bool Form::all_constants_set | ( | ) | const |
Check if all constants associated with the form have been set.
const std::vector< std::pair< std::string, std::shared_ptr< const function::Constant > > > & Form::constants | ( | ) | const |
Access constants.
std::shared_ptr< const function::FunctionSpace > Form::function_space | ( | int | i | ) | const |
Return function space for given argument.
[in] | i | Index of the argument |
std::set< std::string > Form::get_unset_constants | ( | ) | const |
Return names of any constants that have not been set.
std::shared_ptr< const mesh::Mesh > Form::mesh | ( | ) | const |
Extract common mesh from form.
int Form::original_coefficient_position | ( | int | i | ) | const |
Return original coefficient position for each coefficient (0 <= i < n)
int Form::rank | ( | ) | const |
Rank of form (bilinear form = 2, linear form = 1, functional = 0, etc)
void Form::set_cell_domains | ( | const mesh::MeshTags< int > & | cell_domains | ) |
Set cell domains.
[in] | cell_domains | The cell domains |
void Form::set_coefficients | ( | std::map< std::size_t, std::shared_ptr< const function::Function >> | coefficients | ) |
Set coefficient with given number (shared pointer version)
[in] | coefficients | Map from coefficient index to the coefficient |
void Form::set_coefficients | ( | std::map< std::string, std::shared_ptr< const function::Function >> | coefficients | ) |
Set coefficient with given name (shared pointer version)
[in] | coefficients | Map from coefficient name to the coefficient |
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.
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.
void Form::set_exterior_facet_domains | ( | const mesh::MeshTags< int > & | exterior_facet_domains | ) |
Set exterior facet domains.
[in] | exterior_facet_domains | The exterior facet domains |
void Form::set_interior_facet_domains | ( | const mesh::MeshTags< int > & | interior_facet_domains | ) |
Set interior facet domains.
[in] | interior_facet_domains | The interior facet domains |
void Form::set_mesh | ( | std::shared_ptr< const mesh::Mesh > | mesh | ) |
Set mesh, necessary for functionals when there are no function spaces.
[in] | mesh | The mesh |
void Form::set_vertex_domains | ( | const mesh::MeshTags< int > & | vertex_domains | ) |
Set vertex domains.
[in] | vertex_domains | The vertex domains. |