Zoltan2
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"
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
39namespace Zoltan2{
40
41template <typename Adapter>
42class 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 = 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 //serial and parallel local partial distance-2 coloring function
61 template<class ExecutionSpace, typename MemorySpace>
62 void localColoring(const size_t nVtx,
63 Kokkos::View<lno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> adjs_view,
64 Kokkos::View<offset_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> offset_view,
65 Teuchos::RCP<femv_t> femv,
66 Kokkos::View<lno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> vertex_list,
67 size_t vertex_list_size = 0,
68 bool use_vertex_based_coloring = false){
69 using KernelHandle = KokkosKernels::Experimental::KokkosKernelsHandle
70 <offset_t, lno_t, lno_t, ExecutionSpace, MemorySpace, MemorySpace>;
71 KernelHandle kh;
72
73 //Instead of switching between vertex-based and net-based algorithms,
74 //we use only the net-based algorithm, as it is faster than its
75 //vertex-based counterpart.
76 kh.create_distance2_graph_coloring_handle(KokkosGraph::COLORING_D2_NB_BIT);
77
78 //vertex_list_size indicates whether we have provided a list of vertices to recolor
79 //NB_BIT does not currently make use of this.
80 if(vertex_list_size != 0){
81 kh.get_distance2_graph_coloring_handle()->set_vertex_list(vertex_list, vertex_list_size);
82 }
83
84 //the verbose argument should carry through the local coloring
85 kh.get_distance2_graph_coloring_handle()->set_verbose(this->verbose);
86
87 //set initial colors to be the colors from femv
88 auto femvColors = femv->template getLocalView<Kokkos::Device<ExecutionSpace,MemorySpace> >(Tpetra::Access::ReadWrite);
89 auto sv = subview(femvColors,Kokkos::ALL, 0);
90 kh.get_distance2_graph_coloring_handle()->set_vertex_colors(sv);
91
92 //call coloring
93 KokkosGraph::Experimental::bipartite_color_rows(&kh, nVtx, nVtx, offset_view, adjs_view, true);
94
95
96 //output total time
97 if(this->verbose){
98 std::cout<<"\nKokkosKernels Coloring: "
99 <<kh.get_distance2_graph_coloring_handle()->get_overall_coloring_time()
100 <<"\n";
101 }
102 }
103
104 //Entry point for parallel local coloring
105 virtual void colorInterior(const size_t nVtx,
106 Kokkos::View<lno_t*, device_type> adjs_view,
107 Kokkos::View<offset_t*, device_type> offset_view,
108 Teuchos::RCP<femv_t> femv,
109 Kokkos::View<lno_t*, device_type> vertex_list,
110 size_t vertex_list_size=0,
111 bool recolor=false){
112 this->localColoring<execution_space, memory_space>(nVtx,
113 adjs_view,
114 offset_view,
115 femv,
116 vertex_list,
117 vertex_list_size,
118 recolor);
119 }
120 //Entry point for serial local coloring
121 virtual void colorInterior_serial(const size_t nVtx,
122 typename Kokkos::View<lno_t*, device_type >::HostMirror adjs_view,
123 typename Kokkos::View<offset_t*,device_type >::HostMirror offset_view,
124 Teuchos::RCP<femv_t> femv,
125 typename Kokkos::View<lno_t*, device_type>::HostMirror vertex_list,
126 size_t vertex_list_size = 0,
127 bool recolor=false) {
128 this->localColoring<host_exec, host_mem>(nVtx,
129 adjs_view,
130 offset_view,
131 femv,
132 vertex_list,
133 vertex_list_size,
134 recolor);
135
136 }
137 public:
138 //serial and parallel partial distance-2 conflict detection function
139 template <class ExecutionSpace, typename MemorySpace>
140 void detectPD2Conflicts(const size_t n_local,
141 Kokkos::View<offset_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> dist_offsets,
142 Kokkos::View<lno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> dist_adjs,
143 Kokkos::View<int*, Kokkos::Device<ExecutionSpace, MemorySpace>> femv_colors,
144 Kokkos::View<lno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> boundary_verts_view,
145 Kokkos::View<lno_t*,
146 Kokkos::Device<ExecutionSpace, MemorySpace> > verts_to_recolor_view,
147 Kokkos::View<int*,
148 Kokkos::Device<ExecutionSpace, MemorySpace>,
149 Kokkos::MemoryTraits<Kokkos::Atomic> > verts_to_recolor_size_atomic,
150 Kokkos::View<lno_t*,
151 Kokkos::Device<ExecutionSpace, MemorySpace> > verts_to_send_view,
152 Kokkos::View<size_t*,
153 Kokkos::Device<ExecutionSpace, MemorySpace>,
154 Kokkos::MemoryTraits<Kokkos::Atomic> > verts_to_send_size_atomic,
155 Kokkos::View<size_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> recoloringSize,
156 Kokkos::View<int*, Kokkos::Device<ExecutionSpace, MemorySpace>> rand,
157 Kokkos::View<gno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> gid,
158 Kokkos::View<gno_t*, Kokkos::Device<ExecutionSpace, MemorySpace>> ghost_degrees,
159 bool recolor_degrees){
160
161 Kokkos::RangePolicy<ExecutionSpace> policy(0,boundary_verts_view.extent(0));
162 size_t local_recoloring_size;
163 Kokkos::parallel_reduce("PD2 conflict detection",policy, KOKKOS_LAMBDA(const uint64_t& i,size_t& recoloring_size){
164 //we only detect conflicts for vertices in the boundary
165 const size_t curr_lid = boundary_verts_view(i);
166 const int curr_color = femv_colors(curr_lid);
167 const size_t vid_d1_adj_begin = dist_offsets(curr_lid);
168 const size_t vid_d1_adj_end = dist_offsets(curr_lid+1);
169 const size_t curr_degree = vid_d1_adj_end - vid_d1_adj_begin;
170 for(size_t vid_d1_adj = vid_d1_adj_begin; vid_d1_adj < vid_d1_adj_end; vid_d1_adj++){
171 //we skip direct distance-1 conflicts, and only resolve distance-2 conflicts.
172 size_t vid_d1 = dist_adjs(vid_d1_adj);
173 size_t d2_adj_begin = dist_offsets(vid_d1);
174 size_t d2_adj_end = dist_offsets(vid_d1+1);
175
176 //If we find a conflict that uncolors curr_lid, we can safely stop detecting
177 //further conflicts starting at curr_lid. Because this is a nested loop,
178 //we'll use the found boolean to break twice and move to the next vertex.
179 bool found = false;
180 for(size_t vid_d2_adj = d2_adj_begin; vid_d2_adj < d2_adj_end; vid_d2_adj++){
181 const size_t vid_d2 = dist_adjs(vid_d2_adj);
182 size_t vid_d2_degree = 0;
183
184 //calculate the degree for degree-based recoloring
185 if(vid_d2 < n_local){
186 vid_d2_degree = dist_offsets(vid_d2+1) - dist_offsets(vid_d2);
187 } else {
188 vid_d2_degree = ghost_degrees(vid_d2-n_local);
189 }
190 //resolve conflict
191 if(curr_lid != vid_d2 && femv_colors(vid_d2) == curr_color){
192 if(curr_degree < vid_d2_degree && recolor_degrees){
193 found = true;
194 femv_colors(curr_lid) = 0;
195 recoloring_size++;
196 break;//---------------------------------------------------
197 } else if(vid_d2_degree < curr_degree && recolor_degrees){//|
198 femv_colors(vid_d2) = 0; //|
199 recoloring_size++; //|
200 } else if(rand(curr_lid) < rand(vid_d2)){ //|
201 found = true; //|
202 femv_colors(curr_lid) = 0; //|
203 recoloring_size++; //|
204 break;//--------------------------------------------------|
205 } else if(rand(vid_d2) < rand(curr_lid)){ //|
206 femv_colors(vid_d2) = 0; //|
207 recoloring_size++; //|
208 } else { //|
209 if(gid(curr_lid) >= gid(vid_d2)){ //|
210 found = true; //|
211 femv_colors(curr_lid) = 0; //|
212 recoloring_size++; //|
213 break;//------------------------------------------------|
214 } else { //|
215 femv_colors(vid_d2) = 0; //|
216 recoloring_size++;// v
217 }// If we uncolor the vertex whose neighbors we're
218 } // checking, each subsequent conflict check will
219 } // not do anything productive. We need this------
220 } // to completely move on to the next vertex. |
221 if(found) break;//<--------------------------------------------------
222 }
223 },local_recoloring_size);
224 Kokkos::deep_copy(recoloringSize, local_recoloring_size);
225 Kokkos::fence();
226 //update the verts_to_send and verts_to_recolor views
227 Kokkos::parallel_for("rebuild verts_to_send and verts_to_recolor",
228 Kokkos::RangePolicy<ExecutionSpace>(0,femv_colors.size()),
229 KOKKOS_LAMBDA(const uint64_t& i){
230 if(femv_colors(i) == 0){
231 if(i < n_local){
232 //we only send vertices owned by the current process
233 verts_to_send_view(verts_to_send_size_atomic(0)++) = i;
234 }
235 //we need to recolor all vertices for consistency.
236 verts_to_recolor_view(verts_to_recolor_size_atomic(0)++) = i;
237 }
238 });
239 Kokkos::fence();
240
241 }
242 //Entry point for parallel conflict detection
243 virtual void detectConflicts(const size_t n_local,
244 Kokkos::View<offset_t*, device_type > dist_offsets_dev,
245 Kokkos::View<lno_t*, device_type > dist_adjs_dev,
246 Kokkos::View<int*,device_type > femv_colors,
247 Kokkos::View<lno_t*, device_type > boundary_verts_view,
248 Kokkos::View<lno_t*,
249 device_type > verts_to_recolor_view,
250 Kokkos::View<int*,
252 Kokkos::MemoryTraits<Kokkos::Atomic>> verts_to_recolor_size_atomic,
253 Kokkos::View<lno_t*,
254 device_type > verts_to_send_view,
255 Kokkos::View<size_t*,
257 Kokkos::MemoryTraits<Kokkos::Atomic>> verts_to_send_size_atomic,
258 Kokkos::View<size_t*, device_type> recoloringSize,
259 Kokkos::View<int*,
260 device_type> rand,
261 Kokkos::View<gno_t*,
262 device_type> gid,
263 Kokkos::View<gno_t*,
264 device_type> ghost_degrees,
265 bool recolor_degrees){
266
267 this->detectPD2Conflicts<execution_space, memory_space>(n_local,
268 dist_offsets_dev,
269 dist_adjs_dev,
270 femv_colors,
271 boundary_verts_view,
272 verts_to_recolor_view,
273 verts_to_recolor_size_atomic,
274 verts_to_send_view,
275 verts_to_send_size_atomic,
276 recoloringSize,
277 rand,
278 gid,
279 ghost_degrees,
280 recolor_degrees);
281 }
282 //Entry point for serial conflict detection
283 virtual void detectConflicts_serial(const size_t n_local,
284 typename Kokkos::View<offset_t*, device_type >::HostMirror dist_offsets_host,
285 typename Kokkos::View<lno_t*, device_type >::HostMirror dist_adjs_host,
286 typename Kokkos::View<int*,device_type >::HostMirror femv_colors,
287 typename Kokkos::View<lno_t*, device_type >::HostMirror boundary_verts_view,
288 typename Kokkos::View<lno_t*,device_type>::HostMirror verts_to_recolor,
289 typename Kokkos::View<int*,device_type>::HostMirror verts_to_recolor_size,
290 typename Kokkos::View<lno_t*,device_type>::HostMirror verts_to_send,
291 typename Kokkos::View<size_t*,device_type>::HostMirror verts_to_send_size,
292 typename Kokkos::View<size_t*, device_type>::HostMirror recoloringSize,
293 typename Kokkos::View<int*, device_type>::HostMirror rand,
294 typename Kokkos::View<gno_t*,device_type>::HostMirror gid,
295 typename Kokkos::View<gno_t*,device_type>::HostMirror ghost_degrees,
296 bool recolor_degrees) {
297
298 this->detectPD2Conflicts<host_exec, host_mem>(n_local,
299 dist_offsets_host,
300 dist_adjs_host,
301 femv_colors,
302 boundary_verts_view,
303 verts_to_recolor,
304 verts_to_recolor_size,
305 verts_to_send,
306 verts_to_send_size,
307 recoloringSize,
308 rand,
309 gid,
310 ghost_degrees,
311 recolor_degrees);
312 }
313 //Entry point for boundary construction
314 virtual void constructBoundary(const size_t n_local,
315 Kokkos::View<offset_t*, device_type> dist_offsets_dev,
316 Kokkos::View<lno_t*, device_type> dist_adjs_dev,
317 typename Kokkos::View<offset_t*, device_type>::HostMirror dist_offsets_host,
318 typename Kokkos::View<lno_t*, device_type>::HostMirror dist_adjs_host,
319 Kokkos::View<lno_t*, device_type>& boundary_verts,
320 Kokkos::View<lno_t*,
321 device_type > verts_to_send_view,
322 Kokkos::View<size_t*,
324 Kokkos::MemoryTraits<Kokkos::Atomic>> verts_to_send_size_atomic){
325 //count the number of boundary vertices to correctly allocate
326 //the boundary vertex view on device
327 gno_t boundary_size_temp = 0;
328 for(size_t i = 0; i < n_local; i++){
329 for(offset_t j = dist_offsets_host(i); j < dist_offsets_host(i+1); j++){
330 if((size_t)dist_adjs_host(j) >= n_local){
331 boundary_size_temp++;
332 break;
333 }
334 bool found = false;
335 for(offset_t k = dist_offsets_host(dist_adjs_host(j)); k < dist_offsets_host(dist_adjs_host(j)+1); k++){
336 if((size_t)dist_adjs_host(k) >= n_local){
337 boundary_size_temp++;
338 found = true;
339 break;
340 }
341 }
342 if(found) break;
343 }
344 }
345
346 //Initialize the boundary on host
347 boundary_verts = Kokkos::View<lno_t*, device_type>("boundary verts",boundary_size_temp);
348 typename Kokkos::View<lno_t*, device_type>::HostMirror boundary_verts_host = Kokkos::create_mirror_view(boundary_verts);
349
350 //reset the boundary size count to use as an index to construct the view
351 boundary_size_temp = 0;
352
353 for(size_t i = 0; i < n_local; i++){
354 for(offset_t j = dist_offsets_host(i); j < dist_offsets_host(i+1); j++){
355 if((size_t)dist_adjs_host(j) >= n_local){
356 boundary_verts_host(boundary_size_temp++) = i;
357 break;
358 }
359 bool found = false;
360 for(offset_t k = dist_offsets_host(dist_adjs_host(j)); k < dist_offsets_host(dist_adjs_host(j)+1); k++){
361 if((size_t)dist_adjs_host(k) >= n_local){
362 boundary_verts_host(boundary_size_temp++) = i;
363 found = true;
364 break;
365 }
366 }
367 if(found) break;
368 }
369 }
370 //copy boundary over to device views
371 Kokkos::deep_copy(boundary_verts, boundary_verts_host);
372
373 //initialize the verts_to_send views
374 Kokkos::parallel_for(n_local, KOKKOS_LAMBDA(const int& i){
375 for(offset_t j = dist_offsets_dev(i); j < dist_offsets_dev(i+1); j++){
376 if((size_t)dist_adjs_dev(j) >= n_local){
377 verts_to_send_view(verts_to_send_size_atomic(0)++) = i;
378 break;
379 }
380 bool found = false;
381 for(offset_t k = dist_offsets_dev(dist_adjs_dev(j)); k < dist_offsets_dev(dist_adjs_dev(j)+1); k++){
382 if((size_t)dist_adjs_dev(k) >= n_local){
383 verts_to_send_view(verts_to_send_size_atomic(0)++) = i;
384 found = true;
385 break;
386 }
387 }
388 if(found) break;
389 }
390 });
391 Kokkos::fence();
392
393 }
394
395
396 public:
398 const RCP<const base_adapter_t> &adapter_,
399 const RCP<Teuchos::ParameterList> &pl_,
400 const RCP<Environment> &env_,
401 const RCP<const Teuchos::Comm<int> > &comm_)
402 : AlgTwoGhostLayer<Adapter>(adapter_,pl_,env_,comm_){}
403
404
405
406}; //end class
407
408
409
410}//end namespace Zoltan2
411
412#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.
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)
AlgPartialDistance2(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 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)
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)
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)
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.