10#include "FunctionSpace.h"
13#include <basix/mdspan.hpp>
16#include <dolfinx/common/IndexMap.h>
17#include <dolfinx/common/types.h>
18#include <dolfinx/mesh/EntityMap.h>
19#include <dolfinx/mesh/Mesh.h>
32template <dolfinx::scalar T>
34template <dolfinx::scalar T, std::
floating_po
int U>
49template <dolfinx::scalar T, std::
floating_po
int U = scalar_value_t<T>>
57 template <
typename K,
typename V,
typename W>
58 requires std::is_convertible_v<
59 std::remove_cvref_t<K>,
60 std::function<void(T*,
const T*,
const T*,
const U*,
61 const int*,
const uint8_t*,
void*)>>
62 and std::is_convertible_v<std::remove_cvref_t<V>,
63 std::vector<std::int32_t>>
64 and std::is_convertible_v<std::remove_cvref_t<W>,
73 std::function<void(T*,
const T*,
const T*,
const U*,
const int*,
74 const uint8_t*,
void*)>
114template <dolfinx::scalar T, std::
floating_po
int U = dolfinx::scalar_value_t<T>>
151 template <
typename X>
152 requires std::is_convertible_v<
153 std::remove_cvref_t<X>,
154 std::map<std::tuple<IntegralType, int, int>,
165 const std::vector<std::reference_wrapper<const mesh::EntityMap>>&
167 : _function_spaces(V), _integrals(std::forward<X>(integrals)),
172 throw std::runtime_error(
"Form Mesh is null.");
178 auto it = std::ranges::find_if(
182 return ((em.
topology() == mesh0->topology()
188 if (it == entity_maps.end())
190 throw std::runtime_error(
191 "Incompatible mesh. argument entity_maps must be provided.");
199 auto compute_facet_domains
200 = [&](
const auto& int_ents_mesh,
int codim,
const auto& c_to_f,
201 const auto& emap,
bool inverse)
208 std::vector<std::int32_t> entities;
209 entities.reserve(int_ents_mesh.size() / 2);
215 for (std::size_t i = 0; i < int_ents_mesh.size(); i += 2)
216 entities.push_back(int_ents_mesh[i]);
223 for (std::size_t i = 0; i < int_ents_mesh.size(); i += 2)
226 c_to_f->links(int_ents_mesh[i])[int_ents_mesh[i + 1]]);
230 throw std::runtime_error(
"Codimension > 1 not supported.");
234 std::vector<std::int32_t> cells_mesh0
235 = emap.sub_topology_to_topology(entities,
inverse);
241 std::vector<std::int32_t> e = int_ents_mesh;
242 for (std::size_t i = 0; i < cells_mesh0.size(); ++i)
243 e[2 * i] = cells_mesh0[i];
248 for (
auto& space : _function_spaces)
251 std::map<std::tuple<IntegralType, int, int>,
252 std::variant<std::vector<std::int32_t>,
253 std::span<const std::int32_t>>>
256 if (
auto mesh0 = space->mesh(); mesh0 == _mesh)
258 for (
auto& [key, integral] : _integrals)
259 vdata.insert({key, std::span(integral.entities)});
270 for (
auto& [key, itg] : _integrals)
272 auto [type, id, kernel_idx] = key;
273 std::vector<std::int32_t> e;
280 int tdim = topology.
dim();
282 int codim = tdim - mesh0->topology()->dim();
287 e = compute_facet_domains(itg.entities, codim, c_to_f, emap,
291 throw std::runtime_error(
"Integral type not supported.");
293 vdata.insert({key, std::move(e)});
297 _edata.push_back(vdata);
300 for (
auto& [key, integral] : _integrals)
302 auto [type, id, kernel_idx] = key;
303 for (
int c : integral.coeffs)
305 if (
auto mesh0 =
coefficients.at(c)->function_space()->mesh();
308 _cdata.insert({{type, id, c}, std::span(integral.entities)});
316 std::vector<std::int32_t> e;
323 int tdim = topology.
dim();
325 int codim = tdim - mesh0->topology()->dim();
329 e = compute_facet_domains(integral.entities, codim, c_to_f, emap,
333 throw std::runtime_error(
"Integral type not supported.");
335 _cdata.insert({{type, id, c}, std::move(e)});
355 int rank()
const {
return _function_spaces.size(); }
359 std::shared_ptr<const mesh::Mesh<geometry_type>>
mesh()
const
366 const std::vector<std::shared_ptr<const FunctionSpace<geometry_type>>>&
369 return _function_spaces;
382 auto it = _integrals.find({type, id, kernel_idx});
383 if (it == _integrals.end())
384 throw std::runtime_error(
"Requested integral kernel not found.");
385 return it->second.kernel;
392 std::vector<IntegralType> set_data;
393 std::ranges::transform(_integrals, std::back_inserter(set_data),
394 [](
auto& x) {
return std::get<0>(x.first); });
395 return std::set<IntegralType>(set_data.begin(), set_data.end());
410 auto it = std::ranges::find_if(_integrals,
413 auto [t, id_, kernel_idx] = x.first;
414 return t == type and id_ == id;
416 if (it == _integrals.end())
417 throw std::runtime_error(
"Could not find active coefficient list.");
418 return it->second.coeffs;
432 std::vector<int> ids;
433 for (
auto& [key, integral] : _integrals)
435 auto [t, id, kernel_idx] = key;
442 std::sort(ids.begin(), ids.end());
443 auto it = std::unique(ids.begin(), ids.end());
444 ids.erase(it, ids.end());
473 int kernel_idx)
const
475 auto it = _integrals.find({type, id, kernel_idx});
476 if (it == _integrals.end())
477 throw std::runtime_error(
"Requested domain not found.");
478 return it->second.entities;
516 int kernel_idx)
const
518 auto it = _edata.at(
rank).find({type, id, kernel_idx});
519 if (it == _edata.at(
rank).end())
520 throw std::runtime_error(
"Requested domain for argument not found.");
523 return std::get<std::span<const std::int32_t>>(it->second);
525 catch (std::bad_variant_access& e)
527 return std::get<std::vector<std::int32_t>>(it->second);
548 auto it = _cdata.find({type, id, c});
549 if (it == _cdata.end())
550 throw std::runtime_error(
"No domain for requested integral.");
553 return std::get<std::span<const std::int32_t>>(it->second);
555 catch (std::bad_variant_access& e)
557 return std::get<std::vector<std::int32_t>>(it->second);
563 std::shared_ptr<const Function<scalar_type, geometry_type>>>&
566 return _coefficients;
580 std::vector<int> n{0};
581 for (
auto& c : _coefficients)
584 throw std::runtime_error(
"Not all form coefficients have been set.");
585 n.push_back(n.back() + c->function_space()->element()->space_dimension());
591 const std::vector<std::shared_ptr<const Constant<scalar_type>>>&
599 std::vector<std::shared_ptr<const FunctionSpace<geometry_type>>>
603 std::map<std::tuple<IntegralType, int, int>,
608 std::shared_ptr<const mesh::Mesh<geometry_type>> _mesh;
611 std::vector<std::shared_ptr<const Function<scalar_type, geometry_type>>>
615 std::vector<std::shared_ptr<const Constant<scalar_type>>> _constants;
618 bool _needs_facet_permutations;
630 std::vector<std::map<
631 std::tuple<IntegralType, int, int>,
632 std::variant<std::vector<std::int32_t>, std::span<const std::int32_t>>>>
646 std::tuple<IntegralType, int, int>,
647 std::variant<std::vector<std::int32_t>, std::span<const std::int32_t>>>
Constant (in space) value which can be attached to a Form.
Definition Constant.h:22
This class represents a finite element function space defined by a mesh, a finite element,...
Definition FunctionSpace.h:34
A bidirectional map relating entities in one topology to another.
Definition EntityMap.h:20
std::vector< std::int32_t > sub_topology_to_topology(std::span< const std::int32_t > entities, bool inverse) const
Map entities between the sub-topology and the parent topology.
Definition EntityMap.cpp:30
std::shared_ptr< const Topology > sub_topology() const
Get the sub-topology.
Definition EntityMap.cpp:24
std::shared_ptr< const Topology > topology() const
Get the (parent) topology.
Definition EntityMap.cpp:19
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:41
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:832
int dim() const noexcept
Topological dimension of the mesh.
Definition Topology.cpp:785
Finite element method functionality.
Definition assemble_expression_impl.h:23
@ inverse
Inverse.
Definition FiniteElement.h:29
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
Represents integral data, containing the kernel, and a list of entities to integrate over and the ind...
Definition Form.h:51
std::function< void(T *, const T *, const T *, const U *, const int *, const uint8_t *, void *)> kernel
The integration kernel.
Definition Form.h:75
integral_data(K &&kernel, V &&entities, W &&coeffs)
Create a structure to hold integral data.
Definition Form.h:66
std::vector< int > coeffs
Indices of coefficients (from the form) that are in this integral.
Definition Form.h:83
std::vector< std::int32_t > entities
The entities to integrate over for this integral. These are the entities in 'full' mesh.
Definition Form.h:79