amgclSolverBackend.hpp
Go to the documentation of this file.
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
26
27#include <boost/property_tree/ptree.hpp>
28
29#include <amgcl/amg.hpp>
30#include <amgcl/backend/builtin.hpp>
31#include <amgcl/adapter/crs_tuple.hpp>
32#include <amgcl/make_block_solver.hpp>
33#include <amgcl/relaxation/as_preconditioner.hpp>
34#include <amgcl/relaxation/ilu0.hpp>
35#include <amgcl/solver/bicgstab.hpp>
36#include <amgcl/preconditioner/runtime.hpp>
37#include <amgcl/value_type/static_matrix.hpp>
38
39#include <memory>
40#include <mutex>
41#include <vector>
42
43namespace Opm
44{
45namespace Accelerator
46{
47
50template <unsigned int block_size>
51class amgclSolverBackend : public BdaSolver<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::platformID;
61 using Base::deviceID;
62 using Base::maxit;
63 using Base::tolerance;
65
66 typedef amgcl::static_matrix<double, block_size, block_size> dmat_type; // matrix value type in double precision
67 typedef amgcl::static_matrix<double, block_size, 1> dvec_type; // the corresponding vector value type
68 typedef amgcl::backend::builtin<dmat_type> CPU_Backend;
69
70 typedef amgcl::make_solver<amgcl::runtime::preconditioner<CPU_Backend>, amgcl::runtime::solver::wrapper<CPU_Backend> > CPU_Solver;
71
72private:
73
74 // amgcl can use different backends, this lets the user choose
75 enum Amgcl_backend_type {
76 cpu,
77 cuda,
78 vexcl
79 };
80
81 // store matrix in CSR format
82 std::vector<unsigned> A_rows, A_cols;
83 std::vector<double> A_vals, rhs;
84 std::vector<double> x;
85 std::once_flag print_info;
86 Amgcl_backend_type backend_type = cpu;
87
88 boost::property_tree::ptree prm; // amgcl parameters
89 int iters = 0;
90 double error = 0.0;
91
92#if HAVE_CUDA
93 std::once_flag cuda_initialize;
94 void solve_cuda(double *b);
95#endif
96
97#if HAVE_VEXCL
98 std::once_flag vexcl_initialize;
99#endif
103 void initialize(int Nb, int nnzbs);
104
108 void convert_sparsity_pattern(int *rows, int *cols);
109
113 void convert_data(double *vals, int *rows);
114
118 void solve_system(double *b, BdaResult &res);
119
120public:
127 amgclSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, unsigned int platformID, unsigned int deviceID);
128
131
139 SolverStatus solve_system(std::shared_ptr<BlockedMatrix> matrix, double *b,
140 std::shared_ptr<BlockedMatrix> jacMatrix, WellContributions& wellContribs, BdaResult &res) override;
141
144 void get_result(double *x) override;
145
146}; // end class amgclSolverBackend
147
148} // namespace Accelerator
149} // namespace Opm
150
151#endif
152
153
Definition: BdaResult.hpp:31
Definition: BdaSolver.hpp:46
int nnz
Definition: BdaSolver.hpp:63
int verbosity
Definition: BdaSolver.hpp:56
int maxit
Definition: BdaSolver.hpp:58
double tolerance
Definition: BdaSolver.hpp:59
int N
Definition: BdaSolver.hpp:61
int Nb
Definition: BdaSolver.hpp:62
int nnzb
Definition: BdaSolver.hpp:64
bool initialized
Definition: BdaSolver.hpp:69
unsigned int deviceID
Definition: BdaSolver.hpp:67
unsigned int platformID
Definition: BdaSolver.hpp:66
Definition: amgclSolverBackend.hpp:52
~amgclSolverBackend()
Destroy a openclSolver, and free memory.
SolverStatus solve_system(std::shared_ptr< BlockedMatrix > matrix, double *b, std::shared_ptr< BlockedMatrix > jacMatrix, WellContributions &wellContribs, BdaResult &res) override
void get_result(double *x) override
amgclSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, unsigned int platformID, unsigned int deviceID)
Definition: WellContributions.hpp:52
Definition: PropertyTree.hpp:29
SolverStatus
Definition: BdaSolver.hpp:35
Definition: BlackoilPhases.hpp:27