Set hipStream for rocsparseWellcontributions

This commit is contained in:
Tong Dong Qiu
2023-04-14 09:35:15 +02:00
committed by Razvan Nane
parent a49aaf53d2
commit 9bef10a018
3 changed files with 15 additions and 0 deletions

View File

@@ -153,6 +153,12 @@ void rocsparseSolverBackend<block_size>::gpu_pbicgstab([[maybe_unused]] WellCont
Timer t_total, t_prec(false), t_spmv(false), t_rest(false);
// set stream here, the WellContributions object is destroyed every linear solve
// the number of wells can change every linear solve
if(wellContribs.getNumWells() > 0){
static_cast<WellContributionsRocsparse&>(wellContribs).setStream(stream);
}
// HIP_VERSION is defined as (HIP_VERSION_MAJOR * 10000000 + HIP_VERSION_MINOR * 100000 + HIP_VERSION_PATCH)
#if HIP_VERSION >= 50400000
ROCSPARSE_CHECK(rocsparse_dbsrmv_ex(handle, dir, operation,

View File

@@ -163,6 +163,10 @@ void WellContributionsRocsparse::apply(double *d_x, double *d_y){
}
}
void WellContributionsRocsparse::setStream(hipStream_t stream_){
stream = stream_;
}
void WellContributionsRocsparse::APIaddMatrix(MatrixType type,
int* colIndices,
double* values,

View File

@@ -22,6 +22,7 @@
#include <opm/simulators/linalg/bda/WellContributions.hpp>
#include <hip/hip_runtime_api.h>
#include <vector>
@@ -31,10 +32,14 @@ namespace Opm
class WellContributionsRocsparse : public WellContributions
{
private:
hipStream_t stream;
public:
void apply_stdwells(double *d_x, double *d_y);
void apply_mswells(double *d_x, double *d_y);
void apply(double *d_x, double *d_y);
void setStream(hipStream_t stream);
protected:
/// Allocate memory for the StandardWells