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,
const unsigned int *n,
char *which,
32 const 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 dneupd_(
int *rvec,
char *howmany,
int *select,
double *d,
38 double *di,
double *z,
int *ldz,
double *sigmar,
39 double *sigmai,
double *workev,
char *bmat,
const unsigned int *n,
char *which,
40 const unsigned int *nev,
const double *tol,
double *resid,
int *ncv,
41 double *v,
int *ldv,
int *iparam,
int *ipntr,
42 double *workd,
double *workl,
int *lworkl,
int *info);
105 algebraically_largest,
106 algebraically_smallest,
111 largest_imaginary_part,
112 smallest_imaginary_part,
122 const unsigned int number_of_arnoldi_vectors;
125 const unsigned int number_of_arnoldi_vectors = 15,
175 template <
typename VectorType,
typename MatrixType1,
176 typename MatrixType2,
typename INVERSE>
177 void solve (
const MatrixType1 &A,
178 const MatrixType2 &B,
179 const INVERSE &inverse,
180 std::vector<std::complex<double> > &eigenvalues,
181 std::vector<VectorType> &eigenvectors,
182 const unsigned int n_eigenvalues = 0);
203 <<
"Number of wanted eigenvalues " << arg1
204 <<
" is larger that the size of the matrix " << arg2);
207 <<
"Number of Arnoldi vectors " << arg1
208 <<
" is larger that the size of the matrix " << arg2);
211 <<
"Number of Arnoldi vectors " << arg1
212 <<
" is too small to obtain " << arg2
216 <<
" is not supported. Check documentation of ARPACK");
219 <<
" is not supported. Check documentation of ARPACK");
222 <<
"Error with dsaupd, info " << arg1
223 <<
". Check documentation of ARPACK");
226 <<
"Error with dneupd, info " << arg1
227 <<
". Check documentation of ARPACK");
230 <<
"Maximum number " << arg1
231 <<
" of iterations reached.");
234 "No shifts could be applied during implicit" 235 " Arnoldi update, try increasing the number of" 236 " Arnoldi vectors.");
241 ArpackSolver::AdditionalData::
242 AdditionalData (
const unsigned int number_of_arnoldi_vectors,
245 number_of_arnoldi_vectors(number_of_arnoldi_vectors),
246 eigenvalue_of_interest(eigenvalue_of_interest)
260 template <
typename VectorType,
typename MatrixType1,
261 typename MatrixType2,
typename INVERSE>
264 const MatrixType2 &mass_matrix,
265 const INVERSE &inverse,
266 std::vector<std::complex<double> > &eigenvalues,
267 std::vector<VectorType> &eigenvectors,
268 const unsigned int n_eigenvalues)
274 const unsigned int n = eigenvectors[0].size();
275 const unsigned int n_inside_arpack = eigenvectors[0].size();
278 const unsigned int nev = (n_eigenvalues == 0) ? eigenvalues.size() : n_eigenvalues;
284 Assert (n_eigenvalues < n,
285 ExcInvalidNumberofEigenvalues(nev, n));
288 ExcInvalidNumberofArnoldiVectors(
292 ExcSmallNumberofArnoldiVectors(
315 case algebraically_largest:
316 std::strcpy (which,
"LA");
318 case algebraically_smallest:
319 std::strcpy (which,
"SA");
321 case largest_magnitude:
322 std::strcpy (which,
"LM");
324 case smallest_magnitude:
325 std::strcpy (which,
"SM");
327 case largest_real_part:
328 std::strcpy (which,
"LR");
330 case smallest_real_part:
331 std::strcpy (which,
"SR");
333 case largest_imaginary_part:
334 std::strcpy (which,
"LI");
336 case smallest_imaginary_part:
337 std::strcpy (which,
"SI");
340 std::strcpy (which,
"BE");
348 std::vector<double> resid(n, 1.);
355 std::vector<double> v (ldv*ncv, 0.0);
358 std::vector<int> iparam (11, 0);
371 std::vector<int> ipntr (14, 0);
375 workd =
new double[3*n];
377 for (
unsigned int i=0; i<3*n; ++i)
380 int lworkl = 3*ncv*(ncv + 6);
381 std::vector<double> workl (lworkl, 0.);
388 dnaupd_(&ido, bmat, &n_inside_arpack, which, &nev, &tol,
389 &resid[0], &ncv, &v[0], &ldv, &iparam[0], &ipntr[0],
390 workd, &workl[0], &lworkl, &info);
404 VectorType src,dst,tmp;
405 src.reinit(eigenvectors[0]);
411 src(i) = *(workd+ipntr[0]-1+i);
414 mass_matrix.vmult(tmp, src);
416 inverse.vmult(dst,tmp);
419 *(workd+ipntr[1]-1+i) = dst(i);
426 VectorType src,dst,tmp, tmp2;
427 src.reinit(eigenvectors[0]);
434 src(i) = *(workd+ipntr[2]-1+i);
435 tmp(i) = *(workd+ipntr[0]-1+i);
438 inverse.vmult(dst,src);
441 *(workd+ipntr[1]-1+i) = dst(i);
449 src.reinit(eigenvectors[0]);
453 src(i) = *(workd+ipntr[0]-1+i);
456 mass_matrix.vmult(dst, src);
459 *(workd+ipntr[1]-1+i) = dst(i);
465 Assert (
false, ExcArpackIdo(ido));
471 Assert (
false, ExcArpackMode(mode));
478 Assert (
false, ExcArpackInfodsaupd(info));
490 std::vector<int> select (ncv, 1);
494 std::vector<double> z (ldz*ncv, 0.);
500 std::vector<double> workev (lworkev, 0.);
502 std::vector<double> eigenvalues_real (nev, 0.);
503 std::vector<double> eigenvalues_im (nev, 0.);
506 dneupd_(&rvec, &howmany, &select[0], &eigenvalues_real[0],
507 &eigenvalues_im[0], &z[0], &ldz, &sigmar, &sigmai,
508 &workev[0], bmat, &n_inside_arpack, which, &nev, &tol,
509 &resid[0], &ncv, &v[0], &ldv,
510 &iparam[0], &ipntr[0], workd, &workl[0], &lworkl, &info);
518 Assert (
false, ExcArpackNoShifts());
522 Assert (
false, ExcArpackInfodneupd(info));
526 const unsigned int n_eigenvecs = eigenvectors.size();
528 for (
unsigned int j=0; j<n; ++j)
529 eigenvectors[i](j) = v[i*n+j];
536 for (
size_type i=0; i<eigenvalues.size(); ++i)
537 eigenvalues[i] = std::complex<double> (eigenvalues_real[i],
549 DEAL_II_NAMESPACE_CLOSE
SolverControl & solver_control
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
types::global_dof_index size_type
#define DeclException1(Exception1, type1, outsequence)
DeclException2(ExcInvalidNumberofEigenvalues, int, int,<< "Number of wanted eigenvalues "<< arg1<< " is larger that the size of the matrix "<< arg2)
unsigned int global_dof_index
#define Assert(cond, exc)
#define DeclExceptionMsg(Exception, defaulttext)
ArpackSolver(SolverControl &control, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
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
SolverControl & control() const