Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Zoltan2_AlgHybridPD2.hpp
Go to the documentation of this file.
1 #ifndef _ZOLTAN2_PDISTANCE2_HPP_
2 #define _ZOLTAN2_PDISTANCE2_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 AlgPartialDistance2 : 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 = typename femv_t::device_type;
55  using execution_space = typename device_type::execution_space;
56  using memory_space = typename device_type::memory_space;
57  using host_exec = typename femv_t::host_view_type::device_type::execution_space;
58  using host_mem = typename femv_t::host_view_type::device_type::memory_space;
59 
60  private:
61  //serial and parallel local partial distance-2 coloring function
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 use only 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 currently make use of this.
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.get_distance2_graph_coloring_handle()->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::bipartite_color_rows(&kh, nVtx, nVtx, offset_view, adjs_view, true);
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 parallel local 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  this->localColoring<execution_space, memory_space>(nVtx,
114  adjs_view,
115  offset_view,
116  femv,
117  vertex_list,
118  vertex_list_size,
119  recolor);
120  }
121  //Entry point for serial local coloring
122  virtual void colorInterior_serial(const size_t nVtx,
123  typename Kokkos::View<lno_t*, device_type >::HostMirror adjs_view,
124  typename Kokkos::View<offset_t*,device_type >::HostMirror offset_view,
125  Teuchos::RCP<femv_t> femv,
126  typename Kokkos::View<lno_t*, device_type>::HostMirror vertex_list,
127  size_t vertex_list_size = 0,
128  bool recolor=false) {
129  this->localColoring<host_exec, host_mem>(nVtx,
130  adjs_view,
131  offset_view,
132  femv,
133  vertex_list,
134  vertex_list_size,
135  recolor);
136 
137  }
138  public:
139  //serial and parallel partial distance-2 conflict detection function
140  template <class ExecutionSpace, typename MemorySpace>
141  void detectPD2Conflicts(const size_t n_local,
142  Kokkos::View<offset_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> dist_offsets,
143  Kokkos::View<lno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> dist_adjs,
144  Kokkos::View<int*, Kokkos::Device<ExecutionSpace, MemorySpace>> femv_colors,
145  Kokkos::View<lno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> boundary_verts_view,
146  Kokkos::View<lno_t*,
147  Kokkos::Device<ExecutionSpace, MemorySpace> > verts_to_recolor_view,
148  Kokkos::View<int*,
149  Kokkos::Device<ExecutionSpace, MemorySpace>,
150  Kokkos::MemoryTraits<Kokkos::Atomic> > verts_to_recolor_size_atomic,
151  Kokkos::View<lno_t*,
152  Kokkos::Device<ExecutionSpace, MemorySpace> > verts_to_send_view,
153  Kokkos::View<size_t*,
154  Kokkos::Device<ExecutionSpace, MemorySpace>,
155  Kokkos::MemoryTraits<Kokkos::Atomic> > verts_to_send_size_atomic,
156  Kokkos::View<size_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> recoloringSize,
157  Kokkos::View<int*, Kokkos::Device<ExecutionSpace, MemorySpace>> rand,
158  Kokkos::View<gno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> gid,
159  Kokkos::View<gno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> ghost_degrees,
160  bool recolor_degrees){
161 
162  Kokkos::RangePolicy<ExecutionSpace> policy(0,boundary_verts_view.extent(0));
163  size_t local_recoloring_size;
164  Kokkos::parallel_reduce("PD2 conflict detection",policy, KOKKOS_LAMBDA(const uint64_t& i,size_t& recoloring_size){
165  //we only detect conflicts for vertices in the boundary
166  const size_t curr_lid = boundary_verts_view(i);
167  const int curr_color = femv_colors(curr_lid);
168  const size_t vid_d1_adj_begin = dist_offsets(curr_lid);
169  const size_t vid_d1_adj_end = dist_offsets(curr_lid+1);
170  const size_t curr_degree = vid_d1_adj_end - vid_d1_adj_begin;
171  for(size_t vid_d1_adj = vid_d1_adj_begin; vid_d1_adj < vid_d1_adj_end; vid_d1_adj++){
172  //we skip direct distance-1 conflicts, and only resolve distance-2 conflicts.
173  size_t vid_d1 = dist_adjs(vid_d1_adj);
174  size_t d2_adj_begin = dist_offsets(vid_d1);
175  size_t d2_adj_end = dist_offsets(vid_d1+1);
176 
177  //If we find a conflict that uncolors curr_lid, we can safely stop detecting
178  //further conflicts starting at curr_lid. Because this is a nested loop,
179  //we'll use the found boolean to break twice and move to the next vertex.
180  bool found = false;
181  for(size_t vid_d2_adj = d2_adj_begin; vid_d2_adj < d2_adj_end; vid_d2_adj++){
182  const size_t vid_d2 = dist_adjs(vid_d2_adj);
183  size_t vid_d2_degree = 0;
184 
185  //calculate the degree for degree-based recoloring
186  if(vid_d2 < n_local){
187  vid_d2_degree = dist_offsets(vid_d2+1) - dist_offsets(vid_d2);
188  } else {
189  vid_d2_degree = ghost_degrees(vid_d2-n_local);
190  }
191  //resolve conflict
192  if(curr_lid != vid_d2 && femv_colors(vid_d2) == curr_color){
193  if(curr_degree < vid_d2_degree && recolor_degrees){
194  found = true;
195  femv_colors(curr_lid) = 0;
196  recoloring_size++;
197  break;//---------------------------------------------------
198  } else if(vid_d2_degree < curr_degree && recolor_degrees){//|
199  femv_colors(vid_d2) = 0; //|
200  recoloring_size++; //|
201  } else if(rand(curr_lid) < rand(vid_d2)){ //|
202  found = true; //|
203  femv_colors(curr_lid) = 0; //|
204  recoloring_size++; //|
205  break;//--------------------------------------------------|
206  } else if(rand(vid_d2) < rand(curr_lid)){ //|
207  femv_colors(vid_d2) = 0; //|
208  recoloring_size++; //|
209  } else { //|
210  if(gid(curr_lid) >= gid(vid_d2)){ //|
211  found = true; //|
212  femv_colors(curr_lid) = 0; //|
213  recoloring_size++; //|
214  break;//------------------------------------------------|
215  } else { //|
216  femv_colors(vid_d2) = 0; //|
217  recoloring_size++;// v
218  }// If we uncolor the vertex whose neighbors we're
219  } // checking, each subsequent conflict check will
220  } // not do anything productive. We need this------
221  } // to completely move on to the next vertex. |
222  if(found) break;//<--------------------------------------------------
223  }
224  },local_recoloring_size);
225  Kokkos::deep_copy(recoloringSize, local_recoloring_size);
226  Kokkos::fence();
227  //update the verts_to_send and verts_to_recolor views
228  Kokkos::parallel_for("rebuild verts_to_send and verts_to_recolor",
229  Kokkos::RangePolicy<ExecutionSpace>(0,femv_colors.size()),
230  KOKKOS_LAMBDA(const uint64_t& i){
231  if(femv_colors(i) == 0){
232  if(i < n_local){
233  //we only send vertices owned by the current process
234  verts_to_send_view(verts_to_send_size_atomic(0)++) = i;
235  }
236  //we need to recolor all vertices for consistency.
237  verts_to_recolor_view(verts_to_recolor_size_atomic(0)++) = i;
238  }
239  });
240  Kokkos::fence();
241 
242  }
243  //Entry point for parallel conflict detection
244  virtual void detectConflicts(const size_t n_local,
245  Kokkos::View<offset_t*, device_type > dist_offsets_dev,
246  Kokkos::View<lno_t*, device_type > dist_adjs_dev,
247  Kokkos::View<int*,device_type > femv_colors,
248  Kokkos::View<lno_t*, device_type > boundary_verts_view,
249  Kokkos::View<lno_t*,
250  device_type > verts_to_recolor_view,
251  Kokkos::View<int*,
252  device_type,
253  Kokkos::MemoryTraits<Kokkos::Atomic>> verts_to_recolor_size_atomic,
254  Kokkos::View<lno_t*,
255  device_type > verts_to_send_view,
256  Kokkos::View<size_t*,
257  device_type,
258  Kokkos::MemoryTraits<Kokkos::Atomic>> verts_to_send_size_atomic,
259  Kokkos::View<size_t*, device_type> recoloringSize,
260  Kokkos::View<int*,
261  device_type> rand,
262  Kokkos::View<gno_t*,
263  device_type> gid,
264  Kokkos::View<gno_t*,
265  device_type> ghost_degrees,
266  bool recolor_degrees){
267 
268  this->detectPD2Conflicts<execution_space, memory_space>(n_local,
269  dist_offsets_dev,
270  dist_adjs_dev,
271  femv_colors,
272  boundary_verts_view,
273  verts_to_recolor_view,
274  verts_to_recolor_size_atomic,
275  verts_to_send_view,
276  verts_to_send_size_atomic,
277  recoloringSize,
278  rand,
279  gid,
280  ghost_degrees,
281  recolor_degrees);
282  }
283  //Entry point for serial conflict detection
284  virtual void detectConflicts_serial(const size_t n_local,
285  typename Kokkos::View<offset_t*, device_type >::HostMirror dist_offsets_host,
286  typename Kokkos::View<lno_t*, device_type >::HostMirror dist_adjs_host,
287  typename Kokkos::View<int*,device_type >::HostMirror femv_colors,
288  typename Kokkos::View<lno_t*, device_type >::HostMirror boundary_verts_view,
289  typename Kokkos::View<lno_t*,device_type>::HostMirror verts_to_recolor,
290  typename Kokkos::View<int*,device_type>::HostMirror verts_to_recolor_size,
291  typename Kokkos::View<lno_t*,device_type>::HostMirror verts_to_send,
292  typename Kokkos::View<size_t*,device_type>::HostMirror verts_to_send_size,
293  typename Kokkos::View<size_t*, device_type>::HostMirror recoloringSize,
294  typename Kokkos::View<int*, device_type>::HostMirror rand,
295  typename Kokkos::View<gno_t*,device_type>::HostMirror gid,
296  typename Kokkos::View<gno_t*,device_type>::HostMirror ghost_degrees,
297  bool recolor_degrees) {
298 
299  this->detectPD2Conflicts<host_exec, host_mem>(n_local,
300  dist_offsets_host,
301  dist_adjs_host,
302  femv_colors,
303  boundary_verts_view,
304  verts_to_recolor,
305  verts_to_recolor_size,
306  verts_to_send,
307  verts_to_send_size,
308  recoloringSize,
309  rand,
310  gid,
311  ghost_degrees,
312  recolor_degrees);
313  }
314  //Entry point for boundary construction
315  virtual void constructBoundary(const size_t n_local,
316  Kokkos::View<offset_t*, device_type> dist_offsets_dev,
317  Kokkos::View<lno_t*, device_type> dist_adjs_dev,
318  typename Kokkos::View<offset_t*, device_type>::HostMirror dist_offsets_host,
319  typename Kokkos::View<lno_t*, device_type>::HostMirror dist_adjs_host,
320  Kokkos::View<lno_t*, device_type>& boundary_verts,
321  Kokkos::View<lno_t*,
322  device_type > verts_to_send_view,
323  Kokkos::View<size_t*,
324  device_type,
325  Kokkos::MemoryTraits<Kokkos::Atomic>> verts_to_send_size_atomic){
326  //count the number of boundary vertices to correctly allocate
327  //the boundary vertex view on device
328  gno_t boundary_size_temp = 0;
329  for(size_t i = 0; i < n_local; i++){
330  for(offset_t j = dist_offsets_host(i); j < dist_offsets_host(i+1); j++){
331  if((size_t)dist_adjs_host(j) >= n_local){
332  boundary_size_temp++;
333  break;
334  }
335  bool found = false;
336  for(offset_t k = dist_offsets_host(dist_adjs_host(j)); k < dist_offsets_host(dist_adjs_host(j)+1); k++){
337  if((size_t)dist_adjs_host(k) >= n_local){
338  boundary_size_temp++;
339  found = true;
340  break;
341  }
342  }
343  if(found) break;
344  }
345  }
346 
347  //Initialize the boundary on host
348  boundary_verts = Kokkos::View<lno_t*, device_type>("boundary verts",boundary_size_temp);
349  typename Kokkos::View<lno_t*, device_type>::HostMirror boundary_verts_host = Kokkos::create_mirror_view(boundary_verts);
350 
351  //reset the boundary size count to use as an index to construct the view
352  boundary_size_temp = 0;
353 
354  for(size_t i = 0; i < n_local; i++){
355  for(offset_t j = dist_offsets_host(i); j < dist_offsets_host(i+1); j++){
356  if((size_t)dist_adjs_host(j) >= n_local){
357  boundary_verts_host(boundary_size_temp++) = i;
358  break;
359  }
360  bool found = false;
361  for(offset_t k = dist_offsets_host(dist_adjs_host(j)); k < dist_offsets_host(dist_adjs_host(j)+1); k++){
362  if((size_t)dist_adjs_host(k) >= n_local){
363  boundary_verts_host(boundary_size_temp++) = i;
364  found = true;
365  break;
366  }
367  }
368  if(found) break;
369  }
370  }
371  //copy boundary over to device views
372  Kokkos::deep_copy(boundary_verts, boundary_verts_host);
373 
374  //initialize the verts_to_send views
375  Kokkos::parallel_for("init verts to send",
376  Kokkos::RangePolicy<execution_space, int>(0,n_local),
377  KOKKOS_LAMBDA(const int& i){
378  for(offset_t j = dist_offsets_dev(i); j < dist_offsets_dev(i+1); j++){
379  if((size_t)dist_adjs_dev(j) >= n_local){
380  verts_to_send_view(verts_to_send_size_atomic(0)++) = i;
381  break;
382  }
383  bool found = false;
384  for(offset_t k = dist_offsets_dev(dist_adjs_dev(j)); k < dist_offsets_dev(dist_adjs_dev(j)+1); k++){
385  if((size_t)dist_adjs_dev(k) >= n_local){
386  verts_to_send_view(verts_to_send_size_atomic(0)++) = i;
387  found = true;
388  break;
389  }
390  }
391  if(found) break;
392  }
393  });
394  Kokkos::fence();
395 
396  }
397 
398 
399  public:
401  const RCP<const base_adapter_t> &adapter_,
402  const RCP<Teuchos::ParameterList> &pl_,
403  const RCP<Environment> &env_,
404  const RCP<const Teuchos::Comm<int> > &comm_)
405  : AlgTwoGhostLayer<Adapter>(adapter_,pl_,env_,comm_){}
406 
407 
408 
409 }; //end class
410 
411 
412 
413 }//end namespace Zoltan2
414 
415 #endif
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
Tpetra::FEMultiVector< femv_scalar_t, lno_t, gno_t > femv_t
typename Adapter::gno_t gno_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)
typename Adapter::lno_t lno_t
AlgPartialDistance2(const RCP< const base_adapter_t > &adapter_, const RCP< Teuchos::ParameterList > &pl_, const RCP< Environment > &env_, const RCP< const Teuchos::Comm< int > > &comm_)
map_t::global_ordinal_type gno_t
Definition: mapRemotes.cpp:18
typename Adapter::offset_t offset_t
typename femv_t::host_view_type::device_type::memory_space host_mem
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)
map_t::local_ordinal_type lno_t
Definition: mapRemotes.cpp:17
typename femv_t::device_type device_type
Traits class to handle conversions between gno_t/lno_t and TPL data types (e.g., ParMETIS&#39;s idx_t...
typename femv_t::host_view_type::device_type::execution_space host_exec
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)
Defines the ColoringSolution class.
Defines the GraphModel interface.
void detectPD2Conflicts(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)
A gathering of useful namespace methods.
typename Adapter::base_adapter_t base_adapter_t
typename device_type::execution_space execution_space
Tpetra::Map< lno_t, gno_t > map_t
typename Adapter::scalar_t scalar_t
typename device_type::memory_space memory_space
AlltoAll communication methods.
typename Adapter::offset_t offset_t