10#ifndef MUELU_INVERSEAPPROXIMATIONFACTORY_DEF_HPP_
11#define MUELU_INVERSEAPPROXIMATIONFACTORY_DEF_HPP_
13#include <Xpetra_BlockedCrsMatrix.hpp>
14#include <Xpetra_CrsGraph.hpp>
15#include <Xpetra_CrsMatrixWrap.hpp>
16#include <Xpetra_CrsMatrix.hpp>
17#include <Xpetra_VectorFactory.hpp>
18#include <Xpetra_MatrixFactory.hpp>
19#include <Xpetra_Matrix.hpp>
21#include "Kokkos_Sort.hpp"
22#include "KokkosBlas1_set.hpp"
23#include "KokkosBatched_QR_Decl.hpp"
24#include "KokkosBatched_ApplyQ_Decl.hpp"
25#include "KokkosBatched_Trsv_Decl.hpp"
26#include "KokkosBatched_Util.hpp"
30#include "MueLu_Utilities.hpp"
35template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
37 RCP<ParameterList> validParamList = rcp(
new ParameterList());
38 using Magnitude =
typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
40 validParamList->set<RCP<const FactoryBase>>(
"A",
NoFactory::getRCP(),
"Matrix to build the approximate inverse on.\n");
42 validParamList->set<std::string>(
"inverse: approximation type",
"diagonal",
"Method used to approximate the inverse.");
43 validParamList->set<Magnitude>(
"inverse: drop tolerance", 0.0,
"Values below this threshold are dropped from the matrix (or fixed if diagonal fixing is active).");
44 validParamList->set<
bool>(
"inverse: fixing",
false,
"Keep diagonal and fix small entries with 1.0");
46 return validParamList;
49template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
51 Input(currentLevel,
"A");
54template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
58 using STS = Teuchos::ScalarTraits<SC>;
59 const SC one = STS::one();
60 using Magnitude =
typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
62 const ParameterList& pL = GetParameterList();
63 const bool fixing = pL.get<
bool>(
"inverse: fixing");
66 const std::string method = pL.get<std::string>(
"inverse: approximation type");
67 TEUCHOS_TEST_FOR_EXCEPTION(method !=
"diagonal" && method !=
"lumping" && method !=
"sparseapproxinverse",
Exceptions::RuntimeError,
68 "MueLu::InverseApproximationFactory::Build: Approximation type can be 'diagonal' or 'lumping' or "
69 "'sparseapproxinverse'.");
71 RCP<Matrix> A = Get<RCP<Matrix>>(currentLevel,
"A");
72 RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A);
73 const bool isBlocked = (bA == Teuchos::null ? false :
true);
76 if (isBlocked) A = bA->getMatrix(0, 0);
78 const Magnitude tol = pL.get<Magnitude>(
"inverse: drop tolerance");
79 RCP<Matrix> Ainv = Teuchos::null;
81 if (method ==
"diagonal") {
82 const auto diag = VectorFactory::Build(A->getRangeMap(),
true);
83 A->getLocalDiagCopy(*diag);
85 Ainv = MatrixFactory::Build(D);
86 }
else if (method ==
"lumping") {
89 Ainv = MatrixFactory::Build(D);
90 }
else if (method ==
"sparseapproxinverse") {
93 sparsityPattern->computeGlobalConstants();
94 GetOStream(
Statistics1) <<
"NNZ Graph(A): " << A->getCrsGraph()->getGlobalNumEntries() <<
" , NNZ Tresholded Graph(A): " << sparsityPattern->getGlobalNumEntries() << std::endl;
96 RCP<Matrix> pAinv = GetSparseInverse(A, sparsityPattern);
99 rcp_const_cast<CrsGraph>(Ainv->getCrsGraph())->computeGlobalConstants();
100 GetOStream(
Statistics1) <<
"NNZ Ainv: " << pAinv->getGlobalNumEntries() <<
", NNZ Tresholded Ainv (parameter: " << tol <<
"): " << Ainv->getGlobalNumEntries() << std::endl;
104 GetOStream(
Statistics1) <<
"Approximate inverse calculated by: " << method <<
"." << std::endl;
105 GetOStream(
Statistics1) <<
"Ainv has " << Ainv->getGlobalNumRows() <<
"x" << Ainv->getGlobalNumCols() <<
" rows and columns." << std::endl;
107 Set(currentLevel,
"Ainv", Ainv);
110template <
class local_matrix_type>
117 using impl_ATS = KokkosKernels::ArithTraits<impl_scalar_type>;
120 using shared_matrix = Kokkos::View<impl_scalar_type**, typename execution_space::scratch_memory_space, Kokkos::MemoryUnmanaged>;
121 using shared_vector = Kokkos::View<impl_scalar_type*, typename execution_space::scratch_memory_space, Kokkos::MemoryUnmanaged>;
122 using shared_lo_vector = Kokkos::View<local_ordinal_type*, typename execution_space::scratch_memory_space, Kokkos::MemoryUnmanaged>;
137 KOKKOS_INLINE_FUNCTION
138 void operator()(
const typename Kokkos::TeamPolicy<execution_space>::member_type& thread)
const {
139 auto rlid = thread.league_rank();
140 auto rowAinv =
lclAinv.row(rlid);
146 auto i = rowAinv.colidx(ii);
147 auto rowA =
lclA.rowConst(i);
149 auto j = rowA.colidx(jj);
150 column_indices(numColEntries) = j;
160 Kokkos::Experimental::sort_team(thread, Kokkos::subview(column_indices, Kokkos::make_pair(0, numColEntries)));
162 if (numColEntries > 0)
163 ++numUniqeColEntries;
166 if (column_indices(pos) != column_indices(m)) {
167 column_indices(pos + 1) = column_indices(m);
169 ++numUniqeColEntries;
170 if (column_indices(pos) == rlid)
178 KokkosBlas::SerialSet::invoke(impl_ATS::zero(), localA);
182 auto i = rowAinv.colidx(ii);
183 auto rowA =
lclA.rowConst(i);
185 auto j = rowA.colidx(jj);
186 auto v = rowA.value(jj);
190 if (column_indices(m) == j) {
201 ek(i, 0) = (i == diagOffset) ? impl_ATS::one() : impl_ATS::zero();
208 KokkosBatched::SerialQR<KokkosBatched::Algo::QR::Unblocked>::invoke(localA, tau, work);
210 KokkosBatched::SerialApplyQ<KokkosBatched::Side::Left, KokkosBatched::Trans::Transpose, KokkosBatched::Algo::ApplyQ::Unblocked>::invoke(localA, tau, ek, work);
212 auto sub_A = Kokkos::subview(localA, Kokkos::make_pair(0, rowAinv.length), Kokkos::ALL());
213 auto sub_ek = Kokkos::subview(ek, Kokkos::make_pair(0, rowAinv.length), 0);
214 KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Upper, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::NonUnit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(impl_ATS::one(), sub_A, sub_ek);
218 rowAinv.value(i) = sub_ek(i);
223template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
224RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
226 using execution_space =
typename Node::execution_space;
229 RCP<Matrix> Ainv = MatrixFactory::Build(sparsityPattern);
233 RCP<Import> rowImport = ImportFactory::Build(sparsityPattern->getRowMap(), sparsityPattern->getColMap());
234 RCP<Matrix> A = MatrixFactory::Build(Aorg, *rowImport);
236 auto maxRowEntriesAinv = Ainv->getLocalMaxNumRowEntries();
237 auto maxRowEntriesA = A->getLocalMaxNumRowEntries();
238 auto maxUniqueColEntries = maxRowEntriesAinv * maxRowEntriesA;
240 auto lclA = A->getLocalMatrixDevice();
241 auto lclAinv = Ainv->getLocalMatrixDevice();
243 Kokkos::TeamPolicy<execution_space> policy(lclAinv.numRows(), 1);
246 using shared_matrix =
typename spai_functor_type::shared_matrix;
247 using shared_vector =
typename spai_functor_type::shared_vector;
248 using shared_lo_vector =
typename spai_functor_type::shared_lo_vector;
250 int size = shared_matrix::shmem_size(maxUniqueColEntries, maxRowEntriesAinv) + shared_matrix::shmem_size(maxUniqueColEntries, 1) + shared_vector::shmem_size(3 * maxUniqueColEntries) + shared_vector::shmem_size(maxRowEntriesAinv) + shared_lo_vector::shmem_size(maxUniqueColEntries);
252 int scratchLevel = -1;
253 if (size < policy.scratch_size_max((
int)0)) {
254 policy.set_scratch_size((
int)0, Kokkos::PerTeam(size));
256 }
else if (size < policy.scratch_size_max((
int)1)) {
257 policy.set_scratch_size((
int)1, Kokkos::PerTeam(size));
260 throw Exceptions::RuntimeError(
"Neither L0 scratch memory (max size " + std::to_string(policy.scratch_size_max((
int)0)) +
261 "), nor L1 scratch memory (max size " + std::to_string(policy.scratch_size_max((
int)1)) +
262 ") is large enough for requested allocation of size " + std::to_string(size));
264 LocalSPAIFunctor spaiFunctor(lclA, lclAinv, maxUniqueColEntries, scratchLevel);
266 Kokkos::parallel_for(
"MueLu::InverseFactory::LocalSpai", policy, spaiFunctor);
269 Ainv->fillComplete();
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
void Build(Level ¤tLevel) const
Build an object with this factory.
void DeclareInput(Level ¤tLevel) const
Input.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
RCP< Matrix > GetSparseInverse(const RCP< Matrix > &A, const RCP< const CrsGraph > &sparsityPattern) const
Sparse inverse calculation method.
Class that holds all level-specific information.
KokkosKernels::ArithTraits< impl_scalar_type > impl_ATS
Kokkos::View< impl_scalar_type *, typename execution_space::scratch_memory_space, Kokkos::MemoryUnmanaged > shared_vector
Kokkos::View< impl_scalar_type **, typename execution_space::scratch_memory_space, Kokkos::MemoryUnmanaged > shared_matrix
typename local_matrix_type::ordinal_type local_ordinal_type
typename KokkosKernels::ArithTraits< scalar_type >::val_type impl_scalar_type
Kokkos::View< local_ordinal_type *, typename execution_space::scratch_memory_space, Kokkos::MemoryUnmanaged > shared_lo_vector
const local_ordinal_type maxUniqueColEntries
const local_matrix_type lclA
typename local_matrix_type::execution_space execution_space
local_matrix_type lclAinv
LocalSPAIFunctor(const local_matrix_type &lclA_, local_matrix_type &lclAinv_, local_ordinal_type maxUniqueColEntries_, int scratchLevel_)
typename local_matrix_type::value_type scalar_type
KOKKOS_INLINE_FUNCTION void operator()(const typename Kokkos::TeamPolicy< execution_space >::member_type &thread) const
static const RCP< const NoFactory > getRCP()
Static Get() functions.
static Teuchos::RCP< Vector > GetInverse(Teuchos::RCP< const Vector > v, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero())
Return vector containing inverse of input vector.
static RCP< Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > GetThresholdedGraph(const RCP< Matrix > &A, const Magnitude threshold)
Threshold a graph.
static RCP< Matrix > GetThresholdedMatrix(const RCP< Matrix > &Ain, const Magnitude threshold, const bool keepDiagonal=true)
Threshold a matrix.
static Teuchos::RCP< Vector > GetLumpedMatrixDiagonal(Matrix const &A, const bool doReciprocal=false, Magnitude tol=Teuchos::ScalarTraits< Scalar >::magnitude(Teuchos::ScalarTraits< Scalar >::zero()), Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero(), const bool replaceSingleEntryRowWithZero=false, const bool useAverageAbsDiagVal=false)
Extract Matrix Diagonal of lumped matrix.
Namespace for MueLu classes and methods.
@ Statistics1
Print more statistics.