Expression.h

Note

The documentation on this page was automatically extracted from the DOLFIN C++ code and may need to be edited or expanded.

class Expression

Parent class(es)

This class represents a user-defined expression. Expressions can be used as coefficients in variational forms or interpolated into finite element spaces.

An expression is defined by overloading the eval() method. Users may choose to overload either a simple version of eval(), in the case of expressions only depending on the coordinate x, or an optional version for expressions depending on x and mesh data like cell indices or facet normals.

The geometric dimension (the size of x) and the value rank and dimensions of an expression must supplied as arguments to the constructor.

Expression()

Create scalar expression.

explicit Expression(std::size_t dim)

Create vector-valued expression with given dimension.

@param dim (std::size_t)
Dimension of the vector-valued expression.
Expression(std::size_t dim0, std::size_t dim1)

Create matrix-valued expression with given dimensions.

@param dim0 (std::size_t)
Dimension (rows).
@param dim1 (std::size_t)
Dimension (columns).
explicit Expression(std::vector<std::size_t> value_shape)

Create tensor-valued expression with given shape.

@param value_shape (std::vector<std::size_t>)
Shape of expression.
Expression(const Expression &expression)

Copy constructor

@param expression (Expression)
Object to be copied.
void eval(Array<double> &values, const Array<double> &x, const ufc::cell &cell) const override

Evaluate at given point in given cell (deprecated)

@param values (Array<double>)
The values at the point.
@param x (Array<double>)
The coordinates of the point.
@param cell (ufc::cell)
The cell which contains the given point.
void eval(Eigen::Ref<Eigen::VectorXd> values, Eigen::Ref<const Eigen::VectorXd> x, const ufc::cell &cell) const override

Evaluate at given point in given cell

@param values (Eigen::Ref<Eigen::VectorXd>)
The values at the point.
@param x (Eigen::Ref<const Eigen::VectorXd>)
The coordinates of the point.
@param cell (ufc::cell)
The cell which contains the given point.
void eval(Array<double> &values, const Array<double> &x) const override

Evaluate at given point (deprecated)

@param values (Array<double>)
The values at the point.
@param x (Array<double>)
The coordinates of the point.
void eval(Eigen::Ref<Eigen::VectorXd> values, Eigen::Ref<const Eigen::VectorXd> x) const override

Evaluate at given point.

@param values (Eigen::Ref<Eigen::VectorXd>)
The values at the point.
@param x (Eigen::Ref<const Eigen::VectorXd>)
The coordinates of the point.
std::size_t value_rank() const override

Return value rank.

@return std::size_t
The value rank.
std::size_t value_dimension(std::size_t i) const override

Return value dimension for given axis.

@param i (std::size_t)
Integer denoting the axis to use.
@return std::size_t
The value dimension (for the given axis).
std::vector<std::size_t> value_shape() const override

Return value shape

@return std::vector<std::size_t>
The value shape.
void set_property(std::string name, double value)

Property setter for type “double” Used in pybind11 Python interface to attach a value to a python attribute

double get_property(std::string name) const

Property getter for type “double” Used in pybind11 Python interface to get the value of a python attribute

void set_generic_function(std::string name, std::shared_ptr<GenericFunction> f)

Property setter for type “GenericFunction” Used in pybind11 Python interface to attach a value to a python attribute

std::shared_ptr<dolfin::GenericFunction> get_generic_function(std::string name) const

Property getter for type “GenericFunction” Used in pybind11 Python interface to get the value of a python attribute

void restrict(double *w, const FiniteElement &element, const Cell &dolfin_cell, const double *coordinate_dofs, const ufc::cell &ufc_cell) const override

Restrict function to local cell (compute expansion coefficients w).

@param w (list of doubles)
Expansion coefficients.
@param element (FiniteElement)
The element.
@param dolfin_cell (Cell)
The cell.
@param coordinate_dofs (double*)
The coordinates
@param ufc_cell (ufc::cell)
The ufc::cell.
void compute_vertex_values(std::vector<double> &vertex_values, const Mesh &mesh) const override

Compute values at all mesh vertices.

@param vertex_values (Array<double>)
The values at all vertices.
@param mesh (Mesh)
The mesh.
std::shared_ptr<const FunctionSpace> function_space() const override

Return shared pointer to function space (NULL) Expression does not have a FunctionSpace

@return FunctionSpace
Return the shared pointer.