16#ifndef dealii_parpack_solver_h
17#define dealii_parpack_solver_h
31#ifdef DEAL_II_ARPACK_WITH_PARPACK
210template <
typename VectorType>
315 const std::vector<IndexSet> &partitioning);
321 reinit(
const VectorType &distributed_vector);
338 set_shift(
const std::complex<double> sigma);
349 template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
351 solve(
const MatrixType1 & A,
352 const MatrixType2 & B,
353 const INVERSE & inverse,
356 const unsigned int n_eigenvalues);
361 template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
363 solve(
const MatrixType1 & A,
364 const MatrixType2 & B,
365 const INVERSE & inverse,
368 const unsigned int n_eigenvalues);
437 std::vector<double>
v;
460 std::vector<double>
z;
513 << arg1 <<
" eigenpairs were requested, but only " << arg2
519 <<
"Number of wanted eigenvalues " << arg1
520 <<
" is larger that the size of the matrix " << arg2);
525 <<
"Number of wanted eigenvalues " << arg1
526 <<
" is larger that the size of eigenvectors " << arg2);
532 <<
"To store the real and complex parts of " << arg1
533 <<
" eigenvectors in real-valued vectors, their size (currently set to "
534 << arg2 <<
") should be greater than or equal to " << arg1 + 1);
539 <<
"Number of wanted eigenvalues " << arg1
540 <<
" is larger that the size of eigenvalues " << arg2);
545 <<
"Number of Arnoldi vectors " << arg1
546 <<
" is larger that the size of the matrix " << arg2);
551 <<
"Number of Arnoldi vectors " << arg1
552 <<
" is too small to obtain " << arg2 <<
" eigenvalues");
556 <<
"This ido " << arg1
557 <<
" is not supported. Check documentation of ARPACK");
561 <<
"This mode " << arg1
562 <<
" is not supported. Check documentation of ARPACK");
566 <<
"Error with Pdnaupd, info " << arg1
567 <<
". Check documentation of ARPACK");
571 <<
"Error with Pdneupd, info " << arg1
572 <<
". Check documentation of ARPACK");
576 <<
"Maximum number " << arg1 <<
" of iterations reached.");
580 <<
"No shifts could be applied during implicit"
581 <<
" Arnoldi update, try increasing the number of"
582 <<
" Arnoldi vectors.");
587template <
typename VectorType>
592 (workl.size() + workd.size() + v.size() + resid.size() + z.size() +
594 src.memory_consumption() + dst.memory_consumption() +
595 tmp.memory_consumption() +
597 local_indices.size();
602template <
typename VectorType>
604 const unsigned int number_of_arnoldi_vectors,
606 const bool symmetric,
608 : number_of_arnoldi_vectors(number_of_arnoldi_vectors)
609 , eigenvalue_of_interest(eigenvalue_of_interest)
610 , symmetric(symmetric)
619 "'largest real part' can only be used for non-symmetric problems!"));
623 "'smallest real part' can only be used for non-symmetric problems!"));
627 "'largest imaginary part' can only be used for non-symmetric problems!"));
631 "'smallest imaginary part' can only be used for non-symmetric problems!"));
634 ExcMessage(
"Currently, only modes 1, 2 and 3 are supported."));
639template <
typename VectorType>
660template <
typename VectorType>
664 sigmar = sigma.real();
665 sigmai = sigma.imag();
670template <
typename VectorType>
674 initial_vector_provided =
true;
675 Assert(resid.size() == local_indices.size(),
677 vec.extract_subvector_to(local_indices.begin(),
684template <
typename VectorType>
693 ncv = additional_data.number_of_arnoldi_vectors;
699 v.resize(ldv * ncv, 0.0);
701 resid.resize(nloc, 1.0);
704 workd.resize(3 * nloc, 0.0);
707 additional_data.symmetric ? ncv * ncv + 8 * ncv : 3 * ncv * ncv + 6 * ncv;
708 workl.resize(lworkl, 0.);
711 z.resize(ldz * ncv, 0.);
714 lworkev = additional_data.symmetric ? 0
717 workev.resize(lworkev, 0.);
719 select.resize(ncv, 0);
724template <
typename VectorType>
728 internal_reinit(locally_owned_dofs);
731 src.reinit(locally_owned_dofs, mpi_communicator);
732 dst.reinit(locally_owned_dofs, mpi_communicator);
733 tmp.reinit(locally_owned_dofs, mpi_communicator);
738template <
typename VectorType>
742 internal_reinit(distributed_vector.locally_owned_elements());
745 src.reinit(distributed_vector);
746 dst.reinit(distributed_vector);
747 tmp.reinit(distributed_vector);
752template <
typename VectorType>
755 const std::vector<IndexSet> &partitioning)
757 internal_reinit(locally_owned_dofs);
760 src.reinit(partitioning, mpi_communicator);
761 dst.reinit(partitioning, mpi_communicator);
762 tmp.reinit(partitioning, mpi_communicator);
767template <
typename VectorType>
768template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
771 const MatrixType2 & B,
772 const INVERSE & inverse,
775 const unsigned int n_eigenvalues)
777 std::vector<VectorType *> eigenvectors_ptr(
eigenvectors.size());
780 solve(A, B, inverse,
eigenvalues, eigenvectors_ptr, n_eigenvalues);
785template <
typename VectorType>
786template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
789 const MatrixType2 &mass_matrix,
790 const INVERSE & inverse,
793 const unsigned int n_eigenvalues)
795 if (additional_data.symmetric)
798 PArpackExcInvalidEigenvectorSize(n_eigenvalues,
803 PArpackExcInvalidEigenvectorSizeNonsymmetric(n_eigenvalues,
807 PArpackExcInvalidEigenvalueSize(n_eigenvalues,
eigenvalues.size()));
813 PArpackExcInvalidNumberofEigenvalues(n_eigenvalues,
817 PArpackExcInvalidNumberofArnoldiVectors(
818 additional_data.number_of_arnoldi_vectors,
eigenvectors[0]->size()));
820 Assert(additional_data.number_of_arnoldi_vectors > 2 * n_eigenvalues + 1,
821 PArpackExcSmallNumberofArnoldiVectors(
822 additional_data.number_of_arnoldi_vectors, n_eigenvalues));
824 int mode = additional_data.mode;
833 bmat[0] = (mode == 1) ?
'I' :
'G';
847 switch (additional_data.eigenvalue_of_interest)
849 case algebraically_largest:
850 std::strcpy(which,
"LA");
852 case algebraically_smallest:
853 std::strcpy(which,
"SA");
855 case largest_magnitude:
856 std::strcpy(which,
"LM");
858 case smallest_magnitude:
859 std::strcpy(which,
"SM");
861 case largest_real_part:
862 std::strcpy(which,
"LR");
864 case smallest_real_part:
865 std::strcpy(which,
"SR");
867 case largest_imaginary_part:
868 std::strcpy(which,
"LI");
870 case smallest_imaginary_part:
871 std::strcpy(which,
"SI");
874 std::strcpy(which,
"BE");
879 double tol = control().tolerance();
882 std::vector<int> iparam(11, 0);
889 iparam[2] = control().max_steps();
902 std::vector<int> ipntr(14, 0);
910 int info = initial_vector_provided ? 1 : 0;
913 int nev = n_eigenvalues;
914 int n_inside_arpack = nloc;
920 if (additional_data.symmetric)
957 AssertThrow(info == 0, PArpackExcInfoPdnaupd(info));
965 const int shift_x = ipntr[0] - 1;
966 const int shift_y = ipntr[1] - 1;
968 Assert(shift_x + nloc <=
static_cast<int>(workd.size()),
971 Assert(shift_y + nloc <=
static_cast<int>(workd.size()),
977 if ((ido == -1) || (ido == 1 && mode < 3))
980 src.add(nloc, local_indices.data(), workd.data() + shift_x);
986 mass_matrix.vmult(tmp, src);
987 inverse.vmult(dst, tmp);
992 system_matrix.vmult(tmp, src);
994 tmp.extract_subvector_to(local_indices.begin(),
996 workd.data() + shift_x);
997 inverse.vmult(dst, tmp);
1001 system_matrix.vmult(dst, src);
1006 else if (ido == 1 && mode >= 3)
1010 const int shift_b_x = ipntr[2] - 1;
1012 Assert(shift_b_x + nloc <=
static_cast<int>(workd.size()),
1016 src.add(nloc, local_indices.data(), workd.data() + shift_b_x);
1021 inverse.vmult(dst, src);
1026 src.add(nloc, local_indices.data(), workd.data() + shift_x);
1037 mass_matrix.vmult(dst, src);
1046 dst.extract_subvector_to(local_indices.begin(),
1047 local_indices.end(),
1048 workd.data() + shift_y);
1056 char howmany[4] =
"All";
1058 std::vector<double> eigenvalues_real(n_eigenvalues + 1, 0.);
1059 std::vector<double> eigenvalues_im(n_eigenvalues + 1, 0.);
1062 if (additional_data.symmetric)
1063 pdseupd_(&mpi_communicator_fortran,
1067 eigenvalues_real.data(),
1087 pdneupd_(&mpi_communicator_fortran,
1091 eigenvalues_real.data(),
1092 eigenvalues_im.data(),
1116 AssertThrow(
false, PArpackExcInfoMaxIt(control().max_steps()));
1127 for (
int i = 0; i < nev; ++i)
1132 eigenvectors[i]->add(nloc, local_indices.data(), &v[i * nloc]);
1136 for (
size_type i = 0; i < n_eigenvalues; ++i)
1138 std::complex<double>(eigenvalues_real[i], eigenvalues_im[i]);
1141 AssertThrow(iparam[4] >=
static_cast<int>(n_eigenvalues),
1142 PArpackExcConvergedEigenvectors(n_eigenvalues, iparam[4]));
1150 tmp.add(nloc, local_indices.data(), resid.data());
1152 solver_control.check(iparam[2], tmp.l2_norm());
1158template <
typename VectorType>
1162 return solver_control;
size_type n_elements() const
std::vector< size_type > get_index_vector() const
void set_initial_vector(const VectorType &vec)
SolverControl & solver_control
MPI_Comm mpi_communicator
SolverControl & control() const
@ smallest_imaginary_part
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)
PArpackSolver(SolverControl &control, const MPI_Comm mpi_communicator, const AdditionalData &data=AdditionalData())
std::vector< types::global_dof_index > local_indices
std::vector< double > workl
std::vector< int > select
std::size_t memory_consumption() const
MPI_Fint mpi_communicator_fortran
void set_shift(const std::complex< double > sigma)
void reinit(const IndexSet &locally_owned_dofs)
std::vector< double > workd
const AdditionalData additional_data
bool initial_vector_provided
void internal_reinit(const IndexSet &locally_owned_dofs)
std::vector< double > resid
std::vector< double > workev
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & PArpackExcInvalidNumberofEigenvalues(int arg1, int arg2)
static ::ExceptionBase & PArpackExcInvalidEigenvalueSize(int arg1, int arg2)
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & PArpackExcNoShifts(int arg1)
static ::ExceptionBase & PArpackExcInvalidEigenvectorSize(int arg1, int arg2)
static ::ExceptionBase & PArpackExcSmallNumberofArnoldiVectors(int arg1, int arg2)
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & PArpackExcInfoPdneupd(int arg1)
#define AssertIndexRange(index, range)
static ::ExceptionBase & PArpackExcIdo(int arg1)
static ::ExceptionBase & PArpackExcConvergedEigenvectors(int arg1, int arg2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & PArpackExcInfoMaxIt(int arg1)
static ::ExceptionBase & PArpackExcMode(int arg1)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & PArpackExcInfoPdnaupd(int arg1)
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & PArpackExcInvalidEigenvectorSizeNonsymmetric(int arg1, int arg2)
static ::ExceptionBase & PArpackExcInvalidNumberofArnoldiVectors(int arg1, int arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
std::enable_if_t< std::is_fundamental< T >::value, std::size_t > memory_consumption(const T &t)
unsigned int global_dof_index
void pdneupd_(MPI_Fint *comm, int *rvec, char *howmany, int *select, double *d, double *di, double *z, int *ldz, double *sigmar, double *sigmai, double *workev, char *bmat, int *n, char *which, int *nev, double *tol, double *resid, int *ncv, double *v, int *nloc, int *iparam, int *ipntr, double *workd, double *workl, int *lworkl, int *info)
void pdsaupd_(MPI_Fint *comm, int *ido, char *bmat, int *n, char *which, int *nev, double *tol, double *resid, int *ncv, double *v, int *nloc, int *iparam, int *ipntr, double *workd, double *workl, int *lworkl, int *info)
void pdnaupd_(MPI_Fint *comm, int *ido, char *bmat, int *n, char *which, int *nev, double *tol, double *resid, int *ncv, double *v, int *nloc, int *iparam, int *ipntr, double *workd, double *workl, int *lworkl, int *info)
void pdseupd_(MPI_Fint *comm, int *rvec, char *howmany, int *select, double *d, double *z, int *ldz, double *sigmar, char *bmat, int *n, char *which, int *nev, double *tol, double *resid, int *ncv, double *v, int *nloc, int *iparam, int *ipntr, double *workd, double *workl, int *lworkl, int *info)
AdditionalData(const unsigned int number_of_arnoldi_vectors=15, const WhichEigenvalues eigenvalue_of_interest=largest_magnitude, const bool symmetric=false, const int mode=3)
const unsigned int number_of_arnoldi_vectors
const WhichEigenvalues eigenvalue_of_interest
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)