Skip to content
Snippets Groups Projects
Commit 344d756b authored by Timo Koch's avatar Timo Koch
Browse files

[test] Add unit test for sequential linear solver

parent 6f24254d
No related branches found
No related tags found
1 merge request!1850Restructure linearsolver parallel
...@@ -2,6 +2,7 @@ add_subdirectory(common) ...@@ -2,6 +2,7 @@ add_subdirectory(common)
add_subdirectory(geomechanics) add_subdirectory(geomechanics)
add_subdirectory(freeflow) add_subdirectory(freeflow)
add_subdirectory(io) add_subdirectory(io)
add_subdirectory(linear)
add_subdirectory(material) add_subdirectory(material)
add_subdirectory(multidomain) add_subdirectory(multidomain)
add_subdirectory(nonlinear) add_subdirectory(nonlinear)
......
dune_symlink_to_source_files(FILES "params.input")
dumux_add_test(NAME test_linearsolver
SOURCES test_linearsolver.cc
LABELS linear unit)
TestSolver = CG
ProblemSize = 2
[CG.LinearSolver]
Type = cgsolver
Preconditioner.Type = amg
#include <config.h>
#include <iostream>
#include <iomanip>
#include <cmath>
#include <dumux/common/parameters.hh>
#include <dumux/discretization/method.hh>
#include <dumux/linear/linearsolvertraits.hh>
//#include <dumux/linear/istlsolverfactorybackend.hh>
#include <dumux/linear/amgbackend.hh>
#include <dune/common/exceptions.hh>
#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
#include <dune/grid/common/mcmgmapper.hh>
#include <dune/grid/yaspgrid.hh>
#include <dune/istl/bvector.hh>
#include <dune/istl/bcrsmatrix.hh>
#include <dune/istl/test/laplacian.hh>
#include <dune/istl/paamg/test/anisotropic.hh>
namespace Dumux::Test {
struct MockGridGeometry
{
using GridView = Dune::YaspGrid<2>::LeafGridView;
using DofMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
using VertexMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
static constexpr auto discMethod = DiscretizationMethod::box;
};
} // end namespace Dumux::Test
int main(int argc, char* argv[]) try
{
using namespace Dumux;
Dune::MPIHelper::instance(argc, argv);
Parameters::init(argc, argv, "params.input");
using Vector = Dune::BlockVector<Dune::FieldVector<double, 2>>;
using Matrix = Dune::BCRSMatrix<Dune::FieldMatrix<double, 2, 2>>;
// create matrix & vectors (we want to solve Ax=b)
Matrix A; const int numDofs = getParam<int>("ProblemSize");
setupLaplacian(A, numDofs);
Vector x(A.N()); Vector b(A.M());
x = 0; b = 1;
using LinearSolverTraits = Dumux::LinearSolverTraits<Test::MockGridGeometry>;
using LinearSolver = AMGBiCGSTABBackend<LinearSolverTraits>;
const auto testSolverName = getParam<std::string>("TestSolver");
LinearSolver solver(testSolverName);
solver.solve(A, x, b);
if (!solver.result().converged)
DUNE_THROW(Dune::Exception, testSolverName << " did not converge!");
return 0;
}
catch (const Dune::Exception& e)
{
std::cout << e << std::endl;
return 1;
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment