DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
utils.h
Go to the documentation of this file.
1// Copyright (C) 2019-2024 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 "EntityMap.h"
10#include "Mesh.h"
11#include "Topology.h"
12#include "graphbuild.h"
13#include <algorithm>
14#include <basix/mdspan.hpp>
15#include <concepts>
16#include <cstdint>
17#include <dolfinx/graph/AdjacencyList.h>
18#include <dolfinx/graph/ordering.h>
19#include <dolfinx/graph/partition.h>
20#include <functional>
21#include <mpi.h>
22#include <numeric>
23#include <optional>
24#include <span>
25#include <vector>
26
29
30namespace dolfinx::fem
31{
33}
34
35namespace dolfinx::mesh
36{
37enum class CellType : std::int8_t;
38
40enum class GhostMode : std::uint8_t
41{
42 none,
43 shared_facet
44};
45
46namespace impl
47{
53template <typename T>
54void reorder_list(std::span<T> list, std::span<const std::int32_t> nodemap)
55{
56 if (nodemap.empty())
57 return;
58
59 assert(list.size() % nodemap.size() == 0);
60 std::size_t degree = list.size() / nodemap.size();
61 const std::vector<T> orig(list.begin(), list.end());
62 for (std::size_t n = 0; n < nodemap.size(); ++n)
63 {
64 std::span links_old(orig.data() + n * degree, degree);
65 auto links_new = list.subspan(nodemap[n] * degree, degree);
66 std::ranges::copy(links_old, links_new.begin());
67 }
68}
69
83template <std::floating_point T>
84std::tuple<std::vector<std::int32_t>, std::vector<T>, std::vector<std::int32_t>>
86 std::span<const std::int32_t> facets)
87{
88 auto topology = mesh.topology();
89 assert(topology);
90 const int tdim = topology->dim();
91 if (dim == tdim)
92 {
93 throw std::runtime_error(
94 "Cannot use mesh::locate_entities_boundary (boundary) for cells.");
95 }
96
97 // Build set of vertices on boundary and set of boundary entities
98 mesh.topology_mutable()->create_connectivity(tdim - 1, 0);
99 mesh.topology_mutable()->create_connectivity(tdim - 1, dim);
100 std::vector<std::int32_t> vertices, entities;
101 {
102 auto f_to_v = topology->connectivity(tdim - 1, 0);
103 assert(f_to_v);
104 auto f_to_e = topology->connectivity(tdim - 1, dim);
105 assert(f_to_e);
106 for (auto f : facets)
107 {
108 auto v = f_to_v->links(f);
109 vertices.insert(vertices.end(), v.begin(), v.end());
110 auto e = f_to_e->links(f);
111 entities.insert(entities.end(), e.begin(), e.end());
112 }
113
114 // Build vector of boundary vertices
115 {
116 std::ranges::sort(vertices);
117 auto [unique_end, range_end] = std::ranges::unique(vertices);
118 vertices.erase(unique_end, range_end);
119 }
120
121 {
122 std::ranges::sort(entities);
123 auto [unique_end, range_end] = std::ranges::unique(entities);
124 entities.erase(unique_end, range_end);
125 }
126 }
127
128 // Get geometry data
129 auto x_dofmap = mesh.geometry().dofmap();
130 std::span<const T> x_nodes = mesh.geometry().x();
131
132 // Get all vertex 'node' indices
133 mesh.topology_mutable()->create_connectivity(0, tdim);
134 mesh.topology_mutable()->create_connectivity(tdim, 0);
135 auto v_to_c = topology->connectivity(0, tdim);
136 assert(v_to_c);
137 auto c_to_v = topology->connectivity(tdim, 0);
138 assert(c_to_v);
139 std::vector<T> x_vertices(3 * vertices.size(), -1.0);
140 std::vector<std::int32_t> vertex_to_pos(v_to_c->num_nodes(), -1);
141 for (std::size_t i = 0; i < vertices.size(); ++i)
142 {
143 const std::int32_t v = vertices[i];
144
145 // Get first cell and find position
146 const std::int32_t c = v_to_c->links(v).front();
147 auto cell_vertices = c_to_v->links(c);
148 auto it = std::find(cell_vertices.begin(), cell_vertices.end(), v);
149 assert(it != cell_vertices.end());
150 const std::size_t local_pos = std::distance(cell_vertices.begin(), it);
151
152 auto dofs = md::submdspan(x_dofmap, c, md::full_extent);
153 for (std::size_t j = 0; j < 3; ++j)
154 x_vertices[j * vertices.size() + i] = x_nodes[3 * dofs[local_pos] + j];
155 vertex_to_pos[v] = i;
156 }
157
158 return {std::move(entities), std::move(x_vertices), std::move(vertex_to_pos)};
159}
160
161} // namespace impl
162
176std::vector<std::int32_t> exterior_facet_indices(const Topology& topology,
177 int facet_type_idx);
178
190std::vector<std::int32_t> exterior_facet_indices(const Topology& topology);
191
209using CellPartitionFunction = std::function<graph::AdjacencyList<std::int32_t>(
210 MPI_Comm comm, int nparts, const std::vector<CellType>& cell_types,
211 const std::vector<std::span<const std::int64_t>>& cells)>;
212
217using CellReorderFunction = std::function<std::vector<std::int32_t>(
219
228inline auto
230 std::optional<std::int32_t> max_facet_to_cell_links
231 = 2)
232{
244 return [&, max_facet_to_cell_links](
245 const std::vector<CellType>& celltypes,
246 const std::vector<fem::ElementDofLayout>& doflayouts,
247 const std::vector<std::vector<int>>& ghost_owners,
248 std::vector<std::vector<std::int64_t>>& cells,
249 std::vector<std::vector<std::int64_t>>& cells_v,
250 std::vector<std::vector<std::int64_t>>& original_idx)
251 -> std::vector<std::int64_t>
252 {
253 // Build local dual graph for owned cells to (i) get list of vertices
254 // on the process boundary and (ii) apply re-ordering to cells for
255 // locality
256
257 spdlog::info("Build local dual graphs, re-order cells, and compute process "
258 "boundary vertices.");
259
260 std::vector<std::pair<std::vector<std::int64_t>, int>> facets;
261
262 // Build lists of cells (by cell type) that excludes ghosts
263 std::vector<std::span<const std::int64_t>> cells1_v_local;
264 for (std::size_t i = 0; i < celltypes.size(); ++i)
265 {
266 int num_cell_vertices = mesh::num_cell_vertices(celltypes[i]);
267 std::size_t num_owned_cells
268 = cells_v[i].size() / num_cell_vertices - ghost_owners[i].size();
269 cells1_v_local.emplace_back(cells_v[i].data(),
270 num_owned_cells * num_cell_vertices);
271
272 // Build local dual graph for cell type
273 auto [graph, unmatched_facets, max_v, _facet_attached_cells]
274 = build_local_dual_graph(std::vector{celltypes[i]},
275 std::vector{cells1_v_local.back()},
276 max_facet_to_cell_links);
277
278 // Store unmatched_facets for current cell type
279 facets.emplace_back(std::move(unmatched_facets), max_v);
280
281 // Compute re-ordering of graph
282 const std::vector<std::int32_t> remap = reorder_fn(graph);
283
284 // Update 'original' indices
285 const std::vector<std::int64_t>& orig_idx = original_idx[i];
286 std::vector<std::int64_t> _original_idx(orig_idx.size());
287 std::copy_n(orig_idx.rbegin(), ghost_owners[i].size(),
288 _original_idx.rbegin());
289 {
290 for (std::size_t j = 0; j < remap.size(); ++j)
291 _original_idx[remap[j]] = orig_idx[j];
292 }
293 original_idx[i] = _original_idx;
294
295 // Reorder cells
296 impl::reorder_list(
297 std::span(cells_v[i].data(), remap.size() * num_cell_vertices),
298 remap);
299 impl::reorder_list(
300 std::span(cells[i].data(), remap.size() * doflayouts[i].num_dofs()),
301 remap);
302 }
303
304 if (facets.size() == 1) // Optimisation for single cell type
305 {
306 std::vector<std::int64_t>& vertices = facets.front().first;
307
308 // Remove duplicated vertex indices
309 std::ranges::sort(vertices);
310 auto [unique_end, range_end] = std::ranges::unique(vertices);
311 vertices.erase(unique_end, range_end);
312
313 // Remove -1 if it appears as first entity. This can happen in
314 // mixed topology meshes where '-1' is used to pad facet data when
315 // cells facets have differing numbers of vertices.
316 if (!vertices.empty() and vertices.front() == -1)
317 vertices.erase(vertices.begin());
318
319 return vertices;
320 }
321 else
322 {
323 // Pack 'unmatched' facets for all cell types into single array
324 // (facets0)
325 std::vector<std::int64_t> facets0;
326 facets0.reserve(std::accumulate(facets.begin(), facets.end(),
327 std::size_t(0), [](std::size_t x, auto& y)
328 { return x + y.first.size(); }));
329 int max_v = std::ranges::max_element(facets, [](auto& a, auto& b)
330 { return a.second < b.second; })
331 ->second;
332 for (const auto& [v_data, num_v] : facets)
333 {
334 for (auto it = v_data.begin(); it != v_data.end(); it += num_v)
335 {
336 facets0.insert(facets0.end(), it, std::next(it, num_v));
337 facets0.insert(facets0.end(), max_v - num_v, -1);
338 }
339 }
340
341 // Compute row permutation
342 const std::vector<std::int32_t> perm = dolfinx::sort_by_perm(
343 std::span<const std::int64_t>(facets0), max_v);
344
345 // For facets in facets0 that appear only once, store the facet
346 // vertices
347 std::vector<std::int64_t> vertices;
348 // TODO: allocate memory for vertices
349 auto it = perm.begin();
350 while (it != perm.end())
351 {
352 // Find iterator to next facet different from f and trim any -1
353 // padding
354 std::span _f(facets0.data() + (*it) * max_v, max_v);
355 auto end = std::find_if(_f.rbegin(), _f.rend(),
356 [](auto a) { return a >= 0; });
357 auto f = _f.first(std::distance(end, _f.rend()));
358
359 auto it1 = std::find_if_not(
360 it, perm.end(),
361 [f, max_v, it0 = facets0.begin()](auto p) -> bool
362 {
363 return std::equal(f.begin(), f.end(), std::next(it0, p * max_v));
364 });
365
366 // If no repeated facet found, insert f vertices
367 if (std::distance(it, it1) == 1)
368 vertices.insert(vertices.end(), f.begin(), f.end());
369 else if (std::distance(it, it1) > 2)
370 throw std::runtime_error("More than two matching facets found.");
371
372 // Advance iterator
373 it = it1;
374 }
375
376 // Remove duplicate indices
377 std::ranges::sort(vertices);
378 auto [unique_end, range_end] = std::ranges::unique(vertices);
379 vertices.erase(unique_end, range_end);
380
381 return vertices;
382 }
383 };
384}
385
395std::vector<std::int64_t> extract_topology(CellType cell_type,
396 const fem::ElementDofLayout& layout,
397 std::span<const std::int64_t> cells);
398
407template <std::floating_point T>
408std::vector<T> h(const Mesh<T>& mesh, std::span<const std::int32_t> entities,
409 int dim)
410{
411 if (entities.empty())
412 return std::vector<T>();
413 if (dim == 0)
414 return std::vector<T>(entities.size(), 0);
415
416 // Get the geometry dofs for the vertices of each entity
417 const auto [vertex_xdofs, xdof_shape]
418 = entities_to_geometry(mesh, dim, entities, false);
419
420 // Get the geometry coordinate
421 std::span<const T> x = mesh.geometry().x();
422
423 // Function to compute the length of (p0 - p1)
424 auto delta_norm = [](auto&& p0, auto&& p1)
425 {
426 T norm = 0;
427 for (std::size_t i = 0; i < 3; ++i)
428 norm += (p0[i] - p1[i]) * (p0[i] - p1[i]);
429 return std::sqrt(norm);
430 };
431
432 // Compute greatest distance between any to vertices
433 assert(dim > 0);
434 std::vector<T> h(entities.size(), 0);
435 for (std::size_t e = 0; e < entities.size(); ++e)
436 {
437 // Get geometry 'dof' for each vertex of entity e
438 std::span<const std::int32_t> e_vertices(
439 vertex_xdofs.data() + e * xdof_shape[1], xdof_shape[1]);
440
441 // Compute maximum distance between any two vertices
442 for (std::size_t i = 0; i < e_vertices.size(); ++i)
443 {
444 std::span<const T, 3> p0(x.data() + 3 * e_vertices[i], 3);
445 for (std::size_t j = i + 1; j < e_vertices.size(); ++j)
446 {
447 std::span<const T, 3> p1(x.data() + 3 * e_vertices[j], 3);
448 h[e] = std::max(h[e], delta_norm(p0, p1));
449 }
450 }
451 }
452
453 return h;
454}
455
459template <std::floating_point T>
460std::vector<T> cell_normals(const Mesh<T>& mesh, int dim,
461 std::span<const std::int32_t> entities)
462{
463 if (entities.empty())
464 return std::vector<T>();
465
466 auto topology = mesh.topology();
467 assert(topology);
468 if (topology->cell_type() == CellType::prism and dim == 2)
469 {
470 throw std::runtime_error(
471 "Cell normal computation for prism cells not yet supported.");
472 }
473
474 const int gdim = mesh.geometry().dim();
475 const CellType type = cell_entity_type(topology->cell_type(), dim, 0);
476
477 // Find geometry nodes for topology entities
478 std::span<const T> x = mesh.geometry().x();
479 const auto [geometry_entities, eshape]
480 = entities_to_geometry(mesh, dim, entities, false);
481
482 std::vector<T> n(entities.size() * 3);
483 switch (type)
484 {
485 case CellType::interval:
486 {
487 if (gdim > 2)
488 throw std::invalid_argument("Interval cell normal undefined in 3D.");
489 for (std::size_t i = 0; i < entities.size(); ++i)
490 {
491 // Get the two vertices as points
492 std::array vertices{geometry_entities[i * eshape[1]],
493 geometry_entities[i * eshape[1] + 1]};
494 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
495 std::span<const T, 3>(x.data() + 3 * vertices[1], 3)};
496
497 // Define normal by rotating tangent counter-clockwise
498 std::array<T, 3> t;
499 std::ranges::transform(p[1], p[0], t.begin(),
500 [](auto x, auto y) { return x - y; });
501
502 T norm = std::sqrt(t[0] * t[0] + t[1] * t[1]);
503 std::span<T, 3> ni(n.data() + 3 * i, 3);
504 ni[0] = -t[1] / norm;
505 ni[1] = t[0] / norm;
506 ni[2] = 0.0;
507 }
508 return n;
509 }
510 case CellType::triangle:
511 {
512 for (std::size_t i = 0; i < entities.size(); ++i)
513 {
514 // Get the three vertices as points
515 std::array vertices = {geometry_entities[i * eshape[1] + 0],
516 geometry_entities[i * eshape[1] + 1],
517 geometry_entities[i * eshape[1] + 2]};
518 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
519 std::span<const T, 3>(x.data() + 3 * vertices[1], 3),
520 std::span<const T, 3>(x.data() + 3 * vertices[2], 3)};
521
522 // Compute (p1 - p0) and (p2 - p0)
523 std::array<T, 3> dp1, dp2;
524 std::ranges::transform(p[1], p[0], dp1.begin(),
525 [](auto x, auto y) { return x - y; });
526 std::ranges::transform(p[2], p[0], dp2.begin(),
527 [](auto x, auto y) { return x - y; });
528
529 // Define cell normal via cross product of first two edges
530 std::array<T, 3> ni = math::cross(dp1, dp2);
531 T norm = std::sqrt(ni[0] * ni[0] + ni[1] * ni[1] + ni[2] * ni[2]);
532 std::ranges::transform(ni, std::next(n.begin(), 3 * i),
533 [norm](auto x) { return x / norm; });
534 }
535
536 return n;
537 }
538 case CellType::quadrilateral:
539 {
540 // TODO: check
541 for (std::size_t i = 0; i < entities.size(); ++i)
542 {
543 // Get the three vertices as points
544 std::array vertices = {geometry_entities[i * eshape[1] + 0],
545 geometry_entities[i * eshape[1] + 1],
546 geometry_entities[i * eshape[1] + 2]};
547 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
548 std::span<const T, 3>(x.data() + 3 * vertices[1], 3),
549 std::span<const T, 3>(x.data() + 3 * vertices[2], 3)};
550
551 // Compute (p1 - p0) and (p2 - p0)
552 std::array<T, 3> dp1, dp2;
553 std::ranges::transform(p[1], p[0], dp1.begin(),
554 [](auto x, auto y) { return x - y; });
555 std::ranges::transform(p[2], p[0], dp2.begin(),
556 [](auto x, auto y) { return x - y; });
557
558 // Define cell normal via cross product of first two edges
559 std::array<T, 3> ni = math::cross(dp1, dp2);
560 T norm = std::sqrt(ni[0] * ni[0] + ni[1] * ni[1] + ni[2] * ni[2]);
561 std::ranges::transform(ni, std::next(n.begin(), 3 * i),
562 [norm](auto x) { return x / norm; });
563 }
564
565 return n;
566 }
567 default:
568 throw std::invalid_argument(
569 "cell_normal not supported for this cell type.");
570 }
571}
572
576template <std::floating_point T>
577std::vector<T> compute_midpoints(const Mesh<T>& mesh, int dim,
578 std::span<const std::int32_t> entities)
579{
580 if (entities.empty())
581 return std::vector<T>();
582
583 std::span<const T> x = mesh.geometry().x();
584
585 // Build map from entity -> geometry dof
586 const auto [e_to_g, eshape]
587 = entities_to_geometry(mesh, dim, entities, false);
588
589 std::vector<T> x_mid(entities.size() * 3, 0);
590 for (std::size_t e = 0; e < entities.size(); ++e)
591 {
592 std::span<T, 3> p(x_mid.data() + 3 * e, 3);
593 std::span<const std::int32_t> rows(e_to_g.data() + e * eshape[1],
594 eshape[1]);
595 for (auto row : rows)
596 {
597 std::span<const T, 3> xg(x.data() + 3 * row, 3);
598 std::ranges::transform(p, xg, p.begin(),
599 [size = rows.size()](auto x, auto y)
600 { return x + y / size; });
601 }
602 }
603
604 return x_mid;
605}
606
607namespace impl
608{
613template <std::floating_point T>
614std::pair<std::vector<T>, std::array<std::size_t, 2>>
616{
617 auto topology = mesh.topology();
618 assert(topology);
619 const int tdim = topology->dim();
620
621 // Create entities and connectivities
622
623 // Get all vertex 'node' indices
624 const std::int32_t num_vertices = topology->index_map(0)->size_local()
625 + topology->index_map(0)->num_ghosts();
626
627 std::vector<std::int32_t> vertex_to_node(num_vertices);
628 for (int cell_type_idx = 0,
629 num_cell_types = topology->entity_types(tdim).size();
630 cell_type_idx < num_cell_types; ++cell_type_idx)
631 {
632 auto x_dofmap = mesh.geometry().dofmap(cell_type_idx);
633 auto c_to_v = topology->connectivity({tdim, cell_type_idx}, {0, 0});
634 assert(c_to_v);
635 for (int c = 0; c < c_to_v->num_nodes(); ++c)
636 {
637 auto x_dofs = md::submdspan(x_dofmap, c, md::full_extent);
638 auto vertices = c_to_v->links(c);
639 for (std::size_t i = 0; i < vertices.size(); ++i)
640 vertex_to_node[vertices[i]] = x_dofs[i];
641 }
642 }
643
644 // Pack coordinates of vertices
645 std::span<const T> x_nodes = mesh.geometry().x();
646 std::vector<T> x_vertices(3 * vertex_to_node.size(), 0.0);
647 for (std::size_t i = 0; i < vertex_to_node.size(); ++i)
648 {
649 std::int32_t pos = 3 * vertex_to_node[i];
650 for (std::size_t j = 0; j < 3; ++j)
651 x_vertices[j * vertex_to_node.size() + i] = x_nodes[pos + j];
652 }
653
654 return {std::move(x_vertices), {3, vertex_to_node.size()}};
655}
656
657} // namespace impl
658
660template <typename Fn, typename T>
661concept MarkerFn = std::is_invocable_r<
662 std::vector<std::int8_t>, Fn,
663 md::mdspan<const T,
664 md::extents<std::size_t, 3, md::dynamic_extent>>>::value;
665
681template <std::floating_point T, MarkerFn<T> U>
682std::vector<std::int32_t> locate_entities(const Mesh<T>& mesh, int dim,
683 U marker, int entity_type_idx)
684{
685
686 using cmdspan3x_t
687 = md::mdspan<const T, md::extents<std::size_t, 3, md::dynamic_extent>>;
688
689 // Run marker function on vertex coordinates
690 const auto [xdata, xshape] = impl::compute_vertex_coords(mesh);
691
692 cmdspan3x_t x(xdata.data(), xshape);
693 const std::vector<std::int8_t> marked = marker(x);
694 if (marked.size() != x.extent(1))
695 throw std::runtime_error("Length of array of markers is wrong.");
696
697 auto topology = mesh.topology();
698 assert(topology);
699 const int tdim = topology->dim();
700
701 mesh.topology_mutable()->create_entities(dim);
702 if (dim < tdim)
703 mesh.topology_mutable()->create_connectivity(dim, 0);
704
705 // Iterate over entities of dimension 'dim' to build vector of marked
706 // entities
707 auto e_to_v = topology->connectivity({dim, entity_type_idx}, {0, 0});
708 assert(e_to_v);
709 std::vector<std::int32_t> entities;
710 for (int e = 0; e < e_to_v->num_nodes(); ++e)
711 {
712 // Iterate over entity vertices
713 bool all_vertices_marked = true;
714 for (std::int32_t v : e_to_v->links(e))
715 {
716 if (!marked[v])
717 {
718 all_vertices_marked = false;
719 break;
720 }
721 }
722
723 if (all_vertices_marked)
724 entities.push_back(e);
725 }
726
727 return entities;
728}
729
743template <std::floating_point T, MarkerFn<T> U>
744std::vector<std::int32_t> locate_entities(const Mesh<T>& mesh, int dim,
745 U marker)
746{
747 const int num_entity_types = mesh.topology()->entity_types(dim).size();
748 if (num_entity_types > 1)
749 {
750 throw std::runtime_error(
751 "Multiple entity types of this dimension. Specify entity type index");
752 }
753 return locate_entities(mesh, dim, marker, 0);
754}
755
779template <std::floating_point T, MarkerFn<T> U>
780std::vector<std::int32_t> locate_entities_boundary(const Mesh<T>& mesh, int dim,
781 U marker)
782{
783 // TODO Rewrite this function, it should be possible to simplify considerably
784 auto topology = mesh.topology();
785 assert(topology);
786 int tdim = topology->dim();
787 if (dim == tdim)
788 {
789 throw std::runtime_error(
790 "Cannot use mesh::locate_entities_boundary (boundary) for cells.");
791 }
792
793 // Compute list of boundary facets
794 mesh.topology_mutable()->create_entities(tdim - 1);
795 mesh.topology_mutable()->create_connectivity(tdim - 1, tdim);
796 std::vector<std::int32_t> boundary_facets = exterior_facet_indices(*topology);
797
798 using cmdspan3x_t
799 = md::mdspan<const T, md::extents<std::size_t, 3, md::dynamic_extent>>;
800
801 // Run marker function on the vertex coordinates
802 auto [facet_entities, xdata, vertex_to_pos]
803 = impl::compute_vertex_coords_boundary(mesh, dim, boundary_facets);
804 cmdspan3x_t x(xdata.data(), 3, xdata.size() / 3);
805 std::vector<std::int8_t> marked = marker(x);
806 if (marked.size() != x.extent(1))
807 throw std::runtime_error("Length of array of markers is wrong.");
808
809 // Loop over entities and check vertex markers
810 mesh.topology_mutable()->create_entities(dim);
811 auto e_to_v = topology->connectivity(dim, 0);
812 assert(e_to_v);
813 std::vector<std::int32_t> entities;
814 for (auto e : facet_entities)
815 {
816 // Iterate over entity vertices
817 bool all_vertices_marked = true;
818 for (auto v : e_to_v->links(e))
819 {
820 const std::int32_t pos = vertex_to_pos[v];
821 if (!marked[pos])
822 {
823 all_vertices_marked = false;
824 break;
825 }
826 }
827
828 // Mark facet with all vertices marked
829 if (all_vertices_marked)
830 entities.push_back(e);
831 }
832
833 return entities;
834}
835
854template <std::floating_point T>
855std::pair<std::vector<std::int32_t>, std::array<std::size_t, 2>>
857 std::span<const std::int32_t> entities,
858 bool permute = false)
859{
860 auto topology = mesh.topology();
861 assert(topology);
862 CellType cell_type = topology->cell_type();
863 if (cell_type == CellType::prism and dim == 2)
864 {
865 throw std::runtime_error(
866 "mesh::entities_to_geometry for prism cells not yet supported.");
867 }
868
869 const int tdim = topology->dim();
870 const Geometry<T>& geometry = mesh.geometry();
871 auto xdofs = geometry.dofmap();
872
873 // Get the DOF layout and the number of DOFs per entity
874 const fem::CoordinateElement<T>& coord_ele = geometry.cmap();
875 const fem::ElementDofLayout layout = coord_ele.create_dof_layout();
876 const std::size_t num_entity_dofs = layout.num_entity_closure_dofs(dim);
877 std::vector<std::int32_t> entity_xdofs;
878 entity_xdofs.reserve(entities.size() * num_entity_dofs);
879 std::array<std::size_t, 2> eshape{entities.size(), num_entity_dofs};
880
881 // Get the element's closure DOFs
882 const std::vector<std::vector<std::vector<int>>>& closure_dofs_all
883 = layout.entity_closure_dofs_all();
884
885 // Special case when dim == tdim (cells)
886 if (dim == tdim)
887 {
888 for (std::int32_t c : entities)
889 {
890 // Extract degrees of freedom
891 auto x_c = md::submdspan(xdofs, c, md::full_extent);
892 for (std::int32_t entity_dof : closure_dofs_all[tdim][0])
893 entity_xdofs.push_back(x_c[entity_dof]);
894 }
895
896 return {std::move(entity_xdofs), eshape};
897 }
898
899 assert(dim != tdim);
900
901 auto e_to_c = topology->connectivity(dim, tdim);
902 if (!e_to_c)
903 {
904 throw std::runtime_error(
905 "Entity-to-cell connectivity has not been computed. Missing dims "
906 + std::to_string(dim) + "->" + std::to_string(tdim));
907 }
908
909 auto c_to_e = topology->connectivity(tdim, dim);
910 if (!c_to_e)
911 {
912 throw std::runtime_error(
913 "Cell-to-entity connectivity has not been computed. Missing dims "
914 + std::to_string(tdim) + "->" + std::to_string(dim));
915 }
916
917 // Get the cell info, which is needed to permute the closure dofs
918 std::span<const std::uint32_t> cell_info;
919 if (permute)
920 cell_info = std::span(mesh.topology()->get_cell_permutation_info());
921
922 for (std::int32_t e : entities)
923 {
924 // Get a cell connected to the entity
925 assert(!e_to_c->links(e).empty());
926 std::int32_t c = e_to_c->links(e).front();
927
928 // Get the local index of the entity
929 std::span<const std::int32_t> cell_entities = c_to_e->links(c);
930 auto it = std::find(cell_entities.begin(), cell_entities.end(), e);
931 assert(it != cell_entities.end());
932 std::size_t local_entity = std::distance(cell_entities.begin(), it);
933
934 // Cell sub-entities must be permuted so that their local
935 // orientation agrees with their global orientation
936 std::vector<std::int32_t> closure_dofs(closure_dofs_all[dim][local_entity]);
937 if (permute)
938 {
939 mesh::CellType entity_type
940 = mesh::cell_entity_type(cell_type, dim, local_entity);
941 coord_ele.permute_subentity_closure(closure_dofs, cell_info[c],
942 entity_type, local_entity);
943 }
944
945 // Extract degrees of freedom
946 auto x_c = md::submdspan(xdofs, c, md::full_extent);
947 for (std::int32_t entity_dof : closure_dofs)
948 entity_xdofs.push_back(x_c[entity_dof]);
949 }
950
951 return {std::move(entity_xdofs), eshape};
952}
953
965 mesh::GhostMode ghost_mode = mesh::GhostMode::none,
967 std::optional<std::int32_t> max_facet_to_cell_links = 2);
968
976std::vector<std::int32_t>
977compute_incident_entities(const Topology& topology,
978 std::span<const std::int32_t> entities, int d0,
979 int d1);
980
1024template <typename U>
1026 MPI_Comm comm, MPI_Comm commt,
1027 std::vector<std::span<const std::int64_t>> cells,
1028 const std::vector<fem::CoordinateElement<
1029 typename std::remove_reference_t<typename U::value_type>>>& elements,
1030 MPI_Comm commg, const U& x, std::array<std::size_t, 2> xshape,
1031 const CellPartitionFunction& partitioner,
1032 std::optional<std::int32_t> max_facet_to_cell_links,
1033 const CellReorderFunction& reorder_fn = graph::reorder_gps)
1034{
1035 assert(cells.size() == elements.size());
1036 std::vector<CellType> celltypes;
1037 std::ranges::transform(elements, std::back_inserter(celltypes),
1038 [](auto& e) { return e.cell_shape(); });
1039 std::vector<fem::ElementDofLayout> doflayouts;
1040 std::ranges::transform(elements, std::back_inserter(doflayouts),
1041 [](auto& e) { return e.create_dof_layout(); });
1042
1043 // Note: `extract_topology` extracts topology data, i.e. just the
1044 // vertices. For P1 geometry this should just be the identity
1045 // operator. For other elements the filtered lists may have 'gaps',
1046 // i.e. the indices might not be contiguous.
1047 //
1048 // `extract_topology` could be skipped for 'P1 geometry' elements
1049
1050 std::int32_t num_cell_types = cells.size();
1051
1052 // -- Partition topology across ranks of comm
1053 std::vector<std::vector<std::int64_t>> cells1(num_cell_types);
1054 std::vector<std::vector<std::int64_t>> original_idx1(num_cell_types);
1055 std::vector<std::vector<int>> ghost_owners(num_cell_types);
1056 if (partitioner)
1057 {
1058 spdlog::info("Using partitioner with cell data ({} cell types)",
1059 num_cell_types);
1061 if (commt != MPI_COMM_NULL)
1062 {
1063 int size = dolfinx::MPI::size(comm);
1064 std::vector<std::vector<std::int64_t>> t(num_cell_types);
1065 std::vector<std::span<const std::int64_t>> tspan(num_cell_types);
1066 for (std::int32_t i = 0; i < num_cell_types; ++i)
1067 {
1068 t[i] = extract_topology(celltypes[i], doflayouts[i], cells[i]);
1069 tspan[i] = std::span(t[i]);
1070 }
1071 dest = partitioner(commt, size, celltypes, tspan);
1072 }
1073
1074 std::int32_t cell_offset = 0;
1075 for (std::int32_t i = 0; i < num_cell_types; ++i)
1076 {
1077 std::size_t num_cell_nodes = doflayouts[i].num_dofs();
1078 assert(cells[i].size() % num_cell_nodes == 0);
1079 std::size_t num_cells = cells[i].size() / num_cell_nodes;
1080
1081 // Extract destination AdjacencyList for this cell type
1082 std::vector<std::int32_t> offsets_i(
1083 std::next(dest.offsets().begin(), cell_offset),
1084 std::next(dest.offsets().begin(), cell_offset + num_cells + 1));
1085 std::vector<std::int32_t> data_i(
1086 std::next(dest.array().begin(), offsets_i.front()),
1087 std::next(dest.array().begin(), offsets_i.back()));
1088 std::int32_t offset_0 = offsets_i.front();
1089 std::ranges::for_each(offsets_i,
1090 [&offset_0](std::int32_t& j) { j -= offset_0; });
1091 graph::AdjacencyList<std::int32_t> dest_i(data_i, offsets_i);
1092 cell_offset += num_cells;
1093
1094 // Distribute cells (topology, includes higher-order 'nodes') to
1095 // destination rank
1096 std::vector<int> src_ranks;
1097 std::tie(cells1[i], src_ranks, original_idx1[i], ghost_owners[i])
1098 = graph::build::distribute(comm, cells[i],
1099 {num_cells, num_cell_nodes}, dest_i);
1100 spdlog::debug("Got {} cells from distribution", cells1[i].size());
1101 }
1102 }
1103 else
1104 {
1105 // No partitioning, construct a global index
1106 std::int64_t num_owned = 0;
1107 for (std::int32_t i = 0; i < num_cell_types; ++i)
1108 {
1109 cells1[i] = std::vector<std::int64_t>(cells[i].begin(), cells[i].end());
1110 std::int32_t num_cell_nodes = doflayouts[i].num_dofs();
1111 assert(cells1[i].size() % num_cell_nodes == 0);
1112 original_idx1[i].resize(cells1[i].size() / num_cell_nodes);
1113 num_owned += original_idx1[i].size();
1114 }
1115
1116 // Add on global offset
1117 std::int64_t global_offset = 0;
1118 MPI_Exscan(&num_owned, &global_offset, 1, MPI_INT64_T, MPI_SUM, comm);
1119 for (std::int32_t i = 0; i < num_cell_types; ++i)
1120 {
1121 std::iota(original_idx1[i].begin(), original_idx1[i].end(),
1122 global_offset);
1123 global_offset += original_idx1[i].size();
1124 }
1125 }
1126
1127 // Extract cell 'topology', i.e. extract the vertices for each cell
1128 // and discard any 'higher-order' nodes
1129 std::vector<std::vector<std::int64_t>> cells1_v(num_cell_types);
1130 for (std::int32_t i = 0; i < num_cell_types; ++i)
1131 {
1132 cells1_v[i] = extract_topology(celltypes[i], doflayouts[i], cells1[i]);
1133 spdlog::info("Extract basic topology: {}->{}", cells1[i].size(),
1134 cells1_v[i].size());
1135 }
1136
1137 auto boundary_v_fn
1138 = create_boundary_vertices_fn(reorder_fn, max_facet_to_cell_links);
1139 const std::vector<std::int64_t> boundary_v = boundary_v_fn(
1140 celltypes, doflayouts, ghost_owners, cells1, cells1_v, original_idx1);
1141
1142 spdlog::debug("Got {} boundary vertices", boundary_v.size());
1143
1144 // Create Topology
1145 std::vector<std::span<const std::int64_t>> cells1_v_span;
1146 std::ranges::transform(cells1_v, std::back_inserter(cells1_v_span),
1147 [](auto& c) { return std::span(c); });
1148 std::vector<std::span<const std::int64_t>> original_idx1_span;
1149 std::ranges::transform(original_idx1, std::back_inserter(original_idx1_span),
1150 [](auto& c) { return std::span(c); });
1151 std::vector<std::span<const int>> ghost_owners_span;
1152 std::ranges::transform(ghost_owners, std::back_inserter(ghost_owners_span),
1153 [](auto& c) { return std::span(c); });
1154 Topology topology
1155 = create_topology(comm, celltypes, cells1_v_span, original_idx1_span,
1156 ghost_owners_span, boundary_v);
1157
1158 // Create connectivities required higher-order geometries for creating
1159 // a Geometry object
1160 for (int i = 0; i < num_cell_types; ++i)
1161 {
1162 for (int e = 1; e < topology.dim(); ++e)
1163 {
1164 if (doflayouts[i].num_entity_dofs(e) > 0)
1165 topology.create_entities(e);
1166 }
1167
1168 if (elements[i].needs_dof_permutations())
1169 topology.create_entity_permutations();
1170 }
1171
1172 // Build list of unique (global) node indices from cells1 and
1173 // distribute coordinate data
1174 std::vector<std::int64_t> nodes1, nodes2;
1175 for (std::vector<std::int64_t>& c : cells1)
1176 nodes1.insert(nodes1.end(), c.begin(), c.end());
1177 for (std::vector<std::int64_t>& c : cells1)
1178 nodes2.insert(nodes2.end(), c.begin(), c.end());
1179
1180 dolfinx::radix_sort(nodes1);
1181 auto [unique_end, range_end] = std::ranges::unique(nodes1);
1182 nodes1.erase(unique_end, range_end);
1183
1184 std::vector coords
1185 = dolfinx::MPI::distribute_data(comm, nodes1, commg, x, xshape[1]);
1186
1187 // Create geometry object
1189 = create_geometry(topology, elements, nodes1, nodes2, coords, xshape[1]);
1190
1191 return Mesh(comm, std::make_shared<Topology>(std::move(topology)),
1192 std::move(geometry));
1193}
1194
1232template <typename U>
1234 MPI_Comm comm, MPI_Comm commt, std::span<const std::int64_t> cells,
1236 typename std::remove_reference_t<typename U::value_type>>& element,
1237 MPI_Comm commg, const U& x, std::array<std::size_t, 2> xshape,
1238 const CellPartitionFunction& partitioner,
1239 std::optional<std::int32_t> max_facet_to_cell_links = 2,
1240 const CellReorderFunction& reorder_fn = graph::reorder_gps)
1241{
1242 return create_mesh(comm, commt, std::vector{cells}, std::vector{element},
1243 commg, x, xshape, partitioner, max_facet_to_cell_links,
1244 reorder_fn);
1245}
1246
1267template <typename U>
1268Mesh<typename std::remove_reference_t<typename U::value_type>>
1269create_mesh(MPI_Comm comm, std::span<const std::int64_t> cells,
1271 std::remove_reference_t<typename U::value_type>>& elements,
1272 const U& x, std::array<std::size_t, 2> xshape, GhostMode ghost_mode,
1273 std::optional<std::int32_t> max_facet_to_cell_links = 2)
1274{
1275 if (dolfinx::MPI::size(comm) == 1)
1276 {
1277 return create_mesh(comm, comm, std::vector{cells}, std::vector{elements},
1278 comm, x, xshape, nullptr, max_facet_to_cell_links);
1279 }
1280 else
1281 {
1282 return create_mesh(comm, comm, std::vector{cells}, std::vector{elements},
1283 comm, x, xshape, create_cell_partitioner(ghost_mode),
1284 max_facet_to_cell_links);
1285 }
1286}
1287
1301template <std::floating_point T>
1302std::pair<Geometry<T>, std::vector<int32_t>>
1304 std::span<const std::int32_t> subentity_to_entity)
1305{
1306 const Geometry<T>& geometry = mesh.geometry();
1307
1308 // Get the geometry dofs in the sub-geometry based on the entities in
1309 // sub-geometry
1310 const fem::ElementDofLayout layout = geometry.cmap().create_dof_layout();
1311
1312 const std::vector<std::int32_t> x_indices
1313 = entities_to_geometry(mesh, dim, subentity_to_entity, true).first;
1314
1315 std::vector<std::int32_t> sub_x_dofs = x_indices;
1316 std::ranges::sort(sub_x_dofs);
1317 auto [unique_end, range_end] = std::ranges::unique(sub_x_dofs);
1318 sub_x_dofs.erase(unique_end, range_end);
1319
1320 // Get the sub-geometry dofs owned by this process
1321 auto x_index_map = geometry.index_map();
1322 assert(x_index_map);
1323
1324 std::shared_ptr<common::IndexMap> sub_x_dof_index_map;
1325 std::vector<std::int32_t> subx_to_x_dofmap;
1326 {
1327 auto [map, new_to_old] = common::create_sub_index_map(
1328 *x_index_map, sub_x_dofs, common::IndexMapOrder::any, true);
1329 sub_x_dof_index_map = std::make_shared<common::IndexMap>(std::move(map));
1330 subx_to_x_dofmap = std::move(new_to_old);
1331 }
1332
1333 // Create sub-geometry coordinates
1334 std::span<const T> x = geometry.x();
1335 std::int32_t sub_num_x_dofs = subx_to_x_dofmap.size();
1336 std::vector<T> sub_x(3 * sub_num_x_dofs);
1337 for (std::int32_t i = 0; i < sub_num_x_dofs; ++i)
1338 {
1339 std::copy_n(std::next(x.begin(), 3 * subx_to_x_dofmap[i]), 3,
1340 std::next(sub_x.begin(), 3 * i));
1341 }
1342
1343 // Create geometry to sub-geometry map
1344 std::vector<std::int32_t> x_to_subx_dof_map(
1345 x_index_map->size_local() + x_index_map->num_ghosts(), -1);
1346 for (std::size_t i = 0; i < subx_to_x_dofmap.size(); ++i)
1347 x_to_subx_dof_map[subx_to_x_dofmap[i]] = i;
1348
1349 // Create sub-geometry dofmap
1350 std::vector<std::int32_t> sub_x_dofmap;
1351 sub_x_dofmap.reserve(x_indices.size());
1352 std::ranges::transform(x_indices, std::back_inserter(sub_x_dofmap),
1353 [&x_to_subx_dof_map](auto x_dof)
1354 {
1355 assert(x_to_subx_dof_map[x_dof] != -1);
1356 return x_to_subx_dof_map[x_dof];
1357 });
1358
1359 // Sub-geometry coordinate element
1360 CellType sub_xcell = cell_entity_type(geometry.cmap().cell_shape(), dim, 0);
1361
1362 // Special handling of point meshes, as they only support constant
1363 // basis functions
1364 int degree = (sub_xcell == CellType::point) ? 0 : geometry.cmap().degree();
1365 fem::CoordinateElement<T> sub_cmap(sub_xcell, degree,
1366 geometry.cmap().variant());
1367
1368 // Sub-geometry input_global_indices
1369 const std::vector<std::int64_t>& igi = geometry.input_global_indices();
1370 std::vector<std::int64_t> sub_igi;
1371 sub_igi.reserve(subx_to_x_dofmap.size());
1372 std::ranges::transform(subx_to_x_dofmap, std::back_inserter(sub_igi),
1373 [&igi](auto sub_x_dof) { return igi[sub_x_dof]; });
1374
1375 // Create geometry
1376 return {Geometry(
1377 sub_x_dof_index_map,
1378 std::vector<std::vector<std::int32_t>>{std::move(sub_x_dofmap)},
1379 {sub_cmap}, std::move(sub_x), geometry.dim(), std::move(sub_igi)),
1380 std::move(subx_to_x_dofmap)};
1381}
1382
1392template <std::floating_point T>
1393std::tuple<Mesh<T>, EntityMap, EntityMap, std::vector<std::int32_t>>
1395 std::span<const std::int32_t> entities)
1396{
1397 // Create sub-topology
1398 mesh.topology_mutable()->create_connectivity(dim, 0);
1399 auto [topology, subentity_to_entity, subvertex_to_vertex]
1400 = mesh::create_subtopology(*mesh.topology(), dim, entities);
1401
1402 // Create sub-geometry
1403 const int tdim = mesh.topology()->dim();
1404 mesh.topology_mutable()->create_entities(dim);
1405 mesh.topology_mutable()->create_connectivity(dim, tdim);
1406 mesh.topology_mutable()->create_connectivity(tdim, dim);
1407 mesh.topology_mutable()->create_entity_permutations();
1408 auto [geometry, subx_to_x_dofmap]
1409 = mesh::create_subgeometry(mesh, dim, subentity_to_entity);
1410
1411 Mesh<T> submesh
1412 = Mesh(mesh.comm(), std::make_shared<Topology>(std::move(topology)),
1413 std::move(geometry));
1414 EntityMap entity_map(mesh.topology(), submesh.topology(), dim,
1415 subentity_to_entity);
1416 EntityMap vertex_map(mesh.topology(), submesh.topology(), 0,
1417 subvertex_to_vertex);
1418 return {std::move(submesh), std::move(entity_map), std::move(vertex_map),
1419 std::move(subx_to_x_dofmap)};
1420}
1421
1422} // namespace dolfinx::mesh
Definition CoordinateElement.h:38
ElementDofLayout create_dof_layout() const
Compute and return the dof layout.
Definition CoordinateElement.cpp:75
void permute_subentity_closure(std::span< std::int32_t > d, std::uint32_t cell_info, mesh::CellType entity_type, int entity_index) const
Given the closure DOFs of a cell sub-entity in reference ordering, this function computes the permut...
Definition CoordinateElement.cpp:64
Definition ElementDofLayout.h:31
const std::vector< std::vector< std::vector< int > > > & entity_closure_dofs_all() const
Definition ElementDofLayout.cpp:92
int num_entity_closure_dofs(int dim) const
Definition ElementDofLayout.cpp:68
This class provides a static adjacency list data structure.
Definition AdjacencyList.h:38
const std::vector< LinkData > & array() const
Return contiguous array of links for all nodes (const version).
Definition AdjacencyList.h:173
const std::vector< std::int32_t > & offsets() const
Offset for each node in array() (const version).
Definition AdjacencyList.h:179
A bidirectional map relating entities in one topology to another.
Definition EntityMap.h:20
Geometry stores the geometry imposed on a mesh.
Definition Geometry.h:34
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
Requirements on function for geometry marking.
Definition utils.h:661
void reorder_list(std::span< T > list, std::span< const std::int32_t > nodemap)
Re-order the nodes of a fixed-degree adjacency list.
Definition utils.h:54
std::tuple< std::vector< std::int32_t >, std::vector< T >, std::vector< std::int32_t > > compute_vertex_coords_boundary(const mesh::Mesh< T > &mesh, int dim, std::span< const std::int32_t > facets)
Compute the coordinates of 'vertices' for entities of a given dimension that are attached to specifie...
Definition utils.h:85
std::pair< std::vector< T >, std::array< std::size_t, 2 > > compute_vertex_coords(const mesh::Mesh< T > &mesh)
The coordinates for all 'vertices' in the mesh.
Definition utils.h:615
std::vector< typename std::remove_reference_t< typename U::value_type > > distribute_data(MPI_Comm comm0, std::span< const std::int64_t > indices, MPI_Comm comm1, const U &x, int shape1)
Distribute rows of a rectangular data array to ranks where they are required (scalable version).
Definition MPI.h:679
int size(MPI_Comm comm)
Definition MPI.cpp:72
std::pair< IndexMap, std::vector< std::int32_t > > create_sub_index_map(const IndexMap &imap, std::span< const std::int32_t > indices, IndexMapOrder order=IndexMapOrder::any, bool allow_owner_change=false)
Create a new index map from a subset of indices in an existing index map.
Definition IndexMap.cpp:816
@ any
Allow arbitrary ordering of ghost indices in sub-maps.
Definition IndexMap.h:29
Finite element method functionality.
Definition assemble_expression_impl.h:23
Geometry data structures and algorithms.
Definition BoundingBoxTree.h:22
std::tuple< graph::AdjacencyList< std::int64_t >, std::vector< int >, std::vector< std::int64_t >, std::vector< int > > distribute(MPI_Comm comm, const graph::AdjacencyList< std::int64_t > &list, const graph::AdjacencyList< std::int32_t > &destinations)
Distribute adjacency list nodes to destination ranks.
Definition partition.cpp:38
Graph data structures and algorithms.
Definition AdjacencyList.h:20
std::vector< std::int32_t > reorder_gps(const graph::AdjacencyList< std::int32_t > &graph)
Re-order a graph using the Gibbs-Poole-Stockmeyer algorithm.
Definition ordering.cpp:361
std::function< graph::AdjacencyList< std::int32_t >( MPI_Comm, int, const AdjacencyList< std::int64_t > &, bool)> partition_fn
Signature of functions for computing the parallel partitioning of a distributed graph.
Definition partition.h:31
AdjacencyList< std::int32_t > partition_graph(MPI_Comm comm, int nparts, const AdjacencyList< std::int64_t > &local_graph, bool ghosting)
Partition graph across processes using the default graph partitioner.
Definition partition.cpp:21
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32
CellPartitionFunction create_cell_partitioner(mesh::GhostMode ghost_mode=mesh::GhostMode::none, const graph::partition_fn &partfn=&graph::partition_graph, std::optional< std::int32_t > max_facet_to_cell_links=2)
Create a function that computes destination rank for mesh cells on this rank by applying the default ...
Definition utils.cpp:100
Topology create_topology(MPI_Comm comm, const std::vector< CellType > &cell_types, std::vector< std::span< const std::int64_t > > cells, std::vector< std::span< const std::int64_t > > original_cell_index, std::vector< std::span< const int > > ghost_owners, std::span< const std::int64_t > boundary_vertices)
Create a mesh topology.
Definition Topology.cpp:1025
auto create_boundary_vertices_fn(const CellReorderFunction &reorder_fn, std::optional< std::int32_t > max_facet_to_cell_links=2)
Creates the default boundary vertices routine for a given reorder function.
Definition utils.h:229
std::tuple< graph::AdjacencyList< std::int32_t >, std::vector< std::int64_t >, std::size_t, std::vector< std::int32_t > > build_local_dual_graph(std::span< const CellType > celltypes, const std::vector< std::span< const std::int64_t > > &cells, std::optional< std::int32_t > max_facet_to_cell_links=2)
Compute the local part of the dual graph (cell-cell connections via facets) and facets with only one ...
Definition graphbuild.cpp:524
std::vector< T > cell_normals(const Mesh< T > &mesh, int dim, std::span< const std::int32_t > entities)
Compute normal to given cell (viewed as embedded in 3D).
Definition utils.h:460
CellType cell_entity_type(CellType type, int d, int index)
Return type of cell for entity of dimension d at given entity index.
Definition cell_types.cpp:64
std::function< std::vector< std::int32_t >( const graph::AdjacencyList< std::int32_t > &)> CellReorderFunction
Function that reorders (locally) cells that are owned by this process. It takes the local mesh dual g...
Definition utils.h:217
std::tuple< Topology, std::vector< int32_t >, std::vector< int32_t > > create_subtopology(const Topology &topology, int dim, std::span< const std::int32_t > entities)
Create a topology for a subset of entities of a given topological dimension.
Definition Topology.cpp:1313
std::tuple< Mesh< T >, EntityMap, EntityMap, std::vector< std::int32_t > > create_submesh(const Mesh< T > &mesh, int dim, std::span< const std::int32_t > entities)
Create a new mesh consisting of a subset of entities in a mesh.
Definition utils.h:1394
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
std::vector< std::int32_t > locate_entities_boundary(const Mesh< T > &mesh, int dim, U marker)
Compute indices of all mesh entities that are attached to an owned boundary facet and evaluate to tru...
Definition utils.h:780
CellType
Cell type identifier.
Definition cell_types.h:22
int num_cell_vertices(CellType type)
Definition cell_types.cpp:147
std::vector< T > h(const Mesh< T > &mesh, std::span< const std::int32_t > entities, int dim)
Compute greatest distance between any two vertices of the mesh entities (h).
Definition utils.h:408
std::pair< Geometry< T >, std::vector< int32_t > > create_subgeometry(const Mesh< T > &mesh, int dim, std::span< const std::int32_t > subentity_to_entity)
Create a sub-geometry from a mesh and a subset of mesh entities to be included.
Definition utils.h:1303
Mesh< typename std::remove_reference_t< typename U::value_type > > create_mesh(MPI_Comm comm, MPI_Comm commt, std::vector< std::span< const std::int64_t > > cells, const std::vector< fem::CoordinateElement< typename std::remove_reference_t< typename U::value_type > > > &elements, MPI_Comm commg, const U &x, std::array< std::size_t, 2 > xshape, const CellPartitionFunction &partitioner, std::optional< std::int32_t > max_facet_to_cell_links, const CellReorderFunction &reorder_fn=graph::reorder_gps)
Create a distributed mesh::Mesh from mesh data and using the provided graph partitioning function for...
Definition utils.h:1025
Geometry< typename std::remove_reference_t< typename U::value_type > > create_geometry(const Topology &topology, const std::vector< fem::CoordinateElement< std::remove_reference_t< typename U::value_type > > > &elements, std::span< const std::int64_t > nodes, std::span< const std::int64_t > xdofs, const U &x, int dim, const std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> &reorder_fn=nullptr)
Build Geometry from input data.
Definition Geometry.h:237
std::pair< std::vector< std::int32_t >, std::array< std::size_t, 2 > > entities_to_geometry(const Mesh< T > &mesh, int dim, std::span< const std::int32_t > entities, bool permute=false)
Compute the geometry degrees of freedom associated with the closure of a given set of cell entities.
Definition utils.h:856
std::vector< std::int32_t > compute_incident_entities(const Topology &topology, std::span< const std::int32_t > entities, int d0, int d1)
Compute incident entities.
Definition utils.cpp:124
std::vector< std::int64_t > extract_topology(CellType cell_type, const fem::ElementDofLayout &layout, std::span< const std::int64_t > cells)
Extract topology from cell data, i.e. extract cell vertices.
Definition utils.cpp:30
std::vector< std::int32_t > locate_entities(const Mesh< T > &mesh, int dim, U marker, int entity_type_idx)
Compute indices of all mesh entities that evaluate to true for the provided geometric marking functio...
Definition utils.h:682
std::function< graph::AdjacencyList< std::int32_t >( MPI_Comm comm, int nparts, const std::vector< CellType > &cell_types, const std::vector< std::span< const std::int64_t > > &cells)> CellPartitionFunction
Signature for the cell partitioning function. Function that implement this interface compute the dest...
Definition utils.h:209
std::vector< T > compute_midpoints(const Mesh< T > &mesh, int dim, std::span< const std::int32_t > entities)
Compute the midpoints for mesh entities of a given dimension.
Definition utils.h:577
GhostMode
Enum for different partitioning ghost modes.
Definition utils.h:41
std::vector< std::int32_t > sort_by_perm(std::span< const T > x, std::size_t shape1)
Compute the permutation array that sorts a 2D array by row.
Definition sort.h:183
constexpr __radix_sort radix_sort
Radix sort.
Definition sort.h:170