Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/de/d27/namespacedolfinx_1_1geometry.html
DOLFINx 0.7.3
DOLFINx C++ interface
Loading...
Searching...
No Matches
Classes | Functions
dolfinx::geometry Namespace Reference

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>
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.
 

Detailed Description

Geometry data structures and algorithms.

Tools for geometric data structures and operations, e.g. searching.

Function Documentation

◆ compute_closest_entity()

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.

Note
Returns a vector filled with index -1 if the bounding box tree is empty.
Parameters
[in]treeThe bounding box tree for the entities
[in]midpoint_treeA bounding box tree with the midpoints of all the mesh entities. This is used to accelerate the search.
[in]meshThe mesh
[in]pointsThe set of points (shape=(num_points, 3)). Storage is row-major.
Returns
For each point, the index of the closest mesh entity.

◆ compute_colliding_cells()

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.

Note
Uses the GJK algorithm, see geometry::compute_distance_gjk for details.
candidate_cells can for instance be found by using dolfinx::geometry::compute_collisions between a bounding box tree and the set of points.
Parameters
[in]meshThe mesh
[in]candidate_cellsList of candidate colliding cells for the ith point in points
[in]pointsPoints to check for collision (shape=(num_points, 3)). Storage is row-major.
Returns
For each point, the cells that collide with the point.

◆ compute_collisions() [1/2]

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.

Bounding boxes can overlap, therefore points can collide with more than one box.

Parameters
[in]treeThe bounding box tree
[in]pointsThe points (shape=(num_points, 3)). Storage is row-major.
Returns
For each point, the bounding box leaves that collide with the point.

◆ compute_collisions() [2/2]

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.

Parameters
[in]tree0First BoundingBoxTree
[in]tree1Second BoundingBoxTree
Returns
List of pairs of intersecting box indices from each tree, flattened as a vector of size num_intersections*2

◆ compute_distance_gjk()

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.

Uses the Gilbert–Johnson–Keerthi (GJK) distance algorithm.

Parameters
[in]pBody 1 list of points, shape (num_points, 3). Row-major storage.
[in]qBody 2 list of points, shape (num_points, 3). Row-major storage.
Returns
shortest vector between bodies

◆ compute_first_colliding_cell()

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,
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.

Note
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.
Parameters
[in]meshThe mesh
[in]cellsThe candidate cells
[in]pointThe point (shape=(3,))
[in]tolTolerance for accepting a collision (in the squared distance)
Returns
The local cell index, -1 if not found

◆ compute_squared_distance_bbox()

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.

Parameters
[in]bBounding box coordinates
[in]xA point
Returns
The shortest distance between the bounding box b and the point x. Returns zero if x is inside box.

◆ create_midpoint_tree()

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.

Parameters
[in]meshThe mesh
[in]tdimThe topological dimension of the entity
[in]entitiesList of local entity indices
Returns
Bounding box tree for midpoints of entities

◆ determine_point_ownership()

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,
padding 
)

Given a set of points, determine which process is colliding, using the GJK algorithm on cells to determine collisions.

Todo:
This docstring is unclear. Needs fixing.
Parameters
[in]meshThe mesh
[in]pointsPoints to check for collision (shape=(num_points, 3)). Storage is row-major.
[in]paddingAmount 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.
Returns
Tuple (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.
Note
dest_owner is sorted
Returns -1 if no colliding process is found
dest_points is flattened row-major, shape (dest_owner.size(), 3)
Only looks through cells owned by the process
A large padding value can increase the runtime of the function by orders of magnitude, because that for non-colliding cells one has to determine the closest cell among all processes with an intersecting bounding box, which is an expensive operation to perform.

◆ shortest_vector()

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.

Parameters
[in]meshThe mesh
[in]dimTopological dimension of the mesh entity
[in]entitiesList of entities
[in]pointsSet 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.

◆ squared_distance()

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.

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.
Uses a convex hull approximation of linearized geometry
Parameters
[in]meshMesh containing the entities
[in]dimThe topological dimension of the mesh entities
[in]entitiesThe indices of the mesh entities (local to process)
[in]pointsThe 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]