DOLFINx
0.2.0
DOLFINx C++ interface
|
9 #include <dolfinx/common/array2d.h>
10 #include <dolfinx/mesh/Mesh.h>
14 #include <xtensor/xarray.hpp>
15 #include <xtl/xspan.hpp>
48 const std::function<
void(T*,
const T*,
const T*,
const double*)> fn,
63 const std::vector<std::shared_ptr<const fem::Function<T>>>&
74 std::vector<int> n{0};
75 for (
const auto& c : _coefficients)
78 throw std::runtime_error(
"Not all form coefficients have been set.");
79 n.push_back(n.back() + c->function_space()->element()->space_dimension());
90 void eval(
const xtl::span<const std::int32_t>& active_cells, U& values)
const
92 static_assert(std::is_same<T, typename U::value_type>::value,
93 "Expression and array types must be the same");
106 = _mesh->geometry().dofmap();
110 _mesh->topology_mutable().create_entity_permutations();
111 const std::vector<std::uint32_t>& cell_info
112 = _mesh->topology().get_cell_permutation_info();
115 const std::size_t num_dofs_g = x_dofmap.
num_links(0);
116 const xt::xtensor<double, 2>& x_g = _mesh->geometry().x();
119 std::vector<double> coordinate_dofs(3 * num_dofs_g);
123 for (std::size_t c = 0; c < active_cells.size(); ++c)
125 const std::int32_t cell = active_cells[c];
127 auto x_dofs = x_dofmap.
links(cell);
128 for (std::size_t i = 0; i < x_dofs.size(); ++i)
130 std::copy_n(xt::row(x_g, x_dofs[i]).cbegin(), 3,
131 std::next(coordinate_dofs.begin(), 3 * i));
134 auto coeff_cell = coeffs.
row(cell);
135 std::fill(values_e.begin(), values_e.end(), 0.0);
136 fn(values_e.data(), coeff_cell.data(), constant_data.data(),
137 coordinate_dofs.data());
139 for (std::size_t j = 0; j < values_e.size(); ++j)
140 values(c, j) = values_e[j];
146 const std::function<void(T*,
const T*,
const T*,
const double*)>&
156 const std::vector<std::shared_ptr<const fem::Constant<T>>>&
constants()
const
163 std::shared_ptr<const mesh::Mesh>
mesh()
const {
return _mesh; }
182 std::vector<std::shared_ptr<const fem::Function<T>>> _coefficients;
185 std::vector<std::shared_ptr<const fem::Constant<T>>> _constants;
188 std::function<void(T*,
const T*,
const T*,
const double*)> _fn;
194 std::shared_ptr<const mesh::Mesh> _mesh;
197 std::size_t _value_size;
T scalar_type
Scalar type (T).
Definition: Expression.h:178
int num_links(int node) const
Number of connections for given node.
Definition: AdjacencyList.h:120
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: Expression.h:32
const std::vector< std::shared_ptr< const fem::Function< T > > > & coefficients() const
Access coefficients.
Definition: Expression.h:64
std::shared_ptr< const mesh::Mesh > mesh() const
Get mesh.
Definition: Expression.h:163
std::vector< typename U::scalar_type > pack_constants(const U &u)
Pack constants of u of generic type U ready for assembly.
Definition: utils.h:503
This class provides a static adjacency list data structure. It is commonly used to store directed gra...
Definition: AdjacencyList.h:46
const array2d< double > & x() const
Get evaluation points on reference cell.
Definition: Expression.h:167
const std::vector< std::shared_ptr< const fem::Constant< T > > > & constants() const
Access constants.
Definition: Expression.h:156
array2d< typename U::scalar_type > pack_coefficients(const U &u)
Pack coefficients of u of generic type U ready for assembly.
Definition: utils.h:418
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:72
const std::function< void(T *, const T *, const T *, const double *)> & get_tabulate_expression() const
Get function for tabulate_expression.
Definition: Expression.h:147
const std::size_t num_points() const
Get number of points.
Definition: Expression.h:175
const std::size_t value_size() const
Get value size.
Definition: Expression.h:171
Finite element method functionality.
Definition: assemble_matrix_impl.h:22
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 an Expression.
Definition: Expression.h:44
constexpr xtl::span< value_type > row(size_type i)
Access a row in the array.
Definition: array2d.h:114
void eval(const xtl::span< const std::int32_t > &active_cells, U &values) const
Evaluate the expression on cells.
Definition: Expression.h:90
A constant value which can be attached to a Form. Constants may be scalar (rank 0),...
Definition: Constant.h:18
xtl::span< T > links(int node)
Get the links (edges) for given node.
Definition: AdjacencyList.h:130
virtual ~Expression()=default
Destructor.
This class manages coordinate mappings for isoparametric cells.
Definition: CoordinateElement.h:28