19#ifndef OPM_GPUILU0_OPM_Impl_HPP
20#define OPM_GPUILU0_OPM_Impl_HPP
23#include <opm/grid/utility/SparseTable.hpp>
24#include <opm/simulators/linalg/PreconditionerWithUpdate.hpp>
25#include <opm/simulators/linalg/gpuistl/GpuSparseMatrixWrapper.hpp>
26#include <opm/simulators/linalg/gpuistl/GpuVector.hpp>
27#include <opm/simulators/linalg/gpuistl/gpu_resources.hpp>
28#include <opm/simulators/linalg/gpuistl/detail/kernel_enums.hpp>
49template <
class CPUMatrixT,
class X,
class Y,
int l = 1>
85 void pre(X& x, Y&
b)
override;
88 void apply(X&
v,
const Y&
d)
override;
92 void post(X& x)
override;
95 Dune::SolverCategory::Category
category()
const override;
121 virtual bool hasPerfectUpdate()
const override {
133 static constexpr const size_t blocksize_ = CPUMatrixT::block_type::cols;
137 std::vector<int> m_reorderedToNatural;
139 std::vector<int> m_naturalToReordered;
142 std::unique_ptr<GpuMatrix> m_gpuReorderedLU;
144 std::unique_ptr<GpuMatrix> m_gpuMatrixReorderedLower;
145 std::unique_ptr<GpuMatrix> m_gpuMatrixReorderedUpper;
147 std::unique_ptr<FloatMat> m_gpuMatrixReorderedLowerFloat;
148 std::unique_ptr<FloatMat> m_gpuMatrixReorderedUpperFloat;
149 std::optional<GpuVector<float>> m_gpuMatrixReorderedDiagFloat;
151 std::optional<GpuVector<field_type>> m_gpuMatrixReorderedDiag;
153 GpuVector<int> m_gpuNaturalToReorder;
155 GpuVector<int> m_gpuReorderToNatural;
157 GpuVector<field_type> m_gpuDInv;
161 bool m_tuneThreadBlockSizes;
163 const MatrixStorageMPScheme m_mixedPrecisionScheme;
166 int m_upperSolveThreadBlockSize = -1;
167 int m_lowerSolveThreadBlockSize = -1;
168 int m_moveThreadBlockSize = -1;
169 int m_ILU0FactorizationThreadBlockSize = -1;
172 std::map<std::pair<field_type*, const field_type*>,
GPUGraph> m_apply_graphs;
173 std::map<std::pair<field_type*, const field_type*>,
GPUGraphExec> m_executableGraphs;
Interface class adding the update() method to the preconditioner interface.
Definition PreconditionerWithUpdate.hpp:32
Definition BlackoilWellModel.hpp:91
ILU0 preconditioner on the GPU.
Definition OpmGpuILU0.hpp:51
static constexpr bool shouldCallPre()
Definition OpmGpuILU0.hpp:110
Y range_type
The range type of the preconditioner.
Definition OpmGpuILU0.hpp:56
void LUFactorizeMatrix(int factorizationThreadBlockSize)
Compute LU factorization, and update the data of the reordered matrix.
Definition OpmGpuILU0.cpp:326
void pre(X &x, Y &b) override
Prepare the preconditioner.
Definition OpmGpuILU0.cpp:110
void update() final
Updates the matrix data.
Definition OpmGpuILU0.cpp:271
void apply(X &v, const Y &d) override
Apply the preconditoner.
Definition OpmGpuILU0.cpp:116
void reorderAndSplitMatrix(int moveThreadBlockSize)
perform matrix splitting and reordering
Definition OpmGpuILU0.cpp:298
typename X::field_type field_type
The field type of the preconditioner.
Definition OpmGpuILU0.hpp:58
GpuSparseMatrixWrapper< field_type > GpuMatrix
The GPU matrix type.
Definition OpmGpuILU0.hpp:60
Dune::SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
Definition OpmGpuILU0.cpp:264
void post(X &x) override
Post processing.
Definition OpmGpuILU0.cpp:258
void tuneThreadBlockSizes()
function that will experimentally tune the thread block sizes of the important cuda kernels
Definition OpmGpuILU0.cpp:404
static constexpr bool shouldCallPost()
Definition OpmGpuILU0.hpp:116
X domain_type
The domain type of the preconditioner.
Definition OpmGpuILU0.hpp:54
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240