16 #ifndef dealii_arpack_solver_h 17 #define dealii_arpack_solver_h 19 #include <deal.II/base/config.h> 20 #include <deal.II/base/smartpointer.h> 21 #include <deal.II/lac/solver_control.h> 26 #ifdef DEAL_II_WITH_ARPACK 28 DEAL_II_NAMESPACE_OPEN
31 extern "C" void dnaupd_(
int *ido,
char *bmat,
unsigned int *n,
char *which,
32 unsigned int *nev,
const double *tol,
double *resid,
int *ncv,
33 double *v,
int *ldv,
int *iparam,
int *ipntr,
34 double *workd,
double *workl,
int *lworkl,
37 extern "C" void dsaupd_(
int *ido,
char *bmat,
unsigned int *n,
char *which,
38 unsigned int *nev,
double *tol,
double *resid,
int *ncv,
39 double *v,
int *ldv,
int *iparam,
int *ipntr,
40 double *workd,
double *workl,
int *lworkl,
43 extern "C" void dneupd_(
int *rvec,
char *howmany,
int *select,
double *d,
44 double *di,
double *z,
int *ldz,
double *sigmar,
45 double *sigmai,
double *workev,
char *bmat,
unsigned int *n,
char *which,
46 unsigned int *nev,
double *tol,
double *resid,
int *ncv,
47 double *v,
int *ldv,
int *iparam,
int *ipntr,
48 double *workd,
double *workl,
int *lworkl,
int *info);
50 extern "C" void dseupd_(
int *rvec,
char *howmany,
int *select,
double *d,
51 double *z,
int *ldz,
double *sigmar,
52 char *bmat,
unsigned int *n,
char *which,
53 unsigned int *nev,
double *tol,
double *resid,
int *ncv,
54 double *v,
int *ldv,
int *iparam,
int *ipntr,
55 double *workd,
double *workl,
int *lworkl,
int *info);
209 template <
typename VectorType>
217 void set_shift(
const std::complex<double> sigma);
268 template <
typename VectorType,
typename MatrixType1,
269 typename MatrixType2,
typename INVERSE>
270 void solve (
const MatrixType1 &A,
271 const MatrixType2 &B,
272 const INVERSE &inverse,
275 const unsigned int n_eigenvalues = 0);
294 std::vector<double> resid;
313 <<
"Number of wanted eigenvalues " << arg1
314 <<
" is larger that the size of the matrix " << arg2);
317 <<
"Number of wanted eigenvalues " << arg1
318 <<
" is larger that the size of eigenvectors " << arg2);
321 <<
"To store the real and complex parts of " << arg1
322 <<
" eigenvectors in real-valued vectors, their size (currently set to " << arg2
323 <<
") should be greater than or equal to " << arg1+1);
326 <<
"Number of wanted eigenvalues " << arg1
327 <<
" is larger that the size of eigenvalues " << arg2);
330 <<
"Number of Arnoldi vectors " << arg1
331 <<
" is larger that the size of the matrix " << arg2);
334 <<
"Number of Arnoldi vectors " << arg1
335 <<
" is too small to obtain " << arg2
339 <<
" is not supported. Check documentation of ARPACK");
342 <<
" is not supported. Check documentation of ARPACK");
345 <<
"Error with dsaupd, info " << arg1
346 <<
". Check documentation of ARPACK");
349 <<
"Error with dnaupd, info " << arg1
350 <<
". Check documentation of ARPACK");
353 <<
"Error with dseupd, info " << arg1
354 <<
". Check documentation of ARPACK");
357 <<
"Error with dneupd, info " << arg1
358 <<
". Check documentation of ARPACK");
361 <<
"Maximum number " << arg1
362 <<
" of iterations reached.");
365 "No shifts could be applied during implicit" 366 " Arnoldi update, try increasing the number of" 367 " Arnoldi vectors.");
375 const bool symmetric)
377 number_of_arnoldi_vectors(number_of_arnoldi_vectors),
378 eigenvalue_of_interest(eigenvalue_of_interest),
385 ExcMessage(
"'largest real part' can only be used for non-symmetric problems!"));
387 ExcMessage(
"'smallest real part' can only be used for non-symmetric problems!"));
389 ExcMessage(
"'largest imaginary part' can only be used for non-symmetric problems!"));
391 ExcMessage(
"'smallest imaginary part' can only be used for non-symmetric problems!"));
420 template <
typename VectorType>
426 resid.resize(vec.size());
427 for (
size_type i = 0; i < vec.size(); ++i)
432 template <
typename VectorType,
typename MatrixType1,
433 typename MatrixType2,
typename INVERSE>
436 const MatrixType2 &mass_matrix,
437 const INVERSE &inverse,
440 const unsigned int n_eigenvalues)
446 const unsigned int nev_const = (n_eigenvalues == 0) ?
eigenvalues.size() : n_eigenvalues;
448 unsigned int nev = nev_const;
494 std::strcpy (which,
"LA");
497 std::strcpy (which,
"SA");
500 std::strcpy (which,
"LM");
503 std::strcpy (which,
"SM");
506 std::strcpy (which,
"LR");
509 std::strcpy (which,
"SR");
512 std::strcpy (which,
"LI");
515 std::strcpy (which,
"SI");
518 std::strcpy (which,
"BE");
534 std::vector<double> v (ldv*ncv, 0.0);
537 std::vector<int> iparam (11, 0);
548 std::vector<int> ipntr (14, 0);
551 std::vector<double> workd (3*n, 0.);
553 std::vector<double> workl (lworkl, 0.);
562 dsaupd_(&ido, bmat, &n, which, &nev, &tol,
563 resid.data(), &ncv, v.data(), &ldv, iparam.data(), ipntr.data(),
564 workd.data(), workl.data(), &lworkl, &info);
566 dnaupd_(&ido, bmat, &n, which, &nev, &tol,
567 resid.data(), &ncv, v.data(), &ldv, iparam.data(), ipntr.data(),
568 workd.data(), workl.data(), &lworkl, &info);
582 VectorType src,dst,tmp;
589 src(i) = workd[ipntr[0]-1+i];
592 mass_matrix.vmult(tmp, src);
594 inverse.vmult(dst,tmp);
597 workd[ipntr[1]-1+i] = dst(i);
604 VectorType src,dst,tmp, tmp2;
612 src(i) = workd[ipntr[2]-1+i];
613 tmp(i) = workd[ipntr[0]-1+i];
616 inverse.vmult(dst,src);
619 workd[ipntr[1]-1+i] = dst(i);
631 src(i) = workd[ipntr[0]-1+i];
634 mass_matrix.vmult(dst, src);
637 workd[ipntr[1]-1+i] = dst(i);
671 std::vector<int> select (ncv, 1);
675 std::vector<double> eigenvalues_real (nev+1, 0.);
676 std::vector<double> eigenvalues_im (nev+1, 0.);
681 std::vector<double> z (ldz*nev, 0.);
682 dseupd_(&rvec, &howmany, select.data(), eigenvalues_real.data(),
683 z.data(), &ldz, &
sigmar, bmat, &n, which, &nev, &tol,
684 resid.data(), &ncv, v.data(), &ldv,
685 iparam.data(), ipntr.data(), workd.data(), workl.data(), &lworkl, &info);
689 std::vector<double> workev (3*ncv, 0.);
690 dneupd_(&rvec, &howmany, select.data(), eigenvalues_real.data(),
691 eigenvalues_im.data(), v.data(), &ldz, &
sigmar, &
sigmai,
692 workev.data(), bmat, &n, which, &nev, &tol,
693 resid.data(), &ncv, v.data(), &ldv,
694 iparam.data(), ipntr.data(), workd.data(), workl.data(), &lworkl, &info);
715 for (
unsigned int i=0; i<nev; ++i)
716 for (
unsigned int j=0; j<n; ++j)
719 for (
unsigned int i=0; i<nev_const; ++i)
720 eigenvalues[i] = std::complex<double> (eigenvalues_real[i],
732 DEAL_II_NAMESPACE_CLOSE
AdditionalData(const unsigned int number_of_arnoldi_vectors=15, const WhichEigenvalues eigenvalue_of_interest=largest_magnitude, const bool symmetric=false)
SolverControl & solver_control
#define DeclException2(Exception2, type1, type2, outsequence)
std::array< std::pair< Number, Tensor< 1, dim, Number > >, std::integral_constant< int, dim >::value > eigenvectors(const SymmetricTensor< 2, dim, Number > &T, const SymmetricTensorEigenvectorMethod method=SymmetricTensorEigenvectorMethod::ql_implicit_shifts)
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)
void set_shift(const std::complex< double > sigma)
std::array< std::pair< Number, Tensor< 1, dim, Number > >, std::integral_constant< int, dim >::value > eigenvectors(const SymmetricTensor< 2, dim, Number > &T, const SymmetricTensorEigenvectorMethod method=SymmetricTensorEigenvectorMethod::ql_implicit_shifts)
static ::ExceptionBase & ArpackExcInvalidEigenvalueSize(int arg1, int arg2)
static ::ExceptionBase & ArpackExcInvalidEigenvectorSizeNonsymmetric(int arg1, int arg2)
static ::ExceptionBase & ArpackExcArpackMode(int arg1)
types::global_dof_index size_type
static ::ExceptionBase & ExcMessage(std::string arg1)
#define DeclException1(Exception1, type1, outsequence)
unsigned int global_dof_index
#define Assert(cond, exc)
static ::ExceptionBase & ArpackExcInvalidNumberofArnoldiVectors(int arg1, int arg2)
#define DeclExceptionMsg(Exception, defaulttext)
ArpackSolver(SolverControl &control, const AdditionalData &data=AdditionalData())
static ::ExceptionBase & ArpackExcSmallNumberofArnoldiVectors(int arg1, int arg2)
static ::ExceptionBase & ArpackExcInvalidEigenvectorSize(int arg1, int arg2)
const AdditionalData additional_data
static ::ExceptionBase & ArpackExcArpackInfodneupd(int arg1)
static ::ExceptionBase & ArpackExcArpackInfoMaxIt(int arg1)
static ::ExceptionBase & ArpackExcInvalidNumberofEigenvalues(int arg1, int arg2)
void set_initial_vector(const VectorType &vec)
void solve(const MatrixType1 &A, const MatrixType2 &B, const INVERSE &inverse, std::vector< std::complex< double > > &eigenvalues, std::vector< VectorType > &eigenvectors, const unsigned int n_eigenvalues=0)
unsigned int max_steps() const
static ::ExceptionBase & ArpackExcArpackInfodseupd(int arg1)
bool initial_vector_provided
static ::ExceptionBase & ArpackExcArpackInfodsaupd(int arg1)
const unsigned int number_of_arnoldi_vectors
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)
static ::ExceptionBase & ArpackExcArpackIdo(int arg1)
SolverControl & control() const
const WhichEigenvalues eigenvalue_of_interest
static ::ExceptionBase & ArpackExcArpackInfodnaupd(int arg1)
static ::ExceptionBase & ArpackExcArpackNoShifts()