DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
interpolate.h
1// Copyright (C) 2020-2022 Garth N. Wells, Igor A. Baratta, Massimiliano Leoni
2// and Jørgen S.Dokken
3//
4// This file is part of DOLFINx (https://www.fenicsproject.org)
5//
6// SPDX-License-Identifier: LGPL-3.0-or-later
7
8#pragma once
9
10#include "CoordinateElement.h"
11#include "DofMap.h"
12#include "FiniteElement.h"
13#include "FunctionSpace.h"
14#include <basix/mdspan.hpp>
15#include <concepts>
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>
21#include <functional>
22#include <numeric>
23#include <span>
24#include <vector>
25
26namespace dolfinx::fem
27{
28template <dolfinx::scalar T, std::floating_point U>
29class Function;
30
31template <typename T>
32concept MDSpan = requires(T x, std::size_t idx) {
33 x(idx, idx);
34 { x.extent(0) } -> std::integral;
35
36 { x.extent(1) } -> std::integral;
37};
38
49template <std::floating_point T>
50std::vector<T> interpolation_coords(const fem::FiniteElement<T>& element,
51 const mesh::Geometry<T>& geometry,
52 std::span<const std::int32_t> cells)
53{
54 // Get geometry data and the element coordinate map
55 const std::size_t gdim = geometry.dim();
56 auto x_dofmap = geometry.dofmap();
57 std::span<const T> x_g = geometry.x();
58
59 const CoordinateElement<T>& cmap = geometry.cmap();
60 const std::size_t num_dofs_g = cmap.dim();
61
62 // Get the interpolation points on the reference cells
63 const auto [X, Xshape] = element.interpolation_points();
64
65 // Evaluate coordinate element basis at reference points
66 std::array<std::size_t, 4> phi_shape = cmap.tabulate_shape(0, Xshape[0]);
67 std::vector<T> phi_b(
68 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
69 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
70 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>
71 phi_full(phi_b.data(), phi_shape);
72 cmap.tabulate(0, X, Xshape, phi_b);
73 auto phi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
74 phi_full, 0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
75 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
76
77 // Push reference coordinates (X) forward to the physical coordinates
78 // (x) for each cell
79 std::vector<T> coordinate_dofs(num_dofs_g * gdim, 0);
80 std::vector<T> x(3 * (cells.size() * Xshape[0]), 0);
81 for (std::size_t c = 0; c < cells.size(); ++c)
82 {
83 // Get geometry data for current cell
84 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
85 x_dofmap, cells[c], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
86 for (std::size_t i = 0; i < x_dofs.size(); ++i)
87 {
88 std::copy_n(std::next(x_g.begin(), 3 * x_dofs[i]), gdim,
89 std::next(coordinate_dofs.begin(), i * gdim));
90 }
91
92 // Push forward coordinates (X -> x)
93 for (std::size_t p = 0; p < Xshape[0]; ++p)
94 {
95 for (std::size_t j = 0; j < gdim; ++j)
96 {
97 T acc = 0;
98 for (std::size_t k = 0; k < num_dofs_g; ++k)
99 acc += phi(p, k) * coordinate_dofs[k * gdim + j];
100 x[j * (cells.size() * Xshape[0]) + c * Xshape[0] + p] = acc;
101 }
102 }
103 }
104
105 return x;
106}
107
122template <dolfinx::scalar T, std::floating_point 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);
126
127namespace impl
128{
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>>;
133
153template <dolfinx::scalar T>
154void scatter_values(
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>>
159 send_values,
160 std::span<T> recv_values)
161{
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);
165
166 // Build unique set of the sorted src_ranks
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()),
170 out_ranks.end());
171 out_ranks.reserve(out_ranks.size() + 1);
172
173 // Remove negative entries from dest_ranks
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; });
179
180 // Create unique set of sorted in-ranks
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);
184
185 // Create neighborhood communicator
186 MPI_Comm reverse_comm;
187 MPI_Dist_graph_create_adjacent(
188 comm, in_ranks.size(), in_ranks.data(), MPI_UNWEIGHTED, out_ranks.size(),
189 out_ranks.data(), MPI_UNWEIGHTED, MPI_INFO_NULL, false, &reverse_comm);
190
191 std::vector<std::int32_t> comm_to_output;
192 std::vector<std::int32_t> recv_sizes(in_ranks.size());
193 recv_sizes.reserve(1);
194 std::vector<std::int32_t> recv_offsets(in_ranks.size() + 1, 0);
195 {
196 // Build map from parent to neighborhood communicator ranks
197 std::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;
200
201 // Compute receive sizes
202 std::for_each(
203 dest_ranks.begin(), dest_ranks.end(),
204 [&dest_ranks, &rank_to_neighbor, &recv_sizes, block_size](auto rank)
205 {
206 if (rank >= 0)
207 {
208 const int neighbor = rank_to_neighbor[rank];
209 recv_sizes[neighbor] += block_size;
210 }
211 });
212
213 // Compute receiving offsets
214 std::partial_sum(recv_sizes.begin(), recv_sizes.end(),
215 std::next(recv_offsets.begin(), 1));
216
217 // Compute map from receiving values to position in recv_values
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)
221 {
222 if (const std::int32_t rank = dest_ranks[i]; rank >= 0)
223 {
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;
228 }
229 }
230 }
231
232 std::vector<std::int32_t> send_sizes(out_ranks.size());
233 send_sizes.reserve(1);
234 {
235 // Compute map from parent mpi rank to neigbor rank for outgoing data
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;
239
240 // Compute send sizes
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; });
244 }
245
246 // Compute sending offsets
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));
250
251 // Send values to dest ranks
252 std::vector<T> values(recv_offsets.back());
253 values.reserve(1);
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);
259
260 // Insert values received from neighborhood communicator in output span
261 std::fill(recv_values.begin(), recv_values.end(), T(0));
262 for (std::size_t i = 0; i < comm_to_output.size(); i++)
263 {
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);
267 }
268};
269
278template <MDSpan U, MDSpan V, dolfinx::scalar T>
279void interpolation_apply(U&& Pi, V&& data, std::span<T> coeffs, int bs)
280{
281 using X = typename dolfinx::scalar_value_type_t<T>;
282
283 // Compute coefficients = Pi * x (matrix-vector multiply)
284 if (bs == 1)
285 {
286 assert(data.extent(0) * data.extent(1) == Pi.extent(1));
287 for (std::size_t i = 0; i < Pi.extent(0); ++i)
288 {
289 coeffs[i] = 0.0;
290 for (std::size_t k = 0; k < data.extent(1); ++k)
291 for (std::size_t j = 0; j < data.extent(0); ++j)
292 coeffs[i]
293 += static_cast<X>(Pi(i, k * data.extent(0) + j)) * data(j, k);
294 }
295 }
296 else
297 {
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)
302 {
303 for (std::size_t i = 0; i < Pi.extent(0); ++i)
304 {
305 T acc = 0;
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;
309 }
310 }
311 }
312}
313
325template <dolfinx::scalar T, std::floating_point 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)
329{
330 auto V0 = u0.function_space();
331 assert(V0);
332 auto V1 = u1.function_space();
333 assert(V1);
334 auto mesh0 = V0->mesh();
335 assert(mesh0);
336
337 auto mesh1 = V1->mesh();
338 assert(mesh1);
339
340 auto element0 = V0->element();
341 assert(element0);
342 auto element1 = V1->element();
343 assert(element1);
344
345 assert(mesh0->topology()->dim());
346 const int tdim = mesh0->topology()->dim();
347 auto map = mesh0->topology()->index_map(tdim);
348 assert(map);
349 std::span<T> u1_array = u1.x()->mutable_array();
350 std::span<const T> u0_array = u0.x()->array();
351
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())
356 {
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());
361 }
362
363 // Get dofmaps
364 auto dofmap1 = V1->dofmap();
365 auto dofmap0 = V0->dofmap();
366
367 // Get block sizes and dof transformation operators
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>(
375
376 // Create working array
377 std::vector<T> local0(element0->space_dimension());
378 std::vector<T> local1(element1->space_dimension());
379
380 // Create interpolation operator
381 const auto [i_m, im_shape]
382 = element1->create_interpolation_operator(*element0);
383
384 // Iterate over mesh and interpolate on each cell
385 using X = typename dolfinx::scalar_value_type_t<T>;
386 for (std::size_t c = 0; c < cells0.size(); c++)
387 {
388 // Pack and transform cell dofs to reference ordering
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];
393
394 apply_dof_transformation(local0, cell_info0, cells0[c], 1);
395
396 // FIXME: Get compile-time ranges from Basix
397 // Apply interpolation operator
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];
402
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];
408 }
409}
410
422template <dolfinx::scalar T, std::floating_point 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)
426{
427 // Get mesh
428 auto V0 = u0.function_space();
429 assert(V0);
430 auto mesh0 = V0->mesh();
431 assert(mesh0);
432
433 // Mesh dims
434 const int tdim = mesh0->topology()->dim();
435 const int gdim = mesh0->geometry().dim();
436
437 // Get elements
438 auto V1 = u1.function_space();
439 assert(V1);
440 auto mesh1 = V1->mesh();
441 assert(mesh1);
442 auto element0 = V0->element();
443 assert(element0);
444 auto element1 = V1->element();
445 assert(element1);
446
447 std::span<const std::uint32_t> cell_info0;
448 std::span<const std::uint32_t> cell_info1;
449
450 if (element1->needs_dof_transformations()
451 or element0->needs_dof_transformations())
452 {
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());
457 }
458
459 // Get dofmaps
460 auto dofmap0 = V0->dofmap();
461 auto dofmap1 = V1->dofmap();
462
463 const auto [X, Xshape] = element1->interpolation_points();
464
465 // Get block sizes and dof transformation operators
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>(
473
474 // Get sizes of elements
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;
478
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();
483
484 // Evaluate coordinate map basis at reference interpolation points
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);
491
492 // Evaluate v basis functions at reference interpolation points
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);
497
498 // Create working arrays
499 std::vector<T> local1(element1->space_dimension());
500 std::vector<T> coeffs0(element0->space_dimension());
501
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);
504
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);
508
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,
511 V1->value_size());
512
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,
515 V1->value_size());
516
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);
519
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);
526
527 // Get interpolation operator
528 const auto [_Pi_1, pi_shape] = element1->interpolation_operator();
529 impl::mdspan_t<const U, 2> Pi_1(_Pi_1.data(), pi_shape);
530
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>();
540
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>>;
545 auto pull_back_fn1
546 = element1->basix_element().template map_fn<V_t, v_t, K_t, J_t>();
547
548 // Iterate over mesh and interpolate on each cell
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++)
552 {
553 // Get cell geometry (coordinate dofs)
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)
557 {
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];
561 }
562
563 // Compute Jacobians and reference points for current cell
564 std::fill(J_b.begin(), J_b.end(), 0);
565 for (std::size_t p = 0; p < Xshape[0]; ++p)
566 {
567 auto dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
568 phi, std::pair(1, tdim + 1), p,
569 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
570
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);
580 }
581
582 // Copy evaluated basis on reference, apply DOF transformations, and
583 // push forward to physical element
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);
589
590 for (std::size_t p = 0; p < Xshape[0]; ++p)
591 {
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);
596 }
597
598 for (std::size_t i = 0; i < basis0.extent(0); ++i)
599 {
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);
613 }
614
615 // Copy expansion coefficients for v into local array
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];
621
622 // Evaluate v at the interpolation points (physical space values)
623 using X = typename dolfinx::scalar_value_type_t<T>;
624 for (std::size_t p = 0; p < Xshape[0]; ++p)
625 {
626 for (int k = 0; k < bs0; ++k)
627 {
628 for (std::size_t j = 0; j < value_size0; ++j)
629 {
630 T acc = 0;
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;
634 }
635 }
636 }
637
638 // Pull back the physical values to the u reference
639 for (std::size_t i = 0; i < values0.extent(0); ++i)
640 {
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);
654 }
655
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);
661
662 // Copy local coefficients to the correct position in u dof array
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];
668 }
669}
670
671template <dolfinx::scalar T, std::floating_point 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)
678{
679 auto mesh = u.function_space()->mesh();
680 assert(mesh);
681 MPI_Comm comm = mesh->comm();
682
683 {
684 auto mesh_v = v.function_space()->mesh();
685 assert(mesh_v);
686 int result;
687 MPI_Comm_compare(comm, mesh_v->comm(), &result);
688 if (result == MPI_UNEQUAL)
689 {
690 throw std::runtime_error("Interpolation on different meshes is only "
691 "supported with the same communicator.");
692 }
693 }
694
695 assert(mesh->topology());
696 auto cell_map = mesh->topology()->index_map(mesh->topology()->dim());
697 assert(cell_map);
698
699 auto element_u = u.function_space()->element();
700 assert(element_u);
701 const std::size_t value_size = u.function_space()->value_size();
702
703 const auto& [dest_ranks, src_ranks, recv_points, evaluation_cells]
704 = nmm_interpolation_data;
705
706 // Evaluate the interpolating function where possible
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});
710
711 // Send values back to owning process
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));
718
719 // Transpose received data
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);
730
731 // Call local interpolation operator
732 fem::interpolate<T>(u, valuesT_b, {valuesT.extent(0), valuesT.extent(1)},
733 cells);
734}
735//----------------------------------------------------------------------------
736} // namespace impl
737
738template <dolfinx::scalar T, std::floating_point U>
739void interpolate(Function<T, U>& u, std::span<const T> f,
740 std::array<std::size_t, 2> fshape,
741 std::span<const std::int32_t> cells)
742{
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();
752 assert(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)
756 {
757 throw std::runtime_error("Cannot directly interpolate a mixed space. "
758 "Interpolate into subspaces.");
759 }
760
761 // Get mesh
762 assert(u.function_space());
763 auto mesh = u.function_space()->mesh();
764 assert(mesh);
765
766 const int gdim = mesh->geometry().dim();
767 const int tdim = mesh->topology()->dim();
768 const bool symmetric = u.function_space()->symmetric();
769
770 if (fshape[0] != (std::size_t)u.function_space()->value_size())
771 throw std::runtime_error("Interpolation data has the wrong shape/size.");
772
773 std::span<const std::uint32_t> cell_info;
774 if (element->needs_dof_transformations())
775 {
776 mesh->topology_mutable()->create_entity_permutations();
777 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
778 }
779
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);
784
785 // Get dofmap
786 const auto dofmap = u.function_space()->dofmap();
787 assert(dofmap);
788 const int dofmap_bs = dofmap->bs();
789
790 // Loop over cells and compute interpolation dofs
791 const int num_scalar_dofs = element->space_dimension() / element_bs;
792 const int value_size = u.function_space()->value_size() / element_bs;
793
794 std::span<T> coeffs = u.x()->mutable_array();
795 std::vector<T> _coeffs(num_scalar_dofs);
796
797 // This assumes that any element with an identity interpolation matrix
798 // is a point evaluation
799 if (element->map_ident() && element->interpolation_ident())
800 {
801 // Point evaluation element *and* the geometric map is the identity,
802 // e.g. not Piola mapped
803
804 auto apply_inv_transpose_dof_transformation
805 = element->template dof_transformation_fn<T>(
806 doftransform::inverse_transpose, true);
807
808 if (symmetric)
809 {
810 std::size_t matrix_size = 0;
811 while (matrix_size * matrix_size < fshape[0])
812 ++matrix_size;
813 // Loop over cells
814 for (std::size_t c = 0; c < cells.size(); ++c)
815 {
816 // The entries of a symmetric matrix are numbered (for an example 4x4
817 // element):
818 // 0 * * *
819 // 1 2 * *
820 // 3 4 5 *
821 // 6 7 8 9
822 // The loop extracts these elements. In this loop, row is the row of
823 // this matrix, and (k - rowstart) is the column
824 std::size_t row = 0;
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)
829 {
830 if (k - rowstart > row)
831 {
832 ++row;
833 rowstart = k;
834 }
835 // num_scalar_dofs is the number of interpolation points per
836 // cell in this case (interpolation matrix is identity)
837 std::copy_n(
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)
843 {
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];
847 }
848 }
849 }
850 }
851 else
852 {
853 // Loop over cells
854 for (std::size_t c = 0; c < cells.size(); ++c)
855 {
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)
859 {
860 // num_scalar_dofs is the number of interpolation points per
861 // cell in this case (interpolation matrix is identity)
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)
866 {
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];
870 }
871 }
872 }
873 }
874 }
875 else if (element->map_ident())
876 {
877 // Not a point evaluation, but the geometric map is the identity,
878 // e.g. not Piola mapped
879
880 if (symmetric)
881 {
882 throw std::runtime_error(
883 "Interpolation into this element not supported.");
884 }
885
886 const int element_vs = u.function_space()->value_size() / element_bs;
887
888 if (element_vs > 1 && element_bs > 1)
889 {
890 throw std::runtime_error(
891 "Interpolation into this element not supported.");
892 }
893
894 // Get interpolation operator
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);
899
900 auto apply_inv_transpose_dof_transformation
901 = element->template dof_transformation_fn<T>(
902 doftransform::inverse_transpose, true);
903
904 // Loop over cells
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)
911 {
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)
915 {
916 for (int i = 0; i < element_vs; ++i)
917 {
918 std::copy_n(
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));
924 }
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)
928 {
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];
932 }
933 }
934 }
935 }
936 else
937 {
938 if (symmetric)
939 {
940 throw std::runtime_error(
941 "Interpolation into this element not supported.");
942 }
943 // Get the interpolation points on the reference cells
944 const auto [X, Xshape] = element->interpolation_points();
945 if (X.empty())
946 {
947 throw std::runtime_error(
948 "Interpolation into this space is not yet supported.");
949 }
950
951 if (_f.extent(1) != cells.size() * Xshape[0])
952 throw std::runtime_error("Interpolation data has the wrong shape.");
953
954 // Get coordinate map
955 const CoordinateElement<U>& cmap = mesh->geometry().cmap();
956
957 // Get geometry data
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();
961
962 // Create data structures for Jacobian info
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);
969
970 std::vector<U> coord_dofs_b(num_dofs_g * gdim);
971 mdspan2_t coord_dofs(coord_dofs_b.data(), num_dofs_g, gdim);
972
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);
977
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);
982
983 // Tabulate 1st derivative of shape functions at interpolation
984 // coords
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);
989 cmap.tabulate(1, X, Xshape, phi_b);
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);
994
995 const std::function<void(std::span<T>, std::span<const std::uint32_t>,
996 std::int32_t, int)>
997 apply_inverse_transpose_dof_transformation
998 = element->template dof_transformation_fn<T>(
999 doftransform::inverse_transpose);
1000
1001 // Get interpolation operator
1002 const auto [_Pi, pi_shape] = element->interpolation_operator();
1003 cmdspan2_t Pi(_Pi.data(), pi_shape);
1004
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>>;
1013 auto pull_back_fn
1014 = element->basix_element().template map_fn<U_t, u_t, J_t, K_t>();
1015
1016 for (std::size_t c = 0; c < cells.size(); ++c)
1017 {
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)
1022 {
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];
1026 }
1027
1028 // Compute J, detJ and K
1029 std::fill(J_b.begin(), J_b.end(), 0);
1030 for (std::size_t p = 0; p < Xshape[0]; ++p)
1031 {
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);
1038 cmap.compute_jacobian(_dphi, coord_dofs, _J);
1039 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
1040 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
1041 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
1042 cmap.compute_jacobian_inverse(_J, _K);
1043 detJ[p] = cmap.compute_jacobian_determinant(_J, det_scratch);
1044 }
1045
1046 std::span<const std::int32_t> dofs = dofmap->cell_dofs(cell);
1047 for (int k = 0; k < element_bs; ++k)
1048 {
1049 // Extract computed expression values for element block k
1050 for (int m = 0; m < value_size; ++m)
1051 {
1052 for (std::size_t k0 = 0; k0 < Xshape[0]; ++k0)
1053 {
1054 _vals(k0, 0, m)
1055 = f[f_shape1 * (k * value_size + m) + c * Xshape[0] + k0];
1056 }
1057 }
1058
1059 // Get element degrees of freedom for block
1060 for (std::size_t i = 0; i < Xshape[0]; ++i)
1061 {
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);
1075 }
1076
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);
1082
1083 // Copy interpolation dofs into coefficient vector
1084 assert(_coeffs.size() == num_scalar_dofs);
1085 for (int i = 0; i < num_scalar_dofs; ++i)
1086 {
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];
1090 }
1091 }
1092 }
1093 }
1094}
1095
1113template <std::floating_point T>
1114std::tuple<std::vector<std::int32_t>, std::vector<std::int32_t>, std::vector<T>,
1115 std::vector<std::int32_t>>
1117 const mesh::Geometry<T>& geometry0, const FiniteElement<T>& element0,
1118 const mesh::Mesh<T>& mesh1, std::span<const std::int32_t> cells, T padding)
1119{
1120 // Collect all the points at which values are needed to define the
1121 // interpolating function
1122 const std::vector<T> coords
1123 = interpolation_coords(element0, geometry0, cells);
1124
1125 // Transpose interpolation 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];
1131
1132 // Determine ownership of each point
1133 return geometry::determine_point_ownership<T>(mesh1, x, padding);
1134}
1135
1149template <std::floating_point T>
1150std::tuple<std::vector<std::int32_t>, std::vector<std::int32_t>, std::vector<T>,
1151 std::vector<std::int32_t>>
1153 const FiniteElement<T>& element0,
1154 const mesh::Mesh<T>& mesh1,
1155 T padding)
1156{
1157 int tdim = mesh0.topology()->dim();
1158 auto cell_map = mesh0.topology()->index_map(tdim);
1159 assert(cell_map);
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);
1165}
1166
1177template <dolfinx::scalar T, std::floating_point U>
1179 Function<T, U>& u1, const Function<T, U>& u0,
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
1185 = {})
1186{
1187 assert(u1.function_space());
1188 assert(u0.function_space());
1189 auto mesh = u1.function_space()->mesh();
1190 assert(mesh);
1191
1192 auto cell_map0 = mesh->topology()->index_map(mesh->topology()->dim());
1193 assert(cell_map0);
1194 std::size_t num_cells0 = cell_map0->size_local() + cell_map0->num_ghosts();
1195 if (u1.function_space() == u0.function_space()
1196 and cells1.size() == num_cells0)
1197 {
1198 // Same function spaces and on whole mesh
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());
1202 }
1203 else
1204 {
1205 std::vector<std::int32_t> cells0;
1206 cells0.reserve(cells1.size());
1207 // Get mesh and check that functions share the same mesh
1208 if (auto mesh_v = u0.function_space()->mesh(); mesh == mesh_v)
1209 {
1210 cells0.insert(cells0.end(), cells1.begin(), cells1.end());
1211 }
1212 // If meshes are different and input mapping is given
1213 else if (cell_map.size() > 0)
1214 {
1215 std::transform(cells1.begin(), cells1.end(), std::back_inserter(cells0),
1216 [&cell_map](std::int32_t c) { return cell_map[c]; });
1217 }
1218 // Non-matching meshes
1219 if (cells0.empty())
1220 {
1221 impl::interpolate_nonmatching_meshes(u1, u0, cells1,
1222 nmm_interpolation_data);
1223 return;
1224 }
1225
1226 // Get elements and check value shape
1227 auto fs0 = u0.function_space();
1228 auto element0 = fs0->element();
1229 assert(element0);
1230 auto fs1 = u1.function_space();
1231 auto element1 = fs1->element();
1232 assert(element1);
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()))
1236 {
1237 throw std::runtime_error(
1238 "Interpolation: elements have different value dimensions");
1239 }
1240
1241 if (element1 == element0 or *element1 == *element0)
1242 {
1243 // Same element, different dofmaps (or just a subset of cells)
1244
1245 const int tdim = mesh->topology()->dim();
1246 auto cell_map1 = mesh->topology()->index_map(tdim);
1247 assert(cell_map1);
1248
1249 assert(element1->block_size() == element0->block_size());
1250
1251 // Get dofmaps
1252 std::shared_ptr<const DofMap> dofmap0 = u0.function_space()->dofmap();
1253 assert(dofmap0);
1254 std::shared_ptr<const DofMap> dofmap1 = u1.function_space()->dofmap();
1255 assert(dofmap1);
1256
1257 std::span<T> u1_array = u1.x()->mutable_array();
1258 std::span<const T> u0_array = u0.x()->array();
1259
1260 // Iterate over mesh and interpolate on each cell
1261 const int bs0 = dofmap0->bs();
1262 const int bs1 = dofmap1->bs();
1263 for (std::size_t c = 0; c < cells1.size(); ++c)
1264 {
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)
1269 {
1270 for (int k = 0; k < bs0; ++k)
1271 {
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];
1276 }
1277 }
1278 }
1279 }
1280 else if (element1->map_type() == element0->map_type())
1281 {
1282 // Different elements, same basis function map type
1283 impl::interpolate_same_map(u1, u0, cells1, cells0);
1284 }
1285 else
1286 {
1287 // Different elements with different maps for basis functions
1288 impl::interpolate_nonmatching_maps(u1, u0, cells1, cells0);
1289 }
1290 }
1291}
1292} // namespace dolfinx::fem
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
Definition Function.h:45
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