dolfin/function

Documentation for C++ code found in dolfin/function/*.h

Functions

assign

C++ documentation for assign from dolfin/function/assign.h:

void dolfin::assign(std::shared_ptr<Function> receiving_func, std::shared_ptr<const Function> assigning_func)

Assign one function to another. The functions must reside in the same type of FunctionSpace . One or both functions can be sub functions.

Parameters:
  • receiving_func – (std::shared_ptr<Function >) The receiving function
  • assigning_func – (std::shared_ptr<Function >) The assigning function

C++ documentation for assign from dolfin/function/assign.h:

void dolfin::assign(std::shared_ptr<Function> receiving_func, std::vector<std::shared_ptr<const Function>> assigning_funcs)

Assign several functions to sub functions of a mixed receiving function. The number of receiving functions must sum up to the number of sub functions in the assigning mixed function. The sub spaces of the assigning mixed space must be of the same type ans size as the receiving spaces.

Parameters:
  • receiving_func – (std::shared_ptr<Function >) The receiving function
  • assigning_funcs – (std::vector<std::shared_ptr<Function >>) The assigning functions

C++ documentation for assign from dolfin/function/assign.h:

void dolfin::assign(std::vector<std::shared_ptr<Function>> receiving_funcs, std::shared_ptr<const Function> assigning_func)

Assign sub functions of a single mixed function to single receiving functions. The number of sub functions in the assigning mixed function must sum up to the number of receiving functions. The sub spaces of the receiving mixed space must be of the same type ans size as the assigning spaces.

Parameters:
  • receiving_funcs – (std::vector<std::shared_ptr<Function >>) The receiving functions
  • assigning_func – (std::shared_ptr<Function >) The assigning function

Classes

CoefficientAssigner

C++ documentation for CoefficientAssigner from dolfin/function/CoefficientAssigner.h:

class dolfin::CoefficientAssigner

This class is used for assignment of coefficients to forms, which allows magic like a.f = f a.g = g which will insert the coefficients f and g in the correct positions in the list of coefficients for the form.

dolfin::CoefficientAssigner::CoefficientAssigner(Form &form, std::size_t number)

Create coefficient assigner for coefficient with given number.

Parameters:
  • form
  • number
void dolfin::CoefficientAssigner::operator=(std::shared_ptr<const GenericFunction> coefficient)

Assign coefficient.

Parameters:coefficient
dolfin::CoefficientAssigner::~CoefficientAssigner()

Destructor.

Constant

C++ documentation for Constant from dolfin/function/Constant.h:

class dolfin::Constant : public dolfin::Expression

This class represents a constant-valued expression.

dolfin::Constant::Constant(const Constant &constant)

Copy constructor

Parameters:constant – (Constant() ) Object to be copied.
dolfin::Constant::Constant(double value)

Create scalar constant

Constant c(1.0);
Parameters:value – (double) The scalar to create a Constant() object from.
dolfin::Constant::Constant(double value0, double value1)

Create vector constant (dim = 2)

Constant B(0.0, 1.0);
Parameters:
  • value0 – (double) The first vector element.
  • value1 – (double) The second vector element.
dolfin::Constant::Constant(double value0, double value1, double value2)

Create vector constant (dim = 3)

Constant T(0.0, 1.0, 0.0);
Parameters:
  • value0 – (double) The first vector element.
  • value1 – (double) The second vector element.
  • value2 – (double) The third vector element.
dolfin::Constant::Constant(std::vector<double> values)

Create vector-valued constant

Parameters:values – (std::vector<double>) Values to create a vector-valued constant from.
dolfin::Constant::Constant(std::vector<std::size_t> value_shape, std::vector<double> values)

Create tensor-valued constant for flattened array of values

Parameters:
  • value_shape – (std::vector<std::size_t>) Shape of tensor.
  • values – (std::vector<double>) Values to create tensor-valued constant from.
void dolfin::Constant::eval(Array<double> &values, const Array<double> &x) const override

Evaluate at given point (deprecated)

Parameters:
  • values – (Array<double>) The values at the point.
  • x – (Array<double>) The coordinates of the point.
void dolfin::Constant::eval(Eigen::Ref<Eigen::VectorXd> values, Eigen::Ref<const Eigen::VectorXd> x) const override

Evaluate at given point.

Parameters:
  • values – (Eigen::Ref<Eigen::VectorXd>) The values at the point.
  • x – (Eigen::Ref<const Eigen::VectorXd>) The coordinates of the point.
dolfin::Constant::operator double() const

Cast to double (for scalar constants)

Returns:double The scalar value.
const Constant &dolfin::Constant::operator=(const Constant &constant)

Assignment operator

Parameters:constant – (Constant ) Another constant.
const Constant &dolfin::Constant::operator=(double constant)

Assignment operator

Parameters:constant – (double) Another constant.
std::string dolfin::Constant::str(bool verbose) const override

Return informal string representation (pretty-print)

Parameters:verbose
std::vector<double> dolfin::Constant::values() const

Return copy of this Constant ‘s current values

Returns:std::vector<double> The vector of scalar values of the constant.
dolfin::Constant::~Constant()

Destructor.

Expression

C++ documentation for Expression from dolfin/function/Expression.h:

class dolfin::Expression : public dolfin::GenericFunction

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.

dolfin::Expression::Expression()

Create scalar expression.

dolfin::Expression::Expression(const Expression &expression)

Copy constructor

Parameters:expression – (Expression() ) Object to be copied.
dolfin::Expression::Expression(std::size_t dim)

Create vector-valued expression with given dimension.

Parameters:dim – (std::size_t) Dimension of the vector-valued expression.
dolfin::Expression::Expression(std::size_t dim0, std::size_t dim1)

Create matrix-valued expression with given dimensions.

Parameters:
  • dim0 – (std::size_t) Dimension (rows).
  • dim1 – (std::size_t) Dimension (columns).
dolfin::Expression::Expression(std::vector<std::size_t> value_shape)

Create tensor-valued expression with given shape.

Parameters:value_shape – (std::vector<std::size_t>) Shape of expression.
void dolfin::Expression::compute_vertex_values(std::vector<double> &vertex_values, const Mesh &mesh) const override

Compute values at all mesh vertices.

Parameters:
  • vertex_values – (Array<double>) The values at all vertices.
  • mesh – (Mesh ) The mesh.
void dolfin::Expression::eval(Array<double> &values, const Array<double> &x) const override

Evaluate at given point (deprecated)

Parameters:
  • values – (Array<double>) The values at the point.
  • x – (Array<double>) The coordinates of the point.
void dolfin::Expression::eval(Array<double> &values, const Array<double> &x, const ufc::cell &cell) const override

Evaluate at given point in given cell (deprecated)

Parameters:
  • values – (Array<double>) The values at the point.
  • x – (Array<double>) The coordinates of the point.
  • cell – (ufc::cell ) The cell which contains the given point.
void dolfin::Expression::eval(Eigen::Ref<Eigen::VectorXd> values, Eigen::Ref<const Eigen::VectorXd> x) const override

Evaluate at given point.

Parameters:
  • values – (Eigen::Ref<Eigen::VectorXd>) The values at the point.
  • x – (Eigen::Ref<const Eigen::VectorXd>) The coordinates of the point.
void dolfin::Expression::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

Parameters:
  • values – (Eigen::Ref<Eigen::VectorXd>) The values at the point.
  • x – (Eigen::Ref<const Eigen::VectorXd>) The coordinates of the point.
  • cell – (ufc::cell ) The cell which contains the given point.
std::shared_ptr<const FunctionSpace> dolfin::Expression::function_space() const override

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

Returns:FunctionSpace Return the shared pointer.
std::shared_ptr<dolfin::GenericFunction> dolfin::Expression::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

Parameters:name
double dolfin::Expression::get_property(std::string name) const

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

Parameters:name
void dolfin::Expression::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).

Parameters:
  • w – (list of doubles) Expansion coefficients.
  • element – (FiniteElement ) The element.
  • dolfin_cell – (Cell ) The cell.
  • coordinate_dofs – (double*) The coordinates
  • ufc_cell – (ufc::cell ) The ufc::cell .
void dolfin::Expression::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

Parameters:
  • name
  • f
void dolfin::Expression::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

Parameters:
  • name
  • value
std::size_t dolfin::Expression::value_dimension(std::size_t i) const override

Return value dimension for given axis.

Parameters:i – (std::size_t) Integer denoting the axis to use.
Returns:std::size_t The value dimension (for the given axis).
std::size_t dolfin::Expression::value_rank() const override

Return value rank.

Returns:std::size_t The value rank.
std::vector<std::size_t> dolfin::Expression::value_shape() const override

Return value shape

Returns:std::vector<std::size_t> The value shape.
dolfin::Expression::~Expression()

Destructor.

FacetArea

C++ documentation for FacetArea from dolfin/function/SpecialFunctions.h:

class dolfin::FacetArea : public dolfin::Expression

This function represents the area/length of a cell facet on a given mesh.

dolfin::FacetArea::FacetArea(std::shared_ptr<const Mesh> mesh)

Constructor.

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

Evaluate function.

Parameters:
  • values
  • x
  • cell
Event dolfin::FacetArea::not_on_boundary

Function

C++ documentation for Function from dolfin/function/Function.h:

class dolfin::Function : public dolfin::GenericFunction

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

::
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\) .

Friends: FunctionAssigner, FunctionSpace.

dolfin::Function::Function()
dolfin::Function::Function(const Function &v)

Copy constructor Arguments v (Function() ) The object to be copied.

Parameters:v
dolfin::Function::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.

Parameters:
  • v
  • i
dolfin::Function::Function(std::shared_ptr<const FunctionSpace> V)

Create function on given function space (shared data) Arguments V (FunctionSpace ) The function space.

Parameters:V
dolfin::Function::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.

Parameters:
  • V
  • x
dolfin::Function::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.

Parameters:
  • V
  • filename
void dolfin::Function::compute_vertex_values(std::vector<double> &vertex_values)

Compute values at all mesh vertices

Parameters:vertex_values – (Array<double>) The values at all vertices.
void dolfin::Function::compute_vertex_values(std::vector<double> &vertex_values, const Mesh &mesh) const override

Compute values at all mesh vertices

Parameters:
  • vertex_values – (Array<double>) The values at all vertices.
  • mesh – (Mesh ) The mesh.
void dolfin::Function::eval(Array<double> &values, const Array<double> &x) const override

Evaluate function at given coordinates

Parameters:
  • values – (Array<double>) The values.
  • x – (Array<double>) The coordinates.
void dolfin::Function::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

Parameters:
  • values – (Array<double>) The values.
  • x – (Array<double>) The coordinates.
  • dolfin_cell – (Cell ) The cell.
  • ufc_cell – (ufc::cell ) The ufc::cell .
void dolfin::Function::eval(Array<double> &values, const Array<double> &x, const ufc::cell &cell) const override

Evaluate at given point in given cell

Parameters:
  • values – (Array<double>) The values at the point.
  • x – (Array<double>) The coordinates of the point.
  • cell – (ufc::cell ) The cell which contains the given point.
void dolfin::Function::eval(Eigen::Ref<Eigen::VectorXd> values, Eigen::Ref<const Eigen::VectorXd> x) const override

Evaluate function at given coordinates

Parameters:
  • values – (Eigen::Ref<Eigen::VectorXd> values) The values.
  • x – (Eigen::Ref<const Eigen::VectorXd> x) The coordinates.
void dolfin::Function::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

Parameters:
  • values – (Eigen::Ref<Eigen::VectorXd>) The values.
  • x – (Eigen::Ref<const Eigen::VectorXd>) The coordinates.
  • dolfin_cell – (Cell ) The cell.
  • ufc_cell – (ufc::cell ) The ufc::cell .
void dolfin::Function::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

Parameters:
  • values – (Eigen::Ref<Eigen::VectorXd>) The values at the point.
  • x – (Eigen::Ref<const Eigen::VectorXd> The coordinates of the point.
  • cell – (ufc::cell ) The cell which contains the given point.
void dolfin::Function::extrapolate(const Function &v)

Extrapolate function (from a possibly lower-degree function space) Arguments v (Function ) The function to be extrapolated.

Parameters:v
std::shared_ptr<const FunctionSpace> dolfin::Function::function_space() const override

Return shared pointer to function space Returns:cpp:any:FunctionSpace Return the shared pointer.

std::size_t dolfin::Function::geometric_dimension() const

Return geometric dimension Returns std::size_t The geometric dimension.

bool dolfin::Function::get_allow_extrapolation() const

Check if extrapolation is permitted when evaluating the Function

Returns:bool True if extrapolation is permitted, otherwise false
bool dolfin::Function::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.

Parameters:V
void dolfin::Function::init_vector()
void dolfin::Function::interpolate(const GenericFunction &v)

Interpolate function (on possibly non-matching meshes)

Parameters:v – (GenericFunction ) The function to be interpolated.
const Function &dolfin::Function::operator=(const Expression &v)

Assignment from expression using interpolation Arguments v (Expression ) The expression.

Parameters:v
const Function &dolfin::Function::operator=(const Function &v)

Assignment from function Arguments v (Function ) Another function.

Parameters:v
void dolfin::Function::operator=(const FunctionAXPY &axpy)

Assignment from linear combination of function Arguments v (FunctionAXPY ) A linear combination of other Functions

Parameters:axpy
Function &dolfin::Function::operator[](std::size_t i) const

Extract subfunction Arguments i (std::size_t) Index of subfunction. Returns:cpp:any:Function The subfunction.

Parameters:i
void dolfin::Function::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)

Parameters:
  • w – (list of doubles) Expansion coefficients.
  • element – (FiniteElement ) The element.
  • dolfin_cell – (Cell ) The cell.
  • coordinate_dofs – (double *) The coordinates
  • ufc_cell – (ufc::cell ). The ufc::cell .
void dolfin::Function::set_allow_extrapolation(bool allow_extrapolation)

Allow extrapolation when evaluating the Function

Parameters:allow_extrapolation – (bool) Whether or not permit extrapolation.
std::size_t dolfin::Function::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.

Parameters:i
std::size_t dolfin::Function::value_rank() const override

Return value rank Returns std::size_t The value rank.

std::vector<std::size_t> dolfin::Function::value_shape() const override

Return value shape Returns std::vector<std::size_t> The value shape.

std::shared_ptr<GenericVector> dolfin::Function::vector()

Return vector of expansion coefficients (non-const version) Returns:cpp:any:GenericVector The vector of expansion coefficients.

std::shared_ptr<const GenericVector> dolfin::Function::vector() const

Return vector of expansion coefficients (const version) Returns:cpp:any:GenericVector The vector of expansion coefficients (const).

dolfin::Function::~Function()

Destructor.

FunctionAXPY

C++ documentation for FunctionAXPY from dolfin/function/FunctionAXPY.h:

class dolfin::FunctionAXPY

This class represents a linear combination of functions. It is mostly used as an intermediate class for operations such as u = 3*u0 + 4*u1; where the rhs generates an FunctionAXPY .

enum dolfin::FunctionAXPY::Direction

Enum to decide what way AXPY is constructed.

enumerator dolfin::FunctionAXPY::Direction::ADD_ADD = 0
enumerator dolfin::FunctionAXPY::Direction::SUB_ADD = 1
enumerator dolfin::FunctionAXPY::Direction::ADD_SUB = 2
enumerator dolfin::FunctionAXPY::Direction::SUB_SUB = 3
dolfin::FunctionAXPY::FunctionAXPY(const FunctionAXPY &axpy)

Copy constructor.

Parameters:axpy
dolfin::FunctionAXPY::FunctionAXPY(const FunctionAXPY &axpy, double scalar)

Constructor.

Parameters:
  • axpy
  • scalar
dolfin::FunctionAXPY::FunctionAXPY(const FunctionAXPY &axpy, std::shared_ptr<const Function> func, Direction direction)

Constructor.

Parameters:
  • axpy
  • func
  • direction
dolfin::FunctionAXPY::FunctionAXPY(const FunctionAXPY &axpy0, const FunctionAXPY &axpy1, Direction direction)

Constructor.

Parameters:
  • axpy0
  • axpy1
  • direction
dolfin::FunctionAXPY::FunctionAXPY(std::shared_ptr<const Function> func, double scalar)

Constructor.

Parameters:
  • func
  • scalar
dolfin::FunctionAXPY::FunctionAXPY(std::shared_ptr<const Function> func0, std::shared_ptr<const Function> func1, Direction direction)

Constructor.

Parameters:
  • func0
  • func1
  • direction
dolfin::FunctionAXPY::FunctionAXPY(std::vector<std::pair<double, std::shared_ptr<const Function>>> pairs)

Constructor.

Parameters:pairs
FunctionAXPY dolfin::FunctionAXPY::operator*(double scale) const

Scale operator.

Parameters:scale
FunctionAXPY dolfin::FunctionAXPY::operator+(const FunctionAXPY &axpy) const

Addition operator.

Parameters:axpy
FunctionAXPY dolfin::FunctionAXPY::operator+(std::shared_ptr<const Function> func) const

Addition operator.

Parameters:func
FunctionAXPY dolfin::FunctionAXPY::operator-(const FunctionAXPY &axpy) const

Subtraction operator.

Parameters:axpy
FunctionAXPY dolfin::FunctionAXPY::operator-(std::shared_ptr<const Function> func) const

Subtraction operator.

Parameters:func
FunctionAXPY dolfin::FunctionAXPY::operator/(double scale) const

Scale operator.

Parameters:scale
const std::vector<std::pair<double, std::shared_ptr<const Function>>> &dolfin::FunctionAXPY::pairs() const

Return the scalar and Function pairs.

dolfin::FunctionAXPY::~FunctionAXPY()

Destructor.

FunctionAssigner

C++ documentation for FunctionAssigner from dolfin/function/FunctionAssigner.h:

class dolfin::FunctionAssigner

This class facilitate assignments between Function and sub Functions. It builds and caches maps between compatible dofs. These maps are used in the assignment methods which perform the actual assignment. Optionally can a MeshFunction be passed together with a label, facilitating FunctionAssignment over sub domains.

dolfin::FunctionAssigner::FunctionAssigner(std::shared_ptr<const FunctionSpace> receiving_space, std::shared_ptr<const FunctionSpace> assigning_space)

Create a FunctionAssigner() between functions residing in the same type of FunctionSpace . One or both functions can be sub functions. Arguments receiving_space (FunctionSpace ) The function space of the receiving function assigning_space (FunctionSpace ) The function space of the assigning function

Parameters:
  • receiving_space
  • assigning_space
dolfin::FunctionAssigner::FunctionAssigner(std::shared_ptr<const FunctionSpace> receiving_space, std::vector<std::shared_ptr<const FunctionSpace>> assigning_spaces)

Create a FunctionAssigner() between several functions (assigning) and one mixed function (receiving). The number of sub functions in the assigning mixed function must sum up to the number of receiving functions. The sub spaces of the receiving mixed space must be of the same type ans size as the assigning spaces. Arguments receiving_space (std::shared_ptr<FunctionSpace >) The receiving function space assigning_spaces (std::vector<std::shared_ptr<FunctionSpace > >) The assigning function spaces

Parameters:
  • receiving_space
  • assigning_spaces
dolfin::FunctionAssigner::FunctionAssigner(std::vector<std::shared_ptr<const FunctionSpace>> receiving_spaces, std::shared_ptr<const FunctionSpace> assigning_space)

Create a FunctionAssigner() between one mixed function (assigning) and several functions (receiving). The number of receiving functions must sum up to the number of sub functions in the assigning mixed function. The sub spaces of the assigning mixed space must be of the same type ans size as the receiving spaces. Arguments receiving_spaces (std::vector<FunctionSpace >) The receiving function spaces assigning_space (FunctionSpace ) The assigning function space

Parameters:
  • receiving_spaces
  • assigning_space
void dolfin::FunctionAssigner::assign(std::shared_ptr<Function> receiving_func, std::shared_ptr<const Function> assigning_func) const

Assign one function to another Arguments receiving_func (std::shared_ptr<Function >) The receiving function assigning_func (std::shared_ptr<Function >) The assigning function

Parameters:
  • receiving_func
  • assigning_func
void dolfin::FunctionAssigner::assign(std::shared_ptr<Function> receiving_func, std::vector<std::shared_ptr<const Function>> assigning_funcs) const

Assign several functions to sub functions of a mixed receiving function Arguments receiving_func (std::shared_ptr<Function >) The receiving mixed function assigning_funcs (std::vector<std::shared_ptr<Function > >) The assigning functions

Parameters:
  • receiving_func
  • assigning_funcs
void dolfin::FunctionAssigner::assign(std::vector<std::shared_ptr<Function>> receiving_funcs, std::shared_ptr<const Function> assigning_func) const

Assign sub functions of a single mixed function to single receiving functions Arguments receiving_funcs (std::vector<std::shared_ptr<Function > >) The receiving functions assigning_func (std::shared_ptr<Function >) The assigning mixed function

Parameters:
  • receiving_funcs
  • assigning_func
std::size_t dolfin::FunctionAssigner::num_assigning_functions() const

Return the number of assigning functions.

std::size_t dolfin::FunctionAssigner::num_receiving_functions() const

Return the number of receiving functions.

dolfin::FunctionAssigner::~FunctionAssigner()

Destructor.

FunctionSpace

C++ documentation for FunctionSpace from dolfin/function/FunctionSpace.h:

class dolfin::FunctionSpace

This class represents a finite element function space defined by a mesh, a finite element, and a local-to-global mapping of the degrees of freedom (dofmap).

dolfin::FunctionSpace::FunctionSpace(const FunctionSpace &V)

Copy constructor Arguments V (FunctionSpace() ) The object to be copied.

Parameters:V
dolfin::FunctionSpace::FunctionSpace(std::shared_ptr<const Mesh> mesh)

Create empty function space for later initialization. This constructor is intended for use by any sub-classes which need to construct objects before the initialisation of the base class. Data can be attached to the base class using FunctionSpace::attach (...). Arguments mesh (Mesh ) The mesh.

Parameters:mesh
dolfin::FunctionSpace::FunctionSpace(std::shared_ptr<const Mesh> mesh, std::shared_ptr<const FiniteElement> element, std::shared_ptr<const GenericDofMap> dofmap)

Create function space for given mesh, element and dofmap (shared data) Arguments mesh (Mesh ) The mesh. element (FiniteElement ) The element. dofmap (GenericDofMap ) The dofmap.

Parameters:
  • mesh
  • element
  • dofmap
void dolfin::FunctionSpace::attach(std::shared_ptr<const FiniteElement> element, std::shared_ptr<const GenericDofMap> dofmap)

Attach data to an empty function space Arguments element (FiniteElement ) The element. dofmap (GenericDofMap ) The dofmap.

Parameters:
  • element
  • dofmap
std::shared_ptr<FunctionSpace> dolfin::FunctionSpace::collapse() const

Collapse a subspace and return a new function space Returns:cpp:any:FunctionSpace The new function space.

std::shared_ptr<FunctionSpace> dolfin::FunctionSpace::collapse(std::unordered_map<std::size_t, std::size_t> &collapsed_dofs) const

Collapse a subspace and return a new function space and a map from new to old dofs Arguments collapsed_dofs (std::unordered_map<std::size_t, std::size_t>) The map from new to old dofs. Returns:cpp:any:FunctionSpace The new function space.

Parameters:collapsed_dofs
std::vector<std::size_t> dolfin::FunctionSpace::component() const

Return component w.r.t. to root superspace, i.e. W.sub(1).sub(0) == [1, 0]. Returns std::vector<std::size_t> The component (w.r.t to root superspace).

bool dolfin::FunctionSpace::contains(const FunctionSpace &V) const

Check whether V is subspace of this, or this itself Arguments V (FunctionSpace ) The space to be tested for inclusion. Returns bool True if V is contained or equal to this.

Parameters:V
std::size_t dolfin::FunctionSpace::dim() const

Return global dimension of the function space. Equivalent to dofmap() ->global_dimension() Returns std::size_t The dimension of the function space.

std::shared_ptr<const GenericDofMap> dolfin::FunctionSpace::dofmap() const

Return dofmap Returns:cpp:any:GenericDofMap The dofmap.

std::shared_ptr<const FiniteElement> dolfin::FunctionSpace::element() const

Return finite element Returns:cpp:any:FiniteElement The finite element.

std::shared_ptr<FunctionSpace> dolfin::FunctionSpace::extract_sub_space(const std::vector<std::size_t> &component) const

Extract subspace for component Arguments component (std::vector<std::size_t>) The component. Returns:cpp:any:FunctionSpace The subspace.

Parameters:component
bool dolfin::FunctionSpace::has_cell(const Cell &cell) const

Check if function space has given cell Arguments cell (Cell ) The cell. Returns bool True if the function space has the given cell.

Parameters:cell
bool dolfin::FunctionSpace::has_element(const FiniteElement &element) const

Check if function space has given element Arguments element (FiniteElement ) The finite element. Returns bool True if the function space has the given element.

Parameters:element
void dolfin::FunctionSpace::interpolate(GenericVector &expansion_coefficients, const GenericFunction &v) const

Interpolate function v into function space, returning the vector of expansion coefficients Arguments expansion_coefficients (GenericVector ) The expansion coefficients. v (GenericFunction ) The function to be interpolated.

Parameters:
  • expansion_coefficients
  • v
void dolfin::FunctionSpace::interpolate_from_any(GenericVector &expansion_coefficients, const GenericFunction &v) const
Parameters:
  • expansion_coefficients
  • v
void dolfin::FunctionSpace::interpolate_from_parent(GenericVector &expansion_coefficients, const GenericFunction &v) const
Parameters:
  • expansion_coefficients
  • v
std::shared_ptr<const Mesh> dolfin::FunctionSpace::mesh() const

Return mesh Returns:cpp:any:Mesh The mesh.

bool dolfin::FunctionSpace::operator!=(const FunctionSpace &V) const

Inequality operator Arguments V (FunctionSpace ) Another function space.

Parameters:V
const FunctionSpace &dolfin::FunctionSpace::operator=(const FunctionSpace &V)

Assignment operator Arguments V (FunctionSpace ) Another function space.

Parameters:V
bool dolfin::FunctionSpace::operator==(const FunctionSpace &V) const

Equality operator Arguments V (FunctionSpace ) Another function space.

Parameters:V
std::shared_ptr<FunctionSpace> dolfin::FunctionSpace::operator[](std::size_t i) const

Extract subspace for component Arguments i (std::size_t) Index of the subspace. Returns:cpp:any:FunctionSpace The subspace.

Parameters:i
void dolfin::FunctionSpace::print_dofmap() const

Print dofmap (useful for debugging)

void dolfin::FunctionSpace::set_x(GenericVector &x, double value, std::size_t component) const

Set dof entries in vector to value*x[i], where [x][i] is the coordinate of the dof spatial coordinate. Parallel layout of vector must be consistent with dof map range This function is typically used to construct the null space of a matrix operator, e.g. rigid body rotations. Arguments vector (GenericVector ) The vector to set. value (double) The value to multiply to coordinate by. component (std::size_t) The coordinate index. mesh (Mesh ) The mesh.

Parameters:
  • x
  • value
  • component
std::string dolfin::FunctionSpace::str(bool verbose) const

Return informal string representation (pretty-print) Arguments verbose (bool) Flag to turn on additional output. Returns std::string An informal representation of the function space.

Parameters:verbose
std::shared_ptr<FunctionSpace> dolfin::FunctionSpace::sub(const std::vector<std::size_t> &component) const

Extract subspace for component Arguments component (std::vector<std::size_t>) The component. Returns:cpp:any:FunctionSpace The subspace.

Parameters:component
std::shared_ptr<FunctionSpace> dolfin::FunctionSpace::sub(std::size_t component) const

Extract subspace for component Arguments component (std::size_t) Index of the subspace. Returns:cpp:any:FunctionSpace The subspace.

Parameters:component
std::vector<double> dolfin::FunctionSpace::tabulate_dof_coordinates() const

Tabulate the coordinates of all dofs on this process. This function is typically used by preconditioners that require the spatial coordinates of dofs, for example for re-partitioning or nullspace computations. Arguments mesh (Mesh ) The mesh. Returns std::vector<double> The dof coordinates (x0, y0, x1, y1, . . .)

dolfin::FunctionSpace::~FunctionSpace()

Destructor.

GenericFunction

C++ documentation for GenericFunction from dolfin/function/GenericFunction.h:

class dolfin::GenericFunction

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.

dolfin::GenericFunction::GenericFunction()

Constructor.

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

Compute values at all mesh vertices.

Parameters:
  • vertex_values
  • mesh
void dolfin::GenericFunction::eval(Array<double> &values, const Array<double> &x) const

Evaluate at given point (deprecated)

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

Evaluate at given point in given cell (deprecated)

Parameters:
  • values
  • x
  • cell
void dolfin::GenericFunction::eval(Eigen::Ref<Eigen::VectorXd> values, Eigen::Ref<const Eigen::VectorXd> x) const

Evaluate at given point.

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

Evaluate at given point in given cell.

Parameters:
  • values
  • x
  • cell
void dolfin::GenericFunction::evaluate(double *values, const double *coordinates, const ufc::cell &cell) const override

Evaluate function at given point in cell

Parameters:
  • values – (double*)
  • coordinates – (const double*)
  • cell – (ufc::cell &)
std::shared_ptr<const FunctionSpace> dolfin::GenericFunction::function_space() const = 0

Pointer to FunctionSpace , if appropriate, otherwise NULL.

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

Evaluation at given point (vector-valued function)

Parameters:
  • values
  • p
void dolfin::GenericFunction::operator()(Array<double> &values, double x) const

Evaluation at given point (vector-valued function)

Parameters:
  • values
  • x
void dolfin::GenericFunction::operator()(Array<double> &values, double x, double y) const

Evaluation at given point (vector-valued function)

Parameters:
  • values
  • x
  • y
void dolfin::GenericFunction::operator()(Array<double> &values, double x, double y, double z) const

Evaluation at given point (vector-valued function)

Parameters:
  • values
  • x
  • y
  • z
double dolfin::GenericFunction::operator()(const Point &p) const

Evaluation at given point (scalar function)

Parameters:p
double dolfin::GenericFunction::operator()(double x) const

Evaluation at given point (scalar function)

Parameters:x
double dolfin::GenericFunction::operator()(double x, double y) const

Evaluation at given point (scalar function)

Parameters:
  • x
  • y
double dolfin::GenericFunction::operator()(double x, double y, double z) const

Evaluation at given point (scalar function)

Parameters:
  • x
  • y
  • z
void dolfin::GenericFunction::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)

Parameters:
  • w
  • element
  • dolfin_cell
  • coordinate_dofs
  • ufc_cell
void dolfin::GenericFunction::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)

Parameters:
  • w
  • element
  • dolfin_cell
  • coordinate_dofs
  • ufc_cell
void dolfin::GenericFunction::update() const

Update off-process ghost coefficients.

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

Return value dimension for given axis.

Parameters:i
std::size_t dolfin::GenericFunction::value_rank() const = 0

Return value rank.

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

Return value shape.

std::size_t dolfin::GenericFunction::value_size() const

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

dolfin::GenericFunction::~GenericFunction()

Destructor.

LagrangeInterpolator

C++ documentation for LagrangeInterpolator from dolfin/function/LagrangeInterpolator.h:

class dolfin::LagrangeInterpolator

This class interpolates efficiently from a GenericFunction to a Lagrange Function

void dolfin::LagrangeInterpolator::extract_dof_component_map(std::unordered_map<std::size_t, std::size_t> &dof_component_map, const FunctionSpace &V, int *component)
Parameters:
  • dof_component_map
  • V
  • component
bool dolfin::LagrangeInterpolator::in_bounding_box(const std::vector<double> &point, const std::vector<double> &bounding_box, const double tol)
Parameters:
  • point
  • bounding_box
  • tol
void dolfin::LagrangeInterpolator::interpolate(Function &u, const Expression &u0)

Interpolate Expression Arguments u (Function ) The resulting Function u0 (Expression ) The Expression to be interpolated.

Parameters:
  • u
  • u0
void dolfin::LagrangeInterpolator::interpolate(Function &u, const Function &u0)

Interpolate function (on possibly non-matching meshes) Arguments u (Function ) The resulting Function u0 (Function ) The Function to be interpolated.

Parameters:
  • u
  • u0
std::map<std::vector<double>, std::vector<std::size_t>, lt_coordinate> dolfin::LagrangeInterpolator::tabulate_coordinates_to_dofs(const FunctionSpace &V)
Parameters:V

MeshCoordinates

C++ documentation for MeshCoordinates from dolfin/function/SpecialFunctions.h:

class dolfin::MeshCoordinates : public dolfin::Expression

This Function represents the mesh coordinates on a given mesh.

dolfin::MeshCoordinates::MeshCoordinates(std::shared_ptr<const Mesh> mesh)

Constructor.

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

Evaluate function.

Parameters:
  • values
  • x
  • cell

MultiMeshCoefficientAssigner

C++ documentation for MultiMeshCoefficientAssigner from dolfin/function/MultiMeshCoefficientAssigner.h:

class dolfin::MultiMeshCoefficientAssigner

This class is used for assignment of multimesh coefficients to forms, which allows magic like a.f = f a.g = g which will insert the coefficients f and g in the correct positions in the list of coefficients for the form. Note that coefficients can also be assigned manually to the individual parts of a multimesh form. Assigning to a multimesh coefficient assigner will assign the same coefficient to all parts of a form.

dolfin::MultiMeshCoefficientAssigner::MultiMeshCoefficientAssigner(MultiMeshForm &form, std::size_t number)

Create multimesh coefficient assigner for coefficient with given number.

Parameters:
  • form
  • number
void dolfin::MultiMeshCoefficientAssigner::operator=(std::shared_ptr<const GenericFunction> coefficient)

Assign coefficient from GenericFunction .

Parameters:coefficient
void dolfin::MultiMeshCoefficientAssigner::operator=(std::shared_ptr<const MultiMeshFunction> coefficient)

Assign coefficient from MultiMeshFunction .

Parameters:coefficient
dolfin::MultiMeshCoefficientAssigner::~MultiMeshCoefficientAssigner()

Destructor.

MultiMeshFunction

C++ documentation for MultiMeshFunction from dolfin/function/MultiMeshFunction.h:

class dolfin::MultiMeshFunction

This class represents a function on a cut and composite finite element function space (MultiMesh ) defined on one or more possibly intersecting meshes.

dolfin::MultiMeshFunction::MultiMeshFunction()

Constructor.

dolfin::MultiMeshFunction::MultiMeshFunction(std::shared_ptr<const MultiMeshFunctionSpace> V)

Create MultiMesh function on given MultiMesh function space Arguments V (MultiMeshFunctionSpace ) The MultiMesh function space.

Parameters:V
dolfin::MultiMeshFunction::MultiMeshFunction(std::shared_ptr<const MultiMeshFunctionSpace> V, std::shared_ptr<GenericVector> x)

Create MultiMesh function on given MultiMesh function space with a given vector (shared data) Warning: This constructor is intended for internal library use only Arguments V (MultiMeshFunctionSpace ) The multimesh function space. x (GenericVector ) The vector.

Parameters:
  • V
  • x
void dolfin::MultiMeshFunction::assign_part(std::size_t a, const Function &v)

Assign Function to part of a mesh Arguments a (int) Part mesh assigned to V (Function ) The vector

Parameters:
  • a
  • v
void dolfin::MultiMeshFunction::compute_ghost_indices(std::pair<std::size_t, std::size_t> range, std::vector<la_index> &ghost_indices) const
Parameters:
  • range
  • ghost_indices
void dolfin::MultiMeshFunction::eval(Array<double> &values, const Array<double> &x) const

Evaluate at a given point.

Parameters:
  • values
  • x
void dolfin::MultiMeshFunction::eval(Array<double> &values, const Array<double> &x, std::size_t part, const ufc::cell &cell) const

Evaluate at given point in given cell in given part Arguments values (Array <double>) The values at the point. x (Array <double>) The coordinates of the point. cell (ufc::cell ) The cell which contains the given point.

Parameters:
  • values
  • x
  • part
  • cell
std::shared_ptr<const MultiMeshFunctionSpace> dolfin::MultiMeshFunction::function_space() const

Return shared pointer to multi mesh function space Returns:cpp:any:MultiMeshFunctionSpace Return the shared pointer.

void dolfin::MultiMeshFunction::init_vector()
std::shared_ptr<const Function> dolfin::MultiMeshFunction::part(std::size_t i) const

Return function (part) number i Returns:cpp:any:Function Function (part) number i

Parameters:i
std::shared_ptr<const Function> dolfin::MultiMeshFunction::part(std::size_t i, bool deepcopy) const

Return a (deepcopy) function (part) number i Returns:cpp:any:Function Function (part) number i bool deepcopy flag

Parameters:
  • i
  • deepcopy
void dolfin::MultiMeshFunction::restrict(double *w, const FiniteElement &element, std::size_t part, const Cell &dolfin_cell, const double *coordinate_dofs, const ufc::cell &ufc_cell) const

Restrict function to local cell in given part (compute expansion coefficients w) Arguments w (list of doubles) Expansion coefficients. element (FiniteElement ) The element. part (std::size_t) The mesh part dolfin_cell (Cell ) The cell. ufc_cell (ufc::cell ). The ufc::cell .

Parameters:
  • w
  • element
  • part
  • dolfin_cell
  • coordinate_dofs
  • ufc_cell
void dolfin::MultiMeshFunction::restrict_as_ufc_function(double *w, const FiniteElement &element, std::size_t part, const Cell &dolfin_cell, const double *coordinate_dofs, const ufc::cell &ufc_cell) const

Restrict as UFC function (by calling eval)

Parameters:
  • w
  • element
  • part
  • dolfin_cell
  • coordinate_dofs
  • ufc_cell
std::shared_ptr<GenericVector> dolfin::MultiMeshFunction::vector()

Return vector of expansion coefficients (non-const version) Returns:cpp:any:GenericVector The vector of expansion coefficients.

std::shared_ptr<const GenericVector> dolfin::MultiMeshFunction::vector() const

Return vector of expansion coefficients (const version) Returns:cpp:any:GenericVector The vector of expansion coefficients (const).

dolfin::MultiMeshFunction::~MultiMeshFunction()

Destructor.

MultiMeshFunctionSpace

C++ documentation for MultiMeshFunctionSpace from dolfin/function/MultiMeshFunctionSpace.h:

class dolfin::MultiMeshFunctionSpace

This class represents a function space on a multimesh. It may may be created from a set of standard function spaces by repeatedly calling add , followed by a call to build . Note that a multimesh function space is not useful and its data structures are empty until build has been called.

Friends: MultiMeshSubSpace.

dolfin::MultiMeshFunctionSpace::MultiMeshFunctionSpace(std::shared_ptr<const MultiMesh> multimesh)

Create multimesh function space on multimesh (shared pointer version)

Parameters:multimesh
void dolfin::MultiMeshFunctionSpace::add(std::shared_ptr<const FunctionSpace> function_space)

Add function space Arguments function_space (FunctionSpace ) The function space.

Parameters:function_space
void dolfin::MultiMeshFunctionSpace::build()

Build multimesh function space.

void dolfin::MultiMeshFunctionSpace::build(const std::vector<dolfin::la_index> &offsets)

Build multimesh function space. This function uses offsets computed from the full function spaces on each part.

Parameters:offsets
std::size_t dolfin::MultiMeshFunctionSpace::dim() const

Return dimension of the multimesh function space Returns std::size_t The dimension of the multimesh function space.

std::shared_ptr<const MultiMeshDofMap> dolfin::MultiMeshFunctionSpace::dofmap() const

Return multimesh dofmap Returns:cpp:any:MultiMeshDofMap The dofmap.

void dolfin::MultiMeshFunctionSpace::lock_inactive_dofs(GenericMatrix &A, GenericVector &b) const

Lock inactive dofs of a system.

Parameters:
  • A
  • b
std::shared_ptr<const MultiMesh> dolfin::MultiMeshFunctionSpace::multimesh() const

Return multimesh Returns:cpp:any:MultiMesh The multimesh.

std::size_t dolfin::MultiMeshFunctionSpace::num_parts() const

Return the number of function spaces (parts) of the multimesh function space Returns std::size_t The number of function spaces (parts) of the multimesh function space.

std::shared_ptr<const FunctionSpace> dolfin::MultiMeshFunctionSpace::part(std::size_t i) const

Return function space (part) number i Arguments i (std::size_t) The part number Returns:cpp:any:FunctionSpace Function space (part) number i

Parameters:i
std::shared_ptr<const FunctionSpace> dolfin::MultiMeshFunctionSpace::view(std::size_t i) const

Return view of multimesh function space for part number i. This function differs from the part() function in that it does not return the original function space for a part, but rather a view of the common multimesh function space (dofs global to the collection of parts). Arguments i (std::size_t) The part number Returns:cpp:any:FunctionSpace Function space (part) number i

Parameters:i
dolfin::MultiMeshFunctionSpace::~MultiMeshFunctionSpace()

Destructor.

MultiMeshSubSpace

C++ documentation for MultiMeshSubSpace from dolfin/function/MultiMeshSubSpace.h:

class dolfin::MultiMeshSubSpace : public dolfin::MultiMeshFunctionSpace

This class represents a subspace (component) of a multimesh function space. The subspace is specified by an array of indices. For example, the array [3, 0, 2] specifies subspace 2 of subspace 0 of subspace 3. A typical example is the function space W = V x P for Stokes. Here, V = W[0] is the subspace for the velocity component and P = W[1] is the subspace for the pressure component. Furthermore, W[0][0] = V[0] is the first component of the velocity space etc.

dolfin::MultiMeshSubSpace::MultiMeshSubSpace(MultiMeshFunctionSpace &V, const std::vector<std::size_t> &component)

Create subspace for given component (n levels)

Parameters:
  • V
  • component
dolfin::MultiMeshSubSpace::MultiMeshSubSpace(MultiMeshFunctionSpace &V, std::size_t component)

Create subspace for given component (one level)

Parameters:
  • V
  • component
dolfin::MultiMeshSubSpace::MultiMeshSubSpace(MultiMeshFunctionSpace &V, std::size_t component, std::size_t sub_component)

Create subspace for given component (two levels)

Parameters:
  • V
  • component
  • sub_component

SpecialFacetFunction

C++ documentation for SpecialFacetFunction from dolfin/function/SpecialFacetFunction.h:

class dolfin::SpecialFacetFunction : public dolfin::Expression

A SpecialFacetFunction is a representation of a global function that is in P(f) for each Facet f in a Mesh for some FunctionSpace P

dolfin::SpecialFacetFunction::SpecialFacetFunction(std::vector<Function> &f_e)

Create (scalar-valued) SpecialFacetFunction() Arguments f_e (std::vector<Function >) Separate _Function_s for each facet

Parameters:f_e
dolfin::SpecialFacetFunction::SpecialFacetFunction(std::vector<Function> &f_e, std::size_t dim)

Create (vector-valued) SpecialFacetFunction() Arguments f_e (std::vector<Function >) Separate _Function_s for each facet dim (int) The value-dimension of the Functions

Parameters:
  • f_e
  • dim
dolfin::SpecialFacetFunction::SpecialFacetFunction(std::vector<Function> &f_e, std::vector<std::size_t> value_shape)

Create (tensor-valued) SpecialFacetFunction() Arguments f_e (std::vector<Function >) Separate _Function_s for each facet value_shape (std::vector<std::size_t>) The values-shape of the Functions

Parameters:
  • f_e
  • value_shape
void dolfin::SpecialFacetFunction::eval(Array<double> &values, const Array<double> &x, const ufc::cell &cell) const

Evaluate SpecialFacetFunction (cf Expression .eval) Evaluate function for given cell

Parameters:
  • values
  • x
  • cell
Function &dolfin::SpecialFacetFunction::operator[](std::size_t i) const

Extract sub-function i Arguments i (int) component Returns:cpp:any:Function

Parameters:i