11#include "FiniteElement.h"
14#include "FunctionSpace.h"
17#include <basix/mdspan.hpp>
19#include <dolfinx/mesh/Topology.h>
32template <dolfinx::scalar T, std::
floating_po
int U>
38template <dolfinx::scalar T, std::
floating_po
int U>
39std::span<const std::uint32_t>
40get_cell_orientation_info(
const Function<T, U>& coefficient)
42 std::span<const std::uint32_t> cell_info;
43 auto element = coefficient.function_space()->element();
45 if (element->needs_dof_transformations())
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());
56template <
int _bs, dolfinx::scalar T>
58 std::span<const T> v, std::span<const std::uint32_t> cell_info,
59 const DofMap& dofmap,
auto transform)
62 for (std::size_t i = 0; i < dofs.size(); ++i)
64 if constexpr (_bs < 0)
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];
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];
81 transform(coeffs, cell_info,
cell, 1);
94template <dolfinx::scalar T, std::
floating_po
int U>
97 std::span<const std::uint32_t> cell_info,
98 auto cells, std::int32_t offset)
100 static_assert(cells.rank() == 1);
103 std::span<const T> v = u.
x()->array();
107 int space_dim = element->space_dimension();
113 const int bs = dofmap.
bs();
117 for (std::size_t e = 0; e < cells.extent(0); ++e)
119 if (std::int32_t
cell = cells(e);
cell >= 0)
121 auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
128 for (std::size_t e = 0; e < cells.extent(0); ++e)
130 if (std::int32_t
cell = cells(e);
cell >= 0)
132 auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
139 for (std::size_t e = 0; e < cells.extent(0); ++e)
141 if (std::int32_t
cell = cells(e);
cell >= 0)
143 auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
150 for (std::size_t e = 0; e < cells.extent(0); ++e)
152 if (std::int32_t
cell = cells(e);
cell >= 0)
154 auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
170template <dolfinx::scalar T, std::
floating_po
int U>
171std::pair<std::vector<T>,
int>
175 std::size_t num_entities = 0;
177 if (
const std::vector<std::shared_ptr<
const Function<T, U>>>& coefficients
179 !coefficients.empty())
182 cstride = offsets.back();
183 num_entities = form.
domain(integral_type,
id, 0).size();
188 return {std::vector<T>(num_entities * cstride), cstride};
195template <dolfinx::scalar T, std::
floating_po
int U>
196std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>,
int>>
199 std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>,
int>> coeffs;
204 coeffs.emplace_hint(coeffs.end(), std::pair{type, i},
223template <dolfinx::scalar T, std::
floating_po
int U>
225 std::map<std::pair<IntegralType, int>,
226 std::pair<std::vector<T>,
int>>& coeffs)
228 const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients
232 for (
auto& [intergal_data, coeff_data] : coeffs)
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())
239 switch (integral_type)
247 auto mesh = coefficients[coeff]->function_space()->mesh();
254 = form.
mesh()->topology()->dim() -
mesh->topology()->dim();
257 throw std::runtime_error(
"Should not be packing coefficients with "
258 "codim>0 in a cell integral");
261 std::span<const std::int32_t> cells_b
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,
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);
285 std::span<const std::uint32_t> cell_info
286 = impl::get_cell_orientation_info(*coefficients[coeff]);
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,
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]);
310 auto mesh = coefficients[coeff]->function_space()->mesh();
313 std::span<const std::int32_t> entitites_b
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]);
327 throw std::runtime_error(
328 "Could not pack coefficient. Integral type not supported.");
344template <dolfinx::scalar T, std::
floating_po
int U>
348 std::optional<std::reference_wrapper<const dolfinx::mesh::EntityMap>>
354 auto span_to_vector = [](
auto entities)
356 assert(entities.rank() == 1);
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]);
365 if (mesh_c->topology() ==
mesh.topology())
368 if constexpr (entities.rank() == 1)
369 return span_to_vector(entities);
373 return span_to_vector(md::submdspan(entities, md::full_extent, 0));
377 assert(entity_map.has_value());
379 int tdim = topology.
dim();
380 int codim = tdim - mesh_c->topology()->dim();
385 if constexpr (entities.rank() == 1)
391 else if constexpr (entities.rank() == 2)
396 auto cells = md::submdspan(entities, md::full_extent, 0);
405 throw std::runtime_error(
406 "Unsupported mapping. Can only map from submesh to parent mesh.");
409 auto c_to_e = topology.
connectivity(tdim, tdim - codim);
412 throw std::runtime_error(std::format(
413 "Topology connectivity from codim {} to {} not found.", tdim,
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)
421 auto pair = md::submdspan(entities, e, md::full_extent);
422 contiguous_cells.push_back(c_to_e->links(pair[0])[pair[1]]);
445template <dolfinx::scalar T, std::
floating_po
int U>
449 const std::vector<std::reference_wrapper<const dolfinx::mesh::EntityMap>>&
451 std::span<const int> offsets, std::span<T> c)
454 assert(!offsets.empty());
455 const int cstride = offsets.back();
457 if (c.size() < entities.extent(0) * offsets.back())
458 throw std::runtime_error(
"Coefficient packing span is too small.");
461 for (std::size_t coeff = 0; coeff < coeffs.size(); ++coeff)
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())
470 coeffs[coeff].get(),
mesh, entities, std::nullopt);
478 auto it = std::ranges::find_if(
482 return ((em.
topology() == mesh0->topology()
488 if (it == entity_maps.end())
490 throw std::runtime_error(
"Incompatible mesh. argument "
491 "entity_maps must be provided.");
498 coeffs[coeff].get(),
mesh, entities,
499 std::reference_wrapper<const mesh::EntityMap>(emap));
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]);
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(); });
523 std::vector<T> constant_values(size);
524 std::int32_t offset = 0;
525 for (
auto& constant : c)
527 std::ranges::copy(constant.get().value,
528 std::next(constant_values.begin(), offset));
529 offset += constant.get().value.size();
532 return constant_values;
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>>
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),
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
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