24namespace dolfinx::fem::impl
27using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
29 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
32template <dolfinx::scalar T>
35 std::span<
const scalar_value_type_t<T>> x,
36 std::span<const std::int32_t> cells,
37 const std::function<
void(
const std::span<T>&,
38 const std::span<const std::uint32_t>&,
39 std::int32_t,
int)>& dof_transform,
40 mdspan2_t dofmap0,
int bs0,
41 const std::function<
void(
const std::span<T>&,
42 const std::span<const std::uint32_t>&,
43 std::int32_t,
int)>& dof_transform_to_transpose,
44 mdspan2_t dofmap1,
int bs1, std::span<const std::int8_t> bc0,
45 std::span<const std::int8_t> bc1,
FEkernel<T> auto kernel,
46 std::span<const T> coeffs,
int cstride, std::span<const T> constants,
47 std::span<const std::uint32_t> cell_info)
53 const int num_dofs0 = dofmap0.extent(1);
54 const int num_dofs1 = dofmap1.extent(1);
55 const int ndim0 = bs0 * num_dofs0;
56 const int ndim1 = bs1 * num_dofs1;
57 std::vector<T> Ae(ndim0 * ndim1);
59 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
62 for (std::size_t index = 0; index < cells.size(); ++index)
64 std::int32_t c = cells[index];
68 = MDSPAN_IMPL_STANDARD_NAMESPACE::MDSPAN_IMPL_PROPOSED_NAMESPACE::
69 submdspan(x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
70 for (std::size_t i = 0; i < x_dofs.size(); ++i)
72 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
73 std::next(coordinate_dofs.begin(), 3 * i));
77 std::fill(Ae.begin(), Ae.end(), 0);
78 kernel(Ae.data(), coeffs.data() + index * cstride, constants.data(),
79 coordinate_dofs.data(),
nullptr,
nullptr);
81 dof_transform(_Ae, cell_info, c, ndim1);
82 dof_transform_to_transpose(_Ae, cell_info, c, ndim0);
85 auto dofs0 = std::span(dofmap0.data_handle() + c * num_dofs0, num_dofs0);
86 auto dofs1 = std::span(dofmap1.data_handle() + c * num_dofs1, num_dofs1);
90 for (
int i = 0; i < num_dofs0; ++i)
92 for (
int k = 0; k < bs0; ++k)
94 if (bc0[bs0 * dofs0[i] + k])
97 const int row = bs0 * i + k;
98 std::fill_n(std::next(Ae.begin(), ndim1 * row), ndim1, 0.0);
106 for (
int j = 0; j < num_dofs1; ++j)
108 for (
int k = 0; k < bs1; ++k)
110 if (bc1[bs1 * dofs1[j] + k])
113 const int col = bs1 * j + k;
114 for (
int row = 0; row < ndim0; ++row)
115 Ae[row * ndim1 + col] = 0.0;
121 mat_set(dofs0, dofs1, Ae);
126template <dolfinx::scalar T>
127void assemble_exterior_facets(
129 std::span<
const scalar_value_type_t<T>> x,
130 std::span<const std::int32_t> facets,
131 const std::function<
void(
const std::span<T>&,
132 const std::span<const std::uint32_t>&,
133 std::int32_t,
int)>& dof_transform,
134 mdspan2_t dofmap0,
int bs0,
135 const std::function<
void(
const std::span<T>&,
136 const std::span<const std::uint32_t>&,
137 std::int32_t,
int)>& dof_transform_to_transpose,
138 mdspan2_t dofmap1,
int bs1, std::span<const std::int8_t> bc0,
139 std::span<const std::int8_t> bc1,
FEkernel<T> auto kernel,
140 std::span<const T> coeffs,
int cstride, std::span<const T> constants,
141 std::span<const std::uint32_t> cell_info)
147 std::vector<scalar_value_type_t<T>> coordinate_dofs(3 * x_dofmap.extent(1));
148 const int num_dofs0 = dofmap0.extent(1);
149 const int num_dofs1 = dofmap1.extent(1);
150 const int ndim0 = bs0 * num_dofs0;
151 const int ndim1 = bs1 * num_dofs1;
152 std::vector<T> Ae(ndim0 * ndim1);
153 std::span<T> _Ae(Ae);
154 assert(facets.size() % 2 == 0);
155 for (std::size_t index = 0; index < facets.size(); index += 2)
157 std::int32_t
cell = facets[index];
158 std::int32_t local_facet = facets[index + 1];
161 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::
162 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
163 x_dofmap,
cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
164 for (std::size_t i = 0; i < x_dofs.size(); ++i)
166 std::copy_n(std::next(x.begin(), 3 * x_dofs[i]), 3,
167 std::next(coordinate_dofs.begin(), 3 * i));
171 std::fill(Ae.begin(), Ae.end(), 0);
172 kernel(Ae.data(), coeffs.data() + index / 2 * cstride, constants.data(),
173 coordinate_dofs.data(), &local_facet,
nullptr);
175 dof_transform(_Ae, cell_info,
cell, ndim1);
176 dof_transform_to_transpose(_Ae, cell_info,
cell, ndim0);
179 auto dofs0 = std::span(dofmap0.data_handle() +
cell * num_dofs0, num_dofs0);
180 auto dofs1 = std::span(dofmap1.data_handle() +
cell * num_dofs1, num_dofs1);
183 for (
int i = 0; i < num_dofs0; ++i)
185 for (
int k = 0; k < bs0; ++k)
187 if (bc0[bs0 * dofs0[i] + k])
190 const int row = bs0 * i + k;
191 std::fill_n(std::next(Ae.begin(), ndim1 * row), ndim1, 0.0);
198 for (
int j = 0; j < num_dofs1; ++j)
200 for (
int k = 0; k < bs1; ++k)
202 if (bc1[bs1 * dofs1[j] + k])
205 const int col = bs1 * j + k;
206 for (
int row = 0; row < ndim0; ++row)
207 Ae[row * ndim1 + col] = 0.0;
213 mat_set(dofs0, dofs1, Ae);
218template <dolfinx::scalar T>
219void assemble_interior_facets(
221 std::span<
const scalar_value_type_t<T>> x,
int num_cell_facets,
222 std::span<const std::int32_t> facets,
223 const std::function<
void(
const std::span<T>&,
224 const std::span<const std::uint32_t>&,
225 std::int32_t,
int)>& dof_transform,
226 const DofMap& dofmap0,
int bs0,
227 const std::function<
void(
const std::span<T>&,
228 const std::span<const std::uint32_t>&,
229 std::int32_t,
int)>& dof_transform_to_transpose,
230 const DofMap& dofmap1,
int bs1, std::span<const std::int8_t> bc0,
231 std::span<const std::int8_t> bc1,
FEkernel<T> auto kernel,
232 std::span<const T> coeffs,
int cstride, std::span<const int> offsets,
233 std::span<const T> constants, std::span<const std::uint32_t> cell_info,
234 const std::function<std::uint8_t(std::size_t)>& get_perm)
240 using X = scalar_value_type_t<T>;
241 std::vector<X> coordinate_dofs(2 * x_dofmap.extent(1) * 3);
242 std::span<X> cdofs0(coordinate_dofs.data(), x_dofmap.extent(1) * 3);
243 std::span<X> cdofs1(coordinate_dofs.data() + x_dofmap.extent(1) * 3,
244 x_dofmap.extent(1) * 3);
246 std::vector<T> Ae, be;
247 std::vector<T> coeff_array(2 * offsets.back());
248 assert(offsets.back() == cstride);
251 std::vector<std::int32_t> dmapjoint0, dmapjoint1;
252 assert(facets.size() % 4 == 0);
253 for (std::size_t index = 0; index < facets.size(); index += 4)
255 std::array<std::int32_t, 2> cells = {facets[index], facets[index + 2]};
256 std::array<std::int32_t, 2> local_facet
257 = {facets[index + 1], facets[index + 3]};
260 auto x_dofs0 = MDSPAN_IMPL_STANDARD_NAMESPACE::
261 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
262 x_dofmap, cells[0], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
263 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
265 std::copy_n(std::next(x.begin(), 3 * x_dofs0[i]), 3,
266 std::next(cdofs0.begin(), 3 * i));
268 auto x_dofs1 = MDSPAN_IMPL_STANDARD_NAMESPACE::
269 MDSPAN_IMPL_PROPOSED_NAMESPACE::submdspan(
270 x_dofmap, cells[1], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
271 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
273 std::copy_n(std::next(x.begin(), 3 * x_dofs1[i]), 3,
274 std::next(cdofs1.begin(), 3 * i));
278 std::span<const std::int32_t> dmap0_cell0 = dofmap0.cell_dofs(cells[0]);
279 std::span<const std::int32_t> dmap0_cell1 = dofmap0.cell_dofs(cells[1]);
280 dmapjoint0.resize(dmap0_cell0.size() + dmap0_cell1.size());
281 std::copy(dmap0_cell0.begin(), dmap0_cell0.end(), dmapjoint0.begin());
282 std::copy(dmap0_cell1.begin(), dmap0_cell1.end(),
283 std::next(dmapjoint0.begin(), dmap0_cell0.size()));
285 std::span<const std::int32_t> dmap1_cell0 = dofmap1.cell_dofs(cells[0]);
286 std::span<const std::int32_t> dmap1_cell1 = dofmap1.cell_dofs(cells[1]);
287 dmapjoint1.resize(dmap1_cell0.size() + dmap1_cell1.size());
288 std::copy(dmap1_cell0.begin(), dmap1_cell0.end(), dmapjoint1.begin());
289 std::copy(dmap1_cell1.begin(), dmap1_cell1.end(),
290 std::next(dmapjoint1.begin(), dmap1_cell0.size()));
292 const int num_rows = bs0 * dmapjoint0.size();
293 const int num_cols = bs1 * dmapjoint1.size();
296 Ae.resize(num_rows * num_cols);
297 std::fill(Ae.begin(), Ae.end(), 0);
299 const std::array perm{
300 get_perm(cells[0] * num_cell_facets + local_facet[0]),
301 get_perm(cells[1] * num_cell_facets + local_facet[1])};
302 kernel(Ae.data(), coeffs.data() + index / 2 * cstride, constants.data(),
303 coordinate_dofs.data(), local_facet.data(), perm.data());
305 std::span<T> _Ae(Ae);
306 std::span<T> sub_Ae0 = _Ae.subspan(bs0 * dmap0_cell0.size() * num_cols,
307 bs0 * dmap0_cell1.size() * num_cols);
309 = _Ae.subspan(bs1 * dmap1_cell0.size(),
310 num_rows * num_cols - bs1 * dmap1_cell0.size());
317 dof_transform(_Ae, cell_info, cells[0], num_cols);
318 dof_transform(sub_Ae0, cell_info, cells[1], num_cols);
319 dof_transform_to_transpose(_Ae, cell_info, cells[0], num_rows);
320 dof_transform_to_transpose(sub_Ae1, cell_info, cells[1], num_rows);
325 for (std::size_t i = 0; i < dmapjoint0.size(); ++i)
327 for (
int k = 0; k < bs0; ++k)
329 if (bc0[bs0 * dmapjoint0[i] + k])
332 std::fill_n(std::next(Ae.begin(), num_cols * (bs0 * i + k)),
340 for (std::size_t j = 0; j < dmapjoint1.size(); ++j)
342 for (
int k = 0; k < bs1; ++k)
344 if (bc1[bs1 * dmapjoint1[j] + k])
347 for (
int m = 0; m < num_rows; ++m)
348 Ae[m * num_cols + bs1 * j + k] = 0.0;
354 mat_set(dmapjoint0, dmapjoint1, Ae);
363template <dolfinx::scalar T, std::
floating_po
int U>
366 std::span<
const scalar_value_type_t<T>> x, std::span<const T> constants,
367 const std::map<std::pair<IntegralType, int>,
368 std::pair<std::span<const T>,
int>>& coefficients,
369 std::span<const std::int8_t> bc0, std::span<const std::int8_t> bc1)
371 std::shared_ptr<const mesh::Mesh<U>> mesh = a.
mesh();
375 std::shared_ptr<const fem::DofMap> dofmap0
377 std::shared_ptr<const fem::DofMap> dofmap1
381 auto dofs0 = dofmap0->map();
382 const int bs0 = dofmap0->bs();
383 auto dofs1 = dofmap1->map();
384 const int bs1 = dofmap1->bs();
390 const std::function<void(
const std::span<T>&,
391 const std::span<const std::uint32_t>&, std::int32_t,
393 = element0->template get_dof_transformation_function<T>();
394 const std::function<void(
const std::span<T>&,
395 const std::span<const std::uint32_t>&, std::int32_t,
396 int)>& dof_transform_to_transpose
397 = element1->template get_dof_transformation_to_transpose_function<T>();
399 const bool needs_transformation_data
400 = element0->needs_dof_transformations()
401 or element1->needs_dof_transformations()
403 std::span<const std::uint32_t> cell_info;
404 if (needs_transformation_data)
406 mesh->topology_mutable()->create_entity_permutations();
407 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
416 dof_transform, dofs0, bs0, dof_transform_to_transpose,
417 dofs1, bs1, bc0, bc1, fn, coeffs, cstride, constants,
425 auto& [coeffs, cstride]
427 impl::assemble_exterior_facets(
429 dof_transform, dofs0, bs0, dof_transform_to_transpose, dofs1, bs1, bc0,
430 bc1, fn, coeffs, cstride, constants, cell_info);
435 std::function<std::uint8_t(std::size_t)> get_perm;
438 mesh->topology_mutable()->create_entity_permutations();
439 const std::vector<std::uint8_t>& perms
440 = mesh->topology()->get_facet_permutations();
441 get_perm = [&perms](std::size_t i) {
return perms[i]; };
444 get_perm = [](std::size_t) {
return 0; };
446 auto cell_types = mesh->topology()->cell_types();
447 if (cell_types.size() > 1)
448 throw std::runtime_error(
"Multiple cell types in the assembler.");
450 mesh->topology()->dim() - 1);
456 auto& [coeffs, cstride]
458 impl::assemble_interior_facets(
459 mat_set, x_dofmap, x, num_cell_facets,
461 bs0, dof_transform_to_transpose, *dofmap1, bs1, bc0, bc1, fn, coeffs,
462 cstride, c_offsets, constants, cell_info, get_perm);