mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Moved warning prints from cusparseSolver to BdaBridge. Moved solving logic from BdaBridge to cusparseSolver. Added cusparseSolverStatus enum.
This commit is contained in:
parent
b6e13bffd2
commit
a4cb582cfb
@ -20,6 +20,8 @@
|
|||||||
#include <config.h>
|
#include <config.h>
|
||||||
#include <memory>
|
#include <memory>
|
||||||
|
|
||||||
|
#include <opm/common/OpmLog/OpmLog.hpp>
|
||||||
|
|
||||||
#include <opm/simulators/linalg/bda/BdaBridge.hpp>
|
#include <opm/simulators/linalg/bda/BdaBridge.hpp>
|
||||||
#include <opm/simulators/linalg/bda/BdaResult.hpp>
|
#include <opm/simulators/linalg/bda/BdaResult.hpp>
|
||||||
|
|
||||||
@ -91,7 +93,7 @@ int checkZeroDiagonal(BridgeMatrix& mat) {
|
|||||||
// if only_vals, do not convert rowPointers and colIndices
|
// if only_vals, do not convert rowPointers and colIndices
|
||||||
// sparsity pattern should stay the same due to matrix-add-well-contributions
|
// sparsity pattern should stay the same due to matrix-add-well-contributions
|
||||||
template <class BridgeMatrix>
|
template <class BridgeMatrix>
|
||||||
void convertMatrixBsr(BridgeMatrix& mat, std::vector<double> &h_vals, std::vector<int> &h_rows, std::vector<int> &h_cols, int dim, bool only_vals) {
|
void convertMatrixBsr(BridgeMatrix& mat, std::vector<double> &h_vals, std::vector<int> &h_rows, std::vector<int> &h_cols, int dim) {
|
||||||
int sum_nnzs = 0;
|
int sum_nnzs = 0;
|
||||||
int nnz = mat.nonzeroes()*dim*dim;
|
int nnz = mat.nonzeroes()*dim*dim;
|
||||||
|
|
||||||
@ -99,7 +101,7 @@ void convertMatrixBsr(BridgeMatrix& mat, std::vector<double> &h_vals, std::vecto
|
|||||||
memcpy(h_vals.data(), &(mat[0][0][0][0]), sizeof(double)*nnz);
|
memcpy(h_vals.data(), &(mat[0][0][0][0]), sizeof(double)*nnz);
|
||||||
|
|
||||||
// convert colIndices and rowPointers
|
// convert colIndices and rowPointers
|
||||||
if(only_vals == false){
|
if(h_rows.size() == 0){
|
||||||
h_rows.emplace_back(0);
|
h_rows.emplace_back(0);
|
||||||
for(typename BridgeMatrix::const_iterator r = mat.begin(); r != mat.end(); ++r){
|
for(typename BridgeMatrix::const_iterator r = mat.begin(); r != mat.end(); ++r){
|
||||||
int size_row = 0;
|
int size_row = 0;
|
||||||
@ -160,13 +162,11 @@ void BdaBridge::solve_system(BridgeMatrix *mat, BridgeVector &b, InverseOperator
|
|||||||
checkZeroDiagonal(*mat);
|
checkZeroDiagonal(*mat);
|
||||||
#endif
|
#endif
|
||||||
|
|
||||||
bool initialized = backend->isInitialized();
|
|
||||||
|
|
||||||
#if PRINT_TIMERS_BRIDGE
|
#if PRINT_TIMERS_BRIDGE
|
||||||
Dune::Timer t;
|
Dune::Timer t;
|
||||||
#endif
|
#endif
|
||||||
|
|
||||||
convertMatrixBsr(*mat, h_vals, h_rows, h_cols, dim, initialized);
|
convertMatrixBsr(*mat, h_vals, h_rows, h_cols, dim);
|
||||||
convertBlockVectorToArray(b, h_b);
|
convertBlockVectorToArray(b, h_b);
|
||||||
|
|
||||||
#if PRINT_TIMERS_BRIDGE
|
#if PRINT_TIMERS_BRIDGE
|
||||||
@ -176,17 +176,23 @@ void BdaBridge::solve_system(BridgeMatrix *mat, BridgeVector &b, InverseOperator
|
|||||||
/////////////////////////
|
/////////////////////////
|
||||||
// actually solve
|
// actually solve
|
||||||
|
|
||||||
if(initialized == false){
|
typedef cusparseSolverBackend::cusparseSolverStatus cusparseSolverStatus;
|
||||||
backend->initialize(N, nnz, dim);
|
|
||||||
backend->copy_system_to_gpu(h_vals.data(), h_rows.data(), h_cols.data(), h_b.data());
|
cusparseSolverStatus status = backend->solve_system(N, nnz, dim, h_vals.data(), h_rows.data(), h_cols.data(), h_b.data(), result);
|
||||||
backend->analyse_matrix();
|
switch(status){
|
||||||
}else{
|
case cusparseSolverStatus::CUSPARSE_SOLVER_SUCCESS:
|
||||||
backend->update_system_on_gpu(h_vals.data(), h_b.data());
|
//OpmLog::info("cusparseSolver converged");
|
||||||
}
|
break;
|
||||||
backend->reset_prec_on_gpu();
|
case cusparseSolverStatus::CUSPARSE_SOLVER_ANALYSIS_FAILED:
|
||||||
if(backend->create_preconditioner()){
|
OpmLog::warning("cusparseSolver could not analyse level information of matrix, perhaps there is still a 0.0 on the diagonal of a block on the diagonal");
|
||||||
backend->solve_system(result);
|
break;
|
||||||
|
case cusparseSolverStatus::CUSPARSE_SOLVER_CREATE_PRECONDITIONER_FAILED:
|
||||||
|
OpmLog::warning("cusparseSolver could not create preconditioner, perhaps there is still a 0.0 on the diagonal of a block on the diagonal");
|
||||||
|
break;
|
||||||
|
default:
|
||||||
|
OpmLog::warning("cusparseSolver returned unknown status code");
|
||||||
}
|
}
|
||||||
|
|
||||||
res.iterations = result.iterations;
|
res.iterations = result.iterations;
|
||||||
res.reduction = result.reduction;
|
res.reduction = result.reduction;
|
||||||
res.converged = result.converged;
|
res.converged = result.converged;
|
||||||
|
@ -55,8 +55,7 @@ namespace Opm
|
|||||||
finalize();
|
finalize();
|
||||||
}
|
}
|
||||||
|
|
||||||
// return true iff converged
|
void cusparseSolverBackend::gpu_pbicgstab(BdaResult& res){
|
||||||
bool cusparseSolverBackend::gpu_pbicgstab(BdaResult& res){
|
|
||||||
double t_total1, t_total2;
|
double t_total1, t_total2;
|
||||||
int n = N;
|
int n = N;
|
||||||
double rho = 1.0, rhop;
|
double rho = 1.0, rhop;
|
||||||
@ -163,7 +162,6 @@ namespace Opm
|
|||||||
if(verbosity > 0){
|
if(verbosity > 0){
|
||||||
printf("=== converged: %d, conv_rate: %.2f, time: %f, time per iteration: %f, iterations: %.1f\n", res.converged, res.conv_rate, res.elapsed, res.elapsed/it, it);
|
printf("=== converged: %d, conv_rate: %.2f, time: %f, time per iteration: %f, iterations: %.1f\n", res.converged, res.conv_rate, res.elapsed, res.elapsed/it, it);
|
||||||
}
|
}
|
||||||
return res.converged;
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
@ -308,7 +306,7 @@ namespace Opm
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
void cusparseSolverBackend::analyse_matrix(){
|
bool cusparseSolverBackend::analyse_matrix(){
|
||||||
|
|
||||||
int d_bufferSize_M, d_bufferSize_L, d_bufferSize_U, d_bufferSize;
|
int d_bufferSize_M, d_bufferSize_L, d_bufferSize_U, d_bufferSize;
|
||||||
double t1, t2;
|
double t1, t2;
|
||||||
@ -367,9 +365,7 @@ namespace Opm
|
|||||||
int structural_zero;
|
int structural_zero;
|
||||||
cusparseStatus_t status = cusparseXbsrilu02_zeroPivot(cusparseHandle, info_M, &structural_zero);
|
cusparseStatus_t status = cusparseXbsrilu02_zeroPivot(cusparseHandle, info_M, &structural_zero);
|
||||||
if(CUSPARSE_STATUS_ZERO_PIVOT == status){
|
if(CUSPARSE_STATUS_ZERO_PIVOT == status){
|
||||||
fprintf(stderr, "ERROR block U(%d,%d) is not invertible\n", structural_zero, structural_zero);
|
return false;
|
||||||
fprintf(stderr, "cusparse fails when a block has a 0.0 on its diagonal, these should be replaced in BdaBridge::checkZeroDiagonal()\n");
|
|
||||||
exit(1);
|
|
||||||
}
|
}
|
||||||
|
|
||||||
// analysis of ilu apply
|
// analysis of ilu apply
|
||||||
@ -388,6 +384,7 @@ namespace Opm
|
|||||||
printf("cusparseSolver::analyse_matrix(): %f s\n", t2-t1);
|
printf("cusparseSolver::analyse_matrix(): %f s\n", t2-t1);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
return true;
|
||||||
} // end analyse_matrix()
|
} // end analyse_matrix()
|
||||||
|
|
||||||
bool cusparseSolverBackend::create_preconditioner(){
|
bool cusparseSolverBackend::create_preconditioner(){
|
||||||
@ -407,8 +404,6 @@ namespace Opm
|
|||||||
// cusparseXbsrilu02_zeroPivot() calls cudaDeviceSynchronize()
|
// cusparseXbsrilu02_zeroPivot() calls cudaDeviceSynchronize()
|
||||||
cusparseStatus_t status = cusparseXbsrilu02_zeroPivot(cusparseHandle, info_M, &structural_zero);
|
cusparseStatus_t status = cusparseXbsrilu02_zeroPivot(cusparseHandle, info_M, &structural_zero);
|
||||||
if(CUSPARSE_STATUS_ZERO_PIVOT == status){
|
if(CUSPARSE_STATUS_ZERO_PIVOT == status){
|
||||||
fprintf(stderr, "WARNING block U(%d,%d) is not invertible\n", structural_zero, structural_zero);
|
|
||||||
fprintf(stderr, "cusparse fails when a block has a 0.0 on its diagonal, these should be replaced in BdaBridge::checkZeroDiagonal()\n");
|
|
||||||
return false;
|
return false;
|
||||||
}
|
}
|
||||||
|
|
||||||
@ -421,13 +416,11 @@ namespace Opm
|
|||||||
} // end create_preconditioner()
|
} // end create_preconditioner()
|
||||||
|
|
||||||
|
|
||||||
// return true iff converged
|
void cusparseSolverBackend::solve_system(BdaResult &res){
|
||||||
bool cusparseSolverBackend::solve_system(BdaResult &res){
|
|
||||||
// actually solve
|
// actually solve
|
||||||
bool converged = gpu_pbicgstab(res);
|
gpu_pbicgstab(res);
|
||||||
cudaStreamSynchronize(stream);
|
cudaStreamSynchronize(stream);
|
||||||
cudaCheckLastError("Something went wrong during the GPU solve");
|
cudaCheckLastError("Something went wrong during the GPU solve");
|
||||||
return converged;
|
|
||||||
} // end solve_system()
|
} // end solve_system()
|
||||||
|
|
||||||
|
|
||||||
@ -454,9 +447,29 @@ namespace Opm
|
|||||||
} // end post_process()
|
} // end post_process()
|
||||||
|
|
||||||
|
|
||||||
bool cusparseSolverBackend::isInitialized(){
|
typedef cusparseSolverBackend::cusparseSolverStatus cusparseSolverStatus;
|
||||||
return initialized;
|
|
||||||
|
cusparseSolverStatus cusparseSolverBackend::solve_system(int N, int nnz, int dim, double *vals, int *rows, int *cols, double *b, BdaResult &res){
|
||||||
|
if(initialized == false){
|
||||||
|
initialize(N, nnz, dim);
|
||||||
|
copy_system_to_gpu(vals, rows, cols, b);
|
||||||
|
}else{
|
||||||
|
update_system_on_gpu(vals, b);
|
||||||
}
|
}
|
||||||
|
if(analysis_done == false){
|
||||||
|
if(!analyse_matrix()){
|
||||||
|
return cusparseSolverStatus::CUSPARSE_SOLVER_ANALYSIS_FAILED;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
reset_prec_on_gpu();
|
||||||
|
if(create_preconditioner()){
|
||||||
|
solve_system(res);
|
||||||
|
}else{
|
||||||
|
return cusparseSolverStatus::CUSPARSE_SOLVER_CREATE_PRECONDITIONER_FAILED;
|
||||||
|
}
|
||||||
|
return cusparseSolverStatus::CUSPARSE_SOLVER_SUCCESS;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
@ -58,6 +58,7 @@ private:
|
|||||||
int BLOCK_SIZE;
|
int BLOCK_SIZE;
|
||||||
|
|
||||||
bool initialized = false;
|
bool initialized = false;
|
||||||
|
bool analysis_done = false;
|
||||||
|
|
||||||
// verbosity
|
// verbosity
|
||||||
// 0: print nothing during solves, only when initializing
|
// 0: print nothing during solves, only when initializing
|
||||||
@ -67,14 +68,8 @@ private:
|
|||||||
|
|
||||||
int verbosity = 0;
|
int verbosity = 0;
|
||||||
|
|
||||||
public:
|
|
||||||
|
|
||||||
cusparseSolverBackend(int linear_solver_verbosity, int maxit, double tolerance);
|
void gpu_pbicgstab(BdaResult& res);
|
||||||
|
|
||||||
~cusparseSolverBackend();
|
|
||||||
|
|
||||||
// return true iff converged
|
|
||||||
bool gpu_pbicgstab(BdaResult& res);
|
|
||||||
|
|
||||||
void initialize(int N, int nnz, int dim);
|
void initialize(int N, int nnz, int dim);
|
||||||
|
|
||||||
@ -87,17 +82,29 @@ public:
|
|||||||
|
|
||||||
void reset_prec_on_gpu();
|
void reset_prec_on_gpu();
|
||||||
|
|
||||||
void analyse_matrix();
|
bool analyse_matrix();
|
||||||
|
|
||||||
bool create_preconditioner();
|
bool create_preconditioner();
|
||||||
|
|
||||||
// return true iff converged
|
void solve_system(BdaResult &res);
|
||||||
bool solve_system(BdaResult &res);
|
|
||||||
|
public:
|
||||||
|
|
||||||
|
enum class cusparseSolverStatus {
|
||||||
|
CUSPARSE_SOLVER_SUCCESS,
|
||||||
|
CUSPARSE_SOLVER_ANALYSIS_FAILED,
|
||||||
|
CUSPARSE_SOLVER_CREATE_PRECONDITIONER_FAILED,
|
||||||
|
CUSPARSE_SOLVER_UNKNOWN_ERROR
|
||||||
|
};
|
||||||
|
|
||||||
|
cusparseSolverBackend(int linear_solver_verbosity, int maxit, double tolerance);
|
||||||
|
|
||||||
|
~cusparseSolverBackend();
|
||||||
|
|
||||||
|
cusparseSolverStatus solve_system(int N, int nnz, int dim, double *vals, int *rows, int *cols, double *b, BdaResult &res);
|
||||||
|
|
||||||
void post_process(double *x);
|
void post_process(double *x);
|
||||||
|
|
||||||
bool isInitialized();
|
|
||||||
|
|
||||||
}; // end class cusparseSolverBackend
|
}; // end class cusparseSolverBackend
|
||||||
|
|
||||||
}
|
}
|
||||||
|
Loading…
Reference in New Issue
Block a user