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);
168 const unsigned int number_of_arnoldi_vectors;
170 const bool symmetric;
172 const unsigned int number_of_arnoldi_vectors = 15,
174 const bool symmetric =
false);
191 template <
typename VectorType>
199 void set_shift(
const std::complex<double> sigma);
250 template <
typename VectorType,
typename MatrixType1,
251 typename MatrixType2,
typename INVERSE>
252 void solve (
const MatrixType1 &A,
253 const MatrixType2 &B,
254 const INVERSE &inverse,
255 std::vector<std::complex<double> > &eigenvalues,
256 std::vector<VectorType> &eigenvectors,
257 const unsigned int n_eigenvalues = 0);
276 std::vector<double> resid;
295 <<
"Number of wanted eigenvalues " << arg1
296 <<
" is larger that the size of the matrix " << arg2);
299 <<
"Number of wanted eigenvalues " << arg1
300 <<
" is larger that the size of eigenvectors " << arg2);
303 <<
"To store the real and complex parts of " << arg1
304 <<
" eigenvectors in real-valued vectors, their size (currently set to " << arg2
305 <<
") should be greater than or equal to " << arg1+1);
308 <<
"Number of wanted eigenvalues " << arg1
309 <<
" is larger that the size of eigenvalues " << arg2);
312 <<
"Number of Arnoldi vectors " << arg1
313 <<
" is larger that the size of the matrix " << arg2);
316 <<
"Number of Arnoldi vectors " << arg1
317 <<
" is too small to obtain " << arg2
321 <<
" is not supported. Check documentation of ARPACK");
324 <<
" is not supported. Check documentation of ARPACK");
327 <<
"Error with dsaupd, info " << arg1
328 <<
". Check documentation of ARPACK");
331 <<
"Error with dnaupd, info " << arg1
332 <<
". Check documentation of ARPACK");
335 <<
"Error with dseupd, info " << arg1
336 <<
". Check documentation of ARPACK");
339 <<
"Error with dneupd, info " << arg1
340 <<
". Check documentation of ARPACK");
343 <<
"Maximum number " << arg1
344 <<
" of iterations reached.");
347 "No shifts could be applied during implicit" 348 " Arnoldi update, try increasing the number of" 349 " Arnoldi vectors.");
354 ArpackSolver::AdditionalData::
355 AdditionalData (
const unsigned int number_of_arnoldi_vectors,
356 const WhichEigenvalues eigenvalue_of_interest,
357 const bool symmetric)
359 number_of_arnoldi_vectors(number_of_arnoldi_vectors),
360 eigenvalue_of_interest(eigenvalue_of_interest),
367 ExcMessage(
"'largest real part' can only be used for non-symmetric problems!"));
369 ExcMessage(
"'smallest real part' can only be used for non-symmetric problems!"));
371 ExcMessage(
"'largest imaginary part' can only be used for non-symmetric problems!"));
373 ExcMessage(
"'smallest imaginary part' can only be used for non-symmetric problems!"));
382 solver_control (control),
383 additional_data (data),
384 initial_vector_provided(false),
402 template <
typename VectorType>
408 resid.resize(vec.size());
409 for (
size_type i = 0; i < vec.size(); ++i)
414 template <
typename VectorType,
typename MatrixType1,
415 typename MatrixType2,
typename INVERSE>
418 const MatrixType2 &mass_matrix,
419 const INVERSE &inverse,
420 std::vector<std::complex<double> > &eigenvalues,
421 std::vector<VectorType> &eigenvectors,
422 const unsigned int n_eigenvalues)
425 unsigned int n = eigenvectors[0].size();
428 const unsigned int nev_const = (n_eigenvalues == 0) ? eigenvalues.size() : n_eigenvalues;
430 unsigned int nev = nev_const;
435 Assert (nev <= eigenvectors.size(),
439 Assert (nev+1 <= eigenvectors.size(),
442 Assert (nev <= eigenvalues.size(),
482 std::strcpy (which,
"LA");
485 std::strcpy (which,
"SA");
488 std::strcpy (which,
"LM");
491 std::strcpy (which,
"SM");
494 std::strcpy (which,
"LR");
497 std::strcpy (which,
"SR");
500 std::strcpy (which,
"LI");
503 std::strcpy (which,
"SI");
506 std::strcpy (which,
"BE");
522 std::vector<double> v (ldv*ncv, 0.0);
525 std::vector<int> iparam (11, 0);
538 std::vector<int> ipntr (14, 0);
541 std::vector<double> workd (3*n, 0.);
542 int lworkl =
additional_data.symmetric ? ncv*ncv + 8*ncv : 3*ncv*ncv+6*ncv;
543 std::vector<double> workl (lworkl, 0.);
552 dsaupd_(&ido, bmat, &n, which, &nev, &tol,
553 &resid[0], &ncv, &v[0], &ldv, &iparam[0], &ipntr[0],
554 &workd[0], &workl[0], &lworkl, &info);
556 dnaupd_(&ido, bmat, &n, which, &nev, &tol,
557 &resid[0], &ncv, &v[0], &ldv, &iparam[0], &ipntr[0],
558 &workd[0], &workl[0], &lworkl, &info);
572 VectorType src,dst,tmp;
573 src.reinit(eigenvectors[0]);
579 src(i) = workd[ipntr[0]-1+i];
582 mass_matrix.vmult(tmp, src);
584 inverse.vmult(dst,tmp);
587 workd[ipntr[1]-1+i] = dst(i);
594 VectorType src,dst,tmp, tmp2;
595 src.reinit(eigenvectors[0]);
602 src(i) = workd[ipntr[2]-1+i];
603 tmp(i) = workd[ipntr[0]-1+i];
606 inverse.vmult(dst,src);
609 workd[ipntr[1]-1+i] = dst(i);
617 src.reinit(eigenvectors[0]);
621 src(i) = workd[ipntr[0]-1+i];
624 mass_matrix.vmult(dst, src);
627 workd[ipntr[1]-1+i] = dst(i);
663 std::vector<int> select (ncv, 1);
667 std::vector<double> eigenvalues_real (nev+1, 0.);
668 std::vector<double> eigenvalues_im (nev+1, 0.);
673 std::vector<double> z (ldz*nev, 0.);
674 dseupd_(&rvec, &howmany, &select[0], &eigenvalues_real[0],
675 &z[0], &ldz, &
sigmar, bmat, &n, which, &nev, &tol,
676 &resid[0], &ncv, &v[0], &ldv,
677 &iparam[0], &ipntr[0], &workd[0], &workl[0], &lworkl, &info);
681 std::vector<double> workev (3*ncv, 0.);
682 dneupd_(&rvec, &howmany, &select[0], &eigenvalues_real[0],
684 &workev[0], bmat, &n, which, &nev, &tol,
685 &resid[0], &ncv, &v[0], &ldv,
686 &iparam[0], &ipntr[0], &workd[0], &workl[0], &lworkl, &info);
707 for (
unsigned int i=0; i<nev; ++i)
708 for (
unsigned int j=0; j<n; ++j)
709 eigenvectors[i](j) = v[i*n+j];
711 for (
unsigned int i=0; i<nev_const; ++i)
712 eigenvalues[i] = std::complex<double> (eigenvalues_real[i],
724 DEAL_II_NAMESPACE_CLOSE
SolverControl & solver_control
#define DeclException2(Exception2, type1, type2, outsequence)
void set_shift(const std::complex< double > sigma)
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)
static ::ExceptionBase & ArpackExcArpackIdo(int arg1)
SolverControl & control() const
static ::ExceptionBase & ArpackExcArpackInfodnaupd(int arg1)
static ::ExceptionBase & ArpackExcArpackNoShifts()