opm-simulators/opm/simulators/linalg/bda/rocm/rocsparseBILU0.hpp

122 lines
3.9 KiB
C++
Raw Normal View History

/*
Copyright 2024 Equinor ASA
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#ifndef OPM_ROCSPARSEBILU0_HPP
#define OPM_ROCSPARSEBILU0_HPP
#include <opm/simulators/linalg/bda/BlockedMatrix.hpp>
#include <opm/simulators/linalg/bda/rocm/rocsparsePreconditioner.hpp>
#include <rocblas/rocblas.h>
#include <rocsparse/rocsparse.h>
#include <hip/hip_version.h>
namespace Opm::Accelerator {
/// This class implements a Blocked ILU0 preconditioner
/// The decomposition is done on GPU, using exact decomposition, or ChowPatel decomposition
/// The preconditioner is applied via two exact triangular solves
template <class Scalar, unsigned int block_size>
class rocsparseBILU0 : public rocsparsePreconditioner<Scalar, block_size>
{
typedef rocsparsePreconditioner<Scalar, block_size> Base;
using Base::N;
using Base::Nb;
using Base::nnz;
using Base::nnzb;
using Base::verbosity;
private:
rocsparse_mat_descr descr_M, descr_L, descr_U;
rocsparse_mat_info ilu_info;
#if HIP_VERSION >= 50400000
rocsparse_mat_info spmv_info;
#endif
rocsparse_int *d_Mrows, *d_Mcols;
2024-06-05 08:00:47 -05:00
Scalar *d_Mvals, *d_t;
void *d_buffer; // buffer space, used by rocsparse ilu0 analysis
public:
rocsparseBILU0(int verbosity_);
/// Initialize GPU and allocate memory
/// \param[in] matrix matrix A
/// \param[in] jacMatrix matrix for preconditioner
2024-06-05 08:00:47 -05:00
bool initialize(std::shared_ptr<BlockedMatrix<Scalar>> matrix,
std::shared_ptr<BlockedMatrix<Scalar>> jacMatrix,
rocsparse_int *d_Arows,
2024-08-13 04:22:54 -05:00
rocsparse_int *d_Acols) override;
2024-06-05 08:00:47 -05:00
/// Analysis, extract parallelism if specified
/// \param[in] mat matrix A
2024-08-13 04:22:54 -05:00
bool analyze_matrix(BlockedMatrix<Scalar> *mat) override;
2024-06-05 08:00:47 -05:00
/// Analysis, extract parallelism if specified
/// \param[in] mat matrix A
/// \param[in] jacMat matrix for preconditioner, analyze this as well
bool analyze_matrix(BlockedMatrix<Scalar> *mat,
2024-08-13 04:22:54 -05:00
BlockedMatrix<Scalar> *jacMat) override;
2024-06-05 08:00:47 -05:00
/// ILU decomposition
/// \param[in] mat matrix A to decompose
bool create_preconditioner(BlockedMatrix<Scalar> *mat) override;
2024-06-05 08:00:47 -05:00
/// ILU decomposition
/// \param[in] mat matrix A
/// \param[in] jacMat matrix for preconditioner, decompose this one if used
bool create_preconditioner(BlockedMatrix<Scalar> *mat,
BlockedMatrix<Scalar> *jacMat) override;
/// Apply preconditioner, x = prec(y)
/// via Lz = y
/// and Ux = z
/// \param[in] y Input y vector
/// \param[out] x Output x vector
void apply(Scalar& y,
Scalar& x) override;
#if HAVE_OPENCL
// apply preconditioner, x = prec(y)
2024-08-13 04:22:54 -05:00
void apply(const cl::Buffer&, cl::Buffer&) override {}
#endif
2024-06-05 08:00:47 -05:00
/// Copy matrix A values to GPU
/// \param[in] mVals Input values
2024-08-13 04:22:54 -05:00
void copy_system_to_gpu(Scalar *mVals) override;
2024-06-05 08:00:47 -05:00
/// Reassign pointers, in case the addresses of the Dune variables have changed --> TODO: check when/if we need this method
// /// \param[in] vals New values
// /// \param[in] b New b vector
// void update_system(Scalar *vals, Scalar *b);
/// Update GPU values after a new assembly is done
/// \param[in] b New b vector
2024-08-13 04:22:54 -05:00
void update_system_on_gpu(Scalar *b) override;
};
} // namespace Opm
#endif