10 #include "DirichletBC.h"
13 #include "FunctionSpace.h"
15 #include <dolfinx/common/IndexMap.h>
16 #include <dolfinx/common/utils.h>
17 #include <dolfinx/graph/AdjacencyList.h>
18 #include <dolfinx/mesh/Geometry.h>
19 #include <dolfinx/mesh/Mesh.h>
20 #include <dolfinx/mesh/Topology.h>
24 #include <xtensor/xbuilder.hpp>
25 #include <xtensor/xtensor.hpp>
27 namespace dolfinx::fem::impl
39 template <
typename T,
int _bs0 = -1,
int _bs1 = -1>
41 xtl::span<T> b,
const mesh::Geometry& geometry,
42 const std::function<
void(T*,
const T*,
const T*,
const double*,
const int*,
43 const std::uint8_t*)>& kernel,
44 const xtl::span<const std::int32_t>& cells,
45 const std::function<
void(
const xtl::span<T>&,
46 const xtl::span<const std::uint32_t>&,
47 std::int32_t,
int)>& dof_transform,
48 const graph::AdjacencyList<std::int32_t>& dofmap0,
int bs0,
49 const std::function<
void(
const xtl::span<T>&,
50 const xtl::span<const std::uint32_t>&,
51 std::int32_t,
int)>& dof_transform_to_transpose,
52 const graph::AdjacencyList<std::int32_t>& dofmap1,
int bs1,
53 const xtl::span<const T>& constants,
const xtl::span<const T>& coeffs,
54 int cstride,
const xtl::span<const std::uint32_t>& cell_info,
55 const xtl::span<const T>& bc_values1,
56 const xtl::span<const std::int8_t>& bc_markers1,
57 const xtl::span<const T>& x0,
double scale)
59 assert(_bs0 < 0 or _bs0 == bs0);
60 assert(_bs1 < 0 or _bs1 == bs1);
66 const graph::AdjacencyList<std::int32_t>& x_dofmap = geometry.dofmap();
67 const std::size_t num_dofs_g = geometry.cmap().dim();
68 xtl::span<const double> x_g = geometry.x();
71 std::vector<double> coordinate_dofs(3 * num_dofs_g);
72 std::vector<T> Ae, be;
73 for (std::size_t index = 0; index <
cells.size(); ++index)
75 std::int32_t c =
cells[index];
78 auto dmap1 = dofmap1.links(c);
82 for (std::size_t j = 0; j < dmap1.size(); ++j)
84 if constexpr (_bs1 > 0)
86 for (
int k = 0; k < _bs1; ++k)
88 assert(_bs1 * dmap1[j] + k < (
int)bc_markers1.size());
89 if (bc_markers1[_bs1 * dmap1[j] + k])
98 for (
int k = 0; k < bs1; ++k)
100 assert(bs1 * dmap1[j] + k < (
int)bc_markers1.size());
101 if (bc_markers1[bs1 * dmap1[j] + k])
114 auto x_dofs = x_dofmap.links(c);
115 for (std::size_t i = 0; i < x_dofs.size(); ++i)
117 common::impl::copy_N<3>(std::next(x_g.begin(), 3 * x_dofs[i]),
118 std::next(coordinate_dofs.begin(), 3 * i));
122 auto dmap0 = dofmap0.links(c);
124 const int num_rows = bs0 * dmap0.size();
125 const int num_cols = bs1 * dmap1.size();
127 const T* coeff_array = coeffs.data() + index * cstride;
128 Ae.resize(num_rows * num_cols);
129 std::fill(Ae.begin(), Ae.end(), 0);
130 kernel(Ae.data(), coeff_array, constants.data(), coordinate_dofs.data(),
132 dof_transform(Ae, cell_info, c, num_cols);
133 dof_transform_to_transpose(Ae, cell_info, c, num_rows);
137 std::fill(be.begin(), be.end(), 0);
138 for (std::size_t j = 0; j < dmap1.size(); ++j)
140 if constexpr (_bs1 > 0)
142 for (
int k = 0; k < _bs1; ++k)
144 const std::int32_t jj = _bs1 * dmap1[j] + k;
145 assert(jj < (
int)bc_markers1.size());
148 const T bc = bc_values1[jj];
149 const T _x0 = x0.empty() ? 0.0 : x0[jj];
152 for (
int m = 0; m < num_rows; ++m)
153 be[m] -= Ae[m * num_cols + _bs1 * j + k] * scale * (bc - _x0);
159 for (
int k = 0; k < bs1; ++k)
161 const std::int32_t jj = bs1 * dmap1[j] + k;
162 assert(jj < (
int)bc_markers1.size());
165 const T bc = bc_values1[jj];
166 const T _x0 = x0.empty() ? 0.0 : x0[jj];
168 for (
int m = 0; m < num_rows; ++m)
169 be[m] -= Ae[m * num_cols + bs1 * j + k] * scale * (bc - _x0);
175 for (std::size_t i = 0; i < dmap0.size(); ++i)
177 if constexpr (_bs0 > 0)
179 for (
int k = 0; k < _bs0; ++k)
180 b[_bs0 * dmap0[i] + k] += be[_bs0 * i + k];
184 for (
int k = 0; k < bs0; ++k)
185 b[bs0 * dmap0[i] + k] += be[bs0 * i + k];
197 template <
typename T,
int _bs = -1>
198 void _lift_bc_exterior_facets(
199 xtl::span<T> b,
const mesh::Mesh& mesh,
200 const std::function<
void(T*,
const T*,
const T*,
const double*,
const int*,
201 const std::uint8_t*)>& kernel,
202 const xtl::span<
const std::pair<std::int32_t, int>>& facets,
203 const std::function<
void(
const xtl::span<T>&,
204 const xtl::span<const std::uint32_t>&,
205 std::int32_t,
int)>& dof_transform,
206 const graph::AdjacencyList<std::int32_t>& dofmap0,
int bs0,
207 const std::function<
void(
const xtl::span<T>&,
208 const xtl::span<const std::uint32_t>&,
209 std::int32_t,
int)>& dof_transform_to_transpose,
210 const graph::AdjacencyList<std::int32_t>& dofmap1,
int bs1,
211 const xtl::span<const T>& constants,
const xtl::span<const T>& coeffs,
212 int cstride,
const xtl::span<const std::uint32_t>& cell_info,
213 const xtl::span<const T>& bc_values1,
214 const xtl::span<const std::int8_t>& bc_markers1,
215 const xtl::span<const T>& x0,
double scale)
221 const graph::AdjacencyList<std::int32_t>& x_dofmap = mesh.geometry().dofmap();
223 const std::size_t num_dofs_g = mesh.geometry().cmap().dim();
225 xtl::span<const double> x_g = mesh.geometry().x();
228 std::vector<double> coordinate_dofs(3 * num_dofs_g);
229 std::vector<T> Ae, be;
231 for (std::size_t index = 0; index < facets.size(); ++index)
233 std::int32_t
cell = facets[index].first;
234 int local_facet = facets[index].second;
237 auto dmap1 = dofmap1.links(cell);
241 for (std::size_t j = 0; j < dmap1.size(); ++j)
243 for (
int k = 0; k < bs1; ++k)
245 if (bc_markers1[bs1 * dmap1[j] + k])
257 auto x_dofs = x_dofmap.links(cell);
258 for (std::size_t i = 0; i < x_dofs.size(); ++i)
260 common::impl::copy_N<3>(std::next(x_g.begin(), 3 * x_dofs[i]),
261 std::next(coordinate_dofs.begin(), 3 * i));
265 auto dmap0 = dofmap0.links(cell);
267 const int num_rows = bs0 * dmap0.size();
268 const int num_cols = bs1 * dmap1.size();
270 const T* coeff_array = coeffs.data() + index * cstride;
271 Ae.resize(num_rows * num_cols);
272 std::fill(Ae.begin(), Ae.end(), 0);
273 kernel(Ae.data(), coeff_array, constants.data(), coordinate_dofs.data(),
274 &local_facet,
nullptr);
275 dof_transform(Ae, cell_info, cell, num_cols);
276 dof_transform_to_transpose(Ae, cell_info, cell, num_rows);
280 std::fill(be.begin(), be.end(), 0);
281 for (std::size_t j = 0; j < dmap1.size(); ++j)
283 for (
int k = 0; k < bs1; ++k)
285 const std::int32_t jj = bs1 * dmap1[j] + k;
288 const T bc = bc_values1[jj];
289 const T _x0 = x0.empty() ? 0.0 : x0[jj];
291 for (
int m = 0; m < num_rows; ++m)
292 be[m] -= Ae[m * num_cols + bs1 * j + k] * scale * (bc - _x0);
297 for (std::size_t i = 0; i < dmap0.size(); ++i)
298 for (
int k = 0; k < bs0; ++k)
299 b[bs0 * dmap0[i] + k] += be[bs0 * i + k];
309 template <
typename T,
int _bs = -1>
310 void _lift_bc_interior_facets(
311 xtl::span<T> b,
const mesh::Mesh& mesh,
312 const std::function<
void(T*,
const T*,
const T*,
const double*,
const int*,
313 const std::uint8_t*)>& kernel,
314 const xtl::span<
const std::tuple<std::int32_t, int, std::int32_t, int>>&
316 const std::function<
void(
const xtl::span<T>&,
317 const xtl::span<const std::uint32_t>&,
318 std::int32_t,
int)>& dof_transform,
319 const graph::AdjacencyList<std::int32_t>& dofmap0,
int bs0,
320 const std::function<
void(
const xtl::span<T>&,
321 const xtl::span<const std::uint32_t>&,
322 std::int32_t,
int)>& dof_transform_to_transpose,
323 const graph::AdjacencyList<std::int32_t>& dofmap1,
int bs1,
324 const xtl::span<const T>& constants,
const xtl::span<const T>& coeffs,
325 int cstride,
const xtl::span<const std::uint32_t>& cell_info,
326 const std::function<std::uint8_t(std::size_t)>& get_perm,
327 const xtl::span<const T>& bc_values1,
328 const xtl::span<const std::int8_t>& bc_markers1,
329 const xtl::span<const T>& x0,
double scale)
334 const int tdim = mesh.topology().dim();
337 const graph::AdjacencyList<std::int32_t>& x_dofmap = mesh.geometry().dofmap();
338 const std::size_t num_dofs_g = mesh.geometry().cmap().dim();
339 xtl::span<const double> x_g = mesh.geometry().x();
341 const int num_cell_facets
345 xt::xtensor<double, 3> coordinate_dofs({2, num_dofs_g, 3});
346 std::vector<T> Ae, be;
349 std::vector<std::int32_t> dmapjoint0, dmapjoint1;
351 for (std::size_t index = 0; index < facets.size(); ++index)
353 const std::array<std::int32_t, 2>
cells
354 = {std::get<0>(facets[index]), std::get<2>(facets[index])};
355 const std::array<int, 2> local_facet
356 = {std::get<1>(facets[index]), std::get<3>(facets[index])};
359 auto x_dofs0 = x_dofmap.links(cells[0]);
360 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
362 common::impl::copy_N<3>(
363 std::next(x_g.begin(), 3 * x_dofs0[i]),
364 xt::view(coordinate_dofs, 0, i, xt::all()).begin());
366 auto x_dofs1 = x_dofmap.links(cells[1]);
367 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
369 common::impl::copy_N<3>(
370 std::next(x_g.begin(), 3 * x_dofs1[i]),
371 xt::view(coordinate_dofs, 1, i, xt::all()).begin());
375 const xtl::span<const std::int32_t> dmap0_cell0 = dofmap0.links(cells[0]);
376 const xtl::span<const std::int32_t> dmap0_cell1 = dofmap0.links(cells[1]);
377 dmapjoint0.resize(dmap0_cell0.size() + dmap0_cell1.size());
378 std::copy(dmap0_cell0.begin(), dmap0_cell0.end(), dmapjoint0.begin());
379 std::copy(dmap0_cell1.begin(), dmap0_cell1.end(),
380 std::next(dmapjoint0.begin(), dmap0_cell0.size()));
382 const xtl::span<const std::int32_t> dmap1_cell0 = dofmap1.links(cells[0]);
383 const xtl::span<const std::int32_t> dmap1_cell1 = dofmap1.links(cells[1]);
384 dmapjoint1.resize(dmap1_cell0.size() + dmap1_cell1.size());
385 std::copy(dmap1_cell0.begin(), dmap1_cell0.end(), dmapjoint1.begin());
386 std::copy(dmap1_cell1.begin(), dmap1_cell1.end(),
387 std::next(dmapjoint1.begin(), dmap1_cell0.size()));
391 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
393 for (
int k = 0; k < bs1; ++k)
395 if (bc_markers1[bs1 * dmap1_cell0[j] + k])
404 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
406 for (
int k = 0; k < bs1; ++k)
408 if (bc_markers1[bs1 * dmap1_cell1[j] + k])
419 const int num_rows = bs0 * dmapjoint0.size();
420 const int num_cols = bs1 * dmapjoint1.size();
423 Ae.resize(num_rows * num_cols);
424 std::fill(Ae.begin(), Ae.end(), 0);
425 const std::array perm{
426 get_perm(cells[0] * num_cell_facets + local_facet[0]),
427 get_perm(cells[1] * num_cell_facets + local_facet[1])};
428 kernel(Ae.data(), coeffs.data() + index * 2 * cstride, constants.data(),
429 coordinate_dofs.data(), local_facet.data(), perm.data());
431 const xtl::span<T> _Ae(Ae);
433 const xtl::span<T> sub_Ae0
434 = _Ae.subspan(bs0 * dmap0_cell0.size() * num_cols,
435 bs0 * dmap0_cell1.size() * num_cols);
436 const xtl::span<T> sub_Ae1
437 = _Ae.subspan(bs1 * dmap1_cell0.size(),
438 num_rows * num_cols - bs1 * dmap1_cell0.size());
440 dof_transform(_Ae, cell_info, cells[0], num_cols);
441 dof_transform(sub_Ae0, cell_info, cells[1], num_cols);
442 dof_transform_to_transpose(_Ae, cell_info, cells[0], num_rows);
443 dof_transform_to_transpose(sub_Ae1, cell_info, cells[1], num_rows);
446 std::fill(be.begin(), be.end(), 0);
449 for (std::size_t j = 0; j < dmap1_cell0.size(); ++j)
451 for (
int k = 0; k < bs1; ++k)
453 const std::int32_t jj = bs1 * dmap1_cell0[j] + k;
456 const T bc = bc_values1[jj];
457 const T _x0 = x0.empty() ? 0.0 : x0[jj];
459 for (
int m = 0; m < num_rows; ++m)
460 be[m] -= Ae[m * num_cols + bs1 * j + k] * scale * (bc - _x0);
466 const int offset = bs1 * dmap1_cell0.size();
467 for (std::size_t j = 0; j < dmap1_cell1.size(); ++j)
469 for (
int k = 0; k < bs1; ++k)
471 const std::int32_t jj = bs1 * dmap1_cell1[j] + k;
474 const T bc = bc_values1[jj];
475 const T _x0 = x0.empty() ? 0.0 : x0[jj];
477 for (
int m = 0; m < num_rows; ++m)
480 -= Ae[m * num_cols + offset + bs1 * j + k] * scale * (bc - _x0);
486 for (std::size_t i = 0; i < dmap0_cell0.size(); ++i)
487 for (
int k = 0; k < bs0; ++k)
488 b[bs0 * dmap0_cell0[i] + k] += be[bs0 * i + k];
490 const int offset_be = bs0 * dmap0_cell0.size();
491 for (std::size_t i = 0; i < dmap0_cell1.size(); ++i)
492 for (
int k = 0; k < bs0; ++k)
493 b[bs0 * dmap0_cell1[i] + k] += be[offset_be + bs0 * i + k];
502 template <
typename T,
int _bs = -1>
504 const std::function<
void(
const xtl::span<T>&,
505 const xtl::span<const std::uint32_t>&,
506 std::int32_t,
int)>& dof_transform,
507 xtl::span<T> b,
const mesh::Geometry& geometry,
508 const xtl::span<const std::int32_t>& cells,
509 const graph::AdjacencyList<std::int32_t>& dofmap,
int bs,
510 const std::function<
void(T*,
const T*,
const T*,
const double*,
const int*,
511 const std::uint8_t*)>& kernel,
512 const xtl::span<const T>& constants,
const xtl::span<const T>& coeffs,
513 int cstride,
const xtl::span<const std::uint32_t>& cell_info)
515 assert(_bs < 0 or _bs == bs);
521 const graph::AdjacencyList<std::int32_t>& x_dofmap = geometry.dofmap();
522 const std::size_t num_dofs_g = geometry.cmap().dim();
523 xtl::span<const double> x_g = geometry.x();
527 const int num_dofs = dofmap.links(0).size();
528 std::vector<double> coordinate_dofs(3 * num_dofs_g);
529 std::vector<T> be(bs * num_dofs);
530 const xtl::span<T> _be(be);
533 for (std::size_t index = 0; index <
cells.size(); ++index)
535 std::int32_t c =
cells[index];
538 auto x_dofs = x_dofmap.links(c);
539 for (std::size_t i = 0; i < x_dofs.size(); ++i)
541 common::impl::copy_N<3>(std::next(x_g.begin(), 3 * x_dofs[i]),
542 std::next(coordinate_dofs.begin(), 3 * i));
546 std::fill(be.begin(), be.end(), 0);
547 kernel(be.data(), coeffs.data() + index * cstride, constants.data(),
548 coordinate_dofs.data(),
nullptr,
nullptr);
549 dof_transform(_be, cell_info, c, 1);
552 auto dofs = dofmap.links(c);
553 if constexpr (_bs > 0)
555 for (
int i = 0; i < num_dofs; ++i)
556 for (
int k = 0; k < _bs; ++k)
557 b[_bs * dofs[i] + k] += be[_bs * i + k];
561 for (
int i = 0; i < num_dofs; ++i)
562 for (
int k = 0; k < bs; ++k)
563 b[bs * dofs[i] + k] += be[bs * i + k];
574 template <
typename T,
int _bs = -1>
575 void assemble_exterior_facets(
576 const std::function<
void(
const xtl::span<T>&,
577 const xtl::span<const std::uint32_t>&,
578 std::int32_t,
int)>& dof_transform,
579 xtl::span<T> b,
const mesh::Mesh& mesh,
580 const xtl::span<
const std::pair<std::int32_t, int>>& facets,
581 const graph::AdjacencyList<std::int32_t>& dofmap,
int bs,
582 const std::function<
void(T*,
const T*,
const T*,
const double*,
const int*,
583 const std::uint8_t*)>& fn,
584 const xtl::span<const T>& constants,
const xtl::span<const T>& coeffs,
585 int cstride,
const xtl::span<const std::uint32_t>& cell_info)
587 assert(_bs < 0 or _bs == bs);
593 const graph::AdjacencyList<std::int32_t>& x_dofmap = mesh.geometry().dofmap();
594 const std::size_t num_dofs_g = mesh.geometry().cmap().dim();
595 xtl::span<const double> x_g = mesh.geometry().x();
599 const int num_dofs = dofmap.links(0).size();
600 std::vector<double> coordinate_dofs(3 * num_dofs_g);
601 std::vector<T> be(bs * num_dofs);
602 const xtl::span<T> _be(be);
604 for (std::size_t index = 0; index < facets.size(); ++index)
606 std::int32_t
cell = facets[index].first;
607 int local_facet = facets[index].second;
610 auto x_dofs = x_dofmap.links(cell);
611 for (std::size_t i = 0; i < x_dofs.size(); ++i)
613 common::impl::copy_N<3>(std::next(x_g.begin(), 3 * x_dofs[i]),
614 std::next(coordinate_dofs.begin(), 3 * i));
618 std::fill(be.begin(), be.end(), 0);
619 fn(be.data(), coeffs.data() + index * cstride, constants.data(),
620 coordinate_dofs.data(), &local_facet,
nullptr);
622 dof_transform(_be, cell_info, cell, 1);
625 auto dofs = dofmap.links(cell);
626 if constexpr (_bs > 0)
628 for (
int i = 0; i < num_dofs; ++i)
629 for (
int k = 0; k < _bs; ++k)
630 b[_bs * dofs[i] + k] += be[_bs * i + k];
634 for (
int i = 0; i < num_dofs; ++i)
635 for (
int k = 0; k < bs; ++k)
636 b[bs * dofs[i] + k] += be[bs * i + k];
647 template <
typename T,
int _bs = -1>
648 void assemble_interior_facets(
649 const std::function<
void(
const xtl::span<T>&,
650 const xtl::span<const std::uint32_t>&,
651 std::int32_t,
int)>& dof_transform,
652 xtl::span<T> b,
const mesh::Mesh& mesh,
653 const xtl::span<
const std::tuple<std::int32_t, int, std::int32_t, int>>&
655 const fem::DofMap& dofmap,
656 const std::function<
void(T*,
const T*,
const T*,
const double*,
const int*,
657 const std::uint8_t*)>& fn,
658 const xtl::span<const T>& constants,
const xtl::span<const T>& coeffs,
659 int cstride,
const xtl::span<const std::uint32_t>& cell_info,
660 const std::function<std::uint8_t(std::size_t)>& get_perm)
662 const int tdim = mesh.topology().dim();
665 const graph::AdjacencyList<std::int32_t>& x_dofmap = mesh.geometry().dofmap();
666 const std::size_t num_dofs_g = mesh.geometry().cmap().dim();
667 xtl::span<const double> x_g = mesh.geometry().x();
670 xt::xtensor<double, 3> coordinate_dofs({2, num_dofs_g, 3});
673 const int num_cell_facets
676 const int bs = dofmap.bs();
677 assert(_bs < 0 or _bs == bs);
678 for (std::size_t index = 0; index < facets.size(); ++index)
680 const std::array<std::int32_t, 2>
cells
681 = {std::get<0>(facets[index]), std::get<2>(facets[index])};
682 const std::array<int, 2> local_facet
683 = {std::get<1>(facets[index]), std::get<3>(facets[index])};
686 auto x_dofs0 = x_dofmap.links(cells[0]);
687 for (std::size_t i = 0; i < x_dofs0.size(); ++i)
689 common::impl::copy_N<3>(
690 std::next(x_g.begin(), 3 * x_dofs0[i]),
691 xt::view(coordinate_dofs, 0, i, xt::all()).begin());
693 auto x_dofs1 = x_dofmap.links(cells[1]);
694 for (std::size_t i = 0; i < x_dofs1.size(); ++i)
696 common::impl::copy_N<3>(
697 std::next(x_g.begin(), 3 * x_dofs1[i]),
698 xt::view(coordinate_dofs, 1, i, xt::all()).begin());
702 xtl::span<const std::int32_t> dmap0 = dofmap.cell_dofs(cells[0]);
703 xtl::span<const std::int32_t> dmap1 = dofmap.cell_dofs(cells[1]);
706 be.resize(bs * (dmap0.size() + dmap1.size()));
707 std::fill(be.begin(), be.end(), 0);
708 const std::array perm{
709 get_perm(cells[0] * num_cell_facets + local_facet[0]),
710 get_perm(cells[1] * num_cell_facets + local_facet[1])};
711 fn(be.data(), coeffs.data() + index * 2 * cstride, constants.data(),
712 coordinate_dofs.data(), local_facet.data(), perm.data());
714 const xtl::span<T> _be(be);
715 const xtl::span<T> sub_be
716 = _be.subspan(bs * dmap0.size(), bs * dmap1.size());
718 dof_transform(be, cell_info, cells[0], 1);
719 dof_transform(sub_be, cell_info, cells[1], 1);
722 if constexpr (_bs > 0)
724 for (std::size_t i = 0; i < dmap0.size(); ++i)
725 for (
int k = 0; k < _bs; ++k)
726 b[_bs * dmap0[i] + k] += be[_bs * i + k];
727 for (std::size_t i = 0; i < dmap1.size(); ++i)
728 for (
int k = 0; k < _bs; ++k)
729 b[_bs * dmap1[i] + k] += be[_bs * (i + dmap0.size()) + k];
733 for (std::size_t i = 0; i < dmap0.size(); ++i)
734 for (
int k = 0; k < bs; ++k)
735 b[bs * dmap0[i] + k] += be[bs * i + k];
736 for (std::size_t i = 0; i < dmap1.size(); ++i)
737 for (
int k = 0; k < bs; ++k)
738 b[bs * dmap1[i] + k] += be[bs * (i + dmap0.size()) + k];
757 template <
typename T>
758 void lift_bc(xtl::span<T> b,
const Form<T>& a,
759 const xtl::span<const T>& constants,
760 const std::map<std::pair<IntegralType, int>,
761 std::pair<xtl::span<const T>,
int>>& coefficients,
762 const xtl::span<const T>& bc_values1,
763 const xtl::span<const std::int8_t>& bc_markers1,
764 const xtl::span<const T>& x0,
double scale)
766 std::shared_ptr<const mesh::Mesh> mesh = a.mesh();
770 assert(a.function_spaces().at(0));
771 assert(a.function_spaces().at(1));
772 const graph::AdjacencyList<std::int32_t>& dofmap0
773 = a.function_spaces()[0]->dofmap()->list();
774 const int bs0 = a.function_spaces()[0]->dofmap()->bs();
775 std::shared_ptr<const fem::FiniteElement> element0
776 = a.function_spaces()[0]->element();
777 const graph::AdjacencyList<std::int32_t>& dofmap1
778 = a.function_spaces()[1]->dofmap()->list();
779 const int bs1 = a.function_spaces()[1]->dofmap()->bs();
780 std::shared_ptr<const fem::FiniteElement> element1
781 = a.function_spaces()[1]->element();
783 const bool needs_transformation_data
784 = element0->needs_dof_transformations()
785 or element1->needs_dof_transformations()
786 or a.needs_facet_permutations();
788 xtl::span<const std::uint32_t> cell_info;
789 if (needs_transformation_data)
791 mesh->topology_mutable().create_entity_permutations();
792 cell_info = xtl::span(mesh->topology().get_cell_permutation_info());
795 const std::function<void(
const xtl::span<T>&,
796 const xtl::span<const std::uint32_t>&, std::int32_t,
798 dof_transform = element0->get_dof_transformation_function<T>();
799 const std::function<void(
const xtl::span<T>&,
800 const xtl::span<const std::uint32_t>&, std::int32_t,
802 dof_transform_to_transpose
803 = element1->get_dof_transformation_to_transpose_function<T>();
809 const std::vector<std::int32_t>&
cells = a.cell_domains(i);
810 if (bs0 == 1 and bs1 == 1)
812 _lift_bc_cells<T, 1, 1>(b, mesh->geometry(), kernel, cells, dof_transform,
813 dofmap0, bs0, dof_transform_to_transpose, dofmap1,
814 bs1, constants, coeffs, cstride, cell_info,
815 bc_values1, bc_markers1, x0, scale);
817 else if (bs0 == 3 and bs1 == 3)
819 _lift_bc_cells<T, 3, 3>(b, mesh->geometry(), kernel, cells, dof_transform,
820 dofmap0, bs0, dof_transform_to_transpose, dofmap1,
821 bs1, constants, coeffs, cstride, cell_info,
822 bc_values1, bc_markers1, x0, scale);
826 _lift_bc_cells(b, mesh->geometry(), kernel, cells, dof_transform, dofmap0,
827 bs0, dof_transform_to_transpose, dofmap1, bs1, constants,
828 coeffs, cstride, cell_info, bc_values1, bc_markers1, x0,
836 const auto& [coeffs, cstride]
838 const std::vector<std::pair<std::int32_t, int>>& facets
839 = a.exterior_facet_domains(i);
840 _lift_bc_exterior_facets(b, *mesh, kernel, facets, dof_transform, dofmap0,
841 bs0, dof_transform_to_transpose, dofmap1, bs1,
842 constants, coeffs, cstride, cell_info, bc_values1,
843 bc_markers1, x0, scale);
848 std::function<std::uint8_t(std::size_t)> get_perm;
849 if (a.needs_facet_permutations())
851 mesh->topology_mutable().create_entity_permutations();
852 const std::vector<std::uint8_t>& perms
853 = mesh->topology().get_facet_permutations();
854 get_perm = [&perms](std::size_t i) {
return perms[i]; };
857 get_perm = [](std::size_t) {
return 0; };
862 const auto& [coeffs, cstride]
864 const std::vector<std::tuple<std::int32_t, int, std::int32_t, int>>&
866 = a.interior_facet_domains(i);
867 _lift_bc_interior_facets(b, *mesh, kernel, facets, dof_transform, dofmap0,
868 bs0, dof_transform_to_transpose, dofmap1, bs1,
869 constants, coeffs, cstride, cell_info, get_perm,
870 bc_values1, bc_markers1, x0, scale);
894 template <
typename T>
896 xtl::span<T> b,
const std::vector<std::shared_ptr<
const Form<T>>> a,
897 const std::vector<xtl::span<const T>>& constants,
898 const std::vector<std::map<std::pair<IntegralType, int>,
899 std::pair<xtl::span<const T>,
int>>>& coeffs,
900 const std::vector<std::vector<std::shared_ptr<
const DirichletBC<T>>>>& bcs1,
901 const std::vector<xtl::span<const T>>& x0,
double scale)
904 if (!x0.empty() and x0.size() != a.size())
906 throw std::runtime_error(
907 "Mismatch in size between x0 and bilinear form in assembler.");
910 if (a.size() != bcs1.size())
912 throw std::runtime_error(
913 "Mismatch in size between a and bcs in assembler.");
916 for (std::size_t j = 0; j < a.size(); ++j)
918 std::vector<std::int8_t> bc_markers1;
919 std::vector<T> bc_values1;
920 if (a[j] and !bcs1[j].empty())
922 assert(a[j]->function_spaces().at(0));
924 auto V1 = a[j]->function_spaces()[1];
926 auto map1 = V1->dofmap()->index_map;
927 const int bs1 = V1->dofmap()->index_map_bs();
929 const int crange = bs1 * (map1->size_local() + map1->num_ghosts());
930 bc_markers1.assign(crange,
false);
931 bc_values1.assign(crange, 0.0);
932 for (
const std::shared_ptr<
const DirichletBC<T>>& bc : bcs1[j])
934 bc->mark_dofs(bc_markers1);
935 bc->dof_values(bc_values1);
940 lift_bc<T>(b, *a[j], constants[j], coeffs[j], bc_values1, bc_markers1,
945 lift_bc<T>(b, *a[j], constants[j], coeffs[j], bc_values1, bc_markers1,
946 xtl::span<const T>(), scale);
958 template <
typename T>
960 xtl::span<T> b,
const Form<T>& L,
const xtl::span<const T>& constants,
961 const std::map<std::pair<IntegralType, int>,
962 std::pair<xtl::span<const T>,
int>>& coefficients)
964 std::shared_ptr<const mesh::Mesh> mesh = L.mesh();
968 assert(L.function_spaces().at(0));
969 std::shared_ptr<const fem::FiniteElement> element
970 = L.function_spaces().at(0)->element();
971 std::shared_ptr<const fem::DofMap> dofmap
972 = L.function_spaces().at(0)->dofmap();
974 const graph::AdjacencyList<std::int32_t>& dofs = dofmap->list();
975 const int bs = dofmap->bs();
977 const std::function<void(
const xtl::span<T>&,
978 const xtl::span<const std::uint32_t>&, std::int32_t,
980 dof_transform = element->get_dof_transformation_function<T>();
982 const bool needs_transformation_data
983 = element->needs_dof_transformations() or L.needs_facet_permutations();
984 xtl::span<const std::uint32_t> cell_info;
985 if (needs_transformation_data)
987 mesh->topology_mutable().create_entity_permutations();
988 cell_info = xtl::span(mesh->topology().get_cell_permutation_info());
995 const std::vector<std::int32_t>&
cells = L.cell_domains(i);
998 impl::assemble_cells<T, 1>(dof_transform, b, mesh->geometry(), cells,
999 dofs, bs, fn, constants, coeffs, cstride,
1004 impl::assemble_cells<T, 3>(dof_transform, b, mesh->geometry(), cells,
1005 dofs, bs, fn, constants, coeffs, cstride,
1010 impl::assemble_cells(dof_transform, b, mesh->geometry(), cells, dofs, bs,
1011 fn, constants, coeffs, cstride, cell_info);
1018 const auto& [coeffs, cstride]
1020 const std::vector<std::pair<std::int32_t, int>>& facets
1021 = L.exterior_facet_domains(i);
1024 impl::assemble_exterior_facets<T, 1>(dof_transform, b, *mesh, facets,
1025 dofs, bs, fn, constants, coeffs,
1026 cstride, cell_info);
1030 impl::assemble_exterior_facets<T, 3>(dof_transform, b, *mesh, facets,
1031 dofs, bs, fn, constants, coeffs,
1032 cstride, cell_info);
1036 impl::assemble_exterior_facets(dof_transform, b, *mesh, facets, dofs, bs,
1037 fn, constants, coeffs, cstride, cell_info);
1043 std::function<std::uint8_t(std::size_t)> get_perm;
1044 if (L.needs_facet_permutations())
1046 mesh->topology_mutable().create_entity_permutations();
1047 const std::vector<std::uint8_t>& perms
1048 = mesh->topology().get_facet_permutations();
1049 get_perm = [&perms](std::size_t i) {
return perms[i]; };
1052 get_perm = [](std::size_t) {
return 0; };
1057 const auto& [coeffs, cstride]
1059 const std::vector<std::tuple<std::int32_t, int, std::int32_t, int>>&
1061 = L.interior_facet_domains(i);
1064 impl::assemble_interior_facets<T, 1>(dof_transform, b, *mesh, facets,
1065 *dofmap, fn, constants, coeffs,
1066 cstride, cell_info, get_perm);
1070 impl::assemble_interior_facets<T, 3>(dof_transform, b, *mesh, facets,
1071 *dofmap, fn, constants, coeffs,
1072 cstride, cell_info, get_perm);
1076 impl::assemble_interior_facets(dof_transform, b, *mesh, facets, *dofmap,
1077 fn, constants, coeffs, cstride,
1078 cell_info, get_perm);
Degree-of-freedeom map representations ans tools.
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
void apply_lifting(xtl::span< T > b, const std::vector< std::shared_ptr< const Form< T >>> &a, const std::vector< xtl::span< const T >> &constants, const std::vector< std::map< std::pair< IntegralType, int >, std::pair< xtl::span< const T >, int >>> &coeffs, const std::vector< std::vector< std::shared_ptr< const DirichletBC< T >>>> &bcs1, const std::vector< xtl::span< const T >> &x0, double scale)
Modify b such that:
Definition: assembler.h:130
void assemble_vector(xtl::span< T > b, const Form< T > &L, const xtl::span< const T > &constants, const std::map< std::pair< IntegralType, int >, std::pair< xtl::span< const T >, int >> &coefficients)
Assemble linear form into a vector, The caller supplies the form constants and coefficients for this ...
Definition: assembler.h:89
@ 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:185