Dune Core Modules (unstable)

poweriteration.hh
1 // SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2 // SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4 // vi: set et ts=4 sw=2 sts=2:
5 #ifndef DUNE_ISTL_EIGENVALUE_POWERITERATION_HH
6 #define DUNE_ISTL_EIGENVALUE_POWERITERATION_HH
7 
8 #include <cstddef> // provides std::size_t
9 #include <cmath> // provides std::sqrt, std::abs
10 
11 #include <type_traits> // provides std::is_same
12 #include <iostream> // provides std::cout, std::endl
13 #include <limits> // provides std::numeric_limits
14 #include <ios> // provides std::left, std::ios::left
15 #include <iomanip> // provides std::setw, std::resetiosflags
16 #include <memory> // provides std::unique_ptr
17 #include <string> // provides std::string
18 
19 #include <dune/common/exceptions.hh> // provides DUNE_THROW(...)
20 
21 #include <dune/istl/blocklevel.hh> // provides Dune::blockLevel
22 #include <dune/istl/operators.hh> // provides Dune::LinearOperator
23 #include <dune/istl/solvercategory.hh> // provides Dune::SolverCategory::sequential
24 #include <dune/istl/solvertype.hh> // provides Dune::IsDirectSolver
25 #include <dune/istl/operators.hh> // provides Dune::MatrixAdapter
26 #include <dune/istl/istlexception.hh> // provides Dune::ISTLError
27 #include <dune/istl/io.hh> // provides Dune::printvector(...)
28 #include <dune/istl/solvers.hh> // provides Dune::InverseOperatorResult
29 
30 namespace Dune
31 {
32 
37  namespace Impl {
45  template <class X, class Y = X>
46  class ScalingLinearOperator : public Dune::LinearOperator<X,Y>
47  {
48  public:
49  typedef X domain_type;
50  typedef Y range_type;
51  typedef typename X::field_type field_type;
52 
53  ScalingLinearOperator (field_type immutable_scaling,
54  const field_type& mutable_scaling)
55  : immutable_scaling_(immutable_scaling),
56  mutable_scaling_(mutable_scaling)
57  {}
58 
59  virtual void apply (const X& x, Y& y) const
60  {
61  y = x;
62  y *= immutable_scaling_*mutable_scaling_;
63  }
64 
65  virtual void applyscaleadd (field_type alpha, const X& x, Y& y) const
66  {
67  X temp(x);
68  temp *= immutable_scaling_*mutable_scaling_;
69  y.axpy(alpha,temp);
70  }
71 
73  virtual SolverCategory::Category category() const
74  {
76  }
77 
78  protected:
79  const field_type immutable_scaling_;
80  const field_type& mutable_scaling_;
81  };
82 
83 
92  template <class OP1, class OP2>
93  class LinearOperatorSum
94  : public Dune::LinearOperator<typename OP1::domain_type,
95  typename OP1::range_type>
96  {
97  public:
98  typedef typename OP1::domain_type domain_type;
99  typedef typename OP1::range_type range_type;
100  typedef typename domain_type::field_type field_type;
101 
102  LinearOperatorSum (const OP1& op1, const OP2& op2)
103  : op1_(op1), op2_(op2)
104  {
105  static_assert(std::is_same<typename OP2::domain_type,domain_type>::value,
106  "Domain type of both operators doesn't match!");
107  static_assert(std::is_same<typename OP2::range_type,range_type>::value,
108  "Range type of both operators doesn't match!");
109  }
110 
111  virtual void apply (const domain_type& x, range_type& y) const
112  {
113  op1_.apply(x,y);
114  op2_.applyscaleadd(1.0,x,y);
115  }
116 
117  virtual void applyscaleadd (field_type alpha,
118  const domain_type& x, range_type& y) const
119  {
120  range_type temp(y);
121  op1_.apply(x,temp);
122  op2_.applyscaleadd(1.0,x,temp);
123  y.axpy(alpha,temp);
124  }
125 
127  virtual SolverCategory::Category category() const
128  {
130  }
131 
132  protected:
133  const OP1& op1_;
134  const OP2& op2_;
135  };
136  } // end namespace Impl
137 
174  template <typename BCRSMatrix, typename BlockVector>
176  {
177  protected:
178  // Type definitions for type of iteration operator (m_ - mu_*I)
181  typedef Impl::ScalingLinearOperator<BlockVector> ScalingOperator;
182  typedef Impl::LinearOperatorSum<MatrixOperator,ScalingOperator> OperatorSum;
183 
184  public:
186  typedef typename BlockVector::field_type Real;
187 
189  typedef OperatorSum IterationOperator;
190 
191  public:
207  const unsigned int nIterationsMax = 1000,
208  const unsigned int verbosity_level = 0)
209  : m_(m), nIterationsMax_(nIterationsMax),
210  verbosity_level_(verbosity_level),
211  mu_(0.0),
212  matrixOperator_(m_),
213  scalingOperator_(-1.0,mu_),
214  itOperator_(matrixOperator_,scalingOperator_),
215  nIterations_(0),
216  title_(" PowerIteration_Algorithms: "),
217  blank_(title_.length(),' ')
218  {
219  // assert that BCRSMatrix type has blocklevel 2
220  static_assert
221  (blockLevel<BCRSMatrix>() == 2,
222  "Only BCRSMatrices with blocklevel 2 are supported.");
223 
224  // assert that BCRSMatrix type has square blocks
225  static_assert
226  (BCRSMatrix::block_type::rows == BCRSMatrix::block_type::cols,
227  "Only BCRSMatrices with square blocks are supported.");
228 
229  // assert that m_ is square
230  const int nrows = m_.M() * BCRSMatrix::block_type::rows;
231  const int ncols = m_.N() * BCRSMatrix::block_type::cols;
232  if (nrows != ncols)
233  DUNE_THROW(Dune::ISTLError,"Matrix is not square ("
234  << nrows << "x" << ncols << ").");
235  }
236 
241 
247 
260  inline void applyPowerIteration (const Real& epsilon,
261  BlockVector& x, Real& lambda) const
262  {
263  // print verbosity information
264  if (verbosity_level_ > 0)
265  std::cout << title_
266  << "Performing power iteration approximating "
267  << "the dominant eigenvalue." << std::endl;
268 
269  // allocate memory for auxiliary variables
270  BlockVector y(x);
271  BlockVector temp(x);
272 
273  // perform power iteration
274  x *= (1.0 / x.two_norm());
275  m_.mv(x,y);
277  nIterations_ = 0;
278  while (r_norm > epsilon)
279  {
280  // update and check number of iterations
281  if (++nIterations_ > nIterationsMax_)
282  DUNE_THROW(Dune::ISTLError,"Power iteration did not converge "
283  << "in " << nIterationsMax_ << " iterations "
284  << "(║residual║_2 = " << r_norm << ", epsilon = "
285  << epsilon << ").");
286 
287  // do one iteration of the power iteration algorithm
288  // (use that y = m_ * x)
289  x = y;
290  x *= (1.0 / y.two_norm());
291 
292  // get approximated eigenvalue lambda via the Rayleigh quotient
293  m_.mv(x,y);
294  lambda = x * y;
295 
296  // get norm of residual (use that y = m_ * x)
297  temp = y;
298  temp.axpy(-lambda,x);
299  r_norm = temp.two_norm();
300 
301  // print verbosity information
302  if (verbosity_level_ > 1)
303  std::cout << blank_ << std::left
304  << "iteration " << std::setw(3) << nIterations_
305  << " (║residual║_2 = " << std::setw(11) << r_norm
306  << "): λ = " << lambda << std::endl
307  << std::resetiosflags(std::ios::left);
308  }
309 
310  // print verbosity information
311  if (verbosity_level_ > 0)
312  {
313  std::cout << blank_ << "Result ("
314  << "#iterations = " << nIterations_ << ", "
315  << "║residual║_2 = " << r_norm << "): "
316  << "λ = " << lambda << std::endl;
317  if (verbosity_level_ > 2)
318  {
319  // print approximated eigenvector via DUNE-ISTL I/O methods
320  Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
321  }
322  }
323  }
324 
353  template <typename ISTLLinearSolver,
354  bool avoidLinSolverCrime = false>
355  inline void applyInverseIteration (const Real& epsilon,
356  ISTLLinearSolver& solver,
357  BlockVector& x, Real& lambda) const
358  {
359  constexpr Real gamma = 0.0;
360  applyInverseIteration(gamma,epsilon,solver,x,lambda);
361  }
362 
392  template <typename ISTLLinearSolver,
393  bool avoidLinSolverCrime = false>
394  inline void applyInverseIteration (const Real& gamma,
395  const Real& epsilon,
396  ISTLLinearSolver& solver,
397  BlockVector& x, Real& lambda) const
398  {
399  // print verbosity information
400  if (verbosity_level_ > 0)
401  {
402  std::cout << title_;
403  if (gamma == 0.0)
404  std::cout << "Performing inverse iteration approximating "
405  << "the least dominant eigenvalue." << std::endl;
406  else
407  std::cout << "Performing inverse iteration with shift "
408  << "gamma = " << gamma << " approximating the "
409  << "eigenvalue closest to gamma." << std::endl;
410  }
411 
412  // initialize iteration operator,
413  // initialize iteration matrix when needed
414  updateShiftMu(gamma,solver);
415 
416  // allocate memory for linear solver statistics
417  Dune::InverseOperatorResult solver_statistics;
418 
419  // allocate memory for auxiliary variables
420  BlockVector y(x);
421  Real y_norm;
422  BlockVector temp(x);
423 
424  // perform inverse iteration with shift
425  x *= (1.0 / x.two_norm());
427  nIterations_ = 0;
428  while (r_norm > epsilon)
429  {
430  // update and check number of iterations
431  if (++nIterations_ > nIterationsMax_)
432  DUNE_THROW(Dune::ISTLError,"Inverse iteration "
433  << (gamma != 0.0 ? "with shift " : "") << "did not "
434  << "converge in " << nIterationsMax_ << " iterations "
435  << "(║residual║_2 = " << r_norm << ", epsilon = "
436  << epsilon << ").");
437 
438  // do one iteration of the inverse iteration with shift algorithm,
439  // part 1: solve (m_ - gamma*I) * y = x for y
440  // (protect x from being changed)
441  temp = x;
442  solver.apply(y,temp,solver_statistics);
443 
444  // get norm of y
445  y_norm = y.two_norm();
446 
447  // compile time switch between accuracy and efficiency
448  if (avoidLinSolverCrime)
449  {
450  // get approximated eigenvalue lambda via the Rayleigh quotient
451  // (use that x_new = y / y_norm)
452  m_.mv(y,temp);
453  lambda = (y * temp) / (y_norm * y_norm);
454 
455  // get norm of residual
456  // (use that x_new = y / y_norm, additionally use that temp = m_ * y)
457  temp.axpy(-lambda,y);
458  r_norm = temp.two_norm() / y_norm;
459  }
460  else
461  {
462  // get approximated eigenvalue lambda via the Rayleigh quotient
463  // (use that x_new = y / y_norm and use that (m_ - gamma*I) * y = x)
464  lambda = gamma + (y * x) / (y_norm * y_norm);
465 
466  // get norm of residual
467  // (use that x_new = y / y_norm and use that (m_ - gamma*I) * y = x)
468  temp = x; temp.axpy(gamma-lambda,y);
469  r_norm = temp.two_norm() / y_norm;
470  }
471 
472  // do one iteration of the inverse iteration with shift algorithm,
473  // part 2: update x
474  x = y;
475  x *= (1.0 / y_norm);
476 
477  // print verbosity information
478  if (verbosity_level_ > 1)
479  std::cout << blank_ << std::left
480  << "iteration " << std::setw(3) << nIterations_
481  << " (║residual║_2 = " << std::setw(11) << r_norm
482  << "): λ = " << lambda << std::endl
483  << std::resetiosflags(std::ios::left);
484  }
485 
486  // print verbosity information
487  if (verbosity_level_ > 0)
488  {
489  std::cout << blank_ << "Result ("
490  << "#iterations = " << nIterations_ << ", "
491  << "║residual║_2 = " << r_norm << "): "
492  << "λ = " << lambda << std::endl;
493  if (verbosity_level_ > 2)
494  {
495  // print approximated eigenvector via DUNE-ISTL I/O methods
496  Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
497  }
498  }
499  }
500 
531  template <typename ISTLLinearSolver,
532  bool avoidLinSolverCrime = false>
533  inline void applyRayleighQuotientIteration (const Real& epsilon,
534  ISTLLinearSolver& solver,
535  BlockVector& x, Real& lambda) const
536  {
537  // print verbosity information
538  if (verbosity_level_ > 0)
539  std::cout << title_
540  << "Performing Rayleigh quotient iteration for "
541  << "estimated eigenvalue " << lambda << "." << std::endl;
542 
543  // allocate memory for linear solver statistics
544  Dune::InverseOperatorResult solver_statistics;
545 
546  // allocate memory for auxiliary variables
547  BlockVector y(x);
548  Real y_norm;
549  Real lambda_update;
550  BlockVector temp(x);
551 
552  // perform Rayleigh quotient iteration
553  x *= (1.0 / x.two_norm());
555  nIterations_ = 0;
556  while (r_norm > epsilon)
557  {
558  // update and check number of iterations
559  if (++nIterations_ > nIterationsMax_)
560  DUNE_THROW(Dune::ISTLError,"Rayleigh quotient iteration did not "
561  << "converge in " << nIterationsMax_ << " iterations "
562  << "(║residual║_2 = " << r_norm << ", epsilon = "
563  << epsilon << ").");
564 
565  // update iteration operator,
566  // update iteration matrix when needed
567  updateShiftMu(lambda,solver);
568 
569  // do one iteration of the Rayleigh quotient iteration algorithm,
570  // part 1: solve (m_ - lambda*I) * y = x for y
571  // (protect x from being changed)
572  temp = x;
573  solver.apply(y,temp,solver_statistics);
574 
575  // get norm of y
576  y_norm = y.two_norm();
577 
578  // compile time switch between accuracy and efficiency
579  if (avoidLinSolverCrime)
580  {
581  // get approximated eigenvalue lambda via the Rayleigh quotient
582  // (use that x_new = y / y_norm)
583  m_.mv(y,temp);
584  lambda = (y * temp) / (y_norm * y_norm);
585 
586  // get norm of residual
587  // (use that x_new = y / y_norm, additionally use that temp = m_ * y)
588  temp.axpy(-lambda,y);
589  r_norm = temp.two_norm() / y_norm;
590  }
591  else
592  {
593  // get approximated eigenvalue lambda via the Rayleigh quotient
594  // (use that x_new = y / y_norm and use that (m_ - lambda_old*I) * y = x)
595  lambda_update = (y * x) / (y_norm * y_norm);
596  lambda += lambda_update;
597 
598  // get norm of residual
599  // (use that x_new = y / y_norm and use that (m_ - lambda_old*I) * y = x)
600  temp = x; temp.axpy(-lambda_update,y);
601  r_norm = temp.two_norm() / y_norm;
602  }
603 
604  // do one iteration of the Rayleigh quotient iteration algorithm,
605  // part 2: update x
606  x = y;
607  x *= (1.0 / y_norm);
608 
609  // print verbosity information
610  if (verbosity_level_ > 1)
611  std::cout << blank_ << std::left
612  << "iteration " << std::setw(3) << nIterations_
613  << " (║residual║_2 = " << std::setw(11) << r_norm
614  << "): λ = " << lambda << std::endl
615  << std::resetiosflags(std::ios::left);
616  }
617 
618  // print verbosity information
619  if (verbosity_level_ > 0)
620  {
621  std::cout << blank_ << "Result ("
622  << "#iterations = " << nIterations_ << ", "
623  << "║residual║_2 = " << r_norm << "): "
624  << "λ = " << lambda << std::endl;
625  if (verbosity_level_ > 2)
626  {
627  // print approximated eigenvector via DUNE-ISTL I/O methods
628  Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
629  }
630  }
631  }
632 
689  template <typename ISTLLinearSolver,
690  bool avoidLinSolverCrime = false>
691  inline void applyTLIMEIteration (const Real& gamma, const Real& eta,
692  const Real& epsilon,
693  ISTLLinearSolver& solver,
694  const Real& delta, const std::size_t& m,
695  bool& extrnl,
696  BlockVector& x, Real& lambda) const
697  {
698  // use same variable names as in [Szyld, 1988]
699  BlockVector& x_s = x;
700  Real& mu_s = lambda;
701 
702  // print verbosity information
703  if (verbosity_level_ > 0)
704  std::cout << title_
705  << "Performing TLIME iteration for "
706  << "estimated eigenvalue in the "
707  << "interval (" << gamma - eta << ","
708  << gamma + eta << ")." << std::endl;
709 
710  // allocate memory for linear solver statistics
711  Dune::InverseOperatorResult solver_statistics;
712 
713  // allocate memory for auxiliary variables
714  bool doRQI;
715  Real mu;
716  BlockVector y(x_s);
717  Real omega;
718  Real mu_s_old;
719  Real mu_s_update;
720  BlockVector temp(x_s);
721  Real q_norm, r_norm;
722 
723  // perform TLIME iteration
724  x_s *= (1.0 / x_s.two_norm());
725  extrnl = true;
726  doRQI = false;
728  nIterations_ = 0;
729  while (r_norm > epsilon)
730  {
731  // update and check number of iterations
732  if (++nIterations_ > nIterationsMax_)
733  DUNE_THROW(Dune::ISTLError,"TLIME iteration did not "
734  << "converge in " << nIterationsMax_
735  << " iterations (║residual║_2 = " << r_norm
736  << ", epsilon = " << epsilon << ").");
737 
738  // set shift for next iteration according to inverse iteration
739  // with shift (II) resp. Rayleigh quotient iteration (RQI)
740  if (doRQI)
741  mu = mu_s;
742  else
743  mu = gamma;
744 
745  // update II/RQI iteration operator,
746  // update II/RQI iteration matrix when needed
747  updateShiftMu(mu,solver);
748 
749  // do one iteration of the II/RQI algorithm,
750  // part 1: solve (m_ - mu*I) * y = x for y
751  temp = x_s;
752  solver.apply(y,temp,solver_statistics);
753 
754  // do one iteration of the II/RQI algorithm,
755  // part 2: compute omega
756  omega = (1.0 / y.two_norm());
757 
758  // backup the old Rayleigh quotient
759  mu_s_old = mu_s;
760 
761  // compile time switch between accuracy and efficiency
762  if (avoidLinSolverCrime)
763  {
764  // update the Rayleigh quotient mu_s, i.e. the approximated eigenvalue
765  // (use that x_new = y * omega)
766  m_.mv(y,temp);
767  mu_s = (y * temp) * (omega * omega);
768 
769  // get norm of "the residual with respect to the shift used by II",
770  // use normal representation of q
771  // (use that x_new = y * omega, use that temp = m_ * y)
772  temp.axpy(-gamma,y);
773  q_norm = temp.two_norm() * omega;
774 
775  // get norm of "the residual with respect to the Rayleigh quotient"
776  r_norm = q_norm*q_norm - (gamma-mu_s)*(gamma-mu_s);
777  // prevent that truncation errors invalidate the norm
778  // (we don't want to calculate sqrt of a negative number)
779  if (r_norm >= 0)
780  {
781  // use relation between the norms of r and q for efficiency
782  r_norm = std::sqrt(r_norm);
783  }
784  else
785  {
786  // use relation between r and q
787  // (use that x_new = y * omega, use that temp = (m_ - gamma*I) * y = q / omega)
788  temp.axpy(gamma-mu_s,y);
789  r_norm = temp.two_norm() * omega;
790  }
791  }
792  else
793  {
794  // update the Rayleigh quotient mu_s, i.e. the approximated eigenvalue
795  if (!doRQI)
796  {
797  // (use that x_new = y * omega, additionally use that (m_ - gamma*I) * y = x_s)
798  mu_s = gamma + (y * x_s) * (omega * omega);
799  }
800  else
801  {
802  // (use that x_new = y * omega, additionally use that (m_ - mu_s_old*I) * y = x_s)
803  mu_s_update = (y * x_s) * (omega * omega);
804  mu_s += mu_s_update;
805  }
806 
807  // get norm of "the residual with respect to the shift used by II"
808  if (!doRQI)
809  {
810  // use special representation of q in the II case
811  // (use that x_new = y * omega, additionally use that (m_ - gamma*I) * y = x_s)
812  q_norm = omega;
813  }
814  else
815  {
816  // use special representation of q in the RQI case
817  // (use that x_new = y * omega, additionally use that (m_ - mu_s_old*I) * y = x_s)
818  temp = x_s; temp.axpy(mu_s-gamma,y);
819  q_norm = temp.two_norm() * omega;
820  }
821 
822  // get norm of "the residual with respect to the Rayleigh quotient"
823  // don't use efficient relation between the norms of r and q, as
824  // this relation seems to yield a less accurate r_norm in the case
825  // where linear solver crime is admitted
826  if (!doRQI)
827  {
828  // (use that x_new = y * omega and use that (m_ - gamma*I) * y = x_s)
829  temp = x_s; temp.axpy(gamma-lambda,y);
830  r_norm = temp.two_norm() * omega;
831  }
832  else
833  {
834  // (use that x_new = y * omega and use that (m_ - mu_s_old*I) * y = x_s)
835  temp = x_s; temp.axpy(-mu_s_update,y);
836  r_norm = temp.two_norm() * omega;
837  }
838  }
839 
840  // do one iteration of the II/RQI algorithm,
841  // part 3: update x
842  x_s = y; x_s *= omega;
843 
844  // // for relative residual norm mode, scale with mu_s^{-1}
845  // r_norm /= std::abs(mu_s);
846 
847  // print verbosity information
848  if (verbosity_level_ > 1)
849  std::cout << blank_ << "iteration "
850  << std::left << std::setw(3) << nIterations_
851  << " (" << (doRQI ? "RQI," : "II, ")
852  << " " << (doRQI ? "—>" : " ") << " "
853  << "║r║_2 = " << std::setw(11) << r_norm
854  << ", " << (doRQI ? " " : "—>") << " "
855  << "║q║_2 = " << std::setw(11) << q_norm
856  << "): λ = " << lambda << std::endl
857  << std::resetiosflags(std::ios::left);
858 
859  // check if the eigenvalue closest to gamma lies in J
860  if (!doRQI && q_norm < eta)
861  {
862  // J is not free of eigenvalues
863  extrnl = false;
864 
865  // by theory we know now that mu_s also lies in J
866  assert(std::abs(mu_s-gamma) < eta);
867 
868  // switch to RQI
869  doRQI = true;
870  }
871 
872  // revert to II if J is not free of eigenvalues but
873  // at some point mu_s falls back again outside J
874  if (!extrnl && doRQI && std::abs(mu_s-gamma) >= eta)
875  doRQI = false;
876 
877  // if eigenvalue closest to gamma does not lie in J use RQI
878  // solely to accelerate the convergence to this eigenvalue
879  // when II has become stationary
880  if (extrnl && !doRQI)
881  {
882  // switch to RQI if the relative change of the Rayleigh
883  // quotient indicates that II has become stationary
884  if (nIterations_ >= m &&
885  std::abs(mu_s - mu_s_old) / std::abs(mu_s) < delta)
886  doRQI = true;
887  }
888  }
889 
890  // // compute final residual and lambda again (paranoia....)
891  // m_.mv(x_s,temp);
892  // mu_s = x_s * temp;
893  // temp.axpy(-mu_s,x_s);
894  // r_norm = temp.two_norm();
895  // // r_norm /= std::abs(mu_s);
896 
897  // print verbosity information
898  if (verbosity_level_ > 0)
899  {
900  if (extrnl)
901  std::cout << blank_ << "Interval "
902  << "(" << gamma - eta << "," << gamma + eta
903  << ") is free of eigenvalues, approximating "
904  << "the closest eigenvalue." << std::endl;
905  std::cout << blank_ << "Result ("
906  << "#iterations = " << nIterations_ << ", "
907  << "║residual║_2 = " << r_norm << "): "
908  << "λ = " << lambda << std::endl;
909  if (verbosity_level_ > 2)
910  {
911  // print approximated eigenvector via DUNE-ISTL I/O methods
912  Dune::printvector(std::cout,x,blank_+"x",blank_+"row");
913  }
914  }
915  }
916 
926  {
927  // return iteration operator
928  return itOperator_;
929  }
930 
945  inline const BCRSMatrix& getIterationMatrix () const
946  {
947  // create iteration matrix on demand
948  if (!itMatrix_)
949  itMatrix_ = std::make_unique<BCRSMatrix>(m_);
950 
951  // return iteration matrix
952  return *itMatrix_;
953  }
954 
959  inline unsigned int getIterationCount () const
960  {
961  if (nIterations_ == 0)
962  DUNE_THROW(Dune::ISTLError,"No algorithm applied, yet.");
963 
964  return nIterations_;
965  }
966 
967  protected:
982  template <typename ISTLLinearSolver>
983  inline void updateShiftMu (const Real& mu,
984  ISTLLinearSolver& solver) const
985  {
986  // do nothing if new shift equals the old one
987  if (mu == mu_) return;
988 
989  // update shift mu_, i.e. update iteration operator
990  mu_ = mu;
991 
992  // update iteration matrix when needed
993  if (itMatrix_)
994  {
995  // iterate over entries in iteration matrix diagonal
996  constexpr int rowBlockSize = BCRSMatrix::block_type::rows;
997  constexpr int colBlockSize = BCRSMatrix::block_type::cols;
998  for (typename BCRSMatrix::size_type i = 0;
999  i < itMatrix_->M()*rowBlockSize; ++i)
1000  {
1001  // access m_[i,i] where i is the flat index of a row/column
1002  const Real& m_entry = m_
1003  [i/rowBlockSize][i/colBlockSize][i%rowBlockSize][i%colBlockSize];
1004  // access *itMatrix[i,i] where i is the flat index of a row/column
1005  Real& entry = (*itMatrix_)
1006  [i/rowBlockSize][i/colBlockSize][i%rowBlockSize][i%colBlockSize];
1007  // change current entry in iteration matrix diagonal
1008  entry = m_entry - mu_;
1009  }
1010  // notify linear solver about change of the iteration matrix object
1012  (solver,*itMatrix_);
1013  }
1014  }
1015 
1016  protected:
1017  // parameters related to iterative eigenvalue algorithms
1018  const BCRSMatrix& m_;
1019  const unsigned int nIterationsMax_;
1020 
1021  // verbosity setting
1022  const unsigned int verbosity_level_;
1023 
1024  // shift mu_ used by iteration operator/matrix (m_ - mu_*I)
1025  mutable Real mu_;
1026 
1027  // iteration operator (m_ - mu_*I), passing shift mu_ by reference
1028  const MatrixOperator matrixOperator_;
1029  const ScalingOperator scalingOperator_;
1030  OperatorSum itOperator_;
1031 
1032  // iteration matrix (m_ - mu_*I), provided on demand when needed
1033  // (e.g. for preconditioning)
1034  mutable std::unique_ptr<BCRSMatrix> itMatrix_;
1035 
1036  // memory for storing temporary variables (mutable as they shall
1037  // just be effectless auxiliary variables of the const apply*(...)
1038  // methods)
1039  mutable unsigned int nIterations_;
1040 
1041  // constants for printing verbosity information
1042  const std::string title_;
1043  const std::string blank_;
1044  };
1045 
1048 } // namespace Dune
1049 
1050 #endif // DUNE_ISTL_EIGENVALUE_POWERITERATION_HH
Helper functions for determining the vector/matrix block level.
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:466
A::size_type size_type
The type for the index access and the size.
Definition: bcrsmatrix.hh:497
size_type M() const
number of columns (counted in blocks)
Definition: bcrsmatrix.hh:2007
void mv(const X &x, Y &y) const
y = A x
Definition: bcrsmatrix.hh:1641
size_type N() const
number of rows (counted in blocks)
Definition: bcrsmatrix.hh:2001
A vector of blocks with memory management.
Definition: bvector.hh:392
typename Imp::BlockTraits< B >::field_type field_type
export the type representing the field
Definition: bvector.hh:398
derive error class from the base class in common
Definition: istlexception.hh:19
A linear operator.
Definition: operators.hh:67
X::field_type field_type
The field type of the operator.
Definition: operators.hh:74
virtual void applyscaleadd(field_type alpha, const X &x, X &y) const=0
apply operator to x, scale and add:
virtual SolverCategory::Category category() const=0
Category of the linear operator (see SolverCategory::Category)
X range_type
The type of the range of the operator.
Definition: operators.hh:72
virtual void apply(const X &x, X &y) const=0
apply operator to x: The input vector is consistent and the output must also be consistent on the in...
X domain_type
The type of the domain of the operator.
Definition: operators.hh:70
Adapter to turn a matrix into a linear operator.
Definition: operators.hh:137
Iterative eigenvalue algorithms based on power iteration.
Definition: poweriteration.hh:176
PowerIteration_Algorithms(const PowerIteration_Algorithms &)=delete
void applyInverseIteration(const Real &epsilon, ISTLLinearSolver &solver, BlockVector &x, Real &lambda) const
Perform the inverse iteration algorithm to compute an approximation lambda of the least dominant (i....
Definition: poweriteration.hh:355
void applyTLIMEIteration(const Real &gamma, const Real &eta, const Real &epsilon, ISTLLinearSolver &solver, const Real &delta, const std::size_t &m, bool &extrnl, BlockVector &x, Real &lambda) const
Perform the "two-level iterative method for eigenvalue calculations (TLIME)" iteration algorit...
Definition: poweriteration.hh:691
PowerIteration_Algorithms(const BCRSMatrix &m, const unsigned int nIterationsMax=1000, const unsigned int verbosity_level=0)
Construct from required parameters.
Definition: poweriteration.hh:206
void applyPowerIteration(const Real &epsilon, BlockVector &x, Real &lambda) const
Perform the power iteration algorithm to compute an approximation lambda of the dominant (i....
Definition: poweriteration.hh:260
OperatorSum IterationOperator
Type of iteration operator (m_ - mu_*I)
Definition: poweriteration.hh:189
void applyRayleighQuotientIteration(const Real &epsilon, ISTLLinearSolver &solver, BlockVector &x, Real &lambda) const
Perform the Rayleigh quotient iteration algorithm to compute an approximation lambda of an eigenvalue...
Definition: poweriteration.hh:533
void applyInverseIteration(const Real &gamma, const Real &epsilon, ISTLLinearSolver &solver, BlockVector &x, Real &lambda) const
Perform the inverse iteration with shift algorithm to compute an approximation lambda of the eigenval...
Definition: poweriteration.hh:394
void updateShiftMu(const Real &mu, ISTLLinearSolver &solver) const
Update shift mu_, i.e. update iteration operator/matrix (m_ - mu_*I).
Definition: poweriteration.hh:983
const BCRSMatrix & getIterationMatrix() const
Return the iteration matrix (m_ - mu_*I), provided on demand when needed (e.g. for direct solvers or ...
Definition: poweriteration.hh:945
unsigned int getIterationCount() const
Return the number of iterations in last application of an algorithm.
Definition: poweriteration.hh:959
IterationOperator & getIterationOperator()
Return the iteration operator (m_ - mu_*I).
Definition: poweriteration.hh:925
PowerIteration_Algorithms & operator=(const PowerIteration_Algorithms &)=delete
BlockVector::field_type Real
Type of underlying field.
Definition: poweriteration.hh:186
Helper class for notifying a DUNE-ISTL linear solver about a change of the iteration matrix object in...
Definition: solver.hh:522
A few common exception classes.
Some generic functions for pretty printing vectors and matrices.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:484
void printvector(std::ostream &s, const V &v, std::string title, std::string rowtext, int columns=1, int width=10, int precision=2)
Print an ISTL vector.
Definition: io.hh:89
Dune namespace.
Definition: alignedallocator.hh:13
Define general, extensible interface for operators. The available implementation wraps a matrix.
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
Statistics about the application of an inverse operator.
Definition: solver.hh:48
Category
Definition: solvercategory.hh:23
@ sequential
Category for sequential solvers.
Definition: solvercategory.hh:25
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 26, 22:29, 2024)