10#include "DirichletBC.h"
16#include <basix/mdspan.hpp>
18#include <dolfinx/common/IndexMap.h>
19#include <dolfinx/mesh/Geometry.h>
20#include <dolfinx/mesh/Mesh.h>
21#include <dolfinx/mesh/Topology.h>
30template <dolfinx::scalar T, std::
floating_po
int U>
34namespace dolfinx::fem::impl
37using mdspan2_t = md::mdspan<const std::int32_t, md::dextents<std::size_t, 2>>;
79template <
int _bs0 = -1,
int _bs1 = -1,
typename V,
80 dolfinx::scalar T =
typename std::remove_cvref_t<V>::value_type>
81 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
83 V&& b, mdspan2_t x_dofmap,
84 md::mdspan<
const scalar_value_t<T>,
85 md::extents<std::size_t, md::dynamic_extent, 3>>
87 FEkernel<T>
auto kernel, std::span<const std::int32_t> cells,
88 std::tuple<mdspan2_t,
int, std::span<const std::int32_t>> dofmap0,
89 fem::DofTransformKernel<T>
auto P0,
90 std::tuple<mdspan2_t,
int, std::span<const std::int32_t>> dofmap1,
91 fem::DofTransformKernel<T>
auto P1T, std::span<const T> constants,
92 md::mdspan<
const T, md::dextents<std::size_t, 2>> coeffs,
93 std::span<const std::uint32_t> cell_info0,
94 std::span<const std::uint32_t> cell_info1, std::span<const T> bc_values1,
95 std::span<const std::int8_t> bc_markers1, std::span<const T> x0, T alpha)
100 const auto [dmap0, bs0, cells0] = dofmap0;
101 const auto [dmap1, bs1, cells1] = dofmap1;
102 assert(_bs0 < 0 or _bs0 == bs0);
103 assert(_bs1 < 0 or _bs1 == bs1);
105 const int num_rows = bs0 * dmap0.extent(1);
106 const int num_cols = bs1 * dmap1.extent(1);
109 std::vector<scalar_value_t<T>> cdofs(3 * x_dofmap.extent(1));
110 std::vector<T> Ae(num_rows * num_cols), be(num_rows);
111 assert(cells0.size() ==
cells.size());
112 assert(cells1.size() ==
cells.size());
113 for (std::size_t index = 0; index <
cells.size(); ++index)
117 std::int32_t c =
cells[index];
118 std::int32_t c0 = cells0[index];
119 std::int32_t c1 = cells1[index];
122 auto dofs1 = md::submdspan(dmap1, c1, md::full_extent);
126 for (std::size_t j = 0; j < dofs1.size(); ++j)
128 if constexpr (_bs1 > 0)
130 for (
int k = 0; k < _bs1; ++k)
132 assert(_bs1 * dofs1[j] + k < (
int)bc_markers1.size());
133 if (bc_markers1[_bs1 * dofs1[j] + k])
142 for (
int k = 0; k < bs1; ++k)
144 assert(bs1 * dofs1[j] + k < (
int)bc_markers1.size());
145 if (bc_markers1[bs1 * dofs1[j] + k])
158 auto x_dofs = md::submdspan(x_dofmap, c, md::full_extent);
159 for (std::size_t i = 0; i < x_dofs.size(); ++i)
160 std::copy_n(&x(x_dofs[i], 0), 3, std::next(cdofs.begin(), 3 * i));
163 auto dofs0 = md::submdspan(dmap0, c0, md::full_extent);
165 std::ranges::fill(Ae, 0);
166 kernel(Ae.data(), &coeffs(index, 0), constants.data(), cdofs.data(),
167 nullptr,
nullptr,
nullptr);
168 P0(Ae, cell_info0, c0, num_cols);
169 P1T(Ae, cell_info1, c1, num_rows);
172 std::ranges::fill(be, 0);
173 for (std::size_t j = 0; j < dofs1.size(); ++j)
175 if constexpr (_bs1 > 0)
177 for (
int k = 0; k < _bs1; ++k)
179 const std::int32_t jj = _bs1 * dofs1[j] + k;
180 assert(jj < (
int)bc_markers1.size());
183 const T bc = bc_values1[jj];
184 const T _x0 = x0.empty() ? 0 : x0[jj];
187 for (
int m = 0; m < num_rows; ++m)
188 be[m] -= Ae[m * num_cols + _bs1 * j + k] * alpha * (bc - _x0);
194 for (
int k = 0; k < bs1; ++k)
196 const std::int32_t jj = bs1 * dofs1[j] + k;
197 assert(jj < (
int)bc_markers1.size());
200 const T bc = bc_values1[jj];
201 const T _x0 = x0.empty() ? 0 : x0[jj];
203 for (
int m = 0; m < num_rows; ++m)
204 be[m] -= Ae[m * num_cols + bs1 * j + k] * alpha * (bc - _x0);
210 for (std::size_t i = 0; i < dofs0.size(); ++i)
212 if constexpr (_bs0 > 0)
214 for (
int k = 0; k < _bs0; ++k)
215 b[_bs0 * dofs0[i] + k] += be[_bs0 * i + k];
219 for (
int k = 0; k < bs0; ++k)
220 b[bs0 * dofs0[i] + k] += be[bs0 * i + k];
264 dolfinx::scalar T =
typename std::remove_cvref_t<V>::value_type>
265 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
266void _lift_bc_entities(
267 V&& b, mdspan2_t x_dofmap,
268 md::mdspan<
const scalar_value_t<T>,
269 md::extents<std::size_t, md::dynamic_extent, 3>>
271 FEkernel<T>
auto kernel,
272 md::mdspan<
const std::int32_t,
273 md::extents<std::size_t, md::dynamic_extent, 2>>
275 std::tuple<mdspan2_t,
int,
276 md::mdspan<
const std::int32_t,
277 md::extents<std::size_t, md::dynamic_extent, 2>>>
279 fem::DofTransformKernel<T>
auto P0,
280 std::tuple<mdspan2_t,
int,
281 md::mdspan<
const std::int32_t,
282 md::extents<std::size_t, md::dynamic_extent, 2>>>
284 fem::DofTransformKernel<T>
auto P1T, std::span<const T> constants,
285 md::mdspan<
const T, md::dextents<std::size_t, 2>> coeffs,
286 std::span<const std::uint32_t> cell_info0,
287 std::span<const std::uint32_t> cell_info1, std::span<const T> bc_values1,
288 std::span<const std::int8_t> bc_markers1, std::span<const T> x0, T alpha,
289 md::mdspan<
const std::uint8_t, md::dextents<std::size_t, 2>> perms)
291 if (entities.empty())
294 const auto [dmap0, bs0, entities0] = dofmap0;
295 const auto [dmap1, bs1, entities1] = dofmap1;
297 const int num_rows = bs0 * dmap0.extent(1);
298 const int num_cols = bs1 * dmap1.extent(1);
301 std::vector<scalar_value_t<T>> cdofs(3 * x_dofmap.extent(1));
302 std::vector<T> Ae(num_rows * num_cols), be(num_rows);
303 assert(entities0.size() == entities.size());
304 assert(entities1.size() == entities.size());
305 for (std::size_t index = 0; index < entities.extent(0); ++index)
309 std::int32_t
cell = entities(index, 0);
310 std::int32_t cell0 = entities0(index, 0);
311 std::int32_t cell1 = entities1(index, 0);
314 std::int32_t local_entity = entities(index, 1);
317 auto dofs1 = md::submdspan(dmap1, cell1, md::full_extent);
321 for (std::size_t j = 0; j < dofs1.size(); ++j)
323 for (
int k = 0; k < bs1; ++k)
325 if (bc_markers1[bs1 * dofs1[j] + k])
337 auto x_dofs = md::submdspan(x_dofmap,
cell, md::full_extent);
338 for (std::size_t i = 0; i < x_dofs.size(); ++i)
339 std::copy_n(&x(x_dofs[i], 0), 3, std::next(cdofs.begin(), 3 * i));
342 auto dofs0 = md::submdspan(dmap0, cell0, md::full_extent);
345 std::uint8_t perm = perms.empty() ? 0 : perms(
cell, local_entity);
346 std::ranges::fill(Ae, 0);
347 kernel(Ae.data(), &coeffs(index, 0), constants.data(), cdofs.data(),
348 &local_entity, &perm,
nullptr);
349 P0(Ae, cell_info0, cell0, num_cols);
350 P1T(Ae, cell_info1, cell1, num_rows);
353 std::ranges::fill(be, 0);
354 for (std::size_t j = 0; j < dofs1.size(); ++j)
356 for (
int k = 0; k < bs1; ++k)
358 const std::int32_t jj = bs1 * dofs1[j] + k;
361 const T bc = bc_values1[jj];
362 const T _x0 = x0.empty() ? 0 : x0[jj];
364 for (
int m = 0; m < num_rows; ++m)
365 be[m] -= Ae[m * num_cols + bs1 * j + k] * alpha * (bc - _x0);
370 for (std::size_t i = 0; i < dofs0.size(); ++i)
371 for (
int k = 0; k < bs0; ++k)
372 b[bs0 * dofs0[i] + k] += be[bs0 * i + k];
419 dolfinx::scalar T =
typename std::remove_cvref_t<V>::value_type>
420 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
421void _lift_bc_interior_facets(
422 V&& b, mdspan2_t x_dofmap,
423 md::mdspan<
const scalar_value_t<T>,
424 md::extents<std::size_t, md::dynamic_extent, 3>>
426 FEkernel<T>
auto kernel,
427 md::mdspan<
const std::int32_t,
428 md::extents<std::size_t, md::dynamic_extent, 2, 2>>
430 std::tuple<mdspan2_t,
int,
431 md::mdspan<
const std::int32_t,
432 md::extents<std::size_t, md::dynamic_extent, 2, 2>>>
434 fem::DofTransformKernel<T>
auto P0,
435 std::tuple<mdspan2_t,
int,
436 md::mdspan<
const std::int32_t,
437 md::extents<std::size_t, md::dynamic_extent, 2, 2>>>
439 fem::DofTransformKernel<T>
auto P1T, std::span<const T> constants,
440 md::mdspan<
const T, md::extents<std::size_t, md::dynamic_extent, 2,
443 std::span<const std::uint32_t> cell_info0,
444 std::span<const std::uint32_t> cell_info1, std::span<const T> bc_values1,
445 std::span<const std::int8_t> bc_markers1, std::span<const T> x0, T alpha,
446 md::mdspan<
const std::uint8_t, md::dextents<std::size_t, 2>> perms)
451 const auto [dmap0, bs0, facets0] = dofmap0;
452 const auto [dmap1, bs1, facets1] = dofmap1;
454 const int num_dofs0 = dmap0.extent(1);
455 const int num_dofs1 = dmap1.extent(1);
456 const int num_rows = bs0 * 2 * num_dofs0;
457 const int num_cols = bs1 * 2 * num_dofs1;
460 using X = scalar_value_t<T>;
461 std::vector<X> cdofs(2 * x_dofmap.extent(1) * 3);
462 std::span<X> cdofs0(cdofs.data(), x_dofmap.extent(1) * 3);
463 std::span<X> cdofs1(cdofs.data() + x_dofmap.extent(1) * 3,
464 x_dofmap.extent(1) * 3);
465 std::vector<T> Ae(num_rows * num_cols), be(num_rows);
468 std::vector<std::int32_t> dmapjoint0(2 * num_dofs0);
469 std::vector<std::int32_t> dmapjoint1(2 * num_dofs1);
471 assert(facets0.size() == facets.size());
472 assert(facets1.size() == facets.size());
473 for (std::size_t f = 0; f < facets.extent(0); ++f)
477 std::array
cells{facets(f, 0, 0), facets(f, 1, 0)};
478 std::array cells0{facets0(f, 0, 0), facets0(f, 1, 0)};
479 std::array cells1{facets1(f, 0, 0), facets1(f, 1, 0)};
482 std::array local_facet = {facets(f, 0, 1), facets(f, 1, 1)};
485 auto x_dofs0 = md::submdspan(x_dofmap, cells[0], md::full_extent);
486 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
487 std::copy_n(&x(x_dofs0[i], 0), 3, std::next(cdofs0.begin(), 3 * i));
488 auto x_dofs1 = md::submdspan(x_dofmap, cells[1], md::full_extent);
489 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
490 std::copy_n(&x(x_dofs1[i], 0), 3, std::next(cdofs1.begin(), 3 * i));
496 std::span dmap0_cell0
498 ? std::span(dmap0.data_handle() + cells0[0] * num_dofs0,
500 : std::span<const std::int32_t>();
501 std::span dmap0_cell1
503 ? std::span(dmap0.data_handle() + cells0[1] * num_dofs0,
505 : std::span<const std::int32_t>();
507 std::ranges::copy(dmap0_cell0, dmapjoint0.begin());
508 std::ranges::copy(dmap0_cell1, std::next(dmapjoint0.begin(), num_dofs0));
511 std::span<const std::int32_t> dmap1_cell0
513 ? std::span(dmap1.data_handle() + cells1[0] * num_dofs1,
515 : std::span<const std::int32_t>();
516 std::span<const std::int32_t> dmap1_cell1
518 ? std::span(dmap1.data_handle() + cells1[1] * num_dofs1,
520 : std::span<const std::int32_t>();
522 std::ranges::copy(dmap1_cell0, dmapjoint1.begin());
523 std::ranges::copy(dmap1_cell1, std::next(dmapjoint1.begin(), num_dofs1));
527 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
529 for (
int k = 0; k < bs1; ++k)
531 if (bc_markers1[bs1 * dmap1_cell0[j] + k])
540 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
542 for (
int k = 0; k < bs1; ++k)
544 if (bc_markers1[bs1 * dmap1_cell1[j] + k])
556 std::ranges::fill(Ae, 0);
557 std::array perm = perms.empty()
558 ? std::array<std::uint8_t, 2>{0, 0}
559 : std::array{perms(cells[0], local_facet[0]),
560 perms(cells[1], local_facet[1])};
561 kernel(Ae.data(), &coeffs(f, 0, 0), constants.data(), cdofs.data(),
562 local_facet.data(), perm.data(),
nullptr);
565 P0(Ae, cell_info0, cells0[0], num_cols);
568 std::span sub_Ae0(Ae.data() + bs0 * num_dofs0 * num_cols,
569 bs0 * num_dofs1 * num_cols);
570 P0(sub_Ae0, cell_info0, cells0[1], num_cols);
573 P1T(Ae, cell_info1, cells1[0], num_rows);
577 for (
int row = 0; row < num_rows; ++row)
581 std::span sub_Ae1(Ae.data() + row * num_cols + bs1 * num_dofs1,
583 P1T(sub_Ae1, cell_info1, cells1[1], 1);
587 std::ranges::fill(be, 0);
590 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
592 for (
int k = 0; k < bs1; ++k)
594 const std::int32_t jj = bs1 * dmap1_cell0[j] + k;
597 const T bc = bc_values1[jj];
598 const T _x0 = x0.empty() ? 0 : x0[jj];
600 for (
int m = 0; m < num_rows; ++m)
601 be[m] -= Ae[m * num_cols + bs1 * j + k] * alpha * (bc - _x0);
607 const int offset = bs1 * num_dofs1;
608 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
610 for (
int k = 0; k < bs1; ++k)
612 const std::int32_t jj = bs1 * dmap1_cell1[j] + k;
615 const T bc = bc_values1[jj];
616 const T _x0 = x0.empty() ? 0 : x0[jj];
618 for (
int m = 0; m < num_rows; ++m)
621 -= Ae[m * num_cols + offset + bs1 * j + k] * alpha * (bc - _x0);
627 for (std::size_t i = 0; i < dmap0_cell0.size(); ++i)
628 for (
int k = 0; k < bs0; ++k)
629 b[bs0 * dmap0_cell0[i] + k] += be[bs0 * i + k];
631 const int offset_be = bs0 * num_dofs0;
632 for (std::size_t i = 0; i < dmap0_cell1.size(); ++i)
633 for (
int k = 0; k < bs0; ++k)
634 b[bs0 * dmap0_cell1[i] + k] += be[offset_be + bs0 * i + k];
661template <
int _bs = -1,
typename V,
662 dolfinx::scalar T =
typename std::remove_cvref_t<V>::value_type>
663 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
666 fem::DofTransformKernel<T>
auto P0, V&& b, mdspan2_t x_dofmap,
667 md::mdspan<
const scalar_value_t<T>,
668 md::extents<std::size_t, md::dynamic_extent, 3>>
670 std::span<const std::int32_t> cells,
671 std::tuple<mdspan2_t,
int, std::span<const std::int32_t>> dofmap,
672 FEkernel<T>
auto kernel, std::span<const T> constants,
673 md::mdspan<
const T, md::dextents<std::size_t, 2>> coeffs,
674 std::span<const std::uint32_t> cell_info0)
679 const auto [dmap, bs, cells0] = dofmap;
680 assert(_bs < 0 or _bs == bs);
683 std::vector<scalar_value_t<T>> cdofs(3 * x_dofmap.extent(1));
684 std::vector<T> be(bs * dmap.extent(1));
687 for (std::size_t index = 0; index <
cells.size(); ++index)
690 std::int32_t c =
cells[index];
691 std::int32_t c0 = cells0[index];
694 auto x_dofs = md::submdspan(x_dofmap, c, md::full_extent);
695 for (std::size_t i = 0; i < x_dofs.size(); ++i)
696 std::copy_n(&x(x_dofs[i], 0), 3, std::next(cdofs.begin(), 3 * i));
699 std::ranges::fill(be, 0);
700 kernel(be.data(), &coeffs(index, 0), constants.data(), cdofs.data(),
701 nullptr,
nullptr,
nullptr);
702 P0(be, cell_info0, c0, 1);
705 auto dofs = md::submdspan(dmap, c0, md::full_extent);
706 if constexpr (_bs > 0)
708 for (std::size_t i = 0; i < dofs.size(); ++i)
709 for (
int k = 0; k < _bs; ++k)
710 b[_bs * dofs[i] + k] += be[_bs * i + k];
714 for (std::size_t i = 0; i < dofs.size(); ++i)
715 for (
int k = 0; k < bs; ++k)
716 b[bs * dofs[i] + k] += be[bs * i + k];
754template <
int _bs = -1,
typename V,
755 dolfinx::scalar T =
typename std::remove_cvref_t<V>::value_type>
756 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
757void assemble_entities(
758 fem::DofTransformKernel<T>
auto P0, V&& b, mdspan2_t x_dofmap,
759 md::mdspan<
const scalar_value_t<T>,
760 md::extents<std::size_t, md::dynamic_extent, 3>>
762 md::mdspan<
const std::int32_t,
763 std::extents<std::size_t, md::dynamic_extent, 2>>
765 std::tuple<mdspan2_t,
int,
766 md::mdspan<
const std::int32_t,
767 std::extents<std::size_t, md::dynamic_extent, 2>>>
769 FEkernel<T>
auto kernel, std::span<const T> constants,
770 md::mdspan<
const T, md::dextents<std::size_t, 2>> coeffs,
771 std::span<const std::uint32_t> cell_info0,
772 md::mdspan<
const std::uint8_t, md::dextents<std::size_t, 2>> perms)
774 if (entities.empty())
777 const auto [dmap, bs, entities0] = dofmap;
778 assert(_bs < 0 or _bs == bs);
781 const int num_dofs = dmap.extent(1);
782 std::vector<scalar_value_t<T>> cdofs(3 * x_dofmap.extent(1));
783 std::vector<T> be(bs * num_dofs);
784 assert(entities0.size() == entities.size());
785 for (std::size_t f = 0; f < entities.extent(0); ++f)
789 std::int32_t
cell = entities(f, 0);
790 std::int32_t local_entity = entities(f, 1);
791 std::int32_t cell0 = entities0(f, 0);
794 auto x_dofs = md::submdspan(x_dofmap,
cell, md::full_extent);
795 for (std::size_t i = 0; i < x_dofs.size(); ++i)
796 std::copy_n(&x(x_dofs[i], 0), 3, std::next(cdofs.begin(), 3 * i));
799 std::uint8_t perm = perms.empty() ? 0 : perms(
cell, local_entity);
802 std::ranges::fill(be, 0);
803 kernel(be.data(), &coeffs(f, 0), constants.data(), cdofs.data(),
804 &local_entity, &perm,
nullptr);
805 P0(be, cell_info0, cell0, 1);
808 auto dofs = md::submdspan(dmap, cell0, md::full_extent);
809 if constexpr (_bs > 0)
811 for (std::size_t i = 0; i < dofs.size(); ++i)
812 for (
int k = 0; k < _bs; ++k)
813 b[_bs * dofs[i] + k] += be[_bs * i + k];
817 for (std::size_t i = 0; i < dofs.size(); ++i)
818 for (
int k = 0; k < bs; ++k)
819 b[bs * dofs[i] + k] += be[bs * i + k];
849template <
int _bs = -1,
typename V,
850 dolfinx::scalar T =
typename std::remove_cvref_t<V>::value_type>
851 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
852void assemble_interior_facets(
853 fem::DofTransformKernel<T>
auto P0, V&& b, mdspan2_t x_dofmap,
854 md::mdspan<
const scalar_value_t<T>,
855 md::extents<std::size_t, md::dynamic_extent, 3>>
857 md::mdspan<
const std::int32_t,
858 std::extents<std::size_t, md::dynamic_extent, 2, 2>>
860 std::tuple<
const DofMap&,
int,
861 md::mdspan<
const std::int32_t,
862 std::extents<std::size_t, md::dynamic_extent, 2, 2>>>
864 FEkernel<T>
auto kernel, std::span<const T> constants,
865 md::mdspan<
const T, md::extents<std::size_t, md::dynamic_extent, 2,
868 std::span<const std::uint32_t> cell_info0,
869 md::mdspan<
const std::uint8_t, md::dextents<std::size_t, 2>> perms)
871 using X = scalar_value_t<T>;
876 const auto [dmap, bs, facets0] = dofmap;
877 assert(_bs < 0 or _bs == bs);
880 std::vector<X> cdofs(2 * x_dofmap.extent(1) * 3);
881 std::span<X> cdofs0(cdofs.data(), x_dofmap.extent(1) * 3);
882 std::span<X> cdofs1(cdofs.data() + x_dofmap.extent(1) * 3,
883 x_dofmap.extent(1) * 3);
885 const std::size_t dmap_size = dmap.map().extent(1);
886 std::vector<T> be(bs * 2 * dmap_size);
888 assert(facets0.size() == facets.size());
889 for (std::size_t f = 0; f < facets.extent(0); ++f)
892 std::array<std::int32_t, 2>
cells{facets(f, 0, 0), facets(f, 1, 0)};
893 std::array<std::int32_t, 2> cells0{facets0(f, 0, 0), facets0(f, 1, 0)};
896 std::array<std::int32_t, 2> local_facet{facets(f, 0, 1), facets(f, 1, 1)};
899 auto x_dofs0 = md::submdspan(x_dofmap, cells[0], md::full_extent);
900 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
901 std::copy_n(&x(x_dofs0[i], 0), 3, std::next(cdofs0.begin(), 3 * i));
902 auto x_dofs1 = md::submdspan(x_dofmap, cells[1], md::full_extent);
903 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
904 std::copy_n(&x(x_dofs1[i], 0), 3, std::next(cdofs1.begin(), 3 * i));
909 std::span dmap0 = cells0[0] >= 0 ? dmap.cell_dofs(cells0[0])
910 : std::span<const std::int32_t>();
911 std::span dmap1 = cells0[1] >= 0 ? dmap.cell_dofs(cells0[1])
912 : std::span<const std::int32_t>();
915 std::ranges::fill(be, 0);
916 std::array perm = perms.empty()
917 ? std::array<std::uint8_t, 2>{0, 0}
918 : std::array{perms(cells[0], local_facet[0]),
919 perms(cells[1], local_facet[1])};
920 kernel(be.data(), &coeffs(f, 0, 0), constants.data(), cdofs.data(),
921 local_facet.data(), perm.data(),
nullptr);
924 P0(be, cell_info0, cells0[0], 1);
927 std::span sub_be(be.data() + bs * dmap_size, bs * dmap_size);
928 P0(sub_be, cell_info0, cells0[1], 1);
932 if constexpr (_bs > 0)
934 for (std::size_t i = 0; i < dmap0.size(); ++i)
935 for (
int k = 0; k < _bs; ++k)
936 b[_bs * dmap0[i] + k] += be[_bs * i + k];
937 for (std::size_t i = 0; i < dmap1.size(); ++i)
938 for (
int k = 0; k < _bs; ++k)
939 b[_bs * dmap1[i] + k] += be[_bs * (i + dmap_size) + k];
943 for (std::size_t i = 0; i < dmap0.size(); ++i)
944 for (
int k = 0; k < bs; ++k)
945 b[bs * dmap0[i] + k] += be[bs * i + k];
946 for (std::size_t i = 0; i < dmap1.size(); ++i)
947 for (
int k = 0; k < bs; ++k)
948 b[bs * dmap1[i] + k] += be[bs * (i + dmap_size) + k];
969template <
typename V, std::floating_point U,
970 dolfinx::scalar T =
typename std::remove_cvref_t<V>::value_type>
971 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
972void lift_bc(V&& b,
const Form<T, U>& a, mdspan2_t x_dofmap,
973 md::mdspan<
const scalar_value_t<T>,
974 md::extents<std::size_t, md::dynamic_extent, 3>>
976 std::span<const T> constants,
977 const std::map<std::pair<IntegralType, int>,
978 std::pair<std::span<const T>,
int>>& coefficients,
979 std::span<const T> bc_values1,
980 std::span<const std::int8_t> bc_markers1, std::span<const T> x0,
984 std::shared_ptr<const mesh::Mesh<U>> mesh = a.mesh();
988 auto mesh0 = a.function_spaces().at(0)->mesh();
992 auto mesh1 = a.function_spaces().at(1)->mesh();
996 assert(a.function_spaces().at(0));
997 assert(a.function_spaces().at(1));
998 auto dofmap0 = a.function_spaces()[0]->dofmap()->map();
999 const int bs0 = a.function_spaces()[0]->dofmap()->bs();
1000 auto element0 = a.function_spaces()[0]->element();
1001 auto dofmap1 = a.function_spaces()[1]->dofmap()->map();
1002 const int bs1 = a.function_spaces()[1]->dofmap()->bs();
1003 auto element1 = a.function_spaces()[1]->element();
1006 std::span<const std::uint32_t> cell_info0;
1007 std::span<const std::uint32_t> cell_info1;
1009 if (element0->needs_dof_transformations()
1010 or element1->needs_dof_transformations() or a.needs_facet_permutations())
1012 mesh0->topology_mutable()->create_entity_permutations();
1013 mesh1->topology_mutable()->create_entity_permutations();
1014 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
1015 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
1018 fem::DofTransformKernel<T>
auto P0
1020 fem::DofTransformKernel<T>
auto P1T
1021 = element1->template dof_transformation_right_fn<T>(
1032 assert(_coeffs.size() ==
cells.size() * cstride);
1033 auto coeffs = md::mdspan(_coeffs.data(),
cells.size(), cstride);
1034 if (bs0 == 1 and bs1 == 1)
1036 _lift_bc_cells<1, 1>(b, x_dofmap, x, kernel, cells,
1037 {dofmap0, bs0, cells0}, P0, {dofmap1, bs1, cells1},
1038 P1T, constants, coeffs, cell_info0, cell_info1,
1039 bc_values1, bc_markers1, x0, alpha);
1041 else if (bs0 == 3 and bs1 == 3)
1043 _lift_bc_cells<3, 3>(b, x_dofmap, x, kernel, cells,
1044 {dofmap0, bs0, cells0}, P0, {dofmap1, bs1, cells1},
1045 P1T, constants, coeffs, cell_info0, cell_info1,
1046 bc_values1, bc_markers1, x0, alpha);
1050 _lift_bc_cells(b, x_dofmap, x, kernel, cells, {dofmap0, bs0, cells0}, P0,
1051 {dofmap1, bs1, cells1}, P1T, constants, coeffs, cell_info0,
1052 cell_info1, bc_values1, bc_markers1, x0, alpha);
1056 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> facet_perms;
1057 if (a.needs_facet_permutations())
1060 int num_facets_per_cell
1062 mesh->topology_mutable()->create_entity_permutations();
1063 const std::vector<std::uint8_t>& p
1064 = mesh->topology()->get_facet_permutations();
1065 facet_perms = md::mdspan(p.data(), p.size() / num_facets_per_cell,
1066 num_facets_per_cell);
1073 auto& [coeffs, cstride]
1077 = md::mdspan<
const std::int32_t,
1078 md::extents<std::size_t, md::dynamic_extent, 2, 2>>;
1080 = md::mdspan<
const T, md::extents<std::size_t, md::dynamic_extent, 2,
1081 md::dynamic_extent>>;
1083 mdspanx22_t facets(f.data(), f.size() / 4, 2, 2);
1085 mdspanx22_t facets0(f0.data(), f0.size() / 4, 2, 2);
1087 mdspanx22_t facets1(f1.data(), f1.size() / 4, 2, 2);
1088 _lift_bc_interior_facets(
1089 b, x_dofmap, x, kernel, facets, {dofmap0, bs0, facets0}, P0,
1090 {dofmap1, bs1, facets1}, P1T, constants,
1091 mdspanx2x_t(coeffs.data(), facets.extent(0), 2, cstride), cell_info0,
1092 cell_info1, bc_values1, bc_markers1, x0, alpha, facet_perms);
1098 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> perms
1101 : md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>>{};
1103 for (
int i = 0; i < a.num_integrals(itg_type, 0); ++i)
1105 auto kernel = a.kernel(itg_type, i, 0);
1107 auto& [coeffs, cstride] = coefficients.at({itg_type, i});
1110 = md::mdspan<
const std::int32_t,
1111 md::extents<std::size_t, md::dynamic_extent, 2>>;
1112 std::span e = a.domain(itg_type, i, 0);
1113 mdspanx2_t entities(e.data(), e.size() / 2, 2);
1114 std::span e0 = a.domain_arg(itg_type, 0, i, 0);
1115 mdspanx2_t entities0(e0.data(), e0.size() / 2, 2);
1116 std::span e1 = a.domain_arg(itg_type, 1, i, 0);
1117 mdspanx2_t entities1(e1.data(), e1.size() / 2, 2);
1118 assert(coeffs.size() == entities.extent(0) * cstride);
1120 b, x_dofmap, x, kernel, entities, {dofmap0, bs0, entities0}, P0,
1121 {dofmap1, bs1, entities1}, P1T, constants,
1122 md::mdspan(coeffs.data(), entities.extent(0), cstride), cell_info0,
1123 cell_info1, bc_values1, bc_markers1, x0, alpha, perms);
1148template <
typename V, std::floating_point U,
1149 dolfinx::scalar T =
typename std::remove_cvref_t<V>::value_type>
1150 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
1153 std::vector<std::optional<std::reference_wrapper<
const Form<T, U>>>> a,
1154 const std::vector<std::span<const T>>& constants,
1155 const std::vector<std::map<std::pair<IntegralType, int>,
1156 std::pair<std::span<const T>,
int>>>& coeffs,
1158 std::vector<std::reference_wrapper<
const DirichletBC<T, U>>>>& bcs1,
1159 const std::vector<std::span<const T>>& x0, T alpha)
1161 if (!x0.empty() and x0.size() != a.size())
1163 throw std::runtime_error(
1164 "Mismatch in size between x0 and bilinear form in assembler.");
1167 if (a.size() != bcs1.size())
1169 throw std::runtime_error(
1170 "Mismatch in size between a and bcs in assembler.");
1173 for (std::size_t j = 0; j < a.size(); ++j)
1175 std::vector<std::int8_t> bc_markers1;
1176 std::vector<T> bc_values1;
1177 if (a[j] and !bcs1[j].empty())
1180 std::shared_ptr<const mesh::Mesh<U>> mesh = a[j]->get().mesh();
1182 throw std::runtime_error(
"Unable to extract a mesh.");
1183 mdspan2_t x_dofmap = mesh->geometry().dofmap();
1184 std::span _x = mesh->geometry().x();
1185 md::mdspan<const scalar_value_t<T>,
1186 md::extents<std::size_t, md::dynamic_extent, 3>>
1187 x(_x.data(), _x.size() / 3, 3);
1189 assert(a[j]->get().function_spaces().at(0));
1190 auto V1 = a[j]->get().function_spaces()[1];
1192 auto map1 = V1->dofmap()->index_map;
1193 const int bs1 = V1->dofmap()->index_map_bs();
1195 const int crange = bs1 * (map1->size_local() + map1->num_ghosts());
1196 bc_markers1.assign(crange,
false);
1197 bc_values1.assign(crange, 0);
1198 for (
auto& bc : bcs1[j])
1200 bc.get().mark_dofs(bc_markers1);
1201 bc.get().set(bc_values1, std::nullopt, 1);
1206 lift_bc(b, a[j]->get(), x_dofmap, x, constants[j], coeffs[j],
1207 std::span<const T>(bc_values1), bc_markers1, x0[j], alpha);
1211 lift_bc(b, a[j]->get(), x_dofmap, x, constants[j], coeffs[j],
1212 std::span<const T>(bc_values1), bc_markers1,
1213 std::span<const T>(), alpha);
1226template <
typename V, std::floating_point U,
1227 dolfinx::scalar T =
typename std::remove_cvref_t<V>::value_type>
1228 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
1229void assemble_vector(
1230 V&& b,
const Form<T, U>& L,
1231 md::mdspan<
const scalar_value_t<T>,
1232 md::extents<std::size_t, md::dynamic_extent, 3>>
1234 std::span<const T> constants,
1235 const std::map<std::pair<IntegralType, int>,
1236 std::pair<std::span<const T>,
int>>& coefficients)
1239 std::shared_ptr<const mesh::Mesh<U>> mesh = L.mesh();
1243 auto mesh0 = L.function_spaces().at(0)->mesh();
1246 const int num_cell_types = mesh->topology()->cell_types().size();
1247 for (
int cell_type_idx = 0; cell_type_idx < num_cell_types; ++cell_type_idx)
1250 mdspan2_t x_dofmap = mesh->geometry().dofmap(cell_type_idx);
1253 assert(L.function_spaces().at(0));
1254 auto element = L.function_spaces().at(0)->elements(cell_type_idx);
1256 std::shared_ptr<const fem::DofMap> dofmap
1257 = L.function_spaces().at(0)->dofmaps(cell_type_idx);
1259 auto dofs = dofmap->map();
1260 const int bs = dofmap->bs();
1262 fem::DofTransformKernel<T>
auto P0
1265 std::span<const std::uint32_t> cell_info0;
1266 if (element->needs_dof_transformations() or L.needs_facet_permutations())
1268 mesh0->topology_mutable()->create_entity_permutations();
1269 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
1279 assert(
cells.size() * cstride == coeffs.size());
1282 impl::assemble_cells<1>(
1283 P0, b, x_dofmap, x, cells, {dofs, bs, cells0}, fn, constants,
1284 md::mdspan(coeffs.data(),
cells.size(), cstride), cell_info0);
1288 impl::assemble_cells<3>(
1289 P0, b, x_dofmap, x, cells, {dofs, bs, cells0}, fn, constants,
1290 md::mdspan(coeffs.data(),
cells.size(), cstride), cell_info0);
1294 impl::assemble_cells(
1295 P0, b, x_dofmap, x, cells, {dofs, bs, cells0}, fn, constants,
1296 md::mdspan(coeffs.data(),
cells.size(), cstride), cell_info0);
1300 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> facet_perms;
1301 if (L.needs_facet_permutations())
1303 mesh::CellType cell_type = mesh->topology()->cell_types()[cell_type_idx];
1304 int num_facets_per_cell
1306 mesh->topology_mutable()->create_entity_permutations();
1307 const std::vector<std::uint8_t>& p
1308 = mesh->topology()->get_facet_permutations();
1309 facet_perms = md::mdspan(p.data(), p.size() / num_facets_per_cell,
1310 num_facets_per_cell);
1314 = md::mdspan<
const std::int32_t,
1315 md::extents<std::size_t, md::dynamic_extent, 2>>;
1320 = md::mdspan<
const std::int32_t,
1321 md::extents<std::size_t, md::dynamic_extent, 2, 2>>;
1323 = md::mdspan<
const T, md::extents<std::size_t, md::dynamic_extent, 2,
1324 md::dynamic_extent>>;
1328 auto& [coeffs, cstride]
1332 assert((facets.size() / 4) * 2 * cstride == coeffs.size());
1335 impl::assemble_interior_facets<1>(
1337 mdspanx22_t(facets.data(), facets.size() / 4, 2, 2),
1339 mdspanx22_t(facets1.data(), facets1.size() / 4, 2, 2)},
1341 mdspanx2x_t(coeffs.data(), facets.size() / 4, 2, cstride),
1342 cell_info0, facet_perms);
1346 impl::assemble_interior_facets<3>(
1348 mdspanx22_t(facets.data(), facets.size() / 4, 2, 2),
1350 mdspanx22_t(facets1.data(), facets1.size() / 4, 2, 2)},
1352 mdspanx2x_t(coeffs.data(), facets.size() / 4, 2, cstride),
1353 cell_info0, facet_perms);
1357 impl::assemble_interior_facets(
1359 mdspanx22_t(facets.data(), facets.size() / 4, 2, 2),
1361 mdspanx22_t(facets1.data(), facets1.size() / 4, 2, 2)},
1363 mdspanx2x_t(coeffs.data(), facets.size() / 4, 2, cstride),
1364 cell_info0, facet_perms);
1371 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> perms
1374 : md::mdspan<const std::uint8_t,
1375 md::dextents<std::size_t, 2>>{};
1376 for (
int i = 0; i < L.num_integrals(itg_type, 0); ++i)
1378 auto fn = L.kernel(itg_type, i, 0);
1380 auto& [coeffs, cstride] = coefficients.at({itg_type, i});
1381 std::span e = L.domain(itg_type, i, 0);
1382 mdspanx2_t entities(e.data(), e.size() / 2, 2);
1383 std::span e1 = L.domain_arg(itg_type, 0, i, 0);
1384 mdspanx2_t entities1(e1.data(), e1.size() / 2, 2);
1385 assert((entities.size() / 2) * cstride == coeffs.size());
1388 impl::assemble_entities<1>(
1389 P0, b, x_dofmap, x, entities, {dofs, bs, entities1}, fn,
1390 constants, md::mdspan(coeffs.data(), entities.extent(0), cstride),
1395 impl::assemble_entities<3>(
1396 P0, b, x_dofmap, x, entities, {dofs, bs, entities1}, fn,
1398 md::mdspan(coeffs.data(), entities.size() / 2, cstride),
1403 impl::assemble_entities(
1404 P0, b, x_dofmap, x, entities, {dofs, bs, entities1}, fn,
1406 md::mdspan(coeffs.data(), entities.size() / 2, cstride),
1420template <
typename V, std::floating_point U,
1421 dolfinx::scalar T =
typename std::remove_cvref_t<V>::value_type>
1422 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
1423void assemble_vector(
1424 V&& b,
const Form<T, U>& L, std::span<const T> constants,
1425 const std::map<std::pair<IntegralType, int>,
1426 std::pair<std::span<const T>,
int>>& coefficients)
1429 = md::mdspan<const scalar_value_t<T>,
1430 md::extents<std::size_t, md::dynamic_extent, 3>>;
1432 std::shared_ptr<const mesh::Mesh<U>> mesh = L.mesh();
1434 auto x = mesh->geometry().x();
1435 if constexpr (std::is_same_v<U, scalar_value_t<T>>)
1437 impl::assemble_vector(b, L, mdspanx3_t(x.data(), x.size() / 3, 3),
1438 constants, coefficients);
1442 std::vector<scalar_value_t<T>> _x(x.begin(), x.end());
1443 impl::assemble_vector(b, L, mdspanx3_t(_x.data(), _x.size() / 3, 3),
1444 constants, coefficients);
Degree-of-freedom map representations and tools.
Definition DirichletBC.h:255
Functions supporting finite element method operations.
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:16
Finite element method functionality.
Definition assemble_expression_impl.h:23
@ transpose
Transpose.
Definition FiniteElement.h:28
@ standard
Standard.
Definition FiniteElement.h:27
@ vertex
Vertex.
Definition Form.h:43
@ interior_facet
Interior facet.
Definition Form.h:42
@ ridge
Ridge.
Definition Form.h:44
@ cell
Cell.
Definition Form.h:40
@ exterior_facet
Facet.
Definition Form.h:41
CellType
Cell type identifier.
Definition cell_types.h:22
int cell_num_entities(CellType type, int dim)
Definition cell_types.cpp:139