MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_GMRESSolver_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_GMRESSOLVER_DEF_HPP
11#define MUELU_GMRESSOLVER_DEF_HPP
12
13#include <Teuchos_LAPACK.hpp>
14
15#include <Xpetra_MatrixFactory.hpp>
16#include <Xpetra_MatrixFactory2.hpp>
17#include <Xpetra_MatrixMatrix.hpp>
18#include <Xpetra_IO.hpp>
19
21
22#include "MueLu_Constraint.hpp"
23#include "MueLu_Monitor.hpp"
24#include "MueLu_Utilities.hpp"
25
26namespace MueLu {
27
28template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
31
32template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
34 for (int i = 0; i < k; i++) {
35 SC w1 = c[i] * v[i] - s[i] * v[i + 1];
36 SC w2 = s[i] * v[i] + c[i] * v[i + 1];
37 v[i] = w1;
38 v[i + 1] = w2;
39 }
40}
41
42template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
43void GMRESSolver<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Iterate(const Matrix& Aref, const Constraint& C, const Matrix& P0, RCP<Matrix>& finalP) const {
44 PrintMonitor m(*this, "GMRES iterations");
45
46 finalP = MatrixFactory2::BuildCopy(rcpFromRef(P0));
47 if (nIts_ == 0)
48 return;
49
50 TEUCHOS_TEST_FOR_EXCEPTION(nIts_ > 1, Exceptions::RuntimeError,
51 "For now, solving Hessenberg system works only for a single iteration");
52
53 SC one = Teuchos::ScalarTraits<SC>::one(), zero = Teuchos::ScalarTraits<SC>::zero();
54
55 RCP<const Matrix> A = rcpFromRef(Aref);
56 // bool useTpetra = (A->getRowMap()->lib() == Xpetra::UseTpetra);
57
58 // FIXME: Don't know why, but in the MATLAB code we have D = I. Follow that for now.
59#if 0
60 ArrayRCP<const SC> D = Utilities::GetMatrixDiagonal_arcp(*A);
61#else
62 ArrayRCP<const SC> D(A->getLocalNumRows(), one);
63#endif
64
65 Teuchos::FancyOStream& mmfancy = this->GetOStream(Statistics2);
66
67 // Initial P0 would only be used for multiplication
68 RCP<Matrix> X = rcp_const_cast<Matrix>(rcpFromRef(P0)), tmpAP, newV;
69 std::vector<RCP<Matrix> > V(nIts_ + 1);
70
71 // T is used only for projecting onto
72 RCP<CrsMatrix> T_ = CrsMatrixFactory::Build(C.GetPattern());
73 T_->fillComplete(P0.getDomainMap(), P0.getRangeMap());
74 RCP<Matrix> T = rcp(new CrsMatrixWrap(T_));
75
76 SC rho;
77 {
78 // R_0 = -D^{-1}*A*X_0
79 // V_0 = R_0 / ||R_0||_F
80 tmpAP = MatrixMatrix::Multiply(*A, false, *X, false, mmfancy, true /*doFillComplete*/, true /*optimizeStorage*/);
81 C.Apply(*tmpAP, *T);
82
83 V[0] = MatrixFactory2::BuildCopy(T);
84 Utilities::MyOldScaleMatrix(*V[0], D, true /*doInverse*/, true /*doFillComplete*/, false /*doOptimizeStorage*/);
85
86 rho = sqrt(Utilities::Frobenius(*V[0], *V[0]));
87
88 V[0]->scale(-one / rho);
89 }
90
91 std::vector<SC> h((nIts_ + 1) * (nIts_ + 1));
92 std::vector<SC> c(nIts_ + 1, 0.0);
93 std::vector<SC> s(nIts_ + 1, 0.0);
94 std::vector<SC> g(nIts_ + 1, 0.0);
95 g[0] = rho;
96
97#define I(i, j) ((i) + (j) * (nIts_ + 1)) // column ordering
98 for (size_t i = 0; i < nIts_; i++) {
99 // V_{i+1} = D^{-1}*A*V_i
100 tmpAP = MatrixMatrix::Multiply(*A, false, *V[i], false, mmfancy, true /*doFillComplete*/, true /*optimizeStorage*/);
101 C.Apply(*tmpAP, *T);
102
103 V[i + 1] = MatrixFactory2::BuildCopy(T);
104 Utilities::MyOldScaleMatrix(*V[i + 1], D, true /*doInverse*/, true /*doFillComplete*/, false /*doOptimizeStorage*/);
105
106 // Update Hessenberg matrix
107 for (size_t j = 0; j <= i; j++) {
108 h[I(j, i)] = Utilities::Frobenius(*V[i + 1], *V[j]);
109
110 // V_{i+1} = V_{i+1} - h(j,i+1)*V_j
111#ifndef TWO_ARG_MATRIX_ADD
112 newV = Teuchos::null;
113 MatrixMatrix::TwoMatrixAdd(*V[j], false, -h[I(j, i)], *V[i + 1], false, one, newV, mmfancy);
114 newV->fillComplete(V[i + 1]->getDomainMap(), V[i + 1]->getRangeMap());
115 V[i + 1].swap(newV);
116#else
117 // FIXME: this does not work now. Fails with the following exception:
118 // what(): ../../packages/tpetra/core/ext/TpetraExt_MatrixMatrix_def.hpp:408:
119 //
120 // Throw number = 1
121 //
122 // Throw test that evaluated to true: B.isLocallyIndexed()
123 //
124 // TpetraExt::MatrixMatrix::Add(): ERROR, input matrix B must not be locally indexed
125 MatrixMatrix::TwoMatrixAdd(*V[j], false, -h[I(j, i)], *V[i + 1], one);
126#endif
127 }
128 h[I(i + 1, i)] = sqrt(Utilities::Frobenius(*V[i + 1], *V[i + 1]));
129
130 // NOTE: potentially we'll need some reorthogonalization code here
131 // The matching MATLAB code is
132 // normav = norm(v.num(k+1).matrix, 'fro');
133 // normav2 = h(k+1,k);
134 // if (reorth == -1 && normav + .001*normav2 == normav)
135 // for j = 1:k
136 // hr = v(:,j)'*v(:,k+1); % hr=v(:,k+1)'*v(:,j);
137 // h(j,k) = h(j,k)+hr;
138 // v(:,k+1) = v(:,k+1)-hr*v(:,j);
139 // end
140 // h(k+1,k) = norm(v(:,k+1));
141 // end
142
143 // Check for nonsymmetric case
144 if (h[I(i + 1, i)] != zero) {
145 // Normalize V_i
146 V[i + 1]->scale(one / h[I(i + 1, i)]);
147 }
148
149 if (i > 0)
150 givapp(&c[0], &s[0], &h[I(0, i)], i); // Due to column ordering &h[...] is a column
151
152 SC nu = sqrt(h[I(i, i)] * h[I(i, i)] + h[I(i + 1, i)] * h[I(i + 1, i)]);
153 if (nu != zero) {
154 c[i] = h[I(i, i)] / nu;
155 s[i] = -h[I(i + 1, i)] / nu;
156 h[I(i, i)] = c[i] * h[I(i, i)] - s[i] * h[I(i + 1, i)];
157 h[I(i + 1, i)] = zero;
158
159 givapp(&c[i], &s[i], &g[i], 1);
160 }
161 }
162
163 // Solve Hessenberg system
164 // y = solve(H, \rho e_1)
165 std::vector<SC> y(nIts_);
166 if (nIts_ == 1) {
167 y[0] = g[0] / h[I(0, 0)];
168 }
169#undef I
170
171 // Compute final
172 for (size_t i = 0; i < nIts_; i++) {
173#ifndef TWO_ARG_MATRIX_ADD
174 newV = Teuchos::null;
175 MatrixMatrix::TwoMatrixAdd(*V[i], false, y[i], *finalP, false, one, newV, mmfancy);
176 newV->fillComplete(finalP->getDomainMap(), finalP->getRangeMap());
177 finalP.swap(newV);
178#else
179 MatrixMatrix::TwoMatrixAdd(*V[i], false, y[i], *finalP, one);
180#endif
181 }
182}
183
184} // namespace MueLu
185
186#endif // ifndef MUELU_GMRESSOLVER_DECL_HPP
#define I(i, j)
MueLu::DefaultScalar Scalar
Constraint space information for the potential prolongator.
RCP< const CrsGraph > GetPattern() const
void Apply(const Matrix &P, Matrix &Projected) const
Apply constraint.
Exception throws to report errors in the internal logical of the program.
void givapp(SC *c, SC *s, SC *v, int k) const
size_t nIts_
Number of performed iterations.
void Iterate(const Matrix &A, const Constraint &C, const Matrix &P0, RCP< Matrix > &P) const
Iterate.
static Scalar Frobenius(const Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &A, const Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &B)
static void MyOldScaleMatrix(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Teuchos::ArrayRCP< const Scalar > &scalingVector, bool doInverse=true, bool doFillComplete=true, bool doOptimizeStorage=true)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.