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;
127 GpuSparseMatrixWrapper<field_type> m_LU;
129 GpuVector<field_type> m_temporaryStorage;
134 detail::CuSparseResource<bsrsv2Info_t> m_infoL;
135 detail::CuSparseResource<bsrsv2Info_t> m_infoU;
136 detail::CuSparseResource<bsrilu02Info_t> m_infoM;
138 std::unique_ptr<GpuVector<field_type>> m_buffer;
139 detail::CuSparseHandle& m_cuSparseHandle;
141 bool m_analysisDone =
false;
143 void analyzeMatrix();
144 size_t findBufferSize();
148 void updateILUConfiguration();
153template <
class M,
class X,
class Y,
int l>
154template <
typename T,
typename std::enable_if_t<!is_gpu_matrix_v<T>,
int>>
156 : m_underlyingMatrix(A)
159 , m_temporaryStorage(m_LU.N() * m_LU.blockSize())
160 , m_descriptionL(
detail::createLowerDiagonalDescription())
161 , m_descriptionU(
detail::createUpperDiagonalDescription())
162 , m_cuSparseHandle(
detail::CuSparseHandle::getInstance())
165 OPM_ERROR_IF(A.N() != m_LU.N(),
166 fmt::format(
"CuSparse matrix not same size as DUNE matrix. {} vs {}.", m_LU.N(), A.N()));
168 A[0][0].N() != m_LU.blockSize(),
169 fmt::format(
"CuSparse matrix not same blocksize as DUNE matrix. {} vs {}.", m_LU.blockSize(), A[0][0].N()));
171 A.N() * A[0][0].N() != m_LU.dim(),
172 fmt::format(
"CuSparse matrix not same dimension as DUNE matrix. {} vs {}.", m_LU.dim(), A.N() * A[0][0].N()));
173 OPM_ERROR_IF(A.nonzeroes() != m_LU.nonzeroes(),
174 fmt::format(
"CuSparse matrix not same number of non zeroes as DUNE matrix. {} vs {}. ",
179 updateILUConfiguration();
183template <
class M,
class X,
class Y,
int l>
184template <
typename T,
typename std::enable_if_t<is_gpu_matrix_v<T>,
int>>
186 : m_underlyingMatrix(A)
189 , m_temporaryStorage(m_LU.N() * m_LU.blockSize())
190 , m_descriptionL(
detail::createLowerDiagonalDescription())
191 , m_descriptionU(
detail::createUpperDiagonalDescription())
192 , m_cuSparseHandle(
detail::CuSparseHandle::getInstance())
195 OPM_ERROR_IF(A.N() != m_LU.
N(),
196 fmt::format(
"CuSparse matrix not same size as DUNE matrix. {} vs {}.", m_LU.
N(), A.N()));
197 OPM_ERROR_IF(A.nonzeroes() != m_LU.
nonzeroes(),
198 fmt::format(
"CuSparse matrix not same number of non zeroes as DUNE matrix. {} vs {}. ",
203 updateILUConfiguration();
Interface class adding the update() method to the preconditioner interface.
Definition PreconditionerWithUpdate.hpp:32
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:125
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:49
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:112
virtual void apply(X &v, const Y &d) override
Apply the preconditoner.
Definition GpuSeqILU0.cpp:55
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:118
GpuSeqILU0(const M &A, field_type w)
Constructor.
Definition GpuSeqILU0.hpp:155
The GpuSparseMatrixWrapper Checks CUDA/HIP version and dispatches a version either using the old or t...
Definition GpuSparseMatrixWrapper.hpp:61
size_t nonzeroes() const
nonzeroes behaves as the Dune::BCRSMatrix::nonzeros() function and returns the number of non zero blo...
Definition GpuSparseMatrixWrapper.hpp:237
size_t N() const
N returns the number of rows (which is equal to the number of columns).
Definition GpuSparseMatrixWrapper.hpp:228
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