Allow WellContributionsRocsparse to be used

This commit is contained in:
Tong Dong Qiu
2023-03-28 16:31:50 +02:00
committed by Razvan Nane
parent 50ccea0880
commit cb20d041c3
2 changed files with 30 additions and 4 deletions

View File

@@ -32,6 +32,10 @@
#include <opm/simulators/linalg/bda/cuda/cuWellContributions.hpp>
#endif
#ifdef HAVE_ROCSPARSE
#include <opm/simulators/linalg/bda/rocsparseWellContributions.hpp>
#endif
namespace Opm
{
@@ -54,9 +58,14 @@ WellContributions::create(const std::string& accelerator_mode, bool useWellConn)
}
else if(accelerator_mode.compare("rocsparse") == 0){
if (!useWellConn) {
OPM_THROW(std::logic_error, "Error rocsparse requires --matrix-add-well-contributions=true");
#if HAVE_ROCSPARSE
return std::make_unique<WellContributionsRocsparse>();
#else
OPM_THROW(std::runtime_error, "Cannot initialize well contributions: rocsparse is not enabled");
#endif
}
return std::make_unique<WellContributions>();
}
else if(accelerator_mode.compare("amgcl") == 0){
if (!useWellConn) {

View File

@@ -37,6 +37,7 @@
#undef HAVE_CUDA
#include <opm/simulators/linalg/bda/rocsparseSolverBackend.hpp>
#include <opm/simulators/linalg/bda/rocsparseWellContributions.hpp>
#include <opm/simulators/linalg/bda/BdaResult.hpp>
@@ -152,7 +153,7 @@ void rocsparseSolverBackend<block_size>::gpu_pbicgstab([[maybe_unused]] WellCont
double one = 1.0;
double mone = -1.0;
Timer t_total, t_prec(false), t_spmv(false), t_rest(false);
Timer t_total, t_prec(false), t_spmv(false), t_well(false), t_rest(false);
// set stream here, the WellContributions object is destroyed every linear solve
// the number of wells can change every linear solve
@@ -232,10 +233,18 @@ void rocsparseSolverBackend<block_size>::gpu_pbicgstab([[maybe_unused]] WellCont
if (verbosity >= 3) {
HIP_CHECK(hipStreamSynchronize(stream));
t_spmv.stop();
t_rest.start();
t_well.start();
}
// apply wellContributions
if(wellContribs.getNumWells() > 0){
static_cast<WellContributionsRocsparse&>(wellContribs).apply(d_pw, d_v);
}
if (verbosity >= 3) {
HIP_CHECK(hipStreamSynchronize(stream));
t_well.stop();
t_rest.start();
}
ROCBLAS_CHECK(rocblas_ddot(blas_handle, N, d_rw, 1, d_v, 1, &tmp1));
alpha = rho / tmp1;
@@ -285,10 +294,18 @@ void rocsparseSolverBackend<block_size>::gpu_pbicgstab([[maybe_unused]] WellCont
if(verbosity >= 3){
HIP_CHECK(hipStreamSynchronize(stream));
t_spmv.stop();
t_rest.start();
t_well.start();
}
// apply wellContributions
if(wellContribs.getNumWells() > 0){
static_cast<WellContributionsRocsparse&>(wellContribs).apply(d_s, d_t);
}
if (verbosity >= 3) {
HIP_CHECK(hipStreamSynchronize(stream));
t_well.stop();
t_rest.start();
}
ROCBLAS_CHECK(rocblas_ddot(blas_handle, N, d_t, 1, d_r, 1, &tmp1));
ROCBLAS_CHECK(rocblas_ddot(blas_handle, N, d_t, 1, d_t, 1, &tmp2));