Added script front-end for primer-design code
[htsworkflow.git] / htswanalysis / MACS / lib / gsl / gsl-1.11 / doc / multiroots.texi
1 @cindex solving nonlinear systems of equations
2 @cindex nonlinear systems of equations, solution of
3 @cindex systems of equations, nonlinear
4
5 This chapter describes functions for multidimensional root-finding
6 (solving nonlinear systems with @math{n} equations in @math{n}
7 unknowns).  The library provides low level components for a variety of
8 iterative solvers and convergence tests.  These can be combined by the
9 user to achieve the desired solution, with full access to the
10 intermediate steps of the iteration.  Each class of methods uses the
11 same framework, so that you can switch between solvers at runtime
12 without needing to recompile your program.  Each instance of a solver
13 keeps track of its own state, allowing the solvers to be used in
14 multi-threaded programs.  The solvers are based on the original Fortran
15 library @sc{minpack}.
16
17 The header file @file{gsl_multiroots.h} contains prototypes for the
18 multidimensional root finding functions and related declarations.
19
20 @menu
21 * Overview of Multidimensional Root Finding::  
22 * Initializing the Multidimensional Solver::  
23 * Providing the multidimensional system of equations to solve::  
24 * Iteration of the multidimensional solver::  
25 * Search Stopping Parameters for the multidimensional solver::  
26 * Algorithms using Derivatives::  
27 * Algorithms without Derivatives::  
28 * Example programs for Multidimensional Root finding::  
29 * References and Further Reading for Multidimensional Root Finding::  
30 @end menu
31
32 @node Overview of Multidimensional Root Finding
33 @section Overview
34 @cindex multidimensional root finding, overview
35
36 The problem of multidimensional root finding requires the simultaneous
37 solution of @math{n} equations, @math{f_i}, in @math{n} variables,
38 @math{x_i},
39 @tex
40 \beforedisplay
41 $$
42 f_i (x_1, \dots, x_n) = 0 \qquad\hbox{for}~i = 1 \dots n.
43 $$
44 \afterdisplay
45 @end tex
46 @ifinfo
47
48 @example
49 f_i (x_1, ..., x_n) = 0    for i = 1 ... n.
50 @end example
51
52 @end ifinfo
53 @noindent
54 In general there are no bracketing methods available for @math{n}
55 dimensional systems, and no way of knowing whether any solutions
56 exist.  All algorithms proceed from an initial guess using a variant of
57 the Newton iteration,
58 @tex
59 \beforedisplay
60 $$
61 x \to x' = x - J^{-1} f(x)
62 $$
63 \afterdisplay
64 @end tex
65 @ifinfo
66
67 @example
68 x -> x' = x - J^@{-1@} f(x)
69 @end example
70
71 @end ifinfo
72 @noindent
73 where @math{x}, @math{f} are vector quantities and @math{J} is the
74 Jacobian matrix @c{$J_{ij} = \partial f_i / \partial x_j$} 
75 @math{J_@{ij@} = d f_i / d x_j}.
76 Additional strategies can be used to enlarge the region of
77 convergence.  These include requiring a decrease in the norm @math{|f|} on
78 each step proposed by Newton's method, or taking steepest-descent steps in
79 the direction of the negative gradient of @math{|f|}.
80
81 Several root-finding algorithms are available within a single framework.
82 The user provides a high-level driver for the algorithms, and the
83 library provides the individual functions necessary for each of the
84 steps.  There are three main phases of the iteration.  The steps are,
85
86 @itemize @bullet
87 @item
88 initialize solver state, @var{s}, for algorithm @var{T}
89
90 @item
91 update @var{s} using the iteration @var{T}
92
93 @item
94 test @var{s} for convergence, and repeat iteration if necessary
95 @end itemize
96
97 @noindent
98 The evaluation of the Jacobian matrix can be problematic, either because
99 programming the derivatives is intractable or because computation of the
100 @math{n^2} terms of the matrix becomes too expensive.  For these reasons
101 the algorithms provided by the library are divided into two classes according
102 to whether the derivatives are available or not.
103
104 The state for solvers with an analytic Jacobian matrix is held in a
105 @code{gsl_multiroot_fdfsolver} struct.  The updating procedure requires
106 both the function and its derivatives to be supplied by the user.
107
108 The state for solvers which do not use an analytic Jacobian matrix is
109 held in a @code{gsl_multiroot_fsolver} struct.  The updating procedure
110 uses only function evaluations (not derivatives).  The algorithms
111 estimate the matrix @math{J} or @c{$J^{-1}$} 
112 @math{J^@{-1@}} by approximate methods.
113
114 @node Initializing the Multidimensional Solver
115 @section Initializing the Solver
116
117 The following functions initialize a multidimensional solver, either
118 with or without derivatives.  The solver itself depends only on the
119 dimension of the problem and the algorithm and can be reused for
120 different problems.
121
122 @deftypefun {gsl_multiroot_fsolver *} gsl_multiroot_fsolver_alloc (const gsl_multiroot_fsolver_type * @var{T}, size_t @var{n})
123 This function returns a pointer to a newly allocated instance of a
124 solver of type @var{T} for a system of @var{n} dimensions.
125 For example, the following code creates an instance of a hybrid solver, 
126 to solve a 3-dimensional system of equations.
127
128 @example
129 const gsl_multiroot_fsolver_type * T 
130     = gsl_multiroot_fsolver_hybrid;
131 gsl_multiroot_fsolver * s 
132     = gsl_multiroot_fsolver_alloc (T, 3);
133 @end example
134
135 @noindent
136 If there is insufficient memory to create the solver then the function
137 returns a null pointer and the error handler is invoked with an error
138 code of @code{GSL_ENOMEM}.
139 @end deftypefun
140
141 @deftypefun {gsl_multiroot_fdfsolver *} gsl_multiroot_fdfsolver_alloc (const gsl_multiroot_fdfsolver_type * @var{T}, size_t @var{n})
142 This function returns a pointer to a newly allocated instance of a
143 derivative solver of type @var{T} for a system of @var{n} dimensions.
144 For example, the following code creates an instance of a Newton-Raphson solver,
145 for a 2-dimensional system of equations.
146
147 @example
148 const gsl_multiroot_fdfsolver_type * T 
149     = gsl_multiroot_fdfsolver_newton;
150 gsl_multiroot_fdfsolver * s = 
151     gsl_multiroot_fdfsolver_alloc (T, 2);
152 @end example
153
154 @noindent
155 If there is insufficient memory to create the solver then the function
156 returns a null pointer and the error handler is invoked with an error
157 code of @code{GSL_ENOMEM}.
158 @end deftypefun
159
160 @deftypefun int gsl_multiroot_fsolver_set (gsl_multiroot_fsolver * @var{s}, gsl_multiroot_function * @var{f}, const gsl_vector * @var{x})
161 @deftypefunx int gsl_multiroot_fdfsolver_set (gsl_multiroot_fdfsolver * @var{s}, gsl_multiroot_function_fdf * @var{fdf}, const gsl_vector * @var{x})
162 These functions set, or reset, an existing solver @var{s} to use the
163 function @var{f} or function and derivative @var{fdf}, and the initial
164 guess @var{x}.  Note that the initial position is copied from @var{x}, this
165 argument is not modified by subsequent iterations.
166 @end deftypefun
167
168 @deftypefun void gsl_multiroot_fsolver_free (gsl_multiroot_fsolver * @var{s})
169 @deftypefunx void gsl_multiroot_fdfsolver_free (gsl_multiroot_fdfsolver * @var{s})
170 These functions free all the memory associated with the solver @var{s}.
171 @end deftypefun
172
173 @deftypefun {const char *} gsl_multiroot_fsolver_name (const gsl_multiroot_fsolver * @var{s})
174 @deftypefunx {const char *} gsl_multiroot_fdfsolver_name (const gsl_multiroot_fdfsolver * @var{s})
175 These functions return a pointer to the name of the solver.  For example,
176
177 @example
178 printf ("s is a '%s' solver\n", 
179         gsl_multiroot_fdfsolver_name (s));
180 @end example
181
182 @noindent
183 would print something like @code{s is a 'newton' solver}.
184 @end deftypefun
185
186 @node Providing the multidimensional system of equations to solve
187 @section Providing the function to solve
188 @cindex multidimensional root finding, providing a function to solve
189
190 You must provide @math{n} functions of @math{n} variables for the root
191 finders to operate on.  In order to allow for general parameters the
192 functions are defined by the following data types:
193
194 @deftp {Data Type} gsl_multiroot_function 
195 This data type defines a general system of functions with parameters.
196
197 @table @code
198 @item int (* f) (const gsl_vector * @var{x}, void * @var{params}, gsl_vector * @var{f})
199 this function should store the vector result
200 @c{$f(x,\hbox{\it params})$}
201 @math{f(x,params)} in @var{f} for argument @var{x} and parameters @var{params},
202 returning an appropriate error code if the function cannot be computed.
203
204 @item size_t n
205 the dimension of the system, i.e. the number of components of the
206 vectors @var{x} and @var{f}.
207
208 @item void * params
209 a pointer to the parameters of the function.
210 @end table
211 @end deftp
212
213 @noindent
214 Here is an example using Powell's test function,
215 @tex
216 \beforedisplay
217 $$
218 f_1(x) = A x_0 x_1 - 1,
219 f_2(x) = \exp(-x_0) + \exp(-x_1) - (1 + 1/A)
220 $$
221 \afterdisplay
222 @end tex
223 @ifinfo
224
225 @example
226 f_1(x) = A x_0 x_1 - 1,
227 f_2(x) = exp(-x_0) + exp(-x_1) - (1 + 1/A)
228 @end example
229
230 @end ifinfo
231 @noindent
232 with @math{A = 10^4}.  The following code defines a
233 @code{gsl_multiroot_function} system @code{F} which you could pass to a
234 solver:
235
236 @example
237 struct powell_params @{ double A; @};
238
239 int
240 powell (gsl_vector * x, void * p, gsl_vector * f) @{
241    struct powell_params * params 
242      = *(struct powell_params *)p;
243    const double A = (params->A);
244    const double x0 = gsl_vector_get(x,0);
245    const double x1 = gsl_vector_get(x,1);
246
247    gsl_vector_set (f, 0, A * x0 * x1 - 1);
248    gsl_vector_set (f, 1, (exp(-x0) + exp(-x1) 
249                           - (1.0 + 1.0/A)));
250    return GSL_SUCCESS
251 @}
252
253 gsl_multiroot_function F;
254 struct powell_params params = @{ 10000.0 @};
255
256 F.f = &powell;
257 F.n = 2;
258 F.params = &params;
259 @end example
260
261 @deftp {Data Type} gsl_multiroot_function_fdf
262 This data type defines a general system of functions with parameters and
263 the corresponding Jacobian matrix of derivatives,
264
265 @table @code
266 @item int (* f) (const gsl_vector * @var{x}, void * @var{params}, gsl_vector * @var{f})
267 this function should store the vector result
268 @c{$f(x,\hbox{\it params})$}
269 @math{f(x,params)} in @var{f} for argument @var{x} and parameters @var{params},
270 returning an appropriate error code if the function cannot be computed.
271
272 @item int (* df) (const gsl_vector * @var{x}, void * @var{params}, gsl_matrix * @var{J})
273 this function should store the @var{n}-by-@var{n} matrix result
274 @c{$J_{ij} = \partial f_i(x,\hbox{\it params}) / \partial x_j$}
275 @math{J_ij = d f_i(x,params) / d x_j} in @var{J} for argument @var{x} 
276 and parameters @var{params}, returning an appropriate error code if the
277 function cannot be computed.
278
279 @item int (* fdf) (const gsl_vector * @var{x}, void * @var{params}, gsl_vector * @var{f}, gsl_matrix * @var{J})
280 This function should set the values of the @var{f} and @var{J} as above,
281 for arguments @var{x} and parameters @var{params}.  This function
282 provides an optimization of the separate functions for @math{f(x)} and
283 @math{J(x)}---it is always faster to compute the function and its
284 derivative at the same time.
285
286 @item size_t n
287 the dimension of the system, i.e. the number of components of the
288 vectors @var{x} and @var{f}.
289
290 @item void * params
291 a pointer to the parameters of the function.
292 @end table
293 @end deftp
294
295 @noindent
296 The example of Powell's test function defined above can be extended to
297 include analytic derivatives using the following code,
298
299 @example
300 int
301 powell_df (gsl_vector * x, void * p, gsl_matrix * J) 
302 @{
303    struct powell_params * params 
304      = *(struct powell_params *)p;
305    const double A = (params->A);
306    const double x0 = gsl_vector_get(x,0);
307    const double x1 = gsl_vector_get(x,1);
308    gsl_matrix_set (J, 0, 0, A * x1);
309    gsl_matrix_set (J, 0, 1, A * x0);
310    gsl_matrix_set (J, 1, 0, -exp(-x0));
311    gsl_matrix_set (J, 1, 1, -exp(-x1));
312    return GSL_SUCCESS
313 @}
314
315 int
316 powell_fdf (gsl_vector * x, void * p, 
317             gsl_matrix * f, gsl_matrix * J) @{
318    struct powell_params * params 
319      = *(struct powell_params *)p;
320    const double A = (params->A);
321    const double x0 = gsl_vector_get(x,0);
322    const double x1 = gsl_vector_get(x,1);
323
324    const double u0 = exp(-x0);
325    const double u1 = exp(-x1);
326
327    gsl_vector_set (f, 0, A * x0 * x1 - 1);
328    gsl_vector_set (f, 1, u0 + u1 - (1 + 1/A));
329
330    gsl_matrix_set (J, 0, 0, A * x1);
331    gsl_matrix_set (J, 0, 1, A * x0);
332    gsl_matrix_set (J, 1, 0, -u0);
333    gsl_matrix_set (J, 1, 1, -u1);
334    return GSL_SUCCESS
335 @}
336
337 gsl_multiroot_function_fdf FDF;
338
339 FDF.f = &powell_f;
340 FDF.df = &powell_df;
341 FDF.fdf = &powell_fdf;
342 FDF.n = 2;
343 FDF.params = 0;
344 @end example
345
346 @noindent
347 Note that the function @code{powell_fdf} is able to reuse existing terms
348 from the function when calculating the Jacobian, thus saving time.
349
350 @node Iteration of the multidimensional solver
351 @section Iteration
352
353 The following functions drive the iteration of each algorithm.  Each
354 function performs one iteration to update the state of any solver of the
355 corresponding type.  The same functions work for all solvers so that
356 different methods can be substituted at runtime without modifications to
357 the code.
358
359 @deftypefun int gsl_multiroot_fsolver_iterate (gsl_multiroot_fsolver * @var{s})
360 @deftypefunx int gsl_multiroot_fdfsolver_iterate (gsl_multiroot_fdfsolver * @var{s})
361 These functions perform a single iteration of the solver @var{s}.  If the
362 iteration encounters an unexpected problem then an error code will be
363 returned,
364
365 @table @code
366 @item GSL_EBADFUNC
367 the iteration encountered a singular point where the function or its
368 derivative evaluated to @code{Inf} or @code{NaN}.
369
370 @item GSL_ENOPROG
371 the iteration is not making any progress, preventing the algorithm from
372 continuing.
373 @end table
374 @end deftypefun
375
376 The solver maintains a current best estimate of the root @code{s->x}
377 and its function value @code{s->f} at all times.  This information can
378 be accessed with the following auxiliary functions,
379
380 @deftypefun {gsl_vector *} gsl_multiroot_fsolver_root (const gsl_multiroot_fsolver * @var{s})
381 @deftypefunx {gsl_vector *} gsl_multiroot_fdfsolver_root (const gsl_multiroot_fdfsolver * @var{s})
382 These functions return the current estimate of the root for the solver @var{s}, given by @code{s->x}.
383 @end deftypefun
384
385 @deftypefun {gsl_vector *} gsl_multiroot_fsolver_f (const gsl_multiroot_fsolver * @var{s})
386 @deftypefunx {gsl_vector *} gsl_multiroot_fdfsolver_f (const gsl_multiroot_fdfsolver * @var{s})
387 These functions return the function value @math{f(x)} at the current
388 estimate of the root for the solver @var{s}, given by @code{s->f}.
389 @end deftypefun
390
391 @deftypefun {gsl_vector *} gsl_multiroot_fsolver_dx (const gsl_multiroot_fsolver * @var{s})
392 @deftypefunx {gsl_vector *} gsl_multiroot_fdfsolver_dx (const gsl_multiroot_fdfsolver * @var{s})
393 These functions return the last step @math{dx} taken by the solver
394 @var{s}, given by @code{s->dx}.
395 @end deftypefun
396
397 @node Search Stopping Parameters for the multidimensional solver
398 @section Search Stopping Parameters
399 @cindex root finding, stopping parameters
400
401 A root finding procedure should stop when one of the following conditions is
402 true:
403
404 @itemize @bullet
405 @item
406 A multidimensional root has been found to within the user-specified precision.
407
408 @item
409 A user-specified maximum number of iterations has been reached.
410
411 @item
412 An error has occurred.
413 @end itemize
414
415 @noindent
416 The handling of these conditions is under user control.  The functions
417 below allow the user to test the precision of the current result in
418 several standard ways.
419
420 @deftypefun int gsl_multiroot_test_delta (const gsl_vector * @var{dx}, const gsl_vector * @var{x}, double @var{epsabs}, double @var{epsrel})
421
422 This function tests for the convergence of the sequence by comparing the
423 last step @var{dx} with the absolute error @var{epsabs} and relative
424 error @var{epsrel} to the current position @var{x}.  The test returns
425 @code{GSL_SUCCESS} if the following condition is achieved,
426 @tex
427 \beforedisplay
428 $$
429 |dx_i| < \hbox{\it epsabs} + \hbox{\it epsrel\/}\, |x_i|
430 $$
431 \afterdisplay
432 @end tex
433 @ifinfo
434
435 @example
436 |dx_i| < epsabs + epsrel |x_i|
437 @end example
438
439 @end ifinfo
440 @noindent
441 for each component of @var{x} and returns @code{GSL_CONTINUE} otherwise.
442 @end deftypefun
443
444 @cindex residual, in nonlinear systems of equations
445 @deftypefun int gsl_multiroot_test_residual (const gsl_vector * @var{f}, double @var{epsabs})
446 This function tests the residual value @var{f} against the absolute
447 error bound @var{epsabs}.  The test returns @code{GSL_SUCCESS} if the
448 following condition is achieved,
449 @tex
450 \beforedisplay
451 $$
452 \sum_i |f_i| < \hbox{\it epsabs}
453 $$
454 \afterdisplay
455 @end tex
456 @ifinfo
457
458 @example
459 \sum_i |f_i| < epsabs
460 @end example
461
462 @end ifinfo
463 @noindent
464 and returns @code{GSL_CONTINUE} otherwise.  This criterion is suitable
465 for situations where the precise location of the root, @math{x}, is
466 unimportant provided a value can be found where the residual is small
467 enough.
468 @end deftypefun
469
470 @comment ============================================================
471
472 @node Algorithms using Derivatives
473 @section Algorithms using Derivatives
474
475 The root finding algorithms described in this section make use of both
476 the function and its derivative.  They require an initial guess for the
477 location of the root, but there is no absolute guarantee of
478 convergence---the function must be suitable for this technique and the
479 initial guess must be sufficiently close to the root for it to work.
480 When the conditions are satisfied then convergence is quadratic.
481
482
483 @comment ============================================================
484 @cindex HYBRID algorithms for nonlinear systems
485 @deffn {Derivative Solver} gsl_multiroot_fdfsolver_hybridsj
486 @cindex HYBRIDSJ algorithm
487 @cindex MINPACK, minimization algorithms
488 This is a modified version of Powell's Hybrid method as implemented in
489 the @sc{hybrj} algorithm in @sc{minpack}.  Minpack was written by Jorge
490 J. Mor@'e, Burton S. Garbow and Kenneth E. Hillstrom.  The Hybrid
491 algorithm retains the fast convergence of Newton's method but will also
492 reduce the residual when Newton's method is unreliable. 
493
494 The algorithm uses a generalized trust region to keep each step under
495 control.  In order to be accepted a proposed new position @math{x'} must
496 satisfy the condition @math{|D (x' - x)| < \delta}, where @math{D} is a
497 diagonal scaling matrix and @math{\delta} is the size of the trust
498 region.  The components of @math{D} are computed internally, using the
499 column norms of the Jacobian to estimate the sensitivity of the residual
500 to each component of @math{x}.  This improves the behavior of the
501 algorithm for badly scaled functions.
502
503 On each iteration the algorithm first determines the standard Newton
504 step by solving the system @math{J dx = - f}.  If this step falls inside
505 the trust region it is used as a trial step in the next stage.  If not,
506 the algorithm uses the linear combination of the Newton and gradient
507 directions which is predicted to minimize the norm of the function while
508 staying inside the trust region,
509 @tex
510 \beforedisplay
511 $$
512 dx = - \alpha J^{-1} f(x) - \beta \nabla |f(x)|^2.
513 $$
514 \afterdisplay
515 @end tex
516 @ifinfo
517
518 @example
519 dx = - \alpha J^@{-1@} f(x) - \beta \nabla |f(x)|^2.
520 @end example
521
522 @end ifinfo
523 @noindent
524 This combination of Newton and gradient directions is referred to as a
525 @dfn{dogleg step}.
526
527 The proposed step is now tested by evaluating the function at the
528 resulting point, @math{x'}.  If the step reduces the norm of the function
529 sufficiently then it is accepted and size of the trust region is
530 increased.  If the proposed step fails to improve the solution then the
531 size of the trust region is decreased and another trial step is
532 computed.
533
534 The speed of the algorithm is increased by computing the changes to the
535 Jacobian approximately, using a rank-1 update.  If two successive
536 attempts fail to reduce the residual then the full Jacobian is
537 recomputed.  The algorithm also monitors the progress of the solution
538 and returns an error if several steps fail to make any improvement,
539
540 @table @code
541 @item GSL_ENOPROG
542 the iteration is not making any progress, preventing the algorithm from
543 continuing.
544
545 @item GSL_ENOPROGJ
546 re-evaluations of the Jacobian indicate that the iteration is not
547 making any progress, preventing the algorithm from continuing.
548 @end table
549
550 @end deffn
551
552 @deffn {Derivative Solver} gsl_multiroot_fdfsolver_hybridj
553 @cindex HYBRIDJ algorithm
554 This algorithm is an unscaled version of @code{hybridsj}.  The steps are
555 controlled by a spherical trust region @math{|x' - x| < \delta}, instead
556 of a generalized region.  This can be useful if the generalized region
557 estimated by @code{hybridsj} is inappropriate.
558 @end deffn
559
560
561 @deffn {Derivative Solver} gsl_multiroot_fdfsolver_newton
562 @cindex Newton's method for systems of nonlinear equations
563
564 Newton's Method is the standard root-polishing algorithm.  The algorithm
565 begins with an initial guess for the location of the solution.  On each
566 iteration a linear approximation to the function @math{F} is used to
567 estimate the step which will zero all the components of the residual.
568 The iteration is defined by the following sequence,
569 @tex
570 \beforedisplay
571 $$
572 x \to x' = x - J^{-1} f(x)
573 $$
574 \afterdisplay
575 @end tex
576 @ifinfo
577
578 @example
579 x -> x' = x - J^@{-1@} f(x)
580 @end example
581
582 @end ifinfo
583 @noindent
584 where the Jacobian matrix @math{J} is computed from the derivative
585 functions provided by @var{f}.  The step @math{dx} is obtained by solving
586 the linear system,
587 @tex
588 \beforedisplay
589 $$
590 J \,dx = - f(x)
591 $$
592 \afterdisplay
593 @end tex
594 @ifinfo
595
596 @example
597 J dx = - f(x)
598 @end example
599
600 @end ifinfo
601 @noindent
602 using LU decomposition.
603 @end deffn
604
605 @comment ============================================================
606
607 @deffn {Derivative Solver} gsl_multiroot_fdfsolver_gnewton
608 @cindex Modified Newton's method for nonlinear systems
609 @cindex Newton algorithm, globally convergent
610 This is a modified version of Newton's method which attempts to improve
611 global convergence by requiring every step to reduce the Euclidean norm
612 of the residual, @math{|f(x)|}.  If the Newton step leads to an increase
613 in the norm then a reduced step of relative size,
614 @tex
615 \beforedisplay
616 $$
617 t = (\sqrt{1 + 6 r} - 1) / (3 r)
618 $$
619 \afterdisplay
620 @end tex
621 @ifinfo
622
623 @example
624 t = (\sqrt(1 + 6 r) - 1) / (3 r)
625 @end example
626
627 @end ifinfo
628 @noindent
629 is proposed, with @math{r} being the ratio of norms
630 @math{|f(x')|^2/|f(x)|^2}.  This procedure is repeated until a suitable step
631 size is found. 
632 @end deffn
633
634 @comment ============================================================
635
636 @node Algorithms without Derivatives
637 @section Algorithms without Derivatives
638
639 The algorithms described in this section do not require any derivative
640 information to be supplied by the user.  Any derivatives needed are
641 approximated by finite differences.  Note that if the
642 finite-differencing step size chosen by these routines is inappropriate,
643 an explicit user-supplied numerical derivative can always be used with
644 the algorithms described in the previous section.
645
646 @deffn {Solver} gsl_multiroot_fsolver_hybrids
647 @cindex HYBRIDS algorithm, scaled without derivatives
648 This is a version of the Hybrid algorithm which replaces calls to the
649 Jacobian function by its finite difference approximation.  The finite
650 difference approximation is computed using @code{gsl_multiroots_fdjac}
651 with a relative step size of @code{GSL_SQRT_DBL_EPSILON}.  Note that
652 this step size will not be suitable for all problems.
653 @end deffn
654
655 @deffn {Solver} gsl_multiroot_fsolver_hybrid
656 @cindex HYBRID algorithm, unscaled without derivatives
657 This is a finite difference version of the Hybrid algorithm without
658 internal scaling.
659 @end deffn
660
661 @comment ============================================================
662
663 @deffn {Solver} gsl_multiroot_fsolver_dnewton
664
665 @cindex Discrete Newton algorithm for multidimensional roots
666 @cindex Newton algorithm, discrete
667
668 The @dfn{discrete Newton algorithm} is the simplest method of solving a
669 multidimensional system.  It uses the Newton iteration
670 @tex
671 \beforedisplay
672 $$
673 x \to x - J^{-1} f(x)
674 $$
675 \afterdisplay
676 @end tex
677 @ifinfo
678
679 @example
680 x -> x - J^@{-1@} f(x)
681 @end example
682
683 @end ifinfo
684 @noindent
685 where the Jacobian matrix @math{J} is approximated by taking finite
686 differences of the function @var{f}.  The approximation scheme used by
687 this implementation is,
688 @tex
689 \beforedisplay
690 $$
691 J_{ij} = (f_i(x + \delta_j) - f_i(x)) /  \delta_j
692 $$
693 \afterdisplay
694 @end tex
695 @ifinfo
696
697 @example
698 J_@{ij@} = (f_i(x + \delta_j) - f_i(x)) /  \delta_j
699 @end example
700
701 @end ifinfo
702 @noindent
703 where @math{\delta_j} is a step of size @math{\sqrt\epsilon |x_j|} with
704 @math{\epsilon} being the machine precision 
705 (@c{$\epsilon \approx 2.22 \times 10^{-16}$}
706 @math{\epsilon \approx 2.22 \times 10^-16}).
707 The order of convergence of Newton's algorithm is quadratic, but the
708 finite differences require @math{n^2} function evaluations on each
709 iteration.  The algorithm may become unstable if the finite differences
710 are not a good approximation to the true derivatives.
711 @end deffn
712
713 @comment ============================================================
714
715 @deffn {Solver} gsl_multiroot_fsolver_broyden
716 @cindex Broyden algorithm for multidimensional roots
717 @cindex multidimensional root finding, Broyden algorithm
718
719 The @dfn{Broyden algorithm} is a version of the discrete Newton
720 algorithm which attempts to avoids the expensive update of the Jacobian
721 matrix on each iteration.  The changes to the Jacobian are also
722 approximated, using a rank-1 update,
723 @tex
724 \beforedisplay
725 $$
726 J^{-1} \to J^{-1} - (J^{-1} df - dx) dx^T J^{-1} / dx^T J^{-1} df
727 $$
728 \afterdisplay
729 @end tex
730 @ifinfo
731
732 @example
733 J^@{-1@} \to J^@{-1@} - (J^@{-1@} df - dx) dx^T J^@{-1@} / dx^T J^@{-1@} df
734 @end example
735
736 @end ifinfo
737 @noindent
738 where the vectors @math{dx} and @math{df} are the changes in @math{x}
739 and @math{f}.  On the first iteration the inverse Jacobian is estimated
740 using finite differences, as in the discrete Newton algorithm.
741  
742 This approximation gives a fast update but is unreliable if the changes
743 are not small, and the estimate of the inverse Jacobian becomes worse as
744 time passes.  The algorithm has a tendency to become unstable unless it
745 starts close to the root.  The Jacobian is refreshed if this instability
746 is detected (consult the source for details).
747
748 This algorithm is included only for demonstration purposes, and is not
749 recommended for serious use.
750 @end deffn
751
752 @comment ============================================================
753
754
755 @node Example programs for Multidimensional Root finding
756 @section Examples
757
758 The multidimensional solvers are used in a similar way to the
759 one-dimensional root finding algorithms.  This first example
760 demonstrates the @code{hybrids} scaled-hybrid algorithm, which does not
761 require derivatives. The program solves the Rosenbrock system of equations,
762 @tex
763 \beforedisplay
764 $$
765 f_1 (x, y) = a (1 - x),~
766 f_2 (x, y) = b (y - x^2)
767 $$
768 \afterdisplay
769 @end tex
770 @ifinfo
771
772 @example
773 f_1 (x, y) = a (1 - x)
774 f_2 (x, y) = b (y - x^2)
775 @end example
776
777 @end ifinfo
778 @noindent
779 with @math{a = 1, b = 10}. The solution of this system lies at
780 @math{(x,y) = (1,1)} in a narrow valley.
781
782 The first stage of the program is to define the system of equations,
783
784 @example
785 #include <stdlib.h>
786 #include <stdio.h>
787 #include <gsl/gsl_vector.h>
788 #include <gsl/gsl_multiroots.h>
789
790 struct rparams
791   @{
792     double a;
793     double b;
794   @};
795
796 int
797 rosenbrock_f (const gsl_vector * x, void *params, 
798               gsl_vector * f)
799 @{
800   double a = ((struct rparams *) params)->a;
801   double b = ((struct rparams *) params)->b;
802
803   const double x0 = gsl_vector_get (x, 0);
804   const double x1 = gsl_vector_get (x, 1);
805
806   const double y0 = a * (1 - x0);
807   const double y1 = b * (x1 - x0 * x0);
808
809   gsl_vector_set (f, 0, y0);
810   gsl_vector_set (f, 1, y1);
811
812   return GSL_SUCCESS;
813 @}
814 @end example
815
816 @noindent
817 The main program begins by creating the function object @code{f}, with
818 the arguments @code{(x,y)} and parameters @code{(a,b)}. The solver
819 @code{s} is initialized to use this function, with the @code{hybrids}
820 method.
821
822 @example
823 int
824 main (void)
825 @{
826   const gsl_multiroot_fsolver_type *T;
827   gsl_multiroot_fsolver *s;
828
829   int status;
830   size_t i, iter = 0;
831
832   const size_t n = 2;
833   struct rparams p = @{1.0, 10.0@};
834   gsl_multiroot_function f = @{&rosenbrock_f, n, &p@};
835
836   double x_init[2] = @{-10.0, -5.0@};
837   gsl_vector *x = gsl_vector_alloc (n);
838
839   gsl_vector_set (x, 0, x_init[0]);
840   gsl_vector_set (x, 1, x_init[1]);
841
842   T = gsl_multiroot_fsolver_hybrids;
843   s = gsl_multiroot_fsolver_alloc (T, 2);
844   gsl_multiroot_fsolver_set (s, &f, x);
845
846   print_state (iter, s);
847
848   do
849     @{
850       iter++;
851       status = gsl_multiroot_fsolver_iterate (s);
852
853       print_state (iter, s);
854
855       if (status)   /* check if solver is stuck */
856         break;
857
858       status = 
859         gsl_multiroot_test_residual (s->f, 1e-7);
860     @}
861   while (status == GSL_CONTINUE && iter < 1000);
862
863   printf ("status = %s\n", gsl_strerror (status));
864
865   gsl_multiroot_fsolver_free (s);
866   gsl_vector_free (x);
867   return 0;
868 @}
869 @end example
870
871 @noindent
872 Note that it is important to check the return status of each solver
873 step, in case the algorithm becomes stuck.  If an error condition is
874 detected, indicating that the algorithm cannot proceed, then the error
875 can be reported to the user, a new starting point chosen or a different
876 algorithm used.
877
878 The intermediate state of the solution is displayed by the following
879 function.  The solver state contains the vector @code{s->x} which is the
880 current position, and the vector @code{s->f} with corresponding function
881 values.
882
883 @example
884 int
885 print_state (size_t iter, gsl_multiroot_fsolver * s)
886 @{
887   printf ("iter = %3u x = % .3f % .3f "
888           "f(x) = % .3e % .3e\n",
889           iter,
890           gsl_vector_get (s->x, 0), 
891           gsl_vector_get (s->x, 1),
892           gsl_vector_get (s->f, 0), 
893           gsl_vector_get (s->f, 1));
894 @}
895 @end example
896
897 @noindent
898 Here are the results of running the program. The algorithm is started at
899 @math{(-10,-5)} far from the solution.  Since the solution is hidden in
900 a narrow valley the earliest steps follow the gradient of the function
901 downhill, in an attempt to reduce the large value of the residual. Once
902 the root has been approximately located, on iteration 8, the Newton
903 behavior takes over and convergence is very rapid.
904
905 @smallexample
906 iter =  0 x = -10.000  -5.000  f(x) = 1.100e+01 -1.050e+03
907 iter =  1 x = -10.000  -5.000  f(x) = 1.100e+01 -1.050e+03
908 iter =  2 x =  -3.976  24.827  f(x) = 4.976e+00  9.020e+01
909 iter =  3 x =  -3.976  24.827  f(x) = 4.976e+00  9.020e+01
910 iter =  4 x =  -3.976  24.827  f(x) = 4.976e+00  9.020e+01
911 iter =  5 x =  -1.274  -5.680  f(x) = 2.274e+00 -7.302e+01
912 iter =  6 x =  -1.274  -5.680  f(x) = 2.274e+00 -7.302e+01
913 iter =  7 x =   0.249   0.298  f(x) = 7.511e-01  2.359e+00
914 iter =  8 x =   0.249   0.298  f(x) = 7.511e-01  2.359e+00
915 iter =  9 x =   1.000   0.878  f(x) = 1.268e-10 -1.218e+00
916 iter = 10 x =   1.000   0.989  f(x) = 1.124e-11 -1.080e-01
917 iter = 11 x =   1.000   1.000  f(x) = 0.000e+00  0.000e+00
918 status = success
919 @end smallexample
920
921 @noindent
922 Note that the algorithm does not update the location on every
923 iteration. Some iterations are used to adjust the trust-region
924 parameter, after trying a step which was found to be divergent, or to
925 recompute the Jacobian, when poor convergence behavior is detected.
926
927 The next example program adds derivative information, in order to
928 accelerate the solution. There are two derivative functions
929 @code{rosenbrock_df} and @code{rosenbrock_fdf}. The latter computes both
930 the function and its derivative simultaneously. This allows the
931 optimization of any common terms.  For simplicity we substitute calls to
932 the separate @code{f} and @code{df} functions at this point in the code
933 below.
934
935 @example
936 int
937 rosenbrock_df (const gsl_vector * x, void *params, 
938                gsl_matrix * J)
939 @{
940   const double a = ((struct rparams *) params)->a;
941   const double b = ((struct rparams *) params)->b;
942
943   const double x0 = gsl_vector_get (x, 0);
944
945   const double df00 = -a;
946   const double df01 = 0;
947   const double df10 = -2 * b  * x0;
948   const double df11 = b;
949
950   gsl_matrix_set (J, 0, 0, df00);
951   gsl_matrix_set (J, 0, 1, df01);
952   gsl_matrix_set (J, 1, 0, df10);
953   gsl_matrix_set (J, 1, 1, df11);
954
955   return GSL_SUCCESS;
956 @}
957
958 int
959 rosenbrock_fdf (const gsl_vector * x, void *params,
960                 gsl_vector * f, gsl_matrix * J)
961 @{
962   rosenbrock_f (x, params, f);
963   rosenbrock_df (x, params, J);
964
965   return GSL_SUCCESS;
966 @}
967 @end example
968
969 @noindent
970 The main program now makes calls to the corresponding @code{fdfsolver}
971 versions of the functions,
972
973 @example
974 int
975 main (void)
976 @{
977   const gsl_multiroot_fdfsolver_type *T;
978   gsl_multiroot_fdfsolver *s;
979
980   int status;
981   size_t i, iter = 0;
982
983   const size_t n = 2;
984   struct rparams p = @{1.0, 10.0@};
985   gsl_multiroot_function_fdf f = @{&rosenbrock_f, 
986                                   &rosenbrock_df, 
987                                   &rosenbrock_fdf, 
988                                   n, &p@};
989
990   double x_init[2] = @{-10.0, -5.0@};
991   gsl_vector *x = gsl_vector_alloc (n);
992
993   gsl_vector_set (x, 0, x_init[0]);
994   gsl_vector_set (x, 1, x_init[1]);
995
996   T = gsl_multiroot_fdfsolver_gnewton;
997   s = gsl_multiroot_fdfsolver_alloc (T, n);
998   gsl_multiroot_fdfsolver_set (s, &f, x);
999
1000   print_state (iter, s);
1001
1002   do
1003     @{
1004       iter++;
1005
1006       status = gsl_multiroot_fdfsolver_iterate (s);
1007
1008       print_state (iter, s);
1009
1010       if (status)
1011         break;
1012
1013       status = gsl_multiroot_test_residual (s->f, 1e-7);
1014     @}
1015   while (status == GSL_CONTINUE && iter < 1000);
1016
1017   printf ("status = %s\n", gsl_strerror (status));
1018
1019   gsl_multiroot_fdfsolver_free (s);
1020   gsl_vector_free (x);
1021   return 0;
1022 @}
1023 @end example
1024
1025 @noindent
1026 The addition of derivative information to the @code{hybrids} solver does
1027 not make any significant difference to its behavior, since it able to
1028 approximate the Jacobian numerically with sufficient accuracy.  To
1029 illustrate the behavior of a different derivative solver we switch to
1030 @code{gnewton}. This is a traditional Newton solver with the constraint
1031 that it scales back its step if the full step would lead ``uphill''. Here
1032 is the output for the @code{gnewton} algorithm,
1033
1034 @smallexample
1035 iter = 0 x = -10.000  -5.000 f(x) =  1.100e+01 -1.050e+03
1036 iter = 1 x =  -4.231 -65.317 f(x) =  5.231e+00 -8.321e+02
1037 iter = 2 x =   1.000 -26.358 f(x) = -8.882e-16 -2.736e+02
1038 iter = 3 x =   1.000   1.000 f(x) = -2.220e-16 -4.441e-15
1039 status = success
1040 @end smallexample
1041
1042 @noindent
1043 The convergence is much more rapid, but takes a wide excursion out to
1044 the point @math{(-4.23,-65.3)}. This could cause the algorithm to go
1045 astray in a realistic application.  The hybrid algorithm follows the
1046 downhill path to the solution more reliably.
1047
1048 @node References and Further Reading for Multidimensional Root Finding
1049 @section References and Further Reading
1050
1051 The original version of the Hybrid method is described in the following
1052 articles by Powell,
1053
1054 @itemize @asis
1055 @item
1056 M.J.D. Powell, ``A Hybrid Method for Nonlinear Equations'' (Chap 6, p
1057 87--114) and ``A Fortran Subroutine for Solving systems of Nonlinear
1058 Algebraic Equations'' (Chap 7, p 115--161), in @cite{Numerical Methods for
1059 Nonlinear Algebraic Equations}, P. Rabinowitz, editor.  Gordon and
1060 Breach, 1970.
1061 @end itemize
1062
1063 @noindent
1064 The following papers are also relevant to the algorithms described in
1065 this section,
1066
1067 @itemize @asis
1068 @item
1069 J.J. Mor@'e, M.Y. Cosnard, ``Numerical Solution of Nonlinear Equations'',
1070 @cite{ACM Transactions on Mathematical Software}, Vol 5, No 1, (1979), p 64--85
1071
1072 @item
1073 C.G. Broyden, ``A Class of Methods for Solving Nonlinear
1074 Simultaneous Equations'', @cite{Mathematics of Computation}, Vol 19 (1965),
1075 p 577--593
1076
1077 @item 
1078 J.J. Mor@'e, B.S. Garbow, K.E. Hillstrom, ``Testing Unconstrained
1079 Optimization Software'', ACM Transactions on Mathematical Software, Vol
1080 7, No 1 (1981), p 17--41
1081 @end itemize
1082