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 "Mesh.h"
10#include "Topology.h"
11#include "graphbuild.h"
12#include <algorithm>
13#include <basix/mdspan.hpp>
14#include <concepts>
15#include <dolfinx/graph/AdjacencyList.h>
16#include <dolfinx/graph/ordering.h>
17#include <dolfinx/graph/partition.h>
18#include <functional>
19#include <mpi.h>
20#include <numeric>
21#include <span>
22#include <vector>
23
26
27namespace dolfinx::fem
28{
30}
31
32namespace dolfinx::mesh
33{
34enum class CellType;
35
37enum class GhostMode : int
38{
39 none,
40 shared_facet,
41 shared_vertex
42};
43
44namespace impl
45{
51template <typename T>
52void reorder_list(std::span<T> list, std::span<const std::int32_t> nodemap)
53{
54 if (nodemap.empty())
55 return;
56
57 assert(list.size() % nodemap.size() == 0);
58 std::size_t degree = list.size() / nodemap.size();
59 const std::vector<T> orig(list.begin(), list.end());
60 for (std::size_t n = 0; n < nodemap.size(); ++n)
61 {
62 std::span links_old(orig.data() + n * degree, degree);
63 auto links_new = list.subspan(nodemap[n] * degree, degree);
64 std::ranges::copy(links_old, links_new.begin());
65 }
66}
67
81template <std::floating_point T>
82std::tuple<std::vector<std::int32_t>, std::vector<T>, std::vector<std::int32_t>>
84 std::span<const std::int32_t> facets)
85{
86 auto topology = mesh.topology();
87 assert(topology);
88 const int tdim = topology->dim();
89 if (dim == tdim)
90 {
91 throw std::runtime_error(
92 "Cannot use mesh::locate_entities_boundary (boundary) for cells.");
93 }
94
95 // Build set of vertices on boundary and set of boundary entities
96 mesh.topology_mutable()->create_connectivity(tdim - 1, 0);
97 mesh.topology_mutable()->create_connectivity(tdim - 1, dim);
98 std::vector<std::int32_t> vertices, entities;
99 {
100 auto f_to_v = topology->connectivity(tdim - 1, 0);
101 assert(f_to_v);
102 auto f_to_e = topology->connectivity(tdim - 1, dim);
103 assert(f_to_e);
104 for (auto f : facets)
105 {
106 auto v = f_to_v->links(f);
107 vertices.insert(vertices.end(), v.begin(), v.end());
108 auto e = f_to_e->links(f);
109 entities.insert(entities.end(), e.begin(), e.end());
110 }
111
112 // Build vector of boundary vertices
113 {
114 std::ranges::sort(vertices);
115 auto [unique_end, range_end] = std::ranges::unique(vertices);
116 vertices.erase(unique_end, range_end);
117 }
118
119 {
120 std::ranges::sort(entities);
121 auto [unique_end, range_end] = std::ranges::unique(entities);
122 entities.erase(unique_end, range_end);
123 }
124 }
125
126 // Get geometry data
127 auto x_dofmap = mesh.geometry().dofmap();
128 std::span<const T> x_nodes = mesh.geometry().x();
129
130 // Get all vertex 'node' indices
131 mesh.topology_mutable()->create_connectivity(0, tdim);
132 mesh.topology_mutable()->create_connectivity(tdim, 0);
133 auto v_to_c = topology->connectivity(0, tdim);
134 assert(v_to_c);
135 auto c_to_v = topology->connectivity(tdim, 0);
136 assert(c_to_v);
137 std::vector<T> x_vertices(3 * vertices.size(), -1.0);
138 std::vector<std::int32_t> vertex_to_pos(v_to_c->num_nodes(), -1);
139 for (std::size_t i = 0; i < vertices.size(); ++i)
140 {
141 const std::int32_t v = vertices[i];
142
143 // Get first cell and find position
144 const std::int32_t c = v_to_c->links(v).front();
145 auto cell_vertices = c_to_v->links(c);
146 auto it = std::find(cell_vertices.begin(), cell_vertices.end(), v);
147 assert(it != cell_vertices.end());
148 const std::size_t local_pos = std::distance(cell_vertices.begin(), it);
149
150 auto dofs = md::submdspan(x_dofmap, c, md::full_extent);
151 for (std::size_t j = 0; j < 3; ++j)
152 x_vertices[j * vertices.size() + i] = x_nodes[3 * dofs[local_pos] + j];
153 vertex_to_pos[v] = i;
154 }
155
156 return {std::move(entities), std::move(x_vertices), std::move(vertex_to_pos)};
157}
158
159} // namespace impl
160
174std::vector<std::int32_t> exterior_facet_indices(const Topology& topology,
175 int facet_type_idx);
176
188std::vector<std::int32_t> exterior_facet_indices(const Topology& topology);
189
207using CellPartitionFunction = std::function<graph::AdjacencyList<std::int32_t>(
208 MPI_Comm comm, int nparts, const std::vector<CellType>& cell_types,
209 const std::vector<std::span<const std::int64_t>>& cells)>;
210
215using CellReorderFunction = std::function<std::vector<std::int32_t>(
217
225{
237 return [&](const std::vector<CellType>& celltypes,
238 const std::vector<fem::ElementDofLayout>& doflayouts,
239 const std::vector<std::vector<int>>& ghost_owners,
240 std::vector<std::vector<std::int64_t>>& cells,
241 std::vector<std::vector<std::int64_t>>& cells_v,
242 std::vector<std::vector<std::int64_t>>& original_idx)
243 -> std::vector<std::int64_t>
244 {
245 // Build local dual graph for owned cells to (i) get list of vertices
246 // on the process boundary and (ii) apply re-ordering to cells for
247 // locality
248
249 spdlog::info("Build local dual graphs, re-order cells, and compute process "
250 "boundary vertices.");
251
252 std::vector<std::pair<std::vector<std::int64_t>, int>> facets;
253
254 // Build lists of cells (by cell type) that excludes ghosts
255 std::vector<std::span<const std::int64_t>> cells1_v_local;
256 for (std::size_t i = 0; i < celltypes.size(); ++i)
257 {
258 int num_cell_vertices = mesh::num_cell_vertices(celltypes[i]);
259 std::size_t num_owned_cells
260 = cells_v[i].size() / num_cell_vertices - ghost_owners[i].size();
261 cells1_v_local.emplace_back(cells_v[i].data(),
262 num_owned_cells * num_cell_vertices);
263
264 // Build local dual graph for cell type
265 auto [graph, unmatched_facets, max_v, _facet_attached_cells]
266 = build_local_dual_graph(std::vector{celltypes[i]},
267 std::vector{cells1_v_local.back()});
268
269 // Store unmatched_facets for current cell type
270 facets.emplace_back(std::move(unmatched_facets), max_v);
271
272 // Compute re-ordering of graph
273 const std::vector<std::int32_t> remap = reorder_fn(graph);
274
275 // Update 'original' indices
276 const std::vector<std::int64_t>& orig_idx = original_idx[i];
277 std::vector<std::int64_t> _original_idx(orig_idx.size());
278 std::copy_n(orig_idx.rbegin(), ghost_owners[i].size(),
279 _original_idx.rbegin());
280 {
281 for (std::size_t j = 0; j < remap.size(); ++j)
282 _original_idx[remap[j]] = orig_idx[j];
283 }
284 original_idx[i] = _original_idx;
285
286 // Reorder cells
288 std::span(cells_v[i].data(), remap.size() * num_cell_vertices),
289 remap);
291 std::span(cells[i].data(), remap.size() * doflayouts[i].num_dofs()),
292 remap);
293 }
294
295 if (facets.size() == 1) // Optimisation for single cell type
296 {
297 std::vector<std::int64_t>& vertices = facets.front().first;
298
299 // Remove duplicated vertex indices
300 std::ranges::sort(vertices);
301 auto [unique_end, range_end] = std::ranges::unique(vertices);
302 vertices.erase(unique_end, range_end);
303
304 // Remove -1 if it appears as first entity. This can happen in
305 // mixed topology meshes where '-1' is used to pad facet data when
306 // cells facets have differing numbers of vertices.
307 if (!vertices.empty() and vertices.front() == -1)
308 vertices.erase(vertices.begin());
309
310 return vertices;
311 }
312 else
313 {
314 // Pack 'unmatched' facets for all cell types into single array
315 // (facets0)
316 std::vector<std::int64_t> facets0;
317 facets0.reserve(std::accumulate(facets.begin(), facets.end(),
318 std::size_t(0), [](std::size_t x, auto& y)
319 { return x + y.first.size(); }));
320 int max_v = std::ranges::max_element(facets, [](auto& a, auto& b)
321 { return a.second < b.second; })
322 ->second;
323 for (const auto& [v_data, num_v] : facets)
324 {
325 for (auto it = v_data.begin(); it != v_data.end(); it += num_v)
326 {
327 facets0.insert(facets0.end(), it, std::next(it, num_v));
328 facets0.insert(facets0.end(), max_v - num_v, -1);
329 }
330 }
331
332 // Compute row permutation
333 const std::vector<std::int32_t> perm = dolfinx::sort_by_perm(
334 std::span<const std::int64_t>(facets0), max_v);
335
336 // For facets in facets0 that appear only once, store the facet
337 // vertices
338 std::vector<std::int64_t> vertices;
339 // TODO: allocate memory for vertices
340 auto it = perm.begin();
341 while (it != perm.end())
342 {
343 // Find iterator to next facet different from f and trim any -1
344 // padding
345 std::span _f(facets0.data() + (*it) * max_v, max_v);
346 auto end = std::find_if(_f.rbegin(), _f.rend(),
347 [](auto a) { return a >= 0; });
348 auto f = _f.first(std::distance(end, _f.rend()));
349
350 auto it1 = std::find_if_not(
351 it, perm.end(),
352 [f, max_v, it0 = facets0.begin()](auto p) -> bool
353 {
354 return std::equal(f.begin(), f.end(), std::next(it0, p * max_v));
355 });
356
357 // If no repeated facet found, insert f vertices
358 if (std::distance(it, it1) == 1)
359 vertices.insert(vertices.end(), f.begin(), f.end());
360 else if (std::distance(it, it1) > 2)
361 throw std::runtime_error("More than two matching facets found.");
362
363 // Advance iterator
364 it = it1;
365 }
366
367 // Remove duplicate indices
368 std::ranges::sort(vertices);
369 auto [unique_end, range_end] = std::ranges::unique(vertices);
370 vertices.erase(unique_end, range_end);
371
372 return vertices;
373 }
374 };
375}
376
386std::vector<std::int64_t> extract_topology(CellType cell_type,
387 const fem::ElementDofLayout& layout,
388 std::span<const std::int64_t> cells);
389
398template <std::floating_point T>
399std::vector<T> h(const Mesh<T>& mesh, std::span<const std::int32_t> entities,
400 int dim)
401{
402 if (entities.empty())
403 return std::vector<T>();
404 if (dim == 0)
405 return std::vector<T>(entities.size(), 0);
406
407 // Get the geometry dofs for the vertices of each entity
408 const auto [vertex_xdofs, xdof_shape]
409 = entities_to_geometry(mesh, dim, entities, false);
410
411 // Get the geometry coordinate
412 std::span<const T> x = mesh.geometry().x();
413
414 // Function to compute the length of (p0 - p1)
415 auto delta_norm = [](auto&& p0, auto&& p1)
416 {
417 T norm = 0;
418 for (std::size_t i = 0; i < 3; ++i)
419 norm += (p0[i] - p1[i]) * (p0[i] - p1[i]);
420 return std::sqrt(norm);
421 };
422
423 // Compute greatest distance between any to vertices
424 assert(dim > 0);
425 std::vector<T> h(entities.size(), 0);
426 for (std::size_t e = 0; e < entities.size(); ++e)
427 {
428 // Get geometry 'dof' for each vertex of entity e
429 std::span<const std::int32_t> e_vertices(
430 vertex_xdofs.data() + e * xdof_shape[1], xdof_shape[1]);
431
432 // Compute maximum distance between any two vertices
433 for (std::size_t i = 0; i < e_vertices.size(); ++i)
434 {
435 std::span<const T, 3> p0(x.data() + 3 * e_vertices[i], 3);
436 for (std::size_t j = i + 1; j < e_vertices.size(); ++j)
437 {
438 std::span<const T, 3> p1(x.data() + 3 * e_vertices[j], 3);
439 h[e] = std::max(h[e], delta_norm(p0, p1));
440 }
441 }
442 }
443
444 return h;
445}
446
450template <std::floating_point T>
451std::vector<T> cell_normals(const Mesh<T>& mesh, int dim,
452 std::span<const std::int32_t> entities)
453{
454 if (entities.empty())
455 return std::vector<T>();
456
457 auto topology = mesh.topology();
458 assert(topology);
459 if (topology->cell_type() == CellType::prism and dim == 2)
460 {
461 throw std::runtime_error(
462 "Cell normal computation for prism cells not yet supported.");
463 }
464
465 const int gdim = mesh.geometry().dim();
466 const CellType type = cell_entity_type(topology->cell_type(), dim, 0);
467
468 // Find geometry nodes for topology entities
469 std::span<const T> x = mesh.geometry().x();
470 const auto [geometry_entities, eshape]
471 = entities_to_geometry(mesh, dim, entities, false);
472
473 std::vector<T> n(entities.size() * 3);
474 switch (type)
475 {
476 case CellType::interval:
477 {
478 if (gdim > 2)
479 throw std::invalid_argument("Interval cell normal undefined in 3D.");
480 for (std::size_t i = 0; i < entities.size(); ++i)
481 {
482 // Get the two vertices as points
483 std::array vertices{geometry_entities[i * eshape[1]],
484 geometry_entities[i * eshape[1] + 1]};
485 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
486 std::span<const T, 3>(x.data() + 3 * vertices[1], 3)};
487
488 // Define normal by rotating tangent counter-clockwise
489 std::array<T, 3> t;
490 std::ranges::transform(p[1], p[0], t.begin(),
491 [](auto x, auto y) { return x - y; });
492
493 T norm = std::sqrt(t[0] * t[0] + t[1] * t[1]);
494 std::span<T, 3> ni(n.data() + 3 * i, 3);
495 ni[0] = -t[1] / norm;
496 ni[1] = t[0] / norm;
497 ni[2] = 0.0;
498 }
499 return n;
500 }
501 case CellType::triangle:
502 {
503 for (std::size_t i = 0; i < entities.size(); ++i)
504 {
505 // Get the three vertices as points
506 std::array vertices = {geometry_entities[i * eshape[1] + 0],
507 geometry_entities[i * eshape[1] + 1],
508 geometry_entities[i * eshape[1] + 2]};
509 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
510 std::span<const T, 3>(x.data() + 3 * vertices[1], 3),
511 std::span<const T, 3>(x.data() + 3 * vertices[2], 3)};
512
513 // Compute (p1 - p0) and (p2 - p0)
514 std::array<T, 3> dp1, dp2;
515 std::ranges::transform(p[1], p[0], dp1.begin(),
516 [](auto x, auto y) { return x - y; });
517 std::ranges::transform(p[2], p[0], dp2.begin(),
518 [](auto x, auto y) { return x - y; });
519
520 // Define cell normal via cross product of first two edges
521 std::array<T, 3> ni = math::cross(dp1, dp2);
522 T norm = std::sqrt(ni[0] * ni[0] + ni[1] * ni[1] + ni[2] * ni[2]);
523 std::ranges::transform(ni, std::next(n.begin(), 3 * i),
524 [norm](auto x) { return x / norm; });
525 }
526
527 return n;
528 }
529 case CellType::quadrilateral:
530 {
531 // TODO: check
532 for (std::size_t i = 0; i < entities.size(); ++i)
533 {
534 // Get the three vertices as points
535 std::array vertices = {geometry_entities[i * eshape[1] + 0],
536 geometry_entities[i * eshape[1] + 1],
537 geometry_entities[i * eshape[1] + 2]};
538 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
539 std::span<const T, 3>(x.data() + 3 * vertices[1], 3),
540 std::span<const T, 3>(x.data() + 3 * vertices[2], 3)};
541
542 // Compute (p1 - p0) and (p2 - p0)
543 std::array<T, 3> dp1, dp2;
544 std::ranges::transform(p[1], p[0], dp1.begin(),
545 [](auto x, auto y) { return x - y; });
546 std::ranges::transform(p[2], p[0], dp2.begin(),
547 [](auto x, auto y) { return x - y; });
548
549 // Define cell normal via cross product of first two edges
550 std::array<T, 3> ni = math::cross(dp1, dp2);
551 T norm = std::sqrt(ni[0] * ni[0] + ni[1] * ni[1] + ni[2] * ni[2]);
552 std::ranges::transform(ni, std::next(n.begin(), 3 * i),
553 [norm](auto x) { return x / norm; });
554 }
555
556 return n;
557 }
558 default:
559 throw std::invalid_argument(
560 "cell_normal not supported for this cell type.");
561 }
562}
563
567template <std::floating_point T>
568std::vector<T> compute_midpoints(const Mesh<T>& mesh, int dim,
569 std::span<const std::int32_t> entities)
570{
571 if (entities.empty())
572 return std::vector<T>();
573
574 std::span<const T> x = mesh.geometry().x();
575
576 // Build map from entity -> geometry dof
577 const auto [e_to_g, eshape]
578 = entities_to_geometry(mesh, dim, entities, false);
579
580 std::vector<T> x_mid(entities.size() * 3, 0);
581 for (std::size_t e = 0; e < entities.size(); ++e)
582 {
583 std::span<T, 3> p(x_mid.data() + 3 * e, 3);
584 std::span<const std::int32_t> rows(e_to_g.data() + e * eshape[1],
585 eshape[1]);
586 for (auto row : rows)
587 {
588 std::span<const T, 3> xg(x.data() + 3 * row, 3);
589 std::ranges::transform(p, xg, p.begin(),
590 [size = rows.size()](auto x, auto y)
591 { return x + y / size; });
592 }
593 }
594
595 return x_mid;
596}
597
598namespace impl
599{
604template <std::floating_point T>
605std::pair<std::vector<T>, std::array<std::size_t, 2>>
607{
608 auto topology = mesh.topology();
609 assert(topology);
610 const int tdim = topology->dim();
611
612 // Create entities and connectivities
613
614 // Get all vertex 'node' indices
615 const std::int32_t num_vertices = topology->index_map(0)->size_local()
616 + topology->index_map(0)->num_ghosts();
617
618 std::vector<std::int32_t> vertex_to_node(num_vertices);
619 for (int cell_type_idx = 0,
620 num_cell_types = topology->entity_types(tdim).size();
621 cell_type_idx < num_cell_types; ++cell_type_idx)
622 {
623 auto x_dofmap = mesh.geometry().dofmap(cell_type_idx);
624 auto c_to_v = topology->connectivity({tdim, cell_type_idx}, {0, 0});
625 assert(c_to_v);
626 for (int c = 0; c < c_to_v->num_nodes(); ++c)
627 {
628 auto x_dofs = md::submdspan(x_dofmap, c, md::full_extent);
629 auto vertices = c_to_v->links(c);
630 for (std::size_t i = 0; i < vertices.size(); ++i)
631 vertex_to_node[vertices[i]] = x_dofs[i];
632 }
633 }
634
635 // Pack coordinates of vertices
636 std::span<const T> x_nodes = mesh.geometry().x();
637 std::vector<T> x_vertices(3 * vertex_to_node.size(), 0.0);
638 for (std::size_t i = 0; i < vertex_to_node.size(); ++i)
639 {
640 std::int32_t pos = 3 * vertex_to_node[i];
641 for (std::size_t j = 0; j < 3; ++j)
642 x_vertices[j * vertex_to_node.size() + i] = x_nodes[pos + j];
643 }
644
645 return {std::move(x_vertices), {3, vertex_to_node.size()}};
646}
647
648} // namespace impl
649
651template <typename Fn, typename T>
652concept MarkerFn = std::is_invocable_r<
653 std::vector<std::int8_t>, Fn,
654 md::mdspan<const T,
655 md::extents<std::size_t, 3, md::dynamic_extent>>>::value;
656
672template <std::floating_point T, MarkerFn<T> U>
673std::vector<std::int32_t> locate_entities(const Mesh<T>& mesh, int dim,
674 U marker, int entity_type_idx)
675{
676
677 using cmdspan3x_t
678 = md::mdspan<const T, md::extents<std::size_t, 3, md::dynamic_extent>>;
679
680 // Run marker function on vertex coordinates
681 const auto [xdata, xshape] = impl::compute_vertex_coords(mesh);
682
683 cmdspan3x_t x(xdata.data(), xshape);
684 const std::vector<std::int8_t> marked = marker(x);
685 if (marked.size() != x.extent(1))
686 throw std::runtime_error("Length of array of markers is wrong.");
687
688 auto topology = mesh.topology();
689 assert(topology);
690 const int tdim = topology->dim();
691
692 mesh.topology_mutable()->create_entities(dim);
693 if (dim < tdim)
694 mesh.topology_mutable()->create_connectivity(dim, 0);
695
696 // Iterate over entities of dimension 'dim' to build vector of marked
697 // entities
698 auto e_to_v = topology->connectivity({dim, entity_type_idx}, {0, 0});
699 assert(e_to_v);
700 std::vector<std::int32_t> entities;
701 for (int e = 0; e < e_to_v->num_nodes(); ++e)
702 {
703 // Iterate over entity vertices
704 bool all_vertices_marked = true;
705 for (std::int32_t v : e_to_v->links(e))
706 {
707 if (!marked[v])
708 {
709 all_vertices_marked = false;
710 break;
711 }
712 }
713
714 if (all_vertices_marked)
715 entities.push_back(e);
716 }
717
718 return entities;
719}
720
734template <std::floating_point T, MarkerFn<T> U>
735std::vector<std::int32_t> locate_entities(const Mesh<T>& mesh, int dim,
736 U marker)
737{
738 const int num_entity_types = mesh.topology()->entity_types(dim).size();
739 if (num_entity_types > 1)
740 {
741 throw std::runtime_error(
742 "Multiple entity types of this dimension. Specify entity type index");
743 }
744 return locate_entities(mesh, dim, marker, 0);
745}
746
770template <std::floating_point T, MarkerFn<T> U>
771std::vector<std::int32_t> locate_entities_boundary(const Mesh<T>& mesh, int dim,
772 U marker)
773{
774 // TODO Rewrite this function, it should be possible to simplify considerably
775 auto topology = mesh.topology();
776 assert(topology);
777 int tdim = topology->dim();
778 if (dim == tdim)
779 {
780 throw std::runtime_error(
781 "Cannot use mesh::locate_entities_boundary (boundary) for cells.");
782 }
783
784 // Compute list of boundary facets
785 mesh.topology_mutable()->create_entities(tdim - 1);
786 mesh.topology_mutable()->create_connectivity(tdim - 1, tdim);
787 std::vector<std::int32_t> boundary_facets = exterior_facet_indices(*topology);
788
789 using cmdspan3x_t
790 = md::mdspan<const T, md::extents<std::size_t, 3, md::dynamic_extent>>;
791
792 // Run marker function on the vertex coordinates
793 auto [facet_entities, xdata, vertex_to_pos]
794 = impl::compute_vertex_coords_boundary(mesh, dim, boundary_facets);
795 cmdspan3x_t x(xdata.data(), 3, xdata.size() / 3);
796 std::vector<std::int8_t> marked = marker(x);
797 if (marked.size() != x.extent(1))
798 throw std::runtime_error("Length of array of markers is wrong.");
799
800 // Loop over entities and check vertex markers
801 mesh.topology_mutable()->create_entities(dim);
802 auto e_to_v = topology->connectivity(dim, 0);
803 assert(e_to_v);
804 std::vector<std::int32_t> entities;
805 for (auto e : facet_entities)
806 {
807 // Iterate over entity vertices
808 bool all_vertices_marked = true;
809 for (auto v : e_to_v->links(e))
810 {
811 const std::int32_t pos = vertex_to_pos[v];
812 if (!marked[pos])
813 {
814 all_vertices_marked = false;
815 break;
816 }
817 }
818
819 // Mark facet with all vertices marked
820 if (all_vertices_marked)
821 entities.push_back(e);
822 }
823
824 return entities;
825}
826
845template <std::floating_point T>
846std::pair<std::vector<std::int32_t>, std::array<std::size_t, 2>>
848 std::span<const std::int32_t> entities,
849 bool permute = false)
850{
851 auto topology = mesh.topology();
852 assert(topology);
853 CellType cell_type = topology->cell_type();
854 if (cell_type == CellType::prism and dim == 2)
855 {
856 throw std::runtime_error(
857 "mesh::entities_to_geometry for prism cells not yet supported.");
858 }
859
860 const int tdim = topology->dim();
861 const Geometry<T>& geometry = mesh.geometry();
862 auto xdofs = geometry.dofmap();
863
864 // Get the DOF layout and the number of DOFs per entity
865 const fem::CoordinateElement<T>& coord_ele = geometry.cmap();
866 const fem::ElementDofLayout layout = coord_ele.create_dof_layout();
867 const std::size_t num_entity_dofs = layout.num_entity_closure_dofs(dim);
868 std::vector<std::int32_t> entity_xdofs;
869 entity_xdofs.reserve(entities.size() * num_entity_dofs);
870 std::array<std::size_t, 2> eshape{entities.size(), num_entity_dofs};
871
872 // Get the element's closure DOFs
873 const std::vector<std::vector<std::vector<int>>>& closure_dofs_all
874 = layout.entity_closure_dofs_all();
875
876 // Special case when dim == tdim (cells)
877 if (dim == tdim)
878 {
879 for (std::int32_t c : entities)
880 {
881 // Extract degrees of freedom
882 auto x_c = md::submdspan(xdofs, c, md::full_extent);
883 for (std::int32_t entity_dof : closure_dofs_all[tdim][0])
884 entity_xdofs.push_back(x_c[entity_dof]);
885 }
886
887 return {std::move(entity_xdofs), eshape};
888 }
889
890 assert(dim != tdim);
891
892 auto e_to_c = topology->connectivity(dim, tdim);
893 if (!e_to_c)
894 {
895 throw std::runtime_error(
896 "Entity-to-cell connectivity has not been computed. Missing dims "
897 + std::to_string(dim) + "->" + std::to_string(tdim));
898 }
899
900 auto c_to_e = topology->connectivity(tdim, dim);
901 if (!c_to_e)
902 {
903 throw std::runtime_error(
904 "Cell-to-entity connectivity has not been computed. Missing dims "
905 + std::to_string(tdim) + "->" + std::to_string(dim));
906 }
907
908 // Get the cell info, which is needed to permute the closure dofs
909 std::span<const std::uint32_t> cell_info;
910 if (permute)
911 cell_info = std::span(mesh.topology()->get_cell_permutation_info());
912
913 for (std::int32_t e : entities)
914 {
915 // Get a cell connected to the entity
916 assert(!e_to_c->links(e).empty());
917 std::int32_t c = e_to_c->links(e).front();
918
919 // Get the local index of the entity
920 std::span<const std::int32_t> cell_entities = c_to_e->links(c);
921 auto it = std::find(cell_entities.begin(), cell_entities.end(), e);
922 assert(it != cell_entities.end());
923 std::size_t local_entity = std::distance(cell_entities.begin(), it);
924
925 // Cell sub-entities must be permuted so that their local
926 // orientation agrees with their global orientation
927 std::vector<std::int32_t> closure_dofs(closure_dofs_all[dim][local_entity]);
928 if (permute)
929 {
930 mesh::CellType entity_type
931 = mesh::cell_entity_type(cell_type, dim, local_entity);
932 coord_ele.permute_subentity_closure(closure_dofs, cell_info[c],
933 entity_type, local_entity);
934 }
935
936 // Extract degrees of freedom
937 auto x_c = md::submdspan(xdofs, c, md::full_extent);
938 for (std::int32_t entity_dof : closure_dofs)
939 entity_xdofs.push_back(x_c[entity_dof]);
940 }
941
942 return {std::move(entity_xdofs), eshape};
943}
944
950 = mesh::GhostMode::none,
951 const graph::partition_fn& partfn
953
961std::vector<std::int32_t>
962compute_incident_entities(const Topology& topology,
963 std::span<const std::int32_t> entities, int d0,
964 int d1);
965
1007template <typename U>
1009 MPI_Comm comm, MPI_Comm commt,
1010 std::vector<std::span<const std::int64_t>> cells,
1011 const std::vector<fem::CoordinateElement<
1012 typename std::remove_reference_t<typename U::value_type>>>& elements,
1013 MPI_Comm commg, const U& x, std::array<std::size_t, 2> xshape,
1014 const CellPartitionFunction& partitioner,
1015 const CellReorderFunction& reorder_fn = graph::reorder_gps)
1016{
1017 assert(cells.size() == elements.size());
1018 std::vector<CellType> celltypes;
1019 std::ranges::transform(elements, std::back_inserter(celltypes),
1020 [](auto e) { return e.cell_shape(); });
1021 std::vector<fem::ElementDofLayout> doflayouts;
1022 std::ranges::transform(elements, std::back_inserter(doflayouts),
1023 [](auto e) { return e.create_dof_layout(); });
1024
1025 // Note: `extract_topology` extracts topology data, i.e. just the
1026 // vertices. For P1 geometry this should just be the identity
1027 // operator. For other elements the filtered lists may have 'gaps',
1028 // i.e. the indices might not be contiguous.
1029 //
1030 // `extract_topology` could be skipped for 'P1 geometry' elements
1031
1032 std::int32_t num_cell_types = cells.size();
1033
1034 // -- Partition topology across ranks of comm
1035 std::vector<std::vector<std::int64_t>> cells1(num_cell_types);
1036 std::vector<std::vector<std::int64_t>> original_idx1(num_cell_types);
1037 std::vector<std::vector<int>> ghost_owners(num_cell_types);
1038 if (partitioner)
1039 {
1040 spdlog::info("Using partitioner with cell data ({} cell types)",
1041 num_cell_types);
1043 if (commt != MPI_COMM_NULL)
1044 {
1045 int size = dolfinx::MPI::size(comm);
1046 std::vector<std::vector<std::int64_t>> t(num_cell_types);
1047 std::vector<std::span<const std::int64_t>> tspan(num_cell_types);
1048 for (std::int32_t i = 0; i < num_cell_types; ++i)
1049 {
1050 t[i] = extract_topology(celltypes[i], doflayouts[i], cells[i]);
1051 tspan[i] = std::span(t[i]);
1052 }
1053 dest = partitioner(commt, size, celltypes, tspan);
1054 }
1055
1056 std::int32_t cell_offset = 0;
1057 for (std::int32_t i = 0; i < num_cell_types; ++i)
1058 {
1059 std::size_t num_cell_nodes = doflayouts[i].num_dofs();
1060 assert(cells[i].size() % num_cell_nodes == 0);
1061 std::size_t num_cells = cells[i].size() / num_cell_nodes;
1062
1063 // Extract destination AdjacencyList for this cell type
1064 std::vector<std::int32_t> offsets_i(
1065 std::next(dest.offsets().begin(), cell_offset),
1066 std::next(dest.offsets().begin(), cell_offset + num_cells + 1));
1067 std::vector<std::int32_t> data_i(
1068 std::next(dest.array().begin(), offsets_i.front()),
1069 std::next(dest.array().begin(), offsets_i.back()));
1070 std::int32_t offset_0 = offsets_i.front();
1071 std::ranges::for_each(offsets_i,
1072 [&offset_0](std::int32_t& j) { j -= offset_0; });
1073 graph::AdjacencyList<std::int32_t> dest_i(data_i, offsets_i);
1074 cell_offset += num_cells;
1075
1076 // Distribute cells (topology, includes higher-order 'nodes') to
1077 // destination rank
1078 std::vector<int> src_ranks;
1079 std::tie(cells1[i], src_ranks, original_idx1[i], ghost_owners[i])
1080 = graph::build::distribute(comm, cells[i],
1081 {num_cells, num_cell_nodes}, dest_i);
1082 spdlog::debug("Got {} cells from distribution", cells1[i].size());
1083 }
1084 }
1085 else
1086 {
1087 // No partitioning, construct a global index
1088 std::int64_t num_owned = 0;
1089 for (std::int32_t i = 0; i < num_cell_types; ++i)
1090 {
1091 cells1[i] = std::vector<std::int64_t>(cells[i].begin(), cells[i].end());
1092 std::int32_t num_cell_nodes = doflayouts[i].num_dofs();
1093 assert(cells1[i].size() % num_cell_nodes == 0);
1094 original_idx1[i].resize(cells1[i].size() / num_cell_nodes);
1095 num_owned += original_idx1[i].size();
1096 }
1097
1098 // Add on global offset
1099 std::int64_t global_offset = 0;
1100 MPI_Exscan(&num_owned, &global_offset, 1, MPI_INT64_T, MPI_SUM, comm);
1101 for (std::int32_t i = 0; i < num_cell_types; ++i)
1102 {
1103 std::iota(original_idx1[i].begin(), original_idx1[i].end(),
1104 global_offset);
1105 global_offset += original_idx1[i].size();
1106 }
1107 }
1108
1109 // Extract cell 'topology', i.e. extract the vertices for each cell
1110 // and discard any 'higher-order' nodes
1111 std::vector<std::vector<std::int64_t>> cells1_v(num_cell_types);
1112 for (std::int32_t i = 0; i < num_cell_types; ++i)
1113 {
1114 cells1_v[i] = extract_topology(celltypes[i], doflayouts[i], cells1[i]);
1115 spdlog::info("Extract basic topology: {}->{}", cells1[i].size(),
1116 cells1_v[i].size());
1117 }
1118
1119 auto boundary_v_fn = create_boundary_vertices_fn(reorder_fn);
1120 const std::vector<std::int64_t> boundary_v = boundary_v_fn(
1121 celltypes, doflayouts, ghost_owners, cells1, cells1_v, original_idx1);
1122
1123 spdlog::debug("Got {} boundary vertices", boundary_v.size());
1124
1125 // Create Topology
1126 std::vector<std::span<const std::int64_t>> cells1_v_span;
1127 std::ranges::transform(cells1_v, std::back_inserter(cells1_v_span),
1128 [](auto& c) { return std::span(c); });
1129 std::vector<std::span<const std::int64_t>> original_idx1_span;
1130 std::ranges::transform(original_idx1, std::back_inserter(original_idx1_span),
1131 [](auto& c) { return std::span(c); });
1132 std::vector<std::span<const int>> ghost_owners_span;
1133 std::ranges::transform(ghost_owners, std::back_inserter(ghost_owners_span),
1134 [](auto& c) { return std::span(c); });
1135 Topology topology
1136 = create_topology(comm, celltypes, cells1_v_span, original_idx1_span,
1137 ghost_owners_span, boundary_v);
1138
1139 // Create connectivities required higher-order geometries for creating
1140 // a Geometry object
1141 for (int i = 0; i < num_cell_types; ++i)
1142 {
1143 for (int e = 1; e < topology.dim(); ++e)
1144 {
1145 if (doflayouts[i].num_entity_dofs(e) > 0)
1146 topology.create_entities(e);
1147 }
1148
1149 if (elements[i].needs_dof_permutations())
1150 topology.create_entity_permutations();
1151 }
1152
1153 // Build list of unique (global) node indices from cells1 and
1154 // distribute coordinate data
1155 std::vector<std::int64_t> nodes1, nodes2;
1156 for (std::vector<std::int64_t>& c : cells1)
1157 nodes1.insert(nodes1.end(), c.begin(), c.end());
1158 for (std::vector<std::int64_t>& c : cells1)
1159 nodes2.insert(nodes2.end(), c.begin(), c.end());
1160
1161 dolfinx::radix_sort(nodes1);
1162 auto [unique_end, range_end] = std::ranges::unique(nodes1);
1163 nodes1.erase(unique_end, range_end);
1164
1165 std::vector coords
1166 = dolfinx::MPI::distribute_data(comm, nodes1, commg, x, xshape[1]);
1167
1168 // Create geometry object
1170 = create_geometry(topology, elements, nodes1, nodes2, coords, xshape[1]);
1171
1172 return Mesh(comm, std::make_shared<Topology>(std::move(topology)),
1173 std::move(geometry));
1174}
1175
1211template <typename U>
1213 MPI_Comm comm, MPI_Comm commt, std::span<const std::int64_t> cells,
1215 typename std::remove_reference_t<typename U::value_type>>& element,
1216 MPI_Comm commg, const U& x, std::array<std::size_t, 2> xshape,
1217 const CellPartitionFunction& partitioner,
1218 const CellReorderFunction& reorder_fn = graph::reorder_gps)
1219{
1220 return create_mesh(comm, commt, std::vector{cells}, std::vector{element},
1221 commg, x, xshape, partitioner, reorder_fn);
1222}
1223
1242template <typename U>
1243Mesh<typename std::remove_reference_t<typename U::value_type>>
1244create_mesh(MPI_Comm comm, std::span<const std::int64_t> cells,
1246 std::remove_reference_t<typename U::value_type>>& elements,
1247 const U& x, std::array<std::size_t, 2> xshape, GhostMode ghost_mode)
1248{
1249 if (dolfinx::MPI::size(comm) == 1)
1250 return create_mesh(comm, comm, std::vector{cells}, std::vector{elements},
1251 comm, x, xshape, nullptr);
1252 else
1253 {
1254 return create_mesh(comm, comm, std::vector{cells}, std::vector{elements},
1255 comm, x, xshape, create_cell_partitioner(ghost_mode));
1256 }
1257}
1258
1272template <std::floating_point T>
1273std::pair<Geometry<T>, std::vector<int32_t>>
1275 std::span<const std::int32_t> subentity_to_entity)
1276{
1277 const Geometry<T>& geometry = mesh.geometry();
1278
1279 // Get the geometry dofs in the sub-geometry based on the entities in
1280 // sub-geometry
1281 const fem::ElementDofLayout layout = geometry.cmap().create_dof_layout();
1282
1283 const std::vector<std::int32_t> x_indices
1284 = entities_to_geometry(mesh, dim, subentity_to_entity, true).first;
1285
1286 std::vector<std::int32_t> sub_x_dofs = x_indices;
1287 std::ranges::sort(sub_x_dofs);
1288 auto [unique_end, range_end] = std::ranges::unique(sub_x_dofs);
1289 sub_x_dofs.erase(unique_end, range_end);
1290
1291 // Get the sub-geometry dofs owned by this process
1292 auto x_index_map = geometry.index_map();
1293 assert(x_index_map);
1294
1295 std::shared_ptr<common::IndexMap> sub_x_dof_index_map;
1296 std::vector<std::int32_t> subx_to_x_dofmap;
1297 {
1298 auto [map, new_to_old] = common::create_sub_index_map(
1299 *x_index_map, sub_x_dofs, common::IndexMapOrder::any, true);
1300 sub_x_dof_index_map = std::make_shared<common::IndexMap>(std::move(map));
1301 subx_to_x_dofmap = std::move(new_to_old);
1302 }
1303
1304 // Create sub-geometry coordinates
1305 std::span<const T> x = geometry.x();
1306 std::int32_t sub_num_x_dofs = subx_to_x_dofmap.size();
1307 std::vector<T> sub_x(3 * sub_num_x_dofs);
1308 for (std::int32_t i = 0; i < sub_num_x_dofs; ++i)
1309 {
1310 std::copy_n(std::next(x.begin(), 3 * subx_to_x_dofmap[i]), 3,
1311 std::next(sub_x.begin(), 3 * i));
1312 }
1313
1314 // Create geometry to sub-geometry map
1315 std::vector<std::int32_t> x_to_subx_dof_map(
1316 x_index_map->size_local() + x_index_map->num_ghosts(), -1);
1317 for (std::size_t i = 0; i < subx_to_x_dofmap.size(); ++i)
1318 x_to_subx_dof_map[subx_to_x_dofmap[i]] = i;
1319
1320 // Create sub-geometry dofmap
1321 std::vector<std::int32_t> sub_x_dofmap;
1322 sub_x_dofmap.reserve(x_indices.size());
1323 std::ranges::transform(x_indices, std::back_inserter(sub_x_dofmap),
1324 [&x_to_subx_dof_map](auto x_dof)
1325 {
1326 assert(x_to_subx_dof_map[x_dof] != -1);
1327 return x_to_subx_dof_map[x_dof];
1328 });
1329
1330 // Sub-geometry coordinate element
1331 CellType sub_xcell = cell_entity_type(geometry.cmap().cell_shape(), dim, 0);
1332
1333 // Special handling of point meshes, as they only support constant
1334 // basis functions
1335 int degree = (sub_xcell == CellType::point) ? 0 : geometry.cmap().degree();
1336 fem::CoordinateElement<T> sub_cmap(sub_xcell, degree,
1337 geometry.cmap().variant());
1338
1339 // Sub-geometry input_global_indices
1340 const std::vector<std::int64_t>& igi = geometry.input_global_indices();
1341 std::vector<std::int64_t> sub_igi;
1342 sub_igi.reserve(subx_to_x_dofmap.size());
1343 std::ranges::transform(subx_to_x_dofmap, std::back_inserter(sub_igi),
1344 [&igi](auto sub_x_dof) { return igi[sub_x_dof]; });
1345
1346 // Create geometry
1347 return {Geometry(
1348 sub_x_dof_index_map,
1349 std::vector<std::vector<std::int32_t>>{std::move(sub_x_dofmap)},
1350 {sub_cmap}, std::move(sub_x), geometry.dim(), std::move(sub_igi)),
1351 std::move(subx_to_x_dofmap)};
1352}
1353
1363template <std::floating_point T>
1364std::tuple<Mesh<T>, std::vector<std::int32_t>, std::vector<std::int32_t>,
1365 std::vector<std::int32_t>>
1367 std::span<const std::int32_t> entities)
1368{
1369 // Create sub-topology
1370 mesh.topology_mutable()->create_connectivity(dim, 0);
1371 auto [topology, subentity_to_entity, subvertex_to_vertex]
1372 = mesh::create_subtopology(*mesh.topology(), dim, entities);
1373
1374 // Create sub-geometry
1375 const int tdim = mesh.topology()->dim();
1376 mesh.topology_mutable()->create_entities(dim);
1377 mesh.topology_mutable()->create_connectivity(dim, tdim);
1378 mesh.topology_mutable()->create_connectivity(tdim, dim);
1379 mesh.topology_mutable()->create_entity_permutations();
1380 auto [geometry, subx_to_x_dofmap]
1381 = mesh::create_subgeometry(mesh, dim, subentity_to_entity);
1382
1383 return {Mesh(mesh.comm(), std::make_shared<Topology>(std::move(topology)),
1384 std::move(geometry)),
1385 std::move(subentity_to_entity), std::move(subvertex_to_vertex),
1386 std::move(subx_to_x_dofmap)};
1387}
1388
1389} // 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:30
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
Definition AdjacencyList.h:27
const std::vector< T > & array() const
Return contiguous array of links for all nodes (const version)
Definition AdjacencyList.h:128
const std::vector< std::int32_t > & offsets() const
Offset for each node in array() (const version)
Definition AdjacencyList.h:134
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:46
Requirements on function for geometry marking.
Definition utils.h:652
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:52
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:83
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:606
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:680
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:815
@ any
Allow arbitrary ordering of ghost indices in sub-maps.
Definition IndexMap.h:27
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 dofmapbuilder.h:26
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:360
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
GhostMode
Enum for different partitioning ghost modes.
Definition utils.h:38
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:1017
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, std::function< std::vector< int >(const graph::AdjacencyList< std::int32_t > &)> reorder_fn=nullptr)
Build Geometry from input data.
Definition Geometry.h:236
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:451
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:215
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:1303
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:58
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:771
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:399
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:1274
auto create_boundary_vertices_fn(const CellReorderFunction &reorder_fn)
Creates the default boundary vertices routine for a given reorder function.
Definition utils.h:224
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)
Compute the local part of the dual graph (cell-cell connections via facets) and facets with only one ...
Definition graphbuild.cpp:354
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:847
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:123
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:29
CellType
Cell type identifier.
Definition cell_types.h:22
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:673
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:207
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:568
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, 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:1008
std::tuple< Mesh< T >, std::vector< std::int32_t >, std::vector< std::int32_t >, 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:1366
CellPartitionFunction create_cell_partitioner(mesh::GhostMode ghost_mode=mesh::GhostMode::none, const graph::partition_fn &partfn=&graph::partition_graph)
Create a function that computes destination rank for mesh cells on this rank by applying the default ...
Definition utils.cpp:100
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