Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Clean up CUDA solvers and use override specifier. #101

Merged
merged 1 commit into from
Dec 6, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
30 changes: 15 additions & 15 deletions resolve/LinSolverDirectCuSolverGLU.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -59,11 +59,11 @@ namespace ReSolve
n,
nnz,
descr_A_,
A_->getRowData(memory::HOST), //kRowPtr_,
A_->getColData(memory::HOST), //jCol_,
P, /* base-0 */
A_->getRowData(memory::HOST),
A_->getColData(memory::HOST),
P, /* base-0 */
Q, /* base-0 */
M_->getNnz(), /* nnzM */
M_->getNnz(), /* nnzM */
descr_M_,
M_->getRowData(memory::HOST),
M_->getColData(memory::HOST),
Expand All @@ -86,9 +86,9 @@ namespace ReSolve
/* A is original matrix */
nnz,
descr_A_,
A_->getValues( memory::DEVICE), //da_,
A_->getRowData(memory::DEVICE), //kRowPtr_,
A_->getColData(memory::DEVICE), //jCol_,
A_->getValues( memory::DEVICE),
A_->getRowData(memory::DEVICE),
A_->getColData(memory::DEVICE),
info_M_);
error_sum += status_cusolver_;

Expand All @@ -100,7 +100,7 @@ namespace ReSolve

void LinSolverDirectCuSolverGLU::addFactors(matrix::Sparse* L, matrix::Sparse* U)
{
// L and U need to be in CSC format
// L and U need to be in CSC format
index_type n = L->getNumRows();
index_type* Lp = L->getColData(memory::HOST);
index_type* Li = L->getRowData(memory::HOST);
Expand Down Expand Up @@ -162,9 +162,9 @@ namespace ReSolve
/* A is original matrix */
A_->getNnzExpanded(),
descr_A_,
A_->getValues( memory::DEVICE), //da_,
A_->getRowData(memory::DEVICE), //kRowPtr_,
A_->getColData(memory::DEVICE), //jCol_,
A_->getValues( memory::DEVICE),
A_->getRowData(memory::DEVICE),
A_->getColData(memory::DEVICE),
info_M_);
error_sum += status_cusolver_;

Expand All @@ -181,11 +181,11 @@ namespace ReSolve
/* A is original matrix */
A_->getNnz(),
descr_A_,
A_->getValues( memory::DEVICE), //da_,
A_->getRowData(memory::DEVICE), //kRowPtr_,
A_->getColData(memory::DEVICE), //jCol_,
A_->getValues( memory::DEVICE),
A_->getRowData(memory::DEVICE),
A_->getColData(memory::DEVICE),
rhs->getData(memory::DEVICE),/* right hand side */
x->getData(memory::DEVICE),/* left hand side */
x->getData(memory::DEVICE), /* left hand side */
&ite_refine_succ_,
&r_nrminf_,
info_M_,
Expand Down
14 changes: 7 additions & 7 deletions resolve/LinSolverDirectCuSolverGLU.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -28,24 +28,24 @@ namespace ReSolve
LinSolverDirectCuSolverGLU(LinAlgWorkspaceCUDA* workspace);
~LinSolverDirectCuSolverGLU();

int refactorize();
int solve(vector_type* rhs, vector_type* x);
int solve(vector_type* x);
int refactorize() override;
int solve(vector_type* rhs, vector_type* x) override;
int solve(vector_type* x) override;

int setup(matrix::Sparse* A,
matrix::Sparse* L,
matrix::Sparse* U,
index_type* P,
index_type* Q,
vector_type* rhs = nullptr);
vector_type* rhs = nullptr) override;

private:
void addFactors(matrix::Sparse* L, matrix::Sparse* U); //create L+U from sepeate L, U factors
matrix::Sparse* M_;//the matrix that contains added factors
void addFactors(matrix::Sparse* L, matrix::Sparse* U); ///< creates L+U from sepeate L, U factors
matrix::Sparse* M_; ///< the matrix that contains added factors
//note: we need cuSolver handle, we can copy it from the workspace to avoid double allocation
cusparseMatDescr_t descr_M_; //this is NOT sparse matrix descriptor
cusparseMatDescr_t descr_A_; //this is NOT sparse matrix descriptor
LinAlgWorkspaceCUDA* workspace_;// so we can copy cusparse handle
LinAlgWorkspaceCUDA* workspace_; ///< Workspace access so we can copy cusparse handle
cusolverSpHandle_t handle_cusolversp_;
cusolverStatus_t status_cusolver_;
cusparseStatus_t status_cusparse_;
Expand Down
14 changes: 7 additions & 7 deletions resolve/LinSolverDirectCuSolverRf.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -40,9 +40,9 @@ namespace ReSolve
error_sum += status_cusolverrf_;
status_cusolverrf_ = cusolverRfSetupDevice(n,
A_->getNnzExpanded(),
A_->getRowData(memory::DEVICE), //dia_,
A_->getColData(memory::DEVICE), //dja_,
A_->getValues( memory::DEVICE), //da_,
A_->getRowData(memory::DEVICE),
A_->getColData(memory::DEVICE),
A_->getValues( memory::DEVICE),
L->getNnz(),
L->getRowData(memory::DEVICE),
L->getColData(memory::DEVICE),
Expand Down Expand Up @@ -81,9 +81,9 @@ namespace ReSolve
int error_sum = 0;
status_cusolverrf_ = cusolverRfResetValues(A_->getNumRows(),
A_->getNnzExpanded(),
A_->getRowData(memory::DEVICE), //dia_,
A_->getColData(memory::DEVICE), //dja_,
A_->getValues( memory::DEVICE), //da_,
A_->getRowData(memory::DEVICE),
A_->getColData(memory::DEVICE),
A_->getValues( memory::DEVICE),
d_P_,
d_Q_,
handle_cusolverrf_);
Expand Down Expand Up @@ -128,6 +128,6 @@ namespace ReSolve
int LinSolverDirectCuSolverRf::setNumericalProperties(double nzero, double nboost)
{
status_cusolverrf_ = cusolverRfSetNumericProperties(handle_cusolverrf_, nzero, nboost);
return status_cusolverrf_;
return status_cusolverrf_;
}
}// namespace resolve
10 changes: 5 additions & 5 deletions resolve/LinSolverDirectCuSolverRf.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -31,13 +31,13 @@ namespace ReSolve
matrix::Sparse* U,
index_type* P,
index_type* Q,
vector_type* rhs = nullptr);
vector_type* rhs = nullptr) override;

int refactorize() override;
int solve(vector_type* rhs, vector_type* x) override;
int solve(vector_type* rhs) override; // rhs overwritten by solution

void setAlgorithms(cusolverRfFactorization_t fact_alg, cusolverRfTriangularSolve_t solve_alg);

int refactorize();
int solve(vector_type* rhs, vector_type* x);
int solve(vector_type* rhs);// the solutuon is returned IN RHS (rhs is overwritten)
int setNumericalProperties(double nzero, double nboost);//these two NEED TO BE DOUBLE
private:
cusolverRfHandle_t handle_cusolverrf_;
Expand Down
4 changes: 3 additions & 1 deletion resolve/LinSolverDirectCuSparseILU0.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -234,6 +234,7 @@ namespace ReSolve

return error_sum;
}

// solution is returned in RHS
int LinSolverDirectCuSparseILU0::solve(vector_type* rhs)
{
Expand Down Expand Up @@ -309,4 +310,5 @@ namespace ReSolve

return error_sum;
}
}// namespace resolve

} // namespace resolve
28 changes: 14 additions & 14 deletions resolve/LinSolverDirectCuSparseILU0.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -34,40 +34,40 @@ namespace ReSolve
matrix::Sparse* U = nullptr,
index_type* P = nullptr,
index_type* Q = nullptr,
vector_type* rhs = nullptr);
vector_type* rhs = nullptr) override;
// if values of A change, but the nnz pattern does not, redo the analysis only (reuse buffers though)
int reset(matrix::Sparse* A);

int solve(vector_type* rhs, vector_type* x);
int solve(vector_type* rhs);// the solutuon is returned IN RHS (rhs is overwritten)
int solve(vector_type* rhs, vector_type* x) override;
int solve(vector_type* rhs) override;


private:
cusparseStatus_t status_cusparse_;

MemoryHandler mem_; ///< Device memory manager object
LinAlgWorkspaceCUDA* workspace_;
LinAlgWorkspaceCUDA* workspace_{nullptr};

cusparseMatDescr_t descr_A_{nullptr};

cusparseSpMatDescr_t mat_L_;
cusparseSpMatDescr_t mat_U_;
cusparseSpMatDescr_t mat_L_{nullptr};
cusparseSpMatDescr_t mat_U_{nullptr};

cusparseSpSVDescr_t descr_spsv_L_{nullptr};
cusparseSpSVDescr_t descr_spsv_U_{nullptr};
csrilu02Info_t info_A_{nullptr};

void* buffer_;
void* buffer_L_;
void* buffer_U_;
void* buffer_{nullptr};
void* buffer_L_{nullptr};
void* buffer_U_{nullptr};

real_type* d_aux1_;
real_type* d_aux2_;
real_type* d_aux1_{nullptr};
real_type* d_aux2_{nullptr};

cusparseDnVecDescr_t vec_X_;
cusparseDnVecDescr_t vec_Y_;
cusparseDnVecDescr_t vec_X_{nullptr};
cusparseDnVecDescr_t vec_Y_{nullptr};

// since ILU OVERWRITES THE MATRIX values, we need a buffer to keep the values of ILU decomposition.
real_type* d_ILU_vals_;
real_type* d_ILU_vals_{nullptr};
};
}// namespace
Loading