opm-simulators/opm/simulators/linalg/bda/BdaBridge.hpp

88 lines
3.6 KiB
C++
Raw Normal View History

/*
2019-12-05 07:24:37 -06:00
Copyright 2019 Equinor ASA
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/>.
*/
#ifndef BDABRIDGE_HEADER_INCLUDED
#define BDABRIDGE_HEADER_INCLUDED
#include "dune/istl/solver.hh" // for struct InverseOperatorResult
#include "dune/istl/bcrsmatrix.hh"
#include <opm/simulators/linalg/matrixblock.hh>
#include <opm/simulators/linalg/bda/BdaSolver.hpp>
#include <opm/simulators/linalg/bda/ILUReorder.hpp>
#include <opm/simulators/linalg/bda/WellContributions.hpp>
namespace Opm
{
typedef Dune::InverseOperatorResult InverseOperatorResult;
using bda::ILUReorder;
/// BdaBridge acts as interface between opm-simulators with the BdaSolvers
template <class BridgeMatrix, class BridgeVector, int block_size>
class BdaBridge
{
private:
bool use_gpu = false;
std::string gpu_mode;
std::unique_ptr<bda::BdaSolver<block_size> > backend;
public:
2019-12-18 10:05:33 -06:00
/// Construct a BdaBridge
/// \param[in] gpu_mode to select if a gpu solver is used, is passed via command-line: '--gpu-mode=[none|cusparse|opencl]'
/// \param[in] linear_solver_verbosity verbosity of BdaSolver
/// \param[in] maxit maximum number of iterations for BdaSolver
/// \param[in] tolerance required relative tolerance for BdaSolver
/// \param[in] platformID the OpenCL platform ID to be used
/// \param[in] deviceID the device ID to be used by the cusparse- and openclSolvers, too high values could cause runtime errors
/// \param[in] opencl_ilu_reorder select either level_scheduling or graph_coloring, see BILU0.hpp for explanation
BdaBridge(std::string gpu_mode, int linear_solver_verbosity, int maxit, double tolerance, unsigned int platformID, unsigned int deviceID, std::string opencl_ilu_reorder);
2019-12-18 10:05:33 -06:00
/// Solve linear system, A*x = b
/// \warning Values of A might get overwritten!
/// \param[in] mat matrix A, should be of type Dune::BCRSMatrix
/// \param[in] b vector b, should be of type Dune::BlockVector
/// \param[in] wellContribs contains all WellContributions, to apply them separately, instead of adding them to matrix A
/// \param[inout] result summary of solver result
void solve_system(BridgeMatrix *mat, BridgeVector &b, WellContributions& wellContribs, InverseOperatorResult &result);
2019-12-18 10:05:33 -06:00
/// Get the resulting x vector
/// \param[inout] x vector x, should be of type Dune::BlockVector
2019-12-18 08:50:09 -06:00
void get_result(BridgeVector &x);
2020-03-18 07:53:40 -05:00
/// Return whether the BdaBridge will use the GPU or not
/// return whether the BdaBridge will use the GPU or not
bool getUseGpu(){
return use_gpu;
}
/// Initialize the WellContributions object with opencl context and queue
/// those must be set before calling BlackOilWellModel::getWellContributions() in ISTL
/// \param[in] wellContribs container to hold all WellContributions
void initWellContributions(WellContributions& wellContribs);
}; // end class BdaBridge
}
#endif