10#include "DirichletBC.h"
13#include "FunctionSpace.h"
17#include <dolfinx/common/IndexMap.h>
18#include <dolfinx/mesh/Geometry.h>
19#include <dolfinx/mesh/Mesh.h>
20#include <dolfinx/mesh/Topology.h>
26namespace dolfinx::fem::impl
28using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
30 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
43 std::span<T> b, mdspan2_t x_dofmap,
44 std::span<
const scalar_value_type_t<T>> x, FEkernel<T>
auto kernel,
45 std::span<const std::int32_t> cells,
46 const std::function<
void(
const std::span<T>&,
47 const std::span<const std::uint32_t>&,
48 std::int32_t,
int)>& dof_transform,
49 mdspan2_t dofmap0,
int bs0,
50 const std::function<
void(
const std::span<T>&,
51 const std::span<const std::uint32_t>&,
52 std::int32_t,
int)>& dof_transform_to_transpose,
53 mdspan2_t dofmap1,
int bs1, std::span<const T> constants,
54 std::span<const T> coeffs,
int cstride,
55 std::span<const std::uint32_t> cell_info, std::span<const T> bc_values1,
56 std::span<const std::int8_t> bc_markers1, std::span<const T> x0, T scale)
58 assert(_bs0 < 0 or _bs0 == bs0);
59 assert(_bs1 < 0 or _bs1 == bs1);
65 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
66 std::vector<T> Ae, be;
67 for (std::size_t index = 0; index <
cells.size(); ++index)
69 std::int32_t c =
cells[index];
73 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
74 submdspan(dofmap1, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
78 for (std::size_t j = 0; j < dmap1.size(); ++j)
80 if constexpr (_bs1 > 0)
82 for (
int k = 0; k < _bs1; ++k)
84 assert(_bs1 * dmap1[j] + k < (
int)bc_markers1.size());
85 if (bc_markers1[_bs1 * dmap1[j] + k])
94 for (
int k = 0; k < bs1; ++k)
96 assert(bs1 * dmap1[j] + k < (
int)bc_markers1.size());
97 if (bc_markers1[bs1 * dmap1[j] + k])
111 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
112 submdspan(x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
113 for (std::size_t i = 0; i < x_dofs.size(); ++i)
115 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
116 std::next(coordinate_dofs.begin(), 3 * i));
121 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
122 submdspan(dofmap0, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
124 const int num_rows = bs0 * dmap0.size();
125 const int num_cols = bs1 * dmap1.size();
127 const T* coeff_array = coeffs.data() + index * cstride;
128 Ae.resize(num_rows * num_cols);
129 std::fill(Ae.begin(), Ae.end(), 0);
130 kernel(Ae.data(), coeff_array, constants.data(), coordinate_dofs.data(),
132 dof_transform(Ae, cell_info, c, num_cols);
133 dof_transform_to_transpose(Ae, cell_info, c, num_rows);
137 std::fill(be.begin(), be.end(), 0);
138 for (std::size_t j = 0; j < dmap1.size(); ++j)
140 if constexpr (_bs1 > 0)
142 for (
int k = 0; k < _bs1; ++k)
144 const std::int32_t jj = _bs1 * dmap1[j] + k;
145 assert(jj < (
int)bc_markers1.size());
148 const T bc = bc_values1[jj];
149 const T _x0 = x0.empty() ? 0.0 : x0[jj];
152 for (
int m = 0; m < num_rows; ++m)
153 be[m] -= Ae[m * num_cols + _bs1 * j + k] * scale * (bc - _x0);
159 for (
int k = 0; k < bs1; ++k)
161 const std::int32_t jj = bs1 * dmap1[j] + k;
162 assert(jj < (
int)bc_markers1.size());
165 const T bc = bc_values1[jj];
166 const T _x0 = x0.empty() ? 0.0 : x0[jj];
168 for (
int m = 0; m < num_rows; ++m)
169 be[m] -= Ae[m * num_cols + bs1 * j + k] * scale * (bc - _x0);
175 for (std::size_t i = 0; i < dmap0.size(); ++i)
177 if constexpr (_bs0 > 0)
179 for (
int k = 0; k < _bs0; ++k)
180 b[_bs0 * dmap0[i] + k] += be[_bs0 * i + k];
184 for (
int k = 0; k < bs0; ++k)
185 b[bs0 * dmap0[i] + k] += be[bs0 * i + k];
198void _lift_bc_exterior_facets(
199 std::span<T> b, mdspan2_t x_dofmap,
200 std::span<
const scalar_value_type_t<T>> x, FEkernel<T>
auto kernel,
201 std::span<const std::int32_t> facets,
202 const std::function<
void(
const std::span<T>&,
203 const std::span<const std::uint32_t>&,
204 std::int32_t,
int)>& dof_transform,
205 mdspan2_t dofmap0,
int bs0,
206 const std::function<
void(
const std::span<T>&,
207 const std::span<const std::uint32_t>&,
208 std::int32_t,
int)>& dof_transform_to_transpose,
209 mdspan2_t dofmap1,
int bs1, std::span<const T> constants,
210 std::span<const T> coeffs,
int cstride,
211 std::span<const std::uint32_t> cell_info, std::span<const T> bc_values1,
212 std::span<const std::int8_t> bc_markers1, std::span<const T> x0, T scale)
218 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
219 std::vector<T> Ae, be;
220 assert(facets.size() % 2 == 0);
221 for (std::size_t index = 0; index < facets.size(); index += 2)
223 std::int32_t
cell = facets[index];
224 std::int32_t local_facet = facets[index + 1];
227 auto dmap1 = MDSPAN_IMPL_STANDARD_NAMESPACE::
228 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
229 dofmap1,
cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
233 for (std::size_t j = 0; j < dmap1.size(); ++j)
235 for (
int k = 0; k < bs1; ++k)
237 if (bc_markers1[bs1 * dmap1[j] + k])
249 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::
250 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
251 x_dofmap,
cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
252 for (std::size_t i = 0; i < x_dofs.size(); ++i)
254 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
255 std::next(coordinate_dofs.begin(), 3 * i));
259 auto dmap0 = MDSPAN_IMPL_STANDARD_NAMESPACE::
260 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
261 dofmap0,
cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
263 const int num_rows = bs0 * dmap0.size();
264 const int num_cols = bs1 * dmap1.size();
266 const T* coeff_array = coeffs.data() + index / 2 * cstride;
267 Ae.resize(num_rows * num_cols);
268 std::fill(Ae.begin(), Ae.end(), 0);
269 kernel(Ae.data(), coeff_array, constants.data(), coordinate_dofs.data(),
270 &local_facet,
nullptr);
271 dof_transform(Ae, cell_info,
cell, num_cols);
272 dof_transform_to_transpose(Ae, cell_info,
cell, num_rows);
276 std::fill(be.begin(), be.end(), 0);
277 for (std::size_t j = 0; j < dmap1.size(); ++j)
279 for (
int k = 0; k < bs1; ++k)
281 const std::int32_t jj = bs1 * dmap1[j] + k;
284 const T bc = bc_values1[jj];
285 const T _x0 = x0.empty() ? 0.0 : x0[jj];
287 for (
int m = 0; m < num_rows; ++m)
288 be[m] -= Ae[m * num_cols + bs1 * j + k] * scale * (bc - _x0);
293 for (std::size_t i = 0; i < dmap0.size(); ++i)
294 for (
int k = 0; k < bs0; ++k)
295 b[bs0 * dmap0[i] + k] += be[bs0 * i + k];
306void _lift_bc_interior_facets(
307 std::span<T> b, mdspan2_t x_dofmap,
308 std::span<
const scalar_value_type_t<T>> x,
int num_cell_facets,
309 FEkernel<T>
auto kernel, std::span<const std::int32_t> facets,
310 const std::function<
void(
const std::span<T>&,
311 const std::span<const std::uint32_t>&,
312 std::int32_t,
int)>& dof_transform,
313 mdspan2_t dofmap0,
int bs0,
314 const std::function<
void(
const std::span<T>&,
315 const std::span<const std::uint32_t>&,
316 std::int32_t,
int)>& dof_transform_to_transpose,
317 mdspan2_t dofmap1,
int bs1, std::span<const T> constants,
318 std::span<const T> coeffs,
int cstride,
319 std::span<const std::uint32_t> cell_info,
320 const std::function<std::uint8_t(std::size_t)>& get_perm,
321 std::span<const T> bc_values1, std::span<const std::int8_t> bc_markers1,
322 std::span<const T> x0, T scale)
328 using X = scalar_value_type_t<T>;
329 std::vector<X> coordinate_dofs(2 * x_dofmap.extent(1) * 3);
330 std::span<X> cdofs0(coordinate_dofs.data(), x_dofmap.extent(1) * 3);
331 std::span<X> cdofs1(coordinate_dofs.data() + x_dofmap.extent(1) * 3,
332 x_dofmap.extent(1) * 3);
333 std::vector<T> Ae, be;
336 std::vector<std::int32_t> dmapjoint0, dmapjoint1;
337 assert(facets.size() % 4 == 0);
339 const int num_dofs0 = dofmap0.extent(1);
340 const int num_dofs1 = dofmap1.extent(1);
341 for (std::size_t index = 0; index < facets.size(); index += 4)
343 std::array<std::int32_t, 2>
cells = {facets[index], facets[index + 2]};
344 std::array<std::int32_t, 2> local_facet
345 = {facets[index + 1], facets[index + 3]};
348 auto x_dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::
349 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
350 x_dofmap, cells[0], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
351 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
353 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
354 std::next(cdofs0.begin(), 3 * i));
356 auto x_dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::
357 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
358 x_dofmap, cells[1], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
359 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
361 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
362 std::next(cdofs1.begin(), 3 * i));
367 = std::span(dofmap0.data_handle() + cells[0] * num_dofs0, num_dofs0);
369 = std::span(dofmap1.data_handle() + cells[1] * num_dofs1, num_dofs1);
371 dmapjoint0.resize(dmap0_cell0.size() + dmap0_cell1.size());
372 std::copy(dmap0_cell0.begin(), dmap0_cell0.end(), dmapjoint0.begin());
373 std::copy(dmap0_cell1.begin(), dmap0_cell1.end(),
374 std::next(dmapjoint0.begin(), dmap0_cell0.size()));
377 = std::span(dofmap1.data_handle() + cells[0] * num_dofs1, num_dofs1);
379 = std::span(dofmap1.data_handle() + cells[1] * num_dofs1, num_dofs1);
381 dmapjoint1.resize(dmap1_cell0.size() + dmap1_cell1.size());
382 std::copy(dmap1_cell0.begin(), dmap1_cell0.end(), dmapjoint1.begin());
383 std::copy(dmap1_cell1.begin(), dmap1_cell1.end(),
384 std::next(dmapjoint1.begin(), dmap1_cell0.size()));
388 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
390 for (
int k = 0; k < bs1; ++k)
392 if (bc_markers1[bs1 * dmap1_cell0[j] + k])
401 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
403 for (
int k = 0; k < bs1; ++k)
405 if (bc_markers1[bs1 * dmap1_cell1[j] + k])
416 const int num_rows = bs0 * dmapjoint0.size();
417 const int num_cols = bs1 * dmapjoint1.size();
420 Ae.resize(num_rows * num_cols);
421 std::fill(Ae.begin(), Ae.end(), 0);
422 const std::array perm{
423 get_perm(cells[0] * num_cell_facets + local_facet[0]),
424 get_perm(cells[1] * num_cell_facets + local_facet[1])};
425 kernel(Ae.data(), coeffs.data() + index / 2 * cstride, constants.data(),
426 coordinate_dofs.data(), local_facet.data(), perm.data());
428 std::span<T> _Ae(Ae);
429 std::span<T> sub_Ae0 = _Ae.subspan(bs0 * dmap0_cell0.size() * num_cols,
430 bs0 * dmap0_cell1.size() * num_cols);
432 = _Ae.subspan(bs1 * dmap1_cell0.size(),
433 num_rows * num_cols - bs1 * dmap1_cell0.size());
435 dof_transform(_Ae, cell_info, cells[0], num_cols);
436 dof_transform(sub_Ae0, cell_info, cells[1], num_cols);
437 dof_transform_to_transpose(_Ae, cell_info, cells[0], num_rows);
438 dof_transform_to_transpose(sub_Ae1, cell_info, cells[1], num_rows);
441 std::fill(be.begin(), be.end(), 0);
444 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
446 for (
int k = 0; k < bs1; ++k)
448 const std::int32_t jj = bs1 * dmap1_cell0[j] + k;
451 const T bc = bc_values1[jj];
452 const T _x0 = x0.empty() ? 0.0 : x0[jj];
454 for (
int m = 0; m < num_rows; ++m)
455 be[m] -= Ae[m * num_cols + bs1 * j + k] * scale * (bc - _x0);
461 const int offset = bs1 * dmap1_cell0.size();
462 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
464 for (
int k = 0; k < bs1; ++k)
466 const std::int32_t jj = bs1 * dmap1_cell1[j] + k;
469 const T bc = bc_values1[jj];
470 const T _x0 = x0.empty() ? 0.0 : x0[jj];
472 for (
int m = 0; m < num_rows; ++m)
475 -= Ae[m * num_cols + offset + bs1 * j + k] * scale * (bc - _x0);
481 for (std::size_t i = 0; i < dmap0_cell0.size(); ++i)
482 for (
int k = 0; k < bs0; ++k)
483 b[bs0 * dmap0_cell0[i] + k] += be[bs0 * i + k];
485 const int offset_be = bs0 * dmap0_cell0.size();
486 for (std::size_t i = 0; i < dmap0_cell1.size(); ++i)
487 for (
int k = 0; k < bs0; ++k)
488 b[bs0 * dmap0_cell1[i] + k] += be[offset_be + bs0 * i + k];
499 const std::function<
void(
const std::span<T>&,
500 const std::span<const std::uint32_t>&,
501 std::int32_t,
int)>& dof_transform,
502 std::span<T> b, mdspan2_t x_dofmap,
503 std::span<
const scalar_value_type_t<T>> x,
504 std::span<const std::int32_t> cells, mdspan2_t dofmap,
int bs,
505 FEkernel<T>
auto kernel, std::span<const T> constants,
506 std::span<const T> coeffs,
int cstride,
507 std::span<const std::uint32_t> cell_info)
509 assert(_bs < 0 or _bs == bs);
516 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
517 std::vector<T> be(bs * dofmap.extent(1));
518 std::span<T> _be(be);
521 for (std::size_t index = 0; index <
cells.size(); ++index)
523 std::int32_t c =
cells[index];
527 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
528 submdspan(x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
529 for (std::size_t i = 0; i < x_dofs.size(); ++i)
531 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
532 std::next(coordinate_dofs.begin(), 3 * i));
536 std::fill(be.begin(), be.end(), 0);
537 kernel(be.data(), coeffs.data() + index * cstride, constants.data(),
538 coordinate_dofs.data(),
nullptr,
nullptr);
539 dof_transform(_be, cell_info, c, 1);
542 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
543 submdspan(dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
544 if constexpr (_bs > 0)
546 for (std::size_t i = 0; i < dofs.size(); ++i)
547 for (
int k = 0; k < _bs; ++k)
548 b[_bs * dofs[i] + k] += be[_bs * i + k];
552 for (std::size_t i = 0; i < dofs.size(); ++i)
553 for (
int k = 0; k < bs; ++k)
554 b[bs * dofs[i] + k] += be[bs * i + k];
566void assemble_exterior_facets(
567 const std::function<
void(
const std::span<T>&,
568 const std::span<const std::uint32_t>&,
569 std::int32_t,
int)>& dof_transform,
570 std::span<T> b, mdspan2_t x_dofmap,
571 std::span<
const scalar_value_type_t<T>> x,
572 std::span<const std::int32_t> facets, mdspan2_t dofmap,
int bs,
573 FEkernel<T>
auto fn, std::span<const T> constants,
574 std::span<const T> coeffs,
int cstride,
575 std::span<const std::uint32_t> cell_info)
577 assert(_bs < 0 or _bs == bs);
584 const int num_dofs = dofmap.extent(1);
585 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
586 std::vector<T> be(bs * num_dofs);
587 std::span<T> _be(be);
588 assert(facets.size() % 2 == 0);
589 for (std::size_t index = 0; index < facets.size(); index += 2)
591 std::int32_t
cell = facets[index];
592 std::int32_t local_facet = facets[index + 1];
595 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::
596 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
597 x_dofmap,
cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
598 for (std::size_t i = 0; i < x_dofs.size(); ++i)
600 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
601 std::next(coordinate_dofs.begin(), 3 * i));
605 std::fill(be.begin(), be.end(), 0);
606 fn(be.data(), coeffs.data() + index / 2 * cstride, constants.data(),
607 coordinate_dofs.data(), &local_facet,
nullptr);
609 dof_transform(_be, cell_info,
cell, 1);
612 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
613 submdspan(dofmap,
cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
614 if constexpr (_bs > 0)
616 for (std::size_t i = 0; i < dofs.size(); ++i)
617 for (
int k = 0; k < _bs; ++k)
618 b[_bs * dofs[i] + k] += be[_bs * i + k];
622 for (std::size_t i = 0; i < dofs.size(); ++i)
623 for (
int k = 0; k < bs; ++k)
624 b[bs * dofs[i] + k] += be[bs * i + k];
636void assemble_interior_facets(
637 const std::function<
void(
const std::span<T>&,
638 const std::span<const std::uint32_t>&,
639 std::int32_t,
int)>& dof_transform,
640 std::span<T> b, mdspan2_t x_dofmap,
641 std::span<
const scalar_value_type_t<T>> x,
int num_cell_facets,
642 std::span<const std::int32_t> facets,
const fem::DofMap& dofmap,
643 FEkernel<T>
auto fn, std::span<const T> constants,
644 std::span<const T> coeffs,
int cstride,
645 std::span<const std::uint32_t> cell_info,
646 const std::function<std::uint8_t(std::size_t)>& get_perm)
649 using X = scalar_value_type_t<T>;
650 std::vector<X> coordinate_dofs(2 * x_dofmap.extent(1) * 3);
651 std::span<X> cdofs0(coordinate_dofs.data(), x_dofmap.extent(1) * 3);
652 std::span<X> cdofs1(coordinate_dofs.data() + x_dofmap.extent(1) * 3,
653 x_dofmap.extent(1) * 3);
656 const int bs = dofmap.bs();
657 assert(_bs < 0 or _bs == bs);
658 assert(facets.size() % 4 == 0);
659 for (std::size_t index = 0; index < facets.size(); index += 4)
661 std::array<std::int32_t, 2>
cells = {facets[index], facets[index + 2]};
662 std::array<std::int32_t, 2> local_facet
663 = {facets[index + 1], facets[index + 3]};
666 auto x_dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::
667 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
668 x_dofmap, cells[0], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
669 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
671 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
672 std::next(cdofs0.begin(), 3 * i));
674 auto x_dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::
675 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
676 x_dofmap, cells[1], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
677 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
679 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
680 std::next(cdofs1.begin(), 3 * i));
684 std::span<const std::int32_t> dmap0 = dofmap.cell_dofs(cells[0]);
685 std::span<const std::int32_t> dmap1 = dofmap.cell_dofs(cells[1]);
688 be.resize(bs * (dmap0.size() + dmap1.size()));
689 std::fill(be.begin(), be.end(), 0);
690 const std::array perm{
691 get_perm(cells[0] * num_cell_facets + local_facet[0]),
692 get_perm(cells[1] * num_cell_facets + local_facet[1])};
693 fn(be.data(), coeffs.data() + index / 2 * cstride, constants.data(),
694 coordinate_dofs.data(), local_facet.data(), perm.data());
696 std::span<T> _be(be);
697 std::span<T> sub_be = _be.subspan(bs * dmap0.size(), bs * dmap1.size());
699 dof_transform(be, cell_info, cells[0], 1);
700 dof_transform(sub_be, cell_info, cells[1], 1);
703 if constexpr (_bs > 0)
705 for (std::size_t i = 0; i < dmap0.size(); ++i)
706 for (
int k = 0; k < _bs; ++k)
707 b[_bs * dmap0[i] + k] += be[_bs * i + k];
708 for (std::size_t i = 0; i < dmap1.size(); ++i)
709 for (
int k = 0; k < _bs; ++k)
710 b[_bs * dmap1[i] + k] += be[_bs * (i + dmap0.size()) + k];
714 for (std::size_t i = 0; i < dmap0.size(); ++i)
715 for (
int k = 0; k < bs; ++k)
716 b[bs * dmap0[i] + k] += be[bs * i + k];
717 for (std::size_t i = 0; i < dmap1.size(); ++i)
718 for (
int k = 0; k < bs; ++k)
719 b[bs * dmap1[i] + k] += be[bs * (i + dmap0.size()) + k];
740template <dolfinx::scalar T, std::
floating_po
int U>
741void lift_bc(std::span<T> b,
const Form<T, U>& a, mdspan2_t x_dofmap,
742 std::span<
const scalar_value_type_t<T>> x,
743 std::span<const T> constants,
744 const std::map<std::pair<IntegralType, int>,
745 std::pair<std::span<const T>,
int>>& coefficients,
746 std::span<const T> bc_values1,
747 std::span<const std::int8_t> bc_markers1, std::span<const T> x0,
750 std::shared_ptr<const mesh::Mesh<U>> mesh = a.mesh();
754 assert(a.function_spaces().at(0));
755 assert(a.function_spaces().at(1));
756 auto dofmap0 = a.function_spaces()[0]->dofmap()->map();
757 const int bs0 = a.function_spaces()[0]->dofmap()->bs();
758 auto element0 = a.function_spaces()[0]->element();
759 auto dofmap1 = a.function_spaces()[1]->dofmap()->map();
760 const int bs1 = a.function_spaces()[1]->dofmap()->bs();
761 auto element1 = a.function_spaces()[1]->element();
764 const bool needs_transformation_data
765 = element0->needs_dof_transformations()
766 or element1->needs_dof_transformations()
767 or a.needs_facet_permutations();
769 std::span<const std::uint32_t> cell_info;
770 if (needs_transformation_data)
772 mesh->topology_mutable()->create_entity_permutations();
773 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
776 const std::function<void(
const std::span<T>&,
777 const std::span<const std::uint32_t>&, std::int32_t,
779 dof_transform = element0->template get_dof_transformation_function<T>();
780 const std::function<void(
const std::span<T>&,
781 const std::span<const std::uint32_t>&, std::int32_t,
783 dof_transform_to_transpose
784 = element1->template get_dof_transformation_to_transpose_function<T>();
792 if (bs0 == 1 and bs1 == 1)
794 _lift_bc_cells<T, 1, 1>(b, x_dofmap, x, kernel, cells, dof_transform,
795 dofmap0, bs0, dof_transform_to_transpose, dofmap1,
796 bs1, constants, coeffs, cstride, cell_info,
797 bc_values1, bc_markers1, x0, scale);
799 else if (bs0 == 3 and bs1 == 3)
801 _lift_bc_cells<T, 3, 3>(b, x_dofmap, x, kernel, cells, dof_transform,
802 dofmap0, bs0, dof_transform_to_transpose, dofmap1,
803 bs1, constants, coeffs, cstride, cell_info,
804 bc_values1, bc_markers1, x0, scale);
808 _lift_bc_cells(b, x_dofmap, x, kernel, cells, dof_transform, dofmap0, bs0,
809 dof_transform_to_transpose, dofmap1, bs1, constants,
810 coeffs, cstride, cell_info, bc_values1, bc_markers1, x0,
819 auto& [coeffs, cstride]
821 _lift_bc_exterior_facets(
823 dof_transform, dofmap0, bs0, dof_transform_to_transpose, dofmap1, bs1,
824 constants, coeffs, cstride, cell_info, bc_values1, bc_markers1, x0,
830 std::function<std::uint8_t(std::size_t)> get_perm;
831 if (a.needs_facet_permutations())
833 mesh->topology_mutable()->create_entity_permutations();
834 const std::vector<std::uint8_t>& perms
835 = mesh->topology()->get_facet_permutations();
836 get_perm = [&perms](std::size_t i) {
return perms[i]; };
839 get_perm = [](std::size_t) {
return 0; };
841 auto cell_types = mesh->topology()->cell_types();
842 if (cell_types.size() > 1)
843 throw std::runtime_error(
"Multiple cell types in the assembler");
845 mesh->topology()->dim() - 1);
850 auto& [coeffs, cstride]
852 _lift_bc_interior_facets(
853 b, x_dofmap, x, num_cell_facets, kernel,
855 bs0, dof_transform_to_transpose, dofmap1, bs1, constants, coeffs,
856 cstride, cell_info, get_perm, bc_values1, bc_markers1, x0, scale);
882template <dolfinx::scalar T, std::
floating_po
int U>
884 std::span<T> b,
const std::vector<std::shared_ptr<
const Form<T, U>>> a,
885 mdspan2_t x_dofmap, std::span<
const scalar_value_type_t<T>> x,
886 const std::vector<std::span<const T>>& constants,
887 const std::vector<std::map<std::pair<IntegralType, int>,
888 std::pair<std::span<const T>,
int>>>& coeffs,
889 const std::vector<std::vector<std::shared_ptr<
const DirichletBC<T, U>>>>&
891 const std::vector<std::span<const T>>& x0, T scale)
894 if (!x0.empty() and x0.size() != a.size())
896 throw std::runtime_error(
897 "Mismatch in size between x0 and bilinear form in assembler.");
900 if (a.size() != bcs1.size())
902 throw std::runtime_error(
903 "Mismatch in size between a and bcs in assembler.");
906 for (std::size_t j = 0; j < a.size(); ++j)
908 std::vector<std::int8_t> bc_markers1;
909 std::vector<T> bc_values1;
910 if (a[j] and !bcs1[j].empty())
912 assert(a[j]->function_spaces().at(0));
914 auto V1 = a[j]->function_spaces()[1];
916 auto map1 = V1->dofmap()->index_map;
917 const int bs1 = V1->dofmap()->index_map_bs();
919 const int crange = bs1 * (map1->size_local() + map1->num_ghosts());
920 bc_markers1.assign(crange,
false);
921 bc_values1.assign(crange, 0.0);
922 for (
const std::shared_ptr<
const DirichletBC<T, U>>& bc : bcs1[j])
924 bc->mark_dofs(bc_markers1);
925 bc->dof_values(bc_values1);
930 lift_bc<T>(b, *a[j], x_dofmap, x, constants[j], coeffs[j], bc_values1,
931 bc_markers1, x0[j], scale);
935 lift_bc<T>(b, *a[j], x_dofmap, x, constants[j], coeffs[j], bc_values1,
936 bc_markers1, std::span<const T>(), scale);
950template <dolfinx::scalar T, std::
floating_po
int U>
952 std::span<T> b,
const Form<T, U>& L, mdspan2_t x_dofmap,
953 std::span<
const scalar_value_type_t<T>> x, std::span<const T> constants,
954 const std::map<std::pair<IntegralType, int>,
955 std::pair<std::span<const T>,
int>>& coefficients)
957 std::shared_ptr<const mesh::Mesh<U>> mesh = L.mesh();
961 assert(L.function_spaces().at(0));
962 auto element = L.function_spaces().at(0)->element();
964 std::shared_ptr<const fem::DofMap> dofmap
965 = L.function_spaces().at(0)->dofmap();
967 auto dofs = dofmap->map();
968 const int bs = dofmap->bs();
970 const std::function<void(
const std::span<T>&,
971 const std::span<const std::uint32_t>&, std::int32_t,
973 dof_transform = element->template get_dof_transformation_function<T>();
975 const bool needs_transformation_data
976 = element->needs_dof_transformations() or L.needs_facet_permutations();
977 std::span<const std::uint32_t> cell_info;
978 if (needs_transformation_data)
980 mesh->topology_mutable()->create_entity_permutations();
981 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
992 impl::assemble_cells<T, 1>(dof_transform, b, x_dofmap, x, cells, dofs, bs,
993 fn, constants, coeffs, cstride, cell_info);
997 impl::assemble_cells<T, 3>(dof_transform, b, x_dofmap, x, cells, dofs, bs,
998 fn, constants, coeffs, cstride, cell_info);
1002 impl::assemble_cells(dof_transform, b, x_dofmap, x, cells, dofs, bs, fn,
1003 constants, coeffs, cstride, cell_info);
1011 auto& [coeffs, cstride]
1013 std::span<const std::int32_t> facets
1017 impl::assemble_exterior_facets<T, 1>(dof_transform, b, x_dofmap, x,
1018 facets, dofs, bs, fn, constants,
1019 coeffs, cstride, cell_info);
1023 impl::assemble_exterior_facets<T, 3>(dof_transform, b, x_dofmap, x,
1024 facets, dofs, bs, fn, constants,
1025 coeffs, cstride, cell_info);
1029 impl::assemble_exterior_facets(dof_transform, b, x_dofmap, x, facets,
1030 dofs, bs, fn, constants, coeffs, cstride,
1037 std::function<std::uint8_t(std::size_t)> get_perm;
1038 if (L.needs_facet_permutations())
1040 mesh->topology_mutable()->create_entity_permutations();
1041 const std::vector<std::uint8_t>& perms
1042 = mesh->topology()->get_facet_permutations();
1043 get_perm = [&perms](std::size_t i) {
return perms[i]; };
1046 get_perm = [](std::size_t) {
return 0; };
1048 auto cell_types = mesh->topology()->cell_types();
1049 if (cell_types.size() > 1)
1050 throw std::runtime_error(
"Multiple cell types in the assembler");
1052 mesh->topology()->dim() - 1);
1057 auto& [coeffs, cstride]
1059 std::span<const std::int32_t> facets
1063 impl::assemble_interior_facets<T, 1>(
1064 dof_transform, b, x_dofmap, x, num_cell_facets, facets, *dofmap, fn,
1065 constants, coeffs, cstride, cell_info, get_perm);
1069 impl::assemble_interior_facets<T, 3>(
1070 dof_transform, b, x_dofmap, x, num_cell_facets, facets, *dofmap, fn,
1071 constants, coeffs, cstride, cell_info, get_perm);
1075 impl::assemble_interior_facets(
1076 dof_transform, b, x_dofmap, x, num_cell_facets, facets, *dofmap, fn,
1077 constants, coeffs, cstride, cell_info, get_perm);
1089template <dolfinx::scalar T, std::
floating_po
int U>
1090void assemble_vector(
1091 std::span<T> b,
const Form<T, U>& L, std::span<const T> constants,
1092 const std::map<std::pair<IntegralType, int>,
1093 std::pair<std::span<const T>,
int>>& coefficients)
1095 std::shared_ptr<const mesh::Mesh<U>> mesh = L.mesh();
1097 if constexpr (std::is_same_v<U, scalar_value_type_t<T>>)
1098 assemble_vector(b, L, mesh->geometry().dofmap(), mesh->geometry().x(),
1099 constants, coefficients);
1102 auto x = mesh->geometry().x();
1103 std::vector<scalar_value_type_t<T>> _x(x.begin(), x.end());
1104 assemble_vector(b, L, mesh->geometry().dofmap(), _x, constants,
Degree-of-freedom map representations and tools.
This concept is used to constrain the a template type to floating point real or complex types....
Definition types.h:20
void cells(la::SparsityPattern &pattern, std::span< const std::int32_t > cells, std::array< std::reference_wrapper< const DofMap >, 2 > dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition sparsitybuild.cpp:18
IntegralType
Type of integral.
Definition Form.h:33
@ interior_facet
Interior facet.
@ exterior_facet
Exterior facet.
int cell_num_entities(CellType type, int dim)
Number of entities of dimension dim.
Definition cell_types.cpp:139