16 #ifndef dealii_parpack_solver_h 17 #define dealii_parpack_solver_h 31 #ifdef DEAL_II_ARPACK_WITH_PARPACK 212 template <
typename VectorType>
284 const unsigned int number_of_arnoldi_vectors = 15,
286 const bool symmetric =
false,
317 const std::vector<IndexSet> &partitioning);
340 set_shift(
const std::complex<double> sigma);
351 template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
353 solve(
const MatrixType1 &
A,
354 const MatrixType2 & B,
355 const INVERSE & inverse,
358 const unsigned int n_eigenvalues);
363 template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
365 solve(
const MatrixType1 &
A,
366 const MatrixType2 & B,
367 const INVERSE & inverse,
370 const unsigned int n_eigenvalues);
439 std::vector<double>
v;
462 std::vector<double>
z;
515 << arg1 <<
" eigenpairs were requested, but only " << arg2
521 <<
"Number of wanted eigenvalues " << arg1
522 <<
" is larger that the size of the matrix " << arg2);
527 <<
"Number of wanted eigenvalues " << arg1
528 <<
" is larger that the size of eigenvectors " << arg2);
534 <<
"To store the real and complex parts of " << arg1
535 <<
" eigenvectors in real-valued vectors, their size (currently set to " 536 << arg2 <<
") should be greater than or equal to " << arg1 + 1);
541 <<
"Number of wanted eigenvalues " << arg1
542 <<
" is larger that the size of eigenvalues " << arg2);
547 <<
"Number of Arnoldi vectors " << arg1
548 <<
" is larger that the size of the matrix " << arg2);
553 <<
"Number of Arnoldi vectors " << arg1
554 <<
" is too small to obtain " << arg2 <<
" eigenvalues");
558 <<
"This ido " << arg1
559 <<
" is not supported. Check documentation of ARPACK");
563 <<
"This mode " << arg1
564 <<
" is not supported. Check documentation of ARPACK");
568 <<
"Error with Pdnaupd, info " << arg1
569 <<
". Check documentation of ARPACK");
573 <<
"Error with Pdneupd, info " << arg1
574 <<
". Check documentation of ARPACK");
578 <<
"Maximum number " << arg1 <<
" of iterations reached.");
582 <<
"No shifts could be applied during implicit" 583 <<
" Arnoldi update, try increasing the number of" 584 <<
" Arnoldi vectors.");
589 template <
typename VectorType>
596 src.memory_consumption() +
dst.memory_consumption() +
597 tmp.memory_consumption() +
604 template <
typename VectorType>
610 : number_of_arnoldi_vectors(number_of_arnoldi_vectors)
611 , eigenvalue_of_interest(eigenvalue_of_interest)
612 , symmetric(symmetric)
621 "'largest real part' can only be used for non-symmetric problems!"));
625 "'smallest real part' can only be used for non-symmetric problems!"));
629 "'largest imaginary part' can only be used for non-symmetric problems!"));
633 "'smallest imaginary part' can only be used for non-symmetric problems!"));
635 Assert(mode >= 1 && mode <= 3,
636 ExcMessage(
"Currently, only modes 1, 2 and 3 are supported."));
641 template <
typename VectorType>
647 , mpi_communicator(mpi_communicator)
662 template <
typename VectorType>
672 template <
typename VectorType>
686 template <
typename VectorType>
701 v.resize(ldv *
ncv, 0.0);
713 z.resize(
ldz * ncv, 0.);
726 template <
typename VectorType>
740 template <
typename VectorType>
747 src.reinit(distributed_vector);
748 dst.reinit(distributed_vector);
749 tmp.reinit(distributed_vector);
754 template <
typename VectorType>
757 const std::vector<IndexSet> &partitioning)
769 template <
typename VectorType>
770 template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
773 const MatrixType2 & B,
774 const INVERSE & inverse,
777 const unsigned int n_eigenvalues)
779 std::vector<VectorType *> eigenvectors_ptr(eigenvectors.size());
780 for (
unsigned int i = 0; i < eigenvectors.size(); ++i)
781 eigenvectors_ptr[i] = &eigenvectors[i];
787 template <
typename VectorType>
788 template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
792 const INVERSE & inverse,
795 const unsigned int n_eigenvalues)
799 Assert(n_eigenvalues <= eigenvectors.size(),
801 eigenvectors.size()));
804 Assert(n_eigenvalues + 1 <= eigenvectors.size(),
806 eigenvectors.size()));
814 Assert(n_eigenvalues < eigenvectors[0]->size(),
816 eigenvectors[0]->size()));
835 bmat[0] = (mode == 1) ?
'I' :
'G';
852 std::strcpy(which,
"LA");
855 std::strcpy(which,
"SA");
858 std::strcpy(which,
"LM");
861 std::strcpy(which,
"SM");
864 std::strcpy(which,
"LR");
867 std::strcpy(which,
"SR");
870 std::strcpy(which,
"LI");
873 std::strcpy(which,
"SI");
876 std::strcpy(which,
"BE");
884 std::vector<int> iparam(11, 0);
904 std::vector<int> ipntr(14, 0);
915 int nev = n_eigenvalues;
916 int n_inside_arpack =
nloc;
967 const int shift_x = ipntr[0] - 1;
968 const int shift_y = ipntr[1] - 1;
979 if ((ido == -1) || (ido == 1 && mode < 3))
988 mass_matrix.vmult(
tmp, src);
994 system_matrix.vmult(
tmp, src);
998 workd.data() + shift_x);
1003 system_matrix.vmult(
dst, src);
1008 else if (ido == 1 && mode >= 3)
1012 const int shift_b_x = ipntr[2] - 1;
1023 inverse.vmult(
dst, src);
1039 mass_matrix.vmult(
dst, src);
1050 workd.data() + shift_y);
1058 char howmany[4] =
"All";
1060 std::vector<double> eigenvalues_real(n_eigenvalues + 1, 0.);
1061 std::vector<double> eigenvalues_im(n_eigenvalues + 1, 0.);
1069 eigenvalues_real.data(),
1093 eigenvalues_real.data(),
1094 eigenvalues_im.data(),
1129 for (
int i = 0; i < nev; ++i)
1131 (*eigenvectors[i]) = 0.0;
1138 for (
size_type i = 0; i < n_eigenvalues; ++i)
1140 std::complex<double>(eigenvalues_real[i], eigenvalues_im[i]);
1143 AssertThrow(iparam[4] >= static_cast<int>(n_eigenvalues),
1159 template <
typename VectorType>
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)
#define DeclException2(Exception2, type1, type2, outsequence)
static ::ExceptionBase & PArpackExcSmallNumberofArnoldiVectors(int arg1, int arg2)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & PArpackExcInvalidEigenvectorSize(int arg1, int arg2)
virtual State check(const unsigned int step, const double check_value)
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)
const unsigned int number_of_arnoldi_vectors
const WhichEigenvalues eigenvalue_of_interest
#define AssertIndexRange(index, range)
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)
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)
const AdditionalData additional_data
void set_shift(const std::complex< double > sigma)
AdditionalData(const unsigned int number_of_arnoldi_vectors=15, const WhichEigenvalues eigenvalue_of_interest=largest_magnitude, const bool symmetric=false, const int mode=3)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & PArpackExcNoShifts(int arg1)
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)
#define DeclException1(Exception1, type1, outsequence)
void mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
bool initial_vector_provided
SolverControl & control() const
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::vector< int > select
static ::ExceptionBase & PArpackExcConvergedEigenvectors(int arg1, int arg2)
static ::ExceptionBase & PArpackExcInfoPdnaupd(int arg1)
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & PArpackExcInfoMaxIt(int arg1)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
MPI_Comm mpi_communicator
static ::ExceptionBase & PArpackExcInvalidEigenvalueSize(int arg1, int arg2)
void fill_index_vector(std::vector< size_type > &indices) const
std::vector< double > resid
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)
static ::ExceptionBase & PArpackExcIdo(int arg1)
void internal_reinit(const IndexSet &locally_owned_dofs)
#define DEAL_II_NAMESPACE_OPEN
void set_initial_vector(const VectorType &vec)
SolverControl & solver_control
unsigned int max_steps() const
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)
static ::ExceptionBase & ExcNotImplemented()
std::size_t memory_consumption() const
void reinit(const IndexSet &locally_owned_dofs)
static ::ExceptionBase & PArpackExcInvalidNumberofArnoldiVectors(int arg1, int arg2)
std::vector< double > workl
Eigenvalue vector is filled.
size_type n_elements() const
static ::ExceptionBase & PArpackExcInvalidEigenvectorSizeNonsymmetric(int arg1, int arg2)
std::vector< types::global_dof_index > local_indices
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)