template <typename OutputVector>
void
solve (const PETScWrappers::MatrixBase &A,
- std::vector<double> &kr,
- std::vector<OutputVector> &vr,
- const unsigned int n_eigenvectors = 1);
+ std::vector<double> &r_eigenvalues,
+ std::vector<OutputVector> &r_eigenvectors = std::vector<OutputVector> (),
+ const unsigned int n_eigenvectors = 1);
/**
* Same as above, but here a composite method for solving the
void
solve (const PETScWrappers::MatrixBase &A,
const PETScWrappers::MatrixBase &B,
- std::vector<double> &kr,
- std::vector<OutputVector> &vr,
+ std::vector<double> &r_eigenvalues,
+ std::vector<OutputVector> &r_eigenvectors = std::vector<OutputVector> (),
const unsigned int n_eigenvectors = 1);
/**
* \text{n\_converged}-1$.
*/
void
- get_eigenpair (const unsigned int index,
-#ifndef PETSC_USE_COMPLEX
+ get_eigenpair (const unsigned int index,
double &kr,
-#else
- std::complex<double> &kr,
-#endif
PETScWrappers::VectorBase &vr);
/**
std::vector<OutputVector> &vr,
const unsigned int n_eigenvectors)
{
- // Panic if no eigenpairs are wanted.
- AssertThrow (n_eigenvectors != 0, ExcSLEPcWrappersUsageError());
-
- unsigned int n_converged = 0;
+ // Panic if the number of eigenpairs wanted is out of bounds.
+ AssertThrow ((n_eigenvectors > 0) && (n_eigenvectors <= A.m ()),
+ ExcSLEPcWrappersUsageError());
// Set the matrices of the problem
set_matrices (A);
// and solve
+ unsigned int n_converged = 0;
solve (n_eigenvectors, &n_converged);
if (n_converged > n_eigenvectors)
std::vector<OutputVector> &vr,
const unsigned int n_eigenvectors)
{
- // Panic if no eigenpairs are wanted.
- AssertThrow (n_eigenvectors != 0, ExcSLEPcWrappersUsageError());
+ // Guard against incompatible matrix sizes:
+ AssertThrow (A.m() == B.m (), ExcDimensionMismatch(A.m(), B.m()));
+ AssertThrow (A.n() == B.n (), ExcDimensionMismatch(A.n(), B.n()));
- unsigned int n_converged = 0;
+ // Panic if the number of eigenpairs wanted is out of bounds.
+ AssertThrow ((n_eigenvectors > 0) && (n_eigenvectors <= A.m ()),
+ ExcSLEPcWrappersUsageError());
// Set the matrices of the problem
set_matrices (A, B);
// and solve
+ unsigned int n_converged = 0;
solve (n_eigenvectors, &n_converged);
if (n_converged >= n_eigenvectors)