Geometry (dolfinx::geometry
)
-
namespace geometry
Geometry data structures and algorithms.
Tools for geometric data structures and operations, e.g. searching.
Functions
-
std::array<double, 3> compute_distance_gjk(std::span<const double> p, std::span<const double> q)
Calculate the distance between two convex bodies p and q, each defined by a set of points, using the Gilbert–Johnson–Keerthi (GJK) distance algorithm.
- Parameters
p – [in] Body 1 list of points, shape (num_points, 3). Row-major storage.
q – [in] Body 2 list of points, shape (num_points, 3). Row-major storage.
- Returns
shortest vector between bodies
-
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.
- Parameters
mesh – [in] The mesh
tdim – [in] The topological dimension of the entity
entity_indices – [in] List of local entity indices
- Returns
Bounding box tree for midpoints of mesh entities
-
std::vector<std::int32_t> compute_collisions(const BoundingBoxTree &tree0, const BoundingBoxTree &tree1)
Compute all collisions between two BoundingBoxTrees (local to process)
- Parameters
tree0 – [in] First BoundingBoxTree
tree1 – [in] Second BoundingBoxTree
- Returns
List of pairs of intersecting box indices from each tree, flattened as a vector of size num_intersections*2
-
graph::AdjacencyList<std::int32_t> compute_collisions(const BoundingBoxTree &tree, std::span<const double> points)
Compute all collisions between bounding boxes and for a set of points.
- Parameters
tree – [in] The bounding box tree
points – [in] The points (shape=(num_points, 3)). Storage is row-major.
- Returns
An adjacency list where the ith node corresponds to the bounding box leaves that contain the ith point
-
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.
- Parameters
mesh – [in] The mesh
tree – [in] The bounding box tree
point – [in] The point (shape=(3))
- Returns
The local cell index, -1 if not found
-
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.
Note
Returns index -1 if the bounding box tree is empty
- Parameters
tree – [in] The bounding box tree for the entities
midpoint_tree – [in] A bounding box tree with the midpoints of all the mesh entities. This is used to accelerate the search
mesh – [in] The mesh
points – [in] The set of points (shape=(num_points, 3)). Storage is row-major.
- Returns
Index of the closest mesh entity to a point. The ith entry is the index of the closest entity to the ith input point.
-
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.
- Parameters
b – [in] Bounding box coordinates
x – [in] A point
- Returns
The shortest distance between the bounding box
b
and the pointx
. Returns zero ifx
is inside box.
-
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.
- Parameters
mesh – [in] The mesh
dim – [in] The topological dimension of the mesh entity
entities – [in] The list of entities (local to process)
points – [in] The set of points (shape=(num_points, 3)), using row-major storage
- Returns
An array of vectors (shape=(num_points, 3)) where the ith row is the shortest vector between the ith entity and the ith point. Storage is row-major.
-
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 ith input points and the ith input entity.
Note
Uses the GJK algorithm, see geometry::compute_distance_gjk for details.
Note
Uses a convex hull approximation of linearized geometry
- Parameters
mesh – [in] Mesh containing the entities
dim – [in] The topological dimension of the mesh entities
entities – [in] The indices of the mesh entities (local to process)
points – [in] The set points from which to computed the shortest (shape=(num_points, 3)). Storage is row-major.
- Returns
Squared shortest distance from points[i] to entities[i]
-
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.
Note
Uses the GJK algorithm, see geometry::compute_distance_gjk for details
Note
There may be nodes with no entries in the adjacency list
- Parameters
mesh – [in] The mesh
candidate_cells – [in] List of candidate colliding cells for the ith point in
points
points – [in] The points to check for collision (shape=(num_points, 3)). Storage is row-major.
- Returns
Adjacency list where the ith node is the list of entities that collide with the ith point
-
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 cells to determine collisions.
Note
dest_owner is sorted
Note
Returns -1 if no colliding process is found
Note
dest_points is flattened row-major, shape (dest_owner.size(), 3)
Note
Only looks through cells owned by the process
- Parameters
mesh – [in] The mesh
points – [in] The points to check for collision (shape=(num_points, 3)). Storage is row-major.
- Returns
Quadratuplet (src_owner, dest_owner, dest_points, dest_cells), where src_owner is a list of ranks corresponding to the input points. dest_owner is a list of ranks corresponding to dest_points, the points that this process owns. dest_cells contains the corresponding cell for each entry in dest_points.
-
class BoundingBoxTree
- #include <BoundingBoxTree.h>
Axis-Aligned bounding box binary tree. It is used to find entities in a collection (often a mesh::Mesh).
-
std::array<double, 3> compute_distance_gjk(std::span<const double> p, std::span<const double> q)