11#include "FunctionSpace.h"
14#include <basix/mdspan.hpp>
15#include <dolfinx/common/IndexMap.h>
16#include <dolfinx/mesh/Geometry.h>
17#include <dolfinx/mesh/Mesh.h>
18#include <dolfinx/mesh/Topology.h>
22namespace dolfinx::fem::impl
25template <dolfinx::scalar T>
26T assemble_cells(mdspan2_t x_dofmap,
27 md::mdspan<
const scalar_value_t<T>,
28 md::extents<std::size_t, md::dynamic_extent, 3>>
30 std::span<const std::int32_t> cells, FEkernel<T>
auto fn,
31 std::span<const T> constants,
32 md::mdspan<
const T, md::dextents<std::size_t, 2>> coeffs,
33 std::span<scalar_value_t<T>> cdofs_b)
39 assert(cdofs_b.size() >= 3 * x_dofmap.extent(1));
42 for (std::size_t index = 0; index <
cells.size(); ++index)
44 std::int32_t c =
cells[index];
47 auto x_dofs = md::submdspan(x_dofmap, c, md::full_extent);
48 for (std::size_t i = 0; i < x_dofs.size(); ++i)
49 std::copy_n(&x(x_dofs[i], 0), 3, std::next(cdofs_b.begin(), 3 * i));
51 fn(&value, &coeffs(index, 0), constants.data(), cdofs_b.data(),
nullptr,
68template <dolfinx::scalar T>
71 md::mdspan<
const scalar_value_t<T>,
72 md::extents<std::size_t, md::dynamic_extent, 3>>
74 md::mdspan<
const std::int32_t,
75 md::extents<std::size_t, md::dynamic_extent, 2>>
77 FEkernel<T>
auto fn, std::span<const T> constants,
78 md::mdspan<
const T, md::dextents<std::size_t, 2>> coeffs,
79 md::mdspan<
const std::uint8_t, md::dextents<std::size_t, 2>> perms,
80 std::span<scalar_value_t<T>> cdofs_b)
86 assert(cdofs_b.size() >= 3 * x_dofmap.extent(1));
89 for (std::size_t f = 0; f < entities.extent(0); ++f)
91 std::int32_t
cell = entities(f, 0);
92 std::int32_t local_entity = entities(f, 1);
95 auto x_dofs = md::submdspan(x_dofmap,
cell, md::full_extent);
96 for (std::size_t i = 0; i < x_dofs.size(); ++i)
97 std::copy_n(&x(x_dofs[i], 0), 3, std::next(cdofs_b.begin(), 3 * i));
100 std::uint8_t perm = perms.empty() ? 0 : perms(
cell, local_entity);
101 fn(&value, &coeffs(f, 0), constants.data(), cdofs_b.data(), &local_entity,
109template <dolfinx::scalar T>
110T assemble_interior_facets(
112 md::mdspan<
const scalar_value_t<T>,
113 md::extents<std::size_t, md::dynamic_extent, 3>>
115 md::mdspan<
const std::int32_t,
116 md::extents<std::size_t, md::dynamic_extent, 2, 2>>
118 FEkernel<T>
auto fn, std::span<const T> constants,
119 md::mdspan<
const T, md::extents<std::size_t, md::dynamic_extent, 2,
122 md::mdspan<
const std::uint8_t, md::dextents<std::size_t, 2>> perms,
123 std::span<scalar_value_t<T>> cdofs_b)
130 assert(cdofs_b.size() >= 2 * 3 * x_dofmap.extent(1));
131 auto cdofs0 = cdofs_b.first(3 * x_dofmap.extent(1));
132 auto cdofs1 = cdofs_b.last(3 * x_dofmap.extent(1));
135 for (std::size_t f = 0; f < facets.extent(0); ++f)
137 std::array
cells = {facets(f, 0, 0), facets(f, 1, 0)};
138 std::array local_facet = {facets(f, 0, 1), facets(f, 1, 1)};
141 auto x_dofs0 = md::submdspan(x_dofmap, cells[0], md::full_extent);
142 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
143 std::copy_n(&x(x_dofs0[i], 0), 3, std::next(cdofs0.begin(), 3 * i));
144 auto x_dofs1 = md::submdspan(x_dofmap, cells[1], md::full_extent);
145 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
146 std::copy_n(&x(x_dofs1[i], 0), 3, std::next(cdofs1.begin(), 3 * i));
148 std::array perm = perms.empty()
149 ? std::array<std::uint8_t, 2>{0, 0}
150 : std::array{perms(cells[0], local_facet[0]),
151 perms(cells[1], local_facet[1])};
152 fn(&value, &coeffs(f, 0, 0), constants.data(), cdofs_b.data(),
153 local_facet.data(), perm.data(),
nullptr);
160template <dolfinx::scalar T, std::
floating_po
int U>
162 const fem::Form<T, U>& M, mdspan2_t x_dofmap,
163 md::mdspan<
const scalar_value_t<T>,
164 md::extents<std::size_t, md::dynamic_extent, 3>>
166 std::span<const T> constants,
167 const std::map<std::pair<IntegralType, int>,
168 std::pair<std::span<const T>,
int>>& coefficients)
170 std::shared_ptr<const mesh::Mesh<U>> mesh = M.mesh();
173 std::vector<scalar_value_t<T>> cdofs_b(2 * 3 * x_dofmap.extent(1));
182 assert(
cells.size() * cstride == coeffs.size());
183 value += impl::assemble_cells(
184 x_dofmap, x, cells, fn, constants,
185 md::mdspan(coeffs.data(),
cells.size(), cstride), cdofs_b);
189 int num_facets_per_cell
191 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> facet_perms;
192 if (M.needs_facet_permutations())
194 mesh->topology_mutable()->create_entity_permutations();
195 const std::vector<std::uint8_t>& p
196 = mesh->topology()->get_facet_permutations();
197 facet_perms = md::mdspan(p.data(), p.size() / num_facets_per_cell,
198 num_facets_per_cell);
205 auto& [coeffs, cstride]
209 constexpr std::size_t num_adjacent_cells = 2;
211 constexpr std::size_t shape1 = 2 * num_adjacent_cells;
213 assert((facets.size() / shape1) * 2 * cstride == coeffs.size());
214 value += impl::assemble_interior_facets(
216 md::mdspan<
const std::int32_t,
217 md::extents<std::size_t, md::dynamic_extent, 2, 2>>(
218 facets.data(), facets.size() / shape1, 2, 2),
220 md::mdspan<
const T, md::extents<std::size_t, md::dynamic_extent, 2,
221 md::dynamic_extent>>(
222 coeffs.data(), facets.size() / shape1, 2, cstride),
223 facet_perms, cdofs_b);
229 md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>> perms
232 : md::mdspan<const std::uint8_t, md::dextents<std::size_t, 2>>{};
234 for (
int i = 0; i < M.num_integrals(itg_type, 0); ++i)
236 auto fn = M.kernel(itg_type, i, 0);
238 auto& [coeffs, cstride] = coefficients.at({itg_type, i});
240 std::span entities = M.domain(itg_type, i, 0);
243 assert((entities.size() / 2) * cstride == coeffs.size());
244 value += impl::assemble_entities(
246 md::mdspan<
const std::int32_t,
247 md::extents<std::size_t, md::dynamic_extent, 2>>(
248 entities.data(), entities.size() / 2, 2),
250 md::mdspan(coeffs.data(), entities.size() / 2, cstride), perms,
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
@ 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