20 #include "libmesh/libmesh_common.h"
22 #if defined(LIBMESH_HAVE_SLEPC) && defined(LIBMESH_HAVE_PETSC)
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"
44 this->_eigen_solver_type =
ARNOLDI;
45 this->_eigen_problem_type =
NHEP;
65 PetscErrorCode
ierr=0;
67 ierr = LibMeshEPSDestroy(&_eps);
81 PetscErrorCode
ierr=0;
89 ierr = EPSCreate (this->comm().
get(), &_eps);
93 set_slepc_solver_type();
100 std::pair<unsigned int, unsigned int>
105 const unsigned int m_its)
107 LOG_SCOPE(
"solve_standard()",
"SlepcEigenSolver");
117 libmesh_error_msg(
"Error: input matrix to solve_standard() must be a PetscMatrix.");
120 if (this->_close_matrix_before_solve)
123 return _solve_standard_helper(matrix_A->
mat(),
nullptr, nev, ncv, tol, m_its);
127 template <
typename T>
128 std::pair<unsigned int, unsigned int>
133 const unsigned int m_its)
139 PetscErrorCode
ierr=0;
146 ierr = MatCreateShell(this->comm().
get(),
151 const_cast<void *>(static_cast<const void *>(&shell_matrix)),
153 LIBMESH_CHKERR(
ierr);
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);
160 return _solve_standard_helper(mat,
nullptr, nev, ncv, tol, m_its);
163 template <
typename T>
164 std::pair<unsigned int, unsigned int>
170 const unsigned int m_its)
181 return _solve_standard_helper(matrix->
mat(), precond->
mat(), nev, ncv, tol, m_its);
184 template <
typename T>
185 std::pair<unsigned int, unsigned int>
191 const unsigned int m_its)
193 LOG_SCOPE(
"solve_standard()",
"SlepcEigenSolver");
195 PetscErrorCode
ierr=0;
204 PetscReal error, re, im;
211 ierr = EPSSetOperators (_eps, mat, PETSC_NULL);
212 LIBMESH_CHKERR(
ierr);
215 set_slepc_problem_type();
216 set_slepc_position_of_spectrum();
219 #if SLEPC_VERSION_LESS_THAN(3,0,0)
220 ierr = EPSSetDimensions (_eps, nev, ncv);
222 ierr = EPSSetDimensions (_eps, nev, ncv, PETSC_DECIDE);
224 LIBMESH_CHKERR(
ierr);
226 ierr = EPSSetTolerances (_eps, tol, m_its);
227 LIBMESH_CHKERR(
ierr);
234 ierr = EPSSetFromOptions (_eps);
235 LIBMESH_CHKERR(
ierr);
239 ierr = EPSGetST(_eps,&st);LIBMESH_CHKERR(
ierr);
240 ierr = STPrecondSetMatForPC(st,precond);LIBMESH_CHKERR(
ierr);
245 if (this->_solver_configuration)
247 this->_solver_configuration->configure_solver();
251 ierr = EPSSolve (_eps);
252 LIBMESH_CHKERR(
ierr);
255 ierr = EPSGetIterationNumber (_eps, &its);
256 LIBMESH_CHKERR(
ierr);
259 ierr = EPSGetConverged(_eps,&nconv);
260 LIBMESH_CHKERR(
ierr);
269 ierr = PetscPrintf(this->comm().
get(),
270 " k ||Ax-kx||/|kx|\n"
271 " ----------------- -----------------\n" );
272 LIBMESH_CHKERR(
ierr);
274 for (PetscInt i=0; i<nconv; i++ )
276 ierr = EPSGetEigenpair(_eps, i, &kr, &ki, PETSC_NULL, PETSC_NULL);
277 LIBMESH_CHKERR(
ierr);
279 #if SLEPC_VERSION_LESS_THAN(3,6,0)
280 ierr = EPSComputeRelativeError(_eps, i, &error);
282 ierr = EPSComputeError(_eps, i, EPS_ERROR_RELATIVE, &error);
284 LIBMESH_CHKERR(
ierr);
286 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
287 re = PetscRealPart(kr);
288 im = PetscImaginaryPart(kr);
296 ierr = PetscPrintf(this->comm().
get(),
" %9f%+9f i %12f\n", re, im, error);
297 LIBMESH_CHKERR(
ierr);
301 ierr = PetscPrintf(this->comm().
get(),
" %12f %12f\n", re, error);
302 LIBMESH_CHKERR(
ierr);
306 ierr = PetscPrintf(this->comm().
get(),
"\n" );
307 LIBMESH_CHKERR(
ierr);
312 return std::make_pair(nconv, its);
319 template <
typename T>
320 std::pair<unsigned int, unsigned int>
326 const unsigned int m_its)
336 if (!matrix_A || !matrix_B)
337 libmesh_error_msg(
"Error: inputs to solve_generalized() must be of type PetscMatrix.");
340 if (this->_close_matrix_before_solve)
346 return _solve_generalized_helper (matrix_A->
mat(), matrix_B->
mat(),
nullptr, nev, ncv, tol, m_its);
349 template <
typename T>
350 std::pair<unsigned int, unsigned int>
356 const unsigned int m_its)
362 PetscErrorCode
ierr=0;
369 ierr = MatCreateShell(this->comm().
get(),
374 const_cast<void *>(static_cast<const void *>(&shell_matrix_A)),
376 LIBMESH_CHKERR(
ierr);
381 libmesh_error_msg(
"Error: inputs to solve_generalized() must be of type PetscMatrix.");
384 if (this->_close_matrix_before_solve)
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);
392 return _solve_generalized_helper (mat_A, matrix_B->
mat(),
nullptr, nev, ncv, tol, m_its);
395 template <
typename T>
396 std::pair<unsigned int, unsigned int>
402 const unsigned int m_its)
408 PetscErrorCode
ierr=0;
413 libmesh_error_msg(
"Error: inputs to solve_generalized() must be of type PetscMatrix.");
416 if (this->_close_matrix_before_solve)
424 ierr = MatCreateShell(this->comm().
get(),
429 const_cast<void *>(static_cast<const void *>(&shell_matrix_B)),
431 LIBMESH_CHKERR(
ierr);
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);
439 return _solve_generalized_helper (matrix_A->
mat(), mat_B,
nullptr, nev, ncv, tol, m_its);
442 template <
typename T>
443 std::pair<unsigned int, unsigned int>
449 const unsigned int m_its)
455 PetscErrorCode
ierr=0;
462 ierr = MatCreateShell(this->comm().
get(),
467 const_cast<void *>(static_cast<const void *>(&shell_matrix_A)),
469 LIBMESH_CHKERR(
ierr);
472 ierr = MatCreateShell(this->comm().
get(),
477 const_cast<void *>(static_cast<const void *>(&shell_matrix_B)),
479 LIBMESH_CHKERR(
ierr);
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);
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);
491 return _solve_generalized_helper (mat_A, mat_B,
nullptr, nev, ncv, tol, m_its);
494 template <
typename T>
495 std::pair<unsigned int, unsigned int>
502 const unsigned int m_its)
515 return _solve_generalized_helper (matrix_A->
mat(), matrix_B->
mat(), precond->
mat(), nev, ncv, tol, m_its);
518 template <
typename T>
519 std::pair<unsigned int, unsigned int>
526 const unsigned int m_its)
528 LOG_SCOPE(
"solve_generalized()",
"SlepcEigenSolver");
530 PetscErrorCode
ierr=0;
539 PetscReal error, re, im;
546 ierr = EPSSetOperators (_eps, mat_A, mat_B);
547 LIBMESH_CHKERR(
ierr);
550 set_slepc_problem_type();
551 set_slepc_position_of_spectrum();
554 #if SLEPC_VERSION_LESS_THAN(3,0,0)
555 ierr = EPSSetDimensions (_eps, nev, ncv);
557 ierr = EPSSetDimensions (_eps, nev, ncv, PETSC_DECIDE);
559 LIBMESH_CHKERR(
ierr);
563 ierr = EPSSetTolerances (_eps, tol, m_its);
564 LIBMESH_CHKERR(
ierr);
571 ierr = EPSSetFromOptions (_eps);
572 LIBMESH_CHKERR(
ierr);
576 ierr = EPSGetST(_eps,&st);LIBMESH_CHKERR(
ierr);
577 ierr = STPrecondSetMatForPC(st,precond);LIBMESH_CHKERR(
ierr);
582 if (this->_solver_configuration)
584 this->_solver_configuration->configure_solver();
588 ierr = EPSSolve (_eps);
589 LIBMESH_CHKERR(
ierr);
592 ierr = EPSGetIterationNumber (_eps, &its);
593 LIBMESH_CHKERR(
ierr);
596 ierr = EPSGetConverged(_eps,&nconv);
597 LIBMESH_CHKERR(
ierr);
606 ierr = PetscPrintf(this->comm().
get(),
607 " k ||Ax-kx||/|kx|\n"
608 " ----------------- -----------------\n" );
609 LIBMESH_CHKERR(
ierr);
611 for (PetscInt i=0; i<nconv; i++ )
613 ierr = EPSGetEigenpair(_eps, i, &kr, &ki, PETSC_NULL, PETSC_NULL);
614 LIBMESH_CHKERR(
ierr);
616 #if SLEPC_VERSION_LESS_THAN(3,6,0)
617 ierr = EPSComputeRelativeError(_eps, i, &error);
619 ierr = EPSComputeError(_eps, i, EPS_ERROR_RELATIVE, &error);
621 LIBMESH_CHKERR(
ierr);
623 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
624 re = PetscRealPart(kr);
625 im = PetscImaginaryPart(kr);
633 ierr = PetscPrintf(this->comm().
get(),
" %9f%+9f i %12f\n", re, im, error);
634 LIBMESH_CHKERR(
ierr);
638 ierr = PetscPrintf(this->comm().
get(),
" %12f %12f\n", re, error);
639 LIBMESH_CHKERR(
ierr);
643 ierr = PetscPrintf(this->comm().
get(),
"\n" );
644 LIBMESH_CHKERR(
ierr);
649 return std::make_pair(nconv, its);
662 template <
typename T>
665 PetscErrorCode
ierr = 0;
667 switch (this->_eigen_solver_type)
670 ierr = EPSSetType (_eps, EPSPOWER); LIBMESH_CHKERR(
ierr);
return;
672 ierr = EPSSetType (_eps, EPSSUBSPACE); LIBMESH_CHKERR(
ierr);
return;
674 ierr = EPSSetType (_eps, EPSLAPACK); LIBMESH_CHKERR(
ierr);
return;
676 ierr = EPSSetType (_eps, EPSARNOLDI); LIBMESH_CHKERR(
ierr);
return;
678 ierr = EPSSetType (_eps, EPSLANCZOS); LIBMESH_CHKERR(
ierr);
return;
680 ierr = EPSSetType (_eps, EPSKRYLOVSCHUR); LIBMESH_CHKERR(
ierr);
return;
685 libMesh::err <<
"ERROR: Unsupported SLEPc Eigen Solver: "
687 <<
"Continuing with SLEPc defaults" << std::endl;
694 template <
typename T>
697 PetscErrorCode
ierr = 0;
699 switch (this->_eigen_problem_type)
702 ierr = EPSSetProblemType (_eps, EPS_NHEP); LIBMESH_CHKERR(
ierr);
return;
704 ierr = EPSSetProblemType (_eps, EPS_GNHEP); LIBMESH_CHKERR(
ierr);
return;
706 ierr = EPSSetProblemType (_eps, EPS_HEP); LIBMESH_CHKERR(
ierr);
return;
708 ierr = EPSSetProblemType (_eps, EPS_GHEP); LIBMESH_CHKERR(
ierr);
return;
709 #if !SLEPC_VERSION_LESS_THAN(3,3,0)
712 ierr = EPSSetProblemType (_eps, EPS_GHIEP); LIBMESH_CHKERR(
ierr);
return;
716 libMesh::err <<
"ERROR: Unsupported SLEPc Eigen Problem: "
717 << this->_eigen_problem_type << std::endl
718 <<
"Continuing with SLEPc defaults" << std::endl;
724 template <
typename T>
727 PetscErrorCode
ierr = 0;
729 switch (this->_position_of_spectrum)
733 ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_MAGNITUDE);
734 LIBMESH_CHKERR(
ierr);
739 ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_MAGNITUDE);
740 LIBMESH_CHKERR(
ierr);
745 ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_REAL);
746 LIBMESH_CHKERR(
ierr);
751 ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_REAL);
752 LIBMESH_CHKERR(
ierr);
757 ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_IMAGINARY);
758 LIBMESH_CHKERR(
ierr);
763 ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_IMAGINARY);
764 LIBMESH_CHKERR(
ierr);
769 #if !SLEPC_VERSION_LESS_THAN(3,1,0)
772 ierr = EPSSetTarget(_eps, this->_target_val);
773 LIBMESH_CHKERR(
ierr);
774 ierr = EPSSetWhichEigenpairs (_eps, EPS_TARGET_MAGNITUDE);
775 LIBMESH_CHKERR(
ierr);
780 ierr = EPSSetTarget(_eps, this->_target_val);
781 LIBMESH_CHKERR(
ierr);
782 ierr = EPSSetWhichEigenpairs (_eps, EPS_TARGET_REAL);
783 LIBMESH_CHKERR(
ierr);
788 ierr = EPSSetTarget(_eps, this->_target_val);
789 LIBMESH_CHKERR(
ierr);
790 ierr = EPSSetWhichEigenpairs (_eps, EPS_TARGET_IMAGINARY);
791 LIBMESH_CHKERR(
ierr);
797 libmesh_error_msg(
"ERROR: Unsupported SLEPc position of spectrum: " << this->_position_of_spectrum);
806 template <
typename T>
810 PetscErrorCode
ierr=0;
818 libmesh_error_msg(
"Error getting eigenvector: input vector must be a PetscVector.");
825 ierr = EPSGetEigenpair(_eps, i, &kr, &ki, solution->
vec(), PETSC_NULL);
826 LIBMESH_CHKERR(
ierr);
828 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
829 re = PetscRealPart(kr);
830 im = PetscImaginaryPart(kr);
836 return std::make_pair(re, im);
840 template <
typename T>
843 PetscErrorCode
ierr=0;
850 ierr = EPSGetEigenvalue(_eps, i, &kr, &ki);
851 LIBMESH_CHKERR(
ierr);
853 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
854 re = PetscRealPart(kr);
855 im = PetscImaginaryPart(kr);
861 return std::make_pair(re, im);
865 template <
typename T>
868 PetscErrorCode
ierr=0;
871 #if SLEPC_VERSION_LESS_THAN(3,6,0)
872 ierr = EPSComputeRelativeError(_eps, i, &error);
874 ierr = EPSComputeError(_eps, i, EPS_ERROR_RELATIVE, &error);
876 LIBMESH_CHKERR(
ierr);
882 template <
typename T>
887 PetscErrorCode
ierr = 0;
893 if (!deflation_vector_petsc_vec)
894 libmesh_error_msg(
"Error attaching deflation space: input vector must be a PetscVector.");
897 Vec deflation_vector = deflation_vector_petsc_vec->
vec();
899 #if SLEPC_VERSION_LESS_THAN(3,1,0)
900 ierr = EPSAttachDeflationSpace(_eps, 1, &deflation_vector, PETSC_FALSE);
902 ierr = EPSSetDeflationSpace(_eps, 1, &deflation_vector);
904 LIBMESH_CHKERR(
ierr);
907 template <
typename T>
910 #if SLEPC_VERSION_LESS_THAN(3,1,0)
911 libmesh_error_msg(
"SLEPc 3.1 is required to call EigenSolver::set_initial_space()");
915 PetscErrorCode
ierr = 0;
921 if (!initial_space_petsc_vec)
922 libmesh_error_msg(
"Error attaching initial space: input vector must be a PetscVector.");
925 Vec initial_vector = initial_space_petsc_vec->
vec();
927 ierr = EPSSetInitialSpace(_eps, 1, &initial_vector);
928 LIBMESH_CHKERR(
ierr);
932 template <
typename T>
936 PetscErrorCode
ierr=0;
938 ierr = MatShellGetContext(mat,&
ctx);
940 Parallel::communicator comm;
941 PetscObjectGetComm((PetscObject)mat,&comm);
942 CHKERRABORT(comm,
ierr);
957 template <
typename T>
961 PetscErrorCode
ierr=0;
963 ierr = MatShellGetContext(mat,&
ctx);
965 Parallel::communicator comm;
966 PetscObjectGetComm((PetscObject)mat,&comm);
967 CHKERRABORT(comm,
ierr);
989 #endif // #ifdef LIBMESH_HAVE_SLEPC