dolfinx.geometry

Methods for geometric searches and operations.

Functions

bb_tree(mesh, dim[, entities, padding])

Create a bounding box tree for use in collision detection.

compute_colliding_cells(mesh, candidates, x)

From a mesh, find which cells collide with a set of points.

squared_distance(mesh, dim, entities, points)

Compute the squared distance between a point and a mesh entity.

compute_closest_entity(tree, midpoint_tree, ...)

Compute closest mesh entity to a point.

compute_collisions_trees(tree0, tree1)

Compute all collisions between two bounding box trees.

compute_collisions_points(tree, x)

Compute collisions between points and leaf bounding boxes.

compute_distance_gjk(p, q)

Compute the distance between two convex bodies p and q, each defined by a set of points.

create_midpoint_tree(mesh, dim, entities)

Create a bounding box tree for the midpoints of a subset of entities.

Classes

BoundingBoxTree(tree)

Bounding box trees used in collision detection.

PointOwnershipData(ownership_data)

Convenience class for storing data related to the ownership of points.

class dolfinx.geometry.BoundingBoxTree(tree)[source]

Bases: object

Bounding box trees used in collision detection.

Wrap a C++ BoundingBoxTree.

Note

This initializer should not be used in user code. Use

bb_tree.

create_global_tree(comm) BoundingBoxTree[source]
get_bbox(i) ndarray[Any, dtype[floating]][source]

Get lower and upper corners of the ith bounding box.

Parameters:

i – Index of the box.

Returns:

The ‘lower’ and ‘upper’ points of the bounding box. Shape is (2, 3),

property num_bboxes: int

Number of bounding boxes.

class dolfinx.geometry.PointOwnershipData(ownership_data)[source]

Bases: object

Convenience class for storing data related to the ownership of points.

Wrap a C++ PointOwnershipData.

dest_cells() ndarray[Any, dtype[int32]][source]

Cell indices (local to process) where each entry of dest_points is located.

dest_owner() ndarray[Any, dtype[int32]][source]

Ranks that sent dest_points to current process.

dest_points() ndarray[Any, dtype[floating]][source]

Points owned by current rank.

src_owner() ndarray[Any, dtype[int32]][source]

Ranks owning each point sent into ownership determination for current process.

dolfinx.geometry.bb_tree(mesh: Mesh, dim: int, entities: Optional[npt.NDArray[np.int32]] = None, padding: float = 0.0) BoundingBoxTree[source]

Create a bounding box tree for use in collision detection.

Parameters:
  • mesh – The mesh.

  • dim – Dimension of the mesh entities to build bounding box for.

  • entities – List of entity indices (local to process). If not supplied, all owned and ghosted entities are used.

  • padding – Padding for each bounding box.

Returns:

Bounding box tree.

dolfinx.geometry.compute_closest_entity(tree: BoundingBoxTree, midpoint_tree: BoundingBoxTree, mesh: Mesh, points: npt.NDArray[np.floating]) npt.NDArray[np.int32][source]

Compute closest mesh entity to a point.

Parameters:
  • tree – bounding box tree for the entities.

  • midpoint_tree – A bounding box tree with the midpoints of all the mesh entities. This is used to accelerate the search.

  • mesh – The mesh.

  • points – The points to check for collision, shape=(num_points,3).

Returns:

Mesh entity index for each point in points. Returns -1 for a point if the bounding box tree is empty.

dolfinx.geometry.compute_colliding_cells(mesh: Mesh, candidates: AdjacencyList_int32, x: npt.NDArray[np.floating])[source]

From a mesh, find which cells collide with a set of points.

Parameters:
  • mesh – The mesh.

  • candidate_cells – Adjacency list of candidate colliding cells for the ith point in x.

  • points – The points to check for collision shape=(num_points, 3),

Returns:

Adjacency list where the ith node is the list of entities that collide with the ith point.

dolfinx.geometry.compute_collisions_points(tree: BoundingBoxTree, x: ndarray[Any, dtype[floating]]) AdjacencyList_int32[source]

Compute collisions between points and leaf bounding boxes.

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

Parameters:
  • tree – Bounding box tree.

  • x – Points (shape=(num_points, 3)).

Returns:

For each point, the bounding box leaves that collide with the point.

dolfinx.geometry.compute_collisions_trees(tree0: BoundingBoxTree, tree1: BoundingBoxTree) ndarray[Any, dtype[int32]][source]

Compute all collisions between two bounding box trees.

Parameters:
  • tree0 – First bounding box tree.

  • tree1 – Second bounding box tree.

Returns:

List of pairs of intersecting box indices from each tree. Shape is (num_collisions, 2).

dolfinx.geometry.compute_distance_gjk(p: ndarray[Any, dtype[floating]], q: ndarray[Any, dtype[floating]]) ndarray[Any, dtype[floating]][source]

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:
  • p – Body 1 list of points (shape=(num_points, gdim)).

  • q – Body 2 list of points (shape=(num_points, gdim)).

Returns:

Shortest vector between the two bodies.

dolfinx.geometry.create_midpoint_tree(mesh: Mesh, dim: int, entities: npt.NDArray[np.int32]) BoundingBoxTree[source]

Create a bounding box tree for the midpoints of a subset of entities.

Parameters:
  • mesh – The mesh.

  • dim – Topological dimension of the entities.

  • entities – Indices of mesh entities to include.

Returns:

Bounding box tree for midpoints of cell entities.

dolfinx.geometry.squared_distance(mesh: Mesh, dim: int, entities: list[int], points: npt.NDArray[np.floating])[source]

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.

Parameters:
  • mesh – Mesh containing the entities.

  • dim – Topological dimension of the mesh entities.

  • entities – Indices of the mesh entities (local to process).

  • points – Points to compute the shortest distance from (shape=(num_points, 3)).

Returns:

Squared shortest distance from points[i] to entities[i].