DOLFINx 0.7.3
DOLFINx C++ interface
|
Geometry data structures and algorithms. More...
Classes | |
class | BoundingBoxTree |
Axis-Aligned bounding box binary tree. It is used to find entities in a collection (often a mesh::Mesh). More... | |
Functions | |
template<std::floating_point T> | |
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. | |
template<std::floating_point T> | |
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. | |
template<std::floating_point T> | |
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. | |
template<std::floating_point T> | |
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. | |
template<std::floating_point T> | |
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. | |
template<std::floating_point T> | |
std::vector< std::int32_t > | compute_collisions (const BoundingBoxTree< T > &tree0, const BoundingBoxTree< T > &tree1) |
Compute all collisions between two bounding box trees. | |
template<std::floating_point T> | |
graph::AdjacencyList< std::int32_t > | compute_collisions (const BoundingBoxTree< T > &tree, std::span< const T > points) |
Compute collisions between points and leaf bounding boxes. | |
template<std::floating_point T> | |
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. | |
template<std::floating_point T> | |
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. | |
template<std::floating_point T> | |
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. | |
template<std::floating_point T> | |
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 determine collisions. | |
Geometry data structures and algorithms.
Tools for geometric data structures and operations, e.g. searching.
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.
[in] | tree | The bounding box tree for the entities |
[in] | midpoint_tree | A bounding box tree with the midpoints of all the mesh entities. This is used to accelerate the search. |
[in] | mesh | The mesh |
[in] | points | The set of points (shape=(num_points, 3) ). Storage is row-major. |
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.
candidate_cells
can for instance be found by using dolfinx::geometry::compute_collisions
between a bounding box tree and the set of points.[in] | mesh | The mesh |
[in] | candidate_cells | List of candidate colliding cells for the ith point in points |
[in] | points | Points to check for collision (shape=(num_points, 3) ). Storage is row-major. |
graph::AdjacencyList< std::int32_t > compute_collisions | ( | const BoundingBoxTree< T > & | tree, |
std::span< const T > | points | ||
) |
Compute collisions between points and leaf bounding boxes.
Bounding boxes can overlap, therefore points can collide with more than one box.
[in] | tree | The bounding box tree |
[in] | points | The points (shape=(num_points, 3) ). Storage is row-major. |
std::vector< std::int32_t > compute_collisions | ( | const BoundingBoxTree< T > & | tree0, |
const BoundingBoxTree< T > & | tree1 | ||
) |
Compute all collisions between two bounding box trees.
[in] | tree0 | First BoundingBoxTree |
[in] | tree1 | Second BoundingBoxTree |
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.
Uses the Gilbert–Johnson–Keerthi (GJK) distance algorithm.
[in] | p | Body 1 list of points, shape (num_points, 3). Row-major storage. |
[in] | q | Body 2 list of points, shape (num_points, 3). Row-major storage. |
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.
A point can collide with more than one cell. The first cell detected to collide with the point is returned. If no collision is detected, -1 is returned.
cells
can for instance be found by using dolfinx::geometry::compute_collisions
between a bounding box tree for the cells of the mesh and the point.[in] | mesh | The mesh |
[in] | cells | The candidate cells |
[in] | point | The point (shape=(3,) ) |
[in] | tol | Tolerance for accepting a collision (in the squared distance) |
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.
[in] | b | Bounding box coordinates |
[in] | x | A point |
b
and the point x
. Returns zero if x
is inside box. 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.
[in] | mesh | The mesh |
[in] | tdim | The topological dimension of the entity |
[in] | entities | List of local entity indices |
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 determine collisions.
[in] | mesh | The mesh |
[in] | points | Points to check for collision (shape=(num_points, 3) ). Storage is row-major. |
[in] | padding | Amount of absolute padding of bounding boxes of the mesh. Each bounding box of the mesh is padded with this amount, to increase the number of candidates, avoiding rounding errors in determining the owner of a point if the point is on the surface of a cell in the mesh. |
(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. dest_owner
is sorted (dest_owner.size(), 3)
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.
[in] | mesh | The mesh |
[in] | dim | Topological dimension of the mesh entity |
[in] | entities | List of entities |
[in] | points | Set of points (shape=(num_points, 3) ), using row-major storage. |
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.
The distance is computed between the ith input points and the ith input entity.
[in] | mesh | Mesh containing the entities |
[in] | dim | The topological dimension of the mesh entities |
[in] | entities | The indices of the mesh entities (local to process) |
[in] | points | The set points from which to computed the shortest (shape=(num_points, 3)). Storage is row-major. |