DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
pack.h
Go to the documentation of this file.
1// Copyright (C) 2013-2025 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 "Constant.h"
10#include "DofMap.h"
11#include "FiniteElement.h"
12#include "Form.h"
13#include "Function.h"
14#include "FunctionSpace.h"
15#include "traits.h"
16#include <array>
17#include <basix/mdspan.hpp>
18#include <concepts>
19#include <dolfinx/mesh/Topology.h>
20#include <span>
21#include <stdexcept>
22#include <type_traits>
23#include <vector>
24
27
28namespace dolfinx::fem
29{
30template <dolfinx::scalar T, std::floating_point U>
31class Expression;
32
33namespace impl
34{
36template <dolfinx::scalar T, std::floating_point U>
37std::span<const std::uint32_t>
38get_cell_orientation_info(const Function<T, U>& coefficient)
39{
40 std::span<const std::uint32_t> cell_info;
41 auto element = coefficient.function_space()->element();
42 assert(element);
43 if (element->needs_dof_transformations())
44 {
45 auto mesh = coefficient.function_space()->mesh();
46 mesh->topology_mutable()->create_entity_permutations();
47 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
48 }
49
50 return cell_info;
51}
52
54template <int _bs, dolfinx::scalar T>
55void pack_impl(std::span<T> coeffs, std::int32_t cell, int bs,
56 std::span<const T> v, std::span<const std::uint32_t> cell_info,
57 const DofMap& dofmap, auto transform)
58{
59 std::span<const std::int32_t> dofs = dofmap.cell_dofs(cell);
60 for (std::size_t i = 0; i < dofs.size(); ++i)
61 {
62 if constexpr (_bs < 0)
63 {
64 const int pos_c = bs * i;
65 const int pos_v = bs * dofs[i];
66 for (int k = 0; k < bs; ++k)
67 coeffs[pos_c + k] = v[pos_v + k];
68 }
69 else
70 {
71 assert(_bs == bs);
72 const int pos_c = _bs * i;
73 const int pos_v = _bs * dofs[i];
74 for (int k = 0; k < _bs; ++k)
75 coeffs[pos_c + k] = v[pos_v + k];
76 }
77 }
78
79 transform(coeffs, cell_info, cell, 1);
80}
81
92template <dolfinx::scalar T, std::floating_point U>
93void pack_coefficient_entity(std::span<T> c, int cstride,
94 const Function<T, U>& u,
95 std::span<const std::uint32_t> cell_info,
96 auto cells, std::int32_t offset)
97{
98 static_assert(cells.rank() == 1);
99
100 // Read data from coefficient Function u
101 std::span<const T> v = u.x()->array();
102 const DofMap& dofmap = *u.function_space()->dofmap();
103 auto element = u.function_space()->element();
104 assert(element);
105 int space_dim = element->space_dimension();
106
107 // Transformation from conforming degrees-of-freedom to reference
108 // degrees-of-freedom
109 auto transformation
110 = element->template dof_transformation_fn<T>(doftransform::transpose);
111 const int bs = dofmap.bs();
112 switch (bs)
113 {
114 case 1:
115 for (std::size_t e = 0; e < cells.extent(0); ++e)
116 {
117 if (std::int32_t cell = cells(e); cell >= 0)
118 {
119 auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
120 pack_impl<1>(cell_coeff, cell, bs, v, cell_info, dofmap,
121 transformation);
122 }
123 }
124 break;
125 case 2:
126 for (std::size_t e = 0; e < cells.extent(0); ++e)
127 {
128 if (std::int32_t cell = cells(e); cell >= 0)
129 {
130 auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
131 pack_impl<2>(cell_coeff, cell, bs, v, cell_info, dofmap,
132 transformation);
133 }
134 }
135 break;
136 case 3:
137 for (std::size_t e = 0; e < cells.extent(0); ++e)
138 {
139 if (std::int32_t cell = cells(e); cell >= 0)
140 {
141 auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
142 pack_impl<3>(cell_coeff, cell, bs, v, cell_info, dofmap,
143 transformation);
144 }
145 }
146 break;
147 default:
148 for (std::size_t e = 0; e < cells.extent(0); ++e)
149 {
150 if (std::int32_t cell = cells(e); cell >= 0)
151 {
152 auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
153 pack_impl<-1>(cell_coeff, cell, bs, v, cell_info, dofmap,
154 transformation);
155 }
156 }
157 break;
158 }
159}
160} // namespace impl
161
168template <dolfinx::scalar T, std::floating_point U>
169std::pair<std::vector<T>, int>
171 int id)
172{
173 std::size_t num_entities = 0;
174 int cstride = 0;
175 if (const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients
176 = form.coefficients();
177 !coefficients.empty())
178 {
179 const std::vector<int> offsets = form.coefficient_offsets();
180 cstride = offsets.back();
181 num_entities = form.domain(integral_type, id, 0).size();
182 if (integral_type == IntegralType::exterior_facet
183 or integral_type == IntegralType::interior_facet)
184 {
185 num_entities /= 2;
186 }
187 }
188
189 return {std::vector<T>(num_entities * cstride), cstride};
190}
191
196template <dolfinx::scalar T, std::floating_point U>
197std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, int>>
199{
200 std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, int>> coeffs;
201 for (fem::IntegralType type : form.integral_types())
202 {
203 for (int id : form.integral_ids(type))
204 {
205 coeffs.emplace_hint(coeffs.end(), std::pair{type, id},
206 allocate_coefficient_storage(form, type, id));
207 }
208 }
209
210 return coeffs;
211}
212
224template <dolfinx::scalar T, std::floating_point U>
226 std::map<std::pair<IntegralType, int>,
227 std::pair<std::vector<T>, int>>& coeffs)
228{
229 const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients
230 = form.coefficients();
231 const std::vector<int> offsets = form.coefficient_offsets();
232
233 for (auto& [intergal_data, coeff_data] : coeffs)
234 {
235 auto [integral_type, id] = intergal_data;
236 std::vector<T>& c = coeff_data.first;
237 int cstride = coeff_data.second;
238 if (!coefficients.empty())
239 {
240 switch (integral_type)
241 {
243 {
244 // Iterate over coefficients that are active in cell integrals
245 for (int coeff : form.active_coeffs(IntegralType::cell, id))
246 {
247 // Get coefficient mesh
248 auto mesh = coefficients[coeff]->function_space()->mesh();
249 assert(mesh);
250
251 // Other integrals in the form might have coefficients defined
252 // over entities of codim > 0, which don't make sense for cell
253 // integrals, so don't pack them.
254 if (int codim
255 = form.mesh()->topology()->dim() - mesh->topology()->dim();
256 codim > 0)
257 {
258 throw std::runtime_error("Should not be packing coefficients with "
259 "codim>0 in a cell integral");
260 }
261
262 std::span<const std::int32_t> cells_b
263 = form.domain_coeff(IntegralType::cell, id, coeff);
264 md::mdspan cells(cells_b.data(), cells_b.size());
265 std::span<const std::uint32_t> cell_info
266 = impl::get_cell_orientation_info(*coefficients[coeff]);
267 impl::pack_coefficient_entity(std::span(c), cstride,
268 *coefficients[coeff], cell_info, cells,
269 offsets[coeff]);
270 }
271 break;
272 }
274 {
275 // Iterate over coefficients coefficients that are active in
276 // exterior facet integrals
277 for (int coeff : form.active_coeffs(IntegralType::exterior_facet, id))
278 {
279 auto mesh = coefficients[coeff]->function_space()->mesh();
280 std::span<const std::int32_t> facets_b
282 md::mdspan<const std::int32_t,
283 md::extents<std::size_t, md::dynamic_extent, 2>>
284 facets(facets_b.data(), facets_b.size() / 2, 2);
285 auto cells = md::submdspan(facets, md::full_extent, 0);
286
287 std::span<const std::uint32_t> cell_info
288 = impl::get_cell_orientation_info(*coefficients[coeff]);
289 impl::pack_coefficient_entity(std::span(c), cstride,
290 *coefficients[coeff], cell_info, cells,
291 offsets[coeff]);
292 }
293 break;
294 }
296 {
297 // Iterate over coefficients that are active in interior
298 // facet integrals
299 for (int coeff : form.active_coeffs(IntegralType::interior_facet, id))
300 {
301 auto mesh = coefficients[coeff]->function_space()->mesh();
302 std::span<const std::int32_t> facets_b
304 md::mdspan<const std::int32_t,
305 md::extents<std::size_t, md::dynamic_extent, 4>>
306 facets(facets_b.data(), facets_b.size() / 4, 4);
307
308 std::span<const std::uint32_t> cell_info
309 = impl::get_cell_orientation_info(*coefficients[coeff]);
310
311 // Pack coefficient ['+']
312 auto cells0 = md::submdspan(facets, md::full_extent, 0);
313 impl::pack_coefficient_entity(std::span(c), 2 * cstride,
314 *coefficients[coeff], cell_info, cells0,
315 2 * offsets[coeff]);
316
317 // Pack coefficient ['-']
318 auto cells1 = md::submdspan(facets, md::full_extent, 2);
319 impl::pack_coefficient_entity(std::span(c), 2 * cstride,
320 *coefficients[coeff], cell_info, cells1,
321 offsets[coeff] + offsets[coeff + 1]);
322 }
323 break;
324 }
325 default:
326 throw std::runtime_error(
327 "Could not pack coefficient. Integral type not supported.");
328 }
329 }
330 }
331}
332
342template <dolfinx::scalar T, std::floating_point U>
344 std::vector<std::reference_wrapper<const Function<T, U>>> coeffs,
345 std::span<const int> offsets, fem::MDSpan2 auto entities, std::span<T> c)
346{
347 assert(!offsets.empty());
348 const int cstride = offsets.back();
349
350 if (c.size() < entities.extent(0) * offsets.back())
351 throw std::runtime_error("Coefficient packing span is too small.");
352
353 // Iterate over coefficients
354 for (std::size_t coeff = 0; coeff < coeffs.size(); ++coeff)
355 {
356 std::span<const std::uint32_t> cell_info
357 = impl::get_cell_orientation_info(coeffs[coeff].get());
358
359 if constexpr (entities.rank() == 1)
360 {
361 impl::pack_coefficient_entity(std::span(c), cstride, coeffs[coeff].get(),
362 cell_info, entities, offsets[coeff]);
363 }
364 else
365 {
366 auto cells = md::submdspan(entities, md::full_extent, 0);
367 impl::pack_coefficient_entity(std::span(c), cstride, coeffs[coeff].get(),
368 cell_info, cells, offsets[coeff]);
369 }
370 }
371}
372
377template <typename T>
378std::vector<T>
379pack_constants(std::vector<std::reference_wrapper<const fem::Constant<T>>> c)
380{
381 // Calculate size of array needed to store packed constants
382 std::int32_t size = std::accumulate(
383 c.cbegin(), c.cend(), 0, [](std::int32_t sum, auto& constant)
384 { return sum + constant.get().value.size(); });
385
386 // Pack constants
387 std::vector<T> constant_values(size);
388 std::int32_t offset = 0;
389 for (auto& constant : c)
390 {
391 std::ranges::copy(constant.get().value,
392 std::next(constant_values.begin(), offset));
393 offset += constant.get().value.size();
394 }
395
396 return constant_values;
397}
398
403template <typename U>
404 requires std::convertible_to<
406 typename std::decay_t<U>::geometry_type>>
407 or std::convertible_to<
409 typename std::decay_t<U>::geometry_type>>
410std::vector<typename U::scalar_type> pack_constants(const U& u)
411{
412 using T = typename std::decay_t<U>::scalar_type;
413 std::vector<std::reference_wrapper<const Constant<T>>> c;
414 std::ranges::transform(u.constants(), std::back_inserter(c),
415 [](auto c) -> const Constant<T>& { return *c; });
416 return fem::pack_constants(c);
417}
418
419} // namespace dolfinx::fem
Degree-of-freedom map representations and tools.
Constant value which can be attached to a Form.
Definition Constant.h:23
Degree-of-freedom map.
Definition DofMap.h:73
std::span< const std::int32_t > cell_dofs(std::int32_t c) const
Local-to-global mapping of dofs on a cell.
Definition DofMap.h:127
int bs() const noexcept
Return the block size for the dofmap.
Definition DofMap.cpp:168
An Expression represents a mathematical expression evaluated at a pre-defined points on a reference c...
Definition Expression.h:41
A representation of finite element variational forms.
Definition Form.h:175
const std::vector< std::shared_ptr< const Function< scalar_type, geometry_type > > > & coefficients() const
Access coefficients.
Definition Form.h:575
std::vector< int > coefficient_offsets() const
Offset for each coefficient expansion array on a cell.
Definition Form.h:589
std::vector< int > active_coeffs(IntegralType type, int id) const
Indices of coefficients that are active for a given integral (kernel).
Definition Form.h:422
std::shared_ptr< const mesh::Mesh< geometry_type > > mesh() const
Common mesh for the form (the 'integration domain').
Definition Form.h:373
std::span< const std::int32_t > domain_coeff(IntegralType type, int id, int c) const
Coefficient function mesh integration entity indices.
Definition Form.h:556
std::vector< int > integral_ids(IntegralType type) const
Get the IDs for integrals (kernels) for given integral domain type.
Definition Form.h:444
std::set< IntegralType > integral_types() const
Get types of integrals in the form.
Definition Form.h:404
std::span< const std::int32_t > domain(IntegralType type, int id, int kernel_idx) const
Mesh entity indices to integrate over for a given integral (kernel).
Definition Form.h:486
Definition Function.h:47
Concept for mdspan of rank 1 or 2.
Definition traits.h:36
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:170
@ transpose
Transpose.
Definition FiniteElement.h:28
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:225
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:379
IntegralType
Type of integral.
Definition Form.h:36
@ interior_facet
Interior facet.
Definition Form.h:39
@ cell
Cell.
Definition Form.h:37
@ exterior_facet
Exterior facet.
Definition Form.h:38
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32
void pack_impl(std::span< T > coeffs, std::int32_t cell, int bs, std::span< const T > v, std::span< const std::uint32_t > cell_info, const DofMap &dofmap, auto transform)
Pack a single coefficient for a single cell.
Definition pack.h:55
void pack_coefficient_entity(std::span< T > c, int cstride, const Function< T, U > &u, std::span< const std::uint32_t > cell_info, auto cells, std::int32_t offset)
Pack a single coefficient for a set of active entities.
Definition pack.h:93