DOLFIN-X
DOLFIN-X 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 | |
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... | |
Geometry data structures and algorithms.
Tools for geometric data structures and operations, e.g. searching.
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.
[in] | tree | The bounding box tree. It must have been initialised with topological dimension 0. |
[in] | p | The point to compute the distance from |
std::vector< int > dolfinx::geometry::compute_collisions | ( | const BoundingBoxTree & | tree, |
const Eigen::Vector3d & | p | ||
) |
Compute all collisions between bounding boxes and point.
[in] | tree | The bounding box tree |
[in] | p | The point |
std::vector< std::array< int, 2 > > dolfinx::geometry::compute_collisions | ( | const BoundingBoxTree & | tree0, |
const BoundingBoxTree & | tree1 | ||
) |
Compute all collisions between two BoundingBoxTrees.
[in] | tree0 | First BoundingBoxTree |
[in] | tree1 | Second BoundingBoxTree |
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.
[in] | p | Body 1 list of points |
[in] | q | Body 2 list of points |
geometry::BoundingBoxTree dolfinx::geometry::create_midpoint_tree | ( | const mesh::Mesh & | mesh | ) |
Create a boundary box tree for cell midpoints.
[in] | mesh | The mesh build tree of cell midpoints from |
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.
[in] | mesh | Mesh |
[in] | candidate_cells | List of cell indices to test |
[in] | point | Point to check for collision |
[in] | n | Maximum number of positive results to return |
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.
[in] | mesh | Mesh containing the mesh entity |
[in] | dim | The topological dimension of the mesh entity |
[in] | index | The index of the mesh entity |
[in] | p | The point from which to compouted the shortest distance to the mesh to compute the Point |