DOLFINx
0.1.0
DOLFINx C++ interface
|
9 #include <dolfinx/common/array2d.h>
10 #include <dolfinx/fem/evaluate.h>
14 #include <xtl/xspan.hpp>
55 const std::function<
void(T*,
const T*,
const T*,
const double*)> fn,
70 const std::vector<std::shared_ptr<const fem::Function<T>>>&
81 std::vector<int> n{0};
82 for (
const auto& c : _coefficients)
85 throw std::runtime_error(
"Not all form coefficients have been set.");
86 n.push_back(n.back() + c->function_space()->element()->space_dimension());
96 void eval(
const xtl::span<const std::int32_t>& active_cells,
104 const std::function<void(T*,
const T*,
const T*,
const double*)>&
114 const std::vector<std::shared_ptr<const fem::Constant<T>>>&
constants()
const
121 std::shared_ptr<const mesh::Mesh>
mesh()
const {
return _mesh; }
140 std::vector<std::shared_ptr<const fem::Function<T>>> _coefficients;
143 std::vector<std::shared_ptr<const fem::Constant<T>>> _constants;
146 std::function<void(T*,
const T*,
const T*,
const double*)> _fn;
152 std::shared_ptr<const mesh::Mesh> _mesh;
155 std::size_t _value_size;
T scalar_type
Scalar type (T).
Definition: Expression.h:136
void eval(array2d< T > &values, const fem::Expression< T > &e, const xtl::span< const std::int32_t > &active_cells)
Evaluate a UFC expression.
Definition: evaluate.h:28
This class represents a function in a finite element function space , given by.
Definition: Form.h:23
Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell....
Definition: evaluate.h:20
const std::vector< std::shared_ptr< const fem::Function< T > > > & coefficients() const
Access coefficients.
Definition: Expression.h:71
std::shared_ptr< const mesh::Mesh > mesh() const
Get mesh.
Definition: Expression.h:121
void eval(const xtl::span< const std::int32_t > &active_cells, array2d< T > &values) const
Evaluate the expression on cells.
Definition: Expression.h:96
const array2d< double > & x() const
Get evaluation points on reference cell.
Definition: Expression.h:125
const std::vector< std::shared_ptr< const fem::Constant< T > > > & constants() const
Access constants.
Definition: Expression.h:114
std::vector< int > coefficient_offsets() const
Offset for each coefficient expansion array on a cell. Used to pack data for multiple coefficients in...
Definition: Expression.h:79
const std::function< void(T *, const T *, const T *, const double *)> & get_tabulate_expression() const
Get function for tabulate_expression.
Definition: Expression.h:105
const std::size_t num_points() const
Get number of points.
Definition: Expression.h:133
const std::size_t value_size() const
Get value size.
Definition: Expression.h:129
This class provides a dynamic 2-dimensional row-wise array data structure.
Definition: array2d.h:20
Expression(const std::vector< std::shared_ptr< const fem::Function< T >>> &coefficients, const std::vector< std::shared_ptr< const fem::Constant< T >>> &constants, const std::shared_ptr< const mesh::Mesh > &mesh, const array2d< double > &X, const std::function< void(T *, const T *, const T *, const double *)> fn, const std::size_t value_size)
Create Expression.
Definition: Expression.h:51
A constant value which can be attached to a Form. Constants may be scalar (rank 0),...
Definition: Constant.h:18
virtual ~Expression()=default
Destructor.