Removed OpenCL dependencies from rocsparse backend

This commit is contained in:
Justin Chang 2023-09-28 00:22:18 -05:00
parent d1b908421d
commit 54d82d349e
4 changed files with 13 additions and 11 deletions

View File

@ -92,12 +92,12 @@ struct EdgeWeightsMethod {
using type = UndefinedProperty;
};
#if HAVE_OPENCL
#if HAVE_OPENCL || HAVE_ROCSPARSE
template<class TypeTag, class MyTypeTag>
struct NumJacobiBlocks {
using type = UndefinedProperty;
};
#endif // HAVE_OPENCL
#endif // HAVE_OPENCL || HAVE_ROCSPARSE
template<class TypeTag, class MyTypeTag>
struct OwnerCellsFirst {
@ -153,12 +153,12 @@ struct EdgeWeightsMethod<TypeTag, TTag::EclBaseVanguard> {
static constexpr int value = 1;
};
#if HAVE_OPENCL
#if HAVE_OPENCL || HAVE_ROCSPARSE
template<class TypeTag>
struct NumJacobiBlocks<TypeTag, TTag::EclBaseVanguard> {
static constexpr int value = 0;
};
#endif // HAVE_OPENCL
#endif // HAVE_OPENCL || HAVE_ROCSPARSE
template<class TypeTag>
struct OwnerCellsFirst<TypeTag, TTag::EclBaseVanguard> {
@ -248,7 +248,7 @@ public:
EWOMS_REGISTER_PARAM(TypeTag, int, EdgeWeightsMethod,
"Choose edge-weighing strategy: 0=uniform, 1=trans, 2=log(trans).");
#if HAVE_OPENCL
#if HAVE_OPENCL || HAVE_ROCSPARSE
EWOMS_REGISTER_PARAM(TypeTag, int, NumJacobiBlocks,
"Number of blocks to be created for the Block-Jacobi preconditioner.");
#endif
@ -288,7 +288,7 @@ public:
fileName_ = EWOMS_GET_PARAM(TypeTag, std::string, EclDeckFileName);
edgeWeightsMethod_ = Dune::EdgeWeightMethod(EWOMS_GET_PARAM(TypeTag, int, EdgeWeightsMethod));
#if HAVE_OPENCL
#if HAVE_OPENCL || HAVE_ROCSPARSE
numJacobiBlocks_ = EWOMS_GET_PARAM(TypeTag, int, NumJacobiBlocks);
#endif

View File

@ -202,7 +202,7 @@ public:
*/
int numJacobiBlocks() const
{
#if HAVE_OPENCL
#if HAVE_OPENCL || HAVE_ROCSPARSE
return numJacobiBlocks_;
#else
return 0;
@ -283,9 +283,9 @@ protected:
std::string fileName_;
Dune::EdgeWeightMethod edgeWeightsMethod_;
#if HAVE_OPENCL
#if HAVE_OPENCL || HAVE_ROCSPARSE
int numJacobiBlocks_{0};
#endif // HAVE_OPENCL
#endif // HAVE_OPENCL || HAVE_ROCSPARSE
bool ownersFirst_;
#if HAVE_MPI

View File

@ -526,6 +526,9 @@ void rocsparseSolverBackend<block_size>::solve_system(WellContributions &wellCon
Timer t;
// actually solve
gpu_pbicgstab(wellContribs, res);
HIP_CHECK(hipStreamSynchronize(stream));
#if 0
try {
gpu_pbicgstab(wellContribs, res);
} catch (const cl::Error& error) {
@ -538,9 +541,9 @@ void rocsparseSolverBackend<block_size>::solve_system(WellContributions &wellCon
// rethrow exception by OPM_THROW in the try{}, without this, a segfault occurs
throw error;
}
#endif
if (verbosity >= 3) {
HIP_CHECK(hipStreamSynchronize(stream));
std::ostringstream out;
out << "rocsparseSolver::solve_system(): " << t.stop() << " s";
OpmLog::info(out.str());

View File

@ -22,7 +22,6 @@
#include <memory>
#include <opm/simulators/linalg/bda/opencl/opencl.hpp>
#include <opm/simulators/linalg/bda/BdaResult.hpp>
#include <opm/simulators/linalg/bda/BdaSolver.hpp>
#include <opm/simulators/linalg/bda/WellContributions.hpp>