[go: up one dir, main page]

Menu

[d9132d]: / src / Solvers.cpp  Maximize  Restore  History

Download this file

438 lines (394 with data), 14.3 kB

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
#include <vector>
#include "Solvers.h"
#include "math.h"
#include "MatrixMath.h"
#include <iostream>
#include "CoolPropTools.h"
#include <Eigen/Dense>
namespace CoolProp{
/** \brief Calculate the Jacobian using numerical differentiation by column
*/
std::vector<std::vector<double> > FuncWrapperND::Jacobian(const std::vector<double> &x)
{
double epsilon;
std::size_t N = x.size();
std::vector<double> r, xp;
std::vector<std::vector<double> > J(N, std::vector<double>(N, 0));
std::vector<double> r0 = call(x);
// Build the Jacobian by column
for (std::size_t i = 0; i < N; ++i)
{
xp = x;
epsilon = 0.001*x[i];
xp[i] += epsilon;
r = call(xp);
for(std::size_t j = 0; j < N; ++j)
{
J[j][i] = (r[j]-r0[j])/epsilon;
}
}
return J;
}
/**
In this formulation of the Multi-Dimensional Newton-Raphson solver the Jacobian matrix is known.
Therefore, the dx vector can be obtained from
J(x)dx=-f(x)
for a given value of x. The pointer to the class FuncWrapperND that is passed in must implement the call() and Jacobian()
functions, each of which take the vector x. The data is managed using std::vector<double> vectors
@param f A pointer to an subclass of the FuncWrapperND class that implements the call() and Jacobian() functions
@param x0 The initial guess value for the solution
@param tol The root-sum-square of the errors from each of the components
@param maxiter The maximum number of iterations
@param errstring A string with the returned error. If the length of errstring is zero, no errors were found
@returns If no errors are found, the solution. Otherwise, _HUGE, the value for infinity
*/
std::vector<double> NDNewtonRaphson_Jacobian(FuncWrapperND *f, std::vector<double> &x0, double tol, int maxiter, std::string *errstring)
{
int iter=0;
*errstring=std::string("");
std::vector<double> f0,v;
std::vector<std::vector<double> > JJ;
Eigen::VectorXd r(x0.size());
Eigen::Matrix2d J(x0.size(), x0.size());
double error = 999;
while (iter==0 || std::abs(error)>tol){
f0 = f->call(x0);
JJ = f->Jacobian(x0);
for (std::size_t i = 0; i < x0.size(); ++i)
{
r(i) = f0[i];
for (std::size_t j = 0; j < x0.size(); ++j)
{
J(i,j) = JJ[i][j];
}
}
Eigen::Vector2d v = J.colPivHouseholderQr().solve(-r);
// Update the guess
for (std::size_t i = 0; i<x0.size(); i++){ x0[i] += v(i);}
error = root_sum_square(f0);
if (iter>maxiter){
*errstring=std::string("reached maximum number of iterations");
x0[0]=_HUGE;
}
iter++;
}
return x0;
}
/**
In the newton function, a 1-D Newton-Raphson solver is implemented using exact solutions. An initial guess for the solution is provided.
@param f A pointer to an instance of the FuncWrapper1D class that implements the call() function
@param x0 The inital guess for the solution
@param ftol The absolute value of the tolerance accepted for the objective function
@param maxiter Maximum number of iterations
@param errstring A pointer to the std::string that returns the error from Secant. Length is zero if no errors are found
@returns If no errors are found, the solution, otherwise the value _HUGE, the value for infinity
*/
double Newton(FuncWrapper1DWithDeriv* f, double x0, double ftol, int maxiter, std::string &errstring)
{
double x, dx, fval=999;
int iter=1;
errstring.clear();
x = x0;
while (iter < 2 || std::abs(fval) > ftol)
{
fval = f->call(x);
dx = -fval/f->deriv(x);
if (!ValidNumber(fval)){
throw ValueError("Residual function in newton returned invalid number");
};
x += dx;
if (std::abs(dx/x) < 10*DBL_EPSILON)
{
return x;
}
if (iter>maxiter)
{
errstring= "reached maximum number of iterations";
throw SolutionError(format("Newton reached maximum number of iterations"));
}
iter=iter+1;
}
return x;
}
/**
In the Halley's method solver, two derivatives of the input variable are needed, it yields the following method:
\f[
x_{n+1} = x_n - \frac {2 f(x_n) f'(x_n)} {2 {[f'(x_n)]}^2 - f(x_n) f''(x_n)}
\f]
http://en.wikipedia.org/wiki/Halley%27s_method
@param f A pointer to an instance of the FuncWrapper1DWithTwoDerivs class that implements the call() and two derivatives
@param x0 The inital guess for the solution
@param ftol The absolute value of the tolerance accepted for the objective function
@param maxiter Maximum number of iterations
@param errstring A pointer to the std::string that returns the error from Secant. Length is zero if no errors are found
@returns If no errors are found, the solution, otherwise the value _HUGE, the value for infinity
*/
double Halley(FuncWrapper1DWithTwoDerivs* f, double x0, double ftol, int maxiter, std::string &errstring)
{
double x, dx, fval=999, dfdx, d2fdx2;
int iter=1;
errstring.clear();
x = x0;
while (iter < 2 || std::abs(fval) > ftol)
{
fval = f->call(x);
dfdx = f->deriv(x);
d2fdx2 = f->second_deriv(x);
dx = -(2*fval*dfdx)/(2*POW2(dfdx)-fval*d2fdx2);
if (!ValidNumber(fval)){
throw ValueError("Residual function in Halley returned invalid number");
};
if (!ValidNumber(dfdx)){
throw ValueError("Derivative function in Halley returned invalid number");
};
x += dx;
if (std::abs(dx/x) < 10*DBL_EPSILON){
return x;
}
if (iter>maxiter){
errstring= "reached maximum number of iterations";
throw SolutionError(format("Halley reached maximum number of iterations"));
}
iter=iter+1;
}
return x;
}
/**
In the secant function, a 1-D Newton-Raphson solver is implemented. An initial guess for the solution is provided.
@param f A pointer to an instance of the FuncWrapper1D class that implements the call() function
@param x0 The inital guess for the solutionh
@param dx The initial amount that is added to x in order to build the numerical derivative
@param tol The absolute value of the tolerance accepted for the objective function
@param maxiter Maximum number of iterations
@param errstring A pointer to the std::string that returns the error from Secant. Length is zero if no errors are found
@returns If no errors are found, the solution, otherwise the value _HUGE, the value for infinity
*/
double Secant(FuncWrapper1D* f, double x0, double dx, double tol, int maxiter, std::string &errstring)
{
#if defined(COOLPROP_DEEP_DEBUG)
static std::vector<double> xlog, flog;
xlog.clear(); flog.clear();
#endif
double x1=0,x2=0,x3=0,y1=0,y2=0,x,fval=999;
int iter=1;
errstring = "";
if (std::abs(dx)==0){ errstring="dx cannot be zero"; return _HUGE;}
while (iter<=2 || std::abs(fval)>tol)
{
if (iter==1){x1=x0; x=x1;}
if (iter==2){x2=x0+dx; x=x2;}
if (iter>2) {x=x2;}
fval = f->call(x);
#if defined(COOLPROP_DEEP_DEBUG)
xlog.push_back(x);
flog.push_back(fval);
#endif
if (!ValidNumber(fval)){
throw ValueError("Residual function in secant returned invalid number");
};
if (iter==1){y1=fval;}
if (iter>1)
{
double deltax = x2-x1;
if (std::abs(deltax)<1e-14){
return x;
}
y2=fval;
x3=x2-y2/(y2-y1)*(x2-x1);
y1=y2; x1=x2; x2=x3;
}
if (iter>maxiter)
{
errstring=std::string("reached maximum number of iterations");
throw SolutionError(format("Secant reached maximum number of iterations"));
}
iter=iter+1;
}
return x3;
}
/**
In the secant function, a 1-D Newton-Raphson solver is implemented. An initial guess for the solution is provided.
@param f A pointer to an instance of the FuncWrapper1D class that implements the call() function
@param x0 The inital guess for the solution
@param xmax The upper bound for the solution
@param xmin The lower bound for the solution
@param dx The initial amount that is added to x in order to build the numerical derivative
@param tol The absolute value of the tolerance accepted for the objective function
@param maxiter Maximum number of iterations
@param errstring A pointer to the std::string that returns the error from Secant. Length is zero if no errors are found
@returns If no errors are found, the solution, otherwise the value _HUGE, the value for infinity
*/
double BoundedSecant(FuncWrapper1D* f, double x0, double xmin, double xmax, double dx, double tol, int maxiter, std::string &errstring)
{
double x1=0,x2=0,x3=0,y1=0,y2=0,x,fval=999;
int iter=1;
errstring = "";
if (std::abs(dx)==0){ errstring = "dx cannot be zero"; return _HUGE;}
while (iter<=3 || std::abs(fval)>tol)
{
if (iter==1){x1=x0; x=x1;}
else if (iter==2){x2=x0+dx; x=x2;}
else {x=x2;}
fval=f->call(x);
if (iter==1){y1=fval;}
else
{
y2=fval;
x3=x2-y2/(y2-y1)*(x2-x1);
// Check bounds, go half the way to the limit if limit is exceeded
if (x3 < xmin)
{
x3 = (xmin + x2)/2;
}
if (x3 > xmax)
{
x3 = (xmax + x2)/2;
}
y1=y2; x1=x2; x2=x3;
}
if (iter>maxiter)
{
errstring = "reached maximum number of iterations";
throw SolutionError(format("BoundedSecant reached maximum number of iterations"));
}
iter=iter+1;
}
return x3;
}
/**
This function implements a 1-D bounded solver using the algorithm from Brent, R. P., Algorithms for Minimization Without Derivatives.
Englewood Cliffs, NJ: Prentice-Hall, 1973. Ch. 3-4.
a and b must bound the solution of interest and f(a) and f(b) must have opposite signs. If the function is continuous, there must be
at least one solution in the interval [a,b].
@param f A pointer to an instance of the FuncWrapper1D class that must implement the class() function
@param a The minimum bound for the solution of f=0
@param b The maximum bound for the solution of f=0
@param macheps The machine precision
@param t Tolerance (absolute)
@param maxiter Maximum numer of steps allowed. Will throw a SolutionError if the solution cannot be found
@param errstr A pointer to the error string returned. If length is zero, no errors found.
*/
double Brent(FuncWrapper1D* f, double a, double b, double macheps, double t, int maxiter, std::string &errstr)
{
int iter;
errstr.clear();
double fa,fb,c,fc,m,tol,d,e,p,q,s,r;
fa = f->call(a);
fb = f->call(b);
// If one of the boundaries is to within tolerance, just stop
if (std::abs(fb) < t) { return b;}
if (!ValidNumber(fb)){
throw ValueError(format("Brent's method f(b) is NAN for b = %g, other input was a = %g",b,a).c_str());
}
if (std::abs(fa) < t) { return a;}
if (!ValidNumber(fa)){
throw ValueError(format("Brent's method f(a) is NAN for a = %g, other input was b = %g",a,b).c_str());
}
if (fa*fb>0){
throw ValueError(format("Inputs in Brent [%f,%f] do not bracket the root. Function values are [%f,%f]",a,b,fa,fb));
}
c=a;
fc=fa;
iter=1;
if (std::abs(fc)<std::abs(fb)){
// Goto ext: from Brent ALGOL code
a=b;
b=c;
c=a;
fa=fb;
fb=fc;
fc=fa;
}
d=b-a;
e=b-a;
m=0.5*(c-b);
tol=2*macheps*std::abs(b)+t;
while (std::abs(m)>tol && fb!=0){
// See if a bisection is forced
if (std::abs(e)<tol || std::abs(fa) <= std::abs(fb)){
m=0.5*(c-b);
d=e=m;
}
else{
s=fb/fa;
if (a==c){
//Linear interpolation
p=2*m*s;
q=1-s;
}
else{
//Inverse quadratic interpolation
q=fa/fc;
r=fb/fc;
m=0.5*(c-b);
p=s*(2*m*q*(q-r)-(b-a)*(r-1));
q=(q-1)*(r-1)*(s-1);
}
if (p>0){
q=-q;
}
else{
p=-p;
}
s=e;
e=d;
m=0.5*(c-b);
if (2*p<3*m*q-std::abs(tol*q) || p<std::abs(0.5*s*q)){
d=p/q;
}
else{
m=0.5*(c-b);
d=e=m;
}
}
a=b;
fa=fb;
if (std::abs(d)>tol){
b+=d;
}
else if (m>0){
b+=tol;
}
else{
b+=-tol;
}
fb=f->call(b);
if (!ValidNumber(fb)){
throw ValueError(format("Brent's method f(t) is NAN for t = %g",b).c_str());
}
if (std::abs(fb) < macheps){
return b;
}
if (fb*fc>0){
// Goto int: from Brent ALGOL code
c=a;
fc=fa;
d=e=b-a;
}
if (std::abs(fc)<std::abs(fb)){
// Goto ext: from Brent ALGOL code
a=b;
b=c;
c=a;
fa=fb;
fb=fc;
fc=fa;
}
m=0.5*(c-b);
tol=2*macheps*std::abs(b)+t;
iter+=1;
if (!ValidNumber(a)){
throw ValueError(format("Brent's method a is NAN").c_str());}
if (!ValidNumber(b)){
throw ValueError(format("Brent's method b is NAN").c_str());}
if (!ValidNumber(c)){
throw ValueError(format("Brent's method c is NAN").c_str());}
if (iter>maxiter){
throw SolutionError(format("Brent's method reached maximum number of steps of %d ", maxiter));}
if (std::abs(fb)< 2*macheps*std::abs(b)){
return b;
}
}
return b;
}
}; /* namespace CoolProp */