1 #ifndef _COMPADRE_NEIGHBORLISTS_HPP_ 
    2 #define _COMPADRE_NEIGHBORLISTS_HPP_ 
    5 #include <Kokkos_Core.hpp> 
   10 template <
typename view_type>
 
   51     NeighborLists(view_type cr_neighbor_lists, view_type number_of_neighbors_list, 
 
   54                 "cr_neighbor_lists and number_neighbors_list and neighbor_lists_row_offsets must be a 1D Kokkos view.");
 
   78                 && 
"neighbor_lists is not large enough to store all neighbors.");
 
   86     NeighborLists(view_type cr_neighbor_lists, view_type number_of_neighbors_list) {
 
   88                 && 
"cr_neighbor_lists and number_neighbors_list must be a 1D Kokkos view.");
 
  110                 && 
"neighbor_lists is not large enough to store all neighbors.");
 
  121                 && 
"cr_neighbor_lists and number_neighbors_list must be a 1D Kokkos view.");
 
  152     KOKKOS_INLINE_FUNCTION
 
  163         Kokkos::parallel_reduce(
"max number of neighbors", 
 
  165                 KOKKOS_LAMBDA(
const int i, 
int& t_max_num_neighbors) {
 
  166             t_max_num_neighbors = (number_of_neighbors_list(i) > t_max_num_neighbors) ? number_of_neighbors_list(i) : t_max_num_neighbors;
 
  175         Kokkos::parallel_scan(
"number of neighbors offsets", 
 
  178             row_offsets(i) = lsum;
 
  179             lsum += number_of_neighbors_list(i);
 
  202     KOKKOS_INLINE_FUNCTION
 
  213     KOKKOS_INLINE_FUNCTION
 
  224     KOKKOS_INLINE_FUNCTION
 
  233                     && 
"Stale information in host_cr_neighbor_lists. Call CopyDeviceDataToHost() to refresh.");
 
  239     KOKKOS_INLINE_FUNCTION
 
  245     KOKKOS_INLINE_FUNCTION
 
  257     KOKKOS_INLINE_FUNCTION
 
  266 template <
typename view_type>
 
  272 template <
typename view_type>
 
  278 template <
typename view_type_2d, 
typename view_type_1d = Kokkos::View<
int*, 
typename view_type_2d::memory_space, 
typename view_type_2d::memory_traits> >
 
  283     int total_storage_size = 0;
 
  284     Kokkos::parallel_reduce(
"total number of neighbors over all lists", Kokkos::RangePolicy<typename view_type_2d::execution_space>(0, neighbor_lists.extent(0)), 
 
  285             KOKKOS_LAMBDA(
const int i, 
int& t_total_num_neighbors) {
 
  286         t_total_num_neighbors += neighbor_lists(i,0);
 
  287     }, Kokkos::Sum<int>(total_storage_size));
 
  291     view_type_1d new_cr_neighbor_lists(
"compressed row neighbor lists", total_storage_size);
 
  292     view_type_1d new_number_of_neighbors_list(
"number of neighbors list", neighbor_lists.extent(0));
 
  296     auto d_neighbor_lists = create_mirror_view(
typename view_type_1d::execution_space(), neighbor_lists);
 
  297     Kokkos::deep_copy(d_neighbor_lists, neighbor_lists);
 
  299     Kokkos::parallel_for(
"copy number of neighbors to compressed row", 
 
  300             Kokkos::RangePolicy<typename view_type_1d::execution_space>(0, neighbor_lists.extent(0)), 
 
  301             KOKKOS_LAMBDA(
const int i) {
 
  302         new_number_of_neighbors_list(i) = d_neighbor_lists(i,0);
 
  309     auto cr_data = nla.getNeighborLists();
 
  312     if (Kokkos::SpaceAccessibility<device_execution_space, typename view_type_1d::memory_space>::accessible==1) {
 
  313         Kokkos::parallel_for(
"copy neighbor lists to compressed row", Kokkos::RangePolicy<typename view_type_1d::execution_space>(0, neighbor_lists.extent(0)), 
 
  314                 KOKKOS_LAMBDA(
const int i) {
 
  315             for (
int j=0; j<d_neighbor_lists(i,0); ++j) {
 
  316                 cr_data(nla.getRowOffsetDevice(i)+j) = d_neighbor_lists(i,j+1);
 
  319         nla.copyDeviceDataToHost(); 
 
  325         Kokkos::parallel_for(
"copy neighbor lists to compressed row", Kokkos::RangePolicy<host_execution_space>(0, neighbor_lists.extent(0)), 
 
  326                 KOKKOS_LAMBDA(
const int i) {
 
  327             for (
int j=0; j<neighbor_lists(i,0); ++j) {
 
  328                 cr_data(nla.getRowOffsetHost(i)+j) = d_neighbor_lists(i,j+1);
 
view_type::HostMirror _host_number_of_neighbors_list
 
std::size_t global_index_type
 
KOKKOS_INLINE_FUNCTION int getMaxNumNeighbors() const 
Get the maximum number of neighbors of all targets' neighborhoods (host/device) 
 
internal_row_offsets_view_type::HostMirror _host_row_offsets
 
int getNumberOfNeighborsHost(int target_index) const 
Get number of neighbors for a given target (host) 
 
view_type getNeighborLists()
Device view into neighbor lists data (use with caution) 
 
int getNeighborHost(int target_index, int neighbor_num) const 
Offers N(i,j) indexing where N(i,j) is the index of the jth neighbor of i (host) 
 
Kokkos::View< global_index_type *, typename view_type::array_layout, typename view_type::memory_space, typename view_type::memory_traits > internal_row_offsets_view_type
 
void computeMaxNumNeighbors()
Calculate the maximum number of neighbors of all targets' neighborhoods (host) 
 
KOKKOS_INLINE_FUNCTION int getNeighborDevice(int target_index, int neighbor_num) const 
Offers N(i,j) indexing where N(i,j) is the index of the jth neighbor of i (device) ...
 
NeighborLists< view_type > CreateNeighborLists(view_type neighbor_lists, view_type number_of_neighbors_list)
CreateNeighborLists allows for the construction of an object of type NeighborLists with template dedu...
 
#define TO_GLOBAL(variable)
 
KOKKOS_INLINE_FUNCTION int getNumberOfTargets() const 
Get number of total targets having neighborhoods (host/device). 
 
int _max_neighbor_list_row_storage_size
 
view_type internal_view_type
 
KOKKOS_INLINE_FUNCTION int getNumberOfNeighborsDevice(int target_index) const 
Get number of neighbors for a given target (device) 
 
NeighborLists()
Constructor for the purpose of classes who have NeighborLists as a member object. ...
 
KOKKOS_INLINE_FUNCTION void setNeighborDevice(int target_index, int neighbor_num, int new_value)
Setter function for N(i,j) indexing where N(i,j) is the index of the jth neighbor of i...
 
NeighborLists(view_type number_of_neighbors_list)
Constructor for when number_of_neighbors_list is already populated. Will allocate space for compresse...
 
NeighborLists assists in accessing entries of compressed row neighborhood lists. 
 
view_type _cr_neighbor_lists
 
KOKKOS_INLINE_FUNCTION global_index_type getTotalNeighborsOverAllListsDevice() const 
Get the sum of the number of neighbors of all targets' neighborhoods (device) 
 
global_index_type getRowOffsetHost(int target_index) const 
Get offset into compressed row neighbor lists (host) 
 
NeighborLists(view_type cr_neighbor_lists, view_type number_of_neighbors_list)
Constructor for when compressed row cr_neighbor_lists is preallocated/populated, and number_of_neighb...
 
KOKKOS_INLINE_FUNCTION global_index_type getRowOffsetDevice(int target_index) const 
Get offset into compressed row neighbor lists (device) 
 
view_type _number_of_neighbors_list
 
void computeRowOffsets()
Calculate the row offsets for each target's neighborhood (host) 
 
void copyDeviceDataToHost()
Sync the host from the device (copy device data to host) 
 
global_index_type getTotalNeighborsOverAllListsHost() const 
Get the sum of the number of neighbors of all targets' neighborhoods (host) 
 
#define compadre_assert_debug(condition)
compadre_assert_debug is used for assertions that are checked in loops, as these significantly impact...
 
view_type::HostMirror _host_cr_neighbor_lists
 
NeighborLists< view_type_1d > Convert2DToCompressedRowNeighborLists(view_type_2d neighbor_lists)
Converts 2D neighbor lists to compressed row neighbor lists. 
 
NeighborLists(view_type cr_neighbor_lists, view_type number_of_neighbors_list, internal_row_offsets_view_type neighbor_lists_row_offsets, bool compute_max=true)
Constructor for when compressed row cr_neighbor_lists is preallocated/populated, number_of_neighbors_...
 
#define compadre_assert_release(condition)
compadre_assert_release is used for assertions that should always be checked, but generally are not e...
 
internal_row_offsets_view_type _row_offsets
 
#define compadre_kernel_assert_debug(condition)