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/lac/tridiagonal_matrix.h> 27 #include <deal.II/lac/solver_cg.h> 28 #include <deal.II/lac/vector_memory.h> 30 DEAL_II_NAMESPACE_OPEN
34 template <
typename number>
class Vector;
40 template <
typename number>
class Vector;
101 template <
typename MatrixType>
102 void initialize (
const MatrixType &matrix,
108 template<
class VectorType>
109 void vmult (VectorType &,
const VectorType &)
const;
115 template<
class VectorType>
116 void Tvmult (VectorType &,
const VectorType &)
const;
121 template<
class VectorType>
122 void vmult_add (VectorType &,
const VectorType &)
const;
128 template<
class VectorType>
129 void Tvmult_add (VectorType &,
const VectorType &)
const;
144 size_type m ()
const;
153 size_type n ()
const;
224 template <
typename MatrixType>
225 void initialize (
const MatrixType &matrix,
231 template<
class VectorType>
232 void vmult (VectorType &,
const VectorType &)
const;
238 template<
class VectorType>
239 void Tvmult (VectorType &,
const VectorType &)
const;
243 template<
class VectorType>
244 void vmult_add (VectorType &,
const VectorType &)
const;
250 template<
class VectorType>
251 void Tvmult_add (VectorType &,
const VectorType &)
const;
266 size_type m ()
const;
275 size_type n ()
const;
336 template<
typename MatrixType = SparseMatrix<
double>,
class VectorType = Vector<
double> >
343 typedef void ( MatrixType::* function_ptr)(VectorType &,
const VectorType &)
const;
351 const function_ptr method);
357 void vmult (VectorType &dst,
358 const VectorType &src)
const;
382 template<
typename MatrixType = SparseMatrix<
double> >
413 void initialize (
const MatrixType &A,
425 size_type m ()
const;
431 size_type n ()
const;
474 template <
typename MatrixType = SparseMatrix<
double> >
481 template<
class VectorType>
482 void vmult (VectorType &,
const VectorType &)
const;
488 template<
class VectorType>
489 void Tvmult (VectorType &,
const VectorType &)
const;
494 template<
class VectorType>
495 void step (VectorType &x,
const VectorType &rhs)
const;
500 template<
class VectorType>
501 void Tstep (VectorType &x,
const VectorType &rhs)
const;
551 template <
typename MatrixType = SparseMatrix<
double> >
558 template<
class VectorType>
559 void vmult (VectorType &,
const VectorType &)
const;
564 template<
class VectorType>
565 void Tvmult (VectorType &,
const VectorType &)
const;
570 template<
class VectorType>
571 void step (VectorType &x,
const VectorType &rhs)
const;
576 template<
class VectorType>
577 void Tstep (VectorType &x,
const VectorType &rhs)
const;
609 template <
typename MatrixType = SparseMatrix<
double> >
629 void initialize (
const MatrixType &A,
635 template<
class VectorType>
636 void vmult (VectorType &,
const VectorType &)
const;
642 template<
class VectorType>
643 void Tvmult (VectorType &,
const VectorType &)
const;
649 template<
class VectorType>
650 void step (VectorType &x,
const VectorType &rhs)
const;
655 template<
class VectorType>
656 void Tstep (VectorType &x,
const VectorType &rhs)
const;
699 template <
typename MatrixType = SparseMatrix<
double> >
725 const std::vector<size_type> &inverse_permutation,
754 void initialize (
const MatrixType &A,
755 const std::vector<size_type> &permutation,
756 const std::vector<size_type> &inverse_permutation,
770 void initialize (
const MatrixType &A,
776 template<
class VectorType>
777 void vmult (VectorType &,
const VectorType &)
const;
782 template<
class VectorType>
783 void Tvmult (VectorType &,
const VectorType &)
const;
811 template <
typename MatrixType=SparseMatrix<
double>,
class VectorType=Vector<
double> >
830 const double smoothing_range = 0.,
831 const bool nonzero_starting =
false,
832 const unsigned int eig_cg_n_iterations = 8,
833 const double eig_cg_residual = 1e-2,
834 const double max_eigenvalue = 1);
907 void initialize (
const MatrixType &matrix,
914 void vmult (VectorType &dst,
915 const VectorType &src)
const;
921 void Tvmult (VectorType &dst,
922 const VectorType &src)
const;
933 size_type m ()
const;
939 size_type n ()
const;
994 template <
typename MatrixType>
1000 n_columns = matrix.n();
1004 template<
class VectorType>
1013 template<
class VectorType>
1020 template<
class VectorType>
1029 template<
class VectorType>
1039 Assert(n_rows != 0, ExcNotInitialized());
1046 Assert(n_columns != 0, ExcNotInitialized());
1055 relaxation(relaxation)
1067 relaxation=add_data.relaxation;
1081 template <
typename MatrixType>
1084 (
const MatrixType &matrix,
1088 n_rows = matrix.m();
1089 n_columns = matrix.n();
1094 template<
class VectorType>
1098 #ifdef DEAL_II_WITH_CXX11 1100 std::is_same<size_type, typename VectorType::size_type>::value,
1101 "PreconditionRichardson and VectorType must have the same size_type.");
1102 #endif // DEAL_II_WITH_CXX11 1104 dst.equ(relaxation,src);
1109 template<
class VectorType>
1113 #ifdef DEAL_II_WITH_CXX11 1115 std::is_same<size_type, typename VectorType::size_type>::value,
1116 "PreconditionRichardson and VectorType must have the same size_type.");
1117 #endif // DEAL_II_WITH_CXX11 1119 dst.equ(relaxation,src);
1122 template<
class VectorType>
1126 #ifdef DEAL_II_WITH_CXX11 1128 std::is_same<size_type, typename VectorType::size_type>::value,
1129 "PreconditionRichardson and VectorType must have the same size_type.");
1130 #endif // DEAL_II_WITH_CXX11 1132 dst.add(relaxation,src);
1137 template<
class VectorType>
1141 #ifdef DEAL_II_WITH_CXX11 1143 std::is_same<size_type, typename VectorType::size_type>::value,
1144 "PreconditionRichardson and VectorType must have the same size_type.");
1145 #endif // DEAL_II_WITH_CXX11 1147 dst.add(relaxation,src);
1153 Assert(n_rows != 0, ExcNotInitialized());
1160 Assert(n_columns != 0, ExcNotInitialized());
1166 template <
typename MatrixType>
1172 relaxation = parameters.relaxation;
1176 template <
typename MatrixType>
1183 template <
typename MatrixType>
1187 Assert (A!=0, ExcNotInitialized());
1191 template <
typename MatrixType>
1195 Assert (A!=0, ExcNotInitialized());
1201 template <
typename MatrixType>
1202 template<
class VectorType>
1206 #ifdef DEAL_II_WITH_CXX11 1209 "PreconditionJacobi and VectorType must have the same size_type.");
1210 #endif // DEAL_II_WITH_CXX11 1212 Assert (this->A!=0, ExcNotInitialized());
1213 this->A->precondition_Jacobi (dst, src, this->relaxation);
1218 template <
typename MatrixType>
1219 template<
class VectorType>
1223 #ifdef DEAL_II_WITH_CXX11 1226 "PreconditionJacobi and VectorType must have the same size_type.");
1227 #endif // DEAL_II_WITH_CXX11 1229 Assert (this->A!=0, ExcNotInitialized());
1230 this->A->precondition_Jacobi (dst, src, this->relaxation);
1235 template <
typename MatrixType>
1236 template<
class VectorType>
1240 #ifdef DEAL_II_WITH_CXX11 1243 "PreconditionJacobi and VectorType must have the same size_type.");
1244 #endif // DEAL_II_WITH_CXX11 1246 Assert (this->A!=0, ExcNotInitialized());
1247 this->A->Jacobi_step (dst, src, this->relaxation);
1252 template <
typename MatrixType>
1253 template<
class VectorType>
1257 #ifdef DEAL_II_WITH_CXX11 1260 "PreconditionJacobi and VectorType must have the same size_type.");
1261 #endif // DEAL_II_WITH_CXX11 1270 template <
typename MatrixType>
1271 template<
class VectorType>
1275 #ifdef DEAL_II_WITH_CXX11 1278 "PreconditionSOR and VectorType must have the same size_type.");
1279 #endif // DEAL_II_WITH_CXX11 1281 Assert (this->A!=0, ExcNotInitialized());
1282 this->A->precondition_SOR (dst, src, this->relaxation);
1287 template <
typename MatrixType>
1288 template<
class VectorType>
1292 #ifdef DEAL_II_WITH_CXX11 1295 "PreconditionSOR and VectorType must have the same size_type.");
1296 #endif // DEAL_II_WITH_CXX11 1298 Assert (this->A!=0, ExcNotInitialized());
1299 this->A->precondition_TSOR (dst, src, this->relaxation);
1304 template <
typename MatrixType>
1305 template<
class VectorType>
1309 #ifdef DEAL_II_WITH_CXX11 1312 "PreconditionSOR and VectorType must have the same size_type.");
1313 #endif // DEAL_II_WITH_CXX11 1315 Assert (this->A!=0, ExcNotInitialized());
1316 this->A->SOR_step (dst, src, this->relaxation);
1321 template <
typename MatrixType>
1322 template<
class VectorType>
1326 #ifdef DEAL_II_WITH_CXX11 1329 "PreconditionSOR and VectorType must have the same size_type.");
1330 #endif // DEAL_II_WITH_CXX11 1332 Assert (this->A!=0, ExcNotInitialized());
1333 this->A->TSOR_step (dst, src, this->relaxation);
1340 template <
typename MatrixType>
1343 const typename BaseClass::AdditionalData ¶meters)
1355 const size_type
n = this->A->n();
1356 pos_right_of_diagonal.resize(n, static_cast<std::size_t>(-1));
1357 for (size_type row=0; row<
n; ++row)
1364 it = mat->
begin(row)+1;
1365 for ( ; it < mat->
end(row); ++it)
1366 if (it->column() > row)
1368 pos_right_of_diagonal[row] = it - mat->
begin();
1374 template <
typename MatrixType>
1375 template<
class VectorType>
1379 #ifdef DEAL_II_WITH_CXX11 1382 "PreconditionSSOR and VectorType must have the same size_type.");
1383 #endif // DEAL_II_WITH_CXX11 1385 Assert (this->A!=0, ExcNotInitialized());
1386 this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
1391 template <
typename MatrixType>
1392 template<
class VectorType>
1396 #ifdef DEAL_II_WITH_CXX11 1399 "PreconditionSSOR and VectorType must have the same size_type.");
1400 #endif // DEAL_II_WITH_CXX11 1402 Assert (this->A!=0, ExcNotInitialized());
1403 this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
1408 template <
typename MatrixType>
1409 template<
class VectorType>
1413 #ifdef DEAL_II_WITH_CXX11 1416 "PreconditionSSOR and VectorType must have the same size_type.");
1417 #endif // DEAL_II_WITH_CXX11 1419 Assert (this->A!=0, ExcNotInitialized());
1420 this->A->SSOR_step (dst, src, this->relaxation);
1425 template <
typename MatrixType>
1426 template<
class VectorType>
1430 #ifdef DEAL_II_WITH_CXX11 1433 "PreconditionSSOR and VectorType must have the same size_type.");
1434 #endif // DEAL_II_WITH_CXX11 1443 template <
typename MatrixType>
1446 (
const MatrixType &rA,
1447 const std::vector<size_type> &p,
1448 const std::vector<size_type> &ip,
1452 inverse_permutation = &ip;
1457 template <
typename MatrixType>
1463 additional_data.permutation,
1464 additional_data.inverse_permutation,
1465 additional_data.parameters);
1469 template <
typename MatrixType>
1470 template <
typename VectorType>
1474 #ifdef DEAL_II_WITH_CXX11 1477 "PreconditionPSOR and VectorType must have the same size_type.");
1478 #endif // DEAL_II_WITH_CXX11 1480 Assert (this->A!=0, ExcNotInitialized());
1482 this->A->PSOR (dst, *permutation, *inverse_permutation, this->relaxation);
1487 template <
typename MatrixType>
1488 template<
class VectorType>
1492 #ifdef DEAL_II_WITH_CXX11 1495 "PreconditionPSOR and VectorType must have the same size_type.");
1496 #endif // DEAL_II_WITH_CXX11 1498 Assert (this->A!=0, ExcNotInitialized());
1500 this->A->TPSOR (dst, *permutation, *inverse_permutation, this->relaxation);
1503 template <
typename MatrixType>
1505 (
const std::vector<size_type> &permutation,
1506 const std::vector<size_type> &inverse_permutation,
1509 permutation(permutation),
1510 inverse_permutation(inverse_permutation),
1511 parameters(parameters)
1520 template<
typename MatrixType,
class VectorType>
1522 const function_ptr method)
1524 matrix(M), precondition(method)
1529 template<
typename MatrixType,
class VectorType>
1532 const VectorType &src)
const 1534 (matrix.*precondition)(dst, src);
1539 template<
typename MatrixType>
1544 relaxation (relaxation)
1561 template <
typename VectorType>
1564 vector_updates (
const VectorType &src,
1565 const VectorType &matrix_diagonal_inverse,
1566 const bool start_zero,
1567 const double factor1,
1568 const double factor2,
1575 dst.equ (factor2, src);
1576 dst.scale (matrix_diagonal_inverse);
1577 update1.equ(-1.,dst);
1582 update2.scale (matrix_diagonal_inverse);
1584 update1.equ(factor2, update2);
1586 update1.sadd(factor1, factor2, update2);
1594 template <
typename Number>
1595 struct VectorUpdater
1597 VectorUpdater (
const Number *src,
1598 const Number *matrix_diagonal_inverse,
1599 const bool start_zero,
1600 const Number factor1,
1601 const Number factor2,
1607 matrix_diagonal_inverse (matrix_diagonal_inverse),
1608 do_startup (factor1 == Number()),
1609 start_zero (start_zero),
1618 apply_to_subrange (
const std::size_t begin,
1619 const std::size_t end)
const 1625 const Number factor1 = this->factor1;
1626 const Number factor2 = this->factor2;
1630 DEAL_II_OPENMP_SIMD_PRAGMA
1631 for (std::size_t i=begin; i<end; ++i)
1633 dst[i] = factor2 * src[i] * matrix_diagonal_inverse[i];
1634 update1[i] = -dst[i];
1637 DEAL_II_OPENMP_SIMD_PRAGMA
1638 for (std::size_t i=begin; i<end; ++i)
1640 update1[i] = ((update2[i]-src[i]) *
1641 factor2*matrix_diagonal_inverse[i]);
1642 dst[i] -= update1[i];
1646 DEAL_II_OPENMP_SIMD_PRAGMA
1647 for (std::size_t i=begin; i<end; ++i)
1649 const Number update =
1650 factor1 * update1[i] + factor2 *
1651 ((update2[i] - src[i]) * matrix_diagonal_inverse[i]);
1652 update1[i] = update;
1658 const Number *matrix_diagonal_inverse;
1659 const bool do_startup;
1660 const bool start_zero;
1661 const Number factor1;
1662 const Number factor2;
1665 mutable Number *dst;
1668 template<
typename Number>
1671 VectorUpdatesRange(
const VectorUpdater<Number> &updater,
1672 const std::size_t size)
1676 if (size < internal::Vector::minimum_parallel_grain_size)
1677 apply_to_subrange (0, size);
1679 apply_parallel (0, size,
1680 internal::Vector::minimum_parallel_grain_size);
1683 ~VectorUpdatesRange() {}
1686 apply_to_subrange (
const std::size_t begin,
1687 const std::size_t end)
const 1689 updater.apply_to_subrange(begin, end);
1692 const VectorUpdater<Number> &updater;
1696 template <
typename Number>
1699 vector_updates (const ::Vector<Number> &src,
1700 const ::Vector<Number> &matrix_diagonal_inverse,
1701 const bool start_zero,
1702 const double factor1,
1703 const double factor2,
1708 VectorUpdater<Number> upd(src.begin(), matrix_diagonal_inverse.begin(),
1709 start_zero, factor1, factor2,
1711 VectorUpdatesRange<Number>(upd, src.size());
1715 template <
typename Number>
1720 const bool start_zero,
1721 const double factor1,
1722 const double factor2,
1727 VectorUpdater<Number> upd(src.
begin(), matrix_diagonal_inverse.
begin(),
1728 start_zero, factor1, factor2,
1730 VectorUpdatesRange<Number>(upd, src.
local_size());
1733 template <
typename VectorType>
1734 struct DiagonalPreconditioner
1736 DiagonalPreconditioner (
const VectorType &vector)
1738 diagonal_vector(vector)
1741 void vmult (VectorType &dst,
1742 const VectorType &src)
const 1745 dst.scale(diagonal_vector);
1748 const VectorType &diagonal_vector;
1751 struct EigenvalueTracker
1754 void slot(
const std::vector<double> &eigenvalues)
1756 values = eigenvalues;
1759 std::vector<double> values;
1766 template <
typename MatrixType,
class VectorType>
1770 const double smoothing_range,
1771 const bool nonzero_starting,
1772 const unsigned int eig_cg_n_iterations,
1773 const double eig_cg_residual,
1774 const double max_eigenvalue)
1777 smoothing_range (smoothing_range),
1778 nonzero_starting (nonzero_starting),
1779 eig_cg_n_iterations (eig_cg_n_iterations),
1780 eig_cg_residual (eig_cg_residual),
1781 max_eigenvalue (max_eigenvalue)
1786 template <
typename MatrixType,
class VectorType>
1792 #ifdef DEAL_II_WITH_CXX11 1794 std::is_same<size_type, typename VectorType::size_type>::value,
1795 "PreconditionChebyshev and VectorType must have the same size_type.");
1796 #endif // DEAL_II_WITH_CXX11 1801 template <
typename MatrixType,
class VectorType>
1805 (
const MatrixType &matrix,
1809 data = additional_data;
1813 ExcMessage(
"Matrix diagonal vector set but not sized correctly"));
1815 for (
unsigned int i=0; i<matrix.m(); ++i)
1823 double max_eigenvalue, min_eigenvalue;
1827 ExcMessage (
"Need to set at least two iterations to find eigenvalues."));
1832 VectorType *rhs = memory.
alloc();
1833 VectorType *dummy = memory.
alloc();
1842 *rhs = 1./std::sqrt(static_cast<double>(matrix.m()));
1843 if (rhs->locally_owned_elements().is_element(0))
1845 rhs->compress(VectorOperation::insert);
1847 internal::PreconditionChebyshev::EigenvalueTracker eigenvalue_tracker;
1850 &eigenvalue_tracker,
1852 internal::PreconditionChebyshev::DiagonalPreconditioner<VectorType>
1856 solver.
solve(matrix, *dummy, *rhs, preconditioner);
1866 if (eigenvalue_tracker.values.empty())
1867 min_eigenvalue = max_eigenvalue = 1;
1870 min_eigenvalue = eigenvalue_tracker.values.front();
1871 max_eigenvalue = eigenvalue_tracker.values.back();
1876 max_eigenvalue *= 1.2;
1886 std::min(0.9*max_eigenvalue, min_eigenvalue));
1887 delta = (max_eigenvalue-alpha)*0.5;
1888 theta = (max_eigenvalue+alpha)*0.5;
1898 template <
typename MatrixType,
class VectorType>
1902 const VectorType &src)
const 1909 internal::PreconditionChebyshev::vector_updates
1914 internal::PreconditionChebyshev::vector_updates
1921 const double rhokp = 1./(2.*sigma-rhok);
1922 const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/
delta;
1924 internal::PreconditionChebyshev::vector_updates
1932 template <
typename MatrixType,
class VectorType>
1936 const VectorType &src)
const 1943 internal::PreconditionChebyshev::vector_updates
1948 internal::PreconditionChebyshev::vector_updates
1955 const double rhokp = 1./(2.*sigma-rhok);
1956 const double factor1 = rhokp * rhok, factor2 = 2.*rhokp/
delta;
1958 internal::PreconditionChebyshev::vector_updates
1966 template <
typename MatrixType,
typename VectorType>
1978 template <
typename MatrixType,
typename VectorType>
1988 template <
typename MatrixType,
typename VectorType>
1999 DEAL_II_NAMESPACE_CLOSE
void Tvmult(VectorType &, const VectorType &) const
SmartPointer< const MatrixType, PreconditionChebyshev< MatrixType, VectorType > > matrix_ptr
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
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
boost::signals2::connection connect_eigenvalues_slot(const std_cxx11::function< void(const std::vector< double > &)> &slot, const bool every_iteration=false)
void vmult(VectorType &dst, const VectorType &src) const
MatrixType::size_type size_type
void Tvmult_add(VectorType &, const VectorType &) 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)
::ExceptionBase & ExcMessage(std::string arg1)
PreconditionRelaxation< MatrixType >::AdditionalData parameters
void vmult(VectorType &, const VectorType &) const
VectorType matrix_diagonal_inverse
const MatrixType & matrix
AdditionalData(const double relaxation=1.)
void initialize(const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const double relaxation=1.)
types::global_dof_index size_type
types::global_dof_index size_type
const std::vector< size_type > * permutation
const function_ptr precondition
void Tvmult(VectorType &dst, const VectorType &src) const
void vmult(VectorType &, const VectorType &) const
void initialize(const AdditionalData ¶meters)
void Tvmult(VectorType &, const VectorType &) const
types::global_dof_index size_type
void vmult(VectorType &, const VectorType &) const
unsigned int global_dof_index
virtual VectorType * alloc()
const std::vector< size_type > & inverse_permutation
#define Assert(cond, exc)
void step(VectorType &x, const VectorType &rhs) const
unsigned int eig_cg_n_iterations
MatrixType::size_type size_type
virtual void free(const VectorType *const)
void vmult_add(VectorType &, const VectorType &) const
void Tvmult(VectorType &, const VectorType &) const
void initialize(const MatrixType &A, const AdditionalData ¶meters=AdditionalData())
void Tvmult(VectorType &, const VectorType &) const
void initialize(const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
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
AdditionalData(const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const typename PreconditionRelaxation< MatrixType >::AdditionalData ¶meters=typename PreconditionRelaxation< MatrixType >::AdditionalData())
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
const std::vector< size_type > * inverse_permutation
void step(VectorType &x, const VectorType &rhs) const
void vmult(VectorType &dst, const VectorType &src) const
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
size_type local_size() const
std::vector< std::size_t > pos_right_of_diagonal
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())