13 #include "../general/globals.hpp"
56 if (dot_prod_type == 0)
62 double local_dot = (x * y);
65 MPI_Allreduce(&local_dot, &global_dot, 1, MPI_DOUBLE, MPI_SUM, comm);
77 if (dot_prod_type == 0)
84 MPI_Comm_rank(comm, &rank);
156 double r0, nom, nom0, nomold = 1, cf;
172 nom0 = nom =
Dot(
z,
r);
176 nom0 = nom =
Dot(
r,
r);
180 mfem::out <<
" Iteration : " << setw(3) << 0 <<
" (B r, r) = "
219 cf = sqrt(nom/nomold);
221 mfem::out <<
" Iteration : " << setw(3) << i <<
" (B r, r) = "
222 << nom <<
"\tConv. rate: " << cf <<
'\n';
228 mfem::out <<
"Number of SLI iterations: " << i <<
'\n'
229 <<
"Conv. rate: " << cf <<
'\n';
231 mfem::out <<
"(B r_0, r_0) = " << nom0 <<
'\n'
232 <<
"(B r_N, r_N) = " << nom <<
'\n'
233 <<
"Number of SLI iterations: " << i <<
'\n';
247 mfem::err <<
"SLI: No convergence!" <<
'\n';
248 mfem::out <<
"(B r_0, r_0) = " << nom0 <<
'\n'
249 <<
"(B r_N, r_N) = " << nom <<
'\n'
250 <<
"Number of SLI iterations: " <<
final_iter <<
'\n';
254 mfem::out <<
"Average reduction factor = "
261 int print_iter,
int max_num_iter,
262 double RTOLERANCE,
double ATOLERANCE)
274 int print_iter,
int max_num_iter,
275 double RTOLERANCE,
double ATOLERANCE)
298 double r0, den, nom, nom0, betanom,
alpha, beta;
320 nom0 = nom =
Dot(
d,
r);
321 MFEM_ASSERT(
IsFinite(nom),
"nom = " << nom);
325 mfem::out <<
" Iteration : " << setw(3) << 0 <<
" (B r, r) = "
340 MFEM_ASSERT(
IsFinite(den),
"den = " << den);
344 mfem::out <<
"Negative denominator in step 0 of PCG: " << den <<
'\n';
373 MFEM_ASSERT(
IsFinite(betanom),
"betanom = " << betanom);
377 mfem::out <<
" Iteration : " << setw(3) << i <<
" (B r, r) = "
385 mfem::out <<
"Number of PCG iterations: " << i <<
'\n';
389 mfem::out <<
" Iteration : " << setw(3) << i <<
" (B r, r) = "
413 MFEM_ASSERT(
IsFinite(den),
"den = " << den);
417 mfem::out <<
"PCG: The operator is not positive definite. (Ad, d) = "
428 mfem::out <<
" Iteration : " << setw(3) << 0 <<
" (B r, r) = "
434 mfem::out <<
"PCG: No convergence!" <<
'\n';
438 mfem::out <<
"Average reduction factor = "
439 << pow (betanom/nom0, 0.5/
final_iter) <<
'\n';
445 int print_iter,
int max_num_iter,
446 double RTOLERANCE,
double ATOLERANCE)
458 int print_iter,
int max_num_iter,
459 double RTOLERANCE,
double ATOLERANCE)
473 double &cs,
double &sn)
480 else if (fabs(dy) > fabs(dx))
482 double temp = dx / dy;
483 sn = 1.0 / sqrt( 1.0 + temp*temp );
488 double temp = dy / dx;
489 cs = 1.0 / sqrt( 1.0 + temp*temp );
496 double temp = cs * dx + sn * dy;
497 dy = -sn * dx + cs * dy;
507 for (
int i = k; i >= 0; i--)
510 for (
int j = i - 1; j >= 0; j--)
512 y(j) -= h(j,i) * y(i);
516 for (
int j = 0; j <= k; j++)
569 double beta =
Norm(r);
570 MFEM_ASSERT(
IsFinite(beta),
"beta = " << beta);
585 <<
" Iteration : " << setw(3) << 0
586 <<
" ||B r|| = " << beta << (
print_level == 3 ?
" ...\n" :
"\n");
593 if (v[0] == NULL) { v[0] =
new Vector(n); }
594 v[0]->Set(1.0/beta, r);
595 s = 0.0; s(0) = beta;
597 for (i = 0; i <
m && j <=
max_iter; i++, j++)
609 for (k = 0; k <= i; k++)
611 H(k,i) =
Dot(w, *v[k]);
612 w.
Add(-H(k,i), *v[k]);
616 MFEM_ASSERT(
IsFinite(H(i+1,i)),
"Norm(w) = " << H(i+1,i));
617 if (v[i+1] == NULL) { v[i+1] =
new Vector(n); }
618 v[i+1]->Set(1.0/H(i+1,i), w);
620 for (k = 0; k < i; k++)
629 resid = fabs(s(i+1));
630 MFEM_ASSERT(
IsFinite(resid),
"resid = " << resid);
643 mfem::out <<
" Pass : " << setw(2) << (j-1)/
m+1
644 <<
" Iteration : " << setw(3) << j
645 <<
" ||B r|| = " << resid <<
'\n';
667 MFEM_ASSERT(
IsFinite(beta),
"beta = " << beta);
696 for (i = 0; i < v.
Size(); i++)
721 double beta =
Norm(r);
722 MFEM_ASSERT(
IsFinite(beta),
"beta = " << beta);
736 <<
" Iteration : " << setw(3) << 0
737 <<
" || r || = " << beta << endl;
741 for (i= 0; i<=
m; i++)
750 if (v[0] == NULL) { v[0] =
new Vector(b.
Size()); }
752 v[0] ->
Add (1.0/beta, r);
753 s = 0.0; s(0) = beta;
758 if (z[i] == NULL) { z[i] =
new Vector(b.
Size()); }
764 for (k = 0; k <= i; k++)
766 H(k,i) =
Dot( r, *v[k]);
767 r.Add(-H(k,i), (*v[k]));
771 if (v[i+1] == NULL) { v[i+1] =
new Vector(b.
Size()); }
773 v[i+1] ->
Add (1.0/H(i+1,i), r);
775 for (k = 0; k < i; k++)
784 double resid = fabs(s(i+1));
785 MFEM_ASSERT(
IsFinite(resid),
"resid = " << resid);
788 <<
" Iteration : " << setw(3) << i+1
789 <<
" || r || = " << resid << endl;
797 for (i= 0; i<=
m; i++)
799 if (v[i]) {
delete v[i]; }
800 if (z[i]) {
delete z[i]; }
816 MFEM_ASSERT(
IsFinite(beta),
"beta = " << beta);
822 for (i= 0; i<=
m; i++)
824 if (v[i]) {
delete v[i]; }
825 if (z[i]) {
delete z[i]; }
833 for (i = 0; i <=
m; i++)
835 if (v[i]) {
delete v[i]; }
836 if (z[i]) {
delete z[i]; }
845 int &max_iter,
int m,
double &tol,
double atol,
int printit)
862 int print_iter,
int max_num_iter,
int m,
double rtol,
double atol)
864 GMRES(A, x, b, B, max_num_iter, m, rtol, atol, print_iter);
886 double resid, tol_goal;
887 double rho_1, rho_2=1.0,
alpha=1.0, beta, omega=1.0;
902 MFEM_ASSERT(
IsFinite(resid),
"resid = " << resid);
904 mfem::out <<
" Iteration : " << setw(3) << 0
905 <<
" ||r|| = " << resid <<
'\n';
909 if (resid <= tol_goal)
923 mfem::out <<
" Iteration : " << setw(3) << i
924 <<
" ||r|| = " << resid <<
'\n';
936 beta = (rho_1/rho_2) * (
alpha/omega);
952 MFEM_ASSERT(
IsFinite(resid),
"resid = " << resid);
953 if (resid < tol_goal)
957 mfem::out <<
" Iteration : " << setw(3) << i
958 <<
" ||s|| = " << resid <<
'\n';
965 mfem::out <<
" Iteration : " << setw(3) << i
966 <<
" ||s|| = " << resid;
983 MFEM_ASSERT(
IsFinite(resid),
"resid = " << resid);
986 mfem::out <<
" ||r|| = " << resid <<
'\n';
988 if (resid < tol_goal)
1010 int &max_iter,
double &tol,
double atol,
int printit)
1019 bicgstab.
Mult(b, x);
1026 int print_iter,
int max_num_iter,
double rtol,
double atol)
1028 BiCGSTAB(A, x, b, B, max_num_iter, rtol, atol, print_iter);
1054 double beta, eta, gamma0, gamma1, sigma0, sigma1;
1055 double alpha, delta, rho1, rho2, rho3, norm_goal;
1075 eta = beta = sqrt(
Dot(*z,
v1));
1076 MFEM_ASSERT(
IsFinite(eta),
"eta = " << eta);
1077 gamma0 = gamma1 = 1.;
1078 sigma0 = sigma1 = 0.;
1082 if (eta <= norm_goal)
1090 mfem::out <<
"MINRES: iteration " << setw(3) << 0 <<
": ||r||_B = "
1103 MFEM_ASSERT(
IsFinite(alpha),
"alpha = " << alpha);
1110 delta = gamma1*alpha - gamma0*sigma1*beta;
1112 rho2 = sigma1*alpha + gamma0*gamma1*beta;
1122 MFEM_ASSERT(
IsFinite(beta),
"beta = " << beta);
1123 rho1 = hypot(delta, beta);
1127 w0.
Set(1./rho1, *z);
1131 add(1./rho1, *z, -rho2/rho1,
w1,
w0);
1135 add(-rho3/rho1,
w0, -rho2/rho1,
w1,
w0);
1136 w0.
Add(1./rho1, *z);
1140 gamma1 = delta/rho1;
1142 x.
Add(gamma1*eta,
w0);
1148 MFEM_ASSERT(
IsFinite(eta),
"eta = " << eta);
1150 if (fabs(eta) <= norm_goal)
1157 mfem::out <<
"MINRES: iteration " << setw(3) << it <<
": ||r||_B = "
1158 << fabs(eta) <<
'\n';
1193 eta = sqrt(
Dot(*z,
v1));
1194 mfem::out <<
"MINRES: iteration " << setw(3) << it <<
": ||r||_B = "
1195 << eta <<
" (re-computed)" <<
'\n';
1200 mfem::out <<
"MINRES: No convergence!\n";
1205 int max_it,
double rtol,
double atol)
1217 int print_it,
int max_it,
double rtol,
double atol)
1235 MFEM_ASSERT(
height ==
width,
"square Operator is required.");
1243 MFEM_ASSERT(
oper != NULL,
"the Operator is not set (use SetOperator).");
1244 MFEM_ASSERT(
prec != NULL,
"the Solver is not set (use SetSolver).");
1247 double norm0, norm, norm_goal;
1261 norm0 = norm =
Norm(
r);
1267 for (it = 0;
true; it++)
1269 MFEM_ASSERT(
IsFinite(norm),
"norm = " << norm);
1272 mfem::out <<
"Newton iteration " << setw(2) << it
1273 <<
" : ||r|| = " << norm;
1276 mfem::out <<
", ||r||/||r_0|| = " << norm/norm0;
1281 if (norm <= norm_goal)
1303 add(x, -c_scale,
c, x);
1320 int m_max,
int m_min,
int m_step,
double cf,
1321 double &tol,
double &atol,
int printit)
1328 Vector s(m+1), cs(m+1), sn(m+1);
1335 double normb = w.Norml2();
1345 double beta = r.
Norml2();
1347 resid = beta / normb;
1349 if (resid * resid <= tol)
1351 tol = resid * resid;
1358 <<
" Iteration : " << setw(3) << 0
1359 <<
" (r, r) = " << beta*beta <<
'\n';
1361 tol *= (normb*normb);
1362 tol = (atol > tol) ? atol : tol;
1366 for (i= 0; i<=m; i++)
1373 while (j <= max_iter)
1376 v[0] ->
Add (1.0/beta, r);
1377 s = 0.0; s(0) = beta;
1381 for (i = 0; i < m && j <= max_iter; i++)
1386 for (k = 0; k <= i; k++)
1388 H(k,i) = w * (*v[k]);
1389 w.
Add(-H(k,i), (*v[k]));
1392 H(i+1,i) = w.Norml2();
1394 v[i+1] ->
Add (1.0/H(i+1,i), w);
1396 for (k = 0; k < i; k++)
1405 resid = fabs(s(i+1));
1408 <<
" Iteration : " << setw(3) << i+1
1409 <<
" (r, r) = " << resid*resid <<
'\n';
1411 if ( resid*resid < tol)
1414 tol = resid * resid;
1416 for (i= 0; i<=m; i++)
1435 if ( resid*resid < tol)
1437 tol = resid * resid;
1439 for (i= 0; i<=m; i++)
1448 if (m - m_step >= m_min)
1461 tol = resid * resid;
1462 for (i= 0; i<=m; i++)
1484 mfem_error(
"SLBQPOptimizer::SetPreconditioner() : "
1485 "not meaningful for this solver");
1490 mfem_error(
"SLBQPOptimizer::SetOperator() : "
1491 "not meaningful for this solver");
1497 mfem::out <<
"SLBQP iteration " << it <<
": residual = " << r
1498 <<
", lambda = " << l <<
'\n';
1520 const double smin = 0.1;
1527 mfem::out <<
"SLBQP bracketing phase" <<
'\n';
1531 r =
solve(l,xt,x,nclip);
1545 llow = l; rlow = r; l = l + dl;
1548 r =
solve(l,xt,x,nclip);
1551 while ((r < 0) && (nclip <
max_iter))
1555 if (s < smin) { s = smin; }
1560 r =
solve(l,xt,x,nclip);
1568 lupp = l; rupp = r; l = l - dl;
1571 r =
solve(l,xt,x,nclip);
1574 while ((r > 0) && (nclip <
max_iter))
1578 if (s < smin) { s = smin; }
1583 r =
solve(l,xt,x,nclip);
1596 mfem::out <<
"SLBQP secant phase" <<
'\n';
1599 s = 1.0 - rlow/rupp; dl = dl/s; l = lupp - dl;
1602 r =
solve(l,xt,x,nclip);
1605 while ( (fabs(r) > tol) && (nclip <
max_iter) )
1611 lupp = l; rupp = r; s = 1.0 - rlow/rupp;
1612 dl = (lupp - llow)/s; l = lupp - dl;
1617 if (s < smin) { s = smin; }
1619 lnew = 0.75*llow + 0.25*l;
1620 if (lnew < l-dl) { lnew = l-dl; }
1621 lupp = l; rupp = r; l = lnew;
1622 s = (lupp - llow)/(lupp - l);
1630 llow = l; rlow = r; s = 1.0 - rlow/rupp;
1631 dl = (lupp - llow)/s; l = lupp - dl;
1636 if (s < smin) { s = smin; }
1638 lnew = 0.75*lupp + 0.25*l;
1639 if (lnew < l+dl) { lnew = l+dl; }
1640 llow = l; rlow = r; l = lnew;
1641 s = (lupp - llow)/(lupp - l);
1646 r =
solve(l,xt,x,nclip);
1655 mfem::err <<
"SLBQP not converged!" <<
'\n';
1665 mfem::out <<
"SLBQP iterations = " << nclip <<
'\n';
1666 mfem::out <<
"SLBQP lambda = " << l <<
'\n';
1667 mfem::out <<
"SLBQP residual = " << r <<
'\n';
1671 #ifdef MFEM_USE_SUITESPARSE
1698 umfpack_di_free_numeric(&
Numeric);
1702 umfpack_dl_free_numeric(&
Numeric);
1709 MFEM_ABORT(
"not a SparseMatrix");
1719 MFEM_VERIFY(
width ==
height,
"not a square matrix");
1727 int status = umfpack_di_symbolic(
width,
width, Ap, Ai, Ax, &Symbolic,
1732 umfpack_di_report_status(
Control, status);
1734 " umfpack_di_symbolic() failed!");
1737 status = umfpack_di_numeric(Ap, Ai, Ax, Symbolic, &
Numeric,
1742 umfpack_di_report_status(
Control, status);
1744 " umfpack_di_numeric() failed!");
1746 umfpack_di_free_symbolic(&Symbolic);
1750 SuiteSparse_long status;
1754 AI =
new SuiteSparse_long[
width + 1];
1755 AJ =
new SuiteSparse_long[Ap[
width]];
1756 for (
int i = 0; i <=
width; i++)
1758 AI[i] = (SuiteSparse_long)(Ap[i]);
1760 for (
int i = 0; i < Ap[
width]; i++)
1762 AJ[i] = (SuiteSparse_long)(Ai[i]);
1765 status = umfpack_dl_symbolic(
width,
width,
AI,
AJ, Ax, &Symbolic,
1770 umfpack_dl_report_status(
Control, status);
1772 " umfpack_dl_symbolic() failed!");
1775 status = umfpack_dl_numeric(
AI,
AJ, Ax, Symbolic, &
Numeric,
1780 umfpack_dl_report_status(
Control, status);
1782 " umfpack_dl_numeric() failed!");
1784 umfpack_dl_free_symbolic(&Symbolic);
1791 mfem_error(
"UMFPackSolver::Mult : matrix is not set!"
1792 " Call SetOperator first!");
1799 umfpack_di_report_info(Control,
Info);
1802 umfpack_di_report_status(Control, status);
1803 mfem_error(
"UMFPackSolver::Mult : umfpack_di_solve() failed!");
1808 SuiteSparse_long status =
1811 umfpack_dl_report_info(Control,
Info);
1814 umfpack_dl_report_status(Control, status);
1815 mfem_error(
"UMFPackSolver::Mult : umfpack_dl_solve() failed!");
1823 mfem_error(
"UMFPackSolver::MultTranspose : matrix is not set!"
1824 " Call SetOperator first!");
1831 umfpack_di_report_info(Control,
Info);
1834 umfpack_di_report_status(Control, status);
1836 " umfpack_di_solve() failed!");
1841 SuiteSparse_long status =
1844 umfpack_dl_report_info(Control,
Info);
1847 umfpack_dl_report_status(Control, status);
1849 " umfpack_dl_solve() failed!");
1862 umfpack_di_free_numeric(&
Numeric);
1866 umfpack_dl_free_numeric(&
Numeric);
1881 "Had Numeric pointer in KLU, but not Symbolic");
1889 MFEM_VERIFY(
mat != NULL,
"not a SparseMatrix");
1897 MFEM_VERIFY(
width ==
height,
"not a square matrix");
1909 MFEM_VERIFY(
mat != NULL,
1910 "KLUSolver::Mult : matrix is not set! Call SetOperator first!");
1923 MFEM_VERIFY(
mat != NULL,
1924 "KLUSolver::Mult : matrix is not set! Call SetOperator first!");
1943 #endif // MFEM_USE_SUITESPARSE
void ApplyPlaneRotation(double &dx, double &dy, double &cs, double &sn)
int Size() const
Logical size of the array.
Conjugate gradient method.
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
double Info[UMFPACK_INFO]
int GetNumIterations() const
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
void SetSize(int s)
Resize the vector to size s.
virtual Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate ...
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
double Norml2() const
Returns the l2 norm of the vector.
Data type dense matrix using column-major storage.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
int Size() const
Returns the size of the vector.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void SortColumnIndices()
Sort the column indices corresponding to each row.
void add(const Vector &v1, const Vector &v2, Vector &v)
double Norm(const Vector &x) const
double GetFinalNorm() const
int BiCGSTAB(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, double &tol, double atol, int printit)
BiCGSTAB method. (tolerances are squared)
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
bool IsFinite(const double &val)
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
void SetLinearConstraint(const Vector &_w, double _a)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
void SetPrintLevel(int print_lvl)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
double solve(double l, const Vector &xt, Vector &x, int &nclip) const
Solve QP at fixed lambda.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void MINRES(const Operator &A, const Vector &b, Vector &x, int print_it, int max_it, double rtol, double atol)
MINRES method without preconditioner. (tolerances are squared)
void Add(const double c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
void SetMaxIter(int max_it)
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
void print_iteration(int it, double r, double l) const
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
double Dot(const Vector &x, const Vector &y) const
double * GetData() const
Return element data, i.e. array A.
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
int * GetI() const
Return the array I.
Stationary linear iteration: x <- x + B (b - A x)
void Swap(Array< T > &, Array< T > &)
void SetAbsTol(double atol)
void SetRelTol(double rtol)
double Control[UMFPACK_CONTROL]
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
void subtract(const Vector &x, const Vector &y, Vector &z)
void mfem_error(const char *msg)
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
void Update(Vector &x, int k, DenseMatrix &h, Vector &s, Array< Vector * > &v)
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Vector & Set(const double a, const Vector &x)
(*this) = a * x
int height
Dimension of the output / number of rows in the matrix.
void SetBounds(const Vector &_lo, const Vector &_hi)
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const
This method can be overloaded in derived classes to implement line search algorithms.
void GeneratePlaneRotation(double &dx, double &dy, double &cs, double &sn)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void Mult(const Vector &xt, Vector &x) const
Operator application: y=A(x).
virtual void SetOperator(const Operator &op)=0
Set/update the solver for the given operator.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
int aGMRES(const Operator &A, Vector &x, const Vector &b, const Operator &M, int &max_iter, int m_max, int m_min, int m_step, double cf, double &tol, double &atol, int printit)
virtual void SetPreconditioner(Solver &pr)
These are not currently meaningful for this solver and will error out.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
int * GetJ() const
Return the array J.
int width
Dimension of the input / number of columns in the matrix.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
void SLI(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Stationary linear iteration. (tolerances are squared)
int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, int m, double &tol, double atol, int printit)
GMRES method. (tolerances are squared)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.