File: parallelistlbackend.hh

package info (click to toggle)
opm-models 2022.10%2Bds-4
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 6,216 kB
  • sloc: cpp: 37,910; ansic: 1,897; sh: 277; xml: 182; makefile: 10
file content (171 lines) | stat: -rw-r--r-- 6,076 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
  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 2 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/>.

  Consult the COPYING file in the top-level source directory of this
  module for the precise wording of the license and the list of
  copyright holders.
*/
/*!
 * \file
 * \copydoc Opm::Linear::ParallelIstlSolverBackend
 */
#ifndef EWOMS_PARALLEL_ISTL_BACKEND_HH
#define EWOMS_PARALLEL_ISTL_BACKEND_HH

#include "linalgproperties.hh"
#include "parallelbasebackend.hh"
#include "istlsolverwrappers.hh"
#include "istlsparsematrixadapter.hh"

#include <dune/common/version.hh>

namespace Opm::Properties::TTag {

// Create new type tag
struct ParallelIstlLinearSolver { using InheritsFrom = std::tuple<ParallelBaseLinearSolver>; };

} // namespace Opm::Properties::TTag

namespace Opm::Linear {
/*!
 * \ingroup Linear
 *
 * \brief Provides all unmodified linear solvers from dune-istl
 *
 * To set the linear solver, use
 * \code
 * template<class TypeTag>
 * struct LinearSolverWrapper<TypeTag, TTag::YourTypeTag>
 * { using type = Opm::Linear::SolverWrapper$SOLVER<TypeTag>; };
 * \endcode
 *
 * The possible choices for '\c $SOLVER' are:
 * - \c Richardson: A fixpoint solver using the Richardson iteration
 * - \c SteepestDescent: The steepest descent solver
 * - \c ConjugatedGradients: A conjugated gradients solver
 * - \c BiCGStab: A stabilized bi-conjugated gradients solver
 * - \c MinRes: A solver based on the  minimized residual algorithm
 * - \c RestartedGMRes: A restarted GMRES solver
 *
 * Chosing the preconditioner works in an analogous way:
 * \code
 * template<class TypeTag>
 * struct PreconditionerWrapper<TypeTag, TTag::YourTypeTag>
 * { using type = Opm::Linear::PreconditionerWrapper$PRECONDITIONER<TypeTag>; };
 * \endcode
 *
 * Where the choices possible for '\c $PRECONDITIONER' are:
 * - \c Jacobi: A Jacobi preconditioner
 * - \c GaussSeidel: A Gauss-Seidel preconditioner
 * - \c SSOR: A symmetric successive overrelaxation (SSOR) preconditioner
 * - \c SOR: A successive overrelaxation (SOR) preconditioner
 * - \c ILUn: An ILU(n) preconditioner
 * - \c ILU0: A specialized (and optimized) ILU(0) preconditioner
 */
template <class TypeTag>
class ParallelIstlSolverBackend : public ParallelBaseBackend<TypeTag>
{
    using ParentType = ParallelBaseBackend<TypeTag>;

    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using Simulator = GetPropType<TypeTag, Properties::Simulator>;
    using LinearSolverWrapper = GetPropType<TypeTag, Properties::LinearSolverWrapper>;
    using SparseMatrixAdapter = GetPropType<TypeTag, Properties::SparseMatrixAdapter>;

    using ParallelOperator = typename ParentType::ParallelOperator;
    using OverlappingVector = typename ParentType::OverlappingVector;
    using ParallelPreconditioner = typename ParentType::ParallelPreconditioner;
    using ParallelScalarProduct = typename ParentType::ParallelScalarProduct;

    using MatrixBlock = typename SparseMatrixAdapter::MatrixBlock;
    using RawLinearSolver = typename LinearSolverWrapper::RawSolver;

    static_assert(std::is_same<SparseMatrixAdapter, IstlSparseMatrixAdapter<MatrixBlock> >::value,
                  "The ParallelIstlSolverBackend linear solver backend requires the IstlSparseMatrixAdapter");

public:
    ParallelIstlSolverBackend(const Simulator& simulator)
        : ParentType(simulator)
    { }

    /*!
     * \brief Register all run-time parameters for the linear solver.
     */
    static void registerParameters()
    {
        ParentType::registerParameters();

        LinearSolverWrapper::registerParameters();
    }

protected:
    friend ParentType;

    std::shared_ptr<RawLinearSolver> prepareSolver_(ParallelOperator& parOperator,
                                                    ParallelScalarProduct& parScalarProduct,
                                                    ParallelPreconditioner& parPreCond)
    {
        return solverWrapper_.get(parOperator,
                                  parScalarProduct,
                                  parPreCond);
    }

    void cleanupSolver_()
    {
        solverWrapper_.cleanup();
    }

    std::pair<bool, int> runSolver_(std::shared_ptr<RawLinearSolver> solver)
    {
        Dune::InverseOperatorResult result;
        solver->apply(*this->overlappingx_, *this->overlappingb_, result);
        return std::make_pair(result.converged, result.iterations);
    }

    LinearSolverWrapper solverWrapper_;
};

} // namespace Opm::Linear

namespace Opm::Properties {

template<class TypeTag>
struct LinearSolverBackend<TypeTag, TTag::ParallelIstlLinearSolver>
{ using type = Opm::Linear::ParallelIstlSolverBackend<TypeTag>; };

template<class TypeTag>
struct LinearSolverWrapper<TypeTag, TTag::ParallelIstlLinearSolver>
{ using type = Opm::Linear::SolverWrapperBiCGStab<TypeTag>; };

#if DUNE_VERSION_NEWER(DUNE_ISTL, 2,7)
template<class TypeTag>
struct PreconditionerWrapper<TypeTag, TTag::ParallelIstlLinearSolver>
{ using type = Opm::Linear::PreconditionerWrapperILU<TypeTag>; };
#else
template<class TypeTag>
struct PreconditionerWrapper<TypeTag, TTag::ParallelIstlLinearSolver>
{ using type = Opm::Linear::PreconditionerWrapperILU0<TypeTag>; };
#endif

//! set the GMRes restart parameter to 10 by default
template<class TypeTag>
struct GMResRestart<TypeTag, TTag::ParallelIstlLinearSolver> { static constexpr int value = 10; };

} // namespace Opm::Properties

#endif