DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
Linear.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "DNDS/HardEigen.hpp"
4#include "DNDS/Defines.hpp"
5
6namespace DNDS::Linear
7{
8
9 template <class TDATA>
11 {
12 uint32_t kSubspace;
13
14 std::vector<TDATA> Vs;
15 TDATA V_temp;
16 TDATA MLb;
17
18 public:
19 template <class Finit>
21 uint32_t NkSubspace,
22 Finit &&finit = [](TDATA &) {}) : kSubspace(NkSubspace)
23 {
24 Vs.resize(kSubspace + 1);
25 for (auto &i : Vs)
26 finit(i);
27 finit(V_temp);
28 finit(MLb);
29 }
30
31 /**
32 * @brief
33 *
34 * @tparam TFA
35 * @tparam TFML
36 * @tparam TFstop
37 * @param FA void FA(TDATA &x, TDATA &Ax)
38 * @param FML void FML(TDATA &x, TDATA &MLx)
39 * @param b rhs
40 * @param x input and output
41 * @param nRestart
42 * @param FStop bool FStop(iRestart, res, resB)
43 */
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)
46 {
47 FML(b, MLb); // MLb = ML * b
48 real scale_MLb = std::sqrt(fDot(MLb, MLb));
49 MatrixXR h;
50 h.setZero(kSubspace + 1, kSubspace);
51 uint32_t iRestart;
52 for (iRestart = 0; iRestart <= nRestart; iRestart++)
53 {
54 FA(x, V_temp); // V_temp = A * x
55 FML(V_temp, Vs[0]); // Vs[0] = ML * A * x
56 Vs[0].addTo(MLb, -1.0); // Vs[0] = ML * A * x - ML * b = -r
57 real beta = std::sqrt(fDot(Vs[0], Vs[0])); // beta = norm2(r)
58 if (FStop(iRestart, beta, scale_MLb)) // see if converge
59 break;
60 if (std::abs(beta) < verySmallReal || // beta is floating-point negligible
61 (iRestart == nRestart))
62 break;
63 Vs[0] *= -1.0 / beta; // Vs[0] = r/norm2(r)
64
65 uint32_t j = 0;
66 for (; j < kSubspace; j++) // Arnoldi
67 {
68 FA(Vs[j], V_temp); // V_temp = A * Vs[j]
69 FML(V_temp, Vs[j + 1]); // Vs[j + 1] = ML * A * Vs[j]
70 for (uint32_t i = 0; i <= j; i++)
71 {
72 // Gram-Schmidt, calculate projected lengths
73 h(i, j) = fDot(Vs[j + 1], (Vs[i])); // h(i, j) = dot(ML * A * Vs[j],Vs[i])
74 }
75 for (uint32_t i = 0; i <= j; i++)
76 {
77 // Gram-Schmidt, subtract projections
78 Vs[j + 1].addTo(Vs[i], -h(i, j)); // Vs[j + 1] = ML * A * Vs[j] - sum_{i=0,1,...j}(h(i,j) * Vs[i])
79 }
80 h(j + 1, j) = std::sqrt(fDot(Vs[j + 1], Vs[j + 1])); //
81 if (h(j + 1, j) < scale_MLb * 1e-32)
82 {
83 std::cout << "early stop" << std::endl;
84 break;
85 }
86 Vs[j + 1] *= 1.0 / h(j + 1, j); // normalize
87 }
88 // std::cout << beta << std::endl;
89
90 VectorXR eBeta;
91 eBeta.resize(j + 1);
92 eBeta.setZero();
93 eBeta(0) = beta; // eBeta = e_1 * beta
94 if (j < 1)
95 break;
96
97 { // the QR method
98 // auto QR = h.colPivHouseholderQr();
99 // QR.setThreshold(std::sqrt(scale_MLb) * 1e-32);
100 // if (QR.rank() != h.rows())
101 // DNDS_assert_info(false, "GMRES not good");
102 // Eigen::VectorXd y = QR.solve(eBeta);
103 }
104
105 // auto sol = h(Eigen::seq(0, j + 1 - 1), Eigen::seq(0, j - 1)).bdcSvd(Eigen::ComputeThinU | Eigen::ComputeThinV);
106 // Eigen::VectorXd y = sol.solve(eBeta); //! warning: this gmres is dumb and do not lucky stop
107 MatrixXR y;
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");
110 auto rank = HardEigen::EigenLeastSquareSolve(hPart, eBeta, y);
111
112 // int rank;
113 // MPI_Comm_rank(MPI_COMM_WORLD, &rank);
114 // if (rank == 0)
115 // std::cout << h << std::endl;
116
117 for (uint32_t jj = 0; jj < j; jj++) // x = V(:, 0,1,2,...kSubspace-1) * y
118 {
119 x.addTo(Vs[jj], y(jj));
120 // std::cout << iRestart << "::" << Vs[j].transpose() << "::" << y(j) << std::endl;
121 }
122 if (rank < h.cols())
123 {
124
125 break; // do not restart
126 }
127 }
128 return iRestart < nRestart;
129 }
130 };
131
132 template <class TDATA, class TScalar>
134 {
135 bool initialized = false;
136 TDATA r, z, p, Ap;
137 TDATA V_temp;
138
139 TScalar zrDot;
140
141 index pHistorySize{0};
142
143 public:
144 template <class Finit>
146 Finit &&finit = [](TDATA &) {})
147 {
148 finit(r);
149 finit(z);
150 finit(p);
151 finit(Ap);
152 finit(V_temp);
153 }
154
155 void reset() { initialized = false, pHistorySize = 0; }
156
157 [[nodiscard]] index getPHistorySize() const { return pHistorySize; }
158
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)
161 {
162 if (!initialized)
163 {
164 FResPrec(x, z);
165 FM(z, r);
166 p = z;
167 pHistorySize++;
168 FA(p, Ap);
169 }
170 for (int iter = 1; iter <= niter; iter++)
171 {
172 if (initialized)
173 {
174 FResPrec(x, z);
175 FM(z, r);
176 }
177 TScalar zrDotNew = fDot(z, r);
178 TScalar zDotz = fDot(z, z);
179 if (FStop(iter, zrDotNew, zDotz))
180 return true;
181
182 if (initialized)
183 {
184 TScalar beta = zrDotNew / (zrDot + verySmallReal);
185 p *= beta;
186 p += z;
187 pHistorySize++;
188 FA(p, Ap);
189 }
190
191 TScalar alpha = zrDotNew / (fDot(p, Ap) + verySmallReal);
192 x.addTo(p, alpha);
193
194 zrDot = zrDotNew;
195 initialized = true;
196 }
197 return false;
198 }
199 };
200}
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.
Definition Errors.hpp:113
Robust linear-algebra primitives for small / moderately-sized Eigen matrices that are more numericall...
GMRES_LeftPreconditioned(uint32_t NkSubspace, Finit &&finit=[](TDATA &) {})
Definition Linear.hpp:20
bool solve(TFA &&FA, TFML &&FML, TFDot &&fDot, TDATA &b, TDATA &x, uint32_t nRestart, TFstop &&FStop)
Definition Linear.hpp:45
PCG_PreconditionedRes(Finit &&finit=[](TDATA &) {})
Definition Linear.hpp:145
bool solve(TFA &&FA, TFM &&FM, TFResPrec &&FResPrec, TFDot &&fDot, TDATA &x, uint32_t niter, TFstop &&FStop)
Definition Linear.hpp:160
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").
Definition Defines.hpp:194
Eigen::Matrix< real, Eigen::Dynamic, Eigen::Dynamic > MatrixXR
Column-major dynamic Eigen matrix of reals (default layout).
Definition Defines.hpp:205
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
Definition Defines.hpp:107
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Definition Defines.hpp:105
Eigen::Vector< real, Eigen::Dynamic > VectorXR
Dynamic Eigen vector of reals.
Definition Defines.hpp:207
tVec r(NCells)
tVec z(NCells)
tVec x(NCells)
tVec b(NCells)
real alpha