Grid 0.7.0
MixedPrecisionFlexibleGeneralisedMinimalResidual.h
Go to the documentation of this file.
1/*************************************************************************************
2
3Grid physics library, www.github.com/paboyle/Grid
4
5Source file: ./lib/algorithms/iterative/MixedPrecisionFlexibleGeneralisedMinimalResidual.h
6
7Copyright (C) 2015
8
9Author: Daniel Richtmann <daniel.richtmann@ur.de>
10
11This program is free software; you can redistribute it and/or modify
12it under the terms of the GNU General Public License as published by
13the Free Software Foundation; either version 2 of the License, or
14(at your option) any later version.
15
16This program is distributed in the hope that it will be useful,
17but WITHOUT ANY WARRANTY; without even the implied warranty of
18MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19GNU General Public License for more details.
20
21You should have received a copy of the GNU General Public License along
22with this program; if not, write to the Free Software Foundation, Inc.,
2351 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24
25See the full license in the file "LICENSE" in the top level distribution
26directory
27*************************************************************************************/
28/* END LEGAL */
29#ifndef GRID_MIXED_PRECISION_FLEXIBLE_GENERALISED_MINIMAL_RESIDUAL_H
30#define GRID_MIXED_PRECISION_FLEXIBLE_GENERALISED_MINIMAL_RESIDUAL_H
31
32namespace Grid {
33
34template<class FieldD, class FieldF, typename std::enable_if<getPrecision<FieldD>::value == 2, int>::type = 0, typename std::enable_if< getPrecision<FieldF>::value == 1, int>::type = 0>
36 public:
37
38 using OperatorFunction<FieldD>::operator();
39
40 bool ErrorOnNoConverge; // Throw an assert when MPFGMRES fails to converge,
41 // defaults to true
42
44
48 Integer IterationCount; // Number of iterations the MPFGMRES took to finish,
49 // filled in upon completion
50
57
58 Eigen::MatrixXcd H;
59
60 std::vector<ComplexD> y;
61 std::vector<ComplexD> gamma;
62 std::vector<ComplexD> c;
63 std::vector<ComplexD> s;
64
66
68
70 Integer maxit,
71 GridBase * sp_grid,
73 Integer restart_length,
74 bool err_on_no_conv = true)
75 : Tolerance(tol)
76 , MaxIterations(maxit)
77 , RestartLength(restart_length)
79 , ErrorOnNoConverge(err_on_no_conv)
80 , H(Eigen::MatrixXcd::Zero(RestartLength, RestartLength + 1)) // sizes taken from DD-αAMG code base
81 , y(RestartLength + 1, 0.)
82 , gamma(RestartLength + 1, 0.)
83 , c(RestartLength + 1, 0.)
84 , s(RestartLength + 1, 0.)
85 , SinglePrecGrid(sp_grid)
86 , Preconditioner(Prec) {};
87
88 void operator()(LinearOperatorBase<FieldD> &LinOp, const FieldD &src, FieldD &psi) {
89
90 psi.Checkerboard() = src.Checkerboard();
91 conformable(psi, src);
92
93 RealD guess = norm2(psi);
94 assert(std::isnan(guess) == 0);
95
96 RealD cp;
97 RealD ssq = norm2(src);
98 RealD rsq = Tolerance * Tolerance * ssq;
99
100 FieldD r(src.Grid());
101
102 std::cout << std::setprecision(4) << std::scientific;
103 std::cout << GridLogIterative << "MPFGMRES: guess " << guess << std::endl;
104 std::cout << GridLogIterative << "MPFGMRES: src " << ssq << std::endl;
105
106 PrecTimer.Reset();
107 MatrixTimer.Reset();
108 LinalgTimer.Reset();
109 QrTimer.Reset();
110 CompSolutionTimer.Reset();
111 ChangePrecTimer.Reset();
112
113 GridStopWatch SolverTimer;
114 SolverTimer.Start();
115
116 IterationCount = 0;
117
118 for (int k=0; k<MaxNumberOfRestarts; k++) {
119
120 cp = outerLoopBody(LinOp, src, psi, rsq);
121
122 // Stopping condition
123 if (cp <= rsq) {
124
125 SolverTimer.Stop();
126
127 LinOp.Op(psi,r);
128 axpy(r,-1.0,src,r);
129
130 RealD srcnorm = sqrt(ssq);
131 RealD resnorm = sqrt(norm2(r));
132 RealD true_residual = resnorm / srcnorm;
133
134 std::cout << GridLogMessage << "MPFGMRES: Converged on iteration " << IterationCount
135 << " computed residual " << sqrt(cp / ssq)
136 << " true residual " << true_residual
137 << " target " << Tolerance << std::endl;
138
139 std::cout << GridLogMessage << "MPFGMRES Time elapsed: Total " << SolverTimer.Elapsed() << std::endl;
140 std::cout << GridLogMessage << "MPFGMRES Time elapsed: Precon " << PrecTimer.Elapsed() << std::endl;
141 std::cout << GridLogMessage << "MPFGMRES Time elapsed: Matrix " << MatrixTimer.Elapsed() << std::endl;
142 std::cout << GridLogMessage << "MPFGMRES Time elapsed: Linalg " << LinalgTimer.Elapsed() << std::endl;
143 std::cout << GridLogMessage << "MPFGMRES Time elapsed: QR " << QrTimer.Elapsed() << std::endl;
144 std::cout << GridLogMessage << "MPFGMRES Time elapsed: CompSol " << CompSolutionTimer.Elapsed() << std::endl;
145 std::cout << GridLogMessage << "MPFGMRES Time elapsed: PrecChange " << ChangePrecTimer.Elapsed() << std::endl;
146 return;
147 }
148 }
149
150 std::cout << GridLogMessage << "MPFGMRES did NOT converge" << std::endl;
151
153 assert(0);
154 }
155
156 RealD outerLoopBody(LinearOperatorBase<FieldD> &LinOp, const FieldD &src, FieldD &psi, RealD rsq) {
157
158 RealD cp = 0;
159
160 FieldD w(src.Grid());
161 FieldD r(src.Grid());
162
163 // these should probably be made class members so that they are only allocated once, not in every restart
164 std::vector<FieldD> v(RestartLength + 1, src.Grid()); for (auto &elem : v) elem = Zero();
165 std::vector<FieldD> z(RestartLength + 1, src.Grid()); for (auto &elem : z) elem = Zero();
166
167 MatrixTimer.Start();
168 LinOp.Op(psi, w);
169 MatrixTimer.Stop();
170
171 LinalgTimer.Start();
172 r = src - w;
173
174 gamma[0] = sqrt(norm2(r));
175
176 v[0] = (1. / gamma[0]) * r;
177 LinalgTimer.Stop();
178
179 for (int i=0; i<RestartLength; i++) {
180
182
183 arnoldiStep(LinOp, v, z, w, i);
184
185 qrUpdate(i);
186
187 cp = norm(gamma[i+1]);
188
189 std::cout << GridLogIterative << "MPFGMRES: Iteration " << IterationCount
190 << " residual " << cp << " target " << rsq << std::endl;
191
192 if ((i == RestartLength - 1) || (IterationCount == MaxIterations) || (cp <= rsq)) {
193
194 computeSolution(z, psi, i);
195
196 return cp;
197 }
198 }
199
200 assert(0); // Never reached
201 return cp;
202 }
203
204 void arnoldiStep(LinearOperatorBase<FieldD> &LinOp, std::vector<FieldD> &v, std::vector<FieldD> &z, FieldD &w, int iter) {
205
206 FieldF v_f(SinglePrecGrid);
207 FieldF z_f(SinglePrecGrid);
208
209 ChangePrecTimer.Start();
210 precisionChange(v_f, v[iter]);
211 precisionChange(z_f, z[iter]);
212 ChangePrecTimer.Stop();
213
214 PrecTimer.Start();
215 Preconditioner(v_f, z_f);
216 PrecTimer.Stop();
217
218 ChangePrecTimer.Start();
219 precisionChange(z[iter], z_f);
220 ChangePrecTimer.Stop();
221
222 MatrixTimer.Start();
223 LinOp.Op(z[iter], w);
224 MatrixTimer.Stop();
225
226 LinalgTimer.Start();
227 for (int i = 0; i <= iter; ++i) {
228 H(iter, i) = innerProduct(v[i], w);
229 w = w - ComplexD(H(iter, i)) * v[i];
230 }
231
232 H(iter, iter + 1) = sqrt(norm2(w));
233 v[iter + 1] = ComplexD(1. / H(iter, iter + 1)) * w;
234 LinalgTimer.Stop();
235 }
236
237 void qrUpdate(int iter) {
238
239 QrTimer.Start();
240 for (int i = 0; i < iter ; ++i) {
241 auto tmp = -s[i] * ComplexD(H(iter, i)) + c[i] * ComplexD(H(iter, i + 1));
242 H(iter, i) = conjugate(c[i]) * ComplexD(H(iter, i)) + conjugate(s[i]) * ComplexD(H(iter, i + 1));
243 H(iter, i + 1) = tmp;
244 }
245
246 // Compute new Givens Rotation
247 auto nu = sqrt(std::norm(H(iter, iter)) + std::norm(H(iter, iter + 1)));
248 c[iter] = H(iter, iter) / nu;
249 s[iter] = H(iter, iter + 1) / nu;
250
251 // Apply new Givens rotation
252 H(iter, iter) = nu;
253 H(iter, iter + 1) = 0.;
254
255 gamma[iter + 1] = -s[iter] * gamma[iter];
256 gamma[iter] = conjugate(c[iter]) * gamma[iter];
257 QrTimer.Stop();
258 }
259
260 void computeSolution(std::vector<FieldD> const &z, FieldD &psi, int iter) {
261
262 CompSolutionTimer.Start();
263 for (int i = iter; i >= 0; i--) {
264 y[i] = gamma[i];
265 for (int k = i + 1; k <= iter; k++)
266 y[i] = y[i] - ComplexD(H(k, i)) * y[k];
267 y[i] = y[i] / ComplexD(H(i, i));
268 }
269
270 for (int i = 0; i <= iter; i++)
271 psi = psi + z[i] * y[i];
272 CompSolutionTimer.Stop();
273 }
274};
275}
276#endif
accelerator_inline Grid_simd< S, V > sqrt(const Grid_simd< S, V > &r)
void axpy(Lattice< vobj > &ret, sobj a, const Lattice< vobj > &x, const Lattice< vobj > &y)
void conformable(const Lattice< obj1 > &lhs, const Lattice< obj2 > &rhs)
Lattice< vobj > conjugate(const Lattice< vobj > &lhs)
ComplexD innerProduct(const Lattice< vobj > &left, const Lattice< vobj > &right)
RealD norm2(const Lattice< vobj > &arg)
void precisionChange(Lattice< VobjOut > &out, const Lattice< VobjIn > &in, const precisionChangeWorkspace &workspace)
GridLogger GridLogIterative(1, "Iterative", GridLogColours, "BLUE")
GridLogger GridLogMessage(1, "Message", GridLogColours, "NORMAL")
uint32_t Integer
Definition Simd.h:58
std::complex< RealD > ComplexD
Definition Simd.h:79
double RealD
Definition Simd.h:61
void Start(void)
Definition Timer.h:92
GridTime Elapsed(void) const
Definition Timer.h:113
void Stop(void)
Definition Timer.h:99
void operator()(LinearOperatorBase< FieldD > &LinOp, const FieldD &src, FieldD &psi)
void computeSolution(std::vector< FieldD > const &z, FieldD &psi, int iter)
void arnoldiStep(LinearOperatorBase< FieldD > &LinOp, std::vector< FieldD > &v, std::vector< FieldD > &z, FieldD &w, int iter)
RealD outerLoopBody(LinearOperatorBase< FieldD > &LinOp, const FieldD &src, FieldD &psi, RealD rsq)
MixedPrecisionFlexibleGeneralisedMinimalResidual(RealD tol, Integer maxit, GridBase *sp_grid, LinearFunction< FieldF > &Prec, Integer restart_length, bool err_on_no_conv=true)
virtual void Op(const Field &in, Field &out)=0
Definition Simd.h:194