10#include "CoordinateElement.h"
12#include "FiniteElement.h"
13#include "FunctionSpace.h"
15#include <basix/mdspan.hpp>
17#include <dolfinx/common/IndexMap.h>
18#include <dolfinx/common/types.h>
19#include <dolfinx/geometry/utils.h>
20#include <dolfinx/mesh/Mesh.h>
29template <dolfinx::scalar T, std::
floating_po
int U>
33concept MDSpan =
requires(T x, std::size_t idx) {
35 { x.extent(0) } -> std::integral;
36 { x.extent(1) } -> std::integral;
50template <std::
floating_po
int T>
53 std::span<const std::int32_t> cells)
56 const std::size_t gdim = geometry.
dim();
57 auto x_dofmap = geometry.
dofmap();
58 std::span<const T> x_g = geometry.
x();
61 const std::size_t num_dofs_g = cmap.
dim();
64 const auto [X, Xshape] = element.interpolation_points();
67 std::array<std::size_t, 4> phi_shape = cmap.
tabulate_shape(0, Xshape[0]);
69 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
70 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
71 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>
72 phi_full(phi_b.data(), phi_shape);
74 auto phi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
75 phi_full, 0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
76 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
80 std::vector<T> coordinate_dofs(num_dofs_g * gdim, 0);
81 std::vector<T> x(3 * (cells.size() * Xshape[0]), 0);
82 for (std::size_t c = 0; c < cells.size(); ++c)
85 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
86 x_dofmap, cells[c], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
87 for (std::size_t i = 0; i < x_dofs.size(); ++i)
89 std::copy_n(std::next(x_g.begin(), 3 * x_dofs[i]), gdim,
90 std::next(coordinate_dofs.begin(), i * gdim));
94 for (std::size_t p = 0; p < Xshape[0]; ++p)
96 for (std::size_t j = 0; j < gdim; ++j)
99 for (std::size_t k = 0; k < num_dofs_g; ++k)
100 acc += phi(p, k) * coordinate_dofs[k * gdim + j];
101 x[j * (cells.size() * Xshape[0]) + c * Xshape[0] + p] = acc;
124template <dolfinx::scalar T, std::
floating_po
int U>
125void interpolate(Function<T, U>& u, std::span<const T> f,
126 std::array<std::size_t, 2> fshape,
127 std::span<const std::int32_t> cells);
132template <
typename T, std::
size_t D>
133using mdspan_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
134 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, D>>;
155template <dolfinx::scalar T>
156void scatter_values(MPI_Comm comm, std::span<const std::int32_t> src_ranks,
157 std::span<const std::int32_t> dest_ranks,
158 mdspan_t<const T, 2> send_values, std::span<T> recv_values)
160 const std::size_t block_size = send_values.extent(1);
161 assert(src_ranks.size() * block_size == send_values.size());
162 assert(recv_values.size() == dest_ranks.size() * block_size);
165 std::vector<std::int32_t> out_ranks(src_ranks.size());
166 out_ranks.assign(src_ranks.begin(), src_ranks.end());
167 auto [unique_end, range_end] = std::ranges::unique(out_ranks);
168 out_ranks.erase(unique_end, range_end);
169 out_ranks.reserve(out_ranks.size() + 1);
172 std::vector<std::int32_t> in_ranks;
173 in_ranks.reserve(dest_ranks.size());
174 std::copy_if(dest_ranks.begin(), dest_ranks.end(),
175 std::back_inserter(in_ranks),
176 [](
auto rank) { return rank >= 0; });
180 std::ranges::sort(in_ranks);
181 auto [unique_end, range_end] = std::ranges::unique(in_ranks);
182 in_ranks.erase(unique_end, range_end);
184 in_ranks.reserve(in_ranks.size() + 1);
187 MPI_Comm reverse_comm;
188 MPI_Dist_graph_create_adjacent(
189 comm, in_ranks.size(), in_ranks.data(), MPI_UNWEIGHTED, out_ranks.size(),
190 out_ranks.data(), MPI_UNWEIGHTED, MPI_INFO_NULL,
false, &reverse_comm);
192 std::vector<std::int32_t> comm_to_output;
193 std::vector<std::int32_t> recv_sizes(in_ranks.size());
194 recv_sizes.reserve(1);
195 std::vector<std::int32_t> recv_offsets(in_ranks.size() + 1, 0);
198 std::vector<std::pair<std::int32_t, std::int32_t>> rank_to_neighbor;
199 rank_to_neighbor.reserve(in_ranks.size());
200 for (std::size_t i = 0; i < in_ranks.size(); i++)
201 rank_to_neighbor.push_back({in_ranks[i], i});
202 std::ranges::sort(rank_to_neighbor);
205 std::ranges::for_each(
207 [&dest_ranks, &rank_to_neighbor, &recv_sizes, block_size](
auto rank)
211 auto it = std::ranges::lower_bound(rank_to_neighbor, rank,
213 [](
auto e) {
return e.first; });
214 assert(it != rank_to_neighbor.end() and it->first == rank);
215 recv_sizes[it->second] += block_size;
220 std::partial_sum(recv_sizes.begin(), recv_sizes.end(),
221 std::next(recv_offsets.begin(), 1));
224 comm_to_output.resize(recv_offsets.back() / block_size);
225 std::vector<std::int32_t> recv_counter(recv_sizes.size(), 0);
226 for (std::size_t i = 0; i < dest_ranks.size(); ++i)
228 if (
const std::int32_t rank = dest_ranks[i];
rank >= 0)
230 auto it = std::ranges::lower_bound(rank_to_neighbor, rank,
232 [](
auto e) {
return e.first; });
233 assert(it != rank_to_neighbor.end() and it->first == rank);
234 int insert_pos = recv_offsets[it->second] + recv_counter[it->second];
235 comm_to_output[insert_pos / block_size] = i * block_size;
236 recv_counter[it->second] += block_size;
241 std::vector<std::int32_t> send_sizes(out_ranks.size());
242 send_sizes.reserve(1);
247 std::vector<std::pair<std::int32_t, std::int32_t>> rank_to_neighbor;
248 rank_to_neighbor.reserve(out_ranks.size());
249 for (std::size_t i = 0; i < out_ranks.size(); i++)
250 rank_to_neighbor.push_back({out_ranks[i], i});
254 auto start = rank_to_neighbor.begin();
255 std::ranges::for_each(
257 [&rank_to_neighbor, &send_sizes, block_size, &start](
auto rank)
259 auto it = std::ranges::lower_bound(start, rank_to_neighbor.end(),
260 rank, std::ranges::less(),
261 [](
auto e) { return e.first; });
262 assert(it != rank_to_neighbor.end() and it->first == rank);
263 send_sizes[it->second] += block_size;
269 std::vector<std::int32_t> send_offsets(send_sizes.size() + 1, 0);
270 std::partial_sum(send_sizes.begin(), send_sizes.end(),
271 std::next(send_offsets.begin(), 1));
274 std::vector<T> values(recv_offsets.back());
276 MPI_Neighbor_alltoallv(send_values.data_handle(), send_sizes.data(),
278 values.data(), recv_sizes.data(), recv_offsets.data(),
280 MPI_Comm_free(&reverse_comm);
284 std::ranges::fill(recv_values, T(0));
285 for (std::size_t i = 0; i < comm_to_output.size(); i++)
287 auto vals = std::next(recv_values.begin(), comm_to_output[i]);
288 auto vals_from = std::next(values.begin(), i * block_size);
289 std::copy_n(vals_from, block_size, vals);
301template <MDSpan U, MDSpan V, dolfinx::scalar T>
302void interpolation_apply(U&& Pi, V&& data, std::span<T> coeffs,
int bs)
304 using X =
typename dolfinx::scalar_value_type_t<T>;
309 assert(data.extent(0) * data.extent(1) == Pi.extent(1));
310 for (std::size_t i = 0; i < Pi.extent(0); ++i)
313 for (std::size_t k = 0; k < data.extent(1); ++k)
314 for (std::size_t j = 0; j < data.extent(0); ++j)
316 +=
static_cast<X
>(Pi(i, k * data.extent(0) + j)) * data(j, k);
321 const std::size_t cols = Pi.extent(1);
322 assert(data.extent(0) == Pi.extent(1));
323 assert(data.extent(1) == bs);
324 for (
int k = 0; k < bs; ++k)
326 for (std::size_t i = 0; i < Pi.extent(0); ++i)
329 for (std::size_t j = 0; j < cols; ++j)
330 acc +=
static_cast<X
>(Pi(i, j)) * data(j, k);
331 coeffs[bs * i + k] = acc;
356template <dolfinx::scalar T, std::
floating_po
int U>
357void interpolate_same_map(Function<T, U>& u1,
const Function<T, U>& u0,
358 std::span<const std::int32_t> cells1,
359 std::span<const std::int32_t> cells0)
361 auto V0 = u0.function_space();
363 auto V1 = u1.function_space();
365 auto mesh0 = V0->mesh();
368 auto mesh1 = V1->mesh();
371 auto element0 = V0->element();
373 auto element1 = V1->element();
376 assert(mesh0->topology()->dim());
377 const int tdim = mesh0->topology()->dim();
378 auto map = mesh0->topology()->index_map(tdim);
380 std::span<T> u1_array = u1.x()->mutable_array();
381 std::span<const T> u0_array = u0.x()->array();
383 std::span<const std::uint32_t> cell_info0;
384 std::span<const std::uint32_t> cell_info1;
385 if (element1->needs_dof_transformations()
386 or element0->needs_dof_transformations())
388 mesh0->topology_mutable()->create_entity_permutations();
389 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
390 mesh1->topology_mutable()->create_entity_permutations();
391 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
395 auto dofmap1 = V1->dofmap();
396 auto dofmap0 = V0->dofmap();
399 const int bs1 = dofmap1->bs();
400 const int bs0 = dofmap0->bs();
401 auto apply_dof_transformation = element0->template dof_transformation_fn<T>(
403 auto apply_inverse_dof_transform
404 = element1->template dof_transformation_fn<T>(
408 std::vector<T> local0(element0->space_dimension());
409 std::vector<T> local1(element1->space_dimension());
412 const auto [i_m, im_shape]
413 = element1->create_interpolation_operator(*element0);
416 using X =
typename dolfinx::scalar_value_type_t<T>;
417 for (std::size_t c = 0; c < cells0.size(); c++)
420 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(cells0[c]);
421 for (std::size_t i = 0; i < dofs0.size(); ++i)
422 for (
int k = 0; k < bs0; ++k)
423 local0[bs0 * i + k] = u0_array[bs0 * dofs0[i] + k];
425 apply_dof_transformation(local0, cell_info0, cells0[c], 1);
429 std::ranges::fill(local1, 0);
430 for (std::size_t i = 0; i < im_shape[0]; ++i)
431 for (std::size_t j = 0; j < im_shape[1]; ++j)
432 local1[i] +=
static_cast<X
>(i_m[im_shape[1] * i + j]) * local0[j];
434 apply_inverse_dof_transform(local1, cell_info1, cells1[c], 1);
435 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(cells1[c]);
436 for (std::size_t i = 0; i < dofs1.size(); ++i)
437 for (
int k = 0; k < bs1; ++k)
438 u1_array[bs1 * dofs1[i] + k] = local1[bs1 * i + k];
456template <dolfinx::scalar T, std::
floating_po
int U>
457void interpolate_nonmatching_maps(Function<T, U>& u1,
458 std::span<const std::int32_t> cells1,
459 const Function<T, U>& u0,
460 std::span<const std::int32_t> cells0)
463 auto V0 = u0.function_space();
465 auto mesh0 = V0->mesh();
469 const int tdim = mesh0->topology()->dim();
470 const int gdim = mesh0->geometry().dim();
473 auto V1 = u1.function_space();
475 auto mesh1 = V1->mesh();
477 auto element0 = V0->element();
479 auto element1 = V1->element();
482 std::span<const std::uint32_t> cell_info0;
483 std::span<const std::uint32_t> cell_info1;
484 if (element1->needs_dof_transformations()
485 or element0->needs_dof_transformations())
487 mesh0->topology_mutable()->create_entity_permutations();
488 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
489 mesh1->topology_mutable()->create_entity_permutations();
490 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
494 auto dofmap0 = V0->dofmap();
495 auto dofmap1 = V1->dofmap();
497 const auto [X, Xshape] = element1->interpolation_points();
500 const int bs0 = element0->block_size();
501 const int bs1 = element1->block_size();
502 auto apply_dof_transformation0 = element0->template dof_transformation_fn<U>(
504 auto apply_inverse_dof_transform1
505 = element1->template dof_transformation_fn<T>(
509 const std::size_t dim0 = element0->space_dimension() / bs0;
510 const std::size_t value_size_ref0 = element0->reference_value_size();
511 const std::size_t value_size0 = V0->element()->reference_value_size();
513 const CoordinateElement<U>& cmap = mesh0->geometry().cmap();
514 auto x_dofmap = mesh0->geometry().dofmap();
515 const std::size_t num_dofs_g = cmap.dim();
516 std::span<const U> x_g = mesh0->geometry().x();
519 const std::array<std::size_t, 4> phi_shape
520 = cmap.tabulate_shape(1, Xshape[0]);
521 std::vector<U> phi_b(
522 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
523 mdspan_t<const U, 4> phi(phi_b.data(), phi_shape);
524 cmap.tabulate(1, X, Xshape, phi_b);
527 const auto [_basis_derivatives_reference0, b0shape]
528 = element0->tabulate(X, Xshape, 0);
529 mdspan_t<const U, 4> basis_derivatives_reference0(
530 _basis_derivatives_reference0.data(), b0shape);
533 std::vector<T> local1(element1->space_dimension());
534 std::vector<T> coeffs0(element0->space_dimension());
536 std::vector<U> basis0_b(Xshape[0] * dim0 * value_size0);
537 impl::mdspan_t<U, 3> basis0(basis0_b.data(), Xshape[0], dim0, value_size0);
539 std::vector<U> basis_reference0_b(Xshape[0] * dim0 * value_size_ref0);
540 impl::mdspan_t<U, 3> basis_reference0(basis_reference0_b.data(), Xshape[0],
541 dim0, value_size_ref0);
543 std::vector<T> values0_b(Xshape[0] * 1 * V1->element()->value_size());
544 impl::mdspan_t<T, 3> values0(values0_b.data(), Xshape[0], 1,
545 V1->element()->value_size());
547 std::vector<T> mapped_values_b(Xshape[0] * 1 * V1->element()->value_size());
548 impl::mdspan_t<T, 3> mapped_values0(mapped_values_b.data(), Xshape[0], 1,
549 V1->element()->value_size());
551 std::vector<U> coord_dofs_b(num_dofs_g * gdim);
552 impl::mdspan_t<U, 2> coord_dofs(coord_dofs_b.data(), num_dofs_g, gdim);
554 std::vector<U> J_b(Xshape[0] * gdim * tdim);
555 impl::mdspan_t<U, 3> J(J_b.data(), Xshape[0], gdim, tdim);
556 std::vector<U> K_b(Xshape[0] * tdim * gdim);
557 impl::mdspan_t<U, 3> K(K_b.data(), Xshape[0], tdim, gdim);
558 std::vector<U> detJ(Xshape[0]);
559 std::vector<U> det_scratch(2 * gdim * tdim);
562 const auto [_Pi_1, pi_shape] = element1->interpolation_operator();
563 impl::mdspan_t<const U, 2> Pi_1(_Pi_1.data(), pi_shape);
565 using u_t = impl::mdspan_t<U, 2>;
566 using U_t = impl::mdspan_t<const U, 2>;
567 using J_t = impl::mdspan_t<const U, 2>;
568 using K_t = impl::mdspan_t<const U, 2>;
569 auto push_forward_fn0
570 = element0->basix_element().template map_fn<u_t, U_t, J_t, K_t>();
572 using v_t = impl::mdspan_t<const T, 2>;
573 using V_t = impl::mdspan_t<T, 2>;
575 = element1->basix_element().template map_fn<V_t, v_t, K_t, J_t>();
578 std::span<const T> array0 = u0.x()->array();
579 std::span<T> array1 = u1.x()->mutable_array();
580 for (std::size_t c = 0; c < cells0.size(); c++)
583 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
584 x_dofmap, cells0[c], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
585 for (std::size_t i = 0; i < num_dofs_g; ++i)
587 const int pos = 3 * x_dofs[i];
588 for (std::size_t j = 0; j < gdim; ++j)
589 coord_dofs(i, j) = x_g[pos + j];
593 std::ranges::fill(J_b, 0);
594 for (std::size_t p = 0; p < Xshape[0]; ++p)
596 auto dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
597 phi, std::pair(1, tdim + 1), p,
598 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
600 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
601 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
602 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
603 cmap.compute_jacobian(dphi, coord_dofs, _J);
604 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
605 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
606 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
607 cmap.compute_jacobian_inverse(_J, _K);
608 detJ[p] = cmap.compute_jacobian_determinant(_J, det_scratch);
613 for (std::size_t k0 = 0; k0 < basis_reference0.extent(0); ++k0)
614 for (std::size_t k1 = 0; k1 < basis_reference0.extent(1); ++k1)
615 for (std::size_t k2 = 0; k2 < basis_reference0.extent(2); ++k2)
616 basis_reference0(k0, k1, k2)
617 = basis_derivatives_reference0(0, k0, k1, k2);
619 for (std::size_t p = 0; p < Xshape[0]; ++p)
621 apply_dof_transformation0(
622 std::span(basis_reference0_b.data() + p * dim0 * value_size_ref0,
623 dim0 * value_size_ref0),
624 cell_info0, cells0[c], value_size_ref0);
627 for (std::size_t i = 0; i < basis0.extent(0); ++i)
629 auto _u = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
630 basis0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
631 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
632 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
633 basis_reference0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
634 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
635 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
636 K, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
637 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
638 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
639 J, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
640 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
641 push_forward_fn0(_u, _U, _J, detJ[i], _K);
645 const int dof_bs0 = dofmap0->bs();
646 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(cells0[c]);
647 for (std::size_t i = 0; i < dofs0.size(); ++i)
648 for (
int k = 0; k < dof_bs0; ++k)
649 coeffs0[dof_bs0 * i + k] = array0[dof_bs0 * dofs0[i] + k];
652 using X =
typename dolfinx::scalar_value_type_t<T>;
653 for (std::size_t p = 0; p < Xshape[0]; ++p)
655 for (
int k = 0; k < bs0; ++k)
657 for (std::size_t j = 0; j < value_size0; ++j)
660 for (std::size_t i = 0; i < dim0; ++i)
661 acc += coeffs0[bs0 * i + k] *
static_cast<X
>(basis0(p, i, j));
662 values0(p, 0, j * bs0 + k) = acc;
668 for (std::size_t i = 0; i < values0.extent(0); ++i)
670 auto _u = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
671 values0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
672 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
673 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
674 mapped_values0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
675 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
676 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
677 K, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
678 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
679 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
680 J, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
681 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
682 pull_back_fn1(_U, _u, _K, 1.0 / detJ[i], _J);
685 auto values = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
686 mapped_values0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0,
687 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
688 interpolation_apply(Pi_1, values, std::span(local1), bs1);
689 apply_inverse_dof_transform1(local1, cell_info1, cells1[c], 1);
692 const int dof_bs1 = dofmap1->bs();
693 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(c);
694 for (std::size_t i = 0; i < dofs1.size(); ++i)
695 for (
int k = 0; k < dof_bs1; ++k)
696 array1[dof_bs1 * dofs1[i] + k] = local1[dof_bs1 * i + k];
703template <dolfinx::scalar T, std::
floating_po
int U>
705 std::array<std::size_t, 2> fshape,
706 std::span<const std::int32_t> cells)
708 using cmdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
709 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
710 using cmdspan4_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
711 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>;
712 using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
713 U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
714 using mdspan3_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
715 U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>;
716 auto element = u.function_space()->element();
718 const int element_bs = element->block_size();
719 if (
int num_sub = element->num_sub_elements();
720 num_sub > 0 and num_sub != element_bs)
722 throw std::runtime_error(
"Cannot directly interpolate a mixed space. "
723 "Interpolate into subspaces.");
727 assert(u.function_space());
728 auto mesh = u.function_space()->mesh();
731 const int gdim = mesh->geometry().dim();
732 const int tdim = mesh->topology()->dim();
734 if (fshape[0] != (std::size_t)u.function_space()->element()->value_size())
735 throw std::runtime_error(
"Interpolation data has the wrong shape/size.");
737 std::span<const std::uint32_t> cell_info;
738 if (element->needs_dof_transformations())
740 mesh->topology_mutable()->create_entity_permutations();
741 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
744 const std::size_t f_shape1
745 = f.size() / u.function_space()->element()->value_size();
746 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
747 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
748 _f(f.data(), fshape);
751 const auto dofmap = u.function_space()->dofmap();
753 const int dofmap_bs = dofmap->bs();
756 const int num_scalar_dofs = element->space_dimension() / element_bs;
757 const int value_size = u.function_space()->element()->reference_value_size();
759 std::span<T> coeffs = u.x()->mutable_array();
760 std::vector<T> _coeffs(num_scalar_dofs);
764 const bool symmetric = u.function_space()->symmetric();
765 if (element->map_ident() && element->interpolation_ident())
770 auto apply_inv_transpose_dof_transformation
771 = element->template dof_transformation_fn<T>(
776 std::size_t matrix_size = 0;
777 while (matrix_size * matrix_size < fshape[0])
781 for (std::size_t c = 0; c < cells.size(); ++c)
792 std::size_t rowstart = 0;
793 const std::int32_t
cell = cells[c];
794 std::span<const std::int32_t> dofs = dofmap->cell_dofs(
cell);
795 for (
int k = 0; k < element_bs; ++k)
797 if (k - rowstart > row)
806 std::next(f.begin(), (row * matrix_size + k - rowstart) * f_shape1
807 + c * num_scalar_dofs),
808 num_scalar_dofs, _coeffs.begin());
809 apply_inv_transpose_dof_transformation(_coeffs, cell_info,
cell, 1);
810 for (
int i = 0; i < num_scalar_dofs; ++i)
812 const int dof = i * element_bs + k;
813 std::div_t pos = std::div(dof, dofmap_bs);
814 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
822 for (std::size_t c = 0; c < cells.size(); ++c)
824 const std::int32_t
cell = cells[c];
825 std::span<const std::int32_t> dofs = dofmap->cell_dofs(
cell);
826 for (
int k = 0; k < element_bs; ++k)
830 std::copy_n(std::next(f.begin(), k * f_shape1 + c * num_scalar_dofs),
831 num_scalar_dofs, _coeffs.begin());
832 apply_inv_transpose_dof_transformation(_coeffs, cell_info,
cell, 1);
833 for (
int i = 0; i < num_scalar_dofs; ++i)
835 const int dof = i * element_bs + k;
836 std::div_t pos = std::div(dof, dofmap_bs);
837 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
843 else if (element->map_ident())
850 throw std::runtime_error(
851 "Interpolation into this element not supported.");
855 = u.function_space()->element()->reference_value_size();
857 if (element_vs > 1 and element_bs > 1)
859 throw std::runtime_error(
860 "Interpolation into this element not supported.");
864 const auto [_Pi, pi_shape] = element->interpolation_operator();
865 cmdspan2_t Pi(_Pi.data(), pi_shape);
866 const std::size_t num_interp_points = Pi.extent(1);
867 assert(Pi.extent(0) == num_scalar_dofs);
869 auto apply_inv_transpose_dof_transformation
870 = element->template dof_transformation_fn<T>(
874 std::vector<T> ref_data_b(num_interp_points);
875 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
876 T, MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
877 std::size_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent, 1>>
878 ref_data(ref_data_b.data(), num_interp_points, 1);
879 for (std::size_t c = 0; c < cells.size(); ++c)
881 const std::int32_t
cell = cells[c];
882 std::span<const std::int32_t> dofs = dofmap->cell_dofs(
cell);
883 for (
int k = 0; k < element_bs; ++k)
885 for (
int i = 0; i < element_vs; ++i)
888 std::next(f.begin(), (i + k) * f_shape1
889 + c * num_interp_points / element_vs),
890 num_interp_points / element_vs,
891 std::next(ref_data_b.begin(),
892 i * num_interp_points / element_vs));
894 impl::interpolation_apply(Pi, ref_data, std::span(_coeffs), 1);
895 apply_inv_transpose_dof_transformation(_coeffs, cell_info,
cell, 1);
896 for (
int i = 0; i < num_scalar_dofs; ++i)
898 const int dof = i * element_bs + k;
899 std::div_t pos = std::div(dof, dofmap_bs);
900 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
909 throw std::runtime_error(
910 "Interpolation into this element not supported.");
913 const auto [X, Xshape] = element->interpolation_points();
916 throw std::runtime_error(
917 "Interpolation into this space is not yet supported.");
920 if (_f.extent(1) != cells.size() * Xshape[0])
921 throw std::runtime_error(
"Interpolation data has the wrong shape.");
927 auto x_dofmap = mesh->geometry().dofmap();
928 const int num_dofs_g = cmap.
dim();
929 std::span<const U> x_g = mesh->geometry().x();
932 std::vector<U> J_b(Xshape[0] * gdim * tdim);
933 mdspan3_t J(J_b.data(), Xshape[0], gdim, tdim);
934 std::vector<U> K_b(Xshape[0] * tdim * gdim);
935 mdspan3_t K(K_b.data(), Xshape[0], tdim, gdim);
936 std::vector<U> detJ(Xshape[0]);
937 std::vector<U> det_scratch(2 * gdim * tdim);
939 std::vector<U> coord_dofs_b(num_dofs_g * gdim);
940 mdspan2_t coord_dofs(coord_dofs_b.data(), num_dofs_g, gdim);
941 const std::size_t value_size_ref = element->reference_value_size();
942 std::vector<T> ref_data_b(Xshape[0] * 1 * value_size_ref);
943 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
944 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
945 ref_data(ref_data_b.data(), Xshape[0], 1, value_size_ref);
947 std::vector<T> _vals_b(Xshape[0] * 1 * value_size);
948 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
949 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
950 _vals(_vals_b.data(), Xshape[0], 1, value_size);
954 std::array<std::size_t, 4> phi_shape = cmap.
tabulate_shape(1, Xshape[0]);
955 std::vector<U> phi_b(
956 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
957 cmdspan4_t phi(phi_b.data(), phi_shape);
959 auto dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
960 phi, std::pair(1, tdim + 1),
961 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
962 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
964 const std::function<void(std::span<T>, std::span<const std::uint32_t>,
966 apply_inverse_transpose_dof_transformation
967 = element->template dof_transformation_fn<T>(
971 const auto [_Pi, pi_shape] = element->interpolation_operator();
972 cmdspan2_t Pi(_Pi.data(), pi_shape);
974 using u_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
975 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
976 using U_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
977 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
978 using J_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
979 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
980 using K_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
981 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
983 = element->basix_element().template map_fn<U_t, u_t, J_t, K_t>();
985 for (std::size_t c = 0; c < cells.size(); ++c)
987 const std::int32_t
cell = cells[c];
988 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
989 x_dofmap,
cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
990 for (
int i = 0; i < num_dofs_g; ++i)
992 const int pos = 3 * x_dofs[i];
993 for (
int j = 0; j < gdim; ++j)
994 coord_dofs(i, j) = x_g[pos + j];
998 std::ranges::fill(J_b, 0);
999 for (std::size_t p = 0; p < Xshape[0]; ++p)
1001 auto _dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1002 dphi, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, p,
1003 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1004 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1005 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1006 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1008 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1009 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1010 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1015 std::span<const std::int32_t> dofs = dofmap->cell_dofs(
cell);
1016 for (
int k = 0; k < element_bs; ++k)
1019 for (
int m = 0; m < value_size; ++m)
1021 for (std::size_t k0 = 0; k0 < Xshape[0]; ++k0)
1024 = f[f_shape1 * (k * value_size + m) + c * Xshape[0] + k0];
1029 for (std::size_t i = 0; i < Xshape[0]; ++i)
1031 auto _u = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1032 _vals, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1033 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1034 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1035 ref_data, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1036 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1037 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1038 K, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1039 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1040 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1041 J, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1042 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1043 pull_back_fn(_U, _u, _K, 1.0 / detJ[i], _J);
1046 auto ref = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1047 ref_data, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0,
1048 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1049 impl::interpolation_apply(Pi, ref, std::span(_coeffs), element_bs);
1050 apply_inverse_transpose_dof_transformation(_coeffs, cell_info,
cell, 1);
1053 assert(_coeffs.size() == num_scalar_dofs);
1054 for (
int i = 0; i < num_scalar_dofs; ++i)
1056 const int dof = i * element_bs + k;
1057 std::div_t pos = std::div(dof, dofmap_bs);
1058 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
1083template <std::
floating_po
int T>
1086 const mesh::Mesh<T>& mesh1, std::span<const std::int32_t> cells, T padding)
1093 std::vector<T> x(coords.size());
1094 std::size_t num_points = coords.size() / 3;
1095 for (std::size_t i = 0; i < num_points; ++i)
1096 for (std::size_t j = 0; j < 3; ++j)
1097 x[3 * i + j] = coords[i + j * num_points];
1114template <dolfinx::scalar T, std::
floating_po
int U>
1116 std::span<const std::int32_t> cells,
1119 auto mesh = u.function_space()->mesh();
1121 MPI_Comm comm = mesh->comm();
1123 auto mesh_v = v.function_space()->mesh();
1126 MPI_Comm_compare(comm, mesh_v->comm(), &result);
1127 if (result == MPI_UNEQUAL)
1129 throw std::runtime_error(
"Interpolation on different meshes is only "
1130 "supported on the same communicator.");
1134 assert(mesh->topology());
1135 auto cell_map = mesh->topology()->index_map(mesh->topology()->dim());
1137 auto element_u = u.function_space()->element();
1139 const std::size_t value_size = u.function_space()->element()->value_size();
1141 const std::vector<int>& dest_ranks = interpolation_data.
src_owner;
1142 const std::vector<int>& src_ranks = interpolation_data.
dest_owners;
1143 const std::vector<U>& recv_points = interpolation_data.
dest_points;
1144 const std::vector<std::int32_t>& evaluation_cells
1148 std::vector<T> send_values(recv_points.size() / 3 * value_size);
1149 v.eval(recv_points, {recv_points.size() / 3, (std::size_t)3},
1150 evaluation_cells, send_values, {recv_points.size() / 3, value_size});
1152 using dextents2 = MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>;
1155 std::vector<T> values_b(dest_ranks.size() * value_size);
1156 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<const T, dextents2> _send_values(
1157 send_values.data(), src_ranks.size(), value_size);
1158 impl::scatter_values(comm, src_ranks, dest_ranks, _send_values,
1159 std::span(values_b));
1162 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<const T, dextents2> values(
1163 values_b.data(), dest_ranks.size(), value_size);
1164 std::vector<T> valuesT_b(value_size * dest_ranks.size());
1165 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<T, dextents2> valuesT(
1166 valuesT_b.data(), value_size, dest_ranks.size());
1167 for (std::size_t i = 0; i < values.extent(0); ++i)
1168 for (std::size_t j = 0; j < values.extent(1); ++j)
1169 valuesT(j, i) = values(i, j);
1191template <dolfinx::scalar T, std::
floating_po
int U>
1195 if (cells0.size() != cells1.size())
1196 throw std::runtime_error(
"Length of cell lists do not match.");
1199 assert(u0.function_space());
1202 assert(cells0.size() == cells1.size());
1204 auto cell_map0 = mesh->topology()->index_map(mesh->topology()->dim());
1206 std::size_t num_cells0 = cell_map0->size_local() + cell_map0->num_ghosts();
1208 and cells1.size() == num_cells0)
1211 std::span<T> u1_array = u1.
x()->mutable_array();
1212 std::span<const T> u0_array = u0.x()->array();
1213 std::ranges::copy(u0_array, u1_array.begin());
1218 auto fs0 = u0.function_space();
1219 auto element0 = fs0->element();
1222 auto element1 = fs1->element();
1224 if (!std::ranges::equal(fs0->element()->value_shape(),
1225 fs1->element()->value_shape()))
1227 throw std::runtime_error(
1228 "Interpolation: elements have different value dimensions");
1231 if (element1 == element0 or *element1 == *element0)
1234 const int tdim = mesh->topology()->dim();
1235 auto cell_map1 = mesh->topology()->index_map(tdim);
1237 assert(element1->block_size() == element0->block_size());
1240 std::shared_ptr<const DofMap> dofmap0 = u0.function_space()->dofmap();
1242 std::shared_ptr<const DofMap> dofmap1 = u1.
function_space()->dofmap();
1245 std::span<T> u1_array = u1.
x()->mutable_array();
1246 std::span<const T> u0_array = u0.x()->array();
1249 const int bs0 = dofmap0->bs();
1250 const int bs1 = dofmap1->bs();
1251 for (std::size_t c = 0; c < cells1.size(); ++c)
1253 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(cells0[c]);
1254 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(cells1[c]);
1255 assert(bs0 * dofs0.size() == bs1 * dofs1.size());
1256 for (std::size_t i = 0; i < dofs0.size(); ++i)
1258 for (
int k = 0; k < bs0; ++k)
1260 int index = bs0 * i + k;
1261 std::div_t dv1 = std::div(index, bs1);
1262 u1_array[bs1 * dofs1[dv1.quot] + dv1.rem]
1263 = u0_array[bs0 * dofs0[i] + k];
1268 else if (element1->map_type() == element0->map_type())
1271 impl::interpolate_same_map(u1, u0, cells1, cells0);
1276 impl::interpolate_nonmatching_maps(u1, cells1, u0, cells0);
Degree-of-freedom map representations and tools.
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:57
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:33
MPI_Datatype mpi_t
Retrieves the MPI data type associated to the provided type.
Definition MPI.h:282
int rank(MPI_Comm comm)
Return process rank for the communicator.
Definition MPI.cpp:64
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:1084
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:704
@ 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:51
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