|
Packit |
67cb25 |
/* multilarge_nlinear/fdf.c
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* Copyright (C) 2015, 2016 Patrick Alken
|
|
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 |
#include <stdlib.h>
|
|
Packit |
67cb25 |
#include <string.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_math.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_errno.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_blas.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_multilarge_nlinear.h>
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_workspace *
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_alloc (const gsl_multilarge_nlinear_type * T,
|
|
Packit |
67cb25 |
const gsl_multilarge_nlinear_parameters * params,
|
|
Packit |
67cb25 |
const size_t n, const size_t p)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_workspace * w;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (n < p)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR_VAL ("insufficient data points, n < p", GSL_EINVAL, 0);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
w = calloc (1, sizeof (gsl_multilarge_nlinear_workspace));
|
|
Packit |
67cb25 |
if (w == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR_VAL ("failed to allocate space for workspace",
|
|
Packit |
67cb25 |
GSL_ENOMEM, 0);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
w->n = n;
|
|
Packit |
67cb25 |
w->p = p;
|
|
Packit |
67cb25 |
w->type = T;
|
|
Packit |
67cb25 |
w->fdf = NULL;
|
|
Packit |
67cb25 |
w->niter = 0;
|
|
Packit |
67cb25 |
w->params = *params;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* the cgst method uses its own built-in linear solver */
|
|
Packit |
67cb25 |
if (w->params.trs == gsl_multilarge_nlinear_trs_cgst)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
w->params.solver = gsl_multilarge_nlinear_solver_none;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
w->x = gsl_vector_calloc (p);
|
|
Packit |
67cb25 |
if (w->x == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_free (w);
|
|
Packit |
67cb25 |
GSL_ERROR_VAL ("failed to allocate space for x", GSL_ENOMEM, 0);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
w->f = gsl_vector_calloc (n);
|
|
Packit |
67cb25 |
if (w->f == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_free (w);
|
|
Packit |
67cb25 |
GSL_ERROR_VAL ("failed to allocate space for f", GSL_ENOMEM, 0);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
w->dx = gsl_vector_calloc (p);
|
|
Packit |
67cb25 |
if (w->dx == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_free (w);
|
|
Packit |
67cb25 |
GSL_ERROR_VAL ("failed to allocate space for dx", GSL_ENOMEM, 0);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
w->g = gsl_vector_alloc (p);
|
|
Packit |
67cb25 |
if (w->g == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_free (w);
|
|
Packit |
67cb25 |
GSL_ERROR_VAL ("failed to allocate space for g", GSL_ENOMEM, 0);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (w->params.solver == gsl_multilarge_nlinear_solver_cholesky)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
w->JTJ = gsl_matrix_alloc (p, p);
|
|
Packit |
67cb25 |
if (w->JTJ == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_free (w);
|
|
Packit |
67cb25 |
GSL_ERROR_VAL ("failed to allocate space for JTJ", GSL_ENOMEM, 0);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
w->sqrt_wts_work = gsl_vector_calloc (n);
|
|
Packit |
67cb25 |
if (w->sqrt_wts_work == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_free (w);
|
|
Packit |
67cb25 |
GSL_ERROR_VAL ("failed to allocate space for weights", GSL_ENOMEM, 0);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
w->state = (T->alloc)(&(w->params), n, p);
|
|
Packit |
67cb25 |
if (w->state == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_free (w);
|
|
Packit |
67cb25 |
GSL_ERROR_VAL ("failed to allocate space for state", GSL_ENOMEM, 0);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return w;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
void
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_free (gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
RETURN_IF_NULL (w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (w->state)
|
|
Packit |
67cb25 |
(w->type->free) (w->state);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (w->dx)
|
|
Packit |
67cb25 |
gsl_vector_free (w->dx);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (w->x)
|
|
Packit |
67cb25 |
gsl_vector_free (w->x);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (w->f)
|
|
Packit |
67cb25 |
gsl_vector_free (w->f);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (w->sqrt_wts_work)
|
|
Packit |
67cb25 |
gsl_vector_free (w->sqrt_wts_work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (w->g)
|
|
Packit |
67cb25 |
gsl_vector_free (w->g);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (w->JTJ)
|
|
Packit |
67cb25 |
gsl_matrix_free (w->JTJ);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
free (w);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_parameters
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_default_parameters(void)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_parameters params;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
params.trs = gsl_multilarge_nlinear_trs_lm;
|
|
Packit |
67cb25 |
params.scale = gsl_multilarge_nlinear_scale_more;
|
|
Packit |
67cb25 |
params.solver = gsl_multilarge_nlinear_solver_cholesky;
|
|
Packit |
67cb25 |
params.fdtype = GSL_MULTILARGE_NLINEAR_FWDIFF;
|
|
Packit |
67cb25 |
params.factor_up = 3.0;
|
|
Packit |
67cb25 |
params.factor_down = 2.0;
|
|
Packit |
67cb25 |
params.avmax = 0.75;
|
|
Packit |
67cb25 |
params.h_df = GSL_SQRT_DBL_EPSILON;
|
|
Packit |
67cb25 |
params.h_fvv = 0.01;
|
|
Packit |
67cb25 |
params.max_iter = 0;
|
|
Packit |
67cb25 |
params.tol = 1.0e-6;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return params;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_init (const gsl_vector * x,
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_fdf * fdf,
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return gsl_multilarge_nlinear_winit(x, NULL, fdf, w);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_winit (const gsl_vector * x,
|
|
Packit |
67cb25 |
const gsl_vector * wts,
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_fdf * fdf,
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const size_t n = w->f->size;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (n != fdf->n)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("function size does not match workspace", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (w->x->size != x->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("vector length does not match workspace", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (wts != NULL && n != wts->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("weight vector length does not match workspace", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t i;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* initialize counters for function and Jacobian evaluations */
|
|
Packit |
67cb25 |
fdf->nevalf = 0;
|
|
Packit |
67cb25 |
fdf->nevaldfu = 0;
|
|
Packit |
67cb25 |
fdf->nevaldf2 = 0;
|
|
Packit |
67cb25 |
fdf->nevalfvv = 0;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
w->fdf = fdf;
|
|
Packit |
67cb25 |
gsl_vector_memcpy(w->x, x);
|
|
Packit |
67cb25 |
w->niter = 0;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (wts)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
w->sqrt_wts = w->sqrt_wts_work;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 0; i < n; ++i)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
double wi = gsl_vector_get(wts, i);
|
|
Packit |
67cb25 |
gsl_vector_set(w->sqrt_wts, i, sqrt(wi));
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
w->sqrt_wts = NULL;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return (w->type->init) (w->state, w->sqrt_wts, w->fdf,
|
|
Packit |
67cb25 |
w->x, w->f, w->g, w->JTJ);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_iterate (gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status =
|
|
Packit |
67cb25 |
(w->type->iterate) (w->state, w->sqrt_wts, w->fdf,
|
|
Packit |
67cb25 |
w->x, w->f, w->g, w->JTJ, w->dx);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
w->niter++;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_avratio (const gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return (w->type->avratio) (w->state);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_rcond (double * rcond, const gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status = (w->type->rcond) (rcond, w->JTJ, w->state);
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_covar (gsl_matrix * covar, gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
if (covar->size1 != covar->size2)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("covariance matrix must be square", GSL_ENOTSQR);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (covar->size1 != w->p)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("covariance matrix does not match workspace", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status = (w->type->covar) (w->JTJ, covar, w->state);
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_driver()
|
|
Packit |
67cb25 |
Iterate the nonlinear least squares solver until completion
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Inputs: maxiter - maximum iterations to allow
|
|
Packit |
67cb25 |
xtol - tolerance in step x
|
|
Packit |
67cb25 |
gtol - tolerance in gradient
|
|
Packit |
67cb25 |
ftol - tolerance in ||f||
|
|
Packit |
67cb25 |
callback - callback function to call each iteration
|
|
Packit |
67cb25 |
callback_params - parameters to pass to callback function
|
|
Packit |
67cb25 |
info - (output) info flag on why iteration terminated
|
|
Packit |
67cb25 |
1 = stopped due to small step size ||dx|
|
|
Packit |
67cb25 |
2 = stopped due to small gradient
|
|
Packit |
67cb25 |
3 = stopped due to small change in f
|
|
Packit |
67cb25 |
GSL_ETOLX = ||dx|| has converged to within machine
|
|
Packit |
67cb25 |
precision (and xtol is too small)
|
|
Packit |
67cb25 |
GSL_ETOLG = ||g||_inf is smaller than machine
|
|
Packit |
67cb25 |
precision (gtol is too small)
|
|
Packit |
67cb25 |
GSL_ETOLF = change in ||f|| is smaller than machine
|
|
Packit |
67cb25 |
precision (ftol is too small)
|
|
Packit |
67cb25 |
w - workspace
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Return:
|
|
Packit |
67cb25 |
GSL_SUCCESS if converged
|
|
Packit |
67cb25 |
GSL_MAXITER if maxiter exceeded without converging
|
|
Packit |
67cb25 |
GSL_ENOPROG if no accepted step found on first iteration
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_driver (const size_t maxiter,
|
|
Packit |
67cb25 |
const double xtol,
|
|
Packit |
67cb25 |
const double gtol,
|
|
Packit |
67cb25 |
const double ftol,
|
|
Packit |
67cb25 |
void (*callback)(const size_t iter, void *params,
|
|
Packit |
67cb25 |
const gsl_multilarge_nlinear_workspace *w),
|
|
Packit |
67cb25 |
void *callback_params,
|
|
Packit |
67cb25 |
int *info,
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status;
|
|
Packit |
67cb25 |
size_t iter = 0;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* call user callback function prior to any iterations
|
|
Packit |
67cb25 |
* with initial system state */
|
|
Packit |
67cb25 |
if (callback)
|
|
Packit |
67cb25 |
callback(iter, callback_params, w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
do
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
status = gsl_multilarge_nlinear_iterate (w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
* If the solver reports no progress on the first iteration,
|
|
Packit |
67cb25 |
* then it didn't find a single step to reduce the
|
|
Packit |
67cb25 |
* cost function and more iterations won't help so return.
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* If we get a no progress flag on subsequent iterations,
|
|
Packit |
67cb25 |
* it means we did find a good step in a previous iteration,
|
|
Packit |
67cb25 |
* so continue iterating since the solver has now reset
|
|
Packit |
67cb25 |
* mu to its initial value.
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
if (status == GSL_ENOPROG && iter == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
*info = status;
|
|
Packit |
67cb25 |
return GSL_EMAXITER;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
++iter;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (callback)
|
|
Packit |
67cb25 |
callback(iter, callback_params, w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* test for convergence */
|
|
Packit |
67cb25 |
status = gsl_multilarge_nlinear_test(xtol, gtol, ftol, info, w);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
while (status == GSL_CONTINUE && iter < maxiter);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
* the following error codes mean that the solution has converged
|
|
Packit |
67cb25 |
* to within machine precision, so record the error code in info
|
|
Packit |
67cb25 |
* and return success
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
if (status == GSL_ETOLF || status == GSL_ETOLX || status == GSL_ETOLG)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
*info = status;
|
|
Packit |
67cb25 |
status = GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* check if max iterations reached */
|
|
Packit |
67cb25 |
if (iter >= maxiter && status != GSL_SUCCESS)
|
|
Packit |
67cb25 |
status = GSL_EMAXITER;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
} /* gsl_multilarge_nlinear_driver() */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
const char *
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_name (const gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return w->type->name;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector *
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_position (const gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return w->x;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector *
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_residual (const gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return w->f;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector *
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_step (const gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return w->dx;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
size_t
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_niter (const gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return w->niter;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
const char *
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_trs_name (const gsl_multilarge_nlinear_workspace * w)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return w->params.trs->name;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_eval_f()
|
|
Packit |
67cb25 |
Compute residual vector y with user callback function, and apply
|
|
Packit |
67cb25 |
weighting transform if given:
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
y~ = sqrt(W) y
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Inputs: fdf - callback function
|
|
Packit |
67cb25 |
x - model parameters
|
|
Packit |
67cb25 |
swts - weight matrix sqrt(W) = sqrt(diag(w1,w2,...,wn))
|
|
Packit |
67cb25 |
set to NULL for unweighted fit
|
|
Packit |
67cb25 |
y - (output) (weighted) residual vector
|
|
Packit |
67cb25 |
y_i = sqrt(w_i) f_i where f_i is unweighted residual
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_eval_f(gsl_multilarge_nlinear_fdf *fdf,
|
|
Packit |
67cb25 |
const gsl_vector *x,
|
|
Packit |
67cb25 |
const gsl_vector *swts,
|
|
Packit |
67cb25 |
gsl_vector *y)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int s = ((*((fdf)->f)) (x, fdf->params, y));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
++(fdf->nevalf);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* y <- sqrt(W) y */
|
|
Packit |
67cb25 |
if (swts)
|
|
Packit |
67cb25 |
gsl_vector_mul(y, swts);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return s;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_eval_df()
|
|
Packit |
67cb25 |
Compute Jacobian matrix-vector product:
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
v = J * u
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
or
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
v = J^T u
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Inputs: TransJ - use J or J^T
|
|
Packit |
67cb25 |
x - model parameters
|
|
Packit |
67cb25 |
f - residual vector f(x)
|
|
Packit |
67cb25 |
u - input vector u
|
|
Packit |
67cb25 |
swts - weight matrix W = diag(w1,w2,...,wn)
|
|
Packit |
67cb25 |
set to NULL for unweighted fit
|
|
Packit |
67cb25 |
h - finite difference step size
|
|
Packit |
67cb25 |
fdtype - finite difference method
|
|
Packit |
67cb25 |
fdf - callback function
|
|
Packit |
67cb25 |
v - (output) vector v
|
|
Packit |
67cb25 |
JTJ - (output) matrix J^T J
|
|
Packit |
67cb25 |
work - workspace for finite difference, size n
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_eval_df(const CBLAS_TRANSPOSE_t TransJ,
|
|
Packit |
67cb25 |
const gsl_vector *x,
|
|
Packit |
67cb25 |
const gsl_vector *f,
|
|
Packit |
67cb25 |
const gsl_vector *u,
|
|
Packit |
67cb25 |
const gsl_vector *swts,
|
|
Packit |
67cb25 |
const double h,
|
|
Packit |
67cb25 |
const gsl_multilarge_nlinear_fdtype fdtype,
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_fdf *fdf,
|
|
Packit |
67cb25 |
gsl_vector *v,
|
|
Packit |
67cb25 |
gsl_matrix *JTJ,
|
|
Packit |
67cb25 |
gsl_vector *work)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const size_t n = fdf->n;
|
|
Packit |
67cb25 |
const size_t p = fdf->p;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (u != NULL && ((TransJ == CblasNoTrans && u->size != p) ||
|
|
Packit |
67cb25 |
(TransJ == CblasTrans && u->size != n)))
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR("u vector has wrong size", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (v != NULL && ((TransJ == CblasNoTrans && v->size != n) ||
|
|
Packit |
67cb25 |
(TransJ == CblasTrans && v->size != p)))
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR("v vector has wrong size", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (JTJ != NULL && ((JTJ->size1 != p) || (JTJ->size2 != p)))
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR("JTJ matrix has wrong size", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (fdf->df)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
/* call user-supplied function */
|
|
Packit |
67cb25 |
status = ((*((fdf)->df)) (TransJ, x, u, fdf->params, v, JTJ));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (v)
|
|
Packit |
67cb25 |
++(fdf->nevaldfu);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (JTJ)
|
|
Packit |
67cb25 |
++(fdf->nevaldf2);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
#if 0
|
|
Packit |
67cb25 |
/* use finite difference Jacobian approximation */
|
|
Packit |
67cb25 |
status = gsl_multilarge_nlinear_df(h, fdtype, x, swts, fdf, f, df, work);
|
|
Packit |
67cb25 |
#endif
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_eval_fvv()
|
|
Packit |
67cb25 |
Compute second direction derivative vector yvv with user
|
|
Packit |
67cb25 |
callback function, and apply weighting transform if given:
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
yvv~ = sqrt(W) yvv
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Inputs: h - step size for finite difference, if needed
|
|
Packit |
67cb25 |
x - model parameters, size p
|
|
Packit |
67cb25 |
v - unscaled geodesic velocity vector, size p
|
|
Packit |
67cb25 |
f - residual vector f(x), size n
|
|
Packit |
67cb25 |
swts - weight matrix sqrt(W) = sqrt(diag(w1,w2,...,wn))
|
|
Packit |
67cb25 |
set to NULL for unweighted fit
|
|
Packit |
67cb25 |
fdf - callback function
|
|
Packit |
67cb25 |
yvv - (output) (weighted) second directional derivative vector
|
|
Packit |
67cb25 |
yvv_i = sqrt(w_i) fvv_i where f_i is unweighted
|
|
Packit |
67cb25 |
work - workspace, size p
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_eval_fvv(const double h,
|
|
Packit |
67cb25 |
const gsl_vector *x,
|
|
Packit |
67cb25 |
const gsl_vector *v,
|
|
Packit |
67cb25 |
const gsl_vector *f,
|
|
Packit |
67cb25 |
const gsl_vector *swts,
|
|
Packit |
67cb25 |
gsl_multilarge_nlinear_fdf *fdf,
|
|
Packit |
67cb25 |
gsl_vector *yvv,
|
|
Packit |
67cb25 |
gsl_vector *work)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (fdf->fvv != NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
/* call user-supplied function */
|
|
Packit |
67cb25 |
status = ((*((fdf)->fvv)) (x, v, fdf->params, yvv));
|
|
Packit |
67cb25 |
++(fdf->nevalfvv);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
#if 0
|
|
Packit |
67cb25 |
/* use finite difference approximation */
|
|
Packit |
67cb25 |
status = gsl_multilarge_nlinear_fdfvv(h, x, v, f, J,
|
|
Packit |
67cb25 |
swts, fdf, yvv, work);
|
|
Packit |
67cb25 |
#endif
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* yvv <- sqrt(W) yvv */
|
|
Packit |
67cb25 |
if (swts)
|
|
Packit |
67cb25 |
gsl_vector_mul(yvv, swts);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|