ufc

Documentation for C++ code found in ffc/backends/ufc/*.h

UFC, Unified Form-assembly Code, is the interface between the form compiler, FFC, and the DOLFIN problem solving environment. You can, in principle, code your own finite elements and weak forms using the UFC interface to “save” you from using the UFL language to specify forms and elements and FFC to compile the forms to C++. In practice very few people implement code for the UFC interface by hand.

Enumerations

shape

C++ documentation for shape from ufc.h:

enum ufc::shape

Valid cell shapes.

enumerator ufc::shape::interval
enumerator ufc::shape::triangle
enumerator ufc::shape::quadrilateral
enumerator ufc::shape::tetrahedron
enumerator ufc::shape::hexahedron
enumerator ufc::shape::vertex

Classes

cell

C++ documentation for cell from ufc.h:

class ufc::cell

This class defines the data structure for a cell in a mesh.

shape ufc::cell::cell_shape

Shape of the cell.

std::vector<std::vector<std::size_t>> ufc::cell::entity_indices

Array of global indices for the mesh entities of the cell.

std::size_t ufc::cell::geometric_dimension

Geometric dimension of the mesh.

std::size_t ufc::cell::index

Cell index (short-cut for entity_indices[topological_dimension][0])

int ufc::cell::local_facet

Local facet index.

int ufc::cell::mesh_identifier

Unique mesh identifier.

int ufc::cell::orientation

Cell orientation.

std::size_t ufc::cell::topological_dimension

Topological dimension of the mesh.

cell_integral

C++ documentation for cell_integral from ufc.h:

class ufc::cell_integral

This class defines the interface for the tabulation of the cell tensor corresponding to the local contribution to a form from the integral over a cell.

void ufc::cell_integral::tabulate_tensor(double *A, const double *const *w, const double *coordinate_dofs, int cell_orientation) const = 0

Tabulate the tensor for the contribution from a local cell.

Parameters:
  • A
  • w
  • coordinate_dofs
  • cell_orientation
ufc::cell_integral::~cell_integral()

Destructor.

coordinate_mapping

C++ documentation for coordinate_mapping from ufc.h:

class ufc::coordinate_mapping

A representation of a coordinate mapping parameterized by a local finite element basis on each cell.

shape ufc::coordinate_mapping::cell_shape() const = 0

Return cell shape of the coordinate_mapping .

void ufc::coordinate_mapping::compute_geometry(double *x, double *J, double *detJ, double *K, std::size_t num_points, const double *X, const double *coordinate_dofs, int cell_orientation) const = 0

Combined (for convenience) computation of x, J, detJ, K from X and coordinate_dofs on a cell

Parameters:
  • x – Physical coordinates. Dimensions: x[num_points][gdim] [direction=out]
  • J – Jacobian of coordinate field, J = dx/dX. Dimensions: J[num_points][gdim][tdim] [direction=out]
  • detJ – (Pseudo-)Determinant of Jacobian. Dimensions: detJ[num_points] [direction=out]
  • K – (Pseudo-)Inverse of Jacobian of coordinate field. Dimensions: K[num_points][tdim][gdim] [direction=out]
  • num_points – Number of points. [direction=in]
  • X – Reference cell coordinates. Dimensions: X[num_points][tdim] [direction=in]
  • coordinate_dofs – Dofs of the coordinate field on the cell. Dimensions: coordinate_dofs[num_dofs][gdim]. [direction=in]
  • cell_orientation – Orientation of the cell, 1 means flipped w.r.t. reference cell. Only relevant on manifolds (tdim < gdim). [direction=in]
void ufc::coordinate_mapping::compute_jacobian_determinants(double *detJ, std::size_t num_points, const double *J, int cell_orientation) const = 0

Compute determinants of (pseudo-)Jacobians J

Parameters:
  • detJ – (Pseudo-)Determinant of Jacobian. Dimensions: detJ[num_points] [direction=out]
  • num_points – Number of points. [direction=in]
  • J – Jacobian of coordinate field, J = dx/dX. Dimensions: J[num_points][gdim][tdim] [direction=in]
  • cell_orientation – Orientation of the cell, 1 means flipped w.r.t. reference cell. Only relevant on manifolds (tdim < gdim). [direction=in]
void ufc::coordinate_mapping::compute_jacobian_inverses(double *K, std::size_t num_points, const double *J, const double *detJ) const = 0

Compute (pseudo-)inverses K of (pseudo-)Jacobians J

Parameters:
  • K – (Pseudo-)Inverse of Jacobian of coordinate field. Dimensions: K[num_points][tdim][gdim] [direction=out]
  • num_points – Number of points. [direction=in]
  • J – Jacobian of coordinate field, J = dx/dX. Dimensions: J[num_points][gdim][tdim] [direction=in]
  • detJ – (Pseudo-)Determinant of Jacobian. Dimensions: detJ[num_points] [direction=in]
void ufc::coordinate_mapping::compute_jacobians(double *J, std::size_t num_points, const double *X, const double *coordinate_dofs) const = 0

Compute Jacobian of coordinate mapping J = dx/dX at reference coordinates X

Parameters:
  • J – Jacobian of coordinate field, J = dx/dX. Dimensions: J[num_points][gdim][tdim] [direction=out]
  • num_points – Number of points. [direction=in]
  • X – Reference cell coordinates. Dimensions: X[num_points][tdim] [direction=in]
  • coordinate_dofs – Dofs of the coordinate field on the cell. Dimensions: coordinate_dofs[num_dofs][gdim]. [direction=in]
void ufc::coordinate_mapping::compute_midpoint_geometry(double *x, double *J, const double *coordinate_dofs) const = 0

Compute x and J at midpoint of cell

Parameters:
  • x – Physical coordinates. Dimensions: x[gdim] [direction=out]
  • J – Jacobian of coordinate field, J = dx/dX. Dimensions: J[gdim][tdim] [direction=out]
  • coordinate_dofs – Dofs of the coordinate field on the cell. Dimensions: coordinate_dofs[num_dofs][gdim]. [direction=in]
void ufc::coordinate_mapping::compute_physical_coordinates(double *x, std::size_t num_points, const double *X, const double *coordinate_dofs) const = 0

Compute physical coordinates x from reference coordinates X, the inverse of compute_reference_coordinates

Parameters:
  • x – Physical coordinates. Dimensions: x[num_points][gdim] [direction=out]
  • num_points – Number of points. [direction=in]
  • X – Reference cell coordinates. Dimensions: X[num_points][tdim] [direction=in]
  • coordinate_dofs – Dofs of the coordinate field on the cell. Dimensions: coordinate_dofs[num_dofs][gdim]. [direction=in]
void ufc::coordinate_mapping::compute_reference_coordinates(double *X, std::size_t num_points, const double *x, const double *coordinate_dofs, int cell_orientation) const = 0

Compute reference coordinates X from physical coordinates x, the inverse of compute_physical_coordinates

Parameters:
  • X – Reference cell coordinates. Dimensions: X[num_points][tdim] [direction=out]
  • num_points – Number of points. [direction=in]
  • x – Physical coordinates. Dimensions: x[num_points][gdim] [direction=in]
  • coordinate_dofs – Dofs of the coordinate field on the cell. Dimensions: coordinate_dofs[num_dofs][gdim]. [direction=in]
  • cell_orientation – Orientation of the cell, 1 means flipped w.r.t. reference cell. Only relevant on manifolds (tdim < gdim). [direction=in]
void ufc::coordinate_mapping::compute_reference_geometry(double *X, double *J, double *detJ, double *K, std::size_t num_points, const double *x, const double *coordinate_dofs, int cell_orientation) const = 0

Compute X, J, detJ, K from physical coordinates x on a cell

Parameters:
  • X – Reference cell coordinates. Dimensions: X[num_points][tdim] [direction=out]
  • J – Jacobian of coordinate field, J = dx/dX. Dimensions: J[num_points][gdim][tdim] [direction=out]
  • detJ – (Pseudo-)Determinant of Jacobian. Dimensions: detJ[num_points] [direction=out]
  • K – (Pseudo-)Inverse of Jacobian of coordinate field. Dimensions: K[num_points][tdim][gdim] [direction=out]
  • num_points – Number of points. [direction=in]
  • x – Physical coordinates. Dimensions: x[num_points][gdim] [direction=in]
  • coordinate_dofs – Dofs of the coordinate field on the cell. Dimensions: coordinate_dofs[num_dofs][gdim]. [direction=in]
  • cell_orientation – Orientation of the cell, 1 means flipped w.r.t. reference cell. Only relevant on manifolds (tdim < gdim). [direction=in]
coordinate_mapping *ufc::coordinate_mapping::create() const = 0

Create object of the same type.

dofmap *ufc::coordinate_mapping::create_coordinate_dofmap() const = 0

Create dofmap object representing the coordinate parameterization.

finite_element *ufc::coordinate_mapping::create_coordinate_finite_element() const = 0

Create finite_element object representing the coordinate parameterization.

std::size_t ufc::coordinate_mapping::geometric_dimension() const = 0

Return geometric dimension of the coordinate_mapping .

const char *ufc::coordinate_mapping::signature() const = 0

Return coordinate_mapping signature string.

std::size_t ufc::coordinate_mapping::topological_dimension() const = 0

Return topological dimension of the coordinate_mapping .

ufc::coordinate_mapping::~coordinate_mapping()

custom_integral

C++ documentation for custom_integral from ufc.h:

class ufc::custom_integral

This class defines the interface for the tabulation of the tensor corresponding to the local contribution to a form from the integral over a custom domain defined in terms of a set of quadrature points and weights.

ufc::custom_integral::custom_integral()

Constructor.

std::size_t ufc::custom_integral::num_cells() const = 0

Return the number of cells involved in evaluation of the integral.

void ufc::custom_integral::tabulate_tensor(double *A, const double *const *w, const double *coordinate_dofs, std::size_t num_quadrature_points, const double *quadrature_points, const double *quadrature_weights, const double *facet_normals, int cell_orientation) const = 0

Tabulate the tensor for the contribution from a custom domain.

Parameters:
  • A
  • w
  • coordinate_dofs
  • num_quadrature_points
  • quadrature_points
  • quadrature_weights
  • facet_normals
  • cell_orientation
ufc::custom_integral::~custom_integral()

Destructor.

cutcell_integral

C++ documentation for cutcell_integral from ufc.h:

class ufc::cutcell_integral

This class defines the interface for the tabulation of the tensor corresponding to the local contribution to a form from the integral over a cut cell defined in terms of a set of quadrature points and weights.

ufc::cutcell_integral::cutcell_integral()

Constructor.

void ufc::cutcell_integral::tabulate_tensor(double *A, const double *const *w, const double *coordinate_dofs, std::size_t num_quadrature_points, const double *quadrature_points, const double *quadrature_weights, int cell_orientation) const = 0

Tabulate the tensor for the contribution from a cutcell domain.

Parameters:
  • A
  • w
  • coordinate_dofs
  • num_quadrature_points
  • quadrature_points
  • quadrature_weights
  • cell_orientation
ufc::cutcell_integral::~cutcell_integral()

Destructor.

dofmap

C++ documentation for dofmap from ufc.h:

class ufc::dofmap

This class defines the interface for a local-to-global mapping of degrees of freedom (dofs).

dofmap *ufc::dofmap::create() const = 0

Create a new class instance.

dofmap *ufc::dofmap::create_sub_dofmap(std::size_t i) const = 0

Create a new dofmap for sub dofmap i (for a mixed element)

Parameters:i
std::size_t ufc::dofmap::global_dimension(const std::vector<std::size_t> &num_global_mesh_entities) const = 0

Return the dimension of the global finite element function space.

Parameters:num_global_mesh_entities
bool ufc::dofmap::needs_mesh_entities(std::size_t d) const = 0

Return true iff mesh entities of topological dimension d are needed

Parameters:d
std::size_t ufc::dofmap::num_element_dofs() const = 0

Return the dimension of the local finite element function space for a cell (old version including global support dofs)

std::size_t ufc::dofmap::num_element_support_dofs() const = 0

Return the dimension of the local finite element function space for a cell (not including global support dofs)

std::size_t ufc::dofmap::num_entity_closure_dofs(std::size_t d) const = 0

Return the number of dofs associated with the closure of each cell entity dimension d

Parameters:d
std::size_t ufc::dofmap::num_entity_dofs(std::size_t d) const = 0

Return the number of dofs associated with each cell entity of dimension d

Parameters:d
std::size_t ufc::dofmap::num_facet_dofs() const = 0

Return the number of dofs on each cell facet.

std::size_t ufc::dofmap::num_global_support_dofs() const = 0

Return the dimension of the local finite element function space Return the number of dofs with global support (i.e. global constants)

std::size_t ufc::dofmap::num_sub_dofmaps() const = 0

Return the number of sub dofmaps (for a mixed element)

const char *ufc::dofmap::signature() const = 0

Return a string identifying the dofmap.

void ufc::dofmap::tabulate_dofs(std::size_t *dofs, const std::vector<std::size_t> &num_global_entities, const std::vector<std::vector<std::size_t>> &entity_indices) const = 0

Tabulate the local-to-global mapping of dofs on a cell.

Parameters:
  • dofs
  • num_global_entities
  • entity_indices
void ufc::dofmap::tabulate_entity_closure_dofs(std::size_t *dofs, std::size_t d, std::size_t i) const = 0

Tabulate the local-to-local mapping of dofs on the closure of entity (d, i)

Parameters:
  • dofs
  • d
  • i
void ufc::dofmap::tabulate_entity_dofs(std::size_t *dofs, std::size_t d, std::size_t i) const = 0

Tabulate the local-to-local mapping of dofs on entity (d, i)

Parameters:
  • dofs
  • d
  • i
void ufc::dofmap::tabulate_facet_dofs(std::size_t *dofs, std::size_t facet) const = 0

Tabulate the local-to-local mapping from facet dofs to cell dofs.

Parameters:
  • dofs
  • facet
std::size_t ufc::dofmap::topological_dimension() const = 0

Return the topological dimension of the associated cell shape.

ufc::dofmap::~dofmap()

Destructor.

exterior_facet_integral

C++ documentation for exterior_facet_integral from ufc.h:

class ufc::exterior_facet_integral

This class defines the interface for the tabulation of the exterior facet tensor corresponding to the local contribution to a form from the integral over an exterior facet.

void ufc::exterior_facet_integral::tabulate_tensor(double *A, const double *const *w, const double *coordinate_dofs, std::size_t facet, int cell_orientation) const = 0

Tabulate the tensor for the contribution from a local exterior facet.

Parameters:
  • A
  • w
  • coordinate_dofs
  • facet
  • cell_orientation
ufc::exterior_facet_integral::~exterior_facet_integral()

Destructor.

finite_element

C++ documentation for finite_element from ufc.h:

class ufc::finite_element

This class defines the interface for a finite element.

shape ufc::finite_element::cell_shape() const = 0

Return the cell shape.

finite_element *ufc::finite_element::create() const = 0

Create a new class instance.

finite_element *ufc::finite_element::create_sub_element(std::size_t i) const = 0

Create a new finite element for sub element i (for a mixed element)

Parameters:i
std::size_t ufc::finite_element::degree() const = 0

Return the maximum polynomial degree of the finite element function space.

void ufc::finite_element::evaluate_basis(std::size_t i, double *values, const double *x, const double *coordinate_dofs, int cell_orientation, const ufc::coordinate_mapping *cm = nullptr) const = 0

Evaluate basis function i at given point x in cell.

Parameters:
  • i
  • values
  • x
  • coordinate_dofs
  • cell_orientation
  • cm
void ufc::finite_element::evaluate_basis_all(double *values, const double *x, const double *coordinate_dofs, int cell_orientation, const ufc::coordinate_mapping *cm = nullptr) const = 0

Evaluate all basis functions at given point x in cell.

Parameters:
  • values
  • x
  • coordinate_dofs
  • cell_orientation
  • cm
void ufc::finite_element::evaluate_basis_derivatives(std::size_t i, std::size_t n, double *values, const double *x, const double *coordinate_dofs, int cell_orientation, const ufc::coordinate_mapping *cm = nullptr) const = 0

Evaluate order n derivatives of basis function i at given point x in cell.

Parameters:
  • i
  • n
  • values
  • x
  • coordinate_dofs
  • cell_orientation
  • cm
void ufc::finite_element::evaluate_basis_derivatives_all(std::size_t n, double *values, const double *x, const double *coordinate_dofs, int cell_orientation, const ufc::coordinate_mapping *cm = nullptr) const = 0

Evaluate order n derivatives of all basis functions at given point x in cell.

Parameters:
  • n
  • values
  • x
  • coordinate_dofs
  • cell_orientation
  • cm
double ufc::finite_element::evaluate_dof(std::size_t i, const function &f, const double *coordinate_dofs, int cell_orientation, const cell &c, const ufc::coordinate_mapping *cm = nullptr) const = 0

Evaluate linear functional for dof i on the function f The cell argument is only included here so we can pass it to the function.

Parameters:
  • i
  • f
  • coordinate_dofs
  • cell_orientation
  • c
  • cm
void ufc::finite_element::evaluate_dofs(double *values, const function &f, const double *coordinate_dofs, int cell_orientation, const cell &c, const ufc::coordinate_mapping *cm = nullptr) const = 0

Evaluate linear functionals for all dofs on the function f The cell argument is only included here so we can pass it to the function.

Parameters:
  • values
  • f
  • coordinate_dofs
  • cell_orientation
  • c
  • cm
void ufc::finite_element::evaluate_reference_basis(double *reference_values, std::size_t num_points, const double *X) const = 0

Evaluate all basis functions at given point X in reference cell

Parameters:
  • reference_values – Basis function values on reference element. Dimensions: reference_values[num_points][num_dofs][reference_value_size] [direction=out]
  • num_points – Number of points. [direction=in]
  • X – Reference cell coordinates. Dimensions: X[num_points][tdim] [direction=in]
void ufc::finite_element::evaluate_reference_basis_derivatives(double *reference_values, std::size_t order, std::size_t num_points, const double *X) const = 0

Evaluate specific order derivatives of all basis functions at given point X in reference cell

Parameters:
  • reference_values – Basis function derivative values on reference element. Dimensions: reference_values[num_points][num_dofs][num_derivatives][reference_value_size] where num_derivatives = pow(order, tdim). TODO: Document ordering of derivatives for order > 1. [direction=out]
  • order – Derivative order to compute. [direction=in]
  • num_points – Number of points. [direction=in]
  • X – Reference cell coordinates. Dimensions: X[num_points][tdim] [direction=in]
const char *ufc::finite_element::family() const = 0

Return the family of the finite element function space.

std::size_t ufc::finite_element::geometric_dimension() const = 0

Return the geometric dimension of the cell shape.

void ufc::finite_element::interpolate_vertex_values(double *vertex_values, const double *dof_values, const double *coordinate_dofs, int cell_orientation, const ufc::coordinate_mapping *cm = nullptr) const = 0

Interpolate vertex values from dof values.

Parameters:
  • vertex_values
  • dof_values
  • coordinate_dofs
  • cell_orientation
  • cm
std::size_t ufc::finite_element::num_sub_elements() const = 0

Return the number of sub elements (for a mixed element)

std::size_t ufc::finite_element::reference_value_dimension(std::size_t i) const = 0

Return the dimension of the reference value space for axis i.

Parameters:i
std::size_t ufc::finite_element::reference_value_rank() const = 0

Return the rank of the reference value space.

std::size_t ufc::finite_element::reference_value_size() const = 0

Return the number of components of the reference value space.

const char *ufc::finite_element::signature() const = 0

Return a string identifying the finite element.

std::size_t ufc::finite_element::space_dimension() const = 0

Return the dimension of the finite element function space.

void ufc::finite_element::tabulate_dof_coordinates(double *dof_coordinates, const double *coordinate_dofs, const ufc::coordinate_mapping *cm = nullptr) const = 0

Tabulate the coordinates of all dofs on a cell.

Parameters:
  • dof_coordinates
  • coordinate_dofs
  • cm
void ufc::finite_element::tabulate_reference_dof_coordinates(double *reference_dof_coordinates) const = 0

Tabulate the coordinates of all dofs on a reference cell.

Parameters:reference_dof_coordinates
std::size_t ufc::finite_element::topological_dimension() const = 0

Return the topological dimension of the cell shape.

void ufc::finite_element::transform_reference_basis_derivatives(double *values, std::size_t order, std::size_t num_points, const double *reference_values, const double *X, const double *J, const double *detJ, const double *K, int cell_orientation) const = 0

Transform order n derivatives (can be 0) of all basis functions previously evaluated in points X in reference cell with given Jacobian J and its inverse K for each point

Parameters:
  • values – Transformed basis function (derivative) values. Dimensions: values[num_points][num_dofs][num_derivatives][value_size] where num_derivatives = pow(order, tdim). TODO: Document ordering of derivatives for order > 1. [direction=out]
  • order – Derivative order to compute. Can be zero to just apply Piola mappings. [direction=in]
  • num_points – Number of points. [direction=in]
  • reference_values – Basis function derivative values on reference element. Dimensions: reference_values[num_points][num_dofs][num_derivatives][reference_value_size] where num_derivatives = pow(order, tdim). TODO: Document ordering of derivatives for order > 1. [direction=in]
  • X – Reference cell coordinates. Dimensions: X[num_points][tdim] [direction=in]
  • J – Jacobian of coordinate field, J = dx/dX. Dimensions: J[num_points][gdim][tdim] [direction=in]
  • detJ – (Pseudo-)Determinant of Jacobian. Dimensions: detJ[num_points] [direction=in]
  • K – (Pseudo-)Inverse of Jacobian of coordinate field. Dimensions: K[num_points][tdim][gdim] [direction=in]
  • cell_orientation – Orientation of the cell, 1 means flipped w.r.t. reference cell. Only relevant on manifolds (tdim < gdim). [direction=in]
std::size_t ufc::finite_element::value_dimension(std::size_t i) const = 0

Return the dimension of the value space for axis i.

Parameters:i
std::size_t ufc::finite_element::value_rank() const = 0

Return the rank of the value space.

std::size_t ufc::finite_element::value_size() const = 0

Return the number of components of the value space.

ufc::finite_element::~finite_element()

Destructor.

form

C++ documentation for form from ufc.h:

class ufc::form

This class defines the interface for the assembly of the global tensor corresponding to a form with r + n arguments, that is, a mapping

a : V1 x V2 x ... Vr x W1 x W2 x ... x Wn -> R

with arguments v1, v2, ..., vr, w1, w2, ..., wn. The rank r global tensor A is defined by

A = a(V1, V2, ..., Vr, w1, w2, ..., wn),

where each argument Vj represents the application to the sequence of basis functions of Vj and w1, w2, ..., wn are given fixed functions (coefficients).

cell_integral *ufc::form::create_cell_integral(std::size_t subdomain_id) const = 0

Create a new cell integral on sub domain subdomain_id.

Parameters:subdomain_id
dofmap *ufc::form::create_coordinate_dofmap() const = 0

Create a new dofmap for parameterization of coordinates.

finite_element *ufc::form::create_coordinate_finite_element() const = 0

Create a new finite element for parameterization of coordinates.

coordinate_mapping *ufc::form::create_coordinate_mapping() const = 0

Create a new coordinate mapping.

custom_integral *ufc::form::create_custom_integral(std::size_t subdomain_id) const = 0

Create a new custom integral on sub domain subdomain_id.

Parameters:subdomain_id
cutcell_integral *ufc::form::create_cutcell_integral(std::size_t subdomain_id) const = 0

Create a new cutcell integral on sub domain subdomain_id.

Parameters:subdomain_id
cell_integral *ufc::form::create_default_cell_integral() const = 0

Create a new cell integral on everywhere else.

custom_integral *ufc::form::create_default_custom_integral() const = 0

Create a new custom integral on everywhere else.

cutcell_integral *ufc::form::create_default_cutcell_integral() const = 0

Create a new cutcell integral on everywhere else.

exterior_facet_integral *ufc::form::create_default_exterior_facet_integral() const = 0

Create a new exterior facet integral on everywhere else.

interface_integral *ufc::form::create_default_interface_integral() const = 0

Create a new interface integral on everywhere else.

interior_facet_integral *ufc::form::create_default_interior_facet_integral() const = 0

Create a new interior facet integral on everywhere else.

overlap_integral *ufc::form::create_default_overlap_integral() const = 0

Create a new overlap integral on everywhere else.

vertex_integral *ufc::form::create_default_vertex_integral() const = 0

Create a new vertex integral on everywhere else.

dofmap *ufc::form::create_dofmap(std::size_t i) const = 0

Create a new dofmap for argument function 0 <= i < r+n

Parameters:i – Argument number if 0 <= i < r Coefficient number j=i-r if r+j <= i < r+n
exterior_facet_integral *ufc::form::create_exterior_facet_integral(std::size_t subdomain_id) const = 0

Create a new exterior facet integral on sub domain subdomain_id.

Parameters:subdomain_id
finite_element *ufc::form::create_finite_element(std::size_t i) const = 0

Create a new finite element for argument function 0 <= i < r+n

Parameters:i – Argument number if 0 <= i < r Coefficient number j=i-r if r+j <= i < r+n
interface_integral *ufc::form::create_interface_integral(std::size_t subdomain_id) const = 0

Create a new interface integral on sub domain subdomain_id.

Parameters:subdomain_id
interior_facet_integral *ufc::form::create_interior_facet_integral(std::size_t subdomain_id) const = 0

Create a new interior facet integral on sub domain subdomain_id.

Parameters:subdomain_id
overlap_integral *ufc::form::create_overlap_integral(std::size_t subdomain_id) const = 0

Create a new overlap integral on sub domain subdomain_id.

Parameters:subdomain_id
vertex_integral *ufc::form::create_vertex_integral(std::size_t subdomain_id) const = 0

Create a new vertex integral on sub domain subdomain_id.

Parameters:subdomain_id
bool ufc::form::has_cell_integrals() const = 0

Return whether form has any cell integrals.

bool ufc::form::has_custom_integrals() const = 0

Return whether form has any custom integrals.

bool ufc::form::has_cutcell_integrals() const = 0

Return whether form has any cutcell integrals.

bool ufc::form::has_exterior_facet_integrals() const = 0

Return whether form has any exterior facet integrals.

bool ufc::form::has_interface_integrals() const = 0

Return whether form has any interface integrals.

bool ufc::form::has_interior_facet_integrals() const = 0

Return whether form has any interior facet integrals.

bool ufc::form::has_overlap_integrals() const = 0

Return whether form has any overlap integrals.

bool ufc::form::has_vertex_integrals() const = 0

Return whether form has any vertex integrals.

std::size_t ufc::form::max_cell_subdomain_id() const = 0

Return the upper bound on subdomain ids for cell integrals.

std::size_t ufc::form::max_custom_subdomain_id() const = 0

Return the upper bound on subdomain ids for custom integrals.

std::size_t ufc::form::max_cutcell_subdomain_id() const = 0

Return the upper bound on subdomain ids for cutcell integrals.

std::size_t ufc::form::max_exterior_facet_subdomain_id() const = 0

Return the upper bound on subdomain ids for exterior facet integrals.

std::size_t ufc::form::max_interface_subdomain_id() const = 0

Return the upper bound on subdomain ids for interface integrals.

std::size_t ufc::form::max_interior_facet_subdomain_id() const = 0

Return the upper bound on subdomain ids for interior facet integrals.

std::size_t ufc::form::max_overlap_subdomain_id() const = 0

Return the upper bound on subdomain ids for overlap integrals.

std::size_t ufc::form::max_vertex_subdomain_id() const = 0

Return the upper bound on subdomain ids for vertex integrals.

std::size_t ufc::form::num_coefficients() const = 0

Return the number of coefficients (n)

std::size_t ufc::form::original_coefficient_position(std::size_t i) const = 0

Return original coefficient position for each coefficient

Parameters:i – Coefficient number, 0 <= i < n
std::size_t ufc::form::rank() const = 0

Return the rank of the global tensor (r)

const char *ufc::form::signature() const = 0

Return a string identifying the form.

ufc::form::~form()

Destructor.

function

C++ documentation for function from ufc.h:

class ufc::function

This class defines the interface for a general tensor-valued function.

void ufc::function::evaluate(double *values, const double *coordinates, const cell &c) const = 0

Evaluate function at given point in cell.

Parameters:
  • values
  • coordinates
  • c
ufc::function::~function()

Destructor.

integral

C++ documentation for integral from ufc.h:

class ufc::integral

This class defines the shared interface for classes implementing the tabulation of a tensor corresponding to the local contribution to a form from an integral.

const std::vector<bool> &ufc::integral::enabled_coefficients() const = 0

Tabulate which form coefficients are used by this integral.

ufc::integral::~integral()

Destructor.

interface_integral

C++ documentation for interface_integral from ufc.h:

class ufc::interface_integral

This class defines the interface for the tabulation of the tensor corresponding to the local contribution to a form from the integral over a cut cell defined in terms of a set of quadrature points and weights.

ufc::interface_integral::interface_integral()

Constructor.

void ufc::interface_integral::tabulate_tensor(double *A, const double *const *w, const double *coordinate_dofs, std::size_t num_quadrature_points, const double *quadrature_points, const double *quadrature_weights, const double *facet_normals, int cell_orientation) const = 0

Tabulate the tensor for the contribution from an interface domain.

Parameters:
  • A
  • w
  • coordinate_dofs
  • num_quadrature_points
  • quadrature_points
  • quadrature_weights
  • facet_normals
  • cell_orientation
ufc::interface_integral::~interface_integral()

Destructor.

interior_facet_integral

C++ documentation for interior_facet_integral from ufc.h:

class ufc::interior_facet_integral

This class defines the interface for the tabulation of the interior facet tensor corresponding to the local contribution to a form from the integral over an interior facet.

void ufc::interior_facet_integral::tabulate_tensor(double *A, const double *const *w, const double *coordinate_dofs_0, const double *coordinate_dofs_1, std::size_t facet_0, std::size_t facet_1, int cell_orientation_0, int cell_orientation_1) const = 0

Tabulate the tensor for the contribution from a local interior facet.

Parameters:
  • A
  • w
  • coordinate_dofs_0
  • coordinate_dofs_1
  • facet_0
  • facet_1
  • cell_orientation_0
  • cell_orientation_1
ufc::interior_facet_integral::~interior_facet_integral()

Destructor.

overlap_integral

C++ documentation for overlap_integral from ufc.h:

class ufc::overlap_integral

This class defines the interface for the tabulation of the tensor corresponding to the local contribution to a form from the integral over the overlapped portion of a cell defined in terms of a set of quadrature points and weights.

ufc::overlap_integral::overlap_integral()

Constructor.

void ufc::overlap_integral::tabulate_tensor(double *A, const double *const *w, const double *coordinate_dofs, std::size_t num_quadrature_points, const double *quadrature_points, const double *quadrature_weights, int cell_orientation) const = 0

Tabulate the tensor for the contribution from an overlap domain.

Parameters:
  • A
  • w
  • coordinate_dofs
  • num_quadrature_points
  • quadrature_points
  • quadrature_weights
  • cell_orientation
ufc::overlap_integral::~overlap_integral()

Destructor.

vertex_integral

C++ documentation for vertex_integral from ufc.h:

class ufc::vertex_integral

This class defines the interface for the tabulation of an expression evaluated at exactly one point.

void ufc::vertex_integral::tabulate_tensor(double *A, const double *const *w, const double *coordinate_dofs, std::size_t vertex, int cell_orientation) const = 0

Tabulate the tensor for the contribution from the local vertex.

Parameters:
  • A
  • w
  • coordinate_dofs
  • vertex
  • cell_orientation
ufc::vertex_integral::vertex_integral()

Constructor.

ufc::vertex_integral::~vertex_integral()

Destructor.