Blame multifit_nlinear/nielsen.c

Packit 67cb25
/* multifit_nlinear/nielsen.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 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
/*
Packit 67cb25
 * This module contains routines for updating the Levenberg-Marquardt
Packit 67cb25
 * damping parameter on each iteration using Nielsen's method:
Packit 67cb25
 *
Packit 67cb25
 * [1] H. B. Nielsen, K. Madsen, Introduction to Optimization and
Packit 67cb25
 *     Data Fitting, Informatics and Mathematical Modeling,
Packit 67cb25
 *     Technical University of Denmark (DTU), 2010.
Packit 67cb25
 *
Packit 67cb25
 * 5 routines are needed to implement the update procedure:
Packit 67cb25
 *
Packit 67cb25
 * 1. accept - update parameter after a step has been accepted
Packit 67cb25
 * 2. reject - update parameter after a step has been rejected
Packit 67cb25
 * 3. free   - free workspace state
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
#include <config.h>
Packit 67cb25
#include <gsl/gsl_math.h>
Packit 67cb25
#include <gsl/gsl_vector.h>
Packit 67cb25
#include <gsl/gsl_matrix.h>
Packit 67cb25
#include <gsl/gsl_linalg.h>
Packit 67cb25
#include <gsl/gsl_multifit_nlinear.h>
Packit 67cb25
#include <gsl/gsl_errno.h>
Packit 67cb25
#include <gsl/gsl_blas.h>
Packit 67cb25
Packit 67cb25
#define LM_ONE_THIRD         (0.333333333333333)
Packit 67cb25
Packit 67cb25
static int nielsen_init(const gsl_matrix * J, const gsl_vector * diag,
Packit 67cb25
                        double * mu, long * nu);
Packit 67cb25
static int nielsen_accept(const double rho, double * mu, long * nu);
Packit 67cb25
static int nielsen_reject(double * mu, long * nu);
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
nielsen_init(const gsl_matrix * J, const gsl_vector * diag,
Packit 67cb25
             double * mu, long * nu)
Packit 67cb25
{
Packit 67cb25
  const double mu0 = 1.0e-3;
Packit 67cb25
  const size_t p = J->size2;
Packit 67cb25
  size_t j;
Packit 67cb25
  double max = -1.0;
Packit 67cb25
Packit 67cb25
  *nu = 2;
Packit 67cb25
Packit 67cb25
  /* set mu = mu0 * max(diag(J~^T J~)), with J~ = J D^{-1} */
Packit 67cb25
Packit 67cb25
  for (j = 0; j < p; ++j)
Packit 67cb25
    {
Packit 67cb25
      gsl_vector_const_view v = gsl_matrix_const_column(J, j);
Packit 67cb25
      double dj = gsl_vector_get(diag, j);
Packit 67cb25
      double norm = gsl_blas_dnrm2(&v.vector) / dj;
Packit 67cb25
      max = GSL_MAX(max, norm);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  *mu = mu0 * max * max;
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
nielsen_accept(const double rho, double * mu, long * nu)
Packit 67cb25
{
Packit 67cb25
  double  b;
Packit 67cb25
  
Packit 67cb25
  /* reset nu */
Packit 67cb25
  *nu = 2;
Packit 67cb25
Packit 67cb25
  b = 2.0 * rho - 1.0;
Packit 67cb25
  b = 1.0 - b*b*b;
Packit 67cb25
  *mu *= GSL_MAX(LM_ONE_THIRD, b);
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
nielsen_reject(double * mu, long * nu)
Packit 67cb25
{
Packit 67cb25
  *mu *= (double) *nu;
Packit 67cb25
Packit 67cb25
  /* nu := 2*nu */
Packit 67cb25
  *nu <<= 1;
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}