DOLFINx 0.11.0
DOLFINx C++
Loading...
Searching...
No Matches
pack.h
Go to the documentation of this file.
1// Copyright (C) 2013-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 "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 <format>
21#include <ranges>
22#include <span>
23#include <stdexcept>
24#include <type_traits>
25#include <vector>
26
29
30namespace dolfinx::fem
31{
32template <dolfinx::scalar T, std::floating_point U>
33class Expression;
34
35namespace impl
36{
38template <dolfinx::scalar T, std::floating_point U>
39std::span<const std::uint32_t>
40get_cell_orientation_info(const Function<T, U>& coefficient)
41{
42 std::span<const std::uint32_t> cell_info;
43 auto element = coefficient.function_space()->element();
44 assert(element);
45 if (element->needs_dof_transformations())
46 {
47 auto mesh = coefficient.function_space()->mesh();
48 mesh->topology_mutable()->create_entity_permutations();
49 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
50 }
51
52 return cell_info;
53}
54
56template <int _bs, dolfinx::scalar T>
57void pack_impl(std::span<T> coeffs, std::int32_t cell, int bs,
58 std::span<const T> v, std::span<const std::uint32_t> cell_info,
59 const DofMap& dofmap, auto transform)
60{
61 std::span<const std::int32_t> dofs = dofmap.cell_dofs(cell);
62 for (std::size_t i = 0; i < dofs.size(); ++i)
63 {
64 if constexpr (_bs < 0)
65 {
66 const int pos_c = bs * i;
67 const int pos_v = bs * dofs[i];
68 for (int k = 0; k < bs; ++k)
69 coeffs[pos_c + k] = v[pos_v + k];
70 }
71 else
72 {
73 assert(_bs == bs);
74 const int pos_c = _bs * i;
75 const int pos_v = _bs * dofs[i];
76 for (int k = 0; k < _bs; ++k)
77 coeffs[pos_c + k] = v[pos_v + k];
78 }
79 }
80
81 transform(coeffs, cell_info, cell, 1);
82}
83
94template <dolfinx::scalar T, std::floating_point U>
95void pack_coefficient_entity(std::span<T> c, int cstride,
96 const Function<T, U>& u,
97 std::span<const std::uint32_t> cell_info,
98 auto cells, std::int32_t offset)
99{
100 static_assert(cells.rank() == 1);
101
102 // Read data from coefficient Function u
103 std::span<const T> v = u.x()->array();
104 const DofMap& dofmap = *u.function_space()->dofmap();
105 auto element = u.function_space()->element();
106 assert(element);
107 int space_dim = element->space_dimension();
108
109 // Transformation from conforming degrees-of-freedom to reference
110 // degrees-of-freedom
111 auto transformation
112 = element->template dof_transformation_fn<T>(doftransform::transpose);
113 const int bs = dofmap.bs();
114 switch (bs)
115 {
116 case 1:
117 for (std::size_t e = 0; e < cells.extent(0); ++e)
118 {
119 if (std::int32_t cell = cells(e); cell >= 0)
120 {
121 auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
122 pack_impl<1>(cell_coeff, cell, bs, v, cell_info, dofmap,
123 transformation);
124 }
125 }
126 break;
127 case 2:
128 for (std::size_t e = 0; e < cells.extent(0); ++e)
129 {
130 if (std::int32_t cell = cells(e); cell >= 0)
131 {
132 auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
133 pack_impl<2>(cell_coeff, cell, bs, v, cell_info, dofmap,
134 transformation);
135 }
136 }
137 break;
138 case 3:
139 for (std::size_t e = 0; e < cells.extent(0); ++e)
140 {
141 if (std::int32_t cell = cells(e); cell >= 0)
142 {
143 auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
144 pack_impl<3>(cell_coeff, cell, bs, v, cell_info, dofmap,
145 transformation);
146 }
147 }
148 break;
149 default:
150 for (std::size_t e = 0; e < cells.extent(0); ++e)
151 {
152 if (std::int32_t cell = cells(e); cell >= 0)
153 {
154 auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
155 pack_impl<-1>(cell_coeff, cell, bs, v, cell_info, dofmap,
156 transformation);
157 }
158 }
159 break;
160 }
161}
162} // namespace impl
163
170template <dolfinx::scalar T, std::floating_point U>
171std::pair<std::vector<T>, int>
173 int id)
174{
175 std::size_t num_entities = 0;
176 int cstride = 0;
177 if (const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients
178 = form.coefficients();
179 !coefficients.empty())
180 {
181 const std::vector<int> offsets = form.coefficient_offsets();
182 cstride = offsets.back();
183 num_entities = form.domain(integral_type, id, 0).size();
184 if (integral_type != IntegralType::cell)
185 num_entities /= 2;
186 }
187
188 return {std::vector<T>(num_entities * cstride), cstride};
189}
190
195template <dolfinx::scalar T, std::floating_point U>
196std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, int>>
198{
199 std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, int>> coeffs;
200 for (fem::IntegralType type : form.integral_types())
201 {
202 for (int i = 0; i < form.num_integrals(type, 0); ++i)
203 {
204 coeffs.emplace_hint(coeffs.end(), std::pair{type, i},
205 allocate_coefficient_storage(form, type, i));
206 }
207 }
208
209 return coeffs;
210}
211
223template <dolfinx::scalar T, std::floating_point U>
225 std::map<std::pair<IntegralType, int>,
226 std::pair<std::vector<T>, int>>& coeffs)
227{
228 const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients
229 = form.coefficients();
230 const std::vector<int> offsets = form.coefficient_offsets();
231
232 for (auto& [intergal_data, coeff_data] : coeffs)
233 {
234 auto [integral_type, id] = intergal_data;
235 std::vector<T>& c = coeff_data.first;
236 int cstride = coeff_data.second;
237 if (!coefficients.empty())
238 {
239 switch (integral_type)
240 {
242 {
243 // Iterate over coefficients that are active in cell integrals
244 for (int coeff : form.active_coeffs(IntegralType::cell, id))
245 {
246 // Get coefficient mesh
247 auto mesh = coefficients[coeff]->function_space()->mesh();
248 assert(mesh);
249
250 // Other integrals in the form might have coefficients defined
251 // over entities of codim > 0, which don't make sense for cell
252 // integrals, so don't pack them.
253 if (int codim
254 = form.mesh()->topology()->dim() - mesh->topology()->dim();
255 codim > 0)
256 {
257 throw std::runtime_error("Should not be packing coefficients with "
258 "codim>0 in a cell integral");
259 }
260
261 std::span<const std::int32_t> cells_b
262 = form.domain_coeff(IntegralType::cell, id, coeff);
263 md::mdspan cells(cells_b.data(), cells_b.size());
264 std::span<const std::uint32_t> cell_info
265 = impl::get_cell_orientation_info(*coefficients[coeff]);
266 impl::pack_coefficient_entity(std::span(c), cstride,
267 *coefficients[coeff], cell_info, cells,
268 offsets[coeff]);
269 }
270 break;
271 }
273 {
274 // Iterate over coefficients that are active in interior
275 // facet integrals
276 for (int coeff : form.active_coeffs(IntegralType::interior_facet, id))
277 {
278 auto mesh = coefficients[coeff]->function_space()->mesh();
279 std::span<const std::int32_t> facets_b
281 md::mdspan<const std::int32_t,
282 md::extents<std::size_t, md::dynamic_extent, 4>>
283 facets(facets_b.data(), facets_b.size() / 4, 4);
284
285 std::span<const std::uint32_t> cell_info
286 = impl::get_cell_orientation_info(*coefficients[coeff]);
287
288 // Pack coefficient ['+']
289 auto cells0 = md::submdspan(facets, md::full_extent, 0);
290 impl::pack_coefficient_entity(std::span(c), 2 * cstride,
291 *coefficients[coeff], cell_info, cells0,
292 2 * offsets[coeff]);
293
294 // Pack coefficient ['-']
295 auto cells1 = md::submdspan(facets, md::full_extent, 2);
296 impl::pack_coefficient_entity(std::span(c), 2 * cstride,
297 *coefficients[coeff], cell_info, cells1,
298 offsets[coeff] + offsets[coeff + 1]);
299 }
300 break;
301 }
305 {
306 // Iterate over coefficients that are active in vertex integrals
307 for (int coeff : form.active_coeffs(integral_type, id))
308 {
309 // Get coefficient mesh
310 auto mesh = coefficients[coeff]->function_space()->mesh();
311 assert(mesh);
312
313 std::span<const std::int32_t> entitites_b
314 = form.domain_coeff(integral_type, id, coeff);
315 md::mdspan<const std::int32_t,
316 md::extents<std::size_t, md::dynamic_extent, 2>>
317 entities(entitites_b.data(), entitites_b.size() / 2, 2);
318 std::span<const std::uint32_t> cell_info
319 = impl::get_cell_orientation_info(*coefficients[coeff]);
320 impl::pack_coefficient_entity(
321 std::span(c), cstride, *coefficients[coeff], cell_info,
322 md::submdspan(entities, md::full_extent, 0), offsets[coeff]);
323 }
324 break;
325 }
326 default:
327 throw std::runtime_error(
328 "Could not pack coefficient. Integral type not supported.");
329 }
330 }
331 }
332}
333
344template <dolfinx::scalar T, std::floating_point U>
346 const fem::Function<T, U>& coeff, const mesh::Mesh<U>& mesh,
347 fem::MDSpan2 auto entities,
348 std::optional<std::reference_wrapper<const dolfinx::mesh::EntityMap>>
349 entity_map)
350{
351 auto mesh_c = coeff.function_space()->mesh();
352 assert(mesh_c);
353
354 auto span_to_vector = [](auto entities)
355 {
356 assert(entities.rank() == 1);
357
358 std::vector<std::int32_t> vec;
359 vec.reserve(entities.extent(0));
360 for (std::size_t i = 0; i < entities.extent(0); ++i)
361 vec.push_back(entities[i]);
362 return vec;
363 };
364
365 if (mesh_c->topology() == mesh.topology())
366 {
367 // If same mesh no mapping is needed
368 if constexpr (entities.rank() == 1)
369 return span_to_vector(entities);
370
371 else
372 // If (cell, local_index) pairs are given, extract the cells
373 return span_to_vector(md::submdspan(entities, md::full_extent, 0));
374 }
375 else
376 {
377 assert(entity_map.has_value());
378 const mesh::Topology& topology = *mesh.topology();
379 int tdim = topology.dim();
380 int codim = tdim - mesh_c->topology()->dim();
381 const dolfinx::mesh::EntityMap& emap = entity_map.value().get();
382 bool inverse = emap.sub_topology() == mesh_c->topology();
383 // If cells are supplied on the parent mesh, we can directly map them to
384 // cells on the coefficient mesh.
385 if constexpr (entities.rank() == 1)
386 {
387 assert(codim == 0);
388
389 return emap.sub_topology_to_topology(span_to_vector(entities), inverse);
390 }
391 else if constexpr (entities.rank() == 2)
392 {
393 if (codim == 0)
394 {
395 // If codim is zero we extract the cells and map them
396 auto cells = md::submdspan(entities, md::full_extent, 0);
397 return emap.sub_topology_to_topology(span_to_vector(cells), inverse);
398 }
399 else
400 {
401 // Any other codim needs to map (cell, local index) to facets and then
402 // to cells of the submesh
403 if (!inverse)
404 {
405 throw std::runtime_error(
406 "Unsupported mapping. Can only map from submesh to parent mesh.");
407 }
408 assert(codim > 0);
409 auto c_to_e = topology.connectivity(tdim, tdim - codim);
410 if (!c_to_e)
411 {
412 throw std::runtime_error(std::format(
413 "Topology connectivity from codim {} to {} not found.", tdim,
414 tdim - codim));
415 }
416 // Map parent (cell, local_index) to parent facet
417 std::vector<std::int32_t> contiguous_cells;
418 contiguous_cells.reserve(entities.extent(0));
419 for (std::size_t e = 0; e < entities.extent(0); ++e)
420 {
421 auto pair = md::submdspan(entities, e, md::full_extent);
422 contiguous_cells.push_back(c_to_e->links(pair[0])[pair[1]]);
423 }
424 // Map parent facet to submesh cell
425 return emap.sub_topology_to_topology(contiguous_cells, inverse);
426 }
427 }
428 }
429}
430
445template <dolfinx::scalar T, std::floating_point U>
447 std::vector<std::reference_wrapper<const Function<T, U>>> coeffs,
448 const mesh::Mesh<U>& mesh, fem::MDSpan2 auto entities,
449 const std::vector<std::reference_wrapper<const dolfinx::mesh::EntityMap>>&
450 entity_maps,
451 std::span<const int> offsets, std::span<T> c)
452{
453
454 assert(!offsets.empty());
455 const int cstride = offsets.back();
456
457 if (c.size() < entities.extent(0) * offsets.back())
458 throw std::runtime_error("Coefficient packing span is too small.");
459
460 // Iterate over coefficients
461 for (std::size_t coeff = 0; coeff < coeffs.size(); ++coeff)
462 {
463
464 // Get mesh of coefficient and check if entity map is required
465 auto mesh_c = coeffs[coeff].get().function_space()->mesh();
466 std::vector<std::int32_t> coefficient_cells;
467 if (mesh_c->topology() == mesh.topology())
468 {
469 coefficient_cells = extract_coefficient_cells_from_entities(
470 coeffs[coeff].get(), mesh, entities, std::nullopt);
471 }
472 else
473 {
474 // Helper function to get correct entity map
475 auto get_entity_map
476 = [mesh, &entity_maps](auto& mesh0) -> const mesh::EntityMap&
477 {
478 auto it = std::ranges::find_if(
479 entity_maps,
480 [mesh, mesh0](const mesh::EntityMap& em)
481 {
482 return ((em.topology() == mesh0->topology()
483 and em.sub_topology() == mesh.topology()))
484 or ((em.sub_topology() == mesh0->topology()
485 and em.topology() == mesh.topology()));
486 });
487
488 if (it == entity_maps.end())
489 {
490 throw std::runtime_error("Incompatible mesh. argument "
491 "entity_maps must be provided.");
492 }
493 return *it;
494 };
495 // Find correct entity map and determine direction of the map
496 const mesh::EntityMap& emap = get_entity_map(mesh_c);
497 coefficient_cells = extract_coefficient_cells_from_entities(
498 coeffs[coeff].get(), mesh, entities,
499 std::reference_wrapper<const mesh::EntityMap>(emap));
500 }
501
502 std::span<const std::uint32_t> cell_info
503 = impl::get_cell_orientation_info(coeffs[coeff].get());
504 md::mdspan cells(coefficient_cells.data(), coefficient_cells.size());
505 impl::pack_coefficient_entity(std::span(c), cstride, coeffs[coeff].get(),
506 cell_info, cells, offsets[coeff]);
507 }
508}
509
513template <typename T>
514std::vector<T>
515pack_constants(std::vector<std::reference_wrapper<const fem::Constant<T>>> c)
516{
517 // Calculate size of array needed to store packed constants
518 std::int32_t size = std::accumulate(
519 c.cbegin(), c.cend(), 0, [](std::int32_t sum, auto& constant)
520 { return sum + constant.get().value.size(); });
521
522 // Pack constants
523 std::vector<T> constant_values(size);
524 std::int32_t offset = 0;
525 for (auto& constant : c)
526 {
527 std::ranges::copy(constant.get().value,
528 std::next(constant_values.begin(), offset));
529 offset += constant.get().value.size();
530 }
531
532 return constant_values;
533}
534
539template <typename U>
540 requires std::convertible_to<
542 typename std::decay_t<U>::geometry_type>>
543 or std::convertible_to<
545 typename std::decay_t<U>::geometry_type>>
546std::vector<typename U::scalar_type> pack_constants(const U& u)
547{
548 using T = typename std::decay_t<U>::scalar_type;
549 std::vector<std::reference_wrapper<const Constant<T>>> c;
550 std::ranges::transform(u.constants(), std::back_inserter(c),
551 [](auto c) -> const Constant<T>& { return *c; });
552 return fem::pack_constants(c);
553}
554
555} // namespace dolfinx::fem
Degree-of-freedom map representations and tools.
Constant (in space) value which can be attached to a Form.
Definition Constant.h:22
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:170
An Expression represents a mathematical expression evaluated at a pre-defined points on a reference c...
Definition Expression.h:43
A representation of finite element variational forms.
Definition Form.h:117
int num_integrals(IntegralType type, int kernel_idx) const
Get number of integrals (kernels) for a given integral type and kernel index.
Definition Form.h:440
const std::vector< std::shared_ptr< const Function< scalar_type, geometry_type > > > & coefficients() const
Access coefficients.
Definition Form.h:564
std::vector< int > coefficient_offsets() const
Offset for each coefficient expansion array on a cell.
Definition Form.h:578
std::span< const std::int32_t > domain_coeff(IntegralType type, int idx, int c) const
Coefficient function mesh integration entity indices.
Definition Form.h:551
std::vector< int > active_coeffs(IntegralType type, int id) const
Indices of coefficients that are active for a given integral (kernel).
Definition Form.h:414
std::shared_ptr< const mesh::Mesh< geometry_type > > mesh() const
Common mesh for the form (the 'integration domain').
Definition Form.h:362
std::set< IntegralType > integral_types() const
Get types of integrals in the form.
Definition Form.h:396
std::span< const std::int32_t > domain(IntegralType type, int idx, int kernel_idx) const
Mesh entity indices to integrate over for a given integral (kernel).
Definition Form.h:483
Definition Function.h:47
std::shared_ptr< const FunctionSpace< geometry_type > > function_space() const
Access the function space.
Definition Function.h:147
std::shared_ptr< const la::Vector< value_type > > x() const
Underlying vector (const version).
Definition Function.h:153
A bidirectional map relating entities in one topology to another.
Definition EntityMap.h:21
std::vector< std::int32_t > sub_topology_to_topology(CellRange auto &&entities, bool inverse) const
Map entities between the sub-topology and the parent topology.
Definition EntityMap.h:103
std::shared_ptr< const Topology > sub_topology() const
Get the sub-topology.
Definition EntityMap.cpp:23
std::shared_ptr< const Topology > topology() const
Get the (parent) topology.
Definition EntityMap.cpp:18
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition Mesh.h:23
Topology stores the topology of a mesh, consisting of mesh entities and connectivity (incidence relat...
Definition Topology.h:49
std::shared_ptr< const graph::AdjacencyList< std::int32_t > > connectivity(std::array< int, 2 > d0, std::array< int, 2 > d1) const
Get the connectivity from entities of topological dimension d0 to dimension d1.
Definition Topology.cpp:865
int dim() const noexcept
Topological dimension of the mesh.
Definition Topology.cpp:803
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:172
@ transpose
Transpose.
Definition FiniteElement.h:28
@ inverse
Inverse.
Definition FiniteElement.h:29
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
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
std::vector< std::int32_t > extract_coefficient_cells_from_entities(const fem::Function< T, U > &coeff, const mesh::Mesh< U > &mesh, fem::MDSpan2 auto entities, std::optional< std::reference_wrapper< const dolfinx::mesh::EntityMap > > entity_map)
Given a Function and a related mesh and its integration entities, extract the cell indices of the coe...
Definition pack.h:345
IntegralType
Type of integral.
Definition Form.h:39
@ vertex
Vertex.
Definition Form.h:43
@ interior_facet
Interior facet.
Definition Form.h:42
@ ridge
Ridge.
Definition Form.h:44
@ cell
Cell.
Definition Form.h:40
@ exterior_facet
Facet.
Definition Form.h:41
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:57
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:95