DOLFINx 0.9.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
interpolate.h
1// Copyright (C) 2020-2024 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 <span>
24#include <vector>
25
26namespace dolfinx::fem
27{
28template <dolfinx::scalar T, std::floating_point U>
29class Function;
30
31template <typename T>
32concept MDSpan = requires(T x, std::size_t idx) {
33 x(idx, idx);
34 { x.extent(0) } -> std::integral;
35 { x.extent(1) } -> std::integral;
36};
37
49template <std::floating_point T>
50std::vector<T> interpolation_coords(const fem::FiniteElement<T>& element,
51 const mesh::Geometry<T>& geometry,
52 std::span<const std::int32_t> cells)
53{
54 // Get geometry data and the element coordinate map
55 const std::size_t gdim = geometry.dim();
56 auto x_dofmap = geometry.dofmap();
57 std::span<const T> x_g = geometry.x();
58
59 const CoordinateElement<T>& cmap = geometry.cmap();
60 const std::size_t num_dofs_g = cmap.dim();
61
62 // Get the interpolation points on the reference cells
63 const auto [X, Xshape] = element.interpolation_points();
64
65 // Evaluate coordinate element basis at reference points
66 std::array<std::size_t, 4> phi_shape = cmap.tabulate_shape(0, Xshape[0]);
67 std::vector<T> phi_b(
68 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
69 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
70 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>
71 phi_full(phi_b.data(), phi_shape);
72 cmap.tabulate(0, X, Xshape, phi_b);
73 auto phi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
74 phi_full, 0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
75 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
76
77 // Push reference coordinates (X) forward to the physical coordinates
78 // (x) for each cell
79 std::vector<T> coordinate_dofs(num_dofs_g * gdim, 0);
80 std::vector<T> x(3 * (cells.size() * Xshape[0]), 0);
81 for (std::size_t c = 0; c < cells.size(); ++c)
82 {
83 // Get geometry data for current cell
84 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
85 x_dofmap, cells[c], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
86 for (std::size_t i = 0; i < x_dofs.size(); ++i)
87 {
88 std::copy_n(std::next(x_g.begin(), 3 * x_dofs[i]), gdim,
89 std::next(coordinate_dofs.begin(), i * gdim));
90 }
91
92 // Push forward coordinates (X -> x)
93 for (std::size_t p = 0; p < Xshape[0]; ++p)
94 {
95 for (std::size_t j = 0; j < gdim; ++j)
96 {
97 T acc = 0;
98 for (std::size_t k = 0; k < num_dofs_g; ++k)
99 acc += phi(p, k) * coordinate_dofs[k * gdim + j];
100 x[j * (cells.size() * Xshape[0]) + c * Xshape[0] + p] = acc;
101 }
102 }
103 }
104
105 return x;
106}
107
123template <dolfinx::scalar T, std::floating_point U>
124void interpolate(Function<T, U>& u, std::span<const T> f,
125 std::array<std::size_t, 2> fshape,
126 std::span<const std::int32_t> cells);
127
128namespace impl
129{
131template <typename T, std::size_t D>
132using mdspan_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
133 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, D>>;
134
154template <dolfinx::scalar T>
155void scatter_values(MPI_Comm comm, std::span<const std::int32_t> src_ranks,
156 std::span<const std::int32_t> dest_ranks,
157 mdspan_t<const T, 2> send_values, std::span<T> recv_values)
158{
159 const std::size_t block_size = send_values.extent(1);
160 assert(src_ranks.size() * block_size == send_values.size());
161 assert(recv_values.size() == dest_ranks.size() * block_size);
162
163 // Build unique set of the sorted src_ranks
164 std::vector<std::int32_t> out_ranks(src_ranks.size());
165 out_ranks.assign(src_ranks.begin(), src_ranks.end());
166 auto [unique_end, range_end] = std::ranges::unique(out_ranks);
167 out_ranks.erase(unique_end, range_end);
168 out_ranks.reserve(out_ranks.size() + 1);
169
170 // Remove negative entries from dest_ranks
171 std::vector<std::int32_t> in_ranks;
172 in_ranks.reserve(dest_ranks.size());
173 std::copy_if(dest_ranks.begin(), dest_ranks.end(),
174 std::back_inserter(in_ranks),
175 [](auto rank) { return rank >= 0; });
176
177 // Create unique set of sorted in-ranks
178 {
179 std::ranges::sort(in_ranks);
180 auto [unique_end, range_end] = std::ranges::unique(in_ranks);
181 in_ranks.erase(unique_end, range_end);
182 }
183 in_ranks.reserve(in_ranks.size() + 1);
184
185 // Create neighborhood communicator
186 MPI_Comm reverse_comm;
187 MPI_Dist_graph_create_adjacent(
188 comm, in_ranks.size(), in_ranks.data(), MPI_UNWEIGHTED, out_ranks.size(),
189 out_ranks.data(), MPI_UNWEIGHTED, MPI_INFO_NULL, false, &reverse_comm);
190
191 std::vector<std::int32_t> comm_to_output;
192 std::vector<std::int32_t> recv_sizes(in_ranks.size());
193 recv_sizes.reserve(1);
194 std::vector<std::int32_t> recv_offsets(in_ranks.size() + 1, 0);
195 {
196 // Build map from parent to neighborhood communicator ranks
197 std::vector<std::pair<std::int32_t, std::int32_t>> rank_to_neighbor;
198 rank_to_neighbor.reserve(in_ranks.size());
199 for (std::size_t i = 0; i < in_ranks.size(); i++)
200 rank_to_neighbor.push_back({in_ranks[i], i});
201 std::ranges::sort(rank_to_neighbor);
202
203 // Compute receive sizes
204 std::ranges::for_each(
205 dest_ranks,
206 [&dest_ranks, &rank_to_neighbor, &recv_sizes, block_size](auto rank)
207 {
208 if (rank >= 0)
209 {
210 auto it = std::ranges::lower_bound(rank_to_neighbor, rank,
211 std::ranges::less(),
212 [](auto e) { return e.first; });
213 assert(it != rank_to_neighbor.end() and it->first == rank);
214 recv_sizes[it->second] += block_size;
215 }
216 });
217
218 // Compute receiving offsets
219 std::partial_sum(recv_sizes.begin(), recv_sizes.end(),
220 std::next(recv_offsets.begin(), 1));
221
222 // Compute map from receiving values to position in recv_values
223 comm_to_output.resize(recv_offsets.back() / block_size);
224 std::vector<std::int32_t> recv_counter(recv_sizes.size(), 0);
225 for (std::size_t i = 0; i < dest_ranks.size(); ++i)
226 {
227 if (const std::int32_t rank = dest_ranks[i]; rank >= 0)
228 {
229 auto it = std::ranges::lower_bound(rank_to_neighbor, rank,
230 std::ranges::less(),
231 [](auto e) { return e.first; });
232 assert(it != rank_to_neighbor.end() and it->first == rank);
233 int insert_pos = recv_offsets[it->second] + recv_counter[it->second];
234 comm_to_output[insert_pos / block_size] = i * block_size;
235 recv_counter[it->second] += block_size;
236 }
237 }
238 }
239
240 std::vector<std::int32_t> send_sizes(out_ranks.size());
241 send_sizes.reserve(1);
242 {
243 // Compute map from parent MPI rank to neighbor rank for outgoing
244 // data. `out_ranks` is sorted, so rank_to_neighbor will be sorted
245 // too.
246 std::vector<std::pair<std::int32_t, std::int32_t>> rank_to_neighbor;
247 rank_to_neighbor.reserve(out_ranks.size());
248 for (std::size_t i = 0; i < out_ranks.size(); i++)
249 rank_to_neighbor.push_back({out_ranks[i], i});
250
251 // Compute send sizes. As `src_ranks` is sorted, we can move 'start'
252 // in search forward.
253 auto start = rank_to_neighbor.begin();
254 std::ranges::for_each(
255 src_ranks,
256 [&rank_to_neighbor, &send_sizes, block_size, &start](auto rank)
257 {
258 auto it = std::ranges::lower_bound(start, rank_to_neighbor.end(),
259 rank, std::ranges::less(),
260 [](auto e) { return e.first; });
261 assert(it != rank_to_neighbor.end() and it->first == rank);
262 send_sizes[it->second] += block_size;
263 start = it;
264 });
265 }
266
267 // Compute sending offsets
268 std::vector<std::int32_t> send_offsets(send_sizes.size() + 1, 0);
269 std::partial_sum(send_sizes.begin(), send_sizes.end(),
270 std::next(send_offsets.begin(), 1));
271
272 // Send values to dest ranks
273 std::vector<T> values(recv_offsets.back());
274 values.reserve(1);
275 MPI_Neighbor_alltoallv(send_values.data_handle(), send_sizes.data(),
276 send_offsets.data(), dolfinx::MPI::mpi_type<T>(),
277 values.data(), recv_sizes.data(), recv_offsets.data(),
278 dolfinx::MPI::mpi_type<T>(), reverse_comm);
279 MPI_Comm_free(&reverse_comm);
280
281 // Insert values received from neighborhood communicator in output
282 // span
283 std::ranges::fill(recv_values, T(0));
284 for (std::size_t i = 0; i < comm_to_output.size(); i++)
285 {
286 auto vals = std::next(recv_values.begin(), comm_to_output[i]);
287 auto vals_from = std::next(values.begin(), i * block_size);
288 std::copy_n(vals_from, block_size, vals);
289 }
290};
291
300template <MDSpan U, MDSpan V, dolfinx::scalar T>
301void interpolation_apply(U&& Pi, V&& data, std::span<T> coeffs, int bs)
302{
303 using X = typename dolfinx::scalar_value_type_t<T>;
304
305 // Compute coefficients = Pi * x (matrix-vector multiply)
306 if (bs == 1)
307 {
308 assert(data.extent(0) * data.extent(1) == Pi.extent(1));
309 for (std::size_t i = 0; i < Pi.extent(0); ++i)
310 {
311 coeffs[i] = 0.0;
312 for (std::size_t k = 0; k < data.extent(1); ++k)
313 for (std::size_t j = 0; j < data.extent(0); ++j)
314 coeffs[i]
315 += static_cast<X>(Pi(i, k * data.extent(0) + j)) * data(j, k);
316 }
317 }
318 else
319 {
320 const std::size_t cols = Pi.extent(1);
321 assert(data.extent(0) == Pi.extent(1));
322 assert(data.extent(1) == bs);
323 for (int k = 0; k < bs; ++k)
324 {
325 for (std::size_t i = 0; i < Pi.extent(0); ++i)
326 {
327 T acc = 0;
328 for (std::size_t j = 0; j < cols; ++j)
329 acc += static_cast<X>(Pi(i, j)) * data(j, k);
330 coeffs[bs * i + k] = acc;
331 }
332 }
333 }
334}
335
355template <dolfinx::scalar T, std::floating_point U>
356void interpolate_same_map(Function<T, U>& u1, const Function<T, U>& u0,
357 std::span<const std::int32_t> cells1,
358 std::span<const std::int32_t> cells0)
359{
360 auto V0 = u0.function_space();
361 assert(V0);
362 auto V1 = u1.function_space();
363 assert(V1);
364 auto mesh0 = V0->mesh();
365 assert(mesh0);
366
367 auto mesh1 = V1->mesh();
368 assert(mesh1);
369
370 auto element0 = V0->element();
371 assert(element0);
372 auto element1 = V1->element();
373 assert(element1);
374
375 assert(mesh0->topology()->dim());
376 const int tdim = mesh0->topology()->dim();
377 auto map = mesh0->topology()->index_map(tdim);
378 assert(map);
379 std::span<T> u1_array = u1.x()->mutable_array();
380 std::span<const T> u0_array = u0.x()->array();
381
382 std::span<const std::uint32_t> cell_info0;
383 std::span<const std::uint32_t> cell_info1;
384 if (element1->needs_dof_transformations()
385 or element0->needs_dof_transformations())
386 {
387 mesh0->topology_mutable()->create_entity_permutations();
388 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
389 mesh1->topology_mutable()->create_entity_permutations();
390 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
391 }
392
393 // Get dofmaps
394 auto dofmap1 = V1->dofmap();
395 auto dofmap0 = V0->dofmap();
396
397 // Get block sizes and dof transformation operators
398 const int bs1 = dofmap1->bs();
399 const int bs0 = dofmap0->bs();
400 auto apply_dof_transformation = element0->template dof_transformation_fn<T>(
402 auto apply_inverse_dof_transform
403 = element1->template dof_transformation_fn<T>(
405
406 // Create working array
407 std::vector<T> local0(element0->space_dimension());
408 std::vector<T> local1(element1->space_dimension());
409
410 // Create interpolation operator
411 const auto [i_m, im_shape]
412 = element1->create_interpolation_operator(*element0);
413
414 // Iterate over mesh and interpolate on each cell
415 using X = typename dolfinx::scalar_value_type_t<T>;
416 for (std::size_t c = 0; c < cells0.size(); c++)
417 {
418 // Pack and transform cell dofs to reference ordering
419 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(cells0[c]);
420 for (std::size_t i = 0; i < dofs0.size(); ++i)
421 for (int k = 0; k < bs0; ++k)
422 local0[bs0 * i + k] = u0_array[bs0 * dofs0[i] + k];
423
424 apply_dof_transformation(local0, cell_info0, cells0[c], 1);
425
426 // FIXME: Get compile-time ranges from Basix
427 // Apply interpolation operator
428 std::ranges::fill(local1, 0);
429 for (std::size_t i = 0; i < im_shape[0]; ++i)
430 for (std::size_t j = 0; j < im_shape[1]; ++j)
431 local1[i] += static_cast<X>(i_m[im_shape[1] * i + j]) * local0[j];
432
433 apply_inverse_dof_transform(local1, cell_info1, cells1[c], 1);
434 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(cells1[c]);
435 for (std::size_t i = 0; i < dofs1.size(); ++i)
436 for (int k = 0; k < bs1; ++k)
437 u1_array[bs1 * dofs1[i] + k] = local1[bs1 * i + k];
438 }
439}
440
455template <dolfinx::scalar T, std::floating_point U>
456void interpolate_nonmatching_maps(Function<T, U>& u1,
457 std::span<const std::int32_t> cells1,
458 const Function<T, U>& u0,
459 std::span<const std::int32_t> cells0)
460{
461 // Get mesh
462 auto V0 = u0.function_space();
463 assert(V0);
464 auto mesh0 = V0->mesh();
465 assert(mesh0);
466
467 // Mesh dims
468 const int tdim = mesh0->topology()->dim();
469 const int gdim = mesh0->geometry().dim();
470
471 // Get elements
472 auto V1 = u1.function_space();
473 assert(V1);
474 auto mesh1 = V1->mesh();
475 assert(mesh1);
476 auto element0 = V0->element();
477 assert(element0);
478 auto element1 = V1->element();
479 assert(element1);
480
481 std::span<const std::uint32_t> cell_info0;
482 std::span<const std::uint32_t> cell_info1;
483 if (element1->needs_dof_transformations()
484 or element0->needs_dof_transformations())
485 {
486 mesh0->topology_mutable()->create_entity_permutations();
487 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
488 mesh1->topology_mutable()->create_entity_permutations();
489 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
490 }
491
492 // Get dofmaps
493 auto dofmap0 = V0->dofmap();
494 auto dofmap1 = V1->dofmap();
495
496 const auto [X, Xshape] = element1->interpolation_points();
497
498 // Get block sizes and dof transformation operators
499 const int bs0 = element0->block_size();
500 const int bs1 = element1->block_size();
501 auto apply_dof_transformation0 = element0->template dof_transformation_fn<U>(
503 auto apply_inverse_dof_transform1
504 = element1->template dof_transformation_fn<T>(
506
507 // Get sizes of elements
508 const std::size_t dim0 = element0->space_dimension() / bs0;
509 const std::size_t value_size_ref0 = element0->reference_value_size() / bs0;
510 const std::size_t value_size0 = V0->value_size() / bs0;
511
512 const CoordinateElement<U>& cmap = mesh0->geometry().cmap();
513 auto x_dofmap = mesh0->geometry().dofmap();
514 const std::size_t num_dofs_g = cmap.dim();
515 std::span<const U> x_g = mesh0->geometry().x();
516
517 // Evaluate coordinate map basis at reference interpolation points
518 const std::array<std::size_t, 4> phi_shape
519 = cmap.tabulate_shape(1, Xshape[0]);
520 std::vector<U> phi_b(
521 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
522 mdspan_t<const U, 4> phi(phi_b.data(), phi_shape);
523 cmap.tabulate(1, X, Xshape, phi_b);
524
525 // Evaluate v basis functions at reference interpolation points
526 const auto [_basis_derivatives_reference0, b0shape]
527 = element0->tabulate(X, Xshape, 0);
528 mdspan_t<const U, 4> basis_derivatives_reference0(
529 _basis_derivatives_reference0.data(), b0shape);
530
531 // Create working arrays
532 std::vector<T> local1(element1->space_dimension());
533 std::vector<T> coeffs0(element0->space_dimension());
534
535 std::vector<U> basis0_b(Xshape[0] * dim0 * value_size0);
536 impl::mdspan_t<U, 3> basis0(basis0_b.data(), Xshape[0], dim0, value_size0);
537
538 std::vector<U> basis_reference0_b(Xshape[0] * dim0 * value_size_ref0);
539 impl::mdspan_t<U, 3> basis_reference0(basis_reference0_b.data(), Xshape[0],
540 dim0, value_size_ref0);
541
542 std::vector<T> values0_b(Xshape[0] * 1 * V1->value_size());
543 impl::mdspan_t<T, 3> values0(values0_b.data(), Xshape[0], 1,
544 V1->value_size());
545
546 std::vector<T> mapped_values_b(Xshape[0] * 1 * V1->value_size());
547 impl::mdspan_t<T, 3> mapped_values0(mapped_values_b.data(), Xshape[0], 1,
548 V1->value_size());
549
550 std::vector<U> coord_dofs_b(num_dofs_g * gdim);
551 impl::mdspan_t<U, 2> coord_dofs(coord_dofs_b.data(), num_dofs_g, gdim);
552
553 std::vector<U> J_b(Xshape[0] * gdim * tdim);
554 impl::mdspan_t<U, 3> J(J_b.data(), Xshape[0], gdim, tdim);
555 std::vector<U> K_b(Xshape[0] * tdim * gdim);
556 impl::mdspan_t<U, 3> K(K_b.data(), Xshape[0], tdim, gdim);
557 std::vector<U> detJ(Xshape[0]);
558 std::vector<U> det_scratch(2 * gdim * tdim);
559
560 // Get interpolation operator
561 const auto [_Pi_1, pi_shape] = element1->interpolation_operator();
562 impl::mdspan_t<const U, 2> Pi_1(_Pi_1.data(), pi_shape);
563
564 using u_t = impl::mdspan_t<U, 2>;
565 using U_t = impl::mdspan_t<const U, 2>;
566 using J_t = impl::mdspan_t<const U, 2>;
567 using K_t = impl::mdspan_t<const U, 2>;
568 auto push_forward_fn0
569 = element0->basix_element().template map_fn<u_t, U_t, J_t, K_t>();
570
571 using v_t = impl::mdspan_t<const T, 2>;
572 using V_t = impl::mdspan_t<T, 2>;
573 auto pull_back_fn1
574 = element1->basix_element().template map_fn<V_t, v_t, K_t, J_t>();
575
576 // Iterate over mesh and interpolate on each cell
577 std::span<const T> array0 = u0.x()->array();
578 std::span<T> array1 = u1.x()->mutable_array();
579 for (std::size_t c = 0; c < cells0.size(); c++)
580 {
581 // Get cell geometry (coordinate dofs)
582 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
583 x_dofmap, cells0[c], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
584 for (std::size_t i = 0; i < num_dofs_g; ++i)
585 {
586 const int pos = 3 * x_dofs[i];
587 for (std::size_t j = 0; j < gdim; ++j)
588 coord_dofs(i, j) = x_g[pos + j];
589 }
590
591 // Compute Jacobians and reference points for current cell
592 std::ranges::fill(J_b, 0);
593 for (std::size_t p = 0; p < Xshape[0]; ++p)
594 {
595 auto dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
596 phi, std::pair(1, tdim + 1), p,
597 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
598
599 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
600 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
601 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
602 cmap.compute_jacobian(dphi, coord_dofs, _J);
603 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
604 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
605 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
606 cmap.compute_jacobian_inverse(_J, _K);
607 detJ[p] = cmap.compute_jacobian_determinant(_J, det_scratch);
608 }
609
610 // Copy evaluated basis on reference, apply DOF transformations, and
611 // push forward to physical element
612 for (std::size_t k0 = 0; k0 < basis_reference0.extent(0); ++k0)
613 for (std::size_t k1 = 0; k1 < basis_reference0.extent(1); ++k1)
614 for (std::size_t k2 = 0; k2 < basis_reference0.extent(2); ++k2)
615 basis_reference0(k0, k1, k2)
616 = basis_derivatives_reference0(0, k0, k1, k2);
617
618 for (std::size_t p = 0; p < Xshape[0]; ++p)
619 {
620 apply_dof_transformation0(
621 std::span(basis_reference0_b.data() + p * dim0 * value_size_ref0,
622 dim0 * value_size_ref0),
623 cell_info0, cells0[c], value_size_ref0);
624 }
625
626 for (std::size_t i = 0; i < basis0.extent(0); ++i)
627 {
628 auto _u = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
629 basis0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
630 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
631 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
632 basis_reference0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
633 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
634 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
635 K, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
636 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
637 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
638 J, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
639 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
640 push_forward_fn0(_u, _U, _J, detJ[i], _K);
641 }
642
643 // Copy expansion coefficients for v into local array
644 const int dof_bs0 = dofmap0->bs();
645 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(cells0[c]);
646 for (std::size_t i = 0; i < dofs0.size(); ++i)
647 for (int k = 0; k < dof_bs0; ++k)
648 coeffs0[dof_bs0 * i + k] = array0[dof_bs0 * dofs0[i] + k];
649
650 // Evaluate v at the interpolation points (physical space values)
651 using X = typename dolfinx::scalar_value_type_t<T>;
652 for (std::size_t p = 0; p < Xshape[0]; ++p)
653 {
654 for (int k = 0; k < bs0; ++k)
655 {
656 for (std::size_t j = 0; j < value_size0; ++j)
657 {
658 T acc = 0;
659 for (std::size_t i = 0; i < dim0; ++i)
660 acc += coeffs0[bs0 * i + k] * static_cast<X>(basis0(p, i, j));
661 values0(p, 0, j * bs0 + k) = acc;
662 }
663 }
664 }
665
666 // Pull back the physical values to the u reference
667 for (std::size_t i = 0; i < values0.extent(0); ++i)
668 {
669 auto _u = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
670 values0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
671 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
672 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
673 mapped_values0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
674 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
675 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
676 K, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
677 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
678 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
679 J, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
680 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
681 pull_back_fn1(_U, _u, _K, 1.0 / detJ[i], _J);
682 }
683
684 auto values = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
685 mapped_values0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0,
686 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
687 interpolation_apply(Pi_1, values, std::span(local1), bs1);
688 apply_inverse_dof_transform1(local1, cell_info1, cells1[c], 1);
689
690 // Copy local coefficients to the correct position in u dof array
691 const int dof_bs1 = dofmap1->bs();
692 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(c);
693 for (std::size_t i = 0; i < dofs1.size(); ++i)
694 for (int k = 0; k < dof_bs1; ++k)
695 array1[dof_bs1 * dofs1[i] + k] = local1[dof_bs1 * i + k];
696 }
697}
698
699//----------------------------------------------------------------------------
700} // namespace impl
701
702template <dolfinx::scalar T, std::floating_point U>
703void interpolate(Function<T, U>& u, std::span<const T> f,
704 std::array<std::size_t, 2> fshape,
705 std::span<const std::int32_t> cells)
706{
707 using cmdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
708 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
709 using cmdspan4_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
710 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>;
711 using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
712 U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
713 using mdspan3_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
714 U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>;
715 auto element = u.function_space()->element();
716 assert(element);
717 const int element_bs = element->block_size();
718 if (int num_sub = element->num_sub_elements();
719 num_sub > 0 and num_sub != element_bs)
720 {
721 throw std::runtime_error("Cannot directly interpolate a mixed space. "
722 "Interpolate into subspaces.");
723 }
724
725 // Get mesh
726 assert(u.function_space());
727 auto mesh = u.function_space()->mesh();
728 assert(mesh);
729
730 const int gdim = mesh->geometry().dim();
731 const int tdim = mesh->topology()->dim();
732 const bool symmetric = u.function_space()->symmetric();
733
734 if (fshape[0] != (std::size_t)u.function_space()->value_size())
735 throw std::runtime_error("Interpolation data has the wrong shape/size.");
736
737 std::span<const std::uint32_t> cell_info;
738 if (element->needs_dof_transformations())
739 {
740 mesh->topology_mutable()->create_entity_permutations();
741 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
742 }
743
744 const std::size_t f_shape1 = f.size() / u.function_space()->value_size();
745 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
746 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
747 _f(f.data(), fshape);
748
749 // Get dofmap
750 const auto dofmap = u.function_space()->dofmap();
751 assert(dofmap);
752 const int dofmap_bs = dofmap->bs();
753
754 // Loop over cells and compute interpolation dofs
755 const int num_scalar_dofs = element->space_dimension() / element_bs;
756 const int value_size = u.function_space()->value_size() / element_bs;
757
758 std::span<T> coeffs = u.x()->mutable_array();
759 std::vector<T> _coeffs(num_scalar_dofs);
760
761 // This assumes that any element with an identity interpolation matrix
762 // is a point evaluation
763 if (element->map_ident() && element->interpolation_ident())
764 {
765 // Point evaluation element *and* the geometric map is the identity,
766 // e.g. not Piola mapped
767
768 auto apply_inv_transpose_dof_transformation
769 = element->template dof_transformation_fn<T>(
771
772 if (symmetric)
773 {
774 std::size_t matrix_size = 0;
775 while (matrix_size * matrix_size < fshape[0])
776 ++matrix_size;
777 // Loop over cells
778 for (std::size_t c = 0; c < cells.size(); ++c)
779 {
780 // The entries of a symmetric matrix are numbered (for an example 4x4
781 // element):
782 // 0 * * *
783 // 1 2 * *
784 // 3 4 5 *
785 // 6 7 8 9
786 // The loop extracts these elements. In this loop, row is the row of
787 // this matrix, and (k - rowstart) is the column
788 std::size_t row = 0;
789 std::size_t rowstart = 0;
790 const std::int32_t cell = cells[c];
791 std::span<const std::int32_t> dofs = dofmap->cell_dofs(cell);
792 for (int k = 0; k < element_bs; ++k)
793 {
794 if (k - rowstart > row)
795 {
796 ++row;
797 rowstart = k;
798 }
799 // num_scalar_dofs is the number of interpolation points per
800 // cell in this case (interpolation matrix is identity)
801 std::copy_n(
802 std::next(f.begin(), (row * matrix_size + k - rowstart) * f_shape1
803 + c * num_scalar_dofs),
804 num_scalar_dofs, _coeffs.begin());
805 apply_inv_transpose_dof_transformation(_coeffs, cell_info, cell, 1);
806 for (int i = 0; i < num_scalar_dofs; ++i)
807 {
808 const int dof = i * element_bs + k;
809 std::div_t pos = std::div(dof, dofmap_bs);
810 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
811 }
812 }
813 }
814 }
815 else
816 {
817 // Loop over cells
818 for (std::size_t c = 0; c < cells.size(); ++c)
819 {
820 const std::int32_t cell = cells[c];
821 std::span<const std::int32_t> dofs = dofmap->cell_dofs(cell);
822 for (int k = 0; k < element_bs; ++k)
823 {
824 // num_scalar_dofs is the number of interpolation points per
825 // cell in this case (interpolation matrix is identity)
826 std::copy_n(std::next(f.begin(), k * f_shape1 + c * num_scalar_dofs),
827 num_scalar_dofs, _coeffs.begin());
828 apply_inv_transpose_dof_transformation(_coeffs, cell_info, cell, 1);
829 for (int i = 0; i < num_scalar_dofs; ++i)
830 {
831 const int dof = i * element_bs + k;
832 std::div_t pos = std::div(dof, dofmap_bs);
833 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
834 }
835 }
836 }
837 }
838 }
839 else if (element->map_ident())
840 {
841 // Not a point evaluation, but the geometric map is the identity,
842 // e.g. not Piola mapped
843
844 if (symmetric)
845 {
846 throw std::runtime_error(
847 "Interpolation into this element not supported.");
848 }
849
850 const int element_vs = u.function_space()->value_size() / element_bs;
851
852 if (element_vs > 1 && element_bs > 1)
853 {
854 throw std::runtime_error(
855 "Interpolation into this element not supported.");
856 }
857
858 // Get interpolation operator
859 const auto [_Pi, pi_shape] = element->interpolation_operator();
860 cmdspan2_t Pi(_Pi.data(), pi_shape);
861 const std::size_t num_interp_points = Pi.extent(1);
862 assert(Pi.extent(0) == num_scalar_dofs);
863
864 auto apply_inv_transpose_dof_transformation
865 = element->template dof_transformation_fn<T>(
867
868 // Loop over cells
869 std::vector<T> ref_data_b(num_interp_points);
870 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
871 T, MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
872 std::size_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent, 1>>
873 ref_data(ref_data_b.data(), num_interp_points, 1);
874 for (std::size_t c = 0; c < cells.size(); ++c)
875 {
876 const std::int32_t cell = cells[c];
877 std::span<const std::int32_t> dofs = dofmap->cell_dofs(cell);
878 for (int k = 0; k < element_bs; ++k)
879 {
880 for (int i = 0; i < element_vs; ++i)
881 {
882 std::copy_n(
883 std::next(f.begin(), (i + k) * f_shape1
884 + c * num_interp_points / element_vs),
885 num_interp_points / element_vs,
886 std::next(ref_data_b.begin(),
887 i * num_interp_points / element_vs));
888 }
889 impl::interpolation_apply(Pi, ref_data, std::span(_coeffs), 1);
890 apply_inv_transpose_dof_transformation(_coeffs, cell_info, cell, 1);
891 for (int i = 0; i < num_scalar_dofs; ++i)
892 {
893 const int dof = i * element_bs + k;
894 std::div_t pos = std::div(dof, dofmap_bs);
895 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
896 }
897 }
898 }
899 }
900 else
901 {
902 if (symmetric)
903 {
904 throw std::runtime_error(
905 "Interpolation into this element not supported.");
906 }
907 // Get the interpolation points on the reference cells
908 const auto [X, Xshape] = element->interpolation_points();
909 if (X.empty())
910 {
911 throw std::runtime_error(
912 "Interpolation into this space is not yet supported.");
913 }
914
915 if (_f.extent(1) != cells.size() * Xshape[0])
916 throw std::runtime_error("Interpolation data has the wrong shape.");
917
918 // Get coordinate map
919 const CoordinateElement<U>& cmap = mesh->geometry().cmap();
920
921 // Get geometry data
922 auto x_dofmap = mesh->geometry().dofmap();
923 const int num_dofs_g = cmap.dim();
924 std::span<const U> x_g = mesh->geometry().x();
925
926 // Create data structures for Jacobian info
927 std::vector<U> J_b(Xshape[0] * gdim * tdim);
928 mdspan3_t J(J_b.data(), Xshape[0], gdim, tdim);
929 std::vector<U> K_b(Xshape[0] * tdim * gdim);
930 mdspan3_t K(K_b.data(), Xshape[0], tdim, gdim);
931 std::vector<U> detJ(Xshape[0]);
932 std::vector<U> det_scratch(2 * gdim * tdim);
933
934 std::vector<U> coord_dofs_b(num_dofs_g * gdim);
935 mdspan2_t coord_dofs(coord_dofs_b.data(), num_dofs_g, gdim);
936 const std::size_t value_size_ref
937 = element->reference_value_size() / element_bs;
938 std::vector<T> ref_data_b(Xshape[0] * 1 * value_size_ref);
939 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
940 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
941 ref_data(ref_data_b.data(), Xshape[0], 1, value_size_ref);
942
943 std::vector<T> _vals_b(Xshape[0] * 1 * value_size);
944 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
945 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
946 _vals(_vals_b.data(), Xshape[0], 1, value_size);
947
948 // Tabulate 1st derivative of shape functions at interpolation
949 // coords
950 std::array<std::size_t, 4> phi_shape = cmap.tabulate_shape(1, Xshape[0]);
951 std::vector<U> phi_b(
952 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
953 cmdspan4_t phi(phi_b.data(), phi_shape);
954 cmap.tabulate(1, X, Xshape, phi_b);
955 auto dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
956 phi, std::pair(1, tdim + 1),
957 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
958 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
959
960 const std::function<void(std::span<T>, std::span<const std::uint32_t>,
961 std::int32_t, int)>
962 apply_inverse_transpose_dof_transformation
963 = element->template dof_transformation_fn<T>(
965
966 // Get interpolation operator
967 const auto [_Pi, pi_shape] = element->interpolation_operator();
968 cmdspan2_t Pi(_Pi.data(), pi_shape);
969
970 using u_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
971 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
972 using U_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
973 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
974 using J_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
975 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
976 using K_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
977 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
978 auto pull_back_fn
979 = element->basix_element().template map_fn<U_t, u_t, J_t, K_t>();
980
981 for (std::size_t c = 0; c < cells.size(); ++c)
982 {
983 const std::int32_t cell = cells[c];
984 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
985 x_dofmap, cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
986 for (int i = 0; i < num_dofs_g; ++i)
987 {
988 const int pos = 3 * x_dofs[i];
989 for (int j = 0; j < gdim; ++j)
990 coord_dofs(i, j) = x_g[pos + j];
991 }
992
993 // Compute J, detJ and K
994 std::ranges::fill(J_b, 0);
995 for (std::size_t p = 0; p < Xshape[0]; ++p)
996 {
997 auto _dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
998 dphi, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, p,
999 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1000 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1001 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1002 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1003 cmap.compute_jacobian(_dphi, coord_dofs, _J);
1004 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1005 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1006 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1007 cmap.compute_jacobian_inverse(_J, _K);
1008 detJ[p] = cmap.compute_jacobian_determinant(_J, det_scratch);
1009 }
1010
1011 std::span<const std::int32_t> dofs = dofmap->cell_dofs(cell);
1012 for (int k = 0; k < element_bs; ++k)
1013 {
1014 // Extract computed expression values for element block k
1015 for (int m = 0; m < value_size; ++m)
1016 {
1017 for (std::size_t k0 = 0; k0 < Xshape[0]; ++k0)
1018 {
1019 _vals(k0, 0, m)
1020 = f[f_shape1 * (k * value_size + m) + c * Xshape[0] + k0];
1021 }
1022 }
1023
1024 // Get element degrees of freedom for block
1025 for (std::size_t i = 0; i < Xshape[0]; ++i)
1026 {
1027 auto _u = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1028 _vals, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1029 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1030 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1031 ref_data, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1032 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1033 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1034 K, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1035 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1036 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1037 J, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1038 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1039 pull_back_fn(_U, _u, _K, 1.0 / detJ[i], _J);
1040 }
1041
1042 auto ref = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1043 ref_data, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0,
1044 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1045 impl::interpolation_apply(Pi, ref, std::span(_coeffs), element_bs);
1046 apply_inverse_transpose_dof_transformation(_coeffs, cell_info, cell, 1);
1047
1048 // Copy interpolation dofs into coefficient vector
1049 assert(_coeffs.size() == num_scalar_dofs);
1050 for (int i = 0; i < num_scalar_dofs; ++i)
1051 {
1052 const int dof = i * element_bs + k;
1053 std::div_t pos = std::div(dof, dofmap_bs);
1054 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
1055 }
1056 }
1057 }
1058 }
1059}
1060
1079template <std::floating_point T>
1081 const mesh::Geometry<T>& geometry0, const FiniteElement<T>& element0,
1082 const mesh::Mesh<T>& mesh1, std::span<const std::int32_t> cells, T padding)
1083{
1084 // Collect all the points at which values are needed to define the
1085 // interpolating function
1086 std::vector<T> coords = interpolation_coords(element0, geometry0, cells);
1087
1088 // Transpose interpolation coords
1089 std::vector<T> x(coords.size());
1090 std::size_t num_points = coords.size() / 3;
1091 for (std::size_t i = 0; i < num_points; ++i)
1092 for (std::size_t j = 0; j < 3; ++j)
1093 x[3 * i + j] = coords[i + j * num_points];
1094
1095 // Determine ownership of each point
1096 return geometry::determine_point_ownership<T>(mesh1, x, padding);
1097}
1098
1110template <dolfinx::scalar T, std::floating_point U>
1112 std::span<const std::int32_t> cells,
1113 const geometry::PointOwnershipData<U>& interpolation_data)
1114{
1115 auto mesh = u.function_space()->mesh();
1116 assert(mesh);
1117 MPI_Comm comm = mesh->comm();
1118 {
1119 auto mesh_v = v.function_space()->mesh();
1120 assert(mesh_v);
1121 int result;
1122 MPI_Comm_compare(comm, mesh_v->comm(), &result);
1123 if (result == MPI_UNEQUAL)
1124 {
1125 throw std::runtime_error("Interpolation on different meshes is only "
1126 "supported on the same communicator.");
1127 }
1128 }
1129
1130 assert(mesh->topology());
1131 auto cell_map = mesh->topology()->index_map(mesh->topology()->dim());
1132 assert(cell_map);
1133 auto element_u = u.function_space()->element();
1134 assert(element_u);
1135 const std::size_t value_size = u.function_space()->value_size();
1136
1137 const std::vector<int>& dest_ranks = interpolation_data.src_owner;
1138 const std::vector<int>& src_ranks = interpolation_data.dest_owners;
1139 const std::vector<U>& recv_points = interpolation_data.dest_points;
1140 const std::vector<std::int32_t>& evaluation_cells
1141 = interpolation_data.dest_cells;
1142
1143 // Evaluate the interpolating function where possible
1144 std::vector<T> send_values(recv_points.size() / 3 * value_size);
1145 v.eval(recv_points, {recv_points.size() / 3, (std::size_t)3},
1146 evaluation_cells, send_values, {recv_points.size() / 3, value_size});
1147
1148 using dextents2 = MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>;
1149
1150 // Send values back to owning process
1151 std::vector<T> values_b(dest_ranks.size() * value_size);
1152 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<const T, dextents2> _send_values(
1153 send_values.data(), src_ranks.size(), value_size);
1154 impl::scatter_values(comm, src_ranks, dest_ranks, _send_values,
1155 std::span(values_b));
1156
1157 // Transpose received data
1158 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<const T, dextents2> values(
1159 values_b.data(), dest_ranks.size(), value_size);
1160 std::vector<T> valuesT_b(value_size * dest_ranks.size());
1161 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<T, dextents2> valuesT(
1162 valuesT_b.data(), value_size, dest_ranks.size());
1163 for (std::size_t i = 0; i < values.extent(0); ++i)
1164 for (std::size_t j = 0; j < values.extent(1); ++j)
1165 valuesT(j, i) = values(i, j);
1166
1167 // Call local interpolation operator
1168 fem::interpolate<T>(u, valuesT_b, {valuesT.extent(0), valuesT.extent(1)},
1169 cells);
1170}
1171
1187template <dolfinx::scalar T, std::floating_point U>
1188void interpolate(Function<T, U>& u1, std::span<const std::int32_t> cells1,
1189 const Function<T, U>& u0, std::span<const std::int32_t> cells0)
1190{
1191 if (cells0.size() != cells1.size())
1192 throw std::runtime_error("Length of cell lists do not match.");
1193
1194 assert(u1.function_space());
1195 assert(u0.function_space());
1196 auto mesh = u1.function_space()->mesh();
1197 assert(mesh);
1198 assert(cells0.size() == cells1.size());
1199
1200 auto cell_map0 = mesh->topology()->index_map(mesh->topology()->dim());
1201 assert(cell_map0);
1202 std::size_t num_cells0 = cell_map0->size_local() + cell_map0->num_ghosts();
1203 if (u1.function_space() == u0.function_space()
1204 and cells1.size() == num_cells0)
1205 {
1206 // Same function spaces and on whole mesh
1207 std::span<T> u1_array = u1.x()->mutable_array();
1208 std::span<const T> u0_array = u0.x()->array();
1209 std::ranges::copy(u0_array, u1_array.begin());
1210 }
1211 else
1212 {
1213 // Get elements and check value shape
1214 auto fs0 = u0.function_space();
1215 auto element0 = fs0->element();
1216 assert(element0);
1217 auto fs1 = u1.function_space();
1218 auto element1 = fs1->element();
1219 assert(element1);
1220 if (fs0->value_shape().size() != fs1->value_shape().size()
1221 or !std::equal(fs0->value_shape().begin(), fs0->value_shape().end(),
1222 fs1->value_shape().begin()))
1223 {
1224 throw std::runtime_error(
1225 "Interpolation: elements have different value dimensions");
1226 }
1227
1228 if (element1 == element0 or *element1 == *element0)
1229 {
1230 // Same element, different dofmaps (or just a subset of cells)
1231 const int tdim = mesh->topology()->dim();
1232 auto cell_map1 = mesh->topology()->index_map(tdim);
1233 assert(cell_map1);
1234 assert(element1->block_size() == element0->block_size());
1235
1236 // Get dofmaps
1237 std::shared_ptr<const DofMap> dofmap0 = u0.function_space()->dofmap();
1238 assert(dofmap0);
1239 std::shared_ptr<const DofMap> dofmap1 = u1.function_space()->dofmap();
1240 assert(dofmap1);
1241
1242 std::span<T> u1_array = u1.x()->mutable_array();
1243 std::span<const T> u0_array = u0.x()->array();
1244
1245 // Iterate over mesh and interpolate on each cell
1246 const int bs0 = dofmap0->bs();
1247 const int bs1 = dofmap1->bs();
1248 for (std::size_t c = 0; c < cells1.size(); ++c)
1249 {
1250 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(cells0[c]);
1251 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(cells1[c]);
1252 assert(bs0 * dofs0.size() == bs1 * dofs1.size());
1253 for (std::size_t i = 0; i < dofs0.size(); ++i)
1254 {
1255 for (int k = 0; k < bs0; ++k)
1256 {
1257 int index = bs0 * i + k;
1258 std::div_t dv1 = std::div(index, bs1);
1259 u1_array[bs1 * dofs1[dv1.quot] + dv1.rem]
1260 = u0_array[bs0 * dofs0[i] + k];
1261 }
1262 }
1263 }
1264 }
1265 else if (element1->map_type() == element0->map_type())
1266 {
1267 // Different elements, same basis function map type
1268 impl::interpolate_same_map(u1, u0, cells1, cells0);
1269 }
1270 else
1271 {
1272 // Different elements with different maps for basis functions
1273 impl::interpolate_nonmatching_maps(u1, cells1, u0, cells0);
1274 }
1275 }
1276}
1277} // namespace dolfinx::fem
Degree-of-freedom map representations and tools.
Definition XDMFFile.h:27
static void compute_jacobian(const U &dphi, const V &cell_geometry, W &&J)
Definition CoordinateElement.h:126
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
static void compute_jacobian_inverse(const U &J, V &&K)
Compute the inverse of the Jacobian.
Definition CoordinateElement.h:135
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:206
static double compute_jacobian_determinant(const U &J, std::span< typename U::value_type > w)
Compute the determinant of the Jacobian.
Definition CoordinateElement.h:152
Model of a finite element.
Definition FiniteElement.h:38
Definition XDMFFile.h:29
std::shared_ptr< const FunctionSpace< geometry_type > > function_space() const
Access the function space.
Definition Function.h:143
std::shared_ptr< const la::Vector< value_type > > x() const
Underlying vector (const version).
Definition Function.h:149
Geometry stores the geometry imposed on a mesh.
Definition Geometry.h:34
const fem::CoordinateElement< value_type > & cmap() const
The element that describes the geometry map.
Definition Geometry.h:181
int dim() const
Return dimension of the Euclidean coordinate system.
Definition Geometry.h:117
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const std::int32_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents< std::size_t, 2 > > dofmap() const
DofMap for the geometry.
Definition Geometry.h:124
std::span< const value_type > x() const
Access geometry degrees-of-freedom data (const version).
Definition Geometry.h:169
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition Mesh.h:23
Definition interpolate.h:32
int rank(MPI_Comm comm)
Return process rank for the communicator.
Definition MPI.cpp:64
constexpr MPI_Datatype mpi_type()
MPI Type.
Definition MPI.h:273
Finite element method functionality.
Definition assemble_matrix_impl.h:26
geometry::PointOwnershipData< T > create_interpolation_data(const mesh::Geometry< T > &geometry0, const FiniteElement< T > &element0, const mesh::Mesh< T > &mesh1, std::span< const std::int32_t > cells, T padding)
Generate data needed to interpolate finite element Functions across different meshes.
Definition interpolate.h:1080
void interpolate(Function< T, U > &u, std::span< const T > f, std::array< std::size_t, 2 > fshape, std::span< const std::int32_t > cells)
Interpolate an evaluated expression f(x) in a finite element space.
Definition interpolate.h:703
@ inverse_transpose
Transpose inverse.
std::vector< T > interpolation_coords(const fem::FiniteElement< T > &element, const mesh::Geometry< T > &geometry, std::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.h:50
PointOwnershipData< T > determine_point_ownership(const mesh::Mesh< T > &mesh, std::span< const T > points, T padding)
Given a set of points, determine which process is colliding, using the GJK algorithm on cells to dete...
Definition utils.h:686
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