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.6.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 <array>
10#include <cstdint>
11#include <dolfinx/graph/AdjacencyList.h>
12#include <span>
13#include <vector>
14
15namespace dolfinx::mesh
16{
17class Mesh;
18}
19
20namespace dolfinx::geometry
21{
22class BoundingBoxTree;
23
30BoundingBoxTree
31create_midpoint_tree(const mesh::Mesh& mesh, int tdim,
32 std::span<const std::int32_t> entity_indices);
33
40std::vector<std::int32_t> compute_collisions(const BoundingBoxTree& tree0,
41 const BoundingBoxTree& tree1);
42
50graph::AdjacencyList<std::int32_t>
51compute_collisions(const BoundingBoxTree& tree, std::span<const double> points);
52
59int compute_first_colliding_cell(const mesh::Mesh& mesh,
60 const BoundingBoxTree& tree,
61 const std::array<double, 3>& point);
62
73std::vector<std::int32_t>
74compute_closest_entity(const BoundingBoxTree& tree,
75 const BoundingBoxTree& midpoint_tree,
76 const mesh::Mesh& mesh, std::span<const double> points);
77
83double compute_squared_distance_bbox(std::span<const double, 6> b,
84 std::span<const double, 3> x);
85
95std::vector<double> shortest_vector(const mesh::Mesh& mesh, int dim,
96 std::span<const std::int32_t> entities,
97 std::span<const double> points);
98
111std::vector<double> squared_distance(const mesh::Mesh& mesh, int dim,
112 std::span<const std::int32_t> entities,
113 std::span<const double> points);
114
126graph::AdjacencyList<std::int32_t> compute_colliding_cells(
127 const mesh::Mesh& mesh,
128 const graph::AdjacencyList<std::int32_t>& candidate_cells,
129 std::span<const double> points);
130
145std::tuple<std::vector<std::int32_t>, std::vector<std::int32_t>,
146 std::vector<double>, std::vector<std::int32_t>>
147determine_point_ownership(const mesh::Mesh& mesh,
148 std::span<const double> points);
149
150} // namespace dolfinx::geometry
Geometry data structures and algorithms.
Definition: BoundingBoxTree.h:24
std::vector< double > squared_distance(const mesh::Mesh &mesh, int dim, std::span< const std::int32_t > entities, std::span< const double > points)
Compute the squared distance between a point and a mesh entity. The distance is computed between the ...
Definition: utils.cpp:462
std::vector< double > shortest_vector(const mesh::Mesh &mesh, int dim, std::span< const std::int32_t > entities, std::span< const double > points)
Compute the shortest vector from a mesh entity to a point.
Definition: utils.cpp:390
std::vector< std::int32_t > compute_collisions(const BoundingBoxTree &tree0, const BoundingBoxTree &tree1)
Compute all collisions between two BoundingBoxTrees (local to process)
Definition: utils.cpp:276
std::tuple< std::vector< std::int32_t >, std::vector< std::int32_t >, std::vector< double >, std::vector< std::int32_t > > determine_point_ownership(const mesh::Mesh &mesh, std::span< const double > points)
Given a set of points (local on each process) which process is colliding, using the GJK algorithm on ...
Definition: utils.cpp:548
std::vector< std::int32_t > compute_closest_entity(const BoundingBoxTree &tree, const BoundingBoxTree &midpoint_tree, const mesh::Mesh &mesh, std::span< const double > points)
Compute closest mesh entity to a point.
Definition: utils.cpp:316
graph::AdjacencyList< std::int32_t > compute_colliding_cells(const mesh::Mesh &mesh, const graph::AdjacencyList< std::int32_t > &candidate_cells, std::span< const double > points)
From a Mesh, find which cells collide with a set of points.
Definition: utils.cpp:474
BoundingBoxTree create_midpoint_tree(const mesh::Mesh &mesh, int tdim, std::span< const std::int32_t > entity_indices)
Create a bounding box tree for a subset of entities (local to process) based on the entity midpoints.
Definition: utils.cpp:254
int compute_first_colliding_cell(const mesh::Mesh &mesh, const BoundingBoxTree &tree, const std::array< double, 3 > &point)
Compute the first collision between a point and the cells of the mesh.
Definition: utils.cpp:505
double compute_squared_distance_bbox(std::span< const double, 6 > b, std::span< const double, 3 > x)
Compute squared distance between point and bounding box.
Definition: utils.cpp:367
Mesh data structures and algorithms on meshes.
Definition: DofMap.h:31