File: combinedcriterion.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 (242 lines) | stat: -rw-r--r-- 7,874 bytes parent folder | download | duplicates (2)
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
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
// -*- 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::CombinedCriterion
 */
#ifndef EWOMS_COMBINED_CRITERION_HH
#define EWOMS_COMBINED_CRITERION_HH

#include "convergencecriterion.hh"

#include <iostream>

namespace Opm {
namespace Linear {

/*! \addtogroup Linear
 * \{
 */

/*!
 * \brief Convergence criterion which looks at the absolute value of the residual and
 *        fails if the linear solver stagnates.
 *
 * For the CombinedCriterion, the error of the solution is defined as \f[ e^k = \max_i\{
 * \left| r^k_i \right| \}\;, \f]
 *
 * where \f$r^k = \mathbf{A} x^k - b \f$ is the residual for the k-th iterative solution
 * vector \f$x^k\f$.
 *
 * In addition, to the reduction of the maximum residual, the linear solver is aborted
 * early if the residual goes below or above absolute limits.
 */
template <class Vector, class CollectiveCommunication>
class CombinedCriterion : public ConvergenceCriterion<Vector>
{
    using Scalar = typename Vector::field_type;
    using BlockType = typename Vector::block_type;

public:
    CombinedCriterion(const CollectiveCommunication& comm)
        : comm_(comm)
    {}

    CombinedCriterion(const CollectiveCommunication& comm,
                                       Scalar residualReductionTolerance,
                                       Scalar absResidualTolerance = 0.0,
                                       Scalar maxResidual = 0.0)
        : comm_(comm),
          residualReductionTolerance_(residualReductionTolerance),
          absResidualTolerance_(absResidualTolerance),
          maxResidual_(maxResidual)
    { }

    /*!
     * \brief Sets the residual reduction tolerance.
     */
    void setResidualReductionTolerance(Scalar tol)
    { residualReductionTolerance_ = tol; }

    /*!
     * \brief Returns the tolerance of the residual reduction of the solution.
     */
    Scalar residualReductionTolerance() const
    { return residualReductionTolerance_; }

    /*!
     * \brief Returns the reduction of the maximum of the residual compared to the
     *        initial solution.
     */
    Scalar residualReduction() const
    { return residualError_/std::max<Scalar>(1e-20, initialResidualError_); }

    /*!
     * \brief Sets the maximum absolute tolerated residual.
     */
    void setAbsResidualTolerance(Scalar tol)
    { absResidualTolerance_ = tol; }

    /*!
     * \brief Returns the tolerated maximum of the the infinity norm of the absolute
     *        residual.
     */
    Scalar absResidualTolerance() const
    { return absResidualTolerance_; }

    /*!
     * \brief Returns the infinity norm of the absolute residual.
     */
    Scalar absResidual() const
    { return residualError_; }

    /*!
     * \copydoc ConvergenceCriterion::setInitial(const Vector& , const Vector& )
     */
    void setInitial(const Vector& curSol, const Vector& curResid) override
    {
        updateErrors_(curSol, curSol, curResid);
        stagnates_ = false;

        // to avoid divisions by zero, make sure that we don't use an initial error of 0
        residualError_ = std::max<Scalar>(residualError_,
                                          std::numeric_limits<Scalar>::min()*1e10);
        initialResidualError_ = residualError_;
        lastResidualError_ = residualError_;
    }

    /*!
     * \copydoc ConvergenceCriterion::update(const Vector&, const Vector&, const Vector&)
     */
    void update(const Vector& curSol, const Vector& changeIndicator, const Vector& curResid) override
    { updateErrors_(curSol, changeIndicator, curResid);  }

    /*!
     * \copydoc ConvergenceCriterion::converged()
     */
    bool converged() const override
    {
        // we're converged if the solution is better than the tolerance
        // fix-point and residual tolerance.
        return
            residualReduction() <= residualReductionTolerance() ||
            absResidual() <= absResidualTolerance();
    }

    /*!
     * \copydoc ConvergenceCriterion::failed()
     */
    bool failed() const override
    { return !converged() && (stagnates_ || residualError_ > maxResidual_); }

    /*!
     * \copydoc ConvergenceCriterion::accuracy()
     *
     * For the accuracy we only take the residual into account,
     */
    Scalar accuracy() const override
    { return residualError_/initialResidualError_; }

    /*!
     * \copydoc ConvergenceCriterion::printInitial()
     */
    void printInitial(std::ostream& os = std::cout) const override
    {
        os << std::setw(20) << "iteration ";
        os << std::setw(20) << "residual ";
        os << std::setw(20) << "reduction ";
        os << std::setw(20) << "rate ";
        os << std::endl;
    }

    /*!
     * \copydoc ConvergenceCriterion::print()
     */
    void print(Scalar iter, std::ostream& os = std::cout) const override
    {
        const Scalar eps = std::numeric_limits<Scalar>::min()*1e10;

        os << std::setw(20) << iter << " ";
        os << std::setw(20) << absResidual() << " ";
        os << std::setw(20) << accuracy() << " ";
        os << std::setw(20) << lastResidualError_/std::max<Scalar>(residualError_, eps) << " ";
        os << std::endl << std::flush;
    }

private:
    // update the weighted absolute residual
    void updateErrors_(const Vector&, const Vector& changeIndicator,  const Vector& curResid)
    {
        lastResidualError_ = residualError_;
        residualError_ = 0.0;
        stagnates_ = true;
        for (size_t i = 0; i < curResid.size(); ++i) {
            for (unsigned j = 0; j < BlockType::dimension; ++j) {
                residualError_ =
                    std::max<Scalar>(residualError_,
                                     std::abs(curResid[i][j]));

                if (stagnates_ && changeIndicator[i][j] != 0.0)
                    // only stagnation means that we've failed!
                    stagnates_ = false;
            }
        }

        residualError_ = comm_.max(residualError_);

        // the linear solver only stagnates if all processes stagnate
        stagnates_ = comm_.min(stagnates_);
    }

    const CollectiveCommunication& comm_;

    // the infinity norm of the residual of the last iteration
    Scalar lastResidualError_;

    // the infinity norm of the residual of the current iteration
    Scalar residualError_;

    // the infinity norm of the residual of the initial solution
    Scalar initialResidualError_;

    // the minimum reduction of the residual norm where the solution is to be considered
    // converged
    Scalar residualReductionTolerance_;

    // the maximum residual norm for the residual for the solution to be considered to be
    // converged
    Scalar absResidualTolerance_;

    // The maximum error which is tolerated before we fail.
    Scalar maxResidual_;

    // does the linear solver seem to stagnate, i.e. were the last two solutions
    // identical?
    bool stagnates_;
};

//! \} end documentation

}} // end namespace Linear, Opm

#endif