BdaBridge.hpp
Go to the documentation of this file.
1/*
2 Copyright 2019 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 BDABRIDGE_HEADER_INCLUDED
21#define BDABRIDGE_HEADER_INCLUDED
22
23#include "dune/istl/solver.hh" // for struct InverseOperatorResult
24
26
27namespace Opm
28{
29
30template<class Scalar> class WellContributions;
31
33
35template <class BridgeMatrix, class BridgeVector, int block_size>
37{
38private:
39 using Scalar = typename BridgeVector::field_type;
40 int verbosity = 0;
41 bool use_gpu = false;
42 std::string accelerator_mode;
43 std::unique_ptr<Accelerator::BdaSolver<Scalar,block_size>> backend;
44 std::shared_ptr<Accelerator::BlockedMatrix<Scalar>> matrix; // 'stores' matrix, actually points to h_rows, h_cols and the received BridgeMatrix for the nonzeroes
45 std::shared_ptr<Accelerator::BlockedMatrix<Scalar>> jacMatrix; // 'stores' preconditioner matrix, actually points to h_rows, h_cols and the received BridgeMatrix for the nonzeroes
46 std::vector<int> h_rows, h_cols; // store the sparsity pattern of the matrix
47 std::vector<int> h_jacRows, h_jacCols; // store the sparsity pattern of the jacMatrix
48 std::vector<typename BridgeMatrix::size_type> diagIndices; // contains offsets of the diagonal blocks wrt start of the row, used for replaceZeroDiagonal()
49 std::vector<typename BridgeMatrix::size_type> jacDiagIndices; // same but for jacMatrix
50
51public:
61 BdaBridge(std::string accelerator_mode,
62 int linear_solver_verbosity,
63 int maxit,
64 Scalar tolerance,
65 unsigned int platformID,
66 unsigned int deviceID,
67 bool opencl_ilu_parallel,
68 std::string linsolver);
69
70
79 void solve_system(BridgeMatrix* bridgeMat,
80 BridgeMatrix* jacMat,
81 int numJacobiBlocks,
82 BridgeVector& b,
83 WellContributions<Scalar>& wellContribs,
84 InverseOperatorResult &result);
85
88 void get_result(BridgeVector &x);
89
92 bool getUseGpu()
93 {
94 return use_gpu;
95 }
96
101 static void copySparsityPatternFromISTL(const BridgeMatrix& mat,
102 std::vector<int>& h_rows,
103 std::vector<int>& h_cols);
104
109 void initWellContributions(WellContributions<Scalar>& wellContribs, unsigned N);
110
112 std::string getAccleratorName()
113 {
114 return accelerator_mode;
115 }
116}; // end class BdaBridge
117
118}
119
120#endif
BdaBridge acts as interface between opm-simulators with the BdaSolvers.
Definition: BdaBridge.hpp:37
void get_result(BridgeVector &x)
BdaBridge(std::string accelerator_mode, int linear_solver_verbosity, int maxit, Scalar tolerance, unsigned int platformID, unsigned int deviceID, bool opencl_ilu_parallel, std::string linsolver)
void solve_system(BridgeMatrix *bridgeMat, BridgeMatrix *jacMat, int numJacobiBlocks, BridgeVector &b, WellContributions< Scalar > &wellContribs, InverseOperatorResult &result)
bool getUseGpu()
Definition: BdaBridge.hpp:92
std::string getAccleratorName()
Return the selected accelerator mode, this is input via the command-line.
Definition: BdaBridge.hpp:112
void initWellContributions(WellContributions< Scalar > &wellContribs, unsigned N)
static void copySparsityPatternFromISTL(const BridgeMatrix &mat, std::vector< int > &h_rows, std::vector< int > &h_cols)
Definition: WellContributions.hpp:53
Definition: blackoilboundaryratevector.hh:37
Dune::InverseOperatorResult InverseOperatorResult
Definition: BdaBridge.hpp:32