DOLFINx 0.10.0.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 <algorithm>
15#include <cstdint>
16#include <dolfinx/common/types.h>
17#include <memory>
18#include <optional>
19#include <span>
20#include <vector>
21
22namespace dolfinx::fem
23{
24template <dolfinx::scalar T, std::floating_point U>
25class DirichletBC;
26template <dolfinx::scalar T, std::floating_point U>
27class Form;
28template <std::floating_point T>
29class FunctionSpace;
30
31// -- Helper functions -----------------------------------------------------
32
34template <dolfinx::scalar T>
35std::map<std::pair<IntegralType, int>, std::pair<std::span<const T>, int>>
36make_coefficients_span(const std::map<std::pair<IntegralType, int>,
37 std::pair<std::vector<T>, int>>& coeffs)
38{
39 using Key = typename std::remove_reference_t<decltype(coeffs)>::key_type;
40 std::map<Key, std::pair<std::span<const T>, int>> c;
41 std::ranges::transform(
42 coeffs, std::inserter(c, c.end()),
43 [](auto& e) -> typename decltype(c)::value_type
44 { return {e.first, {e.second.first, e.second.second}}; });
45 return c;
46}
47
48// -- Scalar ----------------------------------------------------------------
49
61template <dolfinx::scalar T, std::floating_point U>
62T assemble_scalar(
63 const Form<T, U>& M, std::span<const T> constants,
64 const std::map<std::pair<IntegralType, int>,
65 std::pair<std::span<const T>, int>>& coefficients)
66{
67 std::shared_ptr<const mesh::Mesh<U>> mesh = M.mesh();
68 assert(mesh);
69 if constexpr (std::is_same_v<U, scalar_value_type_t<T>>)
70 {
71 return impl::assemble_scalar(M, mesh->geometry().dofmap(),
72 mesh->geometry().x(), constants, coefficients);
73 }
74 else
75 {
76 auto x = mesh->geometry().x();
77 std::vector<scalar_value_type_t<T>> _x(x.begin(), x.end());
78 return impl::assemble_scalar(M, mesh->geometry().dofmap(), _x, constants,
79 coefficients);
80 }
81}
82
88template <dolfinx::scalar T, std::floating_point U>
89T assemble_scalar(const Form<T, U>& M)
90{
91 const std::vector<T> constants = pack_constants(M);
92 auto coefficients = allocate_coefficient_storage(M);
93 pack_coefficients(M, coefficients);
94 return assemble_scalar(M, std::span(constants),
95 make_coefficients_span(coefficients));
96}
97
98// -- Vectors ----------------------------------------------------------------
99
110template <dolfinx::scalar T, std::floating_point U>
111void assemble_vector(
112 std::span<T> b, const Form<T, U>& L, std::span<const T> constants,
113 const std::map<std::pair<IntegralType, int>,
114 std::pair<std::span<const T>, int>>& coefficients)
115{
116 impl::assemble_vector(b, L, constants, coefficients);
117}
118
123template <dolfinx::scalar T, std::floating_point U>
124void assemble_vector(std::span<T> b, const Form<T, U>& L)
125{
126 auto coefficients = allocate_coefficient_storage(L);
127 pack_coefficients(L, coefficients);
128 const std::vector<T> constants = pack_constants(L);
129 assemble_vector(b, L, std::span(constants),
130 make_coefficients_span(coefficients));
131}
132
133// FIXME: clarify how x0 is used
134// FIXME: if bcs entries are set
135
136// FIXME: need to pass an array of Vec for x0?
137// FIXME: clarify zeroing of vector
138
151template <dolfinx::scalar T, std::floating_point U>
152void apply_lifting(
153 std::span<T> b,
154 std::vector<std::optional<std::reference_wrapper<const Form<T, U>>>> a,
155 const std::vector<std::span<const T>>& constants,
156 const std::vector<std::map<std::pair<IntegralType, int>,
157 std::pair<std::span<const T>, int>>>& coeffs,
158 const std::vector<
159 std::vector<std::reference_wrapper<const DirichletBC<T, U>>>>& bcs1,
160 const std::vector<std::span<const T>>& x0, T alpha)
161{
162 // If all forms are null, there is nothing to do
163 if (std::ranges::all_of(a, [](auto ai) { return !ai; }))
164 return;
165
166 impl::apply_lifting<T>(b, a, constants, coeffs, bcs1, x0, alpha);
167}
168
181template <dolfinx::scalar T, std::floating_point U>
182void apply_lifting(
183 std::span<T> b,
184 std::vector<std::optional<std::reference_wrapper<const Form<T, U>>>> a,
185 const std::vector<
186 std::vector<std::reference_wrapper<const DirichletBC<T, U>>>>& bcs1,
187 const std::vector<std::span<const T>>& x0, T alpha)
188{
189 std::vector<
190 std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, int>>>
191 coeffs;
192 std::vector<std::vector<T>> constants;
193 for (auto _a : a)
194 {
195 if (_a)
196 {
197 auto coefficients = allocate_coefficient_storage(_a->get());
198 pack_coefficients(_a->get(), coefficients);
199 coeffs.push_back(coefficients);
200 constants.push_back(pack_constants(_a->get()));
201 }
202 else
203 {
204 coeffs.emplace_back();
205 constants.emplace_back();
206 }
207 }
208
209 std::vector<std::span<const T>> _constants(constants.begin(),
210 constants.end());
211 std::vector<std::map<std::pair<IntegralType, int>,
212 std::pair<std::span<const T>, int>>>
213 _coeffs;
214 std::ranges::transform(coeffs, std::back_inserter(_coeffs),
215 [](auto& c) { return make_coefficients_span(c); });
216
217 apply_lifting(b, a, _constants, _coeffs, bcs1, x0, alpha);
218}
219
220// -- Matrices ---------------------------------------------------------------
221
234template <dolfinx::scalar T, std::floating_point U>
235void assemble_matrix(
236 la::MatSet<T> auto mat_add, const Form<T, U>& a,
237 std::span<const T> constants,
238 const std::map<std::pair<IntegralType, int>,
239 std::pair<std::span<const T>, int>>& coefficients,
240 std::span<const std::int8_t> dof_marker0,
241 std::span<const std::int8_t> dof_marker1)
242
243{
244 std::shared_ptr<const mesh::Mesh<U>> mesh = a.mesh();
245 assert(mesh);
246 if constexpr (std::is_same_v<U, scalar_value_type_t<T>>)
247 {
248 impl::assemble_matrix(mat_add, a, mesh->geometry().dofmap(),
249 mesh->geometry().x(), constants, coefficients,
250 dof_marker0, dof_marker1);
251 }
252 else
253 {
254 auto x = mesh->geometry().x();
255 std::vector<scalar_value_type_t<T>> _x(x.begin(), x.end());
256 impl::assemble_matrix(mat_add, a, mesh->geometry().dofmap(), _x, constants,
257 coefficients, dof_marker0, dof_marker1);
258 }
259}
260
268template <dolfinx::scalar T, std::floating_point U>
269void assemble_matrix(
270 auto mat_add, const Form<T, U>& a, std::span<const T> constants,
271 const std::map<std::pair<IntegralType, int>,
272 std::pair<std::span<const T>, int>>& coefficients,
273 const std::vector<std::reference_wrapper<const DirichletBC<T, U>>>& bcs)
274{
275 // Index maps for dof ranges
276 auto map0 = a.function_spaces().at(0)->dofmap()->index_map;
277 auto map1 = a.function_spaces().at(1)->dofmap()->index_map;
278 auto bs0 = a.function_spaces().at(0)->dofmap()->index_map_bs();
279 auto bs1 = a.function_spaces().at(1)->dofmap()->index_map_bs();
280
281 // Build dof markers
282 std::vector<std::int8_t> dof_marker0, dof_marker1;
283 assert(map0);
284 std::int32_t dim0 = bs0 * (map0->size_local() + map0->num_ghosts());
285 assert(map1);
286 std::int32_t dim1 = bs1 * (map1->size_local() + map1->num_ghosts());
287 for (std::size_t k = 0; k < bcs.size(); ++k)
288 {
289 assert(bcs[k].get().function_space());
290 if (a.function_spaces().at(0)->contains(*bcs[k].get().function_space()))
291 {
292 dof_marker0.resize(dim0, false);
293 bcs[k].get().mark_dofs(dof_marker0);
294 }
295
296 if (a.function_spaces().at(1)->contains(*bcs[k].get().function_space()))
297 {
298 dof_marker1.resize(dim1, false);
299 bcs[k].get().mark_dofs(dof_marker1);
300 }
301 }
302
303 // Assemble
304 assemble_matrix(mat_add, a, constants, coefficients, dof_marker0,
305 dof_marker1);
306}
307
313template <dolfinx::scalar T, std::floating_point U>
314void assemble_matrix(
315 auto mat_add, const Form<T, U>& a,
316 const std::vector<std::reference_wrapper<const DirichletBC<T, U>>>& bcs)
317{
318 // Prepare constants and coefficients
319 const std::vector<T> constants = pack_constants(a);
320 auto coefficients = allocate_coefficient_storage(a);
321 pack_coefficients(a, coefficients);
322
323 // Assemble
324 assemble_matrix(mat_add, a, std::span(constants),
325 make_coefficients_span(coefficients), bcs);
326}
327
338template <dolfinx::scalar T, std::floating_point U>
339void assemble_matrix(auto mat_add, const Form<T, U>& a,
340 std::span<const std::int8_t> dof_marker0,
341 std::span<const std::int8_t> dof_marker1)
342
343{
344 // Prepare constants and coefficients
345 const std::vector<T> constants = pack_constants(a);
346 auto coefficients = allocate_coefficient_storage(a);
347 pack_coefficients(a, coefficients);
348
349 // Assemble
350 assemble_matrix(mat_add, a, std::span(constants),
351 make_coefficients_span(coefficients), dof_marker0,
352 dof_marker1);
353}
354
366template <dolfinx::scalar T>
367void set_diagonal(auto set_fn, std::span<const std::int32_t> rows,
368 T diagonal = 1.0)
369{
370 for (std::size_t i = 0; i < rows.size(); ++i)
371 {
372 std::span diag_span(&diagonal, 1);
373 set_fn(rows.subspan(i, 1), rows.subspan(i, 1), diag_span);
374 }
375}
376
393template <dolfinx::scalar T, std::floating_point U>
395 auto set_fn, const FunctionSpace<U>& V,
396 const std::vector<std::reference_wrapper<const DirichletBC<T, U>>>& bcs,
397 T diagonal = 1.0)
398{
399 for (auto& bc : bcs)
400 {
401 if (V.contains(*bc.get().function_space()))
402 {
403 const auto [dofs, range] = bc.get().dof_indices();
404 set_diagonal(set_fn, dofs.first(range), diagonal);
405 }
406 }
407}
408} // namespace dolfinx::fem
Definition DirichletBC.h:260
A representation of finite element variational forms.
Definition Form.h:139
This class represents a finite element function space defined by a mesh, a finite element,...
Definition vtk_utils.h:32
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:973
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:367
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:36
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:1314
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:1027