11#include "FiniteElement.h" 
   14#include "FunctionSpace.h" 
   17#include <basix/mdspan.hpp> 
   19#include <dolfinx/mesh/Topology.h> 
   30template <dolfinx::scalar T, std::
floating_po
int U>
 
   36template <dolfinx::scalar T, std::
floating_po
int U>
 
   37std::span<const std::uint32_t>
 
   38get_cell_orientation_info(
const Function<T, U>& coefficient)
 
   40  std::span<const std::uint32_t> cell_info;
 
   41  auto element = coefficient.function_space()->element();
 
   43  if (element->needs_dof_transformations())
 
   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());
 
   54template <
int _bs, dolfinx::scalar T>
 
   56               std::span<const T> v, std::span<const std::uint32_t> cell_info,
 
   57               const DofMap& dofmap, 
auto transform)
 
   60  for (std::size_t i = 0; i < dofs.size(); ++i)
 
   62    if constexpr (_bs < 0)
 
   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];
 
   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];
 
   79  transform(coeffs, cell_info, 
cell, 1);
 
 
   92template <dolfinx::scalar T, std::
floating_po
int U>
 
   95                             std::span<const std::uint32_t> cell_info,
 
   96                             auto cells, std::int32_t offset)
 
   98  static_assert(cells.rank() == 1);
 
  101  std::span<const T> v = u.x()->array();
 
  102  const DofMap& dofmap = *u.function_space()->dofmap();
 
  103  auto element = u.function_space()->element();
 
  105  int space_dim = element->space_dimension();
 
  111  const int bs = dofmap.
bs();
 
  115    for (std::size_t e = 0; e < cells.extent(0); ++e)
 
  117      if (std::int32_t 
cell = cells(e); 
cell >= 0)
 
  119        auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
 
  126    for (std::size_t e = 0; e < cells.extent(0); ++e)
 
  128      if (std::int32_t 
cell = cells(e); 
cell >= 0)
 
  130        auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
 
  137    for (std::size_t e = 0; e < cells.extent(0); ++e)
 
  139      if (std::int32_t 
cell = cells(e); 
cell >= 0)
 
  141        auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
 
  148    for (std::size_t e = 0; e < cells.extent(0); ++e)
 
  150      if (std::int32_t 
cell = cells(e); 
cell >= 0)
 
  152        auto cell_coeff = c.subspan(e * cstride + offset, space_dim);
 
 
  168template <dolfinx::scalar T, std::
floating_po
int U>
 
  169std::pair<std::vector<T>, 
int>
 
  173  std::size_t num_entities = 0;
 
  175  if (
const std::vector<std::shared_ptr<
const Function<T, U>>>& coefficients
 
  177      !coefficients.empty())
 
  180    cstride = offsets.back();
 
  181    num_entities = form.
domain(integral_type, 
id, 0).size();
 
  189  return {std::vector<T>(num_entities * cstride), cstride};
 
 
  196template <dolfinx::scalar T, std::
floating_po
int U>
 
  197std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, 
int>>
 
  200  std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, 
int>> coeffs;
 
  205      coeffs.emplace_hint(coeffs.end(), std::pair{type, i},
 
 
  224template <dolfinx::scalar T, std::
floating_po
int U>
 
  226                       std::map<std::pair<IntegralType, int>,
 
  227                                std::pair<std::vector<T>, 
int>>& coeffs)
 
  229  const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients
 
  233  for (
auto& [intergal_data, coeff_data] : coeffs)
 
  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())
 
  240      switch (integral_type)
 
  248          auto mesh = coefficients[coeff]->function_space()->mesh();
 
  255              = form.
mesh()->topology()->dim() - 
mesh->topology()->dim();
 
  258            throw std::runtime_error(
"Should not be packing coefficients with " 
  259                                     "codim>0 in a cell integral");
 
  262          std::span<const std::int32_t> cells_b
 
  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,
 
  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);
 
  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,
 
  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);
 
  308          std::span<const std::uint32_t> cell_info
 
  309              = impl::get_cell_orientation_info(*coefficients[coeff]);
 
  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,
 
  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]);
 
  331          auto mesh = coefficients[coeff]->function_space()->mesh();
 
  334          std::span<const std::int32_t> vertices_b
 
  336          md::mdspan<
const std::int32_t,
 
  337                     md::extents<std::size_t, md::dynamic_extent, 2>>
 
  338              vertices(vertices_b.data(), vertices_b.size() / 2, 2);
 
  339          std::span<const std::uint32_t> cell_info
 
  340              = impl::get_cell_orientation_info(*coefficients[coeff]);
 
  341          impl::pack_coefficient_entity(
 
  342              std::span(c), cstride, *coefficients[coeff], cell_info,
 
  343              md::submdspan(vertices, md::full_extent, 0), offsets[coeff]);
 
  348        throw std::runtime_error(
 
  349            "Could not pack coefficient. Integral type not supported.");
 
 
  364template <dolfinx::scalar T, std::
floating_po
int U>
 
  367    std::span<const int> offsets, 
fem::MDSpan2 auto entities, std::span<T> c)
 
  369  assert(!offsets.empty());
 
  370  const int cstride = offsets.back();
 
  372  if (c.size() < entities.extent(0) * offsets.back())
 
  373    throw std::runtime_error(
"Coefficient packing span is too small.");
 
  376  for (std::size_t coeff = 0; coeff < coeffs.size(); ++coeff)
 
  378    std::span<const std::uint32_t> cell_info
 
  379        = impl::get_cell_orientation_info(coeffs[coeff].get());
 
  381    if constexpr (entities.rank() == 1)
 
  383      impl::pack_coefficient_entity(std::span(c), cstride, coeffs[coeff].get(),
 
  384                                    cell_info, entities, offsets[coeff]);
 
  388      auto cells = md::submdspan(entities, md::full_extent, 0);
 
  389      impl::pack_coefficient_entity(std::span(c), cstride, coeffs[coeff].get(),
 
  390                                    cell_info, cells, offsets[coeff]);
 
 
  404  std::int32_t size = std::accumulate(
 
  405      c.cbegin(), c.cend(), 0, [](std::int32_t sum, 
auto& constant)
 
  406      { return sum + constant.get().value.size(); });
 
  409  std::vector<T> constant_values(size);
 
  410  std::int32_t offset = 0;
 
  411  for (
auto& constant : c)
 
  413    std::ranges::copy(constant.get().value,
 
  414                      std::next(constant_values.begin(), offset));
 
  415    offset += constant.get().value.size();
 
  418  return constant_values;
 
 
  426  requires std::convertible_to<
 
  428                                  typename std::decay_t<U>::geometry_type>>
 
  429           or std::convertible_to<
 
  431                            typename std::decay_t<U>::geometry_type>>
 
  434  using T = 
typename std::decay_t<U>::scalar_type;
 
  435  std::vector<std::reference_wrapper<const Constant<T>>> c;
 
  436  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:168
 
An Expression represents a mathematical expression evaluated at a pre-defined points on a reference c...
Definition Expression.h:41
 
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:401
 
IntegralType
Type of integral.
Definition Form.h:39
 
@ vertex
Vertex.
Definition Form.h:43
 
@ interior_facet
Interior facet.
Definition Form.h:42
 
@ cell
Cell.
Definition Form.h:40
 
@ exterior_facet
Exterior 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: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