10#ifndef MUELU_DROPPINGCOMMON_HPP
11#define MUELU_DROPPINGCOMMON_HPP
15#include "Kokkos_Core.hpp"
16#include "KokkosKernels_ArithTraits.hpp"
17#include "Tpetra_Access.hpp"
18#include "Xpetra_Matrix.hpp"
19#include "Xpetra_VectorFactory.hpp"
20#include "MueLu_Utilities.hpp"
37template <
class local_ordinal_type>
42 KOKKOS_FORCEINLINE_FUNCTION
51template <
class local_matrix_type>
70 KOKKOS_FORCEINLINE_FUNCTION
72 auto row =
A.rowConst(rlid);
73 const size_t offset =
A.graph.row_map(rlid);
77 auto clid = row.colidx(k);
89template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
92 using local_matrix_type =
typename Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_device_type;
108 A = A_.getLocalMatrixDevice();
109 auto boundaryNodesColumn =
typename boundary_nodes_view::non_const_type(
"boundaryNodesColumn", A_.getColMap()->getLocalNumElements());
110 auto boundaryNodesDomain =
typename boundary_nodes_view::non_const_type(
"boundaryNodesDomain", A_.getDomainMap()->getLocalNumElements());
115 KOKKOS_FORCEINLINE_FUNCTION
117 auto row =
A.rowConst(rlid);
118 const size_t offset =
A.graph.row_map(rlid);
121 auto clid = row.colidx(k);
133template <
class local_matrix_type>
155 KOKKOS_FORCEINLINE_FUNCTION
157 auto row =
A.rowConst(rlid);
158 const size_t offset =
A.graph.row_map(rlid);
162 auto clid = row.colidx(k);
174template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
177 using local_matrix_type =
typename Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_device_type;
198 A = A_.getLocalMatrixDevice();
199 auto boundaryNodesColumn =
typename boundary_nodes_view::non_const_type(
"boundaryNodesColumn", A_.getColMap()->getLocalNumElements());
200 auto boundaryNodesDomain =
typename boundary_nodes_view::non_const_type(
"boundaryNodesDomain", A_.getDomainMap()->getLocalNumElements());
205 KOKKOS_FORCEINLINE_FUNCTION
207 auto row =
A.rowConst(rlid);
208 const size_t offset =
A.graph.row_map(rlid);
211 auto clid = row.colidx(k);
224template <
class local_matrix_type>
240 KOKKOS_FORCEINLINE_FUNCTION
242 auto row =
A.rowConst(rlid);
243 const size_t offset =
A.graph.row_map(rlid);
245 auto clid = row.colidx(k);
258template <
class local_matrix_type>
274 KOKKOS_FORCEINLINE_FUNCTION
276 auto row =
A.rowConst(rlid);
277 const size_t offset =
A.graph.row_map(rlid);
279 auto clid = row.colidx(k);
280 if (clid >=
A.numRows()) {
291template <
class local_matrix_type>
311 KOKKOS_FORCEINLINE_FUNCTION
313 auto row =
A.rowConst(rlid);
314 const size_t offset =
A.graph.row_map(rlid);
316 auto clid = row.colidx(k);
317 if ((
results(offset + k) ==
KEEP) && (rlid != clid))
322 auto clid = row.colidx(k);
335template <
class local_matrix_type>
358 KOKKOS_FORCEINLINE_FUNCTION
360 auto row =
A.rowConst(rlid);
361 const size_t offset =
A.graph.row_map(rlid);
363 auto clid = row.colidx(k);
364 if ((
results(offset + k) ==
KEEP) && (rlid != clid))
370 auto clid = row.colidx(k);
383template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
386 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
405 :
A(A_.getLocalMatrixDevice())
406 ,
point_to_block(point_to_block_.getLocalViewDevice(Tpetra::Access::ReadOnly))
408 auto importer = A_.getCrsGraph()->getImporter();
410 if (!importer.is_null()) {
411 ghosted_point_to_blockMV = Xpetra::VectorFactory<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(),
false);
418 KOKKOS_FORCEINLINE_FUNCTION
420 auto row =
A.rowConst(rlid);
421 const size_t offset =
A.graph.row_map(rlid);
423 auto clid = row.colidx(k);
437template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
440 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
463 :
A(A_.getLocalMatrixDevice())
464 ,
point_to_block(point_to_block_.getLocalViewDevice(Tpetra::Access::ReadOnly))
468 if (!importer.is_null()) {
469 ghosted_point_to_blockMV = Xpetra::VectorFactory<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>::Build(importer->getTargetMap(),
false);
476 KOKKOS_FORCEINLINE_FUNCTION
478 auto row =
A.rowConst(rlid);
479 const size_t offset =
A.graph.row_map(rlid);
482 auto clid = row.colidx(k);
497template <
class local_matrix_type>
515 KOKKOS_FORCEINLINE_FUNCTION
517 auto row =
A.rowConst(rlid);
518 const size_t offset =
A.graph.row_map(rlid);
521 Kokkos::printf(
"No dropping decision was taken for entry (%d, %d)\n", rlid, row.colidx(k));
532template <
class local_matrix_type>
548 KOKKOS_FORCEINLINE_FUNCTION
550 auto row =
A.rowConst(rlid);
551 const size_t offset =
A.graph.row_map(rlid);
554 auto clid = row.colidx(k);
555 if (clid >=
A.numRows())
557 auto row2 =
A.rowConst(clid);
558 const size_t offset2 =
A.graph.row_map(clid);
560 auto clid2 = row2.colidx(k2);
572template <
class view_type,
class comparator_type>
573KOKKOS_INLINE_FUNCTION
void serialHeapSort(view_type& v, comparator_type comparator) {
574 auto N = v.extent(0);
575 size_t start = N / 2;
587 while (2 * root + 1 < end) {
588 size_t child = 2 * root + 1;
589 if ((child + 1 < end) and (comparator(v(child), v(child + 1))))
592 if (comparator(v(root), v(child))) {
638#define MUELU_ETI_SLGN_SoC(CLASSNAME, SC, LO, GO, NO) \
639 template class CLASSNAME<SC, LO, GO, NO, MueLu::Misc::SmoothedAggregationMeasure>; \
640 template class CLASSNAME<SC, LO, GO, NO, MueLu::Misc::SignedRugeStuebenMeasure>; \
641 template class CLASSNAME<SC, LO, GO, NO, MueLu::Misc::SignedSmoothedAggregationMeasure>; \
642 template class CLASSNAME<SC, LO, GO, NO, MueLu::Misc::UnscaledMeasure>;
Functor that drops all entries that are not on the block diagonal.
Teuchos::RCP< block_indices_type > ghosted_point_to_blockMV
typename matrix_type::local_matrix_device_type local_matrix_type
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
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_, results_view &results_)
Functor that drops all entries that are not on the block diagonal.
Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > block_indices_type
id_translation_type col_translation
typename local_matrix_type::memory_space memory_space
Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > importer_type
typename block_indices_type::dual_view_type_const::t_dev local_block_indices_view_type
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
BlockDiagonalizeVectorFunctor(matrix_type &A_, block_indices_type &point_to_block_, const RCP< const importer_type > &importer, results_view &results_, id_translation_type row_translation_, id_translation_type col_translation_)
id_translation_type row_translation
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
local_block_indices_view_type point_to_block
typename local_matrix_type::ordinal_type local_ordinal_type
Teuchos::RCP< block_indices_type > ghosted_point_to_blockMV
Kokkos::View< local_ordinal_type *, memory_space > id_translation_type
typename matrix_type::local_matrix_device_type local_matrix_type
local_block_indices_view_type ghosted_point_to_block
Functor that checks that all entries have been marked.
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
Functor that drops off-rank entries.
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_)
Functor that marks diagonal as kept, unless the are already marked as boundary.
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
Functor that marks singletons (all off-diagonal entries in a row are dropped) as boundary.
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
boundary_nodes_view boundaryNodes
typename local_matrix_type::ordinal_type local_ordinal_type
Kokkos::View< bool *, memory_space > boundary_nodes_view
Functor that marks singletons (all off-diagonal entries in a row are dropped) as boundary.
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
block_indices_view_type point_to_block
boundary_nodes_view boundaryNodes
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Functor that drops boundary nodes for a blockSize == 1 problem.
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
boundary_nodes_view boundaryNodes
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
Functor that drops boundary nodes for a blockSize == 1 problem.
Kokkos::View< DecisionType *, memory_space > results_view
boundary_nodes_view boundaryNodes
typename Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_device_type local_matrix_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
typename local_matrix_type::value_type scalar_type
boundary_nodes_view boundaryNodesCol
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
PointwiseSymmetricDropBoundaryFunctor(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A_, boundary_nodes_view boundaryNodes_, results_view &results_)
Functor that symmetrizes the dropping decisions.
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
Functor that drops boundary nodes for a blockSize > 1 problem.
VectorDropBoundaryFunctor(local_matrix_type &A_, block_indices_view_type point_to_block_, boundary_nodes_view boundaryNodes_, results_view &results_)
block_indices_view_type point_to_block
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
boundary_nodes_view boundaryNodes
typename local_matrix_type::value_type scalar_type
Functor that drops boundary nodes for a blockSize > 1 problem.
boundary_nodes_view boundaryNodes
block_indices_view_type point_to_block
typename local_matrix_type::value_type scalar_type
typename local_matrix_type::memory_space memory_space
typename Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_device_type local_matrix_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(local_ordinal_type rlid) const
Kokkos::View< bool *, memory_space > boundary_nodes_view
Kokkos::View< DecisionType *, memory_space > results_view
boundary_nodes_view boundaryNodesCol
VectorSymmetricDropBoundaryFunctor(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A_, block_indices_view_type point_to_block_, block_indices_view_type ghosted_point_to_block_, boundary_nodes_view boundaryNodes_, results_view &results_)
Kokkos::View< local_ordinal_type *, memory_space > block_indices_view_type
typename local_matrix_type::ordinal_type local_ordinal_type
block_indices_view_type ghosted_point_to_block
static void DetectDirichletColsAndDomains(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayRCP< bool > &dirichletRows, Teuchos::ArrayRCP< bool > dirichletCols, Teuchos::ArrayRCP< bool > dirichletDomain)
Detects Dirichlet columns & domains from a list of Dirichlet rows.
KOKKOS_INLINE_FUNCTION void serialHeapSort(view_type &v, comparator_type comparator)
@ SignedSmoothedAggregationMeasure
@ SignedRugeStuebenMeasure
@ SmoothedAggregationMeasure
Namespace for MueLu classes and methods.