SHOGUN  v2.0.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
brent.cpp
Go to the documentation of this file.
1 # include <cstdlib>
2 # include <iostream>
3 # include <cmath>
4 # include <ctime>
5 
6 using namespace std;
7 
8 # include "brent.h"
9 
10 namespace shogun
11 {
12 
13 //****************************************************************************80
14 
15 double glomin ( double a, double b, double c, double m, double e, double t,
16  func_base& f, double &x )
17 
18 //****************************************************************************80
19 //
20 // Purpose:
21 //
22 // GLOMIN seeks a global minimum of a function F(X) in an interval [A,B].
23 //
24 // Discussion:
25 //
26 // This function assumes that F(X) is twice continuously differentiable
27 // over [A,B] and that F''(X) <= M for all X in [A,B].
28 //
29 // Licensing:
30 //
31 // This code is distributed under the GNU LGPL license.
32 //
33 // Modified:
34 //
35 // 06 May 2012
36 //
37 // Author:
38 //
39 // Original FORTRAN77 version by Richard Brent.
40 // C++ version by John Burkardt.
41 // Modifications by John Denker.
42 //
43 // Reference:
44 //
45 // Richard Brent,
46 // Algorithms for Minimization Without Derivatives,
47 // Dover, 2002,
48 // ISBN: 0-486-41998-3,
49 // LC: QA402.5.B74.
50 //
51 // Parameters:
52 //
53 // Input, double A, B, the endpoints of the interval.
54 // It must be the case that A < B.
55 //
56 // Input, double C, an initial guess for the global
57 // minimizer. If no good guess is known, C = A or B is acceptable.
58 //
59 // Input, double M, the bound on the second derivative.
60 //
61 // Input, double E, a positive tolerance, a bound for the
62 // absolute error in the evaluation of F(X) for any X in [A,B].
63 //
64 // Input, double T, a positive error tolerance.
65 //
66 // Input, func_base& F, a user-supplied c++ functor whose
67 // global minimum is being sought. The input and output
68 // of F() are of type double.
69 //
70 // Output, double &X, the estimated value of the abscissa
71 // for which F attains its global minimum value in [A,B].
72 //
73 // Output, double GLOMIN, the value F(X).
74 //
75 {
76  double a0;
77  double a2;
78  double a3;
79  double d0;
80  double d1;
81  double d2;
82  double h;
83  int k;
84  double m2;
85  double macheps;
86  double p;
87  double q;
88  double qs;
89  double r;
90  double s;
91  double sc;
92  double y;
93  double y0;
94  double y1;
95  double y2;
96  double y3;
97  double yb;
98  double z0;
99  double z1;
100  double z2;
101 
102  a0 = b;
103  x = a0;
104  a2 = a;
105  y0 = f ( b );
106  yb = y0;
107  y2 = f ( a );
108  y = y2;
109 
110  if ( y0 < y )
111  {
112  y = y0;
113  }
114  else
115  {
116  x = a;
117  }
118 
119  if ( m <= 0.0 || b <= a )
120  {
121  return y;
122  }
123 
124  macheps = r8_epsilon ( );
125 
126  m2 = 0.5 * ( 1.0 + 16.0 * macheps ) * m;
127 
128  if ( c <= a || b <= c )
129  {
130  sc = 0.5 * ( a + b );
131  }
132  else
133  {
134  sc = c;
135  }
136 
137  y1 = f ( sc );
138  k = 3;
139  d0 = a2 - sc;
140  h = 9.0 / 11.0;
141 
142  if ( y1 < y )
143  {
144  x = sc;
145  y = y1;
146  }
147 //
148 // Loop.
149 //
150  for ( ; ; )
151  {
152  d1 = a2 - a0;
153  d2 = sc - a0;
154  z2 = b - a2;
155  z0 = y2 - y1;
156  z1 = y2 - y0;
157  r = d1 * d1 * z0 - d0 * d0 * z1;
158  p = r;
159  qs = 2.0 * ( d0 * z1 - d1 * z0 );
160  q = qs;
161 
162  if ( k < 1000000 || y2 <= y )
163  {
164  for ( ; ; )
165  {
166  if ( q * ( r * ( yb - y2 ) + z2 * q * ( ( y2 - y ) + t ) ) <
167  z2 * m2 * r * ( z2 * q - r ) )
168  {
169  a3 = a2 + r / q;
170  y3 = f ( a3 );
171 
172  if ( y3 < y )
173  {
174  x = a3;
175  y = y3;
176  }
177  }
178  k = ( ( 1611 * k ) % 1048576 );
179  q = 1.0;
180  r = ( b - a ) * 0.00001 * ( double ) ( k );
181 
182  if ( z2 <= r )
183  {
184  break;
185  }
186  }
187  }
188  else
189  {
190  k = ( ( 1611 * k ) % 1048576 );
191  q = 1.0;
192  r = ( b - a ) * 0.00001 * ( double ) ( k );
193 
194  while ( r < z2 )
195  {
196  if ( q * ( r * ( yb - y2 ) + z2 * q * ( ( y2 - y ) + t ) ) <
197  z2 * m2 * r * ( z2 * q - r ) )
198  {
199  a3 = a2 + r / q;
200  y3 = f ( a3 );
201 
202  if ( y3 < y )
203  {
204  x = a3;
205  y = y3;
206  }
207  }
208  k = ( ( 1611 * k ) % 1048576 );
209  q = 1.0;
210  r = ( b - a ) * 0.00001 * ( double ) ( k );
211  }
212  }
213 
214  r = m2 * d0 * d1 * d2;
215  s = sqrt ( ( ( y2 - y ) + t ) / m2 );
216  h = 0.5 * ( 1.0 + h );
217  p = h * ( p + 2.0 * r * s );
218  q = q + 0.5 * qs;
219  r = - 0.5 * ( d0 + ( z0 + 2.01 * e ) / ( d0 * m2 ) );
220 
221  if ( r < s || d0 < 0.0 )
222  {
223  r = a2 + s;
224  }
225  else
226  {
227  r = a2 + r;
228  }
229 
230  if ( 0.0 < p * q )
231  {
232  a3 = a2 + p / q;
233  }
234  else
235  {
236  a3 = r;
237  }
238 
239  for ( ; ; )
240  {
241  a3 = r8_max ( a3, r );
242 
243  if ( b <= a3 )
244  {
245  a3 = b;
246  y3 = yb;
247  }
248  else
249  {
250  y3 = f ( a3 );
251  }
252 
253  if ( y3 < y )
254  {
255  x = a3;
256  y = y3;
257  }
258 
259  d0 = a3 - a2;
260 
261  if ( a3 <= r )
262  {
263  break;
264  }
265 
266  p = 2.0 * ( y2 - y3 ) / ( m * d0 );
267 
268  if ( ( 1.0 + 9.0 * macheps ) * d0 <= r8_abs ( p ) )
269  {
270  break;
271  }
272 
273  if ( 0.5 * m2 * ( d0 * d0 + p * p ) <= ( y2 - y ) + ( y3 - y ) + 2.0 * t )
274  {
275  break;
276  }
277  a3 = 0.5 * ( a2 + a3 );
278  h = 0.9 * h;
279  }
280 
281  if ( b <= a3 )
282  {
283  break;
284  }
285 
286  a0 = sc;
287  sc = a2;
288  a2 = a3;
289  y0 = y1;
290  y1 = y2;
291  y2 = y3;
292  }
293 
294  return y;
295 }
296 //****************************************************************************80
297 
298 double local_min ( double a, double b, double t, func_base& f,
299  double &x )
300 
301 //****************************************************************************80
302 //
303 // Purpose:
304 //
305 // LOCAL_MIN seeks a local minimum of a function F(X) in an interval [A,B].
306 //
307 // Discussion:
308 //
309 // The method used is a combination of golden section search and
310 // successive parabolic interpolation. Convergence is never much slower
311 // than that for a Fibonacci search. If F has a continuous second
312 // derivative which is positive at the minimum (which is not at A or
313 // B), then convergence is superlinear, and usually of the order of
314 // about 1.324....
315 //
316 // The values EPS and T define a tolerance TOL = EPS * abs ( X ) + T.
317 // F is never evaluated at two points closer than TOL.
318 //
319 // If F is a unimodal function and the computed values of F are always
320 // unimodal when separated by at least SQEPS * abs ( X ) + (T/3), then
321 // LOCAL_MIN approximates the abscissa of the global minimum of F on the
322 // interval [A,B] with an error less than 3*SQEPS*abs(LOCAL_MIN)+T.
323 //
324 // If F is not unimodal, then LOCAL_MIN may approximate a local, but
325 // perhaps non-global, minimum to the same accuracy.
326 //
327 // Licensing:
328 //
329 // This code is distributed under the GNU LGPL license.
330 //
331 // Modified:
332 //
333 // 17 July 2011
334 //
335 // Author:
336 //
337 // Original FORTRAN77 version by Richard Brent.
338 // C++ version by John Burkardt.
339 // Modifications by John Denker.
340 //
341 // Reference:
342 //
343 // Richard Brent,
344 // Algorithms for Minimization Without Derivatives,
345 // Dover, 2002,
346 // ISBN: 0-486-41998-3,
347 // LC: QA402.5.B74.
348 //
349 // Parameters:
350 //
351 // Input, double A, B, the endpoints of the interval.
352 //
353 // Input, double T, a positive absolute error tolerance.
354 //
355 // Input, func_base& F, a user-supplied c++ functor whose
356 // local minimum is being sought. The input and output
357 // of F() are of type double.
358 //
359 // Output, double &X, the estimated value of an abscissa
360 // for which F attains a local minimum value in [A,B].
361 //
362 // Output, double LOCAL_MIN, the value F(X).
363 //
364 {
365  double c;
366  double d = 0.0;
367  double e;
368  double eps;
369  double fu;
370  double fv;
371  double fw;
372  double fx;
373  double m;
374  double p;
375  double q;
376  double r;
377  double sa;
378  double sb;
379  double t2;
380  double tol;
381  double u;
382  double v;
383  double w;
384 //
385 // C is the square of the inverse of the golden ratio.
386 //
387  c = 0.5 * ( 3.0 - sqrt ( 5.0 ) );
388 
389  eps = sqrt ( r8_epsilon ( ) );
390 
391  sa = a;
392  sb = b;
393  x = sa + c * ( b - a );
394  w = x;
395  v = w;
396  e = 0.0;
397  fx = f ( x );
398  fw = fx;
399  fv = fw;
400 
401  for ( ; ; )
402  {
403  m = 0.5 * ( sa + sb ) ;
404  tol = eps * r8_abs ( x ) + t;
405  t2 = 2.0 * tol;
406 //
407 // Check the stopping criterion.
408 //
409  if ( r8_abs ( x - m ) <= t2 - 0.5 * ( sb - sa ) )
410  {
411  break;
412  }
413 //
414 // Fit a parabola.
415 //
416  r = 0.0;
417  q = r;
418  p = q;
419 
420  if ( tol < r8_abs ( e ) )
421  {
422  r = ( x - w ) * ( fx - fv );
423  q = ( x - v ) * ( fx - fw );
424  p = ( x - v ) * q - ( x - w ) * r;
425  q = 2.0 * ( q - r );
426  if ( 0.0 < q )
427  {
428  p = - p;
429  }
430  q = r8_abs ( q );
431  r = e;
432  e = d;
433  }
434 
435  if ( r8_abs ( p ) < r8_abs ( 0.5 * q * r ) &&
436  q * ( sa - x ) < p &&
437  p < q * ( sb - x ) )
438  {
439 //
440 // Take the parabolic interpolation step.
441 //
442  d = p / q;
443  u = x + d;
444 //
445 // F must not be evaluated too close to A or B.
446 //
447  if ( ( u - sa ) < t2 || ( sb - u ) < t2 )
448  {
449  if ( x < m )
450  {
451  d = tol;
452  }
453  else
454  {
455  d = - tol;
456  }
457  }
458  }
459 //
460 // A golden-section step.
461 //
462  else
463  {
464  if ( x < m )
465  {
466  e = sb - x;
467  }
468  else
469  {
470  e = sa - x;
471  }
472  d = c * e;
473  }
474 //
475 // F must not be evaluated too close to X.
476 //
477  if ( tol <= r8_abs ( d ) )
478  {
479  u = x + d;
480  }
481  else if ( 0.0 < d )
482  {
483  u = x + tol;
484  }
485  else
486  {
487  u = x - tol;
488  }
489 
490  fu = f ( u );
491 //
492 // Update A, B, V, W, and X.
493 //
494  if ( fu <= fx )
495  {
496  if ( u < x )
497  {
498  sb = x;
499  }
500  else
501  {
502  sa = x;
503  }
504  v = w;
505  fv = fw;
506  w = x;
507  fw = fx;
508  x = u;
509  fx = fu;
510  }
511  else
512  {
513  if ( u < x )
514  {
515  sa = u;
516  }
517  else
518  {
519  sb = u;
520  }
521 
522  if ( fu <= fw || w == x )
523  {
524  v = w;
525  fv = fw;
526  w = u;
527  fw = fu;
528  }
529  else if ( fu <= fv || v == x || v== w )
530  {
531  v = u;
532  fv = fu;
533  }
534  }
535  }
536  return fx;
537 }
538 //****************************************************************************80
539 
540 double local_min_rc ( double &a, double &b, int &status, double value )
541 
542 //****************************************************************************80
543 //
544 // Purpose:
545 //
546 // LOCAL_MIN_RC seeks a minimizer of a scalar function of a scalar variable.
547 //
548 // Discussion:
549 //
550 // This routine seeks an approximation to the point where a function
551 // F attains a minimum on the interval (A,B).
552 //
553 // The method used is a combination of golden section search and
554 // successive parabolic interpolation. Convergence is never much
555 // slower than that for a Fibonacci search. If F has a continuous
556 // second derivative which is positive at the minimum (which is not
557 // at A or B), then convergence is superlinear, and usually of the
558 // order of about 1.324...
559 //
560 // The routine is a revised version of the Brent local minimization
561 // algorithm, using reverse communication.
562 //
563 // It is worth stating explicitly that this routine will NOT be
564 // able to detect a minimizer that occurs at either initial endpoint
565 // A or B. If this is a concern to the user, then the user must
566 // either ensure that the initial interval is larger, or to check
567 // the function value at the returned minimizer against the values
568 // at either endpoint.
569 //
570 // Licensing:
571 //
572 // This code is distributed under the GNU LGPL license.
573 //
574 // Modified:
575 //
576 // 17 July 2011
577 //
578 // Author:
579 //
580 // John Burkardt
581 //
582 // Reference:
583 //
584 // Richard Brent,
585 // Algorithms for Minimization Without Derivatives,
586 // Dover, 2002,
587 // ISBN: 0-486-41998-3,
588 // LC: QA402.5.B74.
589 //
590 // David Kahaner, Cleve Moler, Steven Nash,
591 // Numerical Methods and Software,
592 // Prentice Hall, 1989,
593 // ISBN: 0-13-627258-4,
594 // LC: TA345.K34.
595 //
596 // Parameters
597 //
598 // Input/output, double &A, &B. On input, the left and right
599 // endpoints of the initial interval. On output, the lower and upper
600 // bounds for an interval containing the minimizer. It is required
601 // that A < B.
602 //
603 // Input/output, int &STATUS, used to communicate between
604 // the user and the routine. The user only sets STATUS to zero on the first
605 // call, to indicate that this is a startup call. The routine returns STATUS
606 // positive to request that the function be evaluated at ARG, or returns
607 // STATUS as 0, to indicate that the iteration is complete and that
608 // ARG is the estimated minimizer.
609 //
610 // Input, double VALUE, the function value at ARG, as requested
611 // by the routine on the previous call.
612 //
613 // Output, double LOCAL_MIN_RC, the currently considered point.
614 // On return with STATUS positive, the user is requested to evaluate the
615 // function at this point, and return the value in VALUE. On return with
616 // STATUS zero, this is the routine's estimate for the function minimizer.
617 //
618 // Local parameters:
619 //
620 // C is the squared inverse of the golden ratio.
621 //
622 // EPS is the square root of the relative machine precision.
623 //
624 {
625  static double arg;
626  static double c;
627  static double d;
628  static double e;
629  static double eps;
630  static double fu;
631  static double fv;
632  static double fw;
633  static double fx;
634  static double midpoint;
635  static double p;
636  static double q;
637  static double r;
638  static double tol;
639  static double tol1;
640  static double tol2;
641  static double u;
642  static double v;
643  static double w;
644  static double x;
645 //
646 // STATUS (INPUT) = 0, startup.
647 //
648  if ( status == 0 )
649  {
650  if ( b <= a )
651  {
652  cout << "\n";
653  cout << "LOCAL_MIN_RC - Fatal error!\n";
654  cout << " A < B is required, but\n";
655  cout << " A = " << a << "\n";
656  cout << " B = " << b << "\n";
657  status = -1;
658  exit ( 1 );
659  }
660  c = 0.5 * ( 3.0 - sqrt ( 5.0 ) );
661 
662  eps = sqrt ( r8_epsilon ( ) );
663  tol = r8_epsilon ( );
664 
665  v = a + c * ( b - a );
666  w = v;
667  x = v;
668  e = 0.0;
669 
670  status = 1;
671  arg = x;
672 
673  return arg;
674  }
675 //
676 // STATUS (INPUT) = 1, return with initial function value of FX.
677 //
678  else if ( status == 1 )
679  {
680  fx = value;
681  fv = fx;
682  fw = fx;
683  }
684 //
685 // STATUS (INPUT) = 2 or more, update the data.
686 //
687  else if ( 2 <= status )
688  {
689  fu = value;
690 
691  if ( fu <= fx )
692  {
693  if ( x <= u )
694  {
695  a = x;
696  }
697  else
698  {
699  b = x;
700  }
701  v = w;
702  fv = fw;
703  w = x;
704  fw = fx;
705  x = u;
706  fx = fu;
707  }
708  else
709  {
710  if ( u < x )
711  {
712  a = u;
713  }
714  else
715  {
716  b = u;
717  }
718 
719  if ( fu <= fw || w == x )
720  {
721  v = w;
722  fv = fw;
723  w = u;
724  fw = fu;
725  }
726  else if ( fu <= fv || v == x || v == w )
727  {
728  v = u;
729  fv = fu;
730  }
731  }
732  }
733 //
734 // Take the next step.
735 //
736  midpoint = 0.5 * ( a + b );
737  tol1 = eps * r8_abs ( x ) + tol / 3.0;
738  tol2 = 2.0 * tol1;
739 //
740 // If the stopping criterion is satisfied, we can exit.
741 //
742  if ( r8_abs ( x - midpoint ) <= ( tol2 - 0.5 * ( b - a ) ) )
743  {
744  status = 0;
745  return arg;
746  }
747 //
748 // Is golden-section necessary?
749 //
750  if ( r8_abs ( e ) <= tol1 )
751  {
752  if ( midpoint <= x )
753  {
754  e = a - x;
755  }
756  else
757  {
758  e = b - x;
759  }
760  d = c * e;
761  }
762 //
763 // Consider fitting a parabola.
764 //
765  else
766  {
767  r = ( x - w ) * ( fx - fv );
768  q = ( x - v ) * ( fx - fw );
769  p = ( x - v ) * q - ( x - w ) * r;
770  q = 2.0 * ( q - r );
771  if ( 0.0 < q )
772  {
773  p = - p;
774  }
775  q = r8_abs ( q );
776  r = e;
777  e = d;
778 //
779 // Choose a golden-section step if the parabola is not advised.
780 //
781  if (
782  ( r8_abs ( 0.5 * q * r ) <= r8_abs ( p ) ) ||
783  ( p <= q * ( a - x ) ) ||
784  ( q * ( b - x ) <= p ) )
785  {
786  if ( midpoint <= x )
787  {
788  e = a - x;
789  }
790  else
791  {
792  e = b - x;
793  }
794  d = c * e;
795  }
796 //
797 // Choose a parabolic interpolation step.
798 //
799  else
800  {
801  d = p / q;
802  u = x + d;
803 
804  if ( ( u - a ) < tol2 )
805  {
806  d = tol1 * r8_sign ( midpoint - x );
807  }
808 
809  if ( ( b - u ) < tol2 )
810  {
811  d = tol1 * r8_sign ( midpoint - x );
812  }
813  }
814  }
815 //
816 // F must not be evaluated too close to X.
817 //
818  if ( tol1 <= r8_abs ( d ) )
819  {
820  u = x + d;
821  }
822  if ( r8_abs ( d ) < tol1 )
823  {
824  u = x + tol1 * r8_sign ( d );
825  }
826 //
827 // Request value of F(U).
828 //
829  arg = u;
830  status = status + 1;
831 
832  return arg;
833 }
834 //****************************************************************************80
835 
836 double r8_abs ( double x )
837 
838 //****************************************************************************80
839 //
840 // Purpose:
841 //
842 // R8_ABS returns the absolute value of an R8.
843 //
844 // Licensing:
845 //
846 // This code is distributed under the GNU LGPL license.
847 //
848 // Modified:
849 //
850 // 07 May 2006
851 //
852 // Author:
853 //
854 // John Burkardt
855 //
856 // Parameters:
857 //
858 // Input, double X, the quantity whose absolute value is desired.
859 //
860 // Output, double R8_ABS, the absolute value of X.
861 //
862 {
863  double value;
864 
865  if ( 0.0 <= x )
866  {
867  value = x;
868  }
869  else
870  {
871  value = - x;
872  }
873  return value;
874 }
875 //****************************************************************************80
876 
877 double r8_epsilon ( )
878 
879 //****************************************************************************80
880 //
881 // Purpose:
882 //
883 // R8_EPSILON returns the R8 round off unit.
884 //
885 // Discussion:
886 //
887 // R8_EPSILON is a number R which is a power of 2 with the property that,
888 // to the precision of the computer's arithmetic,
889 // 1 < 1 + R
890 // but
891 // 1 = ( 1 + R / 2 )
892 //
893 // Licensing:
894 //
895 // This code is distributed under the GNU LGPL license.
896 //
897 // Modified:
898 //
899 // 08 May 2006
900 //
901 // Author:
902 //
903 // John Burkardt
904 //
905 // Parameters:
906 //
907 // Output, double R8_EPSILON, the double precision round-off unit.
908 //
909 {
910  double r;
911 
912  r = 1.0;
913 
914  while ( 1.0 < ( double ) ( 1.0 + r ) )
915  {
916  r = r / 2.0;
917  }
918 
919  return ( 2.0 * r );
920 }
921 //****************************************************************************80
922 
923 double r8_max ( double x, double y )
924 
925 //****************************************************************************80
926 //
927 // Purpose:
928 //
929 // R8_MAX returns the maximum of two R8's.
930 //
931 // Licensing:
932 //
933 // This code is distributed under the GNU LGPL license.
934 //
935 // Modified:
936 //
937 // 18 August 2004
938 //
939 // Author:
940 //
941 // John Burkardt
942 //
943 // Parameters:
944 //
945 // Input, double X, Y, the quantities to compare.
946 //
947 // Output, double R8_MAX, the maximum of X and Y.
948 //
949 {
950  double value;
951 
952  if ( y < x )
953  {
954  value = x;
955  }
956  else
957  {
958  value = y;
959  }
960  return value;
961 }
962 //****************************************************************************80
963 
964 double r8_sign ( double x )
965 
966 //****************************************************************************80
967 //
968 // Purpose:
969 //
970 // R8_SIGN returns the sign of an R8.
971 //
972 // Licensing:
973 //
974 // This code is distributed under the GNU LGPL license.
975 //
976 // Modified:
977 //
978 // 18 October 2004
979 //
980 // Author:
981 //
982 // John Burkardt
983 //
984 // Parameters:
985 //
986 // Input, double X, the number whose sign is desired.
987 //
988 // Output, double R8_SIGN, the sign of X.
989 //
990 {
991  double value;
992 
993  if ( x < 0.0 )
994  {
995  value = -1.0;
996  }
997  else
998  {
999  value = 1.0;
1000  }
1001  return value;
1002 }
1003 //****************************************************************************80
1004 
1005 void timestamp ( )
1006 
1007 //****************************************************************************80
1008 //
1009 // Purpose:
1010 //
1011 // TIMESTAMP prints the current YMDHMS date as a time stamp.
1012 //
1013 // Example:
1014 //
1015 // 31 May 2001 09:45:54 AM
1016 //
1017 // Licensing:
1018 //
1019 // This code is distributed under the GNU LGPL license.
1020 //
1021 // Modified:
1022 //
1023 // 24 September 2003
1024 //
1025 // Author:
1026 //
1027 // John Burkardt
1028 //
1029 // Parameters:
1030 //
1031 // None
1032 //
1033 {
1034 # define TIME_SIZE 40
1035 
1036  static char time_buffer[TIME_SIZE];
1037  const struct tm *tm;
1038  time_t now;
1039 
1040  now = time ( NULL );
1041  tm = localtime ( &now );
1042 
1043  strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm );
1044 
1045  cout << time_buffer << "\n";
1046 
1047  return;
1048 # undef TIME_SIZE
1049 }
1050 //****************************************************************************80
1051 
1052 double zero ( double a, double b, double t, func_base& f )
1053 
1054 //****************************************************************************80
1055 //
1056 // Purpose:
1057 //
1058 // ZERO seeks the root of a function F(X) in an interval [A,B].
1059 //
1060 // Discussion:
1061 //
1062 // The interval [A,B] must be a change of sign interval for F.
1063 // That is, F(A) and F(B) must be of opposite signs. Then
1064 // assuming that F is continuous implies the existence of at least
1065 // one value C between A and B for which F(C) = 0.
1066 //
1067 // The location of the zero is determined to within an accuracy
1068 // of 6 * MACHEPS * r8_abs ( C ) + 2 * T.
1069 //
1070 // Licensing:
1071 //
1072 // This code is distributed under the GNU LGPL license.
1073 //
1074 // Modified:
1075 //
1076 // 06 May 2012
1077 //
1078 // Author:
1079 //
1080 // Original FORTRAN77 version by Richard Brent.
1081 // C++ version by John Burkardt.
1082 // Modifications by John Denker.
1083 //
1084 // Reference:
1085 //
1086 // Richard Brent,
1087 // Algorithms for Minimization Without Derivatives,
1088 // Dover, 2002,
1089 // ISBN: 0-486-41998-3,
1090 // LC: QA402.5.B74.
1091 //
1092 // Parameters:
1093 //
1094 // Input, double A, B, the endpoints of the change of sign interval.
1095 //
1096 // Input, double T, a positive error tolerance.
1097 //
1098 // Input, func_base& F, the name of a user-supplied c++ functor
1099 // whose zero is being sought. The input and output
1100 // of F() are of type double.
1101 //
1102 // Output, double ZERO, the estimated value of a zero of
1103 // the function F.
1104 //
1105 {
1106  double c;
1107  double d;
1108  double e;
1109  double fa;
1110  double fb;
1111  double fc;
1112  double m;
1113  double macheps;
1114  double p;
1115  double q;
1116  double r;
1117  double s;
1118  double sa;
1119  double sb;
1120  double tol;
1121 //
1122 // Make local copies of A and B.
1123 //
1124  sa = a;
1125  sb = b;
1126  fa = f ( sa );
1127  fb = f ( sb );
1128 
1129  c = sa;
1130  fc = fa;
1131  e = sb - sa;
1132  d = e;
1133 
1134  macheps = r8_epsilon ( );
1135 
1136  for ( ; ; )
1137  {
1138  if ( r8_abs ( fc ) < r8_abs ( fb ) )
1139  {
1140  sa = sb;
1141  sb = c;
1142  c = sa;
1143  fa = fb;
1144  fb = fc;
1145  fc = fa;
1146  }
1147 
1148  tol = 2.0 * macheps * r8_abs ( sb ) + t;
1149  m = 0.5 * ( c - sb );
1150 
1151  if ( r8_abs ( m ) <= tol || fb == 0.0 )
1152  {
1153  break;
1154  }
1155 
1156  if ( r8_abs ( e ) < tol || r8_abs ( fa ) <= r8_abs ( fb ) )
1157  {
1158  e = m;
1159  d = e;
1160  }
1161  else
1162  {
1163  s = fb / fa;
1164 
1165  if ( sa == c )
1166  {
1167  p = 2.0 * m * s;
1168  q = 1.0 - s;
1169  }
1170  else
1171  {
1172  q = fa / fc;
1173  r = fb / fc;
1174  p = s * ( 2.0 * m * a * ( q - r ) - ( sb - sa ) * ( r - 1.0 ) );
1175  q = ( q - 1.0 ) * ( r - 1.0 ) * ( s - 1.0 );
1176  }
1177 
1178  if ( 0.0 < p )
1179  {
1180  q = - q;
1181  }
1182  else
1183  {
1184  p = - p;
1185  }
1186 
1187  s = e;
1188  e = d;
1189 
1190  if ( 2.0 * p < 3.0 * m * q - r8_abs ( tol * q ) &&
1191  p < r8_abs ( 0.5 * s * q ) )
1192  {
1193  d = p / q;
1194  }
1195  else
1196  {
1197  e = m;
1198  d = e;
1199  }
1200  }
1201  sa = sb;
1202  fa = fb;
1203 
1204  if ( tol < r8_abs ( d ) )
1205  {
1206  sb = sb + d;
1207  }
1208  else if ( 0.0 < m )
1209  {
1210  sb = sb + tol;
1211  }
1212  else
1213  {
1214  sb = sb - tol;
1215  }
1216 
1217  fb = f ( sb );
1218 
1219  if ( ( 0.0 < fb && 0.0 < fc ) || ( fb <= 0.0 && fc <= 0.0 ) )
1220  {
1221  c = sa;
1222  fc = fa;
1223  e = sb - sa;
1224  d = e;
1225  }
1226  }
1227  return sb;
1228 }
1229 //****************************************************************************80
1230 
1231 void zero_rc ( double a, double b, double t, double &arg, int &status,
1232  double value )
1233 
1234 //****************************************************************************80
1235 //
1236 // Purpose:
1237 //
1238 // ZERO_RC seeks the root of a function F(X) using reverse communication.
1239 //
1240 // Discussion:
1241 //
1242 // The interval [A,B] must be a change of sign interval for F.
1243 // That is, F(A) and F(B) must be of opposite signs. Then
1244 // assuming that F is continuous implies the existence of at least
1245 // one value C between A and B for which F(C) = 0.
1246 //
1247 // The location of the zero is determined to within an accuracy
1248 // of 6 * MACHEPS * r8_abs ( C ) + 2 * T.
1249 //
1250 // The routine is a revised version of the Brent zero finder
1251 // algorithm, using reverse communication.
1252 //
1253 // Licensing:
1254 //
1255 // This code is distributed under the GNU LGPL license.
1256 //
1257 // Modified:
1258 //
1259 // 17 July 2011
1260 //
1261 // Author:
1262 //
1263 // John Burkardt
1264 //
1265 // Reference:
1266 //
1267 // Richard Brent,
1268 // Algorithms for Minimization Without Derivatives,
1269 // Dover, 2002,
1270 // ISBN: 0-486-41998-3,
1271 // LC: QA402.5.B74.
1272 //
1273 // Parameters:
1274 //
1275 // Input, double A, B, the endpoints of the change of sign interval.
1276 //
1277 // Input, double T, a positive error tolerance.
1278 //
1279 // Output, double &ARG, the currently considered point. The user
1280 // does not need to initialize this value. On return with STATUS positive,
1281 // the user is requested to evaluate the function at ARG, and return
1282 // the value in VALUE. On return with STATUS zero, ARG is the routine's
1283 // estimate for the function's zero.
1284 //
1285 // Input/output, int &STATUS, used to communicate between
1286 // the user and the routine. The user only sets STATUS to zero on the first
1287 // call, to indicate that this is a startup call. The routine returns STATUS
1288 // positive to request that the function be evaluated at ARG, or returns
1289 // STATUS as 0, to indicate that the iteration is complete and that
1290 // ARG is the estimated zero
1291 //
1292 // Input, double VALUE, the function value at ARG, as requested
1293 // by the routine on the previous call.
1294 //
1295 {
1296  static double c;
1297  static double d;
1298  static double e;
1299  static double fa;
1300  static double fb;
1301  static double fc;
1302  double m;
1303  static double macheps;
1304  double p;
1305  double q;
1306  double r;
1307  double s;
1308  static double sa;
1309  static double sb;
1310  double tol;
1311 //
1312 // Input STATUS = 0.
1313 // Initialize, request F(A).
1314 //
1315  if ( status == 0 )
1316  {
1317  macheps = r8_epsilon ( );
1318 
1319  sa = a;
1320  sb = b;
1321  e = sb - sa;
1322  d = e;
1323 
1324  status = 1;
1325  arg = a;
1326  return;
1327  }
1328 //
1329 // Input STATUS = 1.
1330 // Receive F(A), request F(B).
1331 //
1332  else if ( status == 1 )
1333  {
1334  fa = value;
1335  status = 2;
1336  arg = sb;
1337  return;
1338  }
1339 //
1340 // Input STATUS = 2
1341 // Receive F(B).
1342 //
1343  else if ( status == 2 )
1344  {
1345  fb = value;
1346 
1347  if ( 0.0 < fa * fb )
1348  {
1349  status = -1;
1350  return;
1351  }
1352  c = sa;
1353  fc = fa;
1354  }
1355  else
1356  {
1357  fb = value;
1358 
1359  if ( ( 0.0 < fb && 0.0 < fc ) || ( fb <= 0.0 && fc <= 0.0 ) )
1360  {
1361  c = sa;
1362  fc = fa;
1363  e = sb - sa;
1364  d = e;
1365  }
1366  }
1367 //
1368 // Compute the next point at which a function value is requested.
1369 //
1370  if ( r8_abs ( fc ) < r8_abs ( fb ) )
1371  {
1372  sa = sb;
1373  sb = c;
1374  c = sa;
1375  fa = fb;
1376  fb = fc;
1377  fc = fa;
1378  }
1379 
1380  tol = 2.0 * macheps * r8_abs ( sb ) + t;
1381  m = 0.5 * ( c - sb );
1382 
1383  if ( r8_abs ( m ) <= tol || fb == 0.0 )
1384  {
1385  status = 0;
1386  arg = sb;
1387  return;
1388  }
1389 
1390  if ( r8_abs ( e ) < tol || r8_abs ( fa ) <= r8_abs ( fb ) )
1391  {
1392  e = m;
1393  d = e;
1394  }
1395  else
1396  {
1397  s = fb / fa;
1398 
1399  if ( sa == c )
1400  {
1401  p = 2.0 * m * s;
1402  q = 1.0 - s;
1403  }
1404  else
1405  {
1406  q = fa / fc;
1407  r = fb / fc;
1408  p = s * ( 2.0 * m * a * ( q - r ) - ( sb - sa ) * ( r - 1.0 ) );
1409  q = ( q - 1.0 ) * ( r - 1.0 ) * ( s - 1.0 );
1410  }
1411 
1412  if ( 0.0 < p )
1413  {
1414  q = - q;
1415  }
1416  else
1417  {
1418  p = - p;
1419  }
1420  s = e;
1421  e = d;
1422 
1423  if ( 2.0 * p < 3.0 * m * q - r8_abs ( tol * q ) &&
1424  p < r8_abs ( 0.5 * s * q ) )
1425  {
1426  d = p / q;
1427  }
1428  else
1429  {
1430  e = m;
1431  d = e;
1432  }
1433  }
1434 
1435  sa = sb;
1436  fa = fb;
1437 
1438  if ( tol < r8_abs ( d ) )
1439  {
1440  sb = sb + d;
1441  }
1442  else if ( 0.0 < m )
1443  {
1444  sb = sb + tol;
1445  }
1446  else
1447  {
1448  sb = sb - tol;
1449  }
1450 
1451  arg = sb;
1452  status = status + 1;
1453 
1454  return;
1455 }
1456 
1457 // ======================================================================
1458 // === Simple wrapper functions
1459 // === for convenience and/or compatibility.
1460 //
1461 // === The three functions are the same as above,
1462 // === except that they take a plain function F
1463 // === instead of a c++ functor. In all cases, the
1464 // === input and output of F() are of type double.
1465 
1466 typedef double DoubleOfDouble (double);
1467 
1468 class func_wrapper : public func_base {
1469  DoubleOfDouble* func;
1470 public:
1472  func = f;
1473  }
1474  virtual double operator() (double x){
1475  return func(x);
1476  }
1477 };
1478 
1479 //****************************************************************************80
1480 
1481 double glomin ( double a, double b, double c, double m, double e,
1482  double t, double f ( double x ), double &x ){
1483  func_wrapper foo(f);
1484  return glomin(a, b, c, m, e, t, foo, x);
1485 }
1486 
1487 //****************************************************************************80
1488 
1489 double local_min ( double a, double b, double t, double f ( double x ),
1490  double &x ){
1491  func_wrapper foo(f);
1492  return local_min(a, b, t, foo, x);
1493 }
1494 
1495 //****************************************************************************80
1496 
1497 double zero ( double a, double b, double t, double f ( double x ) ){
1498  func_wrapper foo(f);
1499  return zero(a, b, t, foo);
1500 }
1501 
1502 // ======================================================================
1503 // Generally useful functor to evaluate a monic polynomial.
1504 // For details, see class definition in brent.hpp
1505 
1506 double monicPoly::operator()(double x){
1507  double rslt(1);
1508  for (int ii = coeff.size()-1; ii >= 0; ii--){
1509  rslt *= x;
1510  rslt += coeff[ii];
1511  }
1512  return rslt;
1513 }
1514 
1515 // Similarly, evaluate a general polynomial (not necessarily monic):
1516 double Poly::operator()(double x){
1517  double rslt(0);
1518  for (int ii = coeff.size()-1; ii >= 0; ii--){
1519  rslt *= x;
1520  rslt += coeff[ii];
1521  }
1522  return rslt;
1523 }
1524 
1525 }

SHOGUN Machine Learning Toolbox - Documentation