16 #ifndef dealii_parpack_solver_h 17 #define dealii_parpack_solver_h 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> 24 #include <deal.II/lac/solver_control.h> 25 #include <deal.II/lac/vector_operation.h> 31 #ifdef DEAL_II_ARPACK_WITH_PARPACK 33 DEAL_II_NAMESPACE_OPEN
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,
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,
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);
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);
148 template <
typename VectorType>
217 template <
typename MatrixType>
237 void vmult (VectorType &dst,
const VectorType &
src)
const 241 A.vmult_add(dst,
src);
247 void Tvmult (VectorType &dst,
const VectorType &
src)
const 251 A.Tvmult_add(dst,
src);
267 const unsigned int number_of_arnoldi_vectors;
269 const bool symmetric;
272 const unsigned int number_of_arnoldi_vectors = 15,
274 const bool symmetric =
false,
302 const std::vector<IndexSet> &partitioning);
307 void reinit(
const VectorType &distributed_vector);
322 void set_shift(
const std::complex<double> sigma);
333 template <
typename MatrixType1,
334 typename MatrixType2,
typename INVERSE>
336 (
const MatrixType1 &A,
337 const MatrixType2 &B,
338 const INVERSE &inverse,
341 const unsigned int n_eigenvalues);
346 template <
typename MatrixType1,
347 typename MatrixType2,
typename INVERSE>
349 (
const MatrixType1 &A,
350 const MatrixType2 &B,
351 const INVERSE &inverse,
354 const unsigned int n_eigenvalues);
423 std::vector<double>
v;
446 std::vector<double>
z;
497 << arg1 <<
" eigenpairs were requested, but only " 498 << arg2 <<
" converged");
501 <<
"Number of wanted eigenvalues " << arg1
502 <<
" is larger that the size of the matrix " << arg2);
505 <<
"Number of wanted eigenvalues " << arg1
506 <<
" is larger that the size of eigenvectors " << arg2);
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);
514 <<
"Number of wanted eigenvalues " << arg1
515 <<
" is larger that the size of eigenvalues " << arg2);
518 <<
"Number of Arnoldi vectors " << arg1
519 <<
" is larger that the size of the matrix " << arg2);
522 <<
"Number of Arnoldi vectors " << arg1
523 <<
" is too small to obtain " << arg2
527 <<
" is not supported. Check documentation of ARPACK");
530 <<
" is not supported. Check documentation of ARPACK");
533 <<
"Error with Pdnaupd, info " << arg1
534 <<
". Check documentation of ARPACK");
537 <<
"Error with Pdneupd, info " << arg1
538 <<
". Check documentation of ARPACK");
541 <<
"Maximum number " << arg1
542 <<
" of iterations reached.");
545 <<
"No shifts could be applied during implicit" 546 <<
" Arnoldi update, try increasing the number of" 547 <<
" Arnoldi vectors.");
552 template <
typename VectorType>
563 src.memory_consumption() +
564 dst.memory_consumption() +
565 tmp.memory_consumption() +
571 template <
typename VectorType>
574 const WhichEigenvalues eigenvalue_of_interest,
575 const bool symmetric,
578 number_of_arnoldi_vectors(number_of_arnoldi_vectors),
579 eigenvalue_of_interest(eigenvalue_of_interest),
580 symmetric(symmetric),
587 ExcMessage(
"'largest real part' can only be used for non-symmetric problems!"));
589 ExcMessage(
"'smallest real part' can only be used for non-symmetric problems!"));
591 ExcMessage(
"'largest imaginary part' can only be used for non-symmetric problems!"));
593 ExcMessage(
"'smallest imaginary part' can only be used for non-symmetric problems!"));
595 Assert (mode >= 1 && mode <= 3,
596 ExcMessage(
"Currently, only modes 1, 2 and 3 are supported."));
601 template <
typename VectorType>
603 const MPI_Comm &mpi_communicator,
606 solver_control (control),
607 additional_data (data),
608 mpi_communicator( mpi_communicator ),
609 mpi_communicator_fortran ( MPI_Comm_c2f( mpi_communicator ) ),
614 initial_vector_provided(false),
623 template <
typename VectorType>
626 sigmar = sigma.real();
627 sigmai = sigma.imag();
632 template <
typename VectorType>
636 initial_vector_provided =
true;
637 Assert (resid.size() == local_indices.size(),
639 vec.extract_subvector_to (local_indices.begin(),
646 template <
typename VectorType>
655 ncv = additional_data.number_of_arnoldi_vectors;
661 v.resize (ldv*ncv, 0.0);
663 resid.resize(nloc, 1.0);
666 workd.resize(3*nloc,0.0);
668 lworkl = additional_data.symmetric ?
672 workl.resize (lworkl, 0.);
675 z.resize (ldz*ncv, 0.);
678 lworkev = additional_data.symmetric ?
682 workev.resize (lworkev, 0.);
684 select.resize (ncv, 0);
689 template <
typename VectorType>
692 internal_reinit(locally_owned_dofs);
695 src.reinit (locally_owned_dofs,mpi_communicator);
696 dst.reinit (locally_owned_dofs,mpi_communicator);
697 tmp.reinit (locally_owned_dofs,mpi_communicator);
702 template <
typename VectorType>
705 internal_reinit(distributed_vector.locally_owned_elements());
708 src.reinit (distributed_vector);
709 dst.reinit (distributed_vector);
710 tmp.reinit (distributed_vector);
715 template <
typename VectorType>
717 const std::vector<IndexSet> &partitioning)
719 internal_reinit(locally_owned_dofs);
722 src.reinit (partitioning,mpi_communicator);
723 dst.reinit (partitioning,mpi_communicator);
724 tmp.reinit (partitioning,mpi_communicator);
729 template <
typename VectorType>
730 template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
732 (
const MatrixType1 &A,
733 const MatrixType2 &B,
734 const INVERSE &inverse,
737 const unsigned int n_eigenvalues)
739 std::vector<VectorType *> eigenvectors_ptr(
eigenvectors.size());
742 solve(A,B,inverse,
eigenvalues,eigenvectors_ptr,n_eigenvalues);
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,
755 const unsigned int n_eigenvalues)
758 if (additional_data.symmetric)
761 PArpackExcInvalidEigenvectorSize(n_eigenvalues,
eigenvectors.size()));
765 PArpackExcInvalidEigenvectorSizeNonsymmetric(n_eigenvalues,
eigenvectors.size()));
768 PArpackExcInvalidEigenvalueSize(n_eigenvalues,
eigenvalues.size()));
774 PArpackExcInvalidNumberofEigenvalues(n_eigenvalues,
eigenvectors[0]->size()));
777 PArpackExcInvalidNumberofArnoldiVectors(
778 additional_data.number_of_arnoldi_vectors,
eigenvectors[0]->size()));
780 Assert (additional_data.number_of_arnoldi_vectors > 2*n_eigenvalues+1,
781 PArpackExcSmallNumberofArnoldiVectors(
782 additional_data.number_of_arnoldi_vectors, n_eigenvalues));
784 int mode = additional_data.mode;
793 bmat[0] = (mode == 1 ) ?
'I' :
'G';
807 switch (additional_data.eigenvalue_of_interest)
809 case algebraically_largest:
810 std::strcpy (which,
"LA");
812 case algebraically_smallest:
813 std::strcpy (which,
"SA");
815 case largest_magnitude:
816 std::strcpy (which,
"LM");
818 case smallest_magnitude:
819 std::strcpy (which,
"SM");
821 case largest_real_part:
822 std::strcpy (which,
"LR");
824 case smallest_real_part:
825 std::strcpy (which,
"SR");
827 case largest_imaginary_part:
828 std::strcpy (which,
"LI");
830 case smallest_imaginary_part:
831 std::strcpy (which,
"SI");
834 std::strcpy (which,
"BE");
839 double tol = control().tolerance();
842 std::vector<int> iparam (11, 0);
848 iparam[2] = control().max_steps();
861 std::vector<int> ipntr (14, 0);
869 int info = initial_vector_provided? 1 : 0;
872 int nev = n_eigenvalues;
873 int n_inside_arpack = nloc;
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);
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);
888 AssertThrow (info == 0, PArpackExcInfoPdnaupd(info));
896 const int shift_x = ipntr[0]-1;
897 const int shift_y = ipntr[1]-1;
907 (ido == 1 && mode<3))
911 local_indices.data(),
912 workd.data()+shift_x );
918 mass_matrix.vmult(tmp, src);
919 inverse.vmult(dst,tmp);
924 system_matrix.vmult(tmp, src);
926 tmp.extract_subvector_to (local_indices.begin(),
928 workd.data()+shift_x);
929 inverse.vmult(dst,tmp);
933 system_matrix.vmult(dst, src);
939 else if (ido == 1 && mode >= 3)
943 const int shift_b_x = ipntr[2]-1;
949 local_indices.data(),
950 workd.data()+shift_b_x );
955 inverse.vmult(dst,src);
961 local_indices.data(),
962 workd.data()+shift_x );
973 mass_matrix.vmult(dst, src);
982 dst.extract_subvector_to (local_indices.begin(),
984 workd.data()+shift_y);
992 char howmany[4] =
"All";
994 std::vector<double> eigenvalues_real (n_eigenvalues+1, 0.);
995 std::vector<double> eigenvalues_im (n_eigenvalues+1, 0.);
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);
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);
1013 AssertThrow (
false, PArpackExcInfoMaxIt(control().max_steps()));
1024 for (
int i=0; i<nev; ++i)
1030 local_indices.data(),
1035 for (
size_type i=0; i<n_eigenvalues; ++i)
1036 eigenvalues[i] = std::complex<double> (eigenvalues_real[i],
1041 PArpackExcConvergedEigenvectors(n_eigenvalues,iparam[4]));
1051 local_indices.data(),
1053 solver_control.check ( iparam[2], tmp.l2_norm() );
1061 template <
typename VectorType>
1064 return solver_control;
1067 DEAL_II_NAMESPACE_CLOSE
void Tvmult(VectorType &dst, const VectorType &src) const
#define DeclException2(Exception2, type1, type2, outsequence)
static ::ExceptionBase & PArpackExcSmallNumberofArnoldiVectors(int arg1, int arg2)
static ::ExceptionBase & PArpackExcInvalidEigenvectorSize(int arg1, int arg2)
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)
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)
bool initial_vector_provided
unsigned int global_dof_index
SolverControl & control() const
#define Assert(cond, exc)
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)
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
types::global_dof_index size_type
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
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)