MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_CutDrop.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_CUTDROP_HPP
11#define MUELU_CUTDROP_HPP
12
13#include "Kokkos_Core.hpp"
14#include "Kokkos_ArithTraits.hpp"
16#include "MueLu_Utilities.hpp"
17#include "Xpetra_Matrix.hpp"
18#include "Xpetra_MultiVector.hpp"
20
21namespace MueLu::CutDrop {
22
28
33template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
35 public:
36 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
37
38 using local_matrix_type = typename matrix_type::local_matrix_type;
39 using scalar_type = typename local_matrix_type::value_type;
40 using local_ordinal_type = typename local_matrix_type::ordinal_type;
41 using memory_space = typename local_matrix_type::memory_space;
42 using results_view = Kokkos::View<DecisionType*, memory_space>;
43
46
47 private:
48 using ATS = Kokkos::ArithTraits<scalar_type>;
49 using magnitudeType = typename ATS::magnitudeType;
50
51 public:
53 : A(A_.getLocalMatrixDevice())
54 , results(results_) {}
55
56 template <class local_matrix_type2>
57 struct Comparator {
58 private:
59 using scalar_type = typename local_matrix_type2::value_type;
60 using local_ordinal_type = typename local_matrix_type2::ordinal_type;
61 using memory_space = typename local_matrix_type2::memory_space;
62 using results_view = Kokkos::View<DecisionType*, memory_space>;
63
64 using ATS = Kokkos::ArithTraits<scalar_type>;
65 using magnitudeType = typename ATS::magnitudeType;
66
67 const local_matrix_type2 A;
70
71 public:
72 KOKKOS_INLINE_FUNCTION
73 Comparator(const local_matrix_type2& A_, local_ordinal_type rlid_, const results_view& results_)
74 : A(A_)
75 , offset(A_.graph.row_map(rlid_))
76 , results(results_) {}
77
78 KOKKOS_INLINE_FUNCTION
79 magnitudeType get_value(size_t x) const {
80 return ATS::magnitude(A.values(offset + x) * A.values(offset + x));
81 }
82
83 KOKKOS_INLINE_FUNCTION
84 bool operator()(size_t x, size_t y) const {
85 if (results(offset + x) != UNDECIDED) {
86 if (results(offset + y) != UNDECIDED) {
87 // does not matter
88 return (x < y);
89 } else {
90 // sort undecided to the right
91 return true;
92 }
93 } else {
94 if (results(offset + y) != UNDECIDED) {
95 // sort undecided to the right
96 return false;
97 } else {
98 return get_value(x) > get_value(y);
99 }
100 }
101 }
102 };
103
105
106 KOKKOS_INLINE_FUNCTION
110};
111
116template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
118 public:
119 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
120 using local_matrix_type = typename matrix_type::local_matrix_type;
121 using scalar_type = typename local_matrix_type::value_type;
122 using local_ordinal_type = typename local_matrix_type::ordinal_type;
123 using memory_space = typename local_matrix_type::memory_space;
124 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
125 using diag_view_type = typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
126 using results_view = Kokkos::View<DecisionType*, memory_space>;
127
130
131 private:
132 using ATS = Kokkos::ArithTraits<scalar_type>;
133 using magnitudeType = typename ATS::magnitudeType;
134
135 Teuchos::RCP<diag_vec_type> diagVec;
137
138 public:
140 : A(A_.getLocalMatrixDevice())
141 , results(results_) {
143 auto lclDiag2d = diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
144 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
145 }
146
147 template <class local_matrix_type2, class diag_view_type2>
148 struct Comparator {
149 private:
150 using scalar_type = typename local_matrix_type2::value_type;
151 using local_ordinal_type = typename local_matrix_type2::ordinal_type;
152 using memory_space = typename local_matrix_type2::memory_space;
153 using results_view = Kokkos::View<DecisionType*, memory_space>;
154
155 using ATS = Kokkos::ArithTraits<scalar_type>;
156 using magnitudeType = typename ATS::magnitudeType;
157
158 const local_matrix_type2 A;
159 const diag_view_type2 diag;
163
164 public:
165 KOKKOS_INLINE_FUNCTION
166 Comparator(const local_matrix_type2& A_, const diag_view_type2& diag_, const local_ordinal_type rlid_, const results_view& results_)
167 : A(A_)
168 , diag(diag_)
169 , rlid(rlid_)
170 , offset(A_.graph.row_map(rlid_))
171 , results(results_) {}
172
173 KOKKOS_INLINE_FUNCTION
174 magnitudeType get_value(size_t x) const {
175 auto x_aij = ATS::magnitude(A.values(offset + x) * A.values(offset + x));
176 auto x_aiiajj = ATS::magnitude(diag(rlid) * diag(A.graph.entries(offset + x)));
177 return (x_aij / x_aiiajj);
178 }
179
180 KOKKOS_INLINE_FUNCTION
181 bool operator()(size_t x, size_t y) const {
182 if (results(offset + x) != UNDECIDED) {
183 if (results(offset + y) != UNDECIDED) {
184 // does not matter
185 return (x < y);
186 } else {
187 // sort undecided to the right
188 return true;
189 }
190 } else {
191 if (results(offset + y) != UNDECIDED) {
192 // sort undecided to the right
193 return false;
194 } else {
195 return get_value(x) > get_value(y);
196 }
197 }
198 }
199 };
200
202
203 KOKKOS_INLINE_FUNCTION
207};
208
209template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class DistanceFunctorType>
211 public:
212 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
213 using local_matrix_type = typename matrix_type::local_matrix_type;
214 using scalar_type = typename local_matrix_type::value_type;
215 using local_ordinal_type = typename local_matrix_type::ordinal_type;
216 using memory_space = typename local_matrix_type::memory_space;
217 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
218 using diag_view_type = typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
219 using results_view = Kokkos::View<DecisionType*, memory_space>;
220
223
224 private:
225 using ATS = Kokkos::ArithTraits<scalar_type>;
226 using magnitudeType = typename ATS::magnitudeType;
227
228 Teuchos::RCP<diag_vec_type> diagVec;
230 DistanceFunctorType dist2;
231
232 public:
233 UnscaledDistanceLaplacianComparison(matrix_type& A_, DistanceFunctorType& dist2_, results_view& results_)
234 : A(A_.getLocalMatrixDevice())
235 , results(results_)
236 , dist2(dist2_) {
237 // Construct ghosted distance Laplacian diagonal
239 auto lclDiag2d = diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
240 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
241 }
242
243 template <class local_matrix_type2, class DistanceFunctorType2, class diag_view_type2>
244 struct Comparator {
245 private:
246 using scalar_type = typename local_matrix_type2::value_type;
247 using local_ordinal_type = typename local_matrix_type2::ordinal_type;
248 using memory_space = typename local_matrix_type2::memory_space;
249 using results_view = Kokkos::View<DecisionType*, memory_space>;
250
251 using ATS = Kokkos::ArithTraits<scalar_type>;
252 using magnitudeType = typename ATS::magnitudeType;
253
254 const local_matrix_type2 A;
255 const diag_view_type2 diag;
256 const DistanceFunctorType2* dist2;
260
261 const scalar_type one = ATS::one();
262
263 public:
264 KOKKOS_INLINE_FUNCTION
265 Comparator(const local_matrix_type2& A_, const diag_view_type2& diag_, const DistanceFunctorType2* dist2_, local_ordinal_type rlid_, const results_view& results_)
266 : A(A_)
267 , diag(diag_)
268 , dist2(dist2_)
269 , rlid(rlid_)
270 , offset(A_.graph.row_map(rlid_))
271 , results(results_) {}
272
273 KOKKOS_INLINE_FUNCTION
274 magnitudeType get_value(size_t x) const {
275 auto clid = A.graph.entries(offset + x);
276 scalar_type val;
277 if (rlid != clid) {
278 val = one / dist2->distance2(rlid, clid);
279 } else {
280 val = diag(rlid);
281 }
282 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
283 return aij2;
284 }
285
286 KOKKOS_INLINE_FUNCTION
287 bool operator()(size_t x, size_t y) const {
288 if (results(offset + x) != UNDECIDED) {
289 if (results(offset + y) != UNDECIDED) {
290 // does not matter
291 return (x < y);
292 } else {
293 // sort undecided to the right
294 return true;
295 }
296 } else {
297 if (results(offset + y) != UNDECIDED) {
298 // sort undecided to the right
299 return false;
300 } else {
301 return get_value(x) > get_value(y);
302 }
303 }
304 }
305 };
306
308
309 KOKKOS_INLINE_FUNCTION
313};
314
319template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class DistanceFunctorType>
321 public:
322 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
323 using local_matrix_type = typename matrix_type::local_matrix_type;
324 using scalar_type = typename local_matrix_type::value_type;
325 using local_ordinal_type = typename local_matrix_type::ordinal_type;
326 using memory_space = typename local_matrix_type::memory_space;
327 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
328 using diag_view_type = typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
329 using results_view = Kokkos::View<DecisionType*, memory_space>;
330
333
334 private:
335 using ATS = Kokkos::ArithTraits<scalar_type>;
336 using magnitudeType = typename ATS::magnitudeType;
337
338 Teuchos::RCP<diag_vec_type> diagVec;
340 DistanceFunctorType dist2;
341
342 public:
343 ScaledDistanceLaplacianComparison(matrix_type& A_, DistanceFunctorType& dist2_, results_view& results_)
344 : A(A_.getLocalMatrixDevice())
345 , results(results_)
346 , dist2(dist2_) {
347 // Construct ghosted distance Laplacian diagonal
349 auto lclDiag2d = diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
350 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
351 }
352
353 template <class local_matrix_type2, class DistanceFunctorType2, class diag_view_type2>
354 struct Comparator {
355 private:
356 using scalar_type = typename local_matrix_type2::value_type;
357 using local_ordinal_type = typename local_matrix_type2::ordinal_type;
358 using memory_space = typename local_matrix_type2::memory_space;
359 using results_view = Kokkos::View<DecisionType*, memory_space>;
360
361 using ATS = Kokkos::ArithTraits<scalar_type>;
362 using magnitudeType = typename ATS::magnitudeType;
363
364 const local_matrix_type2 A;
365 const diag_view_type2 diag;
366 const DistanceFunctorType2* dist2;
370
371 const scalar_type one = ATS::one();
372
373 public:
374 KOKKOS_INLINE_FUNCTION
375 Comparator(const local_matrix_type2& A_, const diag_view_type2& diag_, const DistanceFunctorType2* dist2_, local_ordinal_type rlid_, const results_view& results_)
376 : A(A_)
377 , diag(diag_)
378 , dist2(dist2_)
379 , rlid(rlid_)
380 , offset(A_.graph.row_map(rlid_))
381 , results(results_) {}
382
383 KOKKOS_INLINE_FUNCTION
384 magnitudeType get_value(size_t x) const {
385 auto clid = A.graph.entries(offset + x);
386 scalar_type val;
387 if (rlid != clid) {
388 val = one / dist2->distance2(rlid, clid);
389 } else {
390 val = diag(rlid);
391 }
392 auto aiiajj = ATS::magnitude(diag(rlid)) * ATS::magnitude(diag(clid)); // |a_ii|*|a_jj|
393 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
394 return (aij2 / aiiajj);
395 }
396
397 KOKKOS_INLINE_FUNCTION
398 bool operator()(size_t x, size_t y) const {
399 if (results(offset + x) != UNDECIDED) {
400 if (results(offset + y) != UNDECIDED) {
401 // does not matter
402 return (x < y);
403 } else {
404 // sort undecided to the right
405 return true;
406 }
407 } else {
408 if (results(offset + y) != UNDECIDED) {
409 // sort undecided to the right
410 return false;
411 } else {
412 return get_value(x) > get_value(y);
413 }
414 }
415 }
416 };
417
419
420 KOKKOS_INLINE_FUNCTION
424};
425
430template <class comparison_type>
432 private:
433 using local_matrix_type = typename comparison_type::local_matrix_type;
434 using scalar_type = typename local_matrix_type::value_type;
435 using local_ordinal_type = typename local_matrix_type::ordinal_type;
436 using memory_space = typename local_matrix_type::memory_space;
437 using results_view = Kokkos::View<DecisionType*, memory_space>;
438
439 using ATS = Kokkos::ArithTraits<scalar_type>;
440 using magnitudeType = typename ATS::magnitudeType;
441 using boundary_nodes_view = Kokkos::View<const bool*, memory_space>;
442
444 comparison_type comparison;
447 Kokkos::View<local_ordinal_type*, memory_space> index;
448
449 public:
450 CutDropFunctor(comparison_type& comparison_, magnitudeType threshold)
451 : A(comparison_.A)
452 , comparison(comparison_)
453 , eps(threshold)
454 , results(comparison_.results) {
455 index = Kokkos::View<local_ordinal_type*, memory_space>("indices", A.nnz());
456 }
457
458 KOKKOS_FORCEINLINE_FUNCTION
459 void operator()(const local_ordinal_type& rlid) const {
460 auto row = A.rowConst(rlid);
461 size_t nnz = row.length;
462
463 auto drop_view = Kokkos::subview(results, Kokkos::make_pair(A.graph.row_map(rlid), A.graph.row_map(rlid + 1)));
464 auto row_permutation = Kokkos::subview(index, Kokkos::make_pair(A.graph.row_map(rlid), A.graph.row_map(rlid + 1)));
465
466 auto comparator = comparison.getComparator(rlid);
467
468 for (size_t i = 0; i < nnz; ++i) {
469 row_permutation(i) = i;
470 }
471 Misc::serialHeapSort(row_permutation, comparator);
472
473 size_t keepStart = 0;
474 size_t dropStart = nnz;
475 // find index where dropping starts
476 for (size_t i = 1; i < nnz; ++i) {
477 auto const& x = row_permutation(i - 1);
478 auto const& y = row_permutation(i);
479 if ((drop_view(x) != UNDECIDED) && (drop_view(y) == UNDECIDED))
480 keepStart = i;
481 if ((drop_view(x) != UNDECIDED) || (drop_view(y) != UNDECIDED))
482 continue;
483 magnitudeType x_aij = comparator.get_value(x);
484 magnitudeType y_aij = comparator.get_value(y);
485 if (eps * eps * x_aij > y_aij) {
486 if (i < dropStart) {
487 dropStart = i;
488 }
489 }
490 }
491
492 // drop everything to the right of where values stop passing threshold
493 for (size_t i = keepStart; i < nnz; ++i) {
494 drop_view(row_permutation(i)) = Kokkos::max(dropStart <= i ? DROP : KEEP, drop_view(row_permutation(i)));
495 }
496 }
497};
498
499} // namespace MueLu::CutDrop
500
501#endif
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< const bool *, memory_space > boundary_nodes_view
KOKKOS_FORCEINLINE_FUNCTION void operator()(const local_ordinal_type &rlid) const
CutDropFunctor(comparison_type &comparison_, magnitudeType threshold)
Kokkos::View< DecisionType *, memory_space > results_view
typename local_matrix_type::memory_space memory_space
Kokkos::ArithTraits< scalar_type > ATS
typename comparison_type::local_matrix_type local_matrix_type
typename ATS::magnitudeType magnitudeType
Kokkos::View< local_ordinal_type *, memory_space > index
typename local_matrix_type::value_type scalar_type
ScaledComparison(matrix_type &A_, results_view &results_)
typename Kokkos::DualView< const scalar_type *, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged >::t_dev diag_view_type
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
Teuchos::RCP< diag_vec_type > diagVec
Kokkos::ArithTraits< scalar_type > ATS
typename matrix_type::local_matrix_type local_matrix_type
Kokkos::View< DecisionType *, memory_space > results_view
typename ATS::magnitudeType magnitudeType
typename local_matrix_type::value_type scalar_type
typename local_matrix_type::memory_space memory_space
Comparator< local_matrix_type, diag_view_type > comparator_type
typename local_matrix_type::ordinal_type local_ordinal_type
KOKKOS_INLINE_FUNCTION comparator_type getComparator(local_ordinal_type rlid) const
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
KOKKOS_INLINE_FUNCTION comparator_type getComparator(local_ordinal_type rlid) const
ScaledDistanceLaplacianComparison(matrix_type &A_, DistanceFunctorType &dist2_, results_view &results_)
typename local_matrix_type::value_type scalar_type
typename Kokkos::DualView< const scalar_type *, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged >::t_dev diag_view_type
typename local_matrix_type::ordinal_type local_ordinal_type
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
typename local_matrix_type::memory_space memory_space
typename matrix_type::local_matrix_type local_matrix_type
Comparator< local_matrix_type, DistanceFunctorType, diag_view_type > comparator_type
Kokkos::View< DecisionType *, memory_space > results_view
typename local_matrix_type::ordinal_type local_ordinal_type
typename ATS::magnitudeType magnitudeType
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
UnscaledComparison(matrix_type &A_, results_view &results_)
typename local_matrix_type::memory_space memory_space
Comparator< local_matrix_type > comparator_type
typename local_matrix_type::value_type scalar_type
Kokkos::ArithTraits< scalar_type > ATS
Kokkos::View< DecisionType *, memory_space > results_view
typename matrix_type::local_matrix_type local_matrix_type
KOKKOS_INLINE_FUNCTION comparator_type getComparator(local_ordinal_type rlid) const
typename matrix_type::local_matrix_type local_matrix_type
typename Kokkos::DualView< const scalar_type *, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged >::t_dev diag_view_type
KOKKOS_INLINE_FUNCTION comparator_type getComparator(local_ordinal_type rlid) const
UnscaledDistanceLaplacianComparison(matrix_type &A_, DistanceFunctorType &dist2_, results_view &results_)
typename local_matrix_type::memory_space memory_space
typename local_matrix_type::value_type scalar_type
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
Kokkos::View< DecisionType *, memory_space > results_view
Comparator< local_matrix_type, DistanceFunctorType, diag_view_type > comparator_type
typename local_matrix_type::ordinal_type local_ordinal_type
Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getDiagonal(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, DistanceFunctorType &distFunctor)
KOKKOS_INLINE_FUNCTION void serialHeapSort(view_type &v, comparator_type comparator)
typename local_matrix_type2::ordinal_type local_ordinal_type
typename local_matrix_type2::memory_space memory_space
KOKKOS_INLINE_FUNCTION magnitudeType get_value(size_t x) const
KOKKOS_INLINE_FUNCTION Comparator(const local_matrix_type2 &A_, const diag_view_type2 &diag_, const local_ordinal_type rlid_, const results_view &results_)
KOKKOS_INLINE_FUNCTION bool operator()(size_t x, size_t y) const
Kokkos::View< DecisionType *, memory_space > results_view
Kokkos::ArithTraits< scalar_type > ATS
typename local_matrix_type2::value_type scalar_type
KOKKOS_INLINE_FUNCTION bool operator()(size_t x, size_t y) const
typename local_matrix_type2::ordinal_type local_ordinal_type
KOKKOS_INLINE_FUNCTION magnitudeType get_value(size_t x) const
typename local_matrix_type2::memory_space memory_space
Kokkos::View< DecisionType *, memory_space > results_view
KOKKOS_INLINE_FUNCTION Comparator(const local_matrix_type2 &A_, const diag_view_type2 &diag_, const DistanceFunctorType2 *dist2_, local_ordinal_type rlid_, const results_view &results_)
KOKKOS_INLINE_FUNCTION magnitudeType get_value(size_t x) const
Kokkos::View< DecisionType *, memory_space > results_view
KOKKOS_INLINE_FUNCTION Comparator(const local_matrix_type2 &A_, local_ordinal_type rlid_, const results_view &results_)
Kokkos::ArithTraits< scalar_type > ATS
typename local_matrix_type2::value_type scalar_type
typename local_matrix_type2::ordinal_type local_ordinal_type
typename local_matrix_type2::memory_space memory_space
KOKKOS_INLINE_FUNCTION bool operator()(size_t x, size_t y) const
KOKKOS_INLINE_FUNCTION magnitudeType get_value(size_t x) const
KOKKOS_INLINE_FUNCTION Comparator(const local_matrix_type2 &A_, const diag_view_type2 &diag_, const DistanceFunctorType2 *dist2_, local_ordinal_type rlid_, const results_view &results_)
Kokkos::View< DecisionType *, memory_space > results_view
KOKKOS_INLINE_FUNCTION bool operator()(size_t x, size_t y) const
typename local_matrix_type2::ordinal_type local_ordinal_type