Reference documentation for deal.II version 8.5.1
arpack_solver.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__arpack_solver_h
17 #define dealii__arpack_solver_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/base/smartpointer.h>
21 #include <deal.II/lac/solver_control.h>
22 
23 #include <cstring>
24 
25 
26 #ifdef DEAL_II_WITH_ARPACK
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 
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,
35  int *info);
36 
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,
41  int *info);
42 
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);
49 
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);
56 
106 class ArpackSolver : public Subscriptor
107 {
108 public:
113 
114 
120  {
160  };
161 
167  {
168  const unsigned int number_of_arnoldi_vectors;
169  const WhichEigenvalues eigenvalue_of_interest;
170  const bool symmetric;
172  const unsigned int number_of_arnoldi_vectors = 15,
173  const WhichEigenvalues eigenvalue_of_interest = largest_magnitude,
174  const bool symmetric = false);
175  };
176 
180  SolverControl &control () const;
181 
186  const AdditionalData &data = AdditionalData());
187 
191  template <typename VectorType>
192  void set_initial_vector(const VectorType &vec);
193 
199  void set_shift(const std::complex<double> sigma);
200 
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);
258 
259 protected:
260 
266 
271 
276  std::vector<double> resid;
277 
281  double sigmar;
282 
286  double sigmai;
287 
288 
289 private:
290 
295  << "Number of wanted eigenvalues " << arg1
296  << " is larger that the size of the matrix " << arg2);
297 
299  << "Number of wanted eigenvalues " << arg1
300  << " is larger that the size of eigenvectors " << arg2);
301 
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);
306 
308  << "Number of wanted eigenvalues " << arg1
309  << " is larger that the size of eigenvalues " << arg2);
310 
312  << "Number of Arnoldi vectors " << arg1
313  << " is larger that the size of the matrix " << arg2);
314 
316  << "Number of Arnoldi vectors " << arg1
317  << " is too small to obtain " << arg2
318  << " eigenvalues");
319 
320  DeclException1 (ArpackExcArpackIdo, int, << "This ido " << arg1
321  << " is not supported. Check documentation of ARPACK");
322 
323  DeclException1 (ArpackExcArpackMode, int, << "This mode " << arg1
324  << " is not supported. Check documentation of ARPACK");
325 
327  << "Error with dsaupd, info " << arg1
328  << ". Check documentation of ARPACK");
329 
331  << "Error with dnaupd, info " << arg1
332  << ". Check documentation of ARPACK");
333 
335  << "Error with dseupd, info " << arg1
336  << ". Check documentation of ARPACK");
337 
339  << "Error with dneupd, info " << arg1
340  << ". Check documentation of ARPACK");
341 
343  << "Maximum number " << arg1
344  << " of iterations reached.");
345 
347  "No shifts could be applied during implicit"
348  " Arnoldi update, try increasing the number of"
349  " Arnoldi vectors.");
350 };
351 
352 
353 inline
354 ArpackSolver::AdditionalData::
355 AdditionalData (const unsigned int number_of_arnoldi_vectors,
356  const WhichEigenvalues eigenvalue_of_interest,
357  const bool symmetric)
358  :
359  number_of_arnoldi_vectors(number_of_arnoldi_vectors),
360  eigenvalue_of_interest(eigenvalue_of_interest),
361  symmetric(symmetric)
362 {
363  //Check for possible options for symmetric problems
364  if (symmetric)
365  {
366  Assert(eigenvalue_of_interest!=largest_real_part,
367  ExcMessage("'largest real part' can only be used for non-symmetric problems!"));
368  Assert(eigenvalue_of_interest!=smallest_real_part,
369  ExcMessage("'smallest real part' can only be used for non-symmetric problems!"));
370  Assert(eigenvalue_of_interest!=largest_imaginary_part,
371  ExcMessage("'largest imaginary part' can only be used for non-symmetric problems!"));
372  Assert(eigenvalue_of_interest!=smallest_imaginary_part,
373  ExcMessage("'smallest imaginary part' can only be used for non-symmetric problems!"));
374  }
375 }
376 
377 
378 inline
380  const AdditionalData &data)
381  :
382  solver_control (control),
383  additional_data (data),
384  initial_vector_provided(false),
385  sigmar(0.0),
386  sigmai(0.0)
387 {}
388 
389 
390 
391 
392 inline
393 void
394 ArpackSolver::set_shift(const std::complex<double> sigma)
395 {
396  sigmar = sigma.real();
397  sigmai = sigma.imag();
398 }
399 
400 
401 
402 template <typename VectorType>
403 inline
404 void ArpackSolver::
405 set_initial_vector(const VectorType &vec)
406 {
408  resid.resize(vec.size());
409  for (size_type i = 0; i < vec.size(); ++i)
410  resid[i] = vec[i];
411 }
412 
413 
414 template <typename VectorType, typename MatrixType1,
415  typename MatrixType2, typename INVERSE>
416 inline
417 void ArpackSolver::solve (const MatrixType1 &/*system_matrix*/,
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)
423 {
424  // Problem size
425  unsigned int n = eigenvectors[0].size();
426 
427  // Number of eigenvalues
428  const unsigned int nev_const = (n_eigenvalues == 0) ? eigenvalues.size() : n_eigenvalues;
429  // nev for arpack, which might change by plus one during dneupd
430  unsigned int nev = nev_const;
431 
432  // check input sizes
433  if (additional_data.symmetric)
434  {
435  Assert (nev <= eigenvectors.size(),
436  ArpackExcInvalidEigenvectorSize(nev, eigenvectors.size()));
437  }
438  else
439  Assert (nev+1 <= eigenvectors.size(),
440  ArpackExcInvalidEigenvectorSizeNonsymmetric(nev, eigenvectors.size()));
441 
442  Assert (nev <= eigenvalues.size(),
443  ArpackExcInvalidEigenvalueSize(nev, eigenvalues.size()));
444 
445  // check large enough problem size
446  Assert (nev < n,
448 
449  Assert (additional_data.number_of_arnoldi_vectors < n,
451  additional_data.number_of_arnoldi_vectors, n));
452 
453  // check whether we have enough Arnoldi vectors
454  Assert (additional_data.number_of_arnoldi_vectors > 2*nev+1,
456  additional_data.number_of_arnoldi_vectors, nev));
457 
458  /* ARPACK mode for dsaupd/dnaupd, here only mode 3,
459  * i.e. shift-invert mode
460  */
461  int mode = 3;
462 
463  // reverse communication parameter
464  int ido = 0;
465 
469  char bmat[2] = "G";
470 
478  char which[3];
479  switch (additional_data.eigenvalue_of_interest)
480  {
482  std::strcpy (which, "LA");
483  break;
485  std::strcpy (which, "SA");
486  break;
487  case largest_magnitude:
488  std::strcpy (which, "LM");
489  break;
490  case smallest_magnitude:
491  std::strcpy (which, "SM");
492  break;
493  case largest_real_part:
494  std::strcpy (which, "LR");
495  break;
496  case smallest_real_part:
497  std::strcpy (which, "SR");
498  break;
500  std::strcpy (which, "LI");
501  break;
503  std::strcpy (which, "SI");
504  break;
505  case both_ends:
506  std::strcpy (which, "BE");
507  break;
508  }
509 
510  // tolerance for ARPACK
511  double tol = control().tolerance();
512 
513  // if the starting vector is used it has to be in resid
514  if (!initial_vector_provided || resid.size() != n)
515  resid.resize(n, 1.);
516 
517  // number of Arnoldi basis vectors specified
518  // in additional_data
519  int ncv = additional_data.number_of_arnoldi_vectors;
520 
521  int ldv = n;
522  std::vector<double> v (ldv*ncv, 0.0);
523 
524  //information to the routines
525  std::vector<int> iparam (11, 0);
526 
527  iparam[0] = 1; // shift strategy
528 
529  // maximum number of iterations
530  iparam[2] = control().max_steps();
531 
537  iparam[6] = mode;
538  std::vector<int> ipntr (14, 0);
539 
540  // work arrays for ARPACK
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.);
544 
545  //information out of the iteration
546  int info = 1;
547 
548  while (ido != 99)
549  {
550  // call of ARPACK dsaupd/dnaupd routine
551  if (additional_data.symmetric)
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);
555  else
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);
559 
560  if (ido == 99)
561  break;
562 
563  switch (mode)
564  {
565  case 3:
566  {
567  switch (ido)
568  {
569  case -1:
570  {
571 
572  VectorType src,dst,tmp;
573  src.reinit(eigenvectors[0]);
574  dst.reinit(src);
575  tmp.reinit(src);
576 
577 
578  for (size_type i=0; i<src.size(); ++i)
579  src(i) = workd[ipntr[0]-1+i];
580 
581  // multiplication with mass matrix M
582  mass_matrix.vmult(tmp, src);
583  // solving linear system
584  inverse.vmult(dst,tmp);
585 
586  for (size_type i=0; i<dst.size(); ++i)
587  workd[ipntr[1]-1+i] = dst(i);
588  }
589  break;
590 
591  case 1:
592  {
593 
594  VectorType src,dst,tmp, tmp2;
595  src.reinit(eigenvectors[0]);
596  dst.reinit(src);
597  tmp.reinit(src);
598  tmp2.reinit(src);
599 
600  for (size_type i=0; i<src.size(); ++i)
601  {
602  src(i) = workd[ipntr[2]-1+i];
603  tmp(i) = workd[ipntr[0]-1+i];
604  }
605  // solving linear system
606  inverse.vmult(dst,src);
607 
608  for (size_type i=0; i<dst.size(); ++i)
609  workd[ipntr[1]-1+i] = dst(i);
610  }
611  break;
612 
613  case 2:
614  {
615 
616  VectorType src,dst;
617  src.reinit(eigenvectors[0]);
618  dst.reinit(src);
619 
620  for (size_type i=0; i<src.size(); ++i)
621  src(i) = workd[ipntr[0]-1+i];
622 
623  // Multiplication with mass matrix M
624  mass_matrix.vmult(dst, src);
625 
626  for (size_type i=0; i<dst.size(); ++i)
627  workd[ipntr[1]-1+i] = dst(i);
628 
629  }
630  break;
631 
632  default:
633  Assert (false, ArpackExcArpackIdo(ido));
634  break;
635  }
636  }
637  break;
638  default:
639  Assert (false, ArpackExcArpackMode(mode));
640  break;
641  }
642  }
643 
644  if (info<0)
645  {
646  if (additional_data.symmetric)
647  {
648  Assert (false, ArpackExcArpackInfodsaupd(info));
649  }
650  else
651  Assert (false, ArpackExcArpackInfodnaupd(info));
652  }
653  else
654  {
658  int rvec = 1;
659 
660  // which eigenvectors
661  char howmany = 'A';
662 
663  std::vector<int> select (ncv, 1);
664 
665  int ldz = n;
666 
667  std::vector<double> eigenvalues_real (nev+1, 0.);
668  std::vector<double> eigenvalues_im (nev+1, 0.);
669 
670  // call of ARPACK dseupd/dneupd routine
671  if (additional_data.symmetric)
672  {
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);
678  }
679  else
680  {
681  std::vector<double> workev (3*ncv, 0.);
682  dneupd_(&rvec, &howmany, &select[0], &eigenvalues_real[0],
683  &eigenvalues_im[0], &v[0], &ldz, &sigmar, &sigmai,
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);
687  }
688 
689  if (info == 1)
690  {
691  Assert (false, ArpackExcArpackInfoMaxIt(control().max_steps()));
692  }
693  else if (info == 3)
694  {
695  Assert (false, ArpackExcArpackNoShifts());
696  }
697  else if (info!=0)
698  {
699  if (additional_data.symmetric)
700  {
701  Assert (false, ArpackExcArpackInfodseupd(info));
702  }
703  else
704  Assert (false, ArpackExcArpackInfodneupd(info));
705  }
706 
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];
710 
711  for (unsigned int i=0; i<nev_const; ++i)
712  eigenvalues[i] = std::complex<double> (eigenvalues_real[i],
713  eigenvalues_im[i]);
714  }
715 }
716 
717 
718 inline
720 {
721  return solver_control;
722 }
723 
724 DEAL_II_NAMESPACE_CLOSE
725 
726 
727 #endif
728 #endif
SolverControl & solver_control
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
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)
Definition: exceptions.h:564
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ArpackExcInvalidNumberofArnoldiVectors(int arg1, int arg2)
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:553
ArpackSolver(SolverControl &control, const AdditionalData &data=AdditionalData())
static ::ExceptionBase & ArpackExcSmallNumberofArnoldiVectors(int arg1, int arg2)
static ::ExceptionBase & ArpackExcInvalidEigenvectorSize(int arg1, int arg2)
double tolerance() const
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()