Compadre
1.5.9
|
PointCloudSearch generates neighbor lists and window sizes for each target site. More...
PointCloudSearch generates neighbor lists and window sizes for each target site.
Search methods can be run in dry-run mode, or not.
neighbors_list
will be populated with number of neighbors found for each target site.
This allows a user to know memory allocation needed before storage of neighbor indices.
neighbors_list_offsets
will be populated with offsets for values (dependent on method) determined by neighbor_list. If a 2D view for neighbors_list
is used, then will store the neighbor of , and will store the number of neighbors for target .
Definition at line 151 of file Compadre_PointCloudSearch.hpp.
#include <Compadre_PointCloudSearch.hpp>
Public Types | |
typedef nanoflann::KDTreeSingleIndexAdaptor < nanoflann::L2_Simple_Adaptor < double, PointCloudSearch < view_type > >, PointCloudSearch < view_type >, 1 > | tree_type_1d |
typedef nanoflann::KDTreeSingleIndexAdaptor < nanoflann::L2_Simple_Adaptor < double, PointCloudSearch < view_type > >, PointCloudSearch < view_type >, 2 > | tree_type_2d |
typedef nanoflann::KDTreeSingleIndexAdaptor < nanoflann::L2_Simple_Adaptor < double, PointCloudSearch < view_type > >, PointCloudSearch < view_type >, 3 > | tree_type_3d |
Public Member Functions | |
PointCloudSearch (view_type src_pts_view, const local_index_type dimension=-1, const local_index_type max_leaf=-1) | |
~PointCloudSearch () | |
template<class BBOX > | |
bool | kdtree_get_bbox (BBOX &bb) const |
Bounding box query method required by Nanoflann. More... | |
int | kdtree_get_point_count () const |
Returns the number of source sites. More... | |
double | kdtree_get_pt (const int idx, int dim) const |
Returns the coordinate value of a point. More... | |
double | kdtree_distance (const double *queryPt, const int idx, long long sz) const |
Returns the distance between a point and a source site, given its index. More... | |
void | generateKDTree () |
template<typename trg_view_type , typename neighbor_lists_view_type , typename epsilons_view_type > | |
size_t | generate2DNeighborListsFromRadiusSearch (bool is_dry_run, trg_view_type trg_pts_view, neighbor_lists_view_type neighbor_lists, epsilons_view_type epsilons, const double uniform_radius=0.0, double max_search_radius=0.0) |
Generates neighbor lists of 2D view by performing a radius search where the radius to be searched is in the epsilons view. If uniform_radius is given, then this overrides the epsilons view radii sizes. Accepts 2D neighbor_lists without number_of_neighbors_list. More... | |
template<typename trg_view_type , typename neighbor_lists_view_type , typename epsilons_view_type > | |
size_t | generateCRNeighborListsFromRadiusSearch (bool is_dry_run, trg_view_type trg_pts_view, neighbor_lists_view_type neighbor_lists, neighbor_lists_view_type number_of_neighbors_list, epsilons_view_type epsilons, const double uniform_radius=0.0, double max_search_radius=0.0) |
Generates compressed row neighbor lists by performing a radius search where the radius to be searched is in the epsilons view. If uniform_radius is given, then this overrides the epsilons view radii sizes. Accepts 1D neighbor_lists with 1D number_of_neighbors_list. More... | |
template<typename trg_view_type , typename neighbor_lists_view_type , typename epsilons_view_type > | |
size_t | generate2DNeighborListsFromKNNSearch (bool is_dry_run, trg_view_type trg_pts_view, neighbor_lists_view_type neighbor_lists, epsilons_view_type epsilons, const int neighbors_needed, const double epsilon_multiplier=1.6, double max_search_radius=0.0) |
Generates neighbor lists as 2D view by performing a k-nearest neighbor search Only accepts 2D neighbor_lists without number_of_neighbors_list. More... | |
template<typename trg_view_type , typename neighbor_lists_view_type , typename epsilons_view_type > | |
size_t | generateCRNeighborListsFromKNNSearch (bool is_dry_run, trg_view_type trg_pts_view, neighbor_lists_view_type neighbor_lists, neighbor_lists_view_type number_of_neighbors_list, epsilons_view_type epsilons, const int neighbors_needed, const double epsilon_multiplier=1.6, double max_search_radius=0.0) |
Generates compressed row neighbor lists by performing a k-nearest neighbor search Only accepts 1D neighbor_lists with 1D number_of_neighbors_list. More... | |
Static Public Member Functions | |
static int | getEstimatedNumberNeighborsUpperBound (int unisolvency_size, const int dimension, const double epsilon_multiplier) |
Returns a liberal estimated upper bound on number of neighbors to be returned by a neighbor search for a given choice of dimension, basis size, and epsilon_multiplier. More... | |
Protected Attributes | |
view_type | _src_pts_view |
source site coordinates More... | |
local_index_type | _dim |
local_index_type | _max_leaf |
std::shared_ptr< tree_type_1d > | _tree_1d |
std::shared_ptr< tree_type_2d > | _tree_2d |
std::shared_ptr< tree_type_3d > | _tree_3d |
typedef nanoflann::KDTreeSingleIndexAdaptor<nanoflann::L2_Simple_Adaptor<double, PointCloudSearch<view_type> >, PointCloudSearch<view_type>, 1> Compadre::PointCloudSearch< view_type >::tree_type_1d |
Definition at line 156 of file Compadre_PointCloudSearch.hpp.
typedef nanoflann::KDTreeSingleIndexAdaptor<nanoflann::L2_Simple_Adaptor<double, PointCloudSearch<view_type> >, PointCloudSearch<view_type>, 2> Compadre::PointCloudSearch< view_type >::tree_type_2d |
Definition at line 158 of file Compadre_PointCloudSearch.hpp.
typedef nanoflann::KDTreeSingleIndexAdaptor<nanoflann::L2_Simple_Adaptor<double, PointCloudSearch<view_type> >, PointCloudSearch<view_type>, 3> Compadre::PointCloudSearch< view_type >::tree_type_3d |
Definition at line 160 of file Compadre_PointCloudSearch.hpp.
|
inline |
Definition at line 175 of file Compadre_PointCloudSearch.hpp.
|
inline |
Definition at line 184 of file Compadre_PointCloudSearch.hpp.
|
inline |
Generates neighbor lists as 2D view by performing a k-nearest neighbor search Only accepts 2D neighbor_lists without number_of_neighbors_list.
is_dry_run | [in] - whether to do a dry-run (find neighbors, but don't store) |
trg_pts_view | [in] - target coordinates from which to seek neighbors |
neighbor_lists | [out] - 2D view of neighbor lists to be populated from search |
epsilons | [in/out] - radius to search, overwritten if uniform_radius != 0 |
neighbors_needed | [in] - k neighbors needed as a minimum |
epsilon_multiplier | [in] - distance to kth neighbor multiplied by epsilon_multiplier for follow-on radius search |
max_search_radius | [in] - largest valid search (useful only for MPI jobs if halo size exists) |
Definition at line 518 of file Compadre_PointCloudSearch.hpp.
|
inline |
Generates neighbor lists of 2D view by performing a radius search where the radius to be searched is in the epsilons view. If uniform_radius is given, then this overrides the epsilons view radii sizes. Accepts 2D neighbor_lists without number_of_neighbors_list.
is_dry_run | [in] - whether to do a dry-run (find neighbors, but don't store) |
trg_pts_view | [in] - target coordinates from which to seek neighbors |
neighbor_lists | [out] - 2D view of neighbor lists to be populated from search |
epsilons | [in/out] - radius to search, overwritten if uniform_radius != 0 |
uniform_radius | [in] - double != 0 determines whether to overwrite all epsilons for uniform search |
max_search_radius | [in] - largest valid search (useful only for MPI jobs if halo size exists) |
Definition at line 240 of file Compadre_PointCloudSearch.hpp.
|
inline |
Generates compressed row neighbor lists by performing a k-nearest neighbor search Only accepts 1D neighbor_lists with 1D number_of_neighbors_list.
is_dry_run | [in] - whether to do a dry-run (find neighbors, but don't store) |
trg_pts_view | [in] - target coordinates from which to seek neighbors |
neighbor_lists | [out] - 1D view of neighbor lists to be populated from search |
number_of_neighbors_list | [in/out] - number of neighbors for each target site |
epsilons | [in/out] - radius to search, overwritten if uniform_radius != 0 |
neighbors_needed | [in] - k neighbors needed as a minimum |
epsilon_multiplier | [in] - distance to kth neighbor multiplied by epsilon_multiplier for follow-on radius search |
max_search_radius | [in] - largest valid search (useful only for MPI jobs if halo size exists) |
Definition at line 660 of file Compadre_PointCloudSearch.hpp.
|
inline |
Generates compressed row neighbor lists by performing a radius search where the radius to be searched is in the epsilons view. If uniform_radius is given, then this overrides the epsilons view radii sizes. Accepts 1D neighbor_lists with 1D number_of_neighbors_list.
is_dry_run | [in] - whether to do a dry-run (find neighbors, but don't store) |
trg_pts_view | [in] - target coordinates from which to seek neighbors |
neighbor_lists | [out] - 1D view of neighbor lists to be populated from search |
number_of_neighbors_list | [in/out] - number of neighbors for each target site |
epsilons | [in/out] - radius to search, overwritten if uniform_radius != 0 |
uniform_radius | [in] - double != 0 determines whether to overwrite all epsilons for uniform search |
max_search_radius | [in] - largest valid search (useful only for MPI jobs if halo size exists) |
Definition at line 374 of file Compadre_PointCloudSearch.hpp.
|
inline |
Definition at line 215 of file Compadre_PointCloudSearch.hpp.
|
inlinestatic |
Returns a liberal estimated upper bound on number of neighbors to be returned by a neighbor search for a given choice of dimension, basis size, and epsilon_multiplier.
Assumes quasiuniform distribution of points. This result can be used to size a preallocated neighbor_lists kokkos view.
Definition at line 189 of file Compadre_PointCloudSearch.hpp.
|
inline |
Returns the distance between a point and a source site, given its index.
Definition at line 205 of file Compadre_PointCloudSearch.hpp.
|
inline |
Bounding box query method required by Nanoflann.
Definition at line 196 of file Compadre_PointCloudSearch.hpp.
|
inline |
Returns the number of source sites.
Definition at line 199 of file Compadre_PointCloudSearch.hpp.
|
inline |
Returns the coordinate value of a point.
Definition at line 202 of file Compadre_PointCloudSearch.hpp.
|
protected |
Definition at line 166 of file Compadre_PointCloudSearch.hpp.
|
protected |
Definition at line 167 of file Compadre_PointCloudSearch.hpp.
|
protected |
source site coordinates
Definition at line 165 of file Compadre_PointCloudSearch.hpp.
|
protected |
Definition at line 169 of file Compadre_PointCloudSearch.hpp.
|
protected |
Definition at line 170 of file Compadre_PointCloudSearch.hpp.
|
protected |
Definition at line 171 of file Compadre_PointCloudSearch.hpp.