16#ifndef dealii_precondition_h
17#define dealii_precondition_h
41template <
typename number>
43template <
typename number>
49 template <
typename,
typename>
112 template <
typename MatrixType>
120 template <
class VectorType>
122 vmult(VectorType &,
const VectorType &)
const;
128 template <
class VectorType>
130 Tvmult(VectorType &,
const VectorType &)
const;
135 template <
class VectorType>
143 template <
class VectorType>
241 template <
typename MatrixType>
248 template <
class VectorType>
250 vmult(VectorType &,
const VectorType &)
const;
256 template <
class VectorType>
258 Tvmult(VectorType &,
const VectorType &)
const;
262 template <
class VectorType>
270 template <
class VectorType>
360template <
typename MatrixType = SparseMatrix<
double>,
361 class VectorType = Vector<
double>>
369 const VectorType &) const;
383 vmult(VectorType &dst,
const VectorType &src)
const;
404template <
typename MatrixType = SparseMatrix<
double>,
405 typename PreconditionerType = IdentityMatrix>
475 template <
class VectorType>
477 vmult(VectorType &,
const VectorType &)
const;
483 template <
class VectorType>
485 Tvmult(VectorType &,
const VectorType &)
const;
490 template <
class VectorType>
492 step(VectorType &x,
const VectorType &rhs)
const;
497 template <
class VectorType>
499 Tstep(VectorType &x,
const VectorType &rhs)
const;
532 template <
typename MatrixType,
typename VectorType>
533 using vmult_functions_t =
decltype(std::declval<MatrixType const>().vmult(
534 std::declval<VectorType &>(),
535 std::declval<const VectorType &>(),
537 const std::function<
void(
const unsigned int,
const unsigned int)> &>(),
539 const std::function<
void(
const unsigned int,
const unsigned int)> &>()));
541 template <
typename MatrixType,
543 typename PreconditionerType>
544 constexpr bool has_vmult_with_std_functions =
545 is_supported_operation<vmult_functions_t, MatrixType, VectorType> &&
546 std::is_same<PreconditionerType,
548 (std::is_same<VectorType,
556 template <
typename MatrixType,
typename VectorType>
557 constexpr bool has_vmult_with_std_functions_for_precondition =
558 is_supported_operation<vmult_functions_t, MatrixType, VectorType>;
562 template <
typename T,
typename VectorType>
563 using Tvmult_t =
decltype(
564 std::declval<T const>().Tvmult(std::declval<VectorType &>(),
565 std::declval<const VectorType &>()));
567 template <
typename T,
typename VectorType>
568 constexpr bool has_Tvmult = is_supported_operation<Tvmult_t, T, VectorType>;
570 template <
typename T,
typename VectorType>
571 using step_t =
decltype(
572 std::declval<T const>().step(std::declval<VectorType &>(),
573 std::declval<const VectorType &>()));
575 template <
typename T,
typename VectorType>
576 constexpr bool has_step = is_supported_operation<step_t, T, VectorType>;
578 template <
typename T,
typename VectorType>
580 decltype(std::declval<T const>().step(std::declval<VectorType &>(),
581 std::declval<const VectorType &>(),
582 std::declval<const double>()));
584 template <
typename T,
typename VectorType>
585 constexpr bool has_step_omega =
586 is_supported_operation<step_omega_t, T, VectorType>;
588 template <
typename T,
typename VectorType>
589 using Tstep_t =
decltype(
590 std::declval<T const>().Tstep(std::declval<VectorType &>(),
591 std::declval<const VectorType &>()));
593 template <
typename T,
typename VectorType>
594 constexpr bool has_Tstep = is_supported_operation<Tstep_t, T, VectorType>;
596 template <
typename T,
typename VectorType>
597 using Tstep_omega_t =
598 decltype(std::declval<T const>().Tstep(std::declval<VectorType &>(),
599 std::declval<const VectorType &>(),
600 std::declval<const double>()));
602 template <
typename T,
typename VectorType>
603 constexpr bool has_Tstep_omega =
604 is_supported_operation<Tstep_omega_t, T, VectorType>;
606 template <
typename T,
typename VectorType>
607 using jacobi_step_t =
decltype(
608 std::declval<T const>().Jacobi_step(std::declval<VectorType &>(),
609 std::declval<const VectorType &>(),
610 std::declval<const double>()));
612 template <
typename T,
typename VectorType>
613 constexpr bool has_jacobi_step =
614 is_supported_operation<jacobi_step_t, T, VectorType>;
616 template <
typename T,
typename VectorType>
617 using SOR_step_t =
decltype(
618 std::declval<T const>().SOR_step(std::declval<VectorType &>(),
619 std::declval<const VectorType &>(),
620 std::declval<const double>()));
622 template <
typename T,
typename VectorType>
623 constexpr bool has_SOR_step =
624 is_supported_operation<SOR_step_t, T, VectorType>;
626 template <
typename T,
typename VectorType>
627 using SSOR_step_t =
decltype(
628 std::declval<T const>().SSOR_step(std::declval<VectorType &>(),
629 std::declval<const VectorType &>(),
630 std::declval<const double>()));
632 template <
typename T,
typename VectorType>
633 constexpr bool has_SSOR_step =
634 is_supported_operation<SSOR_step_t, T, VectorType>;
636 template <
typename MatrixType>
637 class PreconditionJacobiImpl
640 PreconditionJacobiImpl(
const MatrixType &A,
const double relaxation)
642 , relaxation(relaxation)
645 template <
typename VectorType>
647 vmult(VectorType &dst,
const VectorType &src)
const
649 this->A->precondition_Jacobi(dst, src, this->relaxation);
652 template <
typename VectorType>
654 Tvmult(VectorType &dst,
const VectorType &src)
const
657 this->vmult(dst, src);
660 template <
typename VectorType,
661 std::enable_if_t<has_jacobi_step<MatrixType, VectorType>,
662 MatrixType> * =
nullptr>
664 step(VectorType &dst,
const VectorType &src)
const
666 this->A->Jacobi_step(dst, src, this->relaxation);
669 template <
typename VectorType,
670 std::enable_if_t<!has_jacobi_step<MatrixType, VectorType>,
671 MatrixType> * =
nullptr>
673 step(VectorType &,
const VectorType &)
const
677 "Matrix A does not provide a Jacobi_step() function!"));
680 template <
typename VectorType>
682 Tstep(VectorType &dst,
const VectorType &src)
const
685 this->step(dst, src);
690 const double relaxation;
693 template <
typename MatrixType>
694 class PreconditionSORImpl
697 PreconditionSORImpl(
const MatrixType &A,
const double relaxation)
699 , relaxation(relaxation)
702 template <
typename VectorType>
704 vmult(VectorType &dst,
const VectorType &src)
const
706 this->A->precondition_SOR(dst, src, this->relaxation);
709 template <
typename VectorType>
711 Tvmult(VectorType &dst,
const VectorType &src)
const
713 this->A->precondition_TSOR(dst, src, this->relaxation);
716 template <
typename VectorType,
717 std::enable_if_t<has_SOR_step<MatrixType, VectorType>,
718 MatrixType> * =
nullptr>
720 step(VectorType &dst,
const VectorType &src)
const
722 this->A->SOR_step(dst, src, this->relaxation);
725 template <
typename VectorType,
726 std::enable_if_t<!has_SOR_step<MatrixType, VectorType>,
727 MatrixType> * =
nullptr>
729 step(VectorType &,
const VectorType &)
const
733 "Matrix A does not provide a SOR_step() function!"));
736 template <
typename VectorType,
737 std::enable_if_t<has_SOR_step<MatrixType, VectorType>,
738 MatrixType> * =
nullptr>
740 Tstep(VectorType &dst,
const VectorType &src)
const
742 this->A->TSOR_step(dst, src, this->relaxation);
745 template <
typename VectorType,
746 std::enable_if_t<!has_SOR_step<MatrixType, VectorType>,
747 MatrixType> * =
nullptr>
749 Tstep(VectorType &,
const VectorType &)
const
753 "Matrix A does not provide a TSOR_step() function!"));
758 const double relaxation;
761 template <
typename MatrixType>
762 class PreconditionSSORImpl
765 using size_type =
typename MatrixType::size_type;
767 PreconditionSSORImpl(
const MatrixType &A,
const double relaxation)
769 , relaxation(relaxation)
780 const size_type n = this->A->n();
781 pos_right_of_diagonal.resize(n,
static_cast<std::size_t
>(-1));
782 for (size_type row = 0; row < n; ++row)
789 typename MatrixType::value_type>::const_iterator it =
791 for (; it < mat->
end(row); ++it)
792 if (it->column() > row)
794 pos_right_of_diagonal[row] = it - mat->
begin();
799 template <
typename VectorType>
801 vmult(VectorType &dst,
const VectorType &src)
const
803 this->A->precondition_SSOR(dst,
806 pos_right_of_diagonal);
809 template <
typename VectorType>
811 Tvmult(VectorType &dst,
const VectorType &src)
const
813 this->A->precondition_SSOR(dst,
816 pos_right_of_diagonal);
819 template <
typename VectorType,
820 std::enable_if_t<has_SSOR_step<MatrixType, VectorType>,
821 MatrixType> * =
nullptr>
823 step(VectorType &dst,
const VectorType &src)
const
825 this->A->SSOR_step(dst, src, this->relaxation);
828 template <
typename VectorType,
829 std::enable_if_t<!has_SSOR_step<MatrixType, VectorType>,
830 MatrixType> * =
nullptr>
832 step(VectorType &,
const VectorType &)
const
836 "Matrix A does not provide a SSOR_step() function!"));
839 template <
typename VectorType>
841 Tstep(VectorType &dst,
const VectorType &src)
const
844 this->step(dst, src);
849 const double relaxation;
855 std::vector<std::size_t> pos_right_of_diagonal;
858 template <
typename MatrixType>
859 class PreconditionPSORImpl
862 using size_type =
typename MatrixType::size_type;
864 PreconditionPSORImpl(
const MatrixType & A,
865 const double relaxation,
866 const std::vector<size_type> &permutation,
867 const std::vector<size_type> &inverse_permutation)
869 , relaxation(relaxation)
870 , permutation(permutation)
871 , inverse_permutation(inverse_permutation)
874 template <
typename VectorType>
876 vmult(VectorType &dst,
const VectorType &src)
const
879 this->A->PSOR(dst, permutation, inverse_permutation, this->relaxation);
882 template <
typename VectorType>
884 Tvmult(VectorType &dst,
const VectorType &src)
const
887 this->A->TPSOR(dst, permutation, inverse_permutation, this->relaxation);
892 const double relaxation;
894 const std::vector<size_type> &permutation;
895 const std::vector<size_type> &inverse_permutation;
898 template <
typename MatrixType,
899 typename PreconditionerType,
901 std::enable_if_t<has_step_omega<PreconditionerType, VectorType>,
902 PreconditionerType> * =
nullptr>
904 step(
const MatrixType &,
905 const PreconditionerType &preconditioner,
907 const VectorType & src,
908 const double relaxation,
912 preconditioner.step(dst, src, relaxation);
917 typename PreconditionerType,
919 std::enable_if_t<!has_step_omega<PreconditionerType, VectorType> &&
920 has_step<PreconditionerType, VectorType>,
921 PreconditionerType> * =
nullptr>
923 step(
const MatrixType &,
924 const PreconditionerType &preconditioner,
926 const VectorType & src,
927 const double relaxation,
935 preconditioner.step(dst, src);
940 typename PreconditionerType,
942 std::enable_if_t<!has_step_omega<PreconditionerType, VectorType> &&
943 !has_step<PreconditionerType, VectorType>,
944 PreconditionerType> * =
nullptr>
946 step(
const MatrixType & A,
947 const PreconditionerType &preconditioner,
949 const VectorType & src,
950 const double relaxation,
951 VectorType & residual,
954 residual.reinit(dst,
true);
955 tmp.reinit(dst,
true);
957 A.vmult(residual, dst);
958 residual.sadd(-1.0, 1.0, src);
960 preconditioner.vmult(tmp, residual);
961 dst.add(relaxation, tmp);
964 template <
typename MatrixType,
965 typename PreconditionerType,
967 std::enable_if_t<has_Tstep_omega<PreconditionerType, VectorType>,
968 PreconditionerType> * =
nullptr>
970 Tstep(
const MatrixType &,
971 const PreconditionerType &preconditioner,
973 const VectorType & src,
974 const double relaxation,
978 preconditioner.Tstep(dst, src, relaxation);
983 typename PreconditionerType,
985 std::enable_if_t<!has_Tstep_omega<PreconditionerType, VectorType> &&
986 has_Tstep<PreconditionerType, VectorType>,
987 PreconditionerType> * =
nullptr>
989 Tstep(
const MatrixType &,
990 const PreconditionerType &preconditioner,
992 const VectorType & src,
993 const double relaxation,
1001 preconditioner.Tstep(dst, src);
1004 template <
typename MatrixType,
1005 typename VectorType,
1006 std::enable_if_t<has_Tvmult<MatrixType, VectorType>, MatrixType>
1009 Tvmult(
const MatrixType &A, VectorType &dst,
const VectorType &src)
1014 template <
typename MatrixType,
1015 typename VectorType,
1016 std::enable_if_t<!has_Tvmult<MatrixType, VectorType>, MatrixType>
1019 Tvmult(
const MatrixType &, VectorType &,
const VectorType &)
1022 ExcMessage(
"Matrix A does not provide a Tvmult() function!"));
1026 typename MatrixType,
1027 typename PreconditionerType,
1028 typename VectorType,
1029 std::enable_if_t<!has_Tstep_omega<PreconditionerType, VectorType> &&
1030 !has_Tstep<PreconditionerType, VectorType>,
1031 PreconditionerType> * =
nullptr>
1033 Tstep(
const MatrixType & A,
1034 const PreconditionerType &preconditioner,
1036 const VectorType & src,
1037 const double relaxation,
1038 VectorType & residual,
1041 residual.reinit(dst,
true);
1042 tmp.reinit(dst,
true);
1044 Tvmult(A, residual, dst);
1045 residual.sadd(-1.0, 1.0, src);
1047 Tvmult(preconditioner, tmp, residual);
1048 dst.add(relaxation, tmp);
1052 template <
typename MatrixType,
1053 typename PreconditionerType,
1054 typename VectorType,
1055 std::enable_if_t<!has_vmult_with_std_functions_for_precondition<
1060 step_operations(
const MatrixType & A,
1061 const PreconditionerType &preconditioner,
1063 const VectorType & src,
1064 const double relaxation,
1067 const unsigned int i,
1068 const bool transposed)
1073 Tvmult(preconditioner, dst, src);
1075 preconditioner.vmult(dst, src);
1077 if (relaxation != 1.0)
1083 Tstep(A, preconditioner, dst, src, relaxation, tmp1, tmp2);
1085 step(A, preconditioner, dst, src, relaxation, tmp1, tmp2);
1092 typename MatrixType,
1093 typename PreconditionerType,
1094 typename VectorType,
1096 has_vmult_with_std_functions_for_precondition<PreconditionerType,
1098 !has_vmult_with_std_functions_for_precondition<MatrixType,
1102 step_operations(
const MatrixType & A,
1103 const PreconditionerType &preconditioner,
1105 const VectorType & src,
1106 const double relaxation,
1109 const unsigned int i,
1110 const bool transposed)
1113 using Number =
typename VectorType::value_type;
1117 Number * dst_ptr = dst.begin();
1118 const Number *src_ptr = src.begin();
1120 preconditioner.vmult(
1123 [&](
const unsigned int start_range,
const unsigned int end_range) {
1125 if (end_range > start_range)
1126 std::memset(dst.begin() + start_range,
1128 sizeof(Number) * (end_range - start_range));
1130 [&](
const unsigned int start_range,
const unsigned int end_range) {
1131 if (relaxation == 1.0)
1134 const auto src_ptr = src.begin();
1135 const auto dst_ptr = dst.begin();
1138 for (std::size_t i = start_range; i < end_range; ++i)
1139 dst_ptr[i] *= relaxation;
1144 tmp.reinit(src,
true);
1150 preconditioner.vmult(
1153 [&](
const unsigned int start_range,
const unsigned int end_range) {
1154 const auto src_ptr = src.begin();
1155 const auto tmp_ptr = tmp.begin();
1157 if (relaxation == 1.0)
1160 for (std::size_t i = start_range; i < end_range; ++i)
1161 tmp_ptr[i] = src_ptr[i] - tmp_ptr[i];
1169 for (std::size_t i = start_range; i < end_range; ++i)
1170 tmp_ptr[i] = relaxation * (src_ptr[i] - tmp_ptr[i]);
1173 [&](
const unsigned int,
const unsigned int) {
1183 typename MatrixType,
1184 typename PreconditionerType,
1185 typename VectorType,
1187 has_vmult_with_std_functions_for_precondition<PreconditionerType,
1189 has_vmult_with_std_functions_for_precondition<MatrixType, VectorType>,
1192 step_operations(
const MatrixType & A,
1193 const PreconditionerType &preconditioner,
1195 const VectorType & src,
1196 const double relaxation,
1199 const unsigned int i,
1200 const bool transposed)
1203 using Number =
typename VectorType::value_type;
1207 Number * dst_ptr = dst.begin();
1208 const Number *src_ptr = src.begin();
1210 preconditioner.vmult(
1213 [&](
const unsigned int start_range,
const unsigned int end_range) {
1215 if (end_range > start_range)
1216 std::memset(dst.begin() + start_range,
1218 sizeof(Number) * (end_range - start_range));
1220 [&](
const unsigned int start_range,
const unsigned int end_range) {
1221 if (relaxation == 1.0)
1224 const auto src_ptr = src.begin();
1225 const auto dst_ptr = dst.begin();
1228 for (std::size_t i = start_range; i < end_range; ++i)
1229 dst_ptr[i] *= relaxation;
1234 tmp.reinit(src,
true);
1241 [&](
const unsigned int start_range,
const unsigned int end_range) {
1244 if (end_range > start_range)
1245 std::memset(tmp.begin() + start_range,
1247 sizeof(Number) * (end_range - start_range));
1249 [&](
const unsigned int start_range,
const unsigned int end_range) {
1250 const auto src_ptr = src.begin();
1251 const auto tmp_ptr = tmp.begin();
1253 if (relaxation == 1.0)
1256 for (std::size_t i = start_range; i < end_range; ++i)
1257 tmp_ptr[i] = src_ptr[i] - tmp_ptr[i];
1265 for (std::size_t i = start_range; i < end_range; ++i)
1266 tmp_ptr[i] = relaxation * (src_ptr[i] - tmp_ptr[i]);
1270 preconditioner.vmult(dst, tmp, [](
const auto,
const auto) {
1278 template <
typename MatrixType,
1279 typename VectorType,
1280 std::enable_if_t<!IsBlockVector<VectorType>::value &&
1281 !has_vmult_with_std_functions<
1285 VectorType> * =
nullptr>
1287 step_operations(
const MatrixType & A,
1288 const ::DiagonalMatrix<VectorType> &preconditioner,
1290 const VectorType & src,
1291 const double relaxation,
1294 const unsigned int i,
1295 const bool transposed)
1297 using Number =
typename VectorType::value_type;
1301 Number * dst_ptr = dst.begin();
1302 const Number *src_ptr = src.begin();
1303 const Number *diag_ptr = preconditioner.
get_vector().begin();
1305 if (relaxation == 1.0)
1308 for (
unsigned int i = 0; i < dst.locally_owned_size(); ++i)
1309 dst_ptr[i] = src_ptr[i] * diag_ptr[i];
1314 for (
unsigned int i = 0; i < dst.locally_owned_size(); ++i)
1315 dst_ptr[i] = relaxation * src_ptr[i] * diag_ptr[i];
1320 tmp.reinit(src,
true);
1322 Number * dst_ptr = dst.begin();
1323 const Number *src_ptr = src.begin();
1324 const Number *tmp_ptr = tmp.begin();
1325 const Number *diag_ptr = preconditioner.get_vector().begin();
1328 Tvmult(A, tmp, dst);
1332 if (relaxation == 1.0)
1335 for (
unsigned int i = 0; i < dst.locally_owned_size(); ++i)
1336 dst_ptr[i] += (src_ptr[i] - tmp_ptr[i]) * diag_ptr[i];
1341 for (
unsigned int i = 0; i < dst.locally_owned_size(); ++i)
1343 relaxation * (src_ptr[i] - tmp_ptr[i]) * diag_ptr[i];
1350 template <
typename MatrixType,
1351 typename VectorType,
1352 std::enable_if_t<!IsBlockVector<VectorType>::value &&
1353 has_vmult_with_std_functions<
1357 VectorType> * =
nullptr>
1359 step_operations(
const MatrixType & A,
1360 const ::DiagonalMatrix<VectorType> &preconditioner,
1362 const VectorType & src,
1363 const double relaxation,
1366 const unsigned int i,
1367 const bool transposed)
1370 using Number =
typename VectorType::value_type;
1374 Number * dst_ptr = dst.begin();
1375 const Number *src_ptr = src.begin();
1376 const Number *diag_ptr = preconditioner.get_vector().begin();
1378 if (relaxation == 1.0)
1381 for (
unsigned int i = 0; i < dst.locally_owned_size(); ++i)
1382 dst_ptr[i] = src_ptr[i] * diag_ptr[i];
1387 for (
unsigned int i = 0; i < dst.locally_owned_size(); ++i)
1388 dst_ptr[i] = relaxation * src_ptr[i] * diag_ptr[i];
1393 tmp.reinit(src,
true);
1400 [&](
const unsigned int start_range,
const unsigned int end_range) {
1402 if (end_range > start_range)
1403 std::memset(tmp.begin() + start_range,
1405 sizeof(Number) * (end_range - start_range));
1407 [&](
const unsigned int begin,
const unsigned int end) {
1408 const Number *dst_ptr = dst.begin();
1409 const Number *src_ptr = src.begin();
1410 Number * tmp_ptr = tmp.begin();
1411 const Number *diag_ptr = preconditioner.get_vector().begin();
1415 if (relaxation == 1.0)
1418 for (std::size_t i = begin; i < end; ++i)
1420 dst_ptr[i] + (src_ptr[i] - tmp_ptr[i]) * diag_ptr[i];
1425 for (std::size_t i = begin; i < end; ++i)
1426 tmp_ptr[i] = dst_ptr[i] + relaxation *
1427 (src_ptr[i] - tmp_ptr[i]) *
1469template <
typename MatrixType = SparseMatrix<
double>>
1473 internal::PreconditionRelaxation::PreconditionJacobiImpl<MatrixType>>
1476 internal::PreconditionRelaxation::PreconditionJacobiImpl<MatrixType>;
1539template <
typename MatrixType = SparseMatrix<
double>>
1543 internal::PreconditionRelaxation::PreconditionSORImpl<MatrixType>>
1546 internal::PreconditionRelaxation::PreconditionSORImpl<MatrixType>;
1591template <
typename MatrixType = SparseMatrix<
double>>
1595 internal::PreconditionRelaxation::PreconditionSSORImpl<MatrixType>>
1598 internal::PreconditionRelaxation::PreconditionSSORImpl<MatrixType>;
1647template <
typename MatrixType = SparseMatrix<
double>>
1651 internal::PreconditionRelaxation::PreconditionPSORImpl<MatrixType>>
1654 internal::PreconditionRelaxation::PreconditionPSORImpl<MatrixType>;
1711 const std::vector<size_type> & permutation,
1712 const std::vector<size_type> & inverse_permutation,
1929template <
typename MatrixType = SparseMatrix<
double>,
1930 typename VectorType = Vector<
double>,
1931 typename PreconditionerType = DiagonalMatrix<VectorType>>
1989 const unsigned int degree = 1,
2100 vmult(VectorType &dst,
const VectorType &src)
const;
2107 Tvmult(VectorType &dst,
const VectorType &src)
const;
2113 step(VectorType &dst,
const VectorType &src)
const;
2119 Tstep(VectorType &dst,
const VectorType &src)
const;
2187 EigenvalueInformation
2256template <
typename MatrixType>
2266template <
class VectorType>
2275template <
class VectorType>
2282template <
class VectorType>
2291template <
class VectorType>
2323 const double relaxation)
2324 : relaxation(relaxation)
2333 AdditionalData add_data;
2334 relaxation = add_data.relaxation;
2348template <
typename MatrixType>
2351 const MatrixType & matrix,
2361template <
class VectorType>
2366 std::is_same<size_type, typename VectorType::size_type>::value,
2367 "PreconditionRichardson and VectorType must have the same size_type.");
2374template <
class VectorType>
2379 std::is_same<size_type, typename VectorType::size_type>::value,
2380 "PreconditionRichardson and VectorType must have the same size_type.");
2385template <
class VectorType>
2390 std::is_same<size_type, typename VectorType::size_type>::value,
2391 "PreconditionRichardson and VectorType must have the same size_type.");
2398template <
class VectorType>
2403 std::is_same<size_type, typename VectorType::size_type>::value,
2404 "PreconditionRichardson and VectorType must have the same size_type.");
2425template <
typename MatrixType,
typename PreconditionerType>
2428 const MatrixType & rA,
2429 const AdditionalData ¶meters)
2432 relaxation = parameters.relaxation;
2436 preconditioner = parameters.preconditioner;
2437 n_iterations = parameters.n_iterations;
2441template <
typename MatrixType,
typename PreconditionerType>
2446 preconditioner =
nullptr;
2449template <
typename MatrixType,
typename PreconditionerType>
2458template <
typename MatrixType,
typename PreconditionerType>
2467template <
typename MatrixType,
typename PreconditionerType>
2468template <
class VectorType>
2472 const VectorType &src)
const
2477 VectorType tmp1, tmp2;
2479 for (
unsigned int i = 0; i < n_iterations; ++i)
2480 internal::PreconditionRelaxation::step_operations(
2481 *A, *preconditioner, dst, src, relaxation, tmp1, tmp2, i,
false);
2484template <
typename MatrixType,
typename PreconditionerType>
2485template <
class VectorType>
2489 const VectorType &src)
const
2494 VectorType tmp1, tmp2;
2496 for (
unsigned int i = 0; i < n_iterations; ++i)
2497 internal::PreconditionRelaxation::step_operations(
2498 *A, *preconditioner, dst, src, relaxation, tmp1, tmp2, i,
true);
2501template <
typename MatrixType,
typename PreconditionerType>
2502template <
class VectorType>
2506 const VectorType &src)
const
2511 VectorType tmp1, tmp2;
2513 for (
unsigned int i = 1; i <= n_iterations; ++i)
2514 internal::PreconditionRelaxation::step_operations(
2515 *A, *preconditioner, dst, src, relaxation, tmp1, tmp2, i,
false);
2518template <
typename MatrixType,
typename PreconditionerType>
2519template <
class VectorType>
2523 const VectorType &src)
const
2528 VectorType tmp1, tmp2;
2530 for (
unsigned int i = 1; i <= n_iterations; ++i)
2531 internal::PreconditionRelaxation::step_operations(
2532 *A, *preconditioner, dst, src, relaxation, tmp1, tmp2, i,
true);
2537template <
typename MatrixType>
2540 const AdditionalData ¶meters_in)
2544 AdditionalData parameters;
2545 parameters.relaxation = 1.0;
2546 parameters.n_iterations = parameters_in.n_iterations;
2547 parameters.preconditioner =
2548 std::make_shared<PreconditionerType>(A, parameters_in.relaxation);
2550 this->BaseClass::initialize(A, parameters);
2555template <
typename MatrixType>
2558 const AdditionalData ¶meters_in)
2562 AdditionalData parameters;
2563 parameters.relaxation = 1.0;
2564 parameters.n_iterations = parameters_in.n_iterations;
2565 parameters.preconditioner =
2566 std::make_shared<PreconditionerType>(A, parameters_in.relaxation);
2568 this->BaseClass::initialize(A, parameters);
2573template <
typename MatrixType>
2576 const AdditionalData ¶meters_in)
2580 AdditionalData parameters;
2581 parameters.relaxation = 1.0;
2582 parameters.n_iterations = parameters_in.n_iterations;
2583 parameters.preconditioner =
2584 std::make_shared<PreconditionerType>(A, parameters_in.relaxation);
2586 this->BaseClass::initialize(A, parameters);
2593template <
typename MatrixType>
2596 const MatrixType & A,
2597 const std::vector<size_type> & p,
2598 const std::vector<size_type> & ip,
2599 const typename BaseClass::AdditionalData ¶meters_in)
2603 typename BaseClass::AdditionalData parameters;
2604 parameters.relaxation = 1.0;
2605 parameters.n_iterations = parameters_in.n_iterations;
2606 parameters.preconditioner =
2607 std::make_shared<PreconditionerType>(A, parameters_in.relaxation, p, ip);
2609 this->BaseClass::initialize(A, parameters);
2613template <
typename MatrixType>
2616 const AdditionalData &additional_data)
2619 additional_data.permutation,
2620 additional_data.inverse_permutation,
2621 additional_data.parameters);
2624template <
typename MatrixType>
2626 const std::vector<size_type> &permutation,
2627 const std::vector<size_type> &inverse_permutation,
2630 : permutation(permutation)
2631 , inverse_permutation(inverse_permutation)
2632 , parameters(parameters)
2639template <
typename MatrixType,
class VectorType>
2641 const MatrixType & M,
2642 const function_ptr method)
2644 , precondition(method)
2649template <
typename MatrixType,
class VectorType>
2653 const VectorType &src)
const
2655 (
matrix.*precondition)(dst, src);
2660template <
typename MatrixType,
typename PreconditionerType>
2662 AdditionalData(
const double relaxation,
const unsigned int n_iterations)
2663 : relaxation(relaxation)
2664 , n_iterations(n_iterations)
2673 namespace PreconditionChebyshevImplementation
2681 template <
typename VectorType,
typename PreconditionerType>
2683 vector_updates(
const VectorType & rhs,
2684 const PreconditionerType &preconditioner,
2685 const unsigned int iteration_index,
2686 const double factor1,
2687 const double factor2,
2688 VectorType & solution_old,
2689 VectorType & temp_vector1,
2690 VectorType & temp_vector2,
2691 VectorType & solution)
2693 if (iteration_index == 0)
2695 solution.equ(factor2, rhs);
2696 preconditioner.vmult(solution_old, solution);
2698 else if (iteration_index == 1)
2701 temp_vector1.sadd(-1.0, 1.0, rhs);
2702 preconditioner.vmult(solution_old, temp_vector1);
2705 solution_old.sadd(factor2, 1 + factor1, solution);
2710 temp_vector1.sadd(-1.0, 1.0, rhs);
2711 preconditioner.vmult(temp_vector2, temp_vector1);
2714 solution_old.sadd(-factor1, factor2, temp_vector2);
2715 solution_old.add(1 + factor1, solution);
2718 solution.swap(solution_old);
2724 typename PreconditionerType,
2726 !has_vmult_with_std_functions_for_precondition<
2733 const PreconditionerType &preconditioner,
2734 const unsigned int iteration_index,
2735 const double factor1_,
2736 const double factor2_,
2745 const Number factor1 = factor1_;
2746 const Number factor1_plus_1 = 1. + factor1_;
2747 const Number factor2 = factor2_;
2749 if (iteration_index == 0)
2751 const auto solution_old_ptr = solution_old.
begin();
2754 preconditioner.vmult(solution_old, rhs);
2759 solution_old_ptr[i] = solution_old_ptr[i] * factor2;
2761 else if (iteration_index == 1)
2763 const auto solution_ptr = solution.
begin();
2764 const auto solution_old_ptr = solution_old.
begin();
2767 temp_vector1.
sadd(-1.0, 1.0, rhs);
2769 preconditioner.vmult(solution_old, temp_vector1);
2774 solution_old_ptr[i] =
2775 factor1_plus_1 * solution_ptr[i] + solution_old_ptr[i] * factor2;
2779 const auto solution_ptr = solution.
begin();
2780 const auto solution_old_ptr = solution_old.
begin();
2781 const auto temp_vector2_ptr = temp_vector2.
begin();
2784 temp_vector1.
sadd(-1.0, 1.0, rhs);
2786 preconditioner.vmult(temp_vector2, temp_vector1);
2791 solution_old_ptr[i] = factor1_plus_1 * solution_ptr[i] -
2792 factor1 * solution_old_ptr[i] +
2793 temp_vector2_ptr[i] * factor2;
2796 solution.
swap(solution_old);
2801 typename PreconditionerType,
2803 has_vmult_with_std_functions_for_precondition<
2810 const PreconditionerType &preconditioner,
2811 const unsigned int iteration_index,
2812 const double factor1_,
2813 const double factor2_,
2822 const Number factor1 = factor1_;
2823 const Number factor1_plus_1 = 1. + factor1_;
2824 const Number factor2 = factor2_;
2826 const auto rhs_ptr = rhs.
begin();
2827 const auto temp_vector1_ptr = temp_vector1.
begin();
2828 const auto temp_vector2_ptr = temp_vector2.
begin();
2829 const auto solution_ptr = solution.
begin();
2830 const auto solution_old_ptr = solution_old.
begin();
2832 if (iteration_index == 0)
2834 preconditioner.vmult(
2837 [&](
const auto start_range,
const auto end_range) {
2838 if (end_range > start_range)
2839 std::memset(solution.
begin() + start_range,
2841 sizeof(Number) * (end_range - start_range));
2843 [&](
const auto begin,
const auto end) {
2845 for (std::size_t i = begin; i <
end; ++i)
2846 solution_ptr[i] *= factor2;
2851 preconditioner.vmult(
2854 [&](
const auto begin,
const auto end) {
2856 std::memset(temp_vector2.
begin() + begin,
2858 sizeof(Number) * (end - begin));
2861 for (std::size_t i = begin; i <
end; ++i)
2862 temp_vector1_ptr[i] = rhs_ptr[i] - temp_vector1_ptr[i];
2864 [&](
const auto begin,
const auto end) {
2865 if (iteration_index == 1)
2868 for (std::size_t i = begin; i <
end; ++i)
2869 temp_vector2_ptr[i] = factor1_plus_1 * solution_ptr[i] +
2870 factor2 * temp_vector2_ptr[i];
2875 for (std::size_t i = begin; i <
end; ++i)
2876 temp_vector2_ptr[i] = factor1_plus_1 * solution_ptr[i] -
2877 factor1 * solution_old_ptr[i] +
2878 factor2 * temp_vector2_ptr[i];
2883 if (iteration_index > 0)
2885 solution_old.
swap(temp_vector2);
2886 solution_old.
swap(solution);
2893 template <
typename Number>
2894 struct VectorUpdater
2896 VectorUpdater(
const Number * rhs,
2897 const Number * matrix_diagonal_inverse,
2898 const unsigned int iteration_index,
2899 const Number factor1,
2900 const Number factor2,
2901 Number * solution_old,
2902 Number * tmp_vector,
2905 , matrix_diagonal_inverse(matrix_diagonal_inverse)
2906 , iteration_index(iteration_index)
2909 , solution_old(solution_old)
2910 , tmp_vector(tmp_vector)
2911 , solution(solution)
2915 apply_to_subrange(
const std::size_t begin,
const std::size_t end)
const
2921 const Number factor1 = this->factor1;
2922 const Number factor1_plus_1 = 1. + this->factor1;
2923 const Number factor2 = this->factor2;
2924 if (iteration_index == 0)
2927 for (std::size_t i = begin; i <
end; ++i)
2928 solution[i] = factor2 * matrix_diagonal_inverse[i] * rhs[i];
2930 else if (iteration_index == 1)
2934 for (std::size_t i = begin; i <
end; ++i)
2938 factor1_plus_1 * solution[i] +
2939 factor2 * matrix_diagonal_inverse[i] * (rhs[i] - tmp_vector[i]);
2946 for (std::size_t i = begin; i <
end; ++i)
2953 factor1_plus_1 * solution[i] - factor1 * solution_old[i] +
2954 factor2 * matrix_diagonal_inverse[i] * (rhs[i] - tmp_vector[i]);
2959 const Number * matrix_diagonal_inverse;
2960 const unsigned int iteration_index;
2961 const Number factor1;
2962 const Number factor2;
2963 mutable Number * solution_old;
2964 mutable Number * tmp_vector;
2965 mutable Number * solution;
2968 template <
typename Number>
2971 VectorUpdatesRange(
const VectorUpdater<Number> &updater,
2972 const std::size_t size)
2976 VectorUpdatesRange::apply_to_subrange(0, size);
2984 ~VectorUpdatesRange()
override =
default;
2987 apply_to_subrange(
const std::size_t begin,
2988 const std::size_t end)
const override
2990 updater.apply_to_subrange(begin, end);
2993 const VectorUpdater<Number> &updater;
2997 template <
typename Number>
3000 const ::Vector<Number> & rhs,
3002 const unsigned int iteration_index,
3003 const double factor1,
3004 const double factor2,
3010 VectorUpdater<Number> upd(rhs.begin(),
3011 jacobi.get_vector().begin(),
3015 solution_old.
begin(),
3016 temp_vector1.
begin(),
3018 VectorUpdatesRange<Number>(upd, rhs.size());
3021 if (iteration_index == 0)
3028 solution.
swap(temp_vector1);
3029 solution_old.
swap(temp_vector1);
3034 template <
typename Number>
3038 const ::DiagonalMatrix<
3040 const unsigned int iteration_index,
3041 const double factor1,
3042 const double factor2,
3050 VectorUpdater<Number> upd(rhs.
begin(),
3051 jacobi.get_vector().begin(),
3055 solution_old.
begin(),
3056 temp_vector1.
begin(),
3061 if (iteration_index == 0)
3068 solution.
swap(temp_vector1);
3069 solution_old.
swap(temp_vector1);
3078 typename MatrixType,
3079 typename VectorType,
3080 typename PreconditionerType,
3082 !has_vmult_with_std_functions<MatrixType,
3084 PreconditionerType> &&
3085 !(has_vmult_with_std_functions_for_precondition<PreconditionerType,
3087 has_vmult_with_std_functions_for_precondition<MatrixType,
3091 vmult_and_update(
const MatrixType & matrix,
3092 const PreconditionerType &preconditioner,
3093 const VectorType & rhs,
3094 const unsigned int iteration_index,
3095 const double factor1,
3096 const double factor2,
3097 VectorType & solution,
3098 VectorType & solution_old,
3099 VectorType & temp_vector1,
3100 VectorType & temp_vector2)
3102 if (iteration_index > 0)
3103 matrix.vmult(temp_vector1, solution);
3118 typename MatrixType,
3119 typename VectorType,
3120 typename PreconditionerType,
3122 !has_vmult_with_std_functions<MatrixType,
3124 PreconditionerType> &&
3125 (has_vmult_with_std_functions_for_precondition<PreconditionerType,
3127 has_vmult_with_std_functions_for_precondition<MatrixType,
3131 vmult_and_update(
const MatrixType & matrix,
3132 const PreconditionerType &preconditioner,
3133 const VectorType & rhs,
3134 const unsigned int iteration_index,
3135 const double factor1_,
3136 const double factor2_,
3137 VectorType & solution,
3138 VectorType & solution_old,
3139 VectorType & temp_vector1,
3140 VectorType & temp_vector2)
3142 using Number =
typename VectorType::value_type;
3144 const Number factor1 = factor1_;
3145 const Number factor1_plus_1 = 1. + factor1_;
3146 const Number factor2 = factor2_;
3148 if (iteration_index == 0)
3150 preconditioner.vmult(
3153 [&](
const unsigned int start_range,
const unsigned int end_range) {
3155 if (end_range > start_range)
3156 std::memset(solution.begin() + start_range,
3158 sizeof(Number) * (end_range - start_range));
3160 [&](
const unsigned int start_range,
const unsigned int end_range) {
3161 const auto solution_ptr = solution.begin();
3164 for (std::size_t i = start_range; i < end_range; ++i)
3165 solution_ptr[i] *= factor2;
3170 temp_vector1.reinit(rhs,
true);
3171 temp_vector2.reinit(rhs,
true);
3177 [&](
const unsigned int start_range,
const unsigned int end_range) {
3180 if (end_range > start_range)
3181 std::memset(temp_vector1.begin() + start_range,
3183 sizeof(Number) * (end_range - start_range));
3185 [&](
const unsigned int start_range,
const unsigned int end_range) {
3186 const auto rhs_ptr = rhs.begin();
3187 const auto tmp_ptr = temp_vector1.begin();
3190 for (std::size_t i = start_range; i < end_range; ++i)
3191 tmp_ptr[i] = rhs_ptr[i] - tmp_ptr[i];
3195 preconditioner.vmult(
3198 [&](
const unsigned int start_range,
const unsigned int end_range) {
3201 if (end_range > start_range)
3202 std::memset(temp_vector2.begin() + start_range,
3204 sizeof(Number) * (end_range - start_range));
3206 [&](
const unsigned int start_range,
const unsigned int end_range) {
3207 const auto solution_ptr = solution.begin();
3208 const auto solution_old_ptr = solution_old.begin();
3209 const auto tmp_ptr = temp_vector2.begin();
3211 if (iteration_index == 1)
3214 for (std::size_t i = start_range; i < end_range; ++i)
3216 factor1_plus_1 * solution_ptr[i] + factor2 * tmp_ptr[i];
3221 for (std::size_t i = start_range; i < end_range; ++i)
3222 tmp_ptr[i] = factor1_plus_1 * solution_ptr[i] -
3223 factor1 * solution_old_ptr[i] +
3224 factor2 * tmp_ptr[i];
3228 solution.swap(temp_vector2);
3229 solution_old.swap(temp_vector2);
3235 template <
typename MatrixType,
3236 typename VectorType,
3237 typename PreconditionerType,
3238 std::enable_if_t<has_vmult_with_std_functions<MatrixType,
3240 PreconditionerType>,
3243 vmult_and_update(
const MatrixType & matrix,
3244 const PreconditionerType &preconditioner,
3245 const VectorType & rhs,
3246 const unsigned int iteration_index,
3247 const double factor1,
3248 const double factor2,
3249 VectorType & solution,
3250 VectorType & solution_old,
3251 VectorType & temp_vector1,
3254 using Number =
typename VectorType::value_type;
3255 VectorUpdater<Number> updater(rhs.begin(),
3256 preconditioner.get_vector().begin(),
3260 solution_old.begin(),
3261 temp_vector1.begin(),
3263 if (iteration_index > 0)
3267 [&](
const unsigned int start_range,
const unsigned int end_range) {
3270 if (end_range > start_range)
3271 std::memset(temp_vector1.begin() + start_range,
3273 sizeof(Number) * (end_range - start_range));
3275 [&](
const unsigned int start_range,
const unsigned int end_range) {
3276 if (end_range > start_range)
3277 updater.apply_to_subrange(start_range, end_range);
3280 updater.apply_to_subrange(0U, solution.locally_owned_size());
3283 if (iteration_index == 0)
3290 solution.swap(temp_vector1);
3291 solution_old.swap(temp_vector1);
3295 template <
typename MatrixType,
typename PreconditionerType>
3297 initialize_preconditioner(
3298 const MatrixType & matrix,
3299 std::shared_ptr<PreconditionerType> &preconditioner)
3302 (void)preconditioner;
3306 template <
typename MatrixType,
typename VectorType>
3308 initialize_preconditioner(
3309 const MatrixType & matrix,
3312 if (preconditioner.get() ==
nullptr || preconditioner->m() !=
matrix.m())
3314 if (preconditioner.get() ==
nullptr)
3316 std::make_shared<::DiagonalMatrix<VectorType>>();
3319 preconditioner->m() == 0,
3321 "Preconditioner appears to be initialized but not sized correctly"));
3324 if (preconditioner->m() !=
matrix.m())
3326 preconditioner->get_vector().reinit(
matrix.m());
3327 for (
typename VectorType::size_type i = 0; i <
matrix.m(); ++i)
3328 preconditioner->get_vector()(i) = 1. /
matrix.el(i, i);
3333 template <
typename VectorType>
3335 set_initial_guess(VectorType &vector)
3337 vector = 1. /
std::sqrt(
static_cast<double>(vector.size()));
3338 if (vector.locally_owned_elements().is_element(0))
3342 template <
typename Number>
3350 for (
unsigned int i = 0; i < vector.
size(); ++i)
3353 const Number mean_value = vector.
mean_value();
3354 vector.
add(-mean_value);
3357 template <
typename Number>
3362 for (
unsigned int block = 0; block < vector.
n_blocks(); ++block)
3363 set_initial_guess(vector.
block(block));
3366 template <
typename Number,
typename MemorySpace>
3383 Kokkos::RangePolicy<
typename MemorySpace::kokkos_space::execution_space,
3384 Kokkos::IndexType<types::global_dof_index>>
3385 policy(0, n_local_elements);
3386 Kokkos::parallel_for(
3387 "::PreconditionChebyshev::set_initial_guess",
3390 values_ptr[i] = (i + first_local_range) % 11;
3392 const Number mean_value = vector.
mean_value();
3393 vector.
add(-mean_value);
3396 struct EigenvalueTracker
3405 std::vector<double>
values;
3410 template <
typename MatrixType,
3411 typename VectorType,
3412 typename PreconditionerType>
3414 power_iteration(
const MatrixType & matrix,
3415 VectorType & eigenvector,
3416 const PreconditionerType &preconditioner,
3417 const unsigned int n_iterations)
3419 double eigenvalue_estimate = 0.;
3420 eigenvector /= eigenvector.l2_norm();
3421 VectorType vector1, vector2;
3422 vector1.reinit(eigenvector,
true);
3423 if (!std::is_same<PreconditionerType, PreconditionIdentity>::value)
3424 vector2.reinit(eigenvector,
true);
3425 for (
unsigned int i = 0; i < n_iterations; ++i)
3427 if (!std::is_same<PreconditionerType, PreconditionIdentity>::value)
3429 matrix.vmult(vector2, eigenvector);
3430 preconditioner.vmult(vector1, vector2);
3433 matrix.vmult(vector1, eigenvector);
3435 eigenvalue_estimate = eigenvector * vector1;
3437 vector1 /= vector1.l2_norm();
3438 eigenvector.swap(vector1);
3440 return eigenvalue_estimate;
3447template <
typename MatrixType,
class VectorType,
typename PreconditionerType>
3450 const double smoothing_range,
3451 const unsigned int eig_cg_n_iterations,
3452 const double eig_cg_residual,
3453 const double max_eigenvalue,
3454 const EigenvalueAlgorithm eigenvalue_algorithm,
3455 const PolynomialType polynomial_type)
3457 , smoothing_range(smoothing_range)
3458 , eig_cg_n_iterations(eig_cg_n_iterations)
3459 , eig_cg_residual(eig_cg_residual)
3460 , max_eigenvalue(max_eigenvalue)
3461 , eigenvalue_algorithm(eigenvalue_algorithm)
3462 , polynomial_type(polynomial_type)
3467template <
typename MatrixType,
class VectorType,
typename PreconditionerType>
3470 PreconditionerType>::AdditionalData &
3474 degree = other_data.
degree;
3475 smoothing_range = other_data.smoothing_range;
3476 eig_cg_n_iterations = other_data.eig_cg_n_iterations;
3477 eig_cg_residual = other_data.eig_cg_residual;
3478 max_eigenvalue = other_data.max_eigenvalue;
3479 preconditioner = other_data.preconditioner;
3480 eigenvalue_algorithm = other_data.eigenvalue_algorithm;
3481 polynomial_type = other_data.polynomial_type;
3482 constraints.copy_from(other_data.constraints);
3489template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
3494 , eigenvalues_are_initialized(false)
3497 std::is_same<size_type, typename VectorType::size_type>::value,
3498 "PreconditionChebyshev and VectorType must have the same size_type.");
3503template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
3506 const MatrixType & matrix,
3507 const AdditionalData &additional_data)
3510 data = additional_data;
3512 ExcMessage(
"The degree of the Chebyshev method must be positive."));
3513 internal::PreconditionChebyshevImplementation::initialize_preconditioner(
3514 matrix, data.preconditioner);
3515 eigenvalues_are_initialized =
false;
3520template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
3524 eigenvalues_are_initialized =
false;
3525 theta = delta = 1.0;
3526 matrix_ptr =
nullptr;
3528 VectorType empty_vector;
3529 solution_old.reinit(empty_vector);
3530 temp_vector1.reinit(empty_vector);
3531 temp_vector2.reinit(empty_vector);
3533 data.preconditioner.reset();
3538template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
3541 PreconditionerType>::EigenvalueInformation
3549 EigenvalueInformation info{};
3551 solution_old.reinit(src);
3552 temp_vector1.reinit(src,
true);
3554 if (data.eig_cg_n_iterations > 0)
3556 Assert(data.eig_cg_n_iterations > 2,
3558 "Need to set at least two iterations to find eigenvalues."));
3560 internal::PreconditionChebyshevImplementation::EigenvalueTracker
3566 internal::PreconditionChebyshevImplementation::set_initial_guess(
3568 data.constraints.set_zero(temp_vector1);
3570 if (data.eigenvalue_algorithm ==
3571 AdditionalData::EigenvalueAlgorithm::lanczos)
3580 solver.connect_eigenvalues_slot(
3581 [&eigenvalue_tracker](
const std::vector<double> &
eigenvalues) {
3585 solver.solve(*matrix_ptr,
3588 *data.preconditioner);
3590 info.cg_iterations = control.last_step();
3592 else if (data.eigenvalue_algorithm ==
3593 AdditionalData::EigenvalueAlgorithm::power_iteration)
3596 ExcMessage(
"Cannot estimate the minimal eigenvalue with the "
3597 "power iteration"));
3599 eigenvalue_tracker.values.push_back(
3600 internal::PreconditionChebyshevImplementation::power_iteration(
3603 *data.preconditioner,
3604 data.eig_cg_n_iterations));
3610 if (eigenvalue_tracker.values.empty())
3611 info.min_eigenvalue_estimate = info.max_eigenvalue_estimate = 1.;
3614 info.min_eigenvalue_estimate = eigenvalue_tracker.values.front();
3618 info.max_eigenvalue_estimate = 1.2 * eigenvalue_tracker.values.back();
3623 info.max_eigenvalue_estimate = data.max_eigenvalue;
3624 info.min_eigenvalue_estimate = data.max_eigenvalue / data.smoothing_range;
3627 const double alpha = (data.smoothing_range > 1. ?
3628 info.max_eigenvalue_estimate / data.smoothing_range :
3629 std::min(0.9 * info.max_eigenvalue_estimate,
3630 info.min_eigenvalue_estimate));
3639 const double actual_range = info.max_eigenvalue_estimate / alpha;
3640 const double sigma = (1. -
std::sqrt(1. / actual_range)) /
3642 const double eps = data.smoothing_range;
3647 1 +
static_cast<unsigned int>(
3652 info.degree = data.degree;
3657 (data.polynomial_type == AdditionalData::PolynomialType::fourth_kind) ?
3658 (info.max_eigenvalue_estimate) :
3659 ((info.max_eigenvalue_estimate - alpha) * 0.5);
3662 ->theta = (info.max_eigenvalue_estimate + alpha) * 0.5;
3666 using NumberType =
typename VectorType::value_type;
3667 if (std::is_same<PreconditionerType,
3670 ((std::is_same<VectorType,
3674 (
std::is_same<VectorType,
3678 temp_vector2.
reinit(src, true);
3681 VectorType empty_vector;
3682 temp_vector2.reinit(empty_vector);
3687 ->eigenvalues_are_initialized =
true;
3694template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
3697 VectorType & solution,
3698 const VectorType &rhs)
const
3700 std::lock_guard<std::mutex> lock(mutex);
3701 if (eigenvalues_are_initialized ==
false)
3702 estimate_eigenvalues(rhs);
3704 internal::PreconditionChebyshevImplementation::vmult_and_update(
3706 *data.preconditioner,
3710 (data.polynomial_type == AdditionalData::PolynomialType::fourth_kind) ?
3711 (4. / (3. * delta)) :
3720 if (data.degree < 2 ||
std::abs(delta) < 1e-40)
3723 double rhok = delta / theta, sigma = theta / delta;
3724 for (
unsigned int k = 0; k < data.degree - 1; ++k)
3726 double factor1 = 0.0;
3727 double factor2 = 0.0;
3729 if (data.polynomial_type == AdditionalData::PolynomialType::fourth_kind)
3731 factor1 = (2 * k + 1.) / (2 * k + 5.);
3732 factor2 = (8 * k + 12.) / (delta * (2 * k + 5.));
3736 const double rhokp = 1. / (2. * sigma - rhok);
3737 factor1 = rhokp * rhok;
3738 factor2 = 2. * rhokp / delta;
3742 internal::PreconditionChebyshevImplementation::vmult_and_update(
3744 *data.preconditioner,
3758template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
3761 VectorType & solution,
3762 const VectorType &rhs)
const
3764 std::lock_guard<std::mutex> lock(mutex);
3765 if (eigenvalues_are_initialized ==
false)
3766 estimate_eigenvalues(rhs);
3768 internal::PreconditionChebyshevImplementation::vector_updates(
3770 *data.preconditioner,
3773 (data.polynomial_type == AdditionalData::PolynomialType::fourth_kind) ?
3774 (4. / (3. * delta)) :
3781 if (data.degree < 2 ||
std::abs(delta) < 1e-40)
3784 double rhok = delta / theta, sigma = theta / delta;
3785 for (
unsigned int k = 0; k < data.degree - 1; ++k)
3787 double factor1 = 0.0;
3788 double factor2 = 0.0;
3790 if (data.polynomial_type == AdditionalData::PolynomialType::fourth_kind)
3792 factor1 = (2 * k + 1.) / (2 * k + 5.);
3793 factor2 = (8 * k + 12.) / (delta * (2 * k + 5.));
3797 const double rhokp = 1. / (2. * sigma - rhok);
3798 factor1 = rhokp * rhok;
3799 factor2 = 2. * rhokp / delta;
3803 matrix_ptr->Tvmult(temp_vector1, solution);
3804 internal::PreconditionChebyshevImplementation::vector_updates(
3806 *data.preconditioner,
3819template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
3822 VectorType & solution,
3823 const VectorType &rhs)
const
3825 std::lock_guard<std::mutex> lock(mutex);
3826 if (eigenvalues_are_initialized ==
false)
3827 estimate_eigenvalues(rhs);
3829 internal::PreconditionChebyshevImplementation::vmult_and_update(
3831 *data.preconditioner,
3835 (data.polynomial_type == AdditionalData::PolynomialType::fourth_kind) ?
3836 (4. / (3. * delta)) :
3843 if (data.degree < 2 ||
std::abs(delta) < 1e-40)
3846 double rhok = delta / theta, sigma = theta / delta;
3847 for (
unsigned int k = 0; k < data.degree - 1; ++k)
3849 double factor1 = 0.0;
3850 double factor2 = 0.0;
3852 if (data.polynomial_type == AdditionalData::PolynomialType::fourth_kind)
3854 factor1 = (2 * k + 1.) / (2 * k + 5.);
3855 factor2 = (8 * k + 12.) / (delta * (2 * k + 5.));
3859 const double rhokp = 1. / (2. * sigma - rhok);
3860 factor1 = rhokp * rhok;
3861 factor2 = 2. * rhokp / delta;
3865 internal::PreconditionChebyshevImplementation::vmult_and_update(
3867 *data.preconditioner,
3881template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
3884 VectorType & solution,
3885 const VectorType &rhs)
const
3887 std::lock_guard<std::mutex> lock(mutex);
3888 if (eigenvalues_are_initialized ==
false)
3889 estimate_eigenvalues(rhs);
3891 matrix_ptr->Tvmult(temp_vector1, solution);
3892 internal::PreconditionChebyshevImplementation::vector_updates(
3894 *data.preconditioner,
3897 (data.polynomial_type == AdditionalData::PolynomialType::fourth_kind) ?
3898 (4. / (3. * delta)) :
3905 if (data.degree < 2 ||
std::abs(delta) < 1e-40)
3908 double rhok = delta / theta, sigma = theta / delta;
3909 for (
unsigned int k = 0; k < data.degree - 1; ++k)
3911 double factor1 = 0.0;
3912 double factor2 = 0.0;
3914 if (data.polynomial_type == AdditionalData::PolynomialType::fourth_kind)
3916 factor1 = (2 * k + 1.) / (2 * k + 5.);
3917 factor2 = (8 * k + 12.) / (delta * (2 * k + 5.));
3921 const double rhokp = 1. / (2. * sigma - rhok);
3922 factor1 = rhokp * rhok;
3923 factor2 = 2. * rhokp / delta;
3927 matrix_ptr->Tvmult(temp_vector1, solution);
3928 internal::PreconditionChebyshevImplementation::vector_updates(
3930 *data.preconditioner,
3943template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
3946 PreconditionerType>::size_type
3950 return matrix_ptr->m();
3955template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
3958 PreconditionerType>::size_type
3962 return matrix_ptr->n();
unsigned int n_blocks() const
BlockType & block(const unsigned int i)
VectorType & get_vector()
size_type nth_index_in_set(const size_type local_index) const
Number * get_values() const
void swap(Vector< Number, MemorySpace > &v)
size_type locally_owned_size() const
virtual Number mean_value() const override
virtual void sadd(const Number s, const Number a, const VectorSpaceVector< Number > &V) override
virtual void add(const Number a) override
virtual ::IndexSet locally_owned_elements() const override
void Tvmult(VectorType &dst, const VectorType &src) const
void step(VectorType &dst, const VectorType &src) const
EigenvalueInformation estimate_eigenvalues(const VectorType &src) const
void Tstep(VectorType &dst, const VectorType &src) const
SmartPointer< const MatrixType, PreconditionChebyshev< MatrixType, VectorType, PreconditionerType > > matrix_ptr
void vmult(VectorType &dst, const VectorType &src) const
void initialize(const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
bool eigenvalues_are_initialized
void vmult_add(VectorType &, const VectorType &) const
void vmult(VectorType &, const VectorType &) const
void initialize(const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
void Tvmult(VectorType &, const VectorType &) const
void Tvmult_add(VectorType &, const VectorType &) const
typename BaseClass::AdditionalData AdditionalData
internal::PreconditionRelaxation::PreconditionJacobiImpl< MatrixType > PreconditionerType
void initialize(const MatrixType &A, const AdditionalData ¶meters=AdditionalData())
AdditionalData(const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const typename BaseClass::AdditionalData ¶meters=typename BaseClass::AdditionalData())
BaseClass::AdditionalData parameters
const std::vector< size_type > & inverse_permutation
const std::vector< size_type > & permutation
typename BaseClass::size_type size_type
void initialize(const MatrixType &A, const AdditionalData &additional_data)
internal::PreconditionRelaxation::PreconditionPSORImpl< MatrixType > PreconditionerType
void initialize(const MatrixType &A, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const typename BaseClass::AdditionalData ¶meters=typename BaseClass::AdditionalData())
unsigned int n_iterations
std::shared_ptr< PreconditionerType > preconditioner
AdditionalData(const double relaxation=1., const unsigned int n_iterations=1)
SmartPointer< const MatrixType, PreconditionRelaxation< MatrixType > > A
void Tvmult(VectorType &, const VectorType &) const
std::shared_ptr< PreconditionerType > preconditioner
void step(VectorType &x, const VectorType &rhs) const
unsigned int n_iterations
void initialize(const MatrixType &A, const AdditionalData ¶meters=AdditionalData())
void Tstep(VectorType &x, const VectorType &rhs) const
void vmult(VectorType &, const VectorType &) const
types::global_dof_index size_type
AdditionalData(const double relaxation=1.)
types::global_dof_index size_type
void vmult_add(VectorType &, const VectorType &) const
void initialize(const AdditionalData ¶meters)
void initialize(const MatrixType &matrix, const AdditionalData ¶meters)
void vmult(VectorType &, const VectorType &) const
void Tvmult(VectorType &, const VectorType &) const
void Tvmult_add(VectorType &, const VectorType &) const
internal::PreconditionRelaxation::PreconditionSORImpl< MatrixType > PreconditionerType
typename BaseClass::AdditionalData AdditionalData
void initialize(const MatrixType &A, const AdditionalData ¶meters=AdditionalData())
void initialize(const MatrixType &A, const AdditionalData ¶meters=AdditionalData())
typename BaseClass::AdditionalData AdditionalData
internal::PreconditionRelaxation::PreconditionSSORImpl< MatrixType > PreconditionerType
void(MatrixType::*)(VectorType &, const VectorType &) const function_ptr
const function_ptr precondition
const MatrixType & matrix
void vmult(VectorType &dst, const VectorType &src) const
PreconditionUseMatrix(const MatrixType &M, const function_ptr method)
const_iterator end() const
const_iterator begin() const
void add(const std::vector< size_type > &indices, const std::vector< OtherNumber > &values)
Number mean_value() const
virtual void swap(Vector< Number > &v)
#define DEAL_II_OPENMP_SIMD_PRAGMA
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ matrix
Contents is actually a matrix.
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
unsigned int minimum_parallel_grain_size
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > log(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index
PolynomialType polynomial_type
AdditionalData(const unsigned int degree=1, const double smoothing_range=0., const unsigned int eig_cg_n_iterations=8, const double eig_cg_residual=1e-2, const double max_eigenvalue=1, const EigenvalueAlgorithm eigenvalue_algorithm=EigenvalueAlgorithm::lanczos, const PolynomialType polynomial_type=PolynomialType::first_kind)
std::shared_ptr< PreconditionerType > preconditioner
AdditionalData & operator=(const AdditionalData &other_data)
AffineConstraints< double > constraints
unsigned int eig_cg_n_iterations
EigenvalueAlgorithm eigenvalue_algorithm
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)