DOLFINx 0.10.0.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};
45
49template <dolfinx::scalar T, std::floating_point U = scalar_value_t<T>>
51{
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>,
65 std::vector<int>>
67 : kernel(std::forward<K>(kernel)), entities(std::forward<V>(entities)),
68 coeffs(std::forward<W>(coeffs))
69 {
70 }
71
73 std::function<void(T*, const T*, const T*, const U*, const int*,
74 const uint8_t*, void*)>
76
79 std::vector<std::int32_t> entities;
80
83 std::vector<int> coeffs;
84};
85
114template <dolfinx::scalar T, std::floating_point U = dolfinx::scalar_value_t<T>>
115class Form
116{
117public:
119 using scalar_type = T;
120
122 using geometry_type = U;
123
151 template <typename X>
152 requires std::is_convertible_v<
153 std::remove_cvref_t<X>,
154 std::map<std::tuple<IntegralType, int, int>,
157 const std::vector<std::shared_ptr<const FunctionSpace<geometry_type>>>& V,
158 X&& integrals, std::shared_ptr<const mesh::Mesh<geometry_type>> mesh,
159 const std::vector<
160 std::shared_ptr<const Function<scalar_type, geometry_type>>>&
162 const std::vector<std::shared_ptr<const Constant<scalar_type>>>&
163 constants,
165 const std::vector<std::reference_wrapper<const mesh::EntityMap>>&
166 entity_maps)
167 : _function_spaces(V), _integrals(std::forward<X>(integrals)),
168 _mesh(mesh), _coefficients(coefficients), _constants(constants),
169 _needs_facet_permutations(needs_facet_permutations)
170 {
171 if (!_mesh)
172 throw std::runtime_error("Form Mesh is null.");
173
174 // A helper function to find the correct entity map for a given mesh
175 auto get_entity_map
176 = [mesh, &entity_maps](auto& mesh0) -> const mesh::EntityMap&
177 {
178 auto it = std::ranges::find_if(
179 entity_maps,
180 [mesh, mesh0](const mesh::EntityMap& em)
181 {
182 return ((em.topology() == mesh0->topology()
183 and em.sub_topology() == mesh->topology()))
184 or ((em.sub_topology() == mesh0->topology()
185 and em.topology() == mesh->topology()));
186 });
187
188 if (it == entity_maps.end())
189 {
190 throw std::runtime_error(
191 "Incompatible mesh. argument entity_maps must be provided.");
192 }
193 return *it;
194 };
195
196 // A helper function to compute the (cell, local_facet) pairs in the
197 // argument/coefficient domain from the (cell, local_facet) pairs in
198 // `this->mesh()`.
199 auto compute_facet_domains
200 = [&](const auto& int_ents_mesh, int codim, const auto& c_to_f,
201 const auto& emap, bool inverse)
202 {
203 // TODO: This function would be much neater using
204 // `std::views::stride(2)` from C++ 23
205
206 // Get a list of entities to map to the argument/coefficient
207 // domain
208 std::vector<std::int32_t> entities;
209 entities.reserve(int_ents_mesh.size() / 2);
210 if (codim == 0)
211 {
212 // In the codim 0 case, we need to map from cells in
213 // `this->mesh()` to cells in the argument/coefficient mesh, so
214 // here we extract the cells.
215 for (std::size_t i = 0; i < int_ents_mesh.size(); i += 2)
216 entities.push_back(int_ents_mesh[i]);
217 }
218 else if (codim == 1)
219 {
220 // In the codim 1 case, we need to map facets in `this->mesh()`
221 // to cells in the argument/coefficient mesh, so here we extract
222 // the facet index using the cell-to-facet connectivity.
223 for (std::size_t i = 0; i < int_ents_mesh.size(); i += 2)
224 {
225 entities.push_back(
226 c_to_f->links(int_ents_mesh[i])[int_ents_mesh[i + 1]]);
227 }
228 }
229 else
230 throw std::runtime_error("Codimension > 1 not supported.");
231
232 // Map from entity indices in `this->mesh()` to the corresponding
233 // cell indices in the argument/coefficient mesh
234 std::vector<std::int32_t> cells_mesh0
235 = emap.sub_topology_to_topology(entities, inverse);
236
237 // Create a list of (cell, local_facet_index) pairs in the
238 // argument/coefficient domain. Since `create_submesh`preserves
239 // the local facet index (with respect to the cell), we can use
240 // the local facet indices from the input integration entities
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];
244
245 return e;
246 };
247
248 for (auto& space : _function_spaces)
249 {
250 // Working map: [integral type, domain ID, kernel_idx]->entities
251 std::map<std::tuple<IntegralType, int, int>,
252 std::variant<std::vector<std::int32_t>,
253 std::span<const std::int32_t>>>
254 vdata;
255
256 if (auto mesh0 = space->mesh(); mesh0 == _mesh)
257 {
258 for (auto& [key, integral] : _integrals)
259 vdata.insert({key, std::span(integral.entities)});
260 }
261 else
262 {
263 // Find correct entity map
264 const mesh::EntityMap& emap = get_entity_map(mesh0);
265
266 // Determine direction of the map. We need to map from
267 // `this->mesh()` to `mesh0`, so if `emap->sub_topology()` isn't
268 // the source topology, we need the inverse map
269 bool inverse = emap.sub_topology() == mesh0->topology();
270 for (auto& [key, itg] : _integrals)
271 {
272 auto [type, id, kernel_idx] = key;
273 std::vector<std::int32_t> e;
274 if (type == IntegralType::cell)
275 e = emap.sub_topology_to_topology(itg.entities, inverse);
276 else if (type == IntegralType::exterior_facet
278 {
279 const mesh::Topology topology = *_mesh->topology();
280 int tdim = topology.dim();
281 assert(mesh0);
282 int codim = tdim - mesh0->topology()->dim();
283 assert(codim >= 0);
284 auto c_to_f = topology.connectivity(tdim, tdim - 1);
285 assert(c_to_f);
286
287 e = compute_facet_domains(itg.entities, codim, c_to_f, emap,
288 inverse);
289 }
290 else
291 throw std::runtime_error("Integral type not supported.");
292
293 vdata.insert({key, std::move(e)});
294 }
295 }
296
297 _edata.push_back(vdata);
298 }
299
300 for (auto& [key, integral] : _integrals)
301 {
302 auto [type, id, kernel_idx] = key;
303 for (int c : integral.coeffs)
304 {
305 if (auto mesh0 = coefficients.at(c)->function_space()->mesh();
306 mesh0 == _mesh)
307 {
308 _cdata.insert({{type, id, c}, std::span(integral.entities)});
309 }
310 else
311 {
312 // Find correct entity map and determine direction of the map
313 const mesh::EntityMap& emap = get_entity_map(mesh0);
314 bool inverse = emap.sub_topology() == mesh0->topology();
315
316 std::vector<std::int32_t> e;
317 if (type == IntegralType::cell)
318 e = emap.sub_topology_to_topology(integral.entities, inverse);
319 else if (type == IntegralType::exterior_facet
321 {
322 const mesh::Topology topology = *_mesh->topology();
323 int tdim = topology.dim();
324 assert(mesh0);
325 int codim = tdim - mesh0->topology()->dim();
326 auto c_to_f = topology.connectivity(tdim, tdim - 1);
327 assert(c_to_f);
328
329 e = compute_facet_domains(integral.entities, codim, c_to_f, emap,
330 inverse);
331 }
332 else
333 throw std::runtime_error("Integral type not supported.");
334
335 _cdata.insert({{type, id, c}, std::move(e)});
336 }
337 }
338 }
339 }
340
342 Form(const Form& form) = delete;
343
345 Form(Form&& form) = default;
346
348 virtual ~Form() = default;
349
355 int rank() const { return _function_spaces.size(); }
356
359 std::shared_ptr<const mesh::Mesh<geometry_type>> mesh() const
360 {
361 return _mesh;
362 }
363
366 const std::vector<std::shared_ptr<const FunctionSpace<geometry_type>>>&
368 {
369 return _function_spaces;
370 }
371
378 std::function<void(scalar_type*, const scalar_type*, const scalar_type*,
379 const geometry_type*, const int*, const uint8_t*, void*)>
380 kernel(IntegralType type, int id, int kernel_idx) const
381 {
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;
386 }
387
390 std::set<IntegralType> integral_types() const
391 {
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());
396 }
397
408 std::vector<int> active_coeffs(IntegralType type, int id) const
409 {
410 auto it = std::ranges::find_if(_integrals,
411 [type, id](auto& x)
412 {
413 auto [t, id_, kernel_idx] = x.first;
414 return t == type and id_ == id;
415 });
416 if (it == _integrals.end())
417 throw std::runtime_error("Could not find active coefficient list.");
418 return it->second.coeffs;
419 }
420
430 std::vector<int> integral_ids(IntegralType type) const
431 {
432 std::vector<int> ids;
433 for (auto& [key, integral] : _integrals)
434 {
435 auto [t, id, kernel_idx] = key;
436 if (t == type)
437 ids.push_back(id);
438 }
439
440 // IDs may be repeated in mixed-topology meshes, so remove
441 // duplicates
442 std::sort(ids.begin(), ids.end());
443 auto it = std::unique(ids.begin(), ids.end());
444 ids.erase(it, ids.end());
445 return ids;
446 }
447
472 std::span<const std::int32_t> domain(IntegralType type, int id,
473 int kernel_idx) const
474 {
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;
479 }
480
515 std::span<const std::int32_t> domain_arg(IntegralType type, int rank, int id,
516 int kernel_idx) const
517 {
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.");
521 try
522 {
523 return std::get<std::span<const std::int32_t>>(it->second);
524 }
525 catch (std::bad_variant_access& e)
526 {
527 return std::get<std::vector<std::int32_t>>(it->second);
528 }
529 }
530
545 std::span<const std::int32_t> domain_coeff(IntegralType type, int id,
546 int c) const
547 {
548 auto it = _cdata.find({type, id, c});
549 if (it == _cdata.end())
550 throw std::runtime_error("No domain for requested integral.");
551 try
552 {
553 return std::get<std::span<const std::int32_t>>(it->second);
554 }
555 catch (std::bad_variant_access& e)
556 {
557 return std::get<std::vector<std::int32_t>>(it->second);
558 }
559 }
560
562 const std::vector<
563 std::shared_ptr<const Function<scalar_type, geometry_type>>>&
565 {
566 return _coefficients;
567 }
568
572 bool needs_facet_permutations() const { return _needs_facet_permutations; }
573
578 std::vector<int> coefficient_offsets() const
579 {
580 std::vector<int> n{0};
581 for (auto& c : _coefficients)
582 {
583 if (!c)
584 throw std::runtime_error("Not all form coefficients have been set.");
585 n.push_back(n.back() + c->function_space()->element()->space_dimension());
586 }
587 return n;
588 }
589
591 const std::vector<std::shared_ptr<const Constant<scalar_type>>>&
592 constants() const
593 {
594 return _constants;
595 }
596
597private:
598 // Function spaces (one for each argument)
599 std::vector<std::shared_ptr<const FunctionSpace<geometry_type>>>
600 _function_spaces;
601
602 // Integrals (integral type, id, celltype)
603 std::map<std::tuple<IntegralType, int, int>,
605 _integrals;
606
607 // The mesh
608 std::shared_ptr<const mesh::Mesh<geometry_type>> _mesh;
609
610 // Form coefficients
611 std::vector<std::shared_ptr<const Function<scalar_type, geometry_type>>>
612 _coefficients;
613
614 // Constants associated with the Form
615 std::vector<std::shared_ptr<const Constant<scalar_type>>> _constants;
616
617 // True if permutation data needs to be passed into these integrals
618 bool _needs_facet_permutations;
619
620 // Mapped domain index data for argument functions.
621 //
622 // Consider:
623 //
624 // entities = this->domain(IntegralType, integral(id), kernel_idx];
625 // entities0 = _edata[0][IntegralType, integral(id), coefficient_index];
626 //
627 // Then `entities[i]` is a mesh entity index (e.g., cell index) in
628 // `_mesh`, and `entities0[i]` is the index of the same entity but in
629 // the mesh associated with the argument 0 (test function) space.
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>>>>
633 _edata;
634
635 // Mapped domain index data for coefficient functions.
636 //
637 // Consider:
638 //
639 // entities = this->domain(IntegralType, integral(id), kernel_idx];
640 // entities0 = _cdata[IntegralType, integral(id), coefficient_index];
641 //
642 // Then `entities[i]` is a mesh entity index (e.g., cell index) in
643 // `_mesh`, and `entities0[i]` is the index of the same entity but in
644 // the mesh associated with the coefficient Function.
645 std::map<
646 std::tuple<IntegralType, int, int>,
647 std::variant<std::vector<std::int32_t>, std::span<const std::int32_t>>>
648 _cdata;
649};
650} // 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:122
const std::vector< std::shared_ptr< const Function< scalar_type, geometry_type > > > & coefficients() const
Access coefficients.
Definition Form.h:564
bool needs_facet_permutations() const
Get bool indicating whether permutation data needs to be passed into these integrals.
Definition Form.h:572
Form(Form &&form)=default
Move constructor.
int rank() const
Rank of the form.
Definition Form.h:355
Form(const Form &form)=delete
Copy constructor.
virtual ~Form()=default
Destructor.
std::span< const std::int32_t > domain_arg(IntegralType type, int rank, int id, int kernel_idx) const
Argument function mesh integration entity indices.
Definition Form.h:515
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:156
const std::vector< std::shared_ptr< const Constant< scalar_type > > > & constants() const
Access constants.
Definition Form.h:592
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:380
std::vector< int > coefficient_offsets() const
Offset for each coefficient expansion array on a cell.
Definition Form.h:578
std::vector< int > active_coeffs(IntegralType type, int id) const
Indices of coefficients that are active for a given integral (kernel).
Definition Form.h:408
std::shared_ptr< const mesh::Mesh< geometry_type > > mesh() const
Common mesh for the form (the 'integration domain').
Definition Form.h:359
std::span< const std::int32_t > domain_coeff(IntegralType type, int id, int c) const
Coefficient function mesh integration entity indices.
Definition Form.h:545
std::vector< int > integral_ids(IntegralType type) const
Get the IDs for integrals (kernels) for given integral domain type.
Definition Form.h:430
std::set< IntegralType > integral_types() const
Get types of integrals in the form.
Definition Form.h:390
const std::vector< std::shared_ptr< const FunctionSpace< geometry_type > > > & function_spaces() const
Function spaces for all arguments.
Definition Form.h:367
T scalar_type
Scalar type.
Definition Form.h:119
std::span< const std::int32_t > domain(IntegralType type, int id, int kernel_idx) const
Mesh entity indices to integrate over for a given integral (kernel).
Definition Form.h:472
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: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