weightedresidreductioncriterion.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  Copyright (C) 2009-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
25 #ifndef EWOMS_ISTL_WEIGHTED_RESIDUAL_REDUCTION_CRITERION_HH
26 #define EWOMS_ISTL_WEIGHTED_RESIDUAL_REDUCTION_CRITERION_HH
27 
28 #include "convergencecriterion.hh"
29 
30 #include <iostream>
31 
32 namespace Ewoms {
53 template <class Vector, class CollectiveCommunication>
55 {
56  typedef typename Vector::field_type Scalar;
57  typedef typename Vector::block_type BlockType;
58 
59 public:
60  WeightedResidualReductionCriterion(const CollectiveCommunication &comm)
61  : comm_(comm)
62  {}
63 
64  WeightedResidualReductionCriterion(const CollectiveCommunication &comm,
65  const Vector &residWeights,
66  Scalar fixPointTolerance,
68  Scalar absResidualTolerance = 0.0)
69  : comm_(comm),
70  residWeightVec_(residWeights),
71  fixPointTolerance_(fixPointTolerance),
72  residualReductionTolerance_(residualReductionTolerance),
73  absResidualTolerance_(absResidualTolerance)
74  {
75  Scalar minFixPointTolerance = 100*std::numeric_limits<Scalar>::epsilon();
76  fixPointTolerance_ = std::max(fixPointTolerance_, minFixPointTolerance);
77  }
78 
95  void setResidualWeight(const Vector &residWeightVec)
96  { residWeightVec_ = residWeightVec; }
97 
104  Scalar residualWeight(int outerIdx, int innerIdx) const
105  {
106  return (residWeightVec_.size() == 0)
107  ? 1.0
108  : residWeightVec_[outerIdx][innerIdx];
109  }
110 
115  { residualReductionTolerance_ = tol; }
116 
121  { return residualReductionTolerance_; }
122 
126  void setResidualTolerance(Scalar tol)
127  { absResidualTolerance_ = tol; }
128 
132  Scalar absResidualTolerance() const
133  { return absResidualTolerance_; }
134 
139  Scalar residualAccuracy() const
140  { return residualError_/std::max<Scalar>(1e-20, initialResidualError_); }
141 
145  void setFixPointTolerance(Scalar tol)
146  { fixPointTolerance_ = tol; }
147 
153  Scalar fixPointTolerance() const
154  { return fixPointTolerance_; }
155 
160  Scalar fixPointAccuracy() const
161  { return fixPointError_; }
162 
166  void setInitial(const Vector &curSol, const Vector &curResid)
167  {
168  lastResidualError_ = 1e100;
169 
170  lastSolVec_ = curSol;
171  updateErrors_(curSol, curResid);
172  // the fix-point error is not applicable for the initial solution!
173  fixPointError_ = 1e100;
174 
175  // make sure that we don't allow an initial error of 0 to avoid
176  // divisions by zero
177  residualError_ = std::max<Scalar>(residualError_, 1e-20);
178  initialResidualError_ = residualError_;
179  }
180 
184  void update(const Vector &curSol, const Vector &curResid)
185  {
186  lastResidualError_ = residualError_;
187  updateErrors_(curSol, curResid);
188  }
189 
193  bool converged() const
194  {
195  // we're converged if the solution is better than the tolerance
196  // fix-point and residual tolerance.
197  return
198  fixPointError_ <= fixPointTolerance_ ||
200  residualError_ <= absResidualTolerance_;
201  }
202 
208  Scalar accuracy() const
209  { return residualError_; }
210 
214  void printInitial(std::ostream &os = std::cout) const
215  {
216  os << std::setw(20) << " Iter ";
217  os << std::setw(20) << " Delta ";
218  os << std::setw(20) << " Residual ";
219  os << std::setw(20) << " ResidRed ";
220  os << std::setw(20) << " Rate ";
221  os << std::endl;
222 
223  os << std::setw(20) << 0 << " ";
224  os << std::setw(20) << " ";
225  os << std::setw(20) << residualError_ << " ";
226  os << std::setw(20) << 1/residualAccuracy() << " ";
227  os << std::setw(20) << " ";
228  os << std::endl << std::flush;
229  }
230 
234  void print(Scalar iter, std::ostream &os = std::cout) const
235  {
236  os << std::setw(20) << iter << " ";
237  os << std::setw(20) << fixPointAccuracy() << " ";
238  os << std::setw(20) << residualError_ << " ";
239  os << std::setw(20) << 1/residualAccuracy() << " ";
240  os << std::setw(20) << lastResidualError_ / std::max<Scalar>(residualError_, 1e-80) << " ";
241  os << std::endl << std::flush;
242  }
243 
244 private:
245  // update the weighted absolute residual
246  void updateErrors_(const Vector &curSol, const Vector &curResid)
247  {
248  residualError_ = 0.0;
249  fixPointError_ = 0.0;
250  for (size_t i = 0; i < curResid.size(); ++i) {
251  for (size_t j = 0; j < BlockType::dimension; ++j) {
252  residualError_ =
253  std::max<Scalar>(residualError_,
254  residualWeight(i, j)*std::abs(curResid[i][j]));
255  fixPointError_ =
256  std::max<Scalar>(fixPointError_,
257  std::abs(curSol[i][j] - lastSolVec_[i][j])
258  /std::max<Scalar>(1.0, curSol[i][j]));
259  }
260  }
261  lastSolVec_ = curSol;
262 
263  residualError_ = comm_.max(residualError_);
264  fixPointError_ = comm_.max(fixPointError_);
265  }
266 
267  const CollectiveCommunication &comm_;
268 
269  // the weights of the components of the residual
270  Vector residWeightVec_;
271 
272  // solution vector of the last iteration
273  Vector lastSolVec_;
274 
275  // the maximum of the weighted difference between the last two
276  // iterations
277  Scalar fixPointError_;
278 
279  // the maximum allowed relative tolerance for difference of the
280  // solution of two iterations
281  Scalar fixPointTolerance_;
282 
283  // the maximum of the absolute weighted residual of the last
284  // iteration
285  Scalar residualError_;
286 
287  // the maximum of the absolute weighted difference of the last
288  // iteration
289  Scalar lastResidualError_;
290 
291  // the maximum of the absolute weighted residual of the initial
292  // solution
293  Scalar initialResidualError_;
294 
295  // the maximum allowed relative tolerance of the residual for the
296  // solution to be considered converged
297  Scalar residualReductionTolerance_;
298 
299  // the maximum allowed absolute tolerance of the residual for the
300  // solution to be considered converged
301  Scalar absResidualTolerance_;
302 };
303 
305 
306 } // end namespace Ewoms
307 
308 #endif
Scalar accuracy() const
Returns the accuracy of the solution at the last update.
Definition: weightedresidreductioncriterion.hh:208
Scalar residualWeight(int outerIdx, int innerIdx) const
Return the relative weight of a row of the residual.
Definition: weightedresidreductioncriterion.hh:104
Base class for all convergence criteria which only defines an virtual API.
Definition: convergencecriterion.hh:51
Scalar residualAccuracy() const
Returns the reduction of the weighted maximum of the residual compared to the initial solution...
Definition: weightedresidreductioncriterion.hh:139
void setInitial(const Vector &curSol, const Vector &curResid)
Set the initial solution of the linear system of equations.
Definition: weightedresidreductioncriterion.hh:166
bool converged() const
Returns true if and only if the convergence criterion is met.
Definition: weightedresidreductioncriterion.hh:193
WeightedResidualReductionCriterion(const CollectiveCommunication &comm, const Vector &residWeights, Scalar fixPointTolerance, Scalar residualReductionTolerance, Scalar absResidualTolerance=0.0)
Definition: weightedresidreductioncriterion.hh:64
void update(const Vector &curSol, const Vector &curResid)
Update the internal members of the convergence criterion with the current solution.
Definition: weightedresidreductioncriterion.hh:184
Scalar residualReductionTolerance() const
Returns the tolerance of the residual reduction of the solution.
Definition: weightedresidreductioncriterion.hh:120
Convergence criterion which looks at the weighted absolute value of the residual. ...
Definition: weightedresidreductioncriterion.hh:54
void setResidualReductionTolerance(Scalar tol)
Sets the residual reduction tolerance.
Definition: weightedresidreductioncriterion.hh:114
Definition: baseauxiliarymodule.hh:35
void setResidualWeight(const Vector &residWeightVec)
Sets the relative weight of each row of the residual.
Definition: weightedresidreductioncriterion.hh:95
Scalar absResidualTolerance() const
Returns the maximum absolute tolerated residual.
Definition: weightedresidreductioncriterion.hh:132
Base class for all convergence criteria which only defines an virtual API.
Scalar fixPointTolerance() const
Returns the maximum tolerated difference between two iterations to be met before a solution is consid...
Definition: weightedresidreductioncriterion.hh:153
WeightedResidualReductionCriterion(const CollectiveCommunication &comm)
Definition: weightedresidreductioncriterion.hh:60
void setResidualTolerance(Scalar tol)
Sets the maximum absolute tolerated residual.
Definition: weightedresidreductioncriterion.hh:126
void print(Scalar iter, std::ostream &os=std::cout) const
Prints the information about the convergence behaviour for the current iteration. ...
Definition: weightedresidreductioncriterion.hh:234
void setFixPointTolerance(Scalar tol)
Sets the fix-point tolerance.
Definition: weightedresidreductioncriterion.hh:145
Scalar fixPointAccuracy() const
Returns the weighted maximum of the difference between the last two iterative solutions.
Definition: weightedresidreductioncriterion.hh:160
void printInitial(std::ostream &os=std::cout) const
Prints the initial information about the convergence behaviour.
Definition: weightedresidreductioncriterion.hh:214