Thermal-FIST  1.3
Package for hadron resonance gas model applications
Broyden.h
Go to the documentation of this file.
1 /*
2  * Thermal-FIST package
3  *
4  * Copyright (c) 2018-2019 Volodymyr Vovchenko
5  *
6  * GNU General Public License (GPLv3 or later)
7  */
8 
14 #ifndef BROYDEN_H
15 #define BROYDEN_H
16 
17 #include <cstddef>
18 #include <vector>
19 
20 namespace thermalfist {
21 
32  {
33  public:
35  BroydenEquations(void) { }
36 
38  virtual ~BroydenEquations(void) { }
39 
41  virtual int Dimension() const { return m_N; }
42 
48  void SetDimension(int dim) { m_N = dim; }
49 
62  virtual std::vector<double> Equations(const std::vector<double> &x) = 0;
63 
64  protected:
66  int m_N;
67  };
68 
69 
77  {
78  public:
81  static const double EPS;
89  BroydenJacobian(BroydenEquations *eqs = NULL) : m_Equations(eqs), m_dx(EPS) { }
90 
92  virtual ~BroydenJacobian(void) { }
93 
106  virtual std::vector<double> Jacobian(const std::vector<double> &x);
107 
114  void SetDx(double dx) { m_dx = dx; }
115 
119  double CurrentDx() const { return m_dx; }
120 
121  private:
123  BroydenEquations *m_Equations;
125  double m_dx;
126  };
127 
131  class Broyden
132  {
133  public:
145  public:
153  BroydenSolutionCriterium(double maximum_error = TOL) { m_MaximumError = maximum_error; }
154 
160 
172  virtual bool IsSolved(const std::vector<double>& x, const std::vector<double>& f, const std::vector<double>& xdelta = std::vector<double>()) const;
173  protected:
180  };
181  public:
183  static const double TOL;
185  static const int MAX_ITERS;
186 
197  Broyden(BroydenEquations *eqs = NULL, BroydenJacobian *jaco = NULL) : m_Equations(eqs), m_Jacobian(jaco), m_Iterations(0), m_MaxDifference(0.), m_UseNewton(false) { }
198 
203  virtual ~Broyden(void) { }
204 
217  virtual std::vector<double> Solve(const std::vector<double> &x0, BroydenSolutionCriterium *solcrit = NULL, int max_iterations = MAX_ITERS);
218 
229  int Iterations() const { return m_Iterations; }
230 
234  int MaxIterations() const { return m_MaxIterations; }
235 
239  double MaxDifference() const { return m_MaxDifference; }
240 
249  void UseNewton(bool flag) { m_UseNewton = flag; }
250 
255  bool UseNewton() const { return m_UseNewton; }
256 
257  protected:
264  };
265 
266 } // namespace thermalfist
267 
268 #endif
static const double TOL
Default desired solution accuracy.
Definition: Broyden.h:183
Class implementing the Broyden method to solve a system of non-linear equations.
Definition: Broyden.h:131
virtual ~Broyden(void)
Destroy the Broyden object.
Definition: Broyden.h:203
virtual ~BroydenEquations(void)
Destructor.
Definition: Broyden.h:38
int m_N
The number of equations.
Definition: Broyden.h:66
virtual std::vector< double > Equations(const std::vector< double > &x)=0
double m_MaxDifference
Definition: Broyden.h:262
void SetDimension(int dim)
Definition: Broyden.h:48
static const double EPS
Definition: Broyden.h:81
BroydenEquations(void)
Default constructor. Does nothing.
Definition: Broyden.h:35
Abstract class which defines the system of non-linear equations to be solved by the Broyden&#39;s method...
Definition: Broyden.h:31
void UseNewton(bool flag)
Definition: Broyden.h:249
int MaxIterations() const
Definition: Broyden.h:234
Broyden(BroydenEquations *eqs=NULL, BroydenJacobian *jaco=NULL)
Construct a new Broyden object.
Definition: Broyden.h:197
Sub-class where it is determined whether the required accuracy is achieved in the Broyden&#39;s method...
Definition: Broyden.h:144
BroydenSolutionCriterium(double maximum_error=TOL)
Definition: Broyden.h:153
bool UseNewton() const
Definition: Broyden.h:255
double MaxDifference() const
Definition: Broyden.h:239
void SetDx(double dx)
Set the finite variable difference value used for calculating the Jacobian numerically.
Definition: Broyden.h:114
Class which implements calculation of the Jacobian needed for the Broyden&#39;s method.
Definition: Broyden.h:76
BroydenEquations * m_Equations
Definition: Broyden.h:258
virtual ~BroydenSolutionCriterium()
Destroy the BroydenSolutionCriterium object.
Definition: Broyden.h:159
BroydenJacobian * m_Jacobian
Definition: Broyden.h:259
virtual int Dimension() const
Number of equations.
Definition: Broyden.h:41
The main namespace where all classes and functions of the Thermal-FIST library reside.
double CurrentDx() const
Definition: Broyden.h:119
int Iterations() const
Definition: Broyden.h:229
virtual ~BroydenJacobian(void)
Destructor.
Definition: Broyden.h:92
BroydenJacobian(BroydenEquations *eqs=NULL)
Construct a new BroydenJacobian object.
Definition: Broyden.h:89
static const int MAX_ITERS
Maximum number of Broyden iterations before terminating.
Definition: Broyden.h:185