DOLFINx 0.10.0.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{
28using mdspan2_t = md::mdspan<const std::int32_t, md::dextents<std::size_t, 2>>;
29
63template <dolfinx::scalar T>
64void assemble_cells_matrix(
65 la::MatSet<T> auto mat_set, mdspan2_t x_dofmap,
66 md::mdspan<const scalar_value_t<T>,
67 md::extents<std::size_t, md::dynamic_extent, 3>>
68 x,
69 std::span<const std::int32_t> cells,
70 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap0,
71 fem::DofTransformKernel<T> auto P0,
72 std::tuple<mdspan2_t, int, std::span<const std::int32_t>> dofmap1,
73 fem::DofTransformKernel<T> auto P1T, std::span<const std::int8_t> bc0,
74 std::span<const std::int8_t> bc1, FEkernel<T> auto kernel,
75 md::mdspan<const T, md::dextents<std::size_t, 2>> coeffs,
76 std::span<const T> constants, std::span<const std::uint32_t> cell_info0,
77 std::span<const std::uint32_t> cell_info1)
78{
79 if (cells.empty())
80 return;
81
82 const auto [dmap0, bs0, cells0] = dofmap0;
83 const auto [dmap1, bs1, cells1] = dofmap1;
84
85 // Iterate over active cells
86 const int num_dofs0 = dmap0.extent(1);
87 const int num_dofs1 = dmap1.extent(1);
88 const int ndim0 = bs0 * num_dofs0;
89 const int ndim1 = bs1 * num_dofs1;
90 std::vector<T> Ae(ndim0 * ndim1);
91 std::vector<scalar_value_t<T>> cdofs(3 * x_dofmap.extent(1));
92
93 // Iterate over active cells
94 assert(cells0.size() == cells.size());
95 assert(cells1.size() == cells.size());
96 for (std::size_t c = 0; c < cells.size(); ++c)
97 {
98 // Cell index in integration domain mesh (c), test function mesh
99 // (c0) and trial function mesh (c1)
100 std::int32_t cell = cells[c];
101 std::int32_t cell0 = cells0[c];
102 std::int32_t cell1 = cells1[c];
103
104 // Get cell coordinates/geometry
105 auto x_dofs = md::submdspan(x_dofmap, cell, md::full_extent);
106 for (std::size_t i = 0; i < x_dofs.size(); ++i)
107 std::copy_n(&x(x_dofs[i], 0), 3, std::next(cdofs.begin(), 3 * i));
108
109 // Tabulate tensor
110 std::ranges::fill(Ae, 0);
111 kernel(Ae.data(), &coeffs(c, 0), constants.data(), cdofs.data(), nullptr,
112 nullptr, nullptr);
113
114 // Compute A = P_0 \tilde{A} P_1^T (dof transformation)
115 P0(Ae, cell_info0, cell0, ndim1); // B = P0 \tilde{A}
116 P1T(Ae, cell_info1, cell1, ndim0); // A = B P1_T
117
118 // Zero rows/columns for essential bcs
119 std::span dofs0(dmap0.data_handle() + cell0 * num_dofs0, num_dofs0);
120 std::span dofs1(dmap1.data_handle() + cell1 * 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
193template <dolfinx::scalar T>
194void assemble_exterior_facets(
195 la::MatSet<T> auto mat_set, mdspan2_t x_dofmap,
196 md::mdspan<const scalar_value_t<T>,
197 md::extents<std::size_t, md::dynamic_extent, 3>>
198 x,
199 md::mdspan<const std::int32_t,
200 std::extents<std::size_t, md::dynamic_extent, 2>>
201 facets,
202 std::tuple<mdspan2_t, int,
203 md::mdspan<const std::int32_t,
204 std::extents<std::size_t, md::dynamic_extent, 2>>>
205 dofmap0,
206 fem::DofTransformKernel<T> auto P0,
207 std::tuple<mdspan2_t, int,
208 md::mdspan<const std::int32_t,
209 std::extents<std::size_t, md::dynamic_extent, 2>>>
210 dofmap1,
211 fem::DofTransformKernel<T> auto P1T, std::span<const std::int8_t> bc0,
212 std::span<const std::int8_t> bc1, FEkernel<T> auto kernel,
213 md::mdspan<const T, md::dextents<std::size_t, 2>> coeffs,
214 std::span<const T> constants, std::span<const std::uint32_t> cell_info0,
215 std::span<const std::uint32_t> cell_info1,
216 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> perms)
217{
218 if (facets.empty())
219 return;
220
221 const auto [dmap0, bs0, facets0] = dofmap0;
222 const auto [dmap1, bs1, facets1] = dofmap1;
223
224 // Data structures used in assembly
225 std::vector<scalar_value_t<T>> cdofs(3 * x_dofmap.extent(1));
226 const int num_dofs0 = dmap0.extent(1);
227 const int num_dofs1 = dmap1.extent(1);
228 const int ndim0 = bs0 * num_dofs0;
229 const int ndim1 = bs1 * num_dofs1;
230 std::vector<T> Ae(ndim0 * ndim1);
231 assert(facets0.size() == facets.size());
232 assert(facets1.size() == facets.size());
233 for (std::size_t f = 0; f < facets.extent(0); ++f)
234 {
235 // Cell in the integration domain, local facet index relative to the
236 // integration domain cell, and cells in the test and trial function
237 // meshes
238 std::int32_t cell = facets(f, 0);
239 std::int32_t local_facet = facets(f, 1);
240 std::int32_t cell0 = facets0(f, 0);
241 std::int32_t cell1 = facets1(f, 0);
242
243 // Get cell coordinates/geometry
244 auto x_dofs = md::submdspan(x_dofmap, cell, md::full_extent);
245 for (std::size_t i = 0; i < x_dofs.size(); ++i)
246 std::copy_n(&x(x_dofs[i], 0), 3, std::next(cdofs.begin(), 3 * i));
247
248 // Permutations
249 std::uint8_t perm = perms.empty() ? 0 : perms(cell, local_facet);
250
251 // Tabulate tensor
252 std::ranges::fill(Ae, 0);
253 kernel(Ae.data(), &coeffs(f, 0), constants.data(), cdofs.data(),
254 &local_facet, &perm, nullptr);
255 P0(Ae, cell_info0, cell0, ndim1);
256 P1T(Ae, cell_info1, cell1, ndim0);
257
258 // Zero rows/columns for essential bcs
259 std::span dofs0(dmap0.data_handle() + cell0 * num_dofs0, num_dofs0);
260 std::span dofs1(dmap1.data_handle() + cell1 * num_dofs1, num_dofs1);
261 if (!bc0.empty())
262 {
263 for (int i = 0; i < num_dofs0; ++i)
264 {
265 for (int k = 0; k < bs0; ++k)
266 {
267 if (bc0[bs0 * dofs0[i] + k])
268 {
269 // Zero row bs0 * i + k
270 const int row = bs0 * i + k;
271 std::fill_n(std::next(Ae.begin(), ndim1 * row), ndim1, 0);
272 }
273 }
274 }
275 }
276 if (!bc1.empty())
277 {
278 for (int j = 0; j < num_dofs1; ++j)
279 {
280 for (int k = 0; k < bs1; ++k)
281 {
282 if (bc1[bs1 * dofs1[j] + k])
283 {
284 // Zero column bs1 * j + k
285 const int col = bs1 * j + k;
286 for (int row = 0; row < ndim0; ++row)
287 Ae[row * ndim1 + col] = 0;
288 }
289 }
290 }
291 }
292
293 mat_set(dofs0, dofs1, Ae);
294 }
295}
296
333template <dolfinx::scalar T>
334void assemble_interior_facets(
335 la::MatSet<T> auto mat_set, mdspan2_t x_dofmap,
336 md::mdspan<const scalar_value_t<T>,
337 md::extents<std::size_t, md::dynamic_extent, 3>>
338 x,
339 md::mdspan<const std::int32_t,
340 std::extents<std::size_t, md::dynamic_extent, 2, 2>>
341 facets,
342 std::tuple<const DofMap&, int,
343 md::mdspan<const std::int32_t,
344 std::extents<std::size_t, md::dynamic_extent, 2, 2>>>
345 dofmap0,
346 fem::DofTransformKernel<T> auto P0,
347 std::tuple<const DofMap&, int,
348 md::mdspan<const std::int32_t,
349 std::extents<std::size_t, md::dynamic_extent, 2, 2>>>
350 dofmap1,
351 fem::DofTransformKernel<T> auto P1T, std::span<const std::int8_t> bc0,
352 std::span<const std::int8_t> bc1, FEkernel<T> auto kernel,
353 md::mdspan<const T, md::extents<std::size_t, md::dynamic_extent, 2,
354 md::dynamic_extent>>
355 coeffs,
356 std::span<const T> constants, std::span<const std::uint32_t> cell_info0,
357 std::span<const std::uint32_t> cell_info1,
358 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> perms)
359{
360 if (facets.empty())
361 return;
362
363 const auto [dmap0, bs0, facets0] = dofmap0;
364 const auto [dmap1, bs1, facets1] = dofmap1;
365
366 // Data structures used in assembly
367 using X = scalar_value_t<T>;
368 std::vector<X> cdofs(2 * x_dofmap.extent(1) * 3);
369 std::span<X> cdofs0(cdofs.data(), x_dofmap.extent(1) * 3);
370 std::span<X> cdofs1(cdofs.data() + x_dofmap.extent(1) * 3,
371 x_dofmap.extent(1) * 3);
372
373 const std::size_t dmap0_size = dmap0.map().extent(1);
374 const std::size_t dmap1_size = dmap1.map().extent(1);
375 const int num_rows = bs0 * 2 * dmap0_size;
376 const int num_cols = bs1 * 2 * dmap1_size;
377
378 // Temporaries for joint dofmaps
379 std::vector<T> Ae(num_rows * num_cols), be(num_rows);
380 std::vector<std::int32_t> dmapjoint0(2 * dmap0_size);
381 std::vector<std::int32_t> dmapjoint1(2 * dmap1_size);
382 assert(facets0.size() == facets.size());
383 assert(facets1.size() == facets.size());
384 for (std::size_t f = 0; f < facets.extent(0); ++f)
385 {
386 // Cells in integration domain, test function domain and trial
387 // function domain
388 std::array cells{facets(f, 0, 0), facets(f, 1, 0)};
389 std::array cells0{facets0(f, 0, 0), facets0(f, 1, 0)};
390 std::array cells1{facets1(f, 0, 0), facets1(f, 1, 0)};
391
392 // Local facets indices
393 std::array local_facet{facets(f, 0, 1), facets(f, 1, 1)};
394
395 // Get cell geometry
396 auto x_dofs0 = md::submdspan(x_dofmap, cells[0], md::full_extent);
397 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
398 std::copy_n(&x(x_dofs0[i], 0), 3, std::next(cdofs0.begin(), 3 * i));
399 auto x_dofs1 = md::submdspan(x_dofmap, cells[1], md::full_extent);
400 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
401 std::copy_n(&x(x_dofs1[i], 0), 3, std::next(cdofs1.begin(), 3 * i));
402
403 // Get dof maps for cells and pack
404 // When integrating over interfaces between two domains, the test function
405 // might only be defined on one side, so we check which cells exist in the
406 // test function domain
407 std::span<const std::int32_t> dmap0_cell0
408 = cells0[0] >= 0 ? dmap0.cell_dofs(cells0[0])
409 : std::span<const std::int32_t>();
410 std::span<const std::int32_t> dmap0_cell1
411 = cells0[1] >= 0 ? dmap0.cell_dofs(cells0[1])
412 : std::span<const std::int32_t>();
413
414 std::ranges::copy(dmap0_cell0, dmapjoint0.begin());
415 std::ranges::copy(dmap0_cell1, std::next(dmapjoint0.begin(), dmap0_size));
416
417 // Check which cells exist in the trial function domain
418 std::span<const std::int32_t> dmap1_cell0
419 = cells1[0] >= 0 ? dmap1.cell_dofs(cells1[0])
420 : std::span<const std::int32_t>();
421 std::span<const std::int32_t> dmap1_cell1
422 = cells1[1] >= 0 ? dmap1.cell_dofs(cells1[1])
423 : std::span<const std::int32_t>();
424
425 std::ranges::copy(dmap1_cell0, dmapjoint1.begin());
426 std::ranges::copy(dmap1_cell1, std::next(dmapjoint1.begin(), dmap1_size));
427
428 // Tabulate tensor
429 std::ranges::fill(Ae, 0);
430 std::array perm = perms.empty()
431 ? std::array<std::uint8_t, 2>{0, 0}
432 : std::array{perms(cells[0], local_facet[0]),
433 perms(cells[1], local_facet[1])};
434 kernel(Ae.data(), &coeffs(f, 0, 0), constants.data(), cdofs.data(),
435 local_facet.data(), perm.data(), nullptr);
436
437 // Local element layout is a 2x2 block matrix with structure
438 //
439 // cell0cell0 | cell0cell1
440 // cell1cell0 | cell1cell1
441 //
442 // where each block is element tensor of size (dmap0, dmap1).
443
444 // Only apply transformation when cells exist
445 if (cells0[0] >= 0)
446 P0(Ae, cell_info0, cells0[0], num_cols);
447 if (cells0[1] >= 0)
448 {
449 std::span sub_Ae0(Ae.data() + bs0 * dmap0_size * num_cols,
450 bs0 * dmap0_size * num_cols);
451
452 P0(sub_Ae0, cell_info0, cells0[1], num_cols);
453 }
454 if (cells1[0] >= 0)
455 P1T(Ae, cell_info1, cells1[0], num_rows);
456
457 if (cells1[1] >= 0)
458 {
459 for (int row = 0; row < num_rows; ++row)
460 {
461 // DOFs for dmap1 and cell1 are not stored contiguously in the
462 // block matrix, so each row needs a separate span access
463 std::span sub_Ae1(Ae.data() + row * num_cols + bs1 * dmap1_size,
464 bs1 * dmap1_size);
465 P1T(sub_Ae1, cell_info1, cells1[1], 1);
466 }
467 }
468
469 // Zero rows/columns for essential bcs
470 if (!bc0.empty())
471 {
472 for (std::size_t i = 0; i < dmapjoint0.size(); ++i)
473 {
474 for (int k = 0; k < bs0; ++k)
475 {
476 if (bc0[bs0 * dmapjoint0[i] + k])
477 {
478 // Zero row bs0 * i + k
479 std::fill_n(std::next(Ae.begin(), num_cols * (bs0 * i + k)),
480 num_cols, 0);
481 }
482 }
483 }
484 }
485 if (!bc1.empty())
486 {
487 for (std::size_t j = 0; j < dmapjoint1.size(); ++j)
488 {
489 for (int k = 0; k < bs1; ++k)
490 {
491 if (bc1[bs1 * dmapjoint1[j] + k])
492 {
493 // Zero column bs1 * j + k
494 for (int m = 0; m < num_rows; ++m)
495 Ae[m * num_cols + bs1 * j + k] = 0;
496 }
497 }
498 }
499 }
500
501 mat_set(dmapjoint0, dmapjoint1, Ae);
502 }
503}
504
523template <dolfinx::scalar T, std::floating_point U>
524void assemble_matrix(
525 la::MatSet<T> auto mat_set, const Form<T, U>& a,
526 md::mdspan<const scalar_value_t<T>,
527 md::extents<std::size_t, md::dynamic_extent, 3>>
528 x,
529 std::span<const T> constants,
530 const std::map<std::pair<IntegralType, int>,
531 std::pair<std::span<const T>, int>>& coefficients,
532 std::span<const std::int8_t> bc0, std::span<const std::int8_t> bc1)
533{
534 // Integration domain mesh
535 std::shared_ptr<const mesh::Mesh<U>> mesh = a.mesh();
536 assert(mesh);
537
538 // Test function mesh
539 auto mesh0 = a.function_spaces().at(0)->mesh();
540 assert(mesh0);
541
542 // Trial function mesh
543 auto mesh1 = a.function_spaces().at(1)->mesh();
544 assert(mesh1);
545
546 // TODO: Mixed topology with exterior and interior facet integrals.
547 //
548 // NOTE: Can't just loop over cell types for interior facet integrals
549 // because we have a kernel per combination of comparable cell types,
550 // rather than one per cell type. Also, we need the dofmaps for two
551 // different cell types at the same time.
552 const int num_cell_types = mesh->topology()->cell_types().size();
553 for (int cell_type_idx = 0; cell_type_idx < num_cell_types; ++cell_type_idx)
554 {
555 // Geometry dofmap and data
556 mdspan2_t x_dofmap = mesh->geometry().dofmap(cell_type_idx);
557
558 // Get dofmap data
559 std::shared_ptr<const fem::DofMap> dofmap0
560 = a.function_spaces().at(0)->dofmaps(cell_type_idx);
561 std::shared_ptr<const fem::DofMap> dofmap1
562 = a.function_spaces().at(1)->dofmaps(cell_type_idx);
563 assert(dofmap0);
564 assert(dofmap1);
565 auto dofs0 = dofmap0->map();
566 const int bs0 = dofmap0->bs();
567 auto dofs1 = dofmap1->map();
568 const int bs1 = dofmap1->bs();
569
570 auto element0 = a.function_spaces().at(0)->elements(cell_type_idx);
571 assert(element0);
572 auto element1 = a.function_spaces().at(1)->elements(cell_type_idx);
573 assert(element1);
574 fem::DofTransformKernel<T> auto P0
575 = element0->template dof_transformation_fn<T>(doftransform::standard);
576 fem::DofTransformKernel<T> auto P1T
577 = element1->template dof_transformation_right_fn<T>(
579
580 std::span<const std::uint32_t> cell_info0;
581 std::span<const std::uint32_t> cell_info1;
582 if (element0->needs_dof_transformations()
583 or element1->needs_dof_transformations()
584 or a.needs_facet_permutations())
585 {
586 mesh0->topology_mutable()->create_entity_permutations();
587 mesh1->topology_mutable()->create_entity_permutations();
588 cell_info0 = std::span(mesh0->topology()->get_cell_permutation_info());
589 cell_info1 = std::span(mesh1->topology()->get_cell_permutation_info());
590 }
591
592 for (int i = 0; i < a.num_integrals(IntegralType::cell, cell_type_idx); ++i)
593 {
594 auto fn = a.kernel(IntegralType::cell, i, cell_type_idx);
595 assert(fn);
596 std::span cells = a.domain(IntegralType::cell, i, cell_type_idx);
597 std::span cells0 = a.domain_arg(IntegralType::cell, 0, i, cell_type_idx);
598 std::span cells1 = a.domain_arg(IntegralType::cell, 1, i, cell_type_idx);
599 auto& [coeffs, cstride] = coefficients.at({IntegralType::cell, i});
600 assert(cells.size() * cstride == coeffs.size());
601 impl::assemble_cells_matrix(
602 mat_set, x_dofmap, x, cells, {dofs0, bs0, cells0}, P0,
603 {dofs1, bs1, cells1}, P1T, bc0, bc1, fn,
604 md::mdspan(coeffs.data(), cells.size(), cstride), constants,
605 cell_info0, cell_info1);
606 }
607
608 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> perms;
609 if (a.needs_facet_permutations())
610 {
611 mesh::CellType cell_type = mesh->topology()->cell_types()[cell_type_idx];
612 int num_facets_per_cell
613 = mesh::cell_num_entities(cell_type, mesh->topology()->dim() - 1);
614 mesh->topology_mutable()->create_entity_permutations();
615 const std::vector<std::uint8_t>& p
616 = mesh->topology()->get_facet_permutations();
617 perms = md::mdspan(p.data(), p.size() / num_facets_per_cell,
618 num_facets_per_cell);
619 }
620
621 for (int i = 0;
622 i < a.num_integrals(IntegralType::exterior_facet, cell_type_idx); ++i)
623 {
624 if (num_cell_types > 1)
625 {
626 throw std::runtime_error("Exterior facet integrals with mixed "
627 "topology aren't supported yet");
628 }
629
630 using mdspanx2_t
631 = md::mdspan<const std::int32_t,
632 md::extents<std::size_t, md::dynamic_extent, 2>>;
633
634 auto fn = a.kernel(IntegralType::exterior_facet, i, 0);
635 assert(fn);
636 auto& [coeffs, cstride]
637 = coefficients.at({IntegralType::exterior_facet, i});
638
639 std::span f = a.domain(IntegralType::exterior_facet, i, 0);
640 mdspanx2_t facets(f.data(), f.size() / 2, 2);
641 std::span f0 = a.domain_arg(IntegralType::exterior_facet, 0, i, 0);
642 mdspanx2_t facets0(f0.data(), f0.size() / 2, 2);
643 std::span f1 = a.domain_arg(IntegralType::exterior_facet, 1, i, 0);
644 mdspanx2_t facets1(f1.data(), f1.size() / 2, 2);
645 assert((facets.size() / 2) * cstride == coeffs.size());
646 impl::assemble_exterior_facets(
647 mat_set, x_dofmap, x, facets, {dofs0, bs0, facets0}, P0,
648 {dofs1, bs1, facets1}, P1T, bc0, bc1, fn,
649 md::mdspan(coeffs.data(), facets.extent(0), cstride), constants,
650 cell_info0, cell_info1, perms);
651 }
652
653 for (int i = 0;
654 i < a.num_integrals(IntegralType::interior_facet, cell_type_idx); ++i)
655 {
656 if (num_cell_types > 1)
657 {
658 throw std::runtime_error("Interior facet integrals with mixed "
659 "topology aren't supported yet");
660 }
661
662 using mdspanx22_t
663 = md::mdspan<const std::int32_t,
664 md::extents<std::size_t, md::dynamic_extent, 2, 2>>;
665 using mdspanx2x_t
666 = md::mdspan<const T, md::extents<std::size_t, md::dynamic_extent, 2,
667 md::dynamic_extent>>;
668
669 auto fn = a.kernel(IntegralType::interior_facet, i, 0);
670 assert(fn);
671 auto& [coeffs, cstride]
672 = coefficients.at({IntegralType::interior_facet, i});
673
674 std::span facets = a.domain(IntegralType::interior_facet, i, 0);
675 std::span facets0 = a.domain_arg(IntegralType::interior_facet, 0, i, 0);
676 std::span facets1 = a.domain_arg(IntegralType::interior_facet, 1, i, 0);
677 assert((facets.size() / 4) * 2 * cstride == coeffs.size());
678 impl::assemble_interior_facets(
679 mat_set, x_dofmap, x,
680 mdspanx22_t(facets.data(), facets.size() / 4, 2, 2),
681 {*dofmap0, bs0,
682 mdspanx22_t(facets0.data(), facets0.size() / 4, 2, 2)},
683 P0,
684 {*dofmap1, bs1,
685 mdspanx22_t(facets1.data(), facets1.size() / 4, 2, 2)},
686 P1T, bc0, bc1, fn,
687 mdspanx2x_t(coeffs.data(), facets.size() / 4, 2, cstride), constants,
688 cell_info0, cell_info1, perms);
689 }
690 }
691}
692
693} // namespace dolfinx::fem::impl
Degree-of-freedom map representations and tools.
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
@ transpose
Transpose.
Definition FiniteElement.h:28
@ standard
Standard.
Definition FiniteElement.h:27
@ interior_facet
Interior facet.
Definition Form.h:42
@ cell
Cell.
Definition Form.h:40
@ exterior_facet
Exterior facet.
Definition Form.h:41
CellType
Cell type identifier.
Definition cell_types.h:22
int cell_num_entities(CellType type, int dim)
Definition cell_types.cpp:139