Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d0/d47/mesh_2utils_8h_source.html
DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
utils.h
Go to the documentation of this file.
1// Copyright (C) 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 <basix/mdspan.hpp>
13#include <concepts>
14#include <dolfinx/graph/AdjacencyList.h>
15#include <dolfinx/graph/ordering.h>
16#include <dolfinx/graph/partition.h>
17#include <functional>
18#include <mpi.h>
19#include <span>
20
23
24namespace dolfinx::fem
25{
26class ElementDofLayout;
27}
28
29namespace dolfinx::mesh
30{
31enum class CellType;
32
34enum class GhostMode : int
35{
36 none,
37 shared_facet,
38 shared_vertex
39};
40
41namespace impl
42{
44template <typename T>
45void reorder_list(std::span<T> list, std::span<const std::int32_t> nodemap)
46{
47 if (nodemap.empty())
48 return;
49
50 assert(list.size() % nodemap.size() == 0);
51 int degree = list.size() / nodemap.size();
52 const std::vector<T> orig(list.begin(), list.end());
53 for (std::size_t n = 0; n < nodemap.size(); ++n)
54 {
55 auto links_old = std::span(orig.data() + n * degree, degree);
56 auto links_new = list.subspan(nodemap[n] * degree, degree);
57 std::copy(links_old.begin(), links_old.end(), links_new.begin());
58 }
59}
60
74template <std::floating_point T>
75std::tuple<std::vector<std::int32_t>, std::vector<T>, std::vector<std::int32_t>>
77 std::span<const std::int32_t> facets)
78{
79 auto topology = mesh.topology();
80 assert(topology);
81 const int tdim = topology->dim();
82 if (dim == tdim)
83 {
84 throw std::runtime_error(
85 "Cannot use mesh::locate_entities_boundary (boundary) for cells.");
86 }
87
88 // Build set of vertices on boundary and set of boundary entities
89 mesh.topology_mutable()->create_connectivity(tdim - 1, 0);
90 mesh.topology_mutable()->create_connectivity(tdim - 1, dim);
91 std::vector<std::int32_t> vertices, entities;
92 {
93 auto f_to_v = topology->connectivity(tdim - 1, 0);
94 assert(f_to_v);
95 auto f_to_e = topology->connectivity(tdim - 1, dim);
96 assert(f_to_e);
97 for (auto f : facets)
98 {
99 auto v = f_to_v->links(f);
100 vertices.insert(vertices.end(), v.begin(), v.end());
101 auto e = f_to_e->links(f);
102 entities.insert(entities.end(), e.begin(), e.end());
103 }
104
105 // Build vector of boundary vertices
106 std::sort(vertices.begin(), vertices.end());
107 vertices.erase(std::unique(vertices.begin(), vertices.end()),
108 vertices.end());
109 std::sort(entities.begin(), entities.end());
110 entities.erase(std::unique(entities.begin(), entities.end()),
111 entities.end());
112 }
113
114 // Get geometry data
115 auto x_dofmap = mesh.geometry().dofmap();
116 std::span<const T> x_nodes = mesh.geometry().x();
117
118 // Get all vertex 'node' indices
119 mesh.topology_mutable()->create_connectivity(0, tdim);
120 mesh.topology_mutable()->create_connectivity(tdim, 0);
121 auto v_to_c = topology->connectivity(0, tdim);
122 assert(v_to_c);
123 auto c_to_v = topology->connectivity(tdim, 0);
124 assert(c_to_v);
125 std::vector<T> x_vertices(3 * vertices.size(), -1.0);
126 std::vector<std::int32_t> vertex_to_pos(v_to_c->num_nodes(), -1);
127 for (std::size_t i = 0; i < vertices.size(); ++i)
128 {
129 const std::int32_t v = vertices[i];
130
131 // Get first cell and find position
132 const int c = v_to_c->links(v).front();
133 auto cell_vertices = c_to_v->links(c);
134 auto it = std::find(cell_vertices.begin(), cell_vertices.end(), v);
135 assert(it != cell_vertices.end());
136 const int local_pos = std::distance(cell_vertices.begin(), it);
137
138 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
139 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
140 for (std::size_t j = 0; j < 3; ++j)
141 x_vertices[j * vertices.size() + i] = x_nodes[3 * dofs[local_pos] + j];
142 vertex_to_pos[v] = i;
143 }
144
145 return {std::move(entities), std::move(x_vertices), std::move(vertex_to_pos)};
146}
147
148} // namespace impl
149
161std::vector<std::int32_t> exterior_facet_indices(const Topology& topology);
162
178using CellPartitionFunction = std::function<graph::AdjacencyList<std::int32_t>(
179 MPI_Comm comm, int nparts, CellType cell_type,
181
191std::vector<std::int64_t> extract_topology(CellType cell_type,
192 const fem::ElementDofLayout& layout,
193 std::span<const std::int64_t> cells);
194
203template <std::floating_point T>
204std::vector<T> h(const Mesh<T>& mesh, std::span<const std::int32_t> entities,
205 int dim)
206{
207 if (entities.empty())
208 return std::vector<T>();
209 if (dim == 0)
210 return std::vector<T>(entities.size(), 0);
211
212 // Get the geometry dofs for the vertices of each entity
213 const std::vector<std::int32_t> vertex_xdofs
214 = entities_to_geometry(mesh, dim, entities, false);
215 assert(!entities.empty());
216 const std::size_t num_vertices = vertex_xdofs.size() / entities.size();
217
218 // Get the geometry coordinate
219 std::span<const T> x = mesh.geometry().x();
220
221 // Function to compute the length of (p0 - p1)
222 auto delta_norm = [](auto&& p0, auto&& p1)
223 {
224 T norm = 0;
225 for (std::size_t i = 0; i < 3; ++i)
226 norm += (p0[i] - p1[i]) * (p0[i] - p1[i]);
227 return std::sqrt(norm);
228 };
229
230 // Compute greatest distance between any to vertices
231 assert(dim > 0);
232 std::vector<T> h(entities.size(), 0);
233 for (std::size_t e = 0; e < entities.size(); ++e)
234 {
235 // Get geometry 'dof' for each vertex of entity e
236 std::span<const std::int32_t> e_vertices(
237 vertex_xdofs.data() + e * num_vertices, num_vertices);
238
239 // Compute maximum distance between any two vertices
240 for (std::size_t i = 0; i < e_vertices.size(); ++i)
241 {
242 std::span<const T, 3> p0(x.data() + 3 * e_vertices[i], 3);
243 for (std::size_t j = i + 1; j < e_vertices.size(); ++j)
244 {
245 std::span<const T, 3> p1(x.data() + 3 * e_vertices[j], 3);
246 h[e] = std::max(h[e], delta_norm(p0, p1));
247 }
248 }
249 }
250
251 return h;
252}
253
257template <std::floating_point T>
258std::vector<T> cell_normals(const Mesh<T>& mesh, int dim,
259 std::span<const std::int32_t> entities)
260{
261 auto topology = mesh.topology();
262 assert(topology);
263
264 if (entities.empty())
265 return std::vector<T>();
266
267 if (topology->cell_type() == CellType::prism and dim == 2)
268 throw std::runtime_error("More work needed for prism cell");
269
270 const int gdim = mesh.geometry().dim();
271 const CellType type = cell_entity_type(topology->cell_type(), dim, 0);
272
273 // Find geometry nodes for topology entities
274 std::span<const T> x = mesh.geometry().x();
275
276 // Orient cells if they are tetrahedron
277 bool orient = false;
278 if (topology->cell_type() == CellType::tetrahedron)
279 orient = true;
280
281 std::vector<std::int32_t> geometry_entities
282 = entities_to_geometry(mesh, dim, entities, orient);
283
284 const std::size_t shape1 = geometry_entities.size() / entities.size();
285 std::vector<T> n(entities.size() * 3);
286 switch (type)
287 {
288 case CellType::interval:
289 {
290 if (gdim > 2)
291 throw std::invalid_argument("Interval cell normal undefined in 3D");
292 for (std::size_t i = 0; i < entities.size(); ++i)
293 {
294 // Get the two vertices as points
295 std::array vertices{geometry_entities[i * shape1],
296 geometry_entities[i * shape1 + 1]};
297 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
298 std::span<const T, 3>(x.data() + 3 * vertices[1], 3)};
299
300 // Define normal by rotating tangent counter-clockwise
301 std::array<T, 3> t;
302 std::transform(p[1].begin(), p[1].end(), p[0].begin(), t.begin(),
303 [](auto x, auto y) { return x - y; });
304
305 T norm = std::sqrt(t[0] * t[0] + t[1] * t[1]);
306 std::span<T, 3> ni(n.data() + 3 * i, 3);
307 ni[0] = -t[1] / norm;
308 ni[1] = t[0] / norm;
309 ni[2] = 0.0;
310 }
311 return n;
312 }
313 case CellType::triangle:
314 {
315 for (std::size_t i = 0; i < entities.size(); ++i)
316 {
317 // Get the three vertices as points
318 std::array vertices = {geometry_entities[i * shape1 + 0],
319 geometry_entities[i * shape1 + 1],
320 geometry_entities[i * shape1 + 2]};
321 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
322 std::span<const T, 3>(x.data() + 3 * vertices[1], 3),
323 std::span<const T, 3>(x.data() + 3 * vertices[2], 3)};
324
325 // Compute (p1 - p0) and (p2 - p0)
326 std::array<T, 3> dp1, dp2;
327 std::transform(p[1].begin(), p[1].end(), p[0].begin(), dp1.begin(),
328 [](auto x, auto y) { return x - y; });
329 std::transform(p[2].begin(), p[2].end(), p[0].begin(), dp2.begin(),
330 [](auto x, auto y) { return x - y; });
331
332 // Define cell normal via cross product of first two edges
333 std::array<T, 3> ni = math::cross(dp1, dp2);
334 T norm = std::sqrt(ni[0] * ni[0] + ni[1] * ni[1] + ni[2] * ni[2]);
335 std::transform(ni.begin(), ni.end(), std::next(n.begin(), 3 * i),
336 [norm](auto x) { return x / norm; });
337 }
338
339 return n;
340 }
341 case CellType::quadrilateral:
342 {
343 // TODO: check
344 for (std::size_t i = 0; i < entities.size(); ++i)
345 {
346 // Get the three vertices as points
347 std::array vertices = {geometry_entities[i * shape1 + 0],
348 geometry_entities[i * shape1 + 1],
349 geometry_entities[i * shape1 + 2]};
350 std::array p = {std::span<const T, 3>(x.data() + 3 * vertices[0], 3),
351 std::span<const T, 3>(x.data() + 3 * vertices[1], 3),
352 std::span<const T, 3>(x.data() + 3 * vertices[2], 3)};
353
354 // Compute (p1 - p0) and (p2 - p0)
355 std::array<T, 3> dp1, dp2;
356 std::transform(p[1].begin(), p[1].end(), p[0].begin(), dp1.begin(),
357 [](auto x, auto y) { return x - y; });
358 std::transform(p[2].begin(), p[2].end(), p[0].begin(), dp2.begin(),
359 [](auto x, auto y) { return x - y; });
360
361 // Define cell normal via cross product of first two edges
362 std::array<T, 3> ni = math::cross(dp1, dp2);
363 T norm = std::sqrt(ni[0] * ni[0] + ni[1] * ni[1] + ni[2] * ni[2]);
364 std::transform(ni.begin(), ni.end(), std::next(n.begin(), 3 * i),
365 [norm](auto x) { return x / norm; });
366 }
367
368 return n;
369 }
370 default:
371 throw std::invalid_argument(
372 "cell_normal not supported for this cell type.");
373 }
374}
375
379template <std::floating_point T>
380std::vector<T> compute_midpoints(const Mesh<T>& mesh, int dim,
381 std::span<const std::int32_t> entities)
382{
383 if (entities.empty())
384 return std::vector<T>();
385
386 std::span<const T> x = mesh.geometry().x();
387
388 // Build map from entity -> geometry dof
389 // FIXME: This assumes a linear geometry.
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::transform(p.begin(), p.end(), xg.begin(), 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 const 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 const std::vector<std::int32_t> boundary_facets
571 = exterior_facet_indices(*topology);
572
573 using cmdspan3x_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
574 const T,
575 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
576 std::size_t, 3, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>;
577
578 // Run marker function on the vertex coordinates
579 const auto [facet_entities, xdata, vertex_to_pos]
580 = impl::compute_vertex_coords_boundary(mesh, dim, boundary_facets);
581 cmdspan3x_t x(xdata.data(), 3, xdata.size() / 3);
582 const std::vector<std::int8_t> marked = marker(x);
583 if (marked.size() != x.extent(1))
584 throw std::runtime_error("Length of array of markers is wrong.");
585
586 // Loop over entities and check vertex markers
587 mesh.topology_mutable()->create_entities(dim);
588 auto e_to_v = topology->connectivity(dim, 0);
589 assert(e_to_v);
590 std::vector<std::int32_t> entities;
591 for (auto e : facet_entities)
592 {
593 // Iterate over entity vertices
594 bool all_vertices_marked = true;
595 for (auto v : e_to_v->links(e))
596 {
597 const std::int32_t pos = vertex_to_pos[v];
598 if (!marked[pos])
599 {
600 all_vertices_marked = false;
601 break;
602 }
603 }
604
605 // Mark facet with all vertices marked
606 if (all_vertices_marked)
607 entities.push_back(e);
608 }
609
610 return entities;
611}
612
629template <std::floating_point T>
630std::vector<std::int32_t>
631entities_to_geometry(const Mesh<T>& mesh, int dim,
632 std::span<const std::int32_t> entities, bool orient)
633{
634 auto topology = mesh.topology();
635 assert(topology);
636
637 CellType cell_type = topology->cell_type();
638 if (cell_type == CellType::prism and dim == 2)
639 throw std::runtime_error("More work needed for prism cells");
640 if (orient and (cell_type != CellType::tetrahedron or dim != 2))
641 throw std::runtime_error("Can only orient facets of a tetrahedral mesh");
642
643 const Geometry<T>& geometry = mesh.geometry();
644 auto x = geometry.x();
645
646 const int tdim = topology->dim();
647 mesh.topology_mutable()->create_entities(dim);
648 mesh.topology_mutable()->create_connectivity(dim, tdim);
649 mesh.topology_mutable()->create_connectivity(dim, 0);
650 mesh.topology_mutable()->create_connectivity(tdim, 0);
651
652 auto xdofs = geometry.dofmap();
653 auto e_to_c = topology->connectivity(dim, tdim);
654 if (!e_to_c)
655 {
656 throw std::runtime_error(
657 "Entity-to-cell connectivity has not been computed.");
658 }
659
660 auto e_to_v = topology->connectivity(dim, 0);
661 if (!e_to_v)
662 {
663 throw std::runtime_error(
664 "Entity-to-vertex connectivity has not been computed.");
665 }
666
667 auto c_to_v = topology->connectivity(tdim, 0);
668 if (!e_to_v)
669 {
670 throw std::runtime_error(
671 "Cell-to-vertex connectivity has not been computed.");
672 }
673
674 const std::size_t num_vertices
675 = num_cell_vertices(cell_entity_type(cell_type, dim, 0));
676 std::vector<std::int32_t> geometry_idx(entities.size() * num_vertices);
677 for (std::size_t i = 0; i < entities.size(); ++i)
678 {
679 const std::int32_t idx = entities[i];
680 // Always pick the second cell to be consistent with the e_to_v connectivity
681 const std::int32_t cell = e_to_c->links(idx).back();
682 auto ev = e_to_v->links(idx);
683 assert(ev.size() == num_vertices);
684 auto cv = c_to_v->links(cell);
685 std::span<const std::int32_t> xc(
686 xdofs.data_handle() + xdofs.extent(1) * cell, xdofs.extent(1));
687 for (std::size_t j = 0; j < num_vertices; ++j)
688 {
689 int k = std::distance(cv.begin(), std::find(cv.begin(), cv.end(), ev[j]));
690 assert(k < (int)cv.size());
691 geometry_idx[i * num_vertices + j] = xc[k];
692 }
693
694 if (orient)
695 {
696 // Compute cell midpoint
697 std::array<T, 3> midpoint = {0, 0, 0};
698 for (std::int32_t j : xc)
699 for (int k = 0; k < 3; ++k)
700 midpoint[k] += x[3 * j + k];
701 std::transform(midpoint.begin(), midpoint.end(), midpoint.begin(),
702 [size = xc.size()](auto x) { return x / size; });
703
704 // Compute vector triple product of two edges and vector to midpoint
705 std::array<T, 3> p0, p1, p2;
706 std::copy_n(std::next(x.begin(), 3 * geometry_idx[i * num_vertices + 0]),
707 3, p0.begin());
708 std::copy_n(std::next(x.begin(), 3 * geometry_idx[i * num_vertices + 1]),
709 3, p1.begin());
710 std::copy_n(std::next(x.begin(), 3 * geometry_idx[i * num_vertices + 2]),
711 3, p2.begin());
712
713 std::array<T, 9> a;
714 std::transform(midpoint.begin(), midpoint.end(), p0.begin(), a.begin(),
715 [](auto x, auto y) { return x - y; });
716 std::transform(p1.begin(), p1.end(), p0.begin(), std::next(a.begin(), 3),
717 [](auto x, auto y) { return x - y; });
718 std::transform(p2.begin(), p2.end(), p0.begin(), std::next(a.begin(), 6),
719 [](auto x, auto y) { return x - y; });
720
721 // Midpoint direction should be opposite to normal, hence this
722 // should be negative. Switch points if not.
723 if (math::det(a.data(), {3, 3}) > 0.0)
724 {
725 std::swap(geometry_idx[i * num_vertices + 1],
726 geometry_idx[i * num_vertices + 2]);
727 }
728 }
729 }
730
731 return geometry_idx;
732}
733
739 = mesh::GhostMode::none,
740 const graph::partition_fn& partfn
742
750std::vector<std::int32_t>
751compute_incident_entities(const Topology& topology,
752 std::span<const std::int32_t> entities, int d0,
753 int d1);
754
784template <typename U>
786 MPI_Comm comm, MPI_Comm commt, std::span<const std::int64_t> cells,
788 typename std::remove_reference_t<typename U::value_type>>& element,
789 MPI_Comm commg, const U& x, std::array<std::size_t, 2> xshape,
790 const CellPartitionFunction& partitioner)
791{
792 CellType celltype = element.cell_shape();
793 const fem::ElementDofLayout doflayout = element.create_dof_layout();
794
795 const int num_cell_vertices = mesh::num_cell_vertices(element.cell_shape());
796 const int num_cell_nodes = doflayout.num_dofs();
797
798 // Note: `extract_topology` extracts topology data, i.e. just the
799 // vertices. For P1 geometry this should just be the identity
800 // operator. For other elements the filtered lists may have 'gaps',
801 // i.e. the indices might not be contiguous.
802 //
803 // `extract_topology` could be skipped for 'P1 geometry' elements
804
805 // -- Partition topology across ranks of comm
807 // std::vector<std::int64_t> cells1;
808 std::vector<std::int64_t> original_idx1;
809 std::vector<int> ghost_owners;
810 if (partitioner)
811 {
813 if (commt != MPI_COMM_NULL)
814 {
815 int size = dolfinx::MPI::size(comm);
817 extract_topology(element.cell_shape(), doflayout, cells),
819 dest = partitioner(commt, size, celltype, t);
820 }
821
822 // Distribute cells (topology, includes higher-order 'nodes') to
823 // destination rank
824 std::vector<int> src;
825 auto _cells = graph::regular_adjacency_list(
826 std::vector(cells.begin(), cells.end()), num_cell_nodes);
827 std::tie(cells1, src, original_idx1, ghost_owners)
828 = graph::build::distribute(comm, _cells, dest);
829 }
830 else
831 {
833 std::vector(cells.begin(), cells.end()), num_cell_nodes);
834 std::int64_t offset(0), num_owned(cells1.num_nodes());
835 MPI_Exscan(&num_owned, &offset, 1, MPI_INT64_T, MPI_SUM, comm);
836 original_idx1.resize(cells1.num_nodes());
837 std::iota(original_idx1.begin(), original_idx1.end(), offset);
838 }
839
840 // Extract cell 'topology', i.e. extract the vertices for each cell
841 // and discard any 'higher-order' nodes
842 std::vector<std::int64_t> cells1_v
843 = extract_topology(celltype, doflayout, cells1.array());
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 auto [graph, unmatched_facets, max_v, facet_attached_cells]
857 celltype,
858 std::span(cells1_v.data(), num_owned_cells * num_cell_vertices));
859 const std::vector<int> remap = graph::reorder_gps(graph);
860
861 // Create re-ordered cell lists (leaves ghosts unchanged)
862 std::vector<std::int64_t> _original_idx(original_idx1.size());
863 for (std::size_t i = 0; i < remap.size(); ++i)
864 _original_idx[remap[i]] = original_idx1[i];
865 std::copy_n(std::next(original_idx1.cbegin(), num_owned_cells),
866 ghost_owners.size(),
867 std::next(_original_idx.begin(), num_owned_cells));
868 impl::reorder_list(
869 std::span(cells1_v.data(), remap.size() * num_cell_vertices), remap);
870 impl::reorder_list(
871 std::span(cells1.array().data(), remap.size() * num_cell_nodes), remap);
872 original_idx1 = _original_idx;
873
874 // Boundary vertices are marked as 'unknown'
875 boundary_v = unmatched_facets;
876 std::sort(boundary_v.begin(), boundary_v.end());
877 boundary_v.erase(std::unique(boundary_v.begin(), boundary_v.end()),
878 boundary_v.end());
879
880 // Remove -1 if it occurs in boundary vertices (may occur in mixed
881 // topology)
882 if (!boundary_v.empty() > 0 and boundary_v[0] == -1)
883 boundary_v.erase(boundary_v.begin());
884 }
885
886 // Create Topology
887 Topology topology = create_topology(comm, cells1_v, original_idx1,
888 ghost_owners, celltype, boundary_v);
889
890 // Create connectivities required higher-order geometries for creating
891 // a Geometry object
892 for (int e = 1; e < topology.dim(); ++e)
893 if (doflayout.num_entity_dofs(e) > 0)
894 topology.create_entities(e);
895 if (element.needs_dof_permutations())
897
898 // Build list of unique (global) node indices from cells1 and
899 // distribute coordinate data
900 std::vector<std::int64_t> nodes1 = cells1.array();
901 dolfinx::radix_sort(std::span(nodes1));
902 nodes1.erase(std::unique(nodes1.begin(), nodes1.end()), nodes1.end());
903 std::vector coords
904 = dolfinx::MPI::distribute_data(comm, nodes1, commg, x, xshape[1]);
905
906 // Create geometry object
907 Geometry geometry = create_geometry(topology, element, nodes1, cells1.array(),
908 coords, xshape[1]);
909
910 return Mesh(comm, std::make_shared<Topology>(std::move(topology)),
911 std::move(geometry));
912}
913
932template <typename U>
933Mesh<typename std::remove_reference_t<typename U::value_type>>
934create_mesh(MPI_Comm comm, std::span<const std::int64_t> cells,
936 std::remove_reference_t<typename U::value_type>>& elements,
937 const U& x, std::array<std::size_t, 2> xshape, GhostMode ghost_mode)
938{
939 if (dolfinx::MPI::size(comm) == 1)
940 return create_mesh(comm, comm, cells, elements, comm, x, xshape, nullptr);
941 else
942 {
943 return create_mesh(comm, comm, cells, elements, comm, x, xshape,
944 create_cell_partitioner(ghost_mode));
945 }
946}
947
956template <std::floating_point T>
957std::tuple<Mesh<T>, std::vector<std::int32_t>, std::vector<std::int32_t>,
958 std::vector<std::int32_t>>
959create_submesh(const Mesh<T>& mesh, int dim,
960 std::span<const std::int32_t> entities)
961{
962 // Create sub-topology
963 mesh.topology_mutable()->create_connectivity(dim, 0);
964 auto [topology, subentity_to_entity, subvertex_to_vertex]
965 = mesh::create_subtopology(*mesh.topology(), dim, entities);
966
967 // Create sub-geometry
968 const int tdim = mesh.topology()->dim();
969 mesh.topology_mutable()->create_entities(dim);
970 mesh.topology_mutable()->create_connectivity(dim, tdim);
971 mesh.topology_mutable()->create_connectivity(tdim, dim);
972 auto [geometry, subx_to_x_dofmap] = mesh::create_subgeometry(
973 *mesh.topology(), mesh.geometry(), dim, subentity_to_entity);
974
975 return {Mesh(mesh.comm(), std::make_shared<Topology>(std::move(topology)),
976 std::move(geometry)),
977 std::move(subentity_to_entity), std::move(subvertex_to_vertex),
978 std::move(subx_to_x_dofmap)};
979}
980
981} // namespace dolfinx::mesh
Definition CoordinateElement.h:38
Definition ElementDofLayout.h:30
int num_entity_dofs(int dim) const
Definition ElementDofLayout.cpp:63
int num_dofs() const
Definition ElementDofLayout.cpp:61
Definition AdjacencyList.h:28
const std::vector< T > & array() const
Return contiguous array of links for all nodes (const version)
Definition AdjacencyList.h:129
std::int32_t num_nodes() const
Definition AdjacencyList.h:97
Geometry stores the geometry imposed on a mesh.
Definition Geometry.h:33
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:123
std::span< const value_type > x() const
Access geometry degrees-of-freedom data (const version).
Definition Geometry.h:168
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition Mesh.h:23
Topology stores the topology of a mesh, consisting of mesh entities and connectivity (incidence relat...
Definition Topology.h:44
void create_entity_permutations()
Compute entity permutations and reflections.
Definition Topology.cpp:908
std::int32_t create_entities(int dim)
Create entities of given topological dimension.
Definition Topology.cpp:829
int dim() const noexcept
Return the topological dimension of the mesh.
Definition Topology.cpp:792
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:45
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 give dimension that are attached to specified fac...
Definition utils.h:76
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:676
int size(MPI_Comm comm)
Definition MPI.cpp:72
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:32
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
AdjacencyList< typename std::decay_t< U >::value_type > regular_adjacency_list(U &&data, int degree)
Construct a constant degree (valency) adjacency list.
Definition AdjacencyList.h:183
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:1325
GhostMode
Enum for different partitioning ghost modes.
Definition utils.h:35
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:263
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:785
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:258
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:1337
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:204
std::function< graph::AdjacencyList< std::int32_t >( MPI_Comm comm, int nparts, CellType cell_type, const graph::AdjacencyList< std::int64_t > &cells)> CellPartitionFunction
Signature for the cell partitioning function. The function should compute the destination rank for ce...
Definition utils.h:178
std::tuple< graph::AdjacencyList< std::int32_t >, std::vector< std::int64_t >, std::size_t, std::vector< std::int32_t > > build_local_dual_graph(CellType celltype, std::span< const std::int64_t > cells)
Compute the local part of the dual graph (cell-cell connections via facets) and facet with only one a...
Definition graphbuild.cpp:353
std::pair< Geometry< T >, std::vector< int32_t > > create_subgeometry(const Topology &topology, const Geometry< T > &geometry, int dim, std::span< const std::int32_t > subentity_to_entity)
Create a sub-geometry for a subset of entities.
Definition Geometry.h:456
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 > entities_to_geometry(const Mesh< T > &mesh, int dim, std::span< const std::int32_t > entities, bool orient)
Determine the indices in the geometry data for each vertex of the given mesh entities.
Definition utils.h:631
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:109
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< 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:380
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:959
CellPartitionFunction create_cell_partitioner(mesh::GhostMode ghost_mode=mesh::GhostMode::none, const graph::partition_fn &partfn=&graph::partition_graph)
Definition utils.cpp:87
void radix_sort(std::span< T > array)
Definition sort.h:27