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

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. More...
 
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. More...
 
std::vector< std::int32_t > compute_collisions (const BoundingBoxTree &tree0, const BoundingBoxTree &tree1)
 Compute all collisions between two BoundingBoxTrees (local to process) More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 

Detailed Description

Geometry data structures and algorithms.

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

Function Documentation

◆ compute_closest_entity()

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.

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
Index of the closest mesh entity to a point. The ith entry is the index of the closest entity to the ith input point.
Note
Returns index -1 if the bounding box tree is empty

◆ compute_colliding_cells()

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
Parameters
[in]meshThe mesh
[in]candidate_cellsList of candidate colliding cells for the ith point in points
[in]pointsThe 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
Note
There may be nodes with no entries in the adjacency list

◆ compute_collisions() [1/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
[in]treeThe bounding box tree
[in]pointsThe 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

◆ compute_collisions() [2/2]

std::vector< std::int32_t > compute_collisions ( const BoundingBoxTree tree0,
const BoundingBoxTree tree1 
)

Compute all collisions between two BoundingBoxTrees (local to process)

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()

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
[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()

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
[in]meshThe mesh
[in]treeThe bounding box tree
[in]pointThe point (shape=(3))
Returns
The local cell index, -1 if not found

◆ compute_squared_distance_bbox()

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
[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()

geometry::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
[in]meshThe mesh
[in]tdimThe topological dimension of the entity
[in]entity_indicesList of local entity indices
Returns
Bounding box tree for midpoints of mesh entities

◆ determine_point_ownership()

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.

Parameters
[in]meshThe mesh
[in]pointsThe 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.
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

◆ shortest_vector()

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
[in]meshThe mesh
[in]dimThe topological dimension of the mesh entity
[in]entitiesThe list of entities (local to process)
[in]pointsThe 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.

◆ squared_distance()

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