14 std::vector<TDATA> Vs;
19 template <
class Finit>
22 Finit &&finit = [](TDATA &) {}) : kSubspace(NkSubspace)
24 Vs.resize(kSubspace + 1);
44 template <
class TFA,
class TFML,
class TFDot,
class TFstop>
45 bool solve(TFA &&FA, TFML &&FML, TFDot &&fDot, TDATA &
b, TDATA &
x, uint32_t nRestart, TFstop &&FStop)
48 real scale_MLb = std::sqrt(fDot(MLb, MLb));
50 h.setZero(kSubspace + 1, kSubspace);
52 for (iRestart = 0; iRestart <= nRestart; iRestart++)
56 Vs[0].addTo(MLb, -1.0);
57 real beta = std::sqrt(fDot(Vs[0], Vs[0]));
58 if (FStop(iRestart, beta, scale_MLb))
61 (iRestart == nRestart))
66 for (; j < kSubspace; j++)
69 FML(V_temp, Vs[j + 1]);
70 for (uint32_t i = 0; i <= j; i++)
73 h(i, j) = fDot(Vs[j + 1], (Vs[i]));
75 for (uint32_t i = 0; i <= j; i++)
78 Vs[j + 1].addTo(Vs[i], -h(i, j));
80 h(j + 1, j) = std::sqrt(fDot(Vs[j + 1], Vs[j + 1]));
81 if (h(j + 1, j) < scale_MLb * 1e-32)
83 std::cout <<
"early stop" << std::endl;
86 Vs[j + 1] *= 1.0 / h(j + 1, j);
108 MatrixXR hPart = h(Eigen::seq(0, j + 1 - 1), Eigen::seq(0, j - 1));
109 DNDS_assert_info(hPart.allFinite(),
"GMRES_LeftPreconditioned acquired inf or nan coefficient");
117 for (uint32_t jj = 0; jj < j; jj++)
119 x.addTo(Vs[jj], y(jj));
128 return iRestart < nRestart;
132 template <
class TDATA,
class TScalar>
135 bool initialized =
false;
141 index pHistorySize{0};
144 template <
class Finit>
146 Finit &&finit = [](TDATA &) {})
155 void reset() { initialized =
false, pHistorySize = 0; }
159 template <
class TFA,
class TFM,
class TFResPrec,
class TFDot,
class TFstop>
160 bool solve(TFA &&FA, TFM &&FM, TFResPrec &&FResPrec, TFDot &&fDot, TDATA &
x, uint32_t niter, TFstop &&FStop)
170 for (
int iter = 1; iter <= niter; iter++)
177 TScalar zrDotNew = fDot(
z,
r);
178 TScalar zDotz = fDot(
z,
z);
179 if (FStop(iter, zrDotNew, zDotz))
Core type aliases, constants, and metaprogramming utilities for the DNDS framework.
#define DNDS_assert_info(expr, info)
Debug-only assertion with an extra std::string info message.
Robust linear-algebra primitives for small / moderately-sized Eigen matrices that are more numericall...
GMRES_LeftPreconditioned(uint32_t NkSubspace, Finit &&finit=[](TDATA &) {})
bool solve(TFA &&FA, TFML &&FML, TFDot &&fDot, TDATA &b, TDATA &x, uint32_t nRestart, TFstop &&FStop)
PCG_PreconditionedRes(Finit &&finit=[](TDATA &) {})
bool solve(TFA &&FA, TFM &&FM, TFResPrec &&FResPrec, TFDot &&fDot, TDATA &x, uint32_t niter, TFstop &&FStop)
index getPHistorySize() const
Eigen::Index EigenLeastSquareSolve(const Eigen::MatrixXd &A, const Eigen::MatrixXd &B, Eigen::MatrixXd &AIB)
Least-squares solve A * AIB ~= B via a rank-revealing QR-style decomposition; returns the computed ra...
DNDS_CONSTANT const real verySmallReal
Catch-all lower bound ("effectively zero").
Eigen::Matrix< real, Eigen::Dynamic, Eigen::Dynamic > MatrixXR
Column-major dynamic Eigen matrix of reals (default layout).
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Eigen::Vector< real, Eigen::Dynamic > VectorXR
Dynamic Eigen vector of reals.