libMesh
slepc_eigen_solver.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #include "libmesh/libmesh_common.h"
21 
22 #if defined(LIBMESH_HAVE_SLEPC) && defined(LIBMESH_HAVE_PETSC)
23 
24 // Local Includes
25 #include "libmesh/libmesh_logging.h"
26 #include "libmesh/petsc_matrix.h"
27 #include "libmesh/petsc_vector.h"
28 #include "libmesh/slepc_eigen_solver.h"
29 #include "libmesh/shell_matrix.h"
30 #include "libmesh/enum_to_string.h"
31 #include "libmesh/solver_configuration.h"
32 #include "libmesh/enum_eigen_solver_type.h"
33 #include "libmesh/petsc_shell_matrix.h"
34 
35 namespace libMesh
36 {
37 
38 
39 
40 template <typename T>
41 SlepcEigenSolver<T>::SlepcEigenSolver (const Parallel::Communicator & comm_in) :
42  EigenSolver<T>(comm_in)
43 {
44  this->_eigen_solver_type = ARNOLDI;
45  this->_eigen_problem_type = NHEP;
46 }
47 
48 
49 
50 template <typename T>
52 {
53  this->clear ();
54 }
55 
56 
57 
58 template <typename T>
60 {
61  if (this->initialized())
62  {
63  this->_is_initialized = false;
64 
65  PetscErrorCode ierr=0;
66 
67  ierr = LibMeshEPSDestroy(&_eps);
68  LIBMESH_CHKERR(ierr);
69 
70  // SLEPc default eigenproblem solver
71  this->_eigen_solver_type = KRYLOVSCHUR;
72  }
73 }
74 
75 
76 
77 template <typename T>
79 {
80 
81  PetscErrorCode ierr=0;
82 
83  // Initialize the data structures if not done so already.
84  if (!this->initialized())
85  {
86  this->_is_initialized = true;
87 
88  // Create the eigenproblem solver context
89  ierr = EPSCreate (this->comm().get(), &_eps);
90  LIBMESH_CHKERR(ierr);
91 
92  // Set user-specified solver
93  set_slepc_solver_type();
94  }
95 }
96 
97 
98 
99 template <typename T>
100 std::pair<unsigned int, unsigned int>
102  int nev, // number of requested eigenpairs
103  int ncv, // number of basis vectors
104  const double tol, // solver tolerance
105  const unsigned int m_its) // maximum number of iterations
106 {
107  LOG_SCOPE("solve_standard()", "SlepcEigenSolver");
108 
109  this->clear ();
110 
111  this->init ();
112 
113  // Make sure the SparseMatrix passed in is really a PetscMatrix
114  PetscMatrix<T> * matrix_A = dynamic_cast<PetscMatrix<T> *>(&matrix_A_in);
115 
116  if (!matrix_A)
117  libmesh_error_msg("Error: input matrix to solve_standard() must be a PetscMatrix.");
118 
119  // Close the matrix and vectors in case this wasn't already done.
120  if (this->_close_matrix_before_solve)
121  matrix_A->close ();
122 
123  return _solve_standard_helper(matrix_A->mat(), nullptr, nev, ncv, tol, m_its);
124 }
125 
126 
127 template <typename T>
128 std::pair<unsigned int, unsigned int>
130  int nev, // number of requested eigenpairs
131  int ncv, // number of basis vectors
132  const double tol, // solver tolerance
133  const unsigned int m_its) // maximum number of iterations
134 {
135  this->clear ();
136 
137  this->init ();
138 
139  PetscErrorCode ierr=0;
140 
141  // Prepare the matrix. Note that the const_cast is only necessary
142  // because PETSc does not accept a const void *. Inside the member
143  // function _petsc_shell_matrix() below, the pointer is casted back
144  // to a const ShellMatrix<T> *.
145  Mat mat;
146  ierr = MatCreateShell(this->comm().get(),
147  shell_matrix.m(), // Specify the number of local rows
148  shell_matrix.n(), // Specify the number of local columns
149  PETSC_DETERMINE,
150  PETSC_DETERMINE,
151  const_cast<void *>(static_cast<const void *>(&shell_matrix)),
152  &mat);
153  LIBMESH_CHKERR(ierr);
154 
155  ierr = MatShellSetOperation(mat,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult));
156  LIBMESH_CHKERR(ierr);
157  ierr = MatShellSetOperation(mat,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal));
158  LIBMESH_CHKERR(ierr);
159 
160  return _solve_standard_helper(mat, nullptr, nev, ncv, tol, m_its);
161 }
162 
163 template <typename T>
164 std::pair<unsigned int, unsigned int>
166  SparseMatrix<T> & precond_in,
167  int nev, // number of requested eigenpairs
168  int ncv, // number of basis vectors
169  const double tol, // solver tolerance
170  const unsigned int m_its) // maximum number of iterations
171 {
172  this->clear ();
173 
174  this->init ();
175 
176  // Make sure the SparseMatrix passed in is really a PetscMatrix
177  PetscMatrix<T> * precond = dynamic_cast<PetscMatrix<T> *>(&precond_in);
178 
179  PetscShellMatrix<T> * matrix = static_cast<PetscShellMatrix<T> *> (&shell_matrix);
180 
181  return _solve_standard_helper(matrix->mat(), precond->mat(), nev, ncv, tol, m_its);
182 }
183 
184 template <typename T>
185 std::pair<unsigned int, unsigned int>
187  Mat precond,
188  int nev, // number of requested eigenpairs
189  int ncv, // number of basis vectors
190  const double tol, // solver tolerance
191  const unsigned int m_its) // maximum number of iterations
192 {
193  LOG_SCOPE("solve_standard()", "SlepcEigenSolver");
194 
195  PetscErrorCode ierr=0;
196 
197  // converged eigen pairs and number of iterations
198  PetscInt nconv=0;
199  PetscInt its=0;
200  ST st=nullptr;
201 
202 #ifdef DEBUG
203  // The relative error.
204  PetscReal error, re, im;
205 
206  // Pointer to vectors of the real parts, imaginary parts.
207  PetscScalar kr, ki;
208 #endif
209 
210  // Set operators.
211  ierr = EPSSetOperators (_eps, mat, PETSC_NULL);
212  LIBMESH_CHKERR(ierr);
213 
214  //set the problem type and the position of the spectrum
215  set_slepc_problem_type();
216  set_slepc_position_of_spectrum();
217 
218  // Set eigenvalues to be computed.
219 #if SLEPC_VERSION_LESS_THAN(3,0,0)
220  ierr = EPSSetDimensions (_eps, nev, ncv);
221 #else
222  ierr = EPSSetDimensions (_eps, nev, ncv, PETSC_DECIDE);
223 #endif
224  LIBMESH_CHKERR(ierr);
225  // Set the tolerance and maximum iterations.
226  ierr = EPSSetTolerances (_eps, tol, m_its);
227  LIBMESH_CHKERR(ierr);
228 
229  // Set runtime options, e.g.,
230  // -eps_type <type>, -eps_nev <nev>, -eps_ncv <ncv>
231  // Similar to PETSc, these options will override those specified
232  // above as long as EPSSetFromOptions() is called _after_ any
233  // other customization routines.
234  ierr = EPSSetFromOptions (_eps);
235  LIBMESH_CHKERR(ierr);
236 
237  // Set a preconditioning matrix to ST
238  if (precond) {
239  ierr = EPSGetST(_eps,&st);LIBMESH_CHKERR(ierr);
240  ierr = STPrecondSetMatForPC(st,precond);LIBMESH_CHKERR(ierr);
241  }
242 
243  // If the SolverConfiguration object is provided, use it to override
244  // solver options.
245  if (this->_solver_configuration)
246  {
247  this->_solver_configuration->configure_solver();
248  }
249 
250  // Solve the eigenproblem.
251  ierr = EPSSolve (_eps);
252  LIBMESH_CHKERR(ierr);
253 
254  // Get the number of iterations.
255  ierr = EPSGetIterationNumber (_eps, &its);
256  LIBMESH_CHKERR(ierr);
257 
258  // Get number of converged eigenpairs.
259  ierr = EPSGetConverged(_eps,&nconv);
260  LIBMESH_CHKERR(ierr);
261 
262 
263 #ifdef DEBUG
264  // ierr = PetscPrintf(this->comm().get(),
265  // "\n Number of iterations: %d\n"
266  // " Number of converged eigenpairs: %d\n\n", its, nconv);
267 
268  // Display eigenvalues and relative errors.
269  ierr = PetscPrintf(this->comm().get(),
270  " k ||Ax-kx||/|kx|\n"
271  " ----------------- -----------------\n" );
272  LIBMESH_CHKERR(ierr);
273 
274  for (PetscInt i=0; i<nconv; i++ )
275  {
276  ierr = EPSGetEigenpair(_eps, i, &kr, &ki, PETSC_NULL, PETSC_NULL);
277  LIBMESH_CHKERR(ierr);
278 
279 #if SLEPC_VERSION_LESS_THAN(3,6,0)
280  ierr = EPSComputeRelativeError(_eps, i, &error);
281 #else
282  ierr = EPSComputeError(_eps, i, EPS_ERROR_RELATIVE, &error);
283 #endif
284  LIBMESH_CHKERR(ierr);
285 
286 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
287  re = PetscRealPart(kr);
288  im = PetscImaginaryPart(kr);
289 #else
290  re = kr;
291  im = ki;
292 #endif
293 
294  if (im != .0)
295  {
296  ierr = PetscPrintf(this->comm().get()," %9f%+9f i %12f\n", re, im, error);
297  LIBMESH_CHKERR(ierr);
298  }
299  else
300  {
301  ierr = PetscPrintf(this->comm().get()," %12f %12f\n", re, error);
302  LIBMESH_CHKERR(ierr);
303  }
304  }
305 
306  ierr = PetscPrintf(this->comm().get(),"\n" );
307  LIBMESH_CHKERR(ierr);
308 #endif // DEBUG
309 
310  // return the number of converged eigenpairs
311  // and the number of iterations
312  return std::make_pair(nconv, its);
313 }
314 
315 
316 
317 
318 
319 template <typename T>
320 std::pair<unsigned int, unsigned int>
322  SparseMatrix<T> & matrix_B_in,
323  int nev, // number of requested eigenpairs
324  int ncv, // number of basis vectors
325  const double tol, // solver tolerance
326  const unsigned int m_its) // maximum number of iterations
327 {
328  this->clear ();
329 
330  this->init ();
331 
332  // Make sure the data passed in are really of Petsc types
333  PetscMatrix<T> * matrix_A = dynamic_cast<PetscMatrix<T> *>(&matrix_A_in);
334  PetscMatrix<T> * matrix_B = dynamic_cast<PetscMatrix<T> *>(&matrix_B_in);
335 
336  if (!matrix_A || !matrix_B)
337  libmesh_error_msg("Error: inputs to solve_generalized() must be of type PetscMatrix.");
338 
339  // Close the matrix and vectors in case this wasn't already done.
340  if (this->_close_matrix_before_solve)
341  {
342  matrix_A->close ();
343  matrix_B->close ();
344  }
345 
346  return _solve_generalized_helper (matrix_A->mat(), matrix_B->mat(), nullptr, nev, ncv, tol, m_its);
347 }
348 
349 template <typename T>
350 std::pair<unsigned int, unsigned int>
352  SparseMatrix<T> & matrix_B_in,
353  int nev, // number of requested eigenpairs
354  int ncv, // number of basis vectors
355  const double tol, // solver tolerance
356  const unsigned int m_its) // maximum number of iterations
357 {
358  this->clear();
359 
360  this->init ();
361 
362  PetscErrorCode ierr=0;
363 
364  // Prepare the matrix. Note that the const_cast is only necessary
365  // because PETSc does not accept a const void *. Inside the member
366  // function _petsc_shell_matrix() below, the pointer is casted back
367  // to a const ShellMatrix<T> *.
368  Mat mat_A;
369  ierr = MatCreateShell(this->comm().get(),
370  shell_matrix_A.m(), // Specify the number of local rows
371  shell_matrix_A.n(), // Specify the number of local columns
372  PETSC_DETERMINE,
373  PETSC_DETERMINE,
374  const_cast<void *>(static_cast<const void *>(&shell_matrix_A)),
375  &mat_A);
376  LIBMESH_CHKERR(ierr);
377 
378  PetscMatrix<T> * matrix_B = dynamic_cast<PetscMatrix<T> *>(&matrix_B_in);
379 
380  if (!matrix_B)
381  libmesh_error_msg("Error: inputs to solve_generalized() must be of type PetscMatrix.");
382 
383  // Close the matrix and vectors in case this wasn't already done.
384  if (this->_close_matrix_before_solve)
385  matrix_B->close ();
386 
387  ierr = MatShellSetOperation(mat_A,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult));
388  LIBMESH_CHKERR(ierr);
389  ierr = MatShellSetOperation(mat_A,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal));
390  LIBMESH_CHKERR(ierr);
391 
392  return _solve_generalized_helper (mat_A, matrix_B->mat(), nullptr, nev, ncv, tol, m_its);
393 }
394 
395 template <typename T>
396 std::pair<unsigned int, unsigned int>
398  ShellMatrix<T> & shell_matrix_B,
399  int nev, // number of requested eigenpairs
400  int ncv, // number of basis vectors
401  const double tol, // solver tolerance
402  const unsigned int m_its) // maximum number of iterations
403 {
404  this->clear();
405 
406  this->init ();
407 
408  PetscErrorCode ierr=0;
409 
410  PetscMatrix<T> * matrix_A = dynamic_cast<PetscMatrix<T> *>(&matrix_A_in);
411 
412  if (!matrix_A)
413  libmesh_error_msg("Error: inputs to solve_generalized() must be of type PetscMatrix.");
414 
415  // Close the matrix and vectors in case this wasn't already done.
416  if (this->_close_matrix_before_solve)
417  matrix_A->close ();
418 
419  // Prepare the matrix. Note that the const_cast is only necessary
420  // because PETSc does not accept a const void *. Inside the member
421  // function _petsc_shell_matrix() below, the pointer is casted back
422  // to a const ShellMatrix<T> *.
423  Mat mat_B;
424  ierr = MatCreateShell(this->comm().get(),
425  shell_matrix_B.m(), // Specify the number of local rows
426  shell_matrix_B.n(), // Specify the number of local columns
427  PETSC_DETERMINE,
428  PETSC_DETERMINE,
429  const_cast<void *>(static_cast<const void *>(&shell_matrix_B)),
430  &mat_B);
431  LIBMESH_CHKERR(ierr);
432 
433 
434  ierr = MatShellSetOperation(mat_B,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult));
435  LIBMESH_CHKERR(ierr);
436  ierr = MatShellSetOperation(mat_B,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal));
437  LIBMESH_CHKERR(ierr);
438 
439  return _solve_generalized_helper (matrix_A->mat(), mat_B, nullptr, nev, ncv, tol, m_its);
440 }
441 
442 template <typename T>
443 std::pair<unsigned int, unsigned int>
445  ShellMatrix<T> & shell_matrix_B,
446  int nev, // number of requested eigenpairs
447  int ncv, // number of basis vectors
448  const double tol, // solver tolerance
449  const unsigned int m_its) // maximum number of iterations
450 {
451  this->clear();
452 
453  this->init ();
454 
455  PetscErrorCode ierr=0;
456 
457  // Prepare the matrices. Note that the const_casts are only
458  // necessary because PETSc does not accept a const void *. Inside
459  // the member function _petsc_shell_matrix() below, the pointer is
460  // casted back to a const ShellMatrix<T> *.
461  Mat mat_A;
462  ierr = MatCreateShell(this->comm().get(),
463  shell_matrix_A.m(), // Specify the number of local rows
464  shell_matrix_A.n(), // Specify the number of local columns
465  PETSC_DETERMINE,
466  PETSC_DETERMINE,
467  const_cast<void *>(static_cast<const void *>(&shell_matrix_A)),
468  &mat_A);
469  LIBMESH_CHKERR(ierr);
470 
471  Mat mat_B;
472  ierr = MatCreateShell(this->comm().get(),
473  shell_matrix_B.m(), // Specify the number of local rows
474  shell_matrix_B.n(), // Specify the number of local columns
475  PETSC_DETERMINE,
476  PETSC_DETERMINE,
477  const_cast<void *>(static_cast<const void *>(&shell_matrix_B)),
478  &mat_B);
479  LIBMESH_CHKERR(ierr);
480 
481  ierr = MatShellSetOperation(mat_A,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult));
482  LIBMESH_CHKERR(ierr);
483  ierr = MatShellSetOperation(mat_A,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal));
484  LIBMESH_CHKERR(ierr);
485 
486  ierr = MatShellSetOperation(mat_B,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult));
487  LIBMESH_CHKERR(ierr);
488  ierr = MatShellSetOperation(mat_B,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal));
489  LIBMESH_CHKERR(ierr);
490 
491  return _solve_generalized_helper (mat_A, mat_B, nullptr, nev, ncv, tol, m_its);
492 }
493 
494 template <typename T>
495 std::pair<unsigned int, unsigned int>
497  ShellMatrix<T> & shell_matrix_B,
498  SparseMatrix<T> & precond_in,
499  int nev, // number of requested eigenpairs
500  int ncv, // number of basis vectors
501  const double tol, // solver tolerance
502  const unsigned int m_its) // maximum number of iterations
503 {
504  this->clear();
505 
506  this->init ();
507 
508  // Make sure the SparseMatrix passed in is really a PetscMatrix
509  PetscMatrix<T> * precond = static_cast<PetscMatrix<T> *>(&precond_in);
510 
511  PetscShellMatrix<T> * matrix_A = static_cast<PetscShellMatrix<T> *> (&shell_matrix_A);
512 
513  PetscShellMatrix<T> * matrix_B = static_cast<PetscShellMatrix<T> *> (&shell_matrix_B);
514 
515  return _solve_generalized_helper (matrix_A->mat(), matrix_B->mat(), precond->mat(), nev, ncv, tol, m_its);
516 }
517 
518 template <typename T>
519 std::pair<unsigned int, unsigned int>
521  Mat mat_B,
522  Mat precond,
523  int nev, // number of requested eigenpairs
524  int ncv, // number of basis vectors
525  const double tol, // solver tolerance
526  const unsigned int m_its) // maximum number of iterations
527 {
528  LOG_SCOPE("solve_generalized()", "SlepcEigenSolver");
529 
530  PetscErrorCode ierr=0;
531 
532  // converged eigen pairs and number of iterations
533  PetscInt nconv=0;
534  PetscInt its=0;
535  ST st;
536 
537 #ifdef DEBUG
538  // The relative error.
539  PetscReal error, re, im;
540 
541  // Pointer to vectors of the real parts, imaginary parts.
542  PetscScalar kr, ki;
543 #endif
544 
545  // Set operators.
546  ierr = EPSSetOperators (_eps, mat_A, mat_B);
547  LIBMESH_CHKERR(ierr);
548 
549  //set the problem type and the position of the spectrum
550  set_slepc_problem_type();
551  set_slepc_position_of_spectrum();
552 
553  // Set eigenvalues to be computed.
554 #if SLEPC_VERSION_LESS_THAN(3,0,0)
555  ierr = EPSSetDimensions (_eps, nev, ncv);
556 #else
557  ierr = EPSSetDimensions (_eps, nev, ncv, PETSC_DECIDE);
558 #endif
559  LIBMESH_CHKERR(ierr);
560 
561 
562  // Set the tolerance and maximum iterations.
563  ierr = EPSSetTolerances (_eps, tol, m_its);
564  LIBMESH_CHKERR(ierr);
565 
566  // Set runtime options, e.g.,
567  // -eps_type <type>, -eps_nev <nev>, -eps_ncv <ncv>
568  // Similar to PETSc, these options will override those specified
569  // above as long as EPSSetFromOptions() is called _after_ any
570  // other customization routines.
571  ierr = EPSSetFromOptions (_eps);
572  LIBMESH_CHKERR(ierr);
573 
574  // Set a preconditioning matrix to ST
575  if (precond) {
576  ierr = EPSGetST(_eps,&st);LIBMESH_CHKERR(ierr);
577  ierr = STPrecondSetMatForPC(st,precond);LIBMESH_CHKERR(ierr);
578  }
579 
580  // If the SolverConfiguration object is provided, use it to override
581  // solver options.
582  if (this->_solver_configuration)
583  {
584  this->_solver_configuration->configure_solver();
585  }
586 
587  // Solve the eigenproblem.
588  ierr = EPSSolve (_eps);
589  LIBMESH_CHKERR(ierr);
590 
591  // Get the number of iterations.
592  ierr = EPSGetIterationNumber (_eps, &its);
593  LIBMESH_CHKERR(ierr);
594 
595  // Get number of converged eigenpairs.
596  ierr = EPSGetConverged(_eps,&nconv);
597  LIBMESH_CHKERR(ierr);
598 
599 
600 #ifdef DEBUG
601  // ierr = PetscPrintf(this->comm().get(),
602  // "\n Number of iterations: %d\n"
603  // " Number of converged eigenpairs: %d\n\n", its, nconv);
604 
605  // Display eigenvalues and relative errors.
606  ierr = PetscPrintf(this->comm().get(),
607  " k ||Ax-kx||/|kx|\n"
608  " ----------------- -----------------\n" );
609  LIBMESH_CHKERR(ierr);
610 
611  for (PetscInt i=0; i<nconv; i++ )
612  {
613  ierr = EPSGetEigenpair(_eps, i, &kr, &ki, PETSC_NULL, PETSC_NULL);
614  LIBMESH_CHKERR(ierr);
615 
616 #if SLEPC_VERSION_LESS_THAN(3,6,0)
617  ierr = EPSComputeRelativeError(_eps, i, &error);
618 #else
619  ierr = EPSComputeError(_eps, i, EPS_ERROR_RELATIVE, &error);
620 #endif
621  LIBMESH_CHKERR(ierr);
622 
623 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
624  re = PetscRealPart(kr);
625  im = PetscImaginaryPart(kr);
626 #else
627  re = kr;
628  im = ki;
629 #endif
630 
631  if (im != .0)
632  {
633  ierr = PetscPrintf(this->comm().get()," %9f%+9f i %12f\n", re, im, error);
634  LIBMESH_CHKERR(ierr);
635  }
636  else
637  {
638  ierr = PetscPrintf(this->comm().get()," %12f %12f\n", re, error);
639  LIBMESH_CHKERR(ierr);
640  }
641  }
642 
643  ierr = PetscPrintf(this->comm().get(),"\n" );
644  LIBMESH_CHKERR(ierr);
645 #endif // DEBUG
646 
647  // return the number of converged eigenpairs
648  // and the number of iterations
649  return std::make_pair(nconv, its);
650 }
651 
652 
653 
654 
655 
656 
657 
658 
659 
660 
661 
662 template <typename T>
664 {
665  PetscErrorCode ierr = 0;
666 
667  switch (this->_eigen_solver_type)
668  {
669  case POWER:
670  ierr = EPSSetType (_eps, EPSPOWER); LIBMESH_CHKERR(ierr); return;
671  case SUBSPACE:
672  ierr = EPSSetType (_eps, EPSSUBSPACE); LIBMESH_CHKERR(ierr); return;
673  case LAPACK:
674  ierr = EPSSetType (_eps, EPSLAPACK); LIBMESH_CHKERR(ierr); return;
675  case ARNOLDI:
676  ierr = EPSSetType (_eps, EPSARNOLDI); LIBMESH_CHKERR(ierr); return;
677  case LANCZOS:
678  ierr = EPSSetType (_eps, EPSLANCZOS); LIBMESH_CHKERR(ierr); return;
679  case KRYLOVSCHUR:
680  ierr = EPSSetType (_eps, EPSKRYLOVSCHUR); LIBMESH_CHKERR(ierr); return;
681  // case ARPACK:
682  // ierr = EPSSetType (_eps, (char *) EPSARPACK); LIBMESH_CHKERR(ierr); return;
683 
684  default:
685  libMesh::err << "ERROR: Unsupported SLEPc Eigen Solver: "
686  << Utility::enum_to_string(this->_eigen_solver_type) << std::endl
687  << "Continuing with SLEPc defaults" << std::endl;
688  }
689 }
690 
691 
692 
693 
694 template <typename T>
696 {
697  PetscErrorCode ierr = 0;
698 
699  switch (this->_eigen_problem_type)
700  {
701  case NHEP:
702  ierr = EPSSetProblemType (_eps, EPS_NHEP); LIBMESH_CHKERR(ierr); return;
703  case GNHEP:
704  ierr = EPSSetProblemType (_eps, EPS_GNHEP); LIBMESH_CHKERR(ierr); return;
705  case HEP:
706  ierr = EPSSetProblemType (_eps, EPS_HEP); LIBMESH_CHKERR(ierr); return;
707  case GHEP:
708  ierr = EPSSetProblemType (_eps, EPS_GHEP); LIBMESH_CHKERR(ierr); return;
709 #if !SLEPC_VERSION_LESS_THAN(3,3,0)
710  // EPS_GHIEP added in 3.3.0
711  case GHIEP:
712  ierr = EPSSetProblemType (_eps, EPS_GHIEP); LIBMESH_CHKERR(ierr); return;
713 #endif
714 
715  default:
716  libMesh::err << "ERROR: Unsupported SLEPc Eigen Problem: "
717  << this->_eigen_problem_type << std::endl
718  << "Continuing with SLEPc defaults" << std::endl;
719  }
720 }
721 
722 
723 
724 template <typename T>
726 {
727  PetscErrorCode ierr = 0;
728 
729  switch (this->_position_of_spectrum)
730  {
731  case LARGEST_MAGNITUDE:
732  {
733  ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_MAGNITUDE);
734  LIBMESH_CHKERR(ierr);
735  return;
736  }
737  case SMALLEST_MAGNITUDE:
738  {
739  ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_MAGNITUDE);
740  LIBMESH_CHKERR(ierr);
741  return;
742  }
743  case LARGEST_REAL:
744  {
745  ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_REAL);
746  LIBMESH_CHKERR(ierr);
747  return;
748  }
749  case SMALLEST_REAL:
750  {
751  ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_REAL);
752  LIBMESH_CHKERR(ierr);
753  return;
754  }
755  case LARGEST_IMAGINARY:
756  {
757  ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_IMAGINARY);
758  LIBMESH_CHKERR(ierr);
759  return;
760  }
761  case SMALLEST_IMAGINARY:
762  {
763  ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_IMAGINARY);
764  LIBMESH_CHKERR(ierr);
765  return;
766  }
767 
768  // The EPS_TARGET_XXX enums were added in SLEPc 3.1
769 #if !SLEPC_VERSION_LESS_THAN(3,1,0)
770  case TARGET_MAGNITUDE:
771  {
772  ierr = EPSSetTarget(_eps, this->_target_val);
773  LIBMESH_CHKERR(ierr);
774  ierr = EPSSetWhichEigenpairs (_eps, EPS_TARGET_MAGNITUDE);
775  LIBMESH_CHKERR(ierr);
776  return;
777  }
778  case TARGET_REAL:
779  {
780  ierr = EPSSetTarget(_eps, this->_target_val);
781  LIBMESH_CHKERR(ierr);
782  ierr = EPSSetWhichEigenpairs (_eps, EPS_TARGET_REAL);
783  LIBMESH_CHKERR(ierr);
784  return;
785  }
786  case TARGET_IMAGINARY:
787  {
788  ierr = EPSSetTarget(_eps, this->_target_val);
789  LIBMESH_CHKERR(ierr);
790  ierr = EPSSetWhichEigenpairs (_eps, EPS_TARGET_IMAGINARY);
791  LIBMESH_CHKERR(ierr);
792  return;
793  }
794 #endif
795 
796  default:
797  libmesh_error_msg("ERROR: Unsupported SLEPc position of spectrum: " << this->_position_of_spectrum);
798  }
799 }
800 
801 
802 
803 
804 
805 
806 template <typename T>
808  NumericVector<T> & solution_in)
809 {
810  PetscErrorCode ierr=0;
811 
812  PetscReal re, im;
813 
814  // Make sure the NumericVector passed in is really a PetscVector
815  PetscVector<T> * solution = dynamic_cast<PetscVector<T> *>(&solution_in);
816 
817  if (!solution)
818  libmesh_error_msg("Error getting eigenvector: input vector must be a PetscVector.");
819 
820  // real and imaginary part of the ith eigenvalue.
821  PetscScalar kr, ki;
822 
823  solution->close();
824 
825  ierr = EPSGetEigenpair(_eps, i, &kr, &ki, solution->vec(), PETSC_NULL);
826  LIBMESH_CHKERR(ierr);
827 
828 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
829  re = PetscRealPart(kr);
830  im = PetscImaginaryPart(kr);
831 #else
832  re = kr;
833  im = ki;
834 #endif
835 
836  return std::make_pair(re, im);
837 }
838 
839 
840 template <typename T>
842 {
843  PetscErrorCode ierr=0;
844 
845  PetscReal re, im;
846 
847  // real and imaginary part of the ith eigenvalue.
848  PetscScalar kr, ki;
849 
850  ierr = EPSGetEigenvalue(_eps, i, &kr, &ki);
851  LIBMESH_CHKERR(ierr);
852 
853 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
854  re = PetscRealPart(kr);
855  im = PetscImaginaryPart(kr);
856 #else
857  re = kr;
858  im = ki;
859 #endif
860 
861  return std::make_pair(re, im);
862 }
863 
864 
865 template <typename T>
867 {
868  PetscErrorCode ierr=0;
869  PetscReal error;
870 
871 #if SLEPC_VERSION_LESS_THAN(3,6,0)
872  ierr = EPSComputeRelativeError(_eps, i, &error);
873 #else
874  ierr = EPSComputeError(_eps, i, EPS_ERROR_RELATIVE, &error);
875 #endif
876  LIBMESH_CHKERR(ierr);
877 
878  return error;
879 }
880 
881 
882 template <typename T>
884 {
885  this->init();
886 
887  PetscErrorCode ierr = 0;
888 
889  // Make sure the input vector is actually a PetscVector
890  PetscVector<T> * deflation_vector_petsc_vec =
891  dynamic_cast<PetscVector<T> *>(&deflation_vector_in);
892 
893  if (!deflation_vector_petsc_vec)
894  libmesh_error_msg("Error attaching deflation space: input vector must be a PetscVector.");
895 
896  // Get a handle for the underlying Vec.
897  Vec deflation_vector = deflation_vector_petsc_vec->vec();
898 
899 #if SLEPC_VERSION_LESS_THAN(3,1,0)
900  ierr = EPSAttachDeflationSpace(_eps, 1, &deflation_vector, PETSC_FALSE);
901 #else
902  ierr = EPSSetDeflationSpace(_eps, 1, &deflation_vector);
903 #endif
904  LIBMESH_CHKERR(ierr);
905 }
906 
907 template <typename T>
909 {
910 #if SLEPC_VERSION_LESS_THAN(3,1,0)
911  libmesh_error_msg("SLEPc 3.1 is required to call EigenSolver::set_initial_space()");
912 #else
913  this->init();
914 
915  PetscErrorCode ierr = 0;
916 
917  // Make sure the input vector is actually a PetscVector
918  PetscVector<T> * initial_space_petsc_vec =
919  dynamic_cast<PetscVector<T> *>(&initial_space_in);
920 
921  if (!initial_space_petsc_vec)
922  libmesh_error_msg("Error attaching initial space: input vector must be a PetscVector.");
923 
924  // Get a handle for the underlying Vec.
925  Vec initial_vector = initial_space_petsc_vec->vec();
926 
927  ierr = EPSSetInitialSpace(_eps, 1, &initial_vector);
928  LIBMESH_CHKERR(ierr);
929 #endif
930 }
931 
932 template <typename T>
933 PetscErrorCode SlepcEigenSolver<T>::_petsc_shell_matrix_mult(Mat mat, Vec arg, Vec dest)
934 {
935  // Get the matrix context.
936  PetscErrorCode ierr=0;
937  void * ctx;
938  ierr = MatShellGetContext(mat,&ctx);
939 
940  Parallel::communicator comm;
941  PetscObjectGetComm((PetscObject)mat,&comm);
942  CHKERRABORT(comm,ierr);
943 
944  // Get user shell matrix object.
945  const ShellMatrix<T> & shell_matrix = *static_cast<const ShellMatrix<T> *>(ctx);
946 
947  // Make \p NumericVector instances around the vectors.
948  PetscVector<T> arg_global(arg, shell_matrix.comm());
949  PetscVector<T> dest_global(dest, shell_matrix.comm());
950 
951  // Call the user function.
952  shell_matrix.vector_mult(dest_global,arg_global);
953 
954  return ierr;
955 }
956 
957 template <typename T>
959 {
960  // Get the matrix context.
961  PetscErrorCode ierr=0;
962  void * ctx;
963  ierr = MatShellGetContext(mat,&ctx);
964 
965  Parallel::communicator comm;
966  PetscObjectGetComm((PetscObject)mat,&comm);
967  CHKERRABORT(comm,ierr);
968 
969  // Get user shell matrix object.
970  const ShellMatrix<T> & shell_matrix = *static_cast<const ShellMatrix<T> *>(ctx);
971 
972  // Make \p NumericVector instances around the vector.
973  PetscVector<T> dest_global(dest, shell_matrix.comm());
974 
975  // Call the user function.
976  shell_matrix.get_diagonal(dest_global);
977 
978  return ierr;
979 }
980 
981 //------------------------------------------------------------------
982 // Explicit instantiations
983 template class SlepcEigenSolver<Number>;
984 
985 } // namespace libMesh
986 
987 
988 
989 #endif // #ifdef LIBMESH_HAVE_SLEPC
libMesh::ShellMatrix
Generic shell matrix, i.e.
Definition: eigen_preconditioner.h:36
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::TARGET_IMAGINARY
Definition: enum_eigen_solver_type.h:83
libMesh::GHIEP
Definition: enum_eigen_solver_type.h:59
libMesh::HEP
Definition: enum_eigen_solver_type.h:56
libMesh::LAPACK
Definition: enum_eigen_solver_type.h:35
libMesh::POWER
Definition: enum_eigen_solver_type.h:34
libMesh::ShellMatrix::get_diagonal
virtual void get_diagonal(NumericVector< T > &dest) const =0
Copies the diagonal part of the matrix into dest.
libMesh::TARGET_REAL
Definition: enum_eigen_solver_type.h:80
libMesh::SMALLEST_REAL
Definition: enum_eigen_solver_type.h:79
libMesh::LARGEST_REAL
Definition: enum_eigen_solver_type.h:78
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::PetscMatrix::mat
Mat mat()
Definition: petsc_matrix.h:281
libMesh::ShellMatrix::n
virtual numeric_index_type n() const =0
libMesh::ParallelObject::comm
const Parallel::Communicator & comm() const
Definition: parallel_object.h:94
libMesh::SMALLEST_IMAGINARY
Definition: enum_eigen_solver_type.h:82
libMesh::SparseMatrix
Generic sparse matrix.
Definition: vector_fe_ex5.C:45
libMesh::NumericVector
Provides a uniform interface to vector storage schemes for different linear algebra libraries.
Definition: vector_fe_ex5.C:43
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::PetscShellMatrix
This class allows to use a PETSc shell matrix.
Definition: petsc_shell_matrix.h:57
libMesh::PetscVector::close
virtual void close() override
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
Definition: petsc_vector.h:821
libMesh::ShellMatrix::vector_mult
virtual void vector_mult(NumericVector< T > &dest, const NumericVector< T > &arg) const =0
Multiplies the matrix with arg and stores the result in dest.
libMesh::LARGEST_IMAGINARY
Definition: enum_eigen_solver_type.h:81
libMesh::SlepcEigenSolver
This class provides an interface to the SLEPc eigenvalue solver library from http://slepc....
Definition: slepc_eigen_solver.h:50
libMesh::PetscVector
This class provides a nice interface to PETSc's Vec object.
Definition: petsc_vector.h:72
libMesh::Utility::enum_to_string
std::string enum_to_string(const T e)
libMesh::KRYLOVSCHUR
Definition: enum_eigen_solver_type.h:39
libMesh::initialized
bool initialized()
Checks that library initialization has been done.
Definition: libmesh.C:265
libMesh::ARNOLDI
Definition: enum_eigen_solver_type.h:37
libMesh::SMALLEST_MAGNITUDE
Definition: enum_eigen_solver_type.h:76
libMesh::LARGEST_MAGNITUDE
Definition: enum_eigen_solver_type.h:75
libMesh::ReferenceElem::get
const Elem & get(const ElemType type_in)
Definition: reference_elem.C:237
libMesh::NHEP
Definition: enum_eigen_solver_type.h:55
libMesh::libMeshPrivateData::_is_initialized
bool _is_initialized
Flag that tells if init() has been called.
Definition: libmesh.C:246
libMesh::SUBSPACE
Definition: enum_eigen_solver_type.h:36
libMesh::PetscMatrix
This class provides a nice interface to the PETSc C-based data structures for parallel,...
Definition: petsc_matrix.h:87
libMesh::GHEP
Definition: enum_eigen_solver_type.h:58
libMesh::EigenSolver
This class provides an interface to solvers for eigenvalue problems.
Definition: eigen_solver.h:67
libMesh::ctx
void * ctx
Definition: petsc_dm_wrapper.C:71
libMesh::PetscMatrix::close
virtual void close() override
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
Definition: petsc_matrix.C:989
libMesh::TARGET_MAGNITUDE
Definition: enum_eigen_solver_type.h:77
libMesh::PetscVector::vec
Vec vec()
Definition: petsc_vector.h:335
libMesh::err
OStreamProxy err
libMesh::GNHEP
Definition: enum_eigen_solver_type.h:57
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
ierr
PetscErrorCode ierr
Definition: petscdmlibmeshimpl.C:1028
libMesh::SlepcEigenSolver::SlepcEigenSolver
SlepcEigenSolver(const Parallel::Communicator &comm_in)
Constructor.
Definition: slepc_eigen_solver.C:41
libMesh::ShellMatrix::m
virtual numeric_index_type m() const =0
libMesh::PetscShellMatrix::mat
Mat mat()
Definition: petsc_shell_matrix.h:94
libMesh::LANCZOS
Definition: enum_eigen_solver_type.h:38