opm-simulators
Loading...
Searching...
No Matches
Preconditioner.hpp
1/*
2 Copyright 2024 Equinor ASA
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef OPM_PRECONDITIONER_HEADER_INCLUDED
21#define OPM_PRECONDITIONER_HEADER_INCLUDED
22
23#if HAVE_OPENCL
24#include <opm/simulators/linalg/gpubridge/opencl/opencl.hpp>
25#endif
26
27#include <opm/simulators/linalg/gpubridge/WellContributions.hpp>
28
29#include <memory>
30
31namespace Opm::Accelerator {
32
33enum PreconditionerType {
34 BILU0,
35 CPR,
36 BISAI
37};
38
39template<class Scalar> class BlockedMatrix;
40
41template<class Scalar, unsigned int block_size, class ApplyScalar = Scalar>
42class Preconditioner
43{
44protected:
45 int N = 0; // number of rows of the matrix
46 int Nb = 0; // number of blockrows of the matrix
47 int nnz = 0; // number of nonzeroes of the matrix (scalar)
48 int nnzb = 0; // number of blocks of the matrix
49 int verbosity = 0;
50
51 Preconditioner(int verbosity_)
52 : verbosity(verbosity_)
53 {}
54
55public:
56
57 virtual ~Preconditioner() = default;
58
59 static std::unique_ptr<Preconditioner> create(PreconditionerType type,
60 bool opencl_ilu_parallel,
61 int verbosity);
62
63 // apply preconditioner, x = prec(y)
64 virtual void apply(const ApplyScalar& y,
65 ApplyScalar& x,
66 WellContributions<Scalar>& wellContribs) = 0;
67
68 // analyze matrix, e.g. the sparsity pattern
69 // probably only called once
70 // the version with two params can be overloaded, if not, it will default to using the one param version
71 virtual bool analyze_matrix(BlockedMatrix<Scalar>* mat) = 0;
72 virtual bool analyze_matrix(BlockedMatrix<Scalar>* mat,
73 BlockedMatrix<Scalar>* jacMat) = 0;
74
75 // create/update preconditioner, probably used every linear solve
76 // the version with two params can be overloaded, if not, it will default to using the one param version
77 virtual bool create_preconditioner(BlockedMatrix<Scalar>* mat) = 0;
78 virtual bool create_preconditioner(BlockedMatrix<Scalar>* mat,
79 BlockedMatrix<Scalar>* jacMat) = 0;
80};
81
82} // namespace Opm::Accelerator
83
84#endif
This struct resembles a blocked csr matrix, like Dune::BCRSMatrix.
Definition BlockedMatrix.hpp:29
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:51