10#ifndef MUELU_CUTDROP_HPP
11#define MUELU_CUTDROP_HPP
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"
33template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
36 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
48 using ATS = Kokkos::ArithTraits<scalar_type>;
53 :
A(A_.getLocalMatrixDevice())
56 template <
class local_matrix_type2>
64 using ATS = Kokkos::ArithTraits<scalar_type>;
67 const local_matrix_type2
A;
72 KOKKOS_INLINE_FUNCTION
75 ,
offset(A_.graph.row_map(rlid_))
78 KOKKOS_INLINE_FUNCTION
80 return ATS::magnitude(
A.values(
offset + x) *
A.values(
offset + x));
83 KOKKOS_INLINE_FUNCTION
106 KOKKOS_INLINE_FUNCTION
116template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
119 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
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;
132 using ATS = Kokkos::ArithTraits<scalar_type>;
140 :
A(A_.getLocalMatrixDevice())
143 auto lclDiag2d =
diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
144 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
147 template <
class local_matrix_type2,
class diag_view_type2>
155 using ATS = Kokkos::ArithTraits<scalar_type>;
158 const local_matrix_type2
A;
165 KOKKOS_INLINE_FUNCTION
170 ,
offset(A_.graph.row_map(rlid_))
173 KOKKOS_INLINE_FUNCTION
175 auto x_aij = ATS::magnitude(
A.values(
offset + x) *
A.values(
offset + x));
177 return (x_aij / x_aiiajj);
180 KOKKOS_INLINE_FUNCTION
203 KOKKOS_INLINE_FUNCTION
209template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node,
class DistanceFunctorType>
212 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
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;
225 using ATS = Kokkos::ArithTraits<scalar_type>;
234 :
A(A_.getLocalMatrixDevice())
239 auto lclDiag2d =
diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
240 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
243 template <
class local_matrix_type2,
class DistanceFunctorType2,
class diag_view_type2>
251 using ATS = Kokkos::ArithTraits<scalar_type>;
254 const local_matrix_type2
A;
264 KOKKOS_INLINE_FUNCTION
270 ,
offset(A_.graph.row_map(rlid_))
273 KOKKOS_INLINE_FUNCTION
275 auto clid =
A.graph.entries(
offset + x);
282 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val);
286 KOKKOS_INLINE_FUNCTION
309 KOKKOS_INLINE_FUNCTION
319template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node,
class DistanceFunctorType>
322 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
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;
335 using ATS = Kokkos::ArithTraits<scalar_type>;
344 :
A(A_.getLocalMatrixDevice())
349 auto lclDiag2d =
diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
350 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
353 template <
class local_matrix_type2,
class DistanceFunctorType2,
class diag_view_type2>
361 using ATS = Kokkos::ArithTraits<scalar_type>;
364 const local_matrix_type2
A;
374 KOKKOS_INLINE_FUNCTION
380 ,
offset(A_.graph.row_map(rlid_))
383 KOKKOS_INLINE_FUNCTION
385 auto clid =
A.graph.entries(
offset + x);
392 auto aiiajj = ATS::magnitude(
diag(
rlid)) * ATS::magnitude(
diag(clid));
393 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val);
394 return (aij2 / aiiajj);
397 KOKKOS_INLINE_FUNCTION
420 KOKKOS_INLINE_FUNCTION
430template <
class comparison_type>
439 using ATS = Kokkos::ArithTraits<scalar_type>;
447 Kokkos::View<local_ordinal_type*, memory_space>
index;
455 index = Kokkos::View<local_ordinal_type*, memory_space>(
"indices",
A.nnz());
458 KOKKOS_FORCEINLINE_FUNCTION
460 auto row =
A.rowConst(rlid);
461 size_t nnz = row.length;
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)));
466 auto comparator =
comparison.getComparator(rlid);
468 for (
size_t i = 0; i < nnz; ++i) {
469 row_permutation(i) = i;
473 size_t keepStart = 0;
474 size_t dropStart = nnz;
476 for (
size_t i = 1; i < nnz; ++i) {
477 auto const& x = row_permutation(i - 1);
478 auto const& y = row_permutation(i);
485 if (
eps *
eps * x_aij > y_aij) {
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)));
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
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
DistanceFunctorType dist2
Teuchos::RCP< diag_vec_type > diagVec
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 ATS::magnitudeType magnitudeType
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
Kokkos::ArithTraits< scalar_type > ATS
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
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 local_matrix_type::ordinal_type local_ordinal_type
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)
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
const results_view results
typename local_matrix_type2::memory_space memory_space
KOKKOS_INLINE_FUNCTION magnitudeType get_value(size_t x) const
const local_ordinal_type offset
const local_ordinal_type rlid
KOKKOS_INLINE_FUNCTION Comparator(const local_matrix_type2 &A_, const diag_view_type2 &diag_, const local_ordinal_type rlid_, const results_view &results_)
const local_matrix_type A
const diag_view_type diag
KOKKOS_INLINE_FUNCTION bool operator()(size_t x, size_t y) const
typename ATS::magnitudeType magnitudeType
Kokkos::View< DecisionType *, memory_space > results_view
Kokkos::ArithTraits< scalar_type > ATS
typename local_matrix_type2::value_type scalar_type
const DistanceFunctorType * dist2
typename local_matrix_type2::value_type scalar_type
const local_ordinal_type rlid
KOKKOS_INLINE_FUNCTION bool operator()(size_t x, size_t y) const
const diag_view_type diag
typename ATS::magnitudeType magnitudeType
const local_matrix_type A
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
const local_ordinal_type offset
Kokkos::ArithTraits< scalar_type > ATS
const results_view results
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
const local_ordinal_type offset
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
typename ATS::magnitudeType magnitudeType
const local_matrix_type A
const results_view results
KOKKOS_INLINE_FUNCTION bool operator()(size_t x, size_t y) const
const local_ordinal_type offset
KOKKOS_INLINE_FUNCTION magnitudeType get_value(size_t x) const
typename ATS::magnitudeType magnitudeType
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_)
const diag_view_type diag
Kokkos::View< DecisionType *, memory_space > results_view
const local_matrix_type A
Kokkos::ArithTraits< scalar_type > ATS
const results_view results
const DistanceFunctorType * 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