16 #ifndef dealii_parpack_solver_h
17 #define dealii_parpack_solver_h
19 #include <deal.II/base/config.h>
21 #include <deal.II/base/index_set.h>
22 #include <deal.II/base/memory_consumption.h>
23 #include <deal.II/base/smartpointer.h>
25 #include <deal.II/lac/solver_control.h>
26 #include <deal.II/lac/vector_operation.h>
31 #ifdef DEAL_II_ARPACK_WITH_PARPACK
33 DEAL_II_NAMESPACE_OPEN
39 pdnaupd_(MPI_Fint *comm,
59 pdsaupd_(MPI_Fint *comm,
79 pdneupd_(MPI_Fint *comm,
108 pdseupd_(MPI_Fint *comm,
212 template <
typename VectorType>
281 template <
typename MatrixType>
288 Shift(
const MatrixType &A,
const MatrixType &B,
const double sigma)
298 vmult(VectorType &dst,
const VectorType &
src)
const
302 A.vmult_add(dst,
src);
313 A.Tvmult_add(dst,
src);
328 const unsigned int number_of_arnoldi_vectors;
330 const bool symmetric;
333 const unsigned int number_of_arnoldi_vectors = 15,
335 const bool symmetric =
false,
366 const std::vector<IndexSet> &partitioning);
372 reinit(
const VectorType &distributed_vector);
389 set_shift(
const std::complex<double> sigma);
400 template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
402 solve(
const MatrixType1 & A,
403 const MatrixType2 & B,
404 const INVERSE & inverse,
407 const unsigned int n_eigenvalues);
412 template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
414 solve(
const MatrixType1 & A,
415 const MatrixType2 & B,
416 const INVERSE & inverse,
419 const unsigned int n_eigenvalues);
488 std::vector<double>
v;
511 std::vector<double>
z;
564 << arg1 <<
" eigenpairs were requested, but only " << arg2
570 <<
"Number of wanted eigenvalues " << arg1
571 <<
" is larger that the size of the matrix " << arg2);
576 <<
"Number of wanted eigenvalues " << arg1
577 <<
" is larger that the size of eigenvectors " << arg2);
583 <<
"To store the real and complex parts of " << arg1
584 <<
" eigenvectors in real-valued vectors, their size (currently set to "
585 << arg2 <<
") should be greater than or equal to " << arg1 + 1);
590 <<
"Number of wanted eigenvalues " << arg1
591 <<
" is larger that the size of eigenvalues " << arg2);
596 <<
"Number of Arnoldi vectors " << arg1
597 <<
" is larger that the size of the matrix " << arg2);
602 <<
"Number of Arnoldi vectors " << arg1
603 <<
" is too small to obtain " << arg2 <<
" eigenvalues");
607 <<
"This ido " << arg1
608 <<
" is not supported. Check documentation of ARPACK");
612 <<
"This mode " << arg1
613 <<
" is not supported. Check documentation of ARPACK");
617 <<
"Error with Pdnaupd, info " << arg1
618 <<
". Check documentation of ARPACK");
622 <<
"Error with Pdneupd, info " << arg1
623 <<
". Check documentation of ARPACK");
627 <<
"Maximum number " << arg1 <<
" of iterations reached.");
631 <<
"No shifts could be applied during implicit"
632 <<
" Arnoldi update, try increasing the number of"
633 <<
" Arnoldi vectors.");
638 template <
typename VectorType>
643 (workl.size() + workd.size() + v.size() + resid.size() + z.size() +
645 src.memory_consumption() + dst.memory_consumption() +
646 tmp.memory_consumption() +
648 local_indices.size();
653 template <
typename VectorType>
655 const unsigned int number_of_arnoldi_vectors,
656 const WhichEigenvalues eigenvalue_of_interest,
657 const bool symmetric,
659 : number_of_arnoldi_vectors(number_of_arnoldi_vectors)
660 , eigenvalue_of_interest(eigenvalue_of_interest)
661 , symmetric(symmetric)
670 "'largest real part' can only be used for non-symmetric problems!"));
674 "'smallest real part' can only be used for non-symmetric problems!"));
678 "'largest imaginary part' can only be used for non-symmetric problems!"));
682 "'smallest imaginary part' can only be used for non-symmetric problems!"));
684 Assert(mode >= 1 && mode <= 3,
685 ExcMessage(
"Currently, only modes 1, 2 and 3 are supported."));
690 template <
typename VectorType>
692 const MPI_Comm & mpi_communicator,
694 : solver_control(control)
695 , additional_data(data)
696 , mpi_communicator(mpi_communicator)
697 , mpi_communicator_fortran(MPI_Comm_c2f(mpi_communicator))
702 , initial_vector_provided(false)
711 template <
typename VectorType>
715 sigmar = sigma.real();
716 sigmai = sigma.imag();
721 template <
typename VectorType>
725 initial_vector_provided =
true;
726 Assert(resid.size() == local_indices.size(),
728 vec.extract_subvector_to(local_indices.begin(),
735 template <
typename VectorType>
744 ncv = additional_data.number_of_arnoldi_vectors;
750 v.resize(ldv * ncv, 0.0);
752 resid.resize(nloc, 1.0);
755 workd.resize(3 * nloc, 0.0);
758 additional_data.symmetric ? ncv * ncv + 8 * ncv : 3 * ncv * ncv + 6 * ncv;
759 workl.resize(lworkl, 0.);
762 z.resize(ldz * ncv, 0.);
765 lworkev = additional_data.symmetric ? 0
768 workev.resize(lworkev, 0.);
770 select.resize(ncv, 0);
775 template <
typename VectorType>
779 internal_reinit(locally_owned_dofs);
782 src.reinit(locally_owned_dofs, mpi_communicator);
783 dst.reinit(locally_owned_dofs, mpi_communicator);
784 tmp.reinit(locally_owned_dofs, mpi_communicator);
789 template <
typename VectorType>
793 internal_reinit(distributed_vector.locally_owned_elements());
796 src.reinit(distributed_vector);
797 dst.reinit(distributed_vector);
798 tmp.reinit(distributed_vector);
803 template <
typename VectorType>
806 const std::vector<IndexSet> &partitioning)
808 internal_reinit(locally_owned_dofs);
811 src.reinit(partitioning, mpi_communicator);
812 dst.reinit(partitioning, mpi_communicator);
813 tmp.reinit(partitioning, mpi_communicator);
818 template <
typename VectorType>
819 template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
822 const MatrixType2 & B,
823 const INVERSE & inverse,
826 const unsigned int n_eigenvalues)
828 std::vector<VectorType *> eigenvectors_ptr(
eigenvectors.size());
831 solve(A, B, inverse,
eigenvalues, eigenvectors_ptr, n_eigenvalues);
836 template <
typename VectorType>
837 template <
typename MatrixType1,
typename MatrixType2,
typename INVERSE>
840 const MatrixType2 &mass_matrix,
841 const INVERSE & inverse,
844 const unsigned int n_eigenvalues)
846 if (additional_data.symmetric)
849 PArpackExcInvalidEigenvectorSize(n_eigenvalues,
854 PArpackExcInvalidEigenvectorSizeNonsymmetric(n_eigenvalues,
858 PArpackExcInvalidEigenvalueSize(n_eigenvalues,
eigenvalues.size()));
864 PArpackExcInvalidNumberofEigenvalues(n_eigenvalues,
868 PArpackExcInvalidNumberofArnoldiVectors(
869 additional_data.number_of_arnoldi_vectors,
eigenvectors[0]->size()));
871 Assert(additional_data.number_of_arnoldi_vectors > 2 * n_eigenvalues + 1,
872 PArpackExcSmallNumberofArnoldiVectors(
873 additional_data.number_of_arnoldi_vectors, n_eigenvalues));
875 int mode = additional_data.mode;
884 bmat[0] = (mode == 1) ?
'I' :
'G';
898 switch (additional_data.eigenvalue_of_interest)
900 case algebraically_largest:
901 std::strcpy(which,
"LA");
903 case algebraically_smallest:
904 std::strcpy(which,
"SA");
906 case largest_magnitude:
907 std::strcpy(which,
"LM");
909 case smallest_magnitude:
910 std::strcpy(which,
"SM");
912 case largest_real_part:
913 std::strcpy(which,
"LR");
915 case smallest_real_part:
916 std::strcpy(which,
"SR");
918 case largest_imaginary_part:
919 std::strcpy(which,
"LI");
921 case smallest_imaginary_part:
922 std::strcpy(which,
"SI");
925 std::strcpy(which,
"BE");
930 double tol = control().tolerance();
933 std::vector<int> iparam(11, 0);
940 iparam[2] = control().max_steps();
953 std::vector<int> ipntr(14, 0);
961 int info = initial_vector_provided ? 1 : 0;
964 int nev = n_eigenvalues;
965 int n_inside_arpack = nloc;
971 if (additional_data.symmetric)
972 pdsaupd_(&mpi_communicator_fortran,
990 pdnaupd_(&mpi_communicator_fortran,
1008 AssertThrow(info == 0, PArpackExcInfoPdnaupd(info));
1016 const int shift_x = ipntr[0] - 1;
1017 const int shift_y = ipntr[1] - 1;
1019 Assert(shift_x + nloc <= static_cast<int>(workd.size()),
1022 Assert(shift_y + nloc <= static_cast<int>(workd.size()),
1028 if ((ido == -1) || (ido == 1 && mode < 3))
1031 src.add(nloc, local_indices.data(), workd.data() + shift_x);
1037 mass_matrix.vmult(tmp, src);
1038 inverse.vmult(dst, tmp);
1043 system_matrix.vmult(tmp, src);
1045 tmp.extract_subvector_to(local_indices.begin(),
1046 local_indices.end(),
1047 workd.data() + shift_x);
1048 inverse.vmult(dst, tmp);
1052 system_matrix.vmult(dst, src);
1057 else if (ido == 1 && mode >= 3)
1061 const int shift_b_x = ipntr[2] - 1;
1063 Assert(shift_b_x + nloc <= static_cast<int>(workd.size()),
1067 src.add(nloc, local_indices.data(), workd.data() + shift_b_x);
1072 inverse.vmult(dst, src);
1077 src.add(nloc, local_indices.data(), workd.data() + shift_x);
1088 mass_matrix.vmult(dst, src);
1097 dst.extract_subvector_to(local_indices.begin(),
1098 local_indices.end(),
1099 workd.data() + shift_y);
1107 char howmany[4] =
"All";
1109 std::vector<double> eigenvalues_real(n_eigenvalues + 1, 0.);
1110 std::vector<double> eigenvalues_im(n_eigenvalues + 1, 0.);
1113 if (additional_data.symmetric)
1114 pdseupd_(&mpi_communicator_fortran,
1118 eigenvalues_real.data(),
1138 pdneupd_(&mpi_communicator_fortran,
1142 eigenvalues_real.data(),
1143 eigenvalues_im.data(),
1167 AssertThrow(
false, PArpackExcInfoMaxIt(control().max_steps()));
1178 for (
int i = 0; i < nev; ++i)
1183 eigenvectors[i]->add(nloc, local_indices.data(), &v[i * nloc]);
1187 for (
size_type i = 0; i < n_eigenvalues; ++i)
1189 std::complex<double>(eigenvalues_real[i], eigenvalues_im[i]);
1192 AssertThrow(iparam[4] >= static_cast<int>(n_eigenvalues),
1193 PArpackExcConvergedEigenvectors(n_eigenvalues, iparam[4]));
1201 tmp.add(nloc, local_indices.data(), resid.data());
1202 solver_control.check(iparam[2], tmp.l2_norm());
1208 template <
typename VectorType>
1212 return solver_control;
1215 DEAL_II_NAMESPACE_CLOSE