DOLFIN-X
DOLFIN-X C++ interface
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

Eigen::Vector3d compute_distance_gjk (const Eigen::Matrix< double, Eigen::Dynamic, 3, Eigen::RowMajor > &p, const Eigen::Matrix< double, Eigen::Dynamic, 3, Eigen::RowMajor > &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)
 Create a boundary box tree for cell midpoints. More...
 
std::vector< std::array< int, 2 > > compute_collisions (const BoundingBoxTree &tree0, const BoundingBoxTree &tree1)
 Compute all collisions between two BoundingBoxTrees. More...
 
std::vector< int > compute_collisions (const BoundingBoxTree &tree, const Eigen::Vector3d &p)
 Compute all collisions between bounding boxes and point. More...
 
std::vector< int > compute_process_collisions (const BoundingBoxTree &tree, const Eigen::Vector3d &p)
 Compute all collisions between processes and Point returning a list of process ranks.
 
std::pair< int, double > compute_closest_entity (const BoundingBoxTree &tree, const BoundingBoxTree &tree_midpoint, const Eigen::Vector3d &p, const mesh::Mesh &mesh)
 Compute closest mesh entity and distance to the point. The tree must have been initialised with topological co-dimension 0.
 
std::pair< int, double > compute_closest_point (const BoundingBoxTree &tree, const Eigen::Vector3d &p)
 Compute closest point and distance to a given point. More...
 
double compute_squared_distance_bbox (const Eigen::Array< double, 2, 3, Eigen::RowMajor > &b, const Eigen::Vector3d &x)
 Compute squared distance between point and bounding box wih index "node". Returns zero if point is inside box.
 
double squared_distance (const mesh::Mesh &mesh, int dim, std::int32_t index, const Eigen::Vector3d &p)
 Compute squared distance from a given point to the nearest point on a cell (only first order convex cells are supported at this stage) Uses the GJK algorithm, see geometry::compute_distance_gjk for details. More...
 
std::vector< int > select_colliding_cells (const dolfinx::mesh::Mesh &mesh, const std::vector< int > &candidate_cells, const Eigen::Vector3d &point, int n)
 From the given Mesh, select up to n cells from the list which actually collide with point p. n may be zero (selects all valid cells). Less than n cells may be returned. More...
 

Detailed Description

Geometry data structures and algorithms.

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

Function Documentation

◆ compute_closest_point()

std::pair< int, double > dolfinx::geometry::compute_closest_point ( const BoundingBoxTree tree,
const Eigen::Vector3d &  p 
)

Compute closest point and distance to a given point.

Parameters
[in]treeThe bounding box tree. It must have been initialised with topological dimension 0.
[in]pThe point to compute the distance from
Returns
(point index, distance)

◆ compute_collisions() [1/2]

std::vector< int > dolfinx::geometry::compute_collisions ( const BoundingBoxTree tree,
const Eigen::Vector3d &  p 
)

Compute all collisions between bounding boxes and point.

Parameters
[in]treeThe bounding box tree
[in]pThe point
Returns
Bounding box leaves that contain the point

◆ compute_collisions() [2/2]

std::vector< std::array< int, 2 > > dolfinx::geometry::compute_collisions ( const BoundingBoxTree tree0,
const BoundingBoxTree tree1 
)

Compute all collisions between two BoundingBoxTrees.

Parameters
[in]tree0First BoundingBoxTree
[in]tree1Second BoundingBoxTree
Returns
List of pairs of intersecting box indices from each tree

◆ compute_distance_gjk()

Eigen::Vector3d dolfinx::geometry::compute_distance_gjk ( const Eigen::Matrix< double, Eigen::Dynamic, 3, Eigen::RowMajor > &  p,
const Eigen::Matrix< double, Eigen::Dynamic, 3, Eigen::RowMajor > &  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
[in]qBody 2 list of points
Returns
shortest vector between bodies

◆ create_midpoint_tree()

geometry::BoundingBoxTree dolfinx::geometry::create_midpoint_tree ( const mesh::Mesh mesh)

Create a boundary box tree for cell midpoints.

Parameters
[in]meshThe mesh build tree of cell midpoints from
Returns
Bounding box tree for mesh cell midpoints

◆ select_colliding_cells()

std::vector< int > dolfinx::geometry::select_colliding_cells ( const dolfinx::mesh::Mesh mesh,
const std::vector< int > &  candidate_cells,
const Eigen::Vector3d &  point,
int  n 
)

From the given Mesh, select up to n cells from the list which actually collide with point p. n may be zero (selects all valid cells). Less than n cells may be returned.

Parameters
[in]meshMesh
[in]candidate_cellsList of cell indices to test
[in]pointPoint to check for collision
[in]nMaximum number of positive results to return
Returns
List of cells which collide with point

◆ squared_distance()

double dolfinx::geometry::squared_distance ( const mesh::Mesh mesh,
int  dim,
std::int32_t  index,
const Eigen::Vector3d &  p 
)

Compute squared distance from a given point to the nearest point on a cell (only first order convex cells are supported at this stage) Uses the GJK algorithm, see geometry::compute_distance_gjk for details.

Note
Currently a convex hull approximation of linearized geometry.
Parameters
[in]meshMesh containing the mesh entity
[in]dimThe topological dimension of the mesh entity
[in]indexThe index of the mesh entity
[in]pThe point from which to compouted the shortest distance to the mesh to compute the Point
Returns
shortest squared distance from p to entity