mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Update GPU solver and solve_transposed_3x3() tests
This commit is contained in:
parent
374f8276dc
commit
a5ed003418
@ -76,15 +76,16 @@ testCusparseSolver(const boost::property_tree::ptree& prm, const std::string& ma
|
||||
const std::string opencl_ilu_reorder("none"); // unused
|
||||
const int platformID = 0; // unused
|
||||
const int deviceID = 0;
|
||||
const std::string gpu_mode("cusparse");
|
||||
const std::string accelerator_mode("cusparse");
|
||||
const std::string fpga_bitstream("empty"); // unused
|
||||
const std::string linsolver("ilu0");
|
||||
Dune::InverseOperatorResult result;
|
||||
|
||||
Vector x(rhs.size());
|
||||
auto wellContribs = Opm::WellContributions::create("cusparse", false);
|
||||
std::unique_ptr<Opm::BdaBridge<Matrix, Vector, bz> > bridge;
|
||||
try {
|
||||
bridge = std::make_unique<Opm::BdaBridge<Matrix, Vector, bz> >(gpu_mode, fpga_bitstream, linear_solver_verbosity, maxit, tolerance, platformID, deviceID, opencl_ilu_reorder);
|
||||
bridge = std::make_unique<Opm::BdaBridge<Matrix, Vector, bz> >(accelerator_mode, fpga_bitstream, linear_solver_verbosity, maxit, tolerance, platformID, deviceID, opencl_ilu_reorder, linsolver);
|
||||
|
||||
bridge->solve_system(&matrix, rhs, *wellContribs, result);
|
||||
bridge->get_result(x);
|
||||
|
@ -75,15 +75,16 @@ testOpenclSolver(const boost::property_tree::ptree& prm, const std::string& matr
|
||||
const std::string opencl_ilu_reorder("none");
|
||||
const int platformID = 0;
|
||||
const int deviceID = 0;
|
||||
const std::string gpu_mode("opencl");
|
||||
const std::string accelerator_mode("opencl");
|
||||
const std::string fpga_bitstream("empty"); // unused
|
||||
const std::string linsolver("ilu0");
|
||||
Dune::InverseOperatorResult result;
|
||||
|
||||
Vector x(rhs.size());
|
||||
auto wellContribs = Opm::WellContributions::create("opencl", false);
|
||||
std::unique_ptr<Opm::BdaBridge<Matrix, Vector, bz> > bridge;
|
||||
try {
|
||||
bridge = std::make_unique<Opm::BdaBridge<Matrix, Vector, bz> >(gpu_mode, fpga_bitstream, linear_solver_verbosity, maxit, tolerance, platformID, deviceID, opencl_ilu_reorder);
|
||||
bridge = std::make_unique<Opm::BdaBridge<Matrix, Vector, bz> >(accelerator_mode, fpga_bitstream, linear_solver_verbosity, maxit, tolerance, platformID, deviceID, opencl_ilu_reorder, linsolver);
|
||||
} catch (const std::logic_error& error) {
|
||||
BOOST_WARN_MESSAGE(true, error.what());
|
||||
throw PlatformInitException(error.what());
|
||||
@ -100,8 +101,8 @@ void test3(const pt::ptree& prm)
|
||||
{
|
||||
const int bz = 3;
|
||||
auto sol = testOpenclSolver<bz>(prm, "matr33.txt", "rhs3.txt");
|
||||
Dune::BlockVector<Dune::FieldVector<double, bz>> expected {{-1.30307e-2, -3.58263e-6, 1.13836e-9},
|
||||
{-1.25425e-3, -1.4167e-4, -3.2213e-3},
|
||||
Dune::BlockVector<Dune::FieldVector<double, bz>> expected {{-0.0131626, -3.58263e-6, 1.13836e-9},
|
||||
{-1.25425e-3, -1.4167e-4, -0.0029366},
|
||||
{-4.5436e-4, 1.28682e-5, 4.7644e-6}};
|
||||
BOOST_REQUIRE_EQUAL(sol.size(), expected.size());
|
||||
for (size_t i = 0; i < sol.size(); ++i) {
|
||||
|
@ -32,24 +32,48 @@ BOOST_AUTO_TEST_CASE(testsolvetransposed3x3)
|
||||
const unsigned numTests = 3;
|
||||
const unsigned blockSize = 3;
|
||||
|
||||
std::vector<std::vector<double> > A = {{-0.0973322, 1.18758e-08, -0.000185231,
|
||||
0.125114, 1.24919e-11, 0,
|
||||
-11.3854, 1.32995e-06, 0.065447},
|
||||
{-0.00354818, 2.24354e-08, 55.9563,
|
||||
0.284031, 2.83543e-11, 0,
|
||||
-62.1156, 0.000392767, 6350.26},
|
||||
{1, 0, 2,
|
||||
3, 1, 0,
|
||||
0, 0, 3}};
|
||||
std::vector<std::vector<double> > A = {{4, 2, 1,
|
||||
3, 4, 2,
|
||||
2, 4, 3},
|
||||
{1, 2, 4,
|
||||
1, 3, 2,
|
||||
2, 4, 2},
|
||||
{1, 2, 2,
|
||||
1, 3, 5,
|
||||
3, 2, 4}};
|
||||
|
||||
std::vector<double> b = {0, 1, 0};
|
||||
std::vector<double> x(blockSize);
|
||||
std::vector<std::vector<double> > x_expected = {{63886256.67, 66154278, 180813.715},
|
||||
{-290820.58, 556792.09, 2562.61063},
|
||||
{-3, 1, 2}};
|
||||
std::vector<std::vector<double> > b = {{0, 1, 0},
|
||||
{1, 3, 5},
|
||||
{2, 4, 5}};
|
||||
|
||||
std::vector<std::vector<double> > x_expected = {{-0.5, 1, -0.5},
|
||||
{ 1, 1, -0.5},
|
||||
{ 1.3, 0.4, 0.1}};
|
||||
|
||||
for (unsigned testId = 0; testId < numTests; ++testId) {
|
||||
Opm::Accelerator::solve_transposed_3x3(A[testId].data(), b.data(), x.data());
|
||||
std::vector<double> x(blockSize);
|
||||
Opm::Accelerator::solve_transposed_3x3(A[testId].data(), b[testId].data(), x.data());
|
||||
|
||||
for (unsigned i = 0; i < blockSize; ++i) {
|
||||
BOOST_CHECK_CLOSE(x[i], x_expected[testId][i], 1e-6);
|
||||
}
|
||||
}
|
||||
|
||||
// retest cases using Dune methods
|
||||
using Mat3 = Dune::FieldMatrix<double, 3, 3>;
|
||||
using Vec3 = Dune::FieldVector<double, 3>;
|
||||
for (unsigned testId = 0; testId < numTests; ++testId) {
|
||||
Mat3 a3 = {{ A[testId][0], A[testId][1], A[testId][2] },
|
||||
{ A[testId][3], A[testId][4], A[testId][5] },
|
||||
{ A[testId][6], A[testId][7], A[testId][8] } };
|
||||
Vec3 y({b[testId][0], b[testId][1], b[testId][2]});
|
||||
Mat3 b3 = a3;
|
||||
|
||||
// b3 = inv(a3^T)
|
||||
Dune::FMatrixHelp::invertMatrix_retTransposed(a3, b3);
|
||||
|
||||
// x = b3 * y
|
||||
Vec3 x = Dune::FMatrixHelp::mult(b3, y);
|
||||
|
||||
for (unsigned i = 0; i < blockSize; ++i) {
|
||||
BOOST_CHECK_CLOSE(x[i], x_expected[testId][i], 1e-6);
|
||||
|
Loading…
Reference in New Issue
Block a user