MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_InverseApproximationFactory_def.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_INVERSEAPPROXIMATIONFACTORY_DEF_HPP_
11 #define MUELU_INVERSEAPPROXIMATIONFACTORY_DEF_HPP_
12 
14 #include <Xpetra_CrsGraph.hpp>
16 #include <Xpetra_CrsMatrixWrap.hpp>
17 #include <Xpetra_CrsMatrix.hpp>
18 #include <Xpetra_MultiVectorFactory.hpp>
19 #include <Xpetra_VectorFactory.hpp>
20 #include <Xpetra_MatrixFactory.hpp>
21 #include <Xpetra_Matrix.hpp>
22 #include <Xpetra_MatrixMatrix.hpp>
24 
28 
29 #include "MueLu_Level.hpp"
30 #include "MueLu_Monitor.hpp"
31 #include "MueLu_Utilities.hpp"
33 
34 namespace MueLu {
35 
36 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
38  RCP<ParameterList> validParamList = rcp(new ParameterList());
39  using Magnitude = typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
40 
41  validParamList->set<RCP<const FactoryBase>>("A", NoFactory::getRCP(), "Matrix to build the approximate inverse on.\n");
42 
43  validParamList->set<std::string>("inverse: approximation type", "diagonal", "Method used to approximate the inverse.");
44  validParamList->set<Magnitude>("inverse: drop tolerance", 0.0, "Values below this threshold are dropped from the matrix (or fixed if diagonal fixing is active).");
45  validParamList->set<bool>("inverse: fixing", false, "Keep diagonal and fix small entries with 1.0");
46 
47  return validParamList;
48 }
49 
50 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
52  Input(currentLevel, "A");
53 }
54 
55 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
57  FactoryMonitor m(*this, "Build", currentLevel);
58 
59  using STS = Teuchos::ScalarTraits<SC>;
60  const SC one = STS::one();
61  using Magnitude = typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
62 
63  const ParameterList& pL = GetParameterList();
64  const bool fixing = pL.get<bool>("inverse: fixing");
65 
66  // check which approximation type to use
67  const std::string method = pL.get<std::string>("inverse: approximation type");
68  TEUCHOS_TEST_FOR_EXCEPTION(method != "diagonal" && method != "lumping" && method != "sparseapproxinverse", Exceptions::RuntimeError,
69  "MueLu::InverseApproximationFactory::Build: Approximation type can be 'diagonal' or 'lumping' or "
70  "'sparseapproxinverse'.");
71 
72  RCP<Matrix> A = Get<RCP<Matrix>>(currentLevel, "A");
73  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A);
74  const bool isBlocked = (bA == Teuchos::null ? false : true);
75 
76  // if blocked operator is used, defaults to A(0,0)
77  if (isBlocked) A = bA->getMatrix(0, 0);
78 
79  const Magnitude tol = pL.get<Magnitude>("inverse: drop tolerance");
80  RCP<Matrix> Ainv = Teuchos::null;
81 
82  if (method == "diagonal") {
83  const auto diag = VectorFactory::Build(A->getRangeMap(), true);
84  A->getLocalDiagCopy(*diag);
85  const RCP<const Vector> D = (!fixing ? Utilities::GetInverse(diag) : Utilities::GetInverse(diag, tol, one));
86  Ainv = MatrixFactory::Build(D);
87  } else if (method == "lumping") {
88  const auto diag = Utilities::GetLumpedMatrixDiagonal(*A);
89  const RCP<const Vector> D = (!fixing ? Utilities::GetInverse(diag) : Utilities::GetInverse(diag, tol, one));
90  Ainv = MatrixFactory::Build(D);
91  } else if (method == "sparseapproxinverse") {
92  RCP<CrsGraph> sparsityPattern = Utilities::GetThresholdedGraph(A, tol, A->getGlobalMaxNumRowEntries());
93  GetOStream(Statistics1) << "NNZ Graph(A): " << A->getCrsGraph()->getGlobalNumEntries() << " , NNZ Tresholded Graph(A): " << sparsityPattern->getGlobalNumEntries() << std::endl;
94  RCP<Matrix> pAinv = GetSparseInverse(A, sparsityPattern);
95  Ainv = Utilities::GetThresholdedMatrix(pAinv, tol, fixing, pAinv->getGlobalMaxNumRowEntries());
96  GetOStream(Statistics1) << "NNZ Ainv: " << pAinv->getGlobalNumEntries() << ", NNZ Tresholded Ainv (parameter: " << tol << "): " << Ainv->getGlobalNumEntries() << std::endl;
97  }
98 
99  GetOStream(Statistics1) << "Approximate inverse calculated by: " << method << "." << std::endl;
100  GetOStream(Statistics1) << "Ainv has " << Ainv->getGlobalNumRows() << "x" << Ainv->getGlobalNumCols() << " rows and columns." << std::endl;
101 
102  Set(currentLevel, "Ainv", Ainv);
103 }
104 
105 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
108  // construct the inverse matrix with the given sparsity pattern
109  RCP<Matrix> Ainv = MatrixFactory::Build(sparsityPattern);
110  Ainv->resumeFill();
111 
112  // gather missing rows from other procs to generate an overlapping map
113  RCP<Import> rowImport = ImportFactory::Build(sparsityPattern->getRowMap(), sparsityPattern->getColMap());
114  RCP<Matrix> A = MatrixFactory::Build(Aorg, *rowImport);
115 
116  // loop over all rows of the inverse sparsity pattern (this can be done in parallel)
117  for (size_t k = 0; k < sparsityPattern->getLocalNumRows(); k++) {
118  // 1. get column indices Ik of local row k
120  sparsityPattern->getLocalRowView(k, Ik);
121 
122  // 2. get all local A(Ik,:) rows
125  Array<LO> Jk;
126  for (LO i = 0; i < Ik.size(); i++) {
127  A->getLocalRowView(Ik[i], J[i], Ak[i]);
128  for (LO j = 0; j < J[i].size(); j++)
129  Jk.append(J[i][j]);
130  }
131  // set of unique column indices Jk
132  std::sort(Jk.begin(), Jk.end());
133  Jk.erase(std::unique(Jk.begin(), Jk.end()), Jk.end());
134  // create map
135  std::map<LO, LO> G;
136  for (LO i = 0; i < Jk.size(); i++) G.insert(std::pair<LO, LO>(Jk[i], i));
137 
138  // 3. merge rows together
139  Teuchos::SerialDenseMatrix<LO, SC> localA(Jk.size(), Ik.size(), true);
140  for (LO i = 0; i < Ik.size(); i++) {
141  for (LO j = 0; j < J[i].size(); j++) {
142  localA(G.at(J[i][j]), i) = Ak[i][j];
143  }
144  }
145 
146  // 4. get direction-vector
147  // diagonal needs an entry!
149  ek[std::find(Jk.begin(), Jk.end(), k) - Jk.begin()] = Teuchos::ScalarTraits<Scalar>::one();
150  ;
151 
152  // 5. solve linear system for x
155  qrSolver.setMatrix(Teuchos::rcp(&localA, false));
156  qrSolver.setVectors(Teuchos::rcp(&localX, false), Teuchos::rcp(&ek, false));
157  const int err = qrSolver.solve();
159  "MueLu::InverseApproximationFactory::GetSparseInverse: Error in serial QR solve.");
160 
161  // 6. set calculated row into Ainv
162  ArrayView<const SC> Mk(localX.values(), localX.length());
163  Ainv->replaceLocalValues(k, Ik, Mk);
164  }
165  Ainv->fillComplete();
166 
167  return Ainv;
168 }
169 
170 } // namespace MueLu
171 
172 #endif /* MUELU_INVERSEAPPROXIMATIONFACTORY_DEF_HPP_ */
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Array< T > & append(const T &x)
T & get(const std::string &name, T def_value)
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Print more statistics.
static RCP< Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > GetThresholdedGraph(const RCP< Matrix > &A, const Magnitude threshold, const GlobalOrdinal expectedNNZperRow=-1)
Threshold a graph.
size_type size() const
LocalOrdinal LO
iterator erase(iterator position)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
void Build(Level &currentLevel) const
Build an object with this factory.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static RCP< CrsMatrixWrap > GetThresholdedMatrix(const RCP< Matrix > &Ain, const Magnitude threshold, const bool keepDiagonal=true, const GlobalOrdinal expectedNNZperRow=-1)
Threshold a matrix.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
RCP< Matrix > GetSparseInverse(const RCP< Matrix > &A, const RCP< const CrsGraph > &sparsityPattern) const
Sparse inverse calculation method.
void DeclareInput(Level &currentLevel) const
Input.
static const RCP< const NoFactory > getRCP()
Static Get() functions.
iterator end()
size_type size() const
Scalar SC
Exception throws to report errors in the internal logical of the program.
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 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.
iterator begin()