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 xtl::span<const std::int32_t>& cells,
28 const std::function<
void(T*,
const T*,
const T*,
const double*,
29 const int*,
const std::uint8_t*)>& fn,
30 const xtl::span<const T>& constants,
31 const xtl::span<const T>& coeffs,
int cstride)
38 const graph::AdjacencyList<std::int32_t>& x_dofmap = geometry.dofmap();
39 const std::size_t num_dofs_g = geometry.cmap().dim();
40 xtl::span<const double> x_g = geometry.x();
43 std::vector<double> coordinate_dofs(3 * num_dofs_g);
46 for (std::size_t index = 0; index <
cells.size(); ++index)
48 std::int32_t c =
cells[index];
51 auto x_dofs = x_dofmap.links(c);
52 for (std::size_t i = 0; i < x_dofs.size(); ++i)
54 common::impl::copy_N<3>(std::next(x_g.begin(), 3 * x_dofs[i]),
55 std::next(coordinate_dofs.begin(), 3 * i));
58 const T* coeff_cell = coeffs.data() + index * cstride;
59 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(),
nullptr,
68 T assemble_exterior_facets(
69 const mesh::Mesh& mesh,
70 const xtl::span<
const std::pair<std::int32_t, int>>& facets,
71 const std::function<
void(T*,
const T*,
const T*,
const double*,
const int*,
72 const std::uint8_t*)>& fn,
73 const xtl::span<const T>& constants,
const xtl::span<const T>& coeffs,
81 const graph::AdjacencyList<std::int32_t>& x_dofmap = mesh.geometry().dofmap();
82 const std::size_t num_dofs_g = mesh.geometry().cmap().dim();
83 xtl::span<const double> x_g = mesh.geometry().x();
86 std::vector<double> coordinate_dofs(3 * num_dofs_g);
89 for (std::size_t index = 0; index < facets.size(); ++index)
91 std::int32_t
cell = facets[index].first;
92 int local_facet = facets[index].second;
95 auto x_dofs = x_dofmap.links(cell);
96 for (std::size_t i = 0; i < x_dofs.size(); ++i)
98 common::impl::copy_N<3>(std::next(x_g.begin(), 3 * x_dofs[i]),
99 std::next(coordinate_dofs.begin(), 3 * i));
102 const T* coeff_cell = coeffs.data() + index * cstride;
103 fn(&value, coeff_cell, constants.data(), coordinate_dofs.data(),
104 &local_facet,
nullptr);
111 template <
typename T>
112 T assemble_interior_facets(
113 const mesh::Mesh& mesh,
114 const xtl::span<
const std::tuple<std::int32_t, int, std::int32_t, int>>&
116 const std::function<
void(T*,
const T*,
const T*,
const double*,
const int*,
117 const std::uint8_t*)>& fn,
118 const xtl::span<const T>& constants,
const xtl::span<const T>& coeffs,
119 int cstride,
const xtl::span<const int>& offsets,
120 const xtl::span<const std::uint8_t>& perms)
126 const int tdim = mesh.topology().dim();
129 const graph::AdjacencyList<std::int32_t>& x_dofmap = mesh.geometry().dofmap();
130 const std::size_t num_dofs_g = mesh.geometry().cmap().dim();
131 xtl::span<const double> x_g = mesh.geometry().x();
134 xt::xtensor<double, 3> coordinate_dofs({2, num_dofs_g, 3});
135 std::vector<T> coeff_array(2 * offsets.back());
136 assert(offsets.back() == cstride);
138 const int num_cell_facets
142 for (std::size_t index = 0; index < facets.size(); ++index)
144 const std::array<std::int32_t, 2>
cells
145 = {std::get<0>(facets[index]), std::get<2>(facets[index])};
146 const std::array<int, 2> local_facet
147 = {std::get<1>(facets[index]), std::get<3>(facets[index])};
150 auto x_dofs0 = x_dofmap.links(cells[0]);
151 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
153 common::impl::copy_N<3>(
154 std::next(x_g.begin(), 3 * x_dofs0[i]),
155 xt::view(coordinate_dofs, 0, i, xt::all()).begin());
157 auto x_dofs1 = x_dofmap.links(cells[1]);
158 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
160 common::impl::copy_N<3>(
161 std::next(x_g.begin(), 3 * x_dofs1[i]),
162 xt::view(coordinate_dofs, 1, i, xt::all()).begin());
165 const std::array perm{perms[
cells[0] * num_cell_facets + local_facet[0]],
166 perms[
cells[1] * num_cell_facets + local_facet[1]]};
167 fn(&value, coeffs.data() + index * 2 * cstride, constants.data(),
168 coordinate_dofs.data(), local_facet.data(), perm.data());
175 template <
typename T>
177 const fem::Form<T>& M,
const xtl::span<const T>& constants,
178 const std::map<std::pair<IntegralType, int>,
179 std::pair<xtl::span<const T>,
int>>& coefficients)
181 std::shared_ptr<const mesh::Mesh> mesh = M.mesh();
189 const std::vector<std::int32_t>&
cells = M.cell_domains(i);
190 value += impl::assemble_cells(mesh->geometry(), cells, fn, constants,
197 const auto& [coeffs, cstride]
199 const std::vector<std::pair<std::int32_t, int>>& facets
200 = M.exterior_facet_domains(i);
201 value += impl::assemble_exterior_facets(*mesh, facets, fn, constants,
207 mesh->topology_mutable().create_entity_permutations();
209 const std::vector<std::uint8_t>& perms
210 = mesh->topology().get_facet_permutations();
212 const std::vector<int> c_offsets = M.coefficient_offsets();
216 const auto& [coeffs, cstride]
218 const std::vector<std::tuple<std::int32_t, int, std::int32_t, int>>&
220 = 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 fem::DofMap >, 2 > &dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition: sparsitybuild.cpp:18
@ interior_facet
Interior facet.
@ exterior_facet
Exterior facet.
T assemble_scalar(const Form< T > &M, const xtl::span< const T > &constants, const std::map< std::pair< IntegralType, int >, std::pair< xtl::span< const T >, int >> &coefficients)
Assemble functional into scalar. The caller supplies the form constants and coefficients for this ver...
Definition: assembler.h:55
int cell_num_entities(CellType type, int dim)
Number of entities of dimension dim.
Definition: cell_types.cpp:185