2019-12-03 14:10:21 +01:00
/*
2019-12-05 14:24:37 +01:00
Copyright 2019 Equinor ASA
2019-12-03 14:10:21 +01:00
This file is part of the Open Porous Media project ( OPM ) .
OPM is free software : you can redistribute it and / or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation , either version 3 of the License , or
( at your option ) any later version .
OPM is distributed in the hope that it will be useful ,
but WITHOUT ANY WARRANTY ; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE . See the
GNU General Public License for more details .
You should have received a copy of the GNU General Public License
along with OPM . If not , see < http : //www.gnu.org/licenses/>.
*/
# include <config.h>
2021-11-11 13:39:25 +01:00
# include "dune/istl/bcrsmatrix.hh"
# include <opm/simulators/linalg/matrixblock.hh>
2019-12-06 11:05:41 +01:00
# include <opm/common/OpmLog/OpmLog.hpp>
2019-12-18 15:47:35 +01:00
# include <opm/common/ErrorMacros.hpp>
2019-12-06 11:27:17 +01:00
# include <opm/material/common/Unused.hpp>
2019-12-06 11:05:41 +01:00
2019-12-04 16:59:58 +01:00
# include <opm/simulators/linalg/bda/BdaBridge.hpp>
# include <opm/simulators/linalg/bda/BdaResult.hpp>
2019-12-03 14:10:21 +01:00
2021-02-23 12:31:09 +01:00
# if HAVE_CUDA
# include <opm/simulators/linalg/bda/cusparseSolverBackend.hpp>
# endif
# if HAVE_OPENCL
# include <opm/simulators/linalg/bda/openclSolverBackend.hpp>
2021-11-11 14:45:12 +01:00
# include <opm/simulators/linalg/bda/openclWellContributions.hpp>
2021-02-23 12:31:09 +01:00
# endif
# if HAVE_FPGA
# include <opm/simulators/linalg/bda/FPGASolverBackend.hpp>
# endif
2021-06-02 16:19:00 +02:00
# if HAVE_AMGCL
# include <opm/simulators/linalg/bda/amgclSolverBackend.hpp>
# endif
2019-12-03 14:10:21 +01:00
typedef Dune : : InverseOperatorResult InverseOperatorResult ;
namespace Opm
{
2021-10-25 11:08:06 +02:00
using Opm : : Accelerator : : BdaResult ;
using Opm : : Accelerator : : BdaSolver ;
using Opm : : Accelerator : : SolverStatus ;
using Opm : : Accelerator : : ILUReorder ;
2020-06-22 18:26:49 +02:00
2020-06-24 16:46:04 +02:00
template < class BridgeMatrix , class BridgeVector , int block_size >
2021-04-19 09:49:20 +02:00
BdaBridge < BridgeMatrix , BridgeVector , block_size > : : BdaBridge ( std : : string accelerator_mode_ ,
[[maybe_unused]] std : : string fpga_bitstream ,
int linear_solver_verbosity , int maxit ,
2021-05-28 16:05:15 +02:00
double tolerance ,
[[maybe_unused]] unsigned int platformID ,
2021-04-19 09:49:20 +02:00
unsigned int deviceID ,
2021-12-01 11:43:30 +01:00
[[maybe_unused]] std : : string opencl_ilu_reorder ,
[[maybe_unused]] std : : string linsolver )
2021-06-02 16:19:00 +02:00
: verbosity ( linear_solver_verbosity ) , accelerator_mode ( accelerator_mode_ )
2020-02-13 11:04:02 +01:00
{
2020-12-22 12:57:01 +01:00
if ( accelerator_mode . compare ( " cusparse " ) = = 0 ) {
2020-06-23 18:19:33 +02:00
# if HAVE_CUDA
2020-06-22 18:26:49 +02:00
use_gpu = true ;
2021-10-25 11:08:06 +02:00
backend . reset ( new Opm : : Accelerator : : cusparseSolverBackend < block_size > ( linear_solver_verbosity , maxit , tolerance , deviceID ) ) ;
2020-06-23 18:19:33 +02:00
# else
OPM_THROW ( std : : logic_error , " Error cusparseSolver was chosen, but CUDA was not found by CMake " ) ;
# endif
2020-12-22 12:57:01 +01:00
} else if ( accelerator_mode . compare ( " opencl " ) = = 0 ) {
2020-06-22 18:26:49 +02:00
# if HAVE_OPENCL
use_gpu = true ;
2020-12-22 12:57:01 +01:00
ILUReorder ilu_reorder ;
if ( opencl_ilu_reorder = = " " ) {
2021-10-25 11:08:06 +02:00
ilu_reorder = Opm : : Accelerator : : ILUReorder : : GRAPH_COLORING ; // default when not selected by user
2020-12-22 12:57:01 +01:00
} else if ( opencl_ilu_reorder = = " level_scheduling " ) {
2021-10-25 11:08:06 +02:00
ilu_reorder = Opm : : Accelerator : : ILUReorder : : LEVEL_SCHEDULING ;
2020-10-16 15:05:02 +02:00
} else if ( opencl_ilu_reorder = = " graph_coloring " ) {
2021-10-25 11:08:06 +02:00
ilu_reorder = Opm : : Accelerator : : ILUReorder : : GRAPH_COLORING ;
2020-12-17 14:49:59 +01:00
} else if ( opencl_ilu_reorder = = " none " ) {
2021-10-25 11:08:06 +02:00
ilu_reorder = Opm : : Accelerator : : ILUReorder : : NONE ;
2020-10-16 15:05:02 +02:00
} else {
OPM_THROW ( std : : logic_error , " Error invalid argument for --opencl-ilu-reorder, usage: '--opencl-ilu-reorder=[level_scheduling|graph_coloring]' " ) ;
}
2021-12-01 11:43:30 +01:00
backend . reset ( new Opm : : Accelerator : : openclSolverBackend < block_size > ( linear_solver_verbosity , maxit , tolerance , platformID , deviceID , ilu_reorder , linsolver ) ) ;
2020-06-22 18:26:49 +02:00
# else
OPM_THROW ( std : : logic_error , " Error openclSolver was chosen, but OpenCL was not found by CMake " ) ;
# endif
2020-12-22 12:57:01 +01:00
} else if ( accelerator_mode . compare ( " fpga " ) = = 0 ) {
# if HAVE_FPGA
use_fpga = true ;
ILUReorder ilu_reorder ;
if ( opencl_ilu_reorder = = " " ) {
2021-10-25 11:08:06 +02:00
ilu_reorder = Opm : : Accelerator : : ILUReorder : : LEVEL_SCHEDULING ; // default when not selected by user
2020-12-22 12:57:01 +01:00
} else if ( opencl_ilu_reorder = = " level_scheduling " ) {
2021-10-25 11:08:06 +02:00
ilu_reorder = Opm : : Accelerator : : ILUReorder : : LEVEL_SCHEDULING ;
2020-12-22 12:57:01 +01:00
} else if ( opencl_ilu_reorder = = " graph_coloring " ) {
2021-10-25 11:08:06 +02:00
ilu_reorder = Opm : : Accelerator : : ILUReorder : : GRAPH_COLORING ;
2020-12-22 12:57:01 +01:00
} else {
OPM_THROW ( std : : logic_error , " Error invalid argument for --opencl-ilu-reorder, usage: '--opencl-ilu-reorder=[level_scheduling|graph_coloring]' " ) ;
}
2021-10-25 11:08:06 +02:00
backend . reset ( new Opm : : Accelerator : : FpgaSolverBackend < block_size > ( fpga_bitstream , linear_solver_verbosity , maxit , tolerance , ilu_reorder ) ) ;
2020-12-22 12:57:01 +01:00
# else
OPM_THROW ( std : : logic_error , " Error fpgaSolver was chosen, but FPGA was not enabled by CMake " ) ;
2021-06-02 16:19:00 +02:00
# endif
} else if ( accelerator_mode . compare ( " amgcl " ) = = 0 ) {
# if HAVE_AMGCL
use_gpu = true ; // should be replaced by a 'use_bridge' boolean
2021-10-25 11:08:06 +02:00
backend . reset ( new Opm : : Accelerator : : amgclSolverBackend < block_size > ( linear_solver_verbosity , maxit , tolerance , platformID , deviceID ) ) ;
2021-06-02 16:19:00 +02:00
# else
OPM_THROW ( std : : logic_error , " Error amgclSolver was chosen, but amgcl was not found by CMake " ) ;
2020-12-22 12:57:01 +01:00
# endif
} else if ( accelerator_mode . compare ( " none " ) = = 0 ) {
2020-06-23 18:19:33 +02:00
use_gpu = false ;
2020-12-22 12:57:01 +01:00
use_fpga = false ;
2020-06-23 18:19:33 +02:00
} else {
2021-06-02 16:19:00 +02:00
OPM_THROW ( std : : logic_error , " Error unknown value for parameter 'AcceleratorMode', should be passed like '--accelerator-mode=[none|cusparse|opencl|fpga|amgcl] " ) ;
2019-12-18 15:50:09 +01:00
}
2019-12-03 14:10:21 +01:00
}
2019-12-16 10:05:12 +01:00
2019-12-18 15:47:35 +01:00
2019-12-03 14:10:21 +01:00
template < class BridgeMatrix >
int checkZeroDiagonal ( BridgeMatrix & mat ) {
2019-12-18 15:50:09 +01:00
static std : : vector < typename BridgeMatrix : : size_type > diag_indices ; // contains offsets of the diagonal nnzs
int numZeros = 0 ;
2019-12-18 16:37:29 +01:00
const int dim = 3 ; // might be replaced with mat[0][0].N() or BridgeMatrix::block_type::size()
2019-12-18 15:50:09 +01:00
const double zero_replace = 1e-15 ;
2021-08-04 09:09:46 +02:00
if ( diag_indices . empty ( ) ) {
2019-12-18 15:50:09 +01:00
int N = mat . N ( ) ;
diag_indices . reserve ( N ) ;
2019-12-18 15:54:14 +01:00
for ( typename BridgeMatrix : : iterator r = mat . begin ( ) ; r ! = mat . end ( ) ; + + r ) {
2019-12-18 15:50:09 +01:00
auto diag = r - > find ( r . index ( ) ) ; // diag is an iterator
assert ( diag . index ( ) = = r . index ( ) ) ;
2019-12-18 15:54:14 +01:00
for ( int rr = 0 ; rr < dim ; + + rr ) {
2019-12-18 15:50:09 +01:00
auto & val = ( * diag ) [ rr ] [ rr ] ; // reference to easily change the value
2019-12-18 15:54:14 +01:00
if ( val = = 0.0 ) { // could be replaced by '< 1e-30' or similar
2019-12-18 15:50:09 +01:00
val = zero_replace ;
+ + numZeros ;
}
}
diag_indices . emplace_back ( diag . offset ( ) ) ;
}
} else {
2019-12-18 15:54:14 +01:00
for ( typename BridgeMatrix : : iterator r = mat . begin ( ) ; r ! = mat . end ( ) ; + + r ) {
2019-12-18 15:50:09 +01:00
typename BridgeMatrix : : size_type offset = diag_indices [ r . index ( ) ] ;
auto & diag_block = r - > getptr ( ) [ offset ] ; // diag_block is a reference to MatrixBlock, located on column r of row r
2019-12-18 15:54:14 +01:00
for ( int rr = 0 ; rr < dim ; + + rr ) {
2019-12-18 15:50:09 +01:00
auto & val = diag_block [ rr ] [ rr ] ;
2019-12-18 15:54:14 +01:00
if ( val = = 0.0 ) { // could be replaced by '< 1e-30' or similar
2019-12-18 15:50:09 +01:00
val = zero_replace ;
+ + numZeros ;
}
}
}
}
return numZeros ;
2019-12-03 14:10:21 +01:00
}
2019-12-06 14:33:52 +01:00
// iterate sparsity pattern from Matrix and put colIndices and rowPointers in arrays
2021-03-04 11:49:29 +01:00
// sparsity pattern should stay the same
2019-12-18 15:47:35 +01:00
// this could be removed if Dune::BCRSMatrix features an API call that returns colIndices and rowPointers
2021-12-01 11:43:30 +01:00
template < class BridgeMatrix , class BridgeVector , int block_size >
2021-11-24 09:56:54 +01:00
void BdaBridge < BridgeMatrix , BridgeVector , block_size > : : copySparsityPatternFromISTL ( const BridgeMatrix & mat , std : : vector < int > & h_rows , std : : vector < int > & h_cols ) {
2019-12-18 15:50:09 +01:00
2021-12-01 11:43:30 +01:00
h_rows . clear ( ) ;
h_cols . clear ( ) ;
2019-12-18 16:37:29 +01:00
2021-12-01 11:43:30 +01:00
// convert colIndices and rowPointers
h_rows . emplace_back ( 0 ) ;
for ( typename BridgeMatrix : : const_iterator r = mat . begin ( ) ; r ! = mat . end ( ) ; + + r ) {
for ( auto c = r - > begin ( ) ; c ! = r - > end ( ) ; + + c ) {
h_cols . emplace_back ( c . index ( ) ) ;
2019-12-18 15:50:09 +01:00
}
2021-11-22 16:38:04 +01:00
h_rows . emplace_back ( h_cols . size ( ) ) ;
2021-12-01 11:43:30 +01:00
}
// h_rows and h_cols could be changed to 'unsigned int', but cusparse expects 'int'
if ( static_cast < unsigned int > ( h_rows [ mat . N ( ) ] ) ! = mat . nonzeroes ( ) ) {
2021-11-24 09:56:54 +01:00
OPM_THROW ( std : : logic_error , " Error size of rows do not sum to number of nonzeroes in BdaBridge::copySparsityPatternFromISTL() " ) ;
2019-12-18 15:50:09 +01:00
}
2021-11-24 09:56:54 +01:00
}
2019-12-03 14:10:21 +01:00
2020-06-24 16:46:04 +02:00
template < class BridgeMatrix , class BridgeVector , int block_size >
2021-08-02 14:55:41 +02:00
void BdaBridge < BridgeMatrix , BridgeVector , block_size > : : solve_system ( [[maybe_unused]] BridgeMatrix * mat ,
[[maybe_unused]] BridgeVector & b ,
[[maybe_unused]] WellContributions & wellContribs ,
[[maybe_unused]] InverseOperatorResult & res )
2019-12-03 14:10:21 +01:00
{
2020-12-22 12:57:01 +01:00
if ( use_gpu | | use_fpga ) {
2019-12-18 15:50:09 +01:00
BdaResult result ;
result . converged = false ;
static std : : vector < int > h_rows ;
static std : : vector < int > h_cols ;
2019-12-18 16:37:29 +01:00
const int dim = ( * mat ) [ 0 ] [ 0 ] . N ( ) ;
2021-03-04 11:49:29 +01:00
const int Nb = mat - > N ( ) ;
const int N = Nb * dim ;
const int nnzb = ( h_rows . empty ( ) ) ? mat - > nonzeroes ( ) : h_rows . back ( ) ;
const int nnz = nnzb * dim * dim ;
2019-12-18 15:50:09 +01:00
2019-12-18 15:54:14 +01:00
if ( dim ! = 3 ) {
2021-06-02 16:19:00 +02:00
OpmLog : : warning ( " BdaSolver only accepts blocksize = 3 at this time, will use Dune for the remainder of the program " ) ;
use_gpu = use_fpga = false ;
2020-03-19 14:09:42 +01:00
return ;
2019-12-18 15:50:09 +01:00
}
2019-12-18 15:54:14 +01:00
if ( h_rows . capacity ( ) = = 0 ) {
2021-03-04 11:49:29 +01:00
h_rows . reserve ( Nb + 1 ) ;
h_cols . reserve ( nnzb ) ;
2021-11-24 09:56:54 +01:00
copySparsityPatternFromISTL ( * mat , h_rows , h_cols ) ;
2019-12-18 15:50:09 +01:00
}
2019-12-03 14:10:21 +01:00
2019-12-18 15:50:09 +01:00
Dune : : Timer t_zeros ;
int numZeros = checkZeroDiagonal ( * mat ) ;
2021-06-02 16:19:00 +02:00
if ( verbosity > = 2 ) {
std : : ostringstream out ;
out < < " Checking zeros took: " < < t_zeros . stop ( ) < < " s, found " < < numZeros < < " zeros " ;
OpmLog : : info ( out . str ( ) ) ;
}
2019-12-03 14:10:21 +01:00
2019-12-18 15:50:09 +01:00
/////////////////////////
// actually solve
2021-06-02 16:19:00 +02:00
// assume that underlying data (nonzeroes) from mat (Dune::BCRSMatrix) are contiguous, if this is not the case, the chosen BdaSolver is expected to perform undefined behaviour
2020-07-02 12:11:42 +02:00
SolverStatus status = backend - > solve_system ( N , nnz , dim , static_cast < double * > ( & ( ( ( * mat ) [ 0 ] [ 0 ] [ 0 ] [ 0 ] ) ) ) , h_rows . data ( ) , h_cols . data ( ) , static_cast < double * > ( & ( b [ 0 ] [ 0 ] ) ) , wellContribs , result ) ;
2019-12-18 15:54:14 +01:00
switch ( status ) {
2020-07-02 12:11:42 +02:00
case SolverStatus : : BDA_SOLVER_SUCCESS :
2020-06-23 11:30:15 +02:00
//OpmLog::info("BdaSolver converged");
2019-12-18 15:50:09 +01:00
break ;
2020-07-02 12:11:42 +02:00
case SolverStatus : : BDA_SOLVER_ANALYSIS_FAILED :
2020-06-22 18:26:49 +02:00
OpmLog : : warning ( " BdaSolver could not analyse level information of matrix, perhaps there is still a 0.0 on the diagonal of a block on the diagonal " ) ;
2019-12-18 15:50:09 +01:00
break ;
2020-07-02 12:11:42 +02:00
case SolverStatus : : BDA_SOLVER_CREATE_PRECONDITIONER_FAILED :
2020-06-22 18:26:49 +02:00
OpmLog : : warning ( " BdaSolver could not create preconditioner, perhaps there is still a 0.0 on the diagonal of a block on the diagonal " ) ;
2019-12-18 15:50:09 +01:00
break ;
default :
2020-06-22 18:26:49 +02:00
OpmLog : : warning ( " BdaSolver returned unknown status code " ) ;
2019-12-18 15:50:09 +01:00
}
res . iterations = result . iterations ;
res . reduction = result . reduction ;
res . converged = result . converged ;
res . conv_rate = result . conv_rate ;
res . elapsed = result . elapsed ;
2020-12-22 12:57:01 +01:00
} else {
2019-12-18 15:50:09 +01:00
res . converged = false ;
}
2019-12-03 14:10:21 +01:00
}
2020-06-24 16:46:04 +02:00
template < class BridgeMatrix , class BridgeVector , int block_size >
2021-08-02 14:55:41 +02:00
void BdaBridge < BridgeMatrix , BridgeVector , block_size > : : get_result ( [[maybe_unused]] BridgeVector & x ) {
2020-12-22 12:57:01 +01:00
if ( use_gpu | | use_fpga ) {
2020-06-22 18:26:49 +02:00
backend - > get_result ( static_cast < double * > ( & ( x [ 0 ] [ 0 ] ) ) ) ;
2019-12-18 15:50:09 +01:00
}
2019-12-03 14:10:21 +01:00
}
2021-02-23 12:31:09 +01:00
template < class BridgeMatrix , class BridgeVector , int block_size >
2021-05-28 16:05:15 +02:00
void BdaBridge < BridgeMatrix , BridgeVector , block_size > : : initWellContributions ( [[maybe_unused]] WellContributions & wellContribs ) {
2020-12-22 12:57:01 +01:00
if ( accelerator_mode . compare ( " opencl " ) = = 0 ) {
2021-02-23 12:31:09 +01:00
# if HAVE_OPENCL
2021-10-25 11:08:06 +02:00
const auto openclBackend = static_cast < const Opm : : Accelerator : : openclSolverBackend < block_size > * > ( backend . get ( ) ) ;
2021-11-11 14:45:12 +01:00
static_cast < WellContributionsOCL & > ( wellContribs ) . setOpenCLEnv ( openclBackend - > context . get ( ) , openclBackend - > queue . get ( ) ) ;
2021-02-23 12:31:09 +01:00
# else
OPM_THROW ( std : : logic_error , " Error openclSolver was chosen, but OpenCL was not found by CMake " ) ;
# endif
}
}
2021-06-01 13:50:27 +02:00
// the tests use Dune::FieldMatrix, Flow uses Opm::MatrixBlock
# define INSTANTIATE_BDA_FUNCTIONS(n) \
template class BdaBridge < Dune : : BCRSMatrix < Opm : : MatrixBlock < double , n , n > , std : : allocator < Opm : : MatrixBlock < double , n , n > > > , \
Dune : : BlockVector < Dune : : FieldVector < double , n > , std : : allocator < Dune : : FieldVector < double , n > > > , \
n > ; \
\
template class BdaBridge < Dune : : BCRSMatrix < Dune : : FieldMatrix < double , n , n > , std : : allocator < Dune : : FieldMatrix < double , n , n > > > , \
Dune : : BlockVector < Dune : : FieldVector < double , n > , std : : allocator < Dune : : FieldVector < double , n > > > , \
2021-06-01 13:12:52 +02:00
n > ;
2021-02-23 12:31:09 +01:00
2020-05-13 16:29:36 +02:00
INSTANTIATE_BDA_FUNCTIONS ( 1 ) ;
INSTANTIATE_BDA_FUNCTIONS ( 2 ) ;
INSTANTIATE_BDA_FUNCTIONS ( 3 ) ;
INSTANTIATE_BDA_FUNCTIONS ( 4 ) ;
2021-10-06 19:32:35 +02:00
INSTANTIATE_BDA_FUNCTIONS ( 5 ) ;
INSTANTIATE_BDA_FUNCTIONS ( 6 ) ;
2020-05-13 16:29:36 +02:00
# undef INSTANTIATE_BDA_FUNCTIONS
2019-12-03 14:10:21 +01:00
2020-06-24 20:09:14 +02:00
} // namespace Opm