DOLFINx 0.11.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
interpolate.h
1// Copyright (C) 2020-2026 Garth N. Wells, Igor A. Baratta, Massimiliano Leoni
2// and Jørgen S.Dokken
3//
4// This file is part of DOLFINx (https://www.fenicsproject.org)
5//
6// SPDX-License-Identifier: LGPL-3.0-or-later
7
8#pragma once
9
10#include "CoordinateElement.h"
11#include "DofMap.h"
12#include "FiniteElement.h"
13#include "FunctionSpace.h"
14#include <algorithm>
15#include <basix/mdspan.hpp>
16#include <concepts>
17#include <dolfinx/common/IndexMap.h>
18#include <dolfinx/common/types.h>
19#include <dolfinx/geometry/utils.h>
20#include <dolfinx/mesh/Mesh.h>
21#include <functional>
22#include <numeric>
23#include <ranges>
24#include <span>
25#include <vector>
26
27namespace dolfinx::fem
28{
29template <dolfinx::scalar T, std::floating_point U>
30class Function;
31
32template <typename T>
33concept MDSpan = requires(T x, std::size_t idx) {
34 x(idx, idx);
35 { x.extent(0) } -> std::integral;
36 { x.extent(1) } -> std::integral;
37};
38
40template <typename R>
41concept CellRange = std::ranges::input_range<R> and std::ranges::sized_range<R>
42 and std::is_integral_v<
43 std::remove_const_t<std::ranges::range_value_t<R>>>;
44
55template <std::floating_point T>
56std::vector<T> interpolation_coords(const fem::FiniteElement<T>& element,
58 CellRange auto&& cells)
59{
60 // Find CoordinateElement appropriate to element
61 const std::vector<CoordinateElement<T>>& cmaps = geometry.cmaps();
62 mesh::CellType cell_type = element.cell_type();
63 auto it
64 = std::find_if(cmaps.begin(), cmaps.end(), [&cell_type](const auto& cm)
65 { return cell_type == cm.cell_shape(); });
66 if (it == cmaps.end())
67 throw std::runtime_error("Cannot find CoordinateElement for FiniteElement");
68 int index = std::distance(cmaps.begin(), it);
69
70 // Get geometry data and the element coordinate map
71 const std::size_t gdim = geometry.dim();
72 auto x_dofmap = geometry.dofmap(index);
73 std::span<const T> x_g = geometry.x();
74
75 const CoordinateElement<T>& cmap = cmaps.at(index);
76 const std::size_t num_dofs_g = cmap.dim();
77
78 // Get the interpolation points on the reference cells
79 const auto [X, Xshape] = element.interpolation_points();
80
81 // Evaluate coordinate element basis at reference points
82 std::array<std::size_t, 4> phi_shape = cmap.tabulate_shape(0, Xshape[0]);
83 std::vector<T> phi_b(
84 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
85 md::mdspan<const T, md::extents<std::size_t, 1, md::dynamic_extent,
86 md::dynamic_extent, 1>>
87 phi_full(phi_b.data(), phi_shape);
88 cmap.tabulate(0, X, Xshape, phi_b);
89 auto phi = md::submdspan(phi_full, 0, md::full_extent, md::full_extent, 0);
90
91 // Push reference coordinates (X) forward to the physical coordinates
92 // (x) for each cell
93 std::vector<T> coordinate_dofs(num_dofs_g * gdim, 0);
94 std::vector<T> x(3 * (cells.size() * Xshape[0]), 0);
95 for (auto cell_it = cells.begin(); cell_it != cells.end(); ++cell_it)
96 {
97 // Get geometry data for current cell
98 auto x_dofs = md::submdspan(x_dofmap, *cell_it, md::full_extent);
99 for (std::size_t i = 0; i < x_dofs.size(); ++i)
100 {
101 std::copy_n(std::next(x_g.begin(), 3 * x_dofs[i]), gdim,
102 std::next(coordinate_dofs.begin(), i * gdim));
103 }
104
105 // Push forward coordinates (X -> x)
106 std::size_t offset = std::distance(cells.begin(), cell_it);
107 for (std::size_t p = 0; p < Xshape[0]; ++p)
108 {
109 for (std::size_t j = 0; j < gdim; ++j)
110 {
111 T acc = 0;
112 for (std::size_t k = 0; k < num_dofs_g; ++k)
113 acc += phi(p, k) * coordinate_dofs[k * gdim + j];
114 x[j * (cells.size() * Xshape[0]) + offset * Xshape[0] + p] = acc;
115 }
116 }
117 }
118
119 return x;
120}
121
138template <dolfinx::scalar T, std::floating_point U>
139void interpolate(Function<T, U>& u, std::span<const T> f,
140 std::array<std::size_t, 2> fshape, CellRange auto&& cells);
141
142namespace impl
143{
145template <typename T, std::size_t D>
146using mdspan_t = md::mdspan<T, md::dextents<std::size_t, D>>;
147
167template <dolfinx::scalar T>
168void scatter_values(MPI_Comm comm, std::span<const std::int32_t> src_ranks,
169 std::span<const std::int32_t> dest_ranks,
170 mdspan_t<const T, 2> send_values, std::span<T> recv_values)
171{
172 const std::size_t block_size = send_values.extent(1);
173 assert(src_ranks.size() * block_size == send_values.size());
174 assert(recv_values.size() == dest_ranks.size() * block_size);
175
176 // Build unique set of the sorted src_ranks
177 std::vector<std::int32_t> out_ranks(src_ranks.size());
178 out_ranks.assign(src_ranks.begin(), src_ranks.end());
179 auto [unique_end, range_end] = std::ranges::unique(out_ranks);
180 out_ranks.erase(unique_end, range_end);
181 out_ranks.reserve(out_ranks.size() + 1);
182
183 // Remove negative entries from dest_ranks
184 std::vector<std::int32_t> in_ranks;
185 in_ranks.reserve(dest_ranks.size());
186 std::copy_if(dest_ranks.begin(), dest_ranks.end(),
187 std::back_inserter(in_ranks),
188 [](auto rank) { return rank >= 0; });
189
190 // Create unique set of sorted in-ranks
191 {
192 std::ranges::sort(in_ranks);
193 auto [unique_end, range_end] = std::ranges::unique(in_ranks);
194 in_ranks.erase(unique_end, range_end);
195 }
196 in_ranks.reserve(in_ranks.size() + 1);
197
198 // Create neighborhood communicator
199 MPI_Comm reverse_comm;
200 MPI_Dist_graph_create_adjacent(
201 comm, in_ranks.size(), in_ranks.data(), MPI_UNWEIGHTED, out_ranks.size(),
202 out_ranks.data(), MPI_UNWEIGHTED, MPI_INFO_NULL, false, &reverse_comm);
203
204 std::vector<std::int32_t> comm_to_output;
205 std::vector<std::int32_t> recv_sizes(in_ranks.size());
206 recv_sizes.reserve(1);
207 std::vector<std::int32_t> recv_offsets(in_ranks.size() + 1, 0);
208 {
209 // Build map from parent to neighborhood communicator ranks
210 std::vector<std::pair<std::int32_t, std::int32_t>> rank_to_neighbor;
211 rank_to_neighbor.reserve(in_ranks.size());
212 for (std::size_t i = 0; i < in_ranks.size(); i++)
213 rank_to_neighbor.push_back({in_ranks[i], i});
214 std::ranges::sort(rank_to_neighbor);
215
216 // Compute receive sizes
217 std::ranges::for_each(
218 dest_ranks,
219 [&rank_to_neighbor, &recv_sizes, block_size](auto rank)
220 {
221 if (rank >= 0)
222 {
223 auto it = std::ranges::lower_bound(rank_to_neighbor, rank,
224 std::ranges::less(),
225 [](auto e) { return e.first; });
226 assert(it != rank_to_neighbor.end() and it->first == rank);
227 recv_sizes[it->second] += block_size;
228 }
229 });
230
231 // Compute receiving offsets
232 std::partial_sum(recv_sizes.begin(), recv_sizes.end(),
233 std::next(recv_offsets.begin(), 1));
234
235 // Compute map from receiving values to position in recv_values
236 comm_to_output.resize(recv_offsets.back() / block_size);
237 std::vector<std::int32_t> recv_counter(recv_sizes.size(), 0);
238 for (std::size_t i = 0; i < dest_ranks.size(); ++i)
239 {
240 if (const std::int32_t rank = dest_ranks[i]; rank >= 0)
241 {
242 auto it = std::ranges::lower_bound(rank_to_neighbor, rank,
243 std::ranges::less(),
244 [](auto e) { return e.first; });
245 assert(it != rank_to_neighbor.end() and it->first == rank);
246 int insert_pos = recv_offsets[it->second] + recv_counter[it->second];
247 comm_to_output[insert_pos / block_size] = i * block_size;
248 recv_counter[it->second] += block_size;
249 }
250 }
251 }
252
253 std::vector<std::int32_t> send_sizes(out_ranks.size());
254 send_sizes.reserve(1);
255 {
256 // Compute map from parent MPI rank to neighbor rank for outgoing
257 // data. `out_ranks` is sorted, so rank_to_neighbor will be sorted
258 // too.
259 std::vector<std::pair<std::int32_t, std::int32_t>> rank_to_neighbor;
260 rank_to_neighbor.reserve(out_ranks.size());
261 for (std::size_t i = 0; i < out_ranks.size(); i++)
262 rank_to_neighbor.push_back({out_ranks[i], i});
263
264 // Compute send sizes. As `src_ranks` is sorted, we can move 'start'
265 // in search forward.
266 auto start = rank_to_neighbor.begin();
267 std::ranges::for_each(
268 src_ranks,
269 [&rank_to_neighbor, &send_sizes, block_size, &start](auto rank)
270 {
271 auto it = std::ranges::lower_bound(start, rank_to_neighbor.end(),
272 rank, std::ranges::less(),
273 [](auto e) { return e.first; });
274 assert(it != rank_to_neighbor.end() and it->first == rank);
275 send_sizes[it->second] += block_size;
276 start = it;
277 });
278 }
279
280 // Compute sending offsets
281 std::vector<std::int32_t> send_offsets(send_sizes.size() + 1, 0);
282 std::partial_sum(send_sizes.begin(), send_sizes.end(),
283 std::next(send_offsets.begin(), 1));
284
285 // Send values to dest ranks
286 std::vector<T> values(recv_offsets.back());
287 values.reserve(1);
288 MPI_Neighbor_alltoallv(send_values.data_handle(), send_sizes.data(),
289 send_offsets.data(), dolfinx::MPI::mpi_t<T>,
290 values.data(), recv_sizes.data(), recv_offsets.data(),
291 dolfinx::MPI::mpi_t<T>, reverse_comm);
292 MPI_Comm_free(&reverse_comm);
293
294 // Insert values received from neighborhood communicator in output
295 // span
296 std::ranges::fill(recv_values, T(0));
297 for (std::size_t i = 0; i < comm_to_output.size(); i++)
298 {
299 auto vals = std::next(recv_values.begin(), comm_to_output[i]);
300 auto vals_from = std::next(values.begin(), i * block_size);
301 std::copy_n(vals_from, block_size, vals);
302 }
303};
304
313template <MDSpan U, MDSpan V, dolfinx::scalar T>
314void interpolation_apply(U&& Pi, V&& data, std::span<T> coeffs, int bs)
315{
316 using X = typename dolfinx::scalar_value_t<T>;
317
318 // Compute coefficients = Pi * x (matrix-vector multiply)
319 if (bs == 1)
320 {
321 assert(data.extent(0) * data.extent(1) == Pi.extent(1));
322 for (std::size_t i = 0; i < Pi.extent(0); ++i)
323 {
324 coeffs[i] = 0.0;
325 for (std::size_t k = 0; k < data.extent(1); ++k)
326 for (std::size_t j = 0; j < data.extent(0); ++j)
327 coeffs[i]
328 += static_cast<X>(Pi(i, k * data.extent(0) + j)) * data(j, k);
329 }
330 }
331 else
332 {
333 assert(data.extent(0) == Pi.extent(1));
334 assert(static_cast<int>(data.extent(1)) == bs);
335 std::size_t cols = Pi.extent(1);
336 for (int k = 0; k < bs; ++k)
337 {
338 for (std::size_t i = 0; i < Pi.extent(0); ++i)
339 {
340 T acc = 0;
341 for (std::size_t j = 0; j < cols; ++j)
342 acc += static_cast<X>(Pi(i, j)) * data(j, k);
343 coeffs[bs * i + k] = acc;
344 }
345 }
346 }
347}
348
368template <dolfinx::scalar T, std::floating_point U>
369void interpolate_same_map(Function<T, U>& u1, CellRange auto&& cells1,
370 const Function<T, U>& u0, CellRange auto&& cells0)
371{
372 auto V0 = u0.function_space();
373 assert(V0);
374 auto V1 = u1.function_space();
375 assert(V1);
376 auto mesh0 = V0->mesh();
377 assert(mesh0);
378
379 auto mesh1 = V1->mesh();
380 assert(mesh1);
381
382 auto element0 = V0->element();
383 assert(element0);
384 auto element1 = V1->element();
385 assert(element1);
386
387 assert(mesh0->topology()->dim());
388 const int tdim = mesh0->topology()->dim();
389 auto map = mesh0->topology()->index_map(tdim);
390 assert(map);
391 std::span<T> u1_array = u1.x()->array();
392 std::span<const T> u0_array = u0.x()->array();
393
394 std::span<const std::uint32_t> cell_info0;
395 std::span<const std::uint32_t> cell_info1;
396 if (element1->needs_dof_transformations()
397 or element0->needs_dof_transformations())
398 {
399 mesh0->topology_mutable()->create_entity_permutations();
400 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
401 mesh1->topology_mutable()->create_entity_permutations();
402 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
403 }
404
405 // Get dofmaps
406 auto dofmap1 = V1->dofmap();
407 auto dofmap0 = V0->dofmap();
408
409 // Get block sizes and dof transformation operators
410 const int bs1 = dofmap1->bs();
411 const int bs0 = dofmap0->bs();
412 auto apply_dof_transformation = element0->template dof_transformation_fn<T>(
414 auto apply_inverse_dof_transform
415 = element1->template dof_transformation_fn<T>(
417
418 // Create working array
419 std::vector<T> local0(element0->space_dimension());
420 std::vector<T> local1(element1->space_dimension());
421
422 // Create interpolation operator
423 auto [i_m, im_shape] = element1->create_interpolation_operator(*element0);
424
425 // Iterate over mesh and interpolate on each cell
426 using X = typename dolfinx::scalar_value_t<T>;
427 assert(cells0.size() == cells1.size());
428 for (auto cell0_it = cells0.begin(), cell1_it = cells1.begin();
429 cell0_it != cells0.end() and cell1_it != cells1.end();
430 ++cell0_it, ++cell1_it)
431 {
432 // Pack and transform cell dofs to reference ordering
433 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(*cell0_it);
434 for (std::size_t i = 0; i < dofs0.size(); ++i)
435 for (int k = 0; k < bs0; ++k)
436 local0[bs0 * i + k] = u0_array[bs0 * dofs0[i] + k];
437
438 apply_dof_transformation(local0, cell_info0, *cell0_it, 1);
439
440 // FIXME: Get compile-time ranges from Basix
441 // Apply interpolation operator
442 std::ranges::fill(local1, 0);
443 for (std::size_t i = 0; i < im_shape[0]; ++i)
444 for (std::size_t j = 0; j < im_shape[1]; ++j)
445 local1[i] += static_cast<X>(i_m[im_shape[1] * i + j]) * local0[j];
446
447 apply_inverse_dof_transform(local1, cell_info1, *cell1_it, 1);
448 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(*cell1_it);
449 for (std::size_t i = 0; i < dofs1.size(); ++i)
450 for (int k = 0; k < bs1; ++k)
451 u1_array[bs1 * dofs1[i] + k] = local1[bs1 * i + k];
452 }
453}
454
469template <dolfinx::scalar T, std::floating_point U>
470void interpolate_nonmatching_maps(Function<T, U>& u1, CellRange auto&& cells1,
471 const Function<T, U>& u0,
472 CellRange auto&& cells0)
473{
474 // Get mesh
475 auto V0 = u0.function_space();
476 assert(V0);
477 auto mesh0 = V0->mesh();
478 assert(mesh0);
479
480 // Mesh dims
481 const int tdim = mesh0->topology()->dim();
482 const int gdim = mesh0->geometry().dim();
483
484 // Get elements
485 auto V1 = u1.function_space();
486 assert(V1);
487 auto mesh1 = V1->mesh();
488 assert(mesh1);
489 auto element0 = V0->element();
490 assert(element0);
491 auto element1 = V1->element();
492 assert(element1);
493
494 std::span<const std::uint32_t> cell_info0;
495 std::span<const std::uint32_t> cell_info1;
496 if (element1->needs_dof_transformations()
497 or element0->needs_dof_transformations())
498 {
499 mesh0->topology_mutable()->create_entity_permutations();
500 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
501 mesh1->topology_mutable()->create_entity_permutations();
502 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
503 }
504
505 // Get dofmaps
506 auto dofmap0 = V0->dofmap();
507 auto dofmap1 = V1->dofmap();
508
509 const auto [X, Xshape] = element1->interpolation_points();
510
511 // Get block sizes and dof transformation operators
512 const int bs0 = element0->block_size();
513 const int bs1 = element1->block_size();
514 auto apply_dof_transformation0 = element0->template dof_transformation_fn<U>(
516 auto apply_inv_dof_transform1 = element1->template dof_transformation_fn<T>(
518
519 // Get sizes of elements
520 const std::size_t dim0 = element0->space_dimension() / bs0;
521 const std::size_t value_size_ref0 = element0->reference_value_size();
522 const std::size_t value_size0 = V0->element()->reference_value_size();
523
524 const CoordinateElement<U>& cmap = mesh0->geometry().cmap();
525 auto x_dofmap = mesh0->geometry().dofmap();
526 std::span<const U> x_g = mesh0->geometry().x();
527
528 // (0) is derivative index, (1) is the point index, (2) is the basis
529 // function index and (3) is the basis function component.
530
531 // Evaluate coordinate map basis at reference interpolation points
532 const std::array<std::size_t, 4> phi_shape
533 = cmap.tabulate_shape(1, Xshape[0]);
534 std::vector<U> phi_b(
535 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
536 md::mdspan<const U, md::extents<std::size_t, md::dynamic_extent,
537 md::dynamic_extent, md::dynamic_extent, 1>>
538 phi(phi_b.data(), phi_shape);
539 cmap.tabulate(1, X, Xshape, phi_b);
540
541 // Evaluate v basis functions at reference interpolation points
542 const auto [_basis_derivatives_reference0, b0shape]
543 = element0->tabulate(X, Xshape, 0);
544 md::mdspan<const U, std::extents<std::size_t, 1, md::dynamic_extent,
545 md::dynamic_extent, md::dynamic_extent>>
546 basis_derivatives_reference0(_basis_derivatives_reference0.data(),
547 b0shape);
548
549 // Create working arrays
550 std::vector<T> local1(element1->space_dimension());
551 std::vector<T> coeffs0(element0->space_dimension());
552
553 std::vector<U> basis0_b(Xshape[0] * dim0 * value_size0);
554 md::mdspan<U, std::dextents<std::size_t, 3>> basis0(
555 basis0_b.data(), Xshape[0], dim0, value_size0);
556
557 std::vector<U> basis_reference0_b(Xshape[0] * dim0 * value_size_ref0);
558 md::mdspan<U, std::dextents<std::size_t, 3>> basis_reference0(
559 basis_reference0_b.data(), Xshape[0], dim0, value_size_ref0);
560
561 std::vector<T> values0_b(Xshape[0] * 1 * V1->element()->value_size());
562 md::mdspan<
563 T, md::extents<std::size_t, md::dynamic_extent, 1, md::dynamic_extent>>
564 values0(values0_b.data(), Xshape[0], 1, V1->element()->value_size());
565
566 std::vector<T> mapped_values_b(Xshape[0] * 1 * V1->element()->value_size());
567 md::mdspan<
568 T, md::extents<std::size_t, md::dynamic_extent, 1, md::dynamic_extent>>
569 mapped_values0(mapped_values_b.data(), Xshape[0], 1,
570 V1->element()->value_size());
571
572 const std::size_t num_dofs_g = cmap.dim();
573 std::vector<U> coord_dofs_b(num_dofs_g * gdim);
574 md::mdspan<U, std::dextents<std::size_t, 2>> coord_dofs(coord_dofs_b.data(),
575 num_dofs_g, gdim);
576
577 std::vector<U> J_b(Xshape[0] * gdim * tdim);
578 md::mdspan<U, std::dextents<std::size_t, 3>> J(J_b.data(), Xshape[0], gdim,
579 tdim);
580 std::vector<U> K_b(Xshape[0] * tdim * gdim);
581 md::mdspan<U, std::dextents<std::size_t, 3>> K(K_b.data(), Xshape[0], tdim,
582 gdim);
583 std::vector<U> detJ(Xshape[0]);
584 std::vector<U> det_scratch(2 * gdim * tdim);
585
586 // Get interpolation operator
587 const auto [_Pi_1, pi_shape] = element1->interpolation_operator();
588 impl::mdspan_t<const U, 2> Pi_1(_Pi_1.data(), pi_shape);
589
590 using u_t = md::mdspan<U, std::dextents<std::size_t, 2>>;
591 using U_t = md::mdspan<const U, std::dextents<std::size_t, 2>>;
592 using J_t = md::mdspan<const U, std::dextents<std::size_t, 2>>;
593 using K_t = md::mdspan<const U, std::dextents<std::size_t, 2>>;
594 auto push_forward_fn0
595 = element0->basix_element().template map_fn<u_t, U_t, J_t, K_t>();
596
597 using v_t = md::mdspan<const T, std::dextents<std::size_t, 2>>;
598 using V_t = decltype(md::submdspan(mapped_values0, 0, md::full_extent,
599 md::full_extent));
600 auto pull_back_fn1
601 = element1->basix_element().template map_fn<V_t, v_t, K_t, J_t>();
602
603 // Iterate over mesh and interpolate on each cell
604 std::span<const T> array0 = u0.x()->array();
605 std::span<T> array1 = u1.x()->array();
606 assert(cells0.size() == cells1.size());
607 for (auto cell0_it = cells0.begin(), cell1_it = cells1.begin();
608 cell0_it != cells0.end() and cell1_it != cells0.end();
609 ++cell0_it, ++cell1_it)
610 {
611 // Get cell geometry (coordinate dofs)
612 auto x_dofs = md::submdspan(x_dofmap, *cell0_it, md::full_extent);
613 for (std::size_t i = 0; i < num_dofs_g; ++i)
614 {
615 const int pos = 3 * x_dofs[i];
616 for (int j = 0; j < gdim; ++j)
617 coord_dofs(i, j) = x_g[pos + j];
618 }
619
620 // Compute Jacobians and reference points for current cell
621 std::ranges::fill(J_b, 0);
622 for (std::size_t p = 0; p < Xshape[0]; ++p)
623 {
624 auto dphi
625 = md::submdspan(phi, std::pair(1, tdim + 1), p, md::full_extent, 0);
626 auto _J = md::submdspan(J, p, md::full_extent, md::full_extent);
627 cmap.compute_jacobian(dphi, coord_dofs, _J);
628 auto _K = md::submdspan(K, p, md::full_extent, md::full_extent);
629 cmap.compute_jacobian_inverse(_J, _K);
630 detJ[p] = cmap.compute_jacobian_determinant(_J, det_scratch);
631 }
632
633 // Copy evaluated basis on reference, apply DOF transformations, and
634 // push forward to physical element
635 for (std::size_t k0 = 0; k0 < basis_reference0.extent(0); ++k0)
636 for (std::size_t k1 = 0; k1 < basis_reference0.extent(1); ++k1)
637 for (std::size_t k2 = 0; k2 < basis_reference0.extent(2); ++k2)
638 basis_reference0(k0, k1, k2)
639 = basis_derivatives_reference0(0, k0, k1, k2);
640
641 for (std::size_t p = 0; p < Xshape[0]; ++p)
642 {
643 apply_dof_transformation0(
644 std::span(basis_reference0_b.data() + p * dim0 * value_size_ref0,
645 dim0 * value_size_ref0),
646 cell_info0, *cell0_it, value_size_ref0);
647 }
648
649 for (std::size_t i = 0; i < basis0.extent(0); ++i)
650 {
651 auto _u = md::submdspan(basis0, i, md::full_extent, md::full_extent);
652 auto _U = md::submdspan(basis_reference0, i, md::full_extent,
653 md::full_extent);
654 auto _K = md::submdspan(K, i, md::full_extent, md::full_extent);
655 auto _J = md::submdspan(J, i, md::full_extent, md::full_extent);
656 push_forward_fn0(_u, _U, _J, detJ[i], _K);
657 }
658
659 // Copy expansion coefficients for v into local array
660 const int dof_bs0 = dofmap0->bs();
661 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(*cell0_it);
662 for (std::size_t i = 0; i < dofs0.size(); ++i)
663 for (int k = 0; k < dof_bs0; ++k)
664 coeffs0[dof_bs0 * i + k] = array0[dof_bs0 * dofs0[i] + k];
665
666 // Evaluate v at the interpolation points (physical space values)
667 using X = typename dolfinx::scalar_value_t<T>;
668 for (std::size_t p = 0; p < Xshape[0]; ++p)
669 {
670 for (int k = 0; k < bs0; ++k)
671 {
672 for (std::size_t j = 0; j < value_size0; ++j)
673 {
674 T acc = 0;
675 for (std::size_t i = 0; i < dim0; ++i)
676 acc += coeffs0[bs0 * i + k] * static_cast<X>(basis0(p, i, j));
677 values0(p, 0, j * bs0 + k) = acc;
678 }
679 }
680 }
681
682 // Pull back the physical values to the u reference
683 for (std::size_t i = 0; i < values0.extent(0); ++i)
684 {
685 auto _u = md::submdspan(values0, i, md::full_extent, md::full_extent);
686 auto _U
687 = md::submdspan(mapped_values0, i, md::full_extent, md::full_extent);
688 auto _K = md::submdspan(K, i, md::full_extent, md::full_extent);
689 auto _J = md::submdspan(J, i, md::full_extent, md::full_extent);
690 pull_back_fn1(_U, _u, _K, 1.0 / detJ[i], _J);
691 }
692
693 auto values
694 = md::submdspan(mapped_values0, md::full_extent, 0, md::full_extent);
695 interpolation_apply(Pi_1, values, std::span(local1), bs1);
696 apply_inv_dof_transform1(local1, cell_info1, *cell1_it, 1);
697
698 // Copy local coefficients to the correct position in u dof array
699 const int dof_bs1 = dofmap1->bs();
700 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(*cell1_it);
701 for (std::size_t i = 0; i < dofs1.size(); ++i)
702 for (int k = 0; k < dof_bs1; ++k)
703 array1[dof_bs1 * dofs1[i] + k] = local1[dof_bs1 * i + k];
704 }
705}
706
718template <dolfinx::scalar T, std::floating_point U>
719void point_evaluation(const FiniteElement<U>& element, bool symmetric,
720 const DofMap& dofmap, CellRange auto&& cells,
721 std::span<const std::uint32_t> cell_info,
722 std::span<const T> f, std::array<std::size_t, 2> fshape,
723 std::span<T> coeffs)
724{
725 // Point evaluation element *and* the geometric map is the identity,
726 // e.g. not Piola mapped
727
728 const int element_bs = element.block_size();
729 const int num_scalar_dofs = element.space_dimension() / element_bs;
730 const int dofmap_bs = dofmap.bs();
731
732 auto apply_inv_transpose_dof_transformation
733 = element.template dof_transformation_fn<T>(
735 std::vector<T> coeffs_b(num_scalar_dofs);
736 if (symmetric)
737 {
738 std::size_t matrix_size = 0;
739 while (matrix_size * matrix_size < fshape[0])
740 ++matrix_size;
741
742 // Loop over cells
743 for (auto cell_it = cells.begin(); cell_it != cells.end(); ++cell_it)
744 {
745 // The entries of a symmetric matrix are numbered (for an
746 // example 4x4 element):
747 // 0 * * *
748 // 1 2 * *
749 // 3 4 5 *
750 // 6 7 8 9
751 // The loop extracts these elements. In this loop, row is the
752 // row of this matrix, and (k - rowstart) is the column
753 std::size_t row = 0;
754 std::size_t rowstart = 0;
755 std::span<const std::int32_t> dofs = dofmap.cell_dofs(*cell_it);
756 std::size_t offset = std::distance(cells.begin(), cell_it);
757 for (int k = 0; k < element_bs; ++k)
758 {
759 if (k - rowstart > row)
760 {
761 ++row;
762 rowstart = k;
763 }
764
765 // num_scalar_dofs is the number of interpolation points per
766 // cell in this case (interpolation matrix is identity)
767 std::copy_n(
768 std::next(f.begin(), (row * matrix_size + k - rowstart) * fshape[1]
769 + offset * num_scalar_dofs),
770 num_scalar_dofs, coeffs_b.data());
771 apply_inv_transpose_dof_transformation(coeffs_b, cell_info, *cell_it,
772 1);
773 for (int i = 0; i < num_scalar_dofs; ++i)
774 {
775 const int dof = i * element_bs + k;
776 std::div_t pos = std::div(dof, dofmap_bs);
777 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = coeffs_b[i];
778 }
779 }
780 }
781 }
782 else
783 {
784 // Loop over cells
785 for (auto cell_it = cells.begin(); cell_it != cells.end(); ++cell_it)
786 {
787 std::size_t offset = std::distance(cells.begin(), cell_it);
788 std::span<const std::int32_t> dofs = dofmap.cell_dofs(*cell_it);
789 for (int k = 0; k < element_bs; ++k)
790 {
791 // num_scalar_dofs is the number of interpolation points per
792 // cell in this case (interpolation matrix is identity)
793 std::copy_n(
794 std::next(f.begin(), k * fshape[1] + offset * num_scalar_dofs),
795 num_scalar_dofs, coeffs_b.data());
796 apply_inv_transpose_dof_transformation(coeffs_b, cell_info, *cell_it,
797 1);
798 for (int i = 0; i < num_scalar_dofs; ++i)
799 {
800 const int dof = i * element_bs + k;
801 std::div_t pos = std::div(dof, dofmap_bs);
802 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = coeffs_b[i];
803 }
804 }
805 }
806 }
807}
808
820template <dolfinx::scalar T, std::floating_point U>
821void identity_mapped_evaluation(const FiniteElement<U>& element, bool symmetric,
822 const DofMap& dofmap, CellRange auto&& cells,
823 std::span<const std::uint32_t> cell_info,
824 std::span<const T> f,
825 std::array<std::size_t, 2> fshape,
826 std::span<T> coeffs)
827{
828 // Not a point evaluation, but the geometric map is the identity,
829 // e.g. not Piola mapped
830
831 if (symmetric)
832 throw std::runtime_error("Interpolation into this element not supported.");
833
834 const int element_bs = element.block_size();
835 const int num_scalar_dofs = element.space_dimension() / element_bs;
836 const int dofmap_bs = dofmap.bs();
837
838 const int element_vs = element.reference_value_size();
839 if (element_vs > 1 and element_bs > 1)
840 throw std::runtime_error("Interpolation into this element not supported.");
841
842 // Get interpolation operator
843 const auto [_Pi, pi_shape] = element.interpolation_operator();
844 md::mdspan<const U, std::dextents<std::size_t, 2>> Pi(_Pi.data(), pi_shape);
845 const std::size_t num_interp_points = Pi.extent(1);
846 assert(static_cast<int>(Pi.extent(0)) == num_scalar_dofs);
847
848 auto apply_inv_transpose_dof_transformation
849 = element.template dof_transformation_fn<T>(
851
852 // Loop over cells
853 std::vector<T> ref_data_b(num_interp_points);
854 md::mdspan<T, md::extents<std::size_t, md::dynamic_extent, 1>> ref_data(
855 ref_data_b.data(), num_interp_points, 1);
856 std::vector<T> coeffs_b(num_scalar_dofs);
857 for (auto cell_it = cells.begin(); cell_it != cells.end(); ++cell_it)
858 {
859 std::size_t offset = std::distance(cells.begin(), cell_it);
860 std::span<const std::int32_t> dofs = dofmap.cell_dofs(*cell_it);
861 for (int k = 0; k < element_bs; ++k)
862 {
863 for (int i = 0; i < element_vs; ++i)
864 {
865 std::copy_n(
866 std::next(f.begin(), (i + k) * fshape[1]
867 + offset * num_interp_points / element_vs),
868 num_interp_points / element_vs,
869 std::next(ref_data_b.begin(), i * num_interp_points / element_vs));
870 }
871
872 impl::interpolation_apply(Pi, ref_data, std::span(coeffs_b), 1);
873 apply_inv_transpose_dof_transformation(coeffs_b, cell_info, *cell_it, 1);
874 for (int i = 0; i < num_scalar_dofs; ++i)
875 {
876 const int dof = i * element_bs + k;
877 std::div_t pos = std::div(dof, dofmap_bs);
878 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = coeffs_b[i];
879 }
880 }
881 }
882}
883
896template <dolfinx::scalar T, std::floating_point U>
897void piola_mapped_evaluation(const FiniteElement<U>& element, bool symmetric,
898 const DofMap& dofmap, CellRange auto&& cells,
899 std::span<const std::uint32_t> cell_info,
900 std::span<const T> f,
901 std::array<std::size_t, 2> fshape,
902 const mesh::Mesh<U>& mesh, std::span<T> coeffs)
903{
904 if (symmetric)
905 throw std::runtime_error("Interpolation into this element not supported.");
906
907 const int gdim = mesh.geometry().dim();
908 assert(mesh.topology());
909 const int tdim = mesh.topology()->dim();
910
911 const int element_bs = element.block_size();
912 const int num_scalar_dofs = element.space_dimension() / element_bs;
913 const int value_size = element.reference_value_size();
914 const int dofmap_bs = dofmap.bs();
915
916 md::mdspan<const T, md::dextents<std::size_t, 2>> _f(f.data(), fshape);
917
918 // Get the interpolation points on the reference cells
919 const auto [X, Xshape] = element.interpolation_points();
920 if (X.empty())
921 {
922 throw std::runtime_error(
923 "Interpolation into this space is not yet supported.");
924 }
925
926 if (_f.extent(1) != cells.size() * Xshape[0])
927 throw std::runtime_error("Interpolation data has the wrong shape.");
928
929 // Get coordinate map
930 const CoordinateElement<U>& cmap = mesh.geometry().cmap();
931
932 // Get geometry data
933 auto x_dofmap = mesh.geometry().dofmap();
934 const int num_dofs_g = cmap.dim();
935 std::span<const U> x_g = mesh.geometry().x();
936
937 // Create data structures for Jacobian info
938 std::vector<U> J_b(Xshape[0] * gdim * tdim);
939 md::mdspan<U, std::dextents<std::size_t, 3>> J(J_b.data(), Xshape[0], gdim,
940 tdim);
941 std::vector<U> K_b(Xshape[0] * tdim * gdim);
942 md::mdspan<U, std::dextents<std::size_t, 3>> K(K_b.data(), Xshape[0], tdim,
943 gdim);
944 std::vector<U> detJ(Xshape[0]);
945 std::vector<U> det_scratch(2 * gdim * tdim);
946
947 std::vector<U> coord_dofs_b(num_dofs_g * gdim);
948 md::mdspan<U, std::dextents<std::size_t, 2>> coord_dofs(coord_dofs_b.data(),
949 num_dofs_g, gdim);
950 const std::size_t value_size_ref = element.reference_value_size();
951 std::vector<T> ref_data_b(Xshape[0] * 1 * value_size_ref);
952 md::mdspan<
953 T, md::extents<std::size_t, md::dynamic_extent, 1, md::dynamic_extent>>
954 ref_data(ref_data_b.data(), Xshape[0], 1, value_size_ref);
955
956 std::vector<T> _vals_b(Xshape[0] * 1 * value_size);
957 md::mdspan<
958 T, md::extents<std::size_t, md::dynamic_extent, 1, md::dynamic_extent>>
959 _vals(_vals_b.data(), Xshape[0], 1, value_size);
960
961 // Tabulate 1st derivative of shape functions at interpolation
962 // coords
963 std::array<std::size_t, 4> phi_shape = cmap.tabulate_shape(1, Xshape[0]);
964 std::vector<U> phi_b(
965 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
966 md::mdspan<const U, md::extents<std::size_t, md::dynamic_extent,
967 md::dynamic_extent, md::dynamic_extent, 1>>
968 phi(phi_b.data(), phi_shape);
969 cmap.tabulate(1, X, Xshape, phi_b);
970 auto dphi = md::submdspan(phi, std::pair(1, tdim + 1), md::full_extent,
971 md::full_extent, 0);
972
973 std::function<void(std::span<T>, std::span<const std::uint32_t>, std::int32_t,
974 int)>
975 apply_inv_trans_dof_transformation
976 = element.template dof_transformation_fn<T>(
978
979 // Get interpolation operator
980 const auto [_Pi, pi_shape] = element.interpolation_operator();
981 md::mdspan<const U, std::dextents<std::size_t, 2>> Pi(_Pi.data(), pi_shape);
982
983 using u_t = md::mdspan<const T, md::dextents<std::size_t, 2>>;
984 using U_t
985 = decltype(md::submdspan(ref_data, 0, md::full_extent, md::full_extent));
986 using J_t = md::mdspan<const U, md::dextents<std::size_t, 2>>;
987 using K_t = md::mdspan<const U, md::dextents<std::size_t, 2>>;
988 auto pull_back_fn
989 = element.basix_element().template map_fn<U_t, u_t, J_t, K_t>();
990
991 std::vector<T> coeffs_b(num_scalar_dofs);
992 for (auto cell_it = cells.begin(); cell_it != cells.end(); ++cell_it)
993 {
994 auto x_dofs = md::submdspan(x_dofmap, *cell_it, md::full_extent);
995 for (int i = 0; i < num_dofs_g; ++i)
996 {
997 const int pos = 3 * x_dofs[i];
998 for (int j = 0; j < gdim; ++j)
999 coord_dofs(i, j) = x_g[pos + j];
1000 }
1001
1002 // Compute J, detJ and K
1003 std::ranges::fill(J_b, 0);
1004 for (std::size_t p = 0; p < Xshape[0]; ++p)
1005 {
1006 auto _dphi = md::submdspan(dphi, md::full_extent, p, md::full_extent);
1007 auto _J = md::submdspan(J, p, md::full_extent, md::full_extent);
1008 cmap.compute_jacobian(_dphi, coord_dofs, _J);
1009 auto _K = md::submdspan(K, p, md::full_extent, md::full_extent);
1010 cmap.compute_jacobian_inverse(_J, _K);
1011 detJ[p] = cmap.compute_jacobian_determinant(_J, det_scratch);
1012 }
1013
1014 const std::size_t offset = std::distance(cells.begin(), cell_it);
1015 std::span<const std::int32_t> dofs = dofmap.cell_dofs(*cell_it);
1016 for (int k = 0; k < element_bs; ++k)
1017 {
1018 // Extract computed expression values for element block k
1019 for (int m = 0; m < value_size; ++m)
1020 {
1021 for (std::size_t k0 = 0; k0 < Xshape[0]; ++k0)
1022 {
1023 _vals(k0, 0, m)
1024 = f[fshape[1] * (k * value_size + m) + offset * Xshape[0] + k0];
1025 }
1026 }
1027
1028 // Get element degrees of freedom for block
1029 for (std::size_t i = 0; i < Xshape[0]; ++i)
1030 {
1031 auto _u = md::submdspan(_vals, i, md::full_extent, md::full_extent);
1032 auto _U = md::submdspan(ref_data, i, md::full_extent, md::full_extent);
1033 auto _K = md::submdspan(K, i, md::full_extent, md::full_extent);
1034 auto _J = md::submdspan(J, i, md::full_extent, md::full_extent);
1035 pull_back_fn(_U, _u, _K, 1.0 / detJ[i], _J);
1036 }
1037
1038 auto ref = md::submdspan(ref_data, md::full_extent, 0, md::full_extent);
1039 impl::interpolation_apply(Pi, ref, std::span(coeffs_b), element_bs);
1040 apply_inv_trans_dof_transformation(coeffs_b, cell_info, *cell_it, 1);
1041
1042 // Copy interpolation dofs into coefficient vector
1043 assert(coeffs_b.size() == static_cast<std::size_t>(num_scalar_dofs));
1044 for (int i = 0; i < num_scalar_dofs; ++i)
1045 {
1046 const int dof = i * element_bs + k;
1047 std::div_t pos = std::div(dof, dofmap_bs);
1048 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = coeffs_b[i];
1049 }
1050 }
1051 }
1052}
1053
1054//----------------------------------------------------------------------------
1055} // namespace impl
1056
1075template <std::floating_point T>
1077 const mesh::Geometry<T>& geometry0, const FiniteElement<T>& element0,
1078 const mesh::Mesh<T>& mesh1, CellRange auto&& cells, T padding)
1079{
1080 // Collect all the points at which values are needed to define the
1081 // interpolating function
1082 std::vector<T> coords = interpolation_coords(element0, geometry0, cells);
1083
1084 // Transpose interpolation coords
1085 std::vector<T> x(coords.size());
1086 std::size_t num_points = coords.size() / 3;
1087 for (std::size_t i = 0; i < num_points; ++i)
1088 for (std::size_t j = 0; j < 3; ++j)
1089 x[3 * i + j] = coords[i + j * num_points];
1090
1091 // Determine ownership of each point
1092 return geometry::determine_point_ownership<T>(mesh1, x, padding,
1093 std::nullopt);
1094}
1095
1096template <dolfinx::scalar T, std::floating_point U>
1097void interpolate(Function<T, U>& u, std::span<const T> f,
1098 std::array<std::size_t, 2> fshape, CellRange auto&& cells)
1099{
1100 // TODO: Index for mixed-topology, zero for now
1101 const int index = 0;
1102 auto element = u.function_space()->elements(index);
1103 assert(element);
1104 const int element_bs = element->block_size();
1105 if (int num_sub = element->num_sub_elements();
1106 num_sub > 0 and num_sub != element_bs)
1107 {
1108 throw std::runtime_error("Cannot directly interpolate a mixed space. "
1109 "Interpolate into subspaces.");
1110 }
1111
1112 // Get mesh
1113 assert(u.function_space());
1114 auto mesh = u.function_space()->mesh();
1115 assert(mesh);
1116
1117 if (fshape[0]
1118 != (std::size_t)u.function_space()->elements(index)->value_size()
1119 or f.size() != fshape[0] * fshape[1])
1120 {
1121 throw std::runtime_error("Interpolation data has the wrong shape/size.");
1122 }
1123
1124 spdlog::debug("Check for dof transformation");
1125 std::span<const std::uint32_t> cell_info;
1126 if (element->needs_dof_transformations())
1127 {
1128 mesh->topology_mutable()->create_entity_permutations();
1129 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
1130 }
1131
1132 // Get dofmap
1133 spdlog::debug("Interpolate: get dofmap");
1134 const auto dofmap = u.function_space()->dofmaps(index);
1135 assert(dofmap);
1136
1137 // Result will be stored to coeffs
1138 std::span<T> coeffs = u.x()->array();
1139
1140 if (bool symmetric = u.function_space()->symmetric();
1141 element->map_ident() and element->interpolation_ident())
1142 {
1143 // This assumes that any element with an identity interpolation
1144 // matrix is a point evaluation
1145 spdlog::debug("Interpolate: point evaluation");
1146 impl::point_evaluation(*element, symmetric, *dofmap, cells, cell_info, f,
1147 fshape, coeffs);
1148 }
1149 else if (element->map_ident())
1150 {
1151 spdlog::debug("Interpolate: identity-mapped evaluation");
1152 impl::identity_mapped_evaluation(*element, symmetric, *dofmap, cells,
1153 cell_info, f, fshape, coeffs);
1154 }
1155 else
1156 {
1157 spdlog::debug("Interpolate: Piola-mapped evaluation");
1158 impl::piola_mapped_evaluation(*element, symmetric, *dofmap, cells,
1159 cell_info, f, fshape, *mesh, coeffs);
1160 }
1161}
1162
1179template <dolfinx::scalar T, std::floating_point U>
1181 CellRange auto&& cells, double tol, int maxit,
1182 const geometry::PointOwnershipData<U>& interpolation_data)
1183{
1184 auto mesh1 = u1.function_space()->mesh();
1185 assert(mesh1);
1186 MPI_Comm comm = mesh1->comm();
1187 {
1188 assert(u0.function_space());
1189 auto mesh0 = u0.function_space()->mesh();
1190 assert(mesh0);
1191 int result;
1192 MPI_Comm_compare(comm, mesh0->comm(), &result);
1193 if (result == MPI_UNEQUAL)
1194 {
1195 throw std::runtime_error("Interpolation on different meshes is only "
1196 "supported on the same communicator.");
1197 }
1198 }
1199
1200 assert(mesh1->topology());
1201 auto cell_map = mesh1->topology()->index_map(mesh1->topology()->dim());
1202 assert(cell_map);
1203 auto element1 = u1.function_space()->element();
1204 assert(element1);
1205 const std::size_t value_size = element1->value_size();
1206
1207 const std::vector<int>& dest_ranks = interpolation_data.src_owner;
1208 const std::vector<int>& src_ranks = interpolation_data.dest_owners;
1209 const std::vector<U>& recv_points = interpolation_data.dest_points;
1210 const std::vector<std::int32_t>& evaluation_cells
1211 = interpolation_data.dest_cells;
1212
1213 // Evaluate the interpolating function where possible
1214 std::vector<T> send_values(recv_points.size() / 3 * value_size);
1215 u0.eval(recv_points, {recv_points.size() / 3, (std::size_t)3},
1216 evaluation_cells, send_values, {recv_points.size() / 3, value_size},
1217 tol, maxit);
1218
1219 // Send values back to owning process
1220 std::vector<T> values_b(dest_ranks.size() * value_size);
1221 md::mdspan<const T, md::dextents<std::size_t, 2>> _send_values(
1222 send_values.data(), src_ranks.size(), value_size);
1223 impl::scatter_values(comm, src_ranks, dest_ranks, _send_values,
1224 std::span(values_b));
1225
1226 // Transpose received data
1227 md::mdspan<const T, md::dextents<std::size_t, 2>> values(
1228 values_b.data(), dest_ranks.size(), value_size);
1229 std::vector<T> valuesT_b(value_size * dest_ranks.size());
1230 md::mdspan<T, md::dextents<std::size_t, 2>> valuesT(
1231 valuesT_b.data(), value_size, dest_ranks.size());
1232 for (std::size_t i = 0; i < values.extent(0); ++i)
1233 for (std::size_t j = 0; j < values.extent(1); ++j)
1234 valuesT(j, i) = values(i, j);
1235
1236 // Call local interpolation operator
1237 fem::interpolate<T>(u1, valuesT_b, {valuesT.extent(0), valuesT.extent(1)},
1238 cells);
1239}
1240
1257template <dolfinx::scalar T, std::floating_point U>
1258void interpolate(Function<T, U>& u1, CellRange auto&& cells1,
1259 const Function<T, U>& u0, CellRange auto&& cells0)
1260{
1261 if (cells0.size() != cells1.size())
1262 throw std::runtime_error("Length of cell lists do not match.");
1263
1264 auto V1 = u1.function_space();
1265 assert(V1);
1266 auto V0 = u0.function_space();
1267 assert(V0);
1268
1269 // Get elements and check value shape
1270 auto e0 = V0->element();
1271 assert(e0);
1272 auto e1 = V1->element();
1273 assert(e1);
1274 if (!std::ranges::equal(e0->value_shape(), e1->value_shape()))
1275 {
1276 throw std::runtime_error(
1277 "Interpolation: elements have different value dimensions");
1278 }
1279
1280 if (V1->mesh() == V0->mesh() and (e1 == e0 or *e1 == *e0))
1281 {
1282 // Same element and same mesh
1283 if (e1->block_size() != e0->block_size())
1284 throw std::runtime_error("Mismatch in element block size.");
1285
1286 // Get dofmaps
1287 std::shared_ptr<const DofMap> dofmap0 = V0->dofmap();
1288 assert(dofmap0);
1289 std::shared_ptr<const DofMap> dofmap1 = V1->dofmap();
1290 assert(dofmap1);
1291
1292 // Iterate over mesh and interpolate on each cell
1293 const int bs0 = dofmap0->bs();
1294 const int bs1 = dofmap1->bs();
1295 std::span<T> u1_array = u1.x()->array();
1296 std::span<const T> u0_array = u0.x()->array();
1297 assert(cells0.size() == cells1.size());
1298 for (auto cell0_it = cells0.begin(), cell1_it = cells1.begin();
1299 cell0_it != cells0.end() and cell1_it != cells1.end();
1300 ++cell0_it, ++cell1_it)
1301
1302 {
1303 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(*cell0_it);
1304 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(*cell1_it);
1305 assert(bs0 * dofs0.size() == bs1 * dofs1.size());
1306 for (std::size_t i = 0; i < dofs0.size(); ++i)
1307 {
1308 for (int k = 0; k < bs0; ++k)
1309 {
1310 int index = bs0 * i + k;
1311 std::div_t dv1 = std::div(index, bs1);
1312 u1_array[bs1 * dofs1[dv1.quot] + dv1.rem]
1313 = u0_array[bs0 * dofs0[i] + k];
1314 }
1315 }
1316 }
1317 }
1318 else if (e1->map_type() == e0->map_type())
1319 {
1320 // Different elements, same basis function map type
1321 impl::interpolate_same_map(u1, cells1, u0, cells0);
1322 }
1323 else
1324 {
1325 // Different elements with different maps for basis functions
1326 impl::interpolate_nonmatching_maps(u1, cells1, u0, cells0);
1327 }
1328}
1329
1339template <dolfinx::scalar T, std::floating_point U>
1341 std::ranges::input_range auto&& cells)
1342{
1343 assert(u1.function_space());
1344 assert(u0.function_space());
1345 if (u1.function_space()->mesh() == u0.function_space()->mesh())
1346 interpolate<T, U>(u1, cells, u0, cells);
1347 else
1348 throw std::runtime_error("Meshes do no match.");
1349}
1350
1360template <dolfinx::scalar T, std::floating_point U>
1362{
1363 assert(u1.function_space());
1364 assert(u0.function_space());
1365 if (auto V1 = u1.function_space(); V1 == u0.function_space())
1366 std::ranges::copy(u0.x()->array(), u1.x()->array().begin());
1367 else
1368 {
1369 auto mesh = V1->mesh();
1370 assert(mesh);
1371 assert(mesh->topology());
1372 auto map = mesh->topology()->index_map(mesh->topology()->dim());
1373 assert(map);
1374 std::int32_t num_cells = map->size_local() + map->num_ghosts();
1375 interpolate<T, U>(u1, u0, std::ranges::views::iota(0, num_cells));
1376 }
1377}
1378} // namespace dolfinx::fem
Degree-of-freedom map representations and tools.
Definition CoordinateElement.h:38
void tabulate(int nd, std::span< const T > X, std::array< std::size_t, 2 > shape, std::span< T > basis) const
Evaluate basis values and derivatives at set of points.
Definition CoordinateElement.cpp:55
std::array< std::size_t, 4 > tabulate_shape(std::size_t nd, std::size_t num_points) const
Shape of array to fill when calling tabulate.
Definition CoordinateElement.cpp:48
int dim() const
The dimension of the coordinate element space.
Definition CoordinateElement.cpp:205
Model of a finite element.
Definition FiniteElement.h:57
std::pair< std::vector< geometry_type >, std::array< std::size_t, 2 > > interpolation_points() const
Points on the reference cell at which an expression needs to be evaluated in order to interpolate the...
Definition FiniteElement.cpp:464
mesh::CellType cell_type() const noexcept
Cell shape that the element is defined on.
Definition FiniteElement.cpp:279
Definition Function.h:47
std::shared_ptr< const FunctionSpace< geometry_type > > function_space() const
Access the function space.
Definition Function.h:147
void eval(std::span< const geometry_type > x, std::array< std::size_t, 2 > xshape, CellRange auto &&cells, std::span< value_type > u, std::array< std::size_t, 2 > ushape, double tol, int maxit) const
Evaluate the Function at points.
Definition Function.h:457
std::shared_ptr< const la::Vector< value_type > > x() const
Underlying vector (const version).
Definition Function.h:153
Geometry stores the geometry imposed on a mesh.
Definition Geometry.h:34
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition Mesh.h:23
Requirement on range of cell indices.
Definition interpolate.h:41
Definition interpolate.h:33
MPI_Datatype mpi_t
Retrieves the MPI data type associated to the provided type.
Definition MPI.h:280
int rank(MPI_Comm comm)
Return process rank for the communicator.
Definition MPI.cpp:64
void cells(la::SparsityPattern &pattern, const std::pair< R0, R1 > &cells, std::array< std::reference_wrapper< const DofMap >, 2 > dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition sparsitybuild.h:37
Finite element method functionality.
Definition assemble_expression_impl.h:23
std::vector< T > interpolation_coords(const fem::FiniteElement< T > &element, const mesh::Geometry< T > &geometry, CellRange auto &&cells)
Compute the evaluation points in the physical space at which an expression should be computed to inte...
Definition interpolate.h:56
void interpolate(Function< T, U > &u, std::span< const T > f, std::array< std::size_t, 2 > fshape, CellRange auto &&cells)
Interpolate an evaluated expression f(x) in a finite element space.
Definition interpolate.h:1097
@ transpose
Transpose.
Definition FiniteElement.h:28
@ inverse_transpose
Transpose inverse.
Definition FiniteElement.h:30
@ standard
Standard.
Definition FiniteElement.h:27
geometry::PointOwnershipData< T > create_interpolation_data(const mesh::Geometry< T > &geometry0, const FiniteElement< T > &element0, const mesh::Mesh< T > &mesh1, CellRange auto &&cells, T padding)
Generate data needed to interpolate finite element fem::Function's across different meshes.
Definition interpolate.h:1076
Geometry data structures and algorithms.
Definition BoundingBoxTree.h:22
PointOwnershipData< T > determine_point_ownership(const mesh::Mesh< T > &mesh, std::span< const T > points, T padding, std::optional< std::span< const std::int32_t > > cells)
Given a set of points, determine which process is colliding, using the GJK algorithm on cells to dete...
Definition utils.h:683
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32
CellType
Cell type identifier.
Definition cell_types.h:21
Information on the ownership of points distributed across processes.
Definition utils.h:30
std::vector< T > dest_points
Points that are owned by current process.
Definition utils.h:35
std::vector< std::int32_t > dest_cells
Definition utils.h:37
std::vector< int > dest_owners
Ranks that sent dest_points to current process.
Definition utils.h:34
std::vector< int > src_owner
Definition utils.h:31