Reference documentation for deal.II version 8.5.1
parpack_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__parpack_solver_h
17 #define dealii__parpack_solver_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/base/smartpointer.h>
21 #include <deal.II/base/memory_consumption.h>
22 #include <deal.II/lac/solver_control.h>
23 #include <deal.II/lac/vector.h>
24 #include <deal.II/base/index_set.h>
25 
26 #include <cstring>
27 
28 
29 #ifdef DEAL_II_ARPACK_WITH_PARPACK
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 extern "C" {
34 
35  // http://www.mathkeisan.com/usersguide/man/pdnaupd.html
36  void pdnaupd_(MPI_Fint *comm, int *ido, char *bmat, int *n, char *which,
37  int *nev, double *tol, double *resid, int *ncv,
38  double *v, int *nloc, int *iparam, int *ipntr,
39  double *workd, double *workl, int *lworkl,
40  int *info);
41 
42  // http://www.mathkeisan.com/usersguide/man/pdsaupd.html
43  void pdsaupd_(MPI_Fint *comm, int *ido, char *bmat, int *n, char *which,
44  int *nev, double *tol, double *resid, int *ncv,
45  double *v, int *nloc, int *iparam, int *ipntr,
46  double *workd, double *workl, int *lworkl,
47  int *info);
48 
49  // http://www.mathkeisan.com/usersguide/man/pdneupd.html
50  void pdneupd_(MPI_Fint *comm, int *rvec, char *howmany, int *select, double *d,
51  double *di, double *z, int *ldz, double *sigmar,
52  double *sigmai, double *workev, char *bmat, int *n, char *which,
53  int *nev, double *tol, double *resid, int *ncv,
54  double *v, int *nloc, int *iparam, int *ipntr,
55  double *workd, double *workl, int *lworkl, int *info);
56 
57  // http://www.mathkeisan.com/usersguide/man/pdseupd.html
58  void pdseupd_(MPI_Fint *comm, int *rvec, char *howmany, int *select, double *d,
59  double *z, int *ldz, double *sigmar,
60  char *bmat, int *n, char *which,
61  int *nev, double *tol, double *resid, int *ncv,
62  double *v, int *nloc, int *iparam, int *ipntr,
63  double *workd, double *workl, int *lworkl, int *info);
64 
65  // other resources:
66  // http://acts.nersc.gov/superlu/example5/pnslac.c.html
67  // https://github.com/phpisciuneri/tijo/blob/master/dvr_parpack.cpp
68 
69 }
70 
143 template <typename VectorType>
145 {
146 public:
151 
162  {
202  };
203 
207  template <typename MatrixType>
208  class Shift : public ::Subscriptor
209  {
210  public:
211 
215  Shift (const MatrixType &A,
216  const MatrixType &B,
217  const double sigma)
218  :
219  A(A),
220  B(B),
221  sigma(sigma)
222  {}
223 
227  void vmult (VectorType &dst, const VectorType &src) const
228  {
229  B.vmult(dst,src);
230  dst *= (-sigma);
231  A.vmult_add(dst,src);
232  }
233 
237  void Tvmult (VectorType &dst, const VectorType &src) const
238  {
239  B.Tvmult(dst,src);
240  dst *= (-sigma);
241  A.Tvmult_add(dst,src);
242  }
243 
244  private:
245  const MatrixType &A;
246  const MatrixType &B;
247  const double sigma;
248  };
249 
255  {
256  const unsigned int number_of_arnoldi_vectors;
257  const WhichEigenvalues eigenvalue_of_interest;
258  const bool symmetric;
260  const unsigned int number_of_arnoldi_vectors = 15,
261  const WhichEigenvalues eigenvalue_of_interest = largest_magnitude,
262  const bool symmetric = false);
263  };
264 
268  SolverControl &control () const;
269 
274  const MPI_Comm &mpi_communicator,
275  const AdditionalData &data = AdditionalData());
276 
280  void reinit(const IndexSet &locally_owned_dofs );
281 
288  void reinit(const IndexSet &locally_owned_dofs,
289  const std::vector<IndexSet> &partitioning);
290 
294  void reinit(const VectorType &distributed_vector);
295 
299  void set_initial_vector(const VectorType &vec);
300 
306  void set_shift(const std::complex<double> sigma);
307 
313  template <typename MatrixType1,
314  typename MatrixType2, typename INVERSE>
315  void solve
316  (const MatrixType1 &A,
317  const MatrixType2 &B,
318  const INVERSE &inverse,
319  std::vector<std::complex<double> > &eigenvalues,
320  std::vector<VectorType> &eigenvectors,
321  const unsigned int n_eigenvalues);
322 
323  std::size_t memory_consumption() const;
324 
325 protected:
326 
332 
337 
338  // keep MPI communicator non-const as Arpack functions are not const either:
339 
344 
349 
350  // C++98 guarantees that the elements of a vector are stored contiguously
351 
355  int lworkl;
356 
360  std::vector<double> workl;
361 
365  std::vector<double> workd;
366 
370  int nloc;
371 
375  int ncv;
376 
377 
381  int ldv;
382 
387  std::vector<double> v;
388 
393 
398  std::vector<double> resid;
399 
403  int ldz;
404 
410  std::vector<double> z;
411 
415  int lworkev;
416 
420  std::vector<double> workev;
421 
425  std::vector<int> select;
426 
430  VectorType src,dst,tmp;
431 
435  std::vector< types::global_dof_index > local_indices;
436 
440  double sigmar;
441 
445  double sigmai;
446 
447 private:
448 
455  void internal_reinit(const IndexSet &locally_owned_dofs);
456 
461  << arg1 << " eigenpairs were requested, but only "
462  << arg2 << " converged");
463 
465  << "Number of wanted eigenvalues " << arg1
466  << " is larger that the size of the matrix " << arg2);
467 
469  << "Number of wanted eigenvalues " << arg1
470  << " is larger that the size of eigenvectors " << arg2);
471 
473  << "To store the real and complex parts of " << arg1
474  << " eigenvectors in real-valued vectors, their size (currently set to " << arg2
475  << ") should be greater than or equal to " << arg1+1);
476 
478  << "Number of wanted eigenvalues " << arg1
479  << " is larger that the size of eigenvalues " << arg2);
480 
482  << "Number of Arnoldi vectors " << arg1
483  << " is larger that the size of the matrix " << arg2);
484 
486  << "Number of Arnoldi vectors " << arg1
487  << " is too small to obtain " << arg2
488  << " eigenvalues");
489 
490  DeclException1 (PArpackExcIdo, int, << "This ido " << arg1
491  << " is not supported. Check documentation of ARPACK");
492 
493  DeclException1 (PArpackExcMode, int, << "This mode " << arg1
494  << " is not supported. Check documentation of ARPACK");
495 
497  << "Error with Pdnaupd, info " << arg1
498  << ". Check documentation of ARPACK");
499 
501  << "Error with Pdneupd, info " << arg1
502  << ". Check documentation of ARPACK");
503 
505  << "Maximum number " << arg1
506  << " of iterations reached.");
507 
509  << "No shifts could be applied during implicit"
510  << " Arnoldi update, try increasing the number of"
511  << " Arnoldi vectors.");
512 };
513 
514 template <typename VectorType>
515 std::size_t
517 {
518  return MemoryConsumption::memory_consumption (double()) *
519  (workl.size() +
520  workd.size() +
521  v.size() +
522  resid.size() +
523  z.size() +
524  workev.size() ) +
525  src.memory_consumption() +
526  dst.memory_consumption() +
527  tmp.memory_consumption() +
529 }
530 
531 template <typename VectorType>
533 AdditionalData (const unsigned int number_of_arnoldi_vectors,
534  const WhichEigenvalues eigenvalue_of_interest,
535  const bool symmetric)
536  :
537  number_of_arnoldi_vectors(number_of_arnoldi_vectors),
538  eigenvalue_of_interest(eigenvalue_of_interest),
539  symmetric(symmetric)
540 {
541  //Check for possible options for symmetric problems
542  if (symmetric)
543  {
544  Assert(eigenvalue_of_interest!=largest_real_part,
545  ExcMessage("'largest real part' can only be used for non-symmetric problems!"));
546  Assert(eigenvalue_of_interest!=smallest_real_part,
547  ExcMessage("'smallest real part' can only be used for non-symmetric problems!"));
548  Assert(eigenvalue_of_interest!=largest_imaginary_part,
549  ExcMessage("'largest imaginary part' can only be used for non-symmetric problems!"));
550  Assert(eigenvalue_of_interest!=smallest_imaginary_part,
551  ExcMessage("'smallest imaginary part' can only be used for non-symmetric problems!"));
552  }
553 }
554 
555 template <typename VectorType>
557  const MPI_Comm &mpi_communicator,
558  const AdditionalData &data)
559  :
560  solver_control (control),
561  additional_data (data),
562  mpi_communicator( mpi_communicator ),
563  mpi_communicator_fortran ( MPI_Comm_c2f( mpi_communicator ) ),
564  lworkl(0),
565  nloc(0),
566  ncv(0),
567  ldv(0),
568  initial_vector_provided(false),
569  ldz(0),
570  lworkev(0),
571  sigmar(0.0),
572  sigmai(0.0)
573 {}
574 
575 template <typename VectorType>
576 void PArpackSolver<VectorType>::set_shift(const std::complex<double> sigma)
577 {
578  sigmar = sigma.real();
579  sigmai = sigma.imag();
580 }
581 
582 template <typename VectorType>
584 set_initial_vector(const VectorType &vec)
585 {
586  initial_vector_provided = true;
587  Assert (resid.size() == local_indices.size(),
588  ExcDimensionMismatch(resid.size(),local_indices.size()));
589  vec.extract_subvector_to (local_indices.begin(),
590  local_indices.end(),
591  &resid[0]);
592 }
593 
594 
595 template <typename VectorType>
597 internal_reinit(const IndexSet &locally_owned_dofs)
598 {
599  // store local indices to write to vectors
600  locally_owned_dofs.fill_index_vector(local_indices);
601 
602  // scalars
603  nloc = locally_owned_dofs.n_elements ();
604  ncv = additional_data.number_of_arnoldi_vectors;
605 
606  Assert ((int)local_indices.size() == nloc, ExcInternalError() );
607 
608  // vectors
609  ldv = nloc;
610  v.resize (ldv*ncv, 0.0);
611 
612  resid.resize(nloc, 1.0);
613 
614  // work arrays for ARPACK
615  workd.resize(3*nloc,0.0);
616 
617  lworkl = additional_data.symmetric ?
618  ncv*ncv + 8*ncv
619  :
620  3*ncv*ncv+6*ncv;
621  workl.resize (lworkl, 0.);
622 
623  ldz = nloc;
624  z.resize (ldz*ncv, 0.); // TODO we actually need only ldz*nev
625 
626  // WORKEV Double precision work array of dimension 3*NCV.
627  lworkev = additional_data.symmetric ?
628  0 /*not used in symmetric case*/
629  :
630  3*ncv;
631  workev.resize (lworkev, 0.);
632 
633  select.resize (ncv, 0);
634 }
635 
636 template <typename VectorType>
637 void PArpackSolver<VectorType>::reinit(const IndexSet &locally_owned_dofs)
638 {
639  internal_reinit(locally_owned_dofs);
640 
641  // deal.II vectors:
642  src.reinit (locally_owned_dofs,mpi_communicator);
643  dst.reinit (locally_owned_dofs,mpi_communicator);
644  tmp.reinit (locally_owned_dofs,mpi_communicator);
645 }
646 
647 template <typename VectorType>
648 void PArpackSolver<VectorType>::reinit(const VectorType &distributed_vector)
649 {
650  internal_reinit(distributed_vector.locally_owned_elements());
651 
652  // deal.II vectors:
653  src.reinit (distributed_vector);
654  dst.reinit (distributed_vector);
655  tmp.reinit (distributed_vector);
656 }
657 
658 
659 template <typename VectorType>
660 void PArpackSolver<VectorType>::reinit(const IndexSet &locally_owned_dofs,
661  const std::vector<IndexSet> &partitioning)
662 {
663  internal_reinit(locally_owned_dofs);
664 
665  // deal.II vectors:
666  src.reinit (partitioning,mpi_communicator);
667  dst.reinit (partitioning,mpi_communicator);
668  tmp.reinit (partitioning,mpi_communicator);
669 }
670 
671 template <typename VectorType>
672 template <typename MatrixType1,typename MatrixType2, typename INVERSE>
674 (const MatrixType1 &/*system_matrix*/,
675  const MatrixType2 &mass_matrix,
676  const INVERSE &inverse,
677  std::vector<std::complex<double> > &eigenvalues,
678  std::vector<VectorType> &eigenvectors,
679  const unsigned int n_eigenvalues)
680 {
681 
682  if (additional_data.symmetric)
683  {
684  Assert (n_eigenvalues <= eigenvectors.size(),
685  PArpackExcInvalidEigenvectorSize(n_eigenvalues, eigenvectors.size()));
686  }
687  else
688  Assert (n_eigenvalues+1 <= eigenvectors.size(),
689  PArpackExcInvalidEigenvectorSizeNonsymmetric(n_eigenvalues, eigenvectors.size()));
690 
691  Assert (n_eigenvalues <= eigenvalues.size(),
692  PArpackExcInvalidEigenvalueSize(n_eigenvalues, eigenvalues.size()));
693 
694 
695  // use eigenvectors to get the problem size so that it is possible to
696  // employ LinearOperator for mass_matrix.
697  Assert (n_eigenvalues < eigenvectors[0].size(),
698  PArpackExcInvalidNumberofEigenvalues(n_eigenvalues, eigenvectors[0].size()));
699 
700  Assert (additional_data.number_of_arnoldi_vectors < eigenvectors[0].size(),
701  PArpackExcInvalidNumberofArnoldiVectors(
702  additional_data.number_of_arnoldi_vectors, eigenvectors[0].size()));
703 
704  Assert (additional_data.number_of_arnoldi_vectors > 2*n_eigenvalues+1,
705  PArpackExcSmallNumberofArnoldiVectors(
706  additional_data.number_of_arnoldi_vectors, n_eigenvalues));
707  // ARPACK mode for dnaupd, here only
708  // Mode 3: K*x = lambda*M*x, K symmetric, M symmetric positive semi-definite
709  //c ===> OP = (inv[K - sigma*M])*M and B = M.
710  //c ===> Shift-and-Invert mode
711  int mode = 3;
712 
713  // reverse communication parameter
714  // must be zero on the first call to pdnaupd
715  int ido = 0;
716 
717  // 'G' generalized eigenvalue problem
718  // 'I' standard eigenvalue problem
719  char bmat[2] = "G";
720 
721  // Specify the eigenvalues of interest, possible parameters:
722  // "LA" algebraically largest
723  // "SA" algebraically smallest
724  // "LM" largest magnitude
725  // "SM" smallest magnitude
726  // "LR" largest real part
727  // "SR" smallest real part
728  // "LI" largest imaginary part
729  // "SI" smallest imaginary part
730  // "BE" both ends of spectrum simultaneous
731  char which[3];
732  switch (additional_data.eigenvalue_of_interest)
733  {
734  case algebraically_largest:
735  std::strcpy (which, "LA");
736  break;
737  case algebraically_smallest:
738  std::strcpy (which, "SA");
739  break;
740  case largest_magnitude:
741  std::strcpy (which, "LM");
742  break;
743  case smallest_magnitude:
744  std::strcpy (which, "SM");
745  break;
746  case largest_real_part:
747  std::strcpy (which, "LR");
748  break;
749  case smallest_real_part:
750  std::strcpy (which, "SR");
751  break;
752  case largest_imaginary_part:
753  std::strcpy (which, "LI");
754  break;
755  case smallest_imaginary_part:
756  std::strcpy (which, "SI");
757  break;
758  case both_ends:
759  std::strcpy (which, "BE");
760  break;
761  }
762 
763  // tolerance for ARPACK
764  double tol = control().tolerance();
765 
766  //information to the routines
767  std::vector<int> iparam (11, 0);
768 
769  iparam[0] = 1;
770  // shift strategy: exact shifts with respect to the current Hessenberg matrix H.
771 
772  // maximum number of iterations
773  iparam[2] = control().max_steps();
774 
775  // Parpack currently works only for NB = 1
776  iparam[3] = 1;
777 
778  // Sets the mode of dsaupd:
779  // 1 is exact shifting,
780  // 2 is user-supplied shifts,
781  // 3 is shift-invert mode,
782  // 4 is buckling mode,
783  // 5 is Cayley mode.
784 
785  iparam[6] = mode;
786  std::vector<int> ipntr (14, 0);
787 
788  //information out of the iteration
789  // If INFO .EQ. 0, a random initial residual vector is used.
790  // If INFO .NE. 0, RESID contains the initial residual vector,
791  // possibly from a previous run.
792  // Typical choices in this situation might be to use the final value
793  // of the starting vector from the previous eigenvalue calculation
794  int info = initial_vector_provided? 1 : 0;
795 
796  // Number of eigenvalues of OP to be computed. 0 < NEV < N.
797  int nev = n_eigenvalues;
798  int n_inside_arpack = nloc;
799 
800  // IDO = 99: done
801  while (ido != 99)
802  {
803  // call of ARPACK pdnaupd routine
804  if (additional_data.symmetric)
805  pdsaupd_(&mpi_communicator_fortran,&ido, bmat, &n_inside_arpack, which, &nev, &tol,
806  &resid[0], &ncv, &v[0], &ldv, &iparam[0], &ipntr[0],
807  &workd[0], &workl[0], &lworkl, &info);
808  else
809  pdnaupd_(&mpi_communicator_fortran,&ido, bmat, &n_inside_arpack, which, &nev, &tol,
810  &resid[0], &ncv, &v[0], &ldv, &iparam[0], &ipntr[0],
811  &workd[0], &workl[0], &lworkl, &info);
812 
813  if (ido == 99)
814  break;
815 
816  switch (mode)
817  {
818 // OP = (inv[K - sigma*M])*M
819  case 3:
820  {
821  switch (ido)
822  {
823  case -1:
824  // compute Y = OP * X where
825  // IPNTR(1) is the pointer into WORKD for X,
826  // IPNTR(2) is the pointer into WORKD for Y.
827  {
828  const int shift_x = ipntr[0]-1;
829  const int shift_y = ipntr[1]-1;
830  Assert (shift_x>=0, ::ExcInternalError() );
831  Assert (shift_x+nloc <= (int)workd.size(), ::ExcInternalError() );
832  Assert (shift_y>=0, ::ExcInternalError() );
833  Assert (shift_y+nloc <= (int)workd.size(), ::ExcInternalError() );
834 
835  src = 0.0;
836  src.add (nloc,
837  &local_indices[0],
838  &workd[0]+shift_x );
839  src.compress (VectorOperation::add);
840 
841  // multiplication with mass matrix M
842  mass_matrix.vmult(tmp, src);
843  // solving linear system
844  inverse.vmult(dst,tmp);
845 
846  // store the result
847  dst.extract_subvector_to (local_indices.begin(),
848  local_indices.end(),
849  &workd[0]+shift_y );
850  }
851  break;
852 
853  case 1:
854  // compute Y = OP * X where
855  // IPNTR(1) is the pointer into WORKD for X,
856  // IPNTR(2) is the pointer into WORKD for Y.
857  // In mode 3,4 and 5, the vector B * X is already
858  // available in WORKD(ipntr(3)). It does not
859  // need to be recomputed in forming OP * X.
860  {
861  const int shift_x = ipntr[0]-1;
862  const int shift_y = ipntr[1]-1;
863  const int shift_b_x = ipntr[2]-1;
864 
865  Assert (shift_x>=0, ::ExcInternalError() );
866  Assert (shift_x+nloc <= (int)workd.size(), ::ExcInternalError() );
867  Assert (shift_y>=0, ::ExcInternalError() );
868  Assert (shift_y+nloc <= (int)workd.size(), ::ExcInternalError() );
869  Assert (shift_b_x>=0, ::ExcInternalError() );
870  Assert (shift_b_x+nloc <= (int)workd.size(), ::ExcInternalError() );
871  Assert (shift_y>=0, ::ExcInternalError() );
872  Assert (shift_y+nloc <= (int)workd.size(), ::ExcInternalError() );
873 
874  src = 0.0; // B*X
875  src.add (nloc,
876  &local_indices[0],
877  &workd[0]+shift_b_x );
878 
879  tmp = 0.0; // X
880  tmp.add (nloc,
881  &local_indices[0],
882  &workd[0]+shift_x);
883 
884  src.compress (VectorOperation::add);
885  tmp.compress (VectorOperation::add);
886 
887  // solving linear system
888  inverse.vmult(dst,src);
889 
890  // store the result
891  dst.extract_subvector_to (local_indices.begin(),
892  local_indices.end(),
893  &workd[0]+shift_y );
894 
895  }
896  break;
897 
898  case 2:
899  // compute Y = B * X where
900  // IPNTR(1) is the pointer into WORKD for X,
901  // IPNTR(2) is the pointer into WORKD for Y.
902  {
903 
904  const int shift_x = ipntr[0]-1;
905  const int shift_y = ipntr[1]-1;
906  Assert (shift_x>=0, ::ExcInternalError() );
907  Assert (shift_x+nloc <= (int)workd.size(), ::ExcInternalError() );
908  Assert (shift_y>=0, ::ExcInternalError() );
909  Assert (shift_y+nloc <= (int)workd.size(), ::ExcInternalError() );
910 
911  src = 0.0;
912  src.add (nloc,
913  &local_indices[0],
914  &workd[0]+shift_x );
915  src.compress (VectorOperation::add);
916 
917  // Multiplication with mass matrix M
918  mass_matrix.vmult(dst, src);
919 
920  // store the result
921  dst.extract_subvector_to (local_indices.begin(),
922  local_indices.end(),
923  &workd[0]+shift_y);
924 
925  }
926  break;
927 
928  default:
929  AssertThrow (false, PArpackExcIdo(ido));
930  break;
931  }
932  }
933  break;
934  default:
935  AssertThrow (false, PArpackExcMode(mode));
936  break;
937  }
938  }
939 
940  if (info<0)
941  {
942  AssertThrow (false, PArpackExcInfoPdnaupd(info));
943  }
944  else
945  {
946  // 1 - compute eigenvectors,
947  // 0 - only eigenvalues
948  int rvec = 1;
949 
950  // which eigenvectors
951  char howmany[4] = "All";
952 
953  std::vector<double> eigenvalues_real (n_eigenvalues+1, 0.);
954  std::vector<double> eigenvalues_im (n_eigenvalues+1, 0.);
955 
956  // call of ARPACK pdneupd routine
957  if (additional_data.symmetric)
958  pdseupd_(&mpi_communicator_fortran, &rvec, howmany, &select[0], &eigenvalues_real[0],
959  &z[0], &ldz, &sigmar,
960  bmat, &n_inside_arpack, which, &nev, &tol,
961  &resid[0], &ncv, &v[0], &ldv,
962  &iparam[0], &ipntr[0], &workd[0], &workl[0], &lworkl, &info);
963  else
964  pdneupd_(&mpi_communicator_fortran, &rvec, howmany, &select[0], &eigenvalues_real[0],
965  &eigenvalues_im[0], &v[0], &ldz, &sigmar, &sigmai,
966  &workev[0], bmat, &n_inside_arpack, which, &nev, &tol,
967  &resid[0], &ncv, &v[0], &ldv,
968  &iparam[0], &ipntr[0], &workd[0], &workl[0], &lworkl, &info);
969 
970  if (info == 1)
971  {
972  AssertThrow (false, PArpackExcInfoMaxIt(control().max_steps()));
973  }
974  else if (info == 3)
975  {
976  AssertThrow (false, PArpackExcNoShifts(1));
977  }
978  else if (info!=0)
979  {
980  AssertThrow (false, PArpackExcInfoPdneupd(info));
981  }
982 
983  for (int i=0; i<nev; ++i)
984  {
985  eigenvectors[i] = 0.0;
986  Assert (i*nloc + nloc <= (int)v.size(), ::ExcInternalError() );
987 
988  eigenvectors[i].add (nloc,
989  &local_indices[0],
990  &v[i*nloc] );
991  eigenvectors[i].compress (VectorOperation::add);
992  }
993 
994  for (size_type i=0; i<n_eigenvalues; ++i)
995  eigenvalues[i] = std::complex<double> (eigenvalues_real[i],
996  eigenvalues_im[i]);
997  }
998 
999  // Throw an error if the solver did not converge.
1000  AssertThrow (iparam[4] >= (int)n_eigenvalues,
1001  PArpackExcConvergedEigenvectors(n_eigenvalues,iparam[4]));
1002 
1003  // both PDNAUPD and PDSAUPD compute eigenpairs of inv[A - sigma*M]*M
1004  // with respect to a semi-inner product defined by M.
1005 
1006  // resid likely contains residual with respect to M-norm.
1007  {
1008 
1009  tmp = 0.0;
1010  tmp.add (nloc,
1011  &local_indices[0],
1012  &resid[0]);
1013  solver_control.check ( iparam[2], tmp.l2_norm() );
1014  }
1015 
1016 
1017 }
1018 
1019 template <typename VectorType>
1021 {
1022  return solver_control;
1023 }
1024 
1025 DEAL_II_NAMESPACE_CLOSE
1026 
1027 
1028 #endif
1029 #endif
void Tvmult(VectorType &dst, const VectorType &src) const
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
static ::ExceptionBase & PArpackExcSmallNumberofArnoldiVectors(int arg1, int arg2)
static ::ExceptionBase & PArpackExcInvalidEigenvectorSize(int arg1, int arg2)
std::vector< double > z
std::vector< double > workev
static ::ExceptionBase & PArpackExcInvalidNumberofEigenvalues(int arg1, int arg2)
static ::ExceptionBase & PArpackExcMode(int arg1)
PArpackSolver(SolverControl &control, const MPI_Comm &mpi_communicator, const AdditionalData &data=AdditionalData())
MPI_Fint mpi_communicator_fortran
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
const AdditionalData additional_data
void set_shift(const std::complex< double > sigma)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & PArpackExcNoShifts(int arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:564
bool initial_vector_provided
unsigned int global_dof_index
Definition: types.h:88
SolverControl & control() const
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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)
std::vector< int > select
static ::ExceptionBase & PArpackExcConvergedEigenvectors(int arg1, int arg2)
static ::ExceptionBase & PArpackExcInfoPdnaupd(int arg1)
VectorType src
static ::ExceptionBase & PArpackExcInfoMaxIt(int arg1)
MPI_Comm mpi_communicator
static ::ExceptionBase & PArpackExcInvalidEigenvalueSize(int arg1, int arg2)
void fill_index_vector(std::vector< size_type > &indices) const
Definition: index_set.cc:512
types::global_dof_index size_type
std::vector< double > v
std::vector< double > resid
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
void vmult(VectorType &dst, const VectorType &src) const
static ::ExceptionBase & PArpackExcIdo(int arg1)
void internal_reinit(const IndexSet &locally_owned_dofs)
void set_initial_vector(const VectorType &vec)
SolverControl & solver_control
Shift(const MatrixType &A, const MatrixType &B, const double sigma)
void reinit(const IndexSet &locally_owned_dofs)
std::vector< types::global_dof_index > local_indices
static ::ExceptionBase & PArpackExcInvalidNumberofArnoldiVectors(int arg1, int arg2)
std::vector< double > workl
size_type n_elements() const
Definition: index_set.h:1560
static ::ExceptionBase & PArpackExcInvalidEigenvectorSizeNonsymmetric(int arg1, int arg2)
std::vector< double > workd
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & PArpackExcInfoPdneupd(int arg1)