19#ifndef OPM_GPUSEQILU0_HPP
20#define OPM_GPUSEQILU0_HPP
22#include <dune/istl/preconditioner.hh>
23#include <opm/simulators/linalg/PreconditionerWithUpdate.hpp>
24#include <opm/simulators/linalg/gpuistl/GpuSparseMatrixWrapper.hpp>
25#include <opm/simulators/linalg/gpuistl/detail/CuMatrixDescription.hpp>
26#include <opm/simulators/linalg/gpuistl/detail/CuSparseHandle.hpp>
27#include <opm/simulators/linalg/gpuistl/detail/CuSparseResource.hpp>
28#include <opm/simulators/linalg/gpuistl/detail/fix_zero_diagonal.hpp>
29#include <opm/simulators/linalg/is_gpu_operator.hpp>
50template <
class M,
class X,
class Y,
int l = 1>
69 template <
typename T = M,
typename std::enable_if_t<is_gpu_matrix_v<T>,
int> = 0>
79 template <
typename T = M,
typename std::enable_if_t<!is_gpu_matrix_v<T>,
int> = 0>
85 virtual void pre(X& x, Y& b)
override;
88 virtual void apply(X& v,
const Y& d)
override;
92 virtual void post(X& x)
override;
95 virtual Dune::SolverCategory::Category
category()
const override;
98 virtual void update()
override;
113 virtual bool hasPerfectUpdate()
const override {
120 const M& m_underlyingMatrix;
129 GpuSparseMatrixWrapper<field_type, true> m_LU;
131 GpuVector<field_type> m_temporaryStorage;
136 detail::CuSparseResource<bsrsv2Info_t> m_infoL;
137 detail::CuSparseResource<bsrsv2Info_t> m_infoU;
138 detail::CuSparseResource<bsrilu02Info_t> m_infoM;
140 std::unique_ptr<GpuVector<field_type>> m_buffer;
141 detail::CuSparseHandle& m_cuSparseHandle;
143 bool m_analysisDone =
false;
145 void analyzeMatrix();
146 size_t findBufferSize();
150 void updateILUConfiguration();
155template <
class M,
class X,
class Y,
int l>
156template <
typename T,
typename std::enable_if_t<!is_gpu_matrix_v<T>,
int>>
158 : m_underlyingMatrix(A)
161 , m_temporaryStorage(m_LU.N() * m_LU.blockSize())
162 , m_descriptionL(
detail::createLowerDiagonalDescription())
163 , m_descriptionU(
detail::createUpperDiagonalDescription())
164 , m_cuSparseHandle(
detail::CuSparseHandle::getInstance())
167 OPM_ERROR_IF(A.N() != m_LU.N(),
168 fmt::format(
"CuSparse matrix not same size as DUNE matrix. {} vs {}.", m_LU.N(), A.N()));
170 A[0][0].N() != m_LU.blockSize(),
171 fmt::format(
"CuSparse matrix not same blocksize as DUNE matrix. {} vs {}.", m_LU.blockSize(), A[0][0].N()));
173 A.N() * A[0][0].N() != m_LU.dim(),
174 fmt::format(
"CuSparse matrix not same dimension as DUNE matrix. {} vs {}.", m_LU.dim(), A.N() * A[0][0].N()));
175 OPM_ERROR_IF(A.nonzeroes() != m_LU.nonzeroes(),
176 fmt::format(
"CuSparse matrix not same number of non zeroes as DUNE matrix. {} vs {}. ",
181 updateILUConfiguration();
185template <
class M,
class X,
class Y,
int l>
186template <
typename T,
typename std::enable_if_t<is_gpu_matrix_v<T>,
int>>
188 : m_underlyingMatrix(A)
190 , m_LU(A.getNonZeroValues().data(),
191 A.getRowIndices().data(),
192 A.getColumnIndices().data(),
196 , m_temporaryStorage(m_LU.N() * m_LU.blockSize())
197 , m_descriptionL(
detail::createLowerDiagonalDescription())
198 , m_descriptionU(
detail::createUpperDiagonalDescription())
199 , m_cuSparseHandle(
detail::CuSparseHandle::getInstance())
202 OPM_ERROR_IF(A.N() != m_LU.
N(),
203 fmt::format(
"CuSparse matrix not same size as DUNE matrix. {} vs {}.", m_LU.
N(), A.N()));
204 OPM_ERROR_IF(A.nonzeroes() != m_LU.
nonzeroes(),
205 fmt::format(
"CuSparse matrix not same number of non zeroes as DUNE matrix. {} vs {}. ",
210 updateILUConfiguration();
Interface class adding the update() method to the preconditioner interface.
Definition PreconditionerWithUpdate.hpp:34
typename std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition GpuSeqILU0.hpp:55
virtual void update() override
Updates the matrix data.
Definition GpuSeqILU0.cpp:130
Y range_type
The range type of the preconditioner.
Definition GpuSeqILU0.hpp:59
virtual void pre(X &x, Y &b) override
Prepare the preconditioner.
Definition GpuSeqILU0.cpp:54
typename X::field_type field_type
The field type of the preconditioner.
Definition GpuSeqILU0.hpp:61
virtual void post(X &x) override
Post processing.
Definition GpuSeqILU0.cpp:117
virtual void apply(X &v, const Y &d) override
Apply the preconditoner.
Definition GpuSeqILU0.cpp:60
X domain_type
The domain type of the preconditioner.
Definition GpuSeqILU0.hpp:57
static constexpr bool shouldCallPost()
Definition GpuSeqILU0.hpp:108
static constexpr bool shouldCallPre()
Definition GpuSeqILU0.hpp:102
virtual Dune::SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category).
Definition GpuSeqILU0.cpp:123
GpuSeqILU0(const M &A, field_type w)
Constructor.
Definition GpuSeqILU0.hpp:157
The GpuSparseMatrixWrapper Checks CUDA/HIP version and dispatches a version either using the old or t...
Definition GpuSparseMatrixWrapper.hpp:62
std::size_t N() const
N returns the number of rows (which is equal to the number of columns).
Definition GpuSparseMatrixWrapper.hpp:232
std::size_t nonzeroes() const
nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blo...
Definition GpuSparseMatrixWrapper.hpp:241
Contains wrappers to make the CuBLAS library behave as a modern C++ library with function overlading.
Definition autotuner.hpp:30
std::shared_ptr< CuSparseResource< cusparseMatDescr_t > > GpuSparseMatrixDescriptionPtr
Pointer to GpuSparseMatrixDescription holder.
Definition CuMatrixDescription.hpp:35
A small, fixed‑dimension MiniVector class backed by std::array that can be used in both host and CUDA...
Definition AmgxInterface.hpp:38