My Project
CPR.hpp
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/common/version.hh>
26 #if DUNE_VERSION_NEWER(DUNE_ISTL, 2, 7)
27 #include <dune/istl/paamg/matrixhierarchy.hh>
28 #else
29 #include <dune/istl/paamg/hierarchy.hh>
30 #endif
31 #include <dune/istl/umfpack.hh>
32 
33 #include <opm/simulators/linalg/bda/opencl/opencl.hpp>
34 #include <opm/simulators/linalg/bda/opencl/BILU0.hpp>
35 #include <opm/simulators/linalg/bda/Matrix.hpp>
36 #include <opm/simulators/linalg/bda/opencl/OpenclMatrix.hpp>
37 #include <opm/simulators/linalg/bda/ILUReorder.hpp>
38 #include <opm/simulators/linalg/bda/opencl/Preconditioner.hpp>
39 
40 #include <opm/simulators/linalg/bda/opencl/openclSolverBackend.hpp>
41 
42 namespace Opm
43 {
44 namespace Accelerator
45 {
46 
47 class BlockedMatrix;
48 
50 template <unsigned int block_size>
51 class CPR : public Preconditioner<block_size>
52 {
54 
55  using Base::N;
56  using Base::Nb;
57  using Base::nnz;
58  using Base::nnzb;
59  using Base::verbosity;
60  using Base::context;
61  using Base::queue;
62  using Base::events;
63  using Base::err;
64 
65 private:
66  int num_levels;
67  std::vector<double> weights, coarse_vals, coarse_x, coarse_y;
68  std::vector<Matrix> Amatrices, Rmatrices; // scalar matrices that represent the AMG hierarchy
69  std::vector<OpenclMatrix> d_Amatrices, d_Rmatrices; // scalar matrices that represent the AMG hierarchy
70  std::vector<std::vector<int> > PcolIndices; // prolongation does not need a full matrix, only store colIndices
71  std::vector<cl::Buffer> d_PcolIndices;
72  std::vector<std::vector<double> > invDiags; // inverse of diagonal of Amatrices
73  std::vector<cl::Buffer> d_invDiags;
74  std::vector<cl::Buffer> d_t, d_f, d_u; // intermediate vectors used during amg cycle
75  std::unique_ptr<cl::Buffer> d_rs; // use before extracting the pressure
76  std::unique_ptr<cl::Buffer> d_weights; // the quasiimpes weights, used to extract pressure
77  std::unique_ptr<OpenclMatrix> d_mat; // stores blocked matrix
78  std::unique_ptr<cl::Buffer> d_coarse_y, d_coarse_x; // stores the scalar vectors
79  std::once_flag opencl_buffers_allocated; // only allocate OpenCL Buffers once
80 
81  std::unique_ptr<BILU0<block_size> > bilu0; // Blocked ILU0 preconditioner
82  BlockedMatrix *mat = nullptr; // input matrix, blocked
83 
84  using DuneMat = Dune::BCRSMatrix<Dune::FieldMatrix<double, 1, 1> >;
85  using DuneVec = Dune::BlockVector<Dune::FieldVector<double, 1> >;
86  using MatrixOperator = Dune::MatrixAdapter<DuneMat, DuneVec, DuneVec>;
87  using DuneAmg = Dune::Amg::MatrixHierarchy<MatrixOperator, Dune::Amg::SequentialInformation>;
88  std::unique_ptr<DuneAmg> dune_amg;
89  std::unique_ptr<DuneMat> dune_coarse; // extracted pressure matrix, finest level in AMG hierarchy
90  std::shared_ptr<MatrixOperator> dune_op; // operator, input to Dune AMG
91  std::vector<int> level_sizes; // size of each level in the AMG hierarchy
92  std::vector<std::vector<int> > diagIndices; // index of diagonal value for each level
93  Dune::UMFPack<DuneMat> umfpack; // dune/istl/umfpack object used to solve the coarsest level of AMG
94  bool always_recalculate_aggregates = false; // OPM always reuses the aggregates by default
95  bool recalculate_aggregates = true; // only rerecalculate if true
96  const int pressure_idx = 1; // hardcoded to mimic OPM
97  unsigned num_pre_smooth_steps; // number of Jacobi smooth steps before restriction
98  unsigned num_post_smooth_steps; // number of Jacobi smooth steps after prolongation
99 
100  std::unique_ptr<openclSolverBackend<1> > coarse_solver; // coarse solver is scalar
101  ILUReorder opencl_ilu_reorder; // reordering strategy for ILU0 in coarse solver
102 
103  // Analyze the AMG hierarchy build by Dune
104  void analyzeHierarchy();
105 
106  // Analyze the aggregateMaps from the AMG hierarchy
107  // These can be reused, so only use when recalculate_aggregates is true
108  void analyzeAggregateMaps();
109 
110  // Initialize and allocate matrices and vectors
111  void init_opencl_buffers();
112 
113  // Copy matrices and vectors to GPU
114  void opencl_upload();
115 
116  // apply pressure correction to vector
117  void apply_amg(const cl::Buffer& y, cl::Buffer& x);
118 
119  void amg_cycle_gpu(const int level, cl::Buffer &y, cl::Buffer &x);
120 
121  void create_preconditioner_amg(BlockedMatrix *mat);
122 
123 public:
124 
125  CPR(int verbosity, ILUReorder opencl_ilu_reorder);
126 
127  bool analyze_matrix(BlockedMatrix *mat) override;
128 
129  // set own Opencl variables, but also that of the bilu0 preconditioner
130  void setOpencl(std::shared_ptr<cl::Context>& context, std::shared_ptr<cl::CommandQueue>& queue) override;
131 
132  // applies blocked ilu0
133  // also applies amg for pressure component
134  void apply(const cl::Buffer& y, cl::Buffer& x) override;
135 
136  bool create_preconditioner(BlockedMatrix *mat) override;
137 
138  int* getToOrder() override
139  {
140  return bilu0->getToOrder();
141  }
142 
143  int* getFromOrder() override
144  {
145  return bilu0->getFromOrder();
146  }
147 
148  BlockedMatrix* getRMat() override
149  {
150  return bilu0->getRMat();
151  }
152 };
153 
154 // solve A^T * x = b
155 // A should represent a 3x3 matrix
156 // x and b are vectors with 3 elements
157 void solve_transposed_3x3(const double *A, const double *b, double *x);
158 
159 } // namespace Accelerator
160 } // namespace Opm
161 
162 #endif
163 
This struct resembles a blocked csr matrix, like Dune::BCRSMatrix.
Definition: BlockedMatrix.hpp:37
This class implements a Constrained Pressure Residual (CPR) preconditioner.
Definition: CPR.hpp:52
Definition: Preconditioner.hpp:36
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27