Zoltan2
Zoltan2_AlgHybridD2.hpp
Go to the documentation of this file.
1 #ifndef _ZOLTAN2_DISTANCE2_HPP_
2 #define _ZOLTAN2_DISTANCE2_HPP_
3 
4 #include <vector>
5 #include <unordered_map>
6 #include <iostream>
7 #include <queue>
8 #ifdef _WIN32
9 #include <time.h>
10 #else
11 #include <sys/time.h>
12 #endif
13 
14 #include "Zoltan2_Algorithm.hpp"
15 #include "Zoltan2_GraphModel.hpp"
17 #include "Zoltan2_Util.hpp"
18 #include "Zoltan2_TPLTraits.hpp"
19 #include "Zoltan2_AlltoAll.hpp"
20 
21 
22 #include "Tpetra_Core.hpp"
23 #include "Teuchos_RCP.hpp"
24 #include "Tpetra_Import.hpp"
25 #include "Tpetra_FEMultiVector.hpp"
26 
27 #include "Kokkos_Core.hpp"
28 #include "KokkosSparse_CrsMatrix.hpp"
29 #include "KokkosKernels_Handle.hpp"
30 #include "KokkosKernels_IOUtils.hpp"
31 #include "KokkosGraph_Distance2Color.hpp"
32 #include "KokkosGraph_Distance2ColorHandle.hpp"
33 
37 
38 
39 namespace Zoltan2{
40 
41 template <typename Adapter>
42 class AlgDistance2 : public AlgTwoGhostLayer<Adapter> {
43 
44  public:
45 
46  using lno_t = typename Adapter::lno_t;
47  using gno_t = typename Adapter::gno_t;
48  using offset_t = typename Adapter::offset_t;
49  using scalar_t = typename Adapter::scalar_t;
51  using map_t = Tpetra::Map<lno_t,gno_t>;
52  using femv_scalar_t = int;
53  using femv_t = Tpetra::FEMultiVector<femv_scalar_t, lno_t, gno_t>;
54  using device_type = Tpetra::Map<>::device_type;
55  using execution_space = Tpetra::Map<>::execution_space;
56  using memory_space = Tpetra::Map<>::memory_space;
57  using host_exec = typename Kokkos::View<device_type>::HostMirror::execution_space;
58  using host_mem = typename Kokkos::View<device_type>::HostMirror::memory_space;
59  private:
60 
61  //This is both the serial and parallel local coloring.
62  template <class ExecutionSpace, typename MemorySpace>
63  void localColoring(const size_t nVtx,
64  Kokkos::View<lno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> adjs_view,
65  Kokkos::View<offset_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> offset_view,
66  Teuchos::RCP<femv_t> femv,
67  Kokkos::View<lno_t*, Kokkos::Device<ExecutionSpace, MemorySpace> > vertex_list,
68  size_t vertex_list_size = 0,
69  bool use_vertex_based_coloring = false){
70  using KernelHandle = KokkosKernels::Experimental::KokkosKernelsHandle
71  <offset_t, lno_t, lno_t, ExecutionSpace, MemorySpace, MemorySpace>;
72  KernelHandle kh;
73 
74  //Instead of switching between vertex-based and net-based algorithms,
75  //we only use the net-based algorithm, as it is faster than its
76  //vertex-based counterpart.
77  kh.create_distance2_graph_coloring_handle(KokkosGraph::COLORING_D2_NB_BIT);
78 
79  //vertex_list_size indicates whether we have provided a list of vertices to recolor
80  //NB_BIT does not make use of this argument currently.
81  if(vertex_list_size != 0){
82  kh.get_distance2_graph_coloring_handle()->set_vertex_list(vertex_list, vertex_list_size);
83  }
84 
85  //the verbose argument should carry through the local coloring
86  kh.set_verbose(this->verbose);
87 
88  //set initial colors to be the colors from femv
89  auto femvColors = femv->template getLocalView<Kokkos::Device<ExecutionSpace,MemorySpace> >(Tpetra::Access::ReadWrite);
90  auto sv = subview(femvColors, Kokkos::ALL, 0);
91  kh.get_distance2_graph_coloring_handle()->set_vertex_colors(sv);
92 
93  //call coloring
94  KokkosGraph::Experimental::graph_color_distance2(&kh, nVtx, offset_view, adjs_view);
95 
96 
97  //output total time
98  if(this->verbose){
99  std::cout<<"\nKokkosKernels Coloring: "
100  <<kh.get_distance2_graph_coloring_handle()->get_overall_coloring_time()
101  <<"\n";
102  }
103  }
104 
105  //Entry point for device-based coloring
106  virtual void colorInterior(const size_t nVtx,
107  Kokkos::View<lno_t*,device_type> adjs_view,
108  Kokkos::View<offset_t*, device_type> offset_view,
109  Teuchos::RCP<femv_t> femv,
110  Kokkos::View<lno_t*, device_type> vertex_list,
111  size_t vertex_list_size = 0,
112  bool recolor=false){
113 
114  this->localColoring<execution_space, memory_space>(nVtx,
115  adjs_view,
116  offset_view,
117  femv,
118  vertex_list,
119  vertex_list_size,
120  recolor);
121  }
122 
123  //Entry point for serial coloring
124  virtual void colorInterior_serial(const size_t nVtx,
125  typename Kokkos::View<lno_t*, device_type>::HostMirror adjs_view,
126  typename Kokkos::View<offset_t*, device_type>::HostMirror offset_view,
127  Teuchos::RCP<femv_t> femv,
128  typename Kokkos::View<lno_t*, device_type>::HostMirror vertex_list,
129  size_t vertex_list_size = 0,
130  bool recolor=false){
131  this->localColoring<host_exec, host_mem>(nVtx,
132  adjs_view,
133  offset_view,
134  femv,
135  vertex_list,
136  vertex_list_size,
137  recolor);
138  }
139  public:
140  //this function must be public due to Cuda Lambda restrictions.
141  //It is both the serial and parallel conflict detection function.
142  template< class ExecutionSpace, typename MemorySpace>
143  void detectD2Conflicts(const size_t n_local,
144  Kokkos::View<offset_t*, Kokkos::Device<ExecutionSpace,MemorySpace>> dist_offsets,
145  Kokkos::View<lno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> dist_adjs,
146  Kokkos::View<int*, Kokkos::Device<ExecutionSpace, MemorySpace>> femv_colors,
147  Kokkos::View<lno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> boundary_verts_view,
148  Kokkos::View<lno_t*,
149  Kokkos::Device<ExecutionSpace, MemorySpace> > verts_to_recolor_view,
150  Kokkos::View<int*,
151  Kokkos::Device<ExecutionSpace, MemorySpace>,
152  Kokkos::MemoryTraits<Kokkos::Atomic> > verts_to_recolor_size_atomic,
153  Kokkos::View<lno_t*,
154  Kokkos::Device<ExecutionSpace, MemorySpace> > verts_to_send_view,
155  Kokkos::View<size_t*,
156  Kokkos::Device<ExecutionSpace, MemorySpace>,
157  Kokkos::MemoryTraits<Kokkos::Atomic> > verts_to_send_size_atomic,
158  Kokkos::View<size_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> recoloringSize,
159  Kokkos::View<int*, Kokkos::Device<ExecutionSpace, MemorySpace>> rand,
160  Kokkos::View<gno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> gid,
161  Kokkos::View<gno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> ghost_degrees,
162  bool recolor_degrees){
163  Kokkos::RangePolicy<ExecutionSpace> policy(0, boundary_verts_view.extent(0));
164  size_t local_recoloring_size;
165  Kokkos::parallel_reduce("D2 conflict detection",policy, KOKKOS_LAMBDA(const uint64_t& i,size_t& recoloring_size){
166  //we only detect conflicts for vertices in the boundary
167  const size_t curr_lid = boundary_verts_view(i);
168  const int curr_color = femv_colors(curr_lid);
169  const size_t vid_d1_adj_begin = dist_offsets(curr_lid);
170  const size_t vid_d1_adj_end = dist_offsets(curr_lid+1);
171  const size_t curr_degree = vid_d1_adj_end - vid_d1_adj_begin;
172  for(size_t vid_d1_adj = vid_d1_adj_begin; vid_d1_adj < vid_d1_adj_end; vid_d1_adj++){
173  //check all distance-1 neighbors for conflicts
174  size_t vid_d1 = dist_adjs(vid_d1_adj);
175  size_t vid_d1_degree = 0;
176  //calculate the degree for degree-base recoloring
177  if(vid_d1 < n_local){
178  vid_d1_degree = dist_offsets(vid_d1+1) - dist_offsets(vid_d1);
179  } else {
180  vid_d1_degree = ghost_degrees(vid_d1-n_local);
181  }
182  if( vid_d1 != curr_lid && femv_colors(vid_d1) == curr_color){
183  if(curr_degree < vid_d1_degree && recolor_degrees){
184  femv_colors(curr_lid) = 0;
185  recoloring_size++;
186  break;//----------------------------------------------------
187  } else if (vid_d1_degree < curr_degree && recolor_degrees){//|
188  femv_colors(vid_d1) = 0; //|
189  recoloring_size++; //|
190  } else if(rand(curr_lid) < rand(vid_d1)){ //|
191  femv_colors(curr_lid) = 0; //|
192  recoloring_size++; //|
193  break;//---------------------------------------------------|
194  } else if(rand(vid_d1) < rand(curr_lid)){ //|
195  femv_colors(vid_d1) = 0; //|
196  recoloring_size++; //|
197  } else{ //|
198  if(gid(curr_lid) >= gid(vid_d1)){ //|
199  femv_colors(curr_lid) = 0; //|
200  recoloring_size++; //|
201  break;//-------------------------------------------------|
202  } else { // v
203  femv_colors(vid_d1) = 0; //If we uncolor the vertex whose
204  recoloring_size++; //neighbors we're checking, each
205  } //subsquent conflict check will
206  } //not do anything productive.
207  }
208  size_t d2_adj_begin = 0;
209  size_t d2_adj_end = 0;
210  d2_adj_begin = dist_offsets(vid_d1);
211  d2_adj_end = dist_offsets(vid_d1+1);
212 
213  //If we find a conflict that uncolors curr_lid, then we can safely stop
214  //detecting further conflicts. Since this is a nested loop, we need to
215  //break twice, using the found boolean.
216  bool found = false;
217  for(size_t vid_d2_adj = d2_adj_begin; vid_d2_adj < d2_adj_end; vid_d2_adj++){
218  //check all distance-2 neighbors for conflicts
219  const size_t vid_d2 = dist_adjs(vid_d2_adj);
220  size_t vid_d2_degree = 0;
221  //calculate the degree for degree-based recoloring
222  if(vid_d2 < n_local){
223  vid_d2_degree = dist_offsets(vid_d2+1) - dist_offsets(vid_d2);
224  } else {
225  vid_d2_degree = ghost_degrees(vid_d2-n_local);
226  }
227  if(curr_lid != vid_d2 && femv_colors(vid_d2) == curr_color){
228  if(curr_degree < vid_d2_degree && recolor_degrees){
229  found = true;
230  femv_colors(curr_lid) = 0;
231  recoloring_size++;
232  break;//---------------------------------------------------
233  } else if(vid_d2_degree < curr_degree && recolor_degrees){//|
234  femv_colors(vid_d2) = 0; //|
235  recoloring_size++; //|
236  } else if(rand(curr_lid) < rand(vid_d2)){ //|
237  found = true; //|
238  femv_colors(curr_lid) = 0; //|
239  recoloring_size++; //|
240  break;//--------------------------------------------------|
241  } else if(rand(vid_d2) < rand(curr_lid)){ //|
242  femv_colors(vid_d2) = 0; //|
243  recoloring_size++; //|
244  } else { //|
245  if(gid(curr_lid) >= gid(vid_d2)){ //|
246  found = true; //|
247  femv_colors(curr_lid) = 0; //|
248  recoloring_size++; //|
249  break;//------------------------------------------------|
250  } else { //|
251  femv_colors(vid_d2) = 0; //|
252  recoloring_size++;// v
253  }// If we uncolor the vertex whose neighbors we're
254  } // checking, each subsequent conflict check will
255  } // not do anything productive. We need this------
256  } // to completely move on to the next vertex. |
257  if(found) break;//<--------------------------------------------------
258  }
259  },local_recoloring_size);
260  Kokkos::deep_copy(recoloringSize,local_recoloring_size);
261  Kokkos::fence();
262 
263  //update the verts_to_send and verts_to_recolor views.
264  Kokkos::parallel_for("rebuild verts_to_send and verts_to_recolor",
265  Kokkos::RangePolicy<ExecutionSpace>(0,femv_colors.size()),
266  KOKKOS_LAMBDA(const uint64_t& i){
267  if(femv_colors(i) == 0){
268  //we only send vertices owned by the current process
269  if(i < n_local){
270  verts_to_send_view(verts_to_send_size_atomic(0)++) = i;
271  }
272  //we need to recolor all vertices, for consistency.
273  verts_to_recolor_view(verts_to_recolor_size_atomic(0)++) = i;
274  }
275  });
276  Kokkos::fence();
277 
278  }
279 
280  //Entry point for parallel conflict detection
281  virtual void detectConflicts(const size_t n_local,
282  Kokkos::View<offset_t*, device_type > dist_offsets_dev,
283  Kokkos::View<lno_t*, device_type > dist_adjs_dev,
284  Kokkos::View<int*,device_type > femv_colors,
285  Kokkos::View<lno_t*, device_type > boundary_verts_view,
286  Kokkos::View<lno_t*,
287  device_type > verts_to_recolor_view,
288  Kokkos::View<int*,
289  device_type,
290  Kokkos::MemoryTraits<Kokkos::Atomic>> verts_to_recolor_size_atomic,
291  Kokkos::View<lno_t*,
292  device_type > verts_to_send_view,
293  Kokkos::View<size_t*,
294  device_type,
295  Kokkos::MemoryTraits<Kokkos::Atomic>> verts_to_send_size_atomic,
296  Kokkos::View<size_t*, device_type> recoloringSize,
297  Kokkos::View<int*,
298  device_type> rand,
299  Kokkos::View<gno_t*,
300  device_type> gid,
301  Kokkos::View<gno_t*,
302  device_type> ghost_degrees,
303  bool recolor_degrees){
304 
305  this->detectD2Conflicts<execution_space, memory_space>(n_local,
306  dist_offsets_dev,
307  dist_adjs_dev,
308  femv_colors,
309  boundary_verts_view,
310  verts_to_recolor_view,
311  verts_to_recolor_size_atomic,
312  verts_to_send_view,
313  verts_to_send_size_atomic,
314  recoloringSize,
315  rand,
316  gid,
317  ghost_degrees,
318  recolor_degrees);
319  }
320  //Entry point for serial conflict detection
321  virtual void detectConflicts_serial(const size_t n_local,
322  typename Kokkos::View<offset_t*, device_type >::HostMirror dist_offsets_host,
323  typename Kokkos::View<lno_t*, device_type >::HostMirror dist_adjs_host,
324  typename Kokkos::View<int*,device_type >::HostMirror femv_colors,
325  typename Kokkos::View<lno_t*, device_type >::HostMirror boundary_verts_view,
326  typename Kokkos::View<lno_t*,device_type>::HostMirror verts_to_recolor,
327  typename Kokkos::View<int*,device_type>::HostMirror verts_to_recolor_size,
328  typename Kokkos::View<lno_t*,device_type>::HostMirror verts_to_send,
329  typename Kokkos::View<size_t*,device_type>::HostMirror verts_to_send_size,
330  typename Kokkos::View<size_t*, device_type>::HostMirror recoloringSize,
331  typename Kokkos::View<int*, device_type>::HostMirror rand,
332  typename Kokkos::View<gno_t*,device_type>::HostMirror gid,
333  typename Kokkos::View<gno_t*,device_type>::HostMirror ghost_degrees,
334  bool recolor_degrees) {
335  this->detectD2Conflicts<host_exec, host_mem>(n_local,
336  dist_offsets_host,
337  dist_adjs_host,
338  femv_colors,
339  boundary_verts_view,
340  verts_to_recolor,
341  verts_to_recolor_size,
342  verts_to_send,
343  verts_to_send_size,
344  recoloringSize,
345  rand,
346  gid,
347  ghost_degrees,
348  recolor_degrees);
349 
350  }
351 
352  virtual void constructBoundary(const size_t n_local,
353  Kokkos::View<offset_t*, device_type> dist_offsets_dev,
354  Kokkos::View<lno_t*, device_type> dist_adjs_dev,
355  typename Kokkos::View<offset_t*, device_type>::HostMirror dist_offsets_host,
356  typename Kokkos::View<lno_t*, device_type>::HostMirror dist_adjs_host,
357  Kokkos::View<lno_t*, device_type>& boundary_verts,
358  Kokkos::View<lno_t*,
359  device_type > verts_to_send_view,
360  Kokkos::View<size_t*,
361  device_type,
362  Kokkos::MemoryTraits<Kokkos::Atomic>> verts_to_send_size_atomic){
363  //count the number of boundary vertices to correctly allocate
364  //the boundary vertex view on device.
365  gno_t boundary_size_temp = 0;
366  for(size_t i = 0; i < n_local; i++){
367  for(offset_t j = dist_offsets_host(i); j < dist_offsets_host(i+1); j++){
368  if((size_t)dist_adjs_host(j) >= n_local){
369  boundary_size_temp++;
370  break;
371  }
372  bool found = false;
373  for(offset_t k = dist_offsets_host(dist_adjs_host(j)); k < dist_offsets_host(dist_adjs_host(j)+1); k++){
374  if((size_t)dist_adjs_host(k) >= n_local){
375  boundary_size_temp++;
376  found = true;
377  break;
378  }
379  }
380  if(found) break;
381  }
382  }
383 
384  //create a host mirror to fill in the list of boundary vertices
385  boundary_verts = Kokkos::View<lno_t*, device_type>("boundary verts",boundary_size_temp);
386  typename Kokkos::View<lno_t*, device_type>::HostMirror boundary_verts_host = Kokkos::create_mirror_view(boundary_verts);
387 
388  //reset the boundary size count to use as an index to construct the view
389  boundary_size_temp = 0;
390 
391  //a boundary vertex is any vertex within two edges of a ghost vertex.
392  for(size_t i = 0; i < n_local; i++){
393  for(offset_t j = dist_offsets_host(i); j < dist_offsets_host(i+1); j++){
394  if((size_t)dist_adjs_host(j) >= n_local){
395  boundary_verts_host(boundary_size_temp++) = i;
396  break;
397  }
398  bool found = false;
399  for(offset_t k = dist_offsets_host(dist_adjs_host(j)); k < dist_offsets_host(dist_adjs_host(j)+1); k++){
400  if((size_t)dist_adjs_host(k) >= n_local){
401  boundary_verts_host(boundary_size_temp++) = i;
402  found = true;
403  break;
404  }
405  }
406  if(found) break;
407  }
408  }
409  //copy the boundary to the device view
410  Kokkos::deep_copy(boundary_verts, boundary_verts_host);
411 
412  //initialize the list of verts to send
413  Kokkos::parallel_for(n_local, KOKKOS_LAMBDA(const int& i){
414  for(offset_t j = dist_offsets_dev(i); j < dist_offsets_dev(i+1); j++){
415  if((size_t)dist_adjs_dev(j) >= n_local){
416  verts_to_send_view(verts_to_send_size_atomic(0)++) = i;
417  break;
418  }
419  bool found = false;
420  for(offset_t k = dist_offsets_dev(dist_adjs_dev(j)); k < dist_offsets_dev(dist_adjs_dev(j)+1); k++){
421  if((size_t)dist_adjs_dev(k) >= n_local){
422  verts_to_send_view(verts_to_send_size_atomic(0)++) = i;
423  found = true;
424  break;
425  }
426  }
427  if(found) break;
428  }
429  });
430  Kokkos::fence();
431  }
432 
433  public:
435  const RCP<const base_adapter_t> &adapter_,
436  const RCP<Teuchos::ParameterList> &pl_,
437  const RCP<Environment> &env_,
438  const RCP<const Teuchos::Comm<int> > &comm_)
439  : AlgTwoGhostLayer<Adapter>(adapter_,pl_,env_,comm_){}
440 
441 }; //end class
442 
443 
444 }//end namespace Zoltan2
445 
446 #endif
AlltoAll communication methods.
Defines the ColoringSolution class.
Defines the GraphModel interface.
Traits class to handle conversions between gno_t/lno_t and TPL data types (e.g., ParMETIS's idx_t,...
A gathering of useful namespace methods.
virtual void constructBoundary(const size_t n_local, Kokkos::View< offset_t *, device_type > dist_offsets_dev, Kokkos::View< lno_t *, device_type > dist_adjs_dev, typename Kokkos::View< offset_t *, device_type >::HostMirror dist_offsets_host, typename Kokkos::View< lno_t *, device_type >::HostMirror dist_adjs_host, Kokkos::View< lno_t *, device_type > &boundary_verts, Kokkos::View< lno_t *, device_type > verts_to_send_view, Kokkos::View< size_t *, device_type, Kokkos::MemoryTraits< Kokkos::Atomic >> verts_to_send_size_atomic)
typename Adapter::lno_t lno_t
typename Adapter::offset_t offset_t
virtual void detectConflicts(const size_t n_local, Kokkos::View< offset_t *, device_type > dist_offsets_dev, Kokkos::View< lno_t *, device_type > dist_adjs_dev, Kokkos::View< int *, device_type > femv_colors, Kokkos::View< lno_t *, device_type > boundary_verts_view, Kokkos::View< lno_t *, device_type > verts_to_recolor_view, Kokkos::View< int *, device_type, Kokkos::MemoryTraits< Kokkos::Atomic >> verts_to_recolor_size_atomic, Kokkos::View< lno_t *, device_type > verts_to_send_view, Kokkos::View< size_t *, device_type, Kokkos::MemoryTraits< Kokkos::Atomic >> verts_to_send_size_atomic, Kokkos::View< size_t *, device_type > recoloringSize, Kokkos::View< int *, device_type > rand, Kokkos::View< gno_t *, device_type > gid, Kokkos::View< gno_t *, device_type > ghost_degrees, bool recolor_degrees)
AlgDistance2(const RCP< const base_adapter_t > &adapter_, const RCP< Teuchos::ParameterList > &pl_, const RCP< Environment > &env_, const RCP< const Teuchos::Comm< int > > &comm_)
virtual void detectConflicts_serial(const size_t n_local, typename Kokkos::View< offset_t *, device_type >::HostMirror dist_offsets_host, typename Kokkos::View< lno_t *, device_type >::HostMirror dist_adjs_host, typename Kokkos::View< int *, device_type >::HostMirror femv_colors, typename Kokkos::View< lno_t *, device_type >::HostMirror boundary_verts_view, typename Kokkos::View< lno_t *, device_type >::HostMirror verts_to_recolor, typename Kokkos::View< int *, device_type >::HostMirror verts_to_recolor_size, typename Kokkos::View< lno_t *, device_type >::HostMirror verts_to_send, typename Kokkos::View< size_t *, device_type >::HostMirror verts_to_send_size, typename Kokkos::View< size_t *, device_type >::HostMirror recoloringSize, typename Kokkos::View< int *, device_type >::HostMirror rand, typename Kokkos::View< gno_t *, device_type >::HostMirror gid, typename Kokkos::View< gno_t *, device_type >::HostMirror ghost_degrees, bool recolor_degrees)
void detectD2Conflicts(const size_t n_local, Kokkos::View< offset_t *, Kokkos::Device< ExecutionSpace, MemorySpace >> dist_offsets, Kokkos::View< lno_t *, Kokkos::Device< ExecutionSpace, MemorySpace >> dist_adjs, Kokkos::View< int *, Kokkos::Device< ExecutionSpace, MemorySpace >> femv_colors, Kokkos::View< lno_t *, Kokkos::Device< ExecutionSpace, MemorySpace >> boundary_verts_view, Kokkos::View< lno_t *, Kokkos::Device< ExecutionSpace, MemorySpace > > verts_to_recolor_view, Kokkos::View< int *, Kokkos::Device< ExecutionSpace, MemorySpace >, Kokkos::MemoryTraits< Kokkos::Atomic > > verts_to_recolor_size_atomic, Kokkos::View< lno_t *, Kokkos::Device< ExecutionSpace, MemorySpace > > verts_to_send_view, Kokkos::View< size_t *, Kokkos::Device< ExecutionSpace, MemorySpace >, Kokkos::MemoryTraits< Kokkos::Atomic > > verts_to_send_size_atomic, Kokkos::View< size_t *, Kokkos::Device< ExecutionSpace, MemorySpace >> recoloringSize, Kokkos::View< int *, Kokkos::Device< ExecutionSpace, MemorySpace >> rand, Kokkos::View< gno_t *, Kokkos::Device< ExecutionSpace, MemorySpace >> gid, Kokkos::View< gno_t *, Kokkos::Device< ExecutionSpace, MemorySpace >> ghost_degrees, bool recolor_degrees)
Tpetra::Map<>::memory_space memory_space
typename Kokkos::View< device_type >::HostMirror::execution_space host_exec
Tpetra::Map<>::device_type device_type
Tpetra::Map<>::execution_space execution_space
typename Adapter::base_adapter_t base_adapter_t
typename Adapter::gno_t gno_t
typename Adapter::scalar_t scalar_t
typename Adapter::offset_t offset_t
Tpetra::FEMultiVector< femv_scalar_t, lno_t, gno_t > femv_t
typename Kokkos::View< device_type >::HostMirror::memory_space host_mem
typename Adapter::lno_t lno_t
Tpetra::Map< lno_t, gno_t > map_t
map_t::local_ordinal_type lno_t
Definition: mapRemotes.cpp:17
map_t::global_ordinal_type gno_t
Definition: mapRemotes.cpp:18
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
Created by mbenlioglu on Aug 31, 2020.