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 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)
73 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
76 assert(facets.size() % 2 == 0);
77 for (std::size_t index = 0; index < facets.size(); index += 2)
79 std::int32_t
cell = facets[index];
80 std::int32_t local_facet = facets[index + 1];
83 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
84 x_dofmap,
cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
85 for (std::size_t i = 0; i < x_dofs.size(); ++i)
87 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
88 std::next(coordinate_dofs.begin(), 3 * i));
91 const T* coeff_cell = coeffs.data() + index / 2 * cstride;
92 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(),
93 &local_facet,
nullptr);
100template <dolfinx::scalar T>
101T assemble_interior_facets(mdspan2_t x_dofmap,
102 std::span<
const scalar_value_type_t<T>> x,
104 std::span<const std::int32_t> facets,
105 FEkernel<T>
auto fn, std::span<const T> constants,
106 std::span<const T> coeffs,
int cstride,
107 std::span<const int> offsets,
108 std::span<const std::uint8_t> perms)
115 using X = scalar_value_type_t<T>;
116 std::vector<X> coordinate_dofs(2 * x_dofmap.extent(1) * 3);
117 std::span<X> cdofs0(coordinate_dofs.data(), x_dofmap.extent(1) * 3);
118 std::span<X> cdofs1(coordinate_dofs.data() + x_dofmap.extent(1) * 3,
119 x_dofmap.extent(1) * 3);
121 std::vector<T> coeff_array(2 * offsets.back());
122 assert(offsets.back() == cstride);
125 assert(facets.size() % 4 == 0);
126 for (std::size_t index = 0; index < facets.size(); index += 4)
128 std::array<std::int32_t, 2>
cells = {facets[index], facets[index + 2]};
129 std::array<std::int32_t, 2> local_facet
130 = {facets[index + 1], facets[index + 3]};
133 auto x_dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
134 x_dofmap, cells[0], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
135 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
137 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
138 std::next(cdofs0.begin(), 3 * i));
140 auto x_dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
141 x_dofmap, cells[1], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
142 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
144 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
145 std::next(cdofs1.begin(), 3 * i));
148 const std::array perm{perms[
cells[0] * num_cell_facets + local_facet[0]],
149 perms[
cells[1] * num_cell_facets + local_facet[1]]};
150 fn(&value, coeffs.data() + index / 2 * cstride, constants.data(),
151 coordinate_dofs.data(), local_facet.data(), perm.data());
158template <dolfinx::scalar T, std::
floating_po
int U>
160 const fem::Form<T, U>& M, mdspan2_t x_dofmap,
161 std::span<
const scalar_value_type_t<T>> x, std::span<const T> constants,
162 const std::map<std::pair<IntegralType, int>,
163 std::pair<std::span<const T>,
int>>& coefficients)
165 std::shared_ptr<const mesh::Mesh<U>> mesh = M.mesh();
175 value += impl::assemble_cells(x_dofmap, x, cells, fn, constants, coeffs,
183 auto& [coeffs, cstride]
185 value += impl::assemble_exterior_facets(
192 mesh->topology_mutable()->create_entity_permutations();
193 const std::vector<std::uint8_t>& perms
194 = mesh->topology()->get_facet_permutations();
198 const std::vector<int> c_offsets = M.coefficient_offsets();
203 auto& [coeffs, cstride]
205 value += impl::assemble_interior_facets(
206 x_dofmap, x, num_cell_facets,
208 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:15
IntegralType
Type of integral.
Definition Form.h:34
@ 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