Next: , Previous: Algorithms without Derivatives, Up: Multidimensional Root-Finding


34.8 Examples

The multidimensional solvers are used in a similar way to the one-dimensional root finding algorithms. This first example demonstrates the hybrids scaled-hybrid algorithm, which does not require derivatives. The program solves the Rosenbrock system of equations,

     f_1 (x, y) = a (1 - x)
     f_2 (x, y) = b (y - x^2)

with a = 1, b = 10. The solution of this system lies at (x,y) = (1,1) in a narrow valley.

The first stage of the program is to define the system of equations,

     #include <stdlib.h>
     #include <stdio.h>
     #include <gsl/gsl_vector.h>
     #include <gsl/gsl_multiroots.h>
     
     struct rparams
       {
         double a;
         double b;
       };
     
     int
     rosenbrock_f (const gsl_vector * x, void *params,
                   gsl_vector * f)
     {
       double a = ((struct rparams *) params)->a;
       double b = ((struct rparams *) params)->b;
     
       const double x0 = gsl_vector_get (x, 0);
       const double x1 = gsl_vector_get (x, 1);
     
       const double y0 = a * (1 - x0);
       const double y1 = b * (x1 - x0 * x0);
     
       gsl_vector_set (f, 0, y0);
       gsl_vector_set (f, 1, y1);
     
       return GSL_SUCCESS;
     }

The main program begins by creating the function object f, with the arguments (x,y) and parameters (a,b). The solver s is initialized to use this function, with the hybrids method.

     int
     main (void)
     {
       const gsl_multiroot_fsolver_type *T;
       gsl_multiroot_fsolver *s;
     
       int status;
       size_t i, iter = 0;
     
       const size_t n = 2;
       struct rparams p = {1.0, 10.0};
       gsl_multiroot_function f = {&rosenbrock_f, n, &p};
     
       double x_init[2] = {-10.0, -5.0};
       gsl_vector *x = gsl_vector_alloc (n);
     
       gsl_vector_set (x, 0, x_init[0]);
       gsl_vector_set (x, 1, x_init[1]);
     
       T = gsl_multiroot_fsolver_hybrids;
       s = gsl_multiroot_fsolver_alloc (T, 2);
       gsl_multiroot_fsolver_set (s, &f, x);
     
       print_state (iter, s);
     
       do
         {
           iter++;
           status = gsl_multiroot_fsolver_iterate (s);
     
           print_state (iter, s);
     
           if (status)   /* check if solver is stuck */
             break;
     
           status =
             gsl_multiroot_test_residual (s->f, 1e-7);
         }
       while (status == GSL_CONTINUE && iter < 1000);
     
       printf ("status = %s\n", gsl_strerror (status));
     
       gsl_multiroot_fsolver_free (s);
       gsl_vector_free (x);
       return 0;
     }

Note that it is important to check the return status of each solver step, in case the algorithm becomes stuck. If an error condition is detected, indicating that the algorithm cannot proceed, then the error can be reported to the user, a new starting point chosen or a different algorithm used.

The intermediate state of the solution is displayed by the following function. The solver state contains the vector s->x which is the current position, and the vector s->f with corresponding function values.

     int
     print_state (size_t iter, gsl_multiroot_fsolver * s)
     {
       printf ("iter = %3u x = % .3f % .3f "
               "f(x) = % .3e % .3e\n",
               iter,
               gsl_vector_get (s->x, 0),
               gsl_vector_get (s->x, 1),
               gsl_vector_get (s->f, 0),
               gsl_vector_get (s->f, 1));
     }

Here are the results of running the program. The algorithm is started at (-10,-5) far from the solution. Since the solution is hidden in a narrow valley the earliest steps follow the gradient of the function downhill, in an attempt to reduce the large value of the residual. Once the root has been approximately located, on iteration 8, the Newton behavior takes over and convergence is very rapid.

     iter =  0 x = -10.000  -5.000  f(x) = 1.100e+01 -1.050e+03
     iter =  1 x = -10.000  -5.000  f(x) = 1.100e+01 -1.050e+03
     iter =  2 x =  -3.976  24.827  f(x) = 4.976e+00  9.020e+01
     iter =  3 x =  -3.976  24.827  f(x) = 4.976e+00  9.020e+01
     iter =  4 x =  -3.976  24.827  f(x) = 4.976e+00  9.020e+01
     iter =  5 x =  -1.274  -5.680  f(x) = 2.274e+00 -7.302e+01
     iter =  6 x =  -1.274  -5.680  f(x) = 2.274e+00 -7.302e+01
     iter =  7 x =   0.249   0.298  f(x) = 7.511e-01  2.359e+00
     iter =  8 x =   0.249   0.298  f(x) = 7.511e-01  2.359e+00
     iter =  9 x =   1.000   0.878  f(x) = 1.268e-10 -1.218e+00
     iter = 10 x =   1.000   0.989  f(x) = 1.124e-11 -1.080e-01
     iter = 11 x =   1.000   1.000  f(x) = 0.000e+00  0.000e+00
     status = success

Note that the algorithm does not update the location on every iteration. Some iterations are used to adjust the trust-region parameter, after trying a step which was found to be divergent, or to recompute the Jacobian, when poor convergence behavior is detected.

The next example program adds derivative information, in order to accelerate the solution. There are two derivative functions rosenbrock_df and rosenbrock_fdf. The latter computes both the function and its derivative simultaneously. This allows the optimization of any common terms. For simplicity we substitute calls to the separate f and df functions at this point in the code below.

     int
     rosenbrock_df (const gsl_vector * x, void *params,
                    gsl_matrix * J)
     {
       const double a = ((struct rparams *) params)->a;
       const double b = ((struct rparams *) params)->b;
     
       const double x0 = gsl_vector_get (x, 0);
     
       const double df00 = -a;
       const double df01 = 0;
       const double df10 = -2 * b  * x0;
       const double df11 = b;
     
       gsl_matrix_set (J, 0, 0, df00);
       gsl_matrix_set (J, 0, 1, df01);
       gsl_matrix_set (J, 1, 0, df10);
       gsl_matrix_set (J, 1, 1, df11);
     
       return GSL_SUCCESS;
     }
     
     int
     rosenbrock_fdf (const gsl_vector * x, void *params,
                     gsl_vector * f, gsl_matrix * J)
     {
       rosenbrock_f (x, params, f);
       rosenbrock_df (x, params, J);
     
       return GSL_SUCCESS;
     }

The main program now makes calls to the corresponding fdfsolver versions of the functions,

     int
     main (void)
     {
       const gsl_multiroot_fdfsolver_type *T;
       gsl_multiroot_fdfsolver *s;
     
       int status;
       size_t i, iter = 0;
     
       const size_t n = 2;
       struct rparams p = {1.0, 10.0};
       gsl_multiroot_function_fdf f = {&rosenbrock_f,
                                       &rosenbrock_df,
                                       &rosenbrock_fdf,
                                       n, &p};
     
       double x_init[2] = {-10.0, -5.0};
       gsl_vector *x = gsl_vector_alloc (n);
     
       gsl_vector_set (x, 0, x_init[0]);
       gsl_vector_set (x, 1, x_init[1]);
     
       T = gsl_multiroot_fdfsolver_gnewton;
       s = gsl_multiroot_fdfsolver_alloc (T, n);
       gsl_multiroot_fdfsolver_set (s, &f, x);
     
       print_state (iter, s);
     
       do
         {
           iter++;
     
           status = gsl_multiroot_fdfsolver_iterate (s);
     
           print_state (iter, s);
     
           if (status)
             break;
     
           status = gsl_multiroot_test_residual (s->f, 1e-7);
         }
       while (status == GSL_CONTINUE && iter < 1000);
     
       printf ("status = %s\n", gsl_strerror (status));
     
       gsl_multiroot_fdfsolver_free (s);
       gsl_vector_free (x);
       return 0;
     }

The addition of derivative information to the hybrids solver does not make any significant difference to its behavior, since it able to approximate the Jacobian numerically with sufficient accuracy. To illustrate the behavior of a different derivative solver we switch to gnewton. This is a traditional Newton solver with the constraint that it scales back its step if the full step would lead “uphill”. Here is the output for the gnewton algorithm,

     iter = 0 x = -10.000  -5.000 f(x) =  1.100e+01 -1.050e+03
     iter = 1 x =  -4.231 -65.317 f(x) =  5.231e+00 -8.321e+02
     iter = 2 x =   1.000 -26.358 f(x) = -8.882e-16 -2.736e+02
     iter = 3 x =   1.000   1.000 f(x) = -2.220e-16 -4.441e-15
     status = success

The convergence is much more rapid, but takes a wide excursion out to the point (-4.23,-65.3). This could cause the algorithm to go astray in a realistic application. The hybrid algorithm follows the downhill path to the solution more reliably.