13 #include "../general/annotation.hpp"
14 #include "../general/forall.hpp"
15 #include "../general/globals.hpp"
16 #include "../fem/bilinearform.hpp"
53 #endif // MFEM_USE_MPI
60 if (dot_prod_type == 0)
74 int print_level_ = print_lvl;
77 if (dot_prod_type != 0)
80 MPI_Comm_rank(comm, &rank);
99 if (dot_prod_type != 0)
102 MPI_Comm_rank(comm, &rank);
105 derived_print_level = -1;
119 if (comm != MPI_COMM_NULL)
121 MPI_Comm_rank(comm, &rank);
125 switch (print_level_)
142 MFEM_WARNING(
"Unknown print level " << print_level_ <<
143 ". Defaulting to level 0.");
159 else if (print_options_.
summary)
191 const Vector& x,
bool final)
const
202 ess_tdof_list(nullptr),
211 Solver(a.FESpace()->GetTrueVSize()),
214 ess_tdof_list(&ess_tdofs),
233 ess_tdof_list(&ess_tdofs),
266 MFEM_ASSERT(
height ==
width,
"not a square matrix!");
268 ess_tdof_list =
nullptr;
280 const double delta = damping;
281 auto D = diag.
Read();
282 auto DI = dinv.
Write();
283 const bool use_abs_diag_ = use_abs_diag;
288 MFEM_ABORT_KERNEL(
"Zero diagonal entry in OperatorJacobiSmoother");
290 if (!use_abs_diag_) { DI[i] = delta / D[i]; }
291 else { DI[i] = delta / std::abs(D[i]); }
293 if (ess_tdof_list && ess_tdof_list->
Size() > 0)
295 auto I = ess_tdof_list->
Read();
296 MFEM_FORALL(i, ess_tdof_list->
Size(), DI[I[i]] =
delta; );
304 MFEM_ASSERT(x.
Size() ==
Width(),
"invalid input vector");
305 MFEM_ASSERT(y.
Size() ==
Height(),
"invalid output vector");
309 MFEM_VERIFY(oper,
"iterative_mode == true requires the forward operator");
310 oper->
Mult(y, residual);
319 auto DI = dinv.
Read();
320 auto R = residual.
Read();
324 Y[i] += DI[i] * R[i];
331 int order_,
double max_eig_estimate_)
335 max_eig_estimate(max_eig_estimate_),
340 ess_tdof_list(ess_tdofs),
342 oper(&oper_) {
Setup(); }
348 int order_, MPI_Comm comm,
int power_iterations,
double power_tolerance)
353 int order_,
int power_iterations,
double power_tolerance)
361 ess_tdof_list(ess_tdofs),
375 power_iterations, power_tolerance);
383 int order_,
double max_eig_estimate_)
390 int order_, MPI_Comm comm,
int power_iterations,
double power_tolerance)
392 power_iterations, power_tolerance) { }
397 int order_,
int power_iterations,
double power_tolerance)
406 auto D = diag.
Read();
407 auto X = dinv.
Write();
408 MFEM_FORALL(i, N, X[i] = 1.0 / D[i]; );
409 auto I = ess_tdof_list.
Read();
410 MFEM_FORALL(i, ess_tdof_list.
Size(), X[I[i]] = 1.0; );
415 double upper_bound = 1.2 * max_eig_estimate;
416 double lower_bound = 0.3 * max_eig_estimate;
417 double theta = 0.5 * (upper_bound + lower_bound);
418 double delta = 0.5 * (upper_bound - lower_bound);
424 coeffs[0] = 1.0 / theta;
429 double tmp_0 = 1.0/(
pow(delta, 2) - 2*
pow(theta, 2));
430 coeffs[0] = -4*theta*tmp_0;
436 double tmp_0 = 3*
pow(delta, 2);
437 double tmp_1 =
pow(theta, 2);
438 double tmp_2 = 1.0/(-4*
pow(theta, 3) + theta*tmp_0);
439 coeffs[0] = tmp_2*(tmp_0 - 12*tmp_1);
440 coeffs[1] = 12/(tmp_0 - 4*tmp_1);
441 coeffs[2] = -4*tmp_2;
446 double tmp_0 =
pow(delta, 2);
447 double tmp_1 =
pow(theta, 2);
448 double tmp_2 = 8*tmp_0;
449 double tmp_3 = 1.0/(
pow(delta, 4) + 8*
pow(theta, 4) - tmp_1*tmp_2);
450 coeffs[0] = tmp_3*(32*
pow(theta, 3) - 16*theta*tmp_0);
451 coeffs[1] = tmp_3*(-48*tmp_1 + tmp_2);
452 coeffs[2] = 32*theta*tmp_3;
453 coeffs[3] = -8*tmp_3;
458 double tmp_0 = 5*
pow(delta, 4);
459 double tmp_1 =
pow(theta, 4);
460 double tmp_2 =
pow(theta, 2);
461 double tmp_3 =
pow(delta, 2);
462 double tmp_4 = 60*tmp_3;
463 double tmp_5 = 20*tmp_3;
464 double tmp_6 = 1.0/(16*
pow(theta, 5) -
pow(theta, 3)*tmp_5 + theta*tmp_0);
465 double tmp_7 = 160*tmp_2;
466 double tmp_8 = 1.0/(tmp_0 + 16*tmp_1 - tmp_2*tmp_5);
467 coeffs[0] = tmp_6*(tmp_0 + 80*tmp_1 - tmp_2*tmp_4);
468 coeffs[1] = tmp_8*(tmp_4 - tmp_7);
469 coeffs[2] = tmp_6*(-tmp_5 + tmp_7);
470 coeffs[3] = -80*tmp_8;
471 coeffs[4] = 16*tmp_6;
475 MFEM_ABORT(
"Chebyshev smoother not implemented for order = " << order);
483 MFEM_ABORT(
"Chebyshev smoother not implemented for iterative mode");
488 MFEM_ABORT(
"Chebyshev smoother requires operator");
497 for (
int k = 0; k < order; ++k)
502 oper->
Mult(residual, helperVector);
503 residual = helperVector;
508 auto Dinv = dinv.
Read();
510 MFEM_FORALL(i, n, R[i] *= Dinv[i]; );
514 auto C = coeffs.
Read();
515 MFEM_FORALL(i, n, Y[i] += C[k] * R[i]; );
564 double r0, nom, nom0, nomold = 1, cf;
589 mfem::out <<
" Iteration : " << setw(3) << right << 0 <<
" ||Br|| = "
647 mfem::out <<
" Iteration : " << setw(3) << right << (i-1)
648 <<
" ||Br|| = " << setw(11) << left << nom
649 <<
"\tConv. rate: " << cf <<
'\n';
659 <<
"Conv. rate: " << cf <<
'\n'
660 <<
"Average reduction factor: "<< rf <<
'\n';
664 mfem::out <<
"SLI: No convergence!" <<
'\n';
671 int print_iter,
int max_num_iter,
672 double RTOLERANCE,
double ATOLERANCE)
686 int print_iter,
int max_num_iter,
687 double RTOLERANCE,
double ATOLERANCE)
714 double r0, den, nom, nom0, betanom,
alpha, beta;
737 nom0 = nom =
Dot(d,
r);
738 MFEM_ASSERT(
IsFinite(nom),
"nom = " << nom);
741 mfem::out <<
" Iteration : " << setw(3) << 0 <<
" (B r, r) = "
750 mfem::out <<
"PCG: The preconditioner is not positive definite. (Br, r) = "
769 MFEM_ASSERT(
IsFinite(den),
"den = " << den);
774 mfem::out <<
"PCG: The operator is not positive definite. (Ad, d) = "
804 MFEM_ASSERT(
IsFinite(betanom),
"betanom = " << betanom);
809 mfem::out <<
"PCG: The preconditioner is not positive definite. (Br, r) = "
819 mfem::out <<
" Iteration : " << setw(3) << i <<
" (B r, r) = "
848 MFEM_ASSERT(
IsFinite(den),
"den = " << den);
853 mfem::out <<
"PCG: The operator is not positive definite. (Ad, d) = "
877 mfem::out <<
"Average reduction factor = " << arf <<
'\n';
881 mfem::out <<
"PCG: No convergence!" <<
'\n';
890 int print_iter,
int max_num_iter,
891 double RTOLERANCE,
double ATOLERANCE)
905 int print_iter,
int max_num_iter,
906 double RTOLERANCE,
double ATOLERANCE)
922 double &cs,
double &sn)
929 else if (fabs(dy) > fabs(dx))
931 double temp = dx / dy;
932 sn = 1.0 /
sqrt( 1.0 + temp*temp );
937 double temp = dy / dx;
938 cs = 1.0 /
sqrt( 1.0 + temp*temp );
945 double temp = cs * dx + sn * dy;
946 dy = -sn * dx + cs * dy;
956 for (
int i = k; i >= 0; i--)
959 for (
int j = i - 1; j >= 0; j--)
961 y(j) -= h(j,i) * y(i);
965 for (
int j = 0; j <= k; j++)
1018 double beta =
Norm(r);
1019 MFEM_ASSERT(
IsFinite(beta),
"beta = " << beta);
1034 <<
" Iteration : " << setw(3) << 0
1035 <<
" ||B r|| = " << beta
1045 if (v[0] == NULL) { v[0] =
new Vector(n); }
1046 v[0]->Set(1.0/beta, r);
1047 s = 0.0;
s(0) = beta;
1049 for (i = 0; i <
m && j <=
max_iter; i++, j++)
1061 for (k = 0; k <= i; k++)
1063 H(k,i) =
Dot(w, *v[k]);
1064 w.
Add(-H(k,i), *v[k]);
1068 MFEM_ASSERT(
IsFinite(H(i+1,i)),
"Norm(w) = " << H(i+1,i));
1069 if (v[i+1] == NULL) { v[i+1] =
new Vector(n); }
1070 v[i+1]->Set(1.0/H(i+1,i), w);
1072 for (k = 0; k < i; k++)
1081 resid = fabs(
s(i+1));
1082 MFEM_ASSERT(
IsFinite(resid),
"resid = " << resid);
1095 mfem::out <<
" Pass : " << setw(2) << (j-1)/
m+1
1096 <<
" Iteration : " << setw(3) << j
1097 <<
" ||B r|| = " << resid <<
'\n';
1121 MFEM_ASSERT(
IsFinite(beta),
"beta = " << beta);
1138 mfem::out <<
" Pass : " << setw(2) << (j-1)/
m+1
1140 <<
" ||B r|| = " << resid <<
'\n';
1148 mfem::out <<
"GMRES: No convergence!\n";
1153 for (i = 0; i < v.
Size(); i++)
1178 double beta =
Norm(r);
1182 MFEM_ASSERT(
IsFinite(beta),
"beta = " << beta);
1197 <<
" Iteration : " << setw(3) << 0
1198 <<
" || r || = " << beta
1206 for (i= 0; i<=
m; i++)
1215 if (v[0] == NULL) { v[0] =
new Vector(b.
Size()); }
1217 v[0] ->
Add (1.0/beta, r);
1218 s = 0.0;
s(0) = beta;
1220 for (i = 0; i <
m && j <=
max_iter; i++, j++)
1223 if (z[i] == NULL) { z[i] =
new Vector(b.
Size()); }
1236 for (k = 0; k <= i; k++)
1238 H(k,i) =
Dot( r, *v[k]);
1239 r.Add(-H(k,i), (*v[k]));
1243 if (v[i+1] == NULL) { v[i+1] =
new Vector(b.
Size()); }
1245 v[i+1] ->
Add (1.0/H(i+1,i), r);
1247 for (k = 0; k < i; k++)
1256 resid = fabs(
s(i+1));
1257 MFEM_ASSERT(
IsFinite(resid),
"resid = " << resid);
1261 mfem::out <<
" Pass : " << setw(2) << (j-1)/
m+1
1262 <<
" Iteration : " << setw(3) << j
1263 <<
" || r || = " << resid << endl;
1279 for (i= 0; i<=
m; i++)
1281 if (v[i]) {
delete v[i]; }
1282 if (z[i]) {
delete z[i]; }
1298 MFEM_ASSERT(
IsFinite(beta),
"beta = " << beta);
1310 for (i = 0; i <=
m; i++)
1312 if (v[i]) {
delete v[i]; }
1313 if (z[i]) {
delete z[i]; }
1320 mfem::out <<
" Pass : " << setw(2) << (j-1)/
m+1
1321 <<
" Iteration : " << setw(3) << j-1
1322 <<
" || r || = " << resid << endl;
1326 mfem::out <<
"FGMRES: Number of iterations: " << j-1 <<
'\n';
1330 mfem::out <<
"FGMRES: No convergence!\n";
1336 int &max_iter,
int m,
double &tol,
double atol,
int printit)
1355 int print_iter,
int max_num_iter,
int m,
double rtol,
double atol)
1357 GMRES(A, x, b, B, max_num_iter, m, rtol, atol, print_iter);
1379 double resid, tol_goal;
1380 double rho_1, rho_2=1.0,
alpha=1.0, beta,
omega=1.0;
1395 MFEM_ASSERT(
IsFinite(resid),
"resid = " << resid);
1398 mfem::out <<
" Iteration : " << setw(3) << 0
1406 if (resid <= tol_goal)
1421 mfem::out <<
" Iteration : " << setw(3) << i
1422 <<
" ||r|| = " << resid <<
'\n';
1436 mfem::out <<
"BiCGStab: No convergence!\n";
1446 beta = (rho_1/rho_2) * (
alpha/omega);
1462 MFEM_ASSERT(
IsFinite(resid),
"resid = " << resid);
1463 if (resid < tol_goal)
1468 mfem::out <<
" Iteration : " << setw(3) << i
1469 <<
" ||s|| = " << resid <<
'\n';
1482 mfem::out <<
" Iteration : " << setw(3) << i
1483 <<
" ||s|| = " << resid;
1502 MFEM_ASSERT(
IsFinite(resid),
"resid = " << resid);
1505 mfem::out <<
" ||r|| = " << resid <<
'\n';
1508 if (resid < tol_goal)
1515 mfem::out <<
" Iteration : " << setw(3) << i
1516 <<
" ||r|| = " << resid <<
'\n';
1531 mfem::out <<
" Iteration : " << setw(3) << i
1532 <<
" ||r|| = " << resid <<
'\n';
1540 mfem::out <<
"BiCGStab: No convergence!\n";
1553 <<
" ||r|| = " << resid <<
'\n';
1561 mfem::out <<
"BiCGStab: No convergence!\n";
1566 int &max_iter,
double &tol,
double atol,
int printit)
1575 bicgstab.
Mult(b, x);
1582 int print_iter,
int max_num_iter,
double rtol,
double atol)
1584 BiCGSTAB(A, x, b, B, max_num_iter, rtol, atol, print_iter);
1610 double beta, eta, gamma0, gamma1, sigma0, sigma1;
1611 double alpha,
delta, rho1, rho2, rho3, norm_goal;
1632 MFEM_ASSERT(
IsFinite(eta),
"eta = " << eta);
1633 gamma0 = gamma1 = 1.;
1634 sigma0 = sigma1 = 0.;
1638 if (eta <= norm_goal)
1646 mfem::out <<
"MINRES: iteration " << setw(3) << 0 <<
": ||r||_B = "
1660 MFEM_ASSERT(
IsFinite(alpha),
"alpha = " << alpha);
1667 delta = gamma1*alpha - gamma0*sigma1*beta;
1669 rho2 = sigma1*alpha + gamma0*gamma1*beta;
1679 MFEM_ASSERT(
IsFinite(beta),
"beta = " << beta);
1680 rho1 = hypot(delta, beta);
1684 w0.
Set(1./rho1, *z);
1688 add(1./rho1, *z, -rho2/rho1,
w1,
w0);
1692 add(-rho3/rho1,
w0, -rho2/rho1,
w1,
w0);
1693 w0.
Add(1./rho1, *z);
1697 gamma1 = delta/rho1;
1699 x.
Add(gamma1*eta,
w0);
1705 MFEM_ASSERT(
IsFinite(eta),
"eta = " << eta);
1707 if (fabs(eta) <= norm_goal)
1714 mfem::out <<
"MINRES: iteration " << setw(3) << it <<
": ||r||_B = "
1715 << fabs(eta) <<
'\n';
1717 Monitor(it, fabs(eta), *z, x);
1735 mfem::out <<
"MINRES: iteration " << setw(3) << it <<
": ||r||_B = "
1736 << fabs(eta) <<
'\n';
1761 mfem::out <<
"MINRES: No convergence!\n";
1766 int max_it,
double rtol,
double atol)
1780 int print_it,
int max_it,
double rtol,
double atol)
1798 MFEM_ASSERT(
height ==
width,
"square Operator is required.");
1806 MFEM_ASSERT(
oper != NULL,
"the Operator is not set (use SetOperator).");
1807 MFEM_ASSERT(
prec != NULL,
"the Solver is not set (use SetSolver).");
1810 double norm0, norm, norm_goal;
1826 norm0 = norm =
Norm(
r);
1829 mfem::out <<
"Newton iteration " << setw(2) << 0
1830 <<
" : ||r|| = " << norm <<
"...\n";
1837 for (it = 0;
true; it++)
1839 MFEM_ASSERT(
IsFinite(norm),
"norm = " << norm);
1842 mfem::out <<
"Newton iteration " << setw(2) << it
1843 <<
" : ||r|| = " << norm;
1846 mfem::out <<
", ||r||/||r_0|| = " << norm/norm0;
1852 if (norm <= norm_goal)
1885 add(x, -c_scale,
c, x);
1908 mfem::out <<
"Newton: No convergence!\n";
1914 const double rtol_max,
1915 const double alpha_,
1916 const double gamma_)
1921 this->
alpha = alpha_;
1922 this->
gamma = gamma_;
1927 const double fnorm)
const
1935 double sg_threshold = 0.1;
1955 MFEM_ABORT(
"Unknown adaptive linear solver rtol version");
1960 if (sg_eta > sg_threshold) { eta = std::max(eta, sg_eta); }
1964 iterative_solver->SetRelTol(eta);
1968 mfem::out <<
"Eisenstat-Walker rtol = " << eta <<
"\n";
1975 const double fnorm)
const
1994 MFEM_VERIFY(
oper != NULL,
"the Operator is not set (use SetOperator).");
2005 int last_saved_id = -1;
2008 double norm0, norm, norm_goal;
2020 if (have_b) {
r -=
b; }
2024 norm0 = norm =
Norm(
r);
2027 mfem::out <<
"LBFGS iteration " << setw(2) << 0
2028 <<
" : ||r|| = " << norm <<
"...\n";
2031 for (it = 0;
true; it++)
2033 MFEM_ASSERT(
IsFinite(norm),
"norm = " << norm);
2037 <<
" : ||r|| = " << norm;
2040 mfem::out <<
", ||r||/||r_0|| = " << norm/norm0;
2045 if (norm <= norm_goal)
2064 add(x, -c_scale,
c, x);
2076 sk =
c; sk *= -c_scale;
2080 last_saved_id = (last_saved_id ==
m-1) ? 0 : last_saved_id+1;
2085 for (
int i = last_saved_id; i > -1; i--)
2093 for (
int i =
m-1; i > last_saved_id; i--)
2104 for (
int i = last_saved_id+1; i <
m ; i++)
2110 for (
int i = 0; i < last_saved_id+1 ; i++)
2130 mfem::out <<
"LBFGS: No convergence!\n";
2136 int m_max,
int m_min,
int m_step,
double cf,
2137 double &tol,
double &atol,
int printit)
2144 Vector s(m+1), cs(m+1), sn(m+1);
2151 double normb = w.Norml2();
2161 double beta = r.
Norml2();
2163 resid = beta / normb;
2165 if (resid * resid <= tol)
2167 tol = resid * resid;
2175 <<
" Iteration : " << setw(3) << 0
2176 <<
" (r, r) = " << beta*beta <<
'\n';
2179 tol *= (normb*normb);
2180 tol = (atol > tol) ? atol : tol;
2184 for (i= 0; i<=m; i++)
2191 while (j <= max_iter)
2194 v[0] ->
Add (1.0/beta, r);
2195 s = 0.0;
s(0) = beta;
2199 for (i = 0; i < m && j <= max_iter; i++)
2204 for (k = 0; k <= i; k++)
2206 H(k,i) = w * (*v[k]);
2207 w.
Add(-H(k,i), (*v[k]));
2210 H(i+1,i) = w.Norml2();
2212 v[i+1] ->
Add (1.0/H(i+1,i), w);
2214 for (k = 0; k < i; k++)
2223 resid = fabs(
s(i+1));
2227 <<
" Iteration : " << setw(3) << i+1
2228 <<
" (r, r) = " << resid*resid <<
'\n';
2231 if ( resid*resid < tol)
2234 tol = resid * resid;
2236 for (i= 0; i<=m; i++)
2255 if ( resid*resid < tol)
2257 tol = resid * resid;
2259 for (i= 0; i<=m; i++)
2268 if (m - m_step >= m_min)
2281 tol = resid * resid;
2282 for (i= 0; i<=m; i++)
2292 : C(C_), D(D_), c_e(NULL), d_lo(NULL), d_hi(NULL), x_lo(NULL), x_hi(NULL),
2301 MFEM_ASSERT(
C,
"The C operator is unspecified -- can't set constraints.");
2302 MFEM_ASSERT(c.
Size() ==
C->
Height(),
"Wrong size of the constraint.");
2310 MFEM_ASSERT(
D,
"The D operator is unspecified -- can't set constraints.");
2312 "Wrong size of the constraint.");
2320 "Wrong size of the constraint.");
2337 MFEM_WARNING(
"Objective functional is ignored as SLBQP always minimizes"
2338 "the l2 norm of (x - x_target).");
2340 MFEM_ASSERT(prob.
GetC(),
"Linear constraint is not set.");
2341 MFEM_ASSERT(prob.
GetC()->
Height() == 1,
"Solver expects scalar constraint.");
2362 mfem::out <<
"SLBQP iteration " << it <<
": residual = " << r
2363 <<
", lambda = " << l <<
'\n';
2386 const double smin = 0.1;
2393 mfem::out <<
"SLBQP bracketing phase" <<
'\n';
2397 r =
solve(l,xt,x,nclip);
2411 llow = l; rlow = r; l = l + dl;
2414 r =
solve(l,xt,x,nclip);
2417 while ((r < 0) && (nclip <
max_iter))
2421 if (s < smin) { s = smin; }
2426 r =
solve(l,xt,x,nclip);
2434 lupp = l; rupp = r; l = l - dl;
2437 r =
solve(l,xt,x,nclip);
2440 while ((r > 0) && (nclip <
max_iter))
2444 if (s < smin) { s = smin; }
2449 r =
solve(l,xt,x,nclip);
2462 mfem::out <<
"SLBQP secant phase" <<
'\n';
2465 s = 1.0 - rlow/rupp; dl = dl/
s; l = lupp - dl;
2468 r =
solve(l,xt,x,nclip);
2471 while ( (fabs(r) > tol) && (nclip <
max_iter) )
2477 lupp = l; rupp = r; s = 1.0 - rlow/rupp;
2478 dl = (lupp - llow)/s; l = lupp - dl;
2483 if (s < smin) { s = smin; }
2485 lnew = 0.75*llow + 0.25*l;
2486 if (lnew < l-dl) { lnew = l-dl; }
2487 lupp = l; rupp = r; l = lnew;
2488 s = (lupp - llow)/(lupp - l);
2496 llow = l; rlow = r; s = 1.0 - rlow/rupp;
2497 dl = (lupp - llow)/s; l = lupp - dl;
2502 if (s < smin) { s = smin; }
2504 lnew = 0.75*lupp + 0.25*l;
2505 if (lnew < l+dl) { lnew = l+dl; }
2506 llow = l; rlow = r; l = lnew;
2507 s = (lupp - llow)/(lupp - l);
2512 r =
solve(l,xt,x,nclip);
2528 <<
" lambda = " << l <<
'\n'
2533 mfem::out <<
"SLBQP: No convergence!" <<
'\n';
2537 struct WeightMinHeap
2539 const std::vector<double> &w;
2540 std::vector<size_t> c;
2541 std::vector<int> loc;
2543 WeightMinHeap(
const std::vector<double> &w_) : w(w_)
2545 c.reserve(w.size());
2546 loc.resize(w.size());
2547 for (
size_t i=0; i<w.size(); ++i) { push(i); }
2550 size_t percolate_up(
size_t pos,
double val)
2552 for (; pos > 0 && w[c[(pos-1)/2]] > val; pos = (pos-1)/2)
2554 c[pos] = c[(pos-1)/2];
2555 loc[c[(pos-1)/2]] = pos;
2560 size_t percolate_down(
size_t pos,
double val)
2562 while (2*pos+1 < c.size())
2564 size_t left = 2*pos+1;
2565 size_t right = left+1;
2567 if (right < c.size() && w[c[right]] < w[c[left]]) { tgt = right; }
2568 else { tgt = left; }
2569 if (w[c[tgt]] < val)
2587 size_t pos = c.size()-1;
2588 pos = percolate_up(pos, val);
2596 size_t j = c.back();
2600 if (c.empty()) {
return i; }
2603 pos = percolate_down(pos, val);
2609 void update(
size_t i)
2611 size_t pos = loc[i];
2613 pos = percolate_up(pos, val);
2614 pos = percolate_down(pos, val);
2619 bool picked(
size_t i)
2634 for (
int i=0; i<n; ++i)
2636 for (
int j=I[i]; j<I[i+1]; ++j)
2638 V[j] = abs(V[j]/D[i]);
2642 std::vector<double> w(n, 0.0);
2643 for (
int k=0; k<n; ++k)
2646 for (
int ii=I[k]; ii<I[k+1]; ++ii)
2651 for (
int kk=I[i]; kk<I[i+1]; ++kk)
2659 for (
int jj=I[k]; jj<I[k+1]; ++jj)
2662 if (j == k) {
continue; }
2663 double C_kj = V[jj];
2664 bool ij_exists =
false;
2665 for (
int jj2=I[i]; jj2<I[i+1]; ++jj2)
2673 if (!ij_exists) { w[k] +=
pow(C_ik*C_kj,2); }
2679 WeightMinHeap w_heap(w);
2683 for (
int ii=0; ii<n; ++ii)
2685 int pi = w_heap.pop();
2688 for (
int kk=I[pi]; kk<I[pi+1]; ++kk)
2691 if (w_heap.picked(k)) {
continue; }
2695 for (
int ii2=I[k]; ii2<I[k+1]; ++ii2)
2698 if (w_heap.picked(i)) {
continue; }
2701 for (
int kk2=I[i]; kk2<I[i+1]; ++kk2)
2709 for (
int jj=I[k]; jj<I[k+1]; ++jj)
2712 if (j == k || w_heap.picked(j)) {
continue; }
2713 double C_kj = V[jj];
2714 bool ij_exists =
false;
2715 for (
int jj2=I[i]; jj2<I[i+1]; ++jj2)
2723 if (!ij_exists) { w[k] +=
pow(C_ik*C_kj,2); }
2736 block_size(block_size_),
2738 reordering(reordering_)
2745 :
BlockILU(block_size_, reordering_, k_fill_)
2767 MFEM_ABORT(
"BlockILU must be created with a SparseMatrix or HypreParMatrix");
2772 MFEM_ASSERT(A->
Finalized(),
"Matrix must be finalized.");
2773 CreateBlockPattern(*A);
2777 void BlockILU::CreateBlockPattern(
const SparseMatrix &A)
2779 MFEM_VERIFY(k_fill == 0,
"Only block ILU(0) is currently supported.");
2780 if (A.
Height() % block_size != 0)
2782 MFEM_ABORT(
"BlockILU: block size must evenly divide the matrix size");
2786 const int *I = A.
GetI();
2787 const int *J = A.
GetJ();
2788 const double *V = A.
GetData();
2790 int nblockrows = nrows / block_size;
2792 std::vector<std::set<int>> unique_block_cols(nblockrows);
2794 for (
int iblock = 0; iblock < nblockrows; ++iblock)
2796 for (
int bi = 0; bi < block_size; ++bi)
2798 int i = iblock * block_size + bi;
2799 for (
int k = I[i]; k < I[i + 1]; ++k)
2801 unique_block_cols[iblock].insert(J[k] / block_size);
2804 nnz += unique_block_cols[iblock].size();
2809 SparseMatrix C(nblockrows, nblockrows);
2810 for (
int iblock = 0; iblock < nblockrows; ++iblock)
2812 for (
int jblock : unique_block_cols[iblock])
2814 for (
int bi = 0; bi < block_size; ++bi)
2816 int i = iblock * block_size + bi;
2817 for (
int k = I[i]; k < I[i + 1]; ++k)
2820 if (j >= jblock * block_size && j < (jblock + 1) * block_size)
2822 C.Add(iblock, jblock, V[k]*V[k]);
2829 double *CV = C.GetData();
2830 for (
int i=0; i<C.NumNonZeroElems(); ++i)
2832 CV[i] =
sqrt(CV[i]);
2841 MFEM_ABORT(
"BlockILU: unknown reordering")
2848 for (
int i=0; i<nblockrows; ++i)
2856 for (
int i=0; i<nblockrows; ++i)
2862 std::vector<std::vector<int>> unique_block_cols_perminv(nblockrows);
2863 for (
int i=0; i<nblockrows; ++i)
2865 std::vector<int> &cols = unique_block_cols_perminv[i];
2866 for (
int j : unique_block_cols[P[i]])
2868 cols.push_back(Pinv[j]);
2870 std::sort(cols.begin(), cols.end());
2877 AB.
SetSize(block_size, block_size, nnz);
2878 DB.
SetSize(block_size, block_size, nblockrows);
2881 ipiv.
SetSize(block_size*nblockrows);
2884 for (
int iblock = 0; iblock < nblockrows; ++iblock)
2886 int iblock_perm = P[iblock];
2887 for (
int jblock : unique_block_cols_perminv[iblock])
2889 int jblock_perm = P[jblock];
2890 if (iblock == jblock)
2892 ID[iblock] = counter;
2894 JB[counter] = jblock;
2895 for (
int bi = 0; bi < block_size; ++bi)
2897 int i = iblock_perm*block_size + bi;
2898 for (
int k = I[i]; k < I[i + 1]; ++k)
2901 if (j >= jblock_perm*block_size && j < (jblock_perm + 1)*block_size)
2903 int bj = j - jblock_perm*block_size;
2905 AB(bi, bj, counter) = val;
2907 if (iblock == jblock)
2909 DB(bi, bj, iblock) = val;
2916 IB[iblock + 1] = counter;
2920 void BlockILU::Factorize()
2922 int nblockrows =
Height()/block_size;
2925 for (
int i=0; i<nblockrows; ++i)
2927 LUFactors factorization(DB.
GetData(i), &ipiv[i*block_size]);
2928 factorization.Factor(block_size);
2934 DenseMatrix A_ik, A_ij, A_kj;
2936 for (
int i=1; i<nblockrows; ++i)
2939 for (
int kk=IB[i]; kk<IB[i+1]; ++kk)
2943 if (k == i) {
break; }
2946 MFEM_ABORT(
"Matrix must be sorted with nonzero diagonal");
2948 LUFactors A_kk_inv(DB.
GetData(k), &ipiv[k*block_size]);
2949 A_ik.UseExternalData(&AB(0,0,kk), block_size, block_size);
2951 A_kk_inv.RightSolve(block_size, block_size, A_ik.GetData());
2953 for (
int jj=kk+1; jj<IB[i+1]; ++jj)
2956 if (j <= k) {
continue; }
2957 A_ij.UseExternalData(&AB(0,0,jj), block_size, block_size);
2958 for (
int ll=IB[k]; ll<IB[k+1]; ++ll)
2963 A_kj.UseExternalData(&AB(0,0,ll), block_size, block_size);
2970 LUFactors factorization(DB.
GetData(i), &ipiv[i*block_size]);
2971 factorization.Factor(block_size);
2983 MFEM_ASSERT(
height > 0,
"BlockILU(0) preconditioner is not constructed");
2984 int nblockrows =
Height()/block_size;
2993 for (
int i=0; i<nblockrows; ++i)
2996 for (
int ib=0; ib<block_size; ++ib)
2998 yi[ib] = b[ib + P[i]*block_size];
3000 for (
int k=IB[i]; k<ID[i]; ++k)
3010 for (
int i=nblockrows-1; i >= 0; --i)
3013 for (
int ib=0; ib<block_size; ++ib)
3015 xi[ib] = y[ib + i*block_size];
3017 for (
int k=ID[i]+1; k<IB[i+1]; ++k)
3025 LUFactors A_ii_inv(&DB(0,0,i), &ipiv[i*block_size]);
3027 A_ii_inv.
Solve(block_size, 1, xi);
3033 int it,
double norm,
const Vector &r,
bool final)
3037 double bc_norm_squared = 0.0;
3043 bc_norm_squared += r_entry*r_entry;
3048 if (comm != MPI_COMM_NULL)
3050 double glob_bc_norm_squared = 0.0;
3051 MPI_Reduce(&bc_norm_squared, &glob_bc_norm_squared, 1, MPI_DOUBLE,
3053 bc_norm_squared = glob_bc_norm_squared;
3055 MPI_Comm_rank(comm, &rank);
3056 print = (rank == 0);
3059 if ((it == 0 ||
final || bc_norm_squared > 0.0) && print)
3061 mfem::out <<
" ResidualBCMonitor : b.c. residual norm = "
3062 <<
sqrt(bc_norm_squared) << endl;
3067 #ifdef MFEM_USE_SUITESPARSE
3092 umfpack_di_free_numeric(&
Numeric);
3096 umfpack_dl_free_numeric(&
Numeric);
3101 MFEM_VERIFY(
mat,
"not a SparseMatrix");
3110 MFEM_VERIFY(
width ==
height,
"not a square matrix");
3118 int status = umfpack_di_symbolic(
width,
width, Ap, Ai, Ax, &Symbolic,
3123 umfpack_di_report_status(
Control, status);
3125 " umfpack_di_symbolic() failed!");
3128 status = umfpack_di_numeric(Ap, Ai, Ax, Symbolic, &
Numeric,
3133 umfpack_di_report_status(
Control, status);
3135 " umfpack_di_numeric() failed!");
3137 umfpack_di_free_symbolic(&Symbolic);
3141 SuiteSparse_long status;
3145 AI =
new SuiteSparse_long[
width + 1];
3146 AJ =
new SuiteSparse_long[Ap[
width]];
3147 for (
int i = 0; i <=
width; i++)
3149 AI[i] = (SuiteSparse_long)(Ap[i]);
3151 for (
int i = 0; i < Ap[
width]; i++)
3153 AJ[i] = (SuiteSparse_long)(Ai[i]);
3156 status = umfpack_dl_symbolic(
width,
width,
AI,
AJ, Ax, &Symbolic,
3161 umfpack_dl_report_status(
Control, status);
3163 " umfpack_dl_symbolic() failed!");
3166 status = umfpack_dl_numeric(
AI,
AJ, Ax, Symbolic, &
Numeric,
3171 umfpack_dl_report_status(
Control, status);
3173 " umfpack_dl_numeric() failed!");
3175 umfpack_dl_free_symbolic(&Symbolic);
3182 mfem_error(
"UMFPackSolver::Mult : matrix is not set!"
3183 " Call SetOperator first!");
3191 umfpack_di_report_info(Control,
Info);
3194 umfpack_di_report_status(Control, status);
3195 mfem_error(
"UMFPackSolver::Mult : umfpack_di_solve() failed!");
3200 SuiteSparse_long status =
3203 umfpack_dl_report_info(Control,
Info);
3206 umfpack_dl_report_status(Control, status);
3207 mfem_error(
"UMFPackSolver::Mult : umfpack_dl_solve() failed!");
3215 mfem_error(
"UMFPackSolver::MultTranspose : matrix is not set!"
3216 " Call SetOperator first!");
3224 umfpack_di_report_info(Control,
Info);
3227 umfpack_di_report_status(Control, status);
3229 " umfpack_di_solve() failed!");
3234 SuiteSparse_long status =
3237 umfpack_dl_report_info(Control,
Info);
3240 umfpack_dl_report_status(Control, status);
3242 " umfpack_dl_solve() failed!");
3255 umfpack_di_free_numeric(&
Numeric);
3259 umfpack_dl_free_numeric(&
Numeric);
3274 "Had Numeric pointer in KLU, but not Symbolic");
3282 MFEM_VERIFY(
mat != NULL,
"not a SparseMatrix");
3290 MFEM_VERIFY(
width ==
height,
"not a square matrix");
3302 MFEM_VERIFY(
mat != NULL,
3303 "KLUSolver::Mult : matrix is not set! Call SetOperator first!");
3316 MFEM_VERIFY(
mat != NULL,
3317 "KLUSolver::Mult : matrix is not set! Call SetOperator first!");
3336 #endif // MFEM_USE_SUITESPARSE
3344 for (
int i = 0; i < block_dof.
NumRows(); ++i)
3349 block_solvers[i].SetOperator(sub_A);
3358 for (
int i = 0; i < block_dof.
NumRows(); ++i)
3363 block_solvers[i].Mult(sub_rhs, sub_sol);
3377 add(-1.0, z, 1.0, x, z);
3394 add(-1.0, z, 1.0, x, z);
3405 bool op_is_symmetric,
3407 :
Solver(op.NumRows()), aux_map_(aux_map, own_aux_map)
3409 aux_system_.
Reset(
RAP(&op, aux_map));
3412 aux_smoother_.
As<
HypreSmoother>()->SetOperatorSymmetry(op_is_symmetric);
3415 void AuxSpaceSmoother::Mult(
const Vector &x,
Vector &y,
bool transpose)
const
3420 Vector aux_sol(aux_rhs.Size());
3427 aux_smoother_->
Mult(aux_rhs, aux_sol);
3431 aux_map_->
Mult(aux_sol, y);
3433 #endif // MFEM_USE_MPI
const Array< int > * ess_dofs_list
Not owned.
void ApplyPlaneRotation(double &dx, double &dy, double &cs, double &sn)
PowerMethod helper class to estimate the largest eigenvalue of an operator using the iterative power ...
int Size() const
Return the logical size of the array.
int RowSize(const int i) const
Returns the number of elements in row i.
Conjugate gradient method.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
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]
void Mult(const Vector &b, Vector &x) const
Solve the system LUx = b, where L and U are the block ILU factors.
Chebyshev accelerated smoothing with given vector, no matrix necessary.
int GetNumIterations() const
MPI_Comm GetComm() const
Return the associated MPI communicator, or MPI_COMM_NULL if no communicator is set.
Array< Vector * > ykArray
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void SetSize(int s)
Resize the vector to size s.
void Monitor(int it, double norm, const Vector &r, const Vector &x, bool final=false) const
virtual void ProcessNewState(const Vector &x) const
This method can be overloaded in derived classes to perform computations that need knowledge of the n...
void MinimumDiscardedFillOrdering(SparseMatrix &C, Array< int > &p)
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
virtual Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate ...
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
bool summary
A summary of the solver process will be reported after the last iteration to mfem::out.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
bool first_and_last
Information about the first and last iteration will be printed to mfem::out.
double Norml2() const
Returns the l2 norm of the vector.
OperatorJacobiSmoother(const double damping=1.0)
Default constructor: the diagonal will be computed by subsequent calls to SetOperator() using the Ope...
int * GetJ()
Return the array J.
class IterativeSolver * iter_solver
The last IterativeSolver to which this monitor was attached.
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
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.
int * GetI()
Return the array I.
bool errors
If a fatal problem has been detected the failure will be reported to mfem::err.
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 ...
void SetOperator(const Operator &op)
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
const Operator * GetC() const
bool GetConverged() const
Reordering
The reordering method used by the BlockILU factorization.
double * GetData() const
Return a pointer to the beginning of the Vector data.
int print_level
(DEPRECATED) Legacy print level definition, which is left for compatibility with custom iterative sol...
void SortColumnIndices()
Sort the column indices corresponding to each row.
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
void GetSubMatrix(const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm) const
void add(const Vector &v1, const Vector &v2, Vector &v)
double Norm(const Vector &x) const
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
double GetFinalNorm() const
double * GetData()
Return the element data, i.e. the array A.
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)
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
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).
General product operator: x -> (A*B)(x) = A(B(x)).
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
double solve(double l, const Vector &xt, Vector &x, int &nclip) const
Solve QP at fixed lambda.
int GetNumConstraints() const
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
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)
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
void Add(const double c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
Jacobi smoothing for a given bilinear form (no matrix necessary).
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
const double * HostReadData() const
void SetMaxIter(int max_it)
const int * HostReadJ() const
OperatorChebyshevSmoother(const Operator &oper_, const Vector &d, const Array< int > &ess_tdof_list, int order, double max_eig_estimate)
void SetLinearConstraint(const Vector &w_, double a_)
DirectSubBlockSolver(const SparseMatrix &A, const SparseMatrix &block_dof)
void SetKDim(int dim)
Set the number of iteration to perform between restarts, default is 50.
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)
PrintLevel & Iterations()
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
void print_iteration(int it, double r, double l) const
void Setup(const Vector &diag)
Parallel smoothers in hypre.
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 int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
void SetSize(int i, int j, int k, MemoryType mt_=MemoryType::PRESERVE)
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
void Solve(int m, int n, double *X) const
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
Array< Vector * > skArray
PrintLevel & FirstAndLast()
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof...
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
BlockILU(int block_size_, Reordering reordering_=Reordering::MINIMUM_DISCARDED_FILL, int k_fill_=0)
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
bool Finalized() const
Returns whether or not CSR format has been finalized.
void AdaptiveLinRtolPostSolve(const Vector &x, const Vector &b, const int it, const double fnorm) const
Method for the adaptive linear solver rtol invoked after the linear solve.
Stationary linear iteration: x <- x + B (b - A x)
void SetInequalityConstraint(const Vector &dl, const Vector &dh)
void SetEqualityConstraint(const Vector &c)
void Swap(Array< T > &, Array< T > &)
void SetAbsTol(double atol)
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
PrintLevel print_options
Output behavior for the iterative solver.
double p(const Vector &x, double t)
void SetRelTol(double rtol)
double Control[UMFPACK_CONTROL]
Abstract base class for iterative solver.
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
int NumRows() const
Get the number of rows (size of output) of the Operator. Synonym with Height().
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
void subtract(const Vector &x, const Vector &y, Vector &z)
MemoryType
Memory types supported by MFEM.
void SetAdaptiveLinRtol(const int type=2, const double rtol0=0.5, const double rtol_max=0.9, const double alpha=0.5 *(1.0+sqrt(5.0)), const double gamma=1.0)
Enable adaptive linear solver relative tolerance algorithm.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int max_iter
Limit for the number of iterations the solver is allowed to do.
int NumCols() const
Get the number of columns (size of input) of the Operator. Synonym with Width().
void Update(Vector &x, int k, DenseMatrix &h, Vector &s, Array< Vector * > &v)
static int GuessLegacyPrintLevel(PrintLevel)
Use some heuristics to guess a legacy print level corresponding to the given PrintLevel.
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.
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
double InnerProduct(HypreParVector *x, HypreParVector *y)
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
void GetDiag(Vector &d) const
Returns the Diagonal of A.
double rel_tol
Relative tolerance.
virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const
This method can be overloaded in derived classes to implement line search algorithms.
OptimizationProblem(int insize, const Operator *C_, const Operator *D_)
In parallel, insize is the number of the local true dofs.
void SetSolutionBounds(const Vector &xl, const Vector &xh)
virtual void AssembleDiagonal(Vector &diag) const
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operat...
AuxSpaceSmoother(const HypreParMatrix &op, HypreParMatrix *aux_map, bool op_is_symmetric=true, bool own_aux_map=false)
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
FDualNumber< tbase > sqrt(const FDualNumber< tbase > &f)
sqrt([dual number])
void SetOperator(const Operator &op)
Recompute the diagonal using the method AssembleDiagonal of the given new Operator, op.
void GeneratePlaneRotation(double &dx, double &dy, double &cs, double &sn)
MemoryType GetMemoryType(MemoryClass mc)
Return a suitable MemoryType for a given MemoryClass.
void AdaptiveLinRtolPreSolve(const Vector &x, const int it, const double fnorm) const
Method for the adaptive linear solver rtol invoked before the linear solve.
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
PrintLevel FromLegacyPrintLevel(int)
Convert a legacy print level integer to a PrintLevel object.
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
virtual void SetOperator(const Operator &op)=0
Set/update the solver for the given operator.
void MakeRef(T *, int)
Make this Array a reference to a pointer.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
void AddMult_a(double a, const Vector &x, Vector &y) const
y += a * A.x
const int * HostReadI() const
void MonitorResidual(int it, double norm, const Vector &r, bool final) override
Monitor the residual vector r.
virtual void MonitorSolution(int it, double norm, const Vector &x, bool final)
Monitor the solution vector x.
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)
bool warnings
If a non-fatal problem has been detected some context-specific information will be reported to mfem::...
bool iterations
Detailed information about each iteration will be reported to mfem::out.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
const OptimizationProblem * problem
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
double EstimateLargestEigenvalue(Operator &opr, Vector &v0, int numSteps=10, double tolerance=1e-8, int seed=12345)
Returns an estimate of the largest eigenvalue of the operator opr using the iterative power method...
Wrapper for hypre's ParCSR matrix class.
void SetBounds(const Vector &lo_, const Vector &hi_)
void AddMult_a(double alpha, const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A += alpha * B * C.
IterativeSolverMonitor * monitor
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
virtual void MonitorResidual(int it, double norm, const Vector &r, bool final)
Monitor the residual vector r.
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
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.
Settings for the output behavior of the IterativeSolver.
const Operator * C
Not owned, some can remain unused (NULL).
double abs_tol
Absolute tolerance.
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
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.