10#ifndef MUELU_CUTDROP_HPP
11#define MUELU_CUTDROP_HPP
13#include "Kokkos_Core.hpp"
14#if KOKKOS_VERSION >= 40799
15#include "KokkosKernels_ArithTraits.hpp"
17#include "Kokkos_ArithTraits.hpp"
20#include "MueLu_Utilities.hpp"
21#include "Xpetra_Matrix.hpp"
22#include "Xpetra_MultiVector.hpp"
37template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
40 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
52#if KOKKOS_VERSION >= 40799
53 using ATS = KokkosKernels::ArithTraits<scalar_type>;
55 using ATS = Kokkos::ArithTraits<scalar_type>;
58 using values_view = Kokkos::View<magnitudeType*, memory_space>;
63 :
A(A_.getLocalMatrixDevice())
65 ,
values(
"UnscaledComparison::values",
A.nnz()) {}
67 template <
class local_matrix_type2>
75#if KOKKOS_VERSION >= 40799
76 using ATS = KokkosKernels::ArithTraits<scalar_type>;
78 using ATS = Kokkos::ArithTraits<scalar_type>;
81 using values_view = Kokkos::View<magnitudeType*, memory_space>;
83 const local_matrix_type2
A;
89 KOKKOS_INLINE_FUNCTION
92 ,
offset(A_.graph.row_map(rlid_))
94 ,
values(Kokkos::subview(values_, Kokkos::make_pair(
A.graph.row_map(rlid_),
A.graph.row_map(rlid_ + 1)))) {
95 for (
auto i = 0U; i <
values.extent(0); ++i) {
100 KOKKOS_FORCEINLINE_FUNCTION
105 KOKKOS_INLINE_FUNCTION
126 KOKKOS_INLINE_FUNCTION
128 return ATS::magnitude(
A.values(
offset + x) *
A.values(
offset + x));
134 KOKKOS_INLINE_FUNCTION
144template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node, Misc::StrengthMeasure measure>
147 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
152 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
153 using diag_view_type =
typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
160#if KOKKOS_VERSION >= 40799
161 using ATS = KokkosKernels::ArithTraits<scalar_type>;
163 using ATS = Kokkos::ArithTraits<scalar_type>;
174 :
A(A_.getLocalMatrixDevice())
176 ,
values(
"ScaledComparison::values",
A.nnz()) {
179 auto lclDiag2d =
diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
180 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
183 auto lcl2d =
diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
184 diag = Kokkos::subview(lcl2d, Kokkos::ALL(), 0);
188 template <
class local_matrix_type2,
class diag_view_type2>
196#if KOKKOS_VERSION >= 40799
197 using ATS = KokkosKernels::ArithTraits<scalar_type>;
199 using ATS = Kokkos::ArithTraits<scalar_type>;
202#if KOKKOS_VERSION >= 40799
203 using mATS = KokkosKernels::ArithTraits<magnitudeType>;
205 using mATS = Kokkos::ArithTraits<magnitudeType>;
209 const local_matrix_type2
A;
217 KOKKOS_INLINE_FUNCTION
222 ,
offset(A_.graph.row_map(rlid_))
224 ,
values(Kokkos::subview(values_, Kokkos::make_pair(
A.graph.row_map(rlid_),
A.graph.row_map(rlid_ + 1)))) {
225 for (
auto i = 0U; i <
values.extent(0); ++i) {
230 KOKKOS_INLINE_FUNCTION
235 KOKKOS_INLINE_FUNCTION
256 KOKKOS_INLINE_FUNCTION
259 auto x_aij = ATS::magnitude(
A.values(
offset + x) *
A.values(
offset + x));
261 return (x_aij / x_aiiajj);
263 auto val =
A.values(
offset + x);
264 auto neg_aij = -ATS::real(val);
265 auto max_neg_aik = ATS::real(
diag(
rlid));
266 auto v = neg_aij / max_neg_aik;
267 if (ATS::real(v) <= mATS::zero()) {
268 return -ATS::magnitude(v * v);
270 return ATS::magnitude(v * v);
273 auto val =
A.values(
offset + x);
275 const bool is_nonpositive = ATS::real(val) <= mATS::zero();
276 magnitudeType aij2 = ATS::magnitude(val) * ATS::magnitude(val);
279 return (aij2 / x_aiiajj);
286 KOKKOS_INLINE_FUNCTION
293template <Misc::StrengthMeasure measure,
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
305template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node,
class DistanceFunctorType>
308 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
313 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
314 using diag_view_type =
typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
321#if KOKKOS_VERSION >= 40799
322 using ATS = KokkosKernels::ArithTraits<scalar_type>;
324 using ATS = Kokkos::ArithTraits<scalar_type>;
336 :
A(A_.getLocalMatrixDevice())
339 ,
values(
"UnscaledDistanceLaplacianComparison::values",
A.nnz()) {
342 auto lclDiag2d =
diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
343 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
346 template <
class local_matrix_type2,
class DistanceFunctorType2,
class diag_view_type2>
354#if KOKKOS_VERSION >= 40799
355 using ATS = KokkosKernels::ArithTraits<scalar_type>;
357 using ATS = Kokkos::ArithTraits<scalar_type>;
362 const local_matrix_type2
A;
373 KOKKOS_INLINE_FUNCTION
379 ,
offset(A_.graph.row_map(rlid_))
381 ,
values(Kokkos::subview(values_, Kokkos::make_pair(
A.graph.row_map(rlid_),
A.graph.row_map(rlid_ + 1)))) {
382 for (
auto i = 0U; i <
values.extent(0); ++i) {
387 KOKKOS_FORCEINLINE_FUNCTION
392 KOKKOS_INLINE_FUNCTION
413 KOKKOS_INLINE_FUNCTION
415 auto clid =
A.graph.entries(
offset + x);
422 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val);
429 KOKKOS_INLINE_FUNCTION
439template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node,
class DistanceFunctorType, Misc::StrengthMeasure measure>
442 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
447 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
448 using diag_view_type =
typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
455#if KOKKOS_VERSION >= 40799
456 using ATS = KokkosKernels::ArithTraits<scalar_type>;
458 using ATS = Kokkos::ArithTraits<scalar_type>;
471 :
A(A_.getLocalMatrixDevice())
474 ,
values(
"ScaledDistanceLaplacianComparison::values",
A.nnz()) {
478 auto lclDiag2d =
diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
479 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
482 auto lclDiag2d =
diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
483 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
487 template <
class local_matrix_type2,
class DistanceFunctorType2,
class diag_view_type2>
495#if KOKKOS_VERSION >= 40799
496 using ATS = KokkosKernels::ArithTraits<scalar_type>;
498 using ATS = Kokkos::ArithTraits<scalar_type>;
501#if KOKKOS_VERSION >= 40799
502 using mATS = KokkosKernels::ArithTraits<magnitudeType>;
504 using mATS = Kokkos::ArithTraits<magnitudeType>;
509 const local_matrix_type2
A;
522 KOKKOS_INLINE_FUNCTION
528 ,
offset(A_.graph.row_map(rlid_))
530 ,
values(Kokkos::subview(values_, Kokkos::make_pair(
A.graph.row_map(
rlid),
A.graph.row_map(
rlid + 1)))) {
531 for (
auto i = 0U; i <
values.extent(0); ++i) {
536 KOKKOS_FORCEINLINE_FUNCTION
541 KOKKOS_INLINE_FUNCTION
562 KOKKOS_INLINE_FUNCTION
564 auto clid =
A.graph.entries(
offset + x);
573 auto aiiajj = ATS::magnitude(
diag(
rlid)) * ATS::magnitude(
diag(clid));
574 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val);
575 return (aij2 / aiiajj);
577 auto neg_aij = -ATS::real(val);
578 auto max_neg_aik = ATS::real(
diag(
rlid));
579 auto v = ATS::magnitude(neg_aij / max_neg_aik);
580 if (ATS::real(neg_aij) >=
mzero)
585 auto aiiajj = ATS::magnitude(
diag(
rlid)) * ATS::magnitude(
diag(clid));
586 const bool is_nonpositive = ATS::real(val) <= mATS::zero();
587 magnitudeType aij2 = ATS::magnitude(val) * ATS::magnitude(val);
591 return aij2 / aiiajj;
598 KOKKOS_INLINE_FUNCTION
605template <Misc::StrengthMeasure measure,
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node,
class DistanceFunctorType>
607 DistanceFunctorType& dist2_,
622template <
class comparison_type>
631#if KOKKOS_VERSION >= 40799
632 using ATS = KokkosKernels::ArithTraits<scalar_type>;
634 using ATS = Kokkos::ArithTraits<scalar_type>;
643 Kokkos::View<local_ordinal_type*, memory_space>
index;
651 index = Kokkos::View<local_ordinal_type*, memory_space>(
"indices",
A.nnz());
654 KOKKOS_FORCEINLINE_FUNCTION
656 auto row =
A.rowConst(rlid);
657 size_t nnz = row.length;
659 auto drop_view = Kokkos::subview(
results, Kokkos::make_pair(
A.graph.row_map(rlid),
A.graph.row_map(rlid + 1)));
660 auto row_permutation = Kokkos::subview(
index, Kokkos::make_pair(
A.graph.row_map(rlid),
A.graph.row_map(rlid + 1)));
662 auto comparator =
comparison.getComparator(rlid);
664#ifdef MUELU_COALESCE_DROP_DEBUG
666 Kokkos::printf(
"SoC: ");
668 Kokkos::printf(
"%5f ", comparator.get_value(k));
670 Kokkos::printf(
"\n");
674 for (
size_t i = 0; i < nnz; ++i) {
675 row_permutation(i) = i;
679 size_t keepStart = 0;
680 size_t dropStart = nnz;
682 for (
size_t i = 1; i < nnz; ++i) {
683 auto const& x = row_permutation(i - 1);
684 auto const& y = row_permutation(i);
691 if (
eps *
eps * x_aij > y_aij) {
699 for (
size_t i = keepStart; i < nnz; ++i) {
700 drop_view(row_permutation(i)) = Kokkos::max(dropStart <= i ?
DROP :
KEEP, drop_view(row_permutation(i)));
Order each row by a criterion, compare the ratio of values and drop all entries once the ratio is bel...
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
comparison_type comparison
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
Orders entries of row by .
ScaledComparison(matrix_type &A_, results_view &results_)
typename ATS::magnitudeType magnitudeType
typename Kokkos::DualView< const scalar_type *, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged >::t_dev diag_view_type
Kokkos::View< DecisionType *, memory_space > results_view
Comparator< local_matrix_type, diag_view_type > comparator_type
Kokkos::View< magnitudeType *, memory_space > values_view
typename matrix_type::local_matrix_device_type local_matrix_type
Kokkos::ArithTraits< scalar_type > ATS
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
typename local_matrix_type::memory_space memory_space
Teuchos::RCP< diag_vec_type > diagVec
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
typename local_matrix_type::value_type scalar_type
KOKKOS_INLINE_FUNCTION comparator_type getComparator(local_ordinal_type rlid) const
typename local_matrix_type::ordinal_type local_ordinal_type
Orders entries of row by where is the distance Laplacian.
KOKKOS_INLINE_FUNCTION comparator_type getComparator(local_ordinal_type rlid) const
typename local_matrix_type::ordinal_type local_ordinal_type
DistanceFunctorType dist2
Kokkos::ArithTraits< scalar_type > ATS
typename matrix_type::local_matrix_device_type local_matrix_type
typename ATS::magnitudeType magnitudeType
typename Kokkos::DualView< const scalar_type *, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged >::t_dev diag_view_type
Kokkos::View< DecisionType *, memory_space > results_view
ScaledDistanceLaplacianComparison(matrix_type &A_, DistanceFunctorType &dist2_, results_view &results_)
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
Comparator< local_matrix_type, DistanceFunctorType, diag_view_type > comparator_type
Kokkos::View< magnitudeType *, memory_space > values_view
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
typename local_matrix_type::value_type scalar_type
typename local_matrix_type::memory_space memory_space
Teuchos::RCP< diag_vec_type > diagVec
Orders entries of row by .
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< magnitudeType *, memory_space > values_view
Kokkos::View< DecisionType *, memory_space > results_view
KOKKOS_INLINE_FUNCTION comparator_type getComparator(local_ordinal_type rlid) const
typename matrix_type::local_matrix_device_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
Kokkos::ArithTraits< scalar_type > ATS
UnscaledDistanceLaplacianComparison(matrix_type &A_, DistanceFunctorType &dist2_, results_view &results_)
typename ATS::magnitudeType magnitudeType
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
Teuchos::RCP< diag_vec_type > diagVec
DistanceFunctorType dist2
Kokkos::View< DecisionType *, memory_space > results_view
Comparator< local_matrix_type, DistanceFunctorType, diag_view_type > comparator_type
typename matrix_type::local_matrix_device_type local_matrix_type
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< magnitudeType *, memory_space > values_view
static Teuchos::RCP< Vector > GetMatrixMaxMinusOffDiagonal(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Return vector containing: max_{i\not=k}(-a_ik), for each for i in the matrix.
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)
Extract Overlapped Matrix Diagonal.
auto make_comparison_functor(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A_, typename ScaledComparison< Scalar, LocalOrdinal, GlobalOrdinal, Node, measure >::results_view &results_)
auto make_dlap_comparison_functor(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A_, DistanceFunctorType &dist2_, typename ScaledDistanceLaplacianComparison< Scalar, LocalOrdinal, GlobalOrdinal, Node, DistanceFunctorType, measure >::results_view &results_)
Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getDiagonal(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, DistanceFunctorType &distFunctor)
Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getMaxMinusOffDiagonal(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, DistanceFunctorType &distFunctor)
KOKKOS_INLINE_FUNCTION void serialHeapSort(view_type &v, comparator_type comparator)
@ SignedSmoothedAggregationMeasure
@ SignedRugeStuebenMeasure
@ SmoothedAggregationMeasure
KOKKOS_INLINE_FUNCTION magnitudeType get_value(size_t x) const
const local_ordinal_type offset
const results_view results
Kokkos::ArithTraits< magnitudeType > mATS
KOKKOS_INLINE_FUNCTION Comparator(const local_matrix_type2 &A_, const diag_view_type2 &diag_, const local_ordinal_type rlid_, const results_view &results_, values_view &values_)
typename ATS::magnitudeType magnitudeType
Kokkos::View< DecisionType *, memory_space > results_view
Kokkos::View< magnitudeType *, memory_space > values_view
typename local_matrix_type2::ordinal_type local_ordinal_type
const local_ordinal_type rlid
Kokkos::ArithTraits< scalar_type > ATS
const diag_view_type2 diag
KOKKOS_INLINE_FUNCTION bool operator()(size_t x, size_t y) const
typename local_matrix_type2::memory_space memory_space
const local_matrix_type2 A
KOKKOS_INLINE_FUNCTION magnitudeType get_value_impl(size_t x) const
typename local_matrix_type2::value_type scalar_type
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_, values_view &values_)
KOKKOS_FORCEINLINE_FUNCTION magnitudeType get_value(size_t x) const
const DistanceFunctorType2 * dist2
Kokkos::ArithTraits< magnitudeType > mATS
KOKKOS_INLINE_FUNCTION bool operator()(size_t x, size_t y) const
const local_ordinal_type rlid
typename local_matrix_type2::memory_space memory_space
const local_matrix_type2 A
const results_view results
Kokkos::ArithTraits< scalar_type > ATS
Kokkos::View< DecisionType *, memory_space > results_view
typename local_matrix_type2::ordinal_type local_ordinal_type
typename ATS::magnitudeType magnitudeType
const diag_view_type2 diag
const magnitudeType mzero
Kokkos::View< magnitudeType *, memory_space > values_view
const local_ordinal_type offset
KOKKOS_INLINE_FUNCTION magnitudeType get_value_impl(size_t x) const
typename local_matrix_type2::value_type scalar_type
KOKKOS_INLINE_FUNCTION Comparator(const local_matrix_type2 &A_, local_ordinal_type rlid_, const results_view &results_, values_view &values_)
Kokkos::View< DecisionType *, memory_space > results_view
const local_ordinal_type offset
Kokkos::ArithTraits< scalar_type > ATS
typename local_matrix_type2::value_type scalar_type
typename local_matrix_type2::ordinal_type local_ordinal_type
KOKKOS_FORCEINLINE_FUNCTION magnitudeType get_value(size_t x) const
typename local_matrix_type2::memory_space memory_space
typename ATS::magnitudeType magnitudeType
Kokkos::View< magnitudeType *, memory_space > values_view
const local_matrix_type2 A
KOKKOS_INLINE_FUNCTION magnitudeType get_value_impl(size_t x) const
const results_view results
KOKKOS_INLINE_FUNCTION bool operator()(size_t x, size_t y) const
const local_ordinal_type offset
typename ATS::magnitudeType magnitudeType
KOKKOS_INLINE_FUNCTION magnitudeType get_value_impl(size_t x) const
KOKKOS_FORCEINLINE_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_, values_view &values_)
Kokkos::View< magnitudeType *, memory_space > values_view
const diag_view_type2 diag
Kokkos::View< DecisionType *, memory_space > results_view
const local_matrix_type2 A
Kokkos::ArithTraits< scalar_type > ATS
const results_view results
const DistanceFunctorType2 * dist2
KOKKOS_INLINE_FUNCTION bool operator()(size_t x, size_t y) const
const local_ordinal_type rlid
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