9#include "assemble_matrix_impl.h"
10#include "assemble_scalar_impl.h"
11#include "assemble_vector_impl.h"
15#include <dolfinx/common/types.h>
22template <dolfinx::scalar T, std::
floating_po
int U>
24template <dolfinx::scalar T, std::
floating_po
int U>
26template <std::
floating_po
int T>
32template <dolfinx::scalar T>
33std::map<std::pair<IntegralType, int>, std::pair<std::span<const T>,
int>>
35 std::pair<std::vector<T>,
int>>& coeffs)
37 using Key =
typename std::remove_reference_t<
decltype(coeffs)>::key_type;
38 std::map<Key, std::pair<std::span<const T>,
int>> c;
39 std::transform(coeffs.cbegin(), coeffs.cend(), std::inserter(c, c.end()),
40 [](
auto& e) ->
typename decltype(c)::value_type
41 { return {e.first, {e.second.first, e.second.second}}; });
58template <dolfinx::scalar T, std::
floating_po
int U>
60 const Form<T, U>& M, std::span<const T> constants,
61 const std::map<std::pair<IntegralType, int>,
62 std::pair<std::span<const T>,
int>>& coefficients)
64 std::shared_ptr<const mesh::Mesh<U>> mesh = M.mesh();
66 if constexpr (std::is_same_v<U, scalar_value_type_t<T>>)
68 return impl::assemble_scalar(M, mesh->geometry().dofmap(),
69 mesh->geometry().x(), constants, coefficients);
73 auto x = mesh->geometry().x();
74 std::vector<scalar_value_type_t<T>> _x(x.begin(), x.end());
75 return impl::assemble_scalar(M, mesh->geometry().dofmap(), _x, constants,
85template <dolfinx::scalar T, std::
floating_po
int U>
91 return assemble_scalar(M, std::span(constants),
107template <dolfinx::scalar T, std::
floating_po
int U>
109 std::span<T> b,
const Form<T, U>& L, std::span<const T> constants,
110 const std::map<std::pair<IntegralType, int>,
111 std::pair<std::span<const T>,
int>>& coefficients)
113 impl::assemble_vector(b, L, constants, coefficients);
120template <dolfinx::scalar T, std::
floating_po
int U>
126 assemble_vector(b, L, std::span(constants),
148template <dolfinx::scalar T, std::
floating_po
int U>
150 std::span<T> b,
const std::vector<std::shared_ptr<
const Form<T, U>>>& a,
151 const std::vector<std::span<const T>>& constants,
152 const std::vector<std::map<std::pair<IntegralType, int>,
153 std::pair<std::span<const T>,
int>>>& coeffs,
156 const std::vector<std::span<const T>>& x0, T scale)
158 std::shared_ptr<const mesh::Mesh<U>> mesh;
163 if (a_i and mesh and a_i->mesh() != mesh)
164 throw std::runtime_error(
"Mismatch between meshes.");
168 throw std::runtime_error(
"Unable to extract a mesh.");
170 if constexpr (std::is_same_v<U, scalar_value_type_t<T>>)
172 impl::apply_lifting<T>(b, a, mesh->geometry().dofmap(),
173 mesh->geometry().x(), constants, coeffs, bcs1, x0,
178 auto x = mesh->geometry().x();
179 std::vector<scalar_value_type_t<T>> _x(x.begin(), x.end());
180 impl::apply_lifting<T>(b, a, mesh->geometry().dofmap(), _x, constants,
181 coeffs, bcs1, x0, scale);
197template <dolfinx::scalar T, std::
floating_po
int U>
199 std::span<T> b,
const std::vector<std::shared_ptr<
const Form<T, U>>>& a,
202 const std::vector<std::span<const T>>& x0, T scale)
205 std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>,
int>>>
207 std::vector<std::vector<T>> constants;
214 coeffs.push_back(coefficients);
219 coeffs.emplace_back();
220 constants.emplace_back();
224 std::vector<std::span<const T>> _constants(constants.begin(),
226 std::vector<std::map<std::pair<IntegralType, int>,
227 std::pair<std::span<const T>,
int>>>
229 std::transform(coeffs.cbegin(), coeffs.cend(), std::back_inserter(_coeffs),
230 [](
auto& c) { return make_coefficients_span(c); });
232 apply_lifting(b, a, _constants, _coeffs, bcs1, x0, scale);
249template <dolfinx::scalar T, std::
floating_po
int U>
252 std::span<const T> constants,
253 const std::map<std::pair<IntegralType, int>,
254 std::pair<std::span<const T>,
int>>& coefficients,
255 std::span<const std::int8_t> dof_marker0,
256 std::span<const std::int8_t> dof_marker1)
259 std::shared_ptr<const mesh::Mesh<U>> mesh = a.mesh();
261 if constexpr (std::is_same_v<U, scalar_value_type_t<T>>)
263 impl::assemble_matrix(mat_add, a, mesh->geometry().dofmap(),
264 mesh->geometry().x(), constants, coefficients,
265 dof_marker0, dof_marker1);
269 auto x = mesh->geometry().x();
270 std::vector<scalar_value_type_t<T>> _x(x.begin(), x.end());
271 impl::assemble_matrix(mat_add, a, mesh->geometry().dofmap(), _x, constants,
272 coefficients, dof_marker0, dof_marker1);
283template <dolfinx::scalar T, std::
floating_po
int U>
285 auto mat_add,
const Form<T, U>& a, std::span<const T> constants,
286 const std::map<std::pair<IntegralType, int>,
287 std::pair<std::span<const T>,
int>>& coefficients,
291 auto map0 = a.function_spaces().at(0)->dofmap()->index_map;
292 auto map1 = a.function_spaces().at(1)->dofmap()->index_map;
293 auto bs0 = a.function_spaces().at(0)->dofmap()->index_map_bs();
294 auto bs1 = a.function_spaces().at(1)->dofmap()->index_map_bs();
297 std::vector<std::int8_t> dof_marker0, dof_marker1;
299 std::int32_t dim0 = bs0 * (map0->size_local() + map0->num_ghosts());
301 std::int32_t dim1 = bs1 * (map1->size_local() + map1->num_ghosts());
302 for (std::size_t k = 0; k < bcs.size(); ++k)
305 assert(bcs[k]->function_space());
306 if (a.function_spaces().at(0)->contains(*bcs[k]->function_space()))
308 dof_marker0.resize(dim0,
false);
309 bcs[k]->mark_dofs(dof_marker0);
312 if (a.function_spaces().at(1)->contains(*bcs[k]->function_space()))
314 dof_marker1.resize(dim1,
false);
315 bcs[k]->mark_dofs(dof_marker1);
320 assemble_matrix(mat_add, a, constants, coefficients, dof_marker0,
329template <dolfinx::scalar T, std::
floating_po
int U>
340 assemble_matrix(mat_add, a, std::span(constants),
354template <dolfinx::scalar T, std::
floating_po
int U>
356 std::span<const std::int8_t> dof_marker0,
357 std::span<const std::int8_t> dof_marker1)
366 assemble_matrix(mat_add, a, std::span(constants),
382template <dolfinx::scalar T>
386 for (std::size_t i = 0; i < rows.size(); ++i)
388 std::span diag_span(&diagonal, 1);
389 set_fn(rows.subspan(i, 1), rows.subspan(i, 1), diag_span);
409template <dolfinx::scalar T, std::
floating_po
int U>
418 if (V.contains(*bc->function_space()))
420 const auto [dofs, range] = bc->dof_indices();
436template <dolfinx::scalar T, std::
floating_po
int U>
439 std::span<const T> x0, T scale = 1)
441 if (b.size() > x0.size())
442 throw std::runtime_error(
"Size mismatch between b and x0 vectors.");
446 bc->set(b, x0, scale);
453template <dolfinx::scalar T, std::
floating_po
int U>
Definition DirichletBC.h:259
This class represents a finite element function space defined by a mesh, a finite element,...
Definition FunctionSpace.h:34
Matrix accumulate/set concept for functions that can be used in assemblers to accumulate or set value...
Definition utils.h:28
Functions supporting finite element method operations.
Finite element method functionality.
Definition assemble_matrix_impl.h:26
std::pair< std::vector< T >, int > allocate_coefficient_storage(const Form< T, U > &form, IntegralType integral_type, int id)
Allocate storage for coefficients of a pair (integral_type, id) from a Form.
Definition utils.h:965
void set_diagonal(auto set_fn, std::span< const std::int32_t > rows, T diagonal=1.0)
Sets a value to the diagonal of a matrix for specified rows.
Definition assembler.h:383
std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > make_coefficients_span(const std::map< std::pair< IntegralType, int >, std::pair< std::vector< T >, int > > &coeffs)
Create a map of std::spans from a map of std::vectors.
Definition assembler.h:34
std::vector< typename U::scalar_type > pack_constants(const U &u)
Pack constants of u into a single array ready for assembly.
Definition utils.h:1249
void set_bc(std::span< T > b, const std::vector< std::shared_ptr< const DirichletBC< T, U > > > &bcs, std::span< const T > x0, T scale=1)
Definition assembler.h:437
void pack_coefficients(const Form< T, U > &form, IntegralType integral_type, int id, std::span< T > c, int cstride)
Pack coefficients of a Form for a given integral type and domain id.
Definition utils.h:1019