11 #include "FunctionSpace.h"
13 #include <dolfinx/common/IndexMap.h>
14 #include <dolfinx/common/utils.h>
15 #include <dolfinx/mesh/Geometry.h>
16 #include <dolfinx/mesh/Mesh.h>
17 #include <dolfinx/mesh/Topology.h>
21 namespace dolfinx::fem::impl
26 T assemble_cells(
const mesh::Geometry& geometry,
27 const std::span<const std::int32_t>& cells,
28 const std::function<
void(T*,
const T*,
const T*,
29 const scalar_value_type_t<T>*,
30 const int*,
const std::uint8_t*)>& fn,
31 const std::span<const T>& constants,
32 const std::span<const T>& coeffs,
int cstride)
39 const graph::AdjacencyList<std::int32_t>& x_dofmap = geometry.dofmap();
40 const std::size_t num_dofs_g = geometry.cmap().dim();
41 std::span<const double> x_g = geometry.x();
44 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * num_dofs_g);
47 for (std::size_t index = 0; index <
cells.size(); ++index)
49 std::int32_t c =
cells[index];
52 auto x_dofs = x_dofmap.links(c);
53 for (std::size_t i = 0; i < x_dofs.size(); ++i)
55 common::impl::copy_N<3>(std::next(x_g.begin(), 3 * x_dofs[i]),
56 std::next(coordinate_dofs.begin(), 3 * i));
59 const T* coeff_cell = coeffs.data() + index * cstride;
60 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(),
nullptr,
69 T assemble_exterior_facets(
70 const mesh::Mesh& mesh,
const std::span<const std::int32_t>& facets,
71 const std::function<
void(T*,
const T*,
const T*,
72 const scalar_value_type_t<T>*,
const int*,
73 const std::uint8_t*)>& fn,
74 const std::span<const T>& constants,
const std::span<const T>& coeffs,
82 const graph::AdjacencyList<std::int32_t>& x_dofmap = mesh.geometry().dofmap();
83 const std::size_t num_dofs_g = mesh.geometry().cmap().dim();
84 std::span<const double> x_g = mesh.geometry().x();
87 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * num_dofs_g);
90 assert(facets.size() % 2 == 0);
91 for (std::size_t index = 0; index < facets.size(); index += 2)
93 std::int32_t
cell = facets[index];
94 std::int32_t local_facet = facets[index + 1];
97 auto x_dofs = x_dofmap.links(cell);
98 for (std::size_t i = 0; i < x_dofs.size(); ++i)
100 common::impl::copy_N<3>(std::next(x_g.begin(), 3 * x_dofs[i]),
101 std::next(coordinate_dofs.begin(), 3 * i));
104 const T* coeff_cell = coeffs.data() + index / 2 * cstride;
105 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(),
106 &local_facet,
nullptr);
113 template <
typename T>
114 T assemble_interior_facets(
115 const mesh::Mesh& mesh,
const std::span<const std::int32_t>& facets,
116 const std::function<
void(T*,
const T*,
const T*,
117 const scalar_value_type_t<T>*,
const int*,
118 const std::uint8_t*)>& fn,
119 const std::span<const T>& constants,
const std::span<const T>& coeffs,
120 int cstride,
const std::span<const int>& offsets,
121 const std::span<const std::uint8_t>& perms)
127 const int tdim = mesh.topology().dim();
130 const graph::AdjacencyList<std::int32_t>& x_dofmap = mesh.geometry().dofmap();
131 const std::size_t num_dofs_g = mesh.geometry().cmap().dim();
132 std::span<const double> x_g = mesh.geometry().x();
135 using X = scalar_value_type_t<T>;
136 std::vector<X> coordinate_dofs(2 * num_dofs_g * 3);
137 std::span<X> cdofs0(coordinate_dofs.data(), num_dofs_g * 3);
138 std::span<X> cdofs1(coordinate_dofs.data() + num_dofs_g * 3, num_dofs_g * 3);
140 std::vector<T> coeff_array(2 * offsets.back());
141 assert(offsets.back() == cstride);
143 const int num_cell_facets
147 assert(facets.size() % 4 == 0);
148 for (std::size_t index = 0; index < facets.size(); index += 4)
150 std::array<std::int32_t, 2>
cells = {facets[index], facets[index + 2]};
151 std::array<std::int32_t, 2> local_facet
152 = {facets[index + 1], facets[index + 3]};
155 auto x_dofs0 = x_dofmap.links(cells[0]);
156 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
158 common::impl::copy_N<3>(std::next(x_g.begin(), 3 * x_dofs0[i]),
159 std::next(cdofs0.begin(), 3 * i));
161 auto x_dofs1 = x_dofmap.links(cells[1]);
162 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
164 common::impl::copy_N<3>(std::next(x_g.begin(), 3 * x_dofs1[i]),
165 std::next(cdofs1.begin(), 3 * i));
168 const std::array perm{perms[
cells[0] * num_cell_facets + local_facet[0]],
169 perms[
cells[1] * num_cell_facets + local_facet[1]]};
170 fn(&value, coeffs.data() + index / 2 * cstride, constants.data(),
171 coordinate_dofs.data(), local_facet.data(), perm.data());
178 template <
typename T>
180 const fem::Form<T>& M,
const std::span<const T>& constants,
181 const std::map<std::pair<IntegralType, int>,
182 std::pair<std::span<const T>,
int>>& coefficients)
184 std::shared_ptr<const mesh::Mesh> mesh = M.mesh();
192 const std::vector<std::int32_t>&
cells = M.cell_domains(i);
193 value += impl::assemble_cells(mesh->geometry(), cells, fn, constants,
200 const auto& [coeffs, cstride]
202 const std::vector<std::int32_t>& facets = M.exterior_facet_domains(i);
203 value += impl::assemble_exterior_facets(*mesh, facets, fn, constants,
209 mesh->topology_mutable().create_entity_permutations();
211 const std::vector<std::uint8_t>& perms
212 = mesh->topology().get_facet_permutations();
214 const std::vector<int> c_offsets = M.coefficient_offsets();
218 const auto& [coeffs, cstride]
220 const std::vector<std::int32_t>& facets = M.interior_facet_domains(i);
221 value += impl::assemble_interior_facets(
222 *mesh, facets, fn, constants, coeffs, cstride, c_offsets, perms);
void cells(la::SparsityPattern &pattern, const mesh::Topology &topology, const std::array< const std::reference_wrapper< const DofMap >, 2 > &dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition: sparsitybuild.cpp:18
T assemble_scalar(const Form< T > &M, const std::span< const T > &constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int >> &coefficients)
Assemble functional into scalar. The caller supplies the form constants and coefficients for this ver...
Definition: assembler.h:55
@ interior_facet
Interior facet.
@ exterior_facet
Exterior facet.
int cell_num_entities(CellType type, int dim)
Number of entities of dimension dim.
Definition: cell_types.cpp:182