16#ifndef dealii_precondition_h
17#define dealii_precondition_h
40template <
typename number>
42template <
typename number>
48 template <
typename,
typename>
108 template <
typename MatrixType>
116 template <
class VectorType>
118 vmult(VectorType &,
const VectorType &)
const;
124 template <
class VectorType>
126 Tvmult(VectorType &,
const VectorType &)
const;
131 template <
class VectorType>
139 template <
class VectorType>
238 template <
typename MatrixType>
245 template <
class VectorType>
247 vmult(VectorType &,
const VectorType &)
const;
253 template <
class VectorType>
255 Tvmult(VectorType &,
const VectorType &)
const;
259 template <
class VectorType>
267 template <
class VectorType>
357template <
typename MatrixType = SparseMatrix<
double>,
358 class VectorType = Vector<
double>>
366 const VectorType &)
const;
380 vmult(VectorType &dst,
const VectorType &src)
const;
401template <
typename MatrixType = SparseMatrix<
double>>
496template <
typename MatrixType = SparseMatrix<
double>>
509 template <
class VectorType>
511 vmult(VectorType &,
const VectorType &)
const;
517 template <
class VectorType>
519 Tvmult(VectorType &,
const VectorType &)
const;
524 template <
class VectorType>
526 step(VectorType &x,
const VectorType &rhs)
const;
531 template <
class VectorType>
533 Tstep(VectorType &x,
const VectorType &rhs)
const;
582template <
typename MatrixType = SparseMatrix<
double>>
595 template <
class VectorType>
597 vmult(VectorType &,
const VectorType &)
const;
602 template <
class VectorType>
604 Tvmult(VectorType &,
const VectorType &)
const;
609 template <
class VectorType>
611 step(VectorType &x,
const VectorType &rhs)
const;
616 template <
class VectorType>
618 Tstep(VectorType &x,
const VectorType &rhs)
const;
649template <
typename MatrixType = SparseMatrix<
double>>
677 const typename BaseClass::AdditionalData ¶meters =
678 typename BaseClass::AdditionalData());
683 template <
class VectorType>
685 vmult(VectorType &,
const VectorType &)
const;
691 template <
class VectorType>
693 Tvmult(VectorType &,
const VectorType &)
const;
699 template <
class VectorType>
701 step(VectorType &x,
const VectorType &rhs)
const;
706 template <
class VectorType>
708 Tstep(VectorType &x,
const VectorType &rhs)
const;
748template <
typename MatrixType = SparseMatrix<
double>>
829 template <
class VectorType>
831 vmult(VectorType &,
const VectorType &)
const;
836 template <
class VectorType>
838 Tvmult(VectorType &,
const VectorType &)
const;
988template <
typename MatrixType = SparseMatrix<
double>,
989 typename VectorType = Vector<
double>,
990 typename PreconditionerType = DiagonalMatrix<VectorType>>
1105 vmult(VectorType &dst,
const VectorType &src)
const;
1112 Tvmult(VectorType &dst,
const VectorType &src)
const;
1118 step(VectorType &dst,
const VectorType &src)
const;
1124 Tstep(VectorType &dst,
const VectorType &src)
const;
1192 EigenvalueInformation
1261template <
typename MatrixType>
1271template <
class VectorType>
1280template <
class VectorType>
1287template <
class VectorType>
1296template <
class VectorType>
1320 const double relaxation)
1321 : relaxation(relaxation)
1330 AdditionalData add_data;
1331 relaxation = add_data.relaxation;
1345template <
typename MatrixType>
1348 const MatrixType &
matrix,
1358template <
class VectorType>
1363 std::is_same<size_type, typename VectorType::size_type>::value,
1364 "PreconditionRichardson and VectorType must have the same size_type.");
1371template <
class VectorType>
1376 std::is_same<size_type, typename VectorType::size_type>::value,
1377 "PreconditionRichardson and VectorType must have the same size_type.");
1382template <
class VectorType>
1387 std::is_same<size_type, typename VectorType::size_type>::value,
1388 "PreconditionRichardson and VectorType must have the same size_type.");
1395template <
class VectorType>
1400 std::is_same<size_type, typename VectorType::size_type>::value,
1401 "PreconditionRichardson and VectorType must have the same size_type.");
1422template <
typename MatrixType>
1425 const AdditionalData ¶meters)
1432template <
typename MatrixType>
1439template <
typename MatrixType>
1447template <
typename MatrixType>
1457template <
typename MatrixType>
1458template <
class VectorType>
1461 const VectorType &src)
const
1466 "PreconditionJacobi and VectorType must have the same size_type.");
1469 this->
A->precondition_Jacobi(dst, src, this->relaxation);
1474template <
typename MatrixType>
1475template <
class VectorType>
1478 const VectorType &src)
const
1483 "PreconditionJacobi and VectorType must have the same size_type.");
1486 this->
A->precondition_Jacobi(dst, src, this->relaxation);
1491template <
typename MatrixType>
1492template <
class VectorType>
1495 const VectorType &src)
const
1500 "PreconditionJacobi and VectorType must have the same size_type.");
1503 this->
A->Jacobi_step(dst, src, this->relaxation);
1508template <
typename MatrixType>
1509template <
class VectorType>
1512 const VectorType &src)
const
1517 "PreconditionJacobi and VectorType must have the same size_type.");
1526template <
typename MatrixType>
1527template <
class VectorType>
1533 "PreconditionSOR and VectorType must have the same size_type.");
1536 this->
A->precondition_SOR(dst, src, this->relaxation);
1541template <
typename MatrixType>
1542template <
class VectorType>
1545 const VectorType &src)
const
1549 "PreconditionSOR and VectorType must have the same size_type.");
1552 this->
A->precondition_TSOR(dst, src, this->relaxation);
1557template <
typename MatrixType>
1558template <
class VectorType>
1564 "PreconditionSOR and VectorType must have the same size_type.");
1567 this->
A->SOR_step(dst, src, this->relaxation);
1572template <
typename MatrixType>
1573template <
class VectorType>
1579 "PreconditionSOR and VectorType must have the same size_type.");
1582 this->
A->TSOR_step(dst, src, this->relaxation);
1589template <
typename MatrixType>
1592 const MatrixType & rA,
1593 const typename BaseClass::AdditionalData ¶meters)
1607 pos_right_of_diagonal.resize(n,
static_cast<std::size_t
>(-1));
1615 it = mat->
begin(row) + 1;
1616 for (; it < mat->
end(row); ++it)
1617 if (it->column() > row)
1619 pos_right_of_diagonal[row] = it - mat->
begin();
1625template <
typename MatrixType>
1626template <
class VectorType>
1629 const VectorType &src)
const
1634 "PreconditionSSOR and VectorType must have the same size_type.");
1637 this->
A->precondition_SSOR(dst, src, this->relaxation, pos_right_of_diagonal);
1642template <
typename MatrixType>
1643template <
class VectorType>
1646 const VectorType &src)
const
1651 "PreconditionSSOR and VectorType must have the same size_type.");
1654 this->
A->precondition_SSOR(dst, src, this->relaxation, pos_right_of_diagonal);
1659template <
typename MatrixType>
1660template <
class VectorType>
1667 "PreconditionSSOR and VectorType must have the same size_type.");
1670 this->
A->SSOR_step(dst, src, this->relaxation);
1675template <
typename MatrixType>
1676template <
class VectorType>
1679 const VectorType &src)
const
1684 "PreconditionSSOR and VectorType must have the same size_type.");
1693template <
typename MatrixType>
1696 const MatrixType & rA,
1697 const std::vector<size_type> & p,
1698 const std::vector<size_type> & ip,
1702 inverse_permutation = &ip;
1707template <
typename MatrixType>
1710 const AdditionalData &additional_data)
1713 additional_data.permutation,
1714 additional_data.inverse_permutation,
1715 additional_data.parameters);
1719template <
typename MatrixType>
1720template <
typename VectorType>
1723 const VectorType &src)
const
1728 "PreconditionPSOR and VectorType must have the same size_type.");
1732 this->
A->PSOR(dst, *permutation, *inverse_permutation, this->relaxation);
1737template <
typename MatrixType>
1738template <
class VectorType>
1741 const VectorType &src)
const
1746 "PreconditionPSOR and VectorType must have the same size_type.");
1750 this->
A->TPSOR(dst, *permutation, *inverse_permutation, this->relaxation);
1753template <
typename MatrixType>
1755 const std::vector<size_type> &permutation,
1756 const std::vector<size_type> &inverse_permutation,
1758 : permutation(permutation)
1759 , inverse_permutation(inverse_permutation)
1760 , parameters(parameters)
1767template <
typename MatrixType,
class VectorType>
1769 const MatrixType & M,
1770 const function_ptr method)
1772 , precondition(method)
1777template <
typename MatrixType,
class VectorType>
1781 const VectorType &src)
const
1783 (
matrix.*precondition)(dst, src);
1788template <
typename MatrixType>
1790 const double relaxation)
1791 : relaxation(relaxation)
1800 namespace PreconditionChebyshevImplementation
1808 template <
typename VectorType,
typename PreconditionerType>
1810 vector_updates(
const VectorType & rhs,
1811 const PreconditionerType &preconditioner,
1812 const unsigned int iteration_index,
1813 const double factor1,
1814 const double factor2,
1815 VectorType & solution_old,
1816 VectorType & temp_vector1,
1817 VectorType & temp_vector2,
1818 VectorType & solution)
1820 if (iteration_index == 0)
1822 solution.equ(factor2, rhs);
1823 preconditioner.vmult(solution_old, solution);
1825 else if (iteration_index == 1)
1828 temp_vector1.sadd(-1.0, 1.0, rhs);
1829 preconditioner.vmult(solution_old, temp_vector1);
1832 solution_old.sadd(factor2, 1 + factor1, solution);
1837 temp_vector1.sadd(-1.0, 1.0, rhs);
1838 preconditioner.vmult(temp_vector2, temp_vector1);
1841 solution_old.sadd(-factor1, factor2, temp_vector2);
1842 solution_old.add(1 + factor1, solution);
1845 solution.swap(solution_old);
1851 template <
typename Number>
1852 struct VectorUpdater
1854 VectorUpdater(
const Number * rhs,
1855 const Number * matrix_diagonal_inverse,
1856 const unsigned int iteration_index,
1857 const Number factor1,
1858 const Number factor2,
1859 Number * solution_old,
1860 Number * tmp_vector,
1863 , matrix_diagonal_inverse(matrix_diagonal_inverse)
1864 , iteration_index(iteration_index)
1867 , solution_old(solution_old)
1868 , tmp_vector(tmp_vector)
1869 , solution(solution)
1873 apply_to_subrange(
const std::size_t
begin,
const std::size_t
end)
const
1879 const Number factor1 = this->factor1;
1880 const Number factor1_plus_1 = 1. + this->factor1;
1881 const Number factor2 = this->factor2;
1882 if (iteration_index == 0)
1885 for (std::size_t i =
begin; i <
end; ++i)
1886 solution[i] = factor2 * matrix_diagonal_inverse[i] * rhs[i];
1888 else if (iteration_index == 1)
1892 for (std::size_t i =
begin; i <
end; ++i)
1896 factor1_plus_1 * solution[i] +
1897 factor2 * matrix_diagonal_inverse[i] * (rhs[i] - tmp_vector[i]);
1904 for (std::size_t i =
begin; i <
end; ++i)
1906 factor1_plus_1 * solution[i] - factor1 * solution_old[i] +
1907 factor2 * matrix_diagonal_inverse[i] * (rhs[i] - tmp_vector[i]);
1912 const Number * matrix_diagonal_inverse;
1913 const unsigned int iteration_index;
1914 const Number factor1;
1915 const Number factor2;
1916 mutable Number * solution_old;
1917 mutable Number * tmp_vector;
1918 mutable Number * solution;
1921 template <
typename Number>
1924 VectorUpdatesRange(
const VectorUpdater<Number> &updater,
1925 const std::size_t size)
1929 VectorUpdatesRange::apply_to_subrange(0, size);
1937 ~VectorUpdatesRange()
override =
default;
1940 apply_to_subrange(
const std::size_t
begin,
1941 const std::size_t
end)
const override
1943 updater.apply_to_subrange(
begin,
end);
1946 const VectorUpdater<Number> &updater;
1950 template <
typename Number>
1952 vector_updates(const ::Vector<Number> & rhs,
1954 const unsigned int iteration_index,
1955 const double factor1,
1956 const double factor2,
1962 VectorUpdater<Number> upd(rhs.begin(),
1963 jacobi.get_vector().begin(),
1967 solution_old.
begin(),
1968 temp_vector1.
begin(),
1970 VectorUpdatesRange<Number>(upd, rhs.size());
1973 if (iteration_index == 0)
1978 else if (iteration_index == 1)
1980 solution.
swap(temp_vector1);
1981 solution_old.
swap(temp_vector1);
1984 solution.
swap(solution_old);
1988 template <
typename Number>
1994 const unsigned int iteration_index,
1995 const double factor1,
1996 const double factor2,
2004 VectorUpdater<Number> upd(rhs.
begin(),
2005 jacobi.get_vector().begin(),
2009 solution_old.
begin(),
2010 temp_vector1.
begin(),
2015 if (iteration_index == 0)
2020 else if (iteration_index == 1)
2022 solution.
swap(temp_vector1);
2023 solution_old.
swap(temp_vector1);
2026 solution.
swap(solution_old);
2029 template <
typename MatrixType,
typename PreconditionerType>
2031 initialize_preconditioner(
2032 const MatrixType &
matrix,
2033 std::shared_ptr<PreconditionerType> &preconditioner)
2036 (void)preconditioner;
2040 template <
typename MatrixType,
typename VectorType>
2042 initialize_preconditioner(
2043 const MatrixType &
matrix,
2046 if (preconditioner.get() ==
nullptr || preconditioner->m() !=
matrix.m())
2048 if (preconditioner.get() ==
nullptr)
2049 preconditioner = std::make_shared<DiagonalMatrix<VectorType>>();
2052 preconditioner->m() == 0,
2054 "Preconditioner appears to be initialized but not sized correctly"));
2057 if (preconditioner->m() !=
matrix.m())
2059 preconditioner->get_vector().reinit(
matrix.m());
2061 preconditioner->get_vector()(i) = 1. /
matrix.el(i, i);
2066 template <
typename VectorType>
2068 set_initial_guess(VectorType &vector)
2070 vector = 1. /
std::sqrt(
static_cast<double>(vector.size()));
2071 if (vector.locally_owned_elements().is_element(0))
2075 template <
typename Number>
2083 for (
unsigned int i = 0; i < vector.
size(); ++i)
2086 const Number mean_value = vector.
mean_value();
2087 vector.
add(-mean_value);
2090 template <
typename Number>
2108 const Number mean_value = vector.
mean_value();
2109 vector.
add(-mean_value);
2112 template <
typename Number>
2117 for (
unsigned int block = 0; block < vector.
n_blocks(); ++block)
2118 set_initial_guess(vector.
block(block));
2122# ifdef DEAL_II_COMPILER_CUDA_AWARE
2123 template <
typename Number>
2126 const unsigned int locally_owned_size,
2130 const unsigned int index = threadIdx.x + blockDim.x * blockIdx.x;
2131 if (index < locally_owned_size)
2132 values[index] = (index + offset) % 11;
2135 template <
typename Number>
2154 set_initial_guess_kernel<<<n_blocks, CUDAWrappers::block_size>>>(
2155 first_local_range, n_local_elements, vector.
get_values());
2158 const Number mean_value = vector.
mean_value();
2159 vector.
add(-mean_value);
2163 struct EigenvalueTracker
2172 std::vector<double>
values;
2179template <
typename MatrixType,
class VectorType,
typename PreconditionerType>
2182 const double smoothing_range,
2183 const unsigned int eig_cg_n_iterations,
2184 const double eig_cg_residual,
2185 const double max_eigenvalue)
2187 , smoothing_range(smoothing_range)
2188 , eig_cg_n_iterations(eig_cg_n_iterations)
2189 , eig_cg_residual(eig_cg_residual)
2190 , max_eigenvalue(max_eigenvalue)
2195template <
typename MatrixType,
class VectorType,
typename PreconditionerType>
2198 PreconditionerType>::AdditionalData &
2202 degree = other_data.
degree;
2203 smoothing_range = other_data.smoothing_range;
2204 eig_cg_n_iterations = other_data.eig_cg_n_iterations;
2205 eig_cg_residual = other_data.eig_cg_residual;
2206 max_eigenvalue = other_data.max_eigenvalue;
2207 preconditioner = other_data.preconditioner;
2208 constraints.copy_from(other_data.constraints);
2215template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2220 , eigenvalues_are_initialized(false)
2223 std::is_same<size_type, typename VectorType::size_type>::value,
2224 "PreconditionChebyshev and VectorType must have the same size_type.");
2229template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2232 const MatrixType &
matrix,
2233 const AdditionalData &additional_data)
2236 data = additional_data;
2238 ExcMessage(
"The degree of the Chebyshev method must be positive."));
2239 internal::PreconditionChebyshevImplementation::initialize_preconditioner(
2240 matrix, data.preconditioner);
2241 eigenvalues_are_initialized =
false;
2246template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2250 eigenvalues_are_initialized =
false;
2251 theta = delta = 1.0;
2252 matrix_ptr =
nullptr;
2254 VectorType empty_vector;
2255 solution_old.
reinit(empty_vector);
2256 temp_vector1.
reinit(empty_vector);
2257 temp_vector2.reinit(empty_vector);
2259 data.preconditioner.reset();
2264template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2267 PreconditionerType>::EigenvalueInformation
2275 EigenvalueInformation info{};
2277 solution_old.
reinit(src);
2278 temp_vector1.
reinit(src,
true);
2280 if (data.eig_cg_n_iterations > 0)
2282 Assert(data.eig_cg_n_iterations > 2,
2284 "Need to set at least two iterations to find eigenvalues."));
2288 data.eig_cg_n_iterations,
2295 internal::PreconditionChebyshevImplementation::EigenvalueTracker
2298 solver.connect_eigenvalues_slot(
2299 [&eigenvalue_tracker](
const std::vector<double> &
eigenvalues) {
2306 internal::PreconditionChebyshevImplementation::set_initial_guess(
2308 data.constraints.set_zero(temp_vector1);
2312 solver.solve(*matrix_ptr,
2315 *data.preconditioner);
2321 if (eigenvalue_tracker.values.empty())
2322 info.min_eigenvalue_estimate = info.max_eigenvalue_estimate = 1.;
2325 info.min_eigenvalue_estimate = eigenvalue_tracker.values.front();
2329 info.max_eigenvalue_estimate = 1.2 * eigenvalue_tracker.values.back();
2332 info.cg_iterations = control.last_step();
2336 info.max_eigenvalue_estimate = data.max_eigenvalue;
2337 info.min_eigenvalue_estimate = data.max_eigenvalue / data.smoothing_range;
2340 const double alpha = (data.smoothing_range > 1. ?
2341 info.max_eigenvalue_estimate / data.smoothing_range :
2342 std::min(0.9 * info.max_eigenvalue_estimate,
2343 info.min_eigenvalue_estimate));
2352 const double actual_range = info.max_eigenvalue_estimate / alpha;
2353 const double sigma = (1. -
std::sqrt(1. / actual_range)) /
2355 const double eps = data.smoothing_range;
2360 1 +
static_cast<unsigned int>(
2365 info.degree = data.degree;
2369 ->delta = (info.max_eigenvalue_estimate - alpha) * 0.5;
2372 ->theta = (info.max_eigenvalue_estimate + alpha) * 0.5;
2376 using NumberType =
typename VectorType::value_type;
2380 ((std::is_same<VectorType,
2384 (std::is_same<VectorType,
2385 LinearAlgebra::distributed::
2386 Vector<NumberType, MemorySpace::CUDA>>::value ==
2388 temp_vector2.reinit(src,
true);
2391 VectorType empty_vector;
2392 temp_vector2.reinit(empty_vector);
2397 ->eigenvalues_are_initialized =
true;
2404template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2407 VectorType & solution,
2408 const VectorType &rhs)
const
2410 std::lock_guard<std::mutex> lock(mutex);
2411 if (eigenvalues_are_initialized ==
false)
2412 estimate_eigenvalues(rhs);
2414 internal::PreconditionChebyshevImplementation::vector_updates(
2416 *data.preconditioner,
2427 if (data.degree < 2 ||
std::abs(delta) < 1
e-40)
2430 double rhok = delta / theta, sigma = theta / delta;
2431 for (
unsigned int k = 0; k < data.degree - 1; ++k)
2433 matrix_ptr->vmult(temp_vector1, solution);
2434 const double rhokp = 1. / (2. * sigma - rhok);
2435 const double factor1 = rhokp * rhok, factor2 = 2. * rhokp / delta;
2437 internal::PreconditionChebyshevImplementation::vector_updates(
2439 *data.preconditioner,
2452template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2455 VectorType & solution,
2456 const VectorType &rhs)
const
2458 std::lock_guard<std::mutex> lock(mutex);
2459 if (eigenvalues_are_initialized ==
false)
2460 estimate_eigenvalues(rhs);
2462 internal::PreconditionChebyshevImplementation::vector_updates(
2464 *data.preconditioner,
2473 if (data.degree < 2 ||
std::abs(delta) < 1
e-40)
2476 double rhok = delta / theta, sigma = theta / delta;
2477 for (
unsigned int k = 0; k < data.degree - 1; ++k)
2479 matrix_ptr->Tvmult(temp_vector1, solution);
2480 const double rhokp = 1. / (2. * sigma - rhok);
2481 const double factor1 = rhokp * rhok, factor2 = 2. * rhokp / delta;
2483 internal::PreconditionChebyshevImplementation::vector_updates(
2485 *data.preconditioner,
2498template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2501 VectorType & solution,
2502 const VectorType &rhs)
const
2504 std::lock_guard<std::mutex> lock(mutex);
2505 if (eigenvalues_are_initialized ==
false)
2506 estimate_eigenvalues(rhs);
2508 matrix_ptr->vmult(temp_vector1, solution);
2509 internal::PreconditionChebyshevImplementation::vector_updates(
2511 *data.preconditioner,
2520 if (data.degree < 2 ||
std::abs(delta) < 1
e-40)
2523 double rhok = delta / theta, sigma = theta / delta;
2524 for (
unsigned int k = 0; k < data.degree - 1; ++k)
2526 matrix_ptr->vmult(temp_vector1, solution);
2527 const double rhokp = 1. / (2. * sigma - rhok);
2528 const double factor1 = rhokp * rhok, factor2 = 2. * rhokp / delta;
2530 internal::PreconditionChebyshevImplementation::vector_updates(
2532 *data.preconditioner,
2545template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2548 VectorType & solution,
2549 const VectorType &rhs)
const
2551 std::lock_guard<std::mutex> lock(mutex);
2552 if (eigenvalues_are_initialized ==
false)
2553 estimate_eigenvalues(rhs);
2555 matrix_ptr->Tvmult(temp_vector1, solution);
2556 internal::PreconditionChebyshevImplementation::vector_updates(
2558 *data.preconditioner,
2567 if (data.degree < 2 ||
std::abs(delta) < 1
e-40)
2570 double rhok = delta / theta, sigma = theta / delta;
2571 for (
unsigned int k = 0; k < data.degree - 1; ++k)
2573 matrix_ptr->Tvmult(temp_vector1, solution);
2574 const double rhokp = 1. / (2. * sigma - rhok);
2575 const double factor1 = rhokp * rhok, factor2 = 2. * rhokp / delta;
2577 internal::PreconditionChebyshevImplementation::vector_updates(
2579 *data.preconditioner,
2592template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2599 return matrix_ptr->m();
2604template <
typename MatrixType,
typename VectorType,
typename PreconditionerType>
2611 return matrix_ptr->n();
size_type nth_index_in_set(const size_type local_index) const
const_iterator end() const
const_iterator begin() const
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)
#define DEAL_II_OPENMP_SIMD_PRAGMA
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define AssertCudaKernel()
unsigned int n_blocks() const
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
BlockType & block(const unsigned int i)
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
types::global_dof_index size_type
void vmult_add(VectorType &, const VectorType &) const
std::vector< std::size_t > pos_right_of_diagonal
const std::vector< size_type > * inverse_permutation
void Tvmult(VectorType &dst, const VectorType &src) const
void vmult(VectorType &, const VectorType &) const
typename MatrixType::size_type size_type
void Tstep(VectorType &x, const VectorType &rhs) const
SmartPointer< const MatrixType, PreconditionRelaxation< MatrixType > > A
const std::vector< size_type > * permutation
void vmult(VectorType &, const VectorType &) const
std::shared_ptr< PreconditionerType > preconditioner
void(MatrixType::*)(VectorType &, const VectorType &) const function_ptr
void vmult_add(VectorType &, const VectorType &) const
const function_ptr precondition
void step(VectorType &dst, const VectorType &src) const
double min_eigenvalue_estimate
double max_eigenvalue_estimate
typename PreconditionRelaxation< MatrixType >::AdditionalData AdditionalData
AdditionalData(const unsigned int degree=1, const double smoothing_range=0., const unsigned int eig_cg_n_iterations=8, const double eig_cg_residual=1e-2, const double max_eigenvalue=1)
void Tvmult(VectorType &, const VectorType &) const
void step(VectorType &x, const VectorType &rhs) const
void Tvmult(VectorType &, const VectorType &) const
void initialize(const AdditionalData ¶meters)
void vmult(VectorType &, const VectorType &) const
AdditionalData & operator=(const AdditionalData &other_data)
typename MatrixType::size_type size_type
void initialize(const MatrixType &A, const AdditionalData &additional_data)
EigenvalueInformation estimate_eigenvalues(const VectorType &src) const
AffineConstraints< double > constraints
void Tstep(VectorType &dst, const VectorType &src) const
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())
void step(VectorType &x, const VectorType &rhs) const
const std::vector< size_type > & inverse_permutation
typename PreconditionRelaxation< MatrixType >::AdditionalData AdditionalData
void initialize(const MatrixType &A, const AdditionalData ¶meters=AdditionalData())
const MatrixType & matrix
void initialize(const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixType &A, const typename BaseClass::AdditionalData ¶meters=typename BaseClass::AdditionalData())
void initialize(const MatrixType &matrix, const AdditionalData ¶meters)
SmartPointer< const MatrixType, PreconditionChebyshev< MatrixType, VectorType, PreconditionerType > > matrix_ptr
unsigned int eig_cg_n_iterations
unsigned int cg_iterations
AdditionalData(const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const typename PreconditionRelaxation< MatrixType >::AdditionalData ¶meters=typename PreconditionRelaxation< MatrixType >::AdditionalData())
typename PreconditionRelaxation< MatrixType >::AdditionalData AdditionalData
void Tstep(VectorType &x, const VectorType &rhs) const
void Tstep(VectorType &x, const VectorType &rhs) const
AdditionalData(const double relaxation=1.)
PreconditionRelaxation< MatrixType >::AdditionalData parameters
const std::vector< size_type > & permutation
void step(VectorType &x, const VectorType &rhs) const
void vmult(VectorType &, const VectorType &) const
void Tvmult(VectorType &, const VectorType &) const
void Tvmult_add(VectorType &, const VectorType &) const
AdditionalData(const double relaxation=1.)
void vmult(VectorType &, const VectorType &) const
void vmult(VectorType &dst, const VectorType &src) const
typename MatrixType::size_type size_type
void Tvmult(VectorType &, const VectorType &) const
void vmult(VectorType &dst, const VectorType &src) const
void initialize(const MatrixType &matrix, const AdditionalData &additional_data=AdditionalData())
PreconditionUseMatrix(const MatrixType &M, const function_ptr method)
bool eigenvalues_are_initialized
void Tvmult(VectorType &, const VectorType &) const
void Tvmult(VectorType &, const VectorType &) const
void vmult(VectorType &, const VectorType &) const
types::global_dof_index size_type
void Tvmult_add(VectorType &, const VectorType &) const
void add(const std::vector< size_type > &indices, const std::vector< OtherNumber > &values)
Number local_element(const size_type local_index) const
Number * get_values() const
void swap(Vector< Number, MemorySpace > &v)
Number mean_value() const
size_type locally_owned_size() const
virtual Number mean_value() const override
virtual void swap(Vector< Number > &v)
virtual void add(const Number a) override
void reinit(const size_type size, const bool omit_zeroing_entries=false)
virtual ::IndexSet locally_owned_elements() const override
@ matrix
Contents is actually a matrix.
@ eigenvalues
Eigenvalue vector is filled.
types::global_dof_index size_type
std::enable_if< IsBlockVector< VectorType >::value, unsignedint >::type n_blocks(const VectorType &vector)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
std::array< std::pair< Number, Tensor< 1, dim, Number > >, dim > jacobi(::SymmetricTensor< 2, dim, Number > A)
unsigned int minimum_parallel_grain_size
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > log(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index