Reference documentation for deal.II version 8.5.1
slepc_solver.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 2016 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 
17 #ifndef dealii__slepc_solver_h
18 #define dealii__slepc_solver_h
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_SLEPC
23 
24 # include <deal.II/base/std_cxx11/shared_ptr.h>
25 # include <deal.II/lac/exceptions.h>
26 # include <deal.II/lac/solver_control.h>
27 # include <deal.II/lac/petsc_matrix_base.h>
28 # include <deal.II/lac/slepc_spectral_transformation.h>
29 
30 # include <petscconf.h>
31 # include <petscksp.h>
32 # include <slepceps.h>
33 
34 DEAL_II_NAMESPACE_OPEN
35 
130 namespace SLEPcWrappers
131 {
132 
141  {
142  public:
148  const MPI_Comm &mpi_communicator);
149 
153  virtual ~SolverBase ();
154 
173  template <typename OutputVector>
174  void
176  std::vector<PetscScalar> &eigenvalues,
177  std::vector<OutputVector> &eigenvectors,
178  const unsigned int n_eigenpairs = 1);
179 
185  template <typename OutputVector>
186  void
188  const PETScWrappers::MatrixBase &B,
189  std::vector<PetscScalar> &eigenvalues,
190  std::vector<OutputVector> &eigenvectors,
191  const unsigned int n_eigenpairs = 1);
192 
198  template <typename OutputVector>
199  void
201  const PETScWrappers::MatrixBase &B,
202  std::vector<double> &real_eigenvalues,
203  std::vector<double> &imag_eigenvalues,
204  std::vector<OutputVector> &real_eigenvectors,
205  std::vector<OutputVector> &imag_eigenvectors,
206  const unsigned int n_eigenpairs = 1);
207 
211  void
213  (const PETScWrappers::VectorBase &this_initial_vector) DEAL_II_DEPRECATED;
214 
221  template <typename Vector>
222  void
224  (const std::vector<Vector> &initial_space);
225 
229  void
231 
236  void
237  set_target_eigenvalue (const PetscScalar &this_target);
238 
245  void
246  set_which_eigenpairs (EPSWhich set_which);
247 
256  void
257  set_problem_type (EPSProblemType set_problem);
258 
264  void
266 
271 
276  int,
277  << " An error with error number " << arg1
278  << " occurred while calling a SLEPc function");
279 
284  int, int,
285  << " The number of converged eigenvectors is " << arg1
286  << " but " << arg2 << " were requested. ");
287 
291  SolverControl &control () const;
292 
293  protected:
294 
300 
304  const MPI_Comm mpi_communicator;
305 
312  void
313  solve (const unsigned int n_eigenpairs,
314  unsigned int *n_converged);
315 
321  void
322  get_eigenpair (const unsigned int index,
323  PetscScalar &eigenvalues,
324  PETScWrappers::VectorBase &eigenvectors);
325 
331  void
332  get_eigenpair (const unsigned int index,
333  double &real_eigenvalues,
334  double &imag_eigenvalues,
335  PETScWrappers::VectorBase &real_eigenvectors,
336  PETScWrappers::VectorBase &imag_eigenvectors);
337 
342  void
344 
349  void
351  const PETScWrappers::MatrixBase &B);
352 
353  protected:
354 
358  EPS eps;
359 
360  private:
364  EPSConvergedReason reason;
365 
366 
373  static
374  int
375  convergence_test (EPS eps,
376  PetscScalar real_eigenvalue,
377  PetscScalar imag_eigenvalue,
378  PetscReal residual_norm,
379  PetscReal *estimated_error,
380  void *solver_control);
381  };
382 
392  {
393  public:
394 
400  {};
401 
408  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
409  const AdditionalData &data = AdditionalData());
410 
411  protected:
412 
417  };
418 
427  class SolverArnoldi : public SolverBase
428  {
429  public:
435  {
440  AdditionalData (const bool delayed_reorthogonalization = false);
441 
446  };
447 
454  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
455  const AdditionalData &data = AdditionalData());
456 
457  protected:
458 
463  };
464 
473  class SolverLanczos : public SolverBase
474  {
475  public:
481  {
485  EPSLanczosReorthogType reorthog;
486 
491  AdditionalData(const EPSLanczosReorthogType r = EPS_LANCZOS_REORTHOG_FULL);
492  };
493 
500  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
501  const AdditionalData &data = AdditionalData());
502 
503  protected:
504 
509  };
510 
519  class SolverPower : public SolverBase
520  {
521  public:
527  {};
528 
535  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
536  const AdditionalData &data = AdditionalData());
537 
538  protected:
539 
544  };
545 
555  {
556  public:
562  {
567 
571  AdditionalData(bool double_expansion = false);
572  };
573 
580  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
581  const AdditionalData &data = AdditionalData());
582 
583  protected:
584 
589  };
590 
600  {
601  public:
607  {};
608 
615  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
616  const AdditionalData &data = AdditionalData());
617 
618  protected:
619 
624  };
625 
626 
635  class SolverLAPACK : public SolverBase
636  {
637  public:
638 
644  {};
645 
652  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
653  const AdditionalData &data = AdditionalData());
654 
655  protected:
656 
661  };
662 
663  // --------------------------- inline and template functions -----------
668  // todo: The logic of these functions can be simplified without breaking backward compatibility...
669 
670  template <typename OutputVector>
671  void
673  std::vector<PetscScalar> &eigenvalues,
674  std::vector<OutputVector> &eigenvectors,
675  const unsigned int n_eigenpairs)
676  {
677  // Panic if the number of eigenpairs wanted is out of bounds.
678  AssertThrow ((n_eigenpairs > 0) && (n_eigenpairs <= A.m ()),
680 
681  // Set the matrices of the problem
682  set_matrices (A);
683 
684  // and solve
685  unsigned int n_converged = 0;
686  solve (n_eigenpairs, &n_converged);
687 
688  if (n_converged > n_eigenpairs)
689  n_converged = n_eigenpairs;
690  AssertThrow (n_converged == n_eigenpairs,
691  ExcSLEPcEigenvectorConvergenceMismatchError(n_converged, n_eigenpairs));
692 
693  AssertThrow (eigenvectors.size() != 0, ExcSLEPcWrappersUsageError());
694  eigenvectors.resize (n_converged, eigenvectors.front());
695  eigenvalues.resize (n_converged);
696 
697  for (unsigned int index=0; index<n_converged; ++index)
698  get_eigenpair (index, eigenvalues[index], eigenvectors[index]);
699  }
700 
701  template <typename OutputVector>
702  void
704  const PETScWrappers::MatrixBase &B,
705  std::vector<PetscScalar> &eigenvalues,
706  std::vector<OutputVector> &eigenvectors,
707  const unsigned int n_eigenpairs)
708  {
709  // Guard against incompatible matrix sizes:
710  AssertThrow (A.m() == B.m (), ExcDimensionMismatch(A.m(), B.m()));
711  AssertThrow (A.n() == B.n (), ExcDimensionMismatch(A.n(), B.n()));
712 
713  // Panic if the number of eigenpairs wanted is out of bounds.
714  AssertThrow ((n_eigenpairs>0) && (n_eigenpairs<=A.m ()),
716 
717  // Set the matrices of the problem
718  set_matrices (A, B);
719 
720  // and solve
721  unsigned int n_converged = 0;
722  solve (n_eigenpairs, &n_converged);
723 
724  if (n_converged>=n_eigenpairs)
725  n_converged = n_eigenpairs;
726 
727  AssertThrow (n_converged==n_eigenpairs,
728  ExcSLEPcEigenvectorConvergenceMismatchError(n_converged, n_eigenpairs));
729  AssertThrow (eigenvectors.size() != 0, ExcSLEPcWrappersUsageError());
730 
731  eigenvectors.resize (n_converged, eigenvectors.front());
732  eigenvalues.resize (n_converged);
733 
734  for (unsigned int index=0; index<n_converged; ++index)
735  get_eigenpair (index, eigenvalues[index], eigenvectors[index]);
736  }
737 
738  template <typename OutputVector>
739  void
741  const PETScWrappers::MatrixBase &B,
742  std::vector<double> &real_eigenvalues,
743  std::vector<double> &imag_eigenvalues,
744  std::vector<OutputVector> &real_eigenvectors,
745  std::vector<OutputVector> &imag_eigenvectors,
746  const unsigned int n_eigenpairs)
747  {
748  // Guard against incompatible matrix sizes:
749  AssertThrow (A.m() == B.m (), ExcDimensionMismatch(A.m(), B.m()));
750  AssertThrow (A.n() == B.n (), ExcDimensionMismatch(A.n(), B.n()));
751 
752  // and incompatible eigenvalue/eigenvector sizes
753  AssertThrow (real_eigenvalues.size() == imag_eigenvalues.size(),
754  ExcDimensionMismatch(real_eigenvalues.size(), imag_eigenvalues.size()));
755  AssertThrow (real_eigenvectors.size() == imag_eigenvectors.size (),
756  ExcDimensionMismatch(real_eigenvectors.size(), imag_eigenvectors.size()));
757 
758  // Panic if the number of eigenpairs wanted is out of bounds.
759  AssertThrow ((n_eigenpairs>0) && (n_eigenpairs<=A.m ()),
761 
762  // Set the matrices of the problem
763  set_matrices (A, B);
764 
765  // and solve
766  unsigned int n_converged = 0;
767  solve (n_eigenpairs, &n_converged);
768 
769  if (n_converged>=n_eigenpairs)
770  n_converged = n_eigenpairs;
771 
772  AssertThrow (n_converged==n_eigenpairs,
773  ExcSLEPcEigenvectorConvergenceMismatchError(n_converged, n_eigenpairs));
774  AssertThrow ((real_eigenvectors.size()!=0) && (imag_eigenvectors.size()!=0),
776 
777  real_eigenvectors.resize (n_converged, real_eigenvectors.front());
778  imag_eigenvectors.resize (n_converged, imag_eigenvectors.front());
779  real_eigenvalues.resize (n_converged);
780  imag_eigenvalues.resize (n_converged);
781 
782  for (unsigned int index=0; index<n_converged; ++index)
783  get_eigenpair (index,
784  real_eigenvalues[index], imag_eigenvalues[index],
785  real_eigenvectors[index], imag_eigenvectors[index]);
786  }
787 
788  template <typename Vector>
789  void
790  SolverBase::set_initial_space(const std::vector<Vector> &this_initial_space)
791  {
792  std::vector<Vec> vecs(this_initial_space.size());
793 
794  for (unsigned int i = 0; i < this_initial_space.size(); i++)
795  {
796  Assert(this_initial_space[i].l2_norm()>0.0,
797  ExcMessage("Initial vectors should be nonzero."));
798  vecs[i] = this_initial_space[i];
799  }
800 
801  // if the eigensolver supports only a single initial vector, but several
802  // guesses are provided, then all except the first one will be discarded.
803  // One could still build a vector that is rich in the directions of all guesses,
804  // by taking a linear combination of them. (TODO: make function virtual?)
805 
806 #if DEAL_II_PETSC_VERSION_LT(3,1,0)
807  const PetscErrorCode ierr = EPSSetInitialVector (eps, &vecs[0]);
808 #else
809  const PetscErrorCode ierr = EPSSetInitialSpace (eps, vecs.size(), &vecs[0]);
810 #endif
811  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
812  }
813 
814 }
815 
816 DEAL_II_NAMESPACE_CLOSE
817 
818 #endif // DEAL_II_WITH_SLEPC
819 
820 /*---------------------------- slepc_solver.h ---------------------------*/
821 
822 #endif
823 
824 /*---------------------------- slepc_solver.h ---------------------------*/
static ::ExceptionBase & ExcSLEPcWrappersUsageError()
void set_initial_space(const std::vector< Vector > &initial_space)
Definition: slepc_solver.h:790
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
const AdditionalData additional_data
Definition: slepc_solver.h:543
SolverJacobiDavidson(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverPower(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
void solve(const PETScWrappers::MatrixBase &A, std::vector< PetscScalar > &eigenvalues, std::vector< OutputVector > &eigenvectors, const unsigned int n_eigenpairs=1)
Definition: slepc_solver.h:672
const MPI_Comm mpi_communicator
Definition: slepc_solver.h:304
SolverKrylovSchur(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverBase(SolverControl &cn, const MPI_Comm &mpi_communicator)
Definition: slepc_solver.cc:36
AdditionalData(const bool delayed_reorthogonalization=false)
void set_problem_type(EPSProblemType set_problem)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:564
void set_target_eigenvalue(const PetscScalar &this_target)
void set_matrices(const PETScWrappers::MatrixBase &A)
Definition: slepc_solver.cc:75
void set_which_eigenpairs(EPSWhich set_which)
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
SolverGeneralizedDavidson(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: slepc_solver.h:660
#define DeclException0(Exception0)
Definition: exceptions.h:541
EPSConvergedReason reason
Definition: slepc_solver.h:364
const AdditionalData additional_data
Definition: slepc_solver.h:508
const AdditionalData additional_data
Definition: slepc_solver.h:416
void get_eigenpair(const unsigned int index, PetscScalar &eigenvalues, PETScWrappers::VectorBase &eigenvectors)
static int convergence_test(EPS eps, PetscScalar real_eigenvalue, PetscScalar imag_eigenvalue, PetscReal residual_norm, PetscReal *estimated_error, void *solver_control)
SolverLanczos(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
static ::ExceptionBase & ExcSLEPcEigenvectorConvergenceMismatchError(int arg1, int arg2)
SolverLAPACK(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: slepc_solver.h:623
void set_initial_vector(const PETScWrappers::VectorBase &this_initial_vector) 1
AdditionalData(const EPSLanczosReorthogType r=EPS_LANCZOS_REORTHOG_FULL)
const AdditionalData additional_data
Definition: slepc_solver.h:462
static ::ExceptionBase & ExcSLEPcError(int arg1)
void set_transformation(SLEPcWrappers::TransformationBase &this_transformation)
Definition: slepc_solver.cc:92
SolverControl & solver_control
Definition: slepc_solver.h:299
SolverArnoldi(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverControl & control() const
void get_solver_state(const SolverControl::State state)