From: Martin Kronbichler Date: Mon, 23 May 2022 12:13:17 +0000 (+0200) Subject: Add test case X-Git-Tag: v9.4.0-rc1~167^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=359e3cdfac03b60f7be394e49c9b00aa00273ba0;p=dealii.git Add test case --- diff --git a/tests/lac/precondition_chebyshev_06.cc b/tests/lac/precondition_chebyshev_06.cc new file mode 100644 index 0000000000..53cfb9bb55 --- /dev/null +++ b/tests/lac/precondition_chebyshev_06.cc @@ -0,0 +1,89 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2022 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + +// Test PreconditionChebyshev with power iteration rather than SolverCG to +// estimate the eigenvalues (less accurate) + + +#include +#include +#include +#include + +#include "../tests.h" + +#include "../testmatrix.h" + + + +int +main() +{ + initlog(); + deallog << std::setprecision(10); + + for (unsigned int size = 4; size <= 16; size *= 2) + { + unsigned int dim = (size - 1) * (size - 1); + + deallog << "Size " << size << " Unknowns " << dim << std::endl; + + // Make matrix + FDMatrix testproblem(size, size); + SparsityPattern structure(dim, dim, 5); + testproblem.five_point_structure(structure); + structure.compress(); + SparseMatrix A(structure); + testproblem.five_point(A); + + using Chebyshev = PreconditionChebyshev, + Vector, + SparseILU>; + Chebyshev cheby; + Chebyshev::AdditionalData cheby_data; + cheby_data.preconditioner.reset(new SparseILU()); + cheby_data.preconditioner->initialize(A); + cheby_data.degree = 11; + cheby_data.smoothing_range = 40; + cheby_data.eigenvalue_algorithm = + Chebyshev::AdditionalData::EigenvalueAlgorithm::power_iteration; + cheby.initialize(A, cheby_data); + + Vector v(dim); + Vector tmp1(dim), tmp2(dim); + for (unsigned int i = 0; i < 3; ++i) + { + for (unsigned int j = 0; j < dim; ++j) + v(j) = random_value(); + + A.vmult(tmp1, v); + cheby_data.preconditioner->vmult(tmp2, tmp1); + tmp2 -= v; + const double ilu_residual = tmp2.l2_norm(); + + A.vmult(tmp1, v); + cheby.vmult(tmp2, tmp1); + tmp2 -= v; + const double cheby_residual = tmp2.l2_norm(); + + deallog << "Residual step i=" << i << ": " + << " ilu=" << ilu_residual << ", cheby=" << cheby_residual + << std::endl; + } + } + + return 0; +} diff --git a/tests/lac/precondition_chebyshev_06.output b/tests/lac/precondition_chebyshev_06.output new file mode 100644 index 0000000000..0143cd8cf0 --- /dev/null +++ b/tests/lac/precondition_chebyshev_06.output @@ -0,0 +1,13 @@ + +DEAL::Size 4 Unknowns 9 +DEAL::Residual step i=0: ilu=0.3321067622, cheby=0.05099432420 +DEAL::Residual step i=1: ilu=0.3592209788, cheby=0.04600439733 +DEAL::Residual step i=2: ilu=0.1716913688, cheby=0.02710123523 +DEAL::Size 8 Unknowns 49 +DEAL::Residual step i=0: ilu=1.943632434, cheby=0.1724610774 +DEAL::Residual step i=1: ilu=2.148727169, cheby=0.1731956717 +DEAL::Residual step i=2: ilu=1.756318039, cheby=0.1687993258 +DEAL::Size 16 Unknowns 225 +DEAL::Residual step i=0: ilu=6.146258348, cheby=0.4251547823 +DEAL::Residual step i=1: ilu=6.252364591, cheby=0.4316860132 +DEAL::Residual step i=2: ilu=5.653977586, cheby=0.4063348761