DOLFIN-X
DOLFIN-X C++ interface
|
9 #include "BoundingBoxTree.h"
10 #include <Eigen/Dense>
31 std::pair<std::vector<int>, std::vector<int>>
35 std::pair<std::vector<int>, std::vector<int>>
37 const BoundingBoxTree& tree1,
const mesh::Mesh& mesh0,
38 const mesh::Mesh& mesh1);
45 const Eigen::Vector3d& p);
53 const Eigen::Vector3d& p,
54 const mesh::Mesh& mesh);
61 const Eigen::Vector3d& p);
69 const Eigen::Vector3d& p,
70 const mesh::Mesh& mesh);
75 const Eigen::Vector3d& p);
78 bool bbox_in_bbox(
const Eigen::Array<double, 2, 3, Eigen::RowMajor>& a,
79 const Eigen::Array<double, 2, 3, Eigen::RowMajor>& b,
84 std::pair<int, double>
86 const BoundingBoxTree& tree_midpoint,
87 const Eigen::Vector3d& p,
const mesh::Mesh& mesh);
95 const Eigen::Vector3d& p);
103 bool point_in_bbox(
const Eigen::Array<double, 2, 3, Eigen::RowMajor>& b,
104 const Eigen::Vector3d& x,
double rtol = 1e-14);
109 const Eigen::Array<double, 2, 3, Eigen::RowMajor>& b,
110 const Eigen::Vector3d& x);
115 const Eigen::Vector3d& p);
122 const Eigen::Vector3d& a,
123 const Eigen::Vector3d& b,
124 const Eigen::Vector3d& c);
130 const Eigen::Vector3d& a,
131 const Eigen::Vector3d& b);
BoundingBoxTree create_midpoint_tree(const mesh::Mesh &mesh)
Create a boundary box tree for cell midpoints.
Definition: utils.cpp:342
double squared_distance_interval(const Eigen::Vector3d &point, const Eigen::Vector3d &a, const Eigen::Vector3d &b)
Compute squared distance to given point. This version takes the two vertex coordinates as 3D points....
Definition: utils.cpp:749
double squared_distance(const mesh::MeshEntity &entity, const Eigen::Vector3d &p)
Compute squared distance from a given point to the nearest point on a cell (only simplex cells are su...
Definition: utils.cpp:545
std::pair< int, double > compute_closest_point(const BoundingBoxTree &tree, const Eigen::Vector3d &p)
Compute closest point and distance to a given point.
Definition: utils.cpp:510
std::pair< std::vector< int >, std::vector< int > > compute_entity_collisions(const BoundingBoxTree &tree0, const BoundingBoxTree &tree1, const mesh::Mesh &mesh0, const mesh::Mesh &mesh1)
Compute all collisions between entities and BoundingBoxTree.
Definition: utils.cpp:381
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 topol...
Definition: utils.cpp:479
int compute_first_collision(const BoundingBoxTree &tree, const Eigen::Vector3d &p)
Compute first collision between bounding boxes and point.
Definition: utils.cpp:423
bool point_in_bbox(const Eigen::Array< double, 2, 3, Eigen::RowMajor > &b, const Eigen::Vector3d &x, double rtol=1e-14)
Check whether point (x) is in bounding box.
Definition: utils.cpp:536
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 in...
Definition: utils.cpp:469
bool bbox_in_bbox(const Eigen::Array< double, 2, 3, Eigen::RowMajor > &a, const Eigen::Array< double, 2, 3, Eigen::RowMajor > &b, double rtol=1e-14)
Check whether bounding box a collides with bounding box (b)
Definition: utils.cpp:460
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.
Definition: utils.cpp:446
std::pair< std::vector< int >, std::vector< int > > compute_collisions(const BoundingBoxTree &tree0, const BoundingBoxTree &tree1)
Compute all collisions between bounding boxes and BoundingBoxTree.
Definition: utils.cpp:365
int compute_first_entity_collision(const BoundingBoxTree &tree, const Eigen::Vector3d &p, const mesh::Mesh &mesh)
Compute first collision between entities and point.
Definition: utils.cpp:429
double squared_distance_triangle(const Eigen::Vector3d &point, const Eigen::Vector3d &a, const Eigen::Vector3d &b, const Eigen::Vector3d &c)
Compute squared distance to given point. This version takes the three vertex coordinates as 3D points...
Definition: utils.cpp:677