DOLFINx 0.10.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 {
277 sparsitybuild::cells(pattern,
278 {a.domain_arg(type, 0, i, cell_type_idx),
279 a.domain_arg(type, 1, i, cell_type_idx)},
280 {{dofmaps[0], dofmaps[1]}});
281 }
282 break;
284 for (int i = 0; i < a.num_integrals(type, cell_type_idx); ++i)
285 {
287 pattern,
288 {extract_cells(a.domain_arg(type, 0, i, 0)),
289 extract_cells(a.domain_arg(type, 1, i, 0))},
290 {{dofmaps[0], dofmaps[1]}});
291 }
292 break;
296 for (int i = 0; i < a.num_integrals(type, cell_type_idx); ++i)
297 {
298 sparsitybuild::cells(pattern,
299 {extract_cells(a.domain_arg(type, 0, i, 0)),
300 extract_cells(a.domain_arg(type, 1, i, 0))},
301 {{dofmaps[0], dofmaps[1]}});
302 }
303 break;
304 default:
305 throw std::runtime_error("Unsupported integral type");
306 }
307 }
308 }
309
310 t0.stop();
311}
312
314template <std::floating_point T>
316 const std::vector<int>& parent_map
317 = {})
318{
319 // Create subdofmaps and compute offset
320 std::vector<int> offsets(1, 0);
321 std::vector<dolfinx::fem::ElementDofLayout> sub_doflayout;
322 int bs = element.block_size();
323 for (int i = 0; i < element.num_sub_elements(); ++i)
324 {
325 // The ith sub-element. For mixed elements this is subelements()[i]. For
326 // blocked elements, the sub-element will always be the same, so we'll use
327 // sub_elements()[0]
328 std::shared_ptr<const fem::FiniteElement<T>> sub_e
329 = element.sub_elements()[bs > 1 ? 0 : i];
330
331 // In a mixed element DOFs are ordered element by element, so the offset to
332 // the next sub-element is sub_e->space_dimension(). Blocked elements use
333 // xxyyzz ordering, so the offset to the next sub-element is 1
334
335 std::vector<int> parent_map_sub(sub_e->space_dimension(), offsets.back());
336 for (std::size_t j = 0; j < parent_map_sub.size(); ++j)
337 parent_map_sub[j] += bs * j;
338 offsets.push_back(offsets.back() + (bs > 1 ? 1 : sub_e->space_dimension()));
339 sub_doflayout.push_back(
340 dolfinx::fem::create_element_dof_layout(*sub_e, parent_map_sub));
341 }
342
343 return ElementDofLayout(bs, element.entity_dofs(),
344 element.entity_closure_dofs(), parent_map,
345 sub_doflayout);
346}
347
356DofMap
357create_dofmap(MPI_Comm comm, const ElementDofLayout& layout,
358 mesh::Topology& topology,
359 const std::function<void(std::span<std::int32_t>, std::uint32_t)>&
360 permute_inv,
361 const std::function<std::vector<int>(
362 const graph::AdjacencyList<std::int32_t>&)>& reorder_fn);
363
374std::vector<DofMap> create_dofmaps(
375 MPI_Comm comm, const std::vector<ElementDofLayout>& layouts,
376 mesh::Topology& topology,
377 const std::function<void(std::span<std::int32_t>, std::uint32_t)>&
378 permute_inv,
379 const std::function<std::vector<int>(
380 const graph::AdjacencyList<std::int32_t>&)>& reorder_fn);
381
385std::vector<std::string> get_coefficient_names(const ufcx_form& ufcx_form);
386
390std::vector<std::string> get_constant_names(const ufcx_form& ufcx_form);
391
410template <dolfinx::scalar T, std::floating_point U = scalar_value_t<T>>
412 const std::vector<std::reference_wrapper<const ufcx_form>>& ufcx_forms,
413 const std::vector<std::shared_ptr<const FunctionSpace<U>>>& spaces,
414 const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients,
415 const std::vector<std::shared_ptr<const Constant<T>>>& constants,
416 const std::map<
418 std::vector<std::pair<std::int32_t, std::span<const std::int32_t>>>>&
419 subdomains,
420 const std::vector<std::reference_wrapper<const mesh::EntityMap>>&
421 entity_maps,
422 std::shared_ptr<const mesh::Mesh<U>> mesh = nullptr)
423{
424 for (const ufcx_form& ufcx_form : ufcx_forms)
425 {
426 if (ufcx_form.rank != (int)spaces.size())
427 throw std::runtime_error("Wrong number of argument spaces for Form.");
428 if (ufcx_form.num_coefficients != (int)coefficients.size())
429 {
430 throw std::runtime_error("Mismatch between number of expected and "
431 "provided Form coefficients.");
432 }
433
434 // Check Constants for rank and size consistency
435 if (ufcx_form.num_constants != (int)constants.size())
436 {
437 throw std::runtime_error(std::format(
438 "Mismatch between number of expected and "
439 "provided Form Constants. Expected {} constants, but got {}.",
440 ufcx_form.num_constants, constants.size()));
441 }
442 for (std::size_t c = 0; c < constants.size(); ++c)
443 {
444 if (ufcx_form.constant_ranks[c] != (int)constants[c]->shape.size())
445 {
446 throw std::runtime_error(std::format(
447 "Mismatch between expected and actual rank of "
448 "Form Constant. Rank of Constant {} should be {}, but got rank {}.",
449 c, ufcx_form.constant_ranks[c], constants[c]->shape.size()));
450 }
451 if (!std::equal(constants[c]->shape.begin(), constants[c]->shape.end(),
452 ufcx_form.constant_shapes[c]))
453 {
454 throw std::runtime_error(
455 std::format("Mismatch between expected and actual shape of Form "
456 "Constant for Constant {}.",
457 c));
458 }
459 }
460 }
461
462 // Check argument function spaces
463 for (std::size_t form_idx = 0; form_idx < ufcx_forms.size(); ++form_idx)
464 {
465 for (std::size_t i = 0; i < spaces.size(); ++i)
466 {
467 assert(spaces[i]->elements(form_idx));
468 if (auto element_hash
469 = ufcx_forms[form_idx].get().finite_element_hashes[i];
470 element_hash != 0
471 and element_hash
472 != spaces[i]->elements(form_idx)->basix_element().hash())
473 {
474 throw std::runtime_error(
475 "Cannot create form. Elements are different to "
476 "those used to compile the form.");
477 }
478 }
479 }
480
481 // Extract mesh from FunctionSpace, and check they are the same
482 if (!mesh and !spaces.empty())
483 mesh = spaces.front()->mesh();
484 if (!mesh)
485 throw std::runtime_error("No mesh could be associated with the Form.");
486
487 auto topology = mesh->topology();
488 assert(topology);
489 const int tdim = topology->dim();
490
491 // NOTE: This assumes all forms in mixed-topology meshes have the same
492 // integral offsets. Since the UFL forms for each type of cell should be
493 // the same, I think this assumption is OK.
494 const int* integral_offsets = ufcx_forms[0].get().form_integral_offsets;
495 std::array<int, 5> num_integrals_type;
496 for (std::size_t i = 0; i < num_integrals_type.size(); ++i)
497 num_integrals_type[i] = integral_offsets[i + 1] - integral_offsets[i];
498
499 // Create vertices, if required
500 if (num_integrals_type[vertex] > 0)
501 {
502 mesh->topology_mutable()->create_connectivity(0, tdim);
503 mesh->topology_mutable()->create_connectivity(tdim, 0);
504 }
505
506 // Create facets, if required
507 // NOTE: exterior_facet and interior_facet is declared in ufcx.h
508 if (num_integrals_type[exterior_facet] > 0
509 or num_integrals_type[interior_facet] > 0)
510 {
511 mesh->topology_mutable()->create_entities(tdim - 1);
512 mesh->topology_mutable()->create_connectivity(tdim - 1, tdim);
513 mesh->topology_mutable()->create_connectivity(tdim, tdim - 1);
514 }
515
516 // Create ridges, if required
517 if (num_integrals_type[ridge] > 0)
518 {
519 mesh->topology_mutable()->create_entities(tdim - 2);
520 mesh->topology_mutable()->create_connectivity(tdim - 2, tdim);
521 mesh->topology_mutable()->create_connectivity(tdim, tdim - 2);
522 }
523
524 // Get list of integral IDs, and load tabulate tensor into memory for
525 // each
526 std::map<std::tuple<IntegralType, int, int>, integral_data<T, U>> integrals;
527
528 auto check_geometry_hash
529 = [&geo = mesh->geometry()](const ufcx_integral& integral,
530 std::size_t cell_idx)
531 {
532 if (integral.coordinate_element_hash != geo.cmaps().at(cell_idx).hash())
533 {
534 throw std::runtime_error(
535 "Generated integral geometry element does not match mesh geometry: "
536 + std::to_string(integral.coordinate_element_hash) + ", "
537 + std::to_string(geo.cmaps().at(cell_idx).hash()));
538 }
539 };
540
541 // Attach cell kernels
542 bool needs_facet_permutations = false;
543 {
544 std::vector<std::int32_t> default_cells;
545 std::span<const int> ids(ufcx_forms[0].get().form_integral_ids
546 + integral_offsets[cell],
547 num_integrals_type[cell]);
548 auto sd = subdomains.find(IntegralType::cell);
549 for (std::size_t form_idx = 0; form_idx < ufcx_forms.size(); ++form_idx)
550 {
551 const ufcx_form& ufcx_form = ufcx_forms[form_idx];
552 for (int i = 0; i < num_integrals_type[cell]; ++i)
553 {
554 const int id = ids[i];
555 ufcx_integral* integral
556 = ufcx_form.form_integrals[integral_offsets[cell] + i];
557 assert(integral);
558 check_geometry_hash(*integral, form_idx);
559
560 // Build list of active coefficients
561 std::vector<int> active_coeffs;
562 for (int j = 0; j < ufcx_form.num_coefficients; ++j)
563 {
564 if (integral->enabled_coefficients[j])
565 active_coeffs.push_back(j);
566 }
567
568 impl::kernel_t<T, U> k = impl::extract_kernel<T>(integral);
569 if (!k)
570 {
571 throw std::runtime_error(
572 "UFCx kernel function is NULL. Check requested types.");
573 }
574
575 // Build list of entities to assemble over
576 if (id == -1)
577 {
578 // Default kernel, operates on all (owned) cells
579 assert(topology->index_maps(tdim).at(form_idx));
580 default_cells.resize(
581 topology->index_maps(tdim).at(form_idx)->size_local(), 0);
582 std::iota(default_cells.begin(), default_cells.end(), 0);
583 integrals.insert({{IntegralType::cell, i, form_idx},
584 {k, default_cells, active_coeffs}});
585 }
586 else if (sd != subdomains.end())
587 {
588 // NOTE: This requires that pairs are sorted
589 auto it = std::ranges::lower_bound(sd->second, id, std::less<>{},
590 [](auto& a) { return a.first; });
591 if (it != sd->second.end() and it->first == id)
592 {
593 integrals.insert({{IntegralType::cell, i, form_idx},
594 {k,
595 std::vector<std::int32_t>(it->second.begin(),
596 it->second.end()),
597 active_coeffs}});
598 }
599 }
600
601 if (integral->needs_facet_permutations)
602 needs_facet_permutations = true;
603 }
604 }
605 }
606
607 // Attach interior facet kernels
608 {
609 std::vector<std::int32_t> default_facets_int;
610 std::span<const int> ids(ufcx_forms[0].get().form_integral_ids
611 + integral_offsets[interior_facet],
612 num_integrals_type[interior_facet]);
613 auto sd = subdomains.find(IntegralType::interior_facet);
614 for (std::size_t form_idx = 0; form_idx < ufcx_forms.size(); ++form_idx)
615 {
616 const ufcx_form& ufcx_form = ufcx_forms[form_idx];
617
618 // Create indicator for interprocess facets
619 std::vector<std::int8_t> interprocess_marker;
620 if (num_integrals_type[interior_facet] > 0)
621 {
622 assert(topology->index_map(tdim - 1));
623 const std::vector<std::int32_t>& interprocess_facets
624 = topology->interprocess_facets();
625 std::int32_t num_facets = topology->index_map(tdim - 1)->size_local()
626 + topology->index_map(tdim - 1)->num_ghosts();
627 interprocess_marker.resize(num_facets, 0);
628 std::ranges::for_each(interprocess_facets,
629 [&interprocess_marker](auto f)
630 { interprocess_marker[f] = 1; });
631 }
632
633 for (int i = 0; i < num_integrals_type[interior_facet]; ++i)
634 {
635 const int id = ids[i];
636 ufcx_integral* integral
637 = ufcx_form.form_integrals[integral_offsets[interior_facet] + i];
638 assert(integral);
639 check_geometry_hash(*integral, form_idx);
640
641 std::vector<int> active_coeffs;
642 for (int j = 0; j < ufcx_form.num_coefficients; ++j)
643 {
644 if (integral->enabled_coefficients[j])
645 active_coeffs.push_back(j);
646 }
647
648 impl::kernel_t<T, U> k = impl::extract_kernel<T>(integral);
649 assert(k);
650
651 // Build list of entities to assembler over
652 auto f_to_c = topology->connectivity(tdim - 1, tdim);
653 assert(f_to_c);
654 auto c_to_f = topology->connectivity(tdim, tdim - 1);
655 assert(c_to_f);
656 if (id == -1)
657 {
658 // Default kernel, operates on all (owned) interior facets
659 assert(topology->index_map(tdim - 1));
660 std::int32_t num_facets = topology->index_map(tdim - 1)->size_local();
661 default_facets_int.reserve(4 * num_facets);
662 for (std::int32_t f = 0; f < num_facets; ++f)
663 {
664 if (f_to_c->num_links(f) == 2)
665 {
666 std::array<std::int32_t, 4> pairs
667 = impl::get_cell_facet_pairs<2>(f, f_to_c->links(f), *c_to_f);
668 default_facets_int.insert(default_facets_int.end(), pairs.begin(),
669 pairs.end());
670 }
671 else if (interprocess_marker[f])
672 {
673 throw std::runtime_error(
674 "Cannot compute interior facet integral over interprocess "
675 "facet. Please use ghost mode shared facet when creating the "
676 "mesh");
677 }
678 }
679 integrals.insert({{IntegralType::interior_facet, i, form_idx},
680 {k, default_facets_int, active_coeffs}});
681 }
682 else if (sd != subdomains.end())
683 {
684 auto it = std::ranges::lower_bound(sd->second, id, std::less{},
685 [](auto& a) { return a.first; });
686 if (it != sd->second.end() and it->first == id)
687 {
688 integrals.insert({{IntegralType::interior_facet, i, form_idx},
689 {k,
690 std::vector<std::int32_t>(it->second.begin(),
691 it->second.end()),
692 active_coeffs}});
693 }
694 }
695
696 if (integral->needs_facet_permutations)
697 needs_facet_permutations = true;
698 }
699 }
700 }
701
702 // Attach exterior entity integrals
703 {
706 {
707 std::size_t dim;
708 switch (itg_type)
709 {
711 {
712 dim = tdim - 1;
713 break;
714 }
716 {
717 dim = tdim - 2;
718 break;
719 }
721 {
722 dim = 0;
723 break;
724 }
725 default:
726 throw std::runtime_error("Unsupported integral type");
727 }
728
729 const std::function<std::vector<std::int32_t>(const mesh::Topology&,
731 get_default_integration_entities
732 = [dim](const mesh::Topology& topology, IntegralType itg_type)
733 {
734 if (itg_type == IntegralType::exterior_facet)
735 {
736 // Integrate over all owned exterior facets
737 return mesh::exterior_facet_indices(topology);
738 }
739 else
740 {
741 // Integrate over all owned entities
742 std::int32_t num_entities = topology.index_map(dim)->size_local();
743 std::vector<std::int32_t> entities(num_entities);
744 std::iota(entities.begin(), entities.end(), 0);
745 return entities;
746 }
747 };
748
749 std::vector<std::int32_t> default_entities_ext;
750
751 std::span<const int> ids(ufcx_forms[0].get().form_integral_ids
752 + integral_offsets[(std::int8_t)itg_type],
753 num_integrals_type[(std::int8_t)itg_type]);
754 auto sd = subdomains.find(itg_type);
755 for (std::size_t form_idx = 0; form_idx < ufcx_forms.size(); ++form_idx)
756 {
757 const ufcx_form& ufcx_form = ufcx_forms[form_idx];
758 for (int i = 0; i < num_integrals_type[(std::int8_t)itg_type]; ++i)
759 {
760 const int id = ids[i];
761 ufcx_integral* integral
762 = ufcx_form.form_integrals[integral_offsets[(std::int8_t)itg_type]
763 + i];
764 assert(integral);
765 check_geometry_hash(*integral, form_idx);
766
767 std::vector<int> active_coeffs;
768 for (int j = 0; j < ufcx_form.num_coefficients; ++j)
769 {
770 if (integral->enabled_coefficients[j])
771 active_coeffs.push_back(j);
772 }
773
774 impl::kernel_t<T, U> k = impl::extract_kernel<T>(integral);
775
776 // Build list of entities to assembler over
777 auto e_to_c = topology->connectivity(dim, tdim);
778 assert(e_to_c);
779 auto c_to_e = topology->connectivity(tdim, dim);
780 assert(c_to_e);
781 if (id == -1)
782 {
783 std::vector default_entities
784 = get_default_integration_entities(*topology, itg_type);
785 // Default kernel
786 default_entities_ext.reserve(2 * default_entities.size());
787 for (std::int32_t e : default_entities)
788 {
789 // There will only be one pair for an exterior facet integral
790 std::array<std::int32_t, 2> pair = impl::get_cell_entity_pairs<1>(
791 e, e_to_c->links(e), *c_to_e);
792 default_entities_ext.insert(default_entities_ext.end(),
793 pair.begin(), pair.end());
794 }
795 integrals.insert({{itg_type, i, form_idx},
796 {k, default_entities_ext, active_coeffs}});
797 }
798 else if (sd != subdomains.end())
799 {
800 // NOTE: This requires that pairs are sorted
801 auto it = std::ranges::lower_bound(sd->second, id, std::less<>{},
802 [](auto& a) { return a.first; });
803 if (it != sd->second.end() and it->first == id)
804 {
805 integrals.insert({{itg_type, i, form_idx},
806 {k,
807 std::vector<std::int32_t>(it->second.begin(),
808 it->second.end()),
809 active_coeffs}});
810 }
811 }
812
813 if (integral->needs_facet_permutations)
814 needs_facet_permutations = true;
815 }
816 }
817 }
818 }
819
820 return Form<T, U>(spaces, std::move(integrals), mesh, coefficients, constants,
821 needs_facet_permutations, entity_maps);
822}
823
838template <dolfinx::scalar T, std::floating_point U = scalar_value_t<T>>
840 const ufcx_form& ufcx_form,
841 const std::vector<std::shared_ptr<const FunctionSpace<U>>>& spaces,
842 const std::map<std::string, std::shared_ptr<const Function<T, U>>>&
843 coefficients,
844 const std::map<std::string, std::shared_ptr<const Constant<T>>>& constants,
845 const std::map<
847 std::vector<std::pair<std::int32_t, std::span<const std::int32_t>>>>&
848 subdomains,
849 const std::vector<std::reference_wrapper<const mesh::EntityMap>>&
850 entity_maps,
851 std::shared_ptr<const mesh::Mesh<U>> mesh = nullptr)
852{
853 // Place coefficients in appropriate order
854 std::vector<std::shared_ptr<const Function<T, U>>> coeff_map;
855 for (const std::string& name : get_coefficient_names(ufcx_form))
856 {
857 if (auto it = coefficients.find(name); it != coefficients.end())
858 coeff_map.push_back(it->second);
859 else
860 {
861 throw std::runtime_error("Form coefficient \"" + name
862 + "\" not provided.");
863 }
864 }
865
866 // Place constants in appropriate order
867 std::vector<std::shared_ptr<const Constant<T>>> const_map;
868 for (const std::string& name : get_constant_names(ufcx_form))
869 {
870 if (auto it = constants.find(name); it != constants.end())
871 const_map.push_back(it->second);
872 else
873 throw std::runtime_error("Form constant \"" + name + "\" not provided.");
874 }
875
876 return create_form_factory({ufcx_form}, spaces, coeff_map, const_map,
877 subdomains, entity_maps, mesh);
878}
879
898template <dolfinx::scalar T, std::floating_point U = scalar_value_t<T>>
900 ufcx_form* (*fptr)(),
901 const std::vector<std::shared_ptr<const FunctionSpace<U>>>& spaces,
902 const std::map<std::string, std::shared_ptr<const Function<T, U>>>&
903 coefficients,
904 const std::map<std::string, std::shared_ptr<const Constant<T>>>& constants,
905 const std::map<
907 std::vector<std::pair<std::int32_t, std::span<const std::int32_t>>>>&
908 subdomains,
909 const std::vector<std::reference_wrapper<const mesh::EntityMap>>&
910 entity_maps,
911 std::shared_ptr<const mesh::Mesh<U>> mesh = nullptr)
912{
913 ufcx_form* form = fptr();
914 Form<T, U> L = create_form<T, U>(*form, spaces, coefficients, constants,
915 subdomains, entity_maps, mesh);
916 std::free(form);
917 return L;
918}
919
921template <std::floating_point T>
923 std::shared_ptr<mesh::Mesh<T>> mesh,
924 std::shared_ptr<const fem::FiniteElement<T>> e,
925 std::function<std::vector<int>(const graph::AdjacencyList<std::int32_t>&)>
926 reorder_fn
927 = nullptr)
928{
929 // TODO: check cell type of e (need to add method to fem::FiniteElement)
930 assert(e);
931 assert(mesh);
932 assert(mesh->topology());
933 if (e->cell_type() != mesh->topology()->cell_type())
934 throw std::runtime_error("Cell type of element and mesh must match.");
935
936 // Create element dof layout
938
939 // Create a dofmap
940 std::function<void(std::span<std::int32_t>, std::uint32_t)> permute_inv
941 = e->needs_dof_permutations() ? e->dof_permutation_fn(true, true)
942 : nullptr;
943 auto dofmap = std::make_shared<const DofMap>(create_dofmap(
944 mesh->comm(), layout, *mesh->topology(), permute_inv, reorder_fn));
945
946 return FunctionSpace(mesh, e, dofmap);
947}
948
950template <dolfinx::scalar T, std::floating_point U = scalar_value_t<T>>
952 const ufcx_expression& e,
953 const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients,
954 const std::vector<std::shared_ptr<const Constant<T>>>& constants,
955 std::shared_ptr<const FunctionSpace<U>> argument_space = nullptr)
956{
957 if (!coefficients.empty())
958 {
959 assert(coefficients.front());
960 assert(coefficients.front()->function_space());
961 std::shared_ptr<const mesh::Mesh<U>> mesh
962 = coefficients.front()->function_space()->mesh();
963 if (mesh->geometry().cmap().hash() != e.coordinate_element_hash)
964 {
965 throw std::runtime_error(
966 "Expression and mesh geometric maps do not match.");
967 }
968 }
969
970 if (e.rank > 0 and !argument_space)
971 {
972 throw std::runtime_error("Expression has Argument but no Argument "
973 "function space was provided.");
974 }
975
976 std::vector<U> X(e.points, e.points + e.num_points * e.entity_dimension);
977 std::array<std::size_t, 2> Xshape
978 = {static_cast<std::size_t>(e.num_points),
979 static_cast<std::size_t>(e.entity_dimension)};
980 std::vector<std::size_t> value_shape(e.value_shape,
981 e.value_shape + e.num_components);
982 std::function<void(T*, const T*, const T*, const scalar_value_t<T>*,
983 const int*, const std::uint8_t*, void*)>
984 tabulate_tensor = nullptr;
985 if constexpr (std::is_same_v<T, float>)
986 tabulate_tensor = e.tabulate_tensor_float32;
987#ifndef DOLFINX_NO_STDC_COMPLEX_KERNELS
988 else if constexpr (std::is_same_v<T, std::complex<float>>)
989 {
990 tabulate_tensor = reinterpret_cast<void (*)(
991 T*, const T*, const T*, const scalar_value_t<T>*, const int*,
992 const unsigned char*, void*)>(e.tabulate_tensor_complex64);
993 }
994#endif // DOLFINX_NO_STDC_COMPLEX_KERNELS
995 else if constexpr (std::is_same_v<T, double>)
996 tabulate_tensor = e.tabulate_tensor_float64;
997#ifndef DOLFINX_NO_STDC_COMPLEX_KERNELS
998 else if constexpr (std::is_same_v<T, std::complex<double>>)
999 {
1000 tabulate_tensor = reinterpret_cast<void (*)(
1001 T*, const T*, const T*, const scalar_value_t<T>*, const int*,
1002 const unsigned char*, void*)>(e.tabulate_tensor_complex128);
1003 }
1004#endif // DOLFINX_NO_STDC_COMPLEX_KERNELS
1005 else
1006 throw std::runtime_error("Type not supported.");
1007
1008 assert(tabulate_tensor);
1009 return Expression(coefficients, constants, std::span<const U>(X), Xshape,
1010 tabulate_tensor, value_shape, argument_space);
1011}
1012
1015template <dolfinx::scalar T, std::floating_point U = scalar_value_t<T>>
1017 const ufcx_expression& e,
1018 const std::map<std::string, std::shared_ptr<const Function<T, U>>>&
1019 coefficients,
1020 const std::map<std::string, std::shared_ptr<const Constant<T>>>& constants,
1021 std::shared_ptr<const FunctionSpace<U>> argument_space = nullptr)
1022{
1023 // Place coefficients in appropriate order
1024 std::vector<std::shared_ptr<const Function<T, U>>> coeff_map;
1025 std::vector<std::string> coefficient_names;
1026 coefficient_names.reserve(e.num_coefficients);
1027 for (int i = 0; i < e.num_coefficients; ++i)
1028 coefficient_names.push_back(e.coefficient_names[i]);
1029
1030 for (const std::string& name : coefficient_names)
1031 {
1032 if (auto it = coefficients.find(name); it != coefficients.end())
1033 coeff_map.push_back(it->second);
1034 else
1035 {
1036 throw std::runtime_error("Expression coefficient \"" + name
1037 + "\" not provided.");
1038 }
1039 }
1040
1041 // Place constants in appropriate order
1042 std::vector<std::shared_ptr<const Constant<T>>> const_map;
1043 std::vector<std::string> constant_names;
1044 constant_names.reserve(e.num_constants);
1045 for (int i = 0; i < e.num_constants; ++i)
1046 constant_names.push_back(e.constant_names[i]);
1047
1048 for (const std::string& name : constant_names)
1049 {
1050 if (auto it = constants.find(name); it != constants.end())
1051 const_map.push_back(it->second);
1052 else
1053 {
1054 throw std::runtime_error("Expression constant \"" + name
1055 + "\" not provided.");
1056 }
1057 }
1058
1059 return create_expression(e, coeff_map, const_map, argument_space);
1060}
1061
1062} // 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
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:824
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
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:28
void cells(la::SparsityPattern &pattern, std::array< std::span< const std::int32_t >, 2 > cells, std::array< std::reference_wrapper< const DofMap >, 2 > dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition sparsitybuild.cpp:16
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:128
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:922
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:411
std::vector< std::string > get_coefficient_names(const ufcx_form &ufcx_form)
Definition utils.cpp:121
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:951
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:136
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:839
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:315
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:70
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