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