19#ifndef OPM_GPUDILU_HPP
20#define OPM_GPUDILU_HPP
24#include <opm/grid/utility/SparseTable.hpp>
25#include <opm/simulators/linalg/PreconditionerWithUpdate.hpp>
26#include <opm/simulators/linalg/gpuistl/GpuBuffer.hpp>
27#include <opm/simulators/linalg/gpuistl/GpuSparseMatrixWrapper.hpp>
28#include <opm/simulators/linalg/gpuistl/GpuVector.hpp>
29#include <opm/simulators/linalg/gpuistl/detail/kernel_enums.hpp>
30#include <opm/simulators/linalg/gpuistl/detail/preconditionerKernels/ILU_variants_helper_kernels.hpp>
31#include <opm/simulators/linalg/gpuistl/gpu_resources.hpp>
51template <
class CPUMatrixT,
class X,
class Y,
int l = 1>
82 void pre(X& x, Y&
b)
override;
85 void apply(X&
v,
const Y&
d)
override;
89 void post(X& x)
override;
92 Dune::SolverCategory::Category
category()
const override;
119 virtual bool hasPerfectUpdate()
const override
131 static constexpr const size_t blocksize_ = CPUMatrixT::block_type::cols;
135 std::vector<int> m_reorderedToNatural;
137 std::vector<int> m_naturalToReordered;
139 const GPUMatrix& m_gpuMatrix;
141 std::unique_ptr<GPUMatrix> m_gpuMatrixReordered;
143 std::unique_ptr<GPUMatrix> m_gpuMatrixReorderedLower;
144 std::unique_ptr<GPUMatrix> m_gpuMatrixReorderedUpper;
146 std::unique_ptr<GpuVector<field_type>> m_gpuMatrixReorderedDiag;
148 std::unique_ptr<FloatMat> m_gpuMatrixReorderedLowerFloat;
149 std::unique_ptr<FloatMat> m_gpuMatrixReorderedUpperFloat;
150 std::unique_ptr<FloatVec> m_gpuMatrixReorderedDiagFloat;
151 std::unique_ptr<FloatVec> m_gpuDInvFloat;
153 GpuVector<int> m_gpuNaturalToReorder;
155 GpuVector<int> m_gpuReorderToNatural;
159 GpuVector<field_type> m_gpuDInv;
163 bool m_tuneThreadBlockSizes;
165 const MatrixStorageMPScheme m_mixedPrecisionScheme;
168 int m_upperSolveThreadBlockSize = -1;
169 int m_lowerSolveThreadBlockSize = -1;
170 int m_moveThreadBlockSize = -1;
171 int m_DILUFactorizationThreadBlockSize = -1;
179 std::map<std::pair<field_type*, const field_type*>,
GPUGraph> m_apply_graphs;
180 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
DILU preconditioner on the GPU.
Definition GpuDILU.hpp:53
void apply(X &v, const Y &d) override
Apply the preconditoner.
Definition GpuDILU.cpp:145
Y range_type
The range type of the preconditioner.
Definition GpuDILU.hpp:58
void computeDiagonal(int factorizationThreadBlockSize)
Compute the diagonal of the DILU, and update the data of the reordered matrix.
Definition GpuDILU.cpp:401
static constexpr bool shouldCallPre()
Definition GpuDILU.hpp:108
static constexpr bool shouldCallPost()
Definition GpuDILU.hpp:114
void pre(X &x, Y &b) override
Prepare the preconditioner.
Definition GpuDILU.cpp:139
void post(X &x) override
Post processing.
Definition GpuDILU.cpp:329
X domain_type
The domain type of the preconditioner.
Definition GpuDILU.hpp:56
Dune::SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
Definition GpuDILU.cpp:335
void reorderAndSplitMatrix(int moveThreadBlockSize)
perform matrix splitting and reordering
Definition GpuDILU.cpp:373
void tuneThreadBlockSizes()
function that will experimentally tune the thread block sizes of the important cuda kernels
Definition GpuDILU.cpp:504
typename X::field_type field_type
The field type of the preconditioner.
Definition GpuDILU.hpp:60
void update() final
Updates the matrix data.
Definition GpuDILU.cpp:342
Definition gpu_type_detection.hpp:30
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240