DOLFIN
DOLFIN C++ interface
|
#include <Expression.h>
Public Member Functions | |
Expression () | |
Create scalar expression. | |
Expression (std::size_t dim) | |
Expression (std::size_t dim0, std::size_t dim1) | |
Expression (std::vector< std::size_t > value_shape) | |
Expression (const Expression &expression) | |
virtual | ~Expression () |
Destructor. | |
virtual void | eval (Array< double > &values, const Array< double > &x, const ufc::cell &cell) const override |
virtual void | eval (Eigen::Ref< Eigen::VectorXd > values, Eigen::Ref< const Eigen::VectorXd > x, const ufc::cell &cell) const override |
virtual void | eval (Array< double > &values, const Array< double > &x) const override |
virtual void | eval (Eigen::Ref< Eigen::VectorXd > values, Eigen::Ref< const Eigen::VectorXd > x) const override |
virtual std::size_t | value_rank () const override |
virtual std::size_t | value_dimension (std::size_t i) const override |
virtual std::vector< std::size_t > | value_shape () const override |
virtual void | set_property (std::string name, double value) |
virtual double | get_property (std::string name) const |
virtual void | set_generic_function (std::string name, std::shared_ptr< GenericFunction > f) |
virtual std::shared_ptr< dolfin::GenericFunction > | get_generic_function (std::string name) const |
virtual void | restrict (double *w, const FiniteElement &element, const Cell &dolfin_cell, const double *coordinate_dofs, const ufc::cell &ufc_cell) const override |
virtual void | compute_vertex_values (std::vector< double > &vertex_values, const Mesh &mesh) const override |
virtual std::shared_ptr< const FunctionSpace > | function_space () const override |
Public Member Functions inherited from dolfin::GenericFunction | |
GenericFunction () | |
Constructor. | |
virtual | ~GenericFunction () |
Destructor. | |
virtual 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. More... | |
virtual void | evaluate (double *values, const double *coordinates, const ufc::cell &cell) const override |
Public Member Functions inherited from dolfin::Variable | |
Variable () | |
Create unnamed variable. | |
Variable (const std::string name, const std::string label) | |
Create variable with given name and label. | |
Variable (const Variable &variable) | |
Copy constructor. | |
virtual | ~Variable () |
Destructor. | |
const Variable & | operator= (const Variable &variable) |
Assignment operator. | |
void | rename (const std::string name, const std::string label) |
Rename variable. | |
std::string | name () const |
Return name. | |
std::string | label () const |
Return label (description) | |
std::size_t | id () const |
virtual std::string | str (bool verbose) const |
Return informal string representation (pretty-print) | |
Protected Attributes | |
std::vector< std::size_t > | _value_shape |
Additional Inherited Members | |
Public Attributes inherited from dolfin::Variable | |
Parameters | parameters |
Parameters. | |
Protected Member Functions inherited from dolfin::GenericFunction | |
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) | |
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.
|
explicit |
Create vector-valued expression with given dimension.
dim | (std::size_t) Dimension of the vector-valued expression. |
Expression::Expression | ( | std::size_t | dim0, |
std::size_t | dim1 | ||
) |
Create matrix-valued expression with given dimensions.
dim0 | (std::size_t) Dimension (rows). |
dim1 | (std::size_t) Dimension (columns). |
|
explicit |
Create tensor-valued expression with given shape.
value_shape | (std::vector<std::size_t>) Shape of expression. |
Expression::Expression | ( | const Expression & | expression | ) |
Copy constructor
expression | (Expression) Object to be copied. |
|
overridevirtual |
Compute values at all mesh vertices.
vertex_values | (Array<double>) The values at all vertices. |
mesh | (Mesh) The mesh. |
Implements dolfin::GenericFunction.
Reimplemented in dolfin::MeshDisplacement.
|
overridevirtual |
Evaluate at given point in given cell (deprecated)
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. |
Reimplemented from dolfin::GenericFunction.
Reimplemented in dolfin::MeshDisplacement, dolfin::SpecialFacetFunction, dolfin::FacetArea, and dolfin::MeshCoordinates.
|
overridevirtual |
Evaluate at given point in given cell
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. |
Reimplemented from dolfin::GenericFunction.
Evaluate at given point (deprecated)
values | (Array<double>) The values at the point. |
x | (Array<double>) The coordinates of the point. |
Reimplemented from dolfin::GenericFunction.
Reimplemented in dolfin::Constant.
|
overridevirtual |
Evaluate at given point.
values | (Eigen::Ref<Eigen::VectorXd>) The values at the point. |
x | (Eigen::Ref<const Eigen::VectorXd>) The coordinates of the point. |
Reimplemented from dolfin::GenericFunction.
Reimplemented in dolfin::Constant.
|
overridevirtual |
Return shared pointer to function space (NULL) Expression does not have a FunctionSpace
Implements dolfin::GenericFunction.
|
virtual |
Property getter for type "GenericFunction" Used in pybind11 Python interface to get the value of a python attribute
|
virtual |
Property getter for type "double" Used in pybind11 Python interface to get the value of a python attribute
|
overridevirtual |
Restrict function to local cell (compute expansion coefficients w).
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. |
Implements dolfin::GenericFunction.
|
virtual |
Property setter for type "GenericFunction" Used in pybind11 Python interface to attach a value to a python attribute
|
virtual |
Property setter for type "double" Used in pybind11 Python interface to attach a value to a python attribute
|
overridevirtual |
Return value dimension for given axis.
i | (std::size_t) Integer denoting the axis to use. |
Implements dolfin::GenericFunction.
|
overridevirtual |
|
overridevirtual |
Return value shape
Implements dolfin::GenericFunction.