My Project
amgclSolverBackend.hpp
1 /*
2  Copyright 2020 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_AMGCLSOLVER_BACKEND_HEADER_INCLUDED
21 #define OPM_AMGCLSOLVER_BACKEND_HEADER_INCLUDED
22 
23 #include <mutex>
24 
25 #include <opm/simulators/linalg/bda/BdaResult.hpp>
26 #include <opm/simulators/linalg/bda/BdaSolver.hpp>
27 #include <opm/simulators/linalg/bda/WellContributions.hpp>
28 
29 #include <boost/property_tree/ptree.hpp>
30 
31 #include <amgcl/amg.hpp>
32 #include <amgcl/backend/builtin.hpp>
33 #include <amgcl/adapter/crs_tuple.hpp>
34 #include <amgcl/make_block_solver.hpp>
35 #include <amgcl/relaxation/as_preconditioner.hpp>
36 #include <amgcl/relaxation/ilu0.hpp>
37 #include <amgcl/solver/bicgstab.hpp>
38 #include <amgcl/preconditioner/runtime.hpp>
39 #include <amgcl/value_type/static_matrix.hpp>
40 
41 namespace Opm
42 {
43 namespace Accelerator
44 {
45 
48 template <unsigned int block_size>
49 class amgclSolverBackend : public BdaSolver<block_size>
50 {
52 
53  using Base::N;
54  using Base::Nb;
55  using Base::nnz;
56  using Base::nnzb;
57  using Base::verbosity;
58  using Base::platformID;
59  using Base::deviceID;
60  using Base::maxit;
61  using Base::tolerance;
62  using Base::initialized;
63 
64  typedef amgcl::static_matrix<double, block_size, block_size> dmat_type; // matrix value type in double precision
65  typedef amgcl::static_matrix<double, block_size, 1> dvec_type; // the corresponding vector value type
66  typedef amgcl::backend::builtin<dmat_type> CPU_Backend;
67 
68  typedef amgcl::make_solver<amgcl::runtime::preconditioner<CPU_Backend>, amgcl::runtime::solver::wrapper<CPU_Backend> > CPU_Solver;
69 
70 private:
71 
72  // amgcl can use different backends, this lets the user choose
73  enum Amgcl_backend_type {
74  cpu,
75  cuda,
76  vexcl
77  };
78 
79  // store matrix in CSR format
80  std::vector<unsigned> A_rows, A_cols;
81  std::vector<double> A_vals, rhs;
82  std::vector<double> x;
83  std::once_flag print_info;
84  Amgcl_backend_type backend_type = cpu;
85 
86  boost::property_tree::ptree prm; // amgcl parameters
87  int iters = 0;
88  double error = 0.0;
89 
90 #if HAVE_CUDA
91  std::once_flag cuda_initialize;
92  void solve_cuda(double *b);
93 #endif
94 
95 #if HAVE_VEXCL
96  std::once_flag vexcl_initialize;
97 #endif
101  void initialize(int Nb, int nnzbs);
102 
106  void convert_sparsity_pattern(int *rows, int *cols);
107 
111  void convert_data(double *vals, int *rows);
112 
116  void solve_system(double *b, BdaResult &res);
117 
118 public:
125  amgclSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, unsigned int platformID, unsigned int deviceID);
126 
129 
137  SolverStatus solve_system(std::shared_ptr<BlockedMatrix> matrix, double *b,
138  std::shared_ptr<BlockedMatrix> jacMatrix, WellContributions& wellContribs, BdaResult &res) override;
139 
142  void get_result(double *x) override;
143 
144 }; // end class amgclSolverBackend
145 
146 } // namespace Accelerator
147 } // namespace Opm
148 
149 #endif
150 
151 
This class is based on InverseOperatorResult struct from dune/istl/solver.hh It is needed to prevent ...
Definition: BdaResult.hpp:31
This class serves to simplify choosing between different backend solvers, such as cusparseSolver and ...
Definition: BdaSolver.hpp:45
This class does not implement a solver, but converts the BCSR format to normal CSR and uses amgcl for...
Definition: amgclSolverBackend.hpp:50
~amgclSolverBackend()
Destroy a openclSolver, and free memory.
Definition: amgclSolverBackend.cpp:50
SolverStatus solve_system(std::shared_ptr< BlockedMatrix > matrix, double *b, std::shared_ptr< BlockedMatrix > jacMatrix, WellContributions &wellContribs, BdaResult &res) override
Solve linear system, A*x = b, matrix A must be in blocked-CSR format.
void get_result(double *x) override
Get result after linear solve, and peform postprocessing if necessary.
Definition: amgclSolverBackend.cpp:348
amgclSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, unsigned int platformID, unsigned int deviceID)
Construct a openclSolver.
Definition: amgclSolverBackend.cpp:46
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition: WellContributions.hpp:53
Definition: PropertyTree.hpp:29
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27