16 #ifndef dealii__solver_gmres_h 17 #define dealii__solver_gmres_h 21 #include <deal.II/base/config.h> 22 #include <deal.II/base/subscriptor.h> 23 #include <deal.II/base/logstream.h> 24 #include <deal.II/lac/householder.h> 25 #include <deal.II/lac/solver.h> 26 #include <deal.II/lac/solver_control.h> 27 #include <deal.II/lac/full_matrix.h> 28 #include <deal.II/lac/lapack_full_matrix.h> 29 #include <deal.II/lac/vector.h> 35 DEAL_II_NAMESPACE_OPEN
55 template <
typename VectorType>
75 VectorType &
operator[] (
const unsigned int i)
const;
84 const VectorType &temp);
90 unsigned int size()
const;
102 std::vector<VectorType *>
data;
180 template <
class VectorType = Vector<
double> >
270 template<
typename MatrixType,
typename PreconditionerType>
272 solve (
const MatrixType &A,
275 const PreconditionerType &precondition);
283 boost::signals2::connection
285 const bool every_iteration=
false);
293 boost::signals2::connection
295 const std_cxx11::function<
void (
const std::vector<std::complex<double> > &)> &slot,
296 const bool every_iteration=
false);
305 boost::signals2::connection
308 const bool every_iteration=
true);
316 boost::signals2::connection
323 <<
"The number of temporary vectors you gave (" 324 << arg1 <<
") is too small. It should be at least 10 for " 325 <<
"any results, and much more for reasonable ones.");
374 boost::signals2::signal<void (const internal::SolverGMRES::TmpVectors<VectorType> &)>
krylov_space_signal;
400 const unsigned int dim,
401 const unsigned int accumulated_iterations,
404 bool &re_orthogonalize);
416 const unsigned int dim,
417 const boost::signals2::signal<
void (
const std::vector<std::complex<double> > &)> &
eigenvalues_signal,
419 const boost::signals2::signal<
void(
double)> &cond_signal,
420 const bool log_eigenvalues);
461 template <
class VectorType = Vector<
double> >
503 template<
typename MatrixType,
typename PreconditionerType>
505 solve (
const MatrixType &A,
508 const PreconditionerType &precondition);
537 template <
class VectorType>
549 template <
class VectorType>
553 for (
typename std::vector<VectorType *>::iterator v = data.begin();
554 v != data.end(); ++v)
560 template <
class VectorType>
564 Assert (i+offset<data.size(),
568 return *data[i-offset];
572 template <
class VectorType>
575 const VectorType &temp)
577 Assert (i+offset<data.size(),
579 if (data[i-offset] == 0)
581 data[i-offset] = mem.
alloc();
582 data[i-offset]->reinit(temp);
584 return *data[i-offset];
588 template <
class VectorType>
592 return (data.size() > 0 ? data.size()-1 : 0);
598 bool complex_less_pred(
const std::complex<double> &x,
599 const std::complex<double> &y)
601 return x.real() < y.real() || (x.real() == y.real() && x.imag() < y.imag());
608 template <
class VectorType>
612 const bool right_preconditioning,
613 const bool use_default_residual,
614 const bool force_re_orthogonalization)
616 max_n_tmp_vectors(max_n_tmp_vectors),
617 right_preconditioning(right_preconditioning),
618 use_default_residual(use_default_residual),
619 force_re_orthogonalization(force_re_orthogonalization),
620 compute_eigenvalues(false)
625 template <
class VectorType>
629 const bool right_preconditioning,
630 const bool use_default_residual,
631 const bool force_re_orthogonalization,
632 const bool compute_eigenvalues)
634 max_n_tmp_vectors(max_n_tmp_vectors),
635 right_preconditioning(right_preconditioning),
636 use_default_residual(use_default_residual),
637 force_re_orthogonalization(force_re_orthogonalization),
638 compute_eigenvalues(compute_eigenvalues)
643 template <
class VectorType>
646 const AdditionalData &data)
648 Solver<VectorType> (cn,mem),
649 additional_data(data)
654 template <
class VectorType>
656 const AdditionalData &data) :
658 additional_data(data)
663 template <
class VectorType>
672 for (
int i=0 ; i<col ; i++)
674 const double s = si(i);
675 const double c = ci(i);
676 const double dummy = h(i);
677 h(i) = c*dummy + s*h(i+1);
678 h(i+1) = -s*dummy + c*h(i+1);
681 const double r = 1./std::sqrt(h(col)*h(col) + h(col+1)*h(col+1));
682 si(col) = h(col+1) *r;
684 h(col) = ci(col)*h(col) + si(col)*h(col+1);
685 b(col+1)= -si(col)*
b(col);
691 template <
class VectorType>
696 const unsigned int dim,
697 const unsigned int accumulated_iterations,
700 bool &re_orthogonalize)
703 const unsigned int inner_iteration = dim - 1;
706 double norm_vv_start = 0;
707 if (re_orthogonalize ==
false && inner_iteration % 5 == 4)
708 norm_vv_start = vv.l2_norm();
711 h(0) = vv * orthogonal_vectors[0];
712 for (
unsigned int i=1 ; i<dim ; ++i)
713 h(i) = vv.
add_and_dot(-h(i-1), orthogonal_vectors[i-1], orthogonal_vectors[i]);
714 double norm_vv = std::sqrt(vv.add_and_dot(-h(dim-1), orthogonal_vectors[dim-1], vv));
723 if (re_orthogonalize ==
false && inner_iteration % 5 == 4)
725 if (norm_vv > 10. * norm_vv_start *
726 std::sqrt(std::numeric_limits<typename VectorType::value_type>::epsilon()))
731 re_orthogonalize =
true;
732 deallog <<
"Re-orthogonalization enabled at step " 733 << accumulated_iterations << std::endl;
737 if (re_orthogonalize ==
true)
739 double htmp = vv * orthogonal_vectors[0];
741 for (
unsigned int i=1 ; i<dim ; ++i)
743 htmp = vv.
add_and_dot(-htmp, orthogonal_vectors[i-1], orthogonal_vectors[i]);
746 norm_vv = std::sqrt(vv.add_and_dot(-htmp, orthogonal_vectors[dim-1], vv));
754 template<
class VectorType>
758 const unsigned int dim,
759 const boost::signals2::signal<
void (
const std::vector<std::complex<double> > &)> &eigenvalues_signal,
761 const boost::signals2::signal<
void (
double)> &cond_signal,
762 const bool log_eigenvalues)
765 if (!eigenvalues_signal.empty() || !hessenberg_signal.empty()
766 || !cond_signal.empty() || log_eigenvalues )
769 for (
unsigned int i=0; i<dim; ++i)
770 for (
unsigned int j=0; j<dim; ++j)
771 mat(i,j) = H_orig(i,j);
772 hessenberg_signal(H_orig);
774 if (!eigenvalues_signal.empty() || log_eigenvalues )
779 mat_eig.compute_eigenvalues();
780 std::vector<std::complex<double> > eigenvalues(dim);
781 for (
unsigned int i=0; i<mat_eig.n(); ++i)
782 eigenvalues[i] = mat_eig.eigenvalue(i);
784 std::sort(eigenvalues.begin(), eigenvalues.end(),
785 internal::SolverGMRES::complex_less_pred);
786 eigenvalues_signal(eigenvalues);
789 deallog <<
"Eigenvalue estimate: ";
790 for (
unsigned int i=0; i<mat_eig.n(); ++i)
791 deallog <<
' ' << eigenvalues[i];
792 deallog << std::endl;
797 if (!cond_signal.empty() && (mat.n()>1))
800 double condition_number=mat.singular_value(0)/mat.singular_value(mat.n()-1);
801 cond_signal(condition_number);
808 template<
class VectorType>
809 template<
typename MatrixType,
typename PreconditionerType>
814 const PreconditionerType &precondition)
823 deallog.
push(
"GMRES");
832 unsigned int accumulated_iterations = 0;
834 const bool do_eigenvalues=
835 !condition_number_signal.empty()
836 ||!all_condition_numbers_signal.empty()
837 ||!eigenvalues_signal.empty()
838 ||!all_eigenvalues_signal.empty()
839 ||!hessenberg_signal.empty()
840 ||!all_hessenberg_signal.empty()
846 H_orig.
reinit(n_tmp_vectors, n_tmp_vectors-1);
849 H.reinit(n_tmp_vectors, n_tmp_vectors-1);
853 gamma(n_tmp_vectors),
854 ci (n_tmp_vectors-1),
855 si (n_tmp_vectors-1),
859 unsigned int dim = 0;
862 double last_res = -std::numeric_limits<double>::max();
874 VectorType &v = tmp_vectors(0, x);
875 VectorType &p = tmp_vectors(n_tmp_vectors-1, x);
883 if (!use_default_residual)
885 r=this->memory.alloc();
886 x_=this->memory.alloc();
890 gamma_ = new ::Vector<double> (gamma.size());
902 h.
reinit (n_tmp_vectors-1);
904 if (left_precondition)
908 precondition.vmult(v,p);
916 double rho = v.l2_norm();
921 if (use_default_residual)
924 iteration_state = this->iteration_status (accumulated_iterations, rho, x);
931 deallog <<
"default_res=" << rho << std::endl;
933 if (left_precondition)
939 precondition.vmult(*r,v);
941 double res = r->l2_norm();
943 iteration_state = this->iteration_status (accumulated_iterations, res, x);
956 for (
unsigned int inner_iteration=0;
957 ((inner_iteration < n_tmp_vectors-2)
962 ++accumulated_iterations;
964 VectorType &vv = tmp_vectors(inner_iteration+1, x);
966 if (left_precondition)
968 A.vmult(p, tmp_vectors[inner_iteration]);
969 precondition.vmult(vv,p);
973 precondition.vmult(p, tmp_vectors[inner_iteration]);
977 dim = inner_iteration+1;
979 const double s = modified_gram_schmidt(tmp_vectors, dim,
980 accumulated_iterations,
981 vv, h, re_orthogonalize);
982 h(inner_iteration+1) = s;
992 for (
unsigned int i=0; i<dim+1; ++i)
993 H_orig(i,inner_iteration) = h(i);
996 givens_rotation(h,gamma,ci,si,inner_iteration);
999 for (
unsigned int i=0; i<dim; ++i)
1000 H(i,inner_iteration) = h(i);
1003 rho = std::fabs(gamma(dim));
1005 if (use_default_residual)
1008 iteration_state = this->iteration_status (accumulated_iterations, rho, x);
1012 deallog <<
"default_res=" << rho << std::endl;
1019 for (
unsigned int i=0; i<dim+1; ++i)
1020 for (
unsigned int j=0; j<dim; ++j)
1023 H1.backward(h_,*gamma_);
1025 if (left_precondition)
1026 for (
unsigned int i=0 ; i<dim; ++i)
1027 x_->add(h_(i), tmp_vectors[i]);
1031 for (
unsigned int i=0; i<dim; ++i)
1032 p.add(h_(i), tmp_vectors[i]);
1033 precondition.vmult(*r,p);
1039 if (left_precondition)
1041 const double res=r->l2_norm();
1044 iteration_state = this->iteration_status (accumulated_iterations, res, x);
1048 precondition.vmult(*x_, *r);
1049 const double preconditioned_res=x_->l2_norm();
1050 last_res = preconditioned_res;
1052 iteration_state = this->iteration_status (accumulated_iterations,
1053 preconditioned_res, x);
1060 H1.reinit(dim+1,dim);
1062 for (
unsigned int i=0; i<dim+1; ++i)
1063 for (
unsigned int j=0; j<dim; ++j)
1066 compute_eigs_and_cond(H_orig,dim,all_eigenvalues_signal,
1067 all_hessenberg_signal,
1068 all_condition_numbers_signal,
1071 H1.backward(h,gamma);
1073 if (left_precondition)
1074 for (
unsigned int i=0 ; i<dim; ++i)
1075 x.add(h(i), tmp_vectors[i]);
1079 for (
unsigned int i=0; i<dim; ++i)
1080 p.add(h(i), tmp_vectors[i]);
1081 precondition.vmult(v,p);
1089 compute_eigs_and_cond(H_orig,dim,eigenvalues_signal,hessenberg_signal,
1090 condition_number_signal,
1093 if (!krylov_space_signal.empty())
1094 krylov_space_signal(tmp_vectors);
1096 if (!use_default_residual)
1098 this->memory.free(r);
1099 this->memory.free(x_);
1114 template<
class VectorType>
1115 boost::signals2::connection
1117 (
const std_cxx11::function<
void(
double)> &slot,
1118 const bool every_iteration)
1120 if (every_iteration)
1122 return all_condition_numbers_signal.connect(slot);
1126 return condition_number_signal.connect(slot);
1132 template<
class VectorType>
1133 boost::signals2::connection
1135 (
const std_cxx11::function<
void (
const std::vector<std::complex<double> > &)> &slot,
1136 const bool every_iteration)
1138 if (every_iteration)
1140 return all_eigenvalues_signal.connect(slot);
1144 return eigenvalues_signal.connect(slot);
1150 template<
class VectorType>
1151 boost::signals2::connection
1154 const bool every_iteration)
1156 if (every_iteration)
1158 return all_hessenberg_signal.connect(slot);
1162 return hessenberg_signal.connect(slot);
1168 template<
class VectorType>
1169 boost::signals2::connection
1173 return krylov_space_signal.connect(slot);
1178 template<
class VectorType>
1191 template <
class VectorType>
1194 const AdditionalData &data)
1196 Solver<VectorType> (cn, mem),
1197 additional_data(data)
1202 template <
class VectorType>
1204 const AdditionalData &data)
1207 additional_data(data)
1212 template<
class VectorType>
1213 template<
typename MatrixType,
typename PreconditionerType>
1217 const VectorType &b,
1218 const PreconditionerType &precondition)
1220 deallog.
push(
"FGMRES");
1224 const unsigned int basis_size = additional_data.max_basis_size;
1232 unsigned int accumulated_iterations = 0;
1235 H.reinit(basis_size+1, basis_size);
1242 double res = -std::numeric_limits<double>::max();
1244 VectorType *aux = this->memory.alloc();
1249 aux->sadd(-1., 1., b);
1251 double beta = aux->l2_norm();
1253 iteration_state = this->iteration_status(accumulated_iterations, res, x);
1257 H.reinit(basis_size+1, basis_size);
1260 for (
unsigned int j=0; j<basis_size; ++j)
1263 v(j,x).equ(1./a, *aux);
1268 precondition.vmult(z(j,x), v[j]);
1269 A.vmult(*aux, z[j]);
1272 H(0,j) = *aux * v[0];
1273 for (
unsigned int i=1; i<=j; ++i)
1274 H(i,j) = aux->add_and_dot(-H(i-1,j), v[i-1], v[i]);
1275 H(j+1,j) = a = std::sqrt(aux->add_and_dot(-H(j,j), v[j], *aux));
1282 projected_rhs.
reinit(j+1);
1284 projected_rhs(0) = beta;
1292 res = house.least_squares(y, projected_rhs);
1293 iteration_state = this->iteration_status(++accumulated_iterations, res, x);
1300 for (
unsigned int j=0; j<y.
size(); ++j)
1305 this->memory.free(aux);
1316 DEAL_II_NAMESPACE_CLOSE
VectorMemory< VectorType > & mem
SolverFGMRES(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
Number add_and_dot(const Number a, const Vector< Number > &V, const Vector< Number > &W)
virtual VectorType * alloc()=0
boost::signals2::signal< void(const std::vector< std::complex< double > > &)> all_eigenvalues_signal
AdditionalData additional_data
boost::signals2::connection connect_hessenberg_slot(const std_cxx11::function< void(const FullMatrix< double > &)> &slot, const bool every_iteration=true)
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, const bool log_eigenvalues)
static ::ExceptionBase & ExcNotInitialized()
#define AssertThrow(cond, exc)
boost::signals2::signal< void(const internal::SolverGMRES::TmpVectors< VectorType > &)> krylov_space_signal
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
AdditionalData additional_data
virtual void free(const VectorType *const)=0
static ::ExceptionBase & ExcTooFewTmpVectors(int arg1)
bool right_preconditioning
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
VectorType & operator()(const unsigned int i, const VectorType &temp)
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)
#define DeclException1(Exception1, type1, outsequence)
Stop iteration, goal reached.
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
unsigned int max_basis_size
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
boost::signals2::connection connect_krylov_space_slot(const std_cxx11::function< void(const internal::SolverGMRES::TmpVectors< VectorType > &)> &slot)
std::vector< VectorType * > data
boost::signals2::signal< void(const FullMatrix< double > &)> hessenberg_signal
bool force_re_orthogonalization
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
unsigned int size() const
static double modified_gram_schmidt(const internal::SolverGMRES::TmpVectors< VectorType > &orthogonal_vectors, const unsigned int dim, const unsigned int accumulated_iterations, VectorType &vv, Vector< double > &h, bool &re_orthogonalize)
boost::signals2::connection connect_condition_number_slot(const std_cxx11::function< void(double)> &slot, const bool every_iteration=false)
virtual double criterion()
void push(const std::string &text)
boost::signals2::signal< void(double)> condition_number_signal
TmpVectors(const unsigned int max_size, VectorMemory< VectorType > &vmem)
boost::signals2::connection connect_eigenvalues_slot(const std_cxx11::function< void(const std::vector< std::complex< double > > &)> &slot, const bool every_iteration=false)
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
boost::signals2::signal< void(const FullMatrix< double > &)> all_hessenberg_signal
bool use_default_residual
AdditionalData(const unsigned int max_basis_size=30, const bool=true)
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
VectorType & operator[](const unsigned int i) const
boost::signals2::signal< void(const std::vector< std::complex< double > > &)> eigenvalues_signal
unsigned int max_n_tmp_vectors
boost::signals2::signal< void(double)> all_condition_numbers_signal
static ::ExceptionBase & ExcInternalError()