1
0
mirror of https://github.com/paboyle/Grid.git synced 2024-11-14 17:55:38 +00:00
Grid/lib/algorithms/iterative/MinimalResidual.h

217 lines
7.5 KiB
C
Raw Normal View History

2017-07-21 12:39:03 +01:00
/*************************************************************************************
Grid physics library, www.github.com/paboyle/Grid
Source file: ./lib/algorithms/iterative/MinimalResidual.h
Copyright (C) 2015
Author: Azusa Yamaguchi <ayamaguc@staffmail.ed.ac.uk>
Author: Peter Boyle <paboyle@ph.ed.ac.uk>
Author: paboyle <paboyle@ph.ed.ac.uk>
This program is free software; you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation; either version 2 of the License, or
(at your option) any later version.
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License along
with this program; if not, write to the Free Software Foundation, Inc.,
51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
See the full license in the file "LICENSE" in the top level distribution
directory
*************************************************************************************/
/* END LEGAL */
#ifndef GRID_MINIMAL_RESIDUAL_H
#define GRID_MINIMAL_RESIDUAL_H
namespace Grid {
/////////////////////////////////////////////////////////////
// Base classes for iterative processes based on operators
// single input vec, single output vec.
/////////////////////////////////////////////////////////////
2017-10-27 13:09:02 +01:00
template<class Field> class MinimalResidual : public OperatorFunction<Field> {
2017-07-21 12:39:03 +01:00
public:
2017-10-27 13:09:02 +01:00
bool ErrorOnNoConverge; // throw an assert when the MR fails to converge.
// Defaults true.
RealD Tolerance;
2017-07-21 12:39:03 +01:00
Integer MaxIterations;
2017-10-27 13:09:02 +01:00
Integer IterationsToComplete; // Number of iterations the MR took to finish. Filled in upon completion
2017-07-21 12:39:03 +01:00
2017-10-27 13:09:02 +01:00
MinimalResidual(RealD tol, Integer maxit, bool err_on_no_conv = true)
: Tolerance(tol),
MaxIterations(maxit),
ErrorOnNoConverge(err_on_no_conv){};
2017-10-25 09:38:26 +01:00
//! Minimal-residual (MR) algorithm for a generic Linear Operator
/*! \ingroup invert
* This subroutine uses the Minimal Residual (MR) algorithm to determine
2017-10-27 13:09:02 +01:00
* the solution of the set of linear equations. Here we allow M to be
nonhermitian.
2017-10-25 09:38:26 +01:00
*
* M . Psi = src
*
* Algorithm:
*
* Psi[0] Argument
* r[0] := src - M . Psi[0] ; Initial residual
* IF |r[0]| <= RsdCG |src| THEN RETURN; Converged?
* FOR k FROM 1 TO MaxCG DO MR iterations
* a[k-1] := <M.r[k-1],r[k-1]> / <M.r[k-1],M.r[k-1]> ;
* ap[k-1] := MRovpar * a[k] ; Overrelaxtion step
* Psi[k] += ap[k-1] r[k-1] ; New solution vector
* r[k] -= ap[k-1] A . r[k-1] ; New residual
* IF |r[k]| <= RsdCG |src| THEN RETURN; Converged?
* Arguments:
* \param M Linear Operator (Read)
* \param src Source (Read)
* \param psi Solution (Modify)
* \param RsdCG MR residual accuracy (Read)
* \param MRovpar Overrelaxation parameter (Read)
* \param MaxIterations Maximum MR iterations (Read)
* Local Variables:
* r Residual vector
* cp | r[k] |**2
* c | r[k-1] |**2
* k MR iteration counter
* a a[k]
* d < M.r[k], M.r[k] >
* R_Aux Temporary for M.Psi
* Mr Temporary for M.r
* Global Variables:
* MaxIterations Maximum number of MR iterations allowed
* RsdCG Maximum acceptable MR residual (relative to source)
*
* Subroutines:
*
* M Apply matrix to vector
*
* @{
*/
2017-10-27 13:09:02 +01:00
void operator()(LinearOperatorBase<Field> &Linop, const Field &src, Field &psi) {
2017-10-25 09:38:26 +01:00
psi.checkerboard = src.checkerboard;
conformable(psi, src);
Complex a, c;
2017-10-27 13:09:02 +01:00
RealD d;
2017-10-25 09:38:26 +01:00
Field Mr(src);
Field r(src);
// Initial residual computation & set up
RealD guess = norm2(psi);
assert(std::isnan(guess) == 0);
RealD ssq = norm2(src); // flopcount.addSiteFlops(4*Nc*Ns,s); // stands for "source squared"
2017-10-27 13:09:02 +01:00
RealD rsd_sq = Tolerance * Tolerance * ssq; // flopcount.addSiteFlops(4*Nc*Ns,s); //
// stands for "residual squared"
2017-10-25 09:38:26 +01:00
/* r[0] := src - M . Psi[0] */
/* r := M . Psi */
2017-10-27 13:09:02 +01:00
// M(Mr, psi, isign); // flopcount.addFlops(M.nFlops());
Linop.Op(psi, Mr); // flopcount.addFlops(M.nFlops());
2017-10-25 09:38:26 +01:00
r = src - Mr; // flopcount.addSiteFlops(2*Nc*Ns,s);
2017-10-27 13:09:02 +01:00
RealD cp = norm2(r); /* Cp = |r[0]|^2 */
/* 2 Nc Ns flops */ // flopcount.addSiteFlops(4*Nc*Ns, s);
// auto cp = norm2(r); /* Cp = |r[0]|^2 */ /* 2 Nc Ns flops */ //
// flopcount.addSiteFlops(4*Nc*Ns, s);
2017-10-25 09:38:26 +01:00
2017-10-27 13:09:02 +01:00
if(cp <= rsd_sq) { /* IF |r[0]| <= Tolerance|src| THEN RETURN; */
2017-10-25 09:38:26 +01:00
return;
}
std::cout << GridLogIterative << std::setprecision(4)
2017-10-27 13:09:02 +01:00
<< "MinimalResidual: k=0 residual " << cp << " target " << rsd_sq << std::endl;
2017-10-25 09:38:26 +01:00
2017-10-27 13:09:02 +01:00
GridStopWatch LinalgTimer;
GridStopWatch MatrixTimer;
GridStopWatch SolverTimer;
SolverTimer.Start();
2017-10-25 09:38:26 +01:00
auto k = 0;
2017-10-27 13:09:02 +01:00
while((k < MaxIterations) && (cp > rsd_sq)) {
2017-10-25 09:38:26 +01:00
++k;
/* a[k-1] := < M.r[k-1], r[k-1] >/ < M.r[k-1], M.r[k-1] > ; */
2017-10-27 13:09:02 +01:00
MatrixTimer.Start();
// M(Mr, r, isign); /* Mr = M * r */ // flopcount.addFlops(M.nFlops());
Linop.Op(r, Mr); /* Mr = M * r */ // flopcount.addFlops(M.nFlops());
MatrixTimer.Stop();
LinalgTimer.Start();
2017-10-25 09:38:26 +01:00
c = innerProduct(Mr, r); /* c = < M.r, r > */ // flopcount.addSiteFlops(4*Nc*Ns,s);
d = norm2(Mr); /* d = | M.r | ** 2 */ // flopcount.addSiteFlops(4*Nc*Ns,s);
2017-10-27 13:09:02 +01:00
a = c / d;
2017-10-25 09:38:26 +01:00
2017-10-27 13:09:02 +01:00
// a = a * MRovpar; /* a[k-1] *= MRovpar ; */
2017-10-25 09:38:26 +01:00
2017-10-27 13:09:02 +01:00
psi = psi + r * a; /* Psi[k] += a[k-1] r[k-1] ; */ // flopcount.addSiteFlops(4*Nc*Ns,s);
2017-10-25 09:38:26 +01:00
r = r - Mr * a; /* r[k] -= a[k-1] M . r[k-1] ; */ // flopcount.addSiteFlops(4*Nc*Ns,s);
cp = norm2(r); /* cp = | r[k] |**2 */ // flopcount.addSiteFlops(4*Nc*Ns,s);
2017-10-27 13:09:02 +01:00
LinalgTimer.Stop();
std::cout << GridLogIterative << "MinimalResidual: Iteration " << k
<< " residual " << cp << " target " << rsd_sq << std::endl;
2017-10-25 09:38:26 +01:00
}
2017-10-27 13:09:02 +01:00
SolverTimer.Stop();
2017-10-25 09:38:26 +01:00
IterationsToComplete = k;
2017-10-27 13:09:02 +01:00
// res.resid = sqrt(cp);
std::cout << "InvMR: k = " << k << " cp = " << cp << std::endl;
2017-10-25 09:38:26 +01:00
// flopcount.report("invmr", swatch.getTimeInSeconds());
2017-10-27 13:09:02 +01:00
std::cout << GridLogMessage << "MinimalResidual Converged on iteration " << k << std::endl;
std::cout << GridLogMessage << "\tComputed residual " << sqrt(cp / ssq)<<std::endl;
// std::cout << GridLogMessage << "\tTrue residual " << true_residual<<std::endl;
// std::cout << GridLogMessage << "\tTarget " << Tolerance << std::endl;
std::cout << GridLogMessage << "Time breakdown "<<std::endl;
std::cout << GridLogMessage << "\tElapsed " << SolverTimer.Elapsed() <<std::endl;
std::cout << GridLogMessage << "\tMatrix " << MatrixTimer.Elapsed() <<std::endl;
std::cout << GridLogMessage << "\tLinalg " << LinalgTimer.Elapsed() <<std::endl;
2017-10-25 09:38:26 +01:00
// Compute the actual residual
{
2017-10-27 13:09:02 +01:00
// M(Mr, psi, isign);
Linop.Op(psi, Mr);
Field tmp = src - Mr;
// RealD actual_res = norm2(src-Mr);
RealD actual_res = norm2(tmp);
// res.resid = sqrt(actual_res);
2017-10-25 09:38:26 +01:00
}
2017-10-27 13:09:02 +01:00
if(IterationsToComplete == MaxIterations)
std::cerr << "Nonconvergence Warning" << std::endl;
2017-10-25 09:38:26 +01:00
2017-10-27 13:09:02 +01:00
// return res;
2017-10-25 09:38:26 +01:00
}
2017-07-21 12:39:03 +01:00
};
2017-10-27 13:09:02 +01:00
} // namespace Grid
2017-07-21 12:39:03 +01:00
#endif