Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d9/dc0/fem_2utils_8h_source.html
DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
utils.h
Go to the documentation of this file.
1// Copyright (C) 2013-2020 Johan Hake, Jan Blechta and Garth N. Wells
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include "Constant.h"
10#include "CoordinateElement.h"
11#include "DofMap.h"
12#include "ElementDofLayout.h"
13#include "Expression.h"
14#include "Form.h"
15#include "Function.h"
16#include "FunctionSpace.h"
17#include "sparsitybuild.h"
18#include <array>
19#include <concepts>
20#include <dolfinx/common/types.h>
21#include <dolfinx/la/SparsityPattern.h>
22#include <dolfinx/mesh/Topology.h>
23#include <dolfinx/mesh/cell_types.h>
24#include <dolfinx/mesh/utils.h>
25#include <functional>
26#include <memory>
27#include <set>
28#include <span>
29#include <stdexcept>
30#include <string>
31#include <type_traits>
32#include <ufcx.h>
33#include <utility>
34#include <vector>
35
38namespace basix
39{
40template <std::floating_point T>
41class FiniteElement;
42}
43
44namespace dolfinx::common
45{
46class IndexMap;
47}
48
49namespace dolfinx::fem
50{
51
52namespace impl
53{
60template <int num_cells>
61std::array<std::int32_t, 2 * num_cells>
62get_cell_facet_pairs(std::int32_t f, std::span<const std::int32_t> cells,
64{
65 // Loop over cells sharing facet
66 assert(cells.size() == num_cells);
67 std::array<std::int32_t, 2 * num_cells> cell_local_facet_pairs;
68 for (int c = 0; c < num_cells; ++c)
69 {
70 // Get local index of facet with respect to the cell
71 std::int32_t cell = cells[c];
72 auto cell_facets = c_to_f.links(cell);
73 auto facet_it = std::find(cell_facets.begin(), cell_facets.end(), f);
74 assert(facet_it != cell_facets.end());
75 int local_f = std::distance(cell_facets.begin(), facet_it);
76 cell_local_facet_pairs[2 * c] = cell;
77 cell_local_facet_pairs[2 * c + 1] = local_f;
78 }
79
80 return cell_local_facet_pairs;
81}
82
83} // namespace impl
84
110std::vector<std::pair<int, std::vector<std::int32_t>>>
112 const mesh::Topology& topology,
113 std::span<const std::int32_t> entities, int dim,
114 std::span<const int> values);
115
123template <dolfinx::scalar T, std::floating_point U>
124std::vector<std::vector<std::array<std::shared_ptr<const FunctionSpace<U>>, 2>>>
125extract_function_spaces(const std::vector<std::vector<const Form<T, U>*>>& a)
126{
127 std::vector<
128 std::vector<std::array<std::shared_ptr<const FunctionSpace<U>>, 2>>>
129 spaces(
130 a.size(),
131 std::vector<std::array<std::shared_ptr<const FunctionSpace<U>>, 2>>(
132 a.front().size()));
133 for (std::size_t i = 0; i < a.size(); ++i)
134 {
135 for (std::size_t j = 0; j < a[i].size(); ++j)
136 {
137 if (const Form<T, U>* form = a[i][j]; form)
138 spaces[i][j] = {form->function_spaces()[0], form->function_spaces()[1]};
139 }
140 }
141 return spaces;
142}
143
149template <dolfinx::scalar T, std::floating_point U>
151{
152 if (a.rank() != 2)
153 {
154 throw std::runtime_error(
155 "Cannot create sparsity pattern. Form is not a bilinear.");
156 }
157
158 // Get dof maps and mesh
159 std::array<std::reference_wrapper<const DofMap>, 2> dofmaps{
160 *a.function_spaces().at(0)->dofmap(),
161 *a.function_spaces().at(1)->dofmap()};
162 std::shared_ptr mesh = a.mesh();
163 assert(mesh);
164
165 std::shared_ptr mesh0 = a.function_spaces().at(0)->mesh();
166 assert(mesh0);
167 std::shared_ptr mesh1 = a.function_spaces().at(1)->mesh();
168 assert(mesh1);
169
170 const std::set<IntegralType> types = a.integral_types();
171 if (types.find(IntegralType::interior_facet) != types.end()
172 or types.find(IntegralType::exterior_facet) != types.end())
173 {
174 // FIXME: cleanup these calls? Some of the happen internally again.
175 int tdim = mesh->topology()->dim();
176 mesh->topology_mutable()->create_entities(tdim - 1);
177 mesh->topology_mutable()->create_connectivity(tdim - 1, tdim);
178 }
179
180 common::Timer t0("Build sparsity");
181
182 // Get common::IndexMaps for each dimension
183 const std::array index_maps{dofmaps[0].get().index_map,
184 dofmaps[1].get().index_map};
185 const std::array bs
186 = {dofmaps[0].get().index_map_bs(), dofmaps[1].get().index_map_bs()};
187
188 auto extract_cells = [](std::span<const std::int32_t> facets)
189 {
190 assert(facets.size() % 2 == 0);
191 std::vector<std::int32_t> cells;
192 cells.reserve(facets.size() / 2);
193 for (std::size_t i = 0; i < facets.size(); i += 2)
194 cells.push_back(facets[i]);
195 return cells;
196 };
197
198 // Create and build sparsity pattern
199 la::SparsityPattern pattern(mesh->comm(), index_maps, bs);
200 for (auto type : types)
201 {
202 std::vector<int> ids = a.integral_ids(type);
203 switch (type)
204 {
206 for (int id : ids)
207 {
209 pattern, {a.domain(type, id, *mesh0), a.domain(type, id, *mesh1)},
210 {{dofmaps[0], dofmaps[1]}});
211 }
212 break;
214 for (int id : ids)
215 {
217 pattern,
218 {extract_cells(a.domain(type, id, *mesh0)),
219 extract_cells(a.domain(type, id, *mesh1))},
220 {{dofmaps[0], dofmaps[1]}});
221 }
222 break;
224 for (int id : ids)
225 {
226 sparsitybuild::cells(pattern,
227 {extract_cells(a.domain(type, id, *mesh0)),
228 extract_cells(a.domain(type, id, *mesh1))},
229 {{dofmaps[0], dofmaps[1]}});
230 }
231 break;
232 default:
233 throw std::runtime_error("Unsupported integral type");
234 }
235 }
236
237 t0.stop();
238
239 return pattern;
240}
241
243ElementDofLayout create_element_dof_layout(const ufcx_dofmap& dofmap,
244 const mesh::CellType cell_type,
245 const std::vector<int>& parent_map
246 = {});
247
256DofMap create_dofmap(
257 MPI_Comm comm, const ElementDofLayout& layout, mesh::Topology& topology,
258 std::function<void(std::span<std::int32_t>, std::uint32_t)> permute_inv,
259 std::function<std::vector<int>(const graph::AdjacencyList<std::int32_t>&)>
260 reorder_fn);
261
272std::vector<DofMap> create_dofmaps(
273 MPI_Comm comm, const std::vector<ElementDofLayout>& layouts,
274 mesh::Topology& topology,
275 std::function<void(std::span<std::int32_t>, std::uint32_t)> permute_inv,
276 std::function<std::vector<int>(const graph::AdjacencyList<std::int32_t>&)>
277 reorder_fn);
278
282std::vector<std::string> get_coefficient_names(const ufcx_form& ufcx_form);
283
287std::vector<std::string> get_constant_names(const ufcx_form& ufcx_form);
288
306template <dolfinx::scalar T, std::floating_point U = scalar_value_type_t<T>>
308 const ufcx_form& ufcx_form,
309 const std::vector<std::shared_ptr<const FunctionSpace<U>>>& spaces,
310 const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients,
311 const std::vector<std::shared_ptr<const Constant<T>>>& constants,
312 const std::map<
314 std::vector<std::pair<std::int32_t, std::span<const std::int32_t>>>>&
315 subdomains,
316 const std::map<std::shared_ptr<const mesh::Mesh<U>>,
317 std::span<const std::int32_t>>& entity_maps,
318 std::shared_ptr<const mesh::Mesh<U>> mesh = nullptr)
319{
320 if (ufcx_form.rank != (int)spaces.size())
321 throw std::runtime_error("Wrong number of argument spaces for Form.");
322 if (ufcx_form.num_coefficients != (int)coefficients.size())
323 {
324 throw std::runtime_error(
325 "Mismatch between number of expected and provided Form coefficients.");
326 }
327 if (ufcx_form.num_constants != (int)constants.size())
328 {
329 throw std::runtime_error(
330 "Mismatch between number of expected and provided Form constants.");
331 }
332
333 // Check argument function spaces
334#ifndef NDEBUG
335 for (std::size_t i = 0; i < spaces.size(); ++i)
336 {
337 assert(spaces[i]->element());
338 ufcx_finite_element* ufcx_element = ufcx_form.finite_elements[i];
339 assert(ufcx_element);
340 if (std::string(ufcx_element->signature)
341 != spaces[i]->element()->signature())
342 {
343 throw std::runtime_error(
344 "Cannot create form. Wrong type of function space for argument.");
345 }
346 }
347#endif
348
349 // Extract mesh from FunctionSpace, and check they are the same
350 if (!mesh and !spaces.empty())
351 mesh = spaces[0]->mesh();
352 for (auto& V : spaces)
353 {
354 if (mesh != V->mesh() and entity_maps.find(V->mesh()) == entity_maps.end())
355 throw std::runtime_error(
356 "Incompatible mesh. entity_maps must be provided.");
357 }
358 if (!mesh)
359 throw std::runtime_error("No mesh could be associated with the Form.");
360
361 auto topology = mesh->topology();
362 assert(topology);
363 const int tdim = topology->dim();
364
365 const int* integral_offsets = ufcx_form.form_integral_offsets;
366 std::vector<int> num_integrals_type(3);
367 for (int i = 0; i < 3; ++i)
368 num_integrals_type[i] = integral_offsets[i + 1] - integral_offsets[i];
369
370 // Create facets, if required
371 if (num_integrals_type[exterior_facet] > 0
372 or num_integrals_type[interior_facet] > 0)
373 {
374 mesh->topology_mutable()->create_entities(tdim - 1);
375 mesh->topology_mutable()->create_connectivity(tdim - 1, tdim);
376 mesh->topology_mutable()->create_connectivity(tdim, tdim - 1);
377 }
378
379 // Get list of integral IDs, and load tabulate tensor into memory for
380 // each
381 using kern_t = std::function<void(T*, const T*, const T*, const U*,
382 const int*, const std::uint8_t*)>;
383 std::map<IntegralType, std::vector<integral_data<T, U>>> integrals;
384
385 // Attach cell kernels
386 bool needs_facet_permutations = false;
387 std::vector<std::int32_t> default_cells;
388 {
389 std::span<const int> ids(ufcx_form.form_integral_ids
390 + integral_offsets[cell],
391 num_integrals_type[cell]);
392 auto itg = integrals.insert({IntegralType::cell, {}});
393 auto sd = subdomains.find(IntegralType::cell);
394 for (int i = 0; i < num_integrals_type[cell]; ++i)
395 {
396 const int id = ids[i];
397 ufcx_integral* integral
398 = ufcx_form.form_integrals[integral_offsets[cell] + i];
399 assert(integral);
400
401 kern_t k = nullptr;
402 if constexpr (std::is_same_v<T, float>)
403 k = integral->tabulate_tensor_float32;
404 else if constexpr (std::is_same_v<T, std::complex<float>>)
405 {
406 k = reinterpret_cast<void (*)(
407 T*, const T*, const T*,
408 const typename scalar_value_type<T>::value_type*, const int*,
409 const unsigned char*)>(integral->tabulate_tensor_complex64);
410 }
411 else if constexpr (std::is_same_v<T, double>)
412 k = integral->tabulate_tensor_float64;
413 else if constexpr (std::is_same_v<T, std::complex<double>>)
414 {
415 k = reinterpret_cast<void (*)(
416 T*, const T*, const T*,
417 const typename scalar_value_type<T>::value_type*, const int*,
418 const unsigned char*)>(integral->tabulate_tensor_complex128);
419 }
420 if (!k)
421 {
422 throw std::runtime_error(
423 "UFCx kernel function is NULL. Check requested types.");
424 }
425
426 // Build list of entities to assemble over
427 if (id == -1)
428 {
429 // Default kernel, operates on all (owned) cells
430 assert(topology->index_map(tdim));
431 default_cells.resize(topology->index_map(tdim)->size_local(), 0);
432 std::iota(default_cells.begin(), default_cells.end(), 0);
433 itg.first->second.emplace_back(id, k, default_cells);
434 }
435 else if (sd != subdomains.end())
436 {
437 // NOTE: This requires that pairs are sorted
438 auto it = std::lower_bound(sd->second.begin(), sd->second.end(), id,
439 [](auto& pair, auto val)
440 { return pair.first < val; });
441 if (it != sd->second.end() and it->first == id)
442 itg.first->second.emplace_back(id, k, it->second);
443 }
444
445 if (integral->needs_facet_permutations)
446 needs_facet_permutations = true;
447 }
448 }
449
450 // Attach exterior facet kernels
451 std::vector<std::int32_t> default_facets_ext;
452 {
453 std::span<const int> ids(ufcx_form.form_integral_ids
454 + integral_offsets[exterior_facet],
455 num_integrals_type[exterior_facet]);
456 auto itg = integrals.insert({IntegralType::exterior_facet, {}});
457 auto sd = subdomains.find(IntegralType::exterior_facet);
458 for (int i = 0; i < num_integrals_type[exterior_facet]; ++i)
459 {
460 const int id = ids[i];
461 ufcx_integral* integral
462 = ufcx_form.form_integrals[integral_offsets[exterior_facet] + i];
463 assert(integral);
464
465 kern_t k = nullptr;
466 if constexpr (std::is_same_v<T, float>)
467 k = integral->tabulate_tensor_float32;
468 else if constexpr (std::is_same_v<T, std::complex<float>>)
469 {
470 k = reinterpret_cast<void (*)(
471 T*, const T*, const T*,
472 const typename scalar_value_type<T>::value_type*, const int*,
473 const unsigned char*)>(integral->tabulate_tensor_complex64);
474 }
475 else if constexpr (std::is_same_v<T, double>)
476 k = integral->tabulate_tensor_float64;
477 else if constexpr (std::is_same_v<T, std::complex<double>>)
478 {
479 k = reinterpret_cast<void (*)(
480 T*, const T*, const T*,
481 const typename scalar_value_type<T>::value_type*, const int*,
482 const unsigned char*)>(integral->tabulate_tensor_complex128);
483 }
484 assert(k);
485
486 // Build list of entities to assembler over
487 const std::vector bfacets = mesh::exterior_facet_indices(*topology);
488 auto f_to_c = topology->connectivity(tdim - 1, tdim);
489 assert(f_to_c);
490 auto c_to_f = topology->connectivity(tdim, tdim - 1);
491 assert(c_to_f);
492 if (id == -1)
493 {
494 // Default kernel, operates on all (owned) exterior facets
495 default_facets_ext.reserve(2 * bfacets.size());
496 for (std::int32_t f : bfacets)
497 {
498 // There will only be one pair for an exterior facet integral
499 auto pair
500 = impl::get_cell_facet_pairs<1>(f, f_to_c->links(f), *c_to_f);
501 default_facets_ext.insert(default_facets_ext.end(), pair.begin(),
502 pair.end());
503 }
504 itg.first->second.emplace_back(id, k, default_facets_ext);
505 }
506 else if (sd != subdomains.end())
507 {
508 // NOTE: This requires that pairs are sorted
509 auto it = std::lower_bound(sd->second.begin(), sd->second.end(), id,
510 [](auto& pair, auto val)
511 { return pair.first < val; });
512 if (it != sd->second.end() and it->first == id)
513 itg.first->second.emplace_back(id, k, it->second);
514 }
515
516 if (integral->needs_facet_permutations)
517 needs_facet_permutations = true;
518 }
519 }
520
521 // Attach interior facet kernels
522 std::vector<std::int32_t> default_facets_int;
523 {
524 std::span<const int> ids(ufcx_form.form_integral_ids
525 + integral_offsets[interior_facet],
526 num_integrals_type[interior_facet]);
527 auto itg = integrals.insert({IntegralType::interior_facet, {}});
528 auto sd = subdomains.find(IntegralType::interior_facet);
529 for (int i = 0; i < num_integrals_type[interior_facet]; ++i)
530 {
531 const int id = ids[i];
532 ufcx_integral* integral
533 = ufcx_form.form_integrals[integral_offsets[interior_facet] + i];
534 assert(integral);
535
536 kern_t k = nullptr;
537 if constexpr (std::is_same_v<T, float>)
538 k = integral->tabulate_tensor_float32;
539 else if constexpr (std::is_same_v<T, std::complex<float>>)
540 {
541 k = reinterpret_cast<void (*)(
542 T*, const T*, const T*,
543 const typename scalar_value_type<T>::value_type*, const int*,
544 const unsigned char*)>(integral->tabulate_tensor_complex64);
545 }
546 else if constexpr (std::is_same_v<T, double>)
547 k = integral->tabulate_tensor_float64;
548 else if constexpr (std::is_same_v<T, std::complex<double>>)
549 {
550 k = reinterpret_cast<void (*)(
551 T*, const T*, const T*,
552 const typename scalar_value_type<T>::value_type*, const int*,
553 const unsigned char*)>(integral->tabulate_tensor_complex128);
554 }
555 assert(k);
556
557 // Build list of entities to assembler over
558 auto f_to_c = topology->connectivity(tdim - 1, tdim);
559 assert(f_to_c);
560 auto c_to_f = topology->connectivity(tdim, tdim - 1);
561 assert(c_to_f);
562 if (id == -1)
563 {
564 // Default kernel, operates on all (owned) interior facets
565 assert(topology->index_map(tdim - 1));
566 std::int32_t num_facets = topology->index_map(tdim - 1)->size_local();
567 default_facets_int.reserve(4 * num_facets);
568 for (std::int32_t f = 0; f < num_facets; ++f)
569 {
570 if (f_to_c->num_links(f) == 2)
571 {
572 auto pairs
573 = impl::get_cell_facet_pairs<2>(f, f_to_c->links(f), *c_to_f);
574 default_facets_int.insert(default_facets_int.end(), pairs.begin(),
575 pairs.end());
576 }
577 }
578 itg.first->second.emplace_back(id, k, default_facets_int);
579 }
580 else if (sd != subdomains.end())
581 {
582 auto it = std::lower_bound(sd->second.begin(), sd->second.end(), id,
583 [](auto& pair, auto val)
584 { return pair.first < val; });
585 if (it != sd->second.end() and it->first == id)
586 itg.first->second.emplace_back(id, k, it->second);
587 }
588
589 if (integral->needs_facet_permutations)
590 needs_facet_permutations = true;
591 }
592 }
593
594 std::map<IntegralType,
595 std::vector<std::pair<std::int32_t, std::vector<std::int32_t>>>>
596 sd;
597 for (auto& [itg, data] : subdomains)
598 {
599 std::vector<std::pair<std::int32_t, std::vector<std::int32_t>>> x;
600 for (auto& [id, idx] : data)
601 x.emplace_back(id, std::vector(idx.data(), idx.data() + idx.size()));
602 sd.insert({itg, std::move(x)});
603 }
604
605 return Form<T, U>(spaces, integrals, coefficients, constants,
606 needs_facet_permutations, entity_maps, mesh);
607}
608
620template <dolfinx::scalar T, std::floating_point U = scalar_value_type_t<T>>
622 const ufcx_form& ufcx_form,
623 const std::vector<std::shared_ptr<const FunctionSpace<U>>>& spaces,
624 const std::map<std::string, std::shared_ptr<const Function<T, U>>>&
625 coefficients,
626 const std::map<std::string, std::shared_ptr<const Constant<T>>>& constants,
627 const std::map<
629 std::vector<std::pair<std::int32_t, std::span<const std::int32_t>>>>&
630 subdomains,
631 std::shared_ptr<const mesh::Mesh<U>> mesh = nullptr)
632{
633 // Place coefficients in appropriate order
634 std::vector<std::shared_ptr<const Function<T, U>>> coeff_map;
635 for (const std::string& name : get_coefficient_names(ufcx_form))
636 {
637 if (auto it = coefficients.find(name); it != coefficients.end())
638 coeff_map.push_back(it->second);
639 else
640 {
641 throw std::runtime_error("Form coefficient \"" + name
642 + "\" not provided.");
643 }
644 }
645
646 // Place constants in appropriate order
647 std::vector<std::shared_ptr<const Constant<T>>> const_map;
648 for (const std::string& name : get_constant_names(ufcx_form))
649 {
650 if (auto it = constants.find(name); it != constants.end())
651 const_map.push_back(it->second);
652 else
653 throw std::runtime_error("Form constant \"" + name + "\" not provided.");
654 }
655
656 return create_form_factory(ufcx_form, spaces, coeff_map, const_map,
657 subdomains, {}, mesh);
658}
659
675template <dolfinx::scalar T, std::floating_point U = scalar_value_type_t<T>>
677 ufcx_form* (*fptr)(),
678 const std::vector<std::shared_ptr<const FunctionSpace<U>>>& spaces,
679 const std::map<std::string, std::shared_ptr<const Function<T, U>>>&
680 coefficients,
681 const std::map<std::string, std::shared_ptr<const Constant<T>>>& constants,
682 const std::map<
684 std::vector<std::pair<std::int32_t, std::span<const std::int32_t>>>>&
685 subdomains,
686 std::shared_ptr<const mesh::Mesh<U>> mesh = nullptr)
687{
688 ufcx_form* form = fptr();
689 Form<T, U> L = create_form<T, U>(*form, spaces, coefficients, constants,
690 subdomains, mesh);
691 std::free(form);
692 return L;
693}
694
706template <std::floating_point T>
708 std::shared_ptr<mesh::Mesh<T>> mesh, const basix::FiniteElement<T>& e,
709 const std::vector<std::size_t>& value_shape = {},
710 std::function<std::vector<int>(const graph::AdjacencyList<std::int32_t>&)>
711 reorder_fn
712 = nullptr)
713{
714 if (!e.value_shape().empty() and !value_shape.empty())
715 {
716 throw std::runtime_error(
717 "Cannot specify value shape for non-scalar base element.");
718 }
719
720 std::size_t bs = value_shape.empty()
721 ? 1
722 : std::accumulate(value_shape.begin(), value_shape.end(),
723 1, std::multiplies{});
724
725 // Create a DOLFINx element
726 auto _e = std::make_shared<const FiniteElement<T>>(e, bs);
727 assert(_e);
728
729 const std::vector<std::size_t> _value_shape
730 = (value_shape.empty() and !e.value_shape().empty())
731 ? fem::compute_value_shape(_e, mesh->topology()->dim(),
732 mesh->geometry().dim())
733 : value_shape;
734
735 // Create UFC subdofmaps and compute offset
736 const int num_sub_elements = _e->num_sub_elements();
737 std::vector<ElementDofLayout> sub_doflayout;
738 sub_doflayout.reserve(num_sub_elements);
739 for (int i = 0; i < num_sub_elements; ++i)
740 {
741 auto sub_element = _e->extract_sub_element({i});
742 std::vector<int> parent_map_sub(sub_element->space_dimension());
743 for (std::size_t j = 0; j < parent_map_sub.size(); ++j)
744 parent_map_sub[j] = i + _e->block_size() * j;
745 sub_doflayout.emplace_back(1, e.entity_dofs(), e.entity_closure_dofs(),
746 parent_map_sub, std::vector<ElementDofLayout>());
747 }
748
749 // Create a dofmap
750 ElementDofLayout layout(_e->block_size(), e.entity_dofs(),
751 e.entity_closure_dofs(), {}, sub_doflayout);
752 std::function<void(std::span<std::int32_t>, std::uint32_t)> permute_inv
753 = nullptr;
754 if (_e->needs_dof_permutations())
755 permute_inv = _e->dof_permutation_fn(true, true);
756 assert(mesh);
757 assert(mesh->topology());
758 auto dofmap = std::make_shared<const DofMap>(create_dofmap(
759 mesh->comm(), layout, *mesh->topology(), permute_inv, reorder_fn));
760 return FunctionSpace(mesh, _e, dofmap, _value_shape);
761}
762
773template <std::floating_point T>
775 ufcx_function_space* (*fptr)(const char*), const std::string& function_name,
776 std::shared_ptr<mesh::Mesh<T>> mesh,
777 std::function<std::vector<int>(const graph::AdjacencyList<std::int32_t>&)>
778 reorder_fn
779 = nullptr)
780{
781 ufcx_function_space* space = fptr(function_name.c_str());
782 if (!space)
783 {
784 throw std::runtime_error(
785 "Could not create UFC function space with function name "
786 + function_name);
787 }
788
789 ufcx_finite_element* ufcx_element = space->finite_element;
790 assert(ufcx_element);
791 std::vector<std::size_t> value_shape(space->value_shape,
792 space->value_shape + space->value_rank);
793
794 const auto& geometry = mesh->geometry();
795 auto& cmap = geometry.cmap();
796 if (space->geometry_degree != cmap.degree()
797 or static_cast<basix::cell::type>(space->geometry_basix_cell)
798 != mesh::cell_type_to_basix_type(cmap.cell_shape())
799 or static_cast<basix::element::lagrange_variant>(
800 space->geometry_basix_variant)
801 != cmap.variant())
802 {
803 throw std::runtime_error("UFL mesh and CoordinateElement do not match.");
804 }
805
806 auto element = std::make_shared<FiniteElement<T>>(*ufcx_element);
807 assert(element);
808 ufcx_dofmap* ufcx_map = space->dofmap;
809 assert(ufcx_map);
810 const auto topology = mesh->topology();
811 assert(topology);
812 ElementDofLayout layout
813 = create_element_dof_layout(*ufcx_map, topology->cell_type());
814
815 std::function<void(std::span<std::int32_t>, std::uint32_t)> permute_inv;
816 if (element->needs_dof_permutations())
817 permute_inv = element->dof_permutation_fn(true, true);
818 return FunctionSpace(
819 mesh, element,
820 std::make_shared<DofMap>(create_dofmap(mesh->comm(), layout, *topology,
821 permute_inv, reorder_fn)),
822 value_shape);
823}
824
826namespace impl
827{
829template <dolfinx::scalar T, std::floating_point U>
830std::span<const std::uint32_t>
831get_cell_orientation_info(const Function<T, U>& coefficient)
832{
833 std::span<const std::uint32_t> cell_info;
834 auto element = coefficient.function_space()->element();
835 if (element->needs_dof_transformations())
836 {
837 auto mesh = coefficient.function_space()->mesh();
838 mesh->topology_mutable()->create_entity_permutations();
839 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
840 }
841
842 return cell_info;
843}
844
846template <dolfinx::scalar T, int _bs>
847void pack(std::span<T> coeffs, std::int32_t cell, int bs, std::span<const T> v,
848 std::span<const std::uint32_t> cell_info, const DofMap& dofmap,
849 auto transform)
850{
851 auto dofs = dofmap.cell_dofs(cell);
852 for (std::size_t i = 0; i < dofs.size(); ++i)
853 {
854 if constexpr (_bs < 0)
855 {
856 const int pos_c = bs * i;
857 const int pos_v = bs * dofs[i];
858 for (int k = 0; k < bs; ++k)
859 coeffs[pos_c + k] = v[pos_v + k];
860 }
861 else
862 {
863 const int pos_c = _bs * i;
864 const int pos_v = _bs * dofs[i];
865 for (int k = 0; k < _bs; ++k)
866 coeffs[pos_c + k] = v[pos_v + k];
867 }
868 }
869
870 transform(coeffs, cell_info, cell, 1);
871}
872
875template <typename F>
876concept FetchCells = requires(F&& f, std::span<const std::int32_t> v) {
877 requires std::invocable<F, std::span<const std::int32_t>>;
878 { f(v) } -> std::convertible_to<std::int32_t>;
879};
880
894template <dolfinx::scalar T, std::floating_point U>
895void pack_coefficient_entity(std::span<T> c, int cstride,
896 const Function<T, U>& u,
897 std::span<const std::uint32_t> cell_info,
898 std::span<const std::int32_t> entities,
899 std::size_t estride, FetchCells auto&& fetch_cells,
900 std::int32_t offset)
901{
902 // Read data from coefficient Function u
903 std::span<const T> v = u.x()->array();
904 const DofMap& dofmap = *u.function_space()->dofmap();
905 auto element = u.function_space()->element();
906 assert(element);
907 int space_dim = element->space_dimension();
908
909 // Transformation from conforming degrees-of-freedom to reference
910 // degrees-of-freedom
911 auto transformation
912 = element->template dof_transformation_fn<T>(doftransform::transpose);
913 const int bs = dofmap.bs();
914 switch (bs)
915 {
916 case 1:
917 for (std::size_t e = 0; e < entities.size(); e += estride)
918 {
919 auto entity = entities.subspan(e, estride);
920 std::int32_t cell = fetch_cells(entity);
921 auto cell_coeff = c.subspan((e / estride) * cstride + offset, space_dim);
922 pack<T, 1>(cell_coeff, cell, bs, v, cell_info, dofmap, transformation);
923 }
924 break;
925 case 2:
926 for (std::size_t e = 0; e < entities.size(); e += estride)
927 {
928 auto entity = entities.subspan(e, estride);
929 std::int32_t cell = fetch_cells(entity);
930 auto cell_coeff = c.subspan((e / estride) * cstride + offset, space_dim);
931 pack<T, 2>(cell_coeff, cell, bs, v, cell_info, dofmap, transformation);
932 }
933 break;
934 case 3:
935 for (std::size_t e = 0; e < entities.size(); e += estride)
936 {
937 auto entity = entities.subspan(e, estride);
938 std::int32_t cell = fetch_cells(entity);
939 auto cell_coeff = c.subspan(e / estride * cstride + offset, space_dim);
940 pack<T, 3>(cell_coeff, cell, bs, v, cell_info, dofmap, transformation);
941 }
942 break;
943 default:
944 for (std::size_t e = 0; e < entities.size(); e += estride)
945 {
946 auto entity = entities.subspan(e, estride);
947 std::int32_t cell = fetch_cells(entity);
948 auto cell_coeff = c.subspan((e / estride) * cstride + offset, space_dim);
949 pack<T, -1>(cell_coeff, cell, bs, v, cell_info, dofmap, transformation);
950 }
951 break;
952 }
953}
954
955} // namespace impl
956
963template <dolfinx::scalar T, std::floating_point U>
964std::pair<std::vector<T>, int>
966 int id)
967{
968 // Get form coefficient offsets and dofmaps
969 const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients
970 = form.coefficients();
971 const std::vector<int> offsets = form.coefficient_offsets();
972
973 std::size_t num_entities = 0;
974 int cstride = 0;
975 if (!coefficients.empty())
976 {
977 cstride = offsets.back();
978 num_entities = form.domain(integral_type, id).size();
979 if (integral_type == IntegralType::exterior_facet
980 or integral_type == IntegralType::interior_facet)
981 {
982 num_entities /= 2;
983 }
984 }
985
986 return {std::vector<T>(num_entities * cstride), cstride};
987}
988
993template <dolfinx::scalar T, std::floating_point U>
994std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, int>>
996{
997 std::map<std::pair<IntegralType, int>, std::pair<std::vector<T>, int>> coeffs;
998 for (auto integral_type : form.integral_types())
999 {
1000 for (int id : form.integral_ids(integral_type))
1001 {
1002 coeffs.emplace_hint(
1003 coeffs.end(), std::pair(integral_type, id),
1004 allocate_coefficient_storage(form, integral_type, id));
1005 }
1006 }
1007
1008 return coeffs;
1009}
1010
1018template <dolfinx::scalar T, std::floating_point U>
1019void pack_coefficients(const Form<T, U>& form, IntegralType integral_type,
1020 int id, std::span<T> c, int cstride)
1021{
1022 // Get form coefficient offsets and dofmaps
1023 const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients
1024 = form.coefficients();
1025 const std::vector<int> offsets = form.coefficient_offsets();
1026
1027 if (!coefficients.empty())
1028 {
1029 switch (integral_type)
1030 {
1031 case IntegralType::cell:
1032 {
1033 // Iterate over coefficients
1034 for (std::size_t coeff = 0; coeff < coefficients.size(); ++coeff)
1035 {
1036 auto mesh = coefficients[coeff]->function_space()->mesh();
1037 assert(mesh);
1038 std::vector<std::int32_t> cells
1039 = form.domain(IntegralType::cell, id, *mesh);
1040 std::span<const std::uint32_t> cell_info
1041 = impl::get_cell_orientation_info(*coefficients[coeff]);
1042 impl::pack_coefficient_entity(
1043 c, cstride, *coefficients[coeff], cell_info, cells, 1,
1044 [](auto entity) { return entity.front(); }, offsets[coeff]);
1045 }
1046 break;
1047 }
1049 {
1050 // Iterate over coefficients
1051 for (std::size_t coeff = 0; coeff < coefficients.size(); ++coeff)
1052 {
1053 auto mesh = coefficients[coeff]->function_space()->mesh();
1054 std::vector<std::int32_t> facets
1055 = form.domain(IntegralType::exterior_facet, id, *mesh);
1056 std::span<const std::uint32_t> cell_info
1057 = impl::get_cell_orientation_info(*coefficients[coeff]);
1058 impl::pack_coefficient_entity(
1059 c, cstride, *coefficients[coeff], cell_info, facets, 2,
1060 [](auto entity) { return entity.front(); }, offsets[coeff]);
1061 }
1062 break;
1063 }
1065 {
1066 // Iterate over coefficients
1067 for (std::size_t coeff = 0; coeff < coefficients.size(); ++coeff)
1068 {
1069 auto mesh = coefficients[coeff]->function_space()->mesh();
1070 std::vector<std::int32_t> facets
1071 = form.domain(IntegralType::interior_facet, id, *mesh);
1072 std::span<const std::uint32_t> cell_info
1073 = impl::get_cell_orientation_info(*coefficients[coeff]);
1074
1075 // Pack coefficient ['+']
1076 impl::pack_coefficient_entity(
1077 c, 2 * cstride, *coefficients[coeff], cell_info, facets, 4,
1078 [](auto entity) { return entity[0]; }, 2 * offsets[coeff]);
1079 // Pack coefficient ['-']
1080 impl::pack_coefficient_entity(
1081 c, 2 * cstride, *coefficients[coeff], cell_info, facets, 4,
1082 [](auto entity) { return entity[2]; },
1083 offsets[coeff] + offsets[coeff + 1]);
1084 }
1085 break;
1086 }
1087 default:
1088 throw std::runtime_error(
1089 "Could not pack coefficient. Integral type not supported.");
1090 }
1091 }
1092}
1093
1095template <dolfinx::scalar T, std::floating_point U = scalar_value_type_t<T>>
1097 const ufcx_expression& e,
1098 const std::vector<std::shared_ptr<const Function<T, U>>>& coefficients,
1099 const std::vector<std::shared_ptr<const Constant<T>>>& constants,
1100 std::shared_ptr<const FunctionSpace<U>> argument_function_space = nullptr)
1101{
1102 if (e.rank > 0 and !argument_function_space)
1103 {
1104 throw std::runtime_error("Expression has Argument but no Argument "
1105 "function space was provided.");
1106 }
1107
1108 std::vector<U> X(e.points, e.points + e.num_points * e.entity_dimension);
1109 std::array<std::size_t, 2> Xshape
1110 = {static_cast<std::size_t>(e.num_points),
1111 static_cast<std::size_t>(e.entity_dimension)};
1112 std::vector<int> value_shape(e.value_shape, e.value_shape + e.num_components);
1113 std::function<void(T*, const T*, const T*,
1114 const typename scalar_value_type<T>::value_type*,
1115 const int*, const std::uint8_t*)>
1116 tabulate_tensor = nullptr;
1117 if constexpr (std::is_same_v<T, float>)
1118 tabulate_tensor = e.tabulate_tensor_float32;
1119 else if constexpr (std::is_same_v<T, std::complex<float>>)
1120 {
1121 tabulate_tensor = reinterpret_cast<void (*)(
1122 T*, const T*, const T*,
1123 const typename scalar_value_type<T>::value_type*, const int*,
1124 const unsigned char*)>(e.tabulate_tensor_complex64);
1125 }
1126 else if constexpr (std::is_same_v<T, double>)
1127 tabulate_tensor = e.tabulate_tensor_float64;
1128 else if constexpr (std::is_same_v<T, std::complex<double>>)
1129 {
1130 tabulate_tensor = reinterpret_cast<void (*)(
1131 T*, const T*, const T*,
1132 const typename scalar_value_type<T>::value_type*, const int*,
1133 const unsigned char*)>(e.tabulate_tensor_complex128);
1134 }
1135 else
1136 throw std::runtime_error("Type not supported.");
1137
1138 assert(tabulate_tensor);
1139 return Expression(coefficients, constants, std::span<const U>(X), Xshape,
1140 tabulate_tensor, value_shape, argument_function_space);
1141}
1142
1145template <dolfinx::scalar T, std::floating_point U = scalar_value_type_t<T>>
1147 const ufcx_expression& e,
1148 const std::map<std::string, std::shared_ptr<const Function<T, U>>>&
1149 coefficients,
1150 const std::map<std::string, std::shared_ptr<const Constant<T>>>& constants,
1151 std::shared_ptr<const FunctionSpace<U>> argument_function_space = nullptr)
1152{
1153 // Place coefficients in appropriate order
1154 std::vector<std::shared_ptr<const Function<T, U>>> coeff_map;
1155 std::vector<std::string> coefficient_names;
1156 for (int i = 0; i < e.num_coefficients; ++i)
1157 coefficient_names.push_back(e.coefficient_names[i]);
1158
1159 for (const std::string& name : coefficient_names)
1160 {
1161 if (auto it = coefficients.find(name); it != coefficients.end())
1162 coeff_map.push_back(it->second);
1163 else
1164 {
1165 throw std::runtime_error("Expression coefficient \"" + name
1166 + "\" not provided.");
1167 }
1168 }
1169
1170 // Place constants in appropriate order
1171 std::vector<std::shared_ptr<const Constant<T>>> const_map;
1172 std::vector<std::string> constant_names;
1173 for (int i = 0; i < e.num_constants; ++i)
1174 constant_names.push_back(e.constant_names[i]);
1175
1176 for (const std::string& name : constant_names)
1177 {
1178 if (auto it = constants.find(name); it != constants.end())
1179 const_map.push_back(it->second);
1180 else
1181 {
1182 throw std::runtime_error("Expression constant \"" + name
1183 + "\" not provided.");
1184 }
1185 }
1186
1187 return create_expression(e, coeff_map, const_map, argument_function_space);
1188}
1189
1200template <dolfinx::scalar T, std::floating_point U>
1202 std::map<std::pair<IntegralType, int>,
1203 std::pair<std::vector<T>, int>>& coeffs)
1204{
1205 for (auto& [key, val] : coeffs)
1206 pack_coefficients<T>(form, key.first, key.second, val.first, val.second);
1207}
1208
1217template <dolfinx::scalar T, std::floating_point U>
1218std::pair<std::vector<T>, int>
1220 std::span<const std::int32_t> entities, std::size_t estride)
1221{
1222 // Get form coefficient offsets and dofmaps
1223 const std::vector<std::shared_ptr<const Function<T, U>>>& coeffs
1224 = e.coefficients();
1225 const std::vector<int> offsets = e.coefficient_offsets();
1226
1227 // Copy data into coefficient array
1228 const int cstride = offsets.back();
1229 std::vector<T> c(entities.size() / estride * offsets.back());
1230 if (!coeffs.empty())
1231 {
1232 // Iterate over coefficients
1233 for (std::size_t coeff = 0; coeff < coeffs.size(); ++coeff)
1234 {
1235 std::span<const std::uint32_t> cell_info
1236 = impl::get_cell_orientation_info(*coeffs[coeff]);
1237
1238 impl::pack_coefficient_entity(
1239 std::span(c), cstride, *coeffs[coeff], cell_info, entities, estride,
1240 [](auto entity) { return entity[0]; }, offsets[coeff]);
1241 }
1242 }
1243 return {std::move(c), cstride};
1244}
1245
1248template <typename U>
1249std::vector<typename U::scalar_type> pack_constants(const U& u)
1250{
1251 using T = typename U::scalar_type;
1252 const std::vector<std::shared_ptr<const Constant<T>>>& constants
1253 = u.constants();
1254
1255 // Calculate size of array needed to store packed constants
1256 std::int32_t size = std::accumulate(constants.cbegin(), constants.cend(), 0,
1257 [](std::int32_t sum, auto& constant)
1258 { return sum + constant->value.size(); });
1259
1260 // Pack constants
1261 std::vector<T> constant_values(size);
1262 std::int32_t offset = 0;
1263 for (auto& constant : constants)
1264 {
1265 const std::vector<T>& value = constant->value;
1266 std::copy(value.begin(), value.end(),
1267 std::next(constant_values.begin(), offset));
1268 offset += value.size();
1269 }
1270
1271 return constant_values;
1272}
1273
1274} // namespace dolfinx::fem
Degree-of-freedom map representations and tools.
Definition CoordinateElement.h:26
Definition Timer.h:31
double stop()
Definition Timer.cpp:45
Constant value which can be attached to a Form.
Definition Constant.h:23
Degree-of-freedom map.
Definition DofMap.h:76
std::span< const std::int32_t > cell_dofs(std::int32_t c) const
Local-to-global mapping of dofs on a cell.
Definition DofMap.h:130
int bs() const noexcept
Return the block size for the dofmap.
Definition DofMap.cpp:173
Definition ElementDofLayout.h:30
Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell.
Definition Expression.h:41
A representation of finite element variational forms.
Definition Form.h:120
std::vector< int > coefficient_offsets() const
Offset for each coefficient expansion array on a cell.
Definition Form.h:387
std::span< const std::int32_t > domain(IntegralType type, int i) const
Get the list of mesh entity indices for the ith integral (kernel) of a given type.
Definition Form.h:307
std::vector< int > integral_ids(IntegralType type) const
Get the IDs for integrals (kernels) for given integral type.
Definition Form.h:281
std::set< IntegralType > integral_types() const
Get types of integrals in the form.
Definition Form.h:254
const std::vector< std::shared_ptr< const Function< T, U > > > & coefficients() const
Access coefficients.
Definition Form.h:373
This class represents a finite element function space defined by a mesh, a finite element,...
Definition FunctionSpace.h:34
Definition Function.h:45
Definition AdjacencyList.h:28
std::span< T > links(std::size_t node)
Definition AdjacencyList.h:112
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:44
Concepts for function that returns cell index.
Definition utils.h:876
void pack(std::span< T > coeffs, std::int32_t cell, int bs, std::span< const T > v, std::span< const std::uint32_t > cell_info, const DofMap &dofmap, auto transform)
Pack a single coefficient for a single cell.
Definition utils.h:847
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:62
void pack_coefficient_entity(std::span< T > c, int cstride, const Function< T, U > &u, std::span< const std::uint32_t > cell_info, std::span< const std::int32_t > entities, std::size_t estride, FetchCells auto &&fetch_cells, std::int32_t offset)
Pack a single coefficient for a set of active entities.
Definition utils.h:895
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:27
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:15
Finite element method functionality.
Definition assemble_matrix_impl.h:26
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:191
std::pair< std::vector< T >, int > allocate_coefficient_storage(const Form< T, U > &form, IntegralType integral_type, int id)
Allocate storage for coefficients of a pair (integral_type, id) from a Form.
Definition utils.h:965
Form< T, U > create_form_factory(const ufcx_form &ufcx_form, 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::map< std::shared_ptr< const mesh::Mesh< U > >, std::span< const std::int32_t > > &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:307
std::vector< std::pair< int, std::vector< std::int32_t > > > compute_integration_domains(IntegralType integral_type, const mesh::Topology &topology, std::span< const std::int32_t > entities, int dim, std::span< const int > values)
Given an integral type and mesh tag data, compute the entities that should be integrated over.
Definition utils.cpp:199
std::vector< DofMap > create_dofmaps(MPI_Comm comm, const std::vector< ElementDofLayout > &layouts, mesh::Topology &topology, std::function< void(std::span< std::int32_t >, std::uint32_t)> permute_inv, 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:135
DofMap create_dofmap(MPI_Comm comm, const ElementDofLayout &layout, mesh::Topology &topology, std::function< void(std::span< std::int32_t >, std::uint32_t)> permute_inv, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> reorder_fn)
Create a dof map on mesh.
Definition utils.cpp:98
std::vector< std::string > get_coefficient_names(const ufcx_form &ufcx_form)
Definition utils.cpp:184
FunctionSpace(U mesh, V element, W dofmap, X value_shape) -> FunctionSpace< typename std::remove_cvref< typename U::element_type >::type::geometry_type::value_type >
Type deduction.
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:125
std::vector< typename U::scalar_type > pack_constants(const U &u)
Pack constants of u into a single array ready for assembly.
Definition utils.h:1249
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_function_space=nullptr)
Create Expression from UFC.
Definition utils.h:1096
IntegralType
Type of integral.
Definition Form.h:34
@ interior_facet
Interior facet.
@ exterior_facet
Exterior facet.
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, 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:621
FunctionSpace< T > create_functionspace(std::shared_ptr< mesh::Mesh< T > > mesh, const basix::FiniteElement< T > &e, const std::vector< std::size_t > &value_shape={}, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> reorder_fn=nullptr)
Create a function space from a Basix element.
Definition utils.h:707
la::SparsityPattern create_sparsity_pattern(const Form< T, U > &a)
Create a sparsity pattern for a given form.
Definition utils.h:150
std::vector< std::size_t > compute_value_shape(std::shared_ptr< const dolfinx::fem::FiniteElement< T > > element, std::size_t tdim, std::size_t gdim)
Compute the physical value shape of an element for a mesh.
Definition FunctionSpace.h:439
ElementDofLayout create_element_dof_layout(const ufcx_dofmap &dofmap, const mesh::CellType cell_type, const std::vector< int > &parent_map={})
Create an ElementDofLayout from a ufcx_dofmap.
Definition utils.cpp:32
void pack_coefficients(const Form< T, U > &form, IntegralType integral_type, int id, std::span< T > c, int cstride)
Pack coefficients of a Form for a given integral type and domain id.
Definition utils.h:1019
basix::cell::type cell_type_to_basix_type(CellType celltype)
Convert a cell type to a Basix cell type.
Definition cell_types.cpp:212
std::vector< std::int32_t > exterior_facet_indices(const Topology &topology)
Compute the indices of all exterior facets that are owned by the caller.
Definition utils.cpp:58
CellType
Cell type identifier.
Definition cell_types.h:22