DOLFIN-X
DOLFIN-X C++ interface
Function.h
1 // Copyright (C) 2003-2020 Anders Logg and Garth N. Wells
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 "FunctionSpace.h"
10 #include "interpolate.h"
11 #include <dolfinx/common/IndexMap.h>
12 #include <dolfinx/common/UniqueIdGenerator.h>
13 #include <dolfinx/common/array2d.h>
14 #include <dolfinx/fem/DofMap.h>
15 #include <dolfinx/fem/FiniteElement.h>
16 #include <dolfinx/la/PETScVector.h>
17 #include <dolfinx/la/Vector.h>
18 #include <dolfinx/mesh/Geometry.h>
19 #include <dolfinx/mesh/Mesh.h>
20 #include <dolfinx/mesh/Topology.h>
21 #include <functional>
22 #include <memory>
23 #include <numeric>
24 #include <petscvec.h>
25 #include <string>
26 #include <utility>
27 #include <variant>
28 #include <vector>
29 
30 namespace dolfinx::fem
31 {
32 
33 class FunctionSpace;
34 
41 
42 template <typename T>
43 class Function
44 {
45 public:
48  explicit Function(std::shared_ptr<const FunctionSpace> V)
49  : _id(common::UniqueIdGenerator::id()), _function_space(V),
50  _x(std::make_shared<la::Vector<T>>(V->dofmap()->index_map,
51  V->dofmap()->index_map_bs()))
52  {
53  if (!V->component().empty())
54  {
55  throw std::runtime_error("Cannot create Function from subspace. Consider "
56  "collapsing the function space");
57  }
58  }
59 
66  Function(std::shared_ptr<const FunctionSpace> V,
67  std::shared_ptr<la::Vector<T>> x)
68  : _id(common::UniqueIdGenerator::id()), _function_space(V), _x(x)
69  {
70  // We do not check for a subspace since this constructor is used for
71  // creating subfunctions
72 
73  // Assertion uses '<=' to deal with sub-functions
74  assert(V->dofmap());
75  assert(V->dofmap()->index_map->size_global() * V->dofmap()->index_map_bs()
76  <= _x->bs() * _x->map()->size_global());
77  }
78 
79  // Copy constructor
80  Function(const Function& v) = delete;
81 
84  : name(std::move(v.name)), _id(std::move(v._id)),
85  _function_space(std::move(v._function_space)), _x(std::move(v._x)),
86  _petsc_vector(std::exchange(v._petsc_vector, nullptr))
87  {
88  }
89 
91  virtual ~Function()
92  {
93  if (_petsc_vector)
94  VecDestroy(&_petsc_vector);
95  }
96 
98  Function& operator=(Function&& v) noexcept
99  {
100  name = std::move(v.name);
101  _id = std::move(v._id);
102  _function_space = std::move(v._function_space);
103  _x = std::move(v._x);
104  std::swap(_petsc_vector, v._petsc_vector);
105 
106  return *this;
107  }
108 
109  // Assignment
110  Function& operator=(const Function& v) = delete;
111 
115  Function sub(int i) const
116  {
117  auto sub_space = _function_space->sub({i});
118  assert(sub_space);
119  return Function(sub_space, _x);
120  }
121 
126  {
127  // Create new collapsed FunctionSpace
128  const auto [function_space_new, collapsed_map]
129  = _function_space->collapse();
130 
131  // Create new vector
132  assert(function_space_new);
133  auto vector_new = std::make_shared<la::Vector<T>>(
134  function_space_new->dofmap()->index_map,
135  function_space_new->dofmap()->index_map_bs());
136 
137  // Copy values into new vector
138  const std::vector<T>& x_old = _x->array();
139  std::vector<T>& x_new = vector_new->mutable_array();
140  for (std::size_t i = 0; i < collapsed_map.size(); ++i)
141  {
142  assert((int)i < x_new.size());
143  assert(collapsed_map[i] < x_old.size());
144  x_new[i] = x_old[collapsed_map[i]];
145  }
146 
147  return Function(function_space_new, vector_new);
148  }
149 
152  std::shared_ptr<const FunctionSpace> function_space() const
153  {
154  return _function_space;
155  }
156 
160  Vec vector() const
161  {
162  // Check that this is not a sub function
163  assert(_function_space->dofmap());
164  assert(_function_space->dofmap()->index_map);
165  if (_x->bs() * _x->map()->size_global()
166  != _function_space->dofmap()->index_map->size_global()
167  * _function_space->dofmap()->index_map_bs())
168  {
169  throw std::runtime_error(
170  "Cannot access a non-const vector from a subfunction");
171  }
172 
173  // Check that data type is the same as the PETSc build
174  if constexpr (std::is_same<T, PetscScalar>::value)
175  {
176  if (!_petsc_vector)
177  {
178  _petsc_vector = la::create_ghosted_vector(
179  *_function_space->dofmap()->index_map,
180  _function_space->dofmap()->index_map_bs(), _x->mutable_array());
181  }
182 
183  return _petsc_vector;
184  }
185  else
186  {
187  throw std::runtime_error(
188  "Cannot return PETSc vector wrapper. Type mismatch");
189  }
190  }
191 
193  std::shared_ptr<const la::Vector<T>> x() const { return _x; }
194 
196  std::shared_ptr<la::Vector<T>> x() { return _x; }
197 
200  void interpolate(const Function<T>& v) { fem::interpolate(*this, v); }
201 
204  void interpolate(const std::function<std::variant<std::vector<T>, array2d<T>>(
205  const array2d<double>&)>& f)
206  {
207  assert(_function_space);
208  assert(_function_space->element());
209  assert(_function_space->mesh());
210  const int tdim = _function_space->mesh()->topology().dim();
211  auto cell_map = _function_space->mesh()->topology().index_map(tdim);
212  assert(cell_map);
213  const int num_cells = cell_map->size_local() + cell_map->num_ghosts();
214  std::vector<std::int32_t> cells(num_cells, 0);
215  std::iota(cells.begin(), cells.end(), 0);
217  *_function_space->element(), *_function_space->mesh(), cells);
218  fem::interpolate(*this, f, x, cells);
219  }
220 
231  void eval(const array2d<double>& x,
232  const tcb::span<const std::int32_t>& cells, array2d<T>& u) const
233  {
234  // TODO: This could be easily made more efficient by exploiting points
235  // being ordered by the cell to which they belong.
236 
237  if (x.shape[0] != cells.size())
238  {
239  throw std::runtime_error(
240  "Number of points and number of cells must be equal.");
241  }
242  if (x.shape[0] != u.shape[0])
243  {
244  throw std::runtime_error(
245  "Length of array for Function values must be the "
246  "same as the number of points.");
247  }
248 
249  // Get mesh
250  assert(_function_space);
251  std::shared_ptr<const mesh::Mesh> mesh = _function_space->mesh();
252  assert(mesh);
253  const int gdim = mesh->geometry().dim();
254  const int tdim = mesh->topology().dim();
255  auto map = mesh->topology().index_map(tdim);
256 
257  // Get geometry data
258  const graph::AdjacencyList<std::int32_t>& x_dofmap
259  = mesh->geometry().dofmap();
260  // FIXME: Add proper interface for num coordinate dofs
261  const int num_dofs_g = x_dofmap.num_links(0);
262  const array2d<double>& x_g = mesh->geometry().x();
263 
264  // Get coordinate map
265  const fem::CoordinateElement& cmap = mesh->geometry().cmap();
266 
267  // Get element
268  assert(_function_space->element());
269  std::shared_ptr<const fem::FiniteElement> element
270  = _function_space->element();
271  assert(element);
272  const int bs_element = element->block_size();
273  const int reference_value_size
274  = element->reference_value_size() / bs_element;
275  const int value_size = element->value_size() / bs_element;
276  const int space_dimension = element->space_dimension() / bs_element;
277 
278  // If the space has sub elements, concatenate the evaluations on the sub
279  // elements
280  const int num_sub_elements = element->num_sub_elements();
281  if (num_sub_elements > 1 and num_sub_elements != bs_element)
282  {
283  throw std::runtime_error("Function::eval is not supported for mixed "
284  "elements. Extract subspaces.");
285  }
286 
287  // Prepare geometry data structures
288  std::vector<double> J(gdim * tdim);
289  std::array<double, 1> detJ;
290  std::vector<double> K(tdim * gdim);
291  array2d<double> X(1, tdim);
292 
293  // Prepare basis function data structures
294  std::vector<double> basis_reference_values(space_dimension
295  * reference_value_size);
296  std::vector<double> basis_values(space_dimension * value_size);
297 
298  // Create work vector for expansion coefficients
299  std::vector<T> coefficients(space_dimension * bs_element);
300 
301  // Get dofmap
302  std::shared_ptr<const fem::DofMap> dofmap = _function_space->dofmap();
303  assert(dofmap);
304  const int bs_dof = dofmap->bs();
305 
306  mesh->topology_mutable().create_entity_permutations();
307  const std::vector<std::uint32_t>& cell_info
308  = mesh->topology().get_cell_permutation_info();
309  array2d<double> coordinate_dofs(num_dofs_g, gdim);
310 
311  array2d<double> xp(1, gdim);
312 
313  // Loop over points
314  std::fill(u.data(), u.data() + u.size(), 0.0);
315  const std::vector<T>& _v = _x->mutable_array();
316  for (std::size_t p = 0; p < cells.size(); ++p)
317  {
318  const int cell_index = cells[p];
319 
320  // Skip negative cell indices
321  if (cell_index < 0)
322  continue;
323 
324  // Get cell geometry (coordinate dofs)
325  auto x_dofs = x_dofmap.links(cell_index);
326  for (int i = 0; i < num_dofs_g; ++i)
327  for (int j = 0; j < gdim; ++j)
328  coordinate_dofs(i, j) = x_g(x_dofs[i], j);
329 
330  for (int j = 0; j < gdim; ++j)
331  xp(0, j) = x(p, j);
332 
333  // Compute reference coordinates X, and J, detJ and K
334  cmap.compute_reference_geometry(X, J, detJ, K, xp, coordinate_dofs);
335 
336  // Compute basis on reference element
337  element->evaluate_reference_basis(basis_reference_values, X);
338 
339  // Permute the reference values to account for the cell's orientation
340  element->apply_dof_transformation(basis_reference_values.data(),
341  cell_info[cell_index],
342  reference_value_size);
343 
344  // Push basis forward to physical element
345  element->transform_reference_basis(basis_values, basis_reference_values,
346  X, J, detJ, K);
347 
348  // Get degrees of freedom for current cell
349  tcb::span<const std::int32_t> dofs = dofmap->cell_dofs(cell_index);
350  for (std::size_t i = 0; i < dofs.size(); ++i)
351  for (int k = 0; k < bs_dof; ++k)
352  coefficients[bs_dof * i + k] = _v[bs_dof * dofs[i] + k];
353 
354  // Compute expansion
355  auto u_row = u.row(p);
356  for (int k = 0; k < bs_element; ++k)
357  {
358  for (int i = 0; i < space_dimension; ++i)
359  {
360  for (int j = 0; j < value_size; ++j)
361  {
362  u_row[j * bs_element + k] += coefficients[bs_element * i + k]
363  * basis_values[i * value_size + j];
364  }
365  }
366  }
367  }
368  }
369 
373  {
374  assert(_function_space);
375  std::shared_ptr<const mesh::Mesh> mesh = _function_space->mesh();
376  assert(mesh);
377  const int tdim = mesh->topology().dim();
378 
379  // Compute in tensor (one for scalar function, . . .)
380  const int value_size_loc = _function_space->element()->value_size();
381 
382  // Resize Array for holding point values
383  array2d<T> point_values(mesh->geometry().x().shape[0], value_size_loc);
384 
385  // Prepare cell geometry
386  const graph::AdjacencyList<std::int32_t>& x_dofmap
387  = mesh->geometry().dofmap();
388 
389  // FIXME: Add proper interface for num coordinate dofs
390  const int num_dofs_g = x_dofmap.num_links(0);
391  const array2d<double>& x_g = mesh->geometry().x();
392 
393  // Interpolate point values on each cell (using last computed value if
394  // not continuous, e.g. discontinuous Galerkin methods)
395  auto map = mesh->topology().index_map(tdim);
396  assert(map);
397  const std::int32_t num_cells = map->size_local() + map->num_ghosts();
398 
399  std::vector<std::int32_t> cells(x_g.shape[0]);
400  for (std::int32_t c = 0; c < num_cells; ++c)
401  {
402  // Get coordinates for all points in cell
403  tcb::span<const std::int32_t> dofs = x_dofmap.links(c);
404  for (int i = 0; i < num_dofs_g; ++i)
405  cells[dofs[i]] = c;
406  }
407 
408  eval(x_g, cells, point_values);
409 
410  return point_values;
411  }
412 
414  std::string name = "u";
415 
417  std::size_t id() const { return _id; }
418 
419 private:
420  // ID
421  std::size_t _id;
422 
423  // The function space
424  std::shared_ptr<const FunctionSpace> _function_space;
425 
426  // The vector of expansion coefficients (local)
427  std::shared_ptr<la::Vector<T>> _x;
428 
429  // PETSc wrapper of the expansion coefficients
430  mutable Vec _petsc_vector = nullptr;
431 };
432 } // namespace dolfinx::fem
This class provides a dynamic 2-dimensional row-wise array data structure.
Definition: array2d.h:21
std::array< size_type, 2 > shape
The shape of the array.
Definition: array2d.h:157
constexpr tcb::span< value_type > row(size_type i)
Access a row in the array.
Definition: array2d.h:116
constexpr value_type * data() noexcept
Get pointer to the first element of the underlying storage.
Definition: array2d.h:133
constexpr size_type size() const noexcept
Returns the number of elements in the array.
Definition: array2d.h:144
This class manages coordinate mappings for isoparametric cells.
Definition: CoordinateElement.h:25
void compute_reference_geometry(array2d< double > &X, std::vector< double > &J, tcb::span< double > detJ, std::vector< double > &K, const array2d< double > &x, const array2d< double > &cell_geometry) const
Compute reference coordinates X, and J, detJ and K for physical coordinates x.
Definition: CoordinateElement.cpp:89
This class represents a function in a finite element function space , given by.
Definition: Function.h:44
void interpolate(const std::function< std::variant< std::vector< T >, array2d< T >>(const array2d< double > &)> &f)
Interpolate an expression.
Definition: Function.h:204
Vec vector() const
Return vector of expansion coefficients as a PETSc Vec. Throws an exception a PETSc Vec cannot be cre...
Definition: Function.h:160
Function(std::shared_ptr< const FunctionSpace > V, std::shared_ptr< la::Vector< T >> x)
Create function on given function space with a given vector.
Definition: Function.h:66
void interpolate(const Function< T > &v)
Interpolate a Function (on possibly non-matching meshes)
Definition: Function.h:200
Function & operator=(Function &&v) noexcept
Move assignment.
Definition: Function.h:98
Function(Function &&v)
Move constructor.
Definition: Function.h:83
void eval(const array2d< double > &x, const tcb::span< const std::int32_t > &cells, array2d< T > &u) const
Evaluate the Function at points.
Definition: Function.h:231
std::shared_ptr< la::Vector< T > > x()
Underlying vector.
Definition: Function.h:196
std::shared_ptr< const FunctionSpace > function_space() const
Return shared pointer to function space.
Definition: Function.h:152
virtual ~Function()
Destructor.
Definition: Function.h:91
std::size_t id() const
ID.
Definition: Function.h:417
array2d< T > compute_point_values() const
Compute values at all mesh 'nodes'.
Definition: Function.h:372
std::string name
Name.
Definition: Function.h:414
std::shared_ptr< const la::Vector< T > > x() const
Underlying vector.
Definition: Function.h:193
Function collapse() const
Collapse a subfunction (view into the Function) to a stand-alone Function.
Definition: Function.h:125
Function(std::shared_ptr< const FunctionSpace > V)
Create function on given function space.
Definition: Function.h:48
Function sub(int i) const
Extract subfunction (view into the Function)
Definition: Function.h:115
This class provides a static adjacency list data structure. It is commonly used to store directed gra...
Definition: AdjacencyList.h:68
tcb::span< T > links(int node)
Get the links (edges) for given node.
Definition: AdjacencyList.h:151
int num_links(int node) const
Number of connections for given node.
Definition: AdjacencyList.h:141
Distributed vector.
Definition: Vector.h:19
Finite element method functionality.
Definition: assemble_matrix_impl.h:23
array2d< double > interpolation_coords(const fem::FiniteElement &element, const mesh::Mesh &mesh, const tcb::span< const std::int32_t > &cells)
Compute the evaluation points in the physical space at which an expression should be computed to inte...
Definition: interpolate.cpp:15
void interpolate(Function< T > &u, const Function< T > &v)
Interpolate a finite element Function (on possibly non-matching meshes) in another finite element spa...
Definition: interpolate.h:142
Vec create_ghosted_vector(const common::IndexMap &map, int bs, tcb::span< PetscScalar > x)
Create a PETSc Vec that wraps the data in an array.
Definition: PETScVector.cpp:64