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>
50template <dolfinx::scalar T, std::
floating_po
int U = scalar_value_t<T>>
58 template <
typename K,
typename V,
typename W>
59 requires std::is_convertible_v<
60 std::remove_cvref_t<K>,
61 std::function<void(T*,
const T*,
const T*,
const U*,
62 const int*,
const uint8_t*,
void*)>>
63 and std::is_convertible_v<std::remove_cvref_t<V>,
64 std::vector<std::int32_t>>
65 and std::is_convertible_v<std::remove_cvref_t<W>,
74 std::function<void(T*,
const T*,
const T*,
const U*,
const int*,
75 const uint8_t*,
void*)>
115template <dolfinx::scalar T, std::
floating_po
int U = dolfinx::scalar_value_t<T>>
155 template <
typename X>
156 requires std::is_convertible_v<
157 std::remove_cvref_t<X>,
158 std::map<std::tuple<IntegralType, int, int>,
169 const std::vector<std::reference_wrapper<const mesh::EntityMap>>&
171 : _function_spaces(V), _integrals(std::forward<X>(integrals)),
176 throw std::runtime_error(
"Form Mesh is null.");
182 auto it = std::ranges::find_if(
186 return ((em.
topology() == mesh0->topology()
192 if (it == entity_maps.end())
194 throw std::runtime_error(
195 "Incompatible mesh. argument entity_maps must be provided.");
203 auto compute_facet_domains
204 = [&](
const auto& int_ents_mesh,
int codim,
const auto& c_to_f,
205 const auto& emap,
bool inverse)
212 std::vector<std::int32_t> entities;
213 entities.reserve(int_ents_mesh.size() / 2);
219 for (std::size_t i = 0; i < int_ents_mesh.size(); i += 2)
220 entities.push_back(int_ents_mesh[i]);
227 for (std::size_t i = 0; i < int_ents_mesh.size(); i += 2)
230 c_to_f->links(int_ents_mesh[i])[int_ents_mesh[i + 1]]);
234 throw std::runtime_error(
"Codimension > 1 not supported.");
238 std::vector<std::int32_t> cells_mesh0
239 = emap.sub_topology_to_topology(entities,
inverse);
245 std::vector<std::int32_t> e = int_ents_mesh;
246 for (std::size_t i = 0; i < cells_mesh0.size(); ++i)
247 e[2 * i] = cells_mesh0[i];
252 for (
auto& space : _function_spaces)
255 std::map<std::tuple<IntegralType, int, int>,
256 std::variant<std::vector<std::int32_t>,
257 std::span<const std::int32_t>>>
260 if (
auto mesh0 = space->mesh(); mesh0 == _mesh)
262 for (
auto& [key, integral] : _integrals)
263 vdata.insert({key, std::span(integral.entities)});
274 for (
auto& [key, itg] : _integrals)
276 auto [type, idx, kernel_idx] = key;
277 std::vector<std::int32_t> e;
284 int tdim = topology.
dim();
286 int codim = tdim - mesh0->topology()->dim();
291 e = compute_facet_domains(itg.entities, codim, c_to_f, emap,
295 throw std::runtime_error(
"Integral type not supported.");
297 vdata.insert({key, std::move(e)});
301 _edata.push_back(vdata);
304 for (
auto& [key, integral] : _integrals)
306 auto [type, idx, kernel_idx] = key;
307 for (
int c : integral.coeffs)
309 if (
auto mesh0 =
coefficients.at(c)->function_space()->mesh();
312 _cdata.insert({{type, idx, c}, std::span(integral.entities)});
320 std::vector<std::int32_t> e;
327 int tdim = topology.
dim();
329 int codim = tdim - mesh0->topology()->dim();
333 e = compute_facet_domains(integral.entities, codim, c_to_f, emap,
337 throw std::runtime_error(
"Integral type not supported.");
338 _cdata.insert({{type, idx, c}, std::move(e)});
358 int rank()
const {
return _function_spaces.size(); }
362 std::shared_ptr<const mesh::Mesh<geometry_type>>
mesh()
const
369 const std::vector<std::shared_ptr<const FunctionSpace<geometry_type>>>&
372 return _function_spaces;
388 auto it = _integrals.find({type, id, kernel_idx});
389 if (it == _integrals.end())
390 throw std::runtime_error(
"Requested integral kernel not found.");
391 return it->second.kernel;
398 std::vector<IntegralType> set_data;
399 std::ranges::transform(_integrals, std::back_inserter(set_data),
400 [](
auto& x) {
return std::get<0>(x.first); });
401 return std::set<IntegralType>(set_data.begin(), set_data.end());
416 auto it = std::ranges::find_if(_integrals,
419 auto [t, id_, kernel_idx] = x.first;
420 return t == type and id_ == id;
422 if (it == _integrals.end())
423 throw std::runtime_error(
"Could not find active coefficient list.");
424 return it->second.coeffs;
443 int count = std::count_if(_integrals.begin(), _integrals.end(),
444 [type, kernel_idx](
auto& x)
446 auto [t, id, k_idx] = x.first;
447 return t == type and k_idx == kernel_idx;
484 int kernel_idx)
const
486 auto it = _integrals.find({type, idx, kernel_idx});
487 if (it == _integrals.end())
488 throw std::runtime_error(
"Requested domain not found.");
489 return it->second.entities;
527 int kernel_idx)
const
529 auto it = _edata.at(
rank).find({type, idx, kernel_idx});
530 if (it == _edata.at(
rank).end())
531 throw std::runtime_error(
"Requested domain for argument not found.");
534 return std::get<std::span<const std::int32_t>>(it->second);
536 catch (std::bad_variant_access& e)
538 return std::get<std::vector<std::int32_t>>(it->second);
559 auto it = _cdata.find({type, idx, c});
560 if (it == _cdata.end())
561 throw std::runtime_error(
"No domain for requested integral.");
564 return std::get<std::span<const std::int32_t>>(it->second);
566 catch (std::bad_variant_access& e)
568 return std::get<std::vector<std::int32_t>>(it->second);
574 std::shared_ptr<const Function<scalar_type, geometry_type>>>&
577 return _coefficients;
591 std::vector<int> n{0};
592 for (
auto& c : _coefficients)
595 throw std::runtime_error(
"Not all form coefficients have been set.");
596 n.push_back(n.back() + c->function_space()->element()->space_dimension());
602 const std::vector<std::shared_ptr<const Constant<scalar_type>>>&
610 std::vector<std::shared_ptr<const FunctionSpace<geometry_type>>>
614 std::map<std::tuple<IntegralType, int, int>,
619 std::shared_ptr<const mesh::Mesh<geometry_type>> _mesh;
622 std::vector<std::shared_ptr<const Function<scalar_type, geometry_type>>>
626 std::vector<std::shared_ptr<const Constant<scalar_type>>> _constants;
629 bool _needs_facet_permutations;
641 std::vector<std::map<
642 std::tuple<IntegralType, int, int>,
643 std::variant<std::vector<std::int32_t>, std::span<const std::int32_t>>>>
657 std::tuple<IntegralType, int, int>,
658 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:840
int dim() const noexcept
Topological dimension of the mesh.
Definition Topology.cpp:786
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
@ ridge
Ridge.
Definition Form.h:44
@ cell
Cell.
Definition Form.h:40
@ exterior_facet
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:52
std::function< void(T *, const T *, const T *, const U *, const int *, const uint8_t *, void *)> kernel
The integration kernel.
Definition Form.h:76
integral_data(K &&kernel, V &&entities, W &&coeffs)
Create a structure to hold integral data.
Definition Form.h:67
std::vector< int > coeffs
Indices of coefficients (from the form) that are in this integral.
Definition Form.h:84
std::vector< std::int32_t > entities
The entities to integrate over for this integral. These are the entities in 'full' mesh.
Definition Form.h:80