DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
assembler.h
1// Copyright (C) 2018-2022 Garth N. Wells
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include "assemble_matrix_impl.h"
10#include "assemble_scalar_impl.h"
11#include "assemble_vector_impl.h"
12#include "traits.h"
13#include "utils.h"
14#include <cstdint>
15#include <dolfinx/common/types.h>
16#include <memory>
17#include <span>
18#include <vector>
19
20namespace dolfinx::fem
21{
22template <dolfinx::scalar T, std::floating_point U>
23class DirichletBC;
24template <dolfinx::scalar T, std::floating_point U>
25class Form;
26template <std::floating_point T>
27class FunctionSpace;
28
29// -- Helper functions -----------------------------------------------------
30
32template <dolfinx::scalar T>
33std::map<std::pair<IntegralType, int>, std::pair<std::span<const T>, int>>
34make_coefficients_span(const std::map<std::pair<IntegralType, int>,
35 std::pair<std::vector<T>, int>>& coeffs)
36{
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}}; });
42 return c;
43}
44
45// -- Scalar ----------------------------------------------------------------
46
58template <dolfinx::scalar T, std::floating_point U>
59T assemble_scalar(
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)
63{
64 std::shared_ptr<const mesh::Mesh<U>> mesh = M.mesh();
65 assert(mesh);
66 if constexpr (std::is_same_v<U, scalar_value_type_t<T>>)
67 {
68 return impl::assemble_scalar(M, mesh->geometry().dofmap(),
69 mesh->geometry().x(), constants, coefficients);
70 }
71 else
72 {
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,
76 coefficients);
77 }
78}
79
85template <dolfinx::scalar T, std::floating_point U>
86T assemble_scalar(const Form<T, U>& M)
87{
88 const std::vector<T> constants = pack_constants(M);
89 auto coefficients = allocate_coefficient_storage(M);
90 pack_coefficients(M, coefficients);
91 return assemble_scalar(M, std::span(constants),
92 make_coefficients_span(coefficients));
93}
94
95// -- Vectors ----------------------------------------------------------------
96
107template <dolfinx::scalar T, std::floating_point U>
108void assemble_vector(
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)
112{
113 impl::assemble_vector(b, L, constants, coefficients);
114}
115
120template <dolfinx::scalar T, std::floating_point U>
121void assemble_vector(std::span<T> b, const Form<T, U>& L)
122{
123 auto coefficients = allocate_coefficient_storage(L);
124 pack_coefficients(L, coefficients);
125 const std::vector<T> constants = pack_constants(L);
126 assemble_vector(b, L, std::span(constants),
127 make_coefficients_span(coefficients));
128}
129
130// FIXME: clarify how x0 is used
131// FIXME: if bcs entries are set
132
133// FIXME: need to pass an array of Vec for x0?
134// FIXME: clarify zeroing of vector
135
148template <dolfinx::scalar T, std::floating_point U>
149void apply_lifting(
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,
154 const std::vector<std::vector<std::shared_ptr<const DirichletBC<T, U>>>>&
155 bcs1,
156 const std::vector<std::span<const T>>& x0, T scale)
157{
158 std::shared_ptr<const mesh::Mesh<U>> mesh;
159 for (auto& a_i : a)
160 {
161 if (a_i and !mesh)
162 mesh = a_i->mesh();
163 if (a_i and mesh and a_i->mesh() != mesh)
164 throw std::runtime_error("Mismatch between meshes.");
165 }
166
167 if (!mesh)
168 throw std::runtime_error("Unable to extract a mesh.");
169
170 if constexpr (std::is_same_v<U, scalar_value_type_t<T>>)
171 {
172 impl::apply_lifting<T>(b, a, mesh->geometry().dofmap(),
173 mesh->geometry().x(), constants, coeffs, bcs1, x0,
174 scale);
175 }
176 else
177 {
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);
182 }
183}
184
197template <dolfinx::scalar T, std::floating_point U>
198void apply_lifting(
199 std::span<T> b, const std::vector<std::shared_ptr<const Form<T, U>>>& a,
200 const std::vector<std::vector<std::shared_ptr<const DirichletBC<T, U>>>>&
201 bcs1,
202 const std::vector<std::span<const T>>& x0, T scale)
203{
204 std::vector<
205 std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, int>>>
206 coeffs;
207 std::vector<std::vector<T>> constants;
208 for (auto _a : a)
209 {
210 if (_a)
211 {
212 auto coefficients = allocate_coefficient_storage(*_a);
213 pack_coefficients(*_a, coefficients);
214 coeffs.push_back(coefficients);
215 constants.push_back(pack_constants(*_a));
216 }
217 else
218 {
219 coeffs.emplace_back();
220 constants.emplace_back();
221 }
222 }
223
224 std::vector<std::span<const T>> _constants(constants.begin(),
225 constants.end());
226 std::vector<std::map<std::pair<IntegralType, int>,
227 std::pair<std::span<const T>, int>>>
228 _coeffs;
229 std::transform(coeffs.cbegin(), coeffs.cend(), std::back_inserter(_coeffs),
230 [](auto& c) { return make_coefficients_span(c); });
231
232 apply_lifting(b, a, _constants, _coeffs, bcs1, x0, scale);
233}
234
235// -- Matrices ---------------------------------------------------------------
236
249template <dolfinx::scalar T, std::floating_point U>
250void assemble_matrix(
251 la::MatSet<T> auto mat_add, const Form<T, U>& a,
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)
257
258{
259 std::shared_ptr<const mesh::Mesh<U>> mesh = a.mesh();
260 assert(mesh);
261 if constexpr (std::is_same_v<U, scalar_value_type_t<T>>)
262 {
263 impl::assemble_matrix(mat_add, a, mesh->geometry().dofmap(),
264 mesh->geometry().x(), constants, coefficients,
265 dof_marker0, dof_marker1);
266 }
267 else
268 {
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);
273 }
274}
275
283template <dolfinx::scalar T, std::floating_point U>
284void assemble_matrix(
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,
288 const std::vector<std::shared_ptr<const DirichletBC<T, U>>>& bcs)
289{
290 // Index maps for dof ranges
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();
295
296 // Build dof markers
297 std::vector<std::int8_t> dof_marker0, dof_marker1;
298 assert(map0);
299 std::int32_t dim0 = bs0 * (map0->size_local() + map0->num_ghosts());
300 assert(map1);
301 std::int32_t dim1 = bs1 * (map1->size_local() + map1->num_ghosts());
302 for (std::size_t k = 0; k < bcs.size(); ++k)
303 {
304 assert(bcs[k]);
305 assert(bcs[k]->function_space());
306 if (a.function_spaces().at(0)->contains(*bcs[k]->function_space()))
307 {
308 dof_marker0.resize(dim0, false);
309 bcs[k]->mark_dofs(dof_marker0);
310 }
311
312 if (a.function_spaces().at(1)->contains(*bcs[k]->function_space()))
313 {
314 dof_marker1.resize(dim1, false);
315 bcs[k]->mark_dofs(dof_marker1);
316 }
317 }
318
319 // Assemble
320 assemble_matrix(mat_add, a, constants, coefficients, dof_marker0,
321 dof_marker1);
322}
323
329template <dolfinx::scalar T, std::floating_point U>
330void assemble_matrix(
331 auto mat_add, const Form<T, U>& a,
332 const std::vector<std::shared_ptr<const DirichletBC<T, U>>>& bcs)
333{
334 // Prepare constants and coefficients
335 const std::vector<T> constants = pack_constants(a);
336 auto coefficients = allocate_coefficient_storage(a);
337 pack_coefficients(a, coefficients);
338
339 // Assemble
340 assemble_matrix(mat_add, a, std::span(constants),
341 make_coefficients_span(coefficients), bcs);
342}
343
354template <dolfinx::scalar T, std::floating_point U>
355void assemble_matrix(auto mat_add, const Form<T, U>& a,
356 std::span<const std::int8_t> dof_marker0,
357 std::span<const std::int8_t> dof_marker1)
358
359{
360 // Prepare constants and coefficients
361 const std::vector<T> constants = pack_constants(a);
362 auto coefficients = allocate_coefficient_storage(a);
363 pack_coefficients(a, coefficients);
364
365 // Assemble
366 assemble_matrix(mat_add, a, std::span(constants),
367 make_coefficients_span(coefficients), dof_marker0,
368 dof_marker1);
369}
370
382template <dolfinx::scalar T>
383void set_diagonal(auto set_fn, std::span<const std::int32_t> rows,
384 T diagonal = 1.0)
385{
386 for (std::size_t i = 0; i < rows.size(); ++i)
387 {
388 std::span diag_span(&diagonal, 1);
389 set_fn(rows.subspan(i, 1), rows.subspan(i, 1), diag_span);
390 }
391}
392
409template <dolfinx::scalar T, std::floating_point U>
411 auto set_fn, const FunctionSpace<U>& V,
412 const std::vector<std::shared_ptr<const DirichletBC<T, U>>>& bcs,
413 T diagonal = 1.0)
414{
415 for (auto& bc : bcs)
416 {
417 assert(bc);
418 if (V.contains(*bc->function_space()))
419 {
420 const auto [dofs, range] = bc->dof_indices();
421 set_diagonal(set_fn, dofs.first(range), diagonal);
422 }
423 }
424}
425
426// -- Setting bcs ------------------------------------------------------------
427
428// FIXME: Move these function elsewhere?
429
430// FIXME: clarify x0
431// FIXME: clarify what happens with ghosts
432
436template <dolfinx::scalar T, std::floating_point U>
437void set_bc(std::span<T> b,
438 const std::vector<std::shared_ptr<const DirichletBC<T, U>>>& bcs,
439 std::span<const T> x0, T scale = 1)
440{
441 if (b.size() > x0.size())
442 throw std::runtime_error("Size mismatch between b and x0 vectors.");
443 for (auto& bc : bcs)
444 {
445 assert(bc);
446 bc->set(b, x0, scale);
447 }
448}
449
453template <dolfinx::scalar T, std::floating_point U>
454void set_bc(std::span<T> b,
455 const std::vector<std::shared_ptr<const DirichletBC<T, U>>>& bcs,
456 T scale = 1)
457{
458 for (auto& bc : bcs)
459 {
460 assert(bc);
461 bc->set(b, scale);
462 }
463}
464
465} // namespace dolfinx::fem
Definition DirichletBC.h:259
A representation of finite element variational forms.
Definition Form.h:120
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