Reference documentation for deal.II version 8.5.1
slepc_solver.cc
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 #include <deal.II/lac/slepc_solver.h>
17 
18 #ifdef DEAL_II_WITH_SLEPC
19 
20 # include <deal.II/lac/petsc_matrix_base.h>
21 # include <deal.II/lac/petsc_vector_base.h>
22 # include <deal.II/lac/petsc_vector.h>
23 # include <deal.II/lac/slepc_spectral_transformation.h>
24 
25 # include <cmath>
26 # include <vector>
27 
28 # include <petscversion.h>
29 # include <slepcversion.h>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 namespace SLEPcWrappers
34 {
35 
37  const MPI_Comm &mpi_communicator)
38  :
39  solver_control (cn),
40  mpi_communicator (mpi_communicator)
41  {
42  // create eigensolver context
43  PetscErrorCode ierr = EPSCreate (mpi_communicator, &eps);
44  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
45 
46  // hand over the absolute tolerance and the maximum number of
47  // iteration steps to the SLEPc convergence criterion.
48  ierr = EPSSetTolerances(eps, this->solver_control.tolerance(),
49  this->solver_control.max_steps());
50  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
51 
52  // default values:
53  set_which_eigenpairs(EPS_LARGEST_MAGNITUDE);
54  set_problem_type(EPS_GNHEP);
55 
56  // TODO:
57  // By default, EPS initializes the starting vector or the initial subspace randomly.
58  }
59 
61  {
62  if (eps != NULL)
63  {
64  // Destroy the solver object.
65 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
66  const PetscErrorCode ierr = EPSDestroy (eps);
67 #else
68  const PetscErrorCode ierr = EPSDestroy (&eps);
69 #endif
70  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
71  }
72  }
73 
74  void
76  {
77  // standard eigenspectrum problem
78  const PetscErrorCode ierr = EPSSetOperators (eps, A, PETSC_NULL);
79  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
80  }
81 
82  void
85  {
86  // generalized eigenspectrum problem
87  const PetscErrorCode ierr = EPSSetOperators (eps, A, B);
88  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
89  }
90 
91  void
93  {
94  // set transformation type if any
95  // STSetShift is called inside
96  const PetscErrorCode ierr = EPSSetST(eps,transformation.st);
97  AssertThrow (ierr == 0, SolverBase::ExcSLEPcError(ierr));
98  }
99 
100  void
102  {
103  Assert(this_initial_vector.l2_norm()>0.0,
104  ExcMessage("Initial vector should be nonzero."));
105 
106  Vec vec = this_initial_vector;
107 #if DEAL_II_PETSC_VERSION_LT(3,1,0)
108  const PetscErrorCode ierr = EPSSetInitialVector (eps, &vec);
109 #else
110  const PetscErrorCode ierr = EPSSetInitialSpace (eps, 1, &vec);
111 #endif
112  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
113  }
114 
115  void
116  SolverBase::set_target_eigenvalue (const PetscScalar &this_target)
117  {
118  // set target eigenvalues to solve for
119  // in all transformation except STSHIFT there is a direct connection between
120  // the target and the shift, read more on p41 of SLEPc manual.
121  const PetscErrorCode ierr = EPSSetTarget (eps, this_target );
122  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
123  }
124 
125  void
126  SolverBase::set_which_eigenpairs (const EPSWhich eps_which)
127  {
128  // set which portion of the eigenspectrum to solve for
129  const PetscErrorCode ierr = EPSSetWhichEigenpairs (eps, eps_which);
130  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
131  }
132 
133  void
134  SolverBase::set_problem_type (const EPSProblemType eps_problem)
135  {
136  const PetscErrorCode ierr = EPSSetProblemType (eps, eps_problem);
137  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
138  }
139 
140  void
141  SolverBase::solve (const unsigned int n_eigenpairs,
142  unsigned int *n_converged)
143  {
144  // set number of eigenvectors to compute
145  PetscErrorCode ierr = EPSSetDimensions (eps, n_eigenpairs,
146  PETSC_DECIDE, PETSC_DECIDE);
147  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
148 
149  // set the solve options to the eigenvalue problem solver context
150  ierr = EPSSetFromOptions (eps);
151  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
152 
153  // TODO breaks @ref step_36 "step-36"
154  // force Krylov solver to use true residual instead of an estimate.
155  //EPSSetTrueResidual(solver_data->eps, PETSC_TRUE);
156  //AssertThrow (ierr == 0, ExcSLEPcError(ierr));
157 
158  // Set convergence test to be absolute
159  ierr = EPSSetConvergenceTest (eps, EPS_CONV_ABS);
160  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
161 
162  // TODO Set the convergence test function
163  // ierr = EPSSetConvergenceTestFunction (solver_data->eps, &convergence_test,
164  // reinterpret_cast<void *>(&solver_control));
165  // AssertThrow (ierr == 0, ExcSLEPcError(ierr));
166 
167  // solve the eigensystem
168  ierr = EPSSolve (eps);
169  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
170 
171  // get number of converged eigenstates
172  ierr = EPSGetConverged (eps,
173  reinterpret_cast<PetscInt *>(n_converged));
174  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
175 
176  PetscInt n_iterations = 0;
177  PetscReal residual_norm = 0;
178 
179  // @todo Investigate elaborating on some of this to act on the
180  // complete eigenspectrum
181  {
182  // get the number of solver iterations
183  ierr = EPSGetIterationNumber (eps, &n_iterations);
184  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
185 
186  // get the maximum of residual norm among converged eigenvectors.
187  for (unsigned int i = 0; i < *n_converged; i++)
188  {
189  double residual_norm_i = 0.0;
190  // EPSComputeResidualNorm is L2-norm and is not consistent with the stopping criteria
191  // used during the solution process.
192  // Yet, this is the norm which gives error bounds (Saad, 1992, ch3):
193  // | \lambda - \widehat\lambda | <= ||r||_2
194  ierr = EPSComputeResidualNorm (eps, i, &residual_norm_i);
195 
196  // EPSComputeRelativeError may not be consistent with the stopping criteria
197  // used during the solution process. Given EPS_CONV_ABS set above,
198  // this can be either the l2 norm or the mass-matrix induced norm
199  // when EPS_GHEP is set.
200  // ierr = EPSComputeRelativeError (solver_data->eps, i, &residual_norm_i);
201 
202  // EPSGetErrorEstimate is consistent with the residual norm
203  // used during the solution process. However, it is not guaranteed to
204  // be derived from the residual even when EPSSetTrueResidual is set.
205  // ierr = EPSGetErrorEstimate (solver_data->eps, i, &residual_norm_i);
206 
207  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
208  residual_norm = std::max (residual_norm, residual_norm_i);
209  }
210 
211  // check the solver state
212  const SolverControl::State state
213  = solver_control.check (n_iterations, residual_norm);
214 
215  // get the solver state according to SLEPc
216  get_solver_state (state);
217 
218  // as SLEPc uses different stopping criteria, we have to omit this step.
219  // This can be checked only in conjunction with EPSGetErrorEstimate.
220  // and in case of failure: throw exception
221  // if (solver_control.last_check () != SolverControl::success)
222  // AssertThrow(false, SolverControl::NoConvergence (solver_control.last_step(),
223  // solver_control.last_value()));
224  }
225  }
226 
227  void
228  SolverBase::get_eigenpair (const unsigned int index,
229  PetscScalar &eigenvalues,
230  PETScWrappers::VectorBase &eigenvectors)
231  {
232  // get converged eigenpair
233  const PetscErrorCode ierr = EPSGetEigenpair (eps, index,
234  &eigenvalues, PETSC_NULL,
235  eigenvectors, PETSC_NULL);
236  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
237  }
238 
239 
240  void
241  SolverBase::get_eigenpair (const unsigned int index,
242  double &real_eigenvalues,
243  double &imag_eigenvalues,
244  PETScWrappers::VectorBase &real_eigenvectors,
245  PETScWrappers::VectorBase &imag_eigenvectors)
246  {
247 #ifndef PETSC_USE_COMPLEX
248  // get converged eigenpair
249  const PetscErrorCode ierr = EPSGetEigenpair (eps, index,
250  &real_eigenvalues,
251  &imag_eigenvalues,
252  real_eigenvectors,
253  imag_eigenvectors);
254  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
255 #else
256  Assert ((false),
257  ExcMessage ("Your PETSc/SLEPc installation was configured with scalar-type complex "
258  "but this function is not defined for complex types."));
259 
260  // Cast to void to silence compiler warnings
261  (void) index;
262  (void) real_eigenvalues;
263  (void) imag_eigenvalues;
264  (void) real_eigenvectors;
265  (void) imag_eigenvectors;
266 #endif
267  }
268 
269  void
271  {
272  switch (state)
273  {
274  case ::SolverControl::iterate:
275  reason = EPS_CONVERGED_ITERATING;
276  break;
277 
278  case ::SolverControl::success:
279  reason = static_cast<EPSConvergedReason>(1);
280  break;
281 
282  case ::SolverControl::failure:
284  reason = EPS_DIVERGED_ITS;
285  else
286  reason = EPS_DIVERGED_BREAKDOWN;
287  break;
288 
289  default:
290  Assert (false, ExcNotImplemented());
291  }
292  }
293 
294  /* ---------------------- SolverControls ----------------------- */
295  SolverControl &
297  {
298  return solver_control;
299  }
300 
301  int
303  PetscScalar /*real_eigenvalue */,
304  PetscScalar /*imag_eigenvalue */,
305  PetscReal /*residual norm associated to the eigenpair */,
306  PetscReal */*(output) computed error estimate */,
307  void */*solver_control_x*/)
308  {
309  // If the error estimate returned by the convergence test function is less
310  // than the tolerance, then the eigenvalue is accepted as converged.
311  // This function is undefined (future reference only).
312 
313  // return without failure.
314  return 0;
315  }
316 
317  /* ---------------------- SolverKrylovSchur ------------------------ */
319  const MPI_Comm &mpi_communicator,
320  const AdditionalData &data)
321  :
322  SolverBase (cn, mpi_communicator),
323  additional_data (data)
324  {
325  const PetscErrorCode ierr = EPSSetType (eps,
326  const_cast<char *>(EPSKRYLOVSCHUR));
327  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
328  }
329 
330  /* ---------------------- SolverArnoldi ------------------------ */
332  AdditionalData (const bool delayed_reorthogonalization)
333  :
334  delayed_reorthogonalization (delayed_reorthogonalization)
335  {}
336 
338  const MPI_Comm &mpi_communicator,
339  const AdditionalData &data)
340  :
342  additional_data (data)
343  {
344  PetscErrorCode ierr = EPSSetType (eps, const_cast<char *>(EPSARNOLDI));
345  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
346 
347  // if requested, set delayed reorthogonalization in the Arnoldi
348  // iteration.
350  {
351  ierr = EPSArnoldiSetDelayed (eps, PETSC_TRUE);
352  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
353  }
354  }
355 
356 
357  /* ---------------------- Lanczos ------------------------ */
359  AdditionalData(const EPSLanczosReorthogType r)
360  : reorthog(r)
361  {}
362 
364  const MPI_Comm &mpi_communicator,
365  const AdditionalData &data)
366  :
368  additional_data (data)
369  {
370  PetscErrorCode ierr = EPSSetType (eps, const_cast<char *>(EPSLANCZOS));
371  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
372 
373  ierr = EPSLanczosSetReorthog(eps,additional_data.reorthog);
374  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
375  }
376 
377  /* ----------------------- Power ------------------------- */
379  const MPI_Comm &mpi_communicator,
380  const AdditionalData &data)
381  :
382  SolverBase (cn, mpi_communicator),
383  additional_data (data)
384  {
385  PetscErrorCode ierr = EPSSetType (eps, const_cast<char *>(EPSPOWER));
386  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
387  }
388 
389  /* ---------------- Generalized Davidson ----------------- */
391  AdditionalData(bool double_expansion)
392  : double_expansion(double_expansion)
393  {}
394 
396  const MPI_Comm &mpi_communicator,
397  const AdditionalData &data)
398  :
400  additional_data (data)
401  {
402 #if DEAL_II_PETSC_VERSION_GTE(3,1,0)
403  PetscErrorCode ierr = EPSSetType (eps, const_cast<char *>(EPSGD));
404  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
405 
407  {
408  ierr = EPSGDSetDoubleExpansion (eps, PETSC_TRUE);
409  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
410  }
411 #else
412  // PETSc/SLEPc version must be > 3.1.0.
413  Assert ((false),
414  ExcMessage ("Your SLEPc installation does not include a copy of the "
415  "Generalized Davidson solver. A SLEPc version > 3.1.0 is required."));
416 #endif
417  }
418 
419  /* ------------------ Jacobi Davidson -------------------- */
421  const MPI_Comm &mpi_communicator,
422  const AdditionalData &data)
423  :
424  SolverBase (cn, mpi_communicator),
425  additional_data (data)
426  {
427 #if DEAL_II_PETSC_VERSION_GTE(3,1,0)
428  const PetscErrorCode ierr = EPSSetType (eps, const_cast<char *>(EPSJD));
429  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
430 #else
431  // PETSc/SLEPc version must be > 3.1.0.
432  Assert ((false),
433  ExcMessage ("Your SLEPc installation does not include a copy of the "
434  "Jacobi-Davidson solver. A SLEPc version > 3.1.0 is required."));
435 #endif
436  }
437 
438  /* ---------------------- LAPACK ------------------------- */
440  const MPI_Comm &mpi_communicator,
441  const AdditionalData &data)
442  :
443  SolverBase (cn, mpi_communicator),
444  additional_data (data)
445  {
446  // 'Tis overwhelmingly likely that PETSc/SLEPc *always* has
447  // BLAS/LAPACK, but let's be defensive.
448 #if PETSC_HAVE_BLASLAPACK
449  const PetscErrorCode ierr = EPSSetType (eps, const_cast<char *>(EPSLAPACK));
450  AssertThrow (ierr == 0, ExcSLEPcError(ierr));
451 #else
452  Assert ((false),
453  ExcMessage ("Your PETSc/SLEPc installation was not configured with BLAS/LAPACK "
454  "but this is needed to use the LAPACK solver."));
455 #endif
456  }
457 }
458 
459 DEAL_II_NAMESPACE_CLOSE
460 
461 #endif // DEAL_II_WITH_SLEPC
462 
virtual State check(const unsigned int step, const double check_value)
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)
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
SolverGeneralizedDavidson(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
EPSConvergedReason reason
Definition: slepc_solver.h:364
const AdditionalData additional_data
Definition: slepc_solver.h:508
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())
double tolerance() const
unsigned int last_step() const
unsigned int max_steps() const
SolverLAPACK(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
static ::ExceptionBase & ExcNotImplemented()
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)