Reference documentation for deal.II version 9.0.0
arpack_solver.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2017 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 
166  {
172  explicit AdditionalData(
173  const unsigned int number_of_arnoldi_vectors = 15,
175  const bool symmetric = false);
176 
182  const unsigned int number_of_arnoldi_vectors;
183 
188 
192  const bool symmetric;
193  };
194 
198  SolverControl &control () const;
199 
204  const AdditionalData &data = AdditionalData());
205 
209  template <typename VectorType>
210  void set_initial_vector(const VectorType &vec);
211 
217  void set_shift(const std::complex<double> sigma);
218 
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,
273  std::vector<std::complex<double> > &eigenvalues,
274  std::vector<VectorType> &eigenvectors,
275  const unsigned int n_eigenvalues = 0);
276 
277 protected:
278 
284 
289 
294  std::vector<double> resid;
295 
299  double sigmar;
300 
304  double sigmai;
305 
306 
307 private:
308 
313  << "Number of wanted eigenvalues " << arg1
314  << " is larger that the size of the matrix " << arg2);
315 
317  << "Number of wanted eigenvalues " << arg1
318  << " is larger that the size of eigenvectors " << arg2);
319 
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);
324 
326  << "Number of wanted eigenvalues " << arg1
327  << " is larger that the size of eigenvalues " << arg2);
328 
330  << "Number of Arnoldi vectors " << arg1
331  << " is larger that the size of the matrix " << arg2);
332 
334  << "Number of Arnoldi vectors " << arg1
335  << " is too small to obtain " << arg2
336  << " eigenvalues");
337 
338  DeclException1 (ArpackExcArpackIdo, int, << "This ido " << arg1
339  << " is not supported. Check documentation of ARPACK");
340 
341  DeclException1 (ArpackExcArpackMode, int, << "This mode " << arg1
342  << " is not supported. Check documentation of ARPACK");
343 
345  << "Error with dsaupd, info " << arg1
346  << ". Check documentation of ARPACK");
347 
349  << "Error with dnaupd, info " << arg1
350  << ". Check documentation of ARPACK");
351 
353  << "Error with dseupd, info " << arg1
354  << ". Check documentation of ARPACK");
355 
357  << "Error with dneupd, info " << arg1
358  << ". Check documentation of ARPACK");
359 
361  << "Maximum number " << arg1
362  << " of iterations reached.");
363 
365  "No shifts could be applied during implicit"
366  " Arnoldi update, try increasing the number of"
367  " Arnoldi vectors.");
368 };
369 
370 
371 inline
373 AdditionalData (const unsigned int number_of_arnoldi_vectors,
374  const WhichEigenvalues eigenvalue_of_interest,
375  const bool symmetric)
376  :
377  number_of_arnoldi_vectors(number_of_arnoldi_vectors),
378  eigenvalue_of_interest(eigenvalue_of_interest),
379  symmetric(symmetric)
380 {
381  //Check for possible options for symmetric problems
382  if (symmetric)
383  {
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!"));
392  }
393 }
394 
395 
396 inline
398  const AdditionalData &data)
399  :
401  additional_data (data),
403  sigmar(0.0),
404  sigmai(0.0)
405 {}
406 
407 
408 
409 
410 inline
411 void
412 ArpackSolver::set_shift(const std::complex<double> sigma)
413 {
414  sigmar = sigma.real();
415  sigmai = sigma.imag();
416 }
417 
418 
419 
420 template <typename VectorType>
421 inline
422 void ArpackSolver::
423 set_initial_vector(const VectorType &vec)
424 {
426  resid.resize(vec.size());
427  for (size_type i = 0; i < vec.size(); ++i)
428  resid[i] = vec[i];
429 }
430 
431 
432 template <typename VectorType, typename MatrixType1,
433  typename MatrixType2, typename INVERSE>
434 inline
435 void ArpackSolver::solve (const MatrixType1 &/*system_matrix*/,
436  const MatrixType2 &mass_matrix,
437  const INVERSE &inverse,
438  std::vector<std::complex<double> > &eigenvalues,
439  std::vector<VectorType> &eigenvectors,
440  const unsigned int n_eigenvalues)
441 {
442  // Problem size
443  unsigned int n = eigenvectors[0].size();
444 
445  // Number of eigenvalues
446  const unsigned int nev_const = (n_eigenvalues == 0) ? eigenvalues.size() : n_eigenvalues;
447  // nev for arpack, which might change by plus one during dneupd
448  unsigned int nev = nev_const;
449 
450  // check input sizes
452  {
453  Assert (nev <= eigenvectors.size(),
455  }
456  else
457  Assert (nev+1 <= eigenvectors.size(),
459 
460  Assert (nev <= eigenvalues.size(),
462 
463  // check large enough problem size
464  Assert (nev < n,
466 
470 
471  // check whether we have enough Arnoldi vectors
475 
476  // ARPACK mode for dsaupd/dnaupd, here only mode 3, i.e. shift-invert mode
477  int mode = 3;
478 
479  // reverse communication parameter
480  int ido = 0;
481 
482  // 'G' generalized eigenvalue problem 'I' standard eigenvalue problem
483  char bmat[2] = "G";
484 
485  // Specify the eigenvalues of interest, possible parameters "LA" algebraically
486  // largest "SA" algebraically smallest "LM" largest magnitude "SM" smallest
487  // magnitude "LR" largest real part "SR" smallest real part "LI" largest
488  // imaginary part "SI" smallest imaginary part "BE" both ends of spectrum
489  // simultaneous.
490  char which[3];
492  {
494  std::strcpy (which, "LA");
495  break;
497  std::strcpy (which, "SA");
498  break;
499  case largest_magnitude:
500  std::strcpy (which, "LM");
501  break;
502  case smallest_magnitude:
503  std::strcpy (which, "SM");
504  break;
505  case largest_real_part:
506  std::strcpy (which, "LR");
507  break;
508  case smallest_real_part:
509  std::strcpy (which, "SR");
510  break;
512  std::strcpy (which, "LI");
513  break;
515  std::strcpy (which, "SI");
516  break;
517  case both_ends:
518  std::strcpy (which, "BE");
519  break;
520  }
521 
522  // tolerance for ARPACK
523  double tol = control().tolerance();
524 
525  // if the starting vector is used it has to be in resid
526  if (!initial_vector_provided || resid.size() != n)
527  resid.resize(n, 1.);
528 
529  // number of Arnoldi basis vectors specified
530  // in additional_data
532 
533  int ldv = n;
534  std::vector<double> v (ldv*ncv, 0.0);
535 
536  //information to the routines
537  std::vector<int> iparam (11, 0);
538 
539  iparam[0] = 1; // shift strategy
540 
541  // maximum number of iterations
542  iparam[2] = control().max_steps();
543 
544  // Set the mode of dsaupd. 1 is exact shifting, 2 is user-supplied shifts,
545  // 3 is shift-invert mode, 4 is buckling mode, 5 is Cayley mode.
546 
547  iparam[6] = mode;
548  std::vector<int> ipntr (14, 0);
549 
550  // work arrays for ARPACK
551  std::vector<double> workd (3*n, 0.);
552  int lworkl = additional_data.symmetric ? ncv*ncv + 8*ncv : 3*ncv*ncv+6*ncv;
553  std::vector<double> workl (lworkl, 0.);
554 
555  //information out of the iteration
556  int info = 1;
557 
558  while (ido != 99)
559  {
560  // call of ARPACK dsaupd/dnaupd routine
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);
565  else
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);
569 
570  if (ido == 99)
571  break;
572 
573  switch (mode)
574  {
575  case 3:
576  {
577  switch (ido)
578  {
579  case -1:
580  {
581 
582  VectorType src,dst,tmp;
583  src.reinit(eigenvectors[0]);
584  dst.reinit(src);
585  tmp.reinit(src);
586 
587 
588  for (size_type i=0; i<src.size(); ++i)
589  src(i) = workd[ipntr[0]-1+i];
590 
591  // multiplication with mass matrix M
592  mass_matrix.vmult(tmp, src);
593  // solving linear system
594  inverse.vmult(dst,tmp);
595 
596  for (size_type i=0; i<dst.size(); ++i)
597  workd[ipntr[1]-1+i] = dst(i);
598  }
599  break;
600 
601  case 1:
602  {
603 
604  VectorType src,dst,tmp, tmp2;
605  src.reinit(eigenvectors[0]);
606  dst.reinit(src);
607  tmp.reinit(src);
608  tmp2.reinit(src);
609 
610  for (size_type i=0; i<src.size(); ++i)
611  {
612  src(i) = workd[ipntr[2]-1+i];
613  tmp(i) = workd[ipntr[0]-1+i];
614  }
615  // solving linear system
616  inverse.vmult(dst,src);
617 
618  for (size_type i=0; i<dst.size(); ++i)
619  workd[ipntr[1]-1+i] = dst(i);
620  }
621  break;
622 
623  case 2:
624  {
625 
626  VectorType src,dst;
627  src.reinit(eigenvectors[0]);
628  dst.reinit(src);
629 
630  for (size_type i=0; i<src.size(); ++i)
631  src(i) = workd[ipntr[0]-1+i];
632 
633  // Multiplication with mass matrix M
634  mass_matrix.vmult(dst, src);
635 
636  for (size_type i=0; i<dst.size(); ++i)
637  workd[ipntr[1]-1+i] = dst(i);
638 
639  }
640  break;
641 
642  default:
643  Assert (false, ArpackExcArpackIdo(ido));
644  break;
645  }
646  }
647  break;
648  default:
649  Assert (false, ArpackExcArpackMode(mode));
650  break;
651  }
652  }
653 
654  if (info<0)
655  {
657  {
658  Assert (false, ArpackExcArpackInfodsaupd(info));
659  }
660  else
661  Assert (false, ArpackExcArpackInfodnaupd(info));
662  }
663  else
664  {
665  // 1 - compute eigenvectors, 0 - only eigenvalues
666  int rvec = 1;
667 
668  // which eigenvectors
669  char howmany = 'A';
670 
671  std::vector<int> select (ncv, 1);
672 
673  int ldz = n;
674 
675  std::vector<double> eigenvalues_real (nev+1, 0.);
676  std::vector<double> eigenvalues_im (nev+1, 0.);
677 
678  // call of ARPACK dseupd/dneupd routine
680  {
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);
686  }
687  else
688  {
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);
695  }
696 
697  if (info == 1)
698  {
699  Assert (false, ArpackExcArpackInfoMaxIt(control().max_steps()));
700  }
701  else if (info == 3)
702  {
703  Assert (false, ArpackExcArpackNoShifts());
704  }
705  else if (info!=0)
706  {
708  {
709  Assert (false, ArpackExcArpackInfodseupd(info));
710  }
711  else
712  Assert (false, ArpackExcArpackInfodneupd(info));
713  }
714 
715  for (unsigned int i=0; i<nev; ++i)
716  for (unsigned int j=0; j<n; ++j)
717  eigenvectors[i](j) = v[i*n+j];
718 
719  for (unsigned int i=0; i<nev_const; ++i)
720  eigenvalues[i] = std::complex<double> (eigenvalues_real[i],
721  eigenvalues_im[i]);
722  }
723 }
724 
725 
726 inline
728 {
729  return solver_control;
730 }
731 
732 DEAL_II_NAMESPACE_CLOSE
733 
734 
735 #endif
736 #endif
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)
Definition: exceptions.h:358
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)
Definition: exceptions.h:346
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ArpackExcInvalidNumberofArnoldiVectors(int arg1, int arg2)
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:335
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)
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()