DOLFINx 0.10.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Form.h
1// Copyright (C) 2019-2025 Garth N. Wells, Chris Richardson, Joseph P. Dean and
2// Jørgen S. Dokken
3//
4// This file is part of DOLFINx (https://www.fenicsproject.org)
5//
6// SPDX-License-Identifier: LGPL-3.0-or-later
7
8#pragma once
9
10#include "FunctionSpace.h"
11#include "traits.h"
12#include <algorithm>
13#include <basix/mdspan.hpp>
14#include <concepts>
15#include <cstdint>
16#include <dolfinx/common/IndexMap.h>
17#include <dolfinx/common/types.h>
18#include <dolfinx/mesh/EntityMap.h>
19#include <dolfinx/mesh/Mesh.h>
20#include <functional>
21#include <map>
22#include <memory>
23#include <optional>
24#include <ranges>
25#include <span>
26#include <tuple>
27#include <utility>
28#include <vector>
29
30namespace dolfinx::fem
31{
32template <dolfinx::scalar T>
33class Constant;
34template <dolfinx::scalar T, std::floating_point U>
35class Function;
36
38enum class IntegralType : std::int8_t
39{
40 cell = 0,
43 vertex = 3,
44 ridge = 4
45};
46
50template <dolfinx::scalar T, std::floating_point U = scalar_value_t<T>>
52{
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>,
66 std::vector<int>>
68 : kernel(std::forward<K>(kernel)), entities(std::forward<V>(entities)),
69 coeffs(std::forward<W>(coeffs))
70 {
71 }
72
74 std::function<void(T*, const T*, const T*, const U*, const int*,
75 const uint8_t*, void*)>
77
80 std::vector<std::int32_t> entities;
81
84 std::vector<int> coeffs;
85};
86
115template <dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
116class Form
117{
118public:
120 using scalar_type = T;
121
123 using geometry_type = U;
124
155 template <typename X>
156 requires std::is_convertible_v<
157 std::remove_cvref_t<X>,
158 std::map<std::tuple<IntegralType, int, int>,
161 const std::vector<std::shared_ptr<const FunctionSpace<geometry_type>>>& V,
162 X&& integrals, std::shared_ptr<const mesh::Mesh<geometry_type>> mesh,
163 const std::vector<
164 std::shared_ptr<const Function<scalar_type, geometry_type>>>&
166 const std::vector<std::shared_ptr<const Constant<scalar_type>>>&
167 constants,
169 const std::vector<std::reference_wrapper<const mesh::EntityMap>>&
170 entity_maps)
171 : _function_spaces(V), _integrals(std::forward<X>(integrals)),
172 _mesh(mesh), _coefficients(coefficients), _constants(constants),
173 _needs_facet_permutations(needs_facet_permutations)
174 {
175 if (!_mesh)
176 throw std::runtime_error("Form Mesh is null.");
177
178 // A helper function to find the correct entity map for a given mesh
179 auto get_entity_map
180 = [mesh, &entity_maps](auto& mesh0) -> const mesh::EntityMap&
181 {
182 auto it = std::ranges::find_if(
183 entity_maps,
184 [mesh, mesh0](const mesh::EntityMap& em)
185 {
186 return ((em.topology() == mesh0->topology()
187 and em.sub_topology() == mesh->topology()))
188 or ((em.sub_topology() == mesh0->topology()
189 and em.topology() == mesh->topology()));
190 });
191
192 if (it == entity_maps.end())
193 {
194 throw std::runtime_error(
195 "Incompatible mesh. argument entity_maps must be provided.");
196 }
197 return *it;
198 };
199
200 // A helper function to compute the (cell, local_facet) pairs in the
201 // argument/coefficient domain from the (cell, local_facet) pairs in
202 // `this->mesh()`.
203 auto compute_facet_domains
204 = [&](const auto& int_ents_mesh, int codim, const auto& c_to_f,
205 const auto& emap, bool inverse)
206 {
207 // TODO: This function would be much neater using
208 // `std::views::stride(2)` from C++ 23
209
210 // Get a list of entities to map to the argument/coefficient
211 // domain
212 std::vector<std::int32_t> entities;
213 entities.reserve(int_ents_mesh.size() / 2);
214 if (codim == 0)
215 {
216 // In the codim 0 case, we need to map from cells in
217 // `this->mesh()` to cells in the argument/coefficient mesh, so
218 // here we extract the cells.
219 for (std::size_t i = 0; i < int_ents_mesh.size(); i += 2)
220 entities.push_back(int_ents_mesh[i]);
221 }
222 else if (codim == 1)
223 {
224 // In the codim 1 case, we need to map facets in `this->mesh()`
225 // to cells in the argument/coefficient mesh, so here we extract
226 // the facet index using the cell-to-facet connectivity.
227 for (std::size_t i = 0; i < int_ents_mesh.size(); i += 2)
228 {
229 entities.push_back(
230 c_to_f->links(int_ents_mesh[i])[int_ents_mesh[i + 1]]);
231 }
232 }
233 else
234 throw std::runtime_error("Codimension > 1 not supported.");
235
236 // Map from entity indices in `this->mesh()` to the corresponding
237 // cell indices in the argument/coefficient mesh
238 std::vector<std::int32_t> cells_mesh0
239 = emap.sub_topology_to_topology(entities, inverse);
240
241 // Create a list of (cell, local_facet_index) pairs in the
242 // argument/coefficient domain. Since `create_submesh`preserves
243 // the local facet index (with respect to the cell), we can use
244 // the local facet indices from the input integration entities
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];
248
249 return e;
250 };
251
252 for (auto& space : _function_spaces)
253 {
254 // Working map: [integral type, integral_idx, kernel_idx]->entities
255 std::map<std::tuple<IntegralType, int, int>,
256 std::variant<std::vector<std::int32_t>,
257 std::span<const std::int32_t>>>
258 vdata;
259
260 if (auto mesh0 = space->mesh(); mesh0 == _mesh)
261 {
262 for (auto& [key, integral] : _integrals)
263 vdata.insert({key, std::span(integral.entities)});
264 }
265 else
266 {
267 // Find correct entity map
268 const mesh::EntityMap& emap = get_entity_map(mesh0);
269
270 // Determine direction of the map. We need to map from
271 // `this->mesh()` to `mesh0`, so if `emap->sub_topology()` isn't
272 // the source topology, we need the inverse map
273 bool inverse = emap.sub_topology() == mesh0->topology();
274 for (auto& [key, itg] : _integrals)
275 {
276 auto [type, idx, kernel_idx] = key;
277 std::vector<std::int32_t> e;
278 if (type == IntegralType::cell)
279 e = emap.sub_topology_to_topology(itg.entities, inverse);
280 else if (type == IntegralType::exterior_facet
282 {
283 const mesh::Topology topology = *_mesh->topology();
284 int tdim = topology.dim();
285 assert(mesh0);
286 int codim = tdim - mesh0->topology()->dim();
287 assert(codim >= 0);
288 auto c_to_f = topology.connectivity(tdim, tdim - 1);
289 assert(c_to_f);
290
291 e = compute_facet_domains(itg.entities, codim, c_to_f, emap,
292 inverse);
293 }
294 else
295 throw std::runtime_error("Integral type not supported.");
296
297 vdata.insert({key, std::move(e)});
298 }
299 }
300
301 _edata.push_back(vdata);
302 }
303
304 for (auto& [key, integral] : _integrals)
305 {
306 auto [type, idx, kernel_idx] = key;
307 for (int c : integral.coeffs)
308 {
309 if (auto mesh0 = coefficients.at(c)->function_space()->mesh();
310 mesh0 == _mesh)
311 {
312 _cdata.insert({{type, idx, c}, std::span(integral.entities)});
313 }
314 else
315 {
316 // Find correct entity map and determine direction of the map
317 const mesh::EntityMap& emap = get_entity_map(mesh0);
318 bool inverse = emap.sub_topology() == mesh0->topology();
319
320 std::vector<std::int32_t> e;
321 if (type == IntegralType::cell)
322 e = emap.sub_topology_to_topology(integral.entities, inverse);
323 else if (type == IntegralType::exterior_facet
325 {
326 const mesh::Topology topology = *_mesh->topology();
327 int tdim = topology.dim();
328 assert(mesh0);
329 int codim = tdim - mesh0->topology()->dim();
330 auto c_to_f = topology.connectivity(tdim, tdim - 1);
331 assert(c_to_f);
332
333 e = compute_facet_domains(integral.entities, codim, c_to_f, emap,
334 inverse);
335 }
336 else
337 throw std::runtime_error("Integral type not supported.");
338 _cdata.insert({{type, idx, c}, std::move(e)});
339 }
340 }
341 }
342 }
343
345 Form(const Form& form) = delete;
346
348 Form(Form&& form) = default;
349
351 virtual ~Form() = default;
352
358 int rank() const { return _function_spaces.size(); }
359
362 std::shared_ptr<const mesh::Mesh<geometry_type>> mesh() const
363 {
364 return _mesh;
365 }
366
369 const std::vector<std::shared_ptr<const FunctionSpace<geometry_type>>>&
371 {
372 return _function_spaces;
373 }
374
384 std::function<void(scalar_type*, const scalar_type*, const scalar_type*,
385 const geometry_type*, const int*, const uint8_t*, void*)>
386 kernel(IntegralType type, int id, int kernel_idx) const
387 {
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;
392 }
393
396 std::set<IntegralType> integral_types() const
397 {
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());
402 }
403
414 std::vector<int> active_coeffs(IntegralType type, int id) const
415 {
416 auto it = std::ranges::find_if(_integrals,
417 [type, id](auto& x)
418 {
419 auto [t, id_, kernel_idx] = x.first;
420 return t == type and id_ == id;
421 });
422 if (it == _integrals.end())
423 throw std::runtime_error("Could not find active coefficient list.");
424 return it->second.coeffs;
425 }
426
440 int num_integrals(IntegralType type, int kernel_idx) const
441 {
442
443 int count = std::count_if(_integrals.begin(), _integrals.end(),
444 [type, kernel_idx](auto& x)
445 {
446 auto [t, id, k_idx] = x.first;
447 return t == type and k_idx == kernel_idx;
448 });
449 return count;
450 }
451
483 std::span<const std::int32_t> domain(IntegralType type, int idx,
484 int kernel_idx) const
485 {
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;
490 }
491
526 std::span<const std::int32_t> domain_arg(IntegralType type, int rank, int idx,
527 int kernel_idx) const
528 {
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.");
532 try
533 {
534 return std::get<std::span<const std::int32_t>>(it->second);
535 }
536 catch (std::bad_variant_access& e)
537 {
538 return std::get<std::vector<std::int32_t>>(it->second);
539 }
540 }
541
556 std::span<const std::int32_t> domain_coeff(IntegralType type, int idx,
557 int c) const
558 {
559 auto it = _cdata.find({type, idx, c});
560 if (it == _cdata.end())
561 throw std::runtime_error("No domain for requested integral.");
562 try
563 {
564 return std::get<std::span<const std::int32_t>>(it->second);
565 }
566 catch (std::bad_variant_access& e)
567 {
568 return std::get<std::vector<std::int32_t>>(it->second);
569 }
570 }
571
573 const std::vector<
574 std::shared_ptr<const Function<scalar_type, geometry_type>>>&
576 {
577 return _coefficients;
578 }
579
583 bool needs_facet_permutations() const { return _needs_facet_permutations; }
584
589 std::vector<int> coefficient_offsets() const
590 {
591 std::vector<int> n{0};
592 for (auto& c : _coefficients)
593 {
594 if (!c)
595 throw std::runtime_error("Not all form coefficients have been set.");
596 n.push_back(n.back() + c->function_space()->element()->space_dimension());
597 }
598 return n;
599 }
600
602 const std::vector<std::shared_ptr<const Constant<scalar_type>>>&
603 constants() const
604 {
605 return _constants;
606 }
607
608private:
609 // Function spaces (one for each argument)
610 std::vector<std::shared_ptr<const FunctionSpace<geometry_type>>>
611 _function_spaces;
612
613 // Integrals (integral type, id, celltype)
614 std::map<std::tuple<IntegralType, int, int>,
616 _integrals;
617
618 // The mesh
619 std::shared_ptr<const mesh::Mesh<geometry_type>> _mesh;
620
621 // Form coefficients
622 std::vector<std::shared_ptr<const Function<scalar_type, geometry_type>>>
623 _coefficients;
624
625 // Constants associated with the Form
626 std::vector<std::shared_ptr<const Constant<scalar_type>>> _constants;
627
628 // True if permutation data needs to be passed into these integrals
629 bool _needs_facet_permutations;
630
631 // Mapped domain index data for argument functions.
632 //
633 // Consider:
634 //
635 // entities = this->domain(IntegralType, integral(id), kernel_idx];
636 // entities0 = _edata[0][IntegralType, integral(id), coefficient_index];
637 //
638 // Then `entities[i]` is a mesh entity index (e.g., cell index) in
639 // `_mesh`, and `entities0[i]` is the index of the same entity but in
640 // the mesh associated with the argument 0 (test function) space.
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>>>>
644 _edata;
645
646 // Mapped domain index data for coefficient functions.
647 //
648 // Consider:
649 //
650 // entities = this->domain(IntegralType, integral(id), kernel_idx];
651 // entities0 = _cdata[IntegralType, integral(id), coefficient_index];
652 //
653 // Then `entities[i]` is a mesh entity index (e.g., cell index) in
654 // `_mesh`, and `entities0[i]` is the index of the same entity but in
655 // the mesh associated with the coefficient Function.
656 std::map<
657 std::tuple<IntegralType, int, int>,
658 std::variant<std::vector<std::int32_t>, std::span<const std::int32_t>>>
659 _cdata;
660};
661} // namespace dolfinx::fem
Constant (in space) value which can be attached to a Form.
Definition Constant.h:22
U geometry_type
Geometry type.
Definition Form.h:123
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:575
bool needs_facet_permutations() const
Get bool indicating whether permutation data needs to be passed into these integrals.
Definition Form.h:583
Form(Form &&form)=default
Move constructor.
int rank() const
Rank of the form.
Definition Form.h:358
Form(const Form &form)=delete
Copy constructor.
virtual ~Form()=default
Destructor.
Form(const std::vector< std::shared_ptr< const FunctionSpace< geometry_type > > > &V, X &&integrals, std::shared_ptr< const mesh::Mesh< geometry_type > > mesh, const std::vector< std::shared_ptr< const Function< scalar_type, geometry_type > > > &coefficients, const std::vector< std::shared_ptr< const Constant< scalar_type > > > &constants, bool needs_facet_permutations, const std::vector< std::reference_wrapper< const mesh::EntityMap > > &entity_maps)
Create a finite element form.
Definition Form.h:160
const std::vector< std::shared_ptr< const Constant< scalar_type > > > & constants() const
Access constants.
Definition Form.h:603
std::function< void(scalar_type *, const scalar_type *, const scalar_type *, const geometry_type *, const int *, const uint8_t *, void *)> kernel(IntegralType type, int id, int kernel_idx) const
Get the kernel function for an integral.
Definition Form.h:386
std::vector< int > coefficient_offsets() const
Offset for each coefficient expansion array on a cell.
Definition Form.h:589
std::span< const std::int32_t > domain_coeff(IntegralType type, int idx, int c) const
Coefficient function mesh integration entity indices.
Definition Form.h:556
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
const std::vector< std::shared_ptr< const FunctionSpace< geometry_type > > > & function_spaces() const
Function spaces for all arguments.
Definition Form.h:370
std::span< const std::int32_t > domain_arg(IntegralType type, int rank, int idx, int kernel_idx) const
Argument function mesh integration entity indices.
Definition Form.h:526
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
T scalar_type
Scalar type.
Definition Form.h:120
This class represents a finite element function space defined by a mesh, a finite element,...
Definition FunctionSpace.h:34
Definition Function.h:47
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