Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d9/d40/assemble__vector__impl_8h_source.html
DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
assemble_vector_impl.h
1// Copyright (C) 2018-2021 Garth N. Wells
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 "FunctionSpace.h"
14#include "traits.h"
15#include "utils.h"
16#include <algorithm>
17#include <basix/mdspan.hpp>
18#include <cstdint>
19#include <dolfinx/common/IndexMap.h>
20#include <dolfinx/mesh/Geometry.h>
21#include <dolfinx/mesh/Mesh.h>
22#include <dolfinx/mesh/Topology.h>
23#include <functional>
24#include <memory>
25#include <span>
26#include <vector>
27
28namespace dolfinx::fem::impl
29{
30
32using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
33 const std::int32_t,
34 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
36
74template <dolfinx::scalar T, int _bs0 = -1, int _bs1 = -1>
75void _lift_bc_cells(
76 std::span<T> b, mdspan2_t x_dofmap,
77 std::span<const scalar_value_type_t<T>> x, FEkernel<T> auto kernel,
78 std::span<const std::int32_t> cells,
79 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap0,
80 fem::DofTransformKernel<T> auto P0,
81 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap1,
82 fem::DofTransformKernel<T> auto P1T, std::span<const T> constants,
83 std::span<const T> coeffs, int cstride,
84 std::span<const std::uint32_t> cell_info0,
85 std::span<const std::uint32_t> cell_info1, std::span<const T> bc_values1,
86 std::span<const std::int8_t> bc_markers1, std::span<const T> x0, T scale)
87{
88 if (cells.empty())
89 return;
90
91 const auto [dmap0, bs0, cells0] = dofmap0;
92 const auto [dmap1, bs1, cells1] = dofmap1;
93 assert(_bs0 < 0 or _bs0 == bs0);
94 assert(_bs1 < 0 or _bs1 == bs1);
95
96 // Data structures used in bc application
97 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
98 std::vector<T> Ae, be;
99 assert(cells0.size() == cells.size());
100 assert(cells1.size() == cells.size());
101 for (std::size_t index = 0; index < cells.size(); ++index)
102 {
103 // Cell index in integration domain mesh, test function mesh, and trial
104 // function mesh
105 std::int32_t c = cells[index];
106 std::int32_t c0 = cells0[index];
107 std::int32_t c1 = cells1[index];
108
109 // Get dof maps for cell
110 auto dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
111 dmap1, c1, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
112
113 // Check if bc is applied to cell
114 bool has_bc = false;
115 for (std::size_t j = 0; j < dofs1.size(); ++j)
116 {
117 if constexpr (_bs1 > 0)
118 {
119 for (int k = 0; k < _bs1; ++k)
120 {
121 assert(_bs1 * dofs1[j] + k < (int)bc_markers1.size());
122 if (bc_markers1[_bs1 * dofs1[j] + k])
123 {
124 has_bc = true;
125 break;
126 }
127 }
128 }
129 else
130 {
131 for (int k = 0; k < bs1; ++k)
132 {
133 assert(bs1 * dofs1[j] + k < (int)bc_markers1.size());
134 if (bc_markers1[bs1 * dofs1[j] + k])
135 {
136 has_bc = true;
137 break;
138 }
139 }
140 }
141 }
142
143 if (!has_bc)
144 continue;
145
146 // Get cell coordinates/geometry
147 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
148 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
149 for (std::size_t i = 0; i < x_dofs.size(); ++i)
150 {
151 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
152 std::next(coordinate_dofs.begin(), 3 * i));
153 }
154
155 // Size data structure for assembly
156 auto dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
157 dmap0, c0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
158
159 const int num_rows = bs0 * dofs0.size();
160 const int num_cols = bs1 * dofs1.size();
161
162 const T* coeff_array = coeffs.data() + index * cstride;
163 Ae.resize(num_rows * num_cols);
164 std::fill(Ae.begin(), Ae.end(), 0);
165 kernel(Ae.data(), coeff_array, constants.data(), coordinate_dofs.data(),
166 nullptr, nullptr);
167 P0(Ae, cell_info0, c0, num_cols);
168 P1T(Ae, cell_info1, c1, num_rows);
169
170 // Size data structure for assembly
171 be.resize(num_rows);
172 std::fill(be.begin(), be.end(), 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) * scale * (bc - _x0);
187 for (int m = 0; m < num_rows; ++m)
188 be[m] -= Ae[m * num_cols + _bs1 * j + k] * scale * (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) * scale * (bc - _x0);
203 for (int m = 0; m < num_rows; ++m)
204 be[m] -= Ae[m * num_cols + bs1 * j + k] * scale * (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
257template <dolfinx::scalar T, int _bs = -1>
258void _lift_bc_exterior_facets(
259 std::span<T> b, mdspan2_t x_dofmap,
260 std::span<const scalar_value_type_t<T>> x, FEkernel<T> auto kernel,
261 std::span<const std::int32_t> facets,
262 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap0,
263 fem::DofTransformKernel<T> auto P0,
264 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap1,
265 fem::DofTransformKernel<T> auto P1T, std::span<const T> constants,
266 std::span<const T> coeffs, int cstride,
267 std::span<const std::uint32_t> cell_info0,
268 std::span<const std::uint32_t> cell_info1, std::span<const T> bc_values1,
269 std::span<const std::int8_t> bc_markers1, std::span<const T> x0, T scale)
270{
271 if (facets.empty())
272 return;
273
274 const auto [dmap0, bs0, facets0] = dofmap0;
275 const auto [dmap1, bs1, facets1] = dofmap1;
276
277 // Data structures used in bc application
278 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
279 std::vector<T> Ae, be;
280 assert(facets.size() % 2 == 0);
281 assert(facets0.size() == facets.size());
282 assert(facets1.size() == facets.size());
283 for (std::size_t index = 0; index < facets.size(); index += 2)
284 {
285 // Cell in integration domain mesh
286 std::int32_t cell = facets[index];
287 // Cell in test function mesh
288 std::int32_t cell0 = facets0[index];
289 // Cell in trial function mesh
290 std::int32_t cell1 = facets1[index];
291
292 // Local facet index
293 std::int32_t local_facet = facets[index + 1];
294
295 // Get dof maps for cell
296 auto dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
297 dmap1, cell1, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
298
299 // Check if bc is applied to cell
300 bool has_bc = false;
301 for (std::size_t j = 0; j < dofs1.size(); ++j)
302 {
303 for (int k = 0; k < bs1; ++k)
304 {
305 if (bc_markers1[bs1 * dofs1[j] + k])
306 {
307 has_bc = true;
308 break;
309 }
310 }
311 }
312
313 if (!has_bc)
314 continue;
315
316 // Get cell coordinates/geometry
317 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
318 x_dofmap, cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
319 for (std::size_t i = 0; i < x_dofs.size(); ++i)
320 {
321 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
322 std::next(coordinate_dofs.begin(), 3 * i));
323 }
324
325 // Size data structure for assembly
326 auto dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
327 dmap0, cell0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
328
329 const int num_rows = bs0 * dofs0.size();
330 const int num_cols = bs1 * dofs1.size();
331
332 const T* coeff_array = coeffs.data() + index / 2 * cstride;
333 Ae.resize(num_rows * num_cols);
334 std::fill(Ae.begin(), Ae.end(), 0);
335 kernel(Ae.data(), coeff_array, constants.data(), coordinate_dofs.data(),
336 &local_facet, nullptr);
337 P0(Ae, cell_info0, cell0, num_cols);
338 P1T(Ae, cell_info1, cell1, num_rows);
339
340 // Size data structure for assembly
341 be.resize(num_rows);
342 std::fill(be.begin(), be.end(), 0);
343 for (std::size_t j = 0; j < dofs1.size(); ++j)
344 {
345 for (int k = 0; k < bs1; ++k)
346 {
347 const std::int32_t jj = bs1 * dofs1[j] + k;
348 if (bc_markers1[jj])
349 {
350 const T bc = bc_values1[jj];
351 const T _x0 = x0.empty() ? 0 : x0[jj];
352 // be -= Ae.col(bs1 * j + k) * scale * (bc - _x0);
353 for (int m = 0; m < num_rows; ++m)
354 be[m] -= Ae[m * num_cols + bs1 * j + k] * scale * (bc - _x0);
355 }
356 }
357 }
358
359 for (std::size_t i = 0; i < dofs0.size(); ++i)
360 for (int k = 0; k < bs0; ++k)
361 b[bs0 * dofs0[i] + k] += be[bs0 * i + k];
362 }
363}
364
398template <dolfinx::scalar T, int _bs = -1>
399void _lift_bc_interior_facets(
400 std::span<T> b, mdspan2_t x_dofmap,
401 std::span<const scalar_value_type_t<T>> x, int num_cell_facets,
402 FEkernel<T> auto kernel, std::span<const std::int32_t> facets,
403 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap0,
404 fem::DofTransformKernel<T> auto P0,
405 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap1,
406 fem::DofTransformKernel<T> auto P1T, std::span<const T> constants,
407 std::span<const T> coeffs, int cstride,
408 std::span<const std::uint32_t> cell_info0,
409 std::span<const std::uint32_t> cell_info1,
410 const std::function<std::uint8_t(std::size_t)>& get_perm,
411 std::span<const T> bc_values1, std::span<const std::int8_t> bc_markers1,
412 std::span<const T> x0, T scale)
413{
414 if (facets.empty())
415 return;
416
417 const auto [dmap0, bs0, facets0] = dofmap0;
418 const auto [dmap1, bs1, facets1] = dofmap1;
419
420 // Data structures used in assembly
421 using X = scalar_value_type_t<T>;
422 std::vector<X> coordinate_dofs(2 * x_dofmap.extent(1) * 3);
423 std::span<X> cdofs0(coordinate_dofs.data(), x_dofmap.extent(1) * 3);
424 std::span<X> cdofs1(coordinate_dofs.data() + x_dofmap.extent(1) * 3,
425 x_dofmap.extent(1) * 3);
426 std::vector<T> Ae, be;
427
428 // Temporaries for joint dofmaps
429 std::vector<std::int32_t> dmapjoint0, dmapjoint1;
430 assert(facets.size() % 4 == 0);
431
432 const int num_dofs0 = dmap0.extent(1);
433 const int num_dofs1 = dmap1.extent(1);
434 assert(facets0.size() == facets.size());
435 assert(facets1.size() == facets.size());
436 for (std::size_t index = 0; index < facets.size(); index += 4)
437 {
438 // Cells in integration domain, test function domain and trial
439 // function domain meshes
440 std::array cells{facets[index], facets[index + 2]};
441 std::array cells0{facets0[index], facets0[index + 2]};
442 std::array cells1{facets1[index], facets1[index + 2]};
443
444 // Local facet indices
445 std::array<std::int32_t, 2> local_facet
446 = {facets[index + 1], facets[index + 3]};
447
448 // Get cell geometry
449 auto x_dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
450 x_dofmap, cells[0], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
451 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
452 {
453 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
454 std::next(cdofs0.begin(), 3 * i));
455 }
456 auto x_dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
457 x_dofmap, cells[1], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
458 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
459 {
460 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
461 std::next(cdofs1.begin(), 3 * i));
462 }
463
464 // Get dof maps for cells and pack
465 auto dmap0_cell0
466 = std::span(dmap0.data_handle() + cells0[0] * num_dofs0, num_dofs0);
467 auto dmap0_cell1
468 = std::span(dmap0.data_handle() + cells0[1] * num_dofs0, num_dofs0);
469
470 dmapjoint0.resize(dmap0_cell0.size() + dmap0_cell1.size());
471 std::copy(dmap0_cell0.begin(), dmap0_cell0.end(), dmapjoint0.begin());
472 std::copy(dmap0_cell1.begin(), dmap0_cell1.end(),
473 std::next(dmapjoint0.begin(), dmap0_cell0.size()));
474
475 auto dmap1_cell0
476 = std::span(dmap1.data_handle() + cells1[0] * num_dofs1, num_dofs1);
477 auto dmap1_cell1
478 = std::span(dmap1.data_handle() + cells1[1] * num_dofs1, num_dofs1);
479
480 dmapjoint1.resize(dmap1_cell0.size() + dmap1_cell1.size());
481 std::copy(dmap1_cell0.begin(), dmap1_cell0.end(), dmapjoint1.begin());
482 std::copy(dmap1_cell1.begin(), dmap1_cell1.end(),
483 std::next(dmapjoint1.begin(), dmap1_cell0.size()));
484
485 // Check if bc is applied to cell0
486 bool has_bc = false;
487 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
488 {
489 for (int k = 0; k < bs1; ++k)
490 {
491 if (bc_markers1[bs1 * dmap1_cell0[j] + k])
492 {
493 has_bc = true;
494 break;
495 }
496 }
497 }
498
499 // Check if bc is applied to cell1
500 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
501 {
502 for (int k = 0; k < bs1; ++k)
503 {
504 if (bc_markers1[bs1 * dmap1_cell1[j] + k])
505 {
506 has_bc = true;
507 break;
508 }
509 }
510 }
511
512 if (!has_bc)
513 continue;
514
515 const int num_rows = bs0 * dmapjoint0.size();
516 const int num_cols = bs1 * dmapjoint1.size();
517
518 // Tabulate tensor
519 Ae.resize(num_rows * num_cols);
520 std::fill(Ae.begin(), Ae.end(), 0);
521 const std::array perm{
522 get_perm(cells[0] * num_cell_facets + local_facet[0]),
523 get_perm(cells[1] * num_cell_facets + local_facet[1])};
524 kernel(Ae.data(), coeffs.data() + index / 2 * cstride, constants.data(),
525 coordinate_dofs.data(), local_facet.data(), perm.data());
526
527 std::span<T> _Ae(Ae);
528 std::span<T> sub_Ae0 = _Ae.subspan(bs0 * dmap0_cell0.size() * num_cols,
529 bs0 * dmap0_cell1.size() * num_cols);
530
531 P0(_Ae, cell_info0, cells0[0], num_cols);
532 P0(sub_Ae0, cell_info0, cells0[1], num_cols);
533 P1T(_Ae, cell_info1, cells1[0], num_rows);
534
535 for (int row = 0; row < num_rows; ++row)
536 {
537 // DOFs for dmap1 and cell1 are not stored contiguously in
538 // the block matrix, so each row needs a separate span access
539 std::span<T> sub_Ae1 = _Ae.subspan(
540 row * num_cols + bs1 * dmap1_cell0.size(), bs1 * dmap1_cell1.size());
541 P1T(sub_Ae1, cell_info1, cells1[1], 1);
542 }
543
544 be.resize(num_rows);
545 std::fill(be.begin(), be.end(), 0);
546
547 // Compute b = b - A*b for cell0
548 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
549 {
550 for (int k = 0; k < bs1; ++k)
551 {
552 const std::int32_t jj = bs1 * dmap1_cell0[j] + k;
553 if (bc_markers1[jj])
554 {
555 const T bc = bc_values1[jj];
556 const T _x0 = x0.empty() ? 0 : x0[jj];
557 // be -= Ae.col(bs1 * j + k) * scale * (bc - _x0);
558 for (int m = 0; m < num_rows; ++m)
559 be[m] -= Ae[m * num_cols + bs1 * j + k] * scale * (bc - _x0);
560 }
561 }
562 }
563
564 // Compute b = b - A*b for cell1
565 const int offset = bs1 * dmap1_cell0.size();
566 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
567 {
568 for (int k = 0; k < bs1; ++k)
569 {
570 const std::int32_t jj = bs1 * dmap1_cell1[j] + k;
571 if (bc_markers1[jj])
572 {
573 const T bc = bc_values1[jj];
574 const T _x0 = x0.empty() ? 0 : x0[jj];
575 // be -= Ae.col(offset + bs1 * j + k) * scale * (bc - x0[jj]);
576 for (int m = 0; m < num_rows; ++m)
577 {
578 be[m]
579 -= Ae[m * num_cols + offset + bs1 * j + k] * scale * (bc - _x0);
580 }
581 }
582 }
583 }
584
585 for (std::size_t i = 0; i < dmap0_cell0.size(); ++i)
586 for (int k = 0; k < bs0; ++k)
587 b[bs0 * dmap0_cell0[i] + k] += be[bs0 * i + k];
588
589 const int offset_be = bs0 * dmap0_cell0.size();
590 for (std::size_t i = 0; i < dmap0_cell1.size(); ++i)
591 for (int k = 0; k < bs0; ++k)
592 b[bs0 * dmap0_cell1[i] + k] += be[offset_be + bs0 * i + k];
593 }
594}
595
618template <dolfinx::scalar T, int _bs = -1>
619void assemble_cells(
620 fem::DofTransformKernel<T> auto P0, std::span<T> b, mdspan2_t x_dofmap,
621 std::span<const scalar_value_type_t<T>> x,
622 std::span<const std::int32_t> cells,
623 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap,
624 FEkernel<T> auto kernel, std::span<const T> constants,
625 std::span<const T> coeffs, int cstride,
626 std::span<const std::uint32_t> cell_info0)
627{
628 if (cells.empty())
629 return;
630
631 const auto [dmap, bs, cells0] = dofmap;
632 assert(_bs < 0 or _bs == bs);
633
634 // Create data structures used in assembly
635 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
636 std::vector<T> be(bs * dmap.extent(1));
637 std::span<T> _be(be);
638
639 // Iterate over active cells
640 for (std::size_t index = 0; index < cells.size(); ++index)
641 {
642 // Integration domain celland test function cell
643 std::int32_t c = cells[index];
644 std::int32_t c0 = cells0[index];
645
646 // Get cell coordinates/geometry
647 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
648 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
649 for (std::size_t i = 0; i < x_dofs.size(); ++i)
650 {
651 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
652 std::next(coordinate_dofs.begin(), 3 * i));
653 }
654
655 // Tabulate vector for cell
656 std::fill(be.begin(), be.end(), 0);
657 kernel(be.data(), coeffs.data() + index * cstride, constants.data(),
658 coordinate_dofs.data(), nullptr, nullptr);
659 P0(_be, cell_info0, c0, 1);
660
661 // Scatter cell vector to 'global' vector array
662 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
663 dmap, c0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
664 if constexpr (_bs > 0)
665 {
666 for (std::size_t i = 0; i < dofs.size(); ++i)
667 for (int k = 0; k < _bs; ++k)
668 b[_bs * dofs[i] + k] += be[_bs * i + k];
669 }
670 else
671 {
672 for (std::size_t i = 0; i < dofs.size(); ++i)
673 for (int k = 0; k < bs; ++k)
674 b[bs * dofs[i] + k] += be[bs * i + k];
675 }
676 }
677}
678
701template <dolfinx::scalar T, int _bs = -1>
702void assemble_exterior_facets(
703 fem::DofTransformKernel<T> auto P0, std::span<T> b, mdspan2_t x_dofmap,
704 std::span<const scalar_value_type_t<T>> x,
705 std::span<const std::int32_t> facets,
706 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap,
707 FEkernel<T> auto fn, std::span<const T> constants,
708 std::span<const T> coeffs, int cstride,
709 std::span<const std::uint32_t> cell_info0)
710{
711 if (facets.empty())
712 return;
713
714 const auto [dmap, bs, facets0] = dofmap;
715 assert(_bs < 0 or _bs == bs);
716
717 // FIXME: Add proper interface for num_dofs
718 // Create data structures used in assembly
719 const int num_dofs = dmap.extent(1);
720 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
721 std::vector<T> be(bs * num_dofs);
722 std::span<T> _be(be);
723 assert(facets.size() % 2 == 0);
724 assert(facets0.size() == facets.size());
725 for (std::size_t index = 0; index < facets.size(); index += 2)
726 {
727 // Cell in the integration domain, local facet index relative to the
728 // integration domain cell, and cell in the test function mesh
729 std::int32_t cell = facets[index];
730 std::int32_t local_facet = facets[index + 1];
731 std::int32_t cell0 = facets0[index];
732
733 // Get cell coordinates/geometry
734 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
735 x_dofmap, cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
736 for (std::size_t i = 0; i < x_dofs.size(); ++i)
737 {
738 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
739 std::next(coordinate_dofs.begin(), 3 * i));
740 }
741
742 // Tabulate element vector
743 std::fill(be.begin(), be.end(), 0);
744 fn(be.data(), coeffs.data() + index / 2 * cstride, constants.data(),
745 coordinate_dofs.data(), &local_facet, nullptr);
746
747 P0(_be, cell_info0, cell0, 1);
748
749 // Add element vector to global vector
750 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
751 dmap, cell0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
752 if constexpr (_bs > 0)
753 {
754 for (std::size_t i = 0; i < dofs.size(); ++i)
755 for (int k = 0; k < _bs; ++k)
756 b[_bs * dofs[i] + k] += be[_bs * i + k];
757 }
758 else
759 {
760 for (std::size_t i = 0; i < dofs.size(); ++i)
761 for (int k = 0; k < bs; ++k)
762 b[bs * dofs[i] + k] += be[bs * i + k];
763 }
764 }
765}
766
791template <dolfinx::scalar T, int _bs = -1>
792void assemble_interior_facets(
793 fem::DofTransformKernel<T> auto P0, std::span<T> b, mdspan2_t x_dofmap,
794 std::span<const scalar_value_type_t<T>> x, int num_cell_facets,
795 std::span<const std::int32_t> facets,
796 std::tuple<const DofMap&, int, std::span<const std::int32_t>> dofmap,
797 FEkernel<T> auto fn, std::span<const T> constants,
798 std::span<const T> coeffs, int cstride,
799 std::span<const std::uint32_t> cell_info0,
800 const std::function<std::uint8_t(std::size_t)>& get_perm)
801{
802 if (facets.empty())
803 return;
804
805 const auto [dmap, bs, facets0] = dofmap;
806 assert(_bs < 0 or _bs == bs);
807
808 // Create data structures used in assembly
809 using X = scalar_value_type_t<T>;
810 std::vector<X> coordinate_dofs(2 * x_dofmap.extent(1) * 3);
811 std::span<X> cdofs0(coordinate_dofs.data(), x_dofmap.extent(1) * 3);
812 std::span<X> cdofs1(coordinate_dofs.data() + x_dofmap.extent(1) * 3,
813 x_dofmap.extent(1) * 3);
814 std::vector<T> be;
815
816 assert(facets.size() % 4 == 0);
817 assert(facets0.size() == facets.size());
818 for (std::size_t index = 0; index < facets.size(); index += 4)
819 {
820 // Cells in integration domain and test function domain meshes
821 std::array<std::int32_t, 2> cells{facets[index], facets[index + 2]};
822 std::array<std::int32_t, 2> cells0{facets0[index], facets0[index + 2]};
823
824 // Local facet indices
825 std::array<std::int32_t, 2> local_facet{facets[index + 1],
826 facets[index + 3]};
827
828 // Get cell geometry
829 auto x_dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
830 x_dofmap, cells[0], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
831 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
832 {
833 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
834 std::next(cdofs0.begin(), 3 * i));
835 }
836 auto x_dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
837 x_dofmap, cells[1], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
838 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
839 {
840 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
841 std::next(cdofs1.begin(), 3 * i));
842 }
843
844 // Get dofmaps for cells
845 std::span<const std::int32_t> dmap0 = dmap.cell_dofs(cells0[0]);
846 std::span<const std::int32_t> dmap1 = dmap.cell_dofs(cells0[1]);
847
848 // Tabulate element vector
849 be.resize(bs * (dmap0.size() + dmap1.size()));
850 std::fill(be.begin(), be.end(), 0);
851 const std::array perm{
852 get_perm(cells[0] * num_cell_facets + local_facet[0]),
853 get_perm(cells[1] * num_cell_facets + local_facet[1])};
854 fn(be.data(), coeffs.data() + index / 2 * cstride, constants.data(),
855 coordinate_dofs.data(), local_facet.data(), perm.data());
856
857 std::span<T> _be(be);
858 std::span<T> sub_be = _be.subspan(bs * dmap0.size(), bs * dmap1.size());
859
860 P0(be, cell_info0, cells0[0], 1);
861 P0(sub_be, cell_info0, cells0[1], 1);
862
863 // Add element vector to global vector
864 if constexpr (_bs > 0)
865 {
866 for (std::size_t i = 0; i < dmap0.size(); ++i)
867 for (int k = 0; k < _bs; ++k)
868 b[_bs * dmap0[i] + k] += be[_bs * i + k];
869 for (std::size_t i = 0; i < dmap1.size(); ++i)
870 for (int k = 0; k < _bs; ++k)
871 b[_bs * dmap1[i] + k] += be[_bs * (i + dmap0.size()) + k];
872 }
873 else
874 {
875 for (std::size_t i = 0; i < dmap0.size(); ++i)
876 for (int k = 0; k < bs; ++k)
877 b[bs * dmap0[i] + k] += be[bs * i + k];
878 for (std::size_t i = 0; i < dmap1.size(); ++i)
879 for (int k = 0; k < bs; ++k)
880 b[bs * dmap1[i] + k] += be[bs * (i + dmap0.size()) + k];
881 }
882 }
883}
884
901template <dolfinx::scalar T, std::floating_point U>
902void lift_bc(std::span<T> b, const Form<T, U>& a, mdspan2_t x_dofmap,
903 std::span<const scalar_value_type_t<T>> x,
904 std::span<const T> constants,
905 const std::map<std::pair<IntegralType, int>,
906 std::pair<std::span<const T>, int>>& coefficients,
907 std::span<const T> bc_values1,
908 std::span<const std::int8_t> bc_markers1, std::span<const T> x0,
909 T scale)
910{
911 // Integration domain mesh
912 std::shared_ptr<const mesh::Mesh<U>> mesh = a.mesh();
913 assert(mesh);
914 // Test function mesh
915 auto mesh0 = a.function_spaces().at(0)->mesh();
916 assert(mesh0);
917 // Trial function mesh
918 auto mesh1 = a.function_spaces().at(1)->mesh();
919 assert(mesh1);
920
921 // Get dofmap for columns and rows of a
922 assert(a.function_spaces().at(0));
923 assert(a.function_spaces().at(1));
924 auto dofmap0 = a.function_spaces()[0]->dofmap()->map();
925 const int bs0 = a.function_spaces()[0]->dofmap()->bs();
926 auto element0 = a.function_spaces()[0]->element();
927 auto dofmap1 = a.function_spaces()[1]->dofmap()->map();
928 const int bs1 = a.function_spaces()[1]->dofmap()->bs();
929 auto element1 = a.function_spaces()[1]->element();
930 assert(element0);
931
932 std::span<const std::uint32_t> cell_info0;
933 std::span<const std::uint32_t> cell_info1;
934 // TODO: Check for each element instead
935 if (element0->needs_dof_transformations()
936 or element1->needs_dof_transformations() or a.needs_facet_permutations())
937 {
938 mesh0->topology_mutable()->create_entity_permutations();
939 mesh1->topology_mutable()->create_entity_permutations();
940 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
941 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
942 }
943
944 fem::DofTransformKernel<T> auto P0
945 = element0->template dof_transformation_fn<T>(doftransform::standard);
946 fem::DofTransformKernel<T> auto P1T
947 = element1->template dof_transformation_right_fn<T>(
949
950 for (int i : a.integral_ids(IntegralType::cell))
951 {
952 auto kernel = a.kernel(IntegralType::cell, i);
953 assert(kernel);
954 auto& [coeffs, cstride] = coefficients.at({IntegralType::cell, i});
955 std::span<const std::int32_t> cells = a.domain(IntegralType::cell, i);
956 if (bs0 == 1 and bs1 == 1)
957 {
958 _lift_bc_cells<T, 1, 1>(
959 b, x_dofmap, x, kernel, cells,
960 {dofmap0, bs0, a.domain(IntegralType::cell, i, *mesh0)}, P0,
961 {dofmap1, bs1, a.domain(IntegralType::cell, i, *mesh1)}, P1T,
962 constants, coeffs, cstride, cell_info0, cell_info1, bc_values1,
963 bc_markers1, x0, scale);
964 }
965 else if (bs0 == 3 and bs1 == 3)
966 {
967 _lift_bc_cells<T, 3, 3>(
968 b, x_dofmap, x, kernel, cells,
969 {dofmap0, bs0, a.domain(IntegralType::cell, i, *mesh0)}, P0,
970 {dofmap1, bs1, a.domain(IntegralType::cell, i, *mesh1)}, P1T,
971 constants, coeffs, cstride, cell_info0, cell_info1, bc_values1,
972 bc_markers1, x0, scale);
973 }
974 else
975 {
976 _lift_bc_cells(b, x_dofmap, x, kernel, cells,
977 {dofmap0, bs0, a.domain(IntegralType::cell, i, *mesh0)},
978 P0,
979 {dofmap1, bs1, a.domain(IntegralType::cell, i, *mesh1)},
980 P1T, constants, coeffs, cstride, cell_info0, cell_info1,
981 bc_values1, bc_markers1, x0, scale);
982 }
983 }
984
985 for (int i : a.integral_ids(IntegralType::exterior_facet))
986 {
987 auto kernel = a.kernel(IntegralType::exterior_facet, i);
988 assert(kernel);
989 auto& [coeffs, cstride]
990 = coefficients.at({IntegralType::exterior_facet, i});
991 _lift_bc_exterior_facets(
992 b, x_dofmap, x, kernel, a.domain(IntegralType::exterior_facet, i),
993 {dofmap0, bs0, a.domain(IntegralType::exterior_facet, i, *mesh0)}, P0,
994 {dofmap1, bs1, a.domain(IntegralType::exterior_facet, i, *mesh1)}, P1T,
995 constants, coeffs, cstride, cell_info0, cell_info1, bc_values1,
996 bc_markers1, x0, scale);
997 }
998
999 if (a.num_integrals(IntegralType::interior_facet) > 0)
1000 {
1001 std::function<std::uint8_t(std::size_t)> get_perm;
1002 if (a.needs_facet_permutations())
1003 {
1004 mesh->topology_mutable()->create_entity_permutations();
1005 const std::vector<std::uint8_t>& perms
1006 = mesh->topology()->get_facet_permutations();
1007 get_perm = [&perms](std::size_t i) { return perms[i]; };
1008 }
1009 else
1010 get_perm = [](std::size_t) { return 0; };
1011
1012 mesh::CellType cell_type = mesh->topology()->cell_type();
1013 int num_cell_facets
1014 = mesh::cell_num_entities(cell_type, mesh->topology()->dim() - 1);
1015 for (int i : a.integral_ids(IntegralType::interior_facet))
1016 {
1017 auto kernel = a.kernel(IntegralType::interior_facet, i);
1018 assert(kernel);
1019 auto& [coeffs, cstride]
1020 = coefficients.at({IntegralType::interior_facet, i});
1021 _lift_bc_interior_facets(
1022 b, x_dofmap, x, num_cell_facets, kernel,
1023 a.domain(IntegralType::interior_facet, i),
1024 {dofmap0, bs0, a.domain(IntegralType::interior_facet, i, *mesh0)}, P0,
1025 {dofmap1, bs1, a.domain(IntegralType::interior_facet, i, *mesh1)},
1026 P1T, constants, coeffs, cstride, cell_info0, cell_info1, get_perm,
1027 bc_values1, bc_markers1, x0, scale);
1028 }
1029 }
1030}
1031
1053template <dolfinx::scalar T, std::floating_point U>
1054void apply_lifting(
1055 std::span<T> b, const std::vector<std::shared_ptr<const Form<T, U>>> a,
1056 mdspan2_t x_dofmap, std::span<const scalar_value_type_t<T>> x,
1057 const std::vector<std::span<const T>>& constants,
1058 const std::vector<std::map<std::pair<IntegralType, int>,
1059 std::pair<std::span<const T>, int>>>& coeffs,
1060 const std::vector<std::vector<std::shared_ptr<const DirichletBC<T, U>>>>&
1061 bcs1,
1062 const std::vector<std::span<const T>>& x0, T scale)
1063{
1064 // FIXME: make changes to reactivate this check
1065 if (!x0.empty() and x0.size() != a.size())
1066 {
1067 throw std::runtime_error(
1068 "Mismatch in size between x0 and bilinear form in assembler.");
1069 }
1070
1071 if (a.size() != bcs1.size())
1072 {
1073 throw std::runtime_error(
1074 "Mismatch in size between a and bcs in assembler.");
1075 }
1076
1077 for (std::size_t j = 0; j < a.size(); ++j)
1078 {
1079 std::vector<std::int8_t> bc_markers1;
1080 std::vector<T> bc_values1;
1081 if (a[j] and !bcs1[j].empty())
1082 {
1083 assert(a[j]->function_spaces().at(0));
1084
1085 auto V1 = a[j]->function_spaces()[1];
1086 assert(V1);
1087 auto map1 = V1->dofmap()->index_map;
1088 const int bs1 = V1->dofmap()->index_map_bs();
1089 assert(map1);
1090 const int crange = bs1 * (map1->size_local() + map1->num_ghosts());
1091 bc_markers1.assign(crange, false);
1092 bc_values1.assign(crange, 0);
1093 for (const std::shared_ptr<const DirichletBC<T, U>>& bc : bcs1[j])
1094 {
1095 bc->mark_dofs(bc_markers1);
1096 bc->dof_values(bc_values1);
1097 }
1098
1099 if (!x0.empty())
1100 {
1101 lift_bc<T>(b, *a[j], x_dofmap, x, constants[j], coeffs[j], bc_values1,
1102 bc_markers1, x0[j], scale);
1103 }
1104 else
1105 {
1106 lift_bc<T>(b, *a[j], x_dofmap, x, constants[j], coeffs[j], bc_values1,
1107 bc_markers1, std::span<const T>(), scale);
1108 }
1109 }
1110 }
1111}
1112
1121template <dolfinx::scalar T, std::floating_point U>
1122void assemble_vector(
1123 std::span<T> b, const Form<T, U>& L, mdspan2_t x_dofmap,
1124 std::span<const scalar_value_type_t<T>> x, std::span<const T> constants,
1125 const std::map<std::pair<IntegralType, int>,
1126 std::pair<std::span<const T>, int>>& coefficients)
1127{
1128 // Integration domain mesh
1129 std::shared_ptr<const mesh::Mesh<U>> mesh = L.mesh();
1130 assert(mesh);
1131
1132 // Test function mesh
1133 auto mesh0 = L.function_spaces().at(0)->mesh();
1134 assert(mesh0);
1135
1136 // Get dofmap data
1137 assert(L.function_spaces().at(0));
1138 auto element = L.function_spaces().at(0)->element();
1139 assert(element);
1140 std::shared_ptr<const fem::DofMap> dofmap
1141 = L.function_spaces().at(0)->dofmap();
1142 assert(dofmap);
1143 auto dofs = dofmap->map();
1144 const int bs = dofmap->bs();
1145
1146 fem::DofTransformKernel<T> auto P0
1147 = element->template dof_transformation_fn<T>(doftransform::standard);
1148
1149 std::span<const std::uint32_t> cell_info0;
1150 if (element->needs_dof_transformations() or L.needs_facet_permutations())
1151 {
1152 mesh0->topology_mutable()->create_entity_permutations();
1153 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
1154 }
1155
1156 for (int i : L.integral_ids(IntegralType::cell))
1157 {
1158 auto fn = L.kernel(IntegralType::cell, i);
1159 assert(fn);
1160 auto& [coeffs, cstride] = coefficients.at({IntegralType::cell, i});
1161 std::span<const std::int32_t> cells = L.domain(IntegralType::cell, i);
1162 if (bs == 1)
1163 {
1164 impl::assemble_cells<T, 1>(
1165 P0, b, x_dofmap, x, cells,
1166 {dofs, bs, L.domain(IntegralType::cell, i, *mesh0)}, fn, constants,
1167 coeffs, cstride, cell_info0);
1168 }
1169 else if (bs == 3)
1170 {
1171 impl::assemble_cells<T, 3>(
1172 P0, b, x_dofmap, x, cells,
1173 {dofs, bs, L.domain(IntegralType::cell, i, *mesh0)}, fn, constants,
1174 coeffs, cstride, cell_info0);
1175 }
1176 else
1177 {
1178 impl::assemble_cells(P0, b, x_dofmap, x, cells,
1179 {dofs, bs, L.domain(IntegralType::cell, i, *mesh0)},
1180 fn, constants, coeffs, cstride, cell_info0);
1181 }
1182 }
1183
1184 for (int i : L.integral_ids(IntegralType::exterior_facet))
1185 {
1186 auto fn = L.kernel(IntegralType::exterior_facet, i);
1187 assert(fn);
1188 auto& [coeffs, cstride]
1189 = coefficients.at({IntegralType::exterior_facet, i});
1190 std::span<const std::int32_t> facets
1191 = L.domain(IntegralType::exterior_facet, i);
1192 if (bs == 1)
1193 {
1194 impl::assemble_exterior_facets<T, 1>(
1195 P0, b, x_dofmap, x, facets,
1196 {dofs, bs, L.domain(IntegralType::exterior_facet, i, *mesh0)}, fn,
1197 constants, coeffs, cstride, cell_info0);
1198 }
1199 else if (bs == 3)
1200 {
1201 impl::assemble_exterior_facets<T, 3>(
1202 P0, b, x_dofmap, x, facets,
1203 {dofs, bs, L.domain(IntegralType::exterior_facet, i, *mesh0)}, fn,
1204 constants, coeffs, cstride, cell_info0);
1205 }
1206 else
1207 {
1208 impl::assemble_exterior_facets(
1209 P0, b, x_dofmap, x, facets,
1210 {dofs, bs, L.domain(IntegralType::exterior_facet, i, *mesh0)}, fn,
1211 constants, coeffs, cstride, cell_info0);
1212 }
1213 }
1214
1215 if (L.num_integrals(IntegralType::interior_facet) > 0)
1216 {
1217 std::function<std::uint8_t(std::size_t)> get_perm;
1218 if (L.needs_facet_permutations())
1219 {
1220 mesh->topology_mutable()->create_entity_permutations();
1221 const std::vector<std::uint8_t>& perms
1222 = mesh->topology()->get_facet_permutations();
1223 get_perm = [&perms](std::size_t i) { return perms[i]; };
1224 }
1225 else
1226 get_perm = [](std::size_t) { return 0; };
1227
1228 mesh::CellType cell_type = mesh->topology()->cell_type();
1229 int num_cell_facets
1230 = mesh::cell_num_entities(cell_type, mesh->topology()->dim() - 1);
1231 for (int i : L.integral_ids(IntegralType::interior_facet))
1232 {
1233 auto fn = L.kernel(IntegralType::interior_facet, i);
1234 assert(fn);
1235 auto& [coeffs, cstride]
1236 = coefficients.at({IntegralType::interior_facet, i});
1237 std::span<const std::int32_t> facets
1238 = L.domain(IntegralType::interior_facet, i);
1239 if (bs == 1)
1240 {
1241 impl::assemble_interior_facets<T, 1>(
1242 P0, b, x_dofmap, x, num_cell_facets, facets,
1243 {*dofmap, bs, L.domain(IntegralType::interior_facet, i, *mesh0)},
1244 fn, constants, coeffs, cstride, cell_info0, get_perm);
1245 }
1246 else if (bs == 3)
1247 {
1248 impl::assemble_interior_facets<T, 3>(
1249 P0, b, x_dofmap, x, num_cell_facets, facets,
1250 {*dofmap, bs, L.domain(IntegralType::interior_facet, i, *mesh0)},
1251 fn, constants, coeffs, cstride, cell_info0, get_perm);
1252 }
1253 else
1254 {
1255 impl::assemble_interior_facets(
1256 P0, b, x_dofmap, x, num_cell_facets, facets,
1257 {*dofmap, bs, L.domain(IntegralType::interior_facet, i, *mesh0)},
1258 fn, constants, coeffs, cstride, cell_info0, get_perm);
1259 }
1260 }
1261 }
1262}
1263
1270template <dolfinx::scalar T, std::floating_point U>
1271void assemble_vector(
1272 std::span<T> b, const Form<T, U>& L, std::span<const T> constants,
1273 const std::map<std::pair<IntegralType, int>,
1274 std::pair<std::span<const T>, int>>& coefficients)
1275{
1276 std::shared_ptr<const mesh::Mesh<U>> mesh = L.mesh();
1277 assert(mesh);
1278 if constexpr (std::is_same_v<U, scalar_value_type_t<T>>)
1279 assemble_vector(b, L, mesh->geometry().dofmap(), mesh->geometry().x(),
1280 constants, coefficients);
1281 else
1282 {
1283 auto x = mesh->geometry().x();
1284 std::vector<scalar_value_type_t<T>> _x(x.begin(), x.end());
1285 assemble_vector(b, L, mesh->geometry().dofmap(), _x, constants,
1286 coefficients);
1287 }
1288}
1289} // namespace dolfinx::fem::impl
Degree-of-freedom map representations and tools.
Definition types.h:20
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:15
IntegralType
Type of integral.
Definition Form.h:34
@ interior_facet
Interior facet.
@ exterior_facet
Exterior facet.
int cell_num_entities(CellType type, int dim)
Definition cell_types.cpp:139
CellType
Cell type identifier.
Definition cell_types.h:22