78struct BiCGStabLSolverParams
80 using params = BiCGStabLSolverParams;
93 ePreconditionerSideType pside = ePreconditionerSideType::right;
102 double abstol = std::numeric_limits<double>::min();
114 BiCGStabLSolverParams() =
default;
117 : ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, L)
118 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, delta)
119 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, convex)
120 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, pside)
121 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, maxiter)
122 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, tol)
123 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, abstol)
124 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p,
ns_search)
125 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p,
verbose)
127 p.check_params({
"L",
"delta",
"convex",
"pside",
"maxiter",
"tol",
"abstol",
"ns_search",
"verbose" });
130 void get(
PropertyTree& p,
const std::string& path)
const
132 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, L);
133 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, delta);
134 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, convex);
135 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, pside);
136 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, maxiter);
137 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, tol);
138 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, abstol);
139 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path,
ns_search);
140 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path,
verbose);
155 using backend_type = Backend;
156 using BackendType = Backend;
158 typedef typename Backend::vector vector;
159 typedef typename Backend::value_type value_type;
160 typedef typename Backend::params backend_params;
162 typedef typename math::scalar_of<value_type>::type scalar_type;
165 typename math::rhs_of<value_type>::type>::return_type coef_type;
171 const params& prm = params(),
172 const backend_params& backend_prm = backend_params(),
173 const InnerProduct& inner_product = InnerProduct())
176 , Rt(Backend::create_vector(n, backend_prm))
177 , X(Backend::create_vector(n, backend_prm))
178 , B(Backend::create_vector(n, backend_prm))
179 , T(Backend::create_vector(n, backend_prm))
182 , MZa(prm.L + 1, prm.L + 1)
183 , MZb(prm.L + 1, prm.L + 1)
186 , inner_product(inner_product)
188 precondition(prm.L > 0,
"L in BiCGStab(L) should be >=1");
190 for (
int i = 0; i <= prm.L; ++i) {
191 R[i] = Backend::create_vector(n, backend_prm);
192 U[i] = Backend::create_vector(n, backend_prm);
211 template <
class Matrix,
class Precond,
class Vec1,
class Vec2>
212 SolverResult operator()(
const Matrix& A,
const Precond& P,
const Vec1& rhs, Vec2&& x)
const
214 static const coef_type one = math::identity<coef_type>();
215 static const coef_type zero = math::zero<coef_type>();
221 scalar_type norm_rhs = norm(rhs);
224 if (norm_rhs < Alina::detail::eps<scalar_type>(1)) {
226 norm_rhs = math::identity<scalar_type>();
230 return SolverResult(0, norm_rhs);
234 if (prm.pside == ePreconditionerSideType::left) {
235 backend::residual(rhs, A, x, *T);
239 backend::residual(rhs, A, x, *B);
242 scalar_type zeta0 = norm(*B);
243 scalar_type eps = std::max(prm.tol * norm_rhs, prm.abstol);
245 coef_type alpha = zero;
246 coef_type rho0 = one;
247 coef_type omega = one;
250 backend::copy(*B, *R[0]);
251 backend::copy(*B, *Rt);
253 backend::clear(*U[0]);
255 scalar_type zeta = zeta0;
256 scalar_type rnmax_computed = zeta0;
257 scalar_type rnmax_true = zeta0;
260 for (; iter < prm.maxiter && zeta >= eps; iter += L) {
262 rho0 = -omega * rho0;
264 for (
int j = 0; j < L; ++j) {
265 coef_type rho1 = inner_product(*R[j], *Rt);
266 precondition(!math::is_zero(rho1),
"BiCGStab(L) breakdown: diverged (zero rho)");
268 coef_type beta = alpha * (rho1 / rho0);
271 for (
int i = 0; i <= j; ++i)
272 backend::axpby(one, *R[i], -beta, *U[i]);
274 preconditioner_spmv(prm.pside, P, A, *U[j], *U[j + 1], *T);
276 coef_type sigma = inner_product(*U[j + 1], *Rt);
277 precondition(!math::is_zero(sigma),
"BiCGStab(L) breakdown: diverged (zero sigma)");
278 alpha = rho1 / sigma;
280 backend::axpby(alpha, *U[0], one, *X);
282 for (
int i = 0; i <= j; ++i)
283 backend::axpby(-alpha, *U[i + 1], one, *R[i]);
285 preconditioner_spmv(prm.pside, P, A, *R[j], *R[j + 1], *T);
289 rnmax_computed = std::max(zeta, rnmax_computed);
290 rnmax_true = std::max(zeta, rnmax_true);
300 for (
int i = 0; i <= L; ++i) {
301 for (
int j = 0; j <= i; ++j) {
302 MZa(i, j) = inner_product(*R[i], *R[j]);
307 for (
int i = 0; i <= L; ++i) {
308 for (
int j = i + 1; j <= L; ++j) {
309 MZa(i, j) = MZa(j, i) = math::adjoint(MZa(j, i));
313 std::copy(MZa.data(), MZa.data() + MZa.size(), MZb.data());
315 if (prm.convex || L == 1) {
318 qr.solve(L, L, MZa.stride(0), MZa.stride(1),
319 &MZa(1, 1), &MZb(0, 1), &Y0[1]);
324 qr.solve(L - 1, L - 1, MZa.stride(0), MZa.stride(1),
325 &MZa(1, 1), &MZb(0, 1), &Y0[1]);
329 qr.solve(L - 1, L - 1, MZa.stride(0), MZa.stride(1),
330 &MZa(1, 1), &MZb(L, 1), &YL[1],
true);
332 coef_type dot0 = zero;
333 coef_type dot1 = zero;
334 coef_type dotA = zero;
335 for (
int i = 0; i <= L; ++i) {
339 for (
int j = 0; j <= L; ++j) {
340 coef_type M = MZb(i, j);
350 scalar_type kappa0 =
sqrt(std::abs(std::real(dot0)));
351 scalar_type kappa1 =
sqrt(std::abs(std::real(dot1)));
352 scalar_type kappaA = std::real(dotA);
354 if (!math::is_zero(kappa0) && !math::is_zero(kappa1)) {
356 if (kappaA < 0.7 * kappa0 * kappa1) {
357 ghat = (kappaA < 0) ? -0.7 * kappa0 / kappa1 : 0.7 * kappa0 / kappa1;
360 ghat = kappaA / (kappa1 * kappa1);
363 for (
int i = 0; i <= L; ++i)
364 Y0[i] -= ghat * YL[i];
369 for (
int h = L; h > 0 && math::is_zero(omega); --h)
371 precondition(!math::is_zero(omega),
"BiCGStab(L) breakdown: diverged (zero omega)");
373 backend::lin_comb(L, &Y0[1], &R[0], one, *X);
375 for (
int i = 1; i <= L; ++i)
376 Y0[i] = -one * Y0[i];
378 backend::lin_comb(L, &Y0[1], &U[1], one, *U[0]);
379 backend::lin_comb(L, &Y0[1], &R[1], one, *R[0]);
381 for (
int i = 1; i <= L; ++i)
382 Y0[i] = -one * Y0[i];
388 rnmax_computed = std::max(zeta, rnmax_computed);
389 rnmax_true = std::max(zeta, rnmax_true);
391 bool update_x = zeta < prm.delta * zeta0 && zeta0 <= rnmax_computed;
393 if ((zeta < prm.delta * rnmax_true && zeta <= rnmax_true) || update_x) {
394 preconditioner_spmv(prm.pside, P, A, *X, *R[0], *T);
395 backend::axpby(one, *B, -one, *R[0]);
399 if (prm.pside == ePreconditionerSideType::left) {
400 backend::axpby(one, *X, one, x);
403 backend::axpby(one, *T, one, x);
406 backend::copy(*R[0], *B);
408 rnmax_computed = zeta;
412 if (prm.verbose && iter % 5 == 0)
413 std::cout << iter <<
"\t" << std::scientific << zeta / norm_rhs << std::endl;
417 if (prm.pside == ePreconditionerSideType::left) {
418 backend::axpby(one, *X, one, x);
422 backend::axpby(one, *T, one, x);
425 return SolverResult(iter, zeta / norm_rhs);
438 template <
class Precond,
class Vec1,
class Vec2>
441 return (*
this)(P.system_matrix(), P, rhs, x);
448 b += backend::bytes(*Rt);
449 b += backend::bytes(*X);
450 b += backend::bytes(*B);
451 b += backend::bytes(*T);
453 for (
const auto& v : R)
454 b += backend::bytes(*v);
455 for (
const auto& v : U)
456 b += backend::bytes(*v);
458 b += MZa.size() *
sizeof(coef_type);
459 b += MZb.size() *
sizeof(coef_type);
461 b += backend::bytes(Y0);
462 b += backend::bytes(YL);
469 friend std::ostream& operator<<(std::ostream& os,
const BiCGStabLSolver& s)
471 return os <<
"Type: BiCGStab(" << s.prm.L <<
")"
472 <<
"\nUnknowns: " << s.n
473 <<
"\nMemory footprint: " << human_readable_memory(s.
bytes())
485 mutable std::shared_ptr<vector> Rt;
486 mutable std::shared_ptr<vector> X;
487 mutable std::shared_ptr<vector> B;
488 mutable std::shared_ptr<vector> T;
490 mutable std::vector<std::shared_ptr<vector>> R;
491 mutable std::vector<std::shared_ptr<vector>> U;
493 mutable multi_array<coef_type, 2> MZa, MZb;
494 mutable std::vector<coef_type> Y0, YL;
497 InnerProduct inner_product;
500 scalar_type norm(
const Vec& x)
const
502 return sqrt(math::norm(inner_product(x, x)));