16 #ifndef dealii__precondition_h 17 #define dealii__precondition_h 21 #include <deal.II/base/config.h> 22 #include <deal.II/base/smartpointer.h> 23 #include <deal.II/base/utilities.h> 24 #include <deal.II/base/parallel.h> 25 #include <deal.II/base/template_constraints.h> 26 #include <deal.II/base/thread_management.h> 27 #include <deal.II/lac/diagonal_matrix.h> 28 #include <deal.II/lac/solver_cg.h> 29 #include <deal.II/lac/vector_memory.h> 31 DEAL_II_NAMESPACE_OPEN
35 template <
typename number>
class Vector;
41 template <
typename number>
class Vector;
102 template <
typename MatrixType>
104 const AdditionalData &additional_data = AdditionalData());
109 template<
class VectorType>
110 void vmult (VectorType &,
const VectorType &)
const;
116 template<
class VectorType>
117 void Tvmult (VectorType &,
const VectorType &)
const;
122 template<
class VectorType>
123 void vmult_add (VectorType &,
const VectorType &)
const;
129 template<
class VectorType>
130 void Tvmult_add (VectorType &,
const VectorType &)
const;
225 template <
typename MatrixType>
232 template<
class VectorType>
233 void vmult (VectorType &,
const VectorType &)
const;
239 template<
class VectorType>
240 void Tvmult (VectorType &,
const VectorType &)
const;
244 template<
class VectorType>
245 void vmult_add (VectorType &,
const VectorType &)
const;
251 template<
class VectorType>
252 void Tvmult_add (VectorType &,
const VectorType &)
const;
337 template<
typename MatrixType = SparseMatrix<
double>,
class VectorType = Vector<
double> >
344 typedef void ( MatrixType::*
function_ptr)(VectorType &,
const VectorType &)
const;
358 void vmult (VectorType &dst,
359 const VectorType &src)
const;
383 template<
typename MatrixType = SparseMatrix<
double> >
475 template <
typename MatrixType = SparseMatrix<
double> >
487 template<
class VectorType>
488 void vmult (VectorType &,
const VectorType &)
const;
494 template<
class VectorType>
495 void Tvmult (VectorType &,
const VectorType &)
const;
500 template<
class VectorType>
501 void step (VectorType &x,
const VectorType &rhs)
const;
506 template<
class VectorType>
507 void Tstep (VectorType &x,
const VectorType &rhs)
const;
557 template <
typename MatrixType = SparseMatrix<
double> >
569 template<
class VectorType>
570 void vmult (VectorType &,
const VectorType &)
const;
575 template<
class VectorType>
576 void Tvmult (VectorType &,
const VectorType &)
const;
581 template<
class VectorType>
582 void step (VectorType &x,
const VectorType &rhs)
const;
587 template<
class VectorType>
588 void Tstep (VectorType &x,
const VectorType &rhs)
const;
620 template <
typename MatrixType = SparseMatrix<
double> >
651 template<
class VectorType>
652 void vmult (VectorType &,
const VectorType &)
const;
658 template<
class VectorType>
659 void Tvmult (VectorType &,
const VectorType &)
const;
665 template<
class VectorType>
666 void step (VectorType &x,
const VectorType &rhs)
const;
671 template<
class VectorType>
672 void Tstep (VectorType &x,
const VectorType &rhs)
const;
715 template <
typename MatrixType = SparseMatrix<
double> >
792 template<
class VectorType>
793 void vmult (VectorType &,
const VectorType &)
const;
798 template<
class VectorType>
799 void Tvmult (VectorType &,
const VectorType &)
const;
901 template <
typename MatrixType=SparseMatrix<
double>,
typename VectorType=Vector<
double>,
typename PreconditionerType=DiagonalMatrix<VectorType> >
911 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
1002 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
1024 void vmult (VectorType &dst,
1025 const VectorType &src)
const;
1031 void Tvmult (VectorType &dst,
1032 const VectorType &src)
const;
1037 void step (VectorType &dst,
1038 const VectorType &src)
const;
1043 void Tstep (VectorType &dst,
1044 const VectorType &src)
const;
1119 const VectorType &src)
const;
1128 const VectorType &src)
const;
1153 template <
typename MatrixType>
1163 template<
class VectorType>
1172 template<
class VectorType>
1179 template<
class VectorType>
1188 template<
class VectorType>
1214 relaxation(relaxation)
1225 AdditionalData add_data;
1226 relaxation=add_data.relaxation;
1240 template <
typename MatrixType>
1243 (
const MatrixType &matrix,
1253 template<
class VectorType>
1257 #ifdef DEAL_II_WITH_CXX11 1259 std::is_same<size_type, typename VectorType::size_type>::value,
1260 "PreconditionRichardson and VectorType must have the same size_type.");
1261 #endif // DEAL_II_WITH_CXX11 1268 template<
class VectorType>
1272 #ifdef DEAL_II_WITH_CXX11 1274 std::is_same<size_type, typename VectorType::size_type>::value,
1275 "PreconditionRichardson and VectorType must have the same size_type.");
1276 #endif // DEAL_II_WITH_CXX11 1281 template<
class VectorType>
1285 #ifdef DEAL_II_WITH_CXX11 1287 std::is_same<size_type, typename VectorType::size_type>::value,
1288 "PreconditionRichardson and VectorType must have the same size_type.");
1289 #endif // DEAL_II_WITH_CXX11 1296 template<
class VectorType>
1300 #ifdef DEAL_II_WITH_CXX11 1302 std::is_same<size_type, typename VectorType::size_type>::value,
1303 "PreconditionRichardson and VectorType must have the same size_type.");
1304 #endif // DEAL_II_WITH_CXX11 1325 template <
typename MatrixType>
1328 const AdditionalData ¶meters)
1331 relaxation = parameters.relaxation;
1335 template <
typename MatrixType>
1342 template <
typename MatrixType>
1350 template <
typename MatrixType>
1360 template <
typename MatrixType>
1361 template<
class VectorType>
1365 #ifdef DEAL_II_WITH_CXX11 1368 "PreconditionJacobi and VectorType must have the same size_type.");
1369 #endif // DEAL_II_WITH_CXX11 1372 this->A->precondition_Jacobi (dst, src, this->relaxation);
1377 template <
typename MatrixType>
1378 template<
class VectorType>
1382 #ifdef DEAL_II_WITH_CXX11 1385 "PreconditionJacobi and VectorType must have the same size_type.");
1386 #endif // DEAL_II_WITH_CXX11 1389 this->A->precondition_Jacobi (dst, src, this->relaxation);
1394 template <
typename MatrixType>
1395 template<
class VectorType>
1399 #ifdef DEAL_II_WITH_CXX11 1402 "PreconditionJacobi and VectorType must have the same size_type.");
1403 #endif // DEAL_II_WITH_CXX11 1406 this->A->Jacobi_step (dst, src, this->relaxation);
1411 template <
typename MatrixType>
1412 template<
class VectorType>
1416 #ifdef DEAL_II_WITH_CXX11 1419 "PreconditionJacobi and VectorType must have the same size_type.");
1420 #endif // DEAL_II_WITH_CXX11 1429 template <
typename MatrixType>
1430 template<
class VectorType>
1434 #ifdef DEAL_II_WITH_CXX11 1437 "PreconditionSOR and VectorType must have the same size_type.");
1438 #endif // DEAL_II_WITH_CXX11 1441 this->A->precondition_SOR (dst, src, this->relaxation);
1446 template <
typename MatrixType>
1447 template<
class VectorType>
1451 #ifdef DEAL_II_WITH_CXX11 1454 "PreconditionSOR and VectorType must have the same size_type.");
1455 #endif // DEAL_II_WITH_CXX11 1458 this->A->precondition_TSOR (dst, src, this->relaxation);
1463 template <
typename MatrixType>
1464 template<
class VectorType>
1468 #ifdef DEAL_II_WITH_CXX11 1471 "PreconditionSOR and VectorType must have the same size_type.");
1472 #endif // DEAL_II_WITH_CXX11 1475 this->A->SOR_step (dst, src, this->relaxation);
1480 template <
typename MatrixType>
1481 template<
class VectorType>
1485 #ifdef DEAL_II_WITH_CXX11 1488 "PreconditionSOR and VectorType must have the same size_type.");
1489 #endif // DEAL_II_WITH_CXX11 1492 this->A->TSOR_step (dst, src, this->relaxation);
1499 template <
typename MatrixType>
1502 const typename BaseClass::AdditionalData ¶meters)
1515 pos_right_of_diagonal.resize(n, static_cast<std::size_t>(-1));
1516 for (size_type row=0; row<n; ++row)
1523 it = mat->
begin(row)+1;
1524 for ( ; it < mat->
end(row); ++it)
1525 if (it->column() > row)
1527 pos_right_of_diagonal[row] = it - mat->
begin();
1533 template <
typename MatrixType>
1534 template<
class VectorType>
1538 #ifdef DEAL_II_WITH_CXX11 1541 "PreconditionSSOR and VectorType must have the same size_type.");
1542 #endif // DEAL_II_WITH_CXX11 1545 this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
1550 template <
typename MatrixType>
1551 template<
class VectorType>
1555 #ifdef DEAL_II_WITH_CXX11 1558 "PreconditionSSOR and VectorType must have the same size_type.");
1559 #endif // DEAL_II_WITH_CXX11 1562 this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
1567 template <
typename MatrixType>
1568 template<
class VectorType>
1572 #ifdef DEAL_II_WITH_CXX11 1575 "PreconditionSSOR and VectorType must have the same size_type.");
1576 #endif // DEAL_II_WITH_CXX11 1579 this->A->SSOR_step (dst, src, this->relaxation);
1584 template <
typename MatrixType>
1585 template<
class VectorType>
1589 #ifdef DEAL_II_WITH_CXX11 1592 "PreconditionSSOR and VectorType must have the same size_type.");
1593 #endif // DEAL_II_WITH_CXX11 1602 template <
typename MatrixType>
1605 (
const MatrixType &rA,
1606 const std::vector<size_type> &p,
1607 const std::vector<size_type> &ip,
1611 inverse_permutation = &ip;
1616 template <
typename MatrixType>
1619 const AdditionalData &additional_data)
1622 additional_data.permutation,
1623 additional_data.inverse_permutation,
1624 additional_data.parameters);
1628 template <
typename MatrixType>
1629 template <
typename VectorType>
1633 #ifdef DEAL_II_WITH_CXX11 1636 "PreconditionPSOR and VectorType must have the same size_type.");
1637 #endif // DEAL_II_WITH_CXX11 1641 this->A->PSOR (dst, *permutation, *inverse_permutation, this->relaxation);
1646 template <
typename MatrixType>
1647 template<
class VectorType>
1651 #ifdef DEAL_II_WITH_CXX11 1654 "PreconditionPSOR and VectorType must have the same size_type.");
1655 #endif // DEAL_II_WITH_CXX11 1659 this->A->TPSOR (dst, *permutation, *inverse_permutation, this->relaxation);
1662 template <
typename MatrixType>
1664 (
const std::vector<size_type> &permutation,
1665 const std::vector<size_type> &inverse_permutation,
1668 permutation(permutation),
1669 inverse_permutation(inverse_permutation),
1670 parameters(parameters)
1679 template<
typename MatrixType,
class VectorType>
1681 const function_ptr method)
1683 matrix(M), precondition(method)
1688 template<
typename MatrixType,
class VectorType>
1691 const VectorType &src)
const 1693 (matrix.*precondition)(dst, src);
1698 template<
typename MatrixType>
1703 relaxation (relaxation)
1720 template <
typename VectorType,
typename PreconditionerType>
1723 vector_updates (
const VectorType &src,
1724 const PreconditionerType &preconditioner,
1725 const bool start_zero,
1726 const double factor1,
1727 const double factor2,
1728 VectorType &update1,
1729 VectorType &update2,
1730 VectorType &update3,
1735 update1.equ (factor2, src);
1736 preconditioner.
vmult(dst, update1);
1737 update1.equ(-1.,dst);
1742 preconditioner.vmult(update3, update2);
1745 update1.equ(factor2, update2);
1747 update1.sadd(factor1, factor2, update2);
1755 template <
typename Number>
1756 struct VectorUpdater
1758 VectorUpdater (
const Number *src,
1759 const Number *matrix_diagonal_inverse,
1760 const bool start_zero,
1761 const Number factor1,
1762 const Number factor2,
1768 matrix_diagonal_inverse (matrix_diagonal_inverse),
1769 do_startup (factor1 == Number()),
1770 start_zero (start_zero),
1779 apply_to_subrange (
const std::size_t begin,
1780 const std::size_t end)
const 1786 const Number factor1 = this->factor1;
1787 const Number factor2 = this->factor2;
1791 DEAL_II_OPENMP_SIMD_PRAGMA
1792 for (std::size_t i=begin; i<end; ++i)
1794 dst[i] = factor2 * src[i] * matrix_diagonal_inverse[i];
1795 update1[i] = -dst[i];
1798 DEAL_II_OPENMP_SIMD_PRAGMA
1799 for (std::size_t i=begin; i<end; ++i)
1801 update1[i] = ((update2[i]-src[i]) *
1802 factor2*matrix_diagonal_inverse[i]);
1803 dst[i] -= update1[i];
1807 DEAL_II_OPENMP_SIMD_PRAGMA
1808 for (std::size_t i=begin; i<end; ++i)
1810 const Number update =
1811 factor1 * update1[i] + factor2 *
1812 ((update2[i] - src[i]) * matrix_diagonal_inverse[i]);
1813 update1[i] = update;
1819 const Number *matrix_diagonal_inverse;
1820 const bool do_startup;
1821 const bool start_zero;
1822 const Number factor1;
1823 const Number factor2;
1824 mutable Number *update1;
1825 mutable Number *update2;
1826 mutable Number *dst;
1829 template<
typename Number>
1832 VectorUpdatesRange(
const VectorUpdater<Number> &updater,
1833 const std::size_t size)
1837 if (size < internal::Vector::minimum_parallel_grain_size)
1838 apply_to_subrange (0, size);
1840 apply_parallel (0, size,
1841 internal::Vector::minimum_parallel_grain_size);
1844 ~VectorUpdatesRange() {}
1847 apply_to_subrange (
const std::size_t begin,
1848 const std::size_t end)
const 1850 updater.apply_to_subrange(begin, end);
1853 const VectorUpdater<Number> &updater;
1857 template <
typename Number>
1860 vector_updates (const ::Vector<Number> &src,
1862 const bool start_zero,
1863 const double factor1,
1864 const double factor2,
1870 VectorUpdater<Number> upd(src.begin(), jacobi.get_vector().begin(),
1871 start_zero, factor1, factor2,
1873 VectorUpdatesRange<Number>(upd, src.size());
1877 template <
typename Number>
1882 const bool start_zero,
1883 const double factor1,
1884 const double factor2,
1890 VectorUpdater<Number> upd(src.
begin(), jacobi.get_vector().begin(),
1891 start_zero, factor1, factor2,
1893 VectorUpdatesRange<Number>(upd, src.
local_size());
1896 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
1899 initialize_preconditioner(
const MatrixType &matrix,
1900 std_cxx11::shared_ptr<PreconditionerType> &preconditioner,
1904 (void)preconditioner;
1908 template <
typename MatrixType,
typename VectorType>
1911 initialize_preconditioner(
const MatrixType &matrix,
1913 VectorType &diagonal_inverse)
1915 if (preconditioner.get() == NULL ||
1916 preconditioner->m() != matrix.m())
1918 if (preconditioner.get() == NULL)
1921 Assert(preconditioner->m() == 0,
1922 ExcMessage(
"Preconditioner appears to be initialized but not sized correctly"));
1926 preconditioner->reinit(diagonal_inverse);
1928 VectorType empty_vector;
1929 diagonal_inverse.reinit(empty_vector);
1933 if (preconditioner->m() != matrix.m())
1935 preconditioner->get_vector().reinit(matrix.m());
1936 for (
typename VectorType::size_type i=0; i<matrix.m(); ++i)
1937 preconditioner->get_vector()(i) = 1./matrix.el(i,i);
1942 template <
typename VectorType>
1943 void set_initial_guess(VectorType &vector)
1945 vector = 1./std::sqrt(static_cast<double>(vector.size()));
1946 if (vector.locally_owned_elements().is_element(0))
1950 template <
typename Number>
1957 for (
unsigned int i=0; i<vector.
size(); ++i)
1961 template <
typename Number>
1968 for (
unsigned int i=0; i<vector.
local_size(); ++i)
1972 struct EigenvalueTracker
1975 void slot(
const std::vector<double> &eigenvalues)
1977 values = eigenvalues;
1980 std::vector<double> values;
1988 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
1990 template <
typename MatrixType,
class VectorType,
typename PreconditionerType>
1994 const double smoothing_range,
1995 const bool nonzero_starting,
1996 const unsigned int eig_cg_n_iterations,
1997 const double eig_cg_residual,
1998 const double max_eigenvalue)
2001 smoothing_range (smoothing_range),
2002 nonzero_starting (nonzero_starting),
2003 eig_cg_n_iterations (eig_cg_n_iterations),
2004 eig_cg_residual (eig_cg_residual),
2005 max_eigenvalue (max_eigenvalue)
2008 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
2011 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2017 is_initialized (false)
2019 #ifdef DEAL_II_WITH_CXX11 2021 std::is_same<size_type, typename VectorType::size_type>::value,
2022 "PreconditionChebyshev and VectorType must have the same size_type.");
2023 #endif // DEAL_II_WITH_CXX11 2028 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
2030 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2034 (
const MatrixType &matrix,
2035 const AdditionalData &additional_data)
2037 matrix_ptr = &matrix;
2038 data = additional_data;
2039 internal::PreconditionChebyshev::initialize_preconditioner(matrix,
2040 data.preconditioner,
2041 data.matrix_diagonal_inverse);
2042 is_initialized =
false;
2047 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2052 is_initialized =
false;
2053 theta = delta = 1.0;
2056 VectorType empty_vector;
2057 data.matrix_diagonal_inverse.reinit(empty_vector);
2058 update1.
reinit(empty_vector);
2059 update2.
reinit(empty_vector);
2060 update3.reinit(empty_vector);
2062 data.preconditioner.reset();
2065 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
2069 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2073 (
const VectorType &src)
const 2079 update2.
reinit(src,
true);
2084 double max_eigenvalue, min_eigenvalue;
2085 if (data.eig_cg_n_iterations > 0)
2087 Assert (data.eig_cg_n_iterations > 2,
2088 ExcMessage (
"Need to set at least two iterations to find eigenvalues."));
2092 std::sqrt(std::numeric_limits<typename VectorType::value_type>::epsilon()),
2093 1e-10,
false,
false);
2095 internal::PreconditionChebyshev::EigenvalueTracker eigenvalue_tracker;
2097 solver.connect_eigenvalues_slot(std_cxx11::bind(&internal::PreconditionChebyshev::EigenvalueTracker::slot,
2098 &eigenvalue_tracker,
2103 internal::PreconditionChebyshev::set_initial_guess(update2);
2107 solver.solve(*matrix_ptr, update1, update2, *data.preconditioner);
2114 if (eigenvalue_tracker.values.empty())
2115 min_eigenvalue = max_eigenvalue = 1;
2118 min_eigenvalue = eigenvalue_tracker.values.front();
2122 max_eigenvalue = 1.2*eigenvalue_tracker.values.back();
2127 max_eigenvalue = data.max_eigenvalue;
2128 min_eigenvalue = data.max_eigenvalue/data.smoothing_range;
2131 const double alpha = (data.smoothing_range > 1. ?
2132 max_eigenvalue / data.smoothing_range :
2133 std::min(0.9*max_eigenvalue, min_eigenvalue));
2142 const double actual_range = max_eigenvalue / alpha;
2143 const double sigma = (1.-std::sqrt(1./actual_range))/(1.+std::sqrt(1./actual_range));
2144 const double eps = data.smoothing_range;
2146 = 1+std::log(1./eps+std::sqrt(1./eps/eps-1))/std::log(1./sigma);
2159 update3.reinit (src,
true);
2166 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2171 const VectorType &src)
const 2175 if (std::abs(delta) < 1e-40)
2178 double rhok = delta / theta, sigma = theta / delta;
2179 for (
unsigned int k=0; k<data.degree; ++k)
2181 matrix_ptr->vmult (update2, dst);
2182 const double rhokp = 1./(2.*sigma-rhok);
2183 const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/delta;
2185 internal::PreconditionChebyshev::vector_updates
2186 (src, *data.preconditioner,
false, factor1, factor2, update1, update2, update3, dst);
2192 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2197 const VectorType &src)
const 2199 double rhok = delta / theta, sigma = theta / delta;
2200 for (
unsigned int k=0; k<data.degree; ++k)
2202 matrix_ptr->Tvmult (update2, dst);
2203 const double rhokp = 1./(2.*sigma-rhok);
2204 const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/delta;
2206 internal::PreconditionChebyshev::vector_updates
2207 (src, *data.preconditioner,
false, factor1, factor2, update1, update2, update3, dst);
2213 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2218 const VectorType &src)
const 2221 if (is_initialized ==
false)
2222 estimate_eigenvalues(src);
2224 internal::PreconditionChebyshev::vector_updates
2225 (src, *data.preconditioner,
true, 0., 1./theta, update1, update2, update3, dst);
2227 do_chebyshev_loop(dst, src);
2232 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2237 const VectorType &src)
const 2240 if (is_initialized ==
false)
2241 estimate_eigenvalues(src);
2243 internal::PreconditionChebyshev::vector_updates
2244 (src, *data.preconditioner,
true, 0., 1./theta, update1, update2, update3, dst);
2246 do_transpose_chebyshev_loop(dst, src);
2251 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2256 const VectorType &src)
const 2259 if (is_initialized ==
false)
2260 estimate_eigenvalues(src);
2262 if (!dst.all_zero())
2264 matrix_ptr->vmult (update2, dst);
2265 internal::PreconditionChebyshev::vector_updates
2266 (src, *data.preconditioner,
false, 0., 1./theta, update1, update2, update3, dst);
2269 internal::PreconditionChebyshev::vector_updates
2270 (src, *data.preconditioner,
true, 0., 1./theta, update1, update2, update3, dst);
2272 do_chebyshev_loop(dst, src);
2277 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2282 const VectorType &src)
const 2285 if (is_initialized ==
false)
2286 estimate_eigenvalues(src);
2288 if (!dst.all_zero())
2290 matrix_ptr->Tvmult (update2, dst);
2291 internal::PreconditionChebyshev::vector_updates
2292 (src, *data.preconditioner,
false, 0., 1./theta, update1, update2, update3, dst);
2295 internal::PreconditionChebyshev::vector_updates
2296 (src, *data.preconditioner,
true, 0., 1./theta, update1, update2, update3, dst);
2298 do_transpose_chebyshev_loop(dst, src);
2303 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2309 return matrix_ptr->m();
2314 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2320 return matrix_ptr->n();
2325 DEAL_II_NAMESPACE_CLOSE
types::global_dof_index size_type
void reinit(const size_type size, const bool omit_zeroing_entries=false)
void Tvmult(VectorType &, const VectorType &) const
static const unsigned int invalid_unsigned_int
void initialize(const MatrixType &A, const typename BaseClass::AdditionalData ¶meters=typename BaseClass::AdditionalData())
const_iterator end() const
void Tstep(VectorType &x, const VectorType &rhs) const
void vmult(VectorType &dst, const VectorType &src) const
MatrixType::size_type size_type
void Tvmult_add(VectorType &, const VectorType &) const
PreconditionRelaxation< MatrixType >::AdditionalData parameters
void do_chebyshev_loop(VectorType &dst, const VectorType &src) const
void vmult(VectorType &, const VectorType &) const
PreconditionRelaxation< MatrixType >::AdditionalData AdditionalData
const MatrixType & matrix
AdditionalData(const double relaxation=1.)
static ::ExceptionBase & ExcNotInitialized()
void initialize(const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
SmartPointer< const MatrixType, PreconditionChebyshev< MatrixType, VectorType, PreconditionerType > > matrix_ptr
#define AssertThrow(cond, exc)
AdditionalData(const double relaxation=1.)
types::global_dof_index size_type
const std::vector< size_type > * permutation
const function_ptr precondition
void vmult(VectorType &, const VectorType &) const
void initialize(const AdditionalData ¶meters)
void Tvmult(VectorType &, const VectorType &) const
static ::ExceptionBase & ExcMessage(std::string arg1)
void Tvmult(VectorType &dst, const VectorType &src) const
types::global_dof_index size_type
void vmult(VectorType &, const VectorType &) const
unsigned int global_dof_index
std_cxx11::shared_ptr< PreconditionerType > preconditioner
const std::vector< size_type > & inverse_permutation
#define Assert(cond, exc)
void step(VectorType &x, const VectorType &rhs) const
void(MatrixType::* function_ptr)(VectorType &, const VectorType &) const
MatrixType::size_type size_type
types::global_dof_index size_type
void vmult_add(VectorType &, const VectorType &) const
size_type local_size() const
void Tvmult(VectorType &, const VectorType &) const
void initialize(const MatrixType &A, const AdditionalData ¶meters=AdditionalData())
void Tvmult(VectorType &, const VectorType &) const
void estimate_eigenvalues(const VectorType &src) const
AdditionalData(const unsigned int degree=0, const double smoothing_range=0., const bool nonzero_starting=false, const unsigned int eig_cg_n_iterations=8, const double eig_cg_residual=1e-2, const double max_eigenvalue=1)
MatrixType::size_type size_type
PreconditionRelaxation< MatrixType > BaseClass
void Tstep(VectorType &x, const VectorType &rhs) const
void vmult(VectorType &, const VectorType &) const
void Tstep(VectorType &x, const VectorType &rhs) const
void vmult(VectorType &dst, const VectorType &src) const
AdditionalData(const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const typename PreconditionRelaxation< MatrixType >::AdditionalData ¶meters=typename PreconditionRelaxation< MatrixType >::AdditionalData())
VectorType matrix_diagonal_inverse
void Tvmult_add(VectorType &, const VectorType &) const
void Tvmult(VectorType &, const VectorType &) const
void vmult(VectorType &, const VectorType &) const
SmartPointer< const MatrixType, PreconditionRelaxation< MatrixType > > A
void Tvmult(VectorType &, const VectorType &) const
void step(VectorType &x, const VectorType &rhs) const
void initialize(const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
void do_transpose_chebyshev_loop(VectorType &dst, const VectorType &src) const
const std::vector< size_type > * inverse_permutation
PreconditionRelaxation< MatrixType >::AdditionalData AdditionalData
void Tstep(VectorType &dst, const VectorType &src) const
void step(VectorType &dst, const VectorType &src) const
void step(VectorType &x, const VectorType &rhs) const
Number local_element(const size_type local_index) const
unsigned int eig_cg_n_iterations
const_iterator begin() const
void vmult(VectorType &, const VectorType &) const
const std::vector< size_type > & permutation
PreconditionUseMatrix(const MatrixType &M, const function_ptr method)
void vmult_add(VectorType &, const VectorType &) const
PreconditionRelaxation< MatrixType >::AdditionalData AdditionalData
std::vector< std::size_t > pos_right_of_diagonal
static ::ExceptionBase & ExcInternalError()
void initialize(const MatrixType &A, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const typename PreconditionRelaxation< MatrixType >::AdditionalData ¶meters=typename PreconditionRelaxation< MatrixType >::AdditionalData())