BdaBridge acts as interface between opm-simulators with the BdaSolvers.
More...
#include <BdaBridge.hpp>
|
| BdaBridge (std::string accelerator_mode, std::string fpga_bitstream, int linear_solver_verbosity, int maxit, double tolerance, unsigned int platformID, unsigned int deviceID, std::string opencl_ilu_reorder, std::string linsolver) |
| Construct a BdaBridge. More...
|
|
void | solve_system (BridgeMatrix *mat, BridgeVector &b, WellContributions &wellContribs, InverseOperatorResult &result) |
| Solve linear system, A*x = b. More...
|
|
void | get_result (BridgeVector &x) |
| Get the resulting x vector. More...
|
|
bool | getUseGpu () |
| Return whether the BdaBridge will use the GPU or not return whether the BdaBridge will use the GPU or not.
|
|
void | initWellContributions (WellContributions &wellContribs) |
| Initialize the WellContributions object with opencl context and queue those must be set before calling BlackOilWellModel::getWellContributions() in ISTL. More...
|
|
bool | getUseFpga () |
| Return whether the BdaBridge will use the FPGA or not return whether the BdaBridge will use the FPGA or not.
|
|
std::string | getAccleratorName () |
| Return the selected accelerator mode, this is input via the command-line.
|
|
template<class BridgeMatrix, class BridgeVector, int block_size>
class Opm::BdaBridge< BridgeMatrix, BridgeVector, block_size >
BdaBridge acts as interface between opm-simulators with the BdaSolvers.
◆ BdaBridge()
template<class BridgeMatrix , class BridgeVector , int block_size>
Opm::BdaBridge< BridgeMatrix, BridgeVector, block_size >::BdaBridge |
( |
std::string |
accelerator_mode, |
|
|
std::string |
fpga_bitstream, |
|
|
int |
linear_solver_verbosity, |
|
|
int |
maxit, |
|
|
double |
tolerance, |
|
|
unsigned int |
platformID, |
|
|
unsigned int |
deviceID, |
|
|
std::string |
opencl_ilu_reorder, |
|
|
std::string |
linsolver |
|
) |
| |
Construct a BdaBridge.
- Parameters
-
[in] | accelerator_mode | to select if an accelerated solver is used, is passed via command-line: '–accelerator-mode=[none|cusparse|opencl|fpga]' |
[in] | fpga_bitstream | FPGA programming bitstream file name, is passed via command-line: '–fpga-bitstream=[<filename>]' |
[in] | linear_solver_verbosity | verbosity of BdaSolver |
[in] | maxit | maximum number of iterations for BdaSolver |
[in] | tolerance | required relative tolerance for BdaSolver |
[in] | platformID | the OpenCL platform ID to be used |
[in] | deviceID | the device ID to be used by the cusparse- and openclSolvers, too high values could cause runtime errors |
[in] | opencl_ilu_reorder | select either level_scheduling or graph_coloring, see ILUReorder.hpp for explanation |
[in] | linsolver | copy of cmdline argument –linsolver |
◆ copySparsityPatternFromISTL()
template<class BridgeMatrix , class BridgeVector , int block_size>
void Opm::BdaBridge< BridgeMatrix, BridgeVector, block_size >::copySparsityPatternFromISTL |
( |
const BridgeMatrix & |
mat, |
|
|
std::vector< int > & |
h_rows, |
|
|
std::vector< int > & |
h_cols |
|
) |
| |
|
static |
Store sparsity pattern into vectors.
- Parameters
-
[in] | mat | input matrix, probably BCRSMatrix |
[out] | h_rows | rowpointers |
[out] | h_cols | columnindices |
◆ get_result()
template<class BridgeMatrix , class BridgeVector , int block_size>
void Opm::BdaBridge< BridgeMatrix, BridgeVector, block_size >::get_result |
( |
BridgeVector & |
x | ) |
|
Get the resulting x vector.
- Parameters
-
[in,out] | x | vector x, should be of type Dune::BlockVector |
◆ initWellContributions()
template<class BridgeMatrix , class BridgeVector , int block_size>
Initialize the WellContributions object with opencl context and queue those must be set before calling BlackOilWellModel::getWellContributions() in ISTL.
- Parameters
-
◆ solve_system()
template<class BridgeMatrix , class BridgeVector , int block_size>
void Opm::BdaBridge< BridgeMatrix, BridgeVector, block_size >::solve_system |
( |
BridgeMatrix * |
mat, |
|
|
BridgeVector & |
b, |
|
|
WellContributions & |
wellContribs, |
|
|
InverseOperatorResult & |
result |
|
) |
| |
Solve linear system, A*x = b.
- Warning
- Values of A might get overwritten!
- Parameters
-
[in] | mat | matrix A, should be of type Dune::BCRSMatrix |
[in] | b | vector b, should be of type Dune::BlockVector |
[in] | wellContribs | contains all WellContributions, to apply them separately, instead of adding them to matrix A |
[in,out] | result | summary of solver result |
The documentation for this class was generated from the following files:
- opm/simulators/linalg/bda/BdaBridge.hpp
- opm/simulators/linalg/bda/BdaBridge.cpp