11#include "FunctionSpace.h"
14#include <dolfinx/common/IndexMap.h>
15#include <dolfinx/mesh/Geometry.h>
16#include <dolfinx/mesh/Mesh.h>
17#include <dolfinx/mesh/Topology.h>
21namespace dolfinx::fem::impl
25template <dolfinx::scalar T>
26T assemble_cells(mdspan2_t x_dofmap, std::span<
const scalar_value_type_t<T>> x,
27 std::span<const std::int32_t> cells, FEkernel<T>
auto fn,
28 std::span<const T> constants, std::span<const T> coeffs,
36 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
39 for (std::size_t index = 0; index <
cells.size(); ++index)
41 std::int32_t c =
cells[index];
44 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
45 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
46 for (std::size_t i = 0; i < x_dofs.size(); ++i)
48 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
49 std::next(coordinate_dofs.begin(), 3 * i));
52 const T* coeff_cell = coeffs.data() + index * cstride;
53 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(),
nullptr,
61template <dolfinx::scalar T>
62T assemble_exterior_facets(mdspan2_t x_dofmap,
63 std::span<
const scalar_value_type_t<T>> x,
64 int num_facets_per_cell,
65 std::span<const std::int32_t> facets,
66 FEkernel<T>
auto fn, std::span<const T> constants,
67 std::span<const T> coeffs,
int cstride,
68 std::span<const std::uint8_t> perms)
75 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
78 assert(facets.size() % 2 == 0);
79 for (std::size_t index = 0; index < facets.size(); index += 2)
81 std::int32_t
cell = facets[index];
82 std::int32_t local_facet = facets[index + 1];
85 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
86 x_dofmap,
cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
87 for (std::size_t i = 0; i < x_dofs.size(); ++i)
89 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
90 std::next(coordinate_dofs.begin(), 3 * i));
95 = perms.empty() ? 0 : perms[
cell * num_facets_per_cell + local_facet];
96 const T* coeff_cell = coeffs.data() + index / 2 * cstride;
97 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(),
105template <dolfinx::scalar T>
106T assemble_interior_facets(mdspan2_t x_dofmap,
107 std::span<
const scalar_value_type_t<T>> x,
108 int num_facets_per_cell,
109 std::span<const std::int32_t> facets,
110 FEkernel<T>
auto fn, std::span<const T> constants,
111 std::span<const T> coeffs,
int cstride,
112 std::span<const int> offsets,
113 std::span<const std::uint8_t> perms)
120 using X = scalar_value_type_t<T>;
121 std::vector<X> coordinate_dofs(2 * x_dofmap.extent(1) * 3);
122 std::span<X> cdofs0(coordinate_dofs.data(), x_dofmap.extent(1) * 3);
123 std::span<X> cdofs1(coordinate_dofs.data() + x_dofmap.extent(1) * 3,
124 x_dofmap.extent(1) * 3);
126 std::vector<T> coeff_array(2 * offsets.back());
127 assert(offsets.back() == cstride);
130 assert(facets.size() % 4 == 0);
131 for (std::size_t index = 0; index < facets.size(); index += 4)
133 std::array<std::int32_t, 2>
cells = {facets[index], facets[index + 2]};
134 std::array<std::int32_t, 2> local_facet
135 = {facets[index + 1], facets[index + 3]};
138 auto x_dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
139 x_dofmap, cells[0], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
140 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
142 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
143 std::next(cdofs0.begin(), 3 * i));
145 auto x_dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
146 x_dofmap, cells[1], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
147 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
149 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
150 std::next(cdofs1.begin(), 3 * i));
155 ? std::array<std::uint8_t, 2>{0, 0}
157 perms[
cells[0] * num_facets_per_cell + local_facet[0]],
158 perms[
cells[1] * num_facets_per_cell + local_facet[1]]};
159 fn(&value, coeffs.data() + index / 2 * cstride, constants.data(),
160 coordinate_dofs.data(), local_facet.data(), perm.data());
167template <dolfinx::scalar T, std::
floating_po
int U>
169 const fem::Form<T, U>& M, mdspan2_t x_dofmap,
170 std::span<
const scalar_value_type_t<T>> x, std::span<const T> constants,
171 const std::map<std::pair<IntegralType, int>,
172 std::pair<std::span<const T>,
int>>& coefficients)
174 std::shared_ptr<const mesh::Mesh<U>> mesh = M.mesh();
184 value += impl::assemble_cells(x_dofmap, x, cells, fn, constants, coeffs,
188 std::span<const std::uint8_t> perms;
189 if (M.needs_facet_permutations())
191 mesh->topology_mutable()->create_entity_permutations();
192 perms = std::span(mesh->topology()->get_facet_permutations());
196 int num_facets_per_cell
202 auto& [coeffs, cstride]
204 value += impl::assemble_exterior_facets(
205 x_dofmap, x, num_facets_per_cell,
212 const std::vector<int> c_offsets = M.coefficient_offsets();
215 auto& [coeffs, cstride]
217 value += impl::assemble_interior_facets(
218 x_dofmap, x, num_facets_per_cell,
220 cstride, c_offsets, 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
@ 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