Reference documentation for deal.II version 8.5.1
solver_gmres.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__solver_gmres_h
17 #define dealii__solver_gmres_h
18 
19 
20 
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>
30 
31 #include <vector>
32 #include <cmath>
33 #include <algorithm>
34 
35 DEAL_II_NAMESPACE_OPEN
36 
39 
40 namespace internal
41 {
45  namespace SolverGMRES
46  {
55  template <typename VectorType>
56  class TmpVectors
57  {
58  public:
63  TmpVectors(const unsigned int max_size,
65 
69  ~TmpVectors();
70 
75  VectorType &operator[] (const unsigned int i) const;
76 
83  VectorType &operator() (const unsigned int i,
84  const VectorType &temp);
85 
90  unsigned int size() const;
91 
92 
93  private:
98 
102  std::vector<VectorType *> data;
103 
108  unsigned int offset;
109  };
110  }
111 }
112 
180 template <class VectorType = Vector<double> >
181 class SolverGMRES : public Solver<VectorType>
182 {
183 public:
188  {
195  explicit
196  AdditionalData (const unsigned int max_n_tmp_vectors = 30,
197  const bool right_preconditioning = false,
198  const bool use_default_residual = true,
199  const bool force_re_orthogonalization = false);
200 
206  AdditionalData (const unsigned int max_n_tmp_vectors,
207  const bool right_preconditioning,
208  const bool use_default_residual,
209  const bool force_re_orthogonalization,
210  const bool compute_eigenvalues) DEAL_II_DEPRECATED;
211 
217  unsigned int max_n_tmp_vectors;
218 
227 
232 
240 
251  };
252 
258  const AdditionalData &data=AdditionalData());
259 
265  const AdditionalData &data=AdditionalData());
266 
270  template<typename MatrixType, typename PreconditionerType>
271  void
272  solve (const MatrixType &A,
273  VectorType &x,
274  const VectorType &b,
275  const PreconditionerType &precondition);
276 
283  boost::signals2::connection
284  connect_condition_number_slot(const std_cxx11::function<void (double)> &slot,
285  const bool every_iteration=false);
286 
293  boost::signals2::connection
295  const std_cxx11::function<void (const std::vector<std::complex<double> > &)> &slot,
296  const bool every_iteration=false);
297 
305  boost::signals2::connection
307  const std_cxx11::function<void (const FullMatrix<double> &)> &slot,
308  const bool every_iteration=true);
309 
316  boost::signals2::connection
318  const std_cxx11::function<void (const internal::SolverGMRES::TmpVectors<VectorType> &)> &slot);
319 
320 
322  int,
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.");
326 
327 protected:
328 
333 
338  boost::signals2::signal<void (double)> condition_number_signal;
339 
344  boost::signals2::signal<void (double)> all_condition_numbers_signal;
345 
350  boost::signals2::signal<void (const std::vector<std::complex<double> > &)> eigenvalues_signal;
351 
356  boost::signals2::signal<void (const std::vector<std::complex<double> > &)> all_eigenvalues_signal;
357 
362  boost::signals2::signal<void (const FullMatrix<double> &)> hessenberg_signal;
363 
368  boost::signals2::signal<void (const FullMatrix<double> &)> all_hessenberg_signal;
369 
374  boost::signals2::signal<void (const internal::SolverGMRES::TmpVectors<VectorType> &)> krylov_space_signal;
375 
379  virtual double criterion();
380 
387  int col) const;
388 
398  static double
400  const unsigned int dim,
401  const unsigned int accumulated_iterations,
402  VectorType &vv,
403  Vector<double> &h,
404  bool &re_orthogonalize);
405 
413  static void
415  const FullMatrix<double> &H_orig ,
416  const unsigned int dim,
417  const boost::signals2::signal<void (const std::vector<std::complex<double> > &)> &eigenvalues_signal,
418  const boost::signals2::signal<void (const FullMatrix<double> &)> &hessenberg_signal,
419  const boost::signals2::signal<void(double)> &cond_signal,
420  const bool log_eigenvalues);
421 
426 
431 
432 
433 private:
438 };
439 
461 template <class VectorType = Vector<double> >
462 class SolverFGMRES : public Solver<VectorType>
463 {
464 public:
469  {
473  explicit
474  AdditionalData(const unsigned int max_basis_size = 30,
475  const bool /*use_default_residual*/ = true)
476  :
478  {}
479 
483  unsigned int max_basis_size;
484  };
485 
491  const AdditionalData &data=AdditionalData());
492 
498  const AdditionalData &data=AdditionalData());
499 
503  template<typename MatrixType, typename PreconditionerType>
504  void
505  solve (const MatrixType &A,
506  VectorType &x,
507  const VectorType &b,
508  const PreconditionerType &precondition);
509 
510 private:
511 
516 
521 
526 };
527 
529 /* --------------------- Inline and template functions ------------------- */
530 
531 
532 #ifndef DOXYGEN
533 namespace internal
534 {
535  namespace SolverGMRES
536  {
537  template <class VectorType>
538  inline
540  TmpVectors (const unsigned int max_size,
542  :
543  mem(vmem),
544  data (max_size, 0),
545  offset(0)
546  {}
547 
548 
549  template <class VectorType>
550  inline
552  {
553  for (typename std::vector<VectorType *>::iterator v = data.begin();
554  v != data.end(); ++v)
555  if (*v != 0)
556  mem.free(*v);
557  }
558 
559 
560  template <class VectorType>
561  inline VectorType &
562  TmpVectors<VectorType>::operator[] (const unsigned int i) const
563  {
564  Assert (i+offset<data.size(),
565  ExcIndexRange(i, -offset, data.size()-offset));
566 
567  Assert (data[i-offset] != 0, ExcNotInitialized());
568  return *data[i-offset];
569  }
570 
571 
572  template <class VectorType>
573  inline VectorType &
574  TmpVectors<VectorType>::operator() (const unsigned int i,
575  const VectorType &temp)
576  {
577  Assert (i+offset<data.size(),
578  ExcIndexRange(i,-offset, data.size()-offset));
579  if (data[i-offset] == 0)
580  {
581  data[i-offset] = mem.alloc();
582  data[i-offset]->reinit(temp);
583  }
584  return *data[i-offset];
585  }
586 
587 
588  template <class VectorType>
589  unsigned int
591  {
592  return (data.size() > 0 ? data.size()-1 : 0);
593  }
594 
595 
596  // A comparator for better printing eigenvalues
597  inline
598  bool complex_less_pred(const std::complex<double> &x,
599  const std::complex<double> &y)
600  {
601  return x.real() < y.real() || (x.real() == y.real() && x.imag() < y.imag());
602  }
603  }
604 }
605 
606 
607 
608 template <class VectorType>
609 inline
611 AdditionalData (const unsigned int max_n_tmp_vectors,
612  const bool right_preconditioning,
613  const bool use_default_residual,
614  const bool force_re_orthogonalization)
615  :
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)
621 {}
622 
623 
624 
625 template <class VectorType>
626 inline
628 AdditionalData (const unsigned int max_n_tmp_vectors,
629  const bool right_preconditioning,
630  const bool use_default_residual,
631  const bool force_re_orthogonalization,
632  const bool compute_eigenvalues)
633  :
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)
639 {}
640 
641 
642 
643 template <class VectorType>
646  const AdditionalData &data)
647  :
648  Solver<VectorType> (cn,mem),
649  additional_data(data)
650 {}
651 
652 
653 
654 template <class VectorType>
656  const AdditionalData &data) :
657  Solver<VectorType> (cn),
658  additional_data(data)
659 {}
660 
661 
662 
663 template <class VectorType>
664 inline
665 void
667  Vector<double> &b,
668  Vector<double> &ci,
669  Vector<double> &si,
670  int col) const
671 {
672  for (int i=0 ; i<col ; i++)
673  {
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);
679  };
680 
681  const double r = 1./std::sqrt(h(col)*h(col) + h(col+1)*h(col+1));
682  si(col) = h(col+1) *r;
683  ci(col) = h(col) *r;
684  h(col) = ci(col)*h(col) + si(col)*h(col+1);
685  b(col+1)= -si(col)*b(col);
686  b(col) *= ci(col);
687 }
688 
689 
690 
691 template <class VectorType>
692 inline
693 double
695 (const internal::SolverGMRES::TmpVectors<VectorType> &orthogonal_vectors,
696  const unsigned int dim,
697  const unsigned int accumulated_iterations,
698  VectorType &vv,
699  Vector<double> &h,
700  bool &re_orthogonalize)
701 {
702  Assert(dim > 0, ExcInternalError());
703  const unsigned int inner_iteration = dim - 1;
704 
705  // need initial norm for detection of re-orthogonalization, see below
706  double norm_vv_start = 0;
707  if (re_orthogonalize == false && inner_iteration % 5 == 4)
708  norm_vv_start = vv.l2_norm();
709 
710  // Orthogonalization
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));
715 
716  // Re-orthogonalization if loss of orthogonality detected. For the test, use
717  // a strategy discussed in C. T. Kelley, Iterative Methods for Linear and
718  // Nonlinear Equations, SIAM, Philadelphia, 1995: Compare the norm of vv
719  // after orthogonalization with its norm when starting the
720  // orthogonalization. If vv became very small (here: less than the square
721  // root of the machine precision times 10), it is almost in the span of the
722  // previous vectors, which indicates loss of precision.
723  if (re_orthogonalize == false && inner_iteration % 5 == 4)
724  {
725  if (norm_vv > 10. * norm_vv_start *
726  std::sqrt(std::numeric_limits<typename VectorType::value_type>::epsilon()))
727  return norm_vv;
728 
729  else
730  {
731  re_orthogonalize = true;
732  deallog << "Re-orthogonalization enabled at step "
733  << accumulated_iterations << std::endl;
734  }
735  }
736 
737  if (re_orthogonalize == true)
738  {
739  double htmp = vv * orthogonal_vectors[0];
740  h(0) += htmp;
741  for (unsigned int i=1 ; i<dim ; ++i)
742  {
743  htmp = vv.add_and_dot(-htmp, orthogonal_vectors[i-1], orthogonal_vectors[i]);
744  h(i) += htmp;
745  }
746  norm_vv = std::sqrt(vv.add_and_dot(-htmp, orthogonal_vectors[dim-1], vv));
747  }
748 
749  return norm_vv;
750 }
751 
752 
753 
754 template<class VectorType>
755 inline void
757 (const FullMatrix<double> &H_orig,
758  const unsigned int dim,
759  const boost::signals2::signal<void (const std::vector<std::complex<double> > &)> &eigenvalues_signal,
760  const boost::signals2::signal<void (const FullMatrix<double> &)> &hessenberg_signal,
761  const boost::signals2::signal<void (double)> &cond_signal,
762  const bool log_eigenvalues)
763 {
764  //Avoid copying the Hessenberg matrix if it isn't needed.
765  if (!eigenvalues_signal.empty() || !hessenberg_signal.empty()
766  || !cond_signal.empty() || log_eigenvalues )
767  {
768  LAPACKFullMatrix<double> mat(dim,dim);
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);
773  //Avoid computing eigenvalues if they are not needed.
774  if (!eigenvalues_signal.empty() || log_eigenvalues )
775  {
776  //Copy mat so that we can compute svd below. Necessary since
777  //compute_eigenvalues will leave mat in state LAPACKSupport::unusable.
778  LAPACKFullMatrix<double> mat_eig(mat);
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);
783  //Sort eigenvalues for nicer output.
784  std::sort(eigenvalues.begin(), eigenvalues.end(),
785  internal::SolverGMRES::complex_less_pred);
786  eigenvalues_signal(eigenvalues);
787  if (log_eigenvalues)
788  {
789  deallog << "Eigenvalue estimate: ";
790  for (unsigned int i=0; i<mat_eig.n(); ++i)
791  deallog << ' ' << eigenvalues[i];
792  deallog << std::endl;
793  }
794  }
795  //Calculate condition number, avoid calculating the svd if a slot
796  //isn't connected. Need at least a 2-by-2 matrix to do the estimate.
797  if (!cond_signal.empty() && (mat.n()>1))
798  {
799  mat.compute_svd();
800  double condition_number=mat.singular_value(0)/mat.singular_value(mat.n()-1);
801  cond_signal(condition_number);
802  }
803  }
804 }
805 
806 
807 
808 template<class VectorType>
809 template<typename MatrixType, typename PreconditionerType>
810 void
811 SolverGMRES<VectorType>::solve (const MatrixType &A,
812  VectorType &x,
813  const VectorType &b,
814  const PreconditionerType &precondition)
815 {
816  // this code was written a very long time ago by people not associated with
817  // deal.II. we don't make any guarantees to its optimality or that it even
818  // works as expected...
819 
820 //TODO:[?] Check, why there are two different start residuals.
821 //TODO:[GK] Make sure the parameter in the constructor means maximum basis size
822 
823  deallog.push("GMRES");
824  const unsigned int n_tmp_vectors = additional_data.max_n_tmp_vectors;
825 
826  // Generate an object where basis vectors are stored.
827  internal::SolverGMRES::TmpVectors<VectorType> tmp_vectors (n_tmp_vectors, this->memory);
828 
829  // number of the present iteration; this
830  // number is not reset to zero upon a
831  // restart
832  unsigned int accumulated_iterations = 0;
833 
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()
841  ||additional_data.compute_eigenvalues;
842  // for eigenvalue computation, need to collect the Hessenberg matrix (before
843  // applying Givens rotations)
844  FullMatrix<double> H_orig;
845  if (do_eigenvalues)
846  H_orig.reinit(n_tmp_vectors, n_tmp_vectors-1);
847 
848  // matrix used for the orthogonalization process later
849  H.reinit(n_tmp_vectors, n_tmp_vectors-1);
850 
851  // some additional vectors, also used in the orthogonalization
853  gamma(n_tmp_vectors),
854  ci (n_tmp_vectors-1),
855  si (n_tmp_vectors-1),
856  h (n_tmp_vectors-1);
857 
858 
859  unsigned int dim = 0;
860 
862  double last_res = -std::numeric_limits<double>::max();
863 
864  // switch to determine whether we want a left or a right preconditioner. at
865  // present, left is default, but both ways are implemented
866  const bool left_precondition = !additional_data.right_preconditioning;
867 
868  // Per default the left preconditioned GMRes uses the preconditioned
869  // residual and the right preconditioned GMRes uses the unpreconditioned
870  // residual as stopping criterion.
871  const bool use_default_residual = additional_data.use_default_residual;
872 
873  // define two aliases
874  VectorType &v = tmp_vectors(0, x);
875  VectorType &p = tmp_vectors(n_tmp_vectors-1, x);
876 
877  // Following vectors are needed
878  // when not the default residuals
879  // are used as stopping criterion
880  VectorType *r=0;
881  VectorType *x_=0;
882  ::Vector<double> *gamma_=0;
883  if (!use_default_residual)
884  {
885  r=this->memory.alloc();
886  x_=this->memory.alloc();
887  r->reinit(x);
888  x_->reinit(x);
889 
890  gamma_ = new ::Vector<double> (gamma.size());
891  }
892 
893  bool re_orthogonalize = additional_data.force_re_orthogonalization;
894 
896  // outer iteration: loop until we either reach convergence or the maximum
897  // number of iterations is exceeded. each cycle of this loop amounts to one
898  // restart
899  do
900  {
901  // reset this vector to the right size
902  h.reinit (n_tmp_vectors-1);
903 
904  if (left_precondition)
905  {
906  A.vmult(p,x);
907  p.sadd(-1.,1.,b);
908  precondition.vmult(v,p);
909  }
910  else
911  {
912  A.vmult(v,x);
913  v.sadd(-1.,1.,b);
914  };
915 
916  double rho = v.l2_norm();
917 
918  // check the residual here as well since it may be that we got the exact
919  // (or an almost exact) solution vector at the outset. if we wouldn't
920  // check here, the next scaling operation would produce garbage
921  if (use_default_residual)
922  {
923  last_res = rho;
924  iteration_state = this->iteration_status (accumulated_iterations, rho, x);
925 
926  if (iteration_state != SolverControl::iterate)
927  break;
928  }
929  else
930  {
931  deallog << "default_res=" << rho << std::endl;
932 
933  if (left_precondition)
934  {
935  A.vmult(*r,x);
936  r->sadd(-1.,1.,b);
937  }
938  else
939  precondition.vmult(*r,v);
940 
941  double res = r->l2_norm();
942  last_res = res;
943  iteration_state = this->iteration_status (accumulated_iterations, res, x);
944 
945  if (iteration_state != SolverControl::iterate)
946  break;
947  }
948 
949  gamma(0) = rho;
950 
951  v *= 1./rho;
952 
953  // inner iteration doing at most as many steps as there are temporary
954  // vectors. the number of steps actually been done is propagated outside
955  // through the @p dim variable
956  for (unsigned int inner_iteration=0;
957  ((inner_iteration < n_tmp_vectors-2)
958  &&
959  (iteration_state==SolverControl::iterate));
960  ++inner_iteration)
961  {
962  ++accumulated_iterations;
963  // yet another alias
964  VectorType &vv = tmp_vectors(inner_iteration+1, x);
965 
966  if (left_precondition)
967  {
968  A.vmult(p, tmp_vectors[inner_iteration]);
969  precondition.vmult(vv,p);
970  }
971  else
972  {
973  precondition.vmult(p, tmp_vectors[inner_iteration]);
974  A.vmult(vv,p);
975  }
976 
977  dim = inner_iteration+1;
978 
979  const double s = modified_gram_schmidt(tmp_vectors, dim,
980  accumulated_iterations,
981  vv, h, re_orthogonalize);
982  h(inner_iteration+1) = s;
983 
984  //s=0 is a lucky breakdown, the solver will reach convergence,
985  //but we must not divide by zero here.
986  if (s != 0)
987  vv *= 1./s;
988 
989  // for eigenvalues, get the resulting coefficients from the
990  // orthogonalization process
991  if (do_eigenvalues)
992  for (unsigned int i=0; i<dim+1; ++i)
993  H_orig(i,inner_iteration) = h(i);
994 
995  // Transformation into tridiagonal structure
996  givens_rotation(h,gamma,ci,si,inner_iteration);
997 
998  // append vector on matrix
999  for (unsigned int i=0; i<dim; ++i)
1000  H(i,inner_iteration) = h(i);
1001 
1002  // default residual
1003  rho = std::fabs(gamma(dim));
1004 
1005  if (use_default_residual)
1006  {
1007  last_res = rho;
1008  iteration_state = this->iteration_status (accumulated_iterations, rho, x);
1009  }
1010  else
1011  {
1012  deallog << "default_res=" << rho << std::endl;
1013 
1014  ::Vector<double> h_(dim);
1015  *x_=x;
1016  *gamma_=gamma;
1017  H1.reinit(dim+1,dim);
1018 
1019  for (unsigned int i=0; i<dim+1; ++i)
1020  for (unsigned int j=0; j<dim; ++j)
1021  H1(i,j) = H(i,j);
1022 
1023  H1.backward(h_,*gamma_);
1024 
1025  if (left_precondition)
1026  for (unsigned int i=0 ; i<dim; ++i)
1027  x_->add(h_(i), tmp_vectors[i]);
1028  else
1029  {
1030  p = 0.;
1031  for (unsigned int i=0; i<dim; ++i)
1032  p.add(h_(i), tmp_vectors[i]);
1033  precondition.vmult(*r,p);
1034  x_->add(1.,*r);
1035  };
1036  A.vmult(*r,*x_);
1037  r->sadd(-1.,1.,b);
1038  // Now *r contains the unpreconditioned residual!!
1039  if (left_precondition)
1040  {
1041  const double res=r->l2_norm();
1042  last_res = res;
1043 
1044  iteration_state = this->iteration_status (accumulated_iterations, res, x);
1045  }
1046  else
1047  {
1048  precondition.vmult(*x_, *r);
1049  const double preconditioned_res=x_->l2_norm();
1050  last_res = preconditioned_res;
1051 
1052  iteration_state = this->iteration_status (accumulated_iterations,
1053  preconditioned_res, x);
1054  }
1055  }
1056  };
1057  // end of inner iteration. now calculate the solution from the temporary
1058  // vectors
1059  h.reinit(dim);
1060  H1.reinit(dim+1,dim);
1061 
1062  for (unsigned int i=0; i<dim+1; ++i)
1063  for (unsigned int j=0; j<dim; ++j)
1064  H1(i,j) = H(i,j);
1065 
1066  compute_eigs_and_cond(H_orig,dim,all_eigenvalues_signal,
1067  all_hessenberg_signal,
1068  all_condition_numbers_signal,
1069  additional_data.compute_eigenvalues);
1070 
1071  H1.backward(h,gamma);
1072 
1073  if (left_precondition)
1074  for (unsigned int i=0 ; i<dim; ++i)
1075  x.add(h(i), tmp_vectors[i]);
1076  else
1077  {
1078  p = 0.;
1079  for (unsigned int i=0; i<dim; ++i)
1080  p.add(h(i), tmp_vectors[i]);
1081  precondition.vmult(v,p);
1082  x.add(1.,v);
1083  };
1084  // end of outer iteration. restart if no convergence and the number of
1085  // iterations is not exceeded
1086  }
1087  while (iteration_state == SolverControl::iterate);
1088 
1089  compute_eigs_and_cond(H_orig,dim,eigenvalues_signal,hessenberg_signal,
1090  condition_number_signal,
1091  false);
1092 
1093  if (!krylov_space_signal.empty())
1094  krylov_space_signal(tmp_vectors);
1095 
1096  if (!use_default_residual)
1097  {
1098  this->memory.free(r);
1099  this->memory.free(x_);
1100 
1101  delete gamma_;
1102  }
1103 
1104  deallog.pop();
1105 
1106  // in case of failure: throw exception
1107  AssertThrow(iteration_state == SolverControl::success,
1108  SolverControl::NoConvergence (accumulated_iterations,
1109  last_res));
1110 }
1111 
1112 
1113 
1114 template<class VectorType>
1115 boost::signals2::connection
1117 (const std_cxx11::function<void(double)> &slot,
1118  const bool every_iteration)
1119 {
1120  if (every_iteration)
1121  {
1122  return all_condition_numbers_signal.connect(slot);
1123  }
1124  else
1125  {
1126  return condition_number_signal.connect(slot);
1127  }
1128 }
1129 
1130 
1131 
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)
1137 {
1138  if (every_iteration)
1139  {
1140  return all_eigenvalues_signal.connect(slot);
1141  }
1142  else
1143  {
1144  return eigenvalues_signal.connect(slot);
1145  }
1146 }
1147 
1148 
1149 
1150 template<class VectorType>
1151 boost::signals2::connection
1153 (const std_cxx11::function<void (const FullMatrix<double> &)> &slot,
1154  const bool every_iteration)
1155 {
1156  if (every_iteration)
1157  {
1158  return all_hessenberg_signal.connect(slot);
1159  }
1160  else
1161  {
1162  return hessenberg_signal.connect(slot);
1163  }
1164 }
1165 
1166 
1167 
1168 template<class VectorType>
1169 boost::signals2::connection
1171 (const std_cxx11::function<void (const internal::SolverGMRES::TmpVectors<VectorType> &)> &slot)
1172 {
1173  return krylov_space_signal.connect(slot);
1174 }
1175 
1176 
1177 
1178 template<class VectorType>
1179 double
1181 {
1182  // dummy implementation. this function is not needed for the present
1183  // implementation of gmres
1184  Assert (false, ExcInternalError());
1185  return 0;
1186 }
1187 
1188 
1189 //----------------------------------------------------------------------//
1190 
1191 template <class VectorType>
1194  const AdditionalData &data)
1195  :
1196  Solver<VectorType> (cn, mem),
1197  additional_data(data)
1198 {}
1199 
1200 
1201 
1202 template <class VectorType>
1204  const AdditionalData &data)
1205  :
1206  Solver<VectorType> (cn),
1207  additional_data(data)
1208 {}
1209 
1210 
1211 
1212 template<class VectorType>
1213 template<typename MatrixType, typename PreconditionerType>
1214 void
1215 SolverFGMRES<VectorType>::solve (const MatrixType &A,
1216  VectorType &x,
1217  const VectorType &b,
1218  const PreconditionerType &precondition)
1219 {
1220  deallog.push("FGMRES");
1221 
1222  SolverControl::State iteration_state = SolverControl::iterate;
1223 
1224  const unsigned int basis_size = additional_data.max_basis_size;
1225 
1226  // Generate an object where basis vectors are stored.
1227  typename internal::SolverGMRES::TmpVectors<VectorType> v (basis_size, this->memory);
1228  typename internal::SolverGMRES::TmpVectors<VectorType> z (basis_size, this->memory);
1229 
1230  // number of the present iteration; this number is not reset to zero upon a
1231  // restart
1232  unsigned int accumulated_iterations = 0;
1233 
1234  // matrix used for the orthogonalization process later
1235  H.reinit(basis_size+1, basis_size);
1236 
1237  // Vectors for projected system
1238  Vector<double> projected_rhs;
1239  Vector<double> y;
1240 
1241  // Iteration starts here
1242  double res = -std::numeric_limits<double>::max();
1243 
1244  VectorType *aux = this->memory.alloc();
1245  aux->reinit(x);
1246  do
1247  {
1248  A.vmult(*aux, x);
1249  aux->sadd(-1., 1., b);
1250 
1251  double beta = aux->l2_norm();
1252  res = beta;
1253  iteration_state = this->iteration_status(accumulated_iterations, res, x);
1254  if (iteration_state == SolverControl::success)
1255  break;
1256 
1257  H.reinit(basis_size+1, basis_size);
1258  double a = beta;
1259 
1260  for (unsigned int j=0; j<basis_size; ++j)
1261  {
1262  if (a != 0) // treat lucky breakdown
1263  v(j,x).equ(1./a, *aux);
1264  else
1265  v(j,x) = 0.;
1266 
1267 
1268  precondition.vmult(z(j,x), v[j]);
1269  A.vmult(*aux, z[j]);
1270 
1271  // Gram-Schmidt
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));
1276 
1277  // Compute projected solution
1278 
1279  if (j>0)
1280  {
1281  H1.reinit(j+1,j);
1282  projected_rhs.reinit(j+1);
1283  y.reinit(j);
1284  projected_rhs(0) = beta;
1285  H1.fill(H);
1286 
1287  // check convergence. note that the vector 'x' we pass to the
1288  // criterion is not the final solution we compute if we
1289  // decide to jump out of the iteration (we update 'x' again
1290  // right after the current loop)
1291  Householder<double> house(H1);
1292  res = house.least_squares(y, projected_rhs);
1293  iteration_state = this->iteration_status(++accumulated_iterations, res, x);
1294  if (iteration_state != SolverControl::iterate)
1295  break;
1296  }
1297  }
1298 
1299  // Update solution vector
1300  for (unsigned int j=0; j<y.size(); ++j)
1301  x.add(y(j), z[j]);
1302  }
1303  while (iteration_state == SolverControl::iterate);
1304 
1305  this->memory.free(aux);
1306 
1307  deallog.pop();
1308  // in case of failure: throw exception
1309  if (iteration_state != SolverControl::success)
1310  AssertThrow(false, SolverControl::NoConvergence (accumulated_iterations,
1311  res));
1312 }
1313 
1314 #endif // DOXYGEN
1315 
1316 DEAL_II_NAMESPACE_CLOSE
1317 
1318 #endif
VectorMemory< VectorType > & mem
Definition: solver_gmres.h:97
SolverFGMRES(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
void pop()
Definition: logstream.cc:306
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
Definition: solver_gmres.h:356
Continue iteration.
FullMatrix< double > H1
Definition: solver_gmres.h:525
AdditionalData additional_data
Definition: solver_gmres.h:332
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)
Definition: exceptions.h:369
boost::signals2::signal< void(const internal::SolverGMRES::TmpVectors< VectorType > &)> krylov_space_signal
Definition: solver_gmres.h:374
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
AdditionalData additional_data
Definition: solver_gmres.h:515
FullMatrix< double > H1
Definition: solver_gmres.h:430
virtual void free(const VectorType *const)=0
FullMatrix< double > H
Definition: solver_gmres.h:425
static ::ExceptionBase & ExcTooFewTmpVectors(int arg1)
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)
Definition: exceptions.h:564
Stop iteration, goal reached.
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
Definition: exceptions.h:313
std::size_t size() const
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
Definition: solver_gmres.h:102
boost::signals2::signal< void(const FullMatrix< double > &)> hessenberg_signal
Definition: solver_gmres.h:362
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
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)
Definition: solver.h:325
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)
Definition: logstream.cc:294
boost::signals2::signal< void(double)> condition_number_signal
Definition: solver_gmres.h:338
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
Definition: solver_gmres.h:368
FullMatrix< double > H
Definition: solver_gmres.h:520
AdditionalData(const unsigned int max_basis_size=30, const bool=true)
Definition: solver_gmres.h:474
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
Definition: solver_gmres.h:350
boost::signals2::signal< void(double)> all_condition_numbers_signal
Definition: solver_gmres.h:344
static ::ExceptionBase & ExcInternalError()