DOLFINx
0.1.0
DOLFINx C++ interface
|
12 #include <xtl/xspan.hpp>
23 class BoundingBoxTree;
33 const xtl::span<const std::int32_t>& entity_indices);
40 std::vector<std::array<int, 2>>
48 const std::array<double, 3>& p);
58 std::pair<std::int32_t, double>
60 const std::array<double, 3>& p,
const mesh::Mesh& mesh,
67 const std::array<double, 3>& x);
83 const std::array<double, 3>& p);
93 std::vector<std::int32_t>
95 const xtl::span<const std::int32_t>& candidate_cells,
96 const std::array<double, 3>& p,
int n);
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 c...
Definition: utils.cpp:308
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 wit...
Definition: utils.cpp:364
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition: Mesh.h:55
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 ...
Definition: utils.cpp:271
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.
Definition: utils.cpp:209
std::vector< std::array< int, 2 > > compute_collisions(const BoundingBoxTree &tree0, const BoundingBoxTree &tree1)
Compute all collisions between two BoundingBoxTrees (local to process).
Definition: utils.cpp:230
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 in...
Definition: utils.cpp:254