16 #ifndef dealii_solver_gmres_h
17 #define dealii_solver_gmres_h
48 namespace SolverGMRESImplementation
58 template <
typename VectorType>
105 std::vector<typename VectorMemory<VectorType>::Pointer>
data;
177 template <
class VectorType = Vector<
double>>
249 template <
typename MatrixType,
typename PreconditionerType>
251 solve(
const MatrixType &
A,
254 const PreconditionerType &preconditioner);
262 boost::signals2::connection
264 const bool every_iteration =
false);
272 boost::signals2::connection
274 const std::function<
void(
const std::vector<std::complex<double>> &)> &slot,
275 const bool every_iteration =
false);
284 boost::signals2::connection
287 const bool every_iteration =
true);
295 boost::signals2::connection
306 boost::signals2::connection
312 <<
"The number of temporary vectors you gave (" << arg1
313 <<
") is too small. It should be at least 10 for "
314 <<
"any results, and much more for reasonable ones.");
338 boost::signals2::signal<void(
const std::vector<std::complex<double>> &)>
345 boost::signals2::signal<void(
const std::vector<std::complex<double>> &)>
365 boost::signals2::signal<void(
405 & orthogonal_vectors,
406 const unsigned int dim,
407 const unsigned int accumulated_iterations,
410 bool & re_orthogonalize,
412 boost::signals2::signal<
void(
int)>());
423 const unsigned int dim,
424 const boost::signals2::signal<
428 const boost::signals2::signal<
void(
double)> &cond_signal);
463 template <
class VectorType = Vector<
double>>
485 const bool use_default_residual)
488 (void)use_default_residual;
514 template <
typename MatrixType,
typename PreconditionerType>
516 solve(
const MatrixType &
A,
519 const PreconditionerType &preconditioner);
545 namespace SolverGMRESImplementation
547 template <
class VectorType>
556 template <
class VectorType>
568 template <
class VectorType>
574 if (data[i] ==
nullptr)
577 data[i]->reinit(temp);
584 template <
class VectorType>
588 return (data.size() > 0 ? data.size() - 1 : 0);
595 complex_less_pred(
const std::complex<double> &x,
596 const std::complex<double> &y)
598 return x.real() < y.real() ||
599 (x.real() == y.real() && x.imag() < y.imag());
606 template <
class VectorType>
608 const unsigned int max_n_tmp_vectors,
609 const bool right_preconditioning,
610 const bool use_default_residual,
611 const bool force_re_orthogonalization)
612 : max_n_tmp_vectors(max_n_tmp_vectors)
613 , right_preconditioning(right_preconditioning)
614 , use_default_residual(use_default_residual)
615 , force_re_orthogonalization(force_re_orthogonalization)
617 Assert(3 <= max_n_tmp_vectors,
618 ExcMessage(
"SolverGMRES needs at least three "
619 "temporary vectors."));
624 template <
class VectorType>
627 const AdditionalData & data)
629 , additional_data(data)
634 template <
class VectorType>
636 const AdditionalData &data)
638 , additional_data(data)
643 template <
class VectorType>
651 for (
int i = 0; i < col; i++)
653 const double s = si(i);
654 const double c = ci(i);
655 const double dummy = h(i);
656 h(i) = c *
dummy + s * h(i + 1);
657 h(i + 1) = -s *
dummy + c * h(i + 1);
660 const double r = 1. /
std::sqrt(h(col) * h(col) + h(col + 1) * h(col + 1));
661 si(col) = h(col + 1) * r;
662 ci(col) = h(col) * r;
663 h(col) = ci(col) * h(col) + si(col) * h(col + 1);
664 b(col + 1) = -si(col) *
b(col);
670 template <
class VectorType>
674 & orthogonal_vectors,
675 const unsigned int dim,
676 const unsigned int accumulated_iterations,
679 bool & reorthogonalize,
680 const boost::signals2::signal<
void(
int)> &reorthogonalize_signal)
683 const unsigned int inner_iteration = dim - 1;
686 double norm_vv_start = 0;
687 const bool consider_reorthogonalize =
688 (reorthogonalize ==
false) && (inner_iteration % 5 == 4);
689 if (consider_reorthogonalize)
690 norm_vv_start = vv.l2_norm();
693 h(0) = vv * orthogonal_vectors[0];
694 for (
unsigned int i = 1; i < dim; ++i)
695 h(i) = vv.add_and_dot(-h(i - 1),
696 orthogonal_vectors[i - 1],
697 orthogonal_vectors[i]);
699 std::sqrt(vv.add_and_dot(-h(dim - 1), orthogonal_vectors[dim - 1], vv));
708 if (consider_reorthogonalize)
711 10. * norm_vv_start *
718 reorthogonalize =
true;
719 if (!reorthogonalize_signal.empty())
720 reorthogonalize_signal(accumulated_iterations);
724 if (reorthogonalize ==
true)
726 double htmp = vv * orthogonal_vectors[0];
728 for (
unsigned int i = 1; i < dim; ++i)
730 htmp = vv.add_and_dot(-htmp,
731 orthogonal_vectors[i - 1],
732 orthogonal_vectors[i]);
736 std::sqrt(vv.add_and_dot(-htmp, orthogonal_vectors[dim - 1], vv));
744 template <
class VectorType>
748 const unsigned int dim,
749 const boost::signals2::signal<
void(
const std::vector<std::complex<double>> &)>
753 const boost::signals2::signal<
void(
double)> &cond_signal)
756 if ((!eigenvalues_signal.empty() || !hessenberg_signal.empty() ||
757 !cond_signal.empty()) &&
761 for (
unsigned int i = 0; i < dim; ++i)
762 for (
unsigned int j = 0; j < dim; ++j)
763 mat(i, j) = H_orig(i, j);
764 hessenberg_signal(H_orig);
766 if (!eigenvalues_signal.empty())
772 mat_eig.compute_eigenvalues();
773 std::vector<std::complex<double>>
eigenvalues(dim);
774 for (
unsigned int i = 0; i < mat_eig.n(); ++i)
779 internal::SolverGMRESImplementation::complex_less_pred);
784 if (!cond_signal.empty() && (mat.n() > 1))
787 double condition_number =
788 mat.singular_value(0) / mat.singular_value(mat.n() - 1);
789 cond_signal(condition_number);
796 template <
class VectorType>
797 template <
typename MatrixType,
typename PreconditionerType>
802 const PreconditionerType &preconditioner)
812 const unsigned int n_tmp_vectors =
813 std::max(additional_data.max_n_tmp_vectors, 3u);
817 n_tmp_vectors, this->memory);
822 unsigned int accumulated_iterations = 0;
824 const bool do_eigenvalues =
825 !condition_number_signal.empty() || !all_condition_numbers_signal.empty() ||
826 !eigenvalues_signal.empty() || !all_eigenvalues_signal.empty() ||
827 !hessenberg_signal.empty() || !all_hessenberg_signal.empty();
832 H_orig.
reinit(n_tmp_vectors, n_tmp_vectors - 1);
835 H.reinit(n_tmp_vectors, n_tmp_vectors - 1);
839 si(n_tmp_vectors - 1), h(n_tmp_vectors - 1);
842 unsigned int dim = 0;
849 const bool left_precondition = !additional_data.right_preconditioning;
854 const bool use_default_residual = additional_data.use_default_residual;
858 VectorType &p = tmp_vectors(n_tmp_vectors - 1, x);
864 std::unique_ptr<::Vector<double>> gamma_;
865 if (!use_default_residual)
872 gamma_ = std_cxx14::make_unique<::Vector<double>>(
gamma.size());
875 bool re_orthogonalize = additional_data.force_re_orthogonalization;
884 h.
reinit(n_tmp_vectors - 1);
886 if (left_precondition)
890 preconditioner.vmult(v, p);
898 double rho = v.l2_norm();
903 if (use_default_residual)
907 this->iteration_status(accumulated_iterations, rho, x);
914 deallog <<
"default_res=" << rho << std::endl;
916 if (left_precondition)
922 preconditioner.vmult(*r, v);
924 double res = r->l2_norm();
927 this->iteration_status(accumulated_iterations, res, x);
940 for (
unsigned int inner_iteration = 0;
941 ((inner_iteration < n_tmp_vectors - 2) &&
945 ++accumulated_iterations;
947 VectorType &vv = tmp_vectors(inner_iteration + 1, x);
949 if (left_precondition)
951 A.vmult(p, tmp_vectors[inner_iteration]);
952 preconditioner.vmult(vv, p);
956 preconditioner.vmult(p, tmp_vectors[inner_iteration]);
960 dim = inner_iteration + 1;
962 const double s = modified_gram_schmidt(tmp_vectors,
964 accumulated_iterations,
968 re_orthogonalize_signal);
969 h(inner_iteration + 1) = s;
979 for (
unsigned int i = 0; i < dim + 1; ++i)
980 H_orig(i, inner_iteration) = h(i);
986 for (
unsigned int i = 0; i < dim; ++i)
987 H(i, inner_iteration) = h(i);
992 if (use_default_residual)
996 this->iteration_status(accumulated_iterations, rho, x);
1000 deallog <<
"default_res=" << rho << std::endl;
1005 H1.reinit(dim + 1, dim);
1007 for (
unsigned int i = 0; i < dim + 1; ++i)
1008 for (
unsigned int j = 0; j < dim; ++j)
1011 H1.backward(h_, *gamma_);
1013 if (left_precondition)
1014 for (
unsigned int i = 0; i < dim; ++i)
1015 x_->add(h_(i), tmp_vectors[i]);
1019 for (
unsigned int i = 0; i < dim; ++i)
1020 p.add(h_(i), tmp_vectors[i]);
1021 preconditioner.vmult(*r, p);
1025 r->sadd(-1., 1.,
b);
1027 if (left_precondition)
1029 const double res = r->l2_norm();
1033 this->iteration_status(accumulated_iterations, res, x);
1037 preconditioner.vmult(*x_, *r);
1038 const double preconditioned_res = x_->l2_norm();
1039 last_res = preconditioned_res;
1042 this->iteration_status(accumulated_iterations,
1051 H1.reinit(dim + 1, dim);
1053 for (
unsigned int i = 0; i < dim + 1; ++i)
1054 for (
unsigned int j = 0; j < dim; ++j)
1057 compute_eigs_and_cond(H_orig,
1059 all_eigenvalues_signal,
1060 all_hessenberg_signal,
1061 condition_number_signal);
1063 H1.backward(h,
gamma);
1065 if (left_precondition)
1066 for (
unsigned int i = 0; i < dim; ++i)
1067 x.add(h(i), tmp_vectors[i]);
1071 for (
unsigned int i = 0; i < dim; ++i)
1072 p.add(h(i), tmp_vectors[i]);
1073 preconditioner.vmult(v, p);
1081 compute_eigs_and_cond(H_orig,
1085 condition_number_signal);
1087 if (!krylov_space_signal.empty())
1088 krylov_space_signal(tmp_vectors);
1097 template <
class VectorType>
1098 boost::signals2::connection
1100 const std::function<
void(
double)> &slot,
1101 const bool every_iteration)
1103 if (every_iteration)
1105 return all_condition_numbers_signal.connect(slot);
1109 return condition_number_signal.connect(slot);
1115 template <
class VectorType>
1116 boost::signals2::connection
1118 const std::function<
void(
const std::vector<std::complex<double>> &)> &slot,
1119 const bool every_iteration)
1121 if (every_iteration)
1123 return all_eigenvalues_signal.connect(slot);
1127 return eigenvalues_signal.connect(slot);
1133 template <
class VectorType>
1134 boost::signals2::connection
1137 const bool every_iteration)
1139 if (every_iteration)
1141 return all_hessenberg_signal.connect(slot);
1145 return hessenberg_signal.connect(slot);
1151 template <
class VectorType>
1152 boost::signals2::connection
1154 const std::function<
void(
1157 return krylov_space_signal.connect(slot);
1162 template <
class VectorType>
1163 boost::signals2::connection
1165 const std::function<
void(
int)> &slot)
1167 return re_orthogonalize_signal.connect(slot);
1172 template <
class VectorType>
1185 template <
class VectorType>
1188 const AdditionalData & data)
1190 , additional_data(data)
1195 template <
class VectorType>
1197 const AdditionalData &data)
1199 , additional_data(data)
1204 template <
class VectorType>
1205 template <
typename MatrixType,
typename PreconditionerType>
1210 const PreconditionerType &preconditioner)
1216 const unsigned int basis_size = additional_data.max_basis_size;
1220 basis_size, this->memory);
1222 basis_size, this->memory);
1226 unsigned int accumulated_iterations = 0;
1229 H.reinit(basis_size + 1, basis_size);
1243 aux->sadd(-1., 1.,
b);
1245 double beta = aux->l2_norm();
1247 iteration_state = this->iteration_status(accumulated_iterations, res, x);
1251 H.reinit(basis_size + 1, basis_size);
1254 for (
unsigned int j = 0; j < basis_size; ++j)
1257 v(j, x).equ(1. / a, *aux);
1262 preconditioner.vmult(z(j, x), v[j]);
1263 A.vmult(*aux, z[j]);
1266 H(0, j) = *aux * v[0];
1267 for (
unsigned int i = 1; i <= j; ++i)
1268 H(i, j) = aux->add_and_dot(-H(i - 1, j), v[i - 1], v[i]);
1269 H(j + 1, j) = a =
std::sqrt(aux->add_and_dot(-H(j, j), v[j], *aux));
1275 H1.reinit(j + 1, j);
1276 projected_rhs.
reinit(j + 1);
1278 projected_rhs(0) = beta;
1286 res = house.least_squares(y, projected_rhs);
1288 this->iteration_status(++accumulated_iterations, res, x);
1295 for (
unsigned int j = 0; j < y.
size(); ++j)