55 if (dot_prod_type == 0)
61 double local_dot = (x * y);
64 MPI_Allreduce(&local_dot, &global_dot, 1, MPI_DOUBLE, MPI_SUM, comm);
76 if (dot_prod_type == 0)
83 MPI_Comm_rank(comm, &rank);
155 double r0, nom, nom0, nomold = 1, cf;
171 nom0 = nom =
Dot(
z,
r);
175 nom0 = nom =
Dot(
r,
r);
179 cout <<
" Iteration : " << setw(3) << 0 <<
" (B r, r) = "
218 cf = sqrt(nom/nomold);
220 cout <<
" Iteration : " << setw(3) << i <<
" (B r, r) = "
221 << nom <<
"\tConv. rate: " << cf <<
'\n';
227 cout <<
"Number of SLI iterations: " << i <<
'\n'
228 <<
"Conv. rate: " << cf <<
'\n';
230 cout <<
"(B r_0, r_0) = " << nom0 <<
'\n'
231 <<
"(B r_N, r_N) = " << nom <<
'\n'
232 <<
"Number of SLI iterations: " << i <<
'\n';
246 cerr <<
"SLI: No convergence!" <<
'\n';
247 cout <<
"(B r_0, r_0) = " << nom0 <<
'\n'
248 <<
"(B r_N, r_N) = " << nom <<
'\n'
249 <<
"Number of SLI iterations: " <<
final_iter <<
'\n';
253 cout <<
"Average reduction factor = "
260 int print_iter,
int max_num_iter,
261 double RTOLERANCE,
double ATOLERANCE)
273 int print_iter,
int max_num_iter,
274 double RTOLERANCE,
double ATOLERANCE)
297 double r0, den, nom, nom0, betanom, alpha, beta;
319 nom0 = nom =
Dot(
d,
r);
322 cout <<
" Iteration : " << setw(3) << 0 <<
" (B r, r) = "
339 cout <<
"Negative denominator in step 0 of PCG: " << den <<
'\n';
370 cout <<
" Iteration : " << setw(3) << i <<
" (B r, r) = "
377 cout <<
"Number of PCG iterations: " << i <<
'\n';
380 cout <<
"(B r_0, r_0) = " << nom0 <<
'\n'
381 <<
"(B r_N, r_N) = " << betanom <<
'\n'
382 <<
"Number of PCG iterations: " << i <<
'\n';
407 cout <<
"PCG: The operator is not positive definite. (Ad, d) = "
414 cerr <<
"PCG: No convergence!" <<
'\n';
415 cout <<
"(B r_0, r_0) = " << nom0 <<
'\n'
416 <<
"(B r_N, r_N) = " << betanom <<
'\n'
417 <<
"Number of PCG iterations: " <<
final_iter <<
'\n';
421 cout <<
"Average reduction factor = "
422 << pow (betanom/nom0, 0.5/
final_iter) <<
'\n';
428 int print_iter,
int max_num_iter,
429 double RTOLERANCE,
double ATOLERANCE)
441 int print_iter,
int max_num_iter,
442 double RTOLERANCE,
double ATOLERANCE)
456 double &cs,
double &sn)
463 else if (fabs(dy) > fabs(dx))
465 double temp = dx / dy;
466 sn = 1.0 / sqrt( 1.0 + temp*temp );
471 double temp = dy / dx;
472 cs = 1.0 / sqrt( 1.0 + temp*temp );
479 double temp = cs * dx + sn * dy;
480 dy = -sn * dx + cs * dy;
490 for (
int i = k; i >= 0; i--)
493 for (
int j = i - 1; j >= 0; j--)
495 y(j) -= h(j,i) * y(i);
499 for (
int j = 0; j <= k; j++)
551 double beta =
Norm(r);
564 cout <<
" Pass : " << setw(2) << 1
565 <<
" Iteration : " << setw(3) << 0
566 <<
" ||B r|| = " << beta <<
'\n';
573 if (v[0] == NULL) { v[0] =
new Vector(n); }
574 v[0]->Set(1.0/beta, r);
575 s = 0.0; s(0) = beta;
577 for (i = 0; i <
m && j <=
max_iter; i++, j++)
589 for (k = 0; k <= i; k++)
591 H(k,i) =
Dot(w, *v[k]);
592 w.
Add(-H(k,i), *v[k]);
596 if (v[i+1] == NULL) { v[i+1] =
new Vector(n); }
597 v[i+1]->Set(1.0/H(i+1,i), w);
599 for (k = 0; k < i; k++)
608 resid = fabs(s(i+1));
610 cout <<
" Pass : " << setw(2) << (j-1)/
m+1
611 <<
" Iteration : " << setw(3) << j
612 <<
" ||B r|| = " << resid <<
'\n';
619 for (i = 0; i <=
m; i++)
630 cout <<
"Restarting..." <<
'\n';
650 for (i = 0; i <=
m; i++)
661 for (i = 0; i <=
m; i++)
688 double beta =
Norm(r);
701 cout <<
" Pass : " << setw(2) << 1
702 <<
" Iteration : " << setw(3) << 0
703 <<
" || r || = " << beta << endl;
707 for (i= 0; i<=
m; i++)
716 if (v[0] == NULL) { v[0] =
new Vector(b.
Size()); }
718 v[0] ->
Add (1.0/beta, r);
719 s = 0.0; s(0) = beta;
724 if (z[i] == NULL) { z[i] =
new Vector(b.
Size()); }
730 for (k = 0; k <= i; k++)
732 H(k,i) =
Dot( r, *v[k]);
733 r.Add(-H(k,i), (*v[k]));
737 if (v[i+1] == NULL) { v[i+1] =
new Vector(b.
Size()); }
739 v[i+1] ->
Add (1.0/H(i+1,i), r);
741 for (k = 0; k < i; k++)
750 double resid = fabs(s(i+1));
752 cout <<
" Pass : " << setw(2) << j
753 <<
" Iteration : " << setw(3) << i+1
754 <<
" || r || = " << resid << endl;
762 for (i= 0; i<=
m; i++)
764 if (v[i]) {
delete v[i]; }
765 if (z[i]) {
delete z[i]; }
773 cout <<
"Restarting..." << endl;
786 for (i= 0; i<=
m; i++)
788 if (v[i]) {
delete v[i]; }
789 if (z[i]) {
delete z[i]; }
797 for (i = 0; i <=
m; i++)
799 if (v[i]) {
delete v[i]; }
800 if (z[i]) {
delete z[i]; }
809 int &max_iter,
int m,
double &tol,
double atol,
int printit)
826 int print_iter,
int max_num_iter,
int m,
double rtol,
double atol)
828 GMRES(A, x, b, B, max_num_iter, m, rtol, atol, print_iter);
850 double resid, tol_goal;
851 double rho_1, rho_2=1.0, alpha=1.0, beta, omega=1.0;
867 cout <<
" Iteration : " << setw(3) << 0
868 <<
" ||r|| = " << resid <<
'\n';
872 if (resid <= tol_goal)
886 cout <<
" Iteration : " << setw(3) << i
887 <<
" ||r|| = " << resid <<
'\n';
899 beta = (rho_1/rho_2) * (alpha/omega);
915 if (resid < tol_goal)
919 cout <<
" Iteration : " << setw(3) << i
920 <<
" ||s|| = " << resid <<
'\n';
927 cout <<
" Iteration : " << setw(3) << i
928 <<
" ||s|| = " << resid;
947 cout <<
" ||r|| = " << resid <<
'\n';
949 if (resid < tol_goal)
971 int &max_iter,
double &tol,
double atol,
int printit)
987 int print_iter,
int max_num_iter,
double rtol,
double atol)
989 BiCGSTAB(A, x, b, B, max_num_iter, rtol, atol, print_iter);
1015 double beta, eta, gamma0, gamma1, sigma0, sigma1;
1016 double alpha, delta, rho1, rho2, rho3, norm_goal;
1036 eta = beta = sqrt(
Dot(*z,
v1));
1037 gamma0 = gamma1 = 1.;
1038 sigma0 = sigma1 = 0.;
1043 cout <<
"MINRES: iteration " << setw(3) << 0 <<
": ||r||_B = "
1046 if (eta <= norm_goal)
1067 delta = gamma1*alpha - gamma0*sigma1*beta;
1069 rho2 = sigma1*alpha + gamma0*gamma1*beta;
1079 rho1 = hypot(delta, beta);
1083 w0.
Set(1./rho1, *z);
1087 add(1./rho1, *z, -rho2/rho1,
w1,
w0);
1091 add(-rho3/rho1,
w0, -rho2/rho1,
w1,
w0);
1092 w0.
Add(1./rho1, *z);
1096 gamma1 = delta/rho1;
1098 x.
Add(gamma1*eta,
w0);
1106 cout <<
"MINRES: iteration " << setw(3) << it <<
": ||r||_B = "
1107 << fabs(eta) <<
'\n';
1109 if (fabs(eta) <= norm_goal)
1130 cout <<
"MINRES: number of iterations: " << it <<
'\n';
1133 cout <<
"MINRES: iteration " << setw(3) << it <<
": ||r||_B = "
1134 << fabs(eta) <<
'\n';
1144 eta = sqrt(
Dot(*z,
v1));
1145 cout <<
"MINRES: iteration " << setw(3) << it <<
": ||r||_B = "
1146 << eta <<
" (re-computed)" <<
'\n';
1151 cerr <<
"MINRES: No convergence!" <<
'\n';
1156 int max_it,
double rtol,
double atol)
1168 int print_it,
int max_it,
double rtol,
double atol)
1186 MFEM_ASSERT(
height ==
width,
"square Operator is required.");
1194 MFEM_ASSERT(
oper != NULL,
"the Operator is not set (use SetOperator).");
1195 MFEM_ASSERT(
prec != NULL,
"the Solver is not set (use SetSolver).");
1198 double norm, norm_goal;
1218 for (it = 0;
true; it++)
1221 cout <<
"Newton iteration " << setw(2) << it
1222 <<
" : ||r|| = " << norm <<
'\n';
1224 if (norm <= norm_goal)
1257 int m_max,
int m_min,
int m_step,
double cf,
1258 double &tol,
double &atol,
int printit)
1265 Vector s(m+1), cs(m+1), sn(m+1);
1272 double normb = w.Norml2();
1282 double beta = r.
Norml2();
1284 resid = beta / normb;
1286 if (resid * resid <= tol)
1288 tol = resid * resid;
1294 cout <<
" Pass : " << setw(2) << 1
1295 <<
" Iteration : " << setw(3) << 0
1296 <<
" (r, r) = " << beta*beta <<
'\n';
1298 tol *= (normb*normb);
1299 tol = (atol > tol) ? atol : tol;
1303 for (i= 0; i<=m; i++)
1310 while (j <= max_iter)
1313 v[0] ->
Add (1.0/beta, r);
1314 s = 0.0; s(0) = beta;
1318 for (i = 0; i < m && j <= max_iter; i++)
1323 for (k = 0; k <= i; k++)
1325 H(k,i) = w * (*v[k]);
1326 w.
Add(-H(k,i), (*v[k]));
1329 H(i+1,i) = w.Norml2();
1331 v[i+1] ->
Add (1.0/H(i+1,i), w);
1333 for (k = 0; k < i; k++)
1342 resid = fabs(s(i+1));
1344 cout <<
" Pass : " << setw(2) << j
1345 <<
" Iteration : " << setw(3) << i+1
1346 <<
" (r, r) = " << resid*resid <<
'\n';
1348 if ( resid*resid < tol)
1351 tol = resid * resid;
1353 for (i= 0; i<=m; i++)
1363 cout <<
"Restarting..." <<
'\n';
1372 if ( resid*resid < tol)
1374 tol = resid * resid;
1376 for (i= 0; i<=m; i++)
1385 if (m - m_step >= m_min)
1398 tol = resid * resid;
1399 for (i= 0; i<=m; i++)
1421 mfem_error(
"SLBQPOptimizer::SetPreconditioner() : "
1422 "not meaningful for this solver");
1427 mfem_error(
"SLBQPOptimizer::SetOperator() : "
1428 "not meaningful for this solver");
1434 cout <<
"SLBQP iteration " << it <<
": residual = " << r
1435 <<
", lambda = " << l <<
'\n';
1457 const double smin = 0.1;
1464 cout <<
"SLBQP bracketing phase" <<
'\n';
1468 r =
solve(l,xt,x,nclip);
1482 llow = l; rlow = r; l = l + dl;
1485 r =
solve(l,xt,x,nclip);
1488 while ((r < 0) && (nclip <
max_iter))
1492 if (s < smin) { s = smin; }
1497 r =
solve(l,xt,x,nclip);
1505 lupp = l; rupp = r; l = l - dl;
1508 r =
solve(l,xt,x,nclip);
1511 while ((r > 0) && (nclip <
max_iter))
1515 if (s < smin) { s = smin; }
1520 r =
solve(l,xt,x,nclip);
1533 cout <<
"SLBQP secant phase" <<
'\n';
1536 s = 1.0 - rlow/rupp; dl = dl/s; l = lupp - dl;
1539 r =
solve(l,xt,x,nclip);
1542 while ( (fabs(r) > tol) && (nclip <
max_iter) )
1548 lupp = l; rupp = r; s = 1.0 - rlow/rupp;
1549 dl = (lupp - llow)/s; l = lupp - dl;
1554 if (s < smin) { s = smin; }
1556 lnew = 0.75*llow + 0.25*l;
1557 if (lnew < l-dl) { lnew = l-dl; }
1558 lupp = l; rupp = r; l = lnew;
1559 s = (lupp - llow)/(lupp - l);
1567 llow = l; rlow = r; s = 1.0 - rlow/rupp;
1568 dl = (lupp - llow)/s; l = lupp - dl;
1573 if (s < smin) { s = smin; }
1575 lnew = 0.75*lupp + 0.25*l;
1576 if (lnew < l+dl) { lnew = l+dl; }
1577 llow = l; rlow = r; l = lnew;
1578 s = (lupp - llow)/(lupp - l);
1583 r =
solve(l,xt,x,nclip);
1592 cerr <<
"SLBQP not converged!" <<
'\n';
1602 cout <<
"SLBQP iterations = " << nclip <<
'\n';
1603 cout <<
"SLBQP lambda = " << l <<
'\n';
1604 cout <<
"SLBQP residual = " << r <<
'\n';
1608 #ifdef MFEM_USE_SUITESPARSE
1635 umfpack_di_free_numeric(&
Numeric);
1639 umfpack_dl_free_numeric(&
Numeric);
1646 MFEM_ABORT(
"not a SparseMatrix");
1656 MFEM_VERIFY(
width ==
height,
"not a square matrix");
1664 int status = umfpack_di_symbolic(
width,
width, Ap, Ai, Ax, &Symbolic,
1669 umfpack_di_report_status(
Control, status);
1671 " umfpack_di_symbolic() failed!");
1674 status = umfpack_di_numeric(Ap, Ai, Ax, Symbolic, &
Numeric,
1679 umfpack_di_report_status(
Control, status);
1681 " umfpack_di_numeric() failed!");
1683 umfpack_di_free_symbolic(&Symbolic);
1687 SuiteSparse_long status;
1691 AI =
new SuiteSparse_long[
width + 1];
1692 AJ =
new SuiteSparse_long[Ap[
width]];
1693 for (
int i = 0; i <=
width; i++)
1695 AI[i] = (SuiteSparse_long)(Ap[i]);
1697 for (
int i = 0; i < Ap[
width]; i++)
1699 AJ[i] = (SuiteSparse_long)(Ai[i]);
1702 status = umfpack_dl_symbolic(
width,
width,
AI,
AJ, Ax, &Symbolic,
1707 umfpack_dl_report_status(
Control, status);
1709 " umfpack_dl_symbolic() failed!");
1712 status = umfpack_dl_numeric(
AI,
AJ, Ax, Symbolic, &
Numeric,
1717 umfpack_dl_report_status(
Control, status);
1719 " umfpack_dl_numeric() failed!");
1721 umfpack_dl_free_symbolic(&Symbolic);
1728 mfem_error(
"UMFPackSolver::Mult : matrix is not set!"
1729 " Call SetOperator first!");
1736 umfpack_di_report_info(Control,
Info);
1739 umfpack_di_report_status(Control, status);
1740 mfem_error(
"UMFPackSolver::Mult : umfpack_di_solve() failed!");
1745 SuiteSparse_long status =
1748 umfpack_dl_report_info(Control,
Info);
1751 umfpack_dl_report_status(Control, status);
1752 mfem_error(
"UMFPackSolver::Mult : umfpack_dl_solve() failed!");
1760 mfem_error(
"UMFPackSolver::MultTranspose : matrix is not set!"
1761 " Call SetOperator first!");
1768 umfpack_di_report_info(Control,
Info);
1771 umfpack_di_report_status(Control, status);
1773 " umfpack_di_solve() failed!");
1778 SuiteSparse_long status =
1781 umfpack_dl_report_info(Control,
Info);
1784 umfpack_dl_report_status(Control, status);
1786 " umfpack_dl_solve() failed!");
1799 umfpack_di_free_numeric(&
Numeric);
1803 umfpack_dl_free_numeric(&
Numeric);
1808 #endif // MFEM_USE_SUITESPARSE
void ApplyPlaneRotation(double &dx, double &dy, double &cs, double &sn)
Conjugate gradient method.
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator.
double Info[UMFPACK_INFO]
int GetNumIterations() const
void SetSize(int s)
Resizes the vector if the new size is different.
virtual Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x.
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.
int Size() const
Returns the size of the vector.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application.
bool iterative_mode
If true, use the second argument of Mult as an initial guess.
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
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
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
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 &x, Vector &y) const
Operator application.
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.
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
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)
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.
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
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]
void subtract(const Vector &x, const Vector &y, Vector &z)
void mfem_error(const char *msg)
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)
Vector & Set(const double a, const Vector &x)
(*this) = a * x
void SetBounds(const Vector &_lo, const Vector &_hi)
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
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.
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.
int * GetJ() const
Return the array J.
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
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)
Set/update the solver for the given operator.