DOLFINx 0.9.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
utils.h
1// Copyright (C) 2019-2021 Garth N. Wells and Jørgen S. Dokken
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 "BoundingBoxTree.h"
10#include "gjk.h"
11#include <algorithm>
12#include <array>
13#include <concepts>
14#include <cstdint>
15#include <deque>
16#include <dolfinx/graph/AdjacencyList.h>
17#include <dolfinx/mesh/Mesh.h>
18#include <map>
19#include <numeric>
20#include <span>
21#include <vector>
22
23namespace dolfinx::geometry
24{
28template <std::floating_point T>
30{
31 std::vector<int> src_owner;
33 std::vector<int>
35 std::vector<T> dest_points;
36 std::vector<std::int32_t>
39};
40
51template <std::floating_point T>
52std::vector<T> shortest_vector(const mesh::Mesh<T>& mesh, int dim,
53 std::span<const std::int32_t> entities,
54 std::span<const T> points)
55{
56 const int tdim = mesh.topology()->dim();
57 const mesh::Geometry<T>& geometry = mesh.geometry();
58
59 std::span<const T> geom_dofs = geometry.x();
60 auto x_dofmap = geometry.dofmap();
61 std::vector<T> shortest_vectors;
62 shortest_vectors.reserve(3 * entities.size());
63 if (dim == tdim)
64 {
65 for (std::size_t e = 0; e < entities.size(); e++)
66 {
67 // Check that we have sent in valid entities, i.e. that they exist in the
68 // local dofmap. One gets a cryptical memory segfault if entities is -1
69 assert(entities[e] >= 0);
70 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
71 x_dofmap, entities[e], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
72 std::vector<T> nodes(3 * dofs.size());
73 for (std::size_t i = 0; i < dofs.size(); ++i)
74 {
75 const std::int32_t pos = 3 * dofs[i];
76 for (std::size_t j = 0; j < 3; ++j)
77 nodes[3 * i + j] = geom_dofs[pos + j];
78 }
79
80 std::array<T, 3> d
81 = compute_distance_gjk<T>(points.subspan(3 * e, 3), nodes);
82 shortest_vectors.insert(shortest_vectors.end(), d.begin(), d.end());
83 }
84 }
85 else
86 {
87 mesh.topology_mutable()->create_connectivity(dim, tdim);
88 mesh.topology_mutable()->create_connectivity(tdim, dim);
89 auto e_to_c = mesh.topology()->connectivity(dim, tdim);
90 assert(e_to_c);
91 auto c_to_e = mesh.topology_mutable()->connectivity(tdim, dim);
92 assert(c_to_e);
93 for (std::size_t e = 0; e < entities.size(); e++)
94 {
95 const std::int32_t index = entities[e];
96
97 // Find attached cell
98 assert(e_to_c->num_links(index) > 0);
99 const std::int32_t c = e_to_c->links(index)[0];
100
101 // Find local number of entity wrt cell
102 auto cell_entities = c_to_e->links(c);
103 auto it0 = std::find(cell_entities.begin(), cell_entities.end(), index);
104 assert(it0 != cell_entities.end());
105 const int local_cell_entity = std::distance(cell_entities.begin(), it0);
106
107 // Tabulate geometry dofs for the entity
108 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
109 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
110 const std::vector<int> entity_dofs
111 = geometry.cmap().create_dof_layout().entity_closure_dofs(
112 dim, local_cell_entity);
113 std::vector<T> nodes(3 * entity_dofs.size());
114 for (std::size_t i = 0; i < entity_dofs.size(); i++)
115 {
116 const std::int32_t pos = 3 * dofs[entity_dofs[i]];
117 for (std::size_t j = 0; j < 3; ++j)
118 nodes[3 * i + j] = geom_dofs[pos + j];
119 }
120
121 std::array<T, 3> d
122 = compute_distance_gjk<T>(points.subspan(3 * e, 3), nodes);
123 shortest_vectors.insert(shortest_vectors.end(), d.begin(), d.end());
124 }
125 }
126
127 return shortest_vectors;
128}
129
136template <std::floating_point T>
137T compute_squared_distance_bbox(std::span<const T, 6> b,
138 std::span<const T, 3> x)
139{
140 auto b0 = b.template subspan<0, 3>();
141 auto b1 = b.template subspan<3, 3>();
142 return std::transform_reduce(x.begin(), x.end(), b0.begin(), 0.0,
143 std::plus<>{},
144 [](auto x, auto b)
145 {
146 auto dx = x - b;
147 return dx > 0 ? 0 : dx * dx;
148 })
149 + std::transform_reduce(x.begin(), x.end(), b1.begin(), 0.0,
150 std::plus<>{},
151 [](auto x, auto b)
152 {
153 auto dx = x - b;
154 return dx < 0 ? 0 : dx * dx;
155 });
156}
157
173template <std::floating_point T>
174std::vector<T> squared_distance(const mesh::Mesh<T>& mesh, int dim,
175 std::span<const std::int32_t> entities,
176 std::span<const T> points)
177{
178 std::vector<T> v = shortest_vector(mesh, dim, entities, points);
179 std::vector<T> d(v.size() / 3, 0);
180 for (std::size_t i = 0; i < d.size(); ++i)
181 for (std::size_t j = 0; j < 3; ++j)
182 d[i] += v[3 * i + j] * v[3 * i + j];
183 return d;
184}
185
186namespace impl
187{
189constexpr bool is_leaf(std::array<int, 2> bbox)
190{
191 // Leaf nodes are marked by setting child_0 equal to child_1
192 return bbox[0] == bbox[1];
193}
194
199template <std::floating_point T>
200constexpr bool point_in_bbox(const std::array<T, 6>& b, std::span<const T, 3> x)
201{
202 constexpr T rtol = 1e-14;
203 bool in = true;
204 for (std::size_t i = 0; i < 3; i++)
205 {
206 T eps = rtol * (b[i + 3] - b[i]);
207 in &= x[i] >= (b[i] - eps);
208 in &= x[i] <= (b[i + 3] + eps);
209 }
210
211 return in;
212}
213
217template <std::floating_point T>
218constexpr bool bbox_in_bbox(std::span<const T, 6> a, std::span<const T, 6> b)
219{
220 constexpr T rtol = 1e-14;
221 auto a0 = a.template subspan<0, 3>();
222 auto a1 = a.template subspan<3, 3>();
223 auto b0 = b.template subspan<0, 3>();
224 auto b1 = b.template subspan<3, 3>();
225
226 bool in = true;
227 for (std::size_t i = 0; i < 3; i++)
228 {
229 T eps = rtol * (b1[i] - b0[i]);
230 in &= a1[i] >= (b0[i] - eps);
231 in &= a0[i] <= (b1[i] + eps);
232 }
233
234 return in;
235}
236
238template <std::floating_point T>
239std::pair<std::int32_t, T>
240_compute_closest_entity(const geometry::BoundingBoxTree<T>& tree,
241 std::span<const T, 3> point, std::int32_t node,
242 const mesh::Mesh<T>& mesh, std::int32_t closest_entity,
243 T R2)
244{
245 // Get children of current bounding box node (child_1 denotes entity
246 // index for leaves)
247 const std::array<int, 2> bbox = tree.bbox(node);
248 T r2;
249 if (is_leaf(bbox))
250 {
251 // If point cloud tree the exact distance is easy to compute
252 if (tree.tdim() == 0)
253 {
254 std::array<T, 6> diff = tree.get_bbox(node);
255 for (std::size_t k = 0; k < 3; ++k)
256 diff[k] -= point[k];
257 r2 = diff[0] * diff[0] + diff[1] * diff[1] + diff[2] * diff[2];
258 }
259 else
260 {
261 r2 = compute_squared_distance_bbox<T>(tree.get_bbox(node), point);
262
263 // If bounding box closer than previous closest entity, use gjk to
264 // obtain exact distance to the convex hull of the entity
265 if (r2 <= R2)
266 {
267 r2 = squared_distance<T>(mesh, tree.tdim(),
268 std::span(std::next(bbox.begin(), 1), 1),
269 point)
270 .front();
271 }
272 }
273
274 // If entity is closer than best result so far, return it
275 if (r2 <= R2)
276 {
277 closest_entity = bbox.back();
278 R2 = r2;
279 }
280
281 return {closest_entity, R2};
282 }
283 else
284 {
285 // If bounding box is outside radius, then don't search further
286 r2 = compute_squared_distance_bbox<T>(tree.get_bbox(node), point);
287 if (r2 > R2)
288 return {closest_entity, R2};
289
290 // Check both children
291 // We use R2 (as opposed to r2), as a bounding box can be closer
292 // than the actual entity
293 std::pair<int, T> p0 = _compute_closest_entity(tree, point, bbox.front(),
294 mesh, closest_entity, R2);
295 std::pair<int, T> p1 = _compute_closest_entity(tree, point, bbox.back(),
296 mesh, p0.first, p0.second);
297 return p1;
298 }
299}
300
306template <std::floating_point T>
307void _compute_collisions_point(const geometry::BoundingBoxTree<T>& tree,
308 std::span<const T, 3> p,
309 std::vector<std::int32_t>& entities)
310{
311 std::deque<std::int32_t> stack;
312 std::int32_t next = tree.num_bboxes() - 1;
313 while (next != -1)
314 {
315 const std::array<int, 2> bbox = tree.bbox(next);
316 if (is_leaf(bbox) and point_in_bbox(tree.get_bbox(next), p))
317 {
318 // If box is a leaf node then add it to the list of colliding
319 // entities
320 entities.push_back(bbox[1]);
321 next = -1;
322 }
323 else
324 {
325 // Check whether the point collides with child nodes (left and
326 // right)
327 bool left = point_in_bbox(tree.get_bbox(bbox[0]), p);
328 bool right = point_in_bbox(tree.get_bbox(bbox[1]), p);
329 if (left and right)
330 {
331 // If the point collides with both child nodes, add the right
332 // node to the stack (for later visiting) and continue the tree
333 // traversal with the left subtree
334 stack.push_back(bbox[1]);
335 next = bbox[0];
336 }
337 else if (left)
338 {
339 // Traverse the current node's left subtree
340 next = bbox[0];
341 }
342 else if (right)
343 {
344 // Traverse the current node's right subtree
345 next = bbox[1];
346 }
347 else
348 next = -1;
349 }
350
351 // If tree traversal reaches a dead end (box is a leaf node or no
352 // collision detected), check the stack for deferred subtrees
353 if (next == -1 and !stack.empty())
354 {
355 next = stack.back();
356 stack.pop_back();
357 }
358 }
359}
360
361// Compute collisions with tree (recursive)
362template <std::floating_point T>
363void _compute_collisions_tree(const geometry::BoundingBoxTree<T>& A,
364 const geometry::BoundingBoxTree<T>& B,
365 std::int32_t node_A, std::int32_t node_B,
366 std::vector<std::int32_t>& entities)
367{
368 // If bounding boxes don't collide, then don't search further
369 if (!bbox_in_bbox<T>(A.get_bbox(node_A), B.get_bbox(node_B)))
370 return;
371
372 // Get bounding boxes for current nodes
373 const std::array<std::int32_t, 2> bbox_A = A.bbox(node_A);
374 const std::array<std::int32_t, 2> bbox_B = B.bbox(node_B);
375
376 // Check whether we've reached a leaf in A or B
377 const bool is_leaf_A = is_leaf(bbox_A);
378 const bool is_leaf_B = is_leaf(bbox_B);
379 if (is_leaf_A and is_leaf_B)
380 {
381 // If both boxes are leaves (which we know collide), then add them
382 // child_1 denotes entity for leaves
383 entities.push_back(bbox_A[1]);
384 entities.push_back(bbox_B[1]);
385 }
386 else if (is_leaf_A)
387 {
388 // If we reached the leaf in A, then descend B
389 _compute_collisions_tree(A, B, node_A, bbox_B[0], entities);
390 _compute_collisions_tree(A, B, node_A, bbox_B[1], entities);
391 }
392 else if (is_leaf_B)
393 {
394 // If we reached the leaf in B, then descend A
395 _compute_collisions_tree(A, B, bbox_A[0], node_B, entities);
396 _compute_collisions_tree(A, B, bbox_A[1], node_B, entities);
397 }
398 else if (node_A > node_B)
399 {
400 // At this point, we know neither is a leaf so descend the largest
401 // tree first. Note that nodes are added in reverse order with the
402 // top bounding box at the end so the largest tree (the one with the
403 // the most boxes left to traverse) has the largest node number.
404 _compute_collisions_tree(A, B, bbox_A[0], node_B, entities);
405 _compute_collisions_tree(A, B, bbox_A[1], node_B, entities);
406 }
407 else
408 {
409 _compute_collisions_tree(A, B, node_A, bbox_B[0], entities);
410 _compute_collisions_tree(A, B, node_A, bbox_B[1], entities);
411 }
412
413 // Note that cases above can be collected in fewer cases but this way
414 // the logic is easier to follow.
415}
416
417} // namespace impl
418
425template <std::floating_point T>
427 std::span<const std::int32_t> entities)
428{
429 spdlog::info("Building point search tree to accelerate distance queries for "
430 "a given topological dimension and subset of entities.");
431
432 const std::vector<T> midpoints
433 = mesh::compute_midpoints(mesh, tdim, entities);
434 std::vector<std::pair<std::array<T, 3>, std::int32_t>> points(
435 entities.size());
436 for (std::size_t i = 0; i < points.size(); ++i)
437 {
438 for (std::size_t j = 0; j < 3; ++j)
439 points[i].first[j] = midpoints[3 * i + j];
440 points[i].second = entities[i];
441 }
442
443 // Build tree
444 return BoundingBoxTree(points);
445}
446
452template <std::floating_point T>
453std::vector<std::int32_t> compute_collisions(const BoundingBoxTree<T>& tree0,
454 const BoundingBoxTree<T>& tree1)
455{
456 // Call recursive find function
457 std::vector<std::int32_t> entities;
458 if (tree0.num_bboxes() > 0 and tree1.num_bboxes() > 0)
459 {
460 impl::_compute_collisions_tree(tree0, tree1, tree0.num_bboxes() - 1,
461 tree1.num_bboxes() - 1, entities);
462 }
463
464 return entities;
465}
466
477template <std::floating_point T>
479compute_collisions(const BoundingBoxTree<T>& tree, std::span<const T> points)
480{
481 if (tree.num_bboxes() > 0)
482 {
483 std::vector<std::int32_t> entities, offsets(points.size() / 3 + 1, 0);
484 entities.reserve(points.size() / 3);
485 for (std::size_t p = 0; p < points.size() / 3; ++p)
486 {
487 impl::_compute_collisions_point(
488 tree, std::span<const T, 3>(points.data() + 3 * p, 3), entities);
489 offsets[p + 1] = entities.size();
490 }
491
492 return graph::AdjacencyList(std::move(entities), std::move(offsets));
493 }
494 else
495 {
497 std::vector<std::int32_t>(),
498 std::vector<std::int32_t>(points.size() / 3 + 1, 0));
499 }
500}
501
519template <std::floating_point T>
521 std::span<const std::int32_t> cells,
522 std::array<T, 3> point, T tol)
523{
524 if (cells.empty())
525 return -1;
526 else
527 {
528 const mesh::Geometry<T>& geometry = mesh.geometry();
529 std::span<const T> geom_dofs = geometry.x();
530 auto x_dofmap = geometry.dofmap();
531 const std::size_t num_nodes = x_dofmap.extent(1);
532 std::vector<T> coordinate_dofs(num_nodes * 3);
533 for (auto cell : cells)
534 {
535 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
536 x_dofmap, cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
537 for (std::size_t i = 0; i < num_nodes; ++i)
538 {
539 std::copy_n(std::next(geom_dofs.begin(), 3 * dofs[i]), 3,
540 std::next(coordinate_dofs.begin(), 3 * i));
541 }
542
543 std::array<T, 3> shortest_vector
544 = compute_distance_gjk<T>(point, coordinate_dofs);
545 T d2 = std::reduce(shortest_vector.begin(), shortest_vector.end(), T(0),
546 [](auto d, auto e) { return d + e * e; });
547 if (d2 < tol)
548 return cell;
549 }
550
551 return -1;
552 }
553}
554
567template <std::floating_point T>
568std::vector<std::int32_t>
570 const BoundingBoxTree<T>& midpoint_tree,
571 const mesh::Mesh<T>& mesh, std::span<const T> points)
572{
573 if (tree.num_bboxes() == 0)
574 return std::vector<std::int32_t>(points.size() / 3, -1);
575
576 std::vector<std::int32_t> entities;
577 entities.reserve(points.size() / 3);
578 for (std::size_t i = 0; i < points.size() / 3; ++i)
579 {
580 // Use midpoint tree to find initial closest entity to the point.
581 // Start by using a leaf node as the initial guess for the input
582 // entity
583 std::array<int, 2> leaf0 = midpoint_tree.bbox(0);
584 assert(impl::is_leaf(leaf0));
585 std::array<T, 6> diff = midpoint_tree.get_bbox(0);
586 for (std::size_t k = 0; k < 3; ++k)
587 diff[k] -= points[3 * i + k];
588 T R2 = diff[0] * diff[0] + diff[1] * diff[1] + diff[2] * diff[2];
589
590 // Use a recursive search through the bounding box tree
591 // to find determine the entity with the closest midpoint.
592 // As the midpoint tree only consist of points, the distance
593 // queries are lightweight.
594 const auto [m_index, m_distance2] = impl::_compute_closest_entity(
595 midpoint_tree, std::span<const T, 3>(points.data() + 3 * i, 3),
596 midpoint_tree.num_bboxes() - 1, mesh, leaf0[0], R2);
597
598 // Use a recursives search through the bounding box tree to
599 // determine which entity is actually closest.
600 // Uses the entity with the closest midpoint as initial guess, and
601 // the distance from the midpoint to the point of interest as the
602 // initial search radius.
603 const auto [index, distance2] = impl::_compute_closest_entity(
604 tree, std::span<const T, 3>(points.data() + 3 * i, 3),
605 tree.num_bboxes() - 1, mesh, m_index, m_distance2);
606
607 entities.push_back(index);
608 }
609
610 return entities;
611}
612
628template <std::floating_point T>
630 const mesh::Mesh<T>& mesh,
631 const graph::AdjacencyList<std::int32_t>& candidate_cells,
632 std::span<const T> points)
633{
634 std::vector<std::int32_t> offsets = {0};
635 offsets.reserve(candidate_cells.num_nodes() + 1);
636 std::vector<std::int32_t> colliding_cells;
637 constexpr T eps2 = 1e-12;
638 const int tdim = mesh.topology()->dim();
639 for (std::int32_t i = 0; i < candidate_cells.num_nodes(); i++)
640 {
641 auto cells = candidate_cells.links(i);
642 std::vector<T> _point(3 * cells.size());
643 for (std::size_t j = 0; j < cells.size(); ++j)
644 for (std::size_t k = 0; k < 3; ++k)
645 _point[3 * j + k] = points[3 * i + k];
646
647 std::vector distances_sq = squared_distance<T>(mesh, tdim, cells, _point);
648 for (std::size_t j = 0; j < cells.size(); j++)
649 if (distances_sq[j] < eps2)
650 colliding_cells.push_back(cells[j]);
651
652 offsets.push_back(colliding_cells.size());
653 }
654
655 return graph::AdjacencyList(std::move(colliding_cells), std::move(offsets));
656}
657
685template <std::floating_point T>
687 std::span<const T> points,
688 T padding)
689{
690 MPI_Comm comm = mesh.comm();
691
692 // Create a global bounding-box tree to find candidate processes with
693 // cells that could collide with the points
694 const int tdim = mesh.topology()->dim();
695 auto cell_map = mesh.topology()->index_map(tdim);
696 const std::int32_t num_cells = cell_map->size_local();
697 // NOTE: Should we send the cells in as input?
698 std::vector<std::int32_t> cells(num_cells, 0);
699 std::iota(cells.begin(), cells.end(), 0);
700 BoundingBoxTree bb(mesh, tdim, cells, padding);
701 BoundingBoxTree global_bbtree = bb.create_global_tree(comm);
702
703 // Compute collisions:
704 // For each point in `points` get the processes it should be sent to
705 graph::AdjacencyList collisions = compute_collisions(global_bbtree, points);
706
707 // Get unique list of outgoing ranks
708 std::vector<std::int32_t> out_ranks = collisions.array();
709 std::ranges::sort(out_ranks);
710 auto [unique_end, range_end] = std::ranges::unique(out_ranks);
711 out_ranks.erase(unique_end, range_end);
712
713 // Compute incoming edges (source processes)
714 std::vector in_ranks = dolfinx::MPI::compute_graph_edges_nbx(comm, out_ranks);
715 std::ranges::sort(in_ranks);
716
717 // Create neighborhood communicator in forward direction
718 MPI_Comm forward_comm;
719 MPI_Dist_graph_create_adjacent(
720 comm, in_ranks.size(), in_ranks.data(), MPI_UNWEIGHTED, out_ranks.size(),
721 out_ranks.data(), MPI_UNWEIGHTED, MPI_INFO_NULL, false, &forward_comm);
722
723 // Compute map from global mpi rank to neighbor rank, "collisions"
724 // uses global rank
725 std::map<std::int32_t, std::int32_t> rank_to_neighbor;
726 for (std::size_t i = 0; i < out_ranks.size(); i++)
727 rank_to_neighbor[out_ranks[i]] = i;
728
729 // Count the number of points to send per neighbor process
730 std::vector<std::int32_t> send_sizes(out_ranks.size());
731 for (std::size_t i = 0; i < points.size() / 3; ++i)
732 for (auto p : collisions.links(i))
733 send_sizes[rank_to_neighbor[p]] += 3;
734
735 // Compute receive sizes
736 std::vector<std::int32_t> recv_sizes(in_ranks.size());
737 send_sizes.reserve(1);
738 recv_sizes.reserve(1);
739 MPI_Request sizes_request;
740 MPI_Ineighbor_alltoall(send_sizes.data(), 1, MPI_INT, recv_sizes.data(), 1,
741 MPI_INT, forward_comm, &sizes_request);
742
743 // Compute sending offsets
744 std::vector<std::int32_t> send_offsets(send_sizes.size() + 1, 0);
745 std::partial_sum(send_sizes.begin(), send_sizes.end(),
746 std::next(send_offsets.begin(), 1));
747
748 // Pack data to send and store unpack map
749 std::vector<T> send_data(send_offsets.back());
750 std::vector<std::int32_t> counter(send_sizes.size(), 0);
751 // unpack map: [index in adj list][pos in x]
752 std::vector<std::int32_t> unpack_map(send_offsets.back() / 3);
753 for (std::size_t i = 0; i < points.size(); i += 3)
754 {
755 for (auto p : collisions.links(i / 3))
756 {
757 int neighbor = rank_to_neighbor[p];
758 int pos = send_offsets[neighbor] + counter[neighbor];
759 auto it = std::next(send_data.begin(), pos);
760 std::copy_n(std::next(points.begin(), i), 3, it);
761 unpack_map[pos / 3] = i / 3;
762 counter[neighbor] += 3;
763 }
764 }
765
766 MPI_Wait(&sizes_request, MPI_STATUS_IGNORE);
767 std::vector<std::int32_t> recv_offsets(in_ranks.size() + 1, 0);
768 std::partial_sum(recv_sizes.begin(), recv_sizes.end(),
769 std::next(recv_offsets.begin(), 1));
770
771 std::vector<T> received_points((std::size_t)recv_offsets.back());
772 MPI_Neighbor_alltoallv(
773 send_data.data(), send_sizes.data(), send_offsets.data(),
774 dolfinx::MPI::mpi_type<T>(), received_points.data(), recv_sizes.data(),
775 recv_offsets.data(), dolfinx::MPI::mpi_type<T>(), forward_comm);
776
777 // Get mesh geometry for closest entity
778 const mesh::Geometry<T>& geometry = mesh.geometry();
779 std::span<const T> geom_dofs = geometry.x();
780 auto x_dofmap = geometry.dofmap();
781
782 // Compute candidate cells for collisions (and extrapolation)
783 const graph::AdjacencyList<std::int32_t> candidate_collisions
784 = compute_collisions(bb, std::span<const T>(received_points.data(),
785 received_points.size()));
786
787 // Each process checks which points collide with a cell on the process
788 const int rank = dolfinx::MPI::rank(comm);
789 std::vector<std::int32_t> cell_indicator(received_points.size() / 3);
790 std::vector<std::int32_t> closest_cells(received_points.size() / 3);
791 for (std::size_t p = 0; p < received_points.size(); p += 3)
792 {
793 std::array<T, 3> point;
794 std::copy_n(std::next(received_points.begin(), p), 3, point.begin());
795 // Find first colliding cell among the cells with colliding bounding boxes
796 const int colliding_cell = geometry::compute_first_colliding_cell(
797 mesh, candidate_collisions.links(p / 3), point,
798 10 * std::numeric_limits<T>::epsilon());
799 // If a collding cell is found, store the rank of the current process
800 // which will be sent back to the owner of the point
801 cell_indicator[p / 3] = (colliding_cell >= 0) ? rank : -1;
802 // Store the cell index for lookup once the owning processes has determined
803 // the ownership of the point
804 closest_cells[p / 3] = colliding_cell;
805 }
806
807 // Create neighborhood communicator in the reverse direction: send
808 // back col to requesting processes
809 MPI_Comm reverse_comm;
810 MPI_Dist_graph_create_adjacent(
811 comm, out_ranks.size(), out_ranks.data(), MPI_UNWEIGHTED, in_ranks.size(),
812 in_ranks.data(), MPI_UNWEIGHTED, MPI_INFO_NULL, false, &reverse_comm);
813
814 // Reuse sizes and offsets from first communication set
815 // but divide by three
816 {
817 auto rescale = [](auto& x)
818 { std::ranges::transform(x, x.begin(), [](auto e) { return (e / 3); }); };
819 rescale(recv_sizes);
820 rescale(recv_offsets);
821 rescale(send_sizes);
822 rescale(send_offsets);
823
824 // The communication is reversed, so swap recv to send offsets
825 std::swap(recv_sizes, send_sizes);
826 std::swap(recv_offsets, send_offsets);
827 }
828
829 std::vector<std::int32_t> recv_ranks(recv_offsets.back());
830 MPI_Neighbor_alltoallv(cell_indicator.data(), send_sizes.data(),
831 send_offsets.data(), MPI_INT32_T, recv_ranks.data(),
832 recv_sizes.data(), recv_offsets.data(), MPI_INT32_T,
833 reverse_comm);
834
835 std::vector<int> point_owners(points.size() / 3, -1);
836 for (std::size_t i = 0; i < unpack_map.size(); i++)
837 {
838 const std::int32_t pos = unpack_map[i];
839 // Only insert new owner if no owner has previously been found
840 if (recv_ranks[i] >= 0 && point_owners[pos] == -1)
841 point_owners[pos] = recv_ranks[i];
842 }
843
844 // Create extrapolation marker for those points already sent to other
845 // process
846 std::vector<std::uint8_t> send_extrapolate(recv_offsets.back());
847 for (std::int32_t i = 0; i < recv_offsets.back(); i++)
848 {
849 const std::int32_t pos = unpack_map[i];
850 send_extrapolate[i] = point_owners[pos] == -1;
851 }
852
853 // Swap communication direction, to send extrapolation marker to other
854 // processes
855 std::swap(send_sizes, recv_sizes);
856 std::swap(send_offsets, recv_offsets);
857 std::vector<std::uint8_t> dest_extrapolate(recv_offsets.back());
858 MPI_Neighbor_alltoallv(send_extrapolate.data(), send_sizes.data(),
859 send_offsets.data(), MPI_UINT8_T,
860 dest_extrapolate.data(), recv_sizes.data(),
861 recv_offsets.data(), MPI_UINT8_T, forward_comm);
862
863 std::vector<T> squared_distances(received_points.size() / 3, -1);
864
865 for (std::size_t i = 0; i < dest_extrapolate.size(); i++)
866 {
867 if (dest_extrapolate[i] == 1)
868 {
869 assert(closest_cells[i] == -1);
870 std::array<T, 3> point;
871 std::copy_n(std::next(received_points.begin(), 3 * i), 3, point.begin());
872
873 // Find shortest distance among cells with colldiing bounding box
874 T shortest_distance = std::numeric_limits<T>::max();
875 std::int32_t closest_cell = -1;
876 for (auto cell : candidate_collisions.links(i))
877 {
878 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
879 x_dofmap, cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
880 std::vector<T> nodes(3 * dofs.size());
881 for (std::size_t j = 0; j < dofs.size(); ++j)
882 {
883 const int pos = 3 * dofs[j];
884 for (std::size_t k = 0; k < 3; ++k)
885 nodes[3 * j + k] = geom_dofs[pos + k];
886 }
887 const std::array<T, 3> d = compute_distance_gjk<T>(
888 std::span<const T>(point.data(), point.size()), nodes);
889 if (T current_distance = d[0] * d[0] + d[1] * d[1] + d[2] * d[2];
890 current_distance < shortest_distance)
891 {
892 shortest_distance = current_distance;
893 closest_cell = cell;
894 }
895 }
896 closest_cells[i] = closest_cell;
897 squared_distances[i] = shortest_distance;
898 }
899 }
900
901 std::swap(recv_sizes, send_sizes);
902 std::swap(recv_offsets, send_offsets);
903
904 // Get distances from closest entity of points that were on the other process
905 std::vector<T> recv_distances(recv_offsets.back());
906 MPI_Neighbor_alltoallv(
907 squared_distances.data(), send_sizes.data(), send_offsets.data(),
908 dolfinx::MPI::mpi_type<T>(), recv_distances.data(), recv_sizes.data(),
909 recv_offsets.data(), dolfinx::MPI::mpi_type<T>(), reverse_comm);
910
911 // Update point ownership with extrapolation information
912 std::vector<T> closest_distance(point_owners.size(),
913 std::numeric_limits<T>::max());
914 for (std::size_t i = 0; i < out_ranks.size(); i++)
915 {
916 for (std::int32_t j = recv_offsets[i]; j < recv_offsets[i + 1]; j++)
917 {
918 const std::int32_t pos = unpack_map[j];
919 auto current_dist = recv_distances[j];
920 // Update if closer than previous guess and was found
921 if (auto d = closest_distance[pos];
922 (current_dist > 0) and (current_dist < d))
923 {
924 point_owners[pos] = out_ranks[i];
925 closest_distance[pos] = current_dist;
926 }
927 }
928 }
929
930 // Communication is reversed again to send dest ranks to all processes
931 std::swap(send_sizes, recv_sizes);
932 std::swap(send_offsets, recv_offsets);
933
934 // Pack ownership data
935 std::vector<std::int32_t> send_owners(send_offsets.back());
936 std::ranges::fill(counter, 0);
937 for (std::size_t i = 0; i < points.size() / 3; ++i)
938 {
939 for (auto p : collisions.links(i))
940 {
941 int neighbor = rank_to_neighbor[p];
942 send_owners[send_offsets[neighbor] + counter[neighbor]++]
943 = point_owners[i];
944 }
945 }
946
947 // Send ownership info
948 std::vector<std::int32_t> dest_ranks(recv_offsets.back());
949 MPI_Neighbor_alltoallv(send_owners.data(), send_sizes.data(),
950 send_offsets.data(), MPI_INT32_T, dest_ranks.data(),
951 recv_sizes.data(), recv_offsets.data(), MPI_INT32_T,
952 forward_comm);
953
954 // Unpack dest ranks if point owner is this rank
955 std::vector<int> owned_recv_ranks;
956 owned_recv_ranks.reserve(recv_offsets.back());
957 std::vector<T> owned_recv_points;
958 std::vector<std::int32_t> owned_recv_cells;
959 for (std::size_t i = 0; i < in_ranks.size(); i++)
960 {
961 for (std::int32_t j = recv_offsets[i]; j < recv_offsets[i + 1]; j++)
962 {
963 if (rank == dest_ranks[j])
964 {
965 owned_recv_ranks.push_back(in_ranks[i]);
966 owned_recv_points.insert(
967 owned_recv_points.end(), std::next(received_points.cbegin(), 3 * j),
968 std::next(received_points.cbegin(), 3 * (j + 1)));
969 owned_recv_cells.push_back(closest_cells[j]);
970 }
971 }
972 }
973
974 MPI_Comm_free(&forward_comm);
975 MPI_Comm_free(&reverse_comm);
976 return PointOwnershipData<T>{.src_owner = std::move(point_owners),
977 .dest_owners = std::move(owned_recv_ranks),
978 .dest_points = std::move(owned_recv_points),
979 .dest_cells = std::move(owned_recv_cells)};
980}
981
982} // namespace dolfinx::geometry
Definition BoundingBoxTree.h:204
BoundingBoxTree create_global_tree(MPI_Comm comm) const
Definition BoundingBoxTree.h:324
std::int32_t num_bboxes() const
Return number of bounding boxes.
Definition BoundingBoxTree.h:360
std::array< T, 6 > get_bbox(std::size_t node) const
Return bounding box coordinates for a given node in the tree,.
Definition BoundingBoxTree.h:312
std::array< std::int32_t, 2 > bbox(std::size_t node) const
Definition BoundingBoxTree.h:380
Definition topologycomputation.h:24
std::span< T > links(std::size_t node)
Definition AdjacencyList.h:111
const std::vector< T > & array() const
Return contiguous array of links for all nodes (const version)
Definition AdjacencyList.h:128
std::int32_t num_nodes() const
Definition AdjacencyList.h:96
Geometry stores the geometry imposed on a mesh.
Definition Geometry.h:34
const fem::CoordinateElement< value_type > & cmap() const
The element that describes the geometry map.
Definition Geometry.h:181
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
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
std::vector< int > compute_graph_edges_nbx(MPI_Comm comm, std::span< const int > edges)
Determine incoming graph edges using the NBX consensus algorithm.
Definition MPI.cpp:162
int rank(MPI_Comm comm)
Return process rank for the communicator.
Definition MPI.cpp:64
constexpr MPI_Datatype mpi_type()
MPI Type.
Definition MPI.h:273
Geometry data structures and algorithms.
Definition BoundingBoxTree.h:21
PointOwnershipData< T > determine_point_ownership(const mesh::Mesh< T > &mesh, std::span< const T > points, T padding)
Given a set of points, determine which process is colliding, using the GJK algorithm on cells to dete...
Definition utils.h:686
std::vector< T > squared_distance(const mesh::Mesh< T > &mesh, int dim, std::span< const std::int32_t > entities, std::span< const T > points)
Compute the squared distance between a point and a mesh entity.
Definition utils.h:174
std::vector< std::int32_t > compute_collisions(const BoundingBoxTree< T > &tree0, const BoundingBoxTree< T > &tree1)
Compute all collisions between two bounding box trees.
Definition utils.h:453
std::int32_t compute_first_colliding_cell(const mesh::Mesh< T > &mesh, std::span< const std::int32_t > cells, std::array< T, 3 > point, T tol)
Given a set of cells, find the first one that collides with a point.
Definition utils.h:520
T compute_squared_distance_bbox(std::span< const T, 6 > b, std::span< const T, 3 > x)
Compute squared distance between point and bounding box.
Definition utils.h:137
graph::AdjacencyList< std::int32_t > compute_colliding_cells(const mesh::Mesh< T > &mesh, const graph::AdjacencyList< std::int32_t > &candidate_cells, std::span< const T > points)
Compute which cells collide with a point.
Definition utils.h:629
std::vector< T > shortest_vector(const mesh::Mesh< T > &mesh, int dim, std::span< const std::int32_t > entities, std::span< const T > points)
Compute the shortest vector from a mesh entity to a point.
Definition utils.h:52
std::array< T, 3 > compute_distance_gjk(std::span< const T > p, std::span< const T > q)
Compute the distance between two convex bodies p and q, each defined by a set of points.
Definition gjk.h:256
BoundingBoxTree< T > create_midpoint_tree(const mesh::Mesh< T > &mesh, int tdim, std::span< const std::int32_t > entities)
Create a bounding box tree for the midpoints of a subset of entities.
Definition utils.h:426
std::vector< std::int32_t > compute_closest_entity(const BoundingBoxTree< T > &tree, const BoundingBoxTree< T > &midpoint_tree, const mesh::Mesh< T > &mesh, std::span< const T > points)
Compute closest mesh entity to a point.
Definition utils.h:569
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
Information on the ownership of points distributed across processes.
Definition utils.h:30
std::vector< T > dest_points
Points that are owned by current process.
Definition utils.h:35
std::vector< std::int32_t > dest_cells
Definition utils.h:37
std::vector< int > dest_owners
Ranks that sent dest_points to current process.
Definition utils.h:34
std::vector< int > src_owner
Definition utils.h:31