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>
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> >
1022 void vmult (VectorType &dst,
1023 const VectorType &src)
const;
1029 void Tvmult (VectorType &dst,
1030 const VectorType &src)
const;
1035 void step (VectorType &dst,
1036 const VectorType &src)
const;
1041 void Tstep (VectorType &dst,
1042 const VectorType &src)
const;
1117 const VectorType &src)
const;
1126 const VectorType &src)
const;
1151 template <
typename MatrixType>
1161 template <
class VectorType>
1170 template <
class VectorType>
1177 template <
class VectorType>
1186 template <
class VectorType>
1212 relaxation(relaxation)
1223 AdditionalData add_data;
1224 relaxation=add_data.relaxation;
1238 template <
typename MatrixType>
1241 (
const MatrixType &matrix,
1251 template <
class VectorType>
1256 std::is_same<size_type, typename VectorType::size_type>::value,
1257 "PreconditionRichardson and VectorType must have the same size_type.");
1264 template <
class VectorType>
1269 std::is_same<size_type, typename VectorType::size_type>::value,
1270 "PreconditionRichardson and VectorType must have the same size_type.");
1275 template <
class VectorType>
1280 std::is_same<size_type, typename VectorType::size_type>::value,
1281 "PreconditionRichardson and VectorType must have the same size_type.");
1288 template <
class VectorType>
1293 std::is_same<size_type, typename VectorType::size_type>::value,
1294 "PreconditionRichardson and VectorType must have the same size_type.");
1315 template <
typename MatrixType>
1318 const AdditionalData ¶meters)
1321 relaxation = parameters.relaxation;
1325 template <
typename MatrixType>
1332 template <
typename MatrixType>
1340 template <
typename MatrixType>
1350 template <
typename MatrixType>
1351 template <
class VectorType>
1357 "PreconditionJacobi and VectorType must have the same size_type.");
1360 this->A->precondition_Jacobi (dst, src, this->relaxation);
1365 template <
typename MatrixType>
1366 template <
class VectorType>
1372 "PreconditionJacobi and VectorType must have the same size_type.");
1375 this->A->precondition_Jacobi (dst, src, this->relaxation);
1380 template <
typename MatrixType>
1381 template <
class VectorType>
1387 "PreconditionJacobi and VectorType must have the same size_type.");
1390 this->A->Jacobi_step (dst, src, this->relaxation);
1395 template <
typename MatrixType>
1396 template <
class VectorType>
1402 "PreconditionJacobi and VectorType must have the same size_type.");
1411 template <
typename MatrixType>
1412 template <
class VectorType>
1418 "PreconditionSOR and VectorType must have the same size_type.");
1421 this->A->precondition_SOR (dst, src, this->relaxation);
1426 template <
typename MatrixType>
1427 template <
class VectorType>
1433 "PreconditionSOR and VectorType must have the same size_type.");
1436 this->A->precondition_TSOR (dst, src, this->relaxation);
1441 template <
typename MatrixType>
1442 template <
class VectorType>
1448 "PreconditionSOR and VectorType must have the same size_type.");
1451 this->A->SOR_step (dst, src, this->relaxation);
1456 template <
typename MatrixType>
1457 template <
class VectorType>
1463 "PreconditionSOR and VectorType must have the same size_type.");
1466 this->A->TSOR_step (dst, src, this->relaxation);
1473 template <
typename MatrixType>
1476 const typename BaseClass::AdditionalData ¶meters)
1489 pos_right_of_diagonal.resize(n, static_cast<std::size_t>(-1));
1490 for (size_type row=0; row<n; ++row)
1497 it = mat->
begin(row)+1;
1498 for ( ; it < mat->
end(row); ++it)
1499 if (it->column() > row)
1501 pos_right_of_diagonal[row] = it - mat->
begin();
1507 template <
typename MatrixType>
1508 template <
class VectorType>
1514 "PreconditionSSOR and VectorType must have the same size_type.");
1517 this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
1522 template <
typename MatrixType>
1523 template <
class VectorType>
1529 "PreconditionSSOR and VectorType must have the same size_type.");
1532 this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
1537 template <
typename MatrixType>
1538 template <
class VectorType>
1544 "PreconditionSSOR and VectorType must have the same size_type.");
1547 this->A->SSOR_step (dst, src, this->relaxation);
1552 template <
typename MatrixType>
1553 template <
class VectorType>
1559 "PreconditionSSOR and VectorType must have the same size_type.");
1568 template <
typename MatrixType>
1571 (
const MatrixType &rA,
1572 const std::vector<size_type> &p,
1573 const std::vector<size_type> &ip,
1577 inverse_permutation = &ip;
1582 template <
typename MatrixType>
1585 const AdditionalData &additional_data)
1588 additional_data.permutation,
1589 additional_data.inverse_permutation,
1590 additional_data.parameters);
1594 template <
typename MatrixType>
1595 template <
typename VectorType>
1601 "PreconditionPSOR and VectorType must have the same size_type.");
1605 this->A->PSOR (dst, *permutation, *inverse_permutation, this->relaxation);
1610 template <
typename MatrixType>
1611 template <
class VectorType>
1617 "PreconditionPSOR and VectorType must have the same size_type.");
1621 this->A->TPSOR (dst, *permutation, *inverse_permutation, this->relaxation);
1624 template <
typename MatrixType>
1626 (
const std::vector<size_type> &permutation,
1627 const std::vector<size_type> &inverse_permutation,
1630 permutation(permutation),
1631 inverse_permutation(inverse_permutation),
1632 parameters(parameters)
1641 template <
typename MatrixType,
class VectorType>
1643 const function_ptr method)
1645 matrix(M), precondition(method)
1650 template <
typename MatrixType,
class VectorType>
1653 const VectorType &src)
const 1655 (
matrix.*precondition)(dst, src);
1660 template <
typename MatrixType>
1665 relaxation (relaxation)
1674 namespace PreconditionChebyshevImplementation
1682 template <
typename VectorType,
typename PreconditionerType>
1685 vector_updates (
const VectorType &src,
1686 const PreconditionerType &preconditioner,
1687 const bool start_zero,
1688 const double factor1,
1689 const double factor2,
1690 VectorType &update1,
1691 VectorType &update2,
1692 VectorType &update3,
1697 update1.equ (factor2, src);
1698 preconditioner.vmult(dst, update1);
1699 update1.equ(-1.,dst);
1704 preconditioner.vmult(update3, update2);
1707 update1.equ(factor2, update2);
1709 update1.sadd(factor1, factor2, update2);
1717 template <
typename Number>
1718 struct VectorUpdater
1720 VectorUpdater (
const Number *src,
1721 const Number *matrix_diagonal_inverse,
1722 const bool start_zero,
1723 const Number factor1,
1724 const Number factor2,
1730 matrix_diagonal_inverse (matrix_diagonal_inverse),
1731 do_startup (factor1 == Number()),
1732 start_zero (start_zero),
1741 apply_to_subrange (
const std::size_t begin,
1742 const std::size_t end)
const 1748 const Number factor1 = this->factor1;
1749 const Number factor2 = this->factor2;
1753 DEAL_II_OPENMP_SIMD_PRAGMA
1754 for (std::size_t i=begin; i<end; ++i)
1756 dst[i] = factor2 * src[i] * matrix_diagonal_inverse[i];
1757 update1[i] = -dst[i];
1760 DEAL_II_OPENMP_SIMD_PRAGMA
1761 for (std::size_t i=begin; i<end; ++i)
1763 update1[i] = ((update2[i]-src[i]) *
1764 factor2*matrix_diagonal_inverse[i]);
1765 dst[i] -= update1[i];
1769 DEAL_II_OPENMP_SIMD_PRAGMA
1770 for (std::size_t i=begin; i<end; ++i)
1772 const Number update =
1773 factor1 * update1[i] + factor2 *
1774 ((update2[i] - src[i]) * matrix_diagonal_inverse[i]);
1775 update1[i] = update;
1781 const Number *matrix_diagonal_inverse;
1782 const bool do_startup;
1783 const bool start_zero;
1784 const Number factor1;
1785 const Number factor2;
1786 mutable Number *update1;
1787 mutable Number *update2;
1788 mutable Number *dst;
1791 template <
typename Number>
1794 VectorUpdatesRange(
const VectorUpdater<Number> &updater,
1795 const std::size_t size)
1799 if (size < internal::VectorImplementation::minimum_parallel_grain_size)
1800 apply_to_subrange (0, size);
1802 apply_parallel (0, size,
1803 internal::VectorImplementation::minimum_parallel_grain_size);
1806 ~VectorUpdatesRange() =
default;
1809 apply_to_subrange (
const std::size_t begin,
1810 const std::size_t end)
const 1812 updater.apply_to_subrange(begin, end);
1815 const VectorUpdater<Number> &updater;
1819 template <
typename Number>
1822 vector_updates (const ::Vector<Number> &src,
1824 const bool start_zero,
1825 const double factor1,
1826 const double factor2,
1827 ::Vector<Number> &update1,
1828 ::Vector<Number> &update2,
1830 ::Vector<Number> &dst)
1832 VectorUpdater<Number> upd(src.begin(),
jacobi.get_vector().begin(),
1833 start_zero, factor1, factor2,
1834 update1.begin(), update2.begin(), dst.begin());
1835 VectorUpdatesRange<Number>(upd, src.size());
1839 template <
typename Number>
1844 const bool start_zero,
1845 const double factor1,
1846 const double factor2,
1852 VectorUpdater<Number> upd(src.
begin(),
jacobi.get_vector().begin(),
1853 start_zero, factor1, factor2,
1855 VectorUpdatesRange<Number>(upd, src.
local_size());
1858 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
1861 initialize_preconditioner(
const MatrixType &matrix,
1862 std::shared_ptr<PreconditionerType> &preconditioner,
1866 (void)preconditioner;
1870 template <
typename MatrixType,
typename VectorType>
1873 initialize_preconditioner(
const MatrixType &matrix,
1875 VectorType &diagonal_inverse)
1877 if (preconditioner.get() ==
nullptr ||
1878 preconditioner->m() !=
matrix.m())
1880 if (preconditioner.get() ==
nullptr)
1883 Assert(preconditioner->m() == 0,
1884 ExcMessage(
"Preconditioner appears to be initialized but not sized correctly"));
1888 preconditioner->reinit(diagonal_inverse);
1890 VectorType empty_vector;
1891 diagonal_inverse.reinit(empty_vector);
1895 if (preconditioner->m() !=
matrix.m())
1897 preconditioner->get_vector().reinit(
matrix.m());
1898 for (
typename VectorType::size_type i=0; i<
matrix.m(); ++i)
1899 preconditioner->get_vector()(i) = 1./
matrix.el(i,i);
1904 template <
typename VectorType>
1905 void set_initial_guess(VectorType &vector)
1907 vector = 1./std::sqrt(static_cast<double>(vector.size()));
1908 if (vector.locally_owned_elements().is_element(0))
1912 template <
typename Number>
1913 void set_initial_guess(::Vector<Number> &vector)
1919 for (
unsigned int i=0; i<vector.size(); ++i)
1922 const Number mean_value = vector.mean_value();
1923 vector.add(-mean_value);
1926 template <
typename Number>
1938 for (
unsigned int i=0; i<vector.
local_size(); ++i)
1941 const Number mean_value = vector.
mean_value();
1942 vector.
add(-mean_value);
1945 struct EigenvalueTracker
1953 std::vector<double> values;
1962 template <
typename MatrixType,
class VectorType,
typename PreconditionerType>
1966 const double smoothing_range,
1967 const bool nonzero_starting,
1968 const unsigned int eig_cg_n_iterations,
1969 const double eig_cg_residual,
1970 const double max_eigenvalue)
1973 smoothing_range (smoothing_range),
1974 nonzero_starting (nonzero_starting),
1975 eig_cg_n_iterations (eig_cg_n_iterations),
1976 eig_cg_residual (eig_cg_residual),
1977 max_eigenvalue (max_eigenvalue)
1982 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
1988 eigenvalues_are_initialized (false)
1991 std::is_same<size_type, typename VectorType::size_type>::value,
1992 "PreconditionChebyshev and VectorType must have the same size_type.");
1998 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2002 (
const MatrixType &matrix,
2003 const AdditionalData &additional_data)
2006 data = additional_data;
2007 internal::PreconditionChebyshevImplementation::initialize_preconditioner(matrix,
2008 data.preconditioner,
2009 data.matrix_diagonal_inverse);
2010 eigenvalues_are_initialized =
false;
2015 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2020 eigenvalues_are_initialized =
false;
2021 theta = delta = 1.0;
2022 matrix_ptr =
nullptr;
2024 VectorType empty_vector;
2025 data.matrix_diagonal_inverse.reinit(empty_vector);
2026 update1.
reinit(empty_vector);
2027 update2.
reinit(empty_vector);
2028 update3.reinit(empty_vector);
2030 data.preconditioner.reset();
2036 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2040 (
const VectorType &src)
const 2046 update2.
reinit(src,
true);
2051 double max_eigenvalue, min_eigenvalue;
2052 if (data.eig_cg_n_iterations > 0)
2054 Assert (data.eig_cg_n_iterations > 2,
2055 ExcMessage (
"Need to set at least two iterations to find eigenvalues."));
2059 std::sqrt(std::numeric_limits<typename VectorType::value_type>::epsilon()),
2060 1e-10,
false,
false);
2062 internal::PreconditionChebyshevImplementation::EigenvalueTracker eigenvalue_tracker;
2064 solver.connect_eigenvalues_slot(std::bind(&internal::PreconditionChebyshevImplementation::EigenvalueTracker::slot,
2065 &eigenvalue_tracker,
2066 std::placeholders::_1));
2070 internal::PreconditionChebyshevImplementation::set_initial_guess(update2);
2074 solver.solve(*matrix_ptr, update1, update2, *data.preconditioner);
2081 if (eigenvalue_tracker.values.empty())
2082 min_eigenvalue = max_eigenvalue = 1;
2085 min_eigenvalue = eigenvalue_tracker.values.front();
2089 max_eigenvalue = 1.2*eigenvalue_tracker.values.back();
2094 max_eigenvalue = data.max_eigenvalue;
2095 min_eigenvalue = data.max_eigenvalue/data.smoothing_range;
2098 const double alpha = (data.smoothing_range > 1. ?
2099 max_eigenvalue / data.smoothing_range :
2100 std::min(0.9*max_eigenvalue, min_eigenvalue));
2109 const double actual_range = max_eigenvalue / alpha;
2110 const double sigma = (1.-std::sqrt(1./actual_range))/(1.+std::sqrt(1./actual_range));
2111 const double eps = data.smoothing_range;
2113 = 1+std::log(1./eps+std::sqrt(1./eps/eps-1))/std::log(1./sigma);
2126 update3.reinit (src,
true);
2133 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2138 const VectorType &src)
const 2142 if (std::abs(delta) < 1e-40)
2145 double rhok = delta / theta, sigma = theta / delta;
2146 for (
unsigned int k=0; k<data.degree; ++k)
2148 matrix_ptr->vmult (update2, dst);
2149 const double rhokp = 1./(2.*sigma-rhok);
2150 const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/delta;
2152 internal::PreconditionChebyshevImplementation::vector_updates
2153 (src, *data.preconditioner,
false, factor1, factor2, update1, update2, update3, dst);
2159 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2164 const VectorType &src)
const 2166 double rhok = delta / theta, sigma = theta / delta;
2167 for (
unsigned int k=0; k<data.degree; ++k)
2169 matrix_ptr->Tvmult (update2, dst);
2170 const double rhokp = 1./(2.*sigma-rhok);
2171 const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/delta;
2173 internal::PreconditionChebyshevImplementation::vector_updates
2174 (src, *data.preconditioner,
false, factor1, factor2, update1, update2, update3, dst);
2180 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2185 const VectorType &src)
const 2188 if (eigenvalues_are_initialized ==
false)
2189 estimate_eigenvalues(src);
2191 internal::PreconditionChebyshevImplementation::vector_updates
2192 (src, *data.preconditioner,
true, 0., 1./theta, update1, update2, update3, dst);
2194 do_chebyshev_loop(dst, src);
2199 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2204 const VectorType &src)
const 2207 if (eigenvalues_are_initialized ==
false)
2208 estimate_eigenvalues(src);
2210 internal::PreconditionChebyshevImplementation::vector_updates
2211 (src, *data.preconditioner,
true, 0., 1./theta, update1, update2, update3, dst);
2213 do_transpose_chebyshev_loop(dst, src);
2218 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2223 const VectorType &src)
const 2226 if (eigenvalues_are_initialized ==
false)
2227 estimate_eigenvalues(src);
2229 matrix_ptr->vmult (update2, dst);
2230 internal::PreconditionChebyshevImplementation::vector_updates
2231 (src, *data.preconditioner,
false, 0., 1./theta, update1, update2, update3, dst);
2233 do_chebyshev_loop(dst, src);
2238 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2243 const VectorType &src)
const 2246 if (eigenvalues_are_initialized ==
false)
2247 estimate_eigenvalues(src);
2249 matrix_ptr->Tvmult (update2, dst);
2250 internal::PreconditionChebyshevImplementation::vector_updates
2251 (src, *data.preconditioner,
false, 0., 1./theta, update1, update2, update3, dst);
2253 do_transpose_chebyshev_loop(dst, src);
2258 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2264 return matrix_ptr->m();
2269 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2275 return matrix_ptr->n();
2280 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
Contents is actually a matrix.
size_type nth_index_in_set(const unsigned int local_index) const
void Tstep(VectorType &x, const VectorType &rhs) const
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)
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
std::shared_ptr< PreconditionerType > preconditioner
#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
std::array< std::pair< Number, Tensor< 1, dim, Number > >, dim > jacobi(::SymmetricTensor< 2, dim, Number > A)
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
const std::vector< size_type > & inverse_permutation
#define Assert(cond, exc)
virtual void add(const Number a) override
bool eigenvalues_are_initialized
void step(VectorType &x, const VectorType &rhs) const
void(MatrixType::* function_ptr)(VectorType &, const VectorType &) const
virtual ::IndexSet locally_owned_elements() const override
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
virtual Number mean_value() const override
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())