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 bda
42 {
43 
46 template <unsigned int block_size>
47 class amgclSolverBackend : public BdaSolver<block_size>
48 {
50 
51  using Base::N;
52  using Base::Nb;
53  using Base::nnz;
54  using Base::nnzb;
55  using Base::verbosity;
56  using Base::platformID;
57  using Base::deviceID;
58  using Base::maxit;
59  using Base::tolerance;
60  using Base::initialized;
61 
62  typedef amgcl::static_matrix<double, block_size, block_size> dmat_type; // matrix value type in double precision
63  typedef amgcl::static_matrix<double, block_size, 1> dvec_type; // the corresponding vector value type
64  typedef amgcl::backend::builtin<dmat_type> CPU_Backend;
65 
66  typedef amgcl::make_solver<amgcl::runtime::preconditioner<CPU_Backend>, amgcl::runtime::solver::wrapper<CPU_Backend> > CPU_Solver;
67 
68 private:
69 
70  // amgcl can use different backends, this lets the user choose
71  enum Amgcl_backend_type {
72  cpu,
73  cuda,
74  vexcl
75  };
76 
77  // store matrix in CSR format
78  std::vector<unsigned> A_rows, A_cols;
79  std::vector<double> A_vals, rhs;
80  std::vector<double> x;
81  std::once_flag print_info;
82  Amgcl_backend_type backend_type = cpu;
83 
84  boost::property_tree::ptree prm; // amgcl parameters
85  int iters = 0;
86  double error = 0.0;
87 
88 #if HAVE_CUDA
89  std::once_flag cuda_initialize;
90  void solve_cuda(double *b);
91 #endif
92 
93 #if HAVE_VEXCL
94  std::once_flag vexcl_initialize;
95 #endif
100  void initialize(int N, int nnz, int dim);
101 
105  void convert_sparsity_pattern(int *rows, int *cols);
106 
110  void convert_data(double *vals, int *rows);
111 
115  void solve_system(double *b, BdaResult &res);
116 
117 public:
124  amgclSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, unsigned int platformID, unsigned int deviceID);
125 
128 
144  SolverStatus solve_system(int N, int nnz, int dim, double *vals, int *rows, int *cols, double *b, WellContributions& wellContribs, BdaResult &res) override;
145 
148  void get_result(double *x) override;
149 
150 }; // end class amgclSolverBackend
151 
152 } // namespace bda
153 
154 #endif
155 
156 
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition: WellContributions.hpp:61
This class is based on InverseOperatorResult struct from dune/istl/solver.hh It is needed to prevent ...
Definition: BdaResult.hpp:29
This class serves to simplify choosing between different backend solvers, such as cusparseSolver and ...
Definition: BdaSolver.hpp:43
This class does not implement a solver, but converts the BCSR format to normal CSR and uses amgcl for...
Definition: amgclSolverBackend.hpp:48
~amgclSolverBackend()
Destroy a openclSolver, and free memory.
Definition: amgclSolverBackend.cpp:49
amgclSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, unsigned int platformID, unsigned int deviceID)
Construct a openclSolver.
Definition: amgclSolverBackend.cpp:45
void get_result(double *x) override
Get result after linear solve, and peform postprocessing if necessary.
Definition: amgclSolverBackend.cpp:347
Definition: PropertyTree.hpp:29