Blame multiroots/dnewton.c

Packit 67cb25
/* multiroots/dnewton.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Brian Gough
Packit 67cb25
 * 
Packit 67cb25
 * This program is free software; you can redistribute it and/or modify
Packit 67cb25
 * it under the terms of the GNU General Public License as published by
Packit 67cb25
 * the Free Software Foundation; either version 3 of the License, or (at
Packit 67cb25
 * your option) any later version.
Packit 67cb25
 * 
Packit 67cb25
 * This program is distributed in the hope that it will be useful, but
Packit 67cb25
 * WITHOUT ANY WARRANTY; without even the implied warranty of
Packit 67cb25
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
Packit 67cb25
 * General Public License for more details.
Packit 67cb25
 * 
Packit 67cb25
 * You should have received a copy of the GNU General Public License
Packit 67cb25
 * along with this program; if not, write to the Free Software
Packit 67cb25
 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
#include <config.h>
Packit 67cb25
Packit 67cb25
#include <stddef.h>
Packit 67cb25
#include <stdlib.h>
Packit 67cb25
#include <stdio.h>
Packit 67cb25
#include <math.h>
Packit 67cb25
#include <float.h>
Packit 67cb25
Packit 67cb25
#include <gsl/gsl_math.h>
Packit 67cb25
#include <gsl/gsl_errno.h>
Packit 67cb25
#include <gsl/gsl_multiroots.h>
Packit 67cb25
#include <gsl/gsl_linalg.h>
Packit 67cb25
Packit 67cb25
/* Newton method using a finite difference approximation to the jacobian.
Packit 67cb25
   The derivatives are estimated using a step size of 
Packit 67cb25
Packit 67cb25
   h_i = sqrt(DBL_EPSILON) * x_i 
Packit 67cb25
Packit 67cb25
   */
Packit 67cb25
Packit 67cb25
typedef struct
Packit 67cb25
  {
Packit 67cb25
    gsl_matrix * J;
Packit 67cb25
    gsl_matrix * lu;
Packit 67cb25
    gsl_permutation * permutation;
Packit 67cb25
  }
Packit 67cb25
dnewton_state_t;
Packit 67cb25
Packit 67cb25
static int dnewton_alloc (void * vstate, size_t n);
Packit 67cb25
static int dnewton_set (void * vstate, gsl_multiroot_function * function, gsl_vector * x, gsl_vector * f, gsl_vector * dx);
Packit 67cb25
static int dnewton_iterate (void * vstate, gsl_multiroot_function * function, gsl_vector * x, gsl_vector * f, gsl_vector * dx);
Packit 67cb25
static void dnewton_free (void * vstate);
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
dnewton_alloc (void * vstate, size_t n)
Packit 67cb25
{
Packit 67cb25
  dnewton_state_t * state = (dnewton_state_t *) vstate;
Packit 67cb25
  gsl_permutation * p;
Packit 67cb25
  gsl_matrix * m, * J;
Packit 67cb25
Packit 67cb25
  m = gsl_matrix_calloc (n,n);
Packit 67cb25
  
Packit 67cb25
  if (m == 0) 
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("failed to allocate space for lu", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->lu = m ;
Packit 67cb25
Packit 67cb25
  p = gsl_permutation_calloc (n);
Packit 67cb25
Packit 67cb25
  if (p == 0)
Packit 67cb25
    {
Packit 67cb25
      gsl_matrix_free(m);
Packit 67cb25
Packit 67cb25
      GSL_ERROR ("failed to allocate space for permutation", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->permutation = p ;
Packit 67cb25
Packit 67cb25
  J = gsl_matrix_calloc (n,n);
Packit 67cb25
Packit 67cb25
  if (J == 0)
Packit 67cb25
    {
Packit 67cb25
      gsl_permutation_free(p);
Packit 67cb25
      gsl_matrix_free(m);
Packit 67cb25
Packit 67cb25
      GSL_ERROR ("failed to allocate space for d", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->J = J;
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
dnewton_set (void * vstate, gsl_multiroot_function * function, gsl_vector * x, gsl_vector * f, gsl_vector * dx)
Packit 67cb25
{
Packit 67cb25
  dnewton_state_t * state = (dnewton_state_t *) vstate;
Packit 67cb25
  size_t i, n = function->n ;
Packit 67cb25
  int status;
Packit 67cb25
Packit 67cb25
  status = GSL_MULTIROOT_FN_EVAL (function, x, f);
Packit 67cb25
  if (status)
Packit 67cb25
    return status;
Packit 67cb25
Packit 67cb25
  status = gsl_multiroot_fdjacobian (function, x, f, GSL_SQRT_DBL_EPSILON,
Packit 67cb25
      state->J);
Packit 67cb25
  if (status)
Packit 67cb25
    return status;
Packit 67cb25
Packit 67cb25
  for (i = 0; i < n; i++)
Packit 67cb25
    {
Packit 67cb25
      gsl_vector_set (dx, i, 0.0);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
dnewton_iterate (void * vstate, gsl_multiroot_function * function, gsl_vector * x, gsl_vector * f, gsl_vector * dx)
Packit 67cb25
{
Packit 67cb25
  dnewton_state_t * state = (dnewton_state_t *) vstate;
Packit 67cb25
  
Packit 67cb25
  int signum ;
Packit 67cb25
Packit 67cb25
  size_t i;
Packit 67cb25
Packit 67cb25
  size_t n = function->n ;
Packit 67cb25
Packit 67cb25
  gsl_matrix_memcpy (state->lu, state->J);
Packit 67cb25
Packit 67cb25
  {
Packit 67cb25
    int status = gsl_linalg_LU_decomp (state->lu, state->permutation, &signum);
Packit 67cb25
    if (status)
Packit 67cb25
      return status;
Packit 67cb25
  }
Packit 67cb25
  
Packit 67cb25
  {
Packit 67cb25
    int status = gsl_linalg_LU_solve (state->lu, state->permutation, f, dx);
Packit 67cb25
    if (status)
Packit 67cb25
      return status;
Packit 67cb25
  }
Packit 67cb25
Packit 67cb25
  for (i = 0; i < n; i++)
Packit 67cb25
    {
Packit 67cb25
      double e = gsl_vector_get (dx, i);
Packit 67cb25
      double y = gsl_vector_get (x, i);
Packit 67cb25
      gsl_vector_set (dx, i, -e);
Packit 67cb25
      gsl_vector_set (x, i, y - e);
Packit 67cb25
    }
Packit 67cb25
  
Packit 67cb25
  {
Packit 67cb25
    int status = GSL_MULTIROOT_FN_EVAL (function, x, f);
Packit 67cb25
Packit 67cb25
    if (status != GSL_SUCCESS) 
Packit 67cb25
      {
Packit 67cb25
        return GSL_EBADFUNC;
Packit 67cb25
      }
Packit 67cb25
  }
Packit 67cb25
 
Packit 67cb25
  gsl_multiroot_fdjacobian (function, x, f, GSL_SQRT_DBL_EPSILON, state->J);
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
static void
Packit 67cb25
dnewton_free (void * vstate)
Packit 67cb25
{
Packit 67cb25
  dnewton_state_t * state = (dnewton_state_t *) vstate;
Packit 67cb25
Packit 67cb25
  gsl_matrix_free(state->J);
Packit 67cb25
  gsl_matrix_free(state->lu);
Packit 67cb25
  gsl_permutation_free(state->permutation);
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
Packit 67cb25
static const gsl_multiroot_fsolver_type dnewton_type =
Packit 67cb25
{"dnewton",                             /* name */
Packit 67cb25
 sizeof (dnewton_state_t),
Packit 67cb25
 &dnewton_alloc,
Packit 67cb25
 &dnewton_set,
Packit 67cb25
 &dnewton_iterate,
Packit 67cb25
 &dnewton_free};
Packit 67cb25
Packit 67cb25
const gsl_multiroot_fsolver_type  * gsl_multiroot_fsolver_dnewton = &dnewton_type;