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  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_RESID_REDUCTION_CRITERION_HH
26 #define EWOMS_ISTL_RESID_REDUCTION_CRITERION_HH
27 
28 #include "convergencecriterion.hh"
29 
30 #include <dune/istl/scalarproducts.hh>
31 
32 namespace Ewoms {
46 template <class Vector>
48 {
49  typedef typename Vector::field_type Scalar;
50 
51 public:
52  ResidReductionCriterion(Dune::ScalarProduct<Vector> &scalarProduct)
53  : scalarProduct_(scalarProduct)
54  {}
55 
56  ResidReductionCriterion(Dune::ScalarProduct<Vector> &scalarProduct,
57  Scalar reduction)
58  : scalarProduct_(scalarProduct), defectReduction_(reduction)
59  {}
60 
65  void setTolerance(Scalar tol)
66  { defectReduction_ = tol; }
67 
71  Scalar tolerance() const
72  {
73  return initialDefect_*defectReduction_;
74  }
75 
79  void setInitial(const Vector &curSol, const Vector &curResid)
80  {
81  // make sure that we don't allow an initial error of 0 to avoid
82  // divisions by zero
83  initialDefect_ = std::max(scalarProduct_.norm(curResid), 1e-20);
84  curDefect_ = initialDefect_;
85  }
86 
90  void update(const Vector &curSol, const Vector &curResid)
91  { curDefect_ = scalarProduct_.norm(curResid); }
92 
96  bool converged() const
97  { return curDefect_ <= tolerance(); }
98 
102  Scalar accuracy() const
103  { return curDefect_; }
104 
108  void printInitial(std::ostream &os=std::cout) const
109  {
110  os << std::setw(20) << " Iter ";
111  os << std::setw(20) << " Defect ";
112  os << std::setw(20) << " Reduction ";
113  os << std::endl;
114 
115  os << std::setw(20) << 0 << " ";
116  os << std::setw(20) << curDefect_ << " ";
117  os << std::setw(20) << defectReduction_ << " ";
118  os << std::endl << std::flush;
119  }
120 
124  void print(Scalar iter, std::ostream &os=std::cout) const
125  {
126  os << std::setw(20) << iter << " ";
127  os << std::setw(20) << curDefect_ << " ";
128  os << std::setw(20) << defectReduction_ << " ";
129  os << std::endl << std::flush;
130  }
131 
132 private:
133  Dune::ScalarProduct<Vector> &scalarProduct_;
134 
135  Scalar initialDefect_;
136  Scalar curDefect_;
137  Scalar defectReduction_;
138 };
139 
141 
142 } // end namespace Ewoms
143 
144 #endif
void print(Scalar iter, std::ostream &os=std::cout) const
Prints the information about the convergence behaviour for the current iteration. ...
Definition: residreductioncriterion.hh:124
bool converged() const
Returns true if and only if the convergence criterion is met.
Definition: residreductioncriterion.hh:96
Base class for all convergence criteria which only defines an virtual API.
Definition: convergencecriterion.hh:51
ResidReductionCriterion(Dune::ScalarProduct< Vector > &scalarProduct, Scalar reduction)
Definition: residreductioncriterion.hh:56
void setInitial(const Vector &curSol, const Vector &curResid)
Set the initial solution of the linear system of equations.
Definition: residreductioncriterion.hh:79
void update(const Vector &curSol, const Vector &curResid)
Update the internal members of the convergence criterion with the current solution.
Definition: residreductioncriterion.hh:90
Definition: baseauxiliarymodule.hh:35
Provides a convergence criterion which looks at the reduction of the two-norm of the residual for the...
Definition: residreductioncriterion.hh:47
ResidReductionCriterion(Dune::ScalarProduct< Vector > &scalarProduct)
Definition: residreductioncriterion.hh:52
Base class for all convergence criteria which only defines an virtual API.
void printInitial(std::ostream &os=std::cout) const
Prints the initial information about the convergence behaviour.
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
Scalar accuracy() const
Returns the accuracy of the solution at the last update.
Definition: residreductioncriterion.hh:102