16#ifndef dealii_solver_gmres_h
17#define dealii_solver_gmres_h
48 namespace SolverGMRESImplementation
58 template <
typename VectorType>
87 operator()(
const unsigned int i,
const VectorType &temp);
106 std::vector<typename VectorMemory<VectorType>::Pointer>
data;
175template <
class VectorType = Vector<
double>>
247 template <
typename MatrixType,
typename PreconditionerType>
251 const VectorType & b,
252 const PreconditionerType &preconditioner);
260 boost::signals2::connection
262 const bool every_iteration =
false);
270 boost::signals2::connection
272 const std::function<
void(
const std::vector<std::complex<double>> &)> &slot,
273 const bool every_iteration =
false);
282 boost::signals2::connection
285 const bool every_iteration =
true);
293 boost::signals2::connection
304 boost::signals2::connection
310 <<
"The number of temporary vectors you gave (" << arg1
311 <<
") is too small. It should be at least 10 for "
312 <<
"any results, and much more for reasonable ones.");
336 boost::signals2::signal<void(
const std::vector<std::complex<double>> &)>
343 boost::signals2::signal<void(
const std::vector<std::complex<double>> &)>
363 boost::signals2::signal<void(
403 & orthogonal_vectors,
404 const unsigned int dim,
405 const unsigned int accumulated_iterations,
408 bool & re_orthogonalize,
410 boost::signals2::signal<
void(
int)>());
421 const unsigned int dim,
422 const boost::signals2::signal<
426 const boost::signals2::signal<
void(
double)> &cond_signal);
459template <
class VectorType = Vector<
double>>
498 template <
typename MatrixType,
typename PreconditionerType>
502 const VectorType & b,
503 const PreconditionerType &preconditioner);
529 namespace SolverGMRESImplementation
531 template <
class VectorType>
540 template <
class VectorType>
552 template <
class VectorType>
555 const VectorType & temp)
558 if (data[i] ==
nullptr)
561 data[i]->reinit(temp);
568 template <
class VectorType>
572 return (data.size() > 0 ? data.size() - 1 : 0);
579 complex_less_pred(
const std::complex<double> &x,
580 const std::complex<double> &y)
582 return x.real() < y.real() ||
583 (x.real() == y.real() && x.imag() < y.imag());
590template <
class VectorType>
592 const unsigned int max_n_tmp_vectors,
593 const bool right_preconditioning,
594 const bool use_default_residual,
595 const bool force_re_orthogonalization)
596 : max_n_tmp_vectors(max_n_tmp_vectors)
597 , right_preconditioning(right_preconditioning)
598 , use_default_residual(use_default_residual)
599 , force_re_orthogonalization(force_re_orthogonalization)
601 Assert(3 <= max_n_tmp_vectors,
602 ExcMessage(
"SolverGMRES needs at least three "
603 "temporary vectors."));
608template <
class VectorType>
611 const AdditionalData & data)
613 , additional_data(data)
618template <
class VectorType>
620 const AdditionalData &data)
622 , additional_data(data)
627template <
class VectorType>
635 for (
int i = 0; i < col; ++i)
637 const double s = si(i);
638 const double c = ci(i);
639 const double dummy = h(i);
640 h(i) = c * dummy + s * h(i + 1);
641 h(i + 1) = -s * dummy + c * h(i + 1);
644 const double r = 1. /
std::sqrt(h(col) * h(col) + h(col + 1) * h(col + 1));
645 si(col) = h(col + 1) * r;
646 ci(col) = h(col) * r;
647 h(col) = ci(col) * h(col) + si(col) * h(col + 1);
648 b(col + 1) = -si(col) *
b(col);
654template <
class VectorType>
658 & orthogonal_vectors,
659 const unsigned int dim,
660 const unsigned int accumulated_iterations,
663 bool & reorthogonalize,
664 const boost::signals2::signal<
void(
int)> &reorthogonalize_signal)
667 const unsigned int inner_iteration = dim - 1;
670 double norm_vv_start = 0;
671 const bool consider_reorthogonalize =
672 (reorthogonalize ==
false) && (inner_iteration % 5 == 4);
673 if (consider_reorthogonalize)
674 norm_vv_start = vv.l2_norm();
677 h(0) = vv * orthogonal_vectors[0];
678 for (
unsigned int i = 1; i < dim; ++i)
680 orthogonal_vectors[i - 1],
681 orthogonal_vectors[i]);
683 std::sqrt(vv.add_and_dot(-h(dim - 1), orthogonal_vectors[dim - 1], vv));
692 if (consider_reorthogonalize)
695 10. * norm_vv_start *
697 std::numeric_limits<typename VectorType::value_type>::epsilon()))
702 reorthogonalize =
true;
703 if (!reorthogonalize_signal.empty())
704 reorthogonalize_signal(accumulated_iterations);
708 if (reorthogonalize ==
true)
710 double htmp = vv * orthogonal_vectors[0];
712 for (
unsigned int i = 1; i < dim; ++i)
715 orthogonal_vectors[i - 1],
716 orthogonal_vectors[i]);
720 std::sqrt(vv.add_and_dot(-htmp, orthogonal_vectors[dim - 1], vv));
728template <
class VectorType>
732 const unsigned int dim,
733 const boost::signals2::signal<
void(
const std::vector<std::complex<double>> &)>
737 const boost::signals2::signal<
void(
double)> &cond_signal)
740 if ((!eigenvalues_signal.empty() || !hessenberg_signal.empty() ||
741 !cond_signal.empty()) &&
745 for (
unsigned int i = 0; i < dim; ++i)
746 for (
unsigned int j = 0; j < dim; ++j)
747 mat(i, j) = H_orig(i, j);
748 hessenberg_signal(H_orig);
750 if (!eigenvalues_signal.empty())
756 mat_eig.compute_eigenvalues();
757 std::vector<std::complex<double>>
eigenvalues(dim);
758 for (
unsigned int i = 0; i < mat_eig.n(); ++i)
763 internal::SolverGMRESImplementation::complex_less_pred);
768 if (!cond_signal.empty() && (mat.n() > 1))
771 double condition_number =
772 mat.singular_value(0) / mat.singular_value(mat.n() - 1);
773 cond_signal(condition_number);
780template <
class VectorType>
781template <
typename MatrixType,
typename PreconditionerType>
785 const VectorType & b,
786 const PreconditionerType &preconditioner)
796 const unsigned int n_tmp_vectors =
801 n_tmp_vectors, this->memory);
806 unsigned int accumulated_iterations = 0;
808 const bool do_eigenvalues =
809 !condition_number_signal.empty() || !all_condition_numbers_signal.empty() ||
810 !eigenvalues_signal.empty() || !all_eigenvalues_signal.empty() ||
811 !hessenberg_signal.empty() || !all_hessenberg_signal.empty();
816 H_orig.reinit(n_tmp_vectors, n_tmp_vectors - 1);
819 H.reinit(n_tmp_vectors, n_tmp_vectors - 1);
823 si(n_tmp_vectors - 1), h(n_tmp_vectors - 1);
826 unsigned int dim = 0;
829 double last_res = std::numeric_limits<double>::lowest();
841 VectorType &v = tmp_vectors(0, x);
842 VectorType &p = tmp_vectors(n_tmp_vectors - 1, x);
848 std::unique_ptr<::Vector<double>> gamma_;
849 if (!use_default_residual)
856 gamma_ = std::make_unique<::Vector<double>>(
gamma.size());
868 h.
reinit(n_tmp_vectors - 1);
870 if (left_precondition)
874 preconditioner.vmult(v, p);
882 double rho = v.l2_norm();
887 if (use_default_residual)
891 this->iteration_status(accumulated_iterations, rho, x);
898 deallog <<
"default_res=" << rho << std::endl;
900 if (left_precondition)
906 preconditioner.vmult(*r, v);
908 double res = r->l2_norm();
911 this->iteration_status(accumulated_iterations, res, x);
924 for (
unsigned int inner_iteration = 0;
925 ((inner_iteration < n_tmp_vectors - 2) &&
929 ++accumulated_iterations;
931 VectorType &vv = tmp_vectors(inner_iteration + 1, x);
933 if (left_precondition)
935 A.vmult(p, tmp_vectors[inner_iteration]);
936 preconditioner.vmult(vv, p);
940 preconditioner.vmult(p, tmp_vectors[inner_iteration]);
944 dim = inner_iteration + 1;
946 const double s = modified_gram_schmidt(tmp_vectors,
948 accumulated_iterations,
952 re_orthogonalize_signal);
953 h(inner_iteration + 1) = s;
963 for (
unsigned int i = 0; i < dim + 1; ++i)
964 H_orig(i, inner_iteration) = h(i);
970 for (
unsigned int i = 0; i < dim; ++i)
971 H(i, inner_iteration) = h(i);
974 rho = std::fabs(
gamma(dim));
976 if (use_default_residual)
980 this->iteration_status(accumulated_iterations, rho, x);
984 deallog <<
"default_res=" << rho << std::endl;
989 H1.reinit(dim + 1, dim);
991 for (
unsigned int i = 0; i < dim + 1; ++i)
992 for (
unsigned int j = 0; j < dim; ++j)
995 H1.backward(h_, *gamma_);
997 if (left_precondition)
998 for (
unsigned int i = 0; i < dim; ++i)
999 x_->add(h_(i), tmp_vectors[i]);
1003 for (
unsigned int i = 0; i < dim; ++i)
1004 p.add(h_(i), tmp_vectors[i]);
1005 preconditioner.vmult(*r, p);
1009 r->sadd(-1., 1., b);
1011 if (left_precondition)
1013 const double res = r->l2_norm();
1017 this->iteration_status(accumulated_iterations, res, x);
1021 preconditioner.vmult(*x_, *r);
1022 const double preconditioned_res = x_->l2_norm();
1023 last_res = preconditioned_res;
1026 this->iteration_status(accumulated_iterations,
1035 H1.reinit(dim + 1, dim);
1037 for (
unsigned int i = 0; i < dim + 1; ++i)
1038 for (
unsigned int j = 0; j < dim; ++j)
1041 compute_eigs_and_cond(H_orig,
1043 all_eigenvalues_signal,
1044 all_hessenberg_signal,
1045 condition_number_signal);
1047 H1.backward(h, gamma);
1049 if (left_precondition)
1050 for (
unsigned int i = 0; i < dim; ++i)
1051 x.add(h(i), tmp_vectors[i]);
1055 for (
unsigned int i = 0; i < dim; ++i)
1056 p.add(h(i), tmp_vectors[i]);
1057 preconditioner.vmult(v, p);
1065 compute_eigs_and_cond(H_orig,
1069 condition_number_signal);
1071 if (!krylov_space_signal.empty())
1072 krylov_space_signal(tmp_vectors);
1081template <
class VectorType>
1082boost::signals2::connection
1084 const std::function<
void(
double)> &slot,
1085 const bool every_iteration)
1087 if (every_iteration)
1089 return all_condition_numbers_signal.connect(slot);
1093 return condition_number_signal.connect(slot);
1099template <
class VectorType>
1100boost::signals2::connection
1102 const std::function<
void(
const std::vector<std::complex<double>> &)> &slot,
1103 const bool every_iteration)
1105 if (every_iteration)
1107 return all_eigenvalues_signal.connect(slot);
1111 return eigenvalues_signal.connect(slot);
1117template <
class VectorType>
1118boost::signals2::connection
1121 const bool every_iteration)
1123 if (every_iteration)
1125 return all_hessenberg_signal.connect(slot);
1129 return hessenberg_signal.connect(slot);
1135template <
class VectorType>
1136boost::signals2::connection
1138 const std::function<
void(
1141 return krylov_space_signal.connect(slot);
1146template <
class VectorType>
1147boost::signals2::connection
1149 const std::function<
void(
int)> &slot)
1151 return re_orthogonalize_signal.connect(slot);
1156template <
class VectorType>
1169template <
class VectorType>
1172 const AdditionalData & data)
1174 , additional_data(data)
1179template <
class VectorType>
1181 const AdditionalData &data)
1183 , additional_data(data)
1188template <
class VectorType>
1189template <
typename MatrixType,
typename PreconditionerType>
1193 const VectorType & b,
1194 const PreconditionerType &preconditioner)
1200 const unsigned int basis_size = additional_data.max_basis_size;
1204 basis_size, this->memory);
1206 basis_size, this->memory);
1210 unsigned int accumulated_iterations = 0;
1213 H.reinit(basis_size + 1, basis_size);
1220 double res = std::numeric_limits<double>::lowest();
1227 aux->sadd(-1., 1., b);
1229 double beta = aux->l2_norm();
1231 iteration_state = this->iteration_status(accumulated_iterations, res, x);
1235 H.reinit(basis_size + 1, basis_size);
1238 for (
unsigned int j = 0; j < basis_size; ++j)
1241 v(j, x).equ(1. / a, *aux);
1246 preconditioner.vmult(z(j, x), v[j]);
1247 A.vmult(*aux, z[j]);
1250 H(0, j) = *aux * v[0];
1251 for (
unsigned int i = 1; i <= j; ++i)
1252 H(i, j) = aux->add_and_dot(-H(i - 1, j), v[i - 1], v[i]);
1253 H(j + 1, j) = a =
std::sqrt(aux->add_and_dot(-H(j, j), v[j], *aux));
1259 H1.reinit(j + 1, j);
1260 projected_rhs.
reinit(j + 1);
1262 projected_rhs(0) = beta;
1270 res = house.least_squares(y, projected_rhs);
1272 this->iteration_status(++accumulated_iterations, res, x);
1279 for (
unsigned int j = 0; j < y.
size(); ++j)
@ 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())
AdditionalData additional_data
boost::signals2::signal< void(const FullMatrix< double > &)> hessenberg_signal
AdditionalData additional_data
boost::signals2::signal< void(double)> condition_number_signal
SolverGMRES(const SolverGMRES< VectorType > &)=delete
boost::signals2::signal< void(const std::vector< std::complex< double > > &)> eigenvalues_signal
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::signal< void(const internal::SolverGMRESImplementation::TmpVectors< VectorType > &)> krylov_space_signal
boost::signals2::signal< void(const FullMatrix< double > &)> all_hessenberg_signal
boost::signals2::connection connect_re_orthogonalization_slot(const std::function< void(int)> &slot)
SolverGMRES(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
void givens_rotation(Vector< double > &h, Vector< double > &b, Vector< double > &ci, Vector< double > &si, int col) const
boost::signals2::connection connect_hessenberg_slot(const std::function< void(const FullMatrix< double > &)> &slot, const bool every_iteration=true)
boost::signals2::signal< void(const std::vector< std::complex< double > > &)> all_eigenvalues_signal
static void compute_eigs_and_cond(const FullMatrix< double > &H_orig, const unsigned int dim, 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)
static double modified_gram_schmidt(const internal::SolverGMRESImplementation::TmpVectors< VectorType > &orthogonal_vectors, const unsigned int dim, const unsigned int accumulated_iterations, VectorType &vv, Vector< double > &h, bool &re_orthogonalize, const boost::signals2::signal< void(int)> &re_orthogonalize_signal=boost::signals2::signal< void(int)>())
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::signal< void(double)> all_condition_numbers_signal
boost::signals2::connection connect_eigenvalues_slot(const std::function< void(const std::vector< std::complex< double > > &)> &slot, const bool every_iteration=false)
boost::signals2::signal< void(int)> re_orthogonalize_signal
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_NAMESPACE_CLOSE
static ::ExceptionBase & ExcTooFewTmpVectors(int arg1)
#define Assert(cond, exc)
#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)
Number add_and_dot(const Number a, const Vector< Number > &V, const Vector< Number > &W)
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
std::array< NumberType, 3 > givens_rotation(const NumberType &x, const NumberType &y)
long double gamma(const unsigned int n)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
AdditionalData(const unsigned int max_basis_size=30)
unsigned int max_basis_size
AdditionalData(const unsigned int max_n_tmp_vectors=30, const bool right_preconditioning=false, const bool use_default_residual=true, const bool force_re_orthogonalization=false)
bool right_preconditioning
bool force_re_orthogonalization
bool use_default_residual
unsigned int max_n_tmp_vectors
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)