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 <dolfinx::scalar T,
int _bs0 = -1,
int _bs1 = -1>
81 std::span<T> b, mdspan2_t x_dofmap,
82 md::mdspan<
const scalar_value_t<T>,
83 md::extents<std::size_t, md::dynamic_extent, 3>>
85 FEkernel<T>
auto kernel, std::span<const std::int32_t> cells,
86 std::tuple<mdspan2_t,
int, std::span<const std::int32_t>> dofmap0,
87 fem::DofTransformKernel<T>
auto P0,
88 std::tuple<mdspan2_t,
int, std::span<const std::int32_t>> dofmap1,
89 fem::DofTransformKernel<T>
auto P1T, std::span<const T> constants,
90 md::mdspan<
const T, md::dextents<std::size_t, 2>> coeffs,
91 std::span<const std::uint32_t> cell_info0,
92 std::span<const std::uint32_t> cell_info1, std::span<const T> bc_values1,
93 std::span<const std::int8_t> bc_markers1, std::span<const T> x0, T alpha)
98 const auto [dmap0, bs0, cells0] = dofmap0;
99 const auto [dmap1, bs1, cells1] = dofmap1;
100 assert(_bs0 < 0 or _bs0 == bs0);
101 assert(_bs1 < 0 or _bs1 == bs1);
104 std::vector<scalar_value_t<T>> cdofs(3 * x_dofmap.extent(1));
105 std::vector<T> Ae, be;
106 assert(cells0.size() ==
cells.size());
107 assert(cells1.size() ==
cells.size());
108 for (std::size_t index = 0; index <
cells.size(); ++index)
112 std::int32_t c =
cells[index];
113 std::int32_t c0 = cells0[index];
114 std::int32_t c1 = cells1[index];
117 auto dofs1 = md::submdspan(dmap1, c1, md::full_extent);
121 for (std::size_t j = 0; j < dofs1.size(); ++j)
123 if constexpr (_bs1 > 0)
125 for (
int k = 0; k < _bs1; ++k)
127 assert(_bs1 * dofs1[j] + k < (
int)bc_markers1.size());
128 if (bc_markers1[_bs1 * dofs1[j] + k])
137 for (
int k = 0; k < bs1; ++k)
139 assert(bs1 * dofs1[j] + k < (
int)bc_markers1.size());
140 if (bc_markers1[bs1 * dofs1[j] + k])
153 auto x_dofs = md::submdspan(x_dofmap, c, md::full_extent);
154 for (std::size_t i = 0; i < x_dofs.size(); ++i)
155 std::copy_n(&x(x_dofs[i], 0), 3, std::next(cdofs.begin(), 3 * i));
158 auto dofs0 = md::submdspan(dmap0, c0, md::full_extent);
160 const int num_rows = bs0 * dofs0.size();
161 const int num_cols = bs1 * dofs1.size();
163 Ae.resize(num_rows * num_cols);
164 std::ranges::fill(Ae, 0);
165 kernel(Ae.data(), &coeffs(index, 0), constants.data(), cdofs.data(),
166 nullptr,
nullptr,
nullptr);
167 P0(Ae, cell_info0, c0, num_cols);
168 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];
263template <dolfinx::scalar T>
264void _lift_bc_exterior_facets(
265 std::span<T> b, mdspan2_t x_dofmap,
266 md::mdspan<
const scalar_value_t<T>,
267 md::extents<std::size_t, md::dynamic_extent, 3>>
269 FEkernel<T>
auto kernel,
270 md::mdspan<
const std::int32_t,
271 md::extents<std::size_t, md::dynamic_extent, 2>>
273 std::tuple<mdspan2_t,
int,
274 md::mdspan<
const std::int32_t,
275 md::extents<std::size_t, md::dynamic_extent, 2>>>
277 fem::DofTransformKernel<T>
auto P0,
278 std::tuple<mdspan2_t,
int,
279 md::mdspan<
const std::int32_t,
280 md::extents<std::size_t, md::dynamic_extent, 2>>>
282 fem::DofTransformKernel<T>
auto P1T, std::span<const T> constants,
283 md::mdspan<
const T, md::dextents<std::size_t, 2>> coeffs,
284 std::span<const std::uint32_t> cell_info0,
285 std::span<const std::uint32_t> cell_info1, std::span<const T> bc_values1,
286 std::span<const std::int8_t> bc_markers1, std::span<const T> x0, T alpha,
287 md::mdspan<
const std::uint8_t, md::dextents<std::size_t, 2>> perms)
292 const auto [dmap0, bs0, facets0] = dofmap0;
293 const auto [dmap1, bs1, facets1] = dofmap1;
296 std::vector<scalar_value_t<T>> cdofs(3 * x_dofmap.extent(1));
297 std::vector<T> Ae, be;
298 assert(facets0.size() == facets.size());
299 assert(facets1.size() == facets.size());
300 for (std::size_t index = 0; index < facets.extent(0); ++index)
304 std::int32_t
cell = facets(index, 0);
305 std::int32_t cell0 = facets0(index, 0);
306 std::int32_t cell1 = facets1(index, 0);
309 std::int32_t local_facet = facets(index, 1);
312 auto dofs1 = md::submdspan(dmap1, cell1, md::full_extent);
316 for (std::size_t j = 0; j < dofs1.size(); ++j)
318 for (
int k = 0; k < bs1; ++k)
320 if (bc_markers1[bs1 * dofs1[j] + k])
332 auto x_dofs = md::submdspan(x_dofmap,
cell, md::full_extent);
333 for (std::size_t i = 0; i < x_dofs.size(); ++i)
334 std::copy_n(&x(x_dofs[i], 0), 3, std::next(cdofs.begin(), 3 * i));
337 auto dofs0 = md::submdspan(dmap0, cell0, md::full_extent);
339 const int num_rows = bs0 * dofs0.size();
340 const int num_cols = bs1 * dofs1.size();
343 std::uint8_t perm = perms.empty() ? 0 : perms(
cell, local_facet);
345 Ae.resize(num_rows * num_cols);
346 std::ranges::fill(Ae, 0);
347 kernel(Ae.data(), &coeffs(index, 0), constants.data(), cdofs.data(),
348 &local_facet, &perm,
nullptr);
349 P0(Ae, cell_info0, cell0, num_cols);
350 P1T(Ae, cell_info1, cell1, num_rows);
354 std::ranges::fill(be, 0);
355 for (std::size_t j = 0; j < dofs1.size(); ++j)
357 for (
int k = 0; k < bs1; ++k)
359 const std::int32_t jj = bs1 * dofs1[j] + k;
362 const T bc = bc_values1[jj];
363 const T _x0 = x0.empty() ? 0 : x0[jj];
365 for (
int m = 0; m < num_rows; ++m)
366 be[m] -= Ae[m * num_cols + bs1 * j + k] * alpha * (bc - _x0);
371 for (std::size_t i = 0; i < dofs0.size(); ++i)
372 for (
int k = 0; k < bs0; ++k)
373 b[bs0 * dofs0[i] + k] += be[bs0 * i + k];
419template <dolfinx::scalar T>
420void _lift_bc_interior_facets(
421 std::span<T> b, mdspan2_t x_dofmap,
422 md::mdspan<
const scalar_value_t<T>,
423 md::extents<std::size_t, md::dynamic_extent, 3>>
425 FEkernel<T>
auto kernel,
426 md::mdspan<
const std::int32_t,
427 md::extents<std::size_t, md::dynamic_extent, 2, 2>>
429 std::tuple<mdspan2_t,
int,
430 md::mdspan<
const std::int32_t,
431 md::extents<std::size_t, md::dynamic_extent, 2, 2>>>
433 fem::DofTransformKernel<T>
auto P0,
434 std::tuple<mdspan2_t,
int,
435 md::mdspan<
const std::int32_t,
436 md::extents<std::size_t, md::dynamic_extent, 2, 2>>>
438 fem::DofTransformKernel<T>
auto P1T, std::span<const T> constants,
439 md::mdspan<
const T, md::extents<std::size_t, md::dynamic_extent, 2,
442 std::span<const std::uint32_t> cell_info0,
443 std::span<const std::uint32_t> cell_info1, std::span<const T> bc_values1,
444 std::span<const std::int8_t> bc_markers1, std::span<const T> x0, T alpha,
445 md::mdspan<
const std::uint8_t, md::dextents<std::size_t, 2>> perms)
450 const auto [dmap0, bs0, facets0] = dofmap0;
451 const auto [dmap1, bs1, facets1] = dofmap1;
454 using X = scalar_value_t<T>;
455 std::vector<X> cdofs(2 * x_dofmap.extent(1) * 3);
456 std::span<X> cdofs0(cdofs.data(), x_dofmap.extent(1) * 3);
457 std::span<X> cdofs1(cdofs.data() + x_dofmap.extent(1) * 3,
458 x_dofmap.extent(1) * 3);
459 std::vector<T> Ae, be;
462 std::vector<std::int32_t> dmapjoint0, dmapjoint1;
464 const int num_dofs0 = dmap0.extent(1);
465 const int num_dofs1 = dmap1.extent(1);
466 assert(facets0.size() == facets.size());
467 assert(facets1.size() == facets.size());
468 for (std::size_t f = 0; f < facets.extent(0); ++f)
472 std::array
cells{facets(f, 0, 0), facets(f, 1, 0)};
473 std::array cells0{facets0(f, 0, 0), facets0(f, 1, 0)};
474 std::array cells1{facets1(f, 0, 0), facets1(f, 1, 0)};
477 std::array local_facet = {facets(f, 0, 1), facets(f, 1, 1)};
480 auto x_dofs0 = md::submdspan(x_dofmap, cells[0], md::full_extent);
481 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
482 std::copy_n(&x(x_dofs0[i], 0), 3, std::next(cdofs0.begin(), 3 * i));
483 auto x_dofs1 = md::submdspan(x_dofmap, cells[1], md::full_extent);
484 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
485 std::copy_n(&x(x_dofs1[i], 0), 3, std::next(cdofs1.begin(), 3 * i));
491 std::span<const std::int32_t> dmap0_cell0
493 ? std::span(dmap0.data_handle() + cells0[0] * num_dofs0,
495 : std::span<const std::int32_t>();
496 std::span<const std::int32_t> dmap0_cell1
498 ? std::span(dmap0.data_handle() + cells0[1] * num_dofs0,
500 : std::span<const std::int32_t>();
502 dmapjoint0.resize(2 * num_dofs0);
503 std::ranges::copy(dmap0_cell0, dmapjoint0.begin());
504 std::ranges::copy(dmap0_cell1, std::next(dmapjoint0.begin(), num_dofs0));
507 std::span<const std::int32_t> dmap1_cell0
509 ? std::span(dmap1.data_handle() + cells1[0] * num_dofs1,
511 : std::span<const std::int32_t>();
512 std::span<const std::int32_t> dmap1_cell1
514 ? std::span(dmap1.data_handle() + cells1[1] * num_dofs1,
516 : std::span<const std::int32_t>();
518 dmapjoint1.resize(2 * num_dofs1);
519 std::ranges::copy(dmap1_cell0, dmapjoint1.begin());
520 std::ranges::copy(dmap1_cell1, std::next(dmapjoint1.begin(), num_dofs1));
524 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
526 for (
int k = 0; k < bs1; ++k)
528 if (bc_markers1[bs1 * dmap1_cell0[j] + k])
537 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
539 for (
int k = 0; k < bs1; ++k)
541 if (bc_markers1[bs1 * dmap1_cell1[j] + k])
552 const int num_rows = bs0 * 2 * num_dofs0;
553 const int num_cols = bs1 * 2 * num_dofs1;
556 Ae.resize(num_rows * num_cols);
557 std::ranges::fill(Ae, 0);
558 std::array perm = perms.empty()
559 ? std::array<std::uint8_t, 2>{0, 0}
560 : std::array{perms(cells[0], local_facet[0]),
561 perms(cells[1], local_facet[1])};
562 kernel(Ae.data(), &coeffs(f, 0, 0), constants.data(), cdofs.data(),
563 local_facet.data(), perm.data(),
nullptr);
565 std::span<T> _Ae(Ae);
567 = _Ae.subspan(bs0 * num_dofs0 * num_cols, bs0 * num_dofs1 * num_cols);
570 P0(_Ae, cell_info0, cells0[0], num_cols);
572 P0(sub_Ae0, cell_info0, cells0[1], num_cols);
574 P1T(_Ae, cell_info1, cells1[0], num_rows);
578 for (
int row = 0; row < num_rows; ++row)
583 = _Ae.subspan(row * num_cols + bs1 * num_dofs1, bs1 * num_dofs1);
584 P1T(sub_Ae1, cell_info1, cells1[1], 1);
589 std::ranges::fill(be, 0);
592 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
594 for (
int k = 0; k < bs1; ++k)
596 const std::int32_t jj = bs1 * dmap1_cell0[j] + k;
599 const T bc = bc_values1[jj];
600 const T _x0 = x0.empty() ? 0 : x0[jj];
602 for (
int m = 0; m < num_rows; ++m)
603 be[m] -= Ae[m * num_cols + bs1 * j + k] * alpha * (bc - _x0);
609 const int offset = bs1 * num_dofs1;
610 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
612 for (
int k = 0; k < bs1; ++k)
614 const std::int32_t jj = bs1 * dmap1_cell1[j] + k;
617 const T bc = bc_values1[jj];
618 const T _x0 = x0.empty() ? 0 : x0[jj];
620 for (
int m = 0; m < num_rows; ++m)
623 -= Ae[m * num_cols + offset + bs1 * j + k] * alpha * (bc - _x0);
629 for (std::size_t i = 0; i < dmap0_cell0.size(); ++i)
630 for (
int k = 0; k < bs0; ++k)
631 b[bs0 * dmap0_cell0[i] + k] += be[bs0 * i + k];
633 const int offset_be = bs0 * num_dofs0;
634 for (std::size_t i = 0; i < dmap0_cell1.size(); ++i)
635 for (
int k = 0; k < bs0; ++k)
636 b[bs0 * dmap0_cell1[i] + k] += be[offset_be + bs0 * i + k];
663template <dolfinx::scalar T,
int _bs = -1>
665 fem::DofTransformKernel<T>
auto P0, std::span<T> b, mdspan2_t x_dofmap,
666 md::mdspan<
const scalar_value_t<T>,
667 md::extents<std::size_t, md::dynamic_extent, 3>>
669 std::span<const std::int32_t> cells,
670 std::tuple<mdspan2_t,
int, std::span<const std::int32_t>> dofmap,
671 FEkernel<T>
auto kernel, std::span<const T> constants,
672 md::mdspan<
const T, md::dextents<std::size_t, 2>> coeffs,
673 std::span<const std::uint32_t> cell_info0)
678 const auto [dmap, bs, cells0] = dofmap;
679 assert(_bs < 0 or _bs == bs);
682 std::vector<scalar_value_t<T>> cdofs(3 * x_dofmap.extent(1));
683 std::vector<T> be(bs * dmap.extent(1));
684 std::span<T> _be(be);
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];
745template <dolfinx::scalar T,
int _bs = -1>
746void assemble_exterior_facets(
747 fem::DofTransformKernel<T>
auto P0, std::span<T> b, mdspan2_t x_dofmap,
748 md::mdspan<
const scalar_value_t<T>,
749 md::extents<std::size_t, md::dynamic_extent, 3>>
751 md::mdspan<
const std::int32_t,
752 std::extents<std::size_t, md::dynamic_extent, 2>>
754 std::tuple<mdspan2_t,
int,
755 md::mdspan<
const std::int32_t,
756 std::extents<std::size_t, md::dynamic_extent, 2>>>
758 FEkernel<T>
auto fn, std::span<const T> constants,
759 md::mdspan<
const T, md::dextents<std::size_t, 2>> coeffs,
760 std::span<const std::uint32_t> cell_info0,
761 md::mdspan<
const std::uint8_t, md::dextents<std::size_t, 2>> perms)
766 const auto [dmap, bs, facets0] = dofmap;
767 assert(_bs < 0 or _bs == bs);
770 const int num_dofs = dmap.extent(1);
771 std::vector<scalar_value_t<T>> cdofs(3 * x_dofmap.extent(1));
772 std::vector<T> be(bs * num_dofs);
773 std::span<T> _be(be);
774 assert(facets0.size() == facets.size());
775 for (std::size_t f = 0; f < facets.extent(0); ++f)
779 std::int32_t
cell = facets(f, 0);
780 std::int32_t local_facet = facets(f, 1);
781 std::int32_t cell0 = facets0(f, 0);
784 auto x_dofs = md::submdspan(x_dofmap,
cell, md::full_extent);
785 for (std::size_t i = 0; i < x_dofs.size(); ++i)
786 std::copy_n(&x(x_dofs[i], 0), 3, std::next(cdofs.begin(), 3 * i));
789 std::uint8_t perm = perms.empty() ? 0 : perms(
cell, local_facet);
792 std::ranges::fill(be, 0);
793 fn(be.data(), &coeffs(f, 0), constants.data(), cdofs.data(), &local_facet,
796 P0(_be, cell_info0, cell0, 1);
799 auto dofs = md::submdspan(dmap, cell0, md::full_extent);
800 if constexpr (_bs > 0)
802 for (std::size_t i = 0; i < dofs.size(); ++i)
803 for (
int k = 0; k < _bs; ++k)
804 b[_bs * dofs[i] + k] += be[_bs * i + k];
808 for (std::size_t i = 0; i < dofs.size(); ++i)
809 for (
int k = 0; k < bs; ++k)
810 b[bs * dofs[i] + k] += be[bs * i + k];
840template <dolfinx::scalar T,
int _bs = -1>
841void assemble_interior_facets(
842 fem::DofTransformKernel<T>
auto P0, std::span<T> b, mdspan2_t x_dofmap,
843 md::mdspan<
const scalar_value_t<T>,
844 md::extents<std::size_t, md::dynamic_extent, 3>>
846 md::mdspan<
const std::int32_t,
847 std::extents<std::size_t, md::dynamic_extent, 2, 2>>
849 std::tuple<
const DofMap&,
int,
850 md::mdspan<
const std::int32_t,
851 std::extents<std::size_t, md::dynamic_extent, 2, 2>>>
853 FEkernel<T>
auto fn, std::span<const T> constants,
854 md::mdspan<
const T, md::extents<std::size_t, md::dynamic_extent, 2,
857 std::span<const std::uint32_t> cell_info0,
858 md::mdspan<
const std::uint8_t, md::dextents<std::size_t, 2>> perms)
863 const auto [dmap, bs, facets0] = dofmap;
864 assert(_bs < 0 or _bs == bs);
867 using X = scalar_value_t<T>;
868 std::vector<X> cdofs(2 * x_dofmap.extent(1) * 3);
869 std::span<X> cdofs0(cdofs.data(), x_dofmap.extent(1) * 3);
870 std::span<X> cdofs1(cdofs.data() + x_dofmap.extent(1) * 3,
871 x_dofmap.extent(1) * 3);
874 const std::size_t dmap_size = dmap.cell_dofs(0).size();
876 assert(facets0.size() == facets.size());
877 for (std::size_t f = 0; f < facets.extent(0); ++f)
880 std::array<std::int32_t, 2>
cells{facets(f, 0, 0), facets(f, 1, 0)};
881 std::array<std::int32_t, 2> cells0{facets0(f, 0, 0), facets0(f, 1, 0)};
884 std::array<std::int32_t, 2> local_facet{facets(f, 0, 1), facets(f, 1, 1)};
887 auto x_dofs0 = md::submdspan(x_dofmap, cells[0], md::full_extent);
888 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
889 std::copy_n(&x(x_dofs0[i], 0), 3, std::next(cdofs0.begin(), 3 * i));
890 auto x_dofs1 = md::submdspan(x_dofmap, cells[1], md::full_extent);
891 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
892 std::copy_n(&x(x_dofs1[i], 0), 3, std::next(cdofs1.begin(), 3 * i));
898 std::span<const std::int32_t> dmap0 = cells0[0] >= 0
899 ? dmap.cell_dofs(cells0[0])
900 : std::span<const std::int32_t>();
901 std::span<const std::int32_t> dmap1 = cells0[1] >= 0
902 ? dmap.cell_dofs(cells0[1])
903 : std::span<const std::int32_t>();
906 be.resize(bs * 2 * dmap_size);
907 std::ranges::fill(be, 0);
908 std::array perm = perms.empty()
909 ? std::array<std::uint8_t, 2>{0, 0}
910 : std::array{perms(cells[0], local_facet[0]),
911 perms(cells[1], local_facet[1])};
912 fn(be.data(), &coeffs(f, 0, 0), constants.data(), cdofs.data(),
913 local_facet.data(), perm.data(),
nullptr);
915 std::span<T> _be(be);
916 std::span<T> sub_be = _be.subspan(bs * dmap_size, bs * dmap_size);
919 P0(be, cell_info0, cells0[0], 1);
921 P0(sub_be, cell_info0, cells0[1], 1);
924 if constexpr (_bs > 0)
926 for (std::size_t i = 0; i < dmap0.size(); ++i)
927 for (
int k = 0; k < _bs; ++k)
928 b[_bs * dmap0[i] + k] += be[_bs * i + k];
929 for (std::size_t i = 0; i < dmap1.size(); ++i)
930 for (
int k = 0; k < _bs; ++k)
931 b[_bs * dmap1[i] + k] += be[_bs * (i + dmap_size) + k];
935 for (std::size_t i = 0; i < dmap0.size(); ++i)
936 for (
int k = 0; k < bs; ++k)
937 b[bs * dmap0[i] + k] += be[bs * i + k];
938 for (std::size_t i = 0; i < dmap1.size(); ++i)
939 for (
int k = 0; k < bs; ++k)
940 b[bs * dmap1[i] + k] += be[bs * (i + dmap_size) + k];
961template <dolfinx::scalar T, std::
floating_po
int U>
962void lift_bc(std::span<T> b,
const Form<T, U>& a, mdspan2_t x_dofmap,
963 md::mdspan<
const scalar_value_t<T>,
964 md::extents<std::size_t, md::dynamic_extent, 3>>
966 std::span<const T> constants,
967 const std::map<std::pair<IntegralType, int>,
968 std::pair<std::span<const T>,
int>>& coefficients,
969 std::span<const T> bc_values1,
970 std::span<const std::int8_t> bc_markers1, std::span<const T> x0,
974 std::shared_ptr<const mesh::Mesh<U>> mesh = a.mesh();
978 auto mesh0 = a.function_spaces().at(0)->mesh();
982 auto mesh1 = a.function_spaces().at(1)->mesh();
986 assert(a.function_spaces().at(0));
987 assert(a.function_spaces().at(1));
988 auto dofmap0 = a.function_spaces()[0]->dofmap()->map();
989 const int bs0 = a.function_spaces()[0]->dofmap()->bs();
990 auto element0 = a.function_spaces()[0]->element();
991 auto dofmap1 = a.function_spaces()[1]->dofmap()->map();
992 const int bs1 = a.function_spaces()[1]->dofmap()->bs();
993 auto element1 = a.function_spaces()[1]->element();
996 std::span<const std::uint32_t> cell_info0;
997 std::span<const std::uint32_t> cell_info1;
999 if (element0->needs_dof_transformations()
1000 or element1->needs_dof_transformations() or a.needs_facet_permutations())
1002 mesh0->topology_mutable()->create_entity_permutations();
1003 mesh1->topology_mutable()->create_entity_permutations();
1004 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
1005 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
1008 fem::DofTransformKernel<T>
auto P0
1010 fem::DofTransformKernel<T>
auto P1T
1011 = element1->template dof_transformation_right_fn<T>(
1022 assert(_coeffs.size() ==
cells.size() * cstride);
1023 auto coeffs = md::mdspan(_coeffs.data(),
cells.size(), cstride);
1024 if (bs0 == 1 and bs1 == 1)
1026 _lift_bc_cells<T, 1, 1>(
1027 b, x_dofmap, x, kernel, cells, {dofmap0, bs0, cells0}, P0,
1028 {dofmap1, bs1, cells1}, P1T, constants, coeffs, cell_info0,
1029 cell_info1, bc_values1, bc_markers1, x0, alpha);
1031 else if (bs0 == 3 and bs1 == 3)
1033 _lift_bc_cells<T, 3, 3>(
1034 b, x_dofmap, x, kernel, cells, {dofmap0, bs0, cells0}, P0,
1035 {dofmap1, bs1, cells1}, P1T, constants, coeffs, cell_info0,
1036 cell_info1, bc_values1, bc_markers1, x0, alpha);
1040 _lift_bc_cells(b, x_dofmap, x, kernel, cells, {dofmap0, bs0, cells0}, P0,
1041 {dofmap1, bs1, cells1}, P1T, constants, coeffs, cell_info0,
1042 cell_info1, bc_values1, bc_markers1, x0, alpha);
1046 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> perms;
1047 if (a.needs_facet_permutations())
1050 int num_facets_per_cell
1052 mesh->topology_mutable()->create_entity_permutations();
1053 const std::vector<std::uint8_t>& p
1054 = mesh->topology()->get_facet_permutations();
1055 perms = md::mdspan(p.data(), p.size() / num_facets_per_cell,
1056 num_facets_per_cell);
1063 auto& [coeffs, cstride]
1067 = md::mdspan<
const std::int32_t,
1068 md::extents<std::size_t, md::dynamic_extent, 2>>;
1070 mdspanx2_t facets(f.data(), f.size() / 2, 2);
1072 mdspanx2_t facets0(f0.data(), f0.size() / 2, 2);
1074 mdspanx2_t facets1(f1.data(), f1.size() / 2, 2);
1075 assert(coeffs.size() == facets.extent(0) * cstride);
1076 _lift_bc_exterior_facets(
1077 b, x_dofmap, x, kernel, facets, {dofmap0, bs0, facets0}, P0,
1078 {dofmap1, bs1, facets1}, P1T, constants,
1079 md::mdspan(coeffs.data(), facets.extent(0), cstride), cell_info0,
1080 cell_info1, bc_values1, bc_markers1, x0, alpha, perms);
1087 auto& [coeffs, cstride]
1091 = md::mdspan<
const std::int32_t,
1092 md::extents<std::size_t, md::dynamic_extent, 2, 2>>;
1094 = md::mdspan<
const T, md::extents<std::size_t, md::dynamic_extent, 2,
1095 md::dynamic_extent>>;
1097 mdspanx22_t facets(f.data(), f.size() / 4, 2, 2);
1099 mdspanx22_t facets0(f0.data(), f0.size() / 4, 2, 2);
1101 mdspanx22_t facets1(f1.data(), f1.size() / 4, 2, 2);
1102 _lift_bc_interior_facets(
1103 b, x_dofmap, x, kernel, facets, {dofmap0, bs0, facets0}, P0,
1104 {dofmap1, bs1, facets1}, P1T, constants,
1105 mdspanx2x_t(coeffs.data(), facets.extent(0), 2, cstride), cell_info0,
1106 cell_info1, bc_values1, bc_markers1, x0, alpha, perms);
1130template <dolfinx::scalar T, std::
floating_po
int U>
1133 std::vector<std::optional<std::reference_wrapper<
const Form<T, U>>>> a,
1134 const std::vector<std::span<const T>>& constants,
1135 const std::vector<std::map<std::pair<IntegralType, int>,
1136 std::pair<std::span<const T>,
int>>>& coeffs,
1138 std::vector<std::reference_wrapper<
const DirichletBC<T, U>>>>& bcs1,
1139 const std::vector<std::span<const T>>& x0, T alpha)
1141 if (!x0.empty() and x0.size() != a.size())
1143 throw std::runtime_error(
1144 "Mismatch in size between x0 and bilinear form in assembler.");
1147 if (a.size() != bcs1.size())
1149 throw std::runtime_error(
1150 "Mismatch in size between a and bcs in assembler.");
1153 for (std::size_t j = 0; j < a.size(); ++j)
1155 std::vector<std::int8_t> bc_markers1;
1156 std::vector<T> bc_values1;
1157 if (a[j] and !bcs1[j].empty())
1160 std::shared_ptr<const mesh::Mesh<U>> mesh = a[j]->get().mesh();
1162 throw std::runtime_error(
"Unable to extract a mesh.");
1163 mdspan2_t x_dofmap = mesh->geometry().dofmap();
1164 std::span _x = mesh->geometry().x();
1165 md::mdspan<const scalar_value_t<T>,
1166 md::extents<std::size_t, md::dynamic_extent, 3>>
1167 x(_x.data(), _x.size() / 3, 3);
1169 assert(a[j]->get().function_spaces().at(0));
1170 auto V1 = a[j]->get().function_spaces()[1];
1172 auto map1 = V1->dofmap()->index_map;
1173 const int bs1 = V1->dofmap()->index_map_bs();
1175 const int crange = bs1 * (map1->size_local() + map1->num_ghosts());
1176 bc_markers1.assign(crange,
false);
1177 bc_values1.assign(crange, 0);
1178 for (
auto& bc : bcs1[j])
1180 bc.get().mark_dofs(bc_markers1);
1181 bc.get().set(bc_values1, std::nullopt, 1);
1186 lift_bc<T>(b, a[j]->get(), x_dofmap, x, constants[j], coeffs[j],
1187 bc_values1, bc_markers1, x0[j], alpha);
1191 lift_bc<T>(b, a[j]->get(), x_dofmap, x, constants[j], coeffs[j],
1192 bc_values1, bc_markers1, std::span<const T>(), alpha);
1205template <dolfinx::scalar T, std::
floating_po
int U>
1206void assemble_vector(
1207 std::span<T> b,
const Form<T, U>& L,
1208 md::mdspan<
const scalar_value_t<T>,
1209 md::extents<std::size_t, md::dynamic_extent, 3>>
1211 std::span<const T> constants,
1212 const std::map<std::pair<IntegralType, int>,
1213 std::pair<std::span<const T>,
int>>& coefficients)
1216 std::shared_ptr<const mesh::Mesh<U>> mesh = L.mesh();
1220 auto mesh0 = L.function_spaces().at(0)->mesh();
1223 const int num_cell_types = mesh->topology()->cell_types().size();
1224 for (
int cell_type_idx = 0; cell_type_idx < num_cell_types; ++cell_type_idx)
1227 mdspan2_t x_dofmap = mesh->geometry().dofmap(cell_type_idx);
1230 assert(L.function_spaces().at(0));
1231 auto element = L.function_spaces().at(0)->elements(cell_type_idx);
1233 std::shared_ptr<const fem::DofMap> dofmap
1234 = L.function_spaces().at(0)->dofmaps(cell_type_idx);
1236 auto dofs = dofmap->map();
1237 const int bs = dofmap->bs();
1239 fem::DofTransformKernel<T>
auto P0
1242 std::span<const std::uint32_t> cell_info0;
1243 if (element->needs_dof_transformations() or L.needs_facet_permutations())
1245 mesh0->topology_mutable()->create_entity_permutations();
1246 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
1256 assert(
cells.size() * cstride == coeffs.size());
1259 impl::assemble_cells<T, 1>(
1260 P0, b, x_dofmap, x, cells, {dofs, bs, cells0}, fn, constants,
1261 md::mdspan(coeffs.data(),
cells.size(), cstride), cell_info0);
1265 impl::assemble_cells<T, 3>(
1266 P0, b, x_dofmap, x, cells, {dofs, bs, cells0}, fn, constants,
1267 md::mdspan(coeffs.data(),
cells.size(), cstride), cell_info0);
1271 impl::assemble_cells(
1272 P0, b, x_dofmap, x, cells, {dofs, bs, cells0}, fn, constants,
1273 md::mdspan(coeffs.data(),
cells.size(), cstride), cell_info0);
1277 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> perms;
1278 if (L.needs_facet_permutations())
1280 mesh::CellType cell_type = mesh->topology()->cell_types()[cell_type_idx];
1281 int num_facets_per_cell
1283 mesh->topology_mutable()->create_entity_permutations();
1284 const std::vector<std::uint8_t>& p
1285 = mesh->topology()->get_facet_permutations();
1286 perms = md::mdspan(p.data(), p.size() / num_facets_per_cell,
1287 num_facets_per_cell);
1291 = md::mdspan<
const std::int32_t,
1292 md::extents<std::size_t, md::dynamic_extent, 2>>;
1298 auto& [coeffs, cstride]
1301 mdspanx2_t facets(f.data(), f.size() / 2, 2);
1303 mdspanx2_t facets1(f1.data(), f1.size() / 2, 2);
1304 assert((facets.size() / 2) * cstride == coeffs.size());
1307 impl::assemble_exterior_facets<T, 1>(
1308 P0, b, x_dofmap, x, facets, {dofs, bs, facets1}, fn, constants,
1309 md::mdspan(coeffs.data(), facets.extent(0), cstride), cell_info0,
1314 impl::assemble_exterior_facets<T, 3>(
1315 P0, b, x_dofmap, x, facets, {dofs, bs, facets1}, fn, constants,
1316 md::mdspan(coeffs.data(), facets.size() / 2, cstride), cell_info0,
1321 impl::assemble_exterior_facets(
1322 P0, b, x_dofmap, x, facets, {dofs, bs, facets1}, fn, constants,
1323 md::mdspan(coeffs.data(), facets.size() / 2, cstride), cell_info0,
1331 = md::mdspan<
const std::int32_t,
1332 md::extents<std::size_t, md::dynamic_extent, 2, 2>>;
1334 = md::mdspan<
const T, md::extents<std::size_t, md::dynamic_extent, 2,
1335 md::dynamic_extent>>;
1339 auto& [coeffs, cstride]
1343 assert((facets.size() / 4) * 2 * cstride == coeffs.size());
1346 impl::assemble_interior_facets<T, 1>(
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),
1357 impl::assemble_interior_facets<T, 3>(
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),
1368 impl::assemble_interior_facets(
1370 mdspanx22_t(facets.data(), facets.size() / 4, 2, 2),
1372 mdspanx22_t(facets1.data(), facets1.size() / 4, 2, 2)},
1374 mdspanx2x_t(coeffs.data(), facets.size() / 4, 2, cstride),
1387template <dolfinx::scalar T, std::
floating_po
int U>
1388void assemble_vector(
1389 std::span<T> b,
const Form<T, U>& L, std::span<const T> constants,
1390 const std::map<std::pair<IntegralType, int>,
1391 std::pair<std::span<const T>,
int>>& coefficients)
1394 = md::mdspan<const scalar_value_t<T>,
1395 md::extents<std::size_t, md::dynamic_extent, 3>>;
1397 std::shared_ptr<const mesh::Mesh<U>> mesh = L.mesh();
1399 auto x = mesh->geometry().x();
1400 if constexpr (std::is_same_v<U, scalar_value_t<T>>)
1402 assemble_vector(b, L, mdspanx3_t(x.data(), x.size() / 3, 3), constants,
1407 std::vector<scalar_value_t<T>> _x(x.begin(), x.end());
1408 assemble_vector(b, L, mdspanx3_t(_x.data(), _x.size() / 3, 3), constants,
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
@ interior_facet
Interior facet.
Definition Form.h:42
@ cell
Cell.
Definition Form.h:40
@ exterior_facet
Exterior 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