DOLFINx 0.11.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{
247 return [&, max_facet_to_cell_links](
248 const std::vector<CellType>& celltypes,
249 const std::vector<fem::ElementDofLayout>& doflayouts,
250 const std::vector<std::vector<int>>& ghost_owners,
251 std::vector<std::vector<std::int64_t>>& cells,
252 std::vector<std::vector<std::int64_t>>& cells_v,
253 std::vector<std::vector<std::int64_t>>& original_idx)
254 -> std::vector<std::int64_t>
255 {
256 // Build local dual graph for owned cells to (i) get list of vertices
257 // on the process boundary and (ii) apply re-ordering to cells for
258 // locality
259
260 spdlog::info("Build local dual graphs, re-order cells, and compute process "
261 "boundary vertices.");
262
263 std::vector<std::pair<std::vector<std::int64_t>, int>> facets;
264
265 // Build lists of cells (by cell type) that excludes ghosts
266 std::vector<std::span<const std::int64_t>> cells1_v_local;
267 for (std::size_t i = 0; i < celltypes.size(); ++i)
268 {
269 int num_cell_vertices = mesh::num_cell_vertices(celltypes[i]);
270 std::size_t num_owned_cells
271 = cells_v[i].size() / num_cell_vertices - ghost_owners[i].size();
272 cells1_v_local.emplace_back(cells_v[i].data(),
273 num_owned_cells * num_cell_vertices);
274
275 // Build local dual graph for cell type
276 auto [graph, unmatched_facets, max_v, _facet_attached_cells]
277 = build_local_dual_graph(std::vector{celltypes[i]},
278 std::vector{cells1_v_local.back()},
279 max_facet_to_cell_links);
280
281 // Store unmatched_facets for current cell type
282 facets.emplace_back(std::move(unmatched_facets), max_v);
283
284 // Compute re-ordering of graph
285 const std::vector<std::int32_t> remap = reorder_fn(graph);
286
287 // Update 'original' indices
288 const std::vector<std::int64_t>& orig_idx = original_idx[i];
289 std::vector<std::int64_t> _original_idx(orig_idx.size());
290 std::copy_n(orig_idx.rbegin(), ghost_owners[i].size(),
291 _original_idx.rbegin());
292 {
293 for (std::size_t j = 0; j < remap.size(); ++j)
294 _original_idx[remap[j]] = orig_idx[j];
295 }
296 original_idx[i] = _original_idx;
297
298 // Reorder cells
299 impl::reorder_list(
300 std::span(cells_v[i].data(), remap.size() * num_cell_vertices),
301 remap);
302 impl::reorder_list(
303 std::span(cells[i].data(), remap.size() * doflayouts[i].num_dofs()),
304 remap);
305 }
306
307 if (facets.size() == 1) // Optimisation for single cell type
308 {
309 std::vector<std::int64_t>& vertices = facets.front().first;
310
311 // Remove duplicated vertex indices
312 std::ranges::sort(vertices);
313 auto [unique_end, range_end] = std::ranges::unique(vertices);
314 vertices.erase(unique_end, range_end);
315
316 // Remove -1 if it appears as first entity. This can happen in
317 // mixed topology meshes where '-1' is used to pad facet data when
318 // cells facets have differing numbers of vertices.
319 if (!vertices.empty() and vertices.front() == -1)
320 vertices.erase(vertices.begin());
321
322 return vertices;
323 }
324 else
325 {
326 // Pack 'unmatched' facets for all cell types into single array
327 // (facets0)
328 std::vector<std::int64_t> facets0;
329 facets0.reserve(std::accumulate(facets.begin(), facets.end(),
330 std::size_t(0), [](std::size_t x, auto& y)
331 { return x + y.first.size(); }));
332 int max_v = std::ranges::max_element(facets, [](auto& a, auto& b)
333 { return a.second < b.second; })
334 ->second;
335 for (const auto& [v_data, num_v] : facets)
336 {
337 for (auto it = v_data.begin(); it != v_data.end(); it += num_v)
338 {
339 facets0.insert(facets0.end(), it, std::next(it, num_v));
340 facets0.insert(facets0.end(), max_v - num_v, -1);
341 }
342 }
343
344 // Compute row permutation
345 const std::vector<std::int32_t> perm = dolfinx::sort_by_perm(
346 std::span<const std::int64_t>(facets0), max_v);
347
348 // For facets in facets0 that appear only once, store the facet
349 // vertices
350 std::vector<std::int64_t> vertices;
351 // TODO: allocate memory for vertices
352 auto it = perm.begin();
353 while (it != perm.end())
354 {
355 // Find iterator to next facet different from f and trim any -1
356 // padding
357 std::span _f(facets0.data() + (*it) * max_v, max_v);
358 auto end = std::find_if(_f.rbegin(), _f.rend(),
359 [](auto a) { return a >= 0; });
360 auto f = _f.first(std::distance(end, _f.rend()));
361
362 auto it1 = std::find_if_not(
363 it, perm.end(),
364 [f, max_v, it0 = facets0.begin()](auto p) -> bool
365 {
366 return std::equal(f.begin(), f.end(), std::next(it0, p * max_v));
367 });
368
369 // If no repeated facet found, insert f vertices
370 if (std::distance(it, it1) == 1)
371 vertices.insert(vertices.end(), f.begin(), f.end());
372 else if (std::distance(it, it1) > 2)
373 throw std::runtime_error("More than two matching facets found.");
374
375 // Advance iterator
376 it = it1;
377 }
378
379 // Remove duplicate indices
380 std::ranges::sort(vertices);
381 auto [unique_end, range_end] = std::ranges::unique(vertices);
382 vertices.erase(unique_end, range_end);
383
384 return vertices;
385 }
386 };
387}
388
398std::vector<std::int64_t> extract_topology(CellType cell_type,
399 const fem::ElementDofLayout& layout,
400 std::span<const std::int64_t> cells);
401
410template <std::floating_point T>
411std::vector<T> h(const Mesh<T>& mesh, std::span<const std::int32_t> entities,
412 int dim)
413{
414 if (entities.empty())
415 return std::vector<T>();
416 if (dim == 0)
417 return std::vector<T>(entities.size(), 0);
418
419 // Get the geometry dofs for the vertices of each entity
420 const auto [vertex_xdofs, xdof_shape]
421 = entities_to_geometry(mesh, dim, entities, false);
422
423 // Get the geometry coordinate
424 std::span<const T> x = mesh.geometry().x();
425
426 // Function to compute the length of (p0 - p1)
427 auto delta_norm = [](auto&& p0, auto&& p1)
428 {
429 T norm = 0;
430 for (std::size_t i = 0; i < 3; ++i)
431 norm += (p0[i] - p1[i]) * (p0[i] - p1[i]);
432 return std::sqrt(norm);
433 };
434
435 // Compute greatest distance between any to vertices
436 assert(dim > 0);
437 std::vector<T> h(entities.size(), 0);
438 for (std::size_t e = 0; e < entities.size(); ++e)
439 {
440 // Get geometry 'dof' for each vertex of entity e
441 std::span<const std::int32_t> e_vertices(
442 vertex_xdofs.data() + e * xdof_shape[1], xdof_shape[1]);
443
444 // Compute maximum distance between any two vertices
445 for (std::size_t i = 0; i < e_vertices.size(); ++i)
446 {
447 std::span<const T, 3> p0(x.data() + 3 * e_vertices[i], 3);
448 for (std::size_t j = i + 1; j < e_vertices.size(); ++j)
449 {
450 std::span<const T, 3> p1(x.data() + 3 * e_vertices[j], 3);
451 h[e] = std::max(h[e], delta_norm(p0, p1));
452 }
453 }
454 }
455
456 return h;
457}
458
462template <std::floating_point T>
463std::vector<T> cell_normals(const Mesh<T>& mesh, int dim,
464 std::span<const std::int32_t> entities)
465{
466 if (entities.empty())
467 return std::vector<T>();
468
469 auto topology = mesh.topology();
470 assert(topology);
471 if (topology->cell_type() == CellType::prism and dim == 2)
472 {
473 throw std::runtime_error(
474 "Cell normal computation for prism cells not yet supported.");
475 }
476
477 const int gdim = mesh.geometry().dim();
478 const CellType type = cell_entity_type(topology->cell_type(), dim, 0);
479
480 // Find geometry nodes for topology entities
481 std::span<const T> x = mesh.geometry().x();
482 const auto [geometry_entities, eshape]
483 = entities_to_geometry(mesh, dim, entities, false);
484
485 std::vector<T> n(entities.size() * 3);
486 switch (type)
487 {
488 case CellType::interval:
489 {
490 if (gdim > 2)
491 throw std::invalid_argument("Interval cell normal undefined in 3D.");
492 for (std::size_t i = 0; i < entities.size(); ++i)
493 {
494 // Get the two vertices as points
495 std::array vertices{geometry_entities[i * eshape[1]],
496 geometry_entities[i * eshape[1] + 1]};
497 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
498 std::span<const T, 3>(x.data() + 3 * vertices[1], 3)};
499
500 // Define normal by rotating tangent counter-clockwise
501 std::array<T, 3> t;
502 std::ranges::transform(p[1], p[0], t.begin(),
503 [](auto x, auto y) { return x - y; });
504
505 T norm = std::sqrt(t[0] * t[0] + t[1] * t[1]);
506 std::span<T, 3> ni(n.data() + 3 * i, 3);
507 ni[0] = -t[1] / norm;
508 ni[1] = t[0] / norm;
509 ni[2] = 0.0;
510 }
511 return n;
512 }
513 case CellType::triangle:
514 {
515 for (std::size_t i = 0; i < entities.size(); ++i)
516 {
517 // Get the three vertices as points
518 std::array vertices = {geometry_entities[i * eshape[1] + 0],
519 geometry_entities[i * eshape[1] + 1],
520 geometry_entities[i * eshape[1] + 2]};
521 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
522 std::span<const T, 3>(x.data() + 3 * vertices[1], 3),
523 std::span<const T, 3>(x.data() + 3 * vertices[2], 3)};
524
525 // Compute (p1 - p0) and (p2 - p0)
526 std::array<T, 3> dp1, dp2;
527 std::ranges::transform(p[1], p[0], dp1.begin(),
528 [](auto x, auto y) { return x - y; });
529 std::ranges::transform(p[2], p[0], dp2.begin(),
530 [](auto x, auto y) { return x - y; });
531
532 // Define cell normal via cross product of first two edges
533 std::array<T, 3> ni = math::cross(dp1, dp2);
534 T norm = std::sqrt(ni[0] * ni[0] + ni[1] * ni[1] + ni[2] * ni[2]);
535 std::ranges::transform(ni, std::next(n.begin(), 3 * i),
536 [norm](auto x) { return x / norm; });
537 }
538
539 return n;
540 }
541 case CellType::quadrilateral:
542 {
543 // TODO: check
544 for (std::size_t i = 0; i < entities.size(); ++i)
545 {
546 // Get the three vertices as points
547 std::array vertices = {geometry_entities[i * eshape[1] + 0],
548 geometry_entities[i * eshape[1] + 1],
549 geometry_entities[i * eshape[1] + 2]};
550 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
551 std::span<const T, 3>(x.data() + 3 * vertices[1], 3),
552 std::span<const T, 3>(x.data() + 3 * vertices[2], 3)};
553
554 // Compute (p1 - p0) and (p2 - p0)
555 std::array<T, 3> dp1, dp2;
556 std::ranges::transform(p[1], p[0], dp1.begin(),
557 [](auto x, auto y) { return x - y; });
558 std::ranges::transform(p[2], p[0], dp2.begin(),
559 [](auto x, auto y) { return x - y; });
560
561 // Define cell normal via cross product of first two edges
562 std::array<T, 3> ni = math::cross(dp1, dp2);
563 T norm = std::sqrt(ni[0] * ni[0] + ni[1] * ni[1] + ni[2] * ni[2]);
564 std::ranges::transform(ni, std::next(n.begin(), 3 * i),
565 [norm](auto x) { return x / norm; });
566 }
567
568 return n;
569 }
570 default:
571 throw std::invalid_argument(
572 "cell_normal not supported for this cell type.");
573 }
574}
575
579template <std::floating_point T>
580std::vector<T> compute_midpoints(const Mesh<T>& mesh, int dim,
581 std::span<const std::int32_t> entities)
582{
583 if (entities.empty())
584 return std::vector<T>();
585
586 std::span<const T> x = mesh.geometry().x();
587
588 // Build map from entity -> geometry dof
589 const auto [e_to_g, eshape]
590 = entities_to_geometry(mesh, dim, entities, false);
591
592 std::vector<T> x_mid(entities.size() * 3, 0);
593 for (std::size_t e = 0; e < entities.size(); ++e)
594 {
595 std::span<T, 3> p(x_mid.data() + 3 * e, 3);
596 std::span<const std::int32_t> rows(e_to_g.data() + e * eshape[1],
597 eshape[1]);
598 for (auto row : rows)
599 {
600 std::span<const T, 3> xg(x.data() + 3 * row, 3);
601 std::ranges::transform(p, xg, p.begin(),
602 [size = rows.size()](auto x, auto y)
603 { return x + y / size; });
604 }
605 }
606
607 return x_mid;
608}
609
610namespace impl
611{
616template <std::floating_point T>
617std::pair<std::vector<T>, std::array<std::size_t, 2>>
619{
620 auto topology = mesh.topology();
621 assert(topology);
622 const int tdim = topology->dim();
623
624 // Create entities and connectivities
625
626 // Get all vertex 'node' indices
627 const std::int32_t num_vertices = topology->index_map(0)->size_local()
628 + topology->index_map(0)->num_ghosts();
629
630 std::vector<std::int32_t> vertex_to_node(num_vertices);
631 for (int cell_type_idx = 0,
632 num_cell_types = topology->entity_types(tdim).size();
633 cell_type_idx < num_cell_types; ++cell_type_idx)
634 {
635 auto x_dofmap = mesh.geometry().dofmap(cell_type_idx);
636 auto c_to_v = topology->connectivity({tdim, cell_type_idx}, {0, 0});
637 assert(c_to_v);
638 for (int c = 0; c < c_to_v->num_nodes(); ++c)
639 {
640 auto x_dofs = md::submdspan(x_dofmap, c, md::full_extent);
641 auto vertices = c_to_v->links(c);
642 for (std::size_t i = 0; i < vertices.size(); ++i)
643 vertex_to_node[vertices[i]] = x_dofs[i];
644 }
645 }
646
647 // Pack coordinates of vertices
648 std::span<const T> x_nodes = mesh.geometry().x();
649 std::vector<T> x_vertices(3 * vertex_to_node.size(), 0.0);
650 for (std::size_t i = 0; i < vertex_to_node.size(); ++i)
651 {
652 std::int32_t pos = 3 * vertex_to_node[i];
653 for (std::size_t j = 0; j < 3; ++j)
654 x_vertices[j * vertex_to_node.size() + i] = x_nodes[pos + j];
655 }
656
657 return {std::move(x_vertices), {3, vertex_to_node.size()}};
658}
659
660} // namespace impl
661
663template <typename Fn, typename T>
664concept MarkerFn = std::is_invocable_r<
665 std::vector<std::int8_t>, Fn,
666 md::mdspan<const T,
667 md::extents<std::size_t, 3, md::dynamic_extent>>>::value;
668
684template <std::floating_point T, MarkerFn<T> U>
685std::vector<std::int32_t> locate_entities(const Mesh<T>& mesh, int dim,
686 U marker, int entity_type_idx)
687{
688
689 using cmdspan3x_t
690 = md::mdspan<const T, md::extents<std::size_t, 3, md::dynamic_extent>>;
691
692 // Run marker function on vertex coordinates
693 const auto [xdata, xshape] = impl::compute_vertex_coords(mesh);
694
695 cmdspan3x_t x(xdata.data(), xshape);
696 const std::vector<std::int8_t> marked = marker(x);
697 if (marked.size() != x.extent(1))
698 throw std::runtime_error("Length of array of markers is wrong.");
699
700 auto topology = mesh.topology();
701 assert(topology);
702 const int tdim = topology->dim();
703
704 mesh.topology_mutable()->create_entities(dim);
705 if (dim < tdim)
706 mesh.topology_mutable()->create_connectivity(dim, 0);
707
708 // Iterate over entities of dimension 'dim' to build vector of marked
709 // entities
710 auto e_to_v = topology->connectivity({dim, entity_type_idx}, {0, 0});
711 assert(e_to_v);
712 std::vector<std::int32_t> entities;
713 for (int e = 0; e < e_to_v->num_nodes(); ++e)
714 {
715 // Iterate over entity vertices
716 bool all_vertices_marked = true;
717 for (std::int32_t v : e_to_v->links(e))
718 {
719 if (!marked[v])
720 {
721 all_vertices_marked = false;
722 break;
723 }
724 }
725
726 if (all_vertices_marked)
727 entities.push_back(e);
728 }
729
730 return entities;
731}
732
746template <std::floating_point T, MarkerFn<T> U>
747std::vector<std::int32_t> locate_entities(const Mesh<T>& mesh, int dim,
748 U marker)
749{
750 const int num_entity_types = mesh.topology()->entity_types(dim).size();
751 if (num_entity_types > 1)
752 {
753 throw std::runtime_error(
754 "Multiple entity types of this dimension. Specify entity type index");
755 }
756 return locate_entities(mesh, dim, marker, 0);
757}
758
782template <std::floating_point T, MarkerFn<T> U>
783std::vector<std::int32_t> locate_entities_boundary(const Mesh<T>& mesh, int dim,
784 U marker)
785{
786 // TODO Rewrite this function, it should be possible to simplify considerably
787 auto topology = mesh.topology();
788 assert(topology);
789 int tdim = topology->dim();
790 if (dim == tdim)
791 {
792 throw std::runtime_error(
793 "Cannot use mesh::locate_entities_boundary (boundary) for cells.");
794 }
795
796 // Compute list of boundary facets
797 mesh.topology_mutable()->create_entities(tdim - 1);
798 mesh.topology_mutable()->create_connectivity(tdim - 1, tdim);
799 std::vector<std::int32_t> boundary_facets = exterior_facet_indices(*topology);
800
801 using cmdspan3x_t
802 = md::mdspan<const T, md::extents<std::size_t, 3, md::dynamic_extent>>;
803
804 // Run marker function on the vertex coordinates
805 auto [facet_entities, xdata, vertex_to_pos]
806 = impl::compute_vertex_coords_boundary(mesh, dim, boundary_facets);
807 cmdspan3x_t x(xdata.data(), 3, xdata.size() / 3);
808 std::vector<std::int8_t> marked = marker(x);
809 if (marked.size() != x.extent(1))
810 throw std::runtime_error("Length of array of markers is wrong.");
811
812 // Loop over entities and check vertex markers
813 mesh.topology_mutable()->create_entities(dim);
814 auto e_to_v = topology->connectivity(dim, 0);
815 assert(e_to_v);
816 std::vector<std::int32_t> entities;
817 for (auto e : facet_entities)
818 {
819 // Iterate over entity vertices
820 bool all_vertices_marked = true;
821 for (auto v : e_to_v->links(e))
822 {
823 const std::int32_t pos = vertex_to_pos[v];
824 if (!marked[pos])
825 {
826 all_vertices_marked = false;
827 break;
828 }
829 }
830
831 // Mark facet with all vertices marked
832 if (all_vertices_marked)
833 entities.push_back(e);
834 }
835
836 return entities;
837}
838
857template <std::floating_point T>
858std::pair<std::vector<std::int32_t>, std::array<std::size_t, 2>>
860 std::span<const std::int32_t> entities,
861 bool permute = false)
862{
863 auto topology = mesh.topology();
864 assert(topology);
865 CellType cell_type = topology->cell_type();
866 if ((cell_type == CellType::prism or cell_type == CellType::pyramid)
867 and dim == 2)
868 {
869 throw std::runtime_error("mesh::entities_to_geometry for prism/pyramid "
870 "cell facets not yet supported.");
871 }
872
873 const int tdim = topology->dim();
874 const Geometry<T>& geometry = mesh.geometry();
875 auto xdofs = geometry.dofmap();
876
877 // Get the DOF layout and the number of DOFs per entity
878 const fem::CoordinateElement<T>& coord_ele = geometry.cmap();
879 const fem::ElementDofLayout layout = coord_ele.create_dof_layout();
880 const std::size_t num_entity_dofs = layout.entity_closure_dofs(dim, 0).size();
881 std::vector<std::int32_t> entity_xdofs;
882 entity_xdofs.reserve(entities.size() * num_entity_dofs);
883 std::array<std::size_t, 2> eshape{entities.size(), num_entity_dofs};
884
885 // Get the element's closure DOFs
886 const std::vector<std::vector<std::vector<int>>>& closure_dofs_all
887 = layout.entity_closure_dofs_all();
888
889 // Special case when dim == tdim (cells)
890 if (dim == tdim)
891 {
892 for (std::int32_t c : entities)
893 {
894 // Extract degrees of freedom
895 auto x_c = md::submdspan(xdofs, c, md::full_extent);
896 for (std::int32_t entity_dof : closure_dofs_all[tdim][0])
897 entity_xdofs.push_back(x_c[entity_dof]);
898 }
899
900 return {std::move(entity_xdofs), eshape};
901 }
902
903 assert(dim != tdim);
904
905 auto e_to_c = topology->connectivity(dim, tdim);
906 if (!e_to_c)
907 {
908 throw std::runtime_error(
909 "Entity-to-cell connectivity has not been computed. Missing dims "
910 + std::to_string(dim) + "->" + std::to_string(tdim));
911 }
912
913 auto c_to_e = topology->connectivity(tdim, dim);
914 if (!c_to_e)
915 {
916 throw std::runtime_error(
917 "Cell-to-entity connectivity has not been computed. Missing dims "
918 + std::to_string(tdim) + "->" + std::to_string(dim));
919 }
920
921 // Get the cell info, which is needed to permute the closure dofs
922 std::span<const std::uint32_t> cell_info;
923 if (permute)
924 cell_info = std::span(mesh.topology()->get_cell_permutation_info());
925
926 for (std::int32_t e : entities)
927 {
928 // Get a cell connected to the entity
929 assert(!e_to_c->links(e).empty());
930 std::int32_t c = e_to_c->links(e).front();
931
932 // Get the local index of the entity
933 std::span<const std::int32_t> cell_entities = c_to_e->links(c);
934 auto it = std::find(cell_entities.begin(), cell_entities.end(), e);
935 assert(it != cell_entities.end());
936 std::size_t local_entity = std::distance(cell_entities.begin(), it);
937
938 // Cell sub-entities must be permuted so that their local
939 // orientation agrees with their global orientation
940 std::vector<std::int32_t> closure_dofs(closure_dofs_all[dim][local_entity]);
941 if (permute)
942 {
943 mesh::CellType entity_type
944 = mesh::cell_entity_type(cell_type, dim, local_entity);
945 coord_ele.permute_subentity_closure(closure_dofs, cell_info[c],
946 entity_type, local_entity);
947 }
948
949 // Extract degrees of freedom
950 auto x_c = md::submdspan(xdofs, c, md::full_extent);
951 for (std::int32_t entity_dof : closure_dofs)
952 entity_xdofs.push_back(x_c[entity_dof]);
953 }
954
955 return {std::move(entity_xdofs), eshape};
956}
957
971 const graph::partition_fn& partfn,
972 std::optional<std::int32_t> max_facet_to_cell_links);
973
985 std::optional<std::int32_t> max_facet_to_cell_links
986 = 2);
987
995std::vector<std::int32_t>
996compute_incident_entities(const Topology& topology,
997 std::span<const std::int32_t> entities, int d0,
998 int d1);
999
1043template <typename U>
1045 MPI_Comm comm, MPI_Comm commt,
1046 std::vector<std::span<const std::int64_t>> cells,
1047 const std::vector<fem::CoordinateElement<
1048 typename std::remove_reference_t<typename U::value_type>>>& elements,
1049 MPI_Comm commg, const U& x, std::array<std::size_t, 2> xshape,
1050 const CellPartitionFunction& partitioner,
1051 std::optional<std::int32_t> max_facet_to_cell_links,
1052 const CellReorderFunction& reorder_fn = graph::reorder_gps)
1053{
1054 assert(cells.size() == elements.size());
1055 std::vector<CellType> celltypes;
1056 std::ranges::transform(elements, std::back_inserter(celltypes),
1057 [](auto& e) { return e.cell_shape(); });
1058 std::vector<fem::ElementDofLayout> doflayouts;
1059 std::ranges::transform(elements, std::back_inserter(doflayouts),
1060 [](auto& e) { return e.create_dof_layout(); });
1061
1062 // Note: `extract_topology` extracts topology data, i.e. just the
1063 // vertices. For P1 geometry this should just be the identity
1064 // operator. For other elements the filtered lists may have 'gaps',
1065 // i.e. the indices might not be contiguous.
1066 //
1067 // `extract_topology` could be skipped for 'P1 geometry' elements
1068
1069 std::int32_t num_cell_types = cells.size();
1070
1071 // -- Partition topology across ranks of comm
1072 std::vector<std::vector<std::int64_t>> cells1(num_cell_types);
1073 std::vector<std::vector<std::int64_t>> original_idx1(num_cell_types);
1074 std::vector<std::vector<int>> ghost_owners(num_cell_types);
1075 if (partitioner)
1076 {
1077 spdlog::info("Using partitioner with cell data ({} cell types)",
1078 num_cell_types);
1080 if (commt != MPI_COMM_NULL)
1081 {
1082 int size = dolfinx::MPI::size(comm);
1083 std::vector<std::vector<std::int64_t>> t(num_cell_types);
1084 std::vector<std::span<const std::int64_t>> tspan(num_cell_types);
1085 for (std::int32_t i = 0; i < num_cell_types; ++i)
1086 {
1087 t[i] = extract_topology(celltypes[i], doflayouts[i], cells[i]);
1088 tspan[i] = std::span(t[i]);
1089 }
1090 dest = partitioner(commt, size, celltypes, tspan);
1091 }
1092
1093 std::int32_t cell_offset = 0;
1094 for (std::int32_t i = 0; i < num_cell_types; ++i)
1095 {
1096 std::size_t num_cell_nodes = doflayouts[i].num_dofs();
1097 assert(cells[i].size() % num_cell_nodes == 0);
1098 std::size_t num_cells = cells[i].size() / num_cell_nodes;
1099
1100 // Extract destination AdjacencyList for this cell type
1101 std::vector<std::int32_t> offsets_i(
1102 std::next(dest.offsets().begin(), cell_offset),
1103 std::next(dest.offsets().begin(), cell_offset + num_cells + 1));
1104 std::vector<std::int32_t> data_i(
1105 std::next(dest.array().begin(), offsets_i.front()),
1106 std::next(dest.array().begin(), offsets_i.back()));
1107 std::int32_t offset_0 = offsets_i.front();
1108 std::ranges::for_each(offsets_i,
1109 [&offset_0](std::int32_t& j) { j -= offset_0; });
1110 graph::AdjacencyList<std::int32_t> dest_i(data_i, offsets_i);
1111 cell_offset += num_cells;
1112
1113 // Distribute cells (topology, includes higher-order 'nodes') to
1114 // destination rank
1115 std::vector<int> src_ranks;
1116 std::tie(cells1[i], src_ranks, original_idx1[i], ghost_owners[i])
1117 = graph::build::distribute(comm, cells[i],
1118 {num_cells, num_cell_nodes}, dest_i);
1119 spdlog::debug("Got {} cells from distribution", cells1[i].size());
1120 }
1121 }
1122 else
1123 {
1124 // No partitioning, construct a global index
1125 std::int64_t num_owned = 0;
1126 for (std::int32_t i = 0; i < num_cell_types; ++i)
1127 {
1128 cells1[i] = std::vector<std::int64_t>(cells[i].begin(), cells[i].end());
1129 std::int32_t num_cell_nodes = doflayouts[i].num_dofs();
1130 assert(cells1[i].size() % num_cell_nodes == 0);
1131 original_idx1[i].resize(cells1[i].size() / num_cell_nodes);
1132 num_owned += original_idx1[i].size();
1133 }
1134
1135 // Add on global offset
1136 std::int64_t global_offset = 0;
1137 MPI_Exscan(&num_owned, &global_offset, 1, MPI_INT64_T, MPI_SUM, comm);
1138 for (std::int32_t i = 0; i < num_cell_types; ++i)
1139 {
1140 std::iota(original_idx1[i].begin(), original_idx1[i].end(),
1141 global_offset);
1142 global_offset += original_idx1[i].size();
1143 }
1144 }
1145
1146 // Extract cell 'topology', i.e. extract the vertices for each cell
1147 // and discard any 'higher-order' nodes
1148 std::vector<std::vector<std::int64_t>> cells1_v(num_cell_types);
1149 for (std::int32_t i = 0; i < num_cell_types; ++i)
1150 {
1151 cells1_v[i] = extract_topology(celltypes[i], doflayouts[i], cells1[i]);
1152 spdlog::info("Extract basic topology: {}->{}", cells1[i].size(),
1153 cells1_v[i].size());
1154 }
1155
1156 auto boundary_v_fn
1157 = create_boundary_vertices_fn(reorder_fn, max_facet_to_cell_links);
1158 const std::vector<std::int64_t> boundary_v = boundary_v_fn(
1159 celltypes, doflayouts, ghost_owners, cells1, cells1_v, original_idx1);
1160
1161 spdlog::debug("Got {} boundary vertices", boundary_v.size());
1162
1163 // Create Topology
1164 std::vector<std::span<const std::int64_t>> cells1_v_span;
1165 std::ranges::transform(cells1_v, std::back_inserter(cells1_v_span),
1166 [](auto& c) { return std::span(c); });
1167 std::vector<std::span<const std::int64_t>> original_idx1_span;
1168 std::ranges::transform(original_idx1, std::back_inserter(original_idx1_span),
1169 [](auto& c) { return std::span(c); });
1170 std::vector<std::span<const int>> ghost_owners_span;
1171 std::ranges::transform(ghost_owners, std::back_inserter(ghost_owners_span),
1172 [](auto& c) { return std::span(c); });
1173 Topology topology
1174 = create_topology(comm, celltypes, cells1_v_span, original_idx1_span,
1175 ghost_owners_span, boundary_v, 0);
1176
1177 // Create connectivities required higher-order geometries for creating
1178 // a Geometry object
1179 for (int i = 0; i < num_cell_types; ++i)
1180 {
1181 const auto& entity_dofs = doflayouts[i].entity_dofs_all();
1182 for (int dim = 1; dim < topology.dim(); ++dim)
1183 {
1184 // Accumulate count of all dofs on this dimension
1185 int dim_sum
1186 = std::accumulate(entity_dofs[dim].begin(), entity_dofs[dim].end(), 0,
1187 [](int c, auto v) { return c + v.size(); });
1188
1189 spdlog::debug("Counting entity dofs, dim={}: {}", dim, dim_sum);
1190 if (dim_sum > 0)
1191 topology.create_entities(dim);
1192 }
1193
1194 if (elements[i].needs_dof_permutations())
1195 topology.create_entity_permutations();
1196 }
1197
1198 // Build list of unique (global) node indices from cells1 and
1199 // distribute coordinate data
1200 std::vector<std::int64_t> nodes1, nodes2;
1201 for (std::vector<std::int64_t>& c : cells1)
1202 nodes1.insert(nodes1.end(), c.begin(), c.end());
1203 for (std::vector<std::int64_t>& c : cells1)
1204 nodes2.insert(nodes2.end(), c.begin(), c.end());
1205
1206 dolfinx::radix_sort(nodes1);
1207 auto [unique_end, range_end] = std::ranges::unique(nodes1);
1208 nodes1.erase(unique_end, range_end);
1209
1210 std::vector coords
1211 = dolfinx::MPI::distribute_data(comm, nodes1, commg, x, xshape[1]);
1212
1213 // Create geometry object
1215 = create_geometry(topology, elements, nodes1, nodes2, coords, xshape[1]);
1216
1217 return Mesh(comm, std::make_shared<Topology>(std::move(topology)),
1218 std::move(geometry));
1219}
1220
1258template <typename U>
1260 MPI_Comm comm, MPI_Comm commt, std::span<const std::int64_t> cells,
1262 typename std::remove_reference_t<typename U::value_type>>& element,
1263 MPI_Comm commg, const U& x, std::array<std::size_t, 2> xshape,
1264 const CellPartitionFunction& partitioner,
1265 std::optional<std::int32_t> max_facet_to_cell_links,
1266 const CellReorderFunction& reorder_fn = graph::reorder_gps)
1267{
1268 return create_mesh(comm, commt, std::vector{cells}, std::vector{element},
1269 commg, x, xshape, partitioner, max_facet_to_cell_links,
1270 reorder_fn);
1271}
1272
1293template <typename U>
1294Mesh<typename std::remove_reference_t<typename U::value_type>>
1295create_mesh(MPI_Comm comm, std::span<const std::int64_t> cells,
1297 std::remove_reference_t<typename U::value_type>>& elements,
1298 const U& x, std::array<std::size_t, 2> xshape, GhostMode ghost_mode,
1299 std::optional<std::int32_t> max_facet_to_cell_links = 2)
1300{
1301 if (dolfinx::MPI::size(comm) == 1)
1302 {
1303 return create_mesh(comm, comm, std::vector{cells}, std::vector{elements},
1304 comm, x, xshape, nullptr, max_facet_to_cell_links);
1305 }
1306 else
1307 {
1308 return create_mesh(
1309 comm, comm, std::vector{cells}, std::vector{elements}, comm, x, xshape,
1310 create_cell_partitioner(ghost_mode, max_facet_to_cell_links),
1311 max_facet_to_cell_links);
1312 }
1313}
1314
1328template <std::floating_point T>
1329std::pair<Geometry<T>, std::vector<int32_t>>
1331 std::span<const std::int32_t> subentity_to_entity)
1332{
1333 const Geometry<T>& geometry = mesh.geometry();
1334
1335 // Get the geometry dofs in the sub-geometry based on the entities in
1336 // sub-geometry
1337 const fem::ElementDofLayout layout = geometry.cmap().create_dof_layout();
1338
1339 const std::vector<std::int32_t> x_indices
1340 = entities_to_geometry(mesh, dim, subentity_to_entity, true).first;
1341
1342 std::vector<std::int32_t> sub_x_dofs = x_indices;
1343 std::ranges::sort(sub_x_dofs);
1344 auto [unique_end, range_end] = std::ranges::unique(sub_x_dofs);
1345 sub_x_dofs.erase(unique_end, range_end);
1346
1347 // Get the sub-geometry dofs owned by this process
1348 auto x_index_map = geometry.index_map();
1349 assert(x_index_map);
1350
1351 std::shared_ptr<common::IndexMap> sub_x_dof_index_map;
1352 std::vector<std::int32_t> subx_to_x_dofmap;
1353 {
1354 auto [map, new_to_old] = common::create_sub_index_map(
1355 *x_index_map, sub_x_dofs, common::IndexMapOrder::any, true);
1356 sub_x_dof_index_map = std::make_shared<common::IndexMap>(std::move(map));
1357 subx_to_x_dofmap = std::move(new_to_old);
1358 }
1359
1360 // Create sub-geometry coordinates
1361 std::span<const T> x = geometry.x();
1362 std::int32_t sub_num_x_dofs = subx_to_x_dofmap.size();
1363 std::vector<T> sub_x(3 * sub_num_x_dofs);
1364 for (std::int32_t i = 0; i < sub_num_x_dofs; ++i)
1365 {
1366 std::copy_n(std::next(x.begin(), 3 * subx_to_x_dofmap[i]), 3,
1367 std::next(sub_x.begin(), 3 * i));
1368 }
1369
1370 // Create geometry to sub-geometry map
1371 std::vector<std::int32_t> x_to_subx_dof_map(
1372 x_index_map->size_local() + x_index_map->num_ghosts(), -1);
1373 for (std::size_t i = 0; i < subx_to_x_dofmap.size(); ++i)
1374 x_to_subx_dof_map[subx_to_x_dofmap[i]] = i;
1375
1376 // Create sub-geometry dofmap
1377 std::vector<std::int32_t> sub_x_dofmap;
1378 sub_x_dofmap.reserve(x_indices.size());
1379 std::ranges::transform(x_indices, std::back_inserter(sub_x_dofmap),
1380 [&x_to_subx_dof_map](auto x_dof)
1381 {
1382 assert(x_to_subx_dof_map[x_dof] != -1);
1383 return x_to_subx_dof_map[x_dof];
1384 });
1385
1386 // Sub-geometry coordinate element
1387 CellType sub_xcell = cell_entity_type(geometry.cmap().cell_shape(), dim, 0);
1388
1389 // Special handling of point meshes, as they only support constant
1390 // basis functions
1391 int degree = (sub_xcell == CellType::point) ? 0 : geometry.cmap().degree();
1392 fem::CoordinateElement<T> sub_cmap(sub_xcell, degree,
1393 geometry.cmap().variant());
1394
1395 // Sub-geometry input_global_indices
1396 const std::vector<std::int64_t>& igi = geometry.input_global_indices();
1397 std::vector<std::int64_t> sub_igi;
1398 sub_igi.reserve(subx_to_x_dofmap.size());
1399 std::ranges::transform(subx_to_x_dofmap, std::back_inserter(sub_igi),
1400 [&igi](auto sub_x_dof) { return igi[sub_x_dof]; });
1401
1402 // Create geometry
1403 return {Geometry(
1404 sub_x_dof_index_map,
1405 std::vector<std::vector<std::int32_t>>{std::move(sub_x_dofmap)},
1406 {sub_cmap}, std::move(sub_x), geometry.dim(), std::move(sub_igi)),
1407 std::move(subx_to_x_dofmap)};
1408}
1409
1419template <std::floating_point T>
1420std::tuple<Mesh<T>, EntityMap, EntityMap, std::vector<std::int32_t>>
1422 std::span<const std::int32_t> entities)
1423{
1424 // Create sub-topology
1425 mesh.topology_mutable()->create_connectivity(dim, 0);
1426 auto [topology, subentity_to_entity, subvertex_to_vertex]
1427 = mesh::create_subtopology(*mesh.topology(), dim, entities);
1428
1429 // Create sub-geometry
1430 const int tdim = mesh.topology()->dim();
1431 mesh.topology_mutable()->create_entities(dim);
1432 mesh.topology_mutable()->create_connectivity(dim, tdim);
1433 mesh.topology_mutable()->create_connectivity(tdim, dim);
1434 mesh.topology_mutable()->create_entity_permutations();
1435 auto [geometry, subx_to_x_dofmap]
1436 = mesh::create_subgeometry(mesh, dim, subentity_to_entity);
1437
1438 Mesh<T> submesh
1439 = Mesh(mesh.comm(), std::make_shared<Topology>(std::move(topology)),
1440 std::move(geometry));
1441 EntityMap entity_map(mesh.topology(), submesh.topology(), dim,
1442 subentity_to_entity);
1443 EntityMap vertex_map(mesh.topology(), submesh.topology(), 0,
1444 subvertex_to_vertex);
1445 return {std::move(submesh), std::move(entity_map), std::move(vertex_map),
1446 std::move(subx_to_x_dofmap)};
1447}
1448
1449} // 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< int > & entity_closure_dofs(int dim, int entity_index) const
Definition ElementDofLayout.cpp:65
const std::vector< std::vector< std::vector< int > > > & entity_closure_dofs_all() const
Definition ElementDofLayout.cpp:77
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:175
const std::vector< std::int32_t > & offsets() const
Offset for each node in array() (const version).
Definition AdjacencyList.h:181
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
std::shared_ptr< Topology > topology()
Get mesh topology.
Definition Mesh.h:69
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:664
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:618
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:39
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:363
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
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32
CellPartitionFunction create_cell_partitioner(mesh::GhostMode ghost_mode, const graph::partition_fn &partfn, std::optional< std::int32_t > max_facet_to_cell_links)
Create a function that computes destination rank for mesh cells on this rank by applying the default ...
Definition utils.cpp:100
std::tuple< graph::AdjacencyList< std::int32_t >, std::vector< std::int64_t >, int, 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)
Compute the local part of the dual graph (cell-cell connections via facets) and facets with only one ...
Definition graphbuild.cpp:530
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::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 > 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:463
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:1369
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:1421
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:783
CellType
Cell type identifier.
Definition cell_types.h:21
int num_cell_vertices(CellType type)
Number vertices for a cell type.
Definition cell_types.cpp:98
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:411
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:1330
auto create_boundary_vertices_fn(const CellReorderFunction &reorder_fn, std::optional< std::int32_t > max_facet_to_cell_links)
Creates the default boundary vertices routine for a given reorder function.
Definition utils.h:229
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:1044
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:859
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:132
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:685
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:580
GhostMode
Enum for different partitioning ghost modes.
Definition utils.h:41
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, int num_threads)
Create a mesh topology.
Definition Topology.cpp:1051
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.h:110
constexpr void radix_sort(R &&range, P proj={})
Sort a range with radix sorting algorithm. The bucket size is determined by the number of bits to sor...
Definition sort.h:78
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:177