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>
41template <std::
floating_po
int T>
44 std::span<const std::int32_t> cells)
47 const std::size_t gdim = geometry.
dim();
48 auto x_dofmap = geometry.
dofmap();
49 std::span<const T> x_g = geometry.
x();
51 if (geometry.
cmaps().size() > 1)
53 throw std::runtime_error(
"Mixed topology not supported");
57 const std::size_t num_dofs_g = cmap.
dim();
63 std::array<std::size_t, 4> phi_shape = cmap.
tabulate_shape(0, Xshape[0]);
65 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
66 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
67 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>
68 phi_full(phi_b.data(), phi_shape);
70 auto phi = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
71 submdspan(phi_full, 0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
72 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
76 std::vector<T> coordinate_dofs(num_dofs_g * gdim, 0);
77 std::vector<T> x(3 * (cells.size() * Xshape[0]), 0);
78 for (std::size_t c = 0; c < cells.size(); ++c)
81 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::
82 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
83 x_dofmap, cells[c], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
84 for (std::size_t i = 0; i < x_dofs.size(); ++i)
86 std::copy_n(std::next(x_g.begin(), 3 * x_dofs[i]), gdim,
87 std::next(coordinate_dofs.begin(), i * gdim));
91 for (std::size_t p = 0; p < Xshape[0]; ++p)
93 for (std::size_t j = 0; j < gdim; ++j)
96 for (std::size_t k = 0; k < num_dofs_g; ++k)
97 acc += phi(p, k) * coordinate_dofs[k * gdim + j];
98 x[j * (cells.size() * Xshape[0]) + c * Xshape[0] + p] = acc;
120template <dolfinx::scalar T, std::
floating_po
int U>
121void interpolate(Function<T, U>& u, std::span<const T> f,
122 std::array<std::size_t, 2> fshape,
123 std::span<const std::int32_t> cells);
128template <
typename T, std::
size_t D>
129using mdspan_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
130 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, D>>;
151template <dolfinx::scalar T>
153 MPI_Comm comm, std::span<const std::int32_t> src_ranks,
154 std::span<const std::int32_t> dest_ranks,
155 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
156 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
158 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 out_ranks.erase(std::unique(out_ranks.begin(), out_ranks.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; });
179 std::sort(in_ranks.begin(), in_ranks.end());
180 in_ranks.erase(std::unique(in_ranks.begin(), in_ranks.end()), in_ranks.end());
181 in_ranks.reserve(in_ranks.size() + 1);
184 MPI_Comm reverse_comm;
185 MPI_Dist_graph_create_adjacent(
186 comm, in_ranks.size(), in_ranks.data(), MPI_UNWEIGHTED, out_ranks.size(),
187 out_ranks.data(), MPI_UNWEIGHTED, MPI_INFO_NULL,
false, &reverse_comm);
189 std::vector<std::int32_t> comm_to_output;
190 std::vector<std::int32_t> recv_sizes(in_ranks.size());
191 recv_sizes.reserve(1);
192 std::vector<std::int32_t> recv_offsets(in_ranks.size() + 1, 0);
195 std::map<std::int32_t, std::int32_t> rank_to_neighbor;
196 for (std::size_t i = 0; i < in_ranks.size(); i++)
197 rank_to_neighbor[in_ranks[i]] = i;
201 dest_ranks.begin(), dest_ranks.end(),
202 [&dest_ranks, &rank_to_neighbor, &recv_sizes, block_size](
auto rank)
206 const int neighbor = rank_to_neighbor[rank];
207 recv_sizes[neighbor] += block_size;
212 std::partial_sum(recv_sizes.begin(), recv_sizes.end(),
213 std::next(recv_offsets.begin(), 1));
216 comm_to_output.resize(recv_offsets.back() / block_size);
217 std::vector<std::int32_t> recv_counter(recv_sizes.size(), 0);
218 for (std::size_t i = 0; i < dest_ranks.size(); ++i)
220 if (
const std::int32_t rank = dest_ranks[i]; rank >= 0)
222 const int neighbor = rank_to_neighbor[rank];
223 int insert_pos = recv_offsets[neighbor] + recv_counter[neighbor];
224 comm_to_output[insert_pos / block_size] = i * block_size;
225 recv_counter[neighbor] += block_size;
230 std::vector<std::int32_t> send_sizes(out_ranks.size());
231 send_sizes.reserve(1);
234 std::map<std::int32_t, std::int32_t> rank_to_neighbor;
235 for (std::size_t i = 0; i < out_ranks.size(); i++)
236 rank_to_neighbor[out_ranks[i]] = i;
239 std::for_each(src_ranks.begin(), src_ranks.end(),
240 [&rank_to_neighbor, &send_sizes, block_size](
auto rank)
241 { send_sizes[rank_to_neighbor[rank]] += block_size; });
245 std::vector<std::int32_t> send_offsets(send_sizes.size() + 1, 0);
246 std::partial_sum(send_sizes.begin(), send_sizes.end(),
247 std::next(send_offsets.begin(), 1));
250 std::vector<T> values(recv_offsets.back());
252 MPI_Neighbor_alltoallv(send_values.data_handle(), send_sizes.data(),
253 send_offsets.data(), dolfinx::MPI::mpi_type<T>(),
254 values.data(), recv_sizes.data(), recv_offsets.data(),
255 dolfinx::MPI::mpi_type<T>(), reverse_comm);
256 MPI_Comm_free(&reverse_comm);
259 std::fill(recv_values.begin(), recv_values.end(), T(0));
260 for (std::size_t i = 0; i < comm_to_output.size(); i++)
262 auto vals = std::next(recv_values.begin(), comm_to_output[i]);
263 auto vals_from = std::next(values.begin(), i * block_size);
264 std::copy_n(vals_from, block_size, vals);
276template <
typename U,
typename V, dolfinx::scalar T>
278 requires std::convertible_to<
279 U, MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
280 const typename std::decay_t<U>::value_type,
281 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>>;
282 requires std::convertible_to<
283 V, MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
284 const typename std::decay_t<V>::value_type,
285 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>>;
287void interpolation_apply(U&& Pi, V&& data, std::span<T> coeffs,
int bs)
289 using X =
typename dolfinx::scalar_value_type_t<T>;
294 assert(data.extent(0) * data.extent(1) == Pi.extent(1));
295 for (std::size_t i = 0; i < Pi.extent(0); ++i)
298 for (std::size_t k = 0; k < data.extent(1); ++k)
299 for (std::size_t j = 0; j < data.extent(0); ++j)
301 +=
static_cast<X
>(Pi(i, k * data.extent(0) + j)) * data(j, k);
306 const std::size_t cols = Pi.extent(1);
307 assert(data.extent(0) == Pi.extent(1));
308 assert(data.extent(1) == bs);
309 for (
int k = 0; k < bs; ++k)
311 for (std::size_t i = 0; i < Pi.extent(0); ++i)
314 for (std::size_t j = 0; j < cols; ++j)
315 acc +=
static_cast<X
>(Pi(i, j)) * data(j, k);
316 coeffs[bs * i + k] = acc;
332template <dolfinx::scalar T, std::
floating_po
int U>
333void interpolate_same_map(Function<T, U>& u1,
const Function<T, U>& u0,
334 std::span<const std::int32_t> cells)
336 auto V0 = u0.function_space();
338 auto V1 = u1.function_space();
340 auto mesh = V0->mesh();
343 auto element0 = V0->element();
345 auto element1 = V1->element();
348 const int tdim = mesh->topology()->dim();
349 auto map = mesh->topology()->index_map(tdim);
351 std::span<T> u1_array = u1.x()->mutable_array();
352 std::span<const T> u0_array = u0.x()->array();
354 std::span<const std::uint32_t> cell_info;
355 if (element1->needs_dof_transformations()
356 or element0->needs_dof_transformations())
358 mesh->topology_mutable()->create_entity_permutations();
359 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
363 auto dofmap1 = V1->dofmap();
364 auto dofmap0 = V0->dofmap();
367 const int bs1 = dofmap1->bs();
368 const int bs0 = dofmap0->bs();
369 auto apply_dof_transformation
370 = element0->template get_dof_transformation_function<T>(
false,
true,
372 auto apply_inverse_dof_transform
373 = element1->template get_dof_transformation_function<T>(
true,
true,
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>;
388 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(c);
389 for (std::size_t i = 0; i < dofs0.size(); ++i)
390 for (
int k = 0; k < bs0; ++k)
391 local0[bs0 * i + k] = u0_array[bs0 * dofs0[i] + k];
393 apply_dof_transformation(local0, cell_info, c, 1);
397 std::fill(local1.begin(), local1.end(), 0);
398 for (std::size_t i = 0; i < im_shape[0]; ++i)
399 for (std::size_t j = 0; j < im_shape[1]; ++j)
400 local1[i] +=
static_cast<X
>(i_m[im_shape[1] * i + j]) * local0[j];
402 apply_inverse_dof_transform(local1, cell_info, c, 1);
404 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(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];
421template <dolfinx::scalar T, std::
floating_po
int U>
422void interpolate_nonmatching_maps(Function<T, U>& u1,
const Function<T, U>& u0,
423 std::span<const std::int32_t> cells)
426 auto V0 = u0.function_space();
428 auto mesh = V0->mesh();
432 const int tdim = mesh->topology()->dim();
433 const int gdim = mesh->geometry().dim();
436 auto V1 = u1.function_space();
438 auto element0 = V0->element();
440 auto element1 = V1->element();
443 std::span<const std::uint32_t> cell_info;
444 if (element1->needs_dof_transformations()
445 or element0->needs_dof_transformations())
447 mesh->topology_mutable()->create_entity_permutations();
448 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
452 auto dofmap0 = V0->dofmap();
453 auto dofmap1 = V1->dofmap();
455 const auto [X, Xshape] = element1->interpolation_points();
458 const int bs0 = element0->block_size();
459 const int bs1 = element1->block_size();
460 auto apply_dof_transformation0
461 = element0->template get_dof_transformation_function<U>(
false,
false,
463 auto apply_inverse_dof_transform1
464 = element1->template get_dof_transformation_function<T>(
true,
true,
468 const std::size_t dim0 = element0->space_dimension() / bs0;
469 const std::size_t value_size_ref0 = element0->reference_value_size() / bs0;
470 const std::size_t value_size0 = element0->value_size() / bs0;
473 if (mesh->geometry().cmaps().size() > 1)
474 throw std::runtime_error(
"Multiple cmaps");
476 const CoordinateElement<U>& cmap = mesh->geometry().cmaps()[0];
477 auto x_dofmap = mesh->geometry().dofmap();
478 const std::size_t num_dofs_g = cmap.dim();
479 std::span<const U> x_g = mesh->geometry().x();
482 const std::array<std::size_t, 4> phi_shape
483 = cmap.tabulate_shape(1, Xshape[0]);
484 std::vector<U> phi_b(
485 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
486 mdspan_t<const U, 4> phi(phi_b.data(), phi_shape);
487 cmap.tabulate(1, X, Xshape, phi_b);
490 const auto [_basis_derivatives_reference0, b0shape]
491 = element0->tabulate(X, Xshape, 0);
492 mdspan_t<const U, 4> basis_derivatives_reference0(
493 _basis_derivatives_reference0.data(), b0shape);
496 std::vector<T> local1(element1->space_dimension());
497 std::vector<T> coeffs0(element0->space_dimension());
499 std::vector<U> basis0_b(Xshape[0] * dim0 * value_size0);
500 impl::mdspan_t<U, 3> basis0(basis0_b.data(), Xshape[0], dim0, value_size0);
502 std::vector<U> basis_reference0_b(Xshape[0] * dim0 * value_size_ref0);
503 impl::mdspan_t<U, 3> basis_reference0(basis_reference0_b.data(), Xshape[0],
504 dim0, value_size_ref0);
506 std::vector<T> values0_b(Xshape[0] * 1 * element1->value_size());
507 impl::mdspan_t<T, 3> values0(values0_b.data(), Xshape[0], 1,
508 element1->value_size());
510 std::vector<T> mapped_values_b(Xshape[0] * 1 * element1->value_size());
511 impl::mdspan_t<T, 3> mapped_values0(mapped_values_b.data(), Xshape[0], 1,
512 element1->value_size());
514 std::vector<U> coord_dofs_b(num_dofs_g * gdim);
515 impl::mdspan_t<U, 2> coord_dofs(coord_dofs_b.data(), num_dofs_g, gdim);
517 std::vector<U> J_b(Xshape[0] * gdim * tdim);
518 impl::mdspan_t<U, 3> J(J_b.data(), Xshape[0], gdim, tdim);
519 std::vector<U> K_b(Xshape[0] * tdim * gdim);
520 impl::mdspan_t<U, 3> K(K_b.data(), Xshape[0], tdim, gdim);
521 std::vector<U> detJ(Xshape[0]);
522 std::vector<U> det_scratch(2 * gdim * tdim);
525 const auto [_Pi_1, pi_shape] = element1->interpolation_operator();
526 impl::mdspan_t<const U, 2> Pi_1(_Pi_1.data(), pi_shape);
528 using u_t = impl::mdspan_t<U, 2>;
529 using U_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
530 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
531 using J_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
532 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
533 using K_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
534 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
535 auto push_forward_fn0
536 = element0->basix_element().template map_fn<u_t, U_t, J_t, K_t>();
538 using v_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
539 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
540 using V_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
541 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
543 = element1->basix_element().template map_fn<V_t, v_t, K_t, J_t>();
546 std::span<const T> array0 = u0.x()->array();
547 std::span<T> array1 = u1.x()->mutable_array();
552 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE ::
553 submdspan(x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
554 for (std::size_t i = 0; i < num_dofs_g; ++i)
556 const int pos = 3 * x_dofs[i];
557 for (std::size_t j = 0; j < gdim; ++j)
558 coord_dofs(i, j) = x_g[pos + j];
562 std::fill(J_b.begin(), J_b.end(), 0);
563 for (std::size_t p = 0; p < Xshape[0]; ++p)
566 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE ::
567 submdspan(phi, std::pair(1, tdim + 1), p,
568 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
571 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE ::
572 submdspan(J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
573 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
574 cmap.compute_jacobian(dphi, coord_dofs, _J);
576 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE ::
577 submdspan(K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
578 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
579 cmap.compute_jacobian_inverse(_J, _K);
580 detJ[p] = cmap.compute_jacobian_determinant(_J, det_scratch);
585 for (std::size_t k0 = 0; k0 < basis_reference0.extent(0); ++k0)
586 for (std::size_t k1 = 0; k1 < basis_reference0.extent(1); ++k1)
587 for (std::size_t k2 = 0; k2 < basis_reference0.extent(2); ++k2)
588 basis_reference0(k0, k1, k2)
589 = basis_derivatives_reference0(0, k0, k1, k2);
591 for (std::size_t p = 0; p < Xshape[0]; ++p)
593 apply_dof_transformation0(
594 std::span(basis_reference0_b.data() + p * dim0 * value_size_ref0,
595 dim0 * value_size_ref0),
596 cell_info, c, value_size_ref0);
599 for (std::size_t i = 0; i < basis0.extent(0); ++i)
602 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE ::
603 submdspan(basis0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
604 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
606 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE ::
607 submdspan(basis_reference0, i,
608 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
609 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
611 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE ::
612 submdspan(K, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
613 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
615 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE ::
616 submdspan(J, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
617 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
618 push_forward_fn0(_u, _U, _J, detJ[i], _K);
622 const int dof_bs0 = dofmap0->bs();
623 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(c);
624 for (std::size_t i = 0; i < dofs0.size(); ++i)
625 for (
int k = 0; k < dof_bs0; ++k)
626 coeffs0[dof_bs0 * i + k] = array0[dof_bs0 * dofs0[i] + k];
629 using X =
typename dolfinx::scalar_value_type_t<T>;
630 for (std::size_t p = 0; p < Xshape[0]; ++p)
632 for (
int k = 0; k < bs0; ++k)
634 for (std::size_t j = 0; j < value_size0; ++j)
637 for (std::size_t i = 0; i < dim0; ++i)
638 acc += coeffs0[bs0 * i + k] *
static_cast<X
>(basis0(p, i, j));
639 values0(p, 0, j * bs0 + k) = acc;
645 for (std::size_t i = 0; i < values0.extent(0); ++i)
647 auto _u = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
648 submdspan(values0, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
649 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
650 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
651 submdspan(mapped_values0, i,
652 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
653 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
654 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
655 submdspan(K, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
656 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
657 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
658 submdspan(J, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
659 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
660 pull_back_fn1(_U, _u, _K, 1.0 / detJ[i], _J);
664 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE ::
665 submdspan(mapped_values0,
666 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0,
667 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
668 interpolation_apply(Pi_1, values, std::span(local1), bs1);
669 apply_inverse_dof_transform1(local1, cell_info, c, 1);
672 const int dof_bs1 = dofmap1->bs();
673 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(c);
674 for (std::size_t i = 0; i < dofs1.size(); ++i)
675 for (
int k = 0; k < dof_bs1; ++k)
676 array1[dof_bs1 * dofs1[i] + k] = local1[dof_bs1 * i + k];
680template <dolfinx::scalar T, std::
floating_po
int U>
681void interpolate_nonmatching_meshes(
682 Function<T, U>& u,
const Function<T, U>& v,
683 std::span<const std::int32_t> cells,
684 const std::tuple<std::vector<std::int32_t>, std::vector<std::int32_t>,
685 std::vector<U>, std::vector<std::int32_t>>&
686 nmm_interpolation_data)
688 auto mesh = u.function_space()->mesh();
690 MPI_Comm comm = mesh->comm();
693 auto mesh_v = v.function_space()->mesh();
696 MPI_Comm_compare(comm, mesh_v->comm(), &result);
697 if (result == MPI_UNEQUAL)
699 throw std::runtime_error(
"Interpolation on different meshes is only "
700 "supported with the same communicator.");
704 assert(mesh->topology());
705 auto cell_map = mesh->topology()->index_map(mesh->topology()->dim());
708 auto element_u = u.function_space()->element();
710 const std::size_t value_size = element_u->value_size();
712 const auto& [dest_ranks, src_ranks, recv_points, evaluation_cells]
713 = nmm_interpolation_data;
716 std::vector<T> send_values(recv_points.size() / 3 * value_size);
717 v.eval(recv_points, {recv_points.size() / 3, (std::size_t)3},
718 evaluation_cells, send_values, {recv_points.size() / 3, value_size});
721 std::vector<T> values_b(dest_ranks.size() * value_size);
722 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
723 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
724 _send_values(send_values.data(), src_ranks.size(), value_size);
725 impl::scatter_values(comm, src_ranks, dest_ranks, _send_values,
726 std::span(values_b));
729 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
730 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
731 values(values_b.data(), dest_ranks.size(), value_size);
732 std::vector<T> valuesT_b(value_size * dest_ranks.size());
733 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
734 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
735 valuesT(valuesT_b.data(), value_size, dest_ranks.size());
736 for (std::size_t i = 0; i < values.extent(0); ++i)
737 for (std::size_t j = 0; j < values.extent(1); ++j)
738 valuesT(j, i) = values(i, j);
741 fem::interpolate<T>(u, valuesT_b, {valuesT.extent(0), valuesT.extent(1)},
747template <dolfinx::scalar T, std::
floating_po
int U>
749 std::array<std::size_t, 2> fshape,
750 std::span<const std::int32_t> cells)
752 using cmdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
753 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
754 using cmdspan4_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
755 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>;
756 using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
757 U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
758 using mdspan3_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
759 U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>;
762 const int element_bs = element->block_size();
763 if (
int num_sub = element->num_sub_elements();
764 num_sub > 0 and num_sub != element_bs)
766 throw std::runtime_error(
"Cannot directly interpolate a mixed space. "
767 "Interpolate into subspaces.");
770 if (fshape[0] != (std::size_t)element->value_size())
771 throw std::runtime_error(
"Interpolation data has the wrong shape/size.");
778 const int gdim = mesh->geometry().dim();
779 const int tdim = mesh->topology()->dim();
781 std::span<const std::uint32_t> cell_info;
782 if (element->needs_dof_transformations())
784 mesh->topology_mutable()->create_entity_permutations();
785 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
788 const std::size_t f_shape1 = f.size() / element->value_size();
789 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
790 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
791 _f(f.data(), fshape);
796 const int dofmap_bs = dofmap->bs();
799 const int num_scalar_dofs = element->space_dimension() / element_bs;
800 const int value_size = element->value_size() / element_bs;
802 std::span<T> coeffs = u.
x()->mutable_array();
803 std::vector<T> _coeffs(num_scalar_dofs);
807 if (element->map_ident() && element->interpolation_ident())
812 auto apply_inv_transpose_dof_transformation
813 = element->template get_dof_transformation_function<T>(
true,
true,
817 for (std::size_t c = 0; c < cells.size(); ++c)
819 const std::int32_t
cell = cells[c];
820 std::span<const std::int32_t> dofs = dofmap->cell_dofs(
cell);
821 for (
int k = 0; k < element_bs; ++k)
825 std::copy_n(std::next(f.begin(), k * f_shape1 + c * num_scalar_dofs),
826 num_scalar_dofs, _coeffs.begin());
827 apply_inv_transpose_dof_transformation(_coeffs, cell_info,
cell, 1);
828 for (
int i = 0; i < num_scalar_dofs; ++i)
830 const int dof = i * element_bs + k;
831 std::div_t pos = std::div(dof, dofmap_bs);
832 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
837 else if (element->map_ident())
842 const int element_vs = element->value_size() / element_bs;
844 if (element_vs > 1 && element_bs > 1)
846 throw std::runtime_error(
847 "Interpolation into this element not supported.");
851 const auto [_Pi, pi_shape] = element->interpolation_operator();
852 cmdspan2_t Pi(_Pi.data(), pi_shape);
853 const std::size_t num_interp_points = Pi.extent(1);
854 assert(Pi.extent(0) == num_scalar_dofs);
856 auto apply_inv_transpose_dof_transformation
857 = element->template get_dof_transformation_function<T>(
true,
true,
861 std::vector<T> ref_data_b(num_interp_points);
862 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
863 T, MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
864 std::size_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent, 1>>
865 ref_data(ref_data_b.data(), num_interp_points, 1);
866 for (std::size_t c = 0; c < cells.size(); ++c)
868 const std::int32_t
cell = cells[c];
869 std::span<const std::int32_t> dofs = dofmap->cell_dofs(
cell);
870 for (
int k = 0; k < element_bs; ++k)
872 for (
int i = 0; i < element_vs; ++i)
875 std::next(f.begin(), (i + k) * f_shape1
876 + c * num_interp_points / element_vs),
877 num_interp_points / element_vs,
878 std::next(ref_data_b.begin(),
879 i * num_interp_points / element_vs));
881 impl::interpolation_apply(Pi, ref_data, std::span(_coeffs), 1);
882 apply_inv_transpose_dof_transformation(_coeffs, cell_info,
cell, 1);
883 for (
int i = 0; i < num_scalar_dofs; ++i)
885 const int dof = i * element_bs + k;
886 std::div_t pos = std::div(dof, dofmap_bs);
887 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
895 const auto [X, Xshape] = element->interpolation_points();
898 throw std::runtime_error(
899 "Interpolation into this space is not yet supported.");
902 if (_f.extent(1) != cells.size() * Xshape[0])
903 throw std::runtime_error(
"Interpolation data has the wrong shape.");
906 if (mesh->geometry().cmaps().size() > 1)
907 throw std::runtime_error(
"Multiple cmaps");
911 auto x_dofmap = mesh->geometry().dofmap();
912 const int num_dofs_g = cmap.
dim();
913 std::span<const U> x_g = mesh->geometry().x();
916 std::vector<U> J_b(Xshape[0] * gdim * tdim);
917 mdspan3_t J(J_b.data(), Xshape[0], gdim, tdim);
918 std::vector<U> K_b(Xshape[0] * tdim * gdim);
919 mdspan3_t K(K_b.data(), Xshape[0], tdim, gdim);
920 std::vector<U> detJ(Xshape[0]);
921 std::vector<U> det_scratch(2 * gdim * tdim);
923 std::vector<U> coord_dofs_b(num_dofs_g * gdim);
924 mdspan2_t coord_dofs(coord_dofs_b.data(), num_dofs_g, gdim);
926 std::vector<T> ref_data_b(Xshape[0] * 1 * value_size);
927 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
928 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
929 ref_data(ref_data_b.data(), Xshape[0], 1, value_size);
931 std::vector<T> _vals_b(Xshape[0] * 1 * value_size);
932 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
933 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
934 _vals(_vals_b.data(), Xshape[0], 1, value_size);
938 std::array<std::size_t, 4> phi_shape = cmap.
tabulate_shape(1, Xshape[0]);
939 std::vector<U> phi_b(
940 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
941 cmdspan4_t phi(phi_b.data(), phi_shape);
943 auto dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
944 submdspan(phi, std::pair(1, tdim + 1),
945 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
946 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
948 const std::function<void(
const std::span<T>&,
949 const std::span<const std::uint32_t>&,
951 apply_inverse_transpose_dof_transformation
952 = element->template get_dof_transformation_function<T>(
true,
true);
955 const auto [_Pi, pi_shape] = element->interpolation_operator();
956 cmdspan2_t Pi(_Pi.data(), pi_shape);
958 using u_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
959 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
960 using U_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
961 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
962 using J_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
963 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
964 using K_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
965 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
967 = element->basix_element().template map_fn<U_t, u_t, J_t, K_t>();
969 for (std::size_t c = 0; c < cells.size(); ++c)
971 const std::int32_t
cell = cells[c];
972 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::
973 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
974 x_dofmap,
cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
975 for (
int i = 0; i < num_dofs_g; ++i)
977 const int pos = 3 * x_dofs[i];
978 for (
int j = 0; j < gdim; ++j)
979 coord_dofs(i, j) = x_g[pos + j];
983 std::fill(J_b.begin(), J_b.end(), 0);
984 for (std::size_t p = 0; p < Xshape[0]; ++p)
987 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
988 submdspan(dphi, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, p,
989 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
991 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
992 submdspan(J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
993 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
996 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
997 submdspan(K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
998 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1003 std::span<const std::int32_t> dofs = dofmap->cell_dofs(
cell);
1004 for (
int k = 0; k < element_bs; ++k)
1007 for (
int m = 0; m < value_size; ++m)
1009 for (std::size_t k0 = 0; k0 < Xshape[0]; ++k0)
1012 = f[f_shape1 * (k * value_size + m) + c * Xshape[0] + k0];
1017 for (std::size_t i = 0; i < Xshape[0]; ++i)
1020 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
1022 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1023 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1025 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
1026 submdspan(ref_data, i,
1027 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1028 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1030 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
1031 submdspan(K, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1032 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1034 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
1035 submdspan(J, i, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1036 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1037 pull_back_fn(_U, _u, _K, 1.0 / detJ[i], _J);
1041 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
1042 submdspan(ref_data, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1043 0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1044 impl::interpolation_apply(Pi, ref, std::span(_coeffs), element_bs);
1045 apply_inverse_transpose_dof_transformation(_coeffs, cell_info,
cell, 1);
1048 assert(_coeffs.size() == num_scalar_dofs);
1049 for (
int i = 0; i < num_scalar_dofs; ++i)
1051 const int dof = i * element_bs + k;
1052 std::div_t pos = std::div(dof, dofmap_bs);
1053 coeffs[dofmap_bs * dofs[pos.quot] + pos.rem] = _coeffs[i];
1077template <std::
floating_po
int T>
1078std::tuple<std::vector<std::int32_t>, std::vector<std::int32_t>, std::vector<T>,
1079 std::vector<std::int32_t>>
1082 const mesh::Mesh<T>& mesh1, std::span<const std::int32_t> cells, T padding)
1086 const std::vector<T> coords
1090 std::vector<T> x(coords.size());
1091 std::size_t num_points = coords.size() / 3;
1092 for (std::size_t i = 0; i < num_points; ++i)
1093 for (std::size_t j = 0; j < 3; ++j)
1094 x[3 * i + j] = coords[i + j * num_points];
1097 return geometry::determine_point_ownership<T>(mesh1, x, padding);
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>>
1121 int tdim = mesh0.
topology()->dim();
1122 auto cell_map = mesh0.
topology()->index_map(tdim);
1124 std::int32_t num_cells = cell_map->size_local() + cell_map->num_ghosts();
1125 std::vector<std::int32_t> cells(num_cells, 0);
1126 std::iota(cells.begin(), cells.end(), 0);
1128 mesh0.
geometry(), element0, mesh1, cells, padding);
1138template <dolfinx::scalar T, std::
floating_po
int U>
1141 std::span<const std::int32_t> cells,
1142 const std::tuple<std::vector<std::int32_t>, std::vector<std::int32_t>,
1143 std::vector<U>, std::vector<std::int32_t>>&
1144 nmm_interpolation_data
1152 auto cell_map0 = mesh->topology()->index_map(mesh->topology()->dim());
1154 std::size_t num_cells0 = cell_map0->size_local() + cell_map0->num_ghosts();
1158 std::span<T> u1_array = u.
x()->mutable_array();
1159 std::span<const T> u0_array = v.
x()->array();
1160 std::copy(u0_array.begin(), u0_array.end(), u1_array.begin());
1167 impl::interpolate_nonmatching_meshes(u, v, cells, nmm_interpolation_data);
1178 if (element0->value_shape().size() != element1->value_shape().size()
1179 or !std::equal(element0->value_shape().begin(),
1180 element0->value_shape().end(),
1181 element1->value_shape().begin()))
1183 throw std::runtime_error(
1184 "Interpolation: elements have different value dimensions");
1187 if (element1 == element0 or *element1 == *element0)
1191 const int tdim = mesh->topology()->dim();
1192 auto cell_map = mesh->topology()->index_map(tdim);
1195 assert(element1->block_size() == element0->block_size());
1198 std::shared_ptr<const DofMap> dofmap0 = v.
function_space()->dofmap();
1200 std::shared_ptr<const DofMap> dofmap1 = u.
function_space()->dofmap();
1203 std::span<T> u1_array = u.
x()->mutable_array();
1204 std::span<const T> u0_array = v.
x()->array();
1207 const int bs0 = dofmap0->bs();
1208 const int bs1 = dofmap1->bs();
1209 for (
auto c :
cells)
1211 std::span<const std::int32_t> dofs0 = dofmap0->cell_dofs(c);
1212 std::span<const std::int32_t> dofs1 = dofmap1->cell_dofs(c);
1213 assert(bs0 * dofs0.size() == bs1 * dofs1.size());
1214 for (std::size_t i = 0; i < dofs0.size(); ++i)
1216 for (
int k = 0; k < bs0; ++k)
1218 int index = bs0 * i + k;
1219 std::div_t dv1 = std::div(index, bs1);
1220 u1_array[bs1 * dofs1[dv1.quot] + dv1.rem]
1221 = u0_array[bs0 * dofs0[i] + k];
1226 else if (element1->map_type() == element0->map_type())
1229 impl::interpolate_same_map(u, v, cells);
1234 impl::interpolate_nonmatching_maps(u, v, cells);
Degree-of-freedom map representations and tools.
A CoordinateElement manages coordinate mappings for isoparametric cells.
Definition CoordinateElement.h:39
static void compute_jacobian(const U &dphi, const V &cell_geometry, W &&J)
Compute Jacobian for a cell with given geometry using the basis functions and first order derivatives...
Definition CoordinateElement.h:106
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:53
static void compute_jacobian_inverse(const U &J, V &&K)
Compute the inverse of the Jacobian.
Definition CoordinateElement.h:115
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:46
int dim() const
The dimension of the coordinate element space.
Definition CoordinateElement.cpp:193
static double compute_jacobian_determinant(const U &J, std::span< typename U::value_type > w)
Compute the determinant of the Jacobian.
Definition CoordinateElement.h:132
Finite Element, containing the dof layout on a reference element, and various methods for evaluating ...
Definition FiniteElement.h:27
std::pair< std::vector< T >, std::array< std::size_t, 2 > > interpolation_points() const
Points on the reference cell at which an expression needs to be evaluated in order to interpolate the...
Definition FiniteElement.cpp:481
This class represents a function in a finite element function space , given by.
Definition Function.h:45
std::shared_ptr< const FunctionSpace< U > > function_space() const
Access the function space.
Definition Function.h:140
std::shared_ptr< const la::Vector< T > > x() const
Underlying vector.
Definition Function.h:146
Geometry stores the geometry imposed on a mesh.
Definition Geometry.h:33
int dim() const
Return Euclidean dimension of coordinate system.
Definition Geometry.h:88
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const std::int32_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents< std::size_t, 2 > > dofmap() const
DOF map.
Definition Geometry.h:94
std::span< const value_type > x() const
Access geometry degrees-of-freedom data (const version).
Definition Geometry.h:113
const std::vector< fem::CoordinateElement< value_type > > & cmaps() const
The elements that describes the geometry maps.
Definition Geometry.h:125
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
void cells(la::SparsityPattern &pattern, std::span< const std::int32_t > cells, std::array< std::reference_wrapper< const DofMap >, 2 > dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition sparsitybuild.cpp:18
Finite element method functionality.
Definition assemble_matrix_impl.h:25
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: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 expression f(x) in a finite element space.
Definition interpolate.h:748
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:42