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
24template <dolfinx::scalar T>
25T assemble_cells(mdspan2_t x_dofmap, std::span<
const scalar_value_type_t<T>> x,
26 std::span<const std::int32_t> cells, FEkernel<T>
auto fn,
27 std::span<const T> constants, std::span<const T> coeffs,
35 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
38 for (std::size_t index = 0; index <
cells.size(); ++index)
40 std::int32_t c =
cells[index];
43 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
44 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
45 for (std::size_t i = 0; i < x_dofs.size(); ++i)
47 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
48 std::next(coordinate_dofs.begin(), 3 * i));
51 const T* coeff_cell = coeffs.data() + index * cstride;
52 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(),
nullptr,
60template <dolfinx::scalar T>
61T assemble_exterior_facets(mdspan2_t x_dofmap,
62 std::span<
const scalar_value_type_t<T>> x,
63 int num_facets_per_cell,
64 std::span<const std::int32_t> facets,
65 FEkernel<T>
auto fn, std::span<const T> constants,
66 std::span<const T> coeffs,
int cstride,
67 std::span<const std::uint8_t> perms)
74 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
77 assert(facets.size() % 2 == 0);
78 for (std::size_t index = 0; index < facets.size(); index += 2)
80 std::int32_t
cell = facets[index];
81 std::int32_t local_facet = facets[index + 1];
84 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
85 x_dofmap,
cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
86 for (std::size_t i = 0; i < x_dofs.size(); ++i)
88 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
89 std::next(coordinate_dofs.begin(), 3 * i));
94 = perms.empty() ? 0 : perms[
cell * num_facets_per_cell + local_facet];
95 const T* coeff_cell = coeffs.data() + index / 2 * cstride;
96 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(),
104template <dolfinx::scalar T>
105T assemble_interior_facets(mdspan2_t x_dofmap,
106 std::span<
const scalar_value_type_t<T>> x,
107 int num_facets_per_cell,
108 std::span<const std::int32_t> facets,
109 FEkernel<T>
auto fn, std::span<const T> constants,
110 std::span<const T> coeffs,
int cstride,
111 std::span<const int> offsets,
112 std::span<const std::uint8_t> perms)
119 using X = scalar_value_type_t<T>;
120 std::vector<X> coordinate_dofs(2 * x_dofmap.extent(1) * 3);
121 std::span<X> cdofs0(coordinate_dofs.data(), x_dofmap.extent(1) * 3);
122 std::span<X> cdofs1(coordinate_dofs.data() + x_dofmap.extent(1) * 3,
123 x_dofmap.extent(1) * 3);
125 std::vector<T> coeff_array(2 * offsets.back());
126 assert(offsets.back() == cstride);
129 assert(facets.size() % 4 == 0);
130 for (std::size_t index = 0; index < facets.size(); index += 4)
132 std::array<std::int32_t, 2>
cells = {facets[index], facets[index + 2]};
133 std::array<std::int32_t, 2> local_facet
134 = {facets[index + 1], facets[index + 3]};
137 auto x_dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
138 x_dofmap, cells[0], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
139 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
141 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
142 std::next(cdofs0.begin(), 3 * i));
144 auto x_dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
145 x_dofmap, cells[1], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
146 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
148 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
149 std::next(cdofs1.begin(), 3 * i));
154 ? std::array<std::uint8_t, 2>{0, 0}
156 perms[
cells[0] * num_facets_per_cell + local_facet[0]],
157 perms[
cells[1] * num_facets_per_cell + local_facet[1]]};
158 fn(&value, coeffs.data() + index / 2 * cstride, constants.data(),
159 coordinate_dofs.data(), local_facet.data(), perm.data());
166template <dolfinx::scalar T, std::
floating_po
int U>
169 std::span<
const scalar_value_type_t<T>> x, std::span<const T> constants,
170 const std::map<std::pair<IntegralType, int>,
171 std::pair<std::span<const T>,
int>>& coefficients)
173 std::shared_ptr<const mesh::Mesh<U>> mesh = M.mesh();
183 value += impl::assemble_cells(x_dofmap, x, cells, fn, constants, coeffs,
187 std::span<const std::uint8_t> perms;
188 if (M.needs_facet_permutations())
190 mesh->topology_mutable()->create_entity_permutations();
191 perms = std::span(mesh->topology()->get_facet_permutations());
195 int num_facets_per_cell
201 auto& [coeffs, cstride]
203 value += impl::assemble_exterior_facets(
204 x_dofmap, x, num_facets_per_cell,
211 const std::vector<int> c_offsets = M.coefficient_offsets();
214 auto& [coeffs, cstride]
216 value += impl::assemble_interior_facets(
217 x_dofmap, x, num_facets_per_cell,
219 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