16 #ifndef dealii_precondition_h
17 #define dealii_precondition_h
40 template <
typename number>
42 template <
typename number>
48 template <
typename,
typename>
109 template <
typename MatrixType>
117 template <
class VectorType>
119 vmult(VectorType &,
const VectorType &)
const;
125 template <
class VectorType>
127 Tvmult(VectorType &,
const VectorType &)
const;
132 template <
class VectorType>
134 vmult_add(VectorType &,
const VectorType &)
const;
140 template <
class VectorType>
142 Tvmult_add(VectorType &,
const VectorType &)
const;
242 template <
typename MatrixType>
249 template <
class VectorType>
257 template <
class VectorType>
263 template <
class VectorType>
271 template <
class VectorType>
363 template <
typename MatrixType = SparseMatrix<
double>,
364 class VectorType = Vector<
double>>
372 const VectorType &)
const;
386 vmult(VectorType &dst,
const VectorType &src)
const;
410 template <
typename MatrixType = SparseMatrix<
double>>
507 template <
typename MatrixType = SparseMatrix<
double>>
520 template <
class VectorType>
522 vmult(VectorType &,
const VectorType &)
const;
528 template <
class VectorType>
530 Tvmult(VectorType &,
const VectorType &)
const;
535 template <
class VectorType>
537 step(VectorType &x,
const VectorType &rhs)
const;
542 template <
class VectorType>
544 Tstep(VectorType &x,
const VectorType &rhs)
const;
595 template <
typename MatrixType = SparseMatrix<
double>>
608 template <
class VectorType>
615 template <
class VectorType>
622 template <
class VectorType>
629 template <
class VectorType>
664 template <
typename MatrixType = SparseMatrix<
double>>
698 template <
class VectorType>
700 vmult(VectorType &,
const VectorType &)
const;
706 template <
class VectorType>
708 Tvmult(VectorType &,
const VectorType &)
const;
714 template <
class VectorType>
716 step(VectorType &x,
const VectorType &rhs)
const;
721 template <
class VectorType>
723 Tstep(VectorType &x,
const VectorType &rhs)
const;
766 template <
typename MatrixType = SparseMatrix<
double>>
847 template <
class VectorType>
854 template <
class VectorType>
1009 template <
typename MatrixType = SparseMatrix<
double>,
1010 typename VectorType = Vector<
double>,
1011 typename PreconditionerType = DiagonalMatrix<VectorType>>
1213 EigenvalueInformation
1282 template <
typename MatrixType>
1292 template <
class VectorType>
1301 template <
class VectorType>
1308 template <
class VectorType>
1317 template <
class VectorType>
1341 const double relaxation)
1342 : relaxation(relaxation)
1351 AdditionalData add_data;
1352 relaxation = add_data.relaxation;
1366 template <
typename MatrixType>
1369 const MatrixType &
matrix,
1379 template <
class VectorType>
1385 "PreconditionRichardson and VectorType must have the same size_type.");
1392 template <
class VectorType>
1398 "PreconditionRichardson and VectorType must have the same size_type.");
1403 template <
class VectorType>
1409 "PreconditionRichardson and VectorType must have the same size_type.");
1416 template <
class VectorType>
1422 "PreconditionRichardson and VectorType must have the same size_type.");
1443 template <
typename MatrixType>
1446 const AdditionalData ¶meters)
1449 relaxation = parameters.relaxation;
1453 template <
typename MatrixType>
1460 template <
typename MatrixType>
1468 template <
typename MatrixType>
1478 template <
typename MatrixType>
1479 template <
class VectorType>
1487 "PreconditionJacobi and VectorType must have the same size_type.");
1490 this->
A->precondition_Jacobi(dst, src, this->relaxation);
1495 template <
typename MatrixType>
1496 template <
class VectorType>
1504 "PreconditionJacobi and VectorType must have the same size_type.");
1507 this->
A->precondition_Jacobi(dst, src, this->relaxation);
1512 template <
typename MatrixType>
1513 template <
class VectorType>
1521 "PreconditionJacobi and VectorType must have the same size_type.");
1524 this->
A->Jacobi_step(dst, src, this->relaxation);
1529 template <
typename MatrixType>
1530 template <
class VectorType>
1538 "PreconditionJacobi and VectorType must have the same size_type.");
1547 template <
typename MatrixType>
1548 template <
class VectorType>
1554 "PreconditionSOR and VectorType must have the same size_type.");
1557 this->
A->precondition_SOR(dst, src, this->relaxation);
1562 template <
typename MatrixType>
1563 template <
class VectorType>
1570 "PreconditionSOR and VectorType must have the same size_type.");
1573 this->
A->precondition_TSOR(dst, src, this->relaxation);
1578 template <
typename MatrixType>
1579 template <
class VectorType>
1585 "PreconditionSOR and VectorType must have the same size_type.");
1588 this->
A->SOR_step(dst, src, this->relaxation);
1593 template <
typename MatrixType>
1594 template <
class VectorType>
1600 "PreconditionSOR and VectorType must have the same size_type.");
1603 this->
A->TSOR_step(dst, src, this->relaxation);
1610 template <
typename MatrixType>
1613 const MatrixType & rA,
1614 const typename BaseClass::AdditionalData ¶meters)
1628 pos_right_of_diagonal.resize(n,
static_cast<std::size_t
>(-1));
1636 it = mat->
begin(row) + 1;
1637 for (; it < mat->
end(row); ++it)
1638 if (it->column() > row)
1640 pos_right_of_diagonal[row] = it - mat->
begin();
1646 template <
typename MatrixType>
1647 template <
class VectorType>
1655 "PreconditionSSOR and VectorType must have the same size_type.");
1658 this->
A->precondition_SSOR(dst, src, this->relaxation, pos_right_of_diagonal);
1663 template <
typename MatrixType>
1664 template <
class VectorType>
1672 "PreconditionSSOR and VectorType must have the same size_type.");
1675 this->
A->precondition_SSOR(dst, src, this->relaxation, pos_right_of_diagonal);
1680 template <
typename MatrixType>
1681 template <
class VectorType>
1688 "PreconditionSSOR and VectorType must have the same size_type.");
1691 this->
A->SSOR_step(dst, src, this->relaxation);
1696 template <
typename MatrixType>
1697 template <
class VectorType>
1705 "PreconditionSSOR and VectorType must have the same size_type.");
1714 template <
typename MatrixType>
1717 const MatrixType & rA,
1718 const std::vector<size_type> & p,
1719 const std::vector<size_type> & ip,
1723 inverse_permutation = &ip;
1728 template <
typename MatrixType>
1731 const AdditionalData &additional_data)
1734 additional_data.permutation,
1735 additional_data.inverse_permutation,
1736 additional_data.parameters);
1740 template <
typename MatrixType>
1741 template <
typename VectorType>
1749 "PreconditionPSOR and VectorType must have the same size_type.");
1753 this->
A->PSOR(dst, *permutation, *inverse_permutation, this->relaxation);
1758 template <
typename MatrixType>
1759 template <
class VectorType>
1767 "PreconditionPSOR and VectorType must have the same size_type.");
1771 this->
A->TPSOR(dst, *permutation, *inverse_permutation, this->relaxation);
1774 template <
typename MatrixType>
1776 const std::vector<size_type> &permutation,
1777 const std::vector<size_type> &inverse_permutation,
1779 : permutation(permutation)
1780 , inverse_permutation(inverse_permutation)
1781 , parameters(parameters)
1788 template <
typename MatrixType,
class VectorType>
1790 const MatrixType & M,
1791 const function_ptr method)
1793 , precondition(method)
1798 template <
typename MatrixType,
class VectorType>
1804 (
matrix.*precondition)(dst, src);
1809 template <
typename MatrixType>
1811 const double relaxation)
1812 : relaxation(relaxation)
1821 namespace PreconditionChebyshevImplementation
1829 template <
typename VectorType,
typename PreconditionerType>
1832 const PreconditionerType &preconditioner,
1833 const unsigned int iteration_index,
1834 const double factor1,
1835 const double factor2,
1841 if (iteration_index == 0)
1843 solution.equ(factor2, rhs);
1844 preconditioner.vmult(solution_old, solution);
1846 else if (iteration_index == 1)
1849 temp_vector1.sadd(-1.0, 1.0, rhs);
1850 preconditioner.vmult(solution_old, temp_vector1);
1853 solution_old.sadd(factor2, 1 + factor1, solution);
1858 temp_vector1.sadd(-1.0, 1.0, rhs);
1859 preconditioner.vmult(temp_vector2, temp_vector1);
1862 solution_old.sadd(-factor1, factor2, temp_vector2);
1863 solution_old.add(1 + factor1, solution);
1866 solution.swap(solution_old);
1872 template <
typename Number>
1873 struct VectorUpdater
1875 VectorUpdater(
const Number * rhs,
1876 const Number * matrix_diagonal_inverse,
1877 const unsigned int iteration_index,
1878 const Number factor1,
1879 const Number factor2,
1880 Number * solution_old,
1881 Number * tmp_vector,
1884 , matrix_diagonal_inverse(matrix_diagonal_inverse)
1885 , iteration_index(iteration_index)
1888 , solution_old(solution_old)
1889 , tmp_vector(tmp_vector)
1890 , solution(solution)
1894 apply_to_subrange(
const std::size_t
begin,
const std::size_t
end)
const
1900 const Number factor1 = this->factor1;
1901 const Number factor1_plus_1 = 1. + this->factor1;
1902 const Number factor2 = this->factor2;
1903 if (iteration_index == 0)
1906 for (std::size_t i =
begin; i <
end; ++i)
1907 solution[i] = factor2 * matrix_diagonal_inverse[i] * rhs[i];
1909 else if (iteration_index == 1)
1913 for (std::size_t i =
begin; i <
end; ++i)
1917 factor1_plus_1 * solution[i] +
1918 factor2 * matrix_diagonal_inverse[i] * (rhs[i] - tmp_vector[i]);
1925 for (std::size_t i =
begin; i <
end; ++i)
1927 factor1_plus_1 * solution[i] - factor1 * solution_old[i] +
1928 factor2 * matrix_diagonal_inverse[i] * (rhs[i] - tmp_vector[i]);
1933 const Number * matrix_diagonal_inverse;
1934 const unsigned int iteration_index;
1935 const Number factor1;
1936 const Number factor2;
1937 mutable Number * solution_old;
1938 mutable Number * tmp_vector;
1939 mutable Number * solution;
1942 template <
typename Number>
1945 VectorUpdatesRange(
const VectorUpdater<Number> &updater,
1946 const std::size_t size)
1950 VectorUpdatesRange::apply_to_subrange(0, size);
1958 ~VectorUpdatesRange()
override =
default;
1961 apply_to_subrange(
const std::size_t
begin,
1962 const std::size_t
end)
const override
1964 updater.apply_to_subrange(
begin,
end);
1967 const VectorUpdater<Number> &updater;
1971 template <
typename Number>
1973 vector_updates(const ::Vector<Number> & rhs,
1975 const unsigned int iteration_index,
1976 const double factor1,
1977 const double factor2,
1978 ::Vector<Number> &solution_old,
1979 ::Vector<Number> &temp_vector1,
1981 ::Vector<Number> &solution)
1983 VectorUpdater<Number> upd(rhs.begin(),
1984 jacobi.get_vector().begin(),
1988 solution_old.begin(),
1989 temp_vector1.begin(),
1991 VectorUpdatesRange<Number>(upd, rhs.size());
1994 if (iteration_index == 0)
1999 else if (iteration_index == 1)
2001 solution.swap(temp_vector1);
2002 solution_old.swap(temp_vector1);
2005 solution.swap(solution_old);
2009 template <
typename Number>
2015 const unsigned int iteration_index,
2016 const double factor1,
2017 const double factor2,
2025 VectorUpdater<Number> upd(rhs.
begin(),
2026 jacobi.get_vector().begin(),
2030 solution_old.
begin(),
2031 temp_vector1.
begin(),
2033 VectorUpdatesRange<Number>(upd, rhs.
local_size());
2036 if (iteration_index == 0)
2041 else if (iteration_index == 1)
2043 solution.
swap(temp_vector1);
2044 solution_old.
swap(temp_vector1);
2047 solution.
swap(solution_old);
2050 template <
typename MatrixType,
typename PreconditionerType>
2052 initialize_preconditioner(
2053 const MatrixType &
matrix,
2054 std::shared_ptr<PreconditionerType> &preconditioner)
2057 (void)preconditioner;
2061 template <
typename MatrixType,
typename VectorType>
2063 initialize_preconditioner(
2064 const MatrixType &
matrix,
2067 if (preconditioner.get() ==
nullptr || preconditioner->m() !=
matrix.m())
2069 if (preconditioner.get() ==
nullptr)
2070 preconditioner = std::make_shared<DiagonalMatrix<VectorType>>();
2073 preconditioner->m() == 0,
2075 "Preconditioner appears to be initialized but not sized correctly"));
2078 if (preconditioner->m() !=
matrix.m())
2080 preconditioner->get_vector().reinit(
matrix.m());
2082 preconditioner->get_vector()(i) = 1. /
matrix.el(i, i);
2087 template <
typename VectorType>
2091 vector = 1. /
std::sqrt(
static_cast<double>(vector.size()));
2092 if (vector.locally_owned_elements().is_element(0))
2096 template <
typename Number>
2098 set_initial_guess(::Vector<Number> &vector)
2104 for (
unsigned int i = 0; i < vector.size(); ++i)
2107 const Number mean_value = vector.mean_value();
2108 vector.add(-mean_value);
2111 template <
typename Number>
2126 for (
unsigned int i = 0; i < vector.
local_size(); ++i)
2129 const Number mean_value = vector.
mean_value();
2130 vector.
add(-mean_value);
2134 # ifdef DEAL_II_COMPILER_CUDA_AWARE
2135 template <
typename Number>
2138 const unsigned int local_size,
2142 const unsigned int index = threadIdx.x + blockDim.x * blockIdx.x;
2143 if (index < local_size)
2144 values[index] = (index + offset) % 11;
2147 template <
typename Number>
2163 const auto n_local_elements = vector.
local_size();
2166 set_initial_guess_kernel<<<n_blocks, CUDAWrappers::block_size>>>(
2167 first_local_range, n_local_elements, vector.
get_values());
2170 const Number mean_value = vector.
mean_value();
2171 vector.
add(-mean_value);
2173 # endif // DEAL_II_COMPILER_CUDA_AWARE
2175 struct EigenvalueTracker
2184 std::vector<double> values;
2191 template <
typename MatrixType,
class VectorType,
typename PreconditionerType>
2194 const double smoothing_range,
2195 const unsigned int eig_cg_n_iterations,
2196 const double eig_cg_residual,
2197 const double max_eigenvalue)
2199 , smoothing_range(smoothing_range)
2200 , eig_cg_n_iterations(eig_cg_n_iterations)
2201 , eig_cg_residual(eig_cg_residual)
2202 , max_eigenvalue(max_eigenvalue)
2207 template <
typename MatrixType,
class VectorType,
typename PreconditionerType>
2210 PreconditionerType>::AdditionalData &
2214 degree = other_data.
degree;
2215 smoothing_range = other_data.smoothing_range;
2216 eig_cg_n_iterations = other_data.eig_cg_n_iterations;
2217 eig_cg_residual = other_data.eig_cg_residual;
2218 max_eigenvalue = other_data.max_eigenvalue;
2219 preconditioner = other_data.preconditioner;
2220 constraints.copy_from(other_data.constraints);
2227 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2232 , eigenvalues_are_initialized(false)
2236 "PreconditionChebyshev and VectorType must have the same size_type.");
2241 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2244 const MatrixType &
matrix,
2245 const AdditionalData &additional_data)
2248 data = additional_data;
2250 ExcMessage(
"The degree of the Chebyshev method must be positive."));
2251 internal::PreconditionChebyshevImplementation::initialize_preconditioner(
2252 matrix, data.preconditioner);
2253 eigenvalues_are_initialized =
false;
2258 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2262 eigenvalues_are_initialized =
false;
2263 theta = delta = 1.0;
2264 matrix_ptr =
nullptr;
2267 solution_old.
reinit(empty_vector);
2268 temp_vector1.
reinit(empty_vector);
2269 temp_vector2.reinit(empty_vector);
2271 data.preconditioner.reset();
2276 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2279 PreconditionerType>::EigenvalueInformation
2287 EigenvalueInformation info{};
2289 solution_old.
reinit(src);
2290 temp_vector1.
reinit(src,
true);
2295 if (data.eig_cg_n_iterations > 0)
2297 Assert(data.eig_cg_n_iterations > 2,
2299 "Need to set at least two iterations to find eigenvalues."));
2303 data.eig_cg_n_iterations,
2310 internal::PreconditionChebyshevImplementation::EigenvalueTracker
2313 solver.connect_eigenvalues_slot(
2314 [&eigenvalue_tracker](
const std::vector<double> &
eigenvalues) {
2320 internal::PreconditionChebyshevImplementation::set_initial_guess(
2322 data.constraints.set_zero(temp_vector1);
2326 solver.solve(*matrix_ptr,
2329 *data.preconditioner);
2335 if (eigenvalue_tracker.values.empty())
2336 info.min_eigenvalue_estimate = info.max_eigenvalue_estimate = 1.;
2339 info.min_eigenvalue_estimate = eigenvalue_tracker.values.front();
2343 info.max_eigenvalue_estimate = 1.2 * eigenvalue_tracker.values.back();
2346 info.cg_iterations = control.last_step();
2350 info.max_eigenvalue_estimate = data.max_eigenvalue;
2351 info.min_eigenvalue_estimate = data.max_eigenvalue / data.smoothing_range;
2354 const double alpha = (data.smoothing_range > 1. ?
2355 info.max_eigenvalue_estimate / data.smoothing_range :
2356 std::min(0.9 * info.max_eigenvalue_estimate,
2357 info.min_eigenvalue_estimate));
2366 const double actual_range = info.max_eigenvalue_estimate / alpha;
2367 const double sigma = (1. -
std::sqrt(1. / actual_range)) /
2368 (1. + std::sqrt(1. / actual_range));
2369 const double eps = data.smoothing_range;
2374 1 +
static_cast<unsigned int>(
2375 std::log(1. /
eps + std::sqrt(1. /
eps /
eps - 1.)) /
2376 std::log(1. / sigma));
2378 info.degree = data.degree;
2383 ->delta = (info.max_eigenvalue_estimate - alpha) * 0.5;
2386 ->theta = (info.max_eigenvalue_estimate + alpha) * 0.5;
2390 using NumberType =
typename VectorType::value_type;
2396 Vector<NumberType, MemorySpace::Host>>::
value ==
2399 LinearAlgebra::distributed::
2400 Vector<NumberType, MemorySpace::CUDA>>
::value ==
2402 temp_vector2.reinit(src,
true);
2406 temp_vector2.reinit(empty_vector);
2411 ->eigenvalues_are_initialized =
true;
2418 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2424 std::lock_guard<std::mutex> lock(mutex);
2425 if (eigenvalues_are_initialized ==
false)
2426 estimate_eigenvalues(rhs);
2428 internal::PreconditionChebyshevImplementation::vector_updates(
2430 *data.preconditioner,
2441 if (data.degree < 2 || std::abs(delta) < 1
e-40)
2444 double rhok = delta / theta, sigma = theta / delta;
2445 for (
unsigned int k = 0; k < data.degree - 1; ++k)
2447 matrix_ptr->vmult(temp_vector1, solution);
2448 const double rhokp = 1. / (2. * sigma - rhok);
2449 const double factor1 = rhokp * rhok, factor2 = 2. * rhokp / delta;
2451 internal::PreconditionChebyshevImplementation::vector_updates(
2453 *data.preconditioner,
2466 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2472 std::lock_guard<std::mutex> lock(mutex);
2473 if (eigenvalues_are_initialized ==
false)
2474 estimate_eigenvalues(rhs);
2476 internal::PreconditionChebyshevImplementation::vector_updates(
2478 *data.preconditioner,
2487 if (data.degree < 2 || std::abs(delta) < 1
e-40)
2490 double rhok = delta / theta, sigma = theta / delta;
2491 for (
unsigned int k = 0; k < data.degree - 1; ++k)
2493 matrix_ptr->Tvmult(temp_vector1, solution);
2494 const double rhokp = 1. / (2. * sigma - rhok);
2495 const double factor1 = rhokp * rhok, factor2 = 2. * rhokp / delta;
2497 internal::PreconditionChebyshevImplementation::vector_updates(
2499 *data.preconditioner,
2512 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2518 std::lock_guard<std::mutex> lock(mutex);
2519 if (eigenvalues_are_initialized ==
false)
2520 estimate_eigenvalues(rhs);
2522 matrix_ptr->vmult(temp_vector1, solution);
2523 internal::PreconditionChebyshevImplementation::vector_updates(
2525 *data.preconditioner,
2534 if (data.degree < 2 || std::abs(delta) < 1
e-40)
2537 double rhok = delta / theta, sigma = theta / delta;
2538 for (
unsigned int k = 0; k < data.degree - 1; ++k)
2540 matrix_ptr->vmult(temp_vector1, solution);
2541 const double rhokp = 1. / (2. * sigma - rhok);
2542 const double factor1 = rhokp * rhok, factor2 = 2. * rhokp / delta;
2544 internal::PreconditionChebyshevImplementation::vector_updates(
2546 *data.preconditioner,
2559 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2565 std::lock_guard<std::mutex> lock(mutex);
2566 if (eigenvalues_are_initialized ==
false)
2567 estimate_eigenvalues(rhs);
2569 matrix_ptr->Tvmult(temp_vector1, solution);
2570 internal::PreconditionChebyshevImplementation::vector_updates(
2572 *data.preconditioner,
2581 if (data.degree < 2 || std::abs(delta) < 1
e-40)
2584 double rhok = delta / theta, sigma = theta / delta;
2585 for (
unsigned int k = 0; k < data.degree - 1; ++k)
2587 matrix_ptr->Tvmult(temp_vector1, solution);
2588 const double rhokp = 1. / (2. * sigma - rhok);
2589 const double factor1 = rhokp * rhok, factor2 = 2. * rhokp / delta;
2591 internal::PreconditionChebyshevImplementation::vector_updates(
2593 *data.preconditioner,
2606 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2613 return matrix_ptr->m();
2618 template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2625 return matrix_ptr->n();