DOLFIN-X
DOLFIN-X 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 
15 namespace dolfinx
16 {
17 
18 namespace mesh
19 {
20 class Mesh;
21 }
22 
23 namespace fem
24 {
25 template <typename T>
26 class Constant;
27 
36 
37 template <typename T>
39 {
40 public:
51  const std::vector<std::shared_ptr<const fem::Function<T>>>& coefficients,
52  const std::vector<std::shared_ptr<const fem::Constant<T>>>& constants,
53  const std::shared_ptr<const mesh::Mesh>& mesh, const array2d<double>& X,
54  const std::function<void(T*, const T*, const T*, const double*)> fn,
55  const std::size_t value_size)
56  : _coefficients(coefficients), _constants(constants), _mesh(mesh), _x(X),
57  _fn(fn), _value_size(value_size)
58  {
59  // Do nothing
60  }
61 
63  Expression(Expression&& form) = default;
64 
66  virtual ~Expression() = default;
67 
69  const std::vector<std::shared_ptr<const fem::Function<T>>>&
70  coefficients() const
71  {
72  return _coefficients;
73  }
74 
78  std::vector<int> coefficient_offsets() const
79  {
80  std::vector<int> n{0};
81  for (const auto& c : _coefficients)
82  {
83  if (!c)
84  throw std::runtime_error("Not all form coefficients have been set.");
85  n.push_back(n.back() + c->function_space()->element()->space_dimension());
86  }
87  return n;
88  }
89 
95  void eval(const tcb::span<const std::int32_t>& active_cells,
96  array2d<T>& values) const
97  {
98  fem::eval(values, *this, active_cells);
99  }
100 
103  const std::function<void(T*, const T*, const T*, const double*)>&
105  {
106  return _fn;
107  }
108 
113  const std::vector<std::shared_ptr<const fem::Constant<T>>>& constants() const
114  {
115  return _constants;
116  }
117 
120  std::shared_ptr<const mesh::Mesh> mesh() const { return _mesh; }
121 
124  const array2d<double>& x() const { return _x; }
125 
128  const std::size_t value_size() const { return _value_size; }
129 
132  const std::size_t num_points() const { return _x.shape[0]; }
133 
135  using scalar_type = T;
136 
137 private:
138  // Coefficients associated with the Expression
139  std::vector<std::shared_ptr<const fem::Function<T>>> _coefficients;
140 
141  // Constants associated with the Expression
142  std::vector<std::shared_ptr<const fem::Constant<T>>> _constants;
143 
144  // Function to evaluate the Expression
145  std::function<void(T*, const T*, const T*, const double*)> _fn;
146 
147  // Evaluation points on reference cell
148  array2d<double> _x;
149 
150  // The mesh.
151  std::shared_ptr<const mesh::Mesh> _mesh;
152 
153  // Evaluation size
154  std::size_t _value_size;
155 };
156 } // namespace fem
157 } // namespace dolfinx
This class provides a dynamic 2-dimensional row-wise array data structure.
Definition: array2d.h:21
A constant value which can be attached to a Form. Constants may be scalar (rank 0),...
Definition: Constant.h:19
Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell....
Definition: Expression.h:39
const std::size_t value_size() const
Get value size.
Definition: Expression.h:128
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:50
const std::vector< std::shared_ptr< const fem::Constant< T > > > & constants() const
Access constants.
Definition: Expression.h:113
Expression(Expression &&form)=default
Move constructor.
const std::vector< std::shared_ptr< const fem::Function< T > > > & coefficients() const
Access coefficients.
Definition: Expression.h:70
std::shared_ptr< const mesh::Mesh > mesh() const
Get mesh.
Definition: Expression.h:120
const array2d< double > & x() const
Get evaluation points on reference cell.
Definition: Expression.h:124
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:78
virtual ~Expression()=default
Destructor.
const std::size_t num_points() const
Get number of points.
Definition: Expression.h:132
const std::function< void(T *, const T *, const T *, const double *)> & get_tabulate_expression() const
Get function for tabulate_expression.
Definition: Expression.h:104
T scalar_type
Scalar type (T).
Definition: Expression.h:135
void eval(const tcb::span< const std::int32_t > &active_cells, array2d< T > &values) const
Evaluate the expression on cells.
Definition: Expression.h:95
This class represents a function in a finite element function space , given by.
Definition: Function.h:46
void eval(array2d< T > &values, const fem::Expression< T > &e, const tcb::span< const std::int32_t > &active_cells)
Evaluate a UFC expression.
Definition: evaluate.h:28