CprCreation.hpp
Go to the documentation of this file.
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_CPRCREATION_HPP
21#define OPM_CPRCREATION_HPP
22
23
24#include <dune/istl/paamg/matrixhierarchy.hh>
25#include <dune/istl/umfpack.hh>
26
29
30#include <type_traits>
31
32namespace Opm::Accelerator {
33
34template<class Scalar> class BlockedMatrix;
35
37template <class Scalar, unsigned int block_size>
39{
40 int cprN;
41 int cprNb;
42 int cprnnz;
43 int cprnnzb;
44
45public:
47
48protected:
49
51 std::vector<Scalar> weights, coarse_vals, coarse_x, coarse_y;
52 std::vector<Matrix<Scalar>> Amatrices, Rmatrices; // scalar matrices that represent the AMG hierarchy
53 std::vector<std::vector<int> > PcolIndices; // prolongation does not need a full matrix, only store colIndices
54 std::vector<std::vector<Scalar> > invDiags; // inverse of diagonal of Amatrices
55
56 BlockedMatrix<Scalar> *mat = nullptr; // input matrix, blocked
57
58 using DuneMat = Dune::BCRSMatrix<Dune::FieldMatrix<Scalar, 1, 1> >;
59 using DuneVec = Dune::BlockVector<Dune::FieldVector<Scalar, 1> >;
60 using MatrixOperator = Dune::MatrixAdapter<DuneMat, DuneVec, DuneVec>;
61 using DuneAmg = Dune::Amg::MatrixHierarchy<MatrixOperator, Dune::Amg::SequentialInformation>;
62 std::unique_ptr<DuneAmg> dune_amg;
63 std::unique_ptr<DuneMat> dune_coarse; // extracted pressure matrix, finest level in AMG hierarchy
64 std::shared_ptr<MatrixOperator> dune_op; // operator, input to Dune AMG
65 std::vector<int> level_sizes; // size of each level in the AMG hierarchy
66 std::vector<std::vector<int> > diagIndices; // index of diagonal value for each level
67 std::conditional_t<std::is_same_v<Scalar,double>,
68 Dune::UMFPack<DuneMat>, int> umfpack; // dune/istl/umfpack object used to solve the coarsest level of AMG
69 bool always_recalculate_aggregates = false; // OPM always reuses the aggregates by default
70 bool recalculate_aggregates = true; // only rerecalculate if true
71 const int pressure_idx = 1; // hardcoded to mimic OPM
72 unsigned num_pre_smooth_steps; // number of Jacobi smooth steps before restriction
73 unsigned num_post_smooth_steps; // number of Jacobi smooth steps after prolongation
74
75 // Analyze the AMG hierarchy build by Dune
77
78 // Analyze the aggregateMaps from the AMG hierarchy
79 // These can be reused, so only use when recalculate_aggregates is true
81
83};
84
85} // namespace Opm
86
87#endif
88
Definition: MSWellHelpers.hpp:29
Definition: BlockedMatrix.hpp:29
This class implements a Constrained Pressure Residual (CPR) preconditioner.
Definition: CprCreation.hpp:39
BlockedMatrix< Scalar > * mat
Definition: CprCreation.hpp:56
bool always_recalculate_aggregates
Definition: CprCreation.hpp:69
std::shared_ptr< MatrixOperator > dune_op
Definition: CprCreation.hpp:64
int num_levels
Definition: CprCreation.hpp:50
std::vector< Matrix< Scalar > > Rmatrices
Definition: CprCreation.hpp:52
std::unique_ptr< DuneAmg > dune_amg
Definition: CprCreation.hpp:62
Dune::Amg::MatrixHierarchy< MatrixOperator, Dune::Amg::SequentialInformation > DuneAmg
Definition: CprCreation.hpp:61
std::vector< std::vector< int > > diagIndices
Definition: CprCreation.hpp:66
Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, 1, 1 > > DuneMat
Definition: CprCreation.hpp:58
std::vector< Scalar > coarse_y
Definition: CprCreation.hpp:51
std::vector< std::vector< Scalar > > invDiags
Definition: CprCreation.hpp:54
const int pressure_idx
Definition: CprCreation.hpp:71
std::conditional_t< std::is_same_v< Scalar, double >, Dune::UMFPack< DuneMat >, int > umfpack
Definition: CprCreation.hpp:68
std::vector< Scalar > coarse_vals
Definition: CprCreation.hpp:51
std::vector< int > level_sizes
Definition: CprCreation.hpp:65
unsigned num_pre_smooth_steps
Definition: CprCreation.hpp:72
std::vector< Scalar > weights
Definition: CprCreation.hpp:51
std::vector< Scalar > coarse_x
Definition: CprCreation.hpp:51
std::vector< Matrix< Scalar > > Amatrices
Definition: CprCreation.hpp:52
Dune::BlockVector< Dune::FieldVector< Scalar, 1 > > DuneVec
Definition: CprCreation.hpp:59
std::vector< std::vector< int > > PcolIndices
Definition: CprCreation.hpp:53
unsigned num_post_smooth_steps
Definition: CprCreation.hpp:73
std::unique_ptr< DuneMat > dune_coarse
Definition: CprCreation.hpp:63
bool recalculate_aggregates
Definition: CprCreation.hpp:70
Dune::MatrixAdapter< DuneMat, DuneVec, DuneVec > MatrixOperator
Definition: CprCreation.hpp:60
void create_preconditioner_amg(BlockedMatrix< Scalar > *mat)
Definition: amgclSolverBackend.hpp:44