opm-simulators/opm/simulators/linalg/bda/opencl/Preconditioner.hpp

87 lines
2.9 KiB
C++
Raw Normal View History

/*
Copyright 2021 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 OPM_PRECONDITIONER_HEADER_INCLUDED
#define OPM_PRECONDITIONER_HEADER_INCLUDED
2022-02-01 09:51:32 -06:00
#include <opm/simulators/linalg/bda/opencl/opencl.hpp>
#include <memory>
2024-04-15 15:38:04 -05:00
namespace Opm::Accelerator {
2024-04-15 09:13:22 -05:00
template<class Scalar> class BlockedMatrix;
2024-04-15 15:38:04 -05:00
template<class Scalar, unsigned int block_size>
class Preconditioner
{
protected:
int N = 0; // number of rows of the matrix
int Nb = 0; // number of blockrows of the matrix
int nnz = 0; // number of nonzeroes of the matrix (scalar)
int nnzb = 0; // number of blocks of the matrix
int verbosity = 0;
std::shared_ptr<cl::Context> context;
std::shared_ptr<cl::CommandQueue> queue;
std::vector<cl::Event> events;
cl_int err;
Preconditioner(int verbosity_) :
verbosity(verbosity_)
{};
public:
enum class Type {
BILU0,
CPR,
BISAI
};
static std::unique_ptr<Preconditioner> create(Type type,
bool opencl_ilu_parallel,
int verbosity);
virtual ~Preconditioner() = default;
// nested Preconditioners might need to override this
2024-04-15 15:38:04 -05:00
virtual void setOpencl(std::shared_ptr<cl::Context>& context,
std::shared_ptr<cl::CommandQueue>& queue);
// apply preconditioner, x = prec(y)
virtual void apply(const cl::Buffer& y, cl::Buffer& x) = 0;
// analyze matrix, e.g. the sparsity pattern
// probably only called once
// the version with two params can be overloaded, if not, it will default to using the one param version
2024-04-15 15:38:04 -05:00
virtual bool analyze_matrix(BlockedMatrix<Scalar>* mat) = 0;
virtual bool analyze_matrix(BlockedMatrix<Scalar>* mat,
BlockedMatrix<Scalar>* jacMat);
// create/update preconditioner, probably used every linear solve
// the version with two params can be overloaded, if not, it will default to using the one param version
2024-04-15 15:38:04 -05:00
virtual bool create_preconditioner(BlockedMatrix<Scalar>* mat) = 0;
virtual bool create_preconditioner(BlockedMatrix<Scalar>* mat,
BlockedMatrix<Scalar>* jacMat);
};
2024-04-15 15:38:04 -05:00
} // namespace Opm::Accelerator
#endif