DOLFINx 0.11.0
DOLFINx C++
Loading...
Searching...
No Matches
assembler.h
Go to the documentation of this file.
1// Copyright (C) 2018-2026 Garth N. Wells and Jørgen S. Dokken
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 "Function.h"
10#include "FunctionSpace.h"
11#include "assemble_expression_impl.h"
12#include "assemble_matrix_impl.h"
13#include "assemble_scalar_impl.h"
14#include "assemble_vector_impl.h"
15#include "pack.h"
16#include "traits.h"
17#include "utils.h"
18#include <algorithm>
19#include <basix/mdspan.hpp>
20#include <cstdint>
21#include <dolfinx/common/types.h>
22#include <dolfinx/mesh/EntityMap.h>
23#include <memory>
24#include <optional>
25#include <span>
26#include <vector>
27
31
32namespace dolfinx::fem
33{
34template <dolfinx::scalar T, std::floating_point U>
35class DirichletBC;
36template <dolfinx::scalar T, std::floating_point U>
37class Expression;
38template <dolfinx::scalar T, std::floating_point U>
39class Form;
40template <std::floating_point T>
41class FunctionSpace;
42
64template <dolfinx::scalar T, std::floating_point U>
66 std::span<T> values, const fem::Expression<T, U>& e,
67 md::mdspan<const T, md::dextents<std::size_t, 2>> coeffs,
68 std::span<const T> constants, const mesh::Mesh<U>& mesh,
69 fem::MDSpan2 auto entities,
70 std::optional<
71 std::pair<std::reference_wrapper<const FiniteElement<U>>, std::size_t>>
72 element)
73{
74 // Check that domain is the same as mesh of the expression
75 if (e.coordinate_element_hash() != mesh.geometry().cmaps().front().hash())
76 {
77 throw std::runtime_error(
78 "Expression was created on a different mesh. Cannot tabulate.");
79 }
80 auto [X, Xshape] = e.X();
81 impl::tabulate_expression(values, e.kernel(), Xshape, e.value_size(), coeffs,
82 constants, mesh, entities, element);
83}
84
101template <dolfinx::scalar T, std::floating_point U>
102void tabulate_expression(std::span<T> values, const fem::Expression<T, U>& e,
103 const mesh::Mesh<U>& mesh, fem::MDSpan2 auto entities)
104{
105 // Check that domain is the same as mesh of the expression
106 if (e.coordinate_element_hash() != mesh.geometry().cmaps().front().hash())
107 {
108 throw std::runtime_error(
109 "Expression was created on a different mesh. Cannot tabulate.");
110 }
111
112 std::optional<
113 std::pair<std::reference_wrapper<const FiniteElement<U>>, std::size_t>>
114 element = std::nullopt;
115 if (auto V = e.argument_space(); V)
116 {
117 std::size_t num_argument_dofs
118 = V->dofmap()->element_dof_layout().num_dofs() * V->dofmap()->bs();
119 assert(V->element());
120 element = {std::cref(*V->element()), num_argument_dofs};
121 }
122
123 std::vector<int> coffsets = e.coefficient_offsets();
124 const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients
125 = e.coefficients();
126 std::vector<T> coeffs(entities.extent(0) * coffsets.back());
127 int cstride = coffsets.back();
128 {
129 std::vector<std::reference_wrapper<const Function<T, U>>> c;
130 std::ranges::transform(coefficients, std::back_inserter(c),
131 [](auto c) -> const Function<T, U>& { return *c; });
132 fem::pack_coefficients(c, mesh, entities, e.entity_maps(), coffsets,
133 std::span(coeffs));
134 }
135 std::vector<T> constants = fem::pack_constants(e);
136
138 values, e, md::mdspan(coeffs.data(), entities.extent(0), cstride),
139 std::span<const T>(constants), mesh, entities, element);
140}
141
142// -- Helper functions -----------------------------------------------------
143
145template <dolfinx::scalar T>
146std::map<std::pair<IntegralType, int>, std::pair<std::span<const T>, int>>
147make_coefficients_span(const std::map<std::pair<IntegralType, int>,
148 std::pair<std::vector<T>, int>>& coeffs)
149{
150 using Key = typename std::remove_reference_t<decltype(coeffs)>::key_type;
151 std::map<Key, std::pair<std::span<const T>, int>> c;
152 std::ranges::transform(
153 coeffs, std::inserter(c, c.end()),
154 [](auto& e) -> typename decltype(c)::value_type
155 { return {e.first, {e.second.first, e.second.second}}; });
156 return c;
157}
158
159// -- Scalar ----------------------------------------------------------------
160
172template <dolfinx::scalar T, std::floating_point U>
174 const Form<T, U>& M, std::span<const T> constants,
175 const std::map<std::pair<IntegralType, int>,
176 std::pair<std::span<const T>, int>>& coefficients)
177{
178 using mdspanx3_t
179 = md::mdspan<const scalar_value_t<T>,
180 md::extents<std::size_t, md::dynamic_extent, 3>>;
181
182 std::shared_ptr<const mesh::Mesh<U>> mesh = M.mesh();
183 assert(mesh);
184 std::span x = mesh->geometry().x();
185
186 // Accumulate contributions from each cell type
187 const int num_cell_types = mesh->topology()->cell_types().size();
188 T val = 0;
189 for (int cell_type_idx = 0; cell_type_idx < num_cell_types; ++cell_type_idx)
190 {
191 // Geometry dofmap and data
192 md::mdspan<const std::int32_t, md::dextents<std::size_t, 2>> x_dofmap
193 = mesh->geometry().dofmaps().at(cell_type_idx);
194 if constexpr (std::is_same_v<U, scalar_value_t<T>>)
195 {
196 val += impl::assemble_scalar(M, x_dofmap,
197 mdspanx3_t(x.data(), x.size() / 3, 3),
198 constants, coefficients, cell_type_idx);
199 }
200 else
201 {
202 std::vector<scalar_value_t<T>> _x(x.begin(), x.end());
203 val += impl::assemble_scalar(M, x_dofmap,
204 mdspanx3_t(_x.data(), _x.size() / 3, 3),
205 constants, coefficients, cell_type_idx);
206 }
207 }
208 return val;
209}
210
218template <dolfinx::scalar T, std::floating_point U>
220{
221 const std::vector<T> constants = pack_constants(M);
222 auto coefficients = allocate_coefficient_storage(M);
223 pack_coefficients(M, coefficients);
224 return assemble_scalar(M, std::span(constants),
225 make_coefficients_span(coefficients));
226}
227
228// -- Vectors ----------------------------------------------------------------
229
240// template <dolfinx::scalar T, std::floating_point U>
241template <typename V, std::floating_point U,
242 dolfinx::scalar T = typename std::remove_cvref_t<V>::value_type>
243 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
245 V&& b, const Form<T, U>& L, std::span<const T> constants,
246 const std::map<std::pair<IntegralType, int>,
247 std::pair<std::span<const T>, int>>& coefficients)
248{
249 impl::assemble_vector(b, L, constants, coefficients);
250}
251
256// template <dolfinx::scalar T, std::floating_point U>
257// void assemble_vector(std::span<T> b, const Form<T, U>& L)
258template <typename V, std::floating_point U,
259 dolfinx::scalar T = typename std::remove_cvref_t<V>::value_type>
260 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
261void assemble_vector(V&& b, const Form<T, U>& L)
262{
263 auto coefficients = allocate_coefficient_storage(L);
264 pack_coefficients(L, coefficients);
265 const std::vector<T> constants = pack_constants(L);
266 assemble_vector(b, L, std::span(constants),
267 make_coefficients_span(coefficients));
268}
269
347template <typename V,
348 std::floating_point U
349 = scalar_value_t<typename std::remove_cvref_t<V>::value_type>,
350 dolfinx::scalar T = typename std::remove_cvref_t<V>::value_type>
351 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
353 V&& b,
354 std::vector<std::optional<std::reference_wrapper<const Form<T, U>>>> a,
355 const std::vector<std::span<const T>>& constants,
356 const std::vector<std::map<std::pair<IntegralType, int>,
357 std::pair<std::span<const T>, int>>>& coeffs,
358 const std::vector<
359 std::vector<std::reference_wrapper<const DirichletBC<T, U>>>>& bcs1,
360 const std::vector<std::span<const T>>& x0, T alpha)
361{
362 // If all forms are null, there is nothing to do
363 if (std::ranges::all_of(a, [](auto ai) { return !ai; }))
364 return;
365
366 common::Timer t("[Apply lifting]");
367 impl::apply_lifting(b, a, constants, coeffs, bcs1, x0, alpha);
368}
369
398template <typename V,
399 std::floating_point U
400 = scalar_value_t<typename std::remove_cvref_t<V>::value_type>,
401 dolfinx::scalar T = typename std::remove_cvref_t<V>::value_type>
402 requires std::is_same_v<typename std::remove_cvref_t<V>::value_type, T>
404 V&& b,
405 std::vector<std::optional<std::reference_wrapper<const Form<T, U>>>> a,
406 const std::vector<
407 std::vector<std::reference_wrapper<const DirichletBC<T, U>>>>& bcs1,
408 const std::vector<std::span<const T>>& x0, T alpha)
409{
410 std::vector<
411 std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, int>>>
412 coeffs;
413 std::vector<std::vector<T>> constants;
414 for (auto _a : a)
415 {
416 if (_a)
417 {
418 auto coefficients = allocate_coefficient_storage(_a->get());
419 pack_coefficients(_a->get(), coefficients);
420 coeffs.push_back(coefficients);
421 constants.push_back(pack_constants(_a->get()));
422 }
423 else
424 {
425 coeffs.emplace_back();
426 constants.emplace_back();
427 }
428 }
429
430 std::vector<std::span<const T>> _constants(constants.begin(),
431 constants.end());
432 std::vector<std::map<std::pair<IntegralType, int>,
433 std::pair<std::span<const T>, int>>>
434 _coeffs;
435 std::ranges::transform(coeffs, std::back_inserter(_coeffs),
436 [](auto& c) { return make_coefficients_span(c); });
437
438 apply_lifting(b, a, _constants, _coeffs, bcs1, x0, alpha);
439}
440
441// -- Matrices ---------------------------------------------------------------
442
460template <dolfinx::scalar T, std::floating_point U, bool LiftingMode = false>
462 la::MatSet<T> auto mat_add, const Form<T, U>& a,
463 std::span<const T> constants,
464 const std::map<std::pair<IntegralType, int>,
465 std::pair<std::span<const T>, int>>& coefficients,
466 std::span<const std::int8_t> dof_marker0,
467 std::span<const std::int8_t> dof_marker1)
468
469{
470 common::Timer t_assm("[Assemble Matrix]");
471 using mdspanx3_t
472 = md::mdspan<const scalar_value_t<T>,
473 md::extents<std::size_t, md::dynamic_extent, 3>>;
474
475 std::shared_ptr<const mesh::Mesh<U>> mesh = a.mesh();
476 assert(mesh);
477 std::span x = mesh->geometry().x();
478 if constexpr (std::is_same_v<U, scalar_value_t<T>>)
479 {
480 impl::assemble_matrix<T, U, LiftingMode>(
481 mat_add, a, mdspanx3_t(x.data(), x.size() / 3, 3), constants,
482 coefficients, dof_marker0, dof_marker1);
483 }
484 else
485 {
486 std::vector<scalar_value_t<T>> _x(x.begin(), x.end());
487 impl::assemble_matrix<T, U, LiftingMode>(
488 mat_add, a, mdspanx3_t(_x.data(), _x.size() / 3, 3), constants,
489 coefficients, dof_marker0, dof_marker1);
490 }
491}
492
500template <dolfinx::scalar T, std::floating_point U>
502 auto mat_add, const Form<T, U>& a, std::span<const T> constants,
503 const std::map<std::pair<IntegralType, int>,
504 std::pair<std::span<const T>, int>>& coefficients,
505 const std::vector<std::reference_wrapper<const DirichletBC<T, U>>>& bcs)
506{
507 // Index maps for dof ranges
508 // NOTE: For mixed-topology meshes, there will be multiple DOF maps,
509 // but the index maps are the same.
510 auto map0 = a.function_spaces().at(0)->dofmaps().front()->index_map;
511 auto map1 = a.function_spaces().at(1)->dofmaps().front()->index_map;
512 auto bs0 = a.function_spaces().at(0)->dofmaps().front()->index_map_bs();
513 auto bs1 = a.function_spaces().at(1)->dofmaps().front()->index_map_bs();
514
515 // Build dof markers
516 std::vector<std::int8_t> dof_marker0, dof_marker1;
517 assert(map0);
518 std::int32_t dim0 = bs0 * (map0->size_local() + map0->num_ghosts());
519 assert(map1);
520 std::int32_t dim1 = bs1 * (map1->size_local() + map1->num_ghosts());
521 for (std::size_t k = 0; k < bcs.size(); ++k)
522 {
523 assert(bcs[k].get().function_space());
524 if (a.function_spaces().at(0)->contains(*bcs[k].get().function_space()))
525 {
526 dof_marker0.resize(dim0, false);
527 bcs[k].get().mark_dofs(dof_marker0);
528 }
529
530 if (a.function_spaces().at(1)->contains(*bcs[k].get().function_space()))
531 {
532 dof_marker1.resize(dim1, false);
533 bcs[k].get().mark_dofs(dof_marker1);
534 }
535 }
536
537 // Assemble
538 assemble_matrix(mat_add, a, constants, coefficients, dof_marker0,
539 dof_marker1);
540}
541
547template <dolfinx::scalar T, std::floating_point U>
549 auto mat_add, const Form<T, U>& a,
550 const std::vector<std::reference_wrapper<const DirichletBC<T, U>>>& bcs)
551{
552 // Prepare constants and coefficients
553 const std::vector<T> constants = pack_constants(a);
554 auto coefficients = allocate_coefficient_storage(a);
555 pack_coefficients(a, coefficients);
556
557 // Assemble
558 assemble_matrix(mat_add, a, std::span(constants),
559 make_coefficients_span(coefficients), bcs);
560}
561
573template <dolfinx::scalar T, std::floating_point U>
574void assemble_matrix(auto mat_add, const Form<T, U>& a,
575 std::span<const std::int8_t> dof_marker0,
576 std::span<const std::int8_t> dof_marker1)
577
578{
579 // Prepare constants and coefficients
580 const std::vector<T> constants = pack_constants(a);
581 auto coefficients = allocate_coefficient_storage(a);
582 pack_coefficients(a, coefficients);
583
584 // Assemble
585 assemble_matrix(mat_add, a, std::span(constants),
586 make_coefficients_span(coefficients), dof_marker0,
587 dof_marker1);
588}
589
602template <dolfinx::scalar T>
603void set_diagonal(auto set_fn, std::span<const std::int32_t> rows,
604 T diagonal = 1.0)
605{
606 for (std::size_t i = 0; i < rows.size(); ++i)
607 {
608 std::span diag_span(&diagonal, 1);
609 set_fn(rows.subspan(i, 1), rows.subspan(i, 1), diag_span);
610 }
611}
612
629template <dolfinx::scalar T, std::floating_point U>
631 auto set_fn, const FunctionSpace<U>& V,
632 const std::vector<std::reference_wrapper<const DirichletBC<T, U>>>& bcs,
633 T diagonal = 1.0)
634{
635 spdlog::debug("Set diagonal");
636 for (auto& bc : bcs)
637 {
638 if (V.contains(*bc.get().function_space()))
639 {
640 const auto [dofs, range] = bc.get().dof_indices();
641 set_diagonal(set_fn, dofs.first(range), diagonal);
642 }
643 }
644}
645
646} // namespace dolfinx::fem
Timer for measuring and logging elapsed time durations.
Definition Timer.h:40
Definition DirichletBC.h:258
An Expression represents a mathematical expression evaluated at a pre-defined points on a reference c...
Definition Expression.h:43
std::pair< std::vector< geometry_type >, std::array< std::size_t, 2 > > X() const
Evaluation point coordinates on the reference cell.
Definition Expression.h:165
const std::vector< std::shared_ptr< const Function< scalar_type, geometry_type > > > & coefficients() const
Expression coefficients.
Definition Expression.h:114
std::shared_ptr< const FunctionSpace< geometry_type > > argument_space() const
Argument function space.
Definition Expression.h:105
std::uint64_t coordinate_element_hash() const
Hash for coordinate element used to create the expression.
Definition Expression.h:178
const std::function< void(scalar_type *, const scalar_type *, const scalar_type *, const geometry_type *, const int *, const uint8_t *, void *)> & kernel() const
Function for tabulating the Expression.
Definition Expression.h:149
std::vector< int > coefficient_offsets() const
Offset for each coefficient expansion array on a cell.
Definition Expression.h:133
int value_size() const
Value size of the Expression result.
Definition Expression.h:155
const std::vector< std::reference_wrapper< const dolfinx::mesh::EntityMap > > & entity_maps() const
Maps between entities of different meshes.
Definition Expression.h:172
Model of a finite element.
Definition FiniteElement.h:57
A representation of finite element variational forms.
Definition Form.h:117
std::shared_ptr< const mesh::Mesh< geometry_type > > mesh() const
Common mesh for the form (the 'integration domain').
Definition Form.h:362
const std::vector< std::shared_ptr< const FunctionSpace< geometry_type > > > & function_spaces() const
Function spaces for all arguments.
Definition Form.h:370
This class represents a finite element function space defined by a mesh, a finite element,...
Definition FunctionSpace.h:34
bool contains(const FunctionSpace &V) const
Check whether V is subspace of this, or this itself.
Definition FunctionSpace.h:153
Definition Function.h:47
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition Mesh.h:23
Concept for mdspan of rank 1 or 2.
Definition traits.h:36
Matrix accumulate/set concept for functions that can be used in assemblers to accumulate or set value...
Definition utils.h:28
Definition types.h:26
Functions supporting finite element method operations.
Finite element method functionality.
Definition assemble_expression_impl.h:23
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 pack.h:172
T assemble_scalar(const Form< T, U > &M, std::span< const T > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > &coefficients)
Assemble functional into scalar.
Definition assembler.h:173
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:603
void tabulate_expression(std::span< T > values, const fem::Expression< T, U > &e, md::mdspan< const T, md::dextents< std::size_t, 2 > > coeffs, std::span< const T > constants, const mesh::Mesh< U > &mesh, fem::MDSpan2 auto entities, std::optional< std::pair< std::reference_wrapper< const FiniteElement< U > >, std::size_t > > element)
Evaluate an Expression on cells or facets.
Definition assembler.h:65
void apply_lifting(V &&b, std::vector< std::optional< std::reference_wrapper< const Form< T, U > > > > a, const std::vector< std::span< const T > > &constants, const std::vector< std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > > &coeffs, const std::vector< std::vector< std::reference_wrapper< const DirichletBC< T, U > > > > &bcs1, const std::vector< std::span< const T > > &x0, T alpha)
Modify the right-hand side vector to account for constraints (Dirichlet boundary condition constraint...
Definition assembler.h:352
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:147
void pack_coefficients(const Form< T, U > &form, std::map< std::pair< IntegralType, int >, std::pair< std::vector< T >, int > > &coeffs)
Pack coefficients of a Form.
Definition pack.h:224
void assemble_vector(V &&b, const Form< T, U > &L, std::span< const T > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > &coefficients)
Assemble linear form into a vector.
Definition assembler.h:244
std::vector< T > pack_constants(std::vector< std::reference_wrapper< const fem::Constant< T > > > c)
Pack constants of an Expression or Form into a single array ready for assembly.
Definition pack.h:515
void assemble_matrix(la::MatSet< T > auto mat_add, const Form< T, U > &a, std::span< const T > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > &coefficients, std::span< const std::int8_t > dof_marker0, std::span< const std::int8_t > dof_marker1)
Assemble bilinear form into a matrix. Matrix must already be initialised. Does not zero or finalise t...
Definition assembler.h:461
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32
Functions supporting the packing of coefficient data.