Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.1.0/v0.9.0/cpp
DOLFINx  0.1.0
DOLFINx 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

xt::xtensor_fixed< double, xt::xshape< 3 > > compute_distance_gjk (const xt::xtensor< double, 2 > &p, const xt::xtensor< double, 2 > &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, const xtl::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::array< int, 2 > > compute_collisions (const BoundingBoxTree &tree0, const BoundingBoxTree &tree1)
 Compute all collisions between two BoundingBoxTrees (local to process). More...
 
std::vector< int > compute_collisions (const BoundingBoxTree &tree, const std::array< double, 3 > &p)
 Compute all collisions between bounding boxes and point. More...
 
std::pair< std::int32_t, double > compute_closest_entity (const BoundingBoxTree &tree, const std::array< double, 3 > &p, const mesh::Mesh &mesh, double R=-1)
 Compute closest mesh entity (local to process) for the topological distance of the bounding box tree and distance and a point. More...
 
double compute_squared_distance_bbox (const std::array< std::array< double, 3 >, 2 > &b, const std::array< double, 3 > &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 std::array< double, 3 > &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< std::int32_t > select_colliding_cells (const dolfinx::mesh::Mesh &mesh, const xtl::span< const std::int32_t > &candidate_cells, const std::array< double, 3 > &p, int n)
 From the given Mesh, select up to n cells (local to process) 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_entity()

std::pair< int, double > dolfinx::geometry::compute_closest_entity ( const BoundingBoxTree tree,
const std::array< double, 3 > &  p,
const mesh::Mesh mesh,
double  R = -1 
)

Compute closest mesh entity (local to process) for the topological distance of the bounding box tree and distance and a point.

Parameters
[in]treeThe bounding box tree
[in]pThe point
[in]meshThe mesh
[in]RRadius for search. Supplying a negative radius causes the function to estimate an initial search radius.
Returns
The local index of the entity and the distance from the point.

◆ compute_collisions() [1/2]

std::vector< int > dolfinx::geometry::compute_collisions ( const BoundingBoxTree tree,
const std::array< double, 3 > &  p 
)

Compute all collisions between bounding boxes and point.

Parameters
[in]treeThe bounding box tree
[in]pThe point
Returns
Bounding box leaves (local to process) 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 (local to process).

Parameters
[in]tree0First BoundingBoxTree
[in]tree1Second BoundingBoxTree
Returns
List of pairs of intersecting box indices (local to process) from each tree

◆ compute_distance_gjk()

xt::xtensor_fixed< double, xt::xshape< 3 > > dolfinx::geometry::compute_distance_gjk ( const xt::xtensor< double, 2 > &  p,
const xt::xtensor< double, 2 > &  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)
[in]qBody 2 list of points, shape (num_points, 3)
Returns
shortest vector between bodies

◆ create_midpoint_tree()

geometry::BoundingBoxTree dolfinx::geometry::create_midpoint_tree ( const mesh::Mesh mesh,
int  tdim,
const xtl::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

◆ select_colliding_cells()

std::vector< std::int32_t > dolfinx::geometry::select_colliding_cells ( const dolfinx::mesh::Mesh mesh,
const xtl::span< const std::int32_t > &  candidate_cells,
const std::array< double, 3 > &  p,
int  n 
)

From the given Mesh, select up to n cells (local to process) 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]pPoint 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 std::array< double, 3 > &  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 (local to process)
[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