MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_DroppingCommon.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_DROPPINGCOMMON_HPP
11#define MUELU_DROPPINGCOMMON_HPP
12
13#include "Kokkos_Core.hpp"
14#include "Kokkos_ArithTraits.hpp"
15#include "Xpetra_Access.hpp"
16#include "Xpetra_Matrix.hpp"
17
18namespace MueLu {
19
24enum DecisionType : char {
25 UNDECIDED = 0, // no decision has been taken yet, used for initialization
26 KEEP = 1, // keeep the entry
27 DROP = 2, // drop it
28 BOUNDARY = 3 // entry is a boundary
29};
30
31namespace Misc {
32
33template <class local_ordinal_type>
35 public:
37
38 KOKKOS_FORCEINLINE_FUNCTION
39 void operator()(local_ordinal_type rlid) const {
40 }
41};
42
47template <class local_matrix_type>
49 private:
50 using scalar_type = typename local_matrix_type::value_type;
51 using local_ordinal_type = typename local_matrix_type::ordinal_type;
52 using memory_space = typename local_matrix_type::memory_space;
53 using results_view = Kokkos::View<DecisionType*, memory_space>;
54 using boundary_nodes_view = Kokkos::View<const bool*, memory_space>;
55
56 local_matrix_type A;
59
60 public:
61 PointwiseDropBoundaryFunctor(local_matrix_type& A_, boundary_nodes_view boundaryNodes_, results_view& results_)
62 : A(A_)
63 , boundaryNodes(boundaryNodes_)
64 , results(results_) {}
65
66 KOKKOS_FORCEINLINE_FUNCTION
67 void operator()(local_ordinal_type rlid) const {
68 auto row = A.rowConst(rlid);
69 const size_t offset = A.graph.row_map(rlid);
70 const bool isBoundaryRow = boundaryNodes(rlid);
71 if (isBoundaryRow) {
72 for (local_ordinal_type k = 0; k < row.length; ++k) {
73 auto clid = row.colidx(k);
74 results(offset + k) = Kokkos::max(rlid == clid ? KEEP : DROP,
75 results(offset + k));
76 }
77 }
78 }
79};
80
85template <class local_matrix_type>
87 private:
88 using scalar_type = typename local_matrix_type::value_type;
89 using local_ordinal_type = typename local_matrix_type::ordinal_type;
90 using memory_space = typename local_matrix_type::memory_space;
91 using results_view = Kokkos::View<DecisionType*, memory_space>;
92 using boundary_nodes_view = Kokkos::View<const bool*, memory_space>;
93 using block_indices_view_type = Kokkos::View<local_ordinal_type*, memory_space>;
94
95 local_matrix_type A;
99
100 public:
101 VectorDropBoundaryFunctor(local_matrix_type& A_, block_indices_view_type point_to_block_, boundary_nodes_view boundaryNodes_, results_view& results_)
102 : A(A_)
103 , point_to_block(point_to_block_)
104 , boundaryNodes(boundaryNodes_)
105 , results(results_) {}
106
107 KOKKOS_FORCEINLINE_FUNCTION
109 auto row = A.rowConst(rlid);
110 const size_t offset = A.graph.row_map(rlid);
111 const bool isBoundaryRow = boundaryNodes(point_to_block(rlid));
112 if (isBoundaryRow) {
113 for (local_ordinal_type k = 0; k < row.length; ++k) {
114 auto clid = row.colidx(k);
115 results(offset + k) = Kokkos::max(rlid == clid ? KEEP : DROP,
116 results(offset + k));
117 }
118 }
119 }
120};
121
126template <class local_matrix_type>
128 private:
129 using scalar_type = typename local_matrix_type::value_type;
130 using local_ordinal_type = typename local_matrix_type::ordinal_type;
131 using memory_space = typename local_matrix_type::memory_space;
132 using results_view = Kokkos::View<DecisionType*, memory_space>;
133
134 local_matrix_type A;
136
137 public:
138 KeepDiagonalFunctor(local_matrix_type& A_, results_view& results_)
139 : A(A_)
140 , results(results_) {}
141
142 KOKKOS_FORCEINLINE_FUNCTION
144 auto row = A.rowConst(rlid);
145 const size_t offset = A.graph.row_map(rlid);
146 for (local_ordinal_type k = 0; k < row.length; ++k) {
147 auto clid = row.colidx(k);
148 if ((rlid == clid) && (results(offset + k) != BOUNDARY)) {
149 results(offset + k) = KEEP;
150 break;
151 }
152 }
153 }
154};
155
160template <class local_matrix_type>
162 private:
163 using scalar_type = typename local_matrix_type::value_type;
164 using local_ordinal_type = typename local_matrix_type::ordinal_type;
165 using memory_space = typename local_matrix_type::memory_space;
166 using results_view = Kokkos::View<DecisionType*, memory_space>;
167
168 local_matrix_type A;
170
171 public:
172 DropOffRankFunctor(local_matrix_type& A_, results_view& results_)
173 : A(A_)
174 , results(results_) {}
175
176 KOKKOS_FORCEINLINE_FUNCTION
178 auto row = A.rowConst(rlid);
179 const size_t offset = A.graph.row_map(rlid);
180 for (local_ordinal_type k = 0; k < row.length; ++k) {
181 auto clid = row.colidx(k);
182 if (clid >= A.numRows()) {
183 results(offset + k) = Kokkos::max(DROP, results(offset + k));
184 }
185 }
186 }
187};
188
193template <class local_matrix_type>
195 private:
196 using scalar_type = typename local_matrix_type::value_type;
197 using local_ordinal_type = typename local_matrix_type::ordinal_type;
198 using memory_space = typename local_matrix_type::memory_space;
199 using results_view = Kokkos::View<DecisionType*, memory_space>;
200
201 using boundary_nodes_view = Kokkos::View<bool*, memory_space>;
202
203 local_matrix_type A;
206
207 public:
208 MarkSingletonFunctor(local_matrix_type& A_, boundary_nodes_view boundaryNodes_, results_view& results_)
209 : A(A_)
210 , boundaryNodes(boundaryNodes_)
211 , results(results_) {}
212
213 KOKKOS_FORCEINLINE_FUNCTION
215 auto row = A.rowConst(rlid);
216 const size_t offset = A.graph.row_map(rlid);
217 for (local_ordinal_type k = 0; k < row.length; ++k) {
218 auto clid = row.colidx(k);
219 if ((results(offset + k) == KEEP) && (rlid != clid))
220 return;
221 }
222 boundaryNodes(rlid) = true;
223 for (local_ordinal_type k = 0; k < row.length; ++k) {
224 auto clid = row.colidx(k);
225 if (rlid == clid)
226 results(offset + k) = KEEP;
227 else
228 results(offset + k) = BOUNDARY;
229 }
230 }
231};
232
237template <class local_matrix_type>
239 private:
240 using scalar_type = typename local_matrix_type::value_type;
241 using local_ordinal_type = typename local_matrix_type::ordinal_type;
242 using memory_space = typename local_matrix_type::memory_space;
243 using results_view = Kokkos::View<DecisionType*, memory_space>;
244 using block_indices_view_type = Kokkos::View<local_ordinal_type*, memory_space>;
245
246 using boundary_nodes_view = Kokkos::View<bool*, memory_space>;
247
248 local_matrix_type A;
252
253 public:
254 MarkSingletonVectorFunctor(local_matrix_type& A_, block_indices_view_type point_to_block_, boundary_nodes_view boundaryNodes_, results_view& results_)
255 : A(A_)
256 , point_to_block(point_to_block_)
257 , boundaryNodes(boundaryNodes_)
258 , results(results_) {}
259
260 KOKKOS_FORCEINLINE_FUNCTION
262 auto row = A.rowConst(rlid);
263 const size_t offset = A.graph.row_map(rlid);
264 for (local_ordinal_type k = 0; k < row.length; ++k) {
265 auto clid = row.colidx(k);
266 if ((results(offset + k) == KEEP) && (rlid != clid))
267 return;
268 }
269 auto brlid = point_to_block(rlid);
270 boundaryNodes(brlid) = true;
271 for (local_ordinal_type k = 0; k < row.length; ++k) {
272 auto clid = row.colidx(k);
273 if (rlid == clid)
274 results(offset + k) = KEEP;
275 else
276 results(offset + k) = BOUNDARY;
277 }
278 }
279};
280
285template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
287 private:
288 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
289 using local_matrix_type = typename matrix_type::local_matrix_type;
290
291 using scalar_type = typename local_matrix_type::value_type;
292 using local_ordinal_type = typename local_matrix_type::ordinal_type;
293 using memory_space = typename local_matrix_type::memory_space;
294 using results_view = Kokkos::View<DecisionType*, memory_space>;
295
296 using block_indices_type = Xpetra::MultiVector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>;
297 using local_block_indices_view_type = typename block_indices_type::dual_view_type_const::t_dev;
298
303
304 public:
305 BlockDiagonalizeFunctor(matrix_type& A_, block_indices_type& point_to_block_, block_indices_type& ghosted_point_to_block_, results_view& results_)
306 : A(A_.getLocalMatrixDevice())
307 , point_to_block(point_to_block_.getDeviceLocalView(Xpetra::Access::ReadOnly))
308 , ghosted_point_to_block(ghosted_point_to_block_.getDeviceLocalView(Xpetra::Access::ReadOnly))
309 , results(results_) {}
310
311 KOKKOS_FORCEINLINE_FUNCTION
313 auto row = A.rowConst(rlid);
314 const size_t offset = A.graph.row_map(rlid);
315 for (local_ordinal_type k = 0; k < row.length; ++k) {
316 auto clid = row.colidx(k);
317 if (point_to_block(rlid, 0) == ghosted_point_to_block(clid, 0)) {
318 results(offset + k) = Kokkos::max(KEEP, results(offset + k));
319 } else {
320 results(offset + k) = Kokkos::max(DROP, results(offset + k));
321 }
322 }
323 }
324};
325
330template <class local_matrix_type>
332 private:
333 using scalar_type = typename local_matrix_type::value_type;
334 using local_ordinal_type = typename local_matrix_type::ordinal_type;
335 using memory_space = typename local_matrix_type::memory_space;
336 using results_view = Kokkos::View<DecisionType*, memory_space>;
337
338 using boundary_nodes_view = Kokkos::View<bool*, memory_space>;
339
340 local_matrix_type A;
342
343 public:
344 DebugFunctor(local_matrix_type& A_, results_view& results_)
345 : A(A_)
346 , results(results_) {}
347
348 KOKKOS_FORCEINLINE_FUNCTION
350 auto row = A.rowConst(rlid);
351 const size_t offset = A.graph.row_map(rlid);
352 for (local_ordinal_type k = 0; k < row.length; ++k) {
353 if (results(offset + k) == UNDECIDED) {
354 Kokkos::printf("No dropping decision was taken for entry (%d, %d)\n", rlid, row.colidx(k));
355 assert(false);
356 }
357 }
358 }
359};
360
365template <class local_matrix_type>
367 private:
368 using scalar_type = typename local_matrix_type::value_type;
369 using local_ordinal_type = typename local_matrix_type::ordinal_type;
370 using memory_space = typename local_matrix_type::memory_space;
371 using results_view = Kokkos::View<DecisionType*, memory_space>;
372
373 local_matrix_type A;
375
376 public:
377 SymmetrizeFunctor(local_matrix_type& A_, results_view& results_)
378 : A(A_)
379 , results(results_) {}
380
381 KOKKOS_FORCEINLINE_FUNCTION
383 auto row = A.rowConst(rlid);
384 const size_t offset = A.graph.row_map(rlid);
385 for (local_ordinal_type k = 0; k < row.length; ++k) {
386 if (results(offset + k) == KEEP) {
387 auto clid = row.colidx(k);
388 if (clid >= A.numRows())
389 continue;
390 auto row2 = A.rowConst(clid);
391 const size_t offset2 = A.graph.row_map(clid);
392 for (local_ordinal_type k2 = 0; k2 < row2.length; ++k2) {
393 auto clid2 = row2.colidx(k2);
394 if (clid2 == rlid) {
395 if (results(offset2 + k2) == DROP)
396 results(offset2 + k2) = KEEP;
397 break;
398 }
399 }
400 }
401 }
402 }
403};
404
405template <class view_type, class comparator_type>
406KOKKOS_INLINE_FUNCTION void serialHeapSort(view_type& v, comparator_type comparator) {
407 auto N = v.extent(0);
408 size_t start = N / 2;
409 size_t end = N;
410 while (end > 1) {
411 if (start > 0)
412 start = start - 1;
413 else {
414 end = end - 1;
415 auto temp = v(0);
416 v(0) = v(end);
417 v(end) = temp;
418 }
419 size_t root = start;
420 while (2 * root + 1 < end) {
421 size_t child = 2 * root + 1;
422 if ((child + 1 < end) and (comparator(v(child), v(child + 1))))
423 ++child;
424
425 if (comparator(v(root), v(child))) {
426 auto temp = v(root);
427 v(root) = v(child);
428 v(child) = temp;
429 root = child;
430 } else
431 break;
432 }
433 }
434}
435
436} // namespace Misc
437
438} // namespace MueLu
439
440#endif
Xpetra::MultiVector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > block_indices_type
typename local_matrix_type::value_type scalar_type
local_block_indices_view_type point_to_block
typename matrix_type::local_matrix_type local_matrix_type
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
local_block_indices_view_type ghosted_point_to_block
typename block_indices_type::dual_view_type_const::t_dev local_block_indices_view_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename local_matrix_type::memory_space memory_space
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< DecisionType *, memory_space > results_view
BlockDiagonalizeFunctor(matrix_type &A_, block_indices_type &point_to_block_, block_indices_type &ghosted_point_to_block_, results_view &results_)
Kokkos::View< DecisionType *, memory_space > results_view
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
DebugFunctor(local_matrix_type &A_, results_view &results_)
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< bool *, memory_space > boundary_nodes_view
typename local_matrix_type::value_type scalar_type
typename local_matrix_type::memory_space memory_space
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Kokkos::View< DecisionType *, memory_space > results_view
typename local_matrix_type::ordinal_type local_ordinal_type
typename local_matrix_type::memory_space memory_space
typename local_matrix_type::value_type scalar_type
DropOffRankFunctor(local_matrix_type &A_, results_view &results_)
Kokkos::View< DecisionType *, memory_space > results_view
typename local_matrix_type::memory_space memory_space
KeepDiagonalFunctor(local_matrix_type &A_, results_view &results_)
typename local_matrix_type::ordinal_type local_ordinal_type
typename local_matrix_type::value_type scalar_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Kokkos::View< DecisionType *, memory_space > results_view
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename local_matrix_type::value_type scalar_type
MarkSingletonFunctor(local_matrix_type &A_, boundary_nodes_view boundaryNodes_, results_view &results_)
typename local_matrix_type::memory_space memory_space
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< bool *, memory_space > boundary_nodes_view
MarkSingletonVectorFunctor(local_matrix_type &A_, block_indices_view_type point_to_block_, boundary_nodes_view boundaryNodes_, results_view &results_)
Kokkos::View< DecisionType *, memory_space > results_view
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename local_matrix_type::value_type scalar_type
Kokkos::View< local_ordinal_type *, memory_space > block_indices_view_type
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< bool *, memory_space > boundary_nodes_view
typename local_matrix_type::memory_space memory_space
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename local_matrix_type::ordinal_type local_ordinal_type
typename local_matrix_type::memory_space memory_space
Kokkos::View< const bool *, memory_space > boundary_nodes_view
typename local_matrix_type::value_type scalar_type
Kokkos::View< DecisionType *, memory_space > results_view
PointwiseDropBoundaryFunctor(local_matrix_type &A_, boundary_nodes_view boundaryNodes_, results_view &results_)
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename local_matrix_type::value_type scalar_type
SymmetrizeFunctor(local_matrix_type &A_, results_view &results_)
Kokkos::View< DecisionType *, memory_space > results_view
typename local_matrix_type::memory_space memory_space
typename local_matrix_type::ordinal_type local_ordinal_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
VectorDropBoundaryFunctor(local_matrix_type &A_, block_indices_view_type point_to_block_, boundary_nodes_view boundaryNodes_, results_view &results_)
Kokkos::View< local_ordinal_type *, memory_space > block_indices_view_type
Kokkos::View< const bool *, memory_space > boundary_nodes_view
typename local_matrix_type::memory_space memory_space
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< DecisionType *, memory_space > results_view
typename local_matrix_type::value_type scalar_type
KOKKOS_INLINE_FUNCTION void serialHeapSort(view_type &v, comparator_type comparator)
Namespace for MueLu classes and methods.