.\" .de Id .. .de Sp .if n .sp .if t .sp 0.4 .. .TH pminres 4rheolef "rheolef-6.7" "rheolef-6.7" "rheolef-6.7" .\" label: /*D:pminres .SH NAME \fBpminres\fP -- conjugate gradient algorithm. .\" skip: @findex pminres .\" skip: @cindex conjugate gradient algorithm .\" skip: @cindex iterative solver .\" skip: @cindex preconditioner .SH SYNOPSIS .\" begin_example .Sp .nf template int pminres (const Matrix &A, Vector &x, const Vector &b, const Preconditioner &M, int &max_iter, Real &tol, odiststream *p_derr=0); .Sp .fi .\" end_example .SH EXAMPLE The simplest call to 'pminres' has the folling form: .\" begin_example .Sp .nf size_t max_iter = 100; double tol = 1e-7; int status = pminres(a, x, b, EYE, max_iter, tol, &derr); .Sp .fi .\" end_example .SH DESCRIPTION \fBpminres\fP solves the symmetric positive definite linear system Ax=b using the Conjugate Gradient method. .PP The return value indicates convergence within max_iter (input) iterations (0), or no convergence within max_iter iterations (1). Upon successful return, output arguments have the following values: .\" begin table .\" start item .TP .B x approximate solution to Ax = b .\" start item .TP .B max_iter the number of iterations performed before the tolerance was reached .\" start item .TP .B tol the residual after the final iteration .\" end table .SH NOTE \fBpminres\fP follows the algorithm described in "Solution of sparse indefinite systems of linear equations", C. C. Paige and M. A. Saunders, SIAM J. Numer. Anal., 12(4), 1975. For more, see http://www.stanford.edu/group/SOL/software.html and also the PhD "Iterative methods for singular linear equations and least-squares problems", S.-C. T. Choi, Stanford University, 2006, http://www.stanford.edu/group/SOL/dissertations/sou-cheng-choi-thesis.pdf at page 60. The present implementation style is inspired from \fBIML++ 1.2\fP iterative method library, \fBhttp://math.nist.gov/iml++\fP. .\" skip start:AUTHOR: .\" skip start:DATE: .\" skip start:METHODS: .\" END .SH IMPLEMENTATION .\" begin_example .Sp .nf template int pminres(const Matrix &A, Vector &x, const Vector &Mb, const Preconditioner &M, Size &max_iter, Real &tol, odiststream *p_derr = 0, std::string label = "minres") { Vector b = M.solve(Mb); Real norm_b = sqrt(fabs(dot(Mb,b))); if (norm_b == Real(0.)) norm_b = 1; Vector Mr = Mb - A*x; Vector z = M.solve(Mr); Real beta2 = dot(Mr, z); Real norm_r = sqrt(fabs(beta2)); if (p_derr) (*p_derr) << "[" << label << "] # norm_b=" <