fixpointcriterion.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_ISTL_FIXPOINT_CRITERION_HH
28#define EWOMS_ISTL_FIXPOINT_CRITERION_HH
29
31
32namespace Opm {
33namespace Linear {
34
55template <class Vector, class CollectiveCommunication>
57{
58 using Scalar = typename Vector::field_type;
59 using BlockType = typename Vector::block_type;
60
61public:
62 FixPointCriterion(const CollectiveCommunication& comm) : comm_(comm)
63 {}
64
65 FixPointCriterion(const CollectiveCommunication& comm,
66 const Vector& weightVec, Scalar reduction)
67 : comm_(comm), weightVec_(weightVec), tolerance_(reduction)
68 {}
69
86 void setWeight(const Vector& weightVec)
87 { weightVec_ = weightVec; }
88
105 Scalar weight(int outerIdx, int innerIdx) const
106 { return (weightVec_.size() == 0) ? 1.0 : weightVec_[outerIdx][innerIdx]; }
107
116 void setTolerance(Scalar tol)
117 { tolerance_ = tol; }
118
123 Scalar tolerance() const
124 { return tolerance_; }
125
129 void setInitial(const Vector& curSol, const Vector&) override
130 {
131 lastSol_ = curSol;
132 delta_ = 1000 * tolerance_;
133 }
134
138 void update(const Vector& curSol,
139 const Vector&,
140 const Vector&) override
141 {
142 assert(curSol.size() == lastSol_.size());
143
144 delta_ = 0.0;
145 for (size_t i = 0; i < curSol.size(); ++i) {
146 for (size_t j = 0; j < BlockType::dimension; ++j) {
147 delta_ =
148 std::max(delta_, weight(i, j)*std::abs(curSol[i][j] - lastSol_[i][j]));
149 }
150 }
151
152 delta_ = comm_.max(delta_);
153 lastSol_ = curSol;
154 }
155
159 bool converged() const override
160 { return accuracy() < tolerance(); }
161
165 Scalar accuracy() const
166 { return delta_; }
167
168private:
169 const CollectiveCommunication& comm_;
170
171 Vector lastSol_; // solution of the last iteration
172 Vector weightVec_; // solution of the last iteration
173 Scalar delta_; // the maximum of the absolute weighted difference of the
174 // last two iterations
175 Scalar tolerance_; // the maximum allowed delta for the solution to be
176 // considered converged
177};
178
180
181}} // end namespace Linear, Opm
182
183#endif
Base class for all convergence criteria which only defines an virtual API.
Definition: convergencecriterion.hh:56
Provides a convergence criterion for the linear solvers which looks at the weighted maximum of the di...
Definition: fixpointcriterion.hh:57
FixPointCriterion(const CollectiveCommunication &comm)
Definition: fixpointcriterion.hh:62
void update(const Vector &curSol, const Vector &, const Vector &) override
Update the internal members of the convergence criterion with the current solution.
Definition: fixpointcriterion.hh:138
Scalar weight(int outerIdx, int innerIdx) const
Return the relative weight of a primary variable.
Definition: fixpointcriterion.hh:105
void setTolerance(Scalar tol)
Set the maximum allowed weighted maximum difference between two iterations.
Definition: fixpointcriterion.hh:116
Scalar tolerance() const
Return the maximum allowed weighted difference between two iterations for the solution considered to ...
Definition: fixpointcriterion.hh:123
FixPointCriterion(const CollectiveCommunication &comm, const Vector &weightVec, Scalar reduction)
Definition: fixpointcriterion.hh:65
bool converged() const override
Returns true if and only if the convergence criterion is met.
Definition: fixpointcriterion.hh:159
void setWeight(const Vector &weightVec)
Sets the relative weight of a primary variable.
Definition: fixpointcriterion.hh:86
void setInitial(const Vector &curSol, const Vector &) override
Set the initial solution of the linear system of equations.
Definition: fixpointcriterion.hh:129
Scalar accuracy() const
Returns the accuracy of the solution at the last update.
Definition: fixpointcriterion.hh:165
Define some base class for the convergence criteria of the linear solvers of DUNE-ISTL.
Definition: blackoilboundaryratevector.hh:39