From 2da333ae5a09016e6bb6f7a80006454ee32f7024 Mon Sep 17 00:00:00 2001 From: Christoph Lehmann <christoph.lehmann@ufz.de> Date: Tue, 1 Mar 2016 11:25:29 +0100 Subject: [PATCH] [T] adapted to ConfigTree, linear solver constructors --- .../DenseGaussEliminationChecker.cpp | 4 +-- Tests/MathLib/TestDenseGaussAlgorithm.cpp | 21 +++++------ Tests/MathLib/TestLinearSolver.cpp | 36 +++++++++++++------ Tests/MathLib/TestNonlinearNewton.cpp | 4 +-- Tests/MathLib/TestNonlinearPicard.cpp | 4 +-- 5 files changed, 42 insertions(+), 27 deletions(-) diff --git a/SimpleTests/MatrixTests/DenseGaussEliminationChecker.cpp b/SimpleTests/MatrixTests/DenseGaussEliminationChecker.cpp index a3a618b0bc2..57b068c503f 100644 --- a/SimpleTests/MatrixTests/DenseGaussEliminationChecker.cpp +++ b/SimpleTests/MatrixTests/DenseGaussEliminationChecker.cpp @@ -64,8 +64,8 @@ int main(int argc, char *argv[]) b.resize(n_rows); b = mat * x; - MathLib::GaussAlgorithm<MathLib::DenseMatrix<double, std::size_t>> gauss(mat); - gauss.solve(b, true); + MathLib::GaussAlgorithm<MathLib::DenseMatrix<double, std::size_t>> gauss; + gauss.solve(mat, b, true); { std::stringstream stream; diff --git a/Tests/MathLib/TestDenseGaussAlgorithm.cpp b/Tests/MathLib/TestDenseGaussAlgorithm.cpp index 609b54887ea..448672a4442 100644 --- a/Tests/MathLib/TestDenseGaussAlgorithm.cpp +++ b/Tests/MathLib/TestDenseGaussAlgorithm.cpp @@ -58,28 +58,28 @@ TEST(MathLib, DenseGaussAlgorithm) double *x3 (new double[n_cols]); std::generate(x3,x3+n_cols, std::rand); - MathLib::GaussAlgorithm<MathLib::DenseMatrix<double, std::size_t>, double*> gauss(mat); + MathLib::GaussAlgorithm<MathLib::DenseMatrix<double, std::size_t>, double*> gauss; // solve with b0 as right hand side - gauss.solve(b0, true); + gauss.solve(mat, b0, true); for (std::size_t i(0); i<n_rows; i++) { ASSERT_NEAR(b0[i], 0.0, std::numeric_limits<float>::epsilon()); } // solve with b1 as right hand side - gauss.solve(b1, false); + gauss.solve(mat, b1, false); for (std::size_t i(0); i<n_rows; i++) { ASSERT_NEAR(b1[i], 1.0, std::numeric_limits<float>::epsilon()); } // solve with b2 as right hand side - gauss.solve(b2, false); + gauss.solve(mat, b2, false); for (std::size_t i(0); i<n_rows; i++) { ASSERT_NEAR(fabs(b2[i]-x[i])/fabs(x[i]), 0.0, std::numeric_limits<float>::epsilon()); } // solve with b3 as right hand side and x3 as solution vector - gauss.solve(b3, x3, false); + gauss.solve(mat, b3, x3, false); for (std::size_t i(0); i<n_rows; i++) { ASSERT_NEAR(fabs(x3[i]-x[i])/fabs(x[i]), 0.0, std::numeric_limits<float>::epsilon()); } @@ -133,27 +133,28 @@ TEST(MathLib, DenseGaussAlgorithmDenseVector) MathLib::DenseVector<double> x3 (n_cols); std::generate(std::begin(x3),std::end(x3), std::rand); - MathLib::GaussAlgorithm<MathLib::DenseMatrix<double, std::size_t>, MathLib::DenseVector<double>> gauss(mat); + MathLib::GaussAlgorithm<MathLib::DenseMatrix<double, std::size_t>, + MathLib::DenseVector<double>> gauss; // solve with b0 as right hand side - gauss.solve(b0, true); + gauss.solve(mat, b0, true); for (std::size_t i(0); i<n_rows; i++) { ASSERT_NEAR(b0[i], 0.0, 1e5 * std::numeric_limits<float>::epsilon()); } // solve with b1 as right hand side - gauss.solve(b1, false); + gauss.solve(mat, b1, false); for (std::size_t i(0); i<n_rows; i++) { ASSERT_NEAR(b1[i], 1.0, std::numeric_limits<float>::epsilon()); } // solve with b2 as right hand side - gauss.solve(b2, false); + gauss.solve(mat, b2, false); for (std::size_t i(0); i<n_rows; i++) { ASSERT_NEAR(fabs(b2[i]-x[i])/fabs(x[i]), 0.0, std::numeric_limits<float>::epsilon()); } - gauss.solve(b3, x3, false); + gauss.solve(mat, b3, x3, false); for (std::size_t i(0); i<n_rows; i++) { ASSERT_NEAR(fabs(x3[i]-x[i])/fabs(x[i]), 0.0, std::numeric_limits<float>::epsilon()); } diff --git a/Tests/MathLib/TestLinearSolver.cpp b/Tests/MathLib/TestLinearSolver.cpp index 23c8ea8ed5c..beaa6115abf 100644 --- a/Tests/MathLib/TestLinearSolver.cpp +++ b/Tests/MathLib/TestLinearSolver.cpp @@ -34,6 +34,7 @@ #ifdef USE_LIS #include "MathLib/LinAlg/Lis/LisLinearSolver.h" +#include "MathLib/LinAlg/Lis/LisVector.h" #endif #ifdef USE_PETSC @@ -125,8 +126,8 @@ void checkLinearSolverInterface(T_MATRIX &A, BaseLib::ConfigTree const& ls_optio MathLib::finalizeMatrixAssembly(A); // solve - T_LINEAR_SOVLER ls(A, "dummy_name", &ls_option); - ls.solve(rhs, x); + T_LINEAR_SOVLER ls("dummy_name", &ls_option); + ls.solve(A, rhs, x); ASSERT_ARRAY_NEAR(ex1.exH, x, ex1.dim_eqs, 1e-5); @@ -189,8 +190,8 @@ void checkLinearSolverInterface(T_MATRIX& A, T_VECTOR& b, MathLib::finalizeMatrixAssembly(A); // solve - T_LINEAR_SOVLER ls(A, prefix_name, &ls_option); - EXPECT_TRUE(ls.solve(b, x)); + T_LINEAR_SOVLER ls(prefix_name, &ls_option); + EXPECT_TRUE(ls.solve(A, b, x)); EXPECT_GT(ls.getNumberOfIterations(), 0u); @@ -204,7 +205,8 @@ void checkLinearSolverInterface(T_MATRIX& A, T_VECTOR& b, TEST(MathLib, CheckInterface_GaussAlgorithm) { boost::property_tree::ptree t_root; - BaseLib::ConfigTree conf(t_root, ""); + BaseLib::ConfigTree conf(t_root, "", + BaseLib::ConfigTree::onerror, BaseLib::ConfigTree::onwarning); using Example = Example1<std::size_t>; @@ -226,7 +228,8 @@ TEST(Math, CheckInterface_Eigen) t_solver.put("error_tolerance", 1e-15); t_solver.put("max_iteration_step", 1000); t_root.put_child("eigen", t_solver); - BaseLib::ConfigTree conf(t_root, ""); + BaseLib::ConfigTree conf(t_root, "", + BaseLib::ConfigTree::onerror, BaseLib::ConfigTree::onwarning); using IntType = MathLib::EigenMatrix::IndexType; @@ -243,7 +246,8 @@ TEST(Math, CheckInterface_EigenLis) boost::property_tree::ptree t_root; boost::property_tree::ptree t_solver; t_root.put("lis", "-i cg -p none -tol 1e-15 -maxiter 1000"); - BaseLib::ConfigTree conf(t_root, ""); + BaseLib::ConfigTree conf(t_root, "", + BaseLib::ConfigTree::onerror, BaseLib::ConfigTree::onwarning); using IntType = MathLib::LisMatrix::IndexType; @@ -260,7 +264,8 @@ TEST(Math, CheckInterface_Lis) boost::property_tree::ptree t_root; boost::property_tree::ptree t_solver; t_root.put("lis", "-i cg -p none -tol 1e-15 -maxiter 1000"); - BaseLib::ConfigTree conf(t_root, ""); + BaseLib::ConfigTree conf(t_root, "", + BaseLib::ConfigTree::onerror, BaseLib::ConfigTree::onwarning); using IntType = MathLib::LisMatrix::IndexType; @@ -294,7 +299,10 @@ TEST(MPITest_Math, CheckInterface_PETSc_Linear_Solver_basic) checkLinearSolverInterface<MathLib::PETScMatrix, MathLib::PETScVector, MathLib::PETScLinearSolver>( - A, b, "ptest1_", BaseLib::ConfigTree(t_root, "")); + A, b, "ptest1_", + BaseLib::ConfigTree(t_root, "", BaseLib::ConfigTree::onerror, + BaseLib::ConfigTree::onwarning) + ); } TEST(MPITest_Math, CheckInterface_PETSc_Linear_Solver_chebyshev_sor) @@ -320,7 +328,10 @@ TEST(MPITest_Math, CheckInterface_PETSc_Linear_Solver_chebyshev_sor) checkLinearSolverInterface<MathLib::PETScMatrix, MathLib::PETScVector, MathLib::PETScLinearSolver>( - A, b, "ptest2_", BaseLib::ConfigTree(t_root, "")); + A, b, "ptest2_", + BaseLib::ConfigTree(t_root, "", BaseLib::ConfigTree::onerror, + BaseLib::ConfigTree::onwarning) + ); } TEST(MPITest_Math, CheckInterface_PETSc_Linear_Solver_gmres_amg) @@ -348,7 +359,10 @@ TEST(MPITest_Math, CheckInterface_PETSc_Linear_Solver_gmres_amg) checkLinearSolverInterface<MathLib::PETScMatrix, MathLib::PETScVector, MathLib::PETScLinearSolver>( - A, b, "ptest3_", BaseLib::ConfigTree(t_root, "")); + A, b, "ptest3_", + BaseLib::ConfigTree(t_root, "", BaseLib::ConfigTree::onerror, + BaseLib::ConfigTree::onwarning) + ); } #endif diff --git a/Tests/MathLib/TestNonlinearNewton.cpp b/Tests/MathLib/TestNonlinearNewton.cpp index 1cb00f25e41..0f0d3649763 100644 --- a/Tests/MathLib/TestNonlinearNewton.cpp +++ b/Tests/MathLib/TestNonlinearNewton.cpp @@ -51,10 +51,10 @@ public: void operator()(const VectorType &x, const VectorType &r, VectorType &dx) { _f_Jacobian(x, _matJ); - DenseSolverType solver(_matJ); + DenseSolverType solver; VectorType rhs(r); rhs *= -1.; - solver.solve(rhs, dx); + solver.solve(_matJ, rhs, dx); } private: diff --git a/Tests/MathLib/TestNonlinearPicard.cpp b/Tests/MathLib/TestNonlinearPicard.cpp index d225bd26f5d..73401155e27 100644 --- a/Tests/MathLib/TestNonlinearPicard.cpp +++ b/Tests/MathLib/TestNonlinearPicard.cpp @@ -61,8 +61,8 @@ public: b[0] = -2; b[1] = 0.; - DenseSolverType solver(A); - solver.solve(b, x_new); + DenseSolverType solver; + solver.solve(A, b, x_new); } private: MatrixType A; -- GitLab