My Project
Loading...
Searching...
No Matches
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 <opm/simulators/linalg/gpubridge/GpuResult.hpp>
24#include <opm/simulators/linalg/gpubridge/GpuSolver.hpp>
25#include <opm/simulators/linalg/gpubridge/WellContributions.hpp>
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 <type_traits>
42#include <vector>
43
44namespace Opm::Accelerator {
45
48template<class Scalar, unsigned int block_size>
49class amgclSolverBackend : public GpuSolver<Scalar,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 using dmat_type = amgcl::static_matrix<Scalar, block_size, block_size>; // matrix value type in double precision
65 using dvec_type = amgcl::static_matrix<Scalar, block_size, 1>; // the corresponding vector value type
66 using CPU_Backend = std::conditional_t<block_size == 1,
67 amgcl::backend::builtin<Scalar>,
68 amgcl::backend::builtin<dmat_type>>;
69
70 using CPU_Solver = amgcl::make_solver<amgcl::runtime::preconditioner<CPU_Backend>,
71 amgcl::runtime::solver::wrapper<CPU_Backend>>;
72
73private:
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<Scalar> A_vals, rhs;
84 std::vector<Scalar> 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 Scalar error = 0.0;
91
92#if HAVE_CUDA
93 std::once_flag cuda_initialize;
94 void solve_cuda(Scalar* 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(const int *rows, const int *cols);
109
113 void convert_data(const Scalar* vals, const int* rows);
114
118 void solve_system(Scalar* b, GpuResult& res);
119
120public:
128 Scalar tolerance, unsigned int platformID,
129 unsigned int deviceID);
130
132 ~amgclSolverBackend() override;
133
141 SolverStatus solve_system(std::shared_ptr<BlockedMatrix<Scalar>> matrix,
142 Scalar* b,
143 std::shared_ptr<BlockedMatrix<Scalar>> jacMatrix,
145 GpuResult& res) override;
146
149 void get_result(Scalar* x) override;
150
151}; // end class amgclSolverBackend
152
153} // namespace Opm::Accelerator
154
155#endif
This struct resembles a blocked csr matrix, like Dune::BCRSMatrix.
Definition BlockedMatrix.hpp:29
This class is based on InverseOperatorResult struct from dune/istl/solver.hh It is needed to prevent ...
Definition GpuResult.hpp:31
This class serves to simplify choosing between different backend solvers, such as cusparseSolver and ...
Definition GpuSolver.hpp:46
This class does not implement a solver, but converts the BCSR format to normal CSR and uses amgcl for...
Definition amgclSolverBackend.hpp:50
void get_result(Scalar *x) override
Get result after linear solve, and peform postprocessing if necessary.
Definition amgclSolverBackend.cpp:401
~amgclSolverBackend() override
Destroy a openclSolver, and free memory.
Definition amgclSolverBackend.cpp:64
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:51
Definition PropertyTree.hpp:31
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242