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
interpolate.h
1 // Copyright (C) 2020-2021 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 <dolfinx/fem/DofMap.h>
11 #include <dolfinx/fem/FiniteElement.h>
12 #include <dolfinx/mesh/Mesh.h>
13 #include <functional>
14 #include <variant>
15 #include <xtensor/xadapt.hpp>
16 #include <xtensor/xtensor.hpp>
17 #include <xtensor/xview.hpp>
18 #include <xtl/xspan.hpp>
19 
20 namespace dolfinx::fem
21 {
22 
23 template <typename T>
24 class Function;
25 
36 xt::xtensor<double, 2>
37 interpolation_coords(const fem::FiniteElement& element, const mesh::Mesh& mesh,
38  const xtl::span<const std::int32_t>& cells);
39 
44 template <typename T>
45 void interpolate(Function<T>& u, const Function<T>& v);
46 
58 template <typename T>
59 void interpolate(
60  Function<T>& u,
61  const std::function<xt::xarray<T>(const xt::xtensor<double, 2>&)>& f,
62  const xt::xtensor<double, 2>& x,
63  const xtl::span<const std::int32_t>& cells);
64 
81 template <typename T>
82 void interpolate_c(
83  Function<T>& u,
84  const std::function<void(xt::xarray<T>&, const xt::xtensor<double, 2>&)>& f,
85  const xt::xtensor<double, 2>& x,
86  const xtl::span<const std::int32_t>& cells);
87 
88 namespace detail
89 {
90 
91 template <typename T>
92 void interpolate_from_any(Function<T>& u, const Function<T>& v)
93 {
94  assert(v.function_space());
95  const auto element = u.function_space()->element();
96  assert(element);
97  if (v.function_space()->element()->hash() != element->hash())
98  {
99  throw std::runtime_error("Restricting finite elements function in "
100  "different elements not supported.");
101  }
102 
103  const auto mesh = u.function_space()->mesh();
104  assert(mesh);
105  assert(v.function_space()->mesh());
106  if (mesh->id() != v.function_space()->mesh()->id())
107  {
108  throw std::runtime_error(
109  "Interpolation on different meshes not supported (yet).");
110  }
111  const int tdim = mesh->topology().dim();
112 
113  // Get dofmaps
114  assert(v.function_space());
115  std::shared_ptr<const fem::DofMap> dofmap_v = v.function_space()->dofmap();
116  assert(dofmap_v);
117  auto map = mesh->topology().index_map(tdim);
118  assert(map);
119 
120  std::vector<T>& coeffs = u.x()->mutable_array();
121 
122  // Iterate over mesh and interpolate on each cell
123  const auto dofmap_u = u.function_space()->dofmap();
124  const std::vector<T>& v_array = v.x()->array();
125  const int num_cells = map->size_local() + map->num_ghosts();
126  const int bs = dofmap_v->bs();
127  assert(bs == dofmap_u->bs());
128  for (int c = 0; c < num_cells; ++c)
129  {
130  xtl::span<const std::int32_t> dofs_v = dofmap_v->cell_dofs(c);
131  xtl::span<const std::int32_t> cell_dofs = dofmap_u->cell_dofs(c);
132  assert(dofs_v.size() == cell_dofs.size());
133  for (std::size_t i = 0; i < dofs_v.size(); ++i)
134  {
135  for (int k = 0; k < bs; ++k)
136  coeffs[bs * cell_dofs[i] + k] = v_array[bs * dofs_v[i] + k];
137  }
138  }
139 }
140 
141 } // namespace detail
142 
143 //----------------------------------------------------------------------------
144 template <typename T>
146 {
147  assert(u.function_space());
148  const std::shared_ptr<const FiniteElement> element
149  = u.function_space()->element();
150  assert(element);
151 
152  // Check that function ranks match
153  if (int rank_v = v.function_space()->element()->value_rank();
154  element->value_rank() != rank_v)
155  {
156  throw std::runtime_error("Cannot interpolate function into function space. "
157  "Rank of function ("
158  + std::to_string(rank_v)
159  + ") does not match rank of function space ("
160  + std::to_string(element->value_rank()) + ")");
161  }
162 
163  // Check that function dimension match
164  for (int i = 0; i < element->value_rank(); ++i)
165  {
166  if (int v_dim = v.function_space()->element()->value_dimension(i);
167  element->value_dimension(i) != v_dim)
168  {
169  throw std::runtime_error(
170  "Cannot interpolate function into function space. "
171  "Dimension "
172  + std::to_string(i) + " of function (" + std::to_string(v_dim)
173  + ") does not match dimension " + std::to_string(i)
174  + " of function space(" + std::to_string(element->value_dimension(i))
175  + ")");
176  }
177  }
178 
179  detail::interpolate_from_any(u, v);
180 }
181 //----------------------------------------------------------------------------
182 template <typename T>
184  Function<T>& u,
185  const std::function<xt::xarray<T>(const xt::xtensor<double, 2>&)>& f,
186  const xt::xtensor<double, 2>& x, const xtl::span<const std::int32_t>& cells)
187 {
188  const std::shared_ptr<const FiniteElement> element
189  = u.function_space()->element();
190  assert(element);
191  const int element_bs = element->block_size();
192  if (int num_sub = element->num_sub_elements();
193  num_sub > 0 and num_sub != element_bs)
194  {
195  throw std::runtime_error("Cannot directly interpolate a mixed space. "
196  "Interpolate into subspaces.");
197  }
198 
199  // Get mesh
200  assert(u.function_space());
201  auto mesh = u.function_space()->mesh();
202  assert(mesh);
203 
204  const int gdim = mesh->geometry().dim();
205  const int tdim = mesh->topology().dim();
206 
207  // Get the interpolation points on the reference cells
208  const xt::xtensor<double, 2>& X = element->interpolation_points();
209 
210  if (X.shape(0) == 0)
211  throw std::runtime_error(
212  "Interpolation into this space is not yet supported.");
213 
214  mesh->topology_mutable().create_entity_permutations();
215  const std::vector<std::uint32_t>& cell_info
216  = mesh->topology().get_cell_permutation_info();
217 
218  // Evaluate function at physical points. The returned array has a
219  // number of rows equal to the number of components of the function,
220  // and the number of columns is equal to the number of evaluation
221  // points.
222  xt::xarray<T> values = f(x);
223 
224  if (values.dimension() == 1)
225  {
226  if (element->value_size() != 1)
227  throw std::runtime_error("Interpolation data has the wrong shape.");
228  values.reshape({element->value_size(), x.shape(1)});
229  }
230 
231  if (values.shape(0) != element->value_size())
232  throw std::runtime_error("Interpolation data has the wrong shape.");
233 
234  if (values.shape(1) != cells.size() * X.shape(0))
235  throw std::runtime_error("Interpolation data has the wrong shape.");
236 
237  // Get dofmap
238  const auto dofmap = u.function_space()->dofmap();
239  assert(dofmap);
240  const int dofmap_bs = dofmap->bs();
241 
242  // Loop over cells and compute interpolation dofs
243  const int num_scalar_dofs = element->space_dimension() / element_bs;
244  const int value_size = element->value_size() / element_bs;
245 
246  std::vector<T>& coeffs = u.x()->mutable_array();
247  std::vector<T> _coeffs(num_scalar_dofs);
248 
249  // This assumes that any element with an identity interpolation matrix is a
250  // point evaluation
251  if (element->interpolation_ident())
252  {
253  for (std::int32_t c : cells)
254  {
255  xtl::span<const std::int32_t> dofs = dofmap->cell_dofs(c);
256  for (int k = 0; k < element_bs; ++k)
257  {
258  for (int i = 0; i < num_scalar_dofs; ++i)
259  _coeffs[i] = values(k, c * num_scalar_dofs + i);
260  element->apply_inverse_transpose_dof_transformation(
261  tcb::make_span(_coeffs), cell_info[c], 1);
262  for (int i = 0; i < num_scalar_dofs; ++i)
263  {
264  const int dof = i * element_bs + k;
265  std::div_t pos = std::div(dof, dofmap_bs);
266  coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
267  }
268  }
269  }
270  }
271  else
272  {
273  // Get coordinate map
274  const fem::CoordinateElement& cmap = mesh->geometry().cmap();
275 
276  // Get geometry data
277  const graph::AdjacencyList<std::int32_t>& x_dofmap
278  = mesh->geometry().dofmap();
279  // FIXME: Add proper interface for num coordinate dofs
280  const int num_dofs_g = x_dofmap.num_links(0);
281  const xt::xtensor<double, 2>& x_g = mesh->geometry().x();
282 
283  // Create data structures for Jacobian info
284  xt::xtensor<double, 3> J = xt::empty<double>({int(X.shape(0)), gdim, tdim});
285  xt::xtensor<double, 3> K = xt::empty<double>({int(X.shape(0)), tdim, gdim});
286  xt::xtensor<double, 1> detJ = xt::empty<double>({X.shape(0)});
287 
288  xt::xtensor<double, 2> coordinate_dofs
289  = xt::empty<double>({num_dofs_g, gdim});
290 
291  xt::xtensor<T, 3> reference_data({X.shape(0), 1, value_size});
292  xt::xtensor<T, 3> _vals({X.shape(0), 1, value_size});
293 
294  // Tabulate 1st order derivatives of shape functions at interpolation coords
295  xt::xtensor<double, 4> dphi
296  = xt::view(cmap.tabulate(1, X), xt::range(1, tdim + 1), xt::all(),
297  xt::all(), xt::all());
298 
299  for (std::int32_t c : cells)
300  {
301  auto x_dofs = x_dofmap.links(c);
302  for (int i = 0; i < num_dofs_g; ++i)
303  for (int j = 0; j < gdim; ++j)
304  coordinate_dofs(i, j) = x_g(x_dofs[i], j);
305 
306  // Compute J, detJ and K
307  cmap.compute_jacobian(dphi, coordinate_dofs, J);
308  cmap.compute_jacobian_inverse(J, K);
309  cmap.compute_jacobian_determinant(J, detJ);
310 
311  xtl::span<const std::int32_t> dofs = dofmap->cell_dofs(c);
312  for (int k = 0; k < element_bs; ++k)
313  {
314  // Extract computed expression values for element block k
315  for (int m = 0; m < value_size; ++m)
316  {
317  std::copy_n(&values(k * value_size + m, c * X.shape(0)), X.shape(0),
318  xt::view(_vals, xt::all(), 0, m).begin());
319  }
320 
321  // Get element degrees of freedom for block
322  element->map_pull_back(_vals, J, detJ, K, reference_data);
323 
324  xt::xtensor<T, 2> ref_data
325  = xt::transpose(xt::view(reference_data, xt::all(), 0, xt::all()));
326  element->interpolate(ref_data, tcb::make_span(_coeffs));
327  element->apply_inverse_transpose_dof_transformation(
328  tcb::make_span(_coeffs), cell_info[c], 1);
329 
330  assert(_coeffs.size() == num_scalar_dofs);
331 
332  // Copy interpolation dofs into coefficient vector
333  for (int i = 0; i < num_scalar_dofs; ++i)
334  {
335  const int dof = i * element_bs + k;
336  std::div_t pos = std::div(dof, dofmap_bs);
337  coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
338  }
339  }
340  }
341  }
342 }
343 //----------------------------------------------------------------------------
344 template <typename T>
346  Function<T>& u,
347  const std::function<void(xt::xarray<T>&, const xt::xtensor<double, 2>&)>& f,
348  const xt::xtensor<double, 2>& x, const xtl::span<const std::int32_t>& cells)
349 {
350  const std::shared_ptr<const FiniteElement> element
351  = u.function_space()->element();
352  assert(element);
353  std::vector<int> vshape(element->value_rank(), 1);
354  for (std::size_t i = 0; i < vshape.size(); ++i)
355  vshape[i] = element->value_dimension(i);
356  const std::size_t value_size = std::accumulate(
357  std::begin(vshape), std::end(vshape), 1, std::multiplies<>());
358 
359  auto fn = [value_size, &f](const xt::xtensor<double, 2>& x) {
360  xt::xarray<T> values = xt::empty<T>({value_size, x.shape(1)});
361  f(values, x);
362  return values;
363  };
364 
365  interpolate<T>(u, fn, x, cells);
366 }
367 //----------------------------------------------------------------------------
368 
369 } // namespace dolfinx::fem
dolfinx::graph::AdjacencyList::num_links
int num_links(int node) const
Number of connections for given node.
Definition: AdjacencyList.h:120
dolfinx::fem::Function::function_space
std::shared_ptr< const FunctionSpace > function_space() const
Return shared pointer to function space.
Definition: Function.h:156
dolfinx::fem::Function
This class represents a function in a finite element function space , given by.
Definition: Form.h:23
dolfinx::fem::CoordinateElement::compute_jacobian
void compute_jacobian(const xt::xtensor< double, 4 > &dphi, const xt::xtensor< double, 2 > &cell_geometry, xt::xtensor< double, 3 > &J) const
Compute Jacobian for a cell with given geometry using the basis functions and first order derivatives...
Definition: CoordinateElement.cpp:81
dolfinx::fem::CoordinateElement::compute_jacobian_determinant
void compute_jacobian_determinant(const xt::xtensor< double, 3 > &J, xt::xtensor< double, 1 > &detJ) const
Compute the determinant of the Jacobian. If the coordinate element is affine, it computes the determi...
Definition: CoordinateElement.cpp:159
dolfinx::fem::CoordinateElement::compute_jacobian_inverse
void compute_jacobian_inverse(const xt::xtensor< double, 3 > &J, xt::xtensor< double, 3 > &K) const
Compute the inverse of the Jacobian. If the coordinate element is affine, it computes the inverse at ...
Definition: CoordinateElement.cpp:122
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::CoordinateElement::tabulate
xt::xtensor< double, 4 > tabulate(int n, const xt::xtensor< double, 2 > &X) const
Compute basis values and derivatives at set of points.
Definition: CoordinateElement.cpp:76
dolfinx::fem::interpolate_c
void interpolate_c(Function< T > &u, const std::function< void(xt::xarray< T > &, const xt::xtensor< double, 2 > &)> &f, const xt::xtensor< double, 2 > &x, const xtl::span< const std::int32_t > &cells)
Interpolate an expression f(x)
Definition: interpolate.h:345
dolfinx::fem::interpolation_coords
xt::xtensor< double, 2 > interpolation_coords(const fem::FiniteElement &element, const mesh::Mesh &mesh, const xtl::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:18
dolfinx::fem
Finite element method functionality.
Definition: assemble_matrix_impl.h:22
dolfinx::fem::interpolate
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:145
dolfinx::graph::AdjacencyList::links
xtl::span< T > links(int node)
Get the links (edges) for given node.
Definition: AdjacencyList.h:130
dolfinx::fem::Function::x
std::shared_ptr< const la::Vector< T > > x() const
Underlying vector.
Definition: Function.h:197
dolfinx::fem::CoordinateElement
This class manages coordinate mappings for isoparametric cells.
Definition: CoordinateElement.h:30