23 #ifndef __GENERIC_FUNCTION_H 24 #define __GENERIC_FUNCTION_H 27 #include <Eigen/Dense> 29 #include <dolfin/common/Array.h> 30 #include <dolfin/common/Variable.h> 72 virtual std::vector<std::size_t>
value_shape()
const = 0;
76 const ufc::cell& cell)
const;
82 virtual void eval(Eigen::Ref<Eigen::VectorXd> values,
83 Eigen::Ref<const Eigen::VectorXd> x,
84 const ufc::cell& cell)
const;
87 virtual void eval(Eigen::Ref<Eigen::VectorXd> values,
88 Eigen::Ref<const Eigen::VectorXd> x)
const;
93 const Cell& dolfin_cell,
94 const double* coordinate_dofs,
95 const ufc::cell& ufc_cell)
const = 0;
99 const Mesh& mesh)
const = 0;
115 double operator() (
double x,
double y,
double z)
const;
143 virtual void evaluate(
double* values,
144 const double* coordinates,
145 const ufc::cell& cell)
const override;
148 virtual std::shared_ptr<const FunctionSpace>
function_space()
const = 0;
155 const Cell& dolfin_cell,
156 const double* coordinate_dofs,
157 const ufc::cell& ufc_cell)
const;
virtual 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)
Common base class for DOLFIN variables.
Definition: Variable.h:35
virtual ~GenericFunction()
Destructor.
Definition: GenericFunction.cpp:36
virtual void evaluate(double *values, const double *coordinates, const ufc::cell &cell) const override
Definition: GenericFunction.cpp:189
virtual void eval(Array< double > &values, const Array< double > &x, const ufc::cell &cell) const
Evaluate at given point in given cell (deprecated)
Definition: GenericFunction.cpp:41
virtual void compute_vertex_values(std::vector< double > &vertex_values, const Mesh &mesh) const =0
Compute values at all mesh vertices.
std::size_t value_size() const
Evaluation at given point.
Definition: GenericFunction.cpp:181
double operator()(double x) const
Evaluation at given point (scalar function)
Definition: GenericFunction.cpp:71
virtual std::vector< std::size_t > value_shape() const =0
Return value shape.
A Cell is a MeshEntity of topological codimension 0.
Definition: Cell.h:42
GenericFunction()
Constructor.
Definition: GenericFunction.cpp:31
virtual std::size_t value_rank() const =0
Return value rank.
virtual void update() const
Update off-process ghost coefficients.
Definition: GenericFunction.h:104
virtual std::size_t value_dimension(std::size_t i) const =0
Return value dimension for given axis.
Definition: GenericFunction.h:53
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)
Definition: GenericFunction.cpp:205
This is a wrapper for a UFC finite element (ufc::finite_element).
Definition: FiniteElement.h:35
virtual std::shared_ptr< const FunctionSpace > function_space() const =0
Pointer to FunctionSpace, if appropriate, otherwise NULL.