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