Add test with Jacobi matrix for openclSolver

This commit is contained in:
Tong Dong Qiu 2022-10-13 15:58:53 +02:00
parent 09e262bbfd
commit d06b02f88f

View File

@ -93,8 +93,8 @@ getDuneSolution(Matrix<bz>& matrix, Vector<bz>& rhs)
}
template <int bz>
Dune::BlockVector<Dune::FieldVector<double, bz>>
testOpenclSolver(const boost::property_tree::ptree& prm, Matrix<bz>& matrix, Vector<bz>& rhs)
void
createBridge(const boost::property_tree::ptree& prm, std::unique_ptr<Opm::BdaBridge<Matrix<bz>, Vector<bz>, bz> >& bridge)
{
const int linear_solver_verbosity = prm.get<int>("verbosity");
const int maxit = prm.get<int>("maxiter");
@ -105,11 +105,7 @@ testOpenclSolver(const boost::property_tree::ptree& prm, Matrix<bz>& matrix, Vec
const std::string accelerator_mode("opencl");
const std::string fpga_bitstream("empty"); // unused
const std::string linsolver("ilu0");
Dune::InverseOperatorResult result;
Vector<bz> x(rhs.size());
auto wellContribs = Opm::WellContributions::create("opencl", false);
std::unique_ptr<Opm::BdaBridge<Matrix<bz>, Vector<bz>, bz> > bridge;
try {
bridge = std::make_unique<Opm::BdaBridge<Matrix<bz>, Vector<bz>, bz> >(accelerator_mode,
fpga_bitstream,
@ -124,6 +120,15 @@ testOpenclSolver(const boost::property_tree::ptree& prm, Matrix<bz>& matrix, Vec
BOOST_WARN_MESSAGE(true, error.what());
throw PlatformInitException(error.what());
}
}
template <int bz>
Dune::BlockVector<Dune::FieldVector<double, bz>>
testOpenclSolver(std::unique_ptr<Opm::BdaBridge<Matrix<bz>, Vector<bz>, bz> >& bridge, Matrix<bz>& matrix, Vector<bz>& rhs)
{
Dune::InverseOperatorResult result;
Vector<bz> x(rhs.size());
auto wellContribs = Opm::WellContributions::create("opencl", false);
auto mat2 = matrix; // deep copy to make sure nnz values are in contiguous memory
// matrix created by readMatrixMarket() did not have contiguous memory
bridge->solve_system(&mat2, &mat2, /*numJacobiBlocks=*/0, rhs, *wellContribs, result);
@ -132,6 +137,23 @@ testOpenclSolver(const boost::property_tree::ptree& prm, Matrix<bz>& matrix, Vec
return x;
}
template <int bz>
Dune::BlockVector<Dune::FieldVector<double, bz>>
testOpenclSolverJacobi(std::unique_ptr<Opm::BdaBridge<Matrix<bz>, Vector<bz>, bz> >& bridge, Matrix<bz>& matrix, Vector<bz>& rhs)
{
Dune::InverseOperatorResult result;
Vector<bz> x(rhs.size());
auto wellContribs = Opm::WellContributions::create("opencl", false);
auto mat2 = matrix; // deep copy to make sure nnz values are in contiguous memory
// matrix created by readMatrixMarket() did not have contiguous memory
auto mat3 = matrix; // another deep copy, to make sure Jacobi matrix memory is different
// the sparsity pattern and values are actually the same
bridge->solve_system(&mat2, &mat3, /*numJacobiBlocks=*/2, rhs, *wellContribs, result);
bridge->get_result(x);
return x;
}
namespace pt = boost::property_tree;
void test3(const pt::ptree& prm)
@ -139,17 +161,31 @@ void test3(const pt::ptree& prm)
const int bz = 3;
Matrix<bz> matrix;
Vector<bz> rhs;
std::unique_ptr<Opm::BdaBridge<Matrix<bz>, Vector<bz>, bz> > bridge;
readLinearSystem("matr33.txt", "rhs3.txt", matrix, rhs);
Vector<bz> rhs2 = rhs; // deep copy, getDuneSolution() changes values in rhs vector
auto duneSolution = getDuneSolution<bz>(matrix, rhs);
auto sol = testOpenclSolver<bz>(prm, matrix, rhs2);
createBridge(prm, bridge); // create bridge with openclSolver
// should create bridge only once
// test openclSolver without Jacobi matrix
auto sol = testOpenclSolver<bz>(bridge, matrix, rhs2);
BOOST_REQUIRE_EQUAL(sol.size(), duneSolution.size());
for (size_t i = 0; i < sol.size(); ++i) {
for (int row = 0; row < bz; ++row) {
BOOST_CHECK_CLOSE(sol[i][row], duneSolution[i][row], 1e-3);
}
}
// test openclSolver with Jacobi matrix
auto solJacobi = testOpenclSolverJacobi<bz>(bridge, matrix, rhs2);
BOOST_REQUIRE_EQUAL(solJacobi.size(), duneSolution.size());
for (size_t i = 0; i < solJacobi.size(); ++i) {
for (int row = 0; row < bz; ++row) {
BOOST_CHECK_CLOSE(solJacobi[i][row], duneSolution[i][row], 1e-3);
}
}
}