DOLFINx 0.10.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
assemble_vector_impl.h
1// Copyright (C) 2018-2025 Garth N. Wells and Paul T. Kühner
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include "Constant.h"
10#include "DirichletBC.h"
11#include "DofMap.h"
12#include "Form.h"
13#include "traits.h"
14#include "utils.h"
15#include <algorithm>
16#include <basix/mdspan.hpp>
17#include <cstdint>
18#include <dolfinx/common/IndexMap.h>
19#include <dolfinx/mesh/Geometry.h>
20#include <dolfinx/mesh/Mesh.h>
21#include <dolfinx/mesh/Topology.h>
22#include <functional>
23#include <memory>
24#include <optional>
25#include <span>
26#include <vector>
27
28namespace dolfinx::fem
29{
30template <dolfinx::scalar T, std::floating_point U>
31class DirichletBC;
32
33}
34namespace dolfinx::fem::impl
35{
37using mdspan2_t = md::mdspan<const std::int32_t, md::dextents<std::size_t, 2>>;
39
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>
82void _lift_bc_cells(
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>>
86 x,
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)
96{
97 if (cells.empty())
98 return;
99
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);
104
105 const int num_rows = bs0 * dmap0.extent(1);
106 const int num_cols = bs1 * dmap1.extent(1);
107
108 // Data structures used in bc application
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)
114 {
115 // Cell index in integration domain mesh, test function mesh, and trial
116 // function mesh
117 std::int32_t c = cells[index];
118 std::int32_t c0 = cells0[index];
119 std::int32_t c1 = cells1[index];
120
121 // Get dof maps for cell
122 auto dofs1 = md::submdspan(dmap1, c1, md::full_extent);
123
124 // Check if bc is applied to cell
125 bool has_bc = false;
126 for (std::size_t j = 0; j < dofs1.size(); ++j)
127 {
128 if constexpr (_bs1 > 0)
129 {
130 for (int k = 0; k < _bs1; ++k)
131 {
132 assert(_bs1 * dofs1[j] + k < (int)bc_markers1.size());
133 if (bc_markers1[_bs1 * dofs1[j] + k])
134 {
135 has_bc = true;
136 break;
137 }
138 }
139 }
140 else
141 {
142 for (int k = 0; k < bs1; ++k)
143 {
144 assert(bs1 * dofs1[j] + k < (int)bc_markers1.size());
145 if (bc_markers1[bs1 * dofs1[j] + k])
146 {
147 has_bc = true;
148 break;
149 }
150 }
151 }
152 }
153
154 if (!has_bc)
155 continue;
156
157 // Get cell coordinates/geometry
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));
161
162 // Size data structure for assembly
163 auto dofs0 = md::submdspan(dmap0, c0, md::full_extent);
164
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);
170
171 // Size data structure for assembly
172 std::ranges::fill(be, 0);
173 for (std::size_t j = 0; j < dofs1.size(); ++j)
174 {
175 if constexpr (_bs1 > 0)
176 {
177 for (int k = 0; k < _bs1; ++k)
178 {
179 const std::int32_t jj = _bs1 * dofs1[j] + k;
180 assert(jj < (int)bc_markers1.size());
181 if (bc_markers1[jj])
182 {
183 const T bc = bc_values1[jj];
184 const T _x0 = x0.empty() ? 0 : x0[jj];
185 // const T _x0 = 0;
186 // be -= Ae.col(bs1 * j + k) * alpha * (bc - _x0);
187 for (int m = 0; m < num_rows; ++m)
188 be[m] -= Ae[m * num_cols + _bs1 * j + k] * alpha * (bc - _x0);
189 }
190 }
191 }
192 else
193 {
194 for (int k = 0; k < bs1; ++k)
195 {
196 const std::int32_t jj = bs1 * dofs1[j] + k;
197 assert(jj < (int)bc_markers1.size());
198 if (bc_markers1[jj])
199 {
200 const T bc = bc_values1[jj];
201 const T _x0 = x0.empty() ? 0 : x0[jj];
202 // be -= Ae.col(bs1 * j + k) * alpha * (bc - _x0);
203 for (int m = 0; m < num_rows; ++m)
204 be[m] -= Ae[m * num_cols + bs1 * j + k] * alpha * (bc - _x0);
205 }
206 }
207 }
208 }
209
210 for (std::size_t i = 0; i < dofs0.size(); ++i)
211 {
212 if constexpr (_bs0 > 0)
213 {
214 for (int k = 0; k < _bs0; ++k)
215 b[_bs0 * dofs0[i] + k] += be[_bs0 * i + k];
216 }
217 else
218 {
219 for (int k = 0; k < bs0; ++k)
220 b[bs0 * dofs0[i] + k] += be[bs0 * i + k];
221 }
222 }
223 }
224}
225
263template <typename V,
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>>
270 x,
271 FEkernel<T> auto kernel,
272 md::mdspan<const std::int32_t,
273 md::extents<std::size_t, md::dynamic_extent, 2>>
274 entities,
275 std::tuple<mdspan2_t, int,
276 md::mdspan<const std::int32_t,
277 md::extents<std::size_t, md::dynamic_extent, 2>>>
278 dofmap0,
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>>>
283 dofmap1,
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)
290{
291 if (entities.empty())
292 return;
293
294 const auto [dmap0, bs0, entities0] = dofmap0;
295 const auto [dmap1, bs1, entities1] = dofmap1;
296
297 const int num_rows = bs0 * dmap0.extent(1);
298 const int num_cols = bs1 * dmap1.extent(1);
299
300 // Data structures used in bc application
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)
306 {
307 // Cell in integration domain, test function and trial function
308 // meshes
309 std::int32_t cell = entities(index, 0);
310 std::int32_t cell0 = entities0(index, 0);
311 std::int32_t cell1 = entities1(index, 0);
312
313 // Local entity index
314 std::int32_t local_entity = entities(index, 1);
315
316 // Get dof maps for cell
317 auto dofs1 = md::submdspan(dmap1, cell1, md::full_extent);
318
319 // Check if bc is applied to cell
320 bool has_bc = false;
321 for (std::size_t j = 0; j < dofs1.size(); ++j)
322 {
323 for (int k = 0; k < bs1; ++k)
324 {
325 if (bc_markers1[bs1 * dofs1[j] + k])
326 {
327 has_bc = true;
328 break;
329 }
330 }
331 }
332
333 if (!has_bc)
334 continue;
335
336 // Get cell coordinates/geometry
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));
340
341 // Size data structure for assembly
342 auto dofs0 = md::submdspan(dmap0, cell0, md::full_extent);
343
344 // Permutations
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);
351
352 // Size data structure for assembly
353 std::ranges::fill(be, 0);
354 for (std::size_t j = 0; j < dofs1.size(); ++j)
355 {
356 for (int k = 0; k < bs1; ++k)
357 {
358 const std::int32_t jj = bs1 * dofs1[j] + k;
359 if (bc_markers1[jj])
360 {
361 const T bc = bc_values1[jj];
362 const T _x0 = x0.empty() ? 0 : x0[jj];
363 // be -= Ae.col(bs1 * j + k) * alpha * (bc - _x0);
364 for (int m = 0; m < num_rows; ++m)
365 be[m] -= Ae[m * num_cols + bs1 * j + k] * alpha * (bc - _x0);
366 }
367 }
368 }
369
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];
373 }
374}
375
418template <typename V,
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>>
425 x,
426 FEkernel<T> auto kernel,
427 md::mdspan<const std::int32_t,
428 md::extents<std::size_t, md::dynamic_extent, 2, 2>>
429 facets,
430 std::tuple<mdspan2_t, int,
431 md::mdspan<const std::int32_t,
432 md::extents<std::size_t, md::dynamic_extent, 2, 2>>>
433 dofmap0,
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>>>
438 dofmap1,
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,
441 md::dynamic_extent>>
442 coeffs,
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)
447{
448 if (facets.empty())
449 return;
450
451 const auto [dmap0, bs0, facets0] = dofmap0;
452 const auto [dmap1, bs1, facets1] = dofmap1;
453
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;
458
459 // Data structures used in assembly
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);
466
467 // Temporaries for joint dofmaps
468 std::vector<std::int32_t> dmapjoint0(2 * num_dofs0);
469 std::vector<std::int32_t> dmapjoint1(2 * num_dofs1);
470
471 assert(facets0.size() == facets.size());
472 assert(facets1.size() == facets.size());
473 for (std::size_t f = 0; f < facets.extent(0); ++f)
474 {
475 // Cells in integration domain, test function domain and trial
476 // function domain meshes
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)};
480
481 // Local facet indices
482 std::array local_facet = {facets(f, 0, 1), facets(f, 1, 1)};
483
484 // Get cell geometry
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));
491
492 // Get dof maps for cells and pack
493 // When integrating over interfaces between two domains, the test function
494 // might only be defined on one side, so we check which cells exist in the
495 // test function domain
496 std::span dmap0_cell0
497 = cells0[0] >= 0
498 ? std::span(dmap0.data_handle() + cells0[0] * num_dofs0,
499 num_dofs0)
500 : std::span<const std::int32_t>();
501 std::span dmap0_cell1
502 = cells0[1] >= 0
503 ? std::span(dmap0.data_handle() + cells0[1] * num_dofs0,
504 num_dofs0)
505 : std::span<const std::int32_t>();
506
507 std::ranges::copy(dmap0_cell0, dmapjoint0.begin());
508 std::ranges::copy(dmap0_cell1, std::next(dmapjoint0.begin(), num_dofs0));
509
510 // Check which cells exist in the trial function domain
511 std::span<const std::int32_t> dmap1_cell0
512 = cells1[0] >= 0
513 ? std::span(dmap1.data_handle() + cells1[0] * num_dofs1,
514 num_dofs1)
515 : std::span<const std::int32_t>();
516 std::span<const std::int32_t> dmap1_cell1
517 = cells1[1] >= 0
518 ? std::span(dmap1.data_handle() + cells1[1] * num_dofs1,
519 num_dofs1)
520 : std::span<const std::int32_t>();
521
522 std::ranges::copy(dmap1_cell0, dmapjoint1.begin());
523 std::ranges::copy(dmap1_cell1, std::next(dmapjoint1.begin(), num_dofs1));
524
525 // Check if bc is applied to cell0
526 bool has_bc = false;
527 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
528 {
529 for (int k = 0; k < bs1; ++k)
530 {
531 if (bc_markers1[bs1 * dmap1_cell0[j] + k])
532 {
533 has_bc = true;
534 break;
535 }
536 }
537 }
538
539 // Check if bc is applied to cell1
540 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
541 {
542 for (int k = 0; k < bs1; ++k)
543 {
544 if (bc_markers1[bs1 * dmap1_cell1[j] + k])
545 {
546 has_bc = true;
547 break;
548 }
549 }
550 }
551
552 if (!has_bc)
553 continue;
554
555 // Tabulate tensor
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);
563
564 if (cells0[0] >= 0)
565 P0(Ae, cell_info0, cells0[0], num_cols);
566 if (cells0[1] >= 0)
567 {
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);
571 }
572 if (cells1[0] >= 0)
573 P1T(Ae, cell_info1, cells1[0], num_rows);
574
575 if (cells1[1] >= 0)
576 {
577 for (int row = 0; row < num_rows; ++row)
578 {
579 // DOFs for dmap1 and cell1 are not stored contiguously in
580 // the block matrix, so each row needs a separate span access
581 std::span sub_Ae1(Ae.data() + row * num_cols + bs1 * num_dofs1,
582 bs1 * num_dofs1);
583 P1T(sub_Ae1, cell_info1, cells1[1], 1);
584 }
585 }
586
587 std::ranges::fill(be, 0);
588
589 // Compute b = b - A*b for cell0
590 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
591 {
592 for (int k = 0; k < bs1; ++k)
593 {
594 const std::int32_t jj = bs1 * dmap1_cell0[j] + k;
595 if (bc_markers1[jj])
596 {
597 const T bc = bc_values1[jj];
598 const T _x0 = x0.empty() ? 0 : x0[jj];
599 // be -= Ae.col(bs1 * j + k) * alpha * (bc - _x0);
600 for (int m = 0; m < num_rows; ++m)
601 be[m] -= Ae[m * num_cols + bs1 * j + k] * alpha * (bc - _x0);
602 }
603 }
604 }
605
606 // Compute b = b - A*b for cell1
607 const int offset = bs1 * num_dofs1;
608 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
609 {
610 for (int k = 0; k < bs1; ++k)
611 {
612 const std::int32_t jj = bs1 * dmap1_cell1[j] + k;
613 if (bc_markers1[jj])
614 {
615 const T bc = bc_values1[jj];
616 const T _x0 = x0.empty() ? 0 : x0[jj];
617 // be -= Ae.col(offset + bs1 * j + k) * alpha * (bc - x0[jj]);
618 for (int m = 0; m < num_rows; ++m)
619 {
620 be[m]
621 -= Ae[m * num_cols + offset + bs1 * j + k] * alpha * (bc - _x0);
622 }
623 }
624 }
625 }
626
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];
630
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];
635 }
636}
637
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>
664
665void assemble_cells(
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>>
669 x,
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)
675{
676 if (cells.empty())
677 return;
678
679 const auto [dmap, bs, cells0] = dofmap;
680 assert(_bs < 0 or _bs == bs);
681
682 // Create data structures used in assembly
683 std::vector<scalar_value_t<T>> cdofs(3 * x_dofmap.extent(1));
684 std::vector<T> be(bs * dmap.extent(1));
685
686 // Iterate over active cells
687 for (std::size_t index = 0; index < cells.size(); ++index)
688 {
689 // Integration domain celland test function cell
690 std::int32_t c = cells[index];
691 std::int32_t c0 = cells0[index];
692
693 // Get cell coordinates/geometry
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));
697
698 // Tabulate vector for cell
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);
703
704 // Scatter cell vector to 'global' vector array
705 auto dofs = md::submdspan(dmap, c0, md::full_extent);
706 if constexpr (_bs > 0)
707 {
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];
711 }
712 else
713 {
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];
717 }
718 }
719}
720
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>>
761 x,
762 md::mdspan<const std::int32_t,
763 std::extents<std::size_t, md::dynamic_extent, 2>>
764 entities,
765 std::tuple<mdspan2_t, int,
766 md::mdspan<const std::int32_t,
767 std::extents<std::size_t, md::dynamic_extent, 2>>>
768 dofmap,
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)
773{
774 if (entities.empty())
775 return;
776
777 const auto [dmap, bs, entities0] = dofmap;
778 assert(_bs < 0 or _bs == bs);
779
780 // Create data structures used in assembly
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)
786 {
787 // Cell in the integration domain, local facet index relative to the
788 // integration domain cell, and cell in the test function mesh
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);
792
793 // Get cell coordinates/geometry
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));
797
798 // Permutations
799 std::uint8_t perm = perms.empty() ? 0 : perms(cell, local_entity);
800
801 // Tabulate element vector
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);
806
807 // Add element vector to global vector
808 auto dofs = md::submdspan(dmap, cell0, md::full_extent);
809 if constexpr (_bs > 0)
810 {
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];
814 }
815 else
816 {
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];
820 }
821 }
822}
823
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>>
856 x,
857 md::mdspan<const std::int32_t,
858 std::extents<std::size_t, md::dynamic_extent, 2, 2>>
859 facets,
860 std::tuple<const DofMap&, int,
861 md::mdspan<const std::int32_t,
862 std::extents<std::size_t, md::dynamic_extent, 2, 2>>>
863 dofmap,
864 FEkernel<T> auto kernel, std::span<const T> constants,
865 md::mdspan<const T, md::extents<std::size_t, md::dynamic_extent, 2,
866 md::dynamic_extent>>
867 coeffs,
868 std::span<const std::uint32_t> cell_info0,
869 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> perms)
870{
871 using X = scalar_value_t<T>;
872
873 if (facets.empty())
874 return;
875
876 const auto [dmap, bs, facets0] = dofmap;
877 assert(_bs < 0 or _bs == bs);
878
879 // Create data structures used in assembly
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);
884
885 const std::size_t dmap_size = dmap.map().extent(1);
886 std::vector<T> be(bs * 2 * dmap_size);
887
888 assert(facets0.size() == facets.size());
889 for (std::size_t f = 0; f < facets.extent(0); ++f)
890 {
891 // Cells in integration domain and test function domain meshes
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)};
894
895 // Local facet indices
896 std::array<std::int32_t, 2> local_facet{facets(f, 0, 1), facets(f, 1, 1)};
897
898 // Get cell geometry
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));
905
906 // Get dofmaps for cells. When integrating over interfaces between
907 // two domains, the test function might only be defined on one side,
908 // so we check which cells exist in the test function domain.
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>();
913
914 // Tabulate element vector
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);
922
923 if (cells0[0] >= 0)
924 P0(be, cell_info0, cells0[0], 1);
925 if (cells0[1] >= 0)
926 {
927 std::span sub_be(be.data() + bs * dmap_size, bs * dmap_size);
928 P0(sub_be, cell_info0, cells0[1], 1);
929 }
930
931 // Add element vector to global vector
932 if constexpr (_bs > 0)
933 {
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];
940 }
941 else
942 {
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];
949 }
950 }
951}
952
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>>
975 x,
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,
981 T alpha)
982{
983 // Integration domain mesh
984 std::shared_ptr<const mesh::Mesh<U>> mesh = a.mesh();
985 assert(mesh);
986
987 // Test function mesh
988 auto mesh0 = a.function_spaces().at(0)->mesh();
989 assert(mesh0);
990
991 // Trial function mesh
992 auto mesh1 = a.function_spaces().at(1)->mesh();
993 assert(mesh1);
994
995 // Get dofmap for columns and rows of a
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();
1004 assert(element0);
1005
1006 std::span<const std::uint32_t> cell_info0;
1007 std::span<const std::uint32_t> cell_info1;
1008 // TODO: Check for each element instead
1009 if (element0->needs_dof_transformations()
1010 or element1->needs_dof_transformations() or a.needs_facet_permutations())
1011 {
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());
1016 }
1017
1018 fem::DofTransformKernel<T> auto P0
1019 = element0->template dof_transformation_fn<T>(doftransform::standard);
1020 fem::DofTransformKernel<T> auto P1T
1021 = element1->template dof_transformation_right_fn<T>(
1023
1024 for (int i = 0; i < a.num_integrals(IntegralType::cell, 0); ++i)
1025 {
1026 auto kernel = a.kernel(IntegralType::cell, i, 0);
1027 assert(kernel);
1028 auto& [_coeffs, cstride] = coefficients.at({IntegralType::cell, i});
1029 std::span cells = a.domain(IntegralType::cell, i, 0);
1030 std::span cells0 = a.domain_arg(IntegralType::cell, 0, i, 0);
1031 std::span cells1 = a.domain_arg(IntegralType::cell, 1, i, 0);
1032 assert(_coeffs.size() == cells.size() * cstride);
1033 auto coeffs = md::mdspan(_coeffs.data(), cells.size(), cstride);
1034 if (bs0 == 1 and bs1 == 1)
1035 {
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);
1040 }
1041 else if (bs0 == 3 and bs1 == 3)
1042 {
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);
1047 }
1048 else
1049 {
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);
1053 }
1054 }
1055
1056 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> facet_perms;
1057 if (a.needs_facet_permutations())
1058 {
1059 mesh::CellType cell_type = mesh->topology()->cell_type();
1060 int num_facets_per_cell
1061 = mesh::cell_num_entities(cell_type, mesh->topology()->dim() - 1);
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);
1067 }
1068
1069 for (int i = 0; i < a.num_integrals(IntegralType::interior_facet, 0); ++i)
1070 {
1071 auto kernel = a.kernel(IntegralType::interior_facet, i, 0);
1072 assert(kernel);
1073 auto& [coeffs, cstride]
1074 = coefficients.at({IntegralType::interior_facet, i});
1075
1076 using mdspanx22_t
1077 = md::mdspan<const std::int32_t,
1078 md::extents<std::size_t, md::dynamic_extent, 2, 2>>;
1079 using mdspanx2x_t
1080 = md::mdspan<const T, md::extents<std::size_t, md::dynamic_extent, 2,
1081 md::dynamic_extent>>;
1082 std::span f = a.domain(IntegralType::interior_facet, i, 0);
1083 mdspanx22_t facets(f.data(), f.size() / 4, 2, 2);
1084 std::span f0 = a.domain_arg(IntegralType::interior_facet, 0, i, 0);
1085 mdspanx22_t facets0(f0.data(), f0.size() / 4, 2, 2);
1086 std::span f1 = a.domain_arg(IntegralType::interior_facet, 1, i, 0);
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);
1093 }
1094
1095 for (auto itg_type : {fem::IntegralType::exterior_facet,
1097 {
1098 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> perms
1099 = (itg_type == fem::IntegralType::exterior_facet)
1100 ? facet_perms
1101 : md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>>{};
1102
1103 for (int i = 0; i < a.num_integrals(itg_type, 0); ++i)
1104 {
1105 auto kernel = a.kernel(itg_type, i, 0);
1106 assert(kernel);
1107 auto& [coeffs, cstride] = coefficients.at({itg_type, i});
1108
1109 using mdspanx2_t
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);
1119 _lift_bc_entities(
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);
1124 }
1125 }
1126}
1127
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>
1151void apply_lifting(
1152 V&& b,
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,
1157 const std::vector<
1158 std::vector<std::reference_wrapper<const DirichletBC<T, U>>>>& bcs1,
1159 const std::vector<std::span<const T>>& x0, T alpha)
1160{
1161 if (!x0.empty() and x0.size() != a.size())
1162 {
1163 throw std::runtime_error(
1164 "Mismatch in size between x0 and bilinear form in assembler.");
1165 }
1166
1167 if (a.size() != bcs1.size())
1168 {
1169 throw std::runtime_error(
1170 "Mismatch in size between a and bcs in assembler.");
1171 }
1172
1173 for (std::size_t j = 0; j < a.size(); ++j)
1174 {
1175 std::vector<std::int8_t> bc_markers1;
1176 std::vector<T> bc_values1;
1177 if (a[j] and !bcs1[j].empty())
1178 {
1179 // Extract data from mesh
1180 std::shared_ptr<const mesh::Mesh<U>> mesh = a[j]->get().mesh();
1181 if (!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);
1188
1189 assert(a[j]->get().function_spaces().at(0));
1190 auto V1 = a[j]->get().function_spaces()[1];
1191 assert(V1);
1192 auto map1 = V1->dofmap()->index_map;
1193 const int bs1 = V1->dofmap()->index_map_bs();
1194 assert(map1);
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])
1199 {
1200 bc.get().mark_dofs(bc_markers1);
1201 bc.get().set(bc_values1, std::nullopt, 1);
1202 }
1203
1204 if (!x0.empty())
1205 {
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);
1208 }
1209 else
1210 {
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);
1214 }
1215 }
1216 }
1217}
1218
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>>
1233 x,
1234 std::span<const T> constants,
1235 const std::map<std::pair<IntegralType, int>,
1236 std::pair<std::span<const T>, int>>& coefficients)
1237{
1238 // Integration domain mesh
1239 std::shared_ptr<const mesh::Mesh<U>> mesh = L.mesh();
1240 assert(mesh);
1241
1242 // Test function mesh
1243 auto mesh0 = L.function_spaces().at(0)->mesh();
1244 assert(mesh0);
1245
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)
1248 {
1249 // Geometry dofmap and data
1250 mdspan2_t x_dofmap = mesh->geometry().dofmap(cell_type_idx);
1251
1252 // Get dofmap data
1253 assert(L.function_spaces().at(0));
1254 auto element = L.function_spaces().at(0)->elements(cell_type_idx);
1255 assert(element);
1256 std::shared_ptr<const fem::DofMap> dofmap
1257 = L.function_spaces().at(0)->dofmaps(cell_type_idx);
1258 assert(dofmap);
1259 auto dofs = dofmap->map();
1260 const int bs = dofmap->bs();
1261
1262 fem::DofTransformKernel<T> auto P0
1263 = element->template dof_transformation_fn<T>(doftransform::standard);
1264
1265 std::span<const std::uint32_t> cell_info0;
1266 if (element->needs_dof_transformations() or L.needs_facet_permutations())
1267 {
1268 mesh0->topology_mutable()->create_entity_permutations();
1269 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
1270 }
1271
1272 for (int i = 0; i < L.num_integrals(IntegralType::cell, 0); ++i)
1273 {
1274 auto fn = L.kernel(IntegralType::cell, i, cell_type_idx);
1275 assert(fn);
1276 std::span cells = L.domain(IntegralType::cell, i, cell_type_idx);
1277 std::span cells0 = L.domain_arg(IntegralType::cell, 0, i, cell_type_idx);
1278 auto& [coeffs, cstride] = coefficients.at({IntegralType::cell, i});
1279 assert(cells.size() * cstride == coeffs.size());
1280 if (bs == 1)
1281 {
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);
1285 }
1286 else if (bs == 3)
1287 {
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);
1291 }
1292 else
1293 {
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);
1297 }
1298 }
1299
1300 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> facet_perms;
1301 if (L.needs_facet_permutations())
1302 {
1303 mesh::CellType cell_type = mesh->topology()->cell_types()[cell_type_idx];
1304 int num_facets_per_cell
1305 = mesh::cell_num_entities(cell_type, mesh->topology()->dim() - 1);
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);
1311 }
1312
1313 using mdspanx2_t
1314 = md::mdspan<const std::int32_t,
1315 md::extents<std::size_t, md::dynamic_extent, 2>>;
1316
1317 for (int i = 0; i < L.num_integrals(IntegralType::interior_facet, 0); ++i)
1318 {
1319 using mdspanx22_t
1320 = md::mdspan<const std::int32_t,
1321 md::extents<std::size_t, md::dynamic_extent, 2, 2>>;
1322 using mdspanx2x_t
1323 = md::mdspan<const T, md::extents<std::size_t, md::dynamic_extent, 2,
1324 md::dynamic_extent>>;
1325
1326 auto fn = L.kernel(IntegralType::interior_facet, i, 0);
1327 assert(fn);
1328 auto& [coeffs, cstride]
1329 = coefficients.at({IntegralType::interior_facet, i});
1330 std::span facets = L.domain(IntegralType::interior_facet, i, 0);
1331 std::span facets1 = L.domain_arg(IntegralType::interior_facet, 0, i, 0);
1332 assert((facets.size() / 4) * 2 * cstride == coeffs.size());
1333 if (bs == 1)
1334 {
1335 impl::assemble_interior_facets<1>(
1336 P0, b, x_dofmap, x,
1337 mdspanx22_t(facets.data(), facets.size() / 4, 2, 2),
1338 {*dofmap, bs,
1339 mdspanx22_t(facets1.data(), facets1.size() / 4, 2, 2)},
1340 fn, constants,
1341 mdspanx2x_t(coeffs.data(), facets.size() / 4, 2, cstride),
1342 cell_info0, facet_perms);
1343 }
1344 else if (bs == 3)
1345 {
1346 impl::assemble_interior_facets<3>(
1347 P0, b, x_dofmap, x,
1348 mdspanx22_t(facets.data(), facets.size() / 4, 2, 2),
1349 {*dofmap, bs,
1350 mdspanx22_t(facets1.data(), facets1.size() / 4, 2, 2)},
1351 fn, constants,
1352 mdspanx2x_t(coeffs.data(), facets.size() / 4, 2, cstride),
1353 cell_info0, facet_perms);
1354 }
1355 else
1356 {
1357 impl::assemble_interior_facets(
1358 P0, b, x_dofmap, x,
1359 mdspanx22_t(facets.data(), facets.size() / 4, 2, 2),
1360 {*dofmap, bs,
1361 mdspanx22_t(facets1.data(), facets1.size() / 4, 2, 2)},
1362 fn, constants,
1363 mdspanx2x_t(coeffs.data(), facets.size() / 4, 2, cstride),
1364 cell_info0, facet_perms);
1365 }
1366 }
1367
1368 for (auto itg_type : {fem::IntegralType::exterior_facet,
1370 {
1371 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> perms
1372 = (itg_type == fem::IntegralType::exterior_facet)
1373 ? facet_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)
1377 {
1378 auto fn = L.kernel(itg_type, i, 0);
1379 assert(fn);
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());
1386 if (bs == 1)
1387 {
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),
1391 cell_info0, perms);
1392 }
1393 else if (bs == 3)
1394 {
1395 impl::assemble_entities<3>(
1396 P0, b, x_dofmap, x, entities, {dofs, bs, entities1}, fn,
1397 constants,
1398 md::mdspan(coeffs.data(), entities.size() / 2, cstride),
1399 cell_info0, perms);
1400 }
1401 else
1402 {
1403 impl::assemble_entities(
1404 P0, b, x_dofmap, x, entities, {dofs, bs, entities1}, fn,
1405 constants,
1406 md::mdspan(coeffs.data(), entities.size() / 2, cstride),
1407 cell_info0, perms);
1408 }
1409 }
1410 }
1411 }
1412}
1413
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)
1427{
1428 using mdspanx3_t
1429 = md::mdspan<const scalar_value_t<T>,
1430 md::extents<std::size_t, md::dynamic_extent, 3>>;
1431
1432 std::shared_ptr<const mesh::Mesh<U>> mesh = L.mesh();
1433 assert(mesh);
1434 auto x = mesh->geometry().x();
1435 if constexpr (std::is_same_v<U, scalar_value_t<T>>)
1436 {
1437 impl::assemble_vector(b, L, mdspanx3_t(x.data(), x.size() / 3, 3),
1438 constants, coefficients);
1439 }
1440 else
1441 {
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);
1445 }
1446}
1447} // namespace dolfinx::fem::impl
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