Thermal-FIST 1.5
Package for hadron resonance gas model applications
Loading...
Searching...
No Matches
Broyden.cpp
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#include "HRGBase/Broyden.h"
9
10#include <iostream>
11#include <cmath>
12#include <stdexcept>
13
14#include <Eigen/Dense>
15
16
17using namespace Eigen;
18
19using namespace std;
20
21namespace thermalfist {
22
23 const double BroydenJacobian::EPS = 1.0E-6;
24 const double Broyden::TOL = 1.0E-10;
25 const int Broyden::MAX_ITERS = 200;
26
28 if (dim <= 0) {
29 throw std::invalid_argument("Dimension must be a positive integer.");
30 }
31 m_N = dim;
32 }
33
34
35 std::vector<double> BroydenJacobian::Jacobian(const std::vector<double>& x)
36 {
37 if (m_Equations == NULL) {
38 throw std::runtime_error("BroydenJacobian::Jacobian: Equations to solve not specified!");
39 }
40
41 if (m_Equations->Dimension() != static_cast<int>(x.size())) {
42 throw std::runtime_error("**ERROR** BroydenJacobian::Jacobian: Equations dimension does not match that of input!");
43 }
44
45 int N = m_Equations->Dimension();
46
47 std::vector<double> h = x;
48
49 std::vector< std::vector<double> > xh(N);
50
51 for (size_t i = 0; i < x.size(); ++i) {
52 h[i] = m_dx*std::abs(h[i]);
53 if (h[i] == 0.0) h[i] = m_dx;
54 if (h[i] < 1.e-10) h[i] = 1.e-10;
55 //h[i] = max(m_dx, h[i]);
56 }
57
58 for (size_t i = 0; i < x.size(); ++i) {
59 xh[i].resize(N);
60 for (size_t j = 0; j < x.size(); ++j) {
61 if (i != j)
62 xh[i][j] = x[j];
63 else
64 xh[i][j] = x[j] + h[j];
65 }
66 }
67
68 std::vector<double> Jac(N*N);
69
70 std::vector<double> fx = m_Equations->Equations(x);
71
72 for (int j = 0; j < N; ++j) {
73 std::vector<double> dfdxj = m_Equations->Equations(xh[j]);
74 for (size_t i = 0; i < dfdxj.size(); ++i) {
75 dfdxj[i] = (dfdxj[i] - fx[i]) / h[j];
76 Jac[i*N + j] = dfdxj[i];
77 }
78 }
79
80 return Jac;
81 }
82
83 std::vector<double> Broyden::Solve(const std::vector<double> &x0, BroydenSolutionCriterium *solcrit, int max_iterations)
84 {
85 if (m_Equations == NULL) {
86 throw std::runtime_error("Broyden::Solve: Equations to solve not specified!");
87 }
88
89 m_MaxIterations = max_iterations;
90
91 BroydenSolutionCriterium *SolutionCriterium = solcrit;
92 bool UseDefaultSolutionCriterium = false;
93 if (SolutionCriterium == NULL) {
94 SolutionCriterium = new BroydenSolutionCriterium(TOL);
95 UseDefaultSolutionCriterium = true;
96 }
97
98 BroydenJacobian *JacobianInUse = m_Jacobian;
99 bool UseDefaultJacobian = false;
100 if (JacobianInUse == NULL) {
101 JacobianInUse = new BroydenJacobian(m_Equations);
102 UseDefaultJacobian = true;
103 }
104 m_Iterations = 0;
105 double &maxdiff = m_MaxDifference;
106 int N = m_Equations->Dimension();
107 std::vector<double> xcur = x0, tmpvec, xdeltavec = x0;
108 VectorXd xold(N), xnew(N), xdelta(N);
109 VectorXd fold(N), fnew(N), fdelta(N);
110
111 xold = VectorXd::Map(&xcur[0], xcur.size());
112
113 MatrixXd Jac = Eigen::Map< Matrix<double, Dynamic, Dynamic, RowMajor> >(&JacobianInUse->Jacobian(xcur)[0], N, N);
114
115 if (Jac.determinant() == 0.0)
116 {
117 std::cerr << "**WARNING** Singular Jacobian in Broyden::Solve" << std::endl;
118 return xcur;
119 }
120
121 MatrixXd Jinv = Jac.inverse();
122 tmpvec = m_Equations->Equations(xcur);
123 fold = VectorXd::Map(&tmpvec[0], tmpvec.size());
124
125 for (m_Iterations = 1; m_Iterations < max_iterations; ++m_Iterations) {
126 xnew = xold - Jinv * fold;
127
128 VectorXd::Map(&xcur[0], xcur.size()) = xnew;
129
130 tmpvec = m_Equations->Equations(xcur);
131 fnew = VectorXd::Map(&tmpvec[0], tmpvec.size());
132
133
134 maxdiff = 0.;
135 for (size_t i = 0; i < xcur.size(); ++i) {
136 maxdiff = std::max(maxdiff, fabs(fnew[i]));
137 }
138
139 xdelta = xnew - xold;
140 fdelta = fnew - fold;
141
142 VectorXd::Map(&xdeltavec[0], xdeltavec.size()) = xdelta;
143
144 if (SolutionCriterium->IsSolved(xcur, tmpvec, xdeltavec))
145 break;
146 //if (maxdiff < relative_error) break;
147
148 if (!m_UseNewton) // Use Broyden's method
149 {
150
151 double norm = 0.;
152 for (int i = 0; i < N; ++i)
153 for (int j = 0; j < N; ++j)
154 norm += xdelta[i] * Jinv(i, j) * fdelta[j];
155 VectorXd p1(N);
156 p1 = (xdelta - Jinv * fdelta);
157 for (int i = 0; i < N; ++i) p1[i] *= 1. / norm;
158 Jinv = Jinv + (p1 * xdelta.transpose()) * Jinv;
159 }
160 else // Use Newton's method
161 {
162 Jac = Eigen::Map< Matrix<double, Dynamic, Dynamic, RowMajor> >(&JacobianInUse->Jacobian(xcur)[0], N, N);
163 Jinv = Jac.inverse();
164 }
165
166 xold = xnew;
167 fold = fnew;
168 }
169
170 if (UseDefaultSolutionCriterium) {
171 delete SolutionCriterium;
172 SolutionCriterium = NULL;
173 }
174 if (UseDefaultJacobian) {
175 delete JacobianInUse;
176 JacobianInUse = NULL;
177 }
178 return xcur;
179 }
180
181 bool Broyden::BroydenSolutionCriterium::IsSolved(const std::vector<double>& x, const std::vector<double>& f, const std::vector<double>& xdelta) const
182 {
183 double maxdiff = 0., maxdiffxdelta = 0.;
184 for (size_t i = 0; i < x.size(); ++i) {
185 maxdiff = std::max(maxdiff, fabs(f[i]));
186 maxdiffxdelta = std::max(maxdiffxdelta, fabs(xdelta[i]));
187 }
188 return (maxdiff < m_MaximumError) || maxdiffxdelta == 0.;
189 }
190
191} // namespace thermalfist
Implementation of the generic Broyden's method routines.
Sub-class where it is determined whether the required accuracy is achieved in the Broyden's method.
Definition Broyden.h:144
virtual bool IsSolved(const std::vector< double > &x, const std::vector< double > &f, const std::vector< double > &xdelta=std::vector< double >()) const
Definition Broyden.cpp:181
int m_N
The number of equations.
Definition Broyden.h:66
BroydenJacobian * m_Jacobian
Definition Broyden.h:259
BroydenEquations * m_Equations
Definition Broyden.h:258
static const int MAX_ITERS
Maximum number of Broyden iterations before terminating.
Definition Broyden.h:185
virtual std::vector< double > Solve(const std::vector< double > &x0, BroydenSolutionCriterium *solcrit=NULL, int max_iterations=MAX_ITERS)
Definition Broyden.cpp:83
static const double TOL
Default desired solution accuracy.
Definition Broyden.h:183
Class which implements calculation of the Jacobian needed for the Broyden's method.
Definition Broyden.h:77
virtual std::vector< double > Jacobian(const std::vector< double > &x)
Evaluates the Jacobian for given values of the variables.
Definition Broyden.cpp:35
static const double EPS
Definition Broyden.h:81
The main namespace where all classes and functions of the Thermal-FIST library reside.
Definition CosmicEoS.h:9