Compadre  1.5.7
Compadre_NeighborLists.hpp
Go to the documentation of this file.
1 #ifndef _COMPADRE_NEIGHBORLISTS_HPP_
2 #define _COMPADRE_NEIGHBORLISTS_HPP_
3 
4 #include "Compadre_Typedefs.hpp"
5 #include <Kokkos_Core.hpp>
6 
7 namespace Compadre {
8 
9 //! NeighborLists assists in accessing entries of compressed row neighborhood lists
10 template <typename view_type>
11 struct NeighborLists {
12 
13  typedef view_type internal_view_type;
14  typedef Kokkos::View<global_index_type*, typename view_type::array_layout,
15  typename view_type::memory_space, typename view_type::memory_traits> internal_row_offsets_view_type;
16 
21 
23  view_type _cr_neighbor_lists;
25 
26  typename internal_row_offsets_view_type::HostMirror _host_row_offsets;
27  typename view_type::HostMirror _host_cr_neighbor_lists;
28  typename view_type::HostMirror _host_number_of_neighbors_list;
29 
30 /** @name Constructors
31  * Ways to initialize a NeighborLists object
32  */
33 ///@{
34 
35  //! \brief Constructor for the purpose of classes who have NeighborLists as a member object
39  _needs_sync_to_host = true;
41  }
42 
43  /*! \brief Constructor for when compressed row `cr_neighbor_lists` is preallocated/populated,
44  * `number_of_neighbors_list` and `neighbor_lists_row_offsets` have already been populated.
45  */
46  NeighborLists(view_type cr_neighbor_lists, view_type number_of_neighbors_list,
47  internal_row_offsets_view_type neighbor_lists_row_offsets, bool compute_max = true) {
48  compadre_assert_release((view_type::rank==1) &&
49  "cr_neighbor_lists and number_neighbors_list and neighbor_lists_row_offsets must be a 1D Kokkos view.");
50 
51  _number_of_targets = number_of_neighbors_list.extent(0);
52  _number_of_neighbors_list = number_of_neighbors_list;
53  _cr_neighbor_lists = cr_neighbor_lists;
54  _row_offsets = neighbor_lists_row_offsets;
55 
56  _host_cr_neighbor_lists = Kokkos::create_mirror_view(_cr_neighbor_lists);
57  _host_number_of_neighbors_list = Kokkos::create_mirror_view(_number_of_neighbors_list);
58  _host_row_offsets = Kokkos::create_mirror_view(_row_offsets);
59 
60  Kokkos::deep_copy(_host_cr_neighbor_lists, _cr_neighbor_lists);
62  Kokkos::deep_copy(_host_row_offsets, _row_offsets);
63  Kokkos::fence();
64 
65  if (compute_max) {
67  } else {
69  }
71 
72  //check neighbor_lists is large enough
73  compadre_assert_release(((size_t)(this->getTotalNeighborsOverAllListsHost())<=cr_neighbor_lists.extent(0))
74  && "neighbor_lists is not large enough to store all neighbors.");
75 
76  _needs_sync_to_host = false;
77  }
78 
79  /*! \brief Constructor for when compressed row `cr_neighbor_lists` is preallocated/populated,
80  * and `number_of_neighbors_list` is already populated, and row offsets still need to be computed.
81  */
82  NeighborLists(view_type cr_neighbor_lists, view_type number_of_neighbors_list) {
83  compadre_assert_release((view_type::rank==1)
84  && "cr_neighbor_lists and number_neighbors_list must be a 1D Kokkos view.");
85 
86  _number_of_targets = number_of_neighbors_list.extent(0);
87 
88  _row_offsets = internal_row_offsets_view_type("row offsets", number_of_neighbors_list.extent(0));
89  _number_of_neighbors_list = number_of_neighbors_list;
90  _cr_neighbor_lists = cr_neighbor_lists;
91 
92  _host_cr_neighbor_lists = Kokkos::create_mirror_view(_cr_neighbor_lists);
93  _host_number_of_neighbors_list = Kokkos::create_mirror_view(_number_of_neighbors_list);
94  _host_row_offsets = Kokkos::create_mirror_view(_row_offsets);
95 
96  Kokkos::deep_copy(_host_cr_neighbor_lists, _cr_neighbor_lists);
98  Kokkos::deep_copy(_host_row_offsets, _row_offsets);
99  Kokkos::fence();
100 
103 
104  //check neighbor_lists is large enough
105  compadre_assert_release(((size_t)(this->getTotalNeighborsOverAllListsHost())<=cr_neighbor_lists.extent(0))
106  && "neighbor_lists is not large enough to store all neighbors.");
107 
108  _needs_sync_to_host = false;
109  }
110 
111  /*! \brief Constructor for when `number_of_neighbors_list` is already populated.
112  * Will allocate space for compressed row neighbor lists data, and will allocate then
113  * populate information for row offsets.
114  */
115  NeighborLists(view_type number_of_neighbors_list) {
116  compadre_assert_release((view_type::rank==1)
117  && "cr_neighbor_lists and number_neighbors_list must be a 1D Kokkos view.");
118 
119  _number_of_targets = number_of_neighbors_list.extent(0);
120 
121  _row_offsets = internal_row_offsets_view_type("row offsets", number_of_neighbors_list.extent(0));
122  _number_of_neighbors_list = number_of_neighbors_list;
123 
124  _host_number_of_neighbors_list = Kokkos::create_mirror_view(_number_of_neighbors_list);
125  _host_row_offsets = Kokkos::create_mirror_view(_row_offsets);
126 
128  Kokkos::deep_copy(_host_row_offsets, _row_offsets);
129  Kokkos::fence();
130 
133 
134  _cr_neighbor_lists = view_type("compressed row neighbor lists data", this->getTotalNeighborsOverAllListsHost());
135  _host_cr_neighbor_lists = Kokkos::create_mirror_view(_cr_neighbor_lists);
136  Kokkos::deep_copy(_host_cr_neighbor_lists, _cr_neighbor_lists);
137  Kokkos::fence();
138 
139  _needs_sync_to_host = false;
140  }
141 ///@}
142 
143 /** @name Public modifiers
144  */
145 ///@{
146 
147  //! Setter function for N(i,j) indexing where N(i,j) is the index of the jth neighbor of i
148  KOKKOS_INLINE_FUNCTION
149  void setNeighborDevice(int target_index, int neighbor_num, int new_value) {
150  _cr_neighbor_lists(_row_offsets(target_index)+neighbor_num) = new_value;
151  // indicate that host view is now out of sync with device
152  // but only in debug mode (notice the next line is both setting the variable and checking it was set)
154  }
155 
156  //! Calculate the maximum number of neighbors of all targets' neighborhoods (host)
158  if (_number_of_neighbors_list.extent(0)==0) {
160  } else {
161  auto number_of_neighbors_list = _number_of_neighbors_list;
162  Kokkos::parallel_reduce("max number of neighbors",
163  Kokkos::RangePolicy<typename view_type::execution_space>(0, _number_of_neighbors_list.extent(0)),
164  KOKKOS_LAMBDA(const int i, int& t_max_num_neighbors) {
165  t_max_num_neighbors = (number_of_neighbors_list(i) > t_max_num_neighbors) ? number_of_neighbors_list(i) : t_max_num_neighbors;
166  }, Kokkos::Max<int>(_max_neighbor_list_row_storage_size));
167  Kokkos::fence();
168  }
169  }
170 
171  //! Calculate the minimum number of neighbors of all targets' neighborhoods (host)
173  if (_number_of_neighbors_list.extent(0)==0) {
175  } else {
176  auto number_of_neighbors_list = _number_of_neighbors_list;
177  Kokkos::parallel_reduce("min number of neighbors",
178  Kokkos::RangePolicy<typename view_type::execution_space>(0, _number_of_neighbors_list.extent(0)),
179  KOKKOS_LAMBDA(const int i, int& t_min_num_neighbors) {
180  t_min_num_neighbors = (number_of_neighbors_list(i) < t_min_num_neighbors) ? number_of_neighbors_list(i) : t_min_num_neighbors;
181  }, Kokkos::Min<int>(_min_neighbor_list_row_storage_size));
182  Kokkos::fence();
183  }
184  }
185 
186  //! Calculate the row offsets for each target's neighborhood (host)
188  auto number_of_neighbors_list = _number_of_neighbors_list;
189  auto row_offsets = _row_offsets;
190  Kokkos::parallel_scan("number of neighbors offsets",
191  Kokkos::RangePolicy<typename view_type::execution_space>(0, _number_of_neighbors_list.extent(0)),
192  KOKKOS_LAMBDA(const int i, global_index_type& lsum, bool final) {
193  row_offsets(i) = lsum;
194  lsum += number_of_neighbors_list(i);
195  });
196  Kokkos::deep_copy(_host_row_offsets, _row_offsets);
197  Kokkos::fence();
198  }
199 
200  //! Sync the host from the device (copy device data to host)
202  Kokkos::deep_copy(_host_cr_neighbor_lists, _cr_neighbor_lists);
203  Kokkos::fence();
204  _needs_sync_to_host = false;
205  }
206 
207  //! Device view into neighbor lists data (use with caution)
208  view_type getNeighborLists() const {
209  return _cr_neighbor_lists;
210  }
211 
212  //! Device view into number of neighbors list (use with caution)
213  view_type getNumberOfNeighborsList() const {
215  }
216 
217 ///@}
218 /** @name Public accessors
219  */
220 ///@{
221  //! Get number of total targets having neighborhoods (host/device).
222  KOKKOS_INLINE_FUNCTION
223  int getNumberOfTargets() const {
224  return _number_of_targets;
225  }
226 
227  //! Get number of neighbors for a given target (host)
228  int getNumberOfNeighborsHost(int target_index) const {
229  compadre_assert_extreme_debug(target_index < this->getNumberOfTargets());
230  return _host_number_of_neighbors_list(target_index);
231  }
232 
233  //! Get number of neighbors for a given target (device)
234  KOKKOS_INLINE_FUNCTION
235  int getNumberOfNeighborsDevice(int target_index) const {
237  return _number_of_neighbors_list(target_index);
238  }
239 
240  //! Get offset into compressed row neighbor lists (host)
241  global_index_type getRowOffsetHost(int target_index) const {
242  compadre_assert_extreme_debug(target_index < this->getNumberOfTargets());
243  return _host_row_offsets(target_index);
244  }
245 
246  //! Get offset into compressed row neighbor lists (device)
247  KOKKOS_INLINE_FUNCTION
248  global_index_type getRowOffsetDevice(int target_index) const {
250  return _row_offsets(target_index);
251  }
252 
253  //! Offers N(i,j) indexing where N(i,j) is the index of the jth neighbor of i (host)
254  int getNeighborHost(int target_index, int neighbor_num) const {
255  compadre_assert_extreme_debug(target_index < this->getNumberOfTargets());
257  && "Stale information in host_cr_neighbor_lists. Call CopyDeviceDataToHost() to refresh.");
258  compadre_assert_debug((neighbor_num<_host_number_of_neighbors_list(target_index))
259  && "neighor_num exceeds number of neighbors for this target_index.");
260  return _host_cr_neighbor_lists(_host_row_offsets(target_index)+neighbor_num);
261  }
262 
263  //! Offers N(i,j) indexing where N(i,j) is the index of the jth neighbor of i (device)
264  KOKKOS_INLINE_FUNCTION
265  int getNeighborDevice(int target_index, int neighbor_num) const {
267  && "neighor_num exceeds number of neighbors for this target_index.");
269  return _cr_neighbor_lists(_row_offsets(target_index)+neighbor_num);
270  }
271 
272  //! Get the maximum number of neighbors of all targets' neighborhoods (host/device)
273  KOKKOS_INLINE_FUNCTION
274  int getMaxNumNeighbors() const {
275  compadre_kernel_assert_debug((_max_neighbor_list_row_storage_size > -1) && "getMaxNumNeighbors() called but maximum never calculated.");
277  }
278 
279  //! Get the minimum number of neighbors of all targets' neighborhoods (host/device)
280  KOKKOS_INLINE_FUNCTION
281  int getMinNumNeighbors() const {
282  compadre_kernel_assert_debug((_min_neighbor_list_row_storage_size > -1) && "getMinNumNeighbors() called but minimum never calculated.");
284  }
285 
286  //! Get the sum of the number of neighbors of all targets' neighborhoods (host)
288  if (this->getNumberOfTargets()==0) {
289  return 0;
290  } else {
292  }
293  }
294 
295  //! Get the sum of the number of neighbors of all targets' neighborhoods (device)
296  KOKKOS_INLINE_FUNCTION
299  }
300 ///@}
301 
302 }; // NeighborLists
303 
304 //! CreateNeighborLists allows for the construction of an object of type NeighborLists with template deduction
305 template <typename view_type>
306 NeighborLists<view_type> CreateNeighborLists(view_type number_of_neighbors_list) {
307  return NeighborLists<view_type>(number_of_neighbors_list);
308 }
309 
310 //! CreateNeighborLists allows for the construction of an object of type NeighborLists with template deduction
311 template <typename view_type>
312 NeighborLists<view_type> CreateNeighborLists(view_type neighbor_lists, view_type number_of_neighbors_list) {
313  return NeighborLists<view_type>(neighbor_lists, number_of_neighbors_list);
314 }
315 
316 //! CreateNeighborLists allows for the construction of an object of type NeighborLists with template deduction
317 template <typename view_type>
318 NeighborLists<view_type> CreateNeighborLists(view_type neighbor_lists, view_type number_of_neighbors_list, view_type neighbor_lists_row_offsets) {
319  return NeighborLists<view_type>(neighbor_lists, number_of_neighbors_list, neighbor_lists_row_offsets);
320 }
321 
322 //! Converts 2D neighbor lists to compressed row neighbor lists
323 template <typename view_type_2d, typename view_type_1d = Kokkos::View<int*, typename view_type_2d::memory_space, typename view_type_2d::memory_traits> >
325 
326  // gets total number of neighbors over all lists
327  // computes calculation where the data resides (device/host)
328  global_index_type total_storage_size = 0;
329  Kokkos::parallel_reduce("total number of neighbors over all lists", Kokkos::RangePolicy<typename view_type_2d::execution_space>(0, neighbor_lists.extent(0)),
330  KOKKOS_LAMBDA(const int i, global_index_type& t_total_num_neighbors) {
331  t_total_num_neighbors += neighbor_lists(i,0);
332  }, Kokkos::Sum<global_index_type>(total_storage_size));
333  Kokkos::fence();
334 
335  // view_type_1d may be on host or device, and view_type_2d may be either as well (could even be opposite)
336  view_type_1d new_cr_neighbor_lists("compressed row neighbor lists", total_storage_size);
337  view_type_1d new_number_of_neighbors_list("number of neighbors list", neighbor_lists.extent(0));
338 
339  // copy number of neighbors list over to view_type_1d
340  // d_neighbor_lists will be accessible from view_type_1d's execution space
341  auto d_neighbor_lists = create_mirror_view(typename view_type_1d::execution_space(), neighbor_lists);
342  Kokkos::deep_copy(d_neighbor_lists, neighbor_lists);
343  Kokkos::fence();
344  Kokkos::parallel_for("copy number of neighbors to compressed row",
345  Kokkos::RangePolicy<typename view_type_1d::execution_space>(0, neighbor_lists.extent(0)),
346  KOKKOS_LAMBDA(const int i) {
347  new_number_of_neighbors_list(i) = d_neighbor_lists(i,0);
348  });
349  Kokkos::fence();
350 
351 
352  // this will calculate row offsets
353  auto nla = CreateNeighborLists(new_cr_neighbor_lists, new_number_of_neighbors_list);
354  auto cr_data = nla.getNeighborLists();
355 
356  // if device_execution_space can access this view, then write directly into the view
357  if (Kokkos::SpaceAccessibility<device_execution_space, typename view_type_1d::memory_space>::accessible==1) {
358  Kokkos::parallel_for("copy neighbor lists to compressed row", Kokkos::RangePolicy<typename view_type_1d::execution_space>(0, neighbor_lists.extent(0)),
359  KOKKOS_LAMBDA(const int i) {
360  for (int j=0; j<d_neighbor_lists(i,0); ++j) {
361  cr_data(nla.getRowOffsetDevice(i)+j) = d_neighbor_lists(i,j+1);
362  }
363  });
364  Kokkos::fence();
365  nla.copyDeviceDataToHost(); // has a fence at the end
366  }
367  // otherwise we are writing to a view that can't be seen from device (must be host space),
368  // and d_neighbor_lists was already made to be a view_type that is accessible from view_type_1d's execution_space
369  // (which we know is host) so we can do a parallel_for over the host_execution_space
370  else {
371  Kokkos::parallel_for("copy neighbor lists to compressed row", Kokkos::RangePolicy<host_execution_space>(0, neighbor_lists.extent(0)),
372  KOKKOS_LAMBDA(const int i) {
373  for (int j=0; j<neighbor_lists(i,0); ++j) {
374  cr_data(nla.getRowOffsetHost(i)+j) = d_neighbor_lists(i,j+1);
375  }
376  });
377  Kokkos::fence();
378  }
379 
380  return nla;
381 }
382 
383 } // Compadre namespace
384 
385 #endif
386 
#define compadre_kernel_assert_debug(condition)
std::size_t global_index_type
#define compadre_assert_extreme_debug(condition)
compadre_kernel_assert_debug is similar to compadre_assert_debug, but is a call on the device,...
#define compadre_assert_debug(condition)
compadre_assert_debug is used for assertions that are checked in loops, as these significantly impact...
#define TO_GLOBAL(variable)
#define compadre_kernel_assert_extreme_debug(condition)
#define compadre_assert_release(condition)
compadre_assert_release is used for assertions that should always be checked, but generally are not e...
NeighborLists< view_type_1d > Convert2DToCompressedRowNeighborLists(view_type_2d neighbor_lists)
Converts 2D neighbor lists to compressed row neighbor lists.
NeighborLists< view_type > CreateNeighborLists(view_type number_of_neighbors_list)
CreateNeighborLists allows for the construction of an object of type NeighborLists with template dedu...
NeighborLists assists in accessing entries of compressed row neighborhood lists.
KOKKOS_INLINE_FUNCTION global_index_type getTotalNeighborsOverAllListsDevice() const
Get the sum of the number of neighbors of all targets' neighborhoods (device)
view_type::HostMirror _host_number_of_neighbors_list
void computeMinNumNeighbors()
Calculate the minimum number of neighbors of all targets' neighborhoods (host)
KOKKOS_INLINE_FUNCTION int getMinNumNeighbors() const
Get the minimum number of neighbors of all targets' neighborhoods (host/device)
view_type getNeighborLists() const
Device view into neighbor lists data (use with caution)
global_index_type getRowOffsetHost(int target_index) const
Get offset into compressed row neighbor lists (host)
KOKKOS_INLINE_FUNCTION global_index_type getRowOffsetDevice(int target_index) const
Get offset into compressed row neighbor lists (device)
NeighborLists(view_type number_of_neighbors_list)
Constructor for when number_of_neighbors_list is already populated. Will allocate space for compresse...
internal_row_offsets_view_type::HostMirror _host_row_offsets
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
view_type getNumberOfNeighborsList() const
Device view into number of neighbors list (use with caution)
KOKKOS_INLINE_FUNCTION int getNumberOfNeighborsDevice(int target_index) const
Get number of neighbors for a given target (device)
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)
view_type::HostMirror _host_cr_neighbor_lists
global_index_type getTotalNeighborsOverAllListsHost() const
Get the sum of the number of neighbors of all targets' neighborhoods (host)
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_...
void computeRowOffsets()
Calculate the row offsets for each target's neighborhood (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 int getMaxNumNeighbors() const
Get the maximum number of neighbors of all targets' neighborhoods (host/device)
internal_row_offsets_view_type _row_offsets
KOKKOS_INLINE_FUNCTION int getNumberOfTargets() const
Get number of total targets having neighborhoods (host/device).
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()
Constructor for the purpose of classes who have NeighborLists as a member object.
void copyDeviceDataToHost()
Sync the host from the device (copy device data to host)
int getNumberOfNeighborsHost(int target_index) const
Get number of neighbors for a given target (host)
void computeMaxNumNeighbors()
Calculate the maximum number of neighbors of all targets' neighborhoods (host)