My Project
Loading...
Searching...
No Matches
BdaBridge.hpp
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
25#include <opm/simulators/linalg/bda/BdaSolver.hpp>
26
27namespace Opm
28{
29
30class WellContributions;
31
32typedef Dune::InverseOperatorResult InverseOperatorResult;
33
35template <class BridgeMatrix, class BridgeVector, int block_size>
37{
38private:
39 int verbosity = 0;
40 bool use_gpu = false;
41 std::string accelerator_mode;
42 std::unique_ptr<Opm::Accelerator::BdaSolver<block_size> > backend;
43 std::shared_ptr<Opm::Accelerator::BlockedMatrix> matrix; // 'stores' matrix, actually points to h_rows, h_cols and the received BridgeMatrix for the nonzeroes
44 std::shared_ptr<Opm::Accelerator::BlockedMatrix> jacMatrix; // 'stores' preconditioner matrix, actually points to h_rows, h_cols and the received BridgeMatrix for the nonzeroes
45 std::vector<int> h_rows, h_cols; // store the sparsity pattern of the matrix
46 std::vector<int> h_jacRows, h_jacCols; // store the sparsity pattern of the jacMatrix
47 std::vector<typename BridgeMatrix::size_type> diagIndices; // contains offsets of the diagonal blocks wrt start of the row, used for replaceZeroDiagonal()
48 std::vector<typename BridgeMatrix::size_type> jacDiagIndices; // same but for jacMatrix
49
50public:
60 BdaBridge(std::string accelerator_mode, int linear_solver_verbosity, int maxit, double tolerance,
61 unsigned int platformID, unsigned int deviceID, bool opencl_ilu_parallel, std::string linsolver);
62
63
73
76 void get_result(BridgeVector &x);
77
80 bool getUseGpu(){
81 return use_gpu;
82 }
83
88 static void copySparsityPatternFromISTL(const BridgeMatrix& mat, std::vector<int>& h_rows, std::vector<int>& h_cols);
89
95
97 std::string getAccleratorName(){
98 return accelerator_mode;
99 }
100
101}; // end class BdaBridge
102
103}
104
105#endif
Definition AquiferInterface.hpp:35
BdaBridge acts as interface between opm-simulators with the BdaSolvers.
Definition BdaBridge.hpp:37
void get_result(BridgeVector &x)
Get the resulting x vector.
Definition BdaBridge.cpp:291
bool getUseGpu()
Return whether the BdaBridge will use the GPU or not return whether the BdaBridge will use the GPU or...
Definition BdaBridge.hpp:80
void initWellContributions(WellContributions &wellContribs, unsigned N)
Initialize the WellContributions object with opencl context and queue those must be set before callin...
Definition BdaBridge.cpp:298
void solve_system(BridgeMatrix *bridgeMat, BridgeMatrix *jacMat, int numJacobiBlocks, BridgeVector &b, WellContributions &wellContribs, InverseOperatorResult &result)
Solve linear system, A*x = b.
Definition BdaBridge.cpp:204
static void copySparsityPatternFromISTL(const BridgeMatrix &mat, std::vector< int > &h_rows, std::vector< int > &h_cols)
Store sparsity pattern into vectors.
Definition BdaBridge.cpp:159
std::string getAccleratorName()
Return the selected accelerator mode, this is input via the command-line.
Definition BdaBridge.hpp:97
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:52
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27