DOLFINx 0.11.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
utils.h
Go to the documentation of this file.
1// Copyright (C) 2013-2025 Johan Hake, Jan Blechta, Garth N. Wells and Paul T.
2// Kühner
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 "Constant.h"
11#include "CoordinateElement.h"
12#include "DofMap.h"
13#include "ElementDofLayout.h"
14#include "Expression.h"
15#include "Form.h"
16#include "Function.h"
17#include "FunctionSpace.h"
18#include "kernel.h"
19#include "sparsitybuild.h"
20#include <algorithm>
21#include <array>
22#include <concepts>
23#include <cstddef>
24#include <dolfinx/common/defines.h>
25#include <dolfinx/common/types.h>
26#include <dolfinx/la/SparsityPattern.h>
27#include <dolfinx/mesh/EntityMap.h>
28#include <dolfinx/mesh/Topology.h>
29#include <dolfinx/mesh/cell_types.h>
30#include <dolfinx/mesh/utils.h>
31#include <format>
32#include <functional>
33#include <memory>
34#include <optional>
35#include <ranges>
36#include <set>
37#include <span>
38#include <stdexcept>
39#include <string>
40#include <tuple>
41#include <type_traits>
42#include <ufcx.h>
43#include <utility>
44#include <vector>
45
48namespace basix
49{
50template <std::floating_point T>
51class FiniteElement;
52}
53
54namespace dolfinx::common
55{
56class IndexMap;
57}
58
59namespace dolfinx::fem
60{
61template <dolfinx::scalar T, std::floating_point U>
62class Expression;
63
64namespace impl
65{
72template <int num_cells>
73std::array<std::int32_t, 2 * num_cells>
74get_cell_facet_pairs(std::int32_t f, std::span<const std::int32_t> cells,
76{
77 // Loop over cells sharing facet
78 assert(cells.size() == num_cells);
79 std::array<std::int32_t, 2 * num_cells> cell_local_facet_pairs;
80 for (int c = 0; c < num_cells; ++c)
81 {
82 // Get local index of facet with respect to the cell
83 std::int32_t cell = cells[c];
84 auto cell_facets = c_to_f.links(cell);
85 auto facet_it = std::find(cell_facets.begin(), cell_facets.end(), f);
86 assert(facet_it != cell_facets.end());
87 int local_f = std::distance(cell_facets.begin(), facet_it);
88 cell_local_facet_pairs[2 * c] = cell;
89 cell_local_facet_pairs[2 * c + 1] = local_f;
90 }
91
92 return cell_local_facet_pairs;
93}
94
102template <int num_cells>
103std::array<std::int32_t, 2 * num_cells>
104get_cell_entity_pairs(std::int32_t e, std::span<const std::int32_t> cells,
106{
107 static_assert(num_cells == 1); // Patch assembly not supported.
108
109 assert(cells.size() > 0);
110
111 // Use first cell for assembly over by default
112 std::int32_t cell = cells[0];
113
114 // Find local index of entity within cell
115 auto cell_entities = c_to_e.links(cell);
116 auto it = std::ranges::find(cell_entities, e);
117 assert(it != cell_entities.end());
118 std::int32_t local_index = std::distance(cell_entities.begin(), it);
119
120 return {cell, local_index};
121}
122
123} // namespace impl
124
158std::vector<std::int32_t>
160 const mesh::Topology& topology,
161 std::span<const std::int32_t> entities);
162
170template <dolfinx::scalar T, std::floating_point U>
171std::vector<std::vector<std::array<std::shared_ptr<const FunctionSpace<U>>, 2>>>
172extract_function_spaces(const std::vector<std::vector<const Form<T, U>*>>& a)
173{
174 std::vector<
175 std::vector<std::array<std::shared_ptr<const FunctionSpace<U>>, 2>>>
176 spaces(
177 a.size(),
178 std::vector<std::array<std::shared_ptr<const FunctionSpace<U>>, 2>>(
179 a.front().size()));
180 for (std::size_t i = 0; i < a.size(); ++i)
181 {
182 for (std::size_t j = 0; j < a[i].size(); ++j)
183 {
184 if (const Form<T, U>* form = a[i][j]; form)
185 spaces[i][j] = {form->function_spaces()[0], form->function_spaces()[1]};
186 }
187 }
188 return spaces;
189}
190
196template <dolfinx::scalar T, std::floating_point U>
198{
199 std::shared_ptr mesh = a.mesh();
200 assert(mesh);
201
202 // Get index maps and block sizes from the DOF maps. Note that in
203 // mixed-topology meshes, despite there being multiple DOF maps, the
204 // index maps and block sizes are the same.
205 std::array<std::reference_wrapper<const DofMap>, 2> dofmaps{
206 *a.function_spaces().at(0)->dofmaps(0),
207 *a.function_spaces().at(1)->dofmaps(0)};
208
209 const std::array index_maps{dofmaps[0].get().index_map,
210 dofmaps[1].get().index_map};
211 const std::array bs
212 = {dofmaps[0].get().index_map_bs(), dofmaps[1].get().index_map_bs()};
213
214 la::SparsityPattern pattern(mesh->comm(), index_maps, bs);
215 build_sparsity_pattern(pattern, a);
216 return pattern;
217}
218
224template <dolfinx::scalar T, std::floating_point U>
226{
227 if (a.rank() != 2)
228 {
229 throw std::runtime_error(
230 "Cannot create sparsity pattern. Form is not a bilinear.");
231 }
232
233 std::shared_ptr mesh = a.mesh();
234 assert(mesh);
235 std::shared_ptr mesh0 = a.function_spaces().at(0)->mesh();
236 assert(mesh0);
237 std::shared_ptr mesh1 = a.function_spaces().at(1)->mesh();
238 assert(mesh1);
239
240 const std::set<IntegralType> types = a.integral_types();
241 if (types.find(IntegralType::interior_facet) != types.end()
242 or types.find(IntegralType::exterior_facet) != types.end())
243 {
244 // FIXME: cleanup these calls? Some of the happen internally again.
245 int tdim = mesh->topology()->dim();
246 mesh->topology_mutable()->create_entities(tdim - 1);
247 mesh->topology_mutable()->create_connectivity(tdim - 1, tdim);
248 }
249
250 common::Timer t0("Build sparsity");
251
252 auto extract_cells = [](std::span<const std::int32_t> facets)
253 {
254 assert(facets.size() % 2 == 0);
255 std::vector<std::int32_t> cells;
256 cells.reserve(facets.size() / 2);
257 for (std::size_t i = 0; i < facets.size(); i += 2)
258 cells.push_back(facets[i]);
259 return cells;
260 };
261
262 const int num_cell_types = mesh->topology()->cell_types().size();
263 for (int cell_type_idx = 0; cell_type_idx < num_cell_types; ++cell_type_idx)
264 {
265 std::array<std::reference_wrapper<const DofMap>, 2> dofmaps{
266 *a.function_spaces().at(0)->dofmaps(cell_type_idx),
267 *a.function_spaces().at(1)->dofmaps(cell_type_idx)};
268
269 // Create and build sparsity pattern
270 for (auto type : types)
271 {
272 switch (type)
273 {
275 for (int i = 0; i < a.num_integrals(type, cell_type_idx); ++i)
276 {
278 pattern,
279 std::pair{a.domain_arg(type, 0, i, cell_type_idx),
280 a.domain_arg(type, 1, i, cell_type_idx)},
281 {{dofmaps[0], dofmaps[1]}});
282 }
283 break;
285 for (int i = 0; i < a.num_integrals(type, cell_type_idx); ++i)
286 {
288 pattern,
289 {extract_cells(a.domain_arg(type, 0, i, 0)),
290 extract_cells(a.domain_arg(type, 1, i, 0))},
291 {{dofmaps[0], dofmaps[1]}});
292 }
293 break;
297 for (int i = 0; i < a.num_integrals(type, cell_type_idx); ++i)
298 {
300 pattern,
301 std::pair{extract_cells(a.domain_arg(type, 0, i, 0)),
302 extract_cells(a.domain_arg(type, 1, i, 0))},
303 {{dofmaps[0], dofmaps[1]}});
304 }
305 break;
306 default:
307 throw std::runtime_error("Unsupported integral type");
308 }
309 }
310 }
311
312 t0.stop();
313}
314
316template <std::floating_point T>
318 const std::vector<int>& parent_map
319 = {})
320{
321 // Create subdofmaps and compute offset
322 std::vector<int> offsets(1, 0);
323 std::vector<dolfinx::fem::ElementDofLayout> sub_doflayout;
324 int bs = element.block_size();
325 for (int i = 0; i < element.num_sub_elements(); ++i)
326 {
327 // The ith sub-element. For mixed elements this is subelements()[i]. For
328 // blocked elements, the sub-element will always be the same, so we'll use
329 // sub_elements()[0]
330 std::shared_ptr<const fem::FiniteElement<T>> sub_e
331 = element.sub_elements()[bs > 1 ? 0 : i];
332
333 // In a mixed element DOFs are ordered element by element, so the offset to
334 // the next sub-element is sub_e->space_dimension(). Blocked elements use
335 // xxyyzz ordering, so the offset to the next sub-element is 1
336
337 std::vector<int> parent_map_sub(sub_e->space_dimension(), offsets.back());
338 for (std::size_t j = 0; j < parent_map_sub.size(); ++j)
339 parent_map_sub[j] += bs * j;
340 offsets.push_back(offsets.back() + (bs > 1 ? 1 : sub_e->space_dimension()));
341 sub_doflayout.push_back(
342 dolfinx::fem::create_element_dof_layout(*sub_e, parent_map_sub));
343 }
344
345 return ElementDofLayout(bs, element.entity_dofs(),
346 element.entity_closure_dofs(), parent_map,
347 sub_doflayout);
348}
349
358DofMap
359create_dofmap(MPI_Comm comm, const ElementDofLayout& layout,
360 mesh::Topology& topology,
361 const std::function<void(std::span<std::int32_t>, std::uint32_t)>&
362 permute_inv,
363 const std::function<std::vector<int>(
364 const graph::AdjacencyList<std::int32_t>&)>& reorder_fn);
365
376std::vector<DofMap> create_dofmaps(
377 MPI_Comm comm, const std::vector<ElementDofLayout>& layouts,
378 mesh::Topology& topology,
379 const std::function<void(std::span<std::int32_t>, std::uint32_t)>&
380 permute_inv,
381 const std::function<std::vector<int>(
382 const graph::AdjacencyList<std::int32_t>&)>& reorder_fn);
383
387std::vector<std::string> get_coefficient_names(const ufcx_form& ufcx_form);
388
392std::vector<std::string> get_constant_names(const ufcx_form& ufcx_form);
393
412template <dolfinx::scalar T, std::floating_point U = scalar_value_t<T>>
414 const std::vector<std::reference_wrapper<const ufcx_form>>& ufcx_forms,
415 const std::vector<std::shared_ptr<const FunctionSpace<U>>>& spaces,
416 const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients,
417 const std::vector<std::shared_ptr<const Constant<T>>>& constants,
418 const std::map<
420 std::vector<std::pair<std::int32_t, std::span<const std::int32_t>>>>&
421 subdomains,
422 const std::vector<std::reference_wrapper<const mesh::EntityMap>>&
423 entity_maps,
424 std::shared_ptr<const mesh::Mesh<U>> mesh = nullptr)
425{
426 for (const ufcx_form& ufcx_form : ufcx_forms)
427 {
428 if (ufcx_form.rank != (int)spaces.size())
429 throw std::runtime_error("Wrong number of argument spaces for Form.");
430 if (ufcx_form.num_coefficients != (int)coefficients.size())
431 {
432 throw std::runtime_error("Mismatch between number of expected and "
433 "provided Form coefficients.");
434 }
435
436 // Check Constants for rank and size consistency
437 if (ufcx_form.num_constants != (int)constants.size())
438 {
439 throw std::runtime_error(std::format(
440 "Mismatch between number of expected and "
441 "provided Form Constants. Expected {} constants, but got {}.",
442 ufcx_form.num_constants, constants.size()));
443 }
444 for (std::size_t c = 0; c < constants.size(); ++c)
445 {
446 if (ufcx_form.constant_ranks[c] != (int)constants[c]->shape.size())
447 {
448 throw std::runtime_error(std::format(
449 "Mismatch between expected and actual rank of "
450 "Form Constant. Rank of Constant {} should be {}, but got rank {}.",
451 c, ufcx_form.constant_ranks[c], constants[c]->shape.size()));
452 }
453 if (!std::equal(constants[c]->shape.begin(), constants[c]->shape.end(),
454 ufcx_form.constant_shapes[c]))
455 {
456 throw std::runtime_error(
457 std::format("Mismatch between expected and actual shape of Form "
458 "Constant for Constant {}.",
459 c));
460 }
461 }
462 }
463
464 // Check argument function spaces
465 for (std::size_t form_idx = 0; form_idx < ufcx_forms.size(); ++form_idx)
466 {
467 for (std::size_t i = 0; i < spaces.size(); ++i)
468 {
469 assert(spaces[i]->elements(form_idx));
470 if (auto element_hash
471 = ufcx_forms[form_idx].get().finite_element_hashes[i];
472 element_hash != 0
473 and element_hash
474 != spaces[i]->elements(form_idx)->basix_element().hash())
475 {
476 throw std::runtime_error(
477 "Cannot create form. Elements are different to "
478 "those used to compile the form.");
479 }
480 }
481 }
482
483 // Extract mesh from FunctionSpace, and check they are the same
484 if (!mesh and !spaces.empty())
485 mesh = spaces.front()->mesh();
486 if (!mesh)
487 throw std::runtime_error("No mesh could be associated with the Form.");
488
489 auto topology = mesh->topology();
490 assert(topology);
491 const int tdim = topology->dim();
492
493 // NOTE: This assumes all forms in mixed-topology meshes have the same
494 // integral offsets. Since the UFL forms for each type of cell should be
495 // the same, I think this assumption is OK.
496 const int* integral_offsets = ufcx_forms[0].get().form_integral_offsets;
497 std::array<int, 5> num_integrals_type;
498 for (std::size_t i = 0; i < num_integrals_type.size(); ++i)
499 num_integrals_type[i] = integral_offsets[i + 1] - integral_offsets[i];
500
501 // Create vertices, if required
502 if (num_integrals_type[vertex] > 0)
503 {
504 mesh->topology_mutable()->create_connectivity(0, tdim);
505 mesh->topology_mutable()->create_connectivity(tdim, 0);
506 }
507
508 // Create facets, if required
509 // NOTE: exterior_facet and interior_facet is declared in ufcx.h
510 if (num_integrals_type[exterior_facet] > 0
511 or num_integrals_type[interior_facet] > 0)
512 {
513 mesh->topology_mutable()->create_entities(tdim - 1);
514 mesh->topology_mutable()->create_connectivity(tdim - 1, tdim);
515 mesh->topology_mutable()->create_connectivity(tdim, tdim - 1);
516 }
517
518 // Create ridges, if required
519 if (num_integrals_type[ridge] > 0)
520 {
521 mesh->topology_mutable()->create_entities(tdim - 2);
522 mesh->topology_mutable()->create_connectivity(tdim - 2, tdim);
523 mesh->topology_mutable()->create_connectivity(tdim, tdim - 2);
524 }
525
526 // Get list of integral IDs, and load tabulate tensor into memory for
527 // each
528 std::map<std::tuple<IntegralType, int, int>, integral_data<T, U>> integrals;
529
530 auto check_geometry_hash
531 = [&geo = mesh->geometry()](const ufcx_integral& integral,
532 std::size_t cell_idx)
533 {
534 if (integral.coordinate_element_hash != geo.cmaps().at(cell_idx).hash())
535 {
536 throw std::runtime_error(
537 "Generated integral geometry element does not match mesh geometry: "
538 + std::to_string(integral.coordinate_element_hash) + ", "
539 + std::to_string(geo.cmaps().at(cell_idx).hash()));
540 }
541 };
542
543 // Attach cell kernels
544 bool needs_facet_permutations = false;
545 {
546 std::vector<std::int32_t> default_cells;
547 std::span<const int> ids(ufcx_forms[0].get().form_integral_ids
548 + integral_offsets[cell],
549 num_integrals_type[cell]);
550 auto sd = subdomains.find(IntegralType::cell);
551 for (std::size_t form_idx = 0; form_idx < ufcx_forms.size(); ++form_idx)
552 {
553 const ufcx_form& ufcx_form = ufcx_forms[form_idx];
554 for (int i = 0; i < num_integrals_type[cell]; ++i)
555 {
556 const int id = ids[i];
557 ufcx_integral* integral
558 = ufcx_form.form_integrals[integral_offsets[cell] + i];
559 assert(integral);
560 check_geometry_hash(*integral, form_idx);
561
562 // Build list of active coefficients
563 std::vector<int> active_coeffs;
564 for (int j = 0; j < ufcx_form.num_coefficients; ++j)
565 {
566 if (integral->enabled_coefficients[j])
567 active_coeffs.push_back(j);
568 }
569
570 impl::kernel_t<T, U> k = impl::extract_kernel<T>(integral);
571 if (!k)
572 {
573 throw std::runtime_error(
574 "UFCx kernel function is NULL. Check requested types.");
575 }
576
577 // Build list of entities to assemble over
578 if (id == -1)
579 {
580 // Default kernel, operates on all (owned) cells
581 assert(topology->index_maps(tdim).at(form_idx));
582 default_cells.resize(
583 topology->index_maps(tdim).at(form_idx)->size_local(), 0);
584 std::iota(default_cells.begin(), default_cells.end(), 0);
585 integrals.insert({{IntegralType::cell, i, form_idx},
586 {k, default_cells, active_coeffs}});
587 }
588 else if (sd != subdomains.end())
589 {
590 // NOTE: This requires that pairs are sorted
591 auto it = std::ranges::lower_bound(sd->second, id, std::less<>{},
592 [](auto& a) { return a.first; });
593 if (it != sd->second.end() and it->first == id)
594 {
595 integrals.insert({{IntegralType::cell, i, form_idx},
596 {k,
597 std::vector<std::int32_t>(it->second.begin(),
598 it->second.end()),
599 active_coeffs}});
600 }
601 }
602
603 if (integral->needs_facet_permutations)
604 needs_facet_permutations = true;
605 }
606 }
607 }
608
609 // Attach interior facet kernels
610 {
611 std::vector<std::int32_t> default_facets_int;
612 std::span<const int> ids(ufcx_forms[0].get().form_integral_ids
613 + integral_offsets[interior_facet],
614 num_integrals_type[interior_facet]);
615 auto sd = subdomains.find(IntegralType::interior_facet);
616 for (std::size_t form_idx = 0; form_idx < ufcx_forms.size(); ++form_idx)
617 {
618 const ufcx_form& ufcx_form = ufcx_forms[form_idx];
619
620 // Create indicator for interprocess facets
621 std::vector<std::int8_t> interprocess_marker;
622 if (num_integrals_type[interior_facet] > 0)
623 {
624 assert(topology->index_map(tdim - 1));
625 const std::vector<std::int32_t>& interprocess_facets
626 = topology->interprocess_facets();
627 std::int32_t num_facets = topology->index_map(tdim - 1)->size_local()
628 + topology->index_map(tdim - 1)->num_ghosts();
629 interprocess_marker.resize(num_facets, 0);
630 std::ranges::for_each(interprocess_facets,
631 [&interprocess_marker](auto f)
632 { interprocess_marker[f] = 1; });
633 }
634
635 for (int i = 0; i < num_integrals_type[interior_facet]; ++i)
636 {
637 const int id = ids[i];
638 ufcx_integral* integral
639 = ufcx_form.form_integrals[integral_offsets[interior_facet] + i];
640 assert(integral);
641 check_geometry_hash(*integral, form_idx);
642
643 std::vector<int> active_coeffs;
644 for (int j = 0; j < ufcx_form.num_coefficients; ++j)
645 {
646 if (integral->enabled_coefficients[j])
647 active_coeffs.push_back(j);
648 }
649
650 impl::kernel_t<T, U> k = impl::extract_kernel<T>(integral);
651 assert(k);
652
653 // Build list of entities to assembler over
654 auto f_to_c = topology->connectivity(tdim - 1, tdim);
655 assert(f_to_c);
656 auto c_to_f = topology->connectivity(tdim, tdim - 1);
657 assert(c_to_f);
658 if (id == -1)
659 {
660 // Default kernel, operates on all (owned) interior facets
661 assert(topology->index_map(tdim - 1));
662 std::int32_t num_facets = topology->index_map(tdim - 1)->size_local();
663 default_facets_int.reserve(4 * num_facets);
664 for (std::int32_t f = 0; f < num_facets; ++f)
665 {
666 if (f_to_c->num_links(f) == 2)
667 {
668 std::array<std::int32_t, 4> pairs
669 = impl::get_cell_facet_pairs<2>(f, f_to_c->links(f), *c_to_f);
670 default_facets_int.insert(default_facets_int.end(), pairs.begin(),
671 pairs.end());
672 }
673 else if (interprocess_marker[f])
674 {
675 throw std::runtime_error(
676 "Cannot compute interior facet integral over interprocess "
677 "facet. Please use ghost mode shared facet when creating the "
678 "mesh");
679 }
680 }
681 integrals.insert({{IntegralType::interior_facet, i, form_idx},
682 {k, default_facets_int, active_coeffs}});
683 }
684 else if (sd != subdomains.end())
685 {
686 auto it = std::ranges::lower_bound(sd->second, id, std::less{},
687 [](auto& a) { return a.first; });
688 if (it != sd->second.end() and it->first == id)
689 {
690 integrals.insert({{IntegralType::interior_facet, i, form_idx},
691 {k,
692 std::vector<std::int32_t>(it->second.begin(),
693 it->second.end()),
694 active_coeffs}});
695 }
696 }
697
698 if (integral->needs_facet_permutations)
699 needs_facet_permutations = true;
700 }
701 }
702 }
703
704 // Attach exterior entity integrals
705 {
708 {
709 std::size_t dim;
710 switch (itg_type)
711 {
713 {
714 dim = tdim - 1;
715 break;
716 }
718 {
719 dim = tdim - 2;
720 break;
721 }
723 {
724 dim = 0;
725 break;
726 }
727 default:
728 throw std::runtime_error("Unsupported integral type");
729 }
730
731 const std::function<std::vector<std::int32_t>(const mesh::Topology&,
733 get_default_integration_entities
734 = [dim](const mesh::Topology& topology, IntegralType itg_type)
735 {
736 if (itg_type == IntegralType::exterior_facet)
737 {
738 // Integrate over all owned exterior facets
739 return mesh::exterior_facet_indices(topology);
740 }
741 else
742 {
743 // Integrate over all owned entities
744 std::int32_t num_entities = topology.index_map(dim)->size_local();
745 std::vector<std::int32_t> entities(num_entities);
746 std::iota(entities.begin(), entities.end(), 0);
747 return entities;
748 }
749 };
750
751 std::vector<std::int32_t> default_entities_ext;
752
753 std::span<const int> ids(ufcx_forms[0].get().form_integral_ids
754 + integral_offsets[(std::int8_t)itg_type],
755 num_integrals_type[(std::int8_t)itg_type]);
756 auto sd = subdomains.find(itg_type);
757 for (std::size_t form_idx = 0; form_idx < ufcx_forms.size(); ++form_idx)
758 {
759 const ufcx_form& ufcx_form = ufcx_forms[form_idx];
760 for (int i = 0; i < num_integrals_type[(std::int8_t)itg_type]; ++i)
761 {
762 const int id = ids[i];
763 ufcx_integral* integral
764 = ufcx_form.form_integrals[integral_offsets[(std::int8_t)itg_type]
765 + i];
766 assert(integral);
767 check_geometry_hash(*integral, form_idx);
768
769 std::vector<int> active_coeffs;
770 for (int j = 0; j < ufcx_form.num_coefficients; ++j)
771 {
772 if (integral->enabled_coefficients[j])
773 active_coeffs.push_back(j);
774 }
775
776 impl::kernel_t<T, U> k = impl::extract_kernel<T>(integral);
777
778 // Build list of entities to assembler over
779 auto e_to_c = topology->connectivity(dim, tdim);
780 assert(e_to_c);
781 auto c_to_e = topology->connectivity(tdim, dim);
782 assert(c_to_e);
783 if (id == -1)
784 {
785 std::vector default_entities
786 = get_default_integration_entities(*topology, itg_type);
787 // Default kernel
788 default_entities_ext.reserve(2 * default_entities.size());
789 for (std::int32_t e : default_entities)
790 {
791 // There will only be one pair for an exterior facet integral
792 std::array<std::int32_t, 2> pair = impl::get_cell_entity_pairs<1>(
793 e, e_to_c->links(e), *c_to_e);
794 default_entities_ext.insert(default_entities_ext.end(),
795 pair.begin(), pair.end());
796 }
797 integrals.insert({{itg_type, i, form_idx},
798 {k, default_entities_ext, active_coeffs}});
799 }
800 else if (sd != subdomains.end())
801 {
802 // NOTE: This requires that pairs are sorted
803 auto it = std::ranges::lower_bound(sd->second, id, std::less<>{},
804 [](auto& a) { return a.first; });
805 if (it != sd->second.end() and it->first == id)
806 {
807 integrals.insert({{itg_type, i, form_idx},
808 {k,
809 std::vector<std::int32_t>(it->second.begin(),
810 it->second.end()),
811 active_coeffs}});
812 }
813 }
814
815 if (integral->needs_facet_permutations)
816 needs_facet_permutations = true;
817 }
818 }
819 }
820 }
821
822 return Form<T, U>(spaces, std::move(integrals), mesh, coefficients, constants,
823 needs_facet_permutations, entity_maps);
824}
825
840template <dolfinx::scalar T, std::floating_point U = scalar_value_t<T>>
842 const ufcx_form& ufcx_form,
843 const std::vector<std::shared_ptr<const FunctionSpace<U>>>& spaces,
844 const std::map<std::string, std::shared_ptr<const Function<T, U>>>&
845 coefficients,
846 const std::map<std::string, std::shared_ptr<const Constant<T>>>& constants,
847 const std::map<
849 std::vector<std::pair<std::int32_t, std::span<const std::int32_t>>>>&
850 subdomains,
851 const std::vector<std::reference_wrapper<const mesh::EntityMap>>&
852 entity_maps,
853 std::shared_ptr<const mesh::Mesh<U>> mesh = nullptr)
854{
855 // Place coefficients in appropriate order
856 std::vector<std::shared_ptr<const Function<T, U>>> coeff_map;
857 for (const std::string& name : get_coefficient_names(ufcx_form))
858 {
859 if (auto it = coefficients.find(name); it != coefficients.end())
860 coeff_map.push_back(it->second);
861 else
862 {
863 throw std::runtime_error("Form coefficient \"" + name
864 + "\" not provided.");
865 }
866 }
867
868 // Place constants in appropriate order
869 std::vector<std::shared_ptr<const Constant<T>>> const_map;
870 for (const std::string& name : get_constant_names(ufcx_form))
871 {
872 if (auto it = constants.find(name); it != constants.end())
873 const_map.push_back(it->second);
874 else
875 throw std::runtime_error("Form constant \"" + name + "\" not provided.");
876 }
877
878 return create_form_factory({ufcx_form}, spaces, coeff_map, const_map,
879 subdomains, entity_maps, mesh);
880}
881
900template <dolfinx::scalar T, std::floating_point U = scalar_value_t<T>>
902 ufcx_form* (*fptr)(),
903 const std::vector<std::shared_ptr<const FunctionSpace<U>>>& spaces,
904 const std::map<std::string, std::shared_ptr<const Function<T, U>>>&
905 coefficients,
906 const std::map<std::string, std::shared_ptr<const Constant<T>>>& constants,
907 const std::map<
909 std::vector<std::pair<std::int32_t, std::span<const std::int32_t>>>>&
910 subdomains,
911 const std::vector<std::reference_wrapper<const mesh::EntityMap>>&
912 entity_maps,
913 std::shared_ptr<const mesh::Mesh<U>> mesh = nullptr)
914{
915 ufcx_form* form = fptr();
916 Form<T, U> L = create_form<T, U>(*form, spaces, coefficients, constants,
917 subdomains, entity_maps, mesh);
918 std::free(form);
919 return L;
920}
921
923template <std::floating_point T>
925 std::shared_ptr<mesh::Mesh<T>> mesh,
926 std::shared_ptr<const fem::FiniteElement<T>> e,
927 std::function<std::vector<int>(const graph::AdjacencyList<std::int32_t>&)>
928 reorder_fn = nullptr)
929{
930 // TODO: check cell type of e (need to add method to fem::FiniteElement)
931 assert(e);
932 assert(mesh);
933 assert(mesh->topology());
934 if (e->cell_type() != mesh->topology()->cell_type())
935 throw std::runtime_error("Cell type of element and mesh must match.");
936
937 // Create element dof layout
939
940 // Create a dofmap
941 std::function<void(std::span<std::int32_t>, std::uint32_t)> permute_inv
942 = e->needs_dof_permutations() ? e->dof_permutation_fn(true, true)
943 : nullptr;
944 auto dofmap = std::make_shared<const DofMap>(create_dofmap(
945 mesh->comm(), layout, *mesh->topology(), permute_inv, reorder_fn));
946
947 return FunctionSpace(mesh, e, dofmap);
948}
949
951template <dolfinx::scalar T, std::floating_point U = scalar_value_t<T>>
953 const ufcx_expression& e,
954 const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients,
955 const std::vector<std::shared_ptr<const Constant<T>>>& constants,
956 std::shared_ptr<const FunctionSpace<U>> argument_space = nullptr)
957{
958 if (!coefficients.empty())
959 {
960 assert(coefficients.front());
961 assert(coefficients.front()->function_space());
962 std::shared_ptr<const mesh::Mesh<U>> mesh
963 = coefficients.front()->function_space()->mesh();
964 if (mesh->geometry().cmap().hash() != e.coordinate_element_hash)
965 {
966 throw std::runtime_error(
967 "Expression and mesh geometric maps do not match.");
968 }
969 }
970
971 if (e.rank > 0 and !argument_space)
972 {
973 throw std::runtime_error("Expression has Argument but no Argument "
974 "function space was provided.");
975 }
976
977 std::vector<U> X(e.points, e.points + e.num_points * e.entity_dimension);
978 std::array<std::size_t, 2> Xshape
979 = {static_cast<std::size_t>(e.num_points),
980 static_cast<std::size_t>(e.entity_dimension)};
981 std::vector<std::size_t> value_shape(e.value_shape,
982 e.value_shape + e.num_components);
983 std::function<void(T*, const T*, const T*, const scalar_value_t<T>*,
984 const int*, const std::uint8_t*, void*)>
985 tabulate_tensor = nullptr;
986 if constexpr (std::is_same_v<T, float>)
987 tabulate_tensor = e.tabulate_tensor_float32;
988#ifndef DOLFINX_NO_STDC_COMPLEX_KERNELS
989 else if constexpr (std::is_same_v<T, std::complex<float>>)
990 {
991 tabulate_tensor = reinterpret_cast<void (*)(
992 T*, const T*, const T*, const scalar_value_t<T>*, const int*,
993 const unsigned char*, void*)>(e.tabulate_tensor_complex64);
994 }
995#endif // DOLFINX_NO_STDC_COMPLEX_KERNELS
996 else if constexpr (std::is_same_v<T, double>)
997 tabulate_tensor = e.tabulate_tensor_float64;
998#ifndef DOLFINX_NO_STDC_COMPLEX_KERNELS
999 else if constexpr (std::is_same_v<T, std::complex<double>>)
1000 {
1001 tabulate_tensor = reinterpret_cast<void (*)(
1002 T*, const T*, const T*, const scalar_value_t<T>*, const int*,
1003 const unsigned char*, void*)>(e.tabulate_tensor_complex128);
1004 }
1005#endif // DOLFINX_NO_STDC_COMPLEX_KERNELS
1006 else
1007 throw std::runtime_error("Type not supported.");
1008
1009 assert(tabulate_tensor);
1010 return Expression(coefficients, constants, std::span<const U>(X), Xshape,
1011 tabulate_tensor, value_shape, argument_space);
1012}
1013
1016template <dolfinx::scalar T, std::floating_point U = scalar_value_t<T>>
1018 const ufcx_expression& e,
1019 const std::map<std::string, std::shared_ptr<const Function<T, U>>>&
1020 coefficients,
1021 const std::map<std::string, std::shared_ptr<const Constant<T>>>& constants,
1022 std::shared_ptr<const FunctionSpace<U>> argument_space = nullptr)
1023{
1024 // Place coefficients in appropriate order
1025 std::vector<std::shared_ptr<const Function<T, U>>> coeff_map;
1026 std::vector<std::string> coefficient_names;
1027 coefficient_names.reserve(e.num_coefficients);
1028 for (int i = 0; i < e.num_coefficients; ++i)
1029 coefficient_names.push_back(e.coefficient_names[i]);
1030
1031 for (const std::string& name : coefficient_names)
1032 {
1033 if (auto it = coefficients.find(name); it != coefficients.end())
1034 coeff_map.push_back(it->second);
1035 else
1036 {
1037 throw std::runtime_error("Expression coefficient \"" + name
1038 + "\" not provided.");
1039 }
1040 }
1041
1042 // Place constants in appropriate order
1043 std::vector<std::shared_ptr<const Constant<T>>> const_map;
1044 std::vector<std::string> constant_names;
1045 constant_names.reserve(e.num_constants);
1046 for (int i = 0; i < e.num_constants; ++i)
1047 constant_names.push_back(e.constant_names[i]);
1048
1049 for (const std::string& name : constant_names)
1050 {
1051 if (auto it = constants.find(name); it != constants.end())
1052 const_map.push_back(it->second);
1053 else
1054 {
1055 throw std::runtime_error("Expression constant \"" + name
1056 + "\" not provided.");
1057 }
1058 }
1059
1060 return create_expression(e, coeff_map, const_map, argument_space);
1061}
1062
1063} // namespace dolfinx::fem
Degree-of-freedom map representations and tools.
Definition CoordinateElement.h:26
Definition IndexMap.h:97
Timer for measuring and logging elapsed time durations.
Definition Timer.h:40
std::chrono::duration< double, Period > stop()
Stop timer and return elapsed time.
Definition Timer.h:92
Constant (in space) value which can be attached to a Form.
Definition Constant.h:22
Degree-of-freedom map.
Definition DofMap.h:73
Definition ElementDofLayout.h:31
An Expression represents a mathematical expression evaluated at a pre-defined points on a reference c...
Definition Expression.h:41
Model of a finite element.
Definition FiniteElement.h:57
const std::vector< std::shared_ptr< const FiniteElement< geometry_type > > > & sub_elements() const noexcept
Get subelements (if any).
Definition FiniteElement.cpp:396
const std::vector< std::vector< std::vector< int > > > & entity_dofs() const noexcept
Local DOFs associated with each sub-entity of the cell.
Definition FiniteElement.cpp:340
int num_sub_elements() const noexcept
Number of sub elements (for a mixed or blocked element).
Definition FiniteElement.cpp:383
const std::vector< std::vector< std::vector< int > > > & entity_closure_dofs() const noexcept
Local DOFs associated with the closure of each sub-entity of the cell.
Definition FiniteElement.cpp:347
int block_size() const noexcept
Block size of the finite element function space.
Definition FiniteElement.cpp:359
A representation of finite element variational forms.
Definition Form.h:117
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
int rank() const
Rank of the form.
Definition Form.h:358
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
This class represents a finite element function space defined by a mesh, a finite element,...
Definition FunctionSpace.h:34
Definition Function.h:47
This class provides a static adjacency list data structure.
Definition AdjacencyList.h:38
std::span< LinkData > links(std::size_t node)
Get the links (edges) for given node.
Definition AdjacencyList.h:156
Definition SparsityPattern.h:26
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 common::IndexMap > index_map(int dim) const
Get the IndexMap that described the parallel distribution of the mesh entities.
Definition Topology.cpp:843
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:865
std::array< std::int32_t, 2 *num_cells > get_cell_facet_pairs(std::int32_t f, std::span< const std::int32_t > cells, const graph::AdjacencyList< std::int32_t > &c_to_f)
Definition utils.h:74
std::array< std::int32_t, 2 *num_cells > get_cell_entity_pairs(std::int32_t e, std::span< const std::int32_t > cells, const graph::AdjacencyList< std::int32_t > &c_to_e)
Definition utils.h:104
Functions supporting mesh operations.
Miscellaneous classes, functions and types.
Definition dolfinx_common.h:8
void interior_facets(la::SparsityPattern &pattern, std::array< std::span< const std::int32_t >, 2 > cells, std::array< std::reference_wrapper< const DofMap >, 2 > dofmaps)
Iterate over interior facets and insert entries into sparsity pattern.
Definition sparsitybuild.cpp:16
void cells(la::SparsityPattern &pattern, const std::pair< R0, R1 > &cells, std::array< std::reference_wrapper< const DofMap >, 2 > dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition sparsitybuild.h:37
Finite element method functionality.
Definition assemble_expression_impl.h:23
std::vector< std::string > get_constant_names(const ufcx_form &ufcx_form)
Get the name of each constant in a UFC form.
Definition utils.cpp:146
FunctionSpace< T > create_functionspace(std::shared_ptr< mesh::Mesh< T > > mesh, std::shared_ptr< const fem::FiniteElement< T > > e, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> reorder_fn=nullptr)
NEW Create a function space from a fem::FiniteElement.
Definition utils.h:924
Form< T, U > create_form_factory(const std::vector< std::reference_wrapper< const ufcx_form > > &ufcx_forms, const std::vector< std::shared_ptr< const FunctionSpace< U > > > &spaces, const std::vector< std::shared_ptr< const Function< T, U > > > &coefficients, const std::vector< std::shared_ptr< const Constant< T > > > &constants, const std::map< IntegralType, std::vector< std::pair< std::int32_t, std::span< const std::int32_t > > > > &subdomains, const std::vector< std::reference_wrapper< const mesh::EntityMap > > &entity_maps, std::shared_ptr< const mesh::Mesh< U > > mesh=nullptr)
Create a Form from UFCx input with coefficients and constants passed in the required order.
Definition utils.h:413
std::vector< std::string > get_coefficient_names(const ufcx_form &ufcx_form)
Definition utils.cpp:139
std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace< U > >, 2 > > > extract_function_spaces(const std::vector< std::vector< const Form< T, U > * > > &a)
Extract test (0) and trial (1) function spaces pairs for each bilinear form for a rectangular array o...
Definition utils.h:172
Expression< T, U > create_expression(const ufcx_expression &e, const std::vector< std::shared_ptr< const Function< T, U > > > &coefficients, const std::vector< std::shared_ptr< const Constant< T > > > &constants, std::shared_ptr< const FunctionSpace< U > > argument_space=nullptr)
Create Expression from UFC.
Definition utils.h:952
std::vector< std::int32_t > compute_integration_domains(IntegralType integral_type, const mesh::Topology &topology, std::span< const std::int32_t > entities)
Given an integral type and a set of entities, computes and return data for the entities that should b...
Definition utils.cpp:154
Form< T, U > create_form(const ufcx_form &ufcx_form, const std::vector< std::shared_ptr< const FunctionSpace< U > > > &spaces, const std::map< std::string, std::shared_ptr< const Function< T, U > > > &coefficients, const std::map< std::string, std::shared_ptr< const Constant< T > > > &constants, const std::map< IntegralType, std::vector< std::pair< std::int32_t, std::span< const std::int32_t > > > > &subdomains, const std::vector< std::reference_wrapper< const mesh::EntityMap > > &entity_maps, std::shared_ptr< const mesh::Mesh< U > > mesh=nullptr)
Create a Form from UFC input with coefficients and constants resolved by name.
Definition utils.h:841
ElementDofLayout create_element_dof_layout(const fem::FiniteElement< T > &element, const std::vector< int > &parent_map={})
Create an ElementDofLayout from a FiniteElement.
Definition utils.h:317
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
la::SparsityPattern create_sparsity_pattern(const Form< T, U > &a)
Create a sparsity pattern for a given form.
Definition utils.h:197
DofMap create_dofmap(MPI_Comm comm, const ElementDofLayout &layout, mesh::Topology &topology, const std::function< void(std::span< std::int32_t >, std::uint32_t)> &permute_inv, const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &reorder_fn)
Create a dof map on mesh.
Definition utils.cpp:32
std::vector< DofMap > create_dofmaps(MPI_Comm comm, const std::vector< ElementDofLayout > &layouts, mesh::Topology &topology, const std::function< void(std::span< std::int32_t >, std::uint32_t)> &permute_inv, const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &reorder_fn)
Create a set of dofmaps on a given topology.
Definition utils.cpp:78
FunctionSpace(U mesh, V element, W dofmap) -> FunctionSpace< typename std::remove_cvref< typename U::element_type >::type::geometry_type::value_type >
Type deduction.
void build_sparsity_pattern(la::SparsityPattern &pattern, const Form< T, U > &a)
Build a sparsity pattern for a given form.
Definition utils.h:225
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32
std::vector< std::int32_t > exterior_facet_indices(const Topology &topology, int facet_type_idx)
Compute the indices of all exterior facets that are owned by the caller.
Definition utils.cpp:59
Represents integral data, containing the kernel, and a list of entities to integrate over and the ind...
Definition Form.h:52