Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.2.0a0/v0.9.0/cpp
DOLFINx  0.2.0
DOLFINx C++ interface
Expression.h
1 // Copyright (C) 2020 Jack S. Hale
2 //
3 // This file is part of DOLFINx (https://www.fenicsproject.org)
4 //
5 // SPDX-License-Identifier: LGPL-3.0-or-later
6 
7 #pragma once
8 
9 #include <dolfinx/common/array2d.h>
10 #include <dolfinx/mesh/Mesh.h>
11 #include <functional>
12 #include <utility>
13 #include <vector>
14 #include <xtensor/xarray.hpp>
15 #include <xtl/xspan.hpp>
16 
17 namespace dolfinx::fem
18 {
19 template <typename T>
20 class Constant;
21 
30 
31 template <typename T>
33 {
34 public:
45  const std::vector<std::shared_ptr<const fem::Function<T>>>& coefficients,
46  const std::vector<std::shared_ptr<const fem::Constant<T>>>& constants,
47  const std::shared_ptr<const mesh::Mesh>& mesh, const array2d<double>& X,
48  const std::function<void(T*, const T*, const T*, const double*)> fn,
49  const std::size_t value_size)
50  : _coefficients(coefficients), _constants(constants), _mesh(mesh), _x(X),
51  _fn(fn), _value_size(value_size)
52  {
53  // Do nothing
54  }
55 
57  Expression(Expression&& form) = default;
58 
60  virtual ~Expression() = default;
61 
63  const std::vector<std::shared_ptr<const fem::Function<T>>>&
64  coefficients() const
65  {
66  return _coefficients;
67  }
68 
72  std::vector<int> coefficient_offsets() const
73  {
74  std::vector<int> n{0};
75  for (const auto& c : _coefficients)
76  {
77  if (!c)
78  throw std::runtime_error("Not all form coefficients have been set.");
79  n.push_back(n.back() + c->function_space()->element()->space_dimension());
80  }
81  return n;
82  }
83 
89  template <typename U>
90  void eval(const xtl::span<const std::int32_t>& active_cells, U& values) const
91  {
92  static_assert(std::is_same<T, typename U::value_type>::value,
93  "Expression and array types must be the same");
94 
95  // Extract data from Expression
96  assert(_mesh);
97 
98  // Prepare coefficients and constants
99  const array2d<T> coeffs = pack_coefficients(*this);
100  const std::vector<T> constant_data = pack_constants(*this);
101 
102  const auto& fn = this->get_tabulate_expression();
103 
104  // Prepare cell geometry
105  const graph::AdjacencyList<std::int32_t>& x_dofmap
106  = _mesh->geometry().dofmap();
107  const fem::CoordinateElement& cmap = _mesh->geometry().cmap();
108 
109  // Prepate cell permutation info
110  _mesh->topology_mutable().create_entity_permutations();
111  const std::vector<std::uint32_t>& cell_info
112  = _mesh->topology().get_cell_permutation_info();
113 
114  // FIXME: Add proper interface for num coordinate dofs
115  const std::size_t num_dofs_g = x_dofmap.num_links(0);
116  const xt::xtensor<double, 2>& x_g = _mesh->geometry().x();
117 
118  // Create data structures used in evaluation
119  std::vector<double> coordinate_dofs(3 * num_dofs_g);
120 
121  // Iterate over cells and 'assemble' into values
122  std::vector<T> values_e(this->num_points() * this->value_size(), 0);
123  for (std::size_t c = 0; c < active_cells.size(); ++c)
124  {
125  const std::int32_t cell = active_cells[c];
126 
127  auto x_dofs = x_dofmap.links(cell);
128  for (std::size_t i = 0; i < x_dofs.size(); ++i)
129  {
130  std::copy_n(xt::row(x_g, x_dofs[i]).cbegin(), 3,
131  std::next(coordinate_dofs.begin(), 3 * i));
132  }
133 
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());
138 
139  for (std::size_t j = 0; j < values_e.size(); ++j)
140  values(c, j) = values_e[j];
141  }
142  }
143 
146  const std::function<void(T*, const T*, const T*, const double*)>&
148  {
149  return _fn;
150  }
151 
156  const std::vector<std::shared_ptr<const fem::Constant<T>>>& constants() const
157  {
158  return _constants;
159  }
160 
163  std::shared_ptr<const mesh::Mesh> mesh() const { return _mesh; }
164 
167  const array2d<double>& x() const { return _x; }
168 
171  const std::size_t value_size() const { return _value_size; }
172 
175  const std::size_t num_points() const { return _x.shape[0]; }
176 
178  using scalar_type = T;
179 
180 private:
181  // Coefficients associated with the Expression
182  std::vector<std::shared_ptr<const fem::Function<T>>> _coefficients;
183 
184  // Constants associated with the Expression
185  std::vector<std::shared_ptr<const fem::Constant<T>>> _constants;
186 
187  // Function to evaluate the Expression
188  std::function<void(T*, const T*, const T*, const double*)> _fn;
189 
190  // Evaluation points on reference cell
191  array2d<double> _x;
192 
193  // The mesh.
194  std::shared_ptr<const mesh::Mesh> _mesh;
195 
196  // Evaluation size
197  std::size_t _value_size;
198 };
199 } // namespace dolfinx::fem
dolfinx::fem::Expression::scalar_type
T scalar_type
Scalar type (T).
Definition: Expression.h:178
dolfinx::graph::AdjacencyList::num_links
int num_links(int node) const
Number of connections for given node.
Definition: AdjacencyList.h:120
dolfinx::fem::Function
This class represents a function in a finite element function space , given by.
Definition: Form.h:23
dolfinx::fem::Expression
Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell....
Definition: Expression.h:32
dolfinx::fem::Expression::coefficients
const std::vector< std::shared_ptr< const fem::Function< T > > > & coefficients() const
Access coefficients.
Definition: Expression.h:64
dolfinx::fem::Expression::mesh
std::shared_ptr< const mesh::Mesh > mesh() const
Get mesh.
Definition: Expression.h:163
dolfinx::fem::pack_constants
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
dolfinx::graph::AdjacencyList
This class provides a static adjacency list data structure. It is commonly used to store directed gra...
Definition: AdjacencyList.h:46
dolfinx::fem::Expression::x
const array2d< double > & x() const
Get evaluation points on reference cell.
Definition: Expression.h:167
dolfinx::fem::Expression::constants
const std::vector< std::shared_ptr< const fem::Constant< T > > > & constants() const
Access constants.
Definition: Expression.h:156
dolfinx::fem::pack_coefficients
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
dolfinx::fem::Expression::coefficient_offsets
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
dolfinx::fem::Expression::get_tabulate_expression
const std::function< void(T *, const T *, const T *, const double *)> & get_tabulate_expression() const
Get function for tabulate_expression.
Definition: Expression.h:147
dolfinx::fem::Expression::num_points
const std::size_t num_points() const
Get number of points.
Definition: Expression.h:175
dolfinx::fem::Expression::value_size
const std::size_t value_size() const
Get value size.
Definition: Expression.h:171
dolfinx::fem
Finite element method functionality.
Definition: assemble_matrix_impl.h:22
dolfinx::array2d
This class provides a dynamic 2-dimensional row-wise array data structure.
Definition: array2d.h:20
dolfinx::fem::Expression::Expression
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
dolfinx::array2d::row
constexpr xtl::span< value_type > row(size_type i)
Access a row in the array.
Definition: array2d.h:114
dolfinx::fem::Expression::eval
void eval(const xtl::span< const std::int32_t > &active_cells, U &values) const
Evaluate the expression on cells.
Definition: Expression.h:90
dolfinx::fem::Constant
A constant value which can be attached to a Form. Constants may be scalar (rank 0),...
Definition: Constant.h:18
dolfinx::graph::AdjacencyList::links
xtl::span< T > links(int node)
Get the links (edges) for given node.
Definition: AdjacencyList.h:130
dolfinx::fem::Expression::~Expression
virtual ~Expression()=default
Destructor.
dolfinx::fem::CoordinateElement
This class manages coordinate mappings for isoparametric cells.
Definition: CoordinateElement.h:28