10#include "CoordinateElement.h"
12#include "FiniteElement.h"
13#include "FunctionSpace.h"
14#include <basix/mdspan.hpp>
16#include <dolfinx/common/IndexMap.h>
17#include <dolfinx/common/types.h>
18#include <dolfinx/geometry/BoundingBoxTree.h>
19#include <dolfinx/geometry/utils.h>
20#include <dolfinx/mesh/Mesh.h>
28template <dolfinx::scalar T, std::
floating_po
int U>
32concept MDSpan =
requires(T x, std::size_t idx) {
34 { x.extent(0) } -> std::integral;
36 { x.extent(1) } -> std::integral;
49template <std::
floating_po
int T>
52 std::span<const std::int32_t> cells)
55 const std::size_t gdim = geometry.
dim();
56 auto x_dofmap = geometry.
dofmap();
57 std::span<const T> x_g = geometry.
x();
60 const std::size_t num_dofs_g = cmap.
dim();
63 const auto [X, Xshape] = element.interpolation_points();
66 std::array<std::size_t, 4> phi_shape = cmap.
tabulate_shape(0, Xshape[0]);
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);
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);
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)
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)
88 std::copy_n(std::next(x_g.begin(), 3 * x_dofs[i]), gdim,
89 std::next(coordinate_dofs.begin(), i * gdim));
93 for (std::size_t p = 0; p < Xshape[0]; ++p)
95 for (std::size_t j = 0; j < gdim; ++j)
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;
122template <dolfinx::scalar T, std::
floating_po
int U>
123void interpolate(Function<T, U>& u, std::span<const T> f,
124 std::array<std::size_t, 2> fshape,
125 std::span<const std::int32_t> cells);
130template <
typename T, std::
size_t D>
131using mdspan_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
132 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, D>>;
153template <dolfinx::scalar T>
155 MPI_Comm comm, std::span<const std::int32_t> src_ranks,
156 std::span<const std::int32_t> dest_ranks,
157 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
158 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
160 std::span<T> recv_values)
162 const std::size_t block_size = send_values.extent(1);
163 assert(src_ranks.size() * block_size == send_values.size());
164 assert(recv_values.size() == dest_ranks.size() * block_size);
167 std::vector<std::int32_t> out_ranks(src_ranks.size());
168 out_ranks.assign(src_ranks.begin(), src_ranks.end());
169 out_ranks.erase(std::unique(out_ranks.begin(), out_ranks.end()),
171 out_ranks.reserve(out_ranks.size() + 1);
174 std::vector<std::int32_t> in_ranks;
175 in_ranks.reserve(dest_ranks.size());
176 std::copy_if(dest_ranks.begin(), dest_ranks.end(),
177 std::back_inserter(in_ranks),
178 [](
auto rank) { return rank >= 0; });
181 std::sort(in_ranks.begin(), in_ranks.end());
182 in_ranks.erase(std::unique(in_ranks.begin(), in_ranks.end()), in_ranks.end());
183 in_ranks.reserve(in_ranks.size() + 1);
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);
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);
197 std::map<std::int32_t, std::int32_t> rank_to_neighbor;
198 for (std::size_t i = 0; i < in_ranks.size(); i++)
199 rank_to_neighbor[in_ranks[i]] = i;
203 dest_ranks.begin(), dest_ranks.end(),
204 [&dest_ranks, &rank_to_neighbor, &recv_sizes, block_size](
auto rank)
208 const int neighbor = rank_to_neighbor[rank];
209 recv_sizes[neighbor] += block_size;
214 std::partial_sum(recv_sizes.begin(), recv_sizes.end(),
215 std::next(recv_offsets.begin(), 1));
218 comm_to_output.resize(recv_offsets.back() / block_size);
219 std::vector<std::int32_t> recv_counter(recv_sizes.size(), 0);
220 for (std::size_t i = 0; i < dest_ranks.size(); ++i)
222 if (
const std::int32_t rank = dest_ranks[i]; rank >= 0)
224 const int neighbor = rank_to_neighbor[rank];
225 int insert_pos = recv_offsets[neighbor] + recv_counter[neighbor];
226 comm_to_output[insert_pos / block_size] = i * block_size;
227 recv_counter[neighbor] += block_size;
232 std::vector<std::int32_t> send_sizes(out_ranks.size());
233 send_sizes.reserve(1);
236 std::map<std::int32_t, std::int32_t> rank_to_neighbor;
237 for (std::size_t i = 0; i < out_ranks.size(); i++)
238 rank_to_neighbor[out_ranks[i]] = i;
241 std::for_each(src_ranks.begin(), src_ranks.end(),
242 [&rank_to_neighbor, &send_sizes, block_size](
auto rank)
243 { send_sizes[rank_to_neighbor[rank]] += block_size; });
247 std::vector<std::int32_t> send_offsets(send_sizes.size() + 1, 0);
248 std::partial_sum(send_sizes.begin(), send_sizes.end(),
249 std::next(send_offsets.begin(), 1));
252 std::vector<T> values(recv_offsets.back());
254 MPI_Neighbor_alltoallv(send_values.data_handle(), send_sizes.data(),
255 send_offsets.data(), dolfinx::MPI::mpi_type<T>(),
256 values.data(), recv_sizes.data(), recv_offsets.data(),
257 dolfinx::MPI::mpi_type<T>(), reverse_comm);
258 MPI_Comm_free(&reverse_comm);
261 std::fill(recv_values.begin(), recv_values.end(), T(0));
262 for (std::size_t i = 0; i < comm_to_output.size(); i++)
264 auto vals = std::next(recv_values.begin(), comm_to_output[i]);
265 auto vals_from = std::next(values.begin(), i * block_size);
266 std::copy_n(vals_from, block_size, vals);
278template <MDSpan U, MDSpan V, dolfinx::scalar T>
279void interpolation_apply(U&& Pi, V&& data, std::span<T> coeffs,
int bs)
281 using X =
typename dolfinx::scalar_value_type_t<T>;
286 assert(data.extent(0) * data.extent(1) == Pi.extent(1));
287 for (std::size_t i = 0; i < Pi.extent(0); ++i)
290 for (std::size_t k = 0; k < data.extent(1); ++k)
291 for (std::size_t j = 0; j < data.extent(0); ++j)
293 +=
static_cast<X
>(Pi(i, k * data.extent(0) + j)) * data(j, k);
298 const std::size_t cols = Pi.extent(1);
299 assert(data.extent(0) == Pi.extent(1));
300 assert(data.extent(1) == bs);
301 for (
int k = 0; k < bs; ++k)
303 for (std::size_t i = 0; i < Pi.extent(0); ++i)
306 for (std::size_t j = 0; j < cols; ++j)
307 acc +=
static_cast<X
>(Pi(i, j)) * data(j, k);
308 coeffs[bs * i + k] = acc;
325template <dolfinx::scalar T, std::
floating_po
int U>
326void interpolate_same_map(Function<T, U>& u1,
const Function<T, U>& u0,
327 std::span<const std::int32_t> cells1,
328 std::span<const std::int32_t> cells0)
330 auto V0 = u0.function_space();
332 auto V1 = u1.function_space();
334 auto mesh0 = V0->mesh();
337 auto mesh1 = V1->mesh();
340 auto element0 = V0->element();
342 auto element1 = V1->element();
345 assert(mesh0->topology()->dim());
346 const int tdim = mesh0->topology()->dim();
347 auto map = mesh0->topology()->index_map(tdim);
349 std::span<T> u1_array = u1.x()->mutable_array();
350 std::span<const T> u0_array = u0.x()->array();
352 std::span<const std::uint32_t> cell_info0;
353 std::span<const std::uint32_t> cell_info1;
354 if (element1->needs_dof_transformations()
355 or element0->needs_dof_transformations())
357 mesh0->topology_mutable()->create_entity_permutations();
358 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
359 mesh1->topology_mutable()->create_entity_permutations();
360 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
364 auto dofmap1 = V1->dofmap();
365 auto dofmap0 = V0->dofmap();
368 const int bs1 = dofmap1->bs();
369 const int bs0 = dofmap0->bs();
370 auto apply_dof_transformation = element0->template dof_transformation_fn<T>(
372 auto apply_inverse_dof_transform
373 = element1->template dof_transformation_fn<T>(
377 std::vector<T> local0(element0->space_dimension());
378 std::vector<T> local1(element1->space_dimension());
381 const auto [i_m, im_shape]
382 = element1->create_interpolation_operator(*element0);
385 using X =
typename dolfinx::scalar_value_type_t<T>;
386 for (std::size_t c = 0; c < cells0.size(); c++)
389 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(cells0[c]);
390 for (std::size_t i = 0; i < dofs0.size(); ++i)
391 for (
int k = 0; k < bs0; ++k)
392 local0[bs0 * i + k] = u0_array[bs0 * dofs0[i] + k];
394 apply_dof_transformation(local0, cell_info0, cells0[c], 1);
398 std::fill(local1.begin(), local1.end(), 0);
399 for (std::size_t i = 0; i < im_shape[0]; ++i)
400 for (std::size_t j = 0; j < im_shape[1]; ++j)
401 local1[i] +=
static_cast<X
>(i_m[im_shape[1] * i + j]) * local0[j];
403 apply_inverse_dof_transform(local1, cell_info1, cells1[c], 1);
404 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(cells1[c]);
405 for (std::size_t i = 0; i < dofs1.size(); ++i)
406 for (
int k = 0; k < bs1; ++k)
407 u1_array[bs1 * dofs1[i] + k] = local1[bs1 * i + k];
422template <dolfinx::scalar T, std::
floating_po
int U>
423void interpolate_nonmatching_maps(Function<T, U>& u1,
const Function<T, U>& u0,
424 std::span<const std::int32_t> cells1,
425 std::span<const std::int32_t> cells0)
428 auto V0 = u0.function_space();
430 auto mesh0 = V0->mesh();
434 const int tdim = mesh0->topology()->dim();
435 const int gdim = mesh0->geometry().dim();
438 auto V1 = u1.function_space();
440 auto mesh1 = V1->mesh();
442 auto element0 = V0->element();
444 auto element1 = V1->element();
447 std::span<const std::uint32_t> cell_info0;
448 std::span<const std::uint32_t> cell_info1;
450 if (element1->needs_dof_transformations()
451 or element0->needs_dof_transformations())
453 mesh0->topology_mutable()->create_entity_permutations();
454 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
455 mesh1->topology_mutable()->create_entity_permutations();
456 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
460 auto dofmap0 = V0->dofmap();
461 auto dofmap1 = V1->dofmap();
463 const auto [X, Xshape] = element1->interpolation_points();
466 const int bs0 = element0->block_size();
467 const int bs1 = element1->block_size();
468 auto apply_dof_transformation0 = element0->template dof_transformation_fn<U>(
470 auto apply_inverse_dof_transform1
471 = element1->template dof_transformation_fn<T>(
475 const std::size_t dim0 = element0->space_dimension() / bs0;
476 const std::size_t value_size_ref0 = element0->reference_value_size() / bs0;
477 const std::size_t value_size0 = V0->value_size() / bs0;
479 const CoordinateElement<U>& cmap = mesh0->geometry().cmap();
480 auto x_dofmap = mesh0->geometry().dofmap();
481 const std::size_t num_dofs_g = cmap.dim();
482 std::span<const U> x_g = mesh0->geometry().x();
485 const std::array<std::size_t, 4> phi_shape
486 = cmap.tabulate_shape(1, Xshape[0]);
487 std::vector<U> phi_b(
488 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
489 mdspan_t<const U, 4> phi(phi_b.data(), phi_shape);
490 cmap.tabulate(1, X, Xshape, phi_b);
493 const auto [_basis_derivatives_reference0, b0shape]
494 = element0->tabulate(X, Xshape, 0);
495 mdspan_t<const U, 4> basis_derivatives_reference0(
496 _basis_derivatives_reference0.data(), b0shape);
499 std::vector<T> local1(element1->space_dimension());
500 std::vector<T> coeffs0(element0->space_dimension());
502 std::vector<U> basis0_b(Xshape[0] * dim0 * value_size0);
503 impl::mdspan_t<U, 3> basis0(basis0_b.data(), Xshape[0], dim0, value_size0);
505 std::vector<U> basis_reference0_b(Xshape[0] * dim0 * value_size_ref0);
506 impl::mdspan_t<U, 3> basis_reference0(basis_reference0_b.data(), Xshape[0],
507 dim0, value_size_ref0);
509 std::vector<T> values0_b(Xshape[0] * 1 * V1->value_size());
510 impl::mdspan_t<T, 3> values0(values0_b.data(), Xshape[0], 1,
513 std::vector<T> mapped_values_b(Xshape[0] * 1 * V1->value_size());
514 impl::mdspan_t<T, 3> mapped_values0(mapped_values_b.data(), Xshape[0], 1,
517 std::vector<U> coord_dofs_b(num_dofs_g * gdim);
518 impl::mdspan_t<U, 2> coord_dofs(coord_dofs_b.data(), num_dofs_g, gdim);
520 std::vector<U> J_b(Xshape[0] * gdim * tdim);
521 impl::mdspan_t<U, 3> J(J_b.data(), Xshape[0], gdim, tdim);
522 std::vector<U> K_b(Xshape[0] * tdim * gdim);
523 impl::mdspan_t<U, 3> K(K_b.data(), Xshape[0], tdim, gdim);
524 std::vector<U> detJ(Xshape[0]);
525 std::vector<U> det_scratch(2 * gdim * tdim);
528 const auto [_Pi_1, pi_shape] = element1->interpolation_operator();
529 impl::mdspan_t<const U, 2> Pi_1(_Pi_1.data(), pi_shape);
531 using u_t = impl::mdspan_t<U, 2>;
532 using U_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
533 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
534 using J_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
535 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
536 using K_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
537 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
538 auto push_forward_fn0
539 = element0->basix_element().template map_fn<u_t, U_t, J_t, K_t>();
541 using v_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
542 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
543 using V_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
544 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
546 = element1->basix_element().template map_fn<V_t, v_t, K_t, J_t>();
549 std::span<const T> array0 = u0.x()->array();
550 std::span<T> array1 = u1.x()->mutable_array();
551 for (std::size_t c = 0; c < cells0.size(); c++)
554 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
555 x_dofmap, cells0[c], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
556 for (std::size_t i = 0; i < num_dofs_g; ++i)
558 const int pos = 3 * x_dofs[i];
559 for (std::size_t j = 0; j < gdim; ++j)
560 coord_dofs(i, j) = x_g[pos + j];
564 std::fill(J_b.begin(), J_b.end(), 0);
565 for (std::size_t p = 0; p < Xshape[0]; ++p)
567 auto dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
568 phi, std::pair(1, tdim + 1), p,
569 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
571 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
572 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
573 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
574 cmap.compute_jacobian(dphi, coord_dofs, _J);
575 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
576 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
577 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
578 cmap.compute_jacobian_inverse(_J, _K);
579 detJ[p] = cmap.compute_jacobian_determinant(_J, det_scratch);
584 for (std::size_t k0 = 0; k0 < basis_reference0.extent(0); ++k0)
585 for (std::size_t k1 = 0; k1 < basis_reference0.extent(1); ++k1)
586 for (std::size_t k2 = 0; k2 < basis_reference0.extent(2); ++k2)
587 basis_reference0(k0, k1, k2)
588 = basis_derivatives_reference0(0, k0, k1, k2);
590 for (std::size_t p = 0; p < Xshape[0]; ++p)
592 apply_dof_transformation0(
593 std::span(basis_reference0_b.data() + p * dim0 * value_size_ref0,
594 dim0 * value_size_ref0),
595 cell_info0, cells0[c], value_size_ref0);
598 for (std::size_t i = 0; i < basis0.extent(0); ++i)
600 auto _u = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
601 basis0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
602 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
603 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
604 basis_reference0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
605 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
606 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
607 K, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
608 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
609 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
610 J, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
611 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
612 push_forward_fn0(_u, _U, _J, detJ[i], _K);
616 const int dof_bs0 = dofmap0->bs();
617 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(cells0[c]);
618 for (std::size_t i = 0; i < dofs0.size(); ++i)
619 for (
int k = 0; k < dof_bs0; ++k)
620 coeffs0[dof_bs0 * i + k] = array0[dof_bs0 * dofs0[i] + k];
623 using X =
typename dolfinx::scalar_value_type_t<T>;
624 for (std::size_t p = 0; p < Xshape[0]; ++p)
626 for (
int k = 0; k < bs0; ++k)
628 for (std::size_t j = 0; j < value_size0; ++j)
631 for (std::size_t i = 0; i < dim0; ++i)
632 acc += coeffs0[bs0 * i + k] *
static_cast<X
>(basis0(p, i, j));
633 values0(p, 0, j * bs0 + k) = acc;
639 for (std::size_t i = 0; i < values0.extent(0); ++i)
641 auto _u = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
642 values0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
643 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
644 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
645 mapped_values0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
646 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
647 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
648 K, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
649 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
650 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
651 J, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
652 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
653 pull_back_fn1(_U, _u, _K, 1.0 / detJ[i], _J);
656 auto values = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
657 mapped_values0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0,
658 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
659 interpolation_apply(Pi_1, values, std::span(local1), bs1);
660 apply_inverse_dof_transform1(local1, cell_info1, cells1[c], 1);
663 const int dof_bs1 = dofmap1->bs();
664 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(c);
665 for (std::size_t i = 0; i < dofs1.size(); ++i)
666 for (
int k = 0; k < dof_bs1; ++k)
667 array1[dof_bs1 * dofs1[i] + k] = local1[dof_bs1 * i + k];
671template <dolfinx::scalar T, std::
floating_po
int U>
672void interpolate_nonmatching_meshes(
673 Function<T, U>& u,
const Function<T, U>& v,
674 std::span<const std::int32_t> cells,
675 const std::tuple<std::span<const std::int32_t>,
676 std::span<const std::int32_t>, std::span<const U>,
677 std::span<const std::int32_t>>& nmm_interpolation_data)
679 auto mesh = u.function_space()->mesh();
681 MPI_Comm comm = mesh->comm();
684 auto mesh_v = v.function_space()->mesh();
687 MPI_Comm_compare(comm, mesh_v->comm(), &result);
688 if (result == MPI_UNEQUAL)
690 throw std::runtime_error(
"Interpolation on different meshes is only "
691 "supported with the same communicator.");
695 assert(mesh->topology());
696 auto cell_map = mesh->topology()->index_map(mesh->topology()->dim());
699 auto element_u = u.function_space()->element();
701 const std::size_t value_size = u.function_space()->value_size();
703 const auto& [dest_ranks, src_ranks, recv_points, evaluation_cells]
704 = nmm_interpolation_data;
707 std::vector<T> send_values(recv_points.size() / 3 * value_size);
708 v.eval(recv_points, {recv_points.size() / 3, (std::size_t)3},
709 evaluation_cells, send_values, {recv_points.size() / 3, value_size});
712 std::vector<T> values_b(dest_ranks.size() * value_size);
713 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
714 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
715 _send_values(send_values.data(), src_ranks.size(), value_size);
716 impl::scatter_values(comm, src_ranks, dest_ranks, _send_values,
717 std::span(values_b));
720 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
721 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
722 values(values_b.data(), dest_ranks.size(), value_size);
723 std::vector<T> valuesT_b(value_size * dest_ranks.size());
724 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
725 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
726 valuesT(valuesT_b.data(), value_size, dest_ranks.size());
727 for (std::size_t i = 0; i < values.extent(0); ++i)
728 for (std::size_t j = 0; j < values.extent(1); ++j)
729 valuesT(j, i) = values(i, j);
732 fem::interpolate<T>(u, valuesT_b, {valuesT.extent(0), valuesT.extent(1)},
738template <dolfinx::scalar T, std::
floating_po
int U>
740 std::array<std::size_t, 2> fshape,
741 std::span<const std::int32_t> cells)
743 using cmdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
744 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
745 using cmdspan4_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
746 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>;
747 using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
748 U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
749 using mdspan3_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
750 U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>;
751 auto element = u.function_space()->element();
753 const int element_bs = element->block_size();
754 if (
int num_sub = element->num_sub_elements();
755 num_sub > 0 and num_sub != element_bs)
757 throw std::runtime_error(
"Cannot directly interpolate a mixed space. "
758 "Interpolate into subspaces.");
762 assert(u.function_space());
763 auto mesh = u.function_space()->mesh();
766 const int gdim = mesh->geometry().dim();
767 const int tdim = mesh->topology()->dim();
768 const bool symmetric = u.function_space()->symmetric();
770 if (fshape[0] != (std::size_t)u.function_space()->value_size())
771 throw std::runtime_error(
"Interpolation data has the wrong shape/size.");
773 std::span<const std::uint32_t> cell_info;
774 if (element->needs_dof_transformations())
776 mesh->topology_mutable()->create_entity_permutations();
777 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
780 const std::size_t f_shape1 = f.size() / u.function_space()->value_size();
781 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
782 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
783 _f(f.data(), fshape);
786 const auto dofmap = u.function_space()->dofmap();
788 const int dofmap_bs = dofmap->bs();
791 const int num_scalar_dofs = element->space_dimension() / element_bs;
792 const int value_size = u.function_space()->value_size() / element_bs;
794 std::span<T> coeffs = u.x()->mutable_array();
795 std::vector<T> _coeffs(num_scalar_dofs);
799 if (element->map_ident() && element->interpolation_ident())
804 auto apply_inv_transpose_dof_transformation
805 = element->template dof_transformation_fn<T>(
806 doftransform::inverse_transpose,
true);
810 std::size_t matrix_size = 0;
811 while (matrix_size * matrix_size < fshape[0])
814 for (std::size_t c = 0; c < cells.size(); ++c)
825 std::size_t rowstart = 0;
826 const std::int32_t
cell = cells[c];
827 std::span<const std::int32_t> dofs = dofmap->cell_dofs(
cell);
828 for (
int k = 0; k < element_bs; ++k)
830 if (k - rowstart > row)
838 std::next(f.begin(), (row * matrix_size + k - rowstart) * f_shape1
839 + c * num_scalar_dofs),
840 num_scalar_dofs, _coeffs.begin());
841 apply_inv_transpose_dof_transformation(_coeffs, cell_info,
cell, 1);
842 for (
int i = 0; i < num_scalar_dofs; ++i)
844 const int dof = i * element_bs + k;
845 std::div_t pos = std::div(dof, dofmap_bs);
846 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
854 for (std::size_t c = 0; c < cells.size(); ++c)
856 const std::int32_t
cell = cells[c];
857 std::span<const std::int32_t> dofs = dofmap->cell_dofs(
cell);
858 for (
int k = 0; k < element_bs; ++k)
862 std::copy_n(std::next(f.begin(), k * f_shape1 + c * num_scalar_dofs),
863 num_scalar_dofs, _coeffs.begin());
864 apply_inv_transpose_dof_transformation(_coeffs, cell_info,
cell, 1);
865 for (
int i = 0; i < num_scalar_dofs; ++i)
867 const int dof = i * element_bs + k;
868 std::div_t pos = std::div(dof, dofmap_bs);
869 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
875 else if (element->map_ident())
882 throw std::runtime_error(
883 "Interpolation into this element not supported.");
886 const int element_vs = u.function_space()->value_size() / element_bs;
888 if (element_vs > 1 && element_bs > 1)
890 throw std::runtime_error(
891 "Interpolation into this element not supported.");
895 const auto [_Pi, pi_shape] = element->interpolation_operator();
896 cmdspan2_t Pi(_Pi.data(), pi_shape);
897 const std::size_t num_interp_points = Pi.extent(1);
898 assert(Pi.extent(0) == num_scalar_dofs);
900 auto apply_inv_transpose_dof_transformation
901 = element->template dof_transformation_fn<T>(
902 doftransform::inverse_transpose,
true);
905 std::vector<T> ref_data_b(num_interp_points);
906 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
907 T, MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
908 std::size_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent, 1>>
909 ref_data(ref_data_b.data(), num_interp_points, 1);
910 for (std::size_t c = 0; c < cells.size(); ++c)
912 const std::int32_t
cell = cells[c];
913 std::span<const std::int32_t> dofs = dofmap->cell_dofs(
cell);
914 for (
int k = 0; k < element_bs; ++k)
916 for (
int i = 0; i < element_vs; ++i)
919 std::next(f.begin(), (i + k) * f_shape1
920 + c * num_interp_points / element_vs),
921 num_interp_points / element_vs,
922 std::next(ref_data_b.begin(),
923 i * num_interp_points / element_vs));
925 impl::interpolation_apply(Pi, ref_data, std::span(_coeffs), 1);
926 apply_inv_transpose_dof_transformation(_coeffs, cell_info,
cell, 1);
927 for (
int i = 0; i < num_scalar_dofs; ++i)
929 const int dof = i * element_bs + k;
930 std::div_t pos = std::div(dof, dofmap_bs);
931 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
940 throw std::runtime_error(
941 "Interpolation into this element not supported.");
944 const auto [X, Xshape] = element->interpolation_points();
947 throw std::runtime_error(
948 "Interpolation into this space is not yet supported.");
951 if (_f.extent(1) != cells.size() * Xshape[0])
952 throw std::runtime_error(
"Interpolation data has the wrong shape.");
958 auto x_dofmap = mesh->geometry().dofmap();
959 const int num_dofs_g = cmap.
dim();
960 std::span<const U> x_g = mesh->geometry().x();
963 std::vector<U> J_b(Xshape[0] * gdim * tdim);
964 mdspan3_t J(J_b.data(), Xshape[0], gdim, tdim);
965 std::vector<U> K_b(Xshape[0] * tdim * gdim);
966 mdspan3_t K(K_b.data(), Xshape[0], tdim, gdim);
967 std::vector<U> detJ(Xshape[0]);
968 std::vector<U> det_scratch(2 * gdim * tdim);
970 std::vector<U> coord_dofs_b(num_dofs_g * gdim);
971 mdspan2_t coord_dofs(coord_dofs_b.data(), num_dofs_g, gdim);
973 std::vector<T> ref_data_b(Xshape[0] * 1 * value_size);
974 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
975 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
976 ref_data(ref_data_b.data(), Xshape[0], 1, value_size);
978 std::vector<T> _vals_b(Xshape[0] * 1 * value_size);
979 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
980 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
981 _vals(_vals_b.data(), Xshape[0], 1, value_size);
985 std::array<std::size_t, 4> phi_shape = cmap.
tabulate_shape(1, Xshape[0]);
986 std::vector<U> phi_b(
987 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
988 cmdspan4_t phi(phi_b.data(), phi_shape);
990 auto dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
991 phi, std::pair(1, tdim + 1),
992 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
993 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
995 const std::function<void(std::span<T>, std::span<const std::uint32_t>,
997 apply_inverse_transpose_dof_transformation
998 = element->template dof_transformation_fn<T>(
999 doftransform::inverse_transpose);
1002 const auto [_Pi, pi_shape] = element->interpolation_operator();
1003 cmdspan2_t Pi(_Pi.data(), pi_shape);
1005 using u_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
1006 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
1007 using U_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
1008 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
1009 using J_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
1010 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
1011 using K_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
1012 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
1014 = element->basix_element().template map_fn<U_t, u_t, J_t, K_t>();
1016 for (std::size_t c = 0; c < cells.size(); ++c)
1018 const std::int32_t
cell = cells[c];
1019 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1020 x_dofmap,
cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1021 for (
int i = 0; i < num_dofs_g; ++i)
1023 const int pos = 3 * x_dofs[i];
1024 for (
int j = 0; j < gdim; ++j)
1025 coord_dofs(i, j) = x_g[pos + j];
1029 std::fill(J_b.begin(), J_b.end(), 0);
1030 for (std::size_t p = 0; p < Xshape[0]; ++p)
1032 auto _dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1033 dphi, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, p,
1034 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1035 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1036 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1037 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1039 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1040 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1041 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1046 std::span<const std::int32_t> dofs = dofmap->cell_dofs(
cell);
1047 for (
int k = 0; k < element_bs; ++k)
1050 for (
int m = 0; m < value_size; ++m)
1052 for (std::size_t k0 = 0; k0 < Xshape[0]; ++k0)
1055 = f[f_shape1 * (k * value_size + m) + c * Xshape[0] + k0];
1060 for (std::size_t i = 0; i < Xshape[0]; ++i)
1062 auto _u = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1063 _vals, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1064 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1065 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1066 ref_data, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1067 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1068 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1069 K, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1070 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1071 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1072 J, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1073 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1074 pull_back_fn(_U, _u, _K, 1.0 / detJ[i], _J);
1077 auto ref = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1078 ref_data, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0,
1079 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1080 impl::interpolation_apply(Pi, ref, std::span(_coeffs), element_bs);
1081 apply_inverse_transpose_dof_transformation(_coeffs, cell_info,
cell, 1);
1084 assert(_coeffs.size() == num_scalar_dofs);
1085 for (
int i = 0; i < num_scalar_dofs; ++i)
1087 const int dof = i * element_bs + k;
1088 std::div_t pos = std::div(dof, dofmap_bs);
1089 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
1113template <std::
floating_po
int T>
1114std::tuple<std::vector<std::int32_t>, std::vector<std::int32_t>, std::vector<T>,
1115 std::vector<std::int32_t>>
1118 const mesh::Mesh<T>& mesh1, std::span<const std::int32_t> cells, T padding)
1122 const std::vector<T> coords
1126 std::vector<T> x(coords.size());
1127 std::size_t num_points = coords.size() / 3;
1128 for (std::size_t i = 0; i < num_points; ++i)
1129 for (std::size_t j = 0; j < 3; ++j)
1130 x[3 * i + j] = coords[i + j * num_points];
1133 return geometry::determine_point_ownership<T>(mesh1, x, padding);
1149template <std::
floating_po
int T>
1150std::tuple<std::vector<std::int32_t>, std::vector<std::int32_t>, std::vector<T>,
1151 std::vector<std::int32_t>>
1157 int tdim = mesh0.
topology()->dim();
1158 auto cell_map = mesh0.
topology()->index_map(tdim);
1160 std::int32_t num_cells = cell_map->size_local() + cell_map->num_ghosts();
1161 std::vector<std::int32_t> cells(num_cells, 0);
1162 std::iota(cells.begin(), cells.end(), 0);
1164 mesh0.
geometry(), element0, mesh1, cells, padding);
1177template <dolfinx::scalar T, std::
floating_po
int U>
1180 std::span<const std::int32_t> cells1,
1181 std::span<const std::int32_t> cell_map,
1182 const std::tuple<std::span<const std::int32_t>,
1183 std::span<const std::int32_t>, std::span<const U>,
1184 std::span<const std::int32_t>>& nmm_interpolation_data
1192 auto cell_map0 = mesh->topology()->index_map(mesh->topology()->dim());
1194 std::size_t num_cells0 = cell_map0->size_local() + cell_map0->num_ghosts();
1196 and cells1.size() == num_cells0)
1199 std::span<T> u1_array = u1.
x()->mutable_array();
1200 std::span<const T> u0_array = u0.
x()->array();
1201 std::copy(u0_array.begin(), u0_array.end(), u1_array.begin());
1205 std::vector<std::int32_t> cells0;
1206 cells0.reserve(cells1.size());
1210 cells0.insert(cells0.end(), cells1.begin(), cells1.end());
1213 else if (cell_map.size() > 0)
1215 std::transform(cells1.begin(), cells1.end(), std::back_inserter(cells0),
1216 [&cell_map](std::int32_t c) { return cell_map[c]; });
1221 impl::interpolate_nonmatching_meshes(u1, u0, cells1,
1222 nmm_interpolation_data);
1228 auto element0 = fs0->element();
1231 auto element1 = fs1->element();
1233 if (fs0->value_shape().size() != fs1->value_shape().size()
1234 or !std::equal(fs0->value_shape().begin(), fs0->value_shape().end(),
1235 fs1->value_shape().begin()))
1237 throw std::runtime_error(
1238 "Interpolation: elements have different value dimensions");
1241 if (element1 == element0 or *element1 == *element0)
1245 const int tdim = mesh->topology()->dim();
1246 auto cell_map1 = mesh->topology()->index_map(tdim);
1249 assert(element1->block_size() == element0->block_size());
1252 std::shared_ptr<const DofMap> dofmap0 = u0.
function_space()->dofmap();
1254 std::shared_ptr<const DofMap> dofmap1 = u1.
function_space()->dofmap();
1257 std::span<T> u1_array = u1.
x()->mutable_array();
1258 std::span<const T> u0_array = u0.
x()->array();
1261 const int bs0 = dofmap0->bs();
1262 const int bs1 = dofmap1->bs();
1263 for (std::size_t c = 0; c < cells1.size(); ++c)
1265 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(cells0[c]);
1266 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(cells1[c]);
1267 assert(bs0 * dofs0.size() == bs1 * dofs1.size());
1268 for (std::size_t i = 0; i < dofs0.size(); ++i)
1270 for (
int k = 0; k < bs0; ++k)
1272 int index = bs0 * i + k;
1273 std::div_t dv1 = std::div(index, bs1);
1274 u1_array[bs1 * dofs1[dv1.quot] + dv1.rem]
1275 = u0_array[bs0 * dofs0[i] + k];
1280 else if (element1->map_type() == element0->map_type())
1283 impl::interpolate_same_map(u1, u0, cells1, cells0);
1288 impl::interpolate_nonmatching_maps(u1, u0, cells1, cells0);
Degree-of-freedom map representations and tools.
Definition CoordinateElement.h:38
static void compute_jacobian(const U &dphi, const V &cell_geometry, W &&J)
Definition CoordinateElement.h:105
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:54
static void compute_jacobian_inverse(const U &J, V &&K)
Compute the inverse of the Jacobian.
Definition CoordinateElement.h:114
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:47
int dim() const
The dimension of the coordinate element space.
Definition CoordinateElement.cpp:194
static double compute_jacobian_determinant(const U &J, std::span< typename U::value_type > w)
Compute the determinant of the Jacobian.
Definition CoordinateElement.h:131
Model of a finite element.
Definition FiniteElement.h:40
std::shared_ptr< const FunctionSpace< geometry_type > > function_space() const
Access the function space.
Definition Function.h:142
std::shared_ptr< const la::Vector< value_type > > x() const
Underlying vector.
Definition Function.h:148
Geometry stores the geometry imposed on a mesh.
Definition Geometry.h:33
const fem::CoordinateElement< value_type > & cmap() const
The element that describes the geometry map.
Definition Geometry.h:180
int dim() const
Return Euclidean dimension of coordinate system.
Definition Geometry.h:116
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:123
std::span< const value_type > x() const
Access geometry degrees-of-freedom data (const version).
Definition Geometry.h:168
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition Mesh.h:23
std::shared_ptr< Topology > topology()
Get mesh topology.
Definition Mesh.h:64
Geometry< T > & geometry()
Get mesh geometry.
Definition Mesh.h:76
Definition interpolate.h:32
void cells(la::SparsityPattern &pattern, std::array< std::span< const std::int32_t >, 2 > cells, std::array< std::reference_wrapper< const DofMap >, 2 > dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition sparsitybuild.cpp:15
Finite element method functionality.
Definition assemble_matrix_impl.h:26
std::tuple< std::vector< std::int32_t >, std::vector< std::int32_t >, std::vector< T >, std::vector< std::int32_t > > create_nonmatching_meshes_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 discrete functions across different meshes.
Definition interpolate.h:1116
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 expression f(x) in a finite element space.
Definition interpolate.h:739
@ 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