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-2020 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 <span>
21
24
25namespace dolfinx::fem
26{
27class ElementDofLayout;
28}
29
30namespace dolfinx::mesh
31{
32enum class CellType;
33
35enum class GhostMode : int
36{
37 none,
38 shared_facet,
39 shared_vertex
40};
41
42namespace impl
43{
45template <typename T>
46void reorder_list(std::span<T> list, std::span<const std::int32_t> nodemap)
47{
48 if (nodemap.empty())
49 return;
50
51 assert(list.size() % nodemap.size() == 0);
52 int degree = list.size() / nodemap.size();
53 const std::vector<T> orig(list.begin(), list.end());
54 for (std::size_t n = 0; n < nodemap.size(); ++n)
55 {
56 auto links_old = std::span(orig.data() + n * degree, degree);
57 auto links_new = list.subspan(nodemap[n] * degree, degree);
58 std::ranges::copy(links_old, links_new.begin());
59 }
60}
61
75template <std::floating_point T>
76std::tuple<std::vector<std::int32_t>, std::vector<T>, std::vector<std::int32_t>>
78 std::span<const std::int32_t> facets)
79{
80 auto topology = mesh.topology();
81 assert(topology);
82 const int tdim = topology->dim();
83 if (dim == tdim)
84 {
85 throw std::runtime_error(
86 "Cannot use mesh::locate_entities_boundary (boundary) for cells.");
87 }
88
89 // Build set of vertices on boundary and set of boundary entities
90 mesh.topology_mutable()->create_connectivity(tdim - 1, 0);
91 mesh.topology_mutable()->create_connectivity(tdim - 1, dim);
92 std::vector<std::int32_t> vertices, entities;
93 {
94 auto f_to_v = topology->connectivity(tdim - 1, 0);
95 assert(f_to_v);
96 auto f_to_e = topology->connectivity(tdim - 1, dim);
97 assert(f_to_e);
98 for (auto f : facets)
99 {
100 auto v = f_to_v->links(f);
101 vertices.insert(vertices.end(), v.begin(), v.end());
102 auto e = f_to_e->links(f);
103 entities.insert(entities.end(), e.begin(), e.end());
104 }
105
106 // Build vector of boundary vertices
107 {
108 std::ranges::sort(vertices);
109 auto [unique_end, range_end] = std::ranges::unique(vertices);
110 vertices.erase(unique_end, range_end);
111 }
112
113 {
114 std::ranges::sort(entities);
115 auto [unique_end, range_end] = std::ranges::unique(entities);
116 entities.erase(unique_end, range_end);
117 }
118 }
119
120 // Get geometry data
121 auto x_dofmap = mesh.geometry().dofmap();
122 std::span<const T> x_nodes = mesh.geometry().x();
123
124 // Get all vertex 'node' indices
125 mesh.topology_mutable()->create_connectivity(0, tdim);
126 mesh.topology_mutable()->create_connectivity(tdim, 0);
127 auto v_to_c = topology->connectivity(0, tdim);
128 assert(v_to_c);
129 auto c_to_v = topology->connectivity(tdim, 0);
130 assert(c_to_v);
131 std::vector<T> x_vertices(3 * vertices.size(), -1.0);
132 std::vector<std::int32_t> vertex_to_pos(v_to_c->num_nodes(), -1);
133 for (std::size_t i = 0; i < vertices.size(); ++i)
134 {
135 const std::int32_t v = vertices[i];
136
137 // Get first cell and find position
138 const int c = v_to_c->links(v).front();
139 auto cell_vertices = c_to_v->links(c);
140 auto it = std::find(cell_vertices.begin(), cell_vertices.end(), v);
141 assert(it != cell_vertices.end());
142 const int local_pos = std::distance(cell_vertices.begin(), it);
143
144 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
145 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
146 for (std::size_t j = 0; j < 3; ++j)
147 x_vertices[j * vertices.size() + i] = x_nodes[3 * dofs[local_pos] + j];
148 vertex_to_pos[v] = i;
149 }
150
151 return {std::move(entities), std::move(x_vertices), std::move(vertex_to_pos)};
152}
153
154} // namespace impl
155
167std::vector<std::int32_t> exterior_facet_indices(const Topology& topology);
168
185using CellPartitionFunction = std::function<graph::AdjacencyList<std::int32_t>(
186 MPI_Comm comm, int nparts, const std::vector<CellType>& cell_types,
187 const std::vector<std::span<const std::int64_t>>& cells)>;
188
198std::vector<std::int64_t> extract_topology(CellType cell_type,
199 const fem::ElementDofLayout& layout,
200 std::span<const std::int64_t> cells);
201
210template <std::floating_point T>
211std::vector<T> h(const Mesh<T>& mesh, std::span<const std::int32_t> entities,
212 int dim)
213{
214 if (entities.empty())
215 return std::vector<T>();
216 if (dim == 0)
217 return std::vector<T>(entities.size(), 0);
218
219 // Get the geometry dofs for the vertices of each entity
220 const std::vector<std::int32_t> vertex_xdofs
221 = entities_to_geometry(mesh, dim, entities, false);
222 assert(!entities.empty());
223 const std::size_t num_vertices = vertex_xdofs.size() / entities.size();
224
225 // Get the geometry coordinate
226 std::span<const T> x = mesh.geometry().x();
227
228 // Function to compute the length of (p0 - p1)
229 auto delta_norm = [](auto&& p0, auto&& p1)
230 {
231 T norm = 0;
232 for (std::size_t i = 0; i < 3; ++i)
233 norm += (p0[i] - p1[i]) * (p0[i] - p1[i]);
234 return std::sqrt(norm);
235 };
236
237 // Compute greatest distance between any to vertices
238 assert(dim > 0);
239 std::vector<T> h(entities.size(), 0);
240 for (std::size_t e = 0; e < entities.size(); ++e)
241 {
242 // Get geometry 'dof' for each vertex of entity e
243 std::span<const std::int32_t> e_vertices(
244 vertex_xdofs.data() + e * num_vertices, num_vertices);
245
246 // Compute maximum distance between any two vertices
247 for (std::size_t i = 0; i < e_vertices.size(); ++i)
248 {
249 std::span<const T, 3> p0(x.data() + 3 * e_vertices[i], 3);
250 for (std::size_t j = i + 1; j < e_vertices.size(); ++j)
251 {
252 std::span<const T, 3> p1(x.data() + 3 * e_vertices[j], 3);
253 h[e] = std::max(h[e], delta_norm(p0, p1));
254 }
255 }
256 }
257
258 return h;
259}
260
264template <std::floating_point T>
265std::vector<T> cell_normals(const Mesh<T>& mesh, int dim,
266 std::span<const std::int32_t> entities)
267{
268 auto topology = mesh.topology();
269 assert(topology);
270
271 if (entities.empty())
272 return std::vector<T>();
273
274 if (topology->cell_type() == CellType::prism and dim == 2)
275 throw std::runtime_error("More work needed for prism cell");
276
277 const int gdim = mesh.geometry().dim();
278 const CellType type = cell_entity_type(topology->cell_type(), dim, 0);
279
280 // Find geometry nodes for topology entities
281 std::span<const T> x = mesh.geometry().x();
282 std::vector<std::int32_t> geometry_entities
283 = entities_to_geometry(mesh, dim, entities, false);
284
285 const std::size_t shape1 = geometry_entities.size() / entities.size();
286 std::vector<T> n(entities.size() * 3);
287 switch (type)
288 {
289 case CellType::interval:
290 {
291 if (gdim > 2)
292 throw std::invalid_argument("Interval cell normal undefined in 3D");
293 for (std::size_t i = 0; i < entities.size(); ++i)
294 {
295 // Get the two vertices as points
296 std::array vertices{geometry_entities[i * shape1],
297 geometry_entities[i * shape1 + 1]};
298 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
299 std::span<const T, 3>(x.data() + 3 * vertices[1], 3)};
300
301 // Define normal by rotating tangent counter-clockwise
302 std::array<T, 3> t;
303 std::ranges::transform(p[1], p[0], t.begin(),
304 [](auto x, auto y) { return x - y; });
305
306 T norm = std::sqrt(t[0] * t[0] + t[1] * t[1]);
307 std::span<T, 3> ni(n.data() + 3 * i, 3);
308 ni[0] = -t[1] / norm;
309 ni[1] = t[0] / norm;
310 ni[2] = 0.0;
311 }
312 return n;
313 }
314 case CellType::triangle:
315 {
316 for (std::size_t i = 0; i < entities.size(); ++i)
317 {
318 // Get the three vertices as points
319 std::array vertices = {geometry_entities[i * shape1 + 0],
320 geometry_entities[i * shape1 + 1],
321 geometry_entities[i * shape1 + 2]};
322 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
323 std::span<const T, 3>(x.data() + 3 * vertices[1], 3),
324 std::span<const T, 3>(x.data() + 3 * vertices[2], 3)};
325
326 // Compute (p1 - p0) and (p2 - p0)
327 std::array<T, 3> dp1, dp2;
328 std::ranges::transform(p[1], p[0], dp1.begin(),
329 [](auto x, auto y) { return x - y; });
330 std::ranges::transform(p[2], p[0], dp2.begin(),
331 [](auto x, auto y) { return x - y; });
332
333 // Define cell normal via cross product of first two edges
334 std::array<T, 3> ni = math::cross(dp1, dp2);
335 T norm = std::sqrt(ni[0] * ni[0] + ni[1] * ni[1] + ni[2] * ni[2]);
336 std::ranges::transform(ni, std::next(n.begin(), 3 * i),
337 [norm](auto x) { return x / norm; });
338 }
339
340 return n;
341 }
342 case CellType::quadrilateral:
343 {
344 // TODO: check
345 for (std::size_t i = 0; i < entities.size(); ++i)
346 {
347 // Get the three vertices as points
348 std::array vertices = {geometry_entities[i * shape1 + 0],
349 geometry_entities[i * shape1 + 1],
350 geometry_entities[i * shape1 + 2]};
351 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
352 std::span<const T, 3>(x.data() + 3 * vertices[1], 3),
353 std::span<const T, 3>(x.data() + 3 * vertices[2], 3)};
354
355 // Compute (p1 - p0) and (p2 - p0)
356 std::array<T, 3> dp1, dp2;
357 std::ranges::transform(p[1], p[0], dp1.begin(),
358 [](auto x, auto y) { return x - y; });
359 std::ranges::transform(p[2], p[0], dp2.begin(),
360 [](auto x, auto y) { return x - y; });
361
362 // Define cell normal via cross product of first two edges
363 std::array<T, 3> ni = math::cross(dp1, dp2);
364 T norm = std::sqrt(ni[0] * ni[0] + ni[1] * ni[1] + ni[2] * ni[2]);
365 std::ranges::transform(ni, std::next(n.begin(), 3 * i),
366 [norm](auto x) { return x / norm; });
367 }
368
369 return n;
370 }
371 default:
372 throw std::invalid_argument(
373 "cell_normal not supported for this cell type.");
374 }
375}
376
380template <std::floating_point T>
381std::vector<T> compute_midpoints(const Mesh<T>& mesh, int dim,
382 std::span<const std::int32_t> entities)
383{
384 if (entities.empty())
385 return std::vector<T>();
386
387 std::span<const T> x = mesh.geometry().x();
388
389 // Build map from entity -> geometry dof
390 const std::vector<std::int32_t> e_to_g
391 = entities_to_geometry(mesh, dim, entities, false);
392 std::size_t shape1 = e_to_g.size() / entities.size();
393
394 std::vector<T> x_mid(entities.size() * 3, 0);
395 for (std::size_t e = 0; e < entities.size(); ++e)
396 {
397 std::span<T, 3> p(x_mid.data() + 3 * e, 3);
398 std::span<const std::int32_t> rows(e_to_g.data() + e * shape1, shape1);
399 for (auto row : rows)
400 {
401 std::span<const T, 3> xg(x.data() + 3 * row, 3);
402 std::ranges::transform(p, xg, p.begin(),
403 [size = rows.size()](auto x, auto y)
404 { return x + y / size; });
405 }
406 }
407
408 return x_mid;
409}
410
411namespace impl
412{
417template <std::floating_point T>
418std::pair<std::vector<T>, std::array<std::size_t, 2>>
420{
421 auto topology = mesh.topology();
422 assert(topology);
423 const int tdim = topology->dim();
424
425 // Create entities and connectivities
426 mesh.topology_mutable()->create_connectivity(tdim, 0);
427
428 // Get all vertex 'node' indices
429 auto x_dofmap = mesh.geometry().dofmap();
430 const std::int32_t num_vertices = topology->index_map(0)->size_local()
431 + topology->index_map(0)->num_ghosts();
432 auto c_to_v = topology->connectivity(tdim, 0);
433 assert(c_to_v);
434 std::vector<std::int32_t> vertex_to_node(num_vertices);
435 for (int c = 0; c < c_to_v->num_nodes(); ++c)
436 {
437 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
438 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
439 auto vertices = c_to_v->links(c);
440 for (std::size_t i = 0; i < vertices.size(); ++i)
441 vertex_to_node[vertices[i]] = x_dofs[i];
442 }
443
444 // Pack coordinates of vertices
445 std::span<const T> x_nodes = mesh.geometry().x();
446 std::vector<T> x_vertices(3 * vertex_to_node.size(), 0.0);
447 for (std::size_t i = 0; i < vertex_to_node.size(); ++i)
448 {
449 const int pos = 3 * vertex_to_node[i];
450 for (std::size_t j = 0; j < 3; ++j)
451 x_vertices[j * vertex_to_node.size() + i] = x_nodes[pos + j];
452 }
453
454 return {std::move(x_vertices), {3, vertex_to_node.size()}};
455}
456
457} // namespace impl
458
460template <typename Fn, typename T>
461concept MarkerFn = std::is_invocable_r<
462 std::vector<std::int8_t>, Fn,
463 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
464 const T, MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
465 std::size_t, 3,
466 MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>>::value;
467
481template <std::floating_point T, MarkerFn<T> U>
482std::vector<std::int32_t> locate_entities(const Mesh<T>& mesh, int dim,
483 U marker)
484{
485 using cmdspan3x_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
486 const T,
487 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
488 std::size_t, 3, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>;
489
490 // Run marker function on vertex coordinates
491 const auto [xdata, xshape] = impl::compute_vertex_coords(mesh);
492 cmdspan3x_t x(xdata.data(), xshape);
493 const std::vector<std::int8_t> marked = marker(x);
494 if (marked.size() != x.extent(1))
495 throw std::runtime_error("Length of array of markers is wrong.");
496
497 auto topology = mesh.topology();
498 assert(topology);
499 const int tdim = topology->dim();
500
501 mesh.topology_mutable()->create_entities(dim);
502 mesh.topology_mutable()->create_connectivity(tdim, 0);
503 if (dim < tdim)
504 mesh.topology_mutable()->create_connectivity(dim, 0);
505
506 // Iterate over entities of dimension 'dim' to build vector of marked
507 // entities
508 auto e_to_v = topology->connectivity(dim, 0);
509 assert(e_to_v);
510 std::vector<std::int32_t> entities;
511 for (int e = 0; e < e_to_v->num_nodes(); ++e)
512 {
513 // Iterate over entity vertices
514 bool all_vertices_marked = true;
515 for (std::int32_t v : e_to_v->links(e))
516 {
517 if (!marked[v])
518 {
519 all_vertices_marked = false;
520 break;
521 }
522 }
523
524 if (all_vertices_marked)
525 entities.push_back(e);
526 }
527
528 return entities;
529}
530
554template <std::floating_point T, MarkerFn<T> U>
555std::vector<std::int32_t> locate_entities_boundary(const Mesh<T>& mesh, int dim,
556 U marker)
557{
558 auto topology = mesh.topology();
559 assert(topology);
560 int tdim = topology->dim();
561 if (dim == tdim)
562 {
563 throw std::runtime_error(
564 "Cannot use mesh::locate_entities_boundary (boundary) for cells.");
565 }
566
567 // Compute list of boundary facets
568 mesh.topology_mutable()->create_entities(tdim - 1);
569 mesh.topology_mutable()->create_connectivity(tdim - 1, tdim);
570 std::vector<std::int32_t> boundary_facets = exterior_facet_indices(*topology);
571
572 using cmdspan3x_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
573 const T,
574 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
575 std::size_t, 3, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>;
576
577 // Run marker function on the vertex coordinates
578 auto [facet_entities, xdata, vertex_to_pos]
579 = impl::compute_vertex_coords_boundary(mesh, dim, boundary_facets);
580 cmdspan3x_t x(xdata.data(), 3, xdata.size() / 3);
581 std::vector<std::int8_t> marked = marker(x);
582 if (marked.size() != x.extent(1))
583 throw std::runtime_error("Length of array of markers is wrong.");
584
585 // Loop over entities and check vertex markers
586 mesh.topology_mutable()->create_entities(dim);
587 auto e_to_v = topology->connectivity(dim, 0);
588 assert(e_to_v);
589 std::vector<std::int32_t> entities;
590 for (auto e : facet_entities)
591 {
592 // Iterate over entity vertices
593 bool all_vertices_marked = true;
594 for (auto v : e_to_v->links(e))
595 {
596 const std::int32_t pos = vertex_to_pos[v];
597 if (!marked[pos])
598 {
599 all_vertices_marked = false;
600 break;
601 }
602 }
603
604 // Mark facet with all vertices marked
605 if (all_vertices_marked)
606 entities.push_back(e);
607 }
608
609 return entities;
610}
611
630template <std::floating_point T>
631std::vector<std::int32_t>
632entities_to_geometry(const Mesh<T>& mesh, int dim,
633 std::span<const std::int32_t> entities,
634 bool permute = false)
635{
636 auto topology = mesh.topology();
637 assert(topology);
638 CellType cell_type = topology->cell_type();
639 if (cell_type == CellType::prism and dim == 2)
640 throw std::runtime_error("More work needed for prism cells");
641
642 const int tdim = topology->dim();
643 const Geometry<T>& geometry = mesh.geometry();
644 auto xdofs = geometry.dofmap();
645
646 // Get the DOF layout and the number of DOFs per entity
647 const fem::CoordinateElement<T>& coord_ele = geometry.cmap();
648 const fem::ElementDofLayout layout = coord_ele.create_dof_layout();
649 const std::size_t num_entity_dofs = layout.num_entity_closure_dofs(dim);
650 std::vector<std::int32_t> entity_xdofs;
651 entity_xdofs.reserve(entities.size() * num_entity_dofs);
652
653 // Get the element's closure DOFs
654 const std::vector<std::vector<std::vector<int>>>& closure_dofs_all
655 = layout.entity_closure_dofs_all();
656
657 // Special case when dim == tdim (cells)
658 if (dim == tdim)
659 {
660 for (std::size_t i = 0; i < entities.size(); ++i)
661 {
662 const std::int32_t c = entities[i];
663 // Extract degrees of freedom
664 auto x_c = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
665 xdofs, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
666 for (std::int32_t entity_dof : closure_dofs_all[tdim][0])
667 entity_xdofs.push_back(x_c[entity_dof]);
668 }
669 return entity_xdofs;
670 }
671
672 assert(dim != tdim);
673
674 auto e_to_c = topology->connectivity(dim, tdim);
675 if (!e_to_c)
676 {
677 throw std::runtime_error(
678 "Entity-to-cell connectivity has not been computed. Missing dims "
679 + std::to_string(dim) + "->" + std::to_string(tdim));
680 }
681
682 auto c_to_e = topology->connectivity(tdim, dim);
683 if (!c_to_e)
684 {
685 throw std::runtime_error(
686 "Cell-to-entity connectivity has not been computed. Missing dims "
687 + std::to_string(tdim) + "->" + std::to_string(dim));
688 }
689
690 // Get the cell info, which is needed to permute the closure dofs
691 std::span<const std::uint32_t> cell_info;
692 if (permute)
693 cell_info = std::span(mesh.topology()->get_cell_permutation_info());
694
695 for (std::size_t i = 0; i < entities.size(); ++i)
696 {
697 const std::int32_t e = entities[i];
698
699 // Get a cell connected to the entity
700 assert(!e_to_c->links(e).empty());
701 std::int32_t c = e_to_c->links(e).front();
702
703 // Get the local index of the entity
704 std::span<const std::int32_t> cell_entities = c_to_e->links(c);
705 auto it = std::find(cell_entities.begin(), cell_entities.end(), e);
706 assert(it != cell_entities.end());
707 std::size_t local_entity = std::distance(cell_entities.begin(), it);
708
709 std::vector<std::int32_t> closure_dofs(closure_dofs_all[dim][local_entity]);
710
711 // Cell sub-entities must be permuted so that their local orientation agrees
712 // with their global orientation
713 if (permute)
714 {
715 mesh::CellType entity_type
716 = mesh::cell_entity_type(cell_type, dim, local_entity);
717 coord_ele.permute_subentity_closure(closure_dofs, cell_info[c],
718 entity_type, local_entity);
719 }
720
721 // Extract degrees of freedom
722 auto x_c = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
723 xdofs, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
724 for (std::int32_t entity_dof : closure_dofs)
725 entity_xdofs.push_back(x_c[entity_dof]);
726 }
727
728 return entity_xdofs;
729}
730
736 = mesh::GhostMode::none,
737 const graph::partition_fn& partfn
739
747std::vector<std::int32_t>
748compute_incident_entities(const Topology& topology,
749 std::span<const std::int32_t> entities, int d0,
750 int d1);
751
781template <typename U>
783 MPI_Comm comm, MPI_Comm commt, std::span<const std::int64_t> cells,
785 typename std::remove_reference_t<typename U::value_type>>& element,
786 MPI_Comm commg, const U& x, std::array<std::size_t, 2> xshape,
787 const CellPartitionFunction& partitioner)
788{
789 CellType celltype = element.cell_shape();
790 const fem::ElementDofLayout doflayout = element.create_dof_layout();
791
792 const int num_cell_vertices = mesh::num_cell_vertices(element.cell_shape());
793 std::size_t num_cell_nodes = doflayout.num_dofs();
794
795 // Note: `extract_topology` extracts topology data, i.e. just the
796 // vertices. For P1 geometry this should just be the identity
797 // operator. For other elements the filtered lists may have 'gaps',
798 // i.e. the indices might not be contiguous.
799 //
800 // `extract_topology` could be skipped for 'P1 geometry' elements
801
802 // -- Partition topology across ranks of comm
803 std::vector<std::int64_t> cells1;
804 std::vector<std::int64_t> original_idx1;
805 std::vector<int> ghost_owners;
806 if (partitioner)
807 {
808 spdlog::info("Using partitioner with {} cell data", cells.size());
810 if (commt != MPI_COMM_NULL)
811 {
812 int size = dolfinx::MPI::size(comm);
813 auto t = extract_topology(element.cell_shape(), doflayout, cells);
814 dest = partitioner(commt, size, {celltype}, {t});
815 }
816
817 // Distribute cells (topology, includes higher-order 'nodes') to
818 // destination rank
819 assert(cells.size() % num_cell_nodes == 0);
820 std::size_t num_cells = cells.size() / num_cell_nodes;
821 std::vector<int> src_ranks;
822 std::tie(cells1, src_ranks, original_idx1, ghost_owners)
823 = graph::build::distribute(comm, cells, {num_cells, num_cell_nodes},
824 dest);
825 spdlog::debug("Got {} cells from distribution", cells1.size());
826 }
827 else
828 {
829 cells1 = std::vector<std::int64_t>(cells.begin(), cells.end());
830 assert(cells1.size() % num_cell_nodes == 0);
831 std::int64_t offset = 0;
832 std::int64_t num_owned = cells1.size() / num_cell_nodes;
833 MPI_Exscan(&num_owned, &offset, 1, MPI_INT64_T, MPI_SUM, comm);
834 original_idx1.resize(num_owned);
835 std::iota(original_idx1.begin(), original_idx1.end(), offset);
836 }
837
838 // Extract cell 'topology', i.e. extract the vertices for each cell
839 // and discard any 'higher-order' nodes
840 std::vector<std::int64_t> cells1_v
841 = extract_topology(celltype, doflayout, cells1);
842 spdlog::info("Extract basic topology: {}->{}", cells1.size(),
843 cells1_v.size());
844
845 // Build local dual graph for owned cells to (i) get list of vertices
846 // on the process boundary and (ii) apply re-ordering to cells for
847 // locality
848 std::vector<std::int64_t> boundary_v;
849 {
850 std::int32_t num_owned_cells
851 = cells1_v.size() / num_cell_vertices - ghost_owners.size();
852 std::vector<std::int32_t> cell_offsets(num_owned_cells + 1, 0);
853 for (std::size_t i = 1; i < cell_offsets.size(); ++i)
854 cell_offsets[i] = cell_offsets[i - 1] + num_cell_vertices;
855 spdlog::info("Build local dual graph");
856 auto [graph, unmatched_facets, max_v, facet_attached_cells]
858 std::vector{celltype},
859 {std::span(cells1_v.data(), num_owned_cells * num_cell_vertices)});
860 const std::vector<int> remap = graph::reorder_gps(graph);
861
862 // Create re-ordered cell lists (leaves ghosts unchanged)
863 std::vector<std::int64_t> _original_idx(original_idx1.size());
864 for (std::size_t i = 0; i < remap.size(); ++i)
865 _original_idx[remap[i]] = original_idx1[i];
866 std::copy_n(std::next(original_idx1.cbegin(), num_owned_cells),
867 ghost_owners.size(),
868 std::next(_original_idx.begin(), num_owned_cells));
869 impl::reorder_list(
870 std::span(cells1_v.data(), remap.size() * num_cell_vertices), remap);
871 impl::reorder_list(std::span(cells1.data(), remap.size() * num_cell_nodes),
872 remap);
873 original_idx1 = _original_idx;
874
875 // Boundary vertices are marked as 'unknown'
876 boundary_v = unmatched_facets;
877 std::ranges::sort(boundary_v);
878 auto [unique_end, range_end] = std::ranges::unique(boundary_v);
879 boundary_v.erase(unique_end, range_end);
880
881 // Remove -1 if it occurs in boundary vertices (may occur in mixed
882 // topology)
883 if (!boundary_v.empty() > 0 and boundary_v[0] == -1)
884 boundary_v.erase(boundary_v.begin());
885 }
886
887 // Create Topology
888 Topology topology = create_topology(comm, cells1_v, original_idx1,
889 ghost_owners, celltype, boundary_v);
890
891 // Create connectivities required higher-order geometries for creating
892 // a Geometry object
893 for (int e = 1; e < topology.dim(); ++e)
894 if (doflayout.num_entity_dofs(e) > 0)
895 topology.create_entities(e);
896 if (element.needs_dof_permutations())
898
899 // Build list of unique (global) node indices from cells1 and
900 // distribute coordinate data
901 std::vector<std::int64_t> nodes1 = cells1;
902 dolfinx::radix_sort(nodes1);
903 auto [unique_end, range_end] = std::ranges::unique(nodes1);
904 nodes1.erase(unique_end, range_end);
905
906 std::vector coords
907 = dolfinx::MPI::distribute_data(comm, nodes1, commg, x, xshape[1]);
908
909 // Create geometry object
910 Geometry geometry
911 = create_geometry(topology, element, nodes1, cells1, coords, xshape[1]);
912
913 return Mesh(comm, std::make_shared<Topology>(std::move(topology)),
914 std::move(geometry));
915}
916
951template <typename U>
953 MPI_Comm comm, MPI_Comm commt,
954 const std::vector<std::span<const std::int64_t>>& cells,
955 const std::vector<fem::CoordinateElement<
956 typename std::remove_reference_t<typename U::value_type>>>& elements,
957 MPI_Comm commg, const U& x, std::array<std::size_t, 2> xshape,
958 const CellPartitionFunction& partitioner)
959{
960 assert(cells.size() == elements.size());
961 std::int32_t num_cell_types = cells.size();
962 std::vector<CellType> celltypes;
963 std::ranges::transform(elements, std::back_inserter(celltypes),
964 [](auto e) { return e.cell_shape(); });
965 std::vector<fem::ElementDofLayout> doflayouts;
966 std::ranges::transform(elements, std::back_inserter(doflayouts),
967 [](auto e) { return e.create_dof_layout(); });
968
969 // Note: `extract_topology` extracts topology data, i.e. just the
970 // vertices. For P1 geometry this should just be the identity
971 // operator. For other elements the filtered lists may have 'gaps',
972 // i.e. the indices might not be contiguous.
973 //
974 // `extract_topology` could be skipped for 'P1 geometry' elements
975
976 // -- Partition topology across ranks of comm
977 std::vector<std::vector<std::int64_t>> cells1(num_cell_types);
978 std::vector<std::vector<std::int64_t>> original_idx1(num_cell_types);
979 std::vector<std::vector<int>> ghost_owners(num_cell_types);
980 if (partitioner)
981 {
982 spdlog::info("Using partitioner with cell data ({} cell types)",
983 num_cell_types);
985 if (commt != MPI_COMM_NULL)
986 {
987 int size = dolfinx::MPI::size(comm);
988 std::vector<std::vector<std::int64_t>> t(num_cell_types);
989 std::vector<std::span<const std::int64_t>> tspan(num_cell_types);
990 for (std::int32_t i = 0; i < num_cell_types; ++i)
991 {
992 t[i] = extract_topology(celltypes[i], doflayouts[i], cells[i]);
993 tspan[i] = std::span(t[i]);
994 }
995 dest = partitioner(commt, size, celltypes, tspan);
996 }
997
998 std::int32_t cell_offset = 0;
999 for (std::int32_t i = 0; i < num_cell_types; ++i)
1000 {
1001 std::size_t num_cell_nodes = doflayouts[i].num_dofs();
1002 assert(cells[i].size() % num_cell_nodes == 0);
1003 std::size_t num_cells = cells[i].size() / num_cell_nodes;
1004
1005 // Extract destination AdjacencyList for this cell type
1006 std::vector<std::int32_t> offsets_i(
1007 std::next(dest.offsets().begin(), cell_offset),
1008 std::next(dest.offsets().begin(), cell_offset + num_cells + 1));
1009 std::vector<std::int32_t> data_i(
1010 std::next(dest.array().begin(), offsets_i.front()),
1011 std::next(dest.array().begin(), offsets_i.back()));
1012 std::int32_t offset_0 = offsets_i.front();
1013 std::ranges::for_each(offsets_i,
1014 [&offset_0](std::int32_t& j) { j -= offset_0; });
1015 graph::AdjacencyList<std::int32_t> dest_i(data_i, offsets_i);
1016 cell_offset += num_cells;
1017
1018 // Distribute cells (topology, includes higher-order 'nodes') to
1019 // destination rank
1020 std::vector<int> src_ranks;
1021 std::tie(cells1[i], src_ranks, original_idx1[i], ghost_owners[i])
1022 = graph::build::distribute(comm, cells[i],
1023 {num_cells, num_cell_nodes}, dest_i);
1024 spdlog::debug("Got {} cells from distribution", cells1[i].size());
1025 }
1026 }
1027 else
1028 {
1029 // No partitioning, construct a global index
1030 std::int64_t num_owned = 0;
1031 for (std::int32_t i = 0; i < num_cell_types; ++i)
1032 {
1033 cells1[i] = std::vector<std::int64_t>(cells[i].begin(), cells[i].end());
1034 std::int32_t num_cell_nodes = doflayouts[i].num_dofs();
1035 assert(cells1[i].size() % num_cell_nodes == 0);
1036 original_idx1[i].resize(cells1[i].size() / num_cell_nodes);
1037 num_owned += original_idx1[i].size();
1038 }
1039 // Add on global offset
1040 std::int64_t global_offset = 0;
1041 MPI_Exscan(&num_owned, &global_offset, 1, MPI_INT64_T, MPI_SUM, comm);
1042 for (std::int32_t i = 0; i < num_cell_types; ++i)
1043 {
1044 std::iota(original_idx1[i].begin(), original_idx1[i].end(),
1045 global_offset);
1046 global_offset += original_idx1[i].size();
1047 }
1048 }
1049
1050 // Extract cell 'topology', i.e. extract the vertices for each cell
1051 // and discard any 'higher-order' nodes
1052 std::vector<std::vector<std::int64_t>> cells1_v(num_cell_types);
1053 for (std::int32_t i = 0; i < num_cell_types; ++i)
1054 {
1055 cells1_v[i] = extract_topology(celltypes[i], doflayouts[i], cells1[i]);
1056 spdlog::info("Extract basic topology: {}->{}", cells1[i].size(),
1057 cells1_v[i].size());
1058 }
1059
1060 // Build local dual graph for owned cells to (i) get list of vertices
1061 // on the process boundary and (ii) apply re-ordering to cells for
1062 // locality
1063 std::vector<std::int64_t> boundary_v;
1064 {
1065 std::vector<std::span<const std::int64_t>> cells1_v_local_cells;
1066
1067 for (std::int32_t i = 0; i < num_cell_types; ++i)
1068 {
1069 std::int32_t num_cell_vertices = mesh::num_cell_vertices(celltypes[i]);
1070 std::int32_t num_owned_cells
1071 = cells1_v[i].size() / num_cell_vertices - ghost_owners[i].size();
1072 cells1_v_local_cells.push_back(
1073 std::span(cells1_v[i].data(), num_owned_cells * num_cell_vertices));
1074 }
1075 spdlog::info("Build local dual graph");
1076 auto [graph, unmatched_facets, max_v, facet_attached_cells]
1077 = build_local_dual_graph(celltypes, cells1_v_local_cells);
1078
1079 // TODO: in the original create_mesh(), cell reordering is done here.
1080 // This needs reworking for mixed cell topology
1081
1082 // Boundary vertices are marked as 'unknown'
1083 boundary_v = unmatched_facets;
1084 std::ranges::sort(boundary_v);
1085 auto [unique_end, range_end] = std::ranges::unique(boundary_v);
1086 boundary_v.erase(unique_end, range_end);
1087
1088 // Remove -1 if it occurs in boundary vertices (may occur in mixed
1089 // topology)
1090 if (!boundary_v.empty() > 0 and boundary_v[0] == -1)
1091 boundary_v.erase(boundary_v.begin());
1092 }
1093
1094 spdlog::debug("Got {} boundary vertices", boundary_v.size());
1095
1096 // Create Topology
1097
1098 std::vector<std::span<const std::int64_t>> cells1_v_span;
1099 std::ranges::transform(cells1_v, std::back_inserter(cells1_v_span),
1100 [](auto& c) { return std::span(c); });
1101 std::vector<std::span<const std::int64_t>> original_idx1_span;
1102 std::ranges::transform(original_idx1, std::back_inserter(original_idx1_span),
1103 [](auto& c) { return std::span(c); });
1104 std::vector<std::span<const int>> ghost_owners_span;
1105 std::ranges::transform(ghost_owners, std::back_inserter(ghost_owners_span),
1106 [](auto& c) { return std::span(c); });
1107
1108 Topology topology
1109 = create_topology(comm, celltypes, cells1_v_span, original_idx1_span,
1110 ghost_owners_span, boundary_v);
1111
1112 // Create connectivities required higher-order geometries for creating
1113 // a Geometry object
1114 for (int i = 0; i < num_cell_types; ++i)
1115 {
1116 for (int e = 1; e < topology.dim(); ++e)
1117 if (doflayouts[i].num_entity_dofs(e) > 0)
1118 topology.create_entities(e);
1119 if (elements[i].needs_dof_permutations())
1120 topology.create_entity_permutations();
1121 }
1122
1123 // Build list of unique (global) node indices from cells1 and
1124 // distribute coordinate data
1125 std::vector<std::int64_t> nodes1, nodes2;
1126 for (std::vector<std::int64_t>& c : cells1)
1127 nodes1.insert(nodes1.end(), c.begin(), c.end());
1128 for (std::vector<std::int64_t>& c : cells1)
1129 nodes2.insert(nodes2.end(), c.begin(), c.end());
1130
1131 dolfinx::radix_sort(nodes1);
1132 auto [unique_end, range_end] = std::ranges::unique(nodes1);
1133 nodes1.erase(unique_end, range_end);
1134
1135 std::vector coords
1136 = dolfinx::MPI::distribute_data(comm, nodes1, commg, x, xshape[1]);
1137
1138 // Create geometry object
1139 Geometry geometry
1140 = create_geometry(topology, elements, nodes1, nodes2, coords, xshape[1]);
1141
1142 return Mesh(comm, std::make_shared<Topology>(std::move(topology)),
1143 std::move(geometry));
1144}
1145
1164template <typename U>
1165Mesh<typename std::remove_reference_t<typename U::value_type>>
1166create_mesh(MPI_Comm comm, std::span<const std::int64_t> cells,
1168 std::remove_reference_t<typename U::value_type>>& elements,
1169 const U& x, std::array<std::size_t, 2> xshape, GhostMode ghost_mode)
1170{
1171 if (dolfinx::MPI::size(comm) == 1)
1172 return create_mesh(comm, comm, cells, elements, comm, x, xshape, nullptr);
1173 else
1174 {
1175 return create_mesh(comm, comm, cells, elements, comm, x, xshape,
1176 create_cell_partitioner(ghost_mode));
1177 }
1178}
1179
1191template <std::floating_point T>
1192std::pair<Geometry<T>, std::vector<int32_t>>
1193create_subgeometry(const Mesh<T>& mesh, int dim,
1194 std::span<const std::int32_t> subentity_to_entity)
1195{
1196 const Geometry<T>& geometry = mesh.geometry();
1197
1198 // Get the geometry dofs in the sub-geometry based on the entities in
1199 // sub-geometry
1200 const fem::ElementDofLayout layout = geometry.cmap().create_dof_layout();
1201
1202 std::vector<std::int32_t> x_indices
1203 = entities_to_geometry(mesh, dim, subentity_to_entity, true);
1204
1205 std::vector<std::int32_t> sub_x_dofs = x_indices;
1206 std::ranges::sort(sub_x_dofs);
1207 auto [unique_end, range_end] = std::ranges::unique(sub_x_dofs);
1208 sub_x_dofs.erase(unique_end, range_end);
1209
1210 // Get the sub-geometry dofs owned by this process
1211 auto x_index_map = geometry.index_map();
1212 assert(x_index_map);
1213
1214 std::shared_ptr<common::IndexMap> sub_x_dof_index_map;
1215 std::vector<std::int32_t> subx_to_x_dofmap;
1216 {
1217 auto [map, new_to_old] = common::create_sub_index_map(
1218 *x_index_map, sub_x_dofs, common::IndexMapOrder::any, true);
1219 sub_x_dof_index_map = std::make_shared<common::IndexMap>(std::move(map));
1220 subx_to_x_dofmap = std::move(new_to_old);
1221 }
1222
1223 // Create sub-geometry coordinates
1224 std::span<const T> x = geometry.x();
1225 std::int32_t sub_num_x_dofs = subx_to_x_dofmap.size();
1226 std::vector<T> sub_x(3 * sub_num_x_dofs);
1227 for (int i = 0; i < sub_num_x_dofs; ++i)
1228 {
1229 std::copy_n(std::next(x.begin(), 3 * subx_to_x_dofmap[i]), 3,
1230 std::next(sub_x.begin(), 3 * i));
1231 }
1232
1233 // Create geometry to sub-geometry map
1234 std::vector<std::int32_t> x_to_subx_dof_map(
1235 x_index_map->size_local() + x_index_map->num_ghosts(), -1);
1236 for (std::size_t i = 0; i < subx_to_x_dofmap.size(); ++i)
1237 x_to_subx_dof_map[subx_to_x_dofmap[i]] = i;
1238
1239 // Create sub-geometry dofmap
1240 std::vector<std::int32_t> sub_x_dofmap;
1241 sub_x_dofmap.reserve(x_indices.size());
1242 std::ranges::transform(x_indices, std::back_inserter(sub_x_dofmap),
1243 [&x_to_subx_dof_map](auto x_dof)
1244 {
1245 assert(x_to_subx_dof_map[x_dof] != -1);
1246 return x_to_subx_dof_map[x_dof];
1247 });
1248
1249 // Create sub-geometry coordinate element
1250 CellType sub_coord_cell
1251 = cell_entity_type(geometry.cmap().cell_shape(), dim, 0);
1252 // Special handling if point meshes, as they only support constant basis
1253 // functions
1254 int degree = geometry.cmap().degree();
1255 if (sub_coord_cell == CellType::point)
1256 degree = 0;
1257 fem::CoordinateElement<T> sub_cmap(sub_coord_cell, degree,
1258 geometry.cmap().variant());
1259
1260 // Sub-geometry input_global_indices
1261 const std::vector<std::int64_t>& igi = geometry.input_global_indices();
1262 std::vector<std::int64_t> sub_igi;
1263 sub_igi.reserve(subx_to_x_dofmap.size());
1264 std::ranges::transform(subx_to_x_dofmap, std::back_inserter(sub_igi),
1265 [&igi](auto sub_x_dof) { return igi[sub_x_dof]; });
1266
1267 // Create geometry
1268 return {Geometry(sub_x_dof_index_map, std::move(sub_x_dofmap), {sub_cmap},
1269 std::move(sub_x), geometry.dim(), std::move(sub_igi)),
1270 std::move(subx_to_x_dofmap)};
1271}
1272
1281template <std::floating_point T>
1282std::tuple<Mesh<T>, std::vector<std::int32_t>, std::vector<std::int32_t>,
1283 std::vector<std::int32_t>>
1284create_submesh(const Mesh<T>& mesh, int dim,
1285 std::span<const std::int32_t> entities)
1286{
1287 // Create sub-topology
1288 mesh.topology_mutable()->create_connectivity(dim, 0);
1289 auto [topology, subentity_to_entity, subvertex_to_vertex]
1290 = mesh::create_subtopology(*mesh.topology(), dim, entities);
1291
1292 // Create sub-geometry
1293 const int tdim = mesh.topology()->dim();
1294 mesh.topology_mutable()->create_entities(dim);
1295 mesh.topology_mutable()->create_connectivity(dim, tdim);
1296 mesh.topology_mutable()->create_connectivity(tdim, dim);
1297 mesh.topology_mutable()->create_entity_permutations();
1298 auto [geometry, subx_to_x_dofmap]
1299 = mesh::create_subgeometry(mesh, dim, subentity_to_entity);
1300
1301 return {Mesh(mesh.comm(), std::make_shared<Topology>(std::move(topology)),
1302 std::move(geometry)),
1303 std::move(subentity_to_entity), std::move(subvertex_to_vertex),
1304 std::move(subx_to_x_dofmap)};
1305}
1306
1307} // namespace dolfinx::mesh
Definition XDMFFile.h:29
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
int num_entity_dofs(int dim) const
Definition ElementDofLayout.cpp:63
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
int num_dofs() const
Definition ElementDofLayout.cpp:61
Definition topologycomputation.h:24
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
std::shared_ptr< const common::IndexMap > index_map() const
Index map.
Definition Geometry.h:160
const fem::CoordinateElement< value_type > & cmap() const
The element that describes the geometry map.
Definition Geometry.h:181
int dim() const
Return dimension of the Euclidean coordinate system.
Definition Geometry.h:117
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan< const std::int32_t, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents< std::size_t, 2 > > dofmap() const
DofMap for the geometry.
Definition Geometry.h:124
const std::vector< std::int64_t > & input_global_indices() const
Global user indices.
Definition Geometry.h:202
std::span< const value_type > x() const
Access geometry degrees-of-freedom data (const version).
Definition Geometry.h:169
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:45
void create_entity_permutations()
Compute entity permutations and reflections.
Definition Topology.cpp:917
std::int32_t create_entities(int dim)
Create entities of given topological dimension.
Definition Topology.cpp:837
int dim() const noexcept
Return the topological dimension of the mesh.
Definition Topology.cpp:800
Requirements on function for geometry marking.
Definition utils.h:461
void reorder_list(std::span< T > list, std::span< const std::int32_t > nodemap)
Re-order an adjacency list of fixed degree.
Definition utils.h:46
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)
The coordinates of 'vertices' for for entities of a given dimension that are attached to specified fa...
Definition utils.h:77
std::pair< std::vector< T >, std::array< std::size_t, 2 > > compute_vertex_coords(const mesh::Mesh< T > &mesh)
Definition utils.h:419
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:681
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.
Finite element method functionality.
Definition assemble_matrix_impl.h:26
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
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
Topology create_topology(MPI_Comm comm, std::span< const std::int64_t > cells, std::span< const std::int64_t > original_cell_index, std::span< const int > ghost_owners, CellType cell_type, std::span< const std::int64_t > boundary_vertices)
Create a mesh topology.
Definition Topology.cpp:1313
GhostMode
Enum for different partitioning ghost modes.
Definition utils.h:36
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:264
Mesh< typename std::remove_reference_t< typename U::value_type > > create_mesh(MPI_Comm comm, MPI_Comm commt, std::span< const std::int64_t > cells, const fem::CoordinateElement< typename std::remove_reference_t< typename U::value_type > > &element, MPI_Comm commg, const U &x, std::array< std::size_t, 2 > xshape, const CellPartitionFunction &partitioner)
Create a distributed mesh from mesh data using a provided graph partitioning function for determining...
Definition utils.h:782
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:265
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::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:1324
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:555
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:211
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. A sub-geometry is sim...
Definition utils.h:1193
std::vector< std::int32_t > 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:632
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:355
std::vector< std::int32_t > locate_entities(const Mesh< T > &mesh, int dim, U marker)
Compute indices of all mesh entities that evaluate to true for the provided geometric marking functio...
Definition utils.h:482
std::vector< std::int32_t > exterior_facet_indices(const Topology &topology)
Compute the indices of all exterior facets that are owned by the caller.
Definition utils.cpp:58
std::vector< std::int32_t > compute_incident_entities(const Topology &topology, std::span< const std::int32_t > entities, int d0, int d1)
Compute incident indices.
Definition utils.cpp:108
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::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. The function should compute the destination rank for ce...
Definition utils.h:185
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:381
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:1284
CellPartitionFunction create_cell_partitioner(mesh::GhostMode ghost_mode=mesh::GhostMode::none, const graph::partition_fn &partfn=&graph::partition_graph)
Definition utils.cpp:85
constexpr __radix_sort radix_sort
Radix sort.
Definition sort.h:124