Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/dd/dd0/geometry_2utils_8h_source.html
DOLFINx 0.8.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 <array>
12#include <concepts>
13#include <cstdint>
14#include <deque>
15#include <dolfinx/graph/AdjacencyList.h>
16#include <dolfinx/mesh/Mesh.h>
17#include <map>
18#include <numeric>
19#include <span>
20#include <vector>
21
22namespace dolfinx::geometry
23{
24
35template <std::floating_point T>
36std::vector<T> shortest_vector(const mesh::Mesh<T>& mesh, int dim,
37 std::span<const std::int32_t> entities,
38 std::span<const T> points)
39{
40 const int tdim = mesh.topology()->dim();
41 const mesh::Geometry<T>& geometry = mesh.geometry();
42
43 std::span<const T> geom_dofs = geometry.x();
44 auto x_dofmap = geometry.dofmap();
45 std::vector<T> shortest_vectors;
46 shortest_vectors.reserve(3 * entities.size());
47 if (dim == tdim)
48 {
49 for (std::size_t e = 0; e < entities.size(); e++)
50 {
51 // Check that we have sent in valid entities, i.e. that they exist in the
52 // local dofmap. One gets a cryptical memory segfault if entities is -1
53 assert(entities[e] >= 0);
54 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
55 x_dofmap, entities[e], MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
56 std::vector<T> nodes(3 * dofs.size());
57 for (std::size_t i = 0; i < dofs.size(); ++i)
58 {
59 const std::int32_t pos = 3 * dofs[i];
60 for (std::size_t j = 0; j < 3; ++j)
61 nodes[3 * i + j] = geom_dofs[pos + j];
62 }
63
64 std::array<T, 3> d
65 = compute_distance_gjk<T>(points.subspan(3 * e, 3), nodes);
66 shortest_vectors.insert(shortest_vectors.end(), d.begin(), d.end());
67 }
68 }
69 else
70 {
71 mesh.topology_mutable()->create_connectivity(dim, tdim);
72 mesh.topology_mutable()->create_connectivity(tdim, dim);
73 auto e_to_c = mesh.topology()->connectivity(dim, tdim);
74 assert(e_to_c);
75 auto c_to_e = mesh.topology_mutable()->connectivity(tdim, dim);
76 assert(c_to_e);
77 for (std::size_t e = 0; e < entities.size(); e++)
78 {
79 const std::int32_t index = entities[e];
80
81 // Find attached cell
82 assert(e_to_c->num_links(index) > 0);
83 const std::int32_t c = e_to_c->links(index)[0];
84
85 // Find local number of entity wrt cell
86 auto cell_entities = c_to_e->links(c);
87 auto it0 = std::find(cell_entities.begin(), cell_entities.end(), index);
88 assert(it0 != cell_entities.end());
89 const int local_cell_entity = std::distance(cell_entities.begin(), it0);
90
91 // Tabulate geometry dofs for the entity
92 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
93 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
94 const std::vector<int> entity_dofs
95 = geometry.cmap().create_dof_layout().entity_closure_dofs(
96 dim, local_cell_entity);
97 std::vector<T> nodes(3 * entity_dofs.size());
98 for (std::size_t i = 0; i < entity_dofs.size(); i++)
99 {
100 const std::int32_t pos = 3 * dofs[entity_dofs[i]];
101 for (std::size_t j = 0; j < 3; ++j)
102 nodes[3 * i + j] = geom_dofs[pos + j];
103 }
104
105 std::array<T, 3> d
106 = compute_distance_gjk<T>(points.subspan(3 * e, 3), nodes);
107 shortest_vectors.insert(shortest_vectors.end(), d.begin(), d.end());
108 }
109 }
110
111 return shortest_vectors;
112}
113
120template <std::floating_point T>
121T compute_squared_distance_bbox(std::span<const T, 6> b,
122 std::span<const T, 3> x)
123{
124 auto b0 = b.template subspan<0, 3>();
125 auto b1 = b.template subspan<3, 3>();
126 return std::transform_reduce(x.begin(), x.end(), b0.begin(), 0.0,
127 std::plus<>{},
128 [](auto x, auto b)
129 {
130 auto dx = x - b;
131 return dx > 0 ? 0 : dx * dx;
132 })
133 + std::transform_reduce(x.begin(), x.end(), b1.begin(), 0.0,
134 std::plus<>{},
135 [](auto x, auto b)
136 {
137 auto dx = x - b;
138 return dx < 0 ? 0 : dx * dx;
139 });
140}
141
157template <std::floating_point T>
158std::vector<T> squared_distance(const mesh::Mesh<T>& mesh, int dim,
159 std::span<const std::int32_t> entities,
160 std::span<const T> points)
161{
162 std::vector<T> v = shortest_vector(mesh, dim, entities, points);
163 std::vector<T> d(v.size() / 3, 0);
164 for (std::size_t i = 0; i < d.size(); ++i)
165 for (std::size_t j = 0; j < 3; ++j)
166 d[i] += v[3 * i + j] * v[3 * i + j];
167 return d;
168}
169
170namespace impl
171{
173constexpr bool is_leaf(std::array<int, 2> bbox)
174{
175 // Leaf nodes are marked by setting child_0 equal to child_1
176 return bbox[0] == bbox[1];
177}
178
183template <std::floating_point T>
184constexpr bool point_in_bbox(const std::array<T, 6>& b, std::span<const T, 3> x)
185{
186 constexpr T rtol = 1e-14;
187 bool in = true;
188 for (std::size_t i = 0; i < 3; i++)
189 {
190 T eps = rtol * (b[i + 3] - b[i]);
191 in &= x[i] >= (b[i] - eps);
192 in &= x[i] <= (b[i + 3] + eps);
193 }
194
195 return in;
196}
197
201template <std::floating_point T>
202constexpr bool bbox_in_bbox(std::span<const T, 6> a, std::span<const T, 6> b)
203{
204 constexpr T rtol = 1e-14;
205 auto a0 = a.template subspan<0, 3>();
206 auto a1 = a.template subspan<3, 3>();
207 auto b0 = b.template subspan<0, 3>();
208 auto b1 = b.template subspan<3, 3>();
209
210 bool in = true;
211 for (std::size_t i = 0; i < 3; i++)
212 {
213 T eps = rtol * (b1[i] - b0[i]);
214 in &= a1[i] >= (b0[i] - eps);
215 in &= a0[i] <= (b1[i] + eps);
216 }
217
218 return in;
219}
220
222template <std::floating_point T>
223std::pair<std::int32_t, T>
224_compute_closest_entity(const geometry::BoundingBoxTree<T>& tree,
225 std::span<const T, 3> point, std::int32_t node,
226 const mesh::Mesh<T>& mesh, std::int32_t closest_entity,
227 T R2)
228{
229 // Get children of current bounding box node (child_1 denotes entity
230 // index for leaves)
231 const std::array<int, 2> bbox = tree.bbox(node);
232 T r2;
233 if (is_leaf(bbox))
234 {
235 // If point cloud tree the exact distance is easy to compute
236 if (tree.tdim() == 0)
237 {
238 std::array<T, 6> diff = tree.get_bbox(node);
239 for (std::size_t k = 0; k < 3; ++k)
240 diff[k] -= point[k];
241 r2 = diff[0] * diff[0] + diff[1] * diff[1] + diff[2] * diff[2];
242 }
243 else
244 {
245 r2 = compute_squared_distance_bbox<T>(tree.get_bbox(node), point);
246
247 // If bounding box closer than previous closest entity, use gjk to
248 // obtain exact distance to the convex hull of the entity
249 if (r2 <= R2)
250 {
251 r2 = squared_distance<T>(mesh, tree.tdim(),
252 std::span(std::next(bbox.begin(), 1), 1),
253 point)
254 .front();
255 }
256 }
257
258 // If entity is closer than best result so far, return it
259 if (r2 <= R2)
260 {
261 closest_entity = bbox.back();
262 R2 = r2;
263 }
264
265 return {closest_entity, R2};
266 }
267 else
268 {
269 // If bounding box is outside radius, then don't search further
270 r2 = compute_squared_distance_bbox<T>(tree.get_bbox(node), point);
271 if (r2 > R2)
272 return {closest_entity, R2};
273
274 // Check both children
275 // We use R2 (as opposed to r2), as a bounding box can be closer
276 // than the actual entity
277 std::pair<int, T> p0 = _compute_closest_entity(tree, point, bbox.front(),
278 mesh, closest_entity, R2);
279 std::pair<int, T> p1 = _compute_closest_entity(tree, point, bbox.back(),
280 mesh, p0.first, p0.second);
281 return p1;
282 }
283}
284
290template <std::floating_point T>
291void _compute_collisions_point(const geometry::BoundingBoxTree<T>& tree,
292 std::span<const T, 3> p,
293 std::vector<std::int32_t>& entities)
294{
295 std::deque<std::int32_t> stack;
296 std::int32_t next = tree.num_bboxes() - 1;
297 while (next != -1)
298 {
299 const std::array<int, 2> bbox = tree.bbox(next);
300 if (is_leaf(bbox) and point_in_bbox(tree.get_bbox(next), p))
301 {
302 // If box is a leaf node then add it to the list of colliding
303 // entities
304 entities.push_back(bbox[1]);
305 next = -1;
306 }
307 else
308 {
309 // Check whether the point collides with child nodes (left and
310 // right)
311 bool left = point_in_bbox(tree.get_bbox(bbox[0]), p);
312 bool right = point_in_bbox(tree.get_bbox(bbox[1]), p);
313 if (left and right)
314 {
315 // If the point collides with both child nodes, add the right
316 // node to the stack (for later visiting) and continue the tree
317 // traversal with the left subtree
318 stack.push_back(bbox[1]);
319 next = bbox[0];
320 }
321 else if (left)
322 {
323 // Traverse the current node's left subtree
324 next = bbox[0];
325 }
326 else if (right)
327 {
328 // Traverse the current node's right subtree
329 next = bbox[1];
330 }
331 else
332 next = -1;
333 }
334
335 // If tree traversal reaches a dead end (box is a leaf node or no
336 // collision detected), check the stack for deferred subtrees
337 if (next == -1 and !stack.empty())
338 {
339 next = stack.back();
340 stack.pop_back();
341 }
342 }
343}
344
345// Compute collisions with tree (recursive)
346template <std::floating_point T>
347void _compute_collisions_tree(const geometry::BoundingBoxTree<T>& A,
348 const geometry::BoundingBoxTree<T>& B,
349 std::int32_t node_A, std::int32_t node_B,
350 std::vector<std::int32_t>& entities)
351{
352 // If bounding boxes don't collide, then don't search further
353 if (!bbox_in_bbox<T>(A.get_bbox(node_A), B.get_bbox(node_B)))
354 return;
355
356 // Get bounding boxes for current nodes
357 const std::array<std::int32_t, 2> bbox_A = A.bbox(node_A);
358 const std::array<std::int32_t, 2> bbox_B = B.bbox(node_B);
359
360 // Check whether we've reached a leaf in A or B
361 const bool is_leaf_A = is_leaf(bbox_A);
362 const bool is_leaf_B = is_leaf(bbox_B);
363 if (is_leaf_A and is_leaf_B)
364 {
365 // If both boxes are leaves (which we know collide), then add them
366 // child_1 denotes entity for leaves
367 entities.push_back(bbox_A[1]);
368 entities.push_back(bbox_B[1]);
369 }
370 else if (is_leaf_A)
371 {
372 // If we reached the leaf in A, then descend B
373 _compute_collisions_tree(A, B, node_A, bbox_B[0], entities);
374 _compute_collisions_tree(A, B, node_A, bbox_B[1], entities);
375 }
376 else if (is_leaf_B)
377 {
378 // If we reached the leaf in B, then descend A
379 _compute_collisions_tree(A, B, bbox_A[0], node_B, entities);
380 _compute_collisions_tree(A, B, bbox_A[1], node_B, entities);
381 }
382 else if (node_A > node_B)
383 {
384 // At this point, we know neither is a leaf so descend the largest
385 // tree first. Note that nodes are added in reverse order with the
386 // top bounding box at the end so the largest tree (the one with the
387 // the most boxes left to traverse) has the largest node number.
388 _compute_collisions_tree(A, B, bbox_A[0], node_B, entities);
389 _compute_collisions_tree(A, B, bbox_A[1], node_B, entities);
390 }
391 else
392 {
393 _compute_collisions_tree(A, B, node_A, bbox_B[0], entities);
394 _compute_collisions_tree(A, B, node_A, bbox_B[1], entities);
395 }
396
397 // Note that cases above can be collected in fewer cases but this way
398 // the logic is easier to follow.
399}
400
401} // namespace impl
402
409template <std::floating_point T>
411 std::span<const std::int32_t> entities)
412{
413 LOG(INFO) << "Building point search tree to accelerate distance queries for "
414 "a given topological dimension and subset of entities.";
415
416 const std::vector<T> midpoints
417 = mesh::compute_midpoints(mesh, tdim, entities);
418 std::vector<std::pair<std::array<T, 3>, std::int32_t>> points(
419 entities.size());
420 for (std::size_t i = 0; i < points.size(); ++i)
421 {
422 for (std::size_t j = 0; j < 3; ++j)
423 points[i].first[j] = midpoints[3 * i + j];
424 points[i].second = entities[i];
425 }
426
427 // Build tree
428 return BoundingBoxTree(points);
429}
430
436template <std::floating_point T>
437std::vector<std::int32_t> compute_collisions(const BoundingBoxTree<T>& tree0,
438 const BoundingBoxTree<T>& tree1)
439{
440 // Call recursive find function
441 std::vector<std::int32_t> entities;
442 if (tree0.num_bboxes() > 0 and tree1.num_bboxes() > 0)
443 {
444 impl::_compute_collisions_tree(tree0, tree1, tree0.num_bboxes() - 1,
445 tree1.num_bboxes() - 1, entities);
446 }
447
448 return entities;
449}
450
461template <std::floating_point T>
463compute_collisions(const BoundingBoxTree<T>& tree, std::span<const T> points)
464{
465 if (tree.num_bboxes() > 0)
466 {
467 std::vector<std::int32_t> entities, offsets(points.size() / 3 + 1, 0);
468 entities.reserve(points.size() / 3);
469 for (std::size_t p = 0; p < points.size() / 3; ++p)
470 {
471 impl::_compute_collisions_point(
472 tree, std::span<const T, 3>(points.data() + 3 * p, 3), entities);
473 offsets[p + 1] = entities.size();
474 }
475
476 return graph::AdjacencyList(std::move(entities), std::move(offsets));
477 }
478 else
479 {
481 std::vector<std::int32_t>(),
482 std::vector<std::int32_t>(points.size() / 3 + 1, 0));
483 }
484}
485
503template <std::floating_point T>
505 std::span<const std::int32_t> cells,
506 std::array<T, 3> point, T tol)
507{
508 if (cells.empty())
509 return -1;
510 else
511 {
512 const mesh::Geometry<T>& geometry = mesh.geometry();
513 std::span<const T> geom_dofs = geometry.x();
514 auto x_dofmap = geometry.dofmap();
515 const std::size_t num_nodes = x_dofmap.extent(1);
516 std::vector<T> coordinate_dofs(num_nodes * 3);
517 for (auto cell : cells)
518 {
519 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
520 x_dofmap, cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
521 for (std::size_t i = 0; i < num_nodes; ++i)
522 {
523 std::copy_n(std::next(geom_dofs.begin(), 3 * dofs[i]), 3,
524 std::next(coordinate_dofs.begin(), 3 * i));
525 }
526
527 std::array<T, 3> shortest_vector
528 = compute_distance_gjk<T>(point, coordinate_dofs);
529 T d2 = std::reduce(shortest_vector.begin(), shortest_vector.end(), T(0),
530 [](auto d, auto e) { return d + e * e; });
531 if (d2 < tol)
532 return cell;
533 }
534
535 return -1;
536 }
537}
538
551template <std::floating_point T>
552std::vector<std::int32_t>
554 const BoundingBoxTree<T>& midpoint_tree,
555 const mesh::Mesh<T>& mesh, std::span<const T> points)
556{
557 if (tree.num_bboxes() == 0)
558 return std::vector<std::int32_t>(points.size() / 3, -1);
559
560 std::vector<std::int32_t> entities;
561 entities.reserve(points.size() / 3);
562 for (std::size_t i = 0; i < points.size() / 3; ++i)
563 {
564 // Use midpoint tree to find initial closest entity to the point.
565 // Start by using a leaf node as the initial guess for the input
566 // entity
567 std::array<int, 2> leaf0 = midpoint_tree.bbox(0);
568 assert(impl::is_leaf(leaf0));
569 std::array<T, 6> diff = midpoint_tree.get_bbox(0);
570 for (std::size_t k = 0; k < 3; ++k)
571 diff[k] -= points[3 * i + k];
572 T R2 = diff[0] * diff[0] + diff[1] * diff[1] + diff[2] * diff[2];
573
574 // Use a recursive search through the bounding box tree
575 // to find determine the entity with the closest midpoint.
576 // As the midpoint tree only consist of points, the distance
577 // queries are lightweight.
578 const auto [m_index, m_distance2] = impl::_compute_closest_entity(
579 midpoint_tree, std::span<const T, 3>(points.data() + 3 * i, 3),
580 midpoint_tree.num_bboxes() - 1, mesh, leaf0[0], R2);
581
582 // Use a recursives search through the bounding box tree to
583 // determine which entity is actually closest.
584 // Uses the entity with the closest midpoint as initial guess, and
585 // the distance from the midpoint to the point of interest as the
586 // initial search radius.
587 const auto [index, distance2] = impl::_compute_closest_entity(
588 tree, std::span<const T, 3>(points.data() + 3 * i, 3),
589 tree.num_bboxes() - 1, mesh, m_index, m_distance2);
590
591 entities.push_back(index);
592 }
593
594 return entities;
595}
596
612template <std::floating_point T>
614 const mesh::Mesh<T>& mesh,
615 const graph::AdjacencyList<std::int32_t>& candidate_cells,
616 std::span<const T> points)
617{
618 std::vector<std::int32_t> offsets = {0};
619 offsets.reserve(candidate_cells.num_nodes() + 1);
620 std::vector<std::int32_t> colliding_cells;
621 constexpr T eps2 = 1e-12;
622 const int tdim = mesh.topology()->dim();
623 for (std::int32_t i = 0; i < candidate_cells.num_nodes(); i++)
624 {
625 auto cells = candidate_cells.links(i);
626 std::vector<T> _point(3 * cells.size());
627 for (std::size_t j = 0; j < cells.size(); ++j)
628 for (std::size_t k = 0; k < 3; ++k)
629 _point[3 * j + k] = points[3 * i + k];
630
631 std::vector distances_sq = squared_distance<T>(mesh, tdim, cells, _point);
632 for (std::size_t j = 0; j < cells.size(); j++)
633 if (distances_sq[j] < eps2)
634 colliding_cells.push_back(cells[j]);
635
636 offsets.push_back(colliding_cells.size());
637 }
638
639 return graph::AdjacencyList(std::move(colliding_cells), std::move(offsets));
640}
641
659
670template <std::floating_point T>
671std::tuple<std::vector<std::int32_t>, std::vector<std::int32_t>, std::vector<T>,
672 std::vector<std::int32_t>>
673determine_point_ownership(const mesh::Mesh<T>& mesh, std::span<const T> points,
674 T padding)
675{
676 MPI_Comm comm = mesh.comm();
677
678 // Create a global bounding-box tree to find candidate processes with
679 // cells that could collide with the points
680 const int tdim = mesh.topology()->dim();
681 auto cell_map = mesh.topology()->index_map(tdim);
682 const std::int32_t num_cells = cell_map->size_local();
683 // NOTE: Should we send the cells in as input?
684 std::vector<std::int32_t> cells(num_cells, 0);
685 std::iota(cells.begin(), cells.end(), 0);
686 BoundingBoxTree bb(mesh, tdim, cells, padding);
687 BoundingBoxTree global_bbtree = bb.create_global_tree(comm);
688
689 // Compute collisions:
690 // For each point in `points` get the processes it should be sent to
691 graph::AdjacencyList collisions = compute_collisions(global_bbtree, points);
692
693 // Get unique list of outgoing ranks
694 std::vector<std::int32_t> out_ranks = collisions.array();
695 std::sort(out_ranks.begin(), out_ranks.end());
696 out_ranks.erase(std::unique(out_ranks.begin(), out_ranks.end()),
697 out_ranks.end());
698 // Compute incoming edges (source processes)
699 std::vector in_ranks = dolfinx::MPI::compute_graph_edges_nbx(comm, out_ranks);
700 std::sort(in_ranks.begin(), in_ranks.end());
701
702 // Create neighborhood communicator in forward direction
703 MPI_Comm forward_comm;
704 MPI_Dist_graph_create_adjacent(
705 comm, in_ranks.size(), in_ranks.data(), MPI_UNWEIGHTED, out_ranks.size(),
706 out_ranks.data(), MPI_UNWEIGHTED, MPI_INFO_NULL, false, &forward_comm);
707
708 // Compute map from global mpi rank to neighbor rank, "collisions"
709 // uses global rank
710 std::map<std::int32_t, std::int32_t> rank_to_neighbor;
711 for (std::size_t i = 0; i < out_ranks.size(); i++)
712 rank_to_neighbor[out_ranks[i]] = i;
713
714 // Count the number of points to send per neighbor process
715 std::vector<std::int32_t> send_sizes(out_ranks.size());
716 for (std::size_t i = 0; i < points.size() / 3; ++i)
717 for (auto p : collisions.links(i))
718 send_sizes[rank_to_neighbor[p]] += 3;
719
720 // Compute receive sizes
721 std::vector<std::int32_t> recv_sizes(in_ranks.size());
722 send_sizes.reserve(1);
723 recv_sizes.reserve(1);
724 MPI_Request sizes_request;
725 MPI_Ineighbor_alltoall(send_sizes.data(), 1, MPI_INT, recv_sizes.data(), 1,
726 MPI_INT, forward_comm, &sizes_request);
727
728 // Compute sending offsets
729 std::vector<std::int32_t> send_offsets(send_sizes.size() + 1, 0);
730 std::partial_sum(send_sizes.begin(), send_sizes.end(),
731 std::next(send_offsets.begin(), 1));
732
733 // Pack data to send and store unpack map
734 std::vector<T> send_data(send_offsets.back());
735 std::vector<std::int32_t> counter(send_sizes.size(), 0);
736 // unpack map: [index in adj list][pos in x]
737 std::vector<std::int32_t> unpack_map(send_offsets.back() / 3);
738 for (std::size_t i = 0; i < points.size(); i += 3)
739 {
740 for (auto p : collisions.links(i / 3))
741 {
742 int neighbor = rank_to_neighbor[p];
743 int pos = send_offsets[neighbor] + counter[neighbor];
744 auto it = std::next(send_data.begin(), pos);
745 std::copy_n(std::next(points.begin(), i), 3, it);
746 unpack_map[pos / 3] = i / 3;
747 counter[neighbor] += 3;
748 }
749 }
750
751 MPI_Wait(&sizes_request, MPI_STATUS_IGNORE);
752 std::vector<std::int32_t> recv_offsets(in_ranks.size() + 1, 0);
753 std::partial_sum(recv_sizes.begin(), recv_sizes.end(),
754 std::next(recv_offsets.begin(), 1));
755
756 std::vector<T> received_points((std::size_t)recv_offsets.back());
757 MPI_Neighbor_alltoallv(
758 send_data.data(), send_sizes.data(), send_offsets.data(),
759 dolfinx::MPI::mpi_type<T>(), received_points.data(), recv_sizes.data(),
760 recv_offsets.data(), dolfinx::MPI::mpi_type<T>(), forward_comm);
761
762 // Get mesh geometry for closest entity
763 const mesh::Geometry<T>& geometry = mesh.geometry();
764 std::span<const T> geom_dofs = geometry.x();
765 auto x_dofmap = geometry.dofmap();
766
767 // Compute candidate cells for collisions (and extrapolation)
768 const graph::AdjacencyList<std::int32_t> candidate_collisions
769 = compute_collisions(bb, std::span<const T>(received_points.data(),
770 received_points.size()));
771
772 // Each process checks which points collide with a cell on the process
773 const int rank = dolfinx::MPI::rank(comm);
774 std::vector<std::int32_t> cell_indicator(received_points.size() / 3);
775 std::vector<std::int32_t> closest_cells(received_points.size() / 3);
776 for (std::size_t p = 0; p < received_points.size(); p += 3)
777 {
778 std::array<T, 3> point;
779 std::copy_n(std::next(received_points.begin(), p), 3, point.begin());
780 // Find first colliding cell among the cells with colliding bounding boxes
781 const int colliding_cell = geometry::compute_first_colliding_cell(
782 mesh, candidate_collisions.links(p / 3), point,
783 10 * std::numeric_limits<T>::epsilon());
784 // If a collding cell is found, store the rank of the current process
785 // which will be sent back to the owner of the point
786 cell_indicator[p / 3] = (colliding_cell >= 0) ? rank : -1;
787 // Store the cell index for lookup once the owning processes has determined
788 // the ownership of the point
789 closest_cells[p / 3] = colliding_cell;
790 }
791
792 // Create neighborhood communicator in the reverse direction: send
793 // back col to requesting processes
794 MPI_Comm reverse_comm;
795 MPI_Dist_graph_create_adjacent(
796 comm, out_ranks.size(), out_ranks.data(), MPI_UNWEIGHTED, in_ranks.size(),
797 in_ranks.data(), MPI_UNWEIGHTED, MPI_INFO_NULL, false, &reverse_comm);
798
799 // Reuse sizes and offsets from first communication set
800 // but divide by three
801 {
802 auto rescale = [](auto& x)
803 {
804 std::transform(x.cbegin(), x.cend(), x.begin(),
805 [](auto e) { return (e / 3); });
806 };
807 rescale(recv_sizes);
808 rescale(recv_offsets);
809 rescale(send_sizes);
810 rescale(send_offsets);
811
812 // The communication is reversed, so swap recv to send offsets
813 std::swap(recv_sizes, send_sizes);
814 std::swap(recv_offsets, send_offsets);
815 }
816
817 std::vector<std::int32_t> recv_ranks(recv_offsets.back());
818 MPI_Neighbor_alltoallv(cell_indicator.data(), send_sizes.data(),
819 send_offsets.data(), MPI_INT32_T, recv_ranks.data(),
820 recv_sizes.data(), recv_offsets.data(), MPI_INT32_T,
821 reverse_comm);
822
823 std::vector<std::int32_t> point_owners(points.size() / 3, -1);
824 for (std::size_t i = 0; i < unpack_map.size(); i++)
825 {
826 const std::int32_t pos = unpack_map[i];
827 // Only insert new owner if no owner has previously been found
828 if (recv_ranks[i] >= 0 && point_owners[pos] == -1)
829 point_owners[pos] = recv_ranks[i];
830 }
831
832 // Create extrapolation marker for those points already sent to other
833 // process
834 std::vector<std::uint8_t> send_extrapolate(recv_offsets.back());
835 for (std::int32_t i = 0; i < recv_offsets.back(); i++)
836 {
837 const std::int32_t pos = unpack_map[i];
838 send_extrapolate[i] = point_owners[pos] == -1;
839 }
840
841 // Swap communication direction, to send extrapolation marker to other
842 // processes
843 std::swap(send_sizes, recv_sizes);
844 std::swap(send_offsets, recv_offsets);
845 std::vector<std::uint8_t> dest_extrapolate(recv_offsets.back());
846 MPI_Neighbor_alltoallv(send_extrapolate.data(), send_sizes.data(),
847 send_offsets.data(), MPI_UINT8_T,
848 dest_extrapolate.data(), recv_sizes.data(),
849 recv_offsets.data(), MPI_UINT8_T, forward_comm);
850
851 std::vector<T> squared_distances(received_points.size() / 3, -1);
852
853 for (std::size_t i = 0; i < dest_extrapolate.size(); i++)
854 {
855 if (dest_extrapolate[i] == 1)
856 {
857 assert(closest_cells[i] == -1);
858 std::array<T, 3> point;
859 std::copy_n(std::next(received_points.begin(), 3 * i), 3, point.begin());
860
861 // Find shortest distance among cells with colldiing bounding box
862 T shortest_distance = std::numeric_limits<T>::max();
863 std::int32_t closest_cell = -1;
864 for (auto cell : candidate_collisions.links(i))
865 {
866 auto dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
867 x_dofmap, cell, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
868 std::vector<T> nodes(3 * dofs.size());
869 for (std::size_t j = 0; j < dofs.size(); ++j)
870 {
871 const int pos = 3 * dofs[j];
872 for (std::size_t k = 0; k < 3; ++k)
873 nodes[3 * j + k] = geom_dofs[pos + k];
874 }
875 const std::array<T, 3> d = compute_distance_gjk<T>(
876 std::span<const T>(point.data(), point.size()), nodes);
877 if (T current_distance = d[0] * d[0] + d[1] * d[1] + d[2] * d[2];
878 current_distance < shortest_distance)
879 {
880 shortest_distance = current_distance;
881 closest_cell = cell;
882 }
883 }
884 closest_cells[i] = closest_cell;
885 squared_distances[i] = shortest_distance;
886 }
887 }
888
889 std::swap(recv_sizes, send_sizes);
890 std::swap(recv_offsets, send_offsets);
891
892 // Get distances from closest entity of points that were on the other process
893 std::vector<T> recv_distances(recv_offsets.back());
894 MPI_Neighbor_alltoallv(
895 squared_distances.data(), send_sizes.data(), send_offsets.data(),
896 dolfinx::MPI::mpi_type<T>(), recv_distances.data(), recv_sizes.data(),
897 recv_offsets.data(), dolfinx::MPI::mpi_type<T>(), reverse_comm);
898
899 // Update point ownership with extrapolation information
900 std::vector<T> closest_distance(point_owners.size(),
901 std::numeric_limits<T>::max());
902 for (std::size_t i = 0; i < out_ranks.size(); i++)
903 {
904 for (std::int32_t j = recv_offsets[i]; j < recv_offsets[i + 1]; j++)
905 {
906 const std::int32_t pos = unpack_map[j];
907 auto current_dist = recv_distances[j];
908 // Update if closer than previous guess and was found
909 if (auto d = closest_distance[pos];
910 (current_dist > 0) and (current_dist < d))
911 {
912 point_owners[pos] = out_ranks[i];
913 closest_distance[pos] = current_dist;
914 }
915 }
916 }
917
918 // Communication is reversed again to send dest ranks to all processes
919 std::swap(send_sizes, recv_sizes);
920 std::swap(send_offsets, recv_offsets);
921
922 // Pack ownership data
923 std::vector<std::int32_t> send_owners(send_offsets.back());
924 std::fill(counter.begin(), counter.end(), 0);
925 for (std::size_t i = 0; i < points.size() / 3; ++i)
926 {
927 for (auto p : collisions.links(i))
928 {
929 int neighbor = rank_to_neighbor[p];
930 send_owners[send_offsets[neighbor] + counter[neighbor]++]
931 = point_owners[i];
932 }
933 }
934
935 // Send ownership info
936 std::vector<std::int32_t> dest_ranks(recv_offsets.back());
937 MPI_Neighbor_alltoallv(send_owners.data(), send_sizes.data(),
938 send_offsets.data(), MPI_INT32_T, dest_ranks.data(),
939 recv_sizes.data(), recv_offsets.data(), MPI_INT32_T,
940 forward_comm);
941
942 // Unpack dest ranks if point owner is this rank
943 std::vector<std::int32_t> owned_recv_ranks;
944 owned_recv_ranks.reserve(recv_offsets.back());
945 std::vector<T> owned_recv_points;
946 std::vector<std::int32_t> owned_recv_cells;
947 for (std::size_t i = 0; i < in_ranks.size(); i++)
948 {
949 for (std::int32_t j = recv_offsets[i]; j < recv_offsets[i + 1]; j++)
950 {
951 if (rank == dest_ranks[j])
952 {
953 owned_recv_ranks.push_back(in_ranks[i]);
954 owned_recv_points.insert(
955 owned_recv_points.end(), std::next(received_points.cbegin(), 3 * j),
956 std::next(received_points.cbegin(), 3 * (j + 1)));
957 owned_recv_cells.push_back(closest_cells[j]);
958 }
959 }
960 }
961
962 MPI_Comm_free(&forward_comm);
963 MPI_Comm_free(&reverse_comm);
964
965 return std::make_tuple(point_owners, owned_recv_ranks, owned_recv_points,
966 owned_recv_cells);
967}
968
969} // namespace dolfinx::geometry
Definition BoundingBoxTree.h:205
BoundingBoxTree create_global_tree(MPI_Comm comm) const
Definition BoundingBoxTree.h:325
std::int32_t num_bboxes() const
Return number of bounding boxes.
Definition BoundingBoxTree.h:361
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:313
std::array< std::int32_t, 2 > bbox(std::size_t node) const
Definition BoundingBoxTree.h:381
Definition AdjacencyList.h:28
std::span< T > links(std::size_t node)
Definition AdjacencyList.h:112
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
const fem::CoordinateElement< value_type > & cmap() const
The element that describes the geometry map.
Definition Geometry.h:180
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
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
Geometry data structures and algorithms.
Definition BoundingBoxTree.h:21
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:158
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:437
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:504
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:121
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:613
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:36
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:410
std::tuple< std::vector< std::int32_t >, std::vector< std::int32_t >, std::vector< T >, std::vector< std::int32_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:673
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:553
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