Opm::Accelerator::openclSolverBackend< block_size > Class Template Reference

This class implements a opencl-based ilu0-bicgstab solver on GPU. More...

#include <openclSolverBackend.hpp>

Inheritance diagram for Opm::Accelerator::openclSolverBackend< block_size >:
Inheritance graph

Public Member Functions

 openclSolverBackend (int linear_solver_verbosity, int maxit, double tolerance, unsigned int platformID, unsigned int deviceID, bool opencl_ilu_parallel, std::string linsolver)
 
 openclSolverBackend (int linear_solver_verbosity, int maxit, double tolerance, bool opencl_ilu_parallel)
 For the CPR coarse solver. More...
 
SolverStatus solve_system (std::shared_ptr< BlockedMatrix > matrix, double *b, std::shared_ptr< BlockedMatrix > jacMatrix, WellContributions &wellContribs, BdaResult &res) override
 
void get_result (double *x) override
 
void setOpencl (std::shared_ptr< cl::Context > &context, std::shared_ptr< cl::CommandQueue > &queue)
 

Public Attributes

std::shared_ptr< cl::Context > context
 
std::shared_ptr< cl::CommandQueue > queue
 

Detailed Description

template<unsigned int block_size>
class Opm::Accelerator::openclSolverBackend< block_size >

This class implements a opencl-based ilu0-bicgstab solver on GPU.

Constructor & Destructor Documentation

◆ openclSolverBackend() [1/2]

template<unsigned int block_size>
Opm::Accelerator::openclSolverBackend< block_size >::openclSolverBackend ( int  linear_solver_verbosity,
int  maxit,
double  tolerance,
unsigned int  platformID,
unsigned int  deviceID,
bool  opencl_ilu_parallel,
std::string  linsolver 
)

Construct a openclSolver

Parameters
[in]linear_solver_verbosityverbosity of openclSolver
[in]maxitmaximum number of iterations for openclSolver
[in]tolerancerequired relative tolerance for openclSolver
[in]platformIDthe OpenCL platform to be used
[in]deviceIDthe device to be used
[in]opencl_ilu_parallelwhether to parallelize the ILU decomposition and application in OpenCL with level_scheduling
[in]linsolverindicating the preconditioner, equal to the –linear-solver cmdline argument only ilu0, cpr_quasiimpes and isai are supported

◆ openclSolverBackend() [2/2]

template<unsigned int block_size>
Opm::Accelerator::openclSolverBackend< block_size >::openclSolverBackend ( int  linear_solver_verbosity,
int  maxit,
double  tolerance,
bool  opencl_ilu_parallel 
)

For the CPR coarse solver.

Member Function Documentation

◆ get_result()

template<unsigned int block_size>
void Opm::Accelerator::openclSolverBackend< block_size >::get_result ( double *  x)
overridevirtual

Solve scalar linear system, for example a coarse system of an AMG preconditioner Data is already on the GPU Get result after linear solve, and peform postprocessing if necessary

Parameters
[in,out]xresulting x vector, caller must guarantee that x points to a valid array

Implements Opm::Accelerator::BdaSolver< block_size >.

◆ setOpencl()

template<unsigned int block_size>
void Opm::Accelerator::openclSolverBackend< block_size >::setOpencl ( std::shared_ptr< cl::Context > &  context,
std::shared_ptr< cl::CommandQueue > &  queue 
)

Set OpenCL objects This class either creates them based on platformID and deviceID or receives them through this function

Parameters
[in]contextthe opencl context to be used
[in]queuethe opencl queue to be used

◆ solve_system()

template<unsigned int block_size>
SolverStatus Opm::Accelerator::openclSolverBackend< block_size >::solve_system ( std::shared_ptr< BlockedMatrix matrix,
double *  b,
std::shared_ptr< BlockedMatrix jacMatrix,
WellContributions wellContribs,
BdaResult res 
)
overridevirtual

Solve linear system, A*x = b, matrix A must be in blocked-CSR format

Parameters
[in]matrixmatrix A
[in]binput vector, contains N values
[in]jacMatrixmatrix for preconditioner
[in]wellContribsWellContributions, to apply them separately, instead of adding them to matrix A
[in,out]ressummary of solver result
Returns
status code

Implements Opm::Accelerator::BdaSolver< block_size >.

Member Data Documentation

◆ context

template<unsigned int block_size>
std::shared_ptr<cl::Context> Opm::Accelerator::openclSolverBackend< block_size >::context

◆ queue

template<unsigned int block_size>
std::shared_ptr<cl::CommandQueue> Opm::Accelerator::openclSolverBackend< block_size >::queue

The documentation for this class was generated from the following file: