Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.1.0/v0.9.0/cpp
DOLFINx  0.1.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/fem/evaluate.h>
11 #include <functional>
12 #include <utility>
13 #include <vector>
14 #include <xtl/xspan.hpp>
15 
16 namespace dolfinx
17 {
18 
19 namespace mesh
20 {
21 class Mesh;
22 }
23 
24 namespace fem
25 {
26 template <typename T>
27 class Constant;
28 
37 
38 template <typename T>
39 class Expression
40 {
41 public:
52  const std::vector<std::shared_ptr<const fem::Function<T>>>& coefficients,
53  const std::vector<std::shared_ptr<const fem::Constant<T>>>& constants,
54  const std::shared_ptr<const mesh::Mesh>& mesh, const array2d<double>& X,
55  const std::function<void(T*, const T*, const T*, const double*)> fn,
56  const std::size_t value_size)
57  : _coefficients(coefficients), _constants(constants), _mesh(mesh), _x(X),
58  _fn(fn), _value_size(value_size)
59  {
60  // Do nothing
61  }
62 
64  Expression(Expression&& form) = default;
65 
67  virtual ~Expression() = default;
68 
70  const std::vector<std::shared_ptr<const fem::Function<T>>>&
71  coefficients() const
72  {
73  return _coefficients;
74  }
75 
79  std::vector<int> coefficient_offsets() const
80  {
81  std::vector<int> n{0};
82  for (const auto& c : _coefficients)
83  {
84  if (!c)
85  throw std::runtime_error("Not all form coefficients have been set.");
86  n.push_back(n.back() + c->function_space()->element()->space_dimension());
87  }
88  return n;
89  }
90 
96  void eval(const xtl::span<const std::int32_t>& active_cells,
97  array2d<T>& values) const
98  {
99  fem::eval(values, *this, active_cells);
100  }
101 
104  const std::function<void(T*, const T*, const T*, const double*)>&
106  {
107  return _fn;
108  }
109 
114  const std::vector<std::shared_ptr<const fem::Constant<T>>>& constants() const
115  {
116  return _constants;
117  }
118 
121  std::shared_ptr<const mesh::Mesh> mesh() const { return _mesh; }
122 
125  const array2d<double>& x() const { return _x; }
126 
129  const std::size_t value_size() const { return _value_size; }
130 
133  const std::size_t num_points() const { return _x.shape[0]; }
134 
136  using scalar_type = T;
137 
138 private:
139  // Coefficients associated with the Expression
140  std::vector<std::shared_ptr<const fem::Function<T>>> _coefficients;
141 
142  // Constants associated with the Expression
143  std::vector<std::shared_ptr<const fem::Constant<T>>> _constants;
144 
145  // Function to evaluate the Expression
146  std::function<void(T*, const T*, const T*, const double*)> _fn;
147 
148  // Evaluation points on reference cell
149  array2d<double> _x;
150 
151  // The mesh.
152  std::shared_ptr<const mesh::Mesh> _mesh;
153 
154  // Evaluation size
155  std::size_t _value_size;
156 };
157 } // namespace fem
158 } // namespace dolfinx
dolfinx::fem::Expression::scalar_type
T scalar_type
Scalar type (T).
Definition: Expression.h:136
dolfinx::fem::eval
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
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: evaluate.h:20
dolfinx::fem::Expression::coefficients
const std::vector< std::shared_ptr< const fem::Function< T > > > & coefficients() const
Access coefficients.
Definition: Expression.h:71
dolfinx::fem::Expression::mesh
std::shared_ptr< const mesh::Mesh > mesh() const
Get mesh.
Definition: Expression.h:121
dolfinx::fem::Expression::eval
void eval(const xtl::span< const std::int32_t > &active_cells, array2d< T > &values) const
Evaluate the expression on cells.
Definition: Expression.h:96
dolfinx::fem::Expression::x
const array2d< double > & x() const
Get evaluation points on reference cell.
Definition: Expression.h:125
dolfinx::fem::Expression::constants
const std::vector< std::shared_ptr< const fem::Constant< T > > > & constants() const
Access constants.
Definition: Expression.h:114
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:79
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:105
dolfinx::fem::Expression::num_points
const std::size_t num_points() const
Get number of points.
Definition: Expression.h:133
dolfinx::fem::Expression::value_size
const std::size_t value_size() const
Get value size.
Definition: Expression.h:129
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 Expression.
Definition: Expression.h:51
dolfinx::fem::Constant
A constant value which can be attached to a Form. Constants may be scalar (rank 0),...
Definition: Constant.h:18
dolfinx::fem::Expression::~Expression
virtual ~Expression()=default
Destructor.