BdaBridge acts as interface between opm-simulators with the BdaSolvers.
More...
#include <BdaBridge.hpp>
|
| 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) |
|
void | get_result (BridgeVector &x) |
|
bool | getUseGpu () |
|
void | initWellContributions (WellContributions< Scalar > &wellContribs, unsigned N) |
|
std::string | getAccleratorName () |
| Return the selected accelerator mode, this is input via the command-line. More...
|
|
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, |
|
|
int |
linear_solver_verbosity, |
|
|
int |
maxit, |
|
|
Scalar |
tolerance, |
|
|
unsigned int |
platformID, |
|
|
unsigned int |
deviceID, |
|
|
bool |
opencl_ilu_parallel, |
|
|
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|amgcl|rocalution|rocsparse]' |
[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_parallel | whether to parallelize the ILU decomposition and application in OpenCL with level_scheduling |
[in] | linsolver | indicating the preconditioner, equal to the –linear-solver cmdline argument |
◆ copySparsityPatternFromISTL()
template<class BridgeMatrix , class BridgeVector , int block_size>
static 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 |
◆ getAccleratorName()
template<class BridgeMatrix , class BridgeVector , int block_size>
std::string Opm::BdaBridge< BridgeMatrix, BridgeVector, block_size >::getAccleratorName |
( |
| ) |
|
|
inline |
Return the selected accelerator mode, this is input via the command-line.
◆ getUseGpu()
template<class BridgeMatrix , class BridgeVector , int block_size>
bool Opm::BdaBridge< BridgeMatrix, BridgeVector, block_size >::getUseGpu |
( |
| ) |
|
|
inline |
Return whether the BdaBridge will use the GPU or not return whether the BdaBridge will use the GPU or not
◆ 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
-
[in] | wellContribs | container to hold all WellContributions |
[in] | N | number of rows in scalar vector that wellContribs will be applied on |
◆ solve_system()
template<class BridgeMatrix , class BridgeVector , int block_size>
Solve linear system, A*x = b - Warning
- Values of A might get overwritten!
- Parameters
-
[in] | bridgeMat | matrix A, should be of type Dune::BCRSMatrix |
[in] | jacMat | matrix A, but modified for the preconditioner, should be of type Dune::BCRSMatrix |
[in] | numJacobiBlocks | number of jacobi blocks in jacMat |
[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 file:
|