CPR.hpp
Go to the documentation of this file.
1/*
2 Copyright 2021 Equinor ASA
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef OPM_CPR_HPP
21#define OPM_CPR_HPP
22
23#include <mutex>
24
25#include <dune/istl/paamg/matrixhierarchy.hh>
26#include <dune/istl/umfpack.hh>
27
33
35
36namespace Opm
37{
38namespace Accelerator
39{
40
41class BlockedMatrix;
42
44template <unsigned int block_size>
45class CPR : public Preconditioner<block_size>
46{
48
49 using Base::N;
50 using Base::Nb;
51 using Base::nnz;
52 using Base::nnzb;
53 using Base::verbosity;
54 using Base::context;
55 using Base::queue;
56 using Base::events;
57 using Base::err;
58
59private:
60 int num_levels;
61 std::vector<double> weights, coarse_vals, coarse_x, coarse_y;
62 std::vector<Matrix> Amatrices, Rmatrices; // scalar matrices that represent the AMG hierarchy
63 std::vector<OpenclMatrix> d_Amatrices, d_Rmatrices; // scalar matrices that represent the AMG hierarchy
64 std::vector<std::vector<int> > PcolIndices; // prolongation does not need a full matrix, only store colIndices
65 std::vector<cl::Buffer> d_PcolIndices;
66 std::vector<std::vector<double> > invDiags; // inverse of diagonal of Amatrices
67 std::vector<cl::Buffer> d_invDiags;
68 std::vector<cl::Buffer> d_t, d_f, d_u; // intermediate vectors used during amg cycle
69 std::unique_ptr<cl::Buffer> d_rs; // use before extracting the pressure
70 std::unique_ptr<cl::Buffer> d_weights; // the quasiimpes weights, used to extract pressure
71 std::unique_ptr<OpenclMatrix> d_mat; // stores blocked matrix
72 std::unique_ptr<cl::Buffer> d_coarse_y, d_coarse_x; // stores the scalar vectors
73 std::once_flag opencl_buffers_allocated; // only allocate OpenCL Buffers once
74
75 std::unique_ptr<BILU0<block_size> > bilu0; // Blocked ILU0 preconditioner
76 BlockedMatrix *mat = nullptr; // input matrix, blocked
77
78 using DuneMat = Dune::BCRSMatrix<Dune::FieldMatrix<double, 1, 1> >;
79 using DuneVec = Dune::BlockVector<Dune::FieldVector<double, 1> >;
80 using MatrixOperator = Dune::MatrixAdapter<DuneMat, DuneVec, DuneVec>;
81 using DuneAmg = Dune::Amg::MatrixHierarchy<MatrixOperator, Dune::Amg::SequentialInformation>;
82 std::unique_ptr<DuneAmg> dune_amg;
83 std::unique_ptr<DuneMat> dune_coarse; // extracted pressure matrix, finest level in AMG hierarchy
84 std::shared_ptr<MatrixOperator> dune_op; // operator, input to Dune AMG
85 std::vector<int> level_sizes; // size of each level in the AMG hierarchy
86 std::vector<std::vector<int> > diagIndices; // index of diagonal value for each level
87 Dune::UMFPack<DuneMat> umfpack; // dune/istl/umfpack object used to solve the coarsest level of AMG
88 bool always_recalculate_aggregates = false; // OPM always reuses the aggregates by default
89 bool recalculate_aggregates = true; // only rerecalculate if true
90 const int pressure_idx = 1; // hardcoded to mimic OPM
91 unsigned num_pre_smooth_steps; // number of Jacobi smooth steps before restriction
92 unsigned num_post_smooth_steps; // number of Jacobi smooth steps after prolongation
93
94 std::unique_ptr<openclSolverBackend<1> > coarse_solver; // coarse solver is scalar
95 bool opencl_ilu_parallel; // whether ILU0 operation should be parallelized
96
97 // Analyze the AMG hierarchy build by Dune
98 void analyzeHierarchy();
99
100 // Analyze the aggregateMaps from the AMG hierarchy
101 // These can be reused, so only use when recalculate_aggregates is true
102 void analyzeAggregateMaps();
103
104 // Initialize and allocate matrices and vectors
105 void init_opencl_buffers();
106
107 // Copy matrices and vectors to GPU
108 void opencl_upload();
109
110 // apply pressure correction to vector
111 void apply_amg(const cl::Buffer& y, cl::Buffer& x);
112
113 void amg_cycle_gpu(const int level, cl::Buffer &y, cl::Buffer &x);
114
115 void create_preconditioner_amg(BlockedMatrix *mat);
116
117public:
118 CPR(bool opencl_ilu_parallel, int verbosity);
119
120 bool analyze_matrix(BlockedMatrix *mat) override;
121 bool analyze_matrix(BlockedMatrix *mat, BlockedMatrix *jacMat) override;
122
123 // set own Opencl variables, but also that of the bilu0 preconditioner
124 void setOpencl(std::shared_ptr<cl::Context>& context, std::shared_ptr<cl::CommandQueue>& queue) override;
125
126 // applies blocked ilu0
127 // also applies amg for pressure component
128 void apply(const cl::Buffer& y, cl::Buffer& x) override;
129
132};
133
134// solve A^T * x = b
135// A should represent a 3x3 matrix
136// x and b are vectors with 3 elements
137void solve_transposed_3x3(const double *A, const double *b, double *x);
138
139} // namespace Accelerator
140} // namespace Opm
141
142#endif
143
Definition: BlockedMatrix.hpp:31
This class implements a Constrained Pressure Residual (CPR) preconditioner.
Definition: CPR.hpp:46
void setOpencl(std::shared_ptr< cl::Context > &context, std::shared_ptr< cl::CommandQueue > &queue) override
void apply(const cl::Buffer &y, cl::Buffer &x) override
bool create_preconditioner(BlockedMatrix *mat, BlockedMatrix *jacMat) override
bool create_preconditioner(BlockedMatrix *mat) override
CPR(bool opencl_ilu_parallel, int verbosity)
bool analyze_matrix(BlockedMatrix *mat) override
bool analyze_matrix(BlockedMatrix *mat, BlockedMatrix *jacMat) override
Definition: Preconditioner.hpp:36
int nnz
Definition: Preconditioner.hpp:41
std::vector< cl::Event > events
Definition: Preconditioner.hpp:47
std::shared_ptr< cl::CommandQueue > queue
Definition: Preconditioner.hpp:46
int N
Definition: Preconditioner.hpp:39
int verbosity
Definition: Preconditioner.hpp:43
cl_int err
Definition: Preconditioner.hpp:48
int Nb
Definition: Preconditioner.hpp:40
std::shared_ptr< cl::Context > context
Definition: Preconditioner.hpp:45
int nnzb
Definition: Preconditioner.hpp:42
void solve_transposed_3x3(const double *A, const double *b, double *x)
Definition: BlackoilPhases.hpp:27