Function.h

Note

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

class Function

Parent class(es)

This class represents a function \(u_h\) in a finite element function space \(V_h\), given by

\[u_h = \sum_{i=1}^{n} U_i \phi_i\]

where \(\{\phi_i\}_{i=1}^{n}\) is a basis for \(V_h\), and \(U\) is a vector of expansion coefficients for \(u_h\).

explicit Function(std::shared_ptr<const FunctionSpace> V)

Create function on given function space (shared data)

Arguments
V (FunctionSpace)
The function space.
Function(std::shared_ptr<const FunctionSpace> V, std::shared_ptr<GenericVector> x)

Create function on given function space with a given vector (shared data)

Warning: This constructor is intended for internal library use only

Arguments
V (FunctionSpace)
The function space.
x (GenericVector)
The vector.
Function(std::shared_ptr<const FunctionSpace> V, std::string filename)

Create function from vector of dofs stored to file (shared data)

Arguments
V (FunctionSpace)
The function space.
filename_dofdata (std::string)
The name of the file containing the dofmap data.
Function(const Function &v)

Copy constructor

Arguments
v (Function)
The object to be copied.
Function(const Function &v, std::size_t i)

Sub-function constructor with shallow copy of vector (used in Python interface)

Arguments
v (Function)
The function to be copied.
i (std::size_t)
Index of subfunction.
const Function &operator=(const Function &v)

Assignment from function

Arguments
v (Function)
Another function.
const Function &operator=(const Expression &v)

Assignment from expression using interpolation

Arguments
v (Expression)
The expression.
void operator=(const FunctionAXPY &axpy)

Assignment from linear combination of function

Arguments
v (FunctionAXPY)
A linear combination of other Functions
Function &operator[](std::size_t i) const

Extract subfunction

Arguments
i (std::size_t)
Index of subfunction.
Returns
Function
The subfunction.
std::shared_ptr<const FunctionSpace> function_space() const override

Return shared pointer to function space

Returns
FunctionSpace
Return the shared pointer.
std::shared_ptr<GenericVector> vector()

Return vector of expansion coefficients (non-const version)

Returns
GenericVector
The vector of expansion coefficients.
std::shared_ptr<const GenericVector> vector() const

Return vector of expansion coefficients (const version)

Returns
GenericVector
The vector of expansion coefficients (const).
bool in(const FunctionSpace &V) const

Check if function is a member of the given function space

Arguments
V (FunctionSpace)
The function space.
Returns
bool
True if the function is in the function space.
std::size_t geometric_dimension() const

Return geometric dimension

Returns
std::size_t
The geometric dimension.
void eval(Array<double> &values, const Array<double> &x) const override

Evaluate function at given coordinates

@param values (Array<double>)
The values.
@param x (Array<double>)
The coordinates.
void eval(Array<double> &values, const Array<double> &x, const Cell &dolfin_cell, const ufc::cell &ufc_cell) const

Evaluate function at given coordinates in given cell

Arguments @param values (Array<double>)

The values.
@param x (Array<double>)
The coordinates.
@param dolfin_cell (Cell)
The cell.
@param ufc_cell (ufc::cell)
The ufc::cell.
void eval(Eigen::Ref<Eigen::VectorXd> values, Eigen::Ref<const Eigen::VectorXd> x) const override

Evaluate function at given coordinates

@param values (Eigen::Ref<Eigen::VectorXd> values)
The values.
@param x (Eigen::Ref<const Eigen::VectorXd> x)
The coordinates.
void eval(Eigen::Ref<Eigen::VectorXd> values, Eigen::Ref<const Eigen::VectorXd> x, const dolfin::Cell &dolfin_cell, const ufc::cell &ufc_cell) const

Evaluate function at given coordinates in given cell

Arguments @param values (Eigen::Ref<Eigen::VectorXd>)

The values.
@param x (Eigen::Ref<const Eigen::VectorXd>)
The coordinates.
@param dolfin_cell (Cell)
The cell.
@param ufc_cell (ufc::cell)
The ufc::cell.
void interpolate(const GenericFunction &v)

Interpolate function (on possibly non-matching meshes)

@param v (GenericFunction)
The function to be interpolated.
void extrapolate(const Function &v)

Extrapolate function (from a possibly lower-degree function space)

Arguments
v (Function)
The function to be extrapolated.
std::size_t value_rank() const override

Return value rank

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

Return value dimension for given axis

Arguments
i (std::size_t)
The index of the axis.
Returns
std::size_t
The value dimension.
std::vector<std::size_t> value_shape() const override

Return value shape

Returns
std::vector<std::size_t>
The value shape.
void eval(Array<double> &values, const Array<double> &x, const ufc::cell &cell) const override

Evaluate at given point in given cell

@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 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.
void compute_vertex_values(std::vector<double> &vertex_values)

Compute values at all mesh vertices

@param vertex_values (Array<double>)
The values at all vertices.
void set_allow_extrapolation(bool allow_extrapolation)

Allow extrapolation when evaluating the Function

@param allow_extrapolation (bool)
Whether or not permit extrapolation.
bool get_allow_extrapolation() const

Check if extrapolation is permitted when evaluating the Function

@return bool
True if extrapolation is permitted, otherwise false