convergencecriterion.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) 2012-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_CONVERGENCE_CRITERION_HH
26 #define EWOMS_ISTL_CONVERGENCE_CRITERION_HH
27 
28 #include <dune/common/version.hh>
29 #include <dune/common/fvector.hh>
30 
31 #include <cmath>
32 #include <iostream>
33 #include <iomanip>
34 
35 namespace Ewoms {
50 template <class Vector>
52 {
53 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2,4)
54  typedef typename Dune::FieldTraits<typename Vector::field_type>::real_type real_type;
56 #else
57  typedef typename Vector::field_type real_type;
58 #endif
59 
60  typedef real_type Scalar;
61 
62 public:
70  {}
71 
85  virtual void setInitial(const Vector &curSol, const Vector &curResid) = 0;
86 
101  virtual void update(const Vector &curSol, const Vector &curResid) = 0;
102 
107  virtual bool converged() const = 0;
108 
114  virtual Scalar accuracy() const = 0;
115 
126  virtual void printInitial(std::ostream &os = std::cout) const {}
127 
136  virtual void print(Scalar iter, std::ostream &os = std::cout) const {}
137 };
138 
140 
141 } // end namespace Ewoms
142 
143 #endif
virtual Scalar accuracy() const =0
Returns the accuracy of the solution at the last update.
virtual void printInitial(std::ostream &os=std::cout) const
Prints the initial information about the convergence behaviour.
Definition: convergencecriterion.hh:126
Base class for all convergence criteria which only defines an virtual API.
Definition: convergencecriterion.hh:51
virtual void print(Scalar iter, std::ostream &os=std::cout) const
Prints the information about the convergence behaviour for the current iteration. ...
Definition: convergencecriterion.hh:136
virtual void update(const Vector &curSol, const Vector &curResid)=0
Update the internal members of the convergence criterion with the current solution.
Definition: baseauxiliarymodule.hh:35
virtual ~ConvergenceCriterion()
Destructor.
Definition: convergencecriterion.hh:69
virtual void setInitial(const Vector &curSol, const Vector &curResid)=0
Set the initial solution of the linear system of equations.
virtual bool converged() const =0
Returns true if and only if the convergence criterion is met.