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);
95 std::vector<VectorType *>
data;
173 template <
class VectorType = Vector<
double> >
190 const bool right_preconditioning =
false,
191 const bool use_default_residual =
true,
192 const bool force_re_orthogonalization =
false);
200 const bool right_preconditioning,
201 const bool use_default_residual,
202 const bool force_re_orthogonalization,
263 template<
typename MatrixType,
typename PreconditionerType>
265 solve (
const MatrixType &A,
268 const PreconditionerType &precondition);
276 boost::signals2::connection
277 connect_condition_number_slot(
const std_cxx11::function<
void (
double)> &slot,
278 const bool every_iteration=
false);
286 boost::signals2::connection
287 connect_eigenvalues_slot(
288 const std_cxx11::function<
void (
const std::vector<std::complex<double> > &)> &slot,
289 const bool every_iteration=
false);
294 <<
"The number of temporary vectors you gave (" 295 << arg1 <<
") is too small. It should be at least 10 for " 296 <<
"any results, and much more for reasonable ones.");
321 boost::signals2::signal<void (const std::vector<std::complex<double> > &)> eigenvalues_signal;
327 boost::signals2::signal<void (const std::vector<std::complex<double> > &)> all_eigenvalues_signal;
333 virtual double criterion();
354 const unsigned int dim,
355 const unsigned int accumulated_iterations,
358 bool &re_orthogonalize);
368 compute_eigs_and_cond(
370 const unsigned int dim,
371 const boost::signals2::signal<
void (
const std::vector<std::complex<double> > &)> &eigenvalues_signal,
372 const boost::signals2::signal<
void(
double)> &cond_signal,
373 const bool log_eigenvalues);
414 template <
class VectorType = Vector<
double> >
430 max_basis_size(max_basis_size)
456 template<
typename MatrixType,
typename PreconditionerType>
458 solve (
const MatrixType &A,
461 const PreconditionerType &precondition);
490 template <
class VectorType>
502 template <
class VectorType>
506 for (
typename std::vector<VectorType *>::iterator v =
data.begin();
507 v !=
data.end(); ++v)
513 template <
class VectorType>
525 template <
class VectorType>
528 const VectorType &temp)
542 bool complex_less_pred(
const std::complex<double> &x,
543 const std::complex<double> &y)
545 return x.real() < y.real() || (x.real() == y.real() && x.imag() < y.imag());
552 template <
class VectorType>
556 const bool right_preconditioning,
557 const bool use_default_residual,
558 const bool force_re_orthogonalization)
560 max_n_tmp_vectors(max_n_tmp_vectors),
561 right_preconditioning(right_preconditioning),
562 use_default_residual(use_default_residual),
563 force_re_orthogonalization(force_re_orthogonalization),
564 compute_eigenvalues(
false)
569 template <
class VectorType>
573 const bool right_preconditioning,
574 const bool use_default_residual,
575 const bool force_re_orthogonalization,
576 const bool compute_eigenvalues)
578 max_n_tmp_vectors(max_n_tmp_vectors),
579 right_preconditioning(right_preconditioning),
580 use_default_residual(use_default_residual),
581 force_re_orthogonalization(force_re_orthogonalization),
582 compute_eigenvalues(compute_eigenvalues)
587 template <
class VectorType>
593 additional_data(data)
598 template <
class VectorType>
602 additional_data(data)
607 template <
class VectorType>
616 for (
int i=0 ; i<col ; i++)
618 const double s = si(i);
619 const double c = ci(i);
620 const double dummy = h(i);
621 h(i) = c*dummy + s*h(i+1);
622 h(i+1) = -s*dummy + c*h(i+1);
625 const double r = 1./std::sqrt(h(col)*h(col) + h(col+1)*h(col+1));
626 si(col) = h(col+1) *r;
628 h(col) = ci(col)*h(col) + si(col)*h(col+1);
629 b(col+1)= -si(col)*b(col);
635 template <
class VectorType>
640 const unsigned int dim,
641 const unsigned int accumulated_iterations,
644 bool &re_orthogonalize)
646 Assert(dim > 0, ExcInternalError());
647 const unsigned int inner_iteration = dim - 1;
650 double norm_vv_start = 0;
651 if (re_orthogonalize ==
false && inner_iteration % 5 == 4)
652 norm_vv_start = vv.l2_norm();
655 h(0) = vv * orthogonal_vectors[0];
656 for (
unsigned int i=1 ; i<dim ; ++i)
657 h(i) = vv.
add_and_dot(-h(i-1), orthogonal_vectors[i-1], orthogonal_vectors[i]);
658 double norm_vv = std::sqrt(vv.add_and_dot(-h(dim-1), orthogonal_vectors[dim-1], vv));
667 if (re_orthogonalize ==
false && inner_iteration % 5 == 4)
669 if (norm_vv > 10. * norm_vv_start *
670 std::sqrt(std::numeric_limits<typename VectorType::value_type>::epsilon()))
675 re_orthogonalize =
true;
676 deallog <<
"Re-orthogonalization enabled at step " 677 << accumulated_iterations << std::endl;
681 if (re_orthogonalize ==
true)
683 double htmp = vv * orthogonal_vectors[0];
685 for (
unsigned int i=1 ; i<dim ; ++i)
687 htmp = vv.
add_and_dot(-htmp, orthogonal_vectors[i-1], orthogonal_vectors[i]);
690 norm_vv = std::sqrt(vv.add_and_dot(-htmp, orthogonal_vectors[dim-1], vv));
698 template<
class VectorType>
702 const unsigned int dim,
703 const boost::signals2::signal<
void (
const std::vector<std::complex<double> > &)> &eigenvalues_signal,
704 const boost::signals2::signal<
void (
double)> &cond_signal,
705 const bool log_eigenvalues)
708 if (!eigenvalues_signal.empty() || !cond_signal.empty() || log_eigenvalues )
711 for (
unsigned int i=0; i<dim; ++i)
712 for (
unsigned int j=0; j<dim; ++j)
713 mat(i,j) = H_orig(i,j);
715 if (!eigenvalues_signal.empty() || log_eigenvalues )
721 std::vector<std::complex<double> > eigenvalues(dim);
722 for (
unsigned int i=0; i<mat_eig.
n(); ++i)
725 std::sort(eigenvalues.begin(), eigenvalues.end(),
726 internal::SolverGMRES::complex_less_pred);
727 eigenvalues_signal(eigenvalues);
730 deallog <<
"Eigenvalue estimate: ";
731 for (
unsigned int i=0; i<mat_eig.
n(); ++i)
732 deallog <<
' ' << eigenvalues[i];
733 deallog << std::endl;
738 if (!cond_signal.empty() && (mat.
n()>1))
742 cond_signal(condition_number);
749 template<
class VectorType>
750 template<
typename MatrixType,
typename PreconditionerType>
755 const PreconditionerType &precondition)
764 deallog.
push(
"GMRES");
773 unsigned int accumulated_iterations = 0;
775 const bool do_eigenvalues=
776 !condition_number_signal.empty()
777 |!all_condition_numbers_signal.empty()
778 |!eigenvalues_signal.empty()
779 |!all_eigenvalues_signal.empty()
785 H_orig.
reinit(n_tmp_vectors, n_tmp_vectors-1);
788 H.reinit(n_tmp_vectors, n_tmp_vectors-1);
792 gamma(n_tmp_vectors),
793 ci (n_tmp_vectors-1),
794 si (n_tmp_vectors-1),
798 unsigned int dim = 0;
801 double last_res = -std::numeric_limits<double>::max();
813 VectorType &v = tmp_vectors(0, x);
814 VectorType &p = tmp_vectors(n_tmp_vectors-1, x);
822 if (!use_default_residual)
824 r=this->memory.alloc();
825 x_=this->memory.alloc();
829 gamma_ = new ::Vector<double> (gamma.size());
841 h.
reinit (n_tmp_vectors-1);
843 if (left_precondition)
847 precondition.vmult(v,p);
855 double rho = v.l2_norm();
860 if (use_default_residual)
863 iteration_state = this->iteration_status (accumulated_iterations, rho, x);
870 deallog <<
"default_res=" << rho << std::endl;
872 if (left_precondition)
878 precondition.vmult(*r,v);
880 double res = r->l2_norm();
882 iteration_state = this->iteration_status (accumulated_iterations, res, x);
895 for (
unsigned int inner_iteration=0;
896 ((inner_iteration < n_tmp_vectors-2)
901 ++accumulated_iterations;
903 VectorType &vv = tmp_vectors(inner_iteration+1, x);
905 if (left_precondition)
907 A.vmult(p, tmp_vectors[inner_iteration]);
908 precondition.vmult(vv,p);
912 precondition.vmult(p, tmp_vectors[inner_iteration]);
916 dim = inner_iteration+1;
918 const double s = modified_gram_schmidt(tmp_vectors, dim,
919 accumulated_iterations,
920 vv, h, re_orthogonalize);
921 h(inner_iteration+1) = s;
931 for (
unsigned int i=0; i<dim+1; ++i)
932 H_orig(i,inner_iteration) = h(i);
935 givens_rotation(h,gamma,ci,si,inner_iteration);
938 for (
unsigned int i=0; i<dim; ++i)
939 H(i,inner_iteration) = h(i);
942 rho = std::fabs(gamma(dim));
944 if (use_default_residual)
947 iteration_state = this->iteration_status (accumulated_iterations, rho, x);
951 deallog <<
"default_res=" << rho << std::endl;
958 for (
unsigned int i=0; i<dim+1; ++i)
959 for (
unsigned int j=0; j<dim; ++j)
962 H1.backward(h_,*gamma_);
964 if (left_precondition)
965 for (
unsigned int i=0 ; i<dim; ++i)
966 x_->add(h_(i), tmp_vectors[i]);
970 for (
unsigned int i=0; i<dim; ++i)
971 p.add(h_(i), tmp_vectors[i]);
972 precondition.vmult(*r,p);
978 if (left_precondition)
980 const double res=r->l2_norm();
983 iteration_state = this->iteration_status (accumulated_iterations, res, x);
987 precondition.vmult(*x_, *r);
988 const double preconditioned_res=x_->l2_norm();
989 last_res = preconditioned_res;
991 iteration_state = this->iteration_status (accumulated_iterations,
992 preconditioned_res, x);
999 H1.reinit(dim+1,dim);
1001 for (
unsigned int i=0; i<dim+1; ++i)
1002 for (
unsigned int j=0; j<dim; ++j)
1005 compute_eigs_and_cond(H_orig,dim,all_eigenvalues_signal,
1006 all_condition_numbers_signal,
1009 H1.backward(h,gamma);
1011 if (left_precondition)
1012 for (
unsigned int i=0 ; i<dim; ++i)
1013 x.add(h(i), tmp_vectors[i]);
1017 for (
unsigned int i=0; i<dim; ++i)
1018 p.add(h(i), tmp_vectors[i]);
1019 precondition.vmult(v,p);
1027 compute_eigs_and_cond(H_orig,dim,eigenvalues_signal,condition_number_signal,
1029 if (!use_default_residual)
1031 this->memory.free(r);
1032 this->memory.free(x_);
1047 template<
class VectorType>
1048 boost::signals2::connection
1050 (
const std_cxx11::function<
void(
double)> &slot,
1051 const bool every_iteration)
1053 if (every_iteration)
1055 return all_condition_numbers_signal.connect(slot);
1059 return condition_number_signal.connect(slot);
1065 template<
class VectorType>
1066 boost::signals2::connection
1068 (
const std_cxx11::function<
void (
const std::vector<std::complex<double> > &)> &slot,
1069 const bool every_iteration)
1071 if (every_iteration)
1073 return all_eigenvalues_signal.connect(slot);
1077 return eigenvalues_signal.connect(slot);
1083 template<
class VectorType>
1089 Assert (
false, ExcInternalError());
1096 template <
class VectorType>
1102 additional_data(data)
1107 template <
class VectorType>
1112 additional_data(data)
1117 template<
class VectorType>
1118 template<
typename MatrixType,
typename PreconditionerType>
1122 const VectorType &b,
1123 const PreconditionerType &precondition)
1125 deallog.
push(
"FGMRES");
1129 const unsigned int basis_size = additional_data.max_basis_size;
1137 unsigned int accumulated_iterations = 0;
1140 H.reinit(basis_size+1, basis_size);
1147 double res = -std::numeric_limits<double>::max();
1149 VectorType *aux = this->memory.alloc();
1154 aux->sadd(-1., 1., b);
1156 double beta = aux->l2_norm();
1158 iteration_state = this->iteration_status(accumulated_iterations, res, x);
1162 H.reinit(basis_size+1, basis_size);
1165 for (
unsigned int j=0; j<basis_size; ++j)
1168 v(j,x).equ(1./a, *aux);
1173 precondition.vmult(z(j,x), v[j]);
1174 A.vmult(*aux, z[j]);
1177 H(0,j) = *aux * v[0];
1178 for (
unsigned int i=1; i<=j; ++i)
1179 H(i,j) = aux->add_and_dot(-H(i-1,j), v[i-1], v[i]);
1180 H(j+1,j) = a = std::sqrt(aux->add_and_dot(-H(j,j), v[j], *aux));
1187 projected_rhs.
reinit(j+1);
1189 projected_rhs(0) = beta;
1198 iteration_state = this->iteration_status(++accumulated_iterations, res, x);
1205 for (
unsigned int j=0; j<y.
size(); ++j)
1210 this->memory.free(aux);
1221 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
std::complex< number > eigenvalue(const size_type i) const
AdditionalData additional_data
#define AssertThrow(cond, exc)
AdditionalData additional_data
virtual void free(const VectorType *const)=0
bool right_preconditioning
#define DEAL_II_DEPRECATED
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
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(double)> &cond_signal, const bool log_eigenvalues)
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
std::vector< VectorType * > data
bool force_re_orthogonalization
number singular_value(const size_type i) 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)
double least_squares(Vector< number2 > &dst, const Vector< number2 > &src) const
void compute_eigenvalues(const bool right_eigenvectors=false, const bool left_eigenvectors=false)
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
unsigned int max_n_tmp_vectors
boost::signals2::signal< void(double)> all_condition_numbers_signal