residreductioncriterion.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
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 2 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 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef EWOMS_RESID_REDUCTION_CRITERION_HH
28#define EWOMS_RESID_REDUCTION_CRITERION_HH
29
31
32#include <dune/istl/scalarproducts.hh>
33
34namespace Opm {
35namespace Linear {
36
50template <class Vector>
52{
53 using Scalar = typename Vector::field_type;
54
55public:
56 ResidReductionCriterion(Dune::ScalarProduct<Vector>& scalarProduct,
57 Scalar tolerance = 1e-6)
58 : scalarProduct_(scalarProduct), tolerance_(tolerance)
59 {}
60
65 void setTolerance(Scalar tol)
66 { tolerance_ = tol; }
67
71 Scalar tolerance() const
72 { return tolerance_; }
73
77 void setInitial(const Vector&, const Vector& curResid)
78 {
79 static constexpr Scalar eps = std::numeric_limits<Scalar>::min()*1e10;
80
81 // make sure that we don't allow an initial error of 0 to avoid
82 // divisions by zero
83 curDefect_ = scalarProduct_.norm(curResid);
84 lastDefect_ = curDefect_;
85 initialDefect_ = std::max(curDefect_, eps);
86 }
87
91 void update(const Vector&,
92 const Vector&,
93 const Vector& curResid)
94 {
95 lastDefect_ = curDefect_;
96 curDefect_ = scalarProduct_.norm(curResid);
97 }
98
102 bool converged() const
103 { return accuracy() <= tolerance(); }
104
108 Scalar accuracy() const
109 { return curDefect_/initialDefect_; }
110
114 void printInitial(std::ostream& os=std::cout) const
115 {
116 os << std::setw(20) << "iteration ";
117 os << std::setw(20) << "residual ";
118 os << std::setw(20) << "accuracy ";
119 os << std::setw(20) << "rate ";
120 os << std::endl;
121 }
122
126 void print(Scalar iter, std::ostream& os=std::cout) const
127 {
128 static constexpr Scalar eps = std::numeric_limits<Scalar>::min()*1e10;
129
130 os << std::setw(20) << iter << " ";
131 os << std::setw(20) << curDefect_ << " ";
132 os << std::setw(20) << accuracy() << " ";
133 os << std::setw(20) << (lastDefect_/std::max(eps, curDefect_)) << " ";
134 os << std::endl;
135 }
136
137private:
138 Dune::ScalarProduct<Vector>& scalarProduct_;
139
140 Scalar tolerance_;
141 Scalar initialDefect_;
142 Scalar curDefect_;
143 Scalar lastDefect_;
144};
145
147
148}} // end namespace Linear, Opm
149
150#endif
Base class for all convergence criteria which only defines an virtual API.
Definition: convergencecriterion.hh:56
Provides a convergence criterion which looks at the reduction of the two-norm of the residual for the...
Definition: residreductioncriterion.hh:52
ResidReductionCriterion(Dune::ScalarProduct< Vector > &scalarProduct, Scalar tolerance=1e-6)
Definition: residreductioncriterion.hh:56
void printInitial(std::ostream &os=std::cout) const
Prints the initial information about the convergence behaviour.
Definition: residreductioncriterion.hh:114
void update(const Vector &, const Vector &, const Vector &curResid)
Definition: residreductioncriterion.hh:91
void print(Scalar iter, std::ostream &os=std::cout) const
Prints the information about the convergence behaviour for the current iteration.
Definition: residreductioncriterion.hh:126
Scalar accuracy() const
Returns the accuracy of the solution at the last update.
Definition: residreductioncriterion.hh:108
void setTolerance(Scalar tol)
Set the maximum allowed weighted maximum of the reduction of the linear residual.
Definition: residreductioncriterion.hh:65
Scalar tolerance() const
Return the maximum allowed weighted maximum of the reduction of the linear residual.
Definition: residreductioncriterion.hh:71
void setInitial(const Vector &, const Vector &curResid)
Set the initial solution of the linear system of equations.
Definition: residreductioncriterion.hh:77
bool converged() const
Returns true if and only if the convergence criterion is met.
Definition: residreductioncriterion.hh:102
Define some base class for the convergence criteria of the linear solvers of DUNE-ISTL.
Definition: blackoilboundaryratevector.hh:37