Added block_size templates to the BdaBridge. Removed unused gpu_mode variable

This commit is contained in:
T.D. (Tongdong) Qiu 2020-06-24 16:46:04 +02:00
parent 26133c4fd7
commit 3dc368b0b4
3 changed files with 28 additions and 22 deletions

View File

@ -279,7 +279,8 @@ protected:
static const int numEq = Indices::numEq; static const int numEq = Indices::numEq;
#if HAVE_CUDA + HAVE_OPENCL #if HAVE_CUDA + HAVE_OPENCL
std::unique_ptr<BdaBridge> bdaBridge; static const unsigned int block_size = Matrix::block_type::rows;
std::unique_ptr<BdaBridge<Matrix, Vector, block_size>> bdaBridge;
#endif #endif
#if HAVE_MPI #if HAVE_MPI
@ -325,7 +326,7 @@ protected:
const int maxit = EWOMS_GET_PARAM(TypeTag, int, LinearSolverMaxIter); const int maxit = EWOMS_GET_PARAM(TypeTag, int, LinearSolverMaxIter);
const double tolerance = EWOMS_GET_PARAM(TypeTag, double, LinearSolverReduction); const double tolerance = EWOMS_GET_PARAM(TypeTag, double, LinearSolverReduction);
const int linear_solver_verbosity = parameters_.linear_solver_verbosity_; const int linear_solver_verbosity = parameters_.linear_solver_verbosity_;
bdaBridge.reset(new BdaBridge(gpu_mode, linear_solver_verbosity, maxit, tolerance)); bdaBridge.reset(new BdaBridge<Matrix, Vector, block_size>(gpu_mode, linear_solver_verbosity, maxit, tolerance));
#else #else
const std::string gpu_mode = EWOMS_GET_PARAM(TypeTag, std::string, GpuMode); const std::string gpu_mode = EWOMS_GET_PARAM(TypeTag, std::string, GpuMode);
if (gpu_mode.compare("none") != 0) { if (gpu_mode.compare("none") != 0) {

View File

@ -38,10 +38,9 @@ namespace Opm
using bda::BdaResult; using bda::BdaResult;
using bda::BdaSolver; using bda::BdaSolver;
BdaBridge::BdaBridge(std::string gpu_mode_, int linear_solver_verbosity, int maxit, double tolerance) template <class BridgeMatrix, class BridgeVector, int block_size>
: gpu_mode(gpu_mode_) BdaBridge<BridgeMatrix, BridgeVector, block_size>::BdaBridge(std::string gpu_mode, int linear_solver_verbosity, int maxit, double tolerance)
{ {
std::cout << "mode: " << gpu_mode_ << std::endl;
if (gpu_mode.compare("cusparse") == 0) { if (gpu_mode.compare("cusparse") == 0) {
#if HAVE_CUDA #if HAVE_CUDA
use_gpu = true; use_gpu = true;
@ -132,8 +131,8 @@ void getSparsityPattern(BridgeMatrix& mat, std::vector<int> &h_rows, std::vector
} // end getSparsityPattern() } // end getSparsityPattern()
template <class BridgeMatrix, class BridgeVector> template <class BridgeMatrix, class BridgeVector, int block_size>
void BdaBridge::solve_system(BridgeMatrix *mat OPM_UNUSED, BridgeVector &b OPM_UNUSED, WellContributions& wellContribs OPM_UNUSED, InverseOperatorResult &res OPM_UNUSED) void BdaBridge<BridgeMatrix, BridgeVector, block_size>::solve_system(BridgeMatrix *mat OPM_UNUSED, BridgeVector &b OPM_UNUSED, WellContributions& wellContribs OPM_UNUSED, InverseOperatorResult &res OPM_UNUSED)
{ {
if (use_gpu) { if (use_gpu) {
@ -207,20 +206,29 @@ void BdaBridge::solve_system(BridgeMatrix *mat OPM_UNUSED, BridgeVector &b OPM_U
} }
template <class BridgeVector> template <class BridgeMatrix, class BridgeVector, int block_size>
void BdaBridge::get_result(BridgeVector &x OPM_UNUSED) { void BdaBridge<BridgeMatrix, BridgeVector, block_size>::get_result(BridgeVector &x OPM_UNUSED) {
if (use_gpu) { if (use_gpu) {
backend->get_result(static_cast<double*>(&(x[0][0]))); backend->get_result(static_cast<double*>(&(x[0][0])));
} }
} }
#define INSTANTIATE_BDA_FUNCTIONS(n) \ #define INSTANTIATE_BDA_FUNCTIONS(n) \
template void BdaBridge::solve_system \ template BdaBridge<Dune::BCRSMatrix<Opm::MatrixBlock<double, n, n>, std::allocator<Opm::MatrixBlock<double, n, n> > >, \
Dune::BlockVector<Dune::FieldVector<double, n>, std::allocator<Dune::FieldVector<double, n> > >, \
n>::BdaBridge \
(std::string gpu_mode_, int linear_solver_verbosity, int maxit, double tolerance); \
\
template void BdaBridge<Dune::BCRSMatrix<Opm::MatrixBlock<double, n, n>, std::allocator<Opm::MatrixBlock<double, n, n> > >, \
Dune::BlockVector<Dune::FieldVector<double, n>, std::allocator<Dune::FieldVector<double, n> > >, \
n>::solve_system \
(Dune::BCRSMatrix<Opm::MatrixBlock<double, n, n>, std::allocator<Opm::MatrixBlock<double, n, n> > >*, \ (Dune::BCRSMatrix<Opm::MatrixBlock<double, n, n>, std::allocator<Opm::MatrixBlock<double, n, n> > >*, \
Dune::BlockVector<Dune::FieldVector<double, n>, std::allocator<Dune::FieldVector<double, n> > >&, \ Dune::BlockVector<Dune::FieldVector<double, n>, std::allocator<Dune::FieldVector<double, n> > >&, \
WellContributions&, InverseOperatorResult&); \ WellContributions&, InverseOperatorResult&); \
\ \
template void BdaBridge::get_result \ template void BdaBridge<Dune::BCRSMatrix<Opm::MatrixBlock<double, n, n>, std::allocator<Opm::MatrixBlock<double, n, n> > >, \
Dune::BlockVector<Dune::FieldVector<double, n>, std::allocator<Dune::FieldVector<double, n> > >, \
n>::get_result \
(Dune::BlockVector<Dune::FieldVector<double, n>, std::allocator<Dune::FieldVector<double, n> > >&); \ (Dune::BlockVector<Dune::FieldVector<double, n>, std::allocator<Dune::FieldVector<double, n> > >&); \
INSTANTIATE_BDA_FUNCTIONS(1); INSTANTIATE_BDA_FUNCTIONS(1);

View File

@ -42,13 +42,12 @@ namespace Opm
typedef Dune::InverseOperatorResult InverseOperatorResult; typedef Dune::InverseOperatorResult InverseOperatorResult;
/// BdaBridge acts as interface between opm-simulators with the cusparseSolver /// BdaBridge acts as interface between opm-simulators with the BdaSolvers
/// if CUDA was not found during CMake, function bodies of this class are empty template <class BridgeMatrix, class BridgeVector, int block_size>
class BdaBridge class BdaBridge
{ {
private: private:
std::unique_ptr<bda::BdaSolver> backend; std::unique_ptr<bda::BdaSolver> backend;
std::string gpu_mode;
bool use_gpu = false; bool use_gpu = false;
public: public:
@ -66,12 +65,10 @@ public:
/// \param[in] b vector b, should be of type Dune::BlockVector /// \param[in] b vector b, should be of type Dune::BlockVector
/// \param[in] wellContribs contains all WellContributions, to apply them separately, instead of adding them to matrix A /// \param[in] wellContribs contains all WellContributions, to apply them separately, instead of adding them to matrix A
/// \param[inout] result summary of solver result /// \param[inout] result summary of solver result
template <class BridgeMatrix, class BridgeVector>
void solve_system(BridgeMatrix *mat, BridgeVector &b, WellContributions& wellContribs, InverseOperatorResult &result); void solve_system(BridgeMatrix *mat, BridgeVector &b, WellContributions& wellContribs, InverseOperatorResult &result);
/// Get the resulting x vector /// Get the resulting x vector
/// \param[inout] x vector x, should be of type Dune::BlockVector /// \param[inout] x vector x, should be of type Dune::BlockVector
template <class BridgeVector>
void get_result(BridgeVector &x); void get_result(BridgeVector &x);
/// Return whether the BdaBridge will use the GPU or not /// Return whether the BdaBridge will use the GPU or not