|
Packit |
67cb25 |
/* multiroots/gnewton.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 |
#include "enorm.c"
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Simple globally convergent Newton method (rejects uphill steps) */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
typedef struct
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
double phi;
|
|
Packit |
67cb25 |
gsl_vector * x_trial;
|
|
Packit |
67cb25 |
gsl_vector * d;
|
|
Packit |
67cb25 |
gsl_matrix * lu;
|
|
Packit |
67cb25 |
gsl_permutation * permutation;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
gnewton_state_t;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int gnewton_alloc (void * vstate, size_t n);
|
|
Packit |
67cb25 |
static int gnewton_set (void * vstate, gsl_multiroot_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx);
|
|
Packit |
67cb25 |
static int gnewton_iterate (void * vstate, gsl_multiroot_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx);
|
|
Packit |
67cb25 |
static void gnewton_free (void * vstate);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
gnewton_alloc (void * vstate, size_t n)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gnewton_state_t * state = (gnewton_state_t *) vstate;
|
|
Packit |
67cb25 |
gsl_vector * d, * x_trial ;
|
|
Packit |
67cb25 |
gsl_permutation * p;
|
|
Packit |
67cb25 |
gsl_matrix * m;
|
|
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 |
d = gsl_vector_calloc (n);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (d == 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->d = d;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
x_trial = gsl_vector_calloc (n);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (x_trial == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_vector_free(d);
|
|
Packit |
67cb25 |
gsl_permutation_free(p);
|
|
Packit |
67cb25 |
gsl_matrix_free(m);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
GSL_ERROR ("failed to allocate space for x_trial", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->x_trial = x_trial;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
gnewton_set (void * vstate, gsl_multiroot_function_fdf * FDF, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gnewton_state_t * state = (gnewton_state_t *) vstate;
|
|
Packit |
67cb25 |
size_t i, n = FDF->n ;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
GSL_MULTIROOT_FN_EVAL_F_DF (FDF, x, f, J);
|
|
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 |
state->phi = enorm(f);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
gnewton_iterate (void * vstate, gsl_multiroot_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gnewton_state_t * state = (gnewton_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int signum ;
|
|
Packit |
67cb25 |
double t, phi0, phi1;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
size_t i;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
size_t n = fdf->n ;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_matrix_memcpy (state->lu, J);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_linalg_LU_decomp (state->lu, state->permutation, &signum);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status = gsl_linalg_LU_solve (state->lu, state->permutation, f, state->d);
|
|
Packit |
67cb25 |
if (status)
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
t = 1;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
phi0 = state->phi;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
new_step:
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 0; i < n; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
double di = gsl_vector_get (state->d, i);
|
|
Packit |
67cb25 |
double xi = gsl_vector_get (x, i);
|
|
Packit |
67cb25 |
gsl_vector_set (state->x_trial, i, xi - t*di);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status = GSL_MULTIROOT_FN_EVAL_F (fdf, state->x_trial, f);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (status != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return GSL_EBADFUNC;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
phi1 = enorm (f);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (phi1 > phi0 && t > GSL_DBL_EPSILON)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
/* full step goes uphill, take a reduced step instead */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double theta = phi1 / phi0;
|
|
Packit |
67cb25 |
double u = (sqrt(1.0 + 6.0 * theta) - 1.0) / (3.0 * theta);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
t *= u ;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
goto new_step;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* copy x_trial into x */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector_memcpy (x, state->x_trial);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 0; i < n; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
double di = gsl_vector_get (state->d, i);
|
|
Packit |
67cb25 |
gsl_vector_set (dx, i, -t*di);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status = GSL_MULTIROOT_FN_EVAL_DF (fdf, x, J);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (status != GSL_SUCCESS)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return GSL_EBADFUNC;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->phi = phi1;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
gnewton_free (void * vstate)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gnewton_state_t * state = (gnewton_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector_free(state->d);
|
|
Packit |
67cb25 |
gsl_vector_free(state->x_trial);
|
|
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_fdfsolver_type gnewton_type =
|
|
Packit |
67cb25 |
{"gnewton", /* name */
|
|
Packit |
67cb25 |
sizeof (gnewton_state_t),
|
|
Packit |
67cb25 |
&gnewton_alloc,
|
|
Packit |
67cb25 |
&gnewton_set,
|
|
Packit |
67cb25 |
&gnewton_iterate,
|
|
Packit |
67cb25 |
&gnewton_free};
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
const gsl_multiroot_fdfsolver_type * gsl_multiroot_fdfsolver_gnewton = &gnewton_type;
|