17 #ifndef dealii__slepc_solver_h 18 #define dealii__slepc_solver_h 20 #include <deal.II/base/config.h> 22 #ifdef DEAL_II_WITH_SLEPC 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> 30 # include <petscconf.h> 31 # include <petscksp.h> 32 # include <slepceps.h> 34 DEAL_II_NAMESPACE_OPEN
173 template <
typename OutputVector>
176 std::vector<PetscScalar> &eigenvalues,
177 std::vector<OutputVector> &eigenvectors,
178 const unsigned int n_eigenpairs = 1);
185 template <
typename OutputVector>
189 std::vector<PetscScalar> &eigenvalues,
190 std::vector<OutputVector> &eigenvectors,
191 const unsigned int n_eigenpairs = 1);
198 template <
typename OutputVector>
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);
221 template <
typename Vector>
224 (
const std::vector<Vector> &initial_space);
277 <<
" An error with error number " << arg1
278 <<
" occurred while calling a SLEPc function");
285 <<
" The number of converged eigenvectors is " << arg1
286 <<
" but " << arg2 <<
" were requested. ");
313 solve (
const unsigned int n_eigenpairs,
314 unsigned int *n_converged);
323 PetscScalar &eigenvalues,
333 double &real_eigenvalues,
334 double &imag_eigenvalues,
376 PetscScalar real_eigenvalue,
377 PetscScalar imag_eigenvalue,
378 PetscReal residual_norm,
379 PetscReal *estimated_error,
380 void *solver_control);
491 AdditionalData(
const EPSLanczosReorthogType r = EPS_LANCZOS_REORTHOG_FULL);
670 template <
typename OutputVector>
673 std::vector<PetscScalar> &eigenvalues,
674 std::vector<OutputVector> &eigenvectors,
675 const unsigned int n_eigenpairs)
678 AssertThrow ((n_eigenpairs > 0) && (n_eigenpairs <= A.
m ()),
679 ExcSLEPcWrappersUsageError());
685 unsigned int n_converged = 0;
686 solve (n_eigenpairs, &n_converged);
688 if (n_converged > n_eigenpairs)
689 n_converged = n_eigenpairs;
691 ExcSLEPcEigenvectorConvergenceMismatchError(n_converged, n_eigenpairs));
693 AssertThrow (eigenvectors.size() != 0, ExcSLEPcWrappersUsageError());
694 eigenvectors.resize (n_converged, eigenvectors.front());
695 eigenvalues.resize (n_converged);
697 for (
unsigned int index=0; index<n_converged; ++index)
698 get_eigenpair (index, eigenvalues[index], eigenvectors[index]);
701 template <
typename OutputVector>
705 std::vector<PetscScalar> &eigenvalues,
706 std::vector<OutputVector> &eigenvectors,
707 const unsigned int n_eigenpairs)
714 AssertThrow ((n_eigenpairs>0) && (n_eigenpairs<=A.
m ()),
715 ExcSLEPcWrappersUsageError());
721 unsigned int n_converged = 0;
722 solve (n_eigenpairs, &n_converged);
724 if (n_converged>=n_eigenpairs)
725 n_converged = n_eigenpairs;
728 ExcSLEPcEigenvectorConvergenceMismatchError(n_converged, n_eigenpairs));
729 AssertThrow (eigenvectors.size() != 0, ExcSLEPcWrappersUsageError());
731 eigenvectors.resize (n_converged, eigenvectors.front());
732 eigenvalues.resize (n_converged);
734 for (
unsigned int index=0; index<n_converged; ++index)
735 get_eigenpair (index, eigenvalues[index], eigenvectors[index]);
738 template <
typename OutputVector>
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)
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()));
759 AssertThrow ((n_eigenpairs>0) && (n_eigenpairs<=A.
m ()),
760 ExcSLEPcWrappersUsageError());
766 unsigned int n_converged = 0;
767 solve (n_eigenpairs, &n_converged);
769 if (n_converged>=n_eigenpairs)
770 n_converged = n_eigenpairs;
773 ExcSLEPcEigenvectorConvergenceMismatchError(n_converged, n_eigenpairs));
774 AssertThrow ((real_eigenvectors.size()!=0) && (imag_eigenvectors.size()!=0),
775 ExcSLEPcWrappersUsageError());
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);
782 for (
unsigned int index=0; index<n_converged; ++index)
784 real_eigenvalues[index], imag_eigenvalues[index],
785 real_eigenvectors[index], imag_eigenvectors[index]);
788 template <
typename Vector>
793 std::vector<Vec> vecs(this_initial_space.size());
795 for (
unsigned int i = 0; i < this_initial_space.size(); i++)
797 Assert(this_initial_space[i].l2_norm()>0.0,
798 ExcMessage(
"Initial vectors should be nonzero."));
799 vecs[i] = this_initial_space[i];
807 #if DEAL_II_PETSC_VERSION_LT(3,1,0) 808 ierr = EPSSetInitialVector (
eps, &vecs[0]);
810 ierr = EPSSetInitialSpace (
eps, vecs.size(), &vecs[0]);
817 DEAL_II_NAMESPACE_CLOSE
819 #endif // DEAL_II_WITH_SLEPC
void set_initial_space(const std::vector< Vector > &initial_space)
bool delayed_reorthogonalization
const AdditionalData additional_data
void set_initial_vector(const PETScWrappers::VectorBase &this_initial_vector) DEAL_II_DEPRECATED
::ExceptionBase & ExcMessage(std::string arg1)
DeclException0(ExcSLEPcWrappersUsageError)
#define AssertThrow(cond, exc)
void solve(const PETScWrappers::MatrixBase &A, std::vector< PetscScalar > &eigenvalues, std::vector< OutputVector > &eigenvectors, const unsigned int n_eigenpairs=1)
const MPI_Comm mpi_communicator
const AdditionalData additional_data
EPSLanczosReorthogType reorthog
SolverBase(SolverControl &cn, const MPI_Comm &mpi_communicator)
void set_problem_type(EPSProblemType set_problem)
#define DEAL_II_DEPRECATED
DeclException2(ExcSLEPcEigenvectorConvergenceMismatchError, int, int,<< " The number of converged eigenvectors is "<< arg1<< " but "<< arg2<< " were requested. ")
void set_target_eigenvalue(const PetscScalar &this_target)
void set_matrices(const PETScWrappers::MatrixBase &A)
void set_which_eigenpairs(EPSWhich set_which)
#define Assert(cond, exc)
const AdditionalData additional_data
EPSConvergedReason reason
const AdditionalData additional_data
const AdditionalData additional_data
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)
const AdditionalData additional_data
const AdditionalData additional_data
void set_transformation(SLEPcWrappers::TransformationBase &this_transformation)
SolverControl & solver_control
DeclException1(ExcSLEPcError, int,<< " An error with error number "<< arg1<< " occurred while calling a SLEPc function")
SolverControl & control() const
void get_solver_state(const SolverControl::State state)