GenericFunction.h

Note

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

class GenericFunction

Parent class(es)

This is a common base class for functions. Functions can be evaluated at a given point and they can be restricted to a given cell in a finite element mesh. This functionality is implemented by sub-classes that implement the eval() and restrict() functions.

DOLFIN provides two implementations of the GenericFunction interface in the form of the classes Function and Expression.

Sub-classes may optionally implement the update() function that will be called prior to restriction when running in parallel.

GenericFunction()

Constructor

std::size_t value_rank() const = 0

Return value rank

std::size_t value_dimension(std::size_t i) const = 0

Return value dimension for given axis

std::vector<std::size_t> value_shape() const = 0

Return value shape

void eval(Array<double> &values, const Array<double> &x, const ufc::cell &cell) const

Evaluate at given point in given cell (deprecated)

void eval(Array<double> &values, const Array<double> &x) const

Evaluate at given point (deprecated)

void eval(Eigen::Ref<Eigen::VectorXd> values, Eigen::Ref<const Eigen::VectorXd> x, const ufc::cell &cell) const

Evaluate at given point in given cell

void eval(Eigen::Ref<Eigen::VectorXd> values, Eigen::Ref<const Eigen::VectorXd> x) const

Evaluate at given point

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

Restrict function to local cell (compute expansion coefficients w)

void compute_vertex_values(std::vector<double> &vertex_values, const Mesh &mesh) const = 0

Compute values at all mesh vertices

void update() const

Update off-process ghost coefficients

double operator()(double x) const

Evaluation at given point (scalar function)

double operator()(double x, double y) const

Evaluation at given point (scalar function)

double operator()(double x, double y, double z) const

Evaluation at given point (scalar function)

double operator()(const Point &p) const

Evaluation at given point (scalar function)

void operator()(Array<double> &values, double x) const

Evaluation at given point (vector-valued function)

void operator()(Array<double> &values, double x, double y) const

Evaluation at given point (vector-valued function)

void operator()(Array<double> &values, double x, double y, double z) const

Evaluation at given point (vector-valued function)

void operator()(Array<double> &values, const Point &p) const

Evaluation at given point (vector-valued function)

std::size_t value_size() const

Evaluation at given point Return value size (product of value dimensions)

void evaluate(double *values, const double *coordinates, const ufc::cell &cell) const override

Evaluate function at given point in cell @param values (double*) @param coordinates (const double*) @param cell (ufc::cell&)

std::shared_ptr<const FunctionSpace> function_space() const = 0

Pointer to FunctionSpace, if appropriate, otherwise NULL

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

Restrict as UFC function (by calling eval)