DOLFINx 0.9.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
assemble_matrix_impl.h
1// Copyright (C) 2018-2019 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 "DofMap.h"
10#include "Form.h"
11#include "FunctionSpace.h"
12#include "traits.h"
13#include "utils.h"
14#include <algorithm>
15#include <dolfinx/la/utils.h>
16#include <dolfinx/mesh/Geometry.h>
17#include <dolfinx/mesh/Mesh.h>
18#include <dolfinx/mesh/Topology.h>
19#include <functional>
20#include <iterator>
21#include <span>
22#include <tuple>
23#include <vector>
24
25namespace dolfinx::fem::impl
26{
27
28using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
29 const std::int32_t,
30 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
31
60template <dolfinx::scalar T>
61void assemble_cells(
62 la::MatSet<T> auto mat_set, mdspan2_t x_dofmap,
63 std::span<const scalar_value_type_t<T>> x,
64 std::span<const std::int32_t> cells,
65 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap0,
67 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap1,
68 fem::DofTransformKernel<T> auto P1T, std::span<const std::int8_t> bc0,
69 std::span<const std::int8_t> bc1, FEkernel<T> auto kernel,
70 std::span<const T> coeffs, int cstride, std::span<const T> constants,
71 std::span<const std::uint32_t> cell_info0,
72 std::span<const std::uint32_t> cell_info1)
73{
74 if (cells.empty())
75 return;
76
77 const auto [dmap0, bs0, cells0] = dofmap0;
78 const auto [dmap1, bs1, cells1] = dofmap1;
79
80 // Iterate over active cells
81 const int num_dofs0 = dmap0.extent(1);
82 const int num_dofs1 = dmap1.extent(1);
83 const int ndim0 = bs0 * num_dofs0;
84 const int ndim1 = bs1 * num_dofs1;
85 std::vector<T> Ae(ndim0 * ndim1);
86 std::span<T> _Ae(Ae);
87 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
88
89 // Iterate over active cells
90 assert(cells0.size() == cells.size());
91 assert(cells1.size() == cells.size());
92 for (std::size_t index = 0; index < cells.size(); ++index)
93 {
94 // Cell index in integration domain mesh (c), test function mesh
95 // (c0) and trial function mesh (c1)
96 std::int32_t c = cells[index];
97 std::int32_t c0 = cells0[index];
98 std::int32_t c1 = cells1[index];
99
100 // Get cell coordinates/geometry
101 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
102 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
103 for (std::size_t i = 0; i < x_dofs.size(); ++i)
104 {
105 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
106 std::next(coordinate_dofs.begin(), 3 * i));
107 }
108
109 // Tabulate tensor
110 std::ranges::fill(Ae, 0);
111 kernel(Ae.data(), coeffs.data() + index * cstride, constants.data(),
112 coordinate_dofs.data(), nullptr, nullptr);
113
114 // Compute A = P_0 \tilde{A} P_1^T (dof transformation)
115 P0(_Ae, cell_info0, c0, ndim1); // B = P0 \tilde{A}
116 P1T(_Ae, cell_info1, c1, ndim0); // A = B P1_T
117
118 // Zero rows/columns for essential bcs
119 auto dofs0 = std::span(dmap0.data_handle() + c0 * num_dofs0, num_dofs0);
120 auto dofs1 = std::span(dmap1.data_handle() + c1 * num_dofs1, num_dofs1);
121
122 if (!bc0.empty())
123 {
124 for (int i = 0; i < num_dofs0; ++i)
125 {
126 for (int k = 0; k < bs0; ++k)
127 {
128 if (bc0[bs0 * dofs0[i] + k])
129 {
130 // Zero row bs0 * i + k
131 const int row = bs0 * i + k;
132 std::fill_n(std::next(Ae.begin(), ndim1 * row), ndim1, 0);
133 }
134 }
135 }
136 }
137
138 if (!bc1.empty())
139 {
140 for (int j = 0; j < num_dofs1; ++j)
141 {
142 for (int k = 0; k < bs1; ++k)
143 {
144 if (bc1[bs1 * dofs1[j] + k])
145 {
146 // Zero column bs1 * j + k
147 const int col = bs1 * j + k;
148 for (int row = 0; row < ndim0; ++row)
149 Ae[row * ndim1 + col] = 0;
150 }
151 }
152 }
153 }
154
155 mat_set(dofs0, dofs1, Ae);
156 }
157}
158
194template <dolfinx::scalar T>
195void assemble_exterior_facets(
196 la::MatSet<T> auto mat_set, mdspan2_t x_dofmap,
197 std::span<const scalar_value_type_t<T>> x, int num_facets_per_cell,
198 std::span<const std::int32_t> facets,
199 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap0,
201 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap1,
202 fem::DofTransformKernel<T> auto P1T, std::span<const std::int8_t> bc0,
203 std::span<const std::int8_t> bc1, FEkernel<T> auto kernel,
204 std::span<const T> coeffs, int cstride, std::span<const T> constants,
205 std::span<const std::uint32_t> cell_info0,
206 std::span<const std::uint32_t> cell_info1,
207 std::span<const std::uint8_t> perms)
208{
209 if (facets.empty())
210 return;
211
212 const auto [dmap0, bs0, facets0] = dofmap0;
213 const auto [dmap1, bs1, facets1] = dofmap1;
214
215 // Data structures used in assembly
216 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
217 const int num_dofs0 = dmap0.extent(1);
218 const int num_dofs1 = dmap1.extent(1);
219 const int ndim0 = bs0 * num_dofs0;
220 const int ndim1 = bs1 * num_dofs1;
221 std::vector<T> Ae(ndim0 * ndim1);
222 std::span<T> _Ae(Ae);
223 assert(facets.size() % 2 == 0);
224 assert(facets0.size() == facets.size());
225 assert(facets1.size() == facets.size());
226 for (std::size_t index = 0; index < facets.size(); index += 2)
227 {
228 // Cell in the integration domain, local facet index relative to the
229 // integration domain cell, and cells in the test and trial function
230 // meshes
231 std::int32_t cell = facets[index];
232 std::int32_t local_facet = facets[index + 1];
233 std::int32_t cell0 = facets0[index];
234 std::int32_t cell1 = facets1[index];
235
236 // Get cell coordinates/geometry
237 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
238 x_dofmap, cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
239 for (std::size_t i = 0; i < x_dofs.size(); ++i)
240 {
241 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
242 std::next(coordinate_dofs.begin(), 3 * i));
243 }
244
245 // Permutations
246 std::uint8_t perm
247 = perms.empty() ? 0 : perms[cell * num_facets_per_cell + local_facet];
248
249 // Tabulate tensor
250 std::ranges::fill(Ae, 0);
251 kernel(Ae.data(), coeffs.data() + index / 2 * cstride, constants.data(),
252 coordinate_dofs.data(), &local_facet, &perm);
253
254 P0(_Ae, cell_info0, cell0, ndim1);
255 P1T(_Ae, cell_info1, cell1, ndim0);
256
257 // Zero rows/columns for essential bcs
258 auto dofs0 = std::span(dmap0.data_handle() + cell0 * num_dofs0, num_dofs0);
259 auto dofs1 = std::span(dmap1.data_handle() + cell1 * num_dofs1, num_dofs1);
260 if (!bc0.empty())
261 {
262 for (int i = 0; i < num_dofs0; ++i)
263 {
264 for (int k = 0; k < bs0; ++k)
265 {
266 if (bc0[bs0 * dofs0[i] + k])
267 {
268 // Zero row bs0 * i + k
269 const int row = bs0 * i + k;
270 std::fill_n(std::next(Ae.begin(), ndim1 * row), ndim1, 0);
271 }
272 }
273 }
274 }
275 if (!bc1.empty())
276 {
277 for (int j = 0; j < num_dofs1; ++j)
278 {
279 for (int k = 0; k < bs1; ++k)
280 {
281 if (bc1[bs1 * dofs1[j] + k])
282 {
283 // Zero column bs1 * j + k
284 const int col = bs1 * j + k;
285 for (int row = 0; row < ndim0; ++row)
286 Ae[row * ndim1 + col] = 0;
287 }
288 }
289 }
290 }
291
292 mat_set(dofs0, dofs1, Ae);
293 }
294}
295
333template <dolfinx::scalar T>
334void assemble_interior_facets(
335 la::MatSet<T> auto mat_set, mdspan2_t x_dofmap,
336 std::span<const scalar_value_type_t<T>> x, int num_facets_per_cell,
337 std::span<const std::int32_t> facets,
338 std::tuple<const DofMap&, int, std::span<const std::int32_t>> dofmap0,
340 std::tuple<const DofMap&, int, std::span<const std::int32_t>> dofmap1,
341 fem::DofTransformKernel<T> auto P1T, std::span<const std::int8_t> bc0,
342 std::span<const std::int8_t> bc1, FEkernel<T> auto kernel,
343 std::span<const T> coeffs, int cstride, std::span<const int> offsets,
344 std::span<const T> constants, std::span<const std::uint32_t> cell_info0,
345 std::span<const std::uint32_t> cell_info1,
346 std::span<const std::uint8_t> perms)
347{
348 if (facets.empty())
349 return;
350
351 const auto [dmap0, bs0, facets0] = dofmap0;
352 const auto [dmap1, bs1, facets1] = dofmap1;
353
354 // Data structures used in assembly
355 using X = scalar_value_type_t<T>;
356 std::vector<X> coordinate_dofs(2 * x_dofmap.extent(1) * 3);
357 std::span<X> cdofs0(coordinate_dofs.data(), x_dofmap.extent(1) * 3);
358 std::span<X> cdofs1(coordinate_dofs.data() + x_dofmap.extent(1) * 3,
359 x_dofmap.extent(1) * 3);
360
361 std::vector<T> Ae, be;
362 std::vector<T> coeff_array(2 * offsets.back());
363 assert(offsets.back() == cstride);
364
365 // Temporaries for joint dofmaps
366 std::vector<std::int32_t> dmapjoint0, dmapjoint1;
367 assert(facets.size() % 4 == 0);
368 assert(facets0.size() == facets.size());
369 assert(facets1.size() == facets.size());
370 for (std::size_t index = 0; index < facets.size(); index += 4)
371 {
372 // Cells in integration domain, test function domain and trial
373 // function domain
374 std::array cells{facets[index], facets[index + 2]};
375 std::array cells0{facets0[index], facets0[index + 2]};
376 std::array cells1{facets1[index], facets1[index + 2]};
377
378 // Local facets indices
379 std::array local_facet{facets[index + 1], facets[index + 3]};
380
381 // Get cell geometry
382 auto x_dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
383 x_dofmap, cells[0], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
384 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
385 {
386 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
387 std::next(cdofs0.begin(), 3 * i));
388 }
389 auto x_dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
390 x_dofmap, cells[1], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
391 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
392 {
393 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
394 std::next(cdofs1.begin(), 3 * i));
395 }
396
397 // Get dof maps for cells and pack
398 std::span<const std::int32_t> dmap0_cell0 = dmap0.cell_dofs(cells0[0]);
399 std::span<const std::int32_t> dmap0_cell1 = dmap0.cell_dofs(cells0[1]);
400 dmapjoint0.resize(dmap0_cell0.size() + dmap0_cell1.size());
401 std::ranges::copy(dmap0_cell0, dmapjoint0.begin());
402 std::ranges::copy(dmap0_cell1,
403 std::next(dmapjoint0.begin(), dmap0_cell0.size()));
404
405 std::span<const std::int32_t> dmap1_cell0 = dmap1.cell_dofs(cells1[0]);
406 std::span<const std::int32_t> dmap1_cell1 = dmap1.cell_dofs(cells1[1]);
407 dmapjoint1.resize(dmap1_cell0.size() + dmap1_cell1.size());
408 std::ranges::copy(dmap1_cell0, dmapjoint1.begin());
409 std::ranges::copy(dmap1_cell1,
410 std::next(dmapjoint1.begin(), dmap1_cell0.size()));
411
412 const int num_rows = bs0 * dmapjoint0.size();
413 const int num_cols = bs1 * dmapjoint1.size();
414
415 // Tabulate tensor
416 Ae.resize(num_rows * num_cols);
417 std::ranges::fill(Ae, 0);
418
419 std::array perm
420 = perms.empty()
421 ? std::array<std::uint8_t, 2>{0, 0}
422 : std::array{
423 perms[cells[0] * num_facets_per_cell + local_facet[0]],
424 perms[cells[1] * num_facets_per_cell + local_facet[1]]};
425 kernel(Ae.data(), coeffs.data() + index / 2 * cstride, constants.data(),
426 coordinate_dofs.data(), local_facet.data(), perm.data());
427
428 // Local element layout is a 2x2 block matrix with structure
429 //
430 // cell0cell0 | cell0cell1
431 // cell1cell0 | cell1cell1
432 //
433 // where each block is element tensor of size (dmap0, dmap1).
434
435 std::span<T> _Ae(Ae);
436 std::span<T> sub_Ae0 = _Ae.subspan(bs0 * dmap0_cell0.size() * num_cols,
437 bs0 * dmap0_cell1.size() * num_cols);
438
439 P0(_Ae, cell_info0, cells0[0], num_cols);
440 P0(sub_Ae0, cell_info0, cells0[1], num_cols);
441 P1T(_Ae, cell_info1, cells1[0], num_rows);
442
443 for (int row = 0; row < num_rows; ++row)
444 {
445 // DOFs for dmap1 and cell1 are not stored contiguously in
446 // the block matrix, so each row needs a separate span access
447 std::span<T> sub_Ae1 = _Ae.subspan(
448 row * num_cols + bs1 * dmap1_cell0.size(), bs1 * dmap1_cell1.size());
449 P1T(sub_Ae1, cell_info1, cells1[1], 1);
450 }
451
452 // Zero rows/columns for essential bcs
453 if (!bc0.empty())
454 {
455 for (std::size_t i = 0; i < dmapjoint0.size(); ++i)
456 {
457 for (int k = 0; k < bs0; ++k)
458 {
459 if (bc0[bs0 * dmapjoint0[i] + k])
460 {
461 // Zero row bs0 * i + k
462 std::fill_n(std::next(Ae.begin(), num_cols * (bs0 * i + k)),
463 num_cols, 0);
464 }
465 }
466 }
467 }
468 if (!bc1.empty())
469 {
470 for (std::size_t j = 0; j < dmapjoint1.size(); ++j)
471 {
472 for (int k = 0; k < bs1; ++k)
473 {
474 if (bc1[bs1 * dmapjoint1[j] + k])
475 {
476 // Zero column bs1 * j + k
477 for (int m = 0; m < num_rows; ++m)
478 Ae[m * num_cols + bs1 * j + k] = 0;
479 }
480 }
481 }
482 }
483
484 mat_set(dmapjoint0, dmapjoint1, Ae);
485 }
486}
487
493template <dolfinx::scalar T, std::floating_point U>
494void assemble_matrix(
495 la::MatSet<T> auto mat_set, const Form<T, U>& a, mdspan2_t x_dofmap,
496 std::span<const scalar_value_type_t<T>> x, std::span<const T> constants,
497 const std::map<std::pair<IntegralType, int>,
498 std::pair<std::span<const T>, int>>& coefficients,
499 std::span<const std::int8_t> bc0, std::span<const std::int8_t> bc1)
500{
501 // Integration domain mesh
502 std::shared_ptr<const mesh::Mesh<U>> mesh = a.mesh();
503 assert(mesh);
504 // Test function mesh
505 auto mesh0 = a.function_spaces().at(0)->mesh();
506 assert(mesh0);
507 // Trial function mesh
508 auto mesh1 = a.function_spaces().at(1)->mesh();
509 assert(mesh1);
510
511 // Get dofmap data
512 std::shared_ptr<const fem::DofMap> dofmap0
513 = a.function_spaces().at(0)->dofmap();
514 std::shared_ptr<const fem::DofMap> dofmap1
515 = a.function_spaces().at(1)->dofmap();
516 assert(dofmap0);
517 assert(dofmap1);
518 auto dofs0 = dofmap0->map();
519 const int bs0 = dofmap0->bs();
520 auto dofs1 = dofmap1->map();
521 const int bs1 = dofmap1->bs();
522
523 auto element0 = a.function_spaces().at(0)->element();
524 assert(element0);
525 auto element1 = a.function_spaces().at(1)->element();
526 assert(element1);
528 = element0->template dof_transformation_fn<T>(doftransform::standard);
530 = element1->template dof_transformation_right_fn<T>(
532
533 std::span<const std::uint32_t> cell_info0;
534 std::span<const std::uint32_t> cell_info1;
535 if (element0->needs_dof_transformations()
536 or element1->needs_dof_transformations() or a.needs_facet_permutations())
537 {
538 mesh0->topology_mutable()->create_entity_permutations();
539 mesh1->topology_mutable()->create_entity_permutations();
540 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
541 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
542 }
543
544 for (int i : a.integral_ids(IntegralType::cell))
545 {
546 auto fn = a.kernel(IntegralType::cell, i);
547 assert(fn);
548 auto& [coeffs, cstride] = coefficients.at({IntegralType::cell, i});
549 impl::assemble_cells(
550 mat_set, x_dofmap, x, a.domain(IntegralType::cell, i),
551 {dofs0, bs0, a.domain(IntegralType::cell, i, *mesh0)}, P0,
552 {dofs1, bs1, a.domain(IntegralType::cell, i, *mesh1)}, P1T, bc0, bc1,
553 fn, coeffs, cstride, constants, cell_info0, cell_info1);
554 }
555
556 std::span<const std::uint8_t> perms;
557 if (a.needs_facet_permutations())
558 {
559 mesh->topology_mutable()->create_entity_permutations();
560 perms = std::span(mesh->topology()->get_facet_permutations());
561 }
562
563 mesh::CellType cell_type = mesh->topology()->cell_type();
564 int num_facets_per_cell
565 = mesh::cell_num_entities(cell_type, mesh->topology()->dim() - 1);
566 for (int i : a.integral_ids(IntegralType::exterior_facet))
567 {
568 auto fn = a.kernel(IntegralType::exterior_facet, i);
569 assert(fn);
570 auto& [coeffs, cstride]
571 = coefficients.at({IntegralType::exterior_facet, i});
572 impl::assemble_exterior_facets(
573 mat_set, x_dofmap, x, num_facets_per_cell,
574 a.domain(IntegralType::exterior_facet, i),
575 {dofs0, bs0, a.domain(IntegralType::exterior_facet, i, *mesh0)}, P0,
576 {dofs1, bs1, a.domain(IntegralType::exterior_facet, i, *mesh1)}, P1T,
577 bc0, bc1, fn, coeffs, cstride, constants, cell_info0, cell_info1,
578 perms);
579 }
580
581 for (int i : a.integral_ids(IntegralType::interior_facet))
582 {
583 const std::vector<int> c_offsets = a.coefficient_offsets();
584 auto fn = a.kernel(IntegralType::interior_facet, i);
585 assert(fn);
586 auto& [coeffs, cstride]
587 = coefficients.at({IntegralType::interior_facet, i});
588 impl::assemble_interior_facets(
589 mat_set, x_dofmap, x, num_facets_per_cell,
590 a.domain(IntegralType::interior_facet, i),
591 {*dofmap0, bs0, a.domain(IntegralType::interior_facet, i, *mesh0)}, P0,
592 {*dofmap1, bs1, a.domain(IntegralType::interior_facet, i, *mesh1)}, P1T,
593 bc0, bc1, fn, coeffs, cstride, c_offsets, constants, cell_info0,
594 cell_info1, perms);
595 }
596}
597
598} // namespace dolfinx::fem::impl
Degree-of-freedom map representations and tools.
Degree-of-freedom map.
Definition DofMap.h:76
A representation of finite element variational forms.
Definition Form.h:139
DOF transform kernel concept.
Definition traits.h:20
Finite element cell kernel concept.
Definition traits.h:28
Matrix accumulate/set concept for functions that can be used in assemblers to accumulate or set value...
Definition utils.h:28
Functions supporting finite element method operations.
@ 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