DOLFIN
DOLFIN C++ interface
Loading...
Searching...
No Matches
dolfin::Form Class Reference

Base class for UFC code generated by FFC for DOLFIN with option -l. More...

#include <Form.h>

Inheritance diagram for dolfin::Form:
Collaboration diagram for dolfin::Form:

Public Member Functions

 Form (std::size_t rank, std::size_t num_coefficients)
 Form (std::shared_ptr< const ufc::form > ufc_form, std::vector< std::shared_ptr< const FunctionSpace > > function_spaces)
virtual ~Form ()
 Destructor.
std::size_t rank () const
std::size_t num_coefficients () const
std::size_t original_coefficient_position (std::size_t i) const
std::vector< std::size_t > coloring (std::size_t entity_dim) const
void set_mesh (std::shared_ptr< const Mesh > mesh)
std::shared_ptr< const Meshmesh () const
std::shared_ptr< const FunctionSpacefunction_space (std::size_t i) const
std::vector< std::shared_ptr< const FunctionSpace > > function_spaces () const
void set_coefficient (std::size_t i, std::shared_ptr< const GenericFunction > coefficient)
void set_coefficient (std::string name, std::shared_ptr< const GenericFunction > coefficient)
void set_coefficients (std::map< std::string, std::shared_ptr< const GenericFunction > > coefficients)
void set_some_coefficients (std::map< std::string, std::shared_ptr< const GenericFunction > > coefficients)
std::shared_ptr< const GenericFunctioncoefficient (std::size_t i) const
std::shared_ptr< const GenericFunctioncoefficient (std::string name) const
std::vector< std::shared_ptr< const GenericFunction > > coefficients () const
virtual std::size_t coefficient_number (const std::string &name) const
virtual std::string coefficient_name (std::size_t i) const
std::shared_ptr< const MeshFunction< std::size_t > > cell_domains () const
std::shared_ptr< const MeshFunction< std::size_t > > exterior_facet_domains () const
std::shared_ptr< const MeshFunction< std::size_t > > interior_facet_domains () const
std::shared_ptr< const MeshFunction< std::size_t > > vertex_domains () const
void set_cell_domains (std::shared_ptr< const MeshFunction< std::size_t > > cell_domains)
void set_exterior_facet_domains (std::shared_ptr< const MeshFunction< std::size_t > > exterior_facet_domains)
void set_interior_facet_domains (std::shared_ptr< const MeshFunction< std::size_t > > interior_facet_domains)
void set_vertex_domains (std::shared_ptr< const MeshFunction< std::size_t > > vertex_domains)
std::shared_ptr< const ufc::form > ufc_form () const
void check () const
 Check function spaces and coefficients.
Equation operator== (const Form &rhs) const
 Comparison operator, returning equation lhs == rhs.
Equation operator== (int rhs) const
 Comparison operator, returning equation lhs == 0.
Public Member Functions inherited from dolfin::Hierarchical< Form >
 Hierarchical (Form &self)
 Constructor.
virtual ~Hierarchical ()
 Destructor.
std::size_t depth () const
bool has_parent () const
bool has_child () const
Formparent ()
std::shared_ptr< Formparent_shared_ptr ()
Formchild ()
std::shared_ptr< Formchild_shared_ptr ()
Formroot_node ()
std::shared_ptr< Formroot_node_shared_ptr ()
Formleaf_node ()
std::shared_ptr< Formleaf_node_shared_ptr ()
void set_parent (std::shared_ptr< Form > parent)
 Set parent.
void clear_child ()
 Clear child.
void set_child (std::shared_ptr< Form > child)
 Set child.
const Hierarchicaloperator= (const Hierarchical &hierarchical)
 Assignment operator.
void _debug () const
 Function useful for debugging the hierarchy.

Public Attributes

std::shared_ptr< const MeshFunction< std::size_t > > dx
 Domain markers for cells.
std::shared_ptr< const MeshFunction< std::size_t > > ds
 Domain markers for exterior facets.
std::shared_ptr< const MeshFunction< std::size_t > > dS
 Domain markers for interior facets.
std::shared_ptr< const MeshFunction< std::size_t > > dP
 Domain markers for vertices.

Protected Attributes

std::shared_ptr< const ufc::form > _ufc_form
std::vector< std::shared_ptr< const FunctionSpace > > _function_spaces
std::vector< std::shared_ptr< const GenericFunction > > _coefficients
std::shared_ptr< const Mesh_mesh

Detailed Description

Base class for UFC code generated by FFC for DOLFIN with option -l.

embed:rst:leading-slashes
/// 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
///
///     a: V_1 x V_0 -> R
///
/// .. note::
///
///     Figure out how to write this in math mode without it getting
///     messed up in the Python version.
///
/// 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
///
/// .. code-block:: c++
///
///     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 ( std::size_t rank,
std::size_t num_coefficients )

Create form of given rank with given number of coefficients

Parameters
[in]rank(std::size_t) The rank.
[in]num_coefficients(std::size_t) The number of coefficients.

◆ Form() [2/2]

Form::Form ( std::shared_ptr< const ufc::form > ufc_form,
std::vector< std::shared_ptr< const FunctionSpace > > function_spaces )

Create form (shared data)

Parameters
[in]ufc_form(ufc::form) The UFC form.
[in]function_spaces(std::vector<FunctionSpace>) Vector of function spaces.

Member Function Documentation

◆ cell_domains()

std::shared_ptr< const MeshFunction< std::size_t > > Form::cell_domains ( ) const

Return cell domains (zero pointer if no domains have been specified)

Returns
MeshFunction <std::size_t> The cell domains.

◆ coefficient() [1/2]

std::shared_ptr< const GenericFunction > Form::coefficient ( std::size_t i) const

Return coefficient with given number

Parameters
[in]i(std::size_t) Index
Returns
GenericFunction The coefficient.

◆ coefficient() [2/2]

std::shared_ptr< const GenericFunction > Form::coefficient ( std::string name) const

Return coefficient with given name

Parameters
[in]name(std::string) The name.
Returns
GenericFunction The coefficient.

◆ coefficient_name()

std::string Form::coefficient_name ( std::size_t i) const
virtual

Return the name of the coefficient with this number

Parameters
[in]i(std::size_t) The number
Returns
std::string The name of the coefficient with the given number.

◆ coefficient_number()

std::size_t Form::coefficient_number ( const std::string & name) const
virtual

Return the number of the coefficient with this name

Parameters
[in]name(std::string) The name.
Returns
std::size_t The number of the coefficient with the given name.

◆ coefficients()

std::vector< std::shared_ptr< const GenericFunction > > Form::coefficients ( ) const

Return all coefficients

Returns
std::vector<GenericFunction> All coefficients.

◆ coloring()

std::vector< std::size_t > Form::coloring ( std::size_t entity_dim) const

Return coloring type for colored assembly of form over a mesh entity of a given dimension

Parameters
[in]entity_dim(std::size_t) Dimension.
Returns
std::vector<std::size_t> Coloring type.

◆ exterior_facet_domains()

std::shared_ptr< const MeshFunction< std::size_t > > Form::exterior_facet_domains ( ) const

Return exterior facet domains (zero pointer if no domains have been specified)

Returns
std::shared_ptr<MeshFunction <std::size_t>> The exterior facet domains.

◆ function_space()

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

Return function space for given argument

Parameters
i(std::size_t) Index
Returns
FunctionSpace Function space shared pointer.

◆ function_spaces()

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

Return function spaces for arguments

Returns
std::vector<FunctionSpace> Vector of function space shared pointers.

◆ interior_facet_domains()

std::shared_ptr< const MeshFunction< std::size_t > > Form::interior_facet_domains ( ) const

Return interior facet domains (zero pointer if no domains have been specified)

Returns
MeshFunction <std::size_t> The interior facet domains.

◆ mesh()

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

Extract common mesh from form

Returns
Mesh Shared pointer to the mesh.

◆ num_coefficients()

std::size_t Form::num_coefficients ( ) const

Return number of coefficients

Returns
std::size_t The number of coefficients.

◆ original_coefficient_position()

std::size_t Form::original_coefficient_position ( std::size_t i) const

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

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

◆ rank()

std::size_t Form::rank ( ) const

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

Returns
std::size_t The rank of the form.

◆ set_cell_domains()

void Form::set_cell_domains ( std::shared_ptr< const MeshFunction< std::size_t > > cell_domains)

Set cell domains

Parameters
[in]cell_domains(MeshFunction <std::size_t>) The cell domains.

◆ set_coefficient() [1/2]

void Form::set_coefficient ( std::size_t i,
std::shared_ptr< const GenericFunction > coefficient )

Set coefficient with given number (shared pointer version)

Parameters
[in]i(std::size_t) The given number.
[in]coefficient(GenericFunction) The coefficient.

◆ set_coefficient() [2/2]

void Form::set_coefficient ( std::string name,
std::shared_ptr< const GenericFunction > coefficient )

Set coefficient with given name (shared pointer version)

Parameters
[in]name(std::string) The name.
[in]coefficient(GenericFunction) The coefficient.

◆ set_coefficients()

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

Set all coefficients in given map. All coefficients in the given map, which may contain only a subset of the coefficients of the form, will be set.

Parameters
[in]coefficients(std::map<std::string, GenericFunction>) The map of coefficients.

◆ set_exterior_facet_domains()

void Form::set_exterior_facet_domains ( std::shared_ptr< const MeshFunction< std::size_t > > exterior_facet_domains)

Set exterior facet domains

Parameters
[in]exterior_facet_domains(MeshFunction <std::size_t>) The exterior facet domains.

◆ set_interior_facet_domains()

void Form::set_interior_facet_domains ( std::shared_ptr< const MeshFunction< std::size_t > > interior_facet_domains)

Set interior facet domains

Parameters
[in]interior_facet_domains(MeshFunction <std::size_t>) The interior facet domains.

◆ set_mesh()

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

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

Parameters
[in]mesh(Mesh) The mesh.

◆ set_some_coefficients()

void Form::set_some_coefficients ( std::map< std::string, std::shared_ptr< const GenericFunction > > coefficients)

Set some coefficients in given map. Each coefficient in the given map will be set, if the name of the coefficient matches the name of a coefficient in the form.

This is useful when reusing the same coefficient map for several forms, or when some part of the form has been outcommented (for testing) in the UFL file, which means that the coefficient and attaching it to the form does not need to be outcommented in a C++ program using code from the generated UFL file.

Parameters
[in]coefficients(std::map<std::string, GenericFunction>) The map of coefficients.

◆ set_vertex_domains()

void Form::set_vertex_domains ( std::shared_ptr< const MeshFunction< std::size_t > > vertex_domains)

Set vertex domains

Parameters
[in]vertex_domains(MeshFunction <std::size_t>) The vertex domains.

◆ ufc_form()

std::shared_ptr< const ufc::form > Form::ufc_form ( ) const

Return UFC form shared pointer

Returns
ufc::form The UFC form.

◆ vertex_domains()

std::shared_ptr< const MeshFunction< std::size_t > > Form::vertex_domains ( ) const

Return vertex domains (zero pointer if no domains have been specified)

Returns
MeshFunction <std::size_t> The vertex domains.

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