dolfinx.geometry
Methods for geometric searches and operations
Functions
|
From a mesh, find which cells collide with a set of points. |
|
Compute the squared distance between a point and a mesh entity. |
Classes
|
A class for representing bounding box trees used in collision detection. |
- dolfinx.geometry.compute_closest_entity(tree: dolfinx::geometry::BoundingBoxTree, midpoint_tree: dolfinx::geometry::BoundingBoxTree, mesh: dolfinx.cpp.mesh.Mesh, points: numpy.ndarray[numpy.float64]) numpy.ndarray[numpy.int32]
- dolfinx.geometry.compute_colliding_cells(mesh: Mesh, candidates: AdjacencyList_int32, x: numpy.ndarray)[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(*args, **kwargs)
Overloaded function.
compute_collisions(tree: dolfinx::geometry::BoundingBoxTree, points: numpy.ndarray[numpy.float64]) -> dolfinx.cpp.graph.AdjacencyList_int32
compute_collisions(tree0: dolfinx::geometry::BoundingBoxTree, tree1: dolfinx::geometry::BoundingBoxTree) -> numpy.ndarray[numpy.int32]
- dolfinx.geometry.compute_distance_gjk(p: numpy.ndarray[numpy.float64], q: numpy.ndarray[numpy.float64]) numpy.ndarray[numpy.float64]
- dolfinx.geometry.create_midpoint_tree(mesh: dolfinx.cpp.mesh.Mesh, tdim: int, entities: numpy.ndarray[numpy.int32]) dolfinx::geometry::BoundingBoxTree
- dolfinx.geometry.squared_distance(mesh: Mesh, dim: int, entities: List[int], points: numpy.ndarray)[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 – The topological dimension of the mesh entities
entities – Indices of the mesh entities (local to process)
points – Set points from which to computed the shortest distance (
shape=(num_points, 3)
)
- Returns
Squared shortest distance from points[i] to entities[i]