15#ifndef dealii_solver_gmres_h
16#define dealii_solver_gmres_h
49 template <
typename,
typename>
65 namespace SolverGMRESImplementation
74 template <
typename VectorType>
122 std::vector<typename VectorMemory<VectorType>::Pointer>
data;
145 const unsigned int max_basis_size,
146 const bool force_reorthogonalization);
170 template <
typename VectorType>
173 const unsigned int n,
175 const unsigned int accumulated_iterations = 0,
176 const boost::signals2::signal<
void(
int)> &reorthogonalize_signal =
177 boost::signals2::signal<
void(
int)>());
277 std::vector<std::pair<double, double>> &rotations,
350template <
typename VectorType = Vector<
double>>
371 const bool right_preconditioning =
false,
372 const bool use_default_residual =
true,
373 const bool force_re_orthogonalization =
false,
374 const bool batched_mode =
false,
376 orthogonalization_strategy =
378 delayed_classical_gram_schmidt);
459 template <
typename MatrixType,
typename PreconditionerType>
463 void solve(const MatrixType &A,
466 const PreconditionerType &preconditioner);
474 boost::signals2::connection
475 connect_condition_number_slot(const
std::function<
void(
double)> &slot,
476 const
bool every_iteration = false);
484 boost::signals2::connection
485 connect_eigenvalues_slot(
486 const
std::function<
void(const
std::vector<
std::complex<
double>> &)> &slot,
487 const
bool every_iteration = false);
496 boost::signals2::connection
497 connect_hessenberg_slot(
499 const
bool every_iteration = true);
507 boost::signals2::connection
508 connect_krylov_space_slot(
510 void(const
internal::SolverGMRESImplementation::TmpVectors<VectorType> &)>
518 boost::signals2::connection
519 connect_re_orthogonalization_slot(const
std::function<
void(
int)> &slot);
524 << "The number of temporary vectors you gave (" << arg1
525 << ") is too small. It should be at least 10 for "
526 << "any results, and much more for reasonable ones.");
538 boost::signals2::signal<
void(
double)> condition_number_signal;
544 boost::signals2::signal<
void(
double)> all_condition_numbers_signal;
550 boost::signals2::signal<
void(const
std::vector<
std::complex<
double>> &)>
557 boost::signals2::signal<
void(const
std::vector<
std::complex<
double>> &)>
558 all_eigenvalues_signal;
571 all_hessenberg_signal;
577 boost::signals2::signal<
void(
578 const
internal::SolverGMRESImplementation::TmpVectors<VectorType> &)>
585 boost::signals2::signal<
void(
int)> re_orthogonalize_signal;
607 compute_eigs_and_cond(
609 const
unsigned int n,
610 const
boost::signals2::signal<
611 void(const
std::vector<
std::complex<
double>> &)> &eigenvalues_signal,
614 const
boost::signals2::signal<
void(
double)> &cond_signal);
620 internal::SolverGMRESImplementation::ArnoldiProcess arnoldi_process;
645template <typename VectorType =
Vector<
double>>
660 orthogonalization_strategy =
662 delayed_classical_gram_schmidt)
663 : max_basis_size(max_basis_size)
664 , orthogonalization_strategy(orthogonalization_strategy)
695 template <
typename MatrixType,
typename PreconditionerType>
699 void solve(const MatrixType &A,
702 const PreconditionerType &preconditioner);
714 internal::SolverGMRESImplementation::ArnoldiProcess arnoldi_process;
723template <
typename VectorType>
726 const unsigned int max_basis_size,
727 const bool right_preconditioning,
728 const bool use_default_residual,
729 const bool force_re_orthogonalization,
730 const bool batched_mode,
732 : max_n_tmp_vectors(0)
733 , max_basis_size(max_basis_size)
734 , right_preconditioning(right_preconditioning)
735 , use_default_residual(use_default_residual)
736 , force_re_orthogonalization(force_re_orthogonalization)
737 , batched_mode(batched_mode)
738 , orthogonalization_strategy(orthogonalization_strategy)
740 Assert(max_basis_size >= 1,
741 ExcMessage(
"SolverGMRES needs at least one vector in the "
747template <
typename VectorType>
751 const AdditionalData &data)
753 , additional_data(data)
759template <
typename VectorType>
762 const AdditionalData &data)
764 , additional_data(data)
772 namespace SolverGMRESImplementation
774 template <
typename VectorType>
783 template <
typename VectorType>
785 TmpVectors<VectorType>::operator[](
const unsigned int i)
const
795 template <
typename VectorType>
797 TmpVectors<VectorType>::operator()(
const unsigned int i,
798 const VectorType &temp)
801 if (data[i] ==
nullptr)
804 data[i]->reinit(temp,
true);
811 template <
typename VectorType>
813 TmpVectors<VectorType>::size()
const
815 return (data.size() > 0 ? data.size() - 1 : 0);
820 template <
typename VectorType,
typename Enable =
void>
821 struct is_dealii_compatible_distributed_vector;
823 template <
typename VectorType>
824 struct is_dealii_compatible_distributed_vector<
826 std::enable_if_t<!internal::is_block_vector<VectorType>>>
828 static constexpr bool value = std::is_same_v<
836 template <
typename VectorType>
837 struct is_dealii_compatible_distributed_vector<
839 std::enable_if_t<internal::is_block_vector<VectorType>>>
841 static constexpr bool value = std::is_same_v<
842 typename VectorType::BlockType,
849 template <
typename VectorType,
850 std::enable_if_t<!IsBlockVector<VectorType>::value, VectorType>
860 template <
typename VectorType,
861 std::enable_if_t<IsBlockVector<VectorType>::value, VectorType> * =
866 return vector.n_blocks();
871 template <
typename VectorType,
872 std::enable_if_t<!IsBlockVector<VectorType>::value, VectorType>
875 block(VectorType &vector,
const unsigned int b)
884 template <
typename VectorType,
885 std::enable_if_t<!IsBlockVector<VectorType>::value, VectorType>
888 block(
const VectorType &vector,
const unsigned int b)
897 template <
typename VectorType,
898 std::enable_if_t<IsBlockVector<VectorType>::value, VectorType> * =
900 typename VectorType::BlockType &
901 block(VectorType &vector,
const unsigned int b)
903 return vector.block(b);
908 template <
typename VectorType,
909 std::enable_if_t<IsBlockVector<VectorType>::value, VectorType> * =
911 const typename VectorType::BlockType &
912 block(
const VectorType &vector,
const unsigned int b)
914 return vector.block(b);
919 template <
bool delayed_reorthogonalization,
922 !is_dealii_compatible_distributed_vector<VectorType>::value,
923 VectorType> * =
nullptr>
925 Tvmult_add(
const unsigned int n,
926 const VectorType &vv,
927 const TmpVectors<VectorType> &orthogonal_vectors,
930 for (
unsigned int i = 0; i < n; ++i)
932 h(i) += vv * orthogonal_vectors[i];
933 if (delayed_reorthogonalization)
934 h(n + i) += orthogonal_vectors[i] * orthogonal_vectors[n - 1];
936 if (delayed_reorthogonalization)
942 template <
bool delayed_reorthogonalization,
945 is_dealii_compatible_distributed_vector<VectorType>::value,
946 VectorType> * =
nullptr>
948 Tvmult_add(
const unsigned int n,
949 const VectorType &vv,
950 const TmpVectors<VectorType> &orthogonal_vectors,
953 for (
unsigned int b = 0;
b <
n_blocks(vv); ++
b)
960 static constexpr unsigned int n_lanes =
964 for (
unsigned int i = 0; i < n; ++i)
967 correct[delayed_reorthogonalization ? 129 : 1];
968 if (delayed_reorthogonalization)
969 for (
unsigned int i = 0; i < n + 1; ++i)
974 constexpr unsigned int inner_batch_size =
975 delayed_reorthogonalization ? 6 : 12;
977 for (; c < block(vv, b).locally_owned_size() / n_lanes /
979 ++c, j += n_lanes * inner_batch_size)
982 for (
unsigned int k = 0; k < inner_batch_size; ++k)
983 vvec[k].load(block(vv, b).begin() + j + k * n_lanes);
985 for (
unsigned int k = 0; k < inner_batch_size; ++k)
987 block(orthogonal_vectors[n - 1], b).begin() + j +
992 last_vector[0] * vvec[0];
994 last_vector[0] * last_vector[0];
996 for (
unsigned int k = 1; k < inner_batch_size; ++k)
998 local_sum_0 += last_vector[k] * vvec[k];
999 if (delayed_reorthogonalization)
1001 local_sum_1 += last_vector[k] * last_vector[k];
1002 local_sum_2 += vvec[k] * vvec[k];
1005 hs[n - 1] += local_sum_0;
1006 if (delayed_reorthogonalization)
1008 correct[n - 1] += local_sum_1;
1009 correct[n] += local_sum_2;
1013 for (
unsigned int i = 0; i < n - 1; ++i)
1019 temp.
load(block(orthogonal_vectors[i], b).begin() + j);
1022 delayed_reorthogonalization ? temp * last_vector[0] :
1024 for (
unsigned int k = 1; k < inner_batch_size; ++k)
1026 temp.
load(block(orthogonal_vectors[i], b).begin() +
1028 local_sum_0 += temp * vvec[k];
1029 if (delayed_reorthogonalization)
1030 local_sum_1 += temp * last_vector[k];
1032 hs[i] += local_sum_0;
1033 if (delayed_reorthogonalization)
1034 correct[i] += local_sum_1;
1038 c *= inner_batch_size;
1039 for (; c < block(vv, b).locally_owned_size() / n_lanes;
1043 vvec.
load(block(vv, b).begin() + j);
1044 last_vector.
load(block(orthogonal_vectors[n - 1], b).begin() +
1046 hs[n - 1] += last_vector * vvec;
1047 if (delayed_reorthogonalization)
1049 correct[n - 1] += last_vector * last_vector;
1050 correct[n] += vvec * vvec;
1053 for (
unsigned int i = 0; i < n - 1; ++i)
1056 temp.
load(block(orthogonal_vectors[i], b).begin() + j);
1057 hs[i] += temp * vvec;
1058 if (delayed_reorthogonalization)
1059 correct[i] += temp * last_vector;
1063 for (
unsigned int i = 0; i < n; ++i)
1065 h(i) += hs[i].
sum();
1066 if (delayed_reorthogonalization)
1067 h(i + n) += correct[i].
sum();
1069 if (delayed_reorthogonalization)
1070 h(n + n) += correct[n].
sum();
1075 for (; j < block(vv, b).locally_owned_size(); ++j)
1077 const double vvec = block(vv, b).local_element(j);
1078 const double last_vector =
1079 block(orthogonal_vectors[n - 1], b).local_element(j);
1080 h(n - 1) += last_vector * vvec;
1081 if (delayed_reorthogonalization)
1083 h(n + n - 1) += last_vector * last_vector;
1084 h(n + n) += vvec * vvec;
1086 for (
unsigned int i = 0; i < n - 1; ++i)
1089 block(orthogonal_vectors[i], b).local_element(j);
1090 h(i) += temp * vvec;
1091 if (delayed_reorthogonalization)
1092 h(n + i) += temp * last_vector;
1102 template <
bool delayed_reorthogonalization,
1103 typename VectorType,
1105 !is_dealii_compatible_distributed_vector<VectorType>::value,
1106 VectorType> * =
nullptr>
1108 subtract_and_norm(
const unsigned int n,
1109 const TmpVectors<VectorType> &orthogonal_vectors,
1115 VectorType &last_vector =
1116 const_cast<VectorType &
>(orthogonal_vectors[n - 1]);
1117 for (
unsigned int i = 0; i < n - 1; ++i)
1119 if (delayed_reorthogonalization && i + 2 < n)
1120 last_vector.add(-h(n + i), orthogonal_vectors[i]);
1121 vv.add(-h(i), orthogonal_vectors[i]);
1124 if (delayed_reorthogonalization)
1127 last_vector.sadd(1. / h(n + n - 1),
1128 -h(n + n - 2) / h(n + n - 1),
1129 orthogonal_vectors[n - 2]);
1132 const double scaling_factor_vv = h(n + n) > 0.0 ?
1133 1. / (h(n + n - 1) * h(n + n)) :
1134 1. / (h(n + n - 1) * h(n + n - 1));
1135 vv.sadd(scaling_factor_vv,
1136 -h(n - 1) * scaling_factor_vv,
1141 return std::numeric_limits<double>::signaling_NaN();
1145 vv.add_and_dot(-h(n - 1), orthogonal_vectors[n - 1], vv));
1150 template <
bool delayed_reorthogonalization,
1151 typename VectorType,
1153 is_dealii_compatible_distributed_vector<VectorType>::value,
1154 VectorType> * =
nullptr>
1156 subtract_and_norm(
const unsigned int n,
1157 const TmpVectors<VectorType> &orthogonal_vectors,
1163 double norm_vv_temp = 0.0;
1164 VectorType &last_vector =
1165 const_cast<VectorType &
>(orthogonal_vectors[n - 1]);
1166 const double inverse_norm_previous =
1167 delayed_reorthogonalization ? 1. / h(n + n - 1) : 0.;
1168 const double scaling_factor_vv =
1169 delayed_reorthogonalization ?
1170 (h(n + n) > 0.0 ? inverse_norm_previous / h(n + n) :
1171 inverse_norm_previous / h(n + n - 1)) :
1174 for (
unsigned int b = 0;
b <
n_blocks(vv); ++
b)
1178 constexpr unsigned int inner_batch_size =
1179 delayed_reorthogonalization ? 6 : 12;
1184 block(vv, b).locally_owned_size() / n_lanes / inner_batch_size;
1185 ++c, j += n_lanes * inner_batch_size)
1190 const double last_factor = h(n - 1);
1191 for (
unsigned int k = 0; k < inner_batch_size; ++k)
1193 temp[k].
load(block(vv, b).
begin() + j + k * n_lanes);
1194 last_vec[k].
load(block(last_vector, b).
begin() + j +
1196 if (!delayed_reorthogonalization)
1197 temp[k] -= last_factor * last_vec[k];
1200 for (
unsigned int i = 0; i < n - 1; ++i)
1202 const double factor = h(i);
1203 const double correction_factor =
1204 (delayed_reorthogonalization ? h(n + i) : 0.0);
1205 for (
unsigned int k = 0; k < inner_batch_size; ++k)
1208 vec.
load(block(orthogonal_vectors[i], b).
begin() + j +
1210 temp[k] -= factor * vec;
1211 if (delayed_reorthogonalization)
1212 last_vec[k] -= correction_factor * vec;
1216 if (delayed_reorthogonalization)
1217 for (
unsigned int k = 0; k < inner_batch_size; ++k)
1219 last_vec[k] = last_vec[k] * inverse_norm_previous;
1220 last_vec[k].
store(block(last_vector, b).
begin() + j +
1222 temp[k] -= last_factor * last_vec[k];
1223 temp[k] = temp[k] * scaling_factor_vv;
1224 temp[k].
store(block(vv, b).
begin() + j + k * n_lanes);
1227 for (
unsigned int k = 0; k < inner_batch_size; ++k)
1229 temp[k].
store(block(vv, b).
begin() + j + k * n_lanes);
1230 norm_vv_temp_vectorized += temp[k] * temp[k];
1234 c *= inner_batch_size;
1235 for (; c < block(vv, b).locally_owned_size() / n_lanes;
1239 temp.
load(block(vv, b).begin() + j);
1240 last_vec.
load(block(last_vector, b).begin() + j);
1241 if (!delayed_reorthogonalization)
1242 temp -= h(n - 1) * last_vec;
1244 for (
unsigned int i = 0; i < n - 1; ++i)
1247 vec.
load(block(orthogonal_vectors[i], b).
begin() + j);
1249 if (delayed_reorthogonalization)
1250 last_vec -= h(n + i) * vec;
1253 if (delayed_reorthogonalization)
1255 last_vec = last_vec * inverse_norm_previous;
1256 last_vec.
store(block(last_vector, b).
begin() + j);
1257 temp -= h(n - 1) * last_vec;
1258 temp = temp * scaling_factor_vv;
1264 norm_vv_temp_vectorized += temp * temp;
1268 if (!delayed_reorthogonalization)
1269 norm_vv_temp += norm_vv_temp_vectorized.
sum();
1271 for (; j < block(vv, b).locally_owned_size(); ++j)
1273 double temp = block(vv, b).local_element(j);
1274 double last_vec = block(last_vector, b).local_element(j);
1275 if (delayed_reorthogonalization)
1277 for (
unsigned int i = 0; i < n - 1; ++i)
1280 block(orthogonal_vectors[i], b).local_element(j);
1282 last_vec -= h(n + i) * vec;
1284 last_vec *= inverse_norm_previous;
1285 block(last_vector, b).local_element(j) = last_vec;
1286 temp -= h(n - 1) * last_vec;
1287 temp *= scaling_factor_vv;
1291 temp -= h(n - 1) * last_vec;
1292 for (
unsigned int i = 0; i < n - 1; ++i)
1294 h(i) * block(orthogonal_vectors[i], b).local_element(j);
1295 norm_vv_temp += temp * temp;
1297 block(vv, b).local_element(j) = temp;
1307 template <
typename VectorType,
1309 !is_dealii_compatible_distributed_vector<VectorType>::value,
1310 VectorType> * =
nullptr>
1313 const unsigned int n,
1315 const TmpVectors<VectorType> &tmp_vectors,
1316 const bool zero_out)
1319 p.equ(h(0), tmp_vectors[0]);
1321 p.add(h(0), tmp_vectors[0]);
1323 for (
unsigned int i = 1; i < n; ++i)
1324 p.add(h(i), tmp_vectors[i]);
1329 template <
typename VectorType,
1331 is_dealii_compatible_distributed_vector<VectorType>::value,
1332 VectorType> * =
nullptr>
1335 const unsigned int n,
1337 const TmpVectors<VectorType> &tmp_vectors,
1338 const bool zero_out)
1340 for (
unsigned int b = 0;
b <
n_blocks(p); ++
b)
1341 for (
unsigned int j = 0; j < block(p, b).locally_owned_size(); ++j)
1343 double temp = zero_out ? 0 : block(p, b).local_element(j);
1344 for (
unsigned int i = 0; i < n; ++i)
1345 temp += block(tmp_vectors[i], b).local_element(j) * h(i);
1346 block(p, b).local_element(j) = temp;
1353 ArnoldiProcess::initialize(
1355 const unsigned int basis_size,
1356 const bool force_reorthogonalization)
1358 this->orthogonalization_strategy = orthogonalization_strategy;
1359 this->do_reorthogonalization = force_reorthogonalization;
1361 hessenberg_matrix.reinit(basis_size + 1, basis_size);
1362 triangular_matrix.reinit(basis_size + 1, basis_size,
true);
1365 projected_rhs.reinit(basis_size + 1,
true);
1366 givens_rotations.reserve(basis_size);
1368 if (orthogonalization_strategy ==
1371 h.
reinit(2 * basis_size + 3);
1373 h.
reinit(basis_size + 1);
1378 template <
typename VectorType>
1380 ArnoldiProcess::orthonormalize_nth_vector(
1381 const unsigned int n,
1382 TmpVectors<VectorType> &orthogonal_vectors,
1383 const unsigned int accumulated_iterations,
1384 const boost::signals2::signal<
void(
int)> &reorthogonalize_signal)
1389 VectorType &vv = orthogonal_vectors[n];
1391 double residual_estimate = std::numeric_limits<double>::signaling_NaN();
1394 givens_rotations.clear();
1395 residual_estimate = vv.l2_norm();
1396 if (residual_estimate != 0.)
1397 vv /= residual_estimate;
1398 projected_rhs(0) = residual_estimate;
1400 else if (orthogonalization_strategy ==
1410 const double previous_scaling = n > 0 ? h(n + n - 2) : 1.;
1413 h.reinit(n + n + 1);
1416 Tvmult_add<true>(n, vv, orthogonal_vectors, h);
1420 for (
unsigned int i = 0; i < n - 1; ++i)
1421 tmp += h(n + i) * h(n + i);
1422 const double alpha_j = h(n + n - 1) > tmp ?
1425 h(n + n - 1) = alpha_j;
1428 for (
unsigned int i = 0; i < n - 1; ++i)
1429 tmp += h(i) * h(n + i);
1430 h(n - 1) = (h(n - 1) - tmp) / alpha_j;
1435 for (
unsigned int i = 0; i < n - 1; ++i)
1436 hessenberg_matrix(i, n - 2) += h(n + i) * previous_scaling;
1437 hessenberg_matrix(n - 1, n - 2) = alpha_j * previous_scaling;
1439 for (
unsigned int i = 0; i < n; ++i)
1442 for (
unsigned int j = (i == 0 ? 0 : i - 1); j < n - 1; ++j)
1443 sum += hessenberg_matrix(i, j) * h(n + j);
1444 hessenberg_matrix(i, n - 1) = (h(i) -
sum) / alpha_j;
1450 for (
unsigned int i = 0; i < n - 1; ++i)
1452 sum += (2. - 1.) * h(n - 1) * h(n - 1);
1453 hessenberg_matrix(n, n - 1) =
1460 h(n + n) = hessenberg_matrix(n, n - 1);
1461 subtract_and_norm<true>(n, orthogonal_vectors, h, vv);
1465 residual_estimate = do_givens_rotation(
1466 true, n - 2, triangular_matrix, givens_rotations, projected_rhs);
1471 double norm_vv = 0.0;
1472 double norm_vv_start = 0;
1473 const bool consider_reorthogonalize =
1474 (do_reorthogonalization ==
false) && (n % 5 == 0);
1475 if (consider_reorthogonalize)
1476 norm_vv_start = vv.l2_norm();
1482 for (
unsigned int c = 0; c < 2; ++c)
1485 if (orthogonalization_strategy ==
1489 double htmp = vv * orthogonal_vectors[0];
1491 for (
unsigned int i = 1; i < n; ++i)
1493 htmp = vv.add_and_dot(-htmp,
1494 orthogonal_vectors[i - 1],
1495 orthogonal_vectors[i]);
1500 vv.add_and_dot(-htmp, orthogonal_vectors[n - 1], vv));
1502 else if (orthogonalization_strategy ==
1506 Tvmult_add<false>(n, vv, orthogonal_vectors, h);
1508 subtract_and_norm<false>(n, orthogonal_vectors, h, vv);
1526 if (consider_reorthogonalize)
1529 10. * norm_vv_start *
1531 typename VectorType::value_type>::epsilon()))
1536 do_reorthogonalization =
true;
1537 if (!reorthogonalize_signal.empty())
1538 reorthogonalize_signal(accumulated_iterations);
1542 if (do_reorthogonalization ==
false)
1546 for (
unsigned int i = 0; i < n; ++i)
1547 hessenberg_matrix(i, n - 1) = h(i);
1548 hessenberg_matrix(n, n - 1) = norm_vv;
1555 residual_estimate = do_givens_rotation(
1556 false, n - 1, triangular_matrix, givens_rotations, projected_rhs);
1559 return residual_estimate;
1565 ArnoldiProcess::do_givens_rotation(
1566 const bool delayed_reorthogonalization,
1569 std::vector<std::pair<double, double>> &rotations,
1577 if (delayed_reorthogonalization)
1582 matrix(0, col) = hessenberg_matrix(0, col);
1584 double H_next = hessenberg_matrix(0, col + 1);
1585 for (
int i = 0; i < col; ++i)
1587 const double c = rotations[i].first;
1588 const double s = rotations[i].second;
1589 const double Hi =
matrix(i, col);
1590 const double Hi1 = hessenberg_matrix(i + 1, col);
1591 H_next = -s * H_next + c * hessenberg_matrix(i + 1, col + 1);
1592 matrix(i, col) = c * Hi + s * Hi1;
1593 matrix(i + 1, col) = -s * Hi + c * Hi1;
1598 const double H_col1 = hessenberg_matrix(col + 1, col);
1599 const double H_col =
matrix(col, col);
1600 const double r = 1. /
std::sqrt(H_col * H_col + H_col1 * H_col1);
1601 rotations.emplace_back(H_col * r, H_col1 * r);
1603 rotations[col].first * H_col + rotations[col].second * H_col1;
1605 rhs(col + 1) = -rotations[col].second * rhs(col);
1606 rhs(col) *= rotations[col].first;
1609 -rotations[col].second * H_next +
1610 rotations[col].first * hessenberg_matrix(col + 1, col + 1);
1613 const double H_last = hessenberg_matrix(col + 2, col + 1);
1614 const double r = 1. /
std::sqrt(H_next * H_next + H_last * H_last);
1615 return std::abs(H_last * r * rhs(col + 1));
1621 matrix(0, col) = hessenberg_matrix(0, col);
1622 for (
int i = 0; i < col; ++i)
1624 const double c = rotations[i].first;
1625 const double s = rotations[i].second;
1626 const double Hi =
matrix(i, col);
1627 const double Hi1 = hessenberg_matrix(i + 1, col);
1628 matrix(i, col) = c * Hi + s * Hi1;
1629 matrix(i + 1, col) = -s * Hi + c * Hi1;
1632 const double Hi =
matrix(col, col);
1633 const double Hi1 = hessenberg_matrix(col + 1, col);
1634 const double r = 1. /
std::sqrt(Hi * Hi + Hi1 * Hi1);
1635 rotations.emplace_back(Hi * r, Hi1 * r);
1637 rotations[col].first * Hi + rotations[col].second * Hi1;
1639 rhs(col + 1) = -rotations[col].second * rhs(col);
1640 rhs(col) *= rotations[col].first;
1649 ArnoldiProcess::solve_projected_system(
1650 const bool orthogonalization_finished)
1656 unsigned int n = givens_rotations.
size();
1665 if (orthogonalization_strategy ==
1670 if (!orthogonalization_finished)
1672 tmp_triangular_matrix = triangular_matrix;
1673 tmp_rhs = projected_rhs;
1674 std::vector<std::pair<double, double>> tmp_givens_rotations(
1676 do_givens_rotation(
false,
1677 givens_rotations.size(),
1678 tmp_triangular_matrix,
1679 tmp_givens_rotations,
1681 matrix = &tmp_triangular_matrix;
1685 do_givens_rotation(
false,
1686 givens_rotations.size(),
1693 projected_solution.reinit(n);
1694 for (
int i = n - 1; i >= 0; --i)
1696 double s = (*rhs)(i);
1697 for (
unsigned int j = i + 1; j < n; ++j)
1698 s -= projected_solution(j) * (*matrix)(i, j);
1699 projected_solution(i) = s / (*matrix)(i, i);
1703 return projected_solution;
1709 ArnoldiProcess::get_hessenberg_matrix()
const
1711 return hessenberg_matrix;
1718 complex_less_pred(
const std::complex<double> &x,
1719 const std::complex<double> &y)
1721 return x.real() < y.real() ||
1722 (x.real() == y.real() && x.imag() < y.imag());
1729template <
typename VectorType>
1733 const unsigned int n,
1734 const boost::signals2::signal<
void(
const std::vector<std::complex<double>> &)>
1735 &eigenvalues_signal,
1738 const boost::signals2::signal<
void(
double)> &cond_signal)
1741 if ((!eigenvalues_signal.empty() || !hessenberg_signal.empty() ||
1742 !cond_signal.empty()) &&
1746 for (
unsigned int i = 0; i < n; ++i)
1747 for (
unsigned int j = 0; j < n; ++j)
1748 mat(i, j) = H_orig(i, j);
1749 hessenberg_signal(H_orig);
1751 if (!eigenvalues_signal.empty())
1757 mat_eig.compute_eigenvalues();
1759 for (
unsigned int i = 0; i < mat_eig.n(); ++i)
1764 internal::SolverGMRESImplementation::complex_less_pred);
1769 if (!cond_signal.empty() && (mat.n() > 1))
1772 double condition_number =
1773 mat.singular_value(0) / mat.singular_value(mat.n() - 1);
1774 cond_signal(condition_number);
1781template <
typename VectorType>
1783template <
typename MatrixType,
typename PreconditionerType>
1789 const VectorType &b,
1790 const PreconditionerType &preconditioner)
1792 std::unique_ptr<LogStream::Prefix> prefix;
1793 if (!additional_data.batched_mode)
1794 prefix = std::make_unique<LogStream::Prefix>(
"GMRES");
1798 const unsigned int basis_size =
1801 std::max(additional_data.max_n_tmp_vectors, 3u) - 2);
1805 basis_size + 2, this->memory);
1809 unsigned int accumulated_iterations = 0;
1811 const bool do_eigenvalues =
1812 !additional_data.batched_mode &&
1813 (!condition_number_signal.empty() ||
1814 !all_condition_numbers_signal.empty() || !eigenvalues_signal.empty() ||
1815 !all_eigenvalues_signal.empty() || !hessenberg_signal.empty() ||
1816 !all_hessenberg_signal.empty());
1819 double res = std::numeric_limits<double>::lowest();
1823 const bool left_precondition = !additional_data.right_preconditioning;
1828 const bool use_default_residual = additional_data.use_default_residual;
1831 VectorType &p = basis_vectors(basis_size + 1, x);
1837 if (!use_default_residual)
1847 additional_data.force_re_orthogonalization);
1855 VectorType &v = basis_vectors(0, x);
1864 if (left_precondition)
1866 if (accumulated_iterations == 0 && x.all_zero())
1867 preconditioner.vmult(v, b);
1872 preconditioner.vmult(v, p);
1877 if (accumulated_iterations == 0 && x.all_zero())
1886 const double norm_v =
1887 arnoldi_process.orthonormalize_nth_vector(0,
1889 accumulated_iterations);
1894 if (use_default_residual)
1897 if (additional_data.batched_mode)
1898 iteration_state = solver_control.
check(accumulated_iterations, res);
1901 this->iteration_status(accumulated_iterations, res, x);
1908 deallog <<
"default_res=" << norm_v << std::endl;
1910 if (left_precondition)
1913 r->sadd(-1., 1., b);
1916 preconditioner.vmult(*r, v);
1919 if (additional_data.batched_mode)
1920 iteration_state = solver_control.
check(accumulated_iterations, res);
1923 this->iteration_status(accumulated_iterations, res, x);
1931 unsigned int inner_iteration = 0;
1932 for (; (inner_iteration < basis_size &&
1936 ++accumulated_iterations;
1938 VectorType &vv = basis_vectors(inner_iteration + 1, x);
1940 if (left_precondition)
1942 A.vmult(p, basis_vectors[inner_iteration]);
1943 preconditioner.vmult(vv, p);
1947 preconditioner.vmult(p, basis_vectors[inner_iteration]);
1952 arnoldi_process.orthonormalize_nth_vector(inner_iteration + 1,
1954 accumulated_iterations,
1955 re_orthogonalize_signal);
1957 if (use_default_residual)
1959 if (additional_data.batched_mode)
1961 solver_control.
check(accumulated_iterations, res);
1964 this->iteration_status(accumulated_iterations, res, x);
1968 if (!additional_data.batched_mode)
1969 deallog <<
"default_res=" << res << std::endl;
1973 arnoldi_process.solve_projected_system(
false);
1975 if (left_precondition)
1976 for (
unsigned int i = 0; i < inner_iteration + 1; ++i)
1977 x_->add(projected_solution(i), basis_vectors[i]);
1981 for (
unsigned int i = 0; i < inner_iteration + 1; ++i)
1982 p.add(projected_solution(i), basis_vectors[i]);
1983 preconditioner.vmult(*r, p);
1987 r->sadd(-1., 1., b);
1990 if (left_precondition)
1994 this->iteration_status(accumulated_iterations, res, x);
1998 preconditioner.vmult(*x_, *r);
1999 res = x_->l2_norm();
2001 if (additional_data.batched_mode)
2003 solver_control.
check(accumulated_iterations, res);
2006 this->iteration_status(accumulated_iterations, res, x);
2014 arnoldi_process.solve_projected_system(
true);
2017 compute_eigs_and_cond(arnoldi_process.get_hessenberg_matrix(),
2019 all_eigenvalues_signal,
2020 all_hessenberg_signal,
2021 condition_number_signal);
2023 if (left_precondition)
2024 ::internal::SolverGMRESImplementation::add(
2025 x, inner_iteration, projected_solution, basis_vectors,
false);
2028 ::internal::SolverGMRESImplementation::add(
2029 p, inner_iteration, projected_solution, basis_vectors,
true);
2030 preconditioner.vmult(v, p);
2038 compute_eigs_and_cond(arnoldi_process.get_hessenberg_matrix(),
2042 condition_number_signal);
2044 if (!additional_data.batched_mode && !krylov_space_signal.empty())
2045 krylov_space_signal(basis_vectors);
2060template <
typename VectorType>
2062boost::signals2::connection
2064 const std::function<
void(
double)> &slot,
2065 const bool every_iteration)
2067 if (every_iteration)
2069 return all_condition_numbers_signal.connect(slot);
2073 return condition_number_signal.connect(slot);
2079template <
typename VectorType>
2082 const std::function<
void(
const std::vector<std::complex<double>> &)> &slot,
2083 const bool every_iteration)
2085 if (every_iteration)
2087 return all_eigenvalues_signal.connect(slot);
2091 return eigenvalues_signal.connect(slot);
2097template <
typename VectorType>
2101 const bool every_iteration)
2103 if (every_iteration)
2105 return all_hessenberg_signal.connect(slot);
2109 return hessenberg_signal.connect(slot);
2115template <
typename VectorType>
2118 const std::function<
void(
2121 return krylov_space_signal.connect(slot);
2126template <
typename VectorType>
2128boost::signals2::connection
2130 const std::function<
void(
int)> &slot)
2132 return re_orthogonalize_signal.connect(slot);
2137template <
typename VectorType>
2150template <
typename VectorType>
2154 const AdditionalData &data)
2156 , additional_data(data)
2161template <
typename VectorType>
2164 const AdditionalData &data)
2166 , additional_data(data)
2171template <
typename VectorType>
2173template <
typename MatrixType,
typename PreconditionerType>
2179 const VectorType &b,
2180 const PreconditionerType &preconditioner)
2190 basis_size + 1, this->memory);
2192 basis_size, this->memory);
2196 unsigned int accumulated_iterations = 0;
2204 double res = std::numeric_limits<double>::lowest();
2214 if (accumulated_iterations == 0 && x.all_zero())
2218 A.vmult(v(0, x), x);
2219 v[0].sadd(-1., 1., b);
2222 res = arnoldi_process.orthonormalize_nth_vector(0, v);
2223 iteration_state = this->iteration_status(accumulated_iterations, res, x);
2227 unsigned int inner_iteration = 0;
2228 for (; (inner_iteration < basis_size &&
2232 preconditioner.vmult(z(inner_iteration, x), v[inner_iteration]);
2233 A.vmult(v(inner_iteration + 1, x), z[inner_iteration]);
2236 arnoldi_process.orthonormalize_nth_vector(inner_iteration + 1, v);
2243 this->iteration_status(++accumulated_iterations, res, x);
2249 arnoldi_process.solve_projected_system(
true);
2250 ::internal::SolverGMRESImplementation::add(
2251 x, inner_iteration, projected_solution, z,
false);
virtual State check(const unsigned int step, const double check_value)
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
SolverFGMRES(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
SolverFGMRES(SolverControl &cn, const AdditionalData &data=AdditionalData())
SolverGMRES(const SolverGMRES< VectorType > &)=delete
boost::signals2::connection connect_condition_number_slot(const std::function< void(double)> &slot, const bool every_iteration=false)
boost::signals2::connection connect_krylov_space_slot(const std::function< void(const internal::SolverGMRESImplementation::TmpVectors< VectorType > &)> &slot)
boost::signals2::connection connect_re_orthogonalization_slot(const std::function< void(int)> &slot)
static void compute_eigs_and_cond(const FullMatrix< double > &H_orig, const unsigned int n, const boost::signals2::signal< void(const std::vector< std::complex< double > > &)> &eigenvalues_signal, const boost::signals2::signal< void(const FullMatrix< double > &)> &hessenberg_signal, const boost::signals2::signal< void(double)> &cond_signal)
SolverGMRES(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
boost::signals2::connection connect_hessenberg_slot(const std::function< void(const FullMatrix< double > &)> &slot, const bool every_iteration=true)
SolverGMRES(SolverControl &cn, const AdditionalData &data=AdditionalData())
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
virtual double criterion()
boost::signals2::connection connect_eigenvalues_slot(const std::function< void(const std::vector< std::complex< double > > &)> &slot, const bool every_iteration=false)
virtual size_type size() const override
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
static constexpr std::size_t size()
void store(OtherNumber *ptr) const
void load(const OtherNumber *ptr)
void initialize(const LinearAlgebra::OrthogonalizationStrategy orthogonalization_strategy, const unsigned int max_basis_size, const bool force_reorthogonalization)
double orthonormalize_nth_vector(const unsigned int n, TmpVectors< VectorType > &orthogonal_vectors, const unsigned int accumulated_iterations=0, const boost::signals2::signal< void(int)> &reorthogonalize_signal=boost::signals2::signal< void(int)>())
Vector< double > projected_rhs
FullMatrix< double > hessenberg_matrix
Vector< double > projected_solution
LinearAlgebra::OrthogonalizationStrategy orthogonalization_strategy
std::vector< std::pair< double, double > > givens_rotations
bool do_reorthogonalization
const Vector< double > & solve_projected_system(const bool orthogonalization_finished)
FullMatrix< double > triangular_matrix
double do_givens_rotation(const bool delayed_reorthogonalization, const int col, FullMatrix< double > &matrix, std::vector< std::pair< double, double > > &rotations, Vector< double > &rhs)
const FullMatrix< double > & get_hessenberg_matrix() const
VectorMemory< VectorType > & mem
std::vector< typename VectorMemory< VectorType >::Pointer > data
TmpVectors(const unsigned int max_size, VectorMemory< VectorType > &vmem)
unsigned int size() const
VectorType & operator[](const unsigned int i) const
VectorType & operator()(const unsigned int i, const VectorType &temp)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertIsFinite(number)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNotInitialized()
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
#define DEAL_II_ASSERT_UNREACHABLE()
@ matrix
Contents is actually a matrix.
OrthogonalizationStrategy
@ delayed_classical_gram_schmidt
std::enable_if_t< IsBlockVector< VectorType >::value, unsigned int > n_blocks(const VectorType &vector)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
VectorType::value_type * begin(VectorType &V)
T sum(const T &t, const MPI_Comm mpi_communicator)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
AdditionalData(const unsigned int max_basis_size=30, const LinearAlgebra::OrthogonalizationStrategy orthogonalization_strategy=LinearAlgebra::OrthogonalizationStrategy::delayed_classical_gram_schmidt)
unsigned int max_basis_size
LinearAlgebra::OrthogonalizationStrategy orthogonalization_strategy
AdditionalData(const unsigned int max_basis_size=30, const bool right_preconditioning=false, const bool use_default_residual=true, const bool force_re_orthogonalization=false, const bool batched_mode=false, const LinearAlgebra::OrthogonalizationStrategy orthogonalization_strategy=LinearAlgebra::OrthogonalizationStrategy::delayed_classical_gram_schmidt)
bool right_preconditioning
bool force_re_orthogonalization
LinearAlgebra::OrthogonalizationStrategy orthogonalization_strategy
bool use_default_residual
unsigned int max_basis_size
unsigned int max_n_tmp_vectors
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)