Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d4/dec/assembler_8h_source.html
DOLFINx 0.6.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 <cstdint>
13#include <memory>
14#include <span>
15#include <vector>
16
17namespace dolfinx::fem
18{
19template <typename T>
20class DirichletBC;
21template <typename T>
22class Form;
23class FunctionSpace;
24
25// -- Helper functions -----------------------------------------------------
26
28template <typename T>
29std::map<std::pair<fem::IntegralType, int>, std::pair<std::span<const T>, int>>
30make_coefficients_span(const std::map<std::pair<IntegralType, int>,
31 std::pair<std::vector<T>, int>>& coeffs)
32{
33 using Key = typename std::remove_reference_t<decltype(coeffs)>::key_type;
34 std::map<Key, std::pair<std::span<const T>, int>> c;
35 std::transform(coeffs.cbegin(), coeffs.cend(), std::inserter(c, c.end()),
36 [](auto& e) -> typename decltype(c)::value_type {
37 return {e.first, {e.second.first, e.second.second}};
38 });
39 return c;
40}
41
42// -- Scalar ----------------------------------------------------------------
43
55template <typename T>
56T assemble_scalar(
57 const Form<T>& M, std::span<const T> constants,
58 const std::map<std::pair<IntegralType, int>,
59 std::pair<std::span<const T>, int>>& coefficients)
60{
61 return impl::assemble_scalar(M, constants, coefficients);
62}
63
69template <typename T>
70T assemble_scalar(const Form<T>& M)
71{
72 const std::vector<T> constants = pack_constants(M);
73 auto coefficients = allocate_coefficient_storage(M);
74 pack_coefficients(M, coefficients);
75 return assemble_scalar(M, std::span(constants),
76 make_coefficients_span(coefficients));
77}
78
79// -- Vectors ----------------------------------------------------------------
80
91template <typename T>
92void assemble_vector(
93 std::span<T> b, const Form<T>& L, std::span<const T> constants,
94 const std::map<std::pair<IntegralType, int>,
95 std::pair<std::span<const T>, int>>& coefficients)
96{
97 impl::assemble_vector(b, L, constants, coefficients);
98}
99
104template <typename T>
105void assemble_vector(std::span<T> b, const Form<T>& L)
106{
107 auto coefficients = allocate_coefficient_storage(L);
108 pack_coefficients(L, coefficients);
109 const std::vector<T> constants = pack_constants(L);
110 assemble_vector(b, L, std::span(constants),
111 make_coefficients_span(coefficients));
112}
113
114// FIXME: clarify how x0 is used
115// FIXME: if bcs entries are set
116
117// FIXME: need to pass an array of Vec for x0?
118// FIXME: clarify zeroing of vector
119
132template <typename T>
133void apply_lifting(
134 std::span<T> b, const std::vector<std::shared_ptr<const Form<T>>>& a,
135 const std::vector<std::span<const T>>& constants,
136 const std::vector<std::map<std::pair<IntegralType, int>,
137 std::pair<std::span<const T>, int>>>& coeffs,
138 const std::vector<std::vector<std::shared_ptr<const DirichletBC<T>>>>& bcs1,
139 const std::vector<std::span<const T>>& x0, double scale)
140{
141 impl::apply_lifting(b, a, constants, coeffs, bcs1, x0, scale);
142}
143
156template <typename T>
157void apply_lifting(
158 std::span<T> b, const std::vector<std::shared_ptr<const Form<T>>>& a,
159 const std::vector<std::vector<std::shared_ptr<const DirichletBC<T>>>>& bcs1,
160 const std::vector<std::span<const T>>& x0, double scale)
161{
162 std::vector<
163 std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, int>>>
164 coeffs;
165 std::vector<std::vector<T>> constants;
166 for (auto _a : a)
167 {
168 if (_a)
169 {
170 auto coefficients = allocate_coefficient_storage(*_a);
171 pack_coefficients(*_a, coefficients);
172 coeffs.push_back(coefficients);
173 constants.push_back(pack_constants(*_a));
174 }
175 else
176 {
177 coeffs.emplace_back();
178 constants.emplace_back();
179 }
180 }
181
182 std::vector<std::span<const T>> _constants(constants.begin(),
183 constants.end());
184 std::vector<std::map<std::pair<IntegralType, int>,
185 std::pair<std::span<const T>, int>>>
186 _coeffs;
187 std::transform(coeffs.cbegin(), coeffs.cend(), std::back_inserter(_coeffs),
188 [](auto& c) { return make_coefficients_span(c); });
189 apply_lifting(b, a, _constants, _coeffs, bcs1, x0, scale);
190}
191
192// -- Matrices ---------------------------------------------------------------
193
201template <typename T>
202void assemble_matrix(
203 auto mat_add, const Form<T>& a, std::span<const T> constants,
204 const std::map<std::pair<IntegralType, int>,
205 std::pair<std::span<const T>, int>>& coefficients,
206 const std::vector<std::shared_ptr<const DirichletBC<T>>>& bcs)
207{
208 // Index maps for dof ranges
209 auto map0 = a.function_spaces().at(0)->dofmap()->index_map;
210 auto map1 = a.function_spaces().at(1)->dofmap()->index_map;
211 auto bs0 = a.function_spaces().at(0)->dofmap()->index_map_bs();
212 auto bs1 = a.function_spaces().at(1)->dofmap()->index_map_bs();
213
214 // Build dof markers
215 std::vector<std::int8_t> dof_marker0, dof_marker1;
216 assert(map0);
217 std::int32_t dim0 = bs0 * (map0->size_local() + map0->num_ghosts());
218 assert(map1);
219 std::int32_t dim1 = bs1 * (map1->size_local() + map1->num_ghosts());
220 for (std::size_t k = 0; k < bcs.size(); ++k)
221 {
222 assert(bcs[k]);
223 assert(bcs[k]->function_space());
224 if (a.function_spaces().at(0)->contains(*bcs[k]->function_space()))
225 {
226 dof_marker0.resize(dim0, false);
227 bcs[k]->mark_dofs(dof_marker0);
228 }
229
230 if (a.function_spaces().at(1)->contains(*bcs[k]->function_space()))
231 {
232 dof_marker1.resize(dim1, false);
233 bcs[k]->mark_dofs(dof_marker1);
234 }
235 }
236
237 // Assemble
238 impl::assemble_matrix(mat_add, a, constants, coefficients, dof_marker0,
239 dof_marker1);
240}
241
247template <typename T>
248void assemble_matrix(
249 auto mat_add, const Form<T>& a,
250 const std::vector<std::shared_ptr<const DirichletBC<T>>>& bcs)
251{
252 // Prepare constants and coefficients
253 const std::vector<T> constants = pack_constants(a);
254 auto coefficients = allocate_coefficient_storage(a);
255 pack_coefficients(a, coefficients);
256
257 // Assemble
258 assemble_matrix(mat_add, a, std::span(constants),
259 make_coefficients_span(coefficients), bcs);
260}
261
274template <typename T>
275void assemble_matrix(
276 auto mat_add, const Form<T>& a, std::span<const T> constants,
277 const std::map<std::pair<IntegralType, int>,
278 std::pair<std::span<const T>, int>>& coefficients,
279 std::span<const std::int8_t> dof_marker0,
280 std::span<const std::int8_t> dof_marker1)
281
282{
283 impl::assemble_matrix(mat_add, a, constants, coefficients, dof_marker0,
284 dof_marker1);
285}
286
297template <typename T>
298void assemble_matrix(auto mat_add, const Form<T>& a,
299 std::span<const std::int8_t> dof_marker0,
300 std::span<const std::int8_t> dof_marker1)
301
302{
303 // Prepare constants and coefficients
304 const std::vector<T> constants = pack_constants(a);
305 auto coefficients = allocate_coefficient_storage(a);
306 pack_coefficients(a, coefficients);
307
308 // Assemble
309 assemble_matrix(mat_add, a, std::span(constants),
310 make_coefficients_span(coefficients), dof_marker0,
311 dof_marker1);
312}
313
325template <typename T>
326void set_diagonal(auto set_fn, std::span<const std::int32_t> rows,
327 T diagonal = 1.0)
328{
329 for (std::size_t i = 0; i < rows.size(); ++i)
330 {
331 std::span diag_span(&diagonal, 1);
332 set_fn(rows.subspan(i, 1), rows.subspan(i, 1), diag_span);
333 }
334}
335
352template <typename T>
353void set_diagonal(auto set_fn, const fem::FunctionSpace& V,
354 const std::vector<std::shared_ptr<const DirichletBC<T>>>& bcs,
355 T diagonal = 1.0)
356{
357 for (const auto& bc : bcs)
358 {
359 assert(bc);
360 if (V.contains(*bc->function_space()))
361 {
362 const auto [dofs, range] = bc->dof_indices();
363 set_diagonal(set_fn, dofs.first(range), diagonal);
364 }
365 }
366}
367
368// -- Setting bcs ------------------------------------------------------------
369
370// FIXME: Move these function elsewhere?
371
372// FIXME: clarify x0
373// FIXME: clarify what happens with ghosts
374
378template <typename T>
379void set_bc(std::span<T> b,
380 const std::vector<std::shared_ptr<const DirichletBC<T>>>& bcs,
381 std::span<const T> x0, double scale = 1.0)
382{
383 if (b.size() > x0.size())
384 throw std::runtime_error("Size mismatch between b and x0 vectors.");
385 for (const auto& bc : bcs)
386 {
387 assert(bc);
388 bc->set(b, x0, scale);
389 }
390}
391
395template <typename T>
396void set_bc(std::span<T> b,
397 const std::vector<std::shared_ptr<const DirichletBC<T>>>& bcs,
398 double scale = 1.0)
399{
400 for (const auto& bc : bcs)
401 {
402 assert(bc);
403 bc->set(b, scale);
404 }
405}
406
407} // namespace dolfinx::fem
Object for setting (strong) Dirichlet boundary conditions.
Definition: DirichletBC.h:133
A representation of finite element variational forms.
Definition: Form.h:64
const std::vector< std::shared_ptr< const FunctionSpace > > & function_spaces() const
Return function spaces for all arguments.
Definition: Form.h:182
This class represents a finite element function space defined by a mesh, a finite element,...
Definition: FunctionSpace.h:31
bool contains(const FunctionSpace &V) const
Check whether V is subspace of this, or this itself.
Definition: FunctionSpace.cpp:69
Finite element method functionality.
Definition: assemble_matrix_impl.h:25
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:326
std::map< std::pair< fem::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:30
void set_bc(std::span< T > b, const std::vector< std::shared_ptr< const DirichletBC< T > > > &bcs, std::span< const T > x0, double scale=1.0)
Set bc values in owned (local) part of the vector, multiplied by 'scale'. The vectors b and x0 must h...
Definition: assembler.h:379
void pack_coefficients(const Form< T > &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:755
std::vector< typename U::scalar_type > pack_constants(const U &u)
Pack constants of u of generic type U ready for assembly.
Definition: utils.h:985
std::pair< std::vector< T >, int > allocate_coefficient_storage(const Form< T > &form, IntegralType integral_type, int id)
Allocate storage for coefficients of a pair (integral_type, id) from a fem::Form form.
Definition: utils.h:692