|
Packit |
67cb25 |
/* multifit_nlinear/cholesky.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 |
/*
|
|
Packit |
67cb25 |
* This module calculates the solution of the normal equations least squares
|
|
Packit |
67cb25 |
* system:
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* [ J~^T J~ + mu D~^T D~ ] p~ = -J~^T f
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* using the modified Cholesky decomposition. Quantities are scaled
|
|
Packit |
67cb25 |
* according to:
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* J~ = J S
|
|
Packit |
67cb25 |
* D~ = D S
|
|
Packit |
67cb25 |
* p~ = S^{-1} p
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* where S is a diagonal matrix and S_jj = || J_j || and J_j is column
|
|
Packit |
67cb25 |
* j of the Jacobian. This balancing transformation seems to be more
|
|
Packit |
67cb25 |
* numerically stable for some Jacobians.
|
|
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 |
#include <gsl/gsl_permutation.h>
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
#include "common.c"
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
typedef struct
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_matrix *JTJ; /* J^T J */
|
|
Packit |
67cb25 |
gsl_matrix *work_JTJ; /* copy of J^T J */
|
|
Packit |
67cb25 |
gsl_vector *rhs; /* -J^T f, size p */
|
|
Packit |
67cb25 |
gsl_permutation *perm; /* permutation matrix for modified Cholesky */
|
|
Packit |
67cb25 |
gsl_vector *work3p; /* workspace, size 3*p */
|
|
Packit |
67cb25 |
double mu; /* current regularization parameter */
|
|
Packit |
67cb25 |
} cholesky_state_t;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void *cholesky_alloc (const size_t n, const size_t p);
|
|
Packit |
67cb25 |
static int cholesky_init(const void * vtrust_state, void * vstate);
|
|
Packit |
67cb25 |
static int cholesky_presolve(const double mu, const void * vtrust_state, void * vstate);
|
|
Packit |
67cb25 |
static int cholesky_solve(const gsl_vector * f, gsl_vector *x,
|
|
Packit |
67cb25 |
const void * vtrust_state, void *vstate);
|
|
Packit |
67cb25 |
static int cholesky_solve_rhs(const gsl_vector * b, gsl_vector *x, cholesky_state_t *state);
|
|
Packit |
67cb25 |
static int cholesky_regularize(const double mu, const gsl_vector * diag, gsl_matrix * A,
|
|
Packit |
67cb25 |
cholesky_state_t * state);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void *
|
|
Packit |
67cb25 |
cholesky_alloc (const size_t n, const size_t p)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
cholesky_state_t *state;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
(void)n;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state = calloc(1, sizeof(cholesky_state_t));
|
|
Packit |
67cb25 |
if (state == NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate cholesky state", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->JTJ = gsl_matrix_alloc(p, p);
|
|
Packit |
67cb25 |
if (state->JTJ == NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for JTJ", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->work_JTJ = gsl_matrix_alloc(p, p);
|
|
Packit |
67cb25 |
if (state->work_JTJ == NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for JTJ workspace",
|
|
Packit |
67cb25 |
GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->rhs = gsl_vector_alloc(p);
|
|
Packit |
67cb25 |
if (state->rhs == NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for rhs", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->perm = gsl_permutation_alloc(p);
|
|
Packit |
67cb25 |
if (state->perm == NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for perm", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->work3p = gsl_vector_alloc(3 * p);
|
|
Packit |
67cb25 |
if (state->work3p == NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR_NULL ("failed to allocate space for work3p", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->mu = -1.0;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return state;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
cholesky_free(void *vstate)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
cholesky_state_t *state = (cholesky_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->JTJ)
|
|
Packit |
67cb25 |
gsl_matrix_free(state->JTJ);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->work_JTJ)
|
|
Packit |
67cb25 |
gsl_matrix_free(state->work_JTJ);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->rhs)
|
|
Packit |
67cb25 |
gsl_vector_free(state->rhs);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->perm)
|
|
Packit |
67cb25 |
gsl_permutation_free(state->perm);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->work3p)
|
|
Packit |
67cb25 |
gsl_vector_free(state->work3p);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
free(state);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
cholesky_init(const void * vtrust_state, void * vstate)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const gsl_multifit_nlinear_trust_state *trust_state =
|
|
Packit |
67cb25 |
(const gsl_multifit_nlinear_trust_state *) vtrust_state;
|
|
Packit |
67cb25 |
cholesky_state_t *state = (cholesky_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* compute J^T J */
|
|
Packit |
67cb25 |
gsl_blas_dsyrk(CblasLower, CblasTrans, 1.0, trust_state->J, 0.0, state->JTJ);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
cholesky_presolve()
|
|
Packit |
67cb25 |
Compute the modified Cholesky decomposition of J^T J + mu D^T D.
|
|
Packit |
67cb25 |
Modified Cholesky is used in case mu = 0 and there are rounding
|
|
Packit |
67cb25 |
errors in forming J^T J which could lead to an indefinite matrix.
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Inputs: mu - LM parameter
|
|
Packit |
67cb25 |
vstate - workspace
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Notes:
|
|
Packit |
67cb25 |
1) On output, state->work_JTJ contains the Cholesky decomposition of
|
|
Packit |
67cb25 |
J^T J + mu D^T D
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
cholesky_presolve(const double mu, const void * vtrust_state, void * vstate)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const gsl_multifit_nlinear_trust_state *trust_state =
|
|
Packit |
67cb25 |
(const gsl_multifit_nlinear_trust_state *) vtrust_state;
|
|
Packit |
67cb25 |
cholesky_state_t *state = (cholesky_state_t *) vstate;
|
|
Packit |
67cb25 |
gsl_matrix *JTJ = state->work_JTJ;
|
|
Packit |
67cb25 |
const gsl_vector *diag = trust_state->diag;
|
|
Packit |
67cb25 |
int status;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* copy lower triangle of A to workspace */
|
|
Packit |
67cb25 |
gsl_matrix_tricpy('L', 1, JTJ, state->JTJ);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* augment normal equations: A -> A + mu D^T D */
|
|
Packit |
67cb25 |
status = cholesky_regularize(mu, diag, JTJ, state);
|
|
Packit |
67cb25 |
if (status)
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* compute modified Cholesky decomposition */
|
|
Packit |
67cb25 |
status = gsl_linalg_mcholesky_decomp(JTJ, state->perm, NULL);
|
|
Packit |
67cb25 |
if (status)
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->mu = mu;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
cholesky_solve()
|
|
Packit |
67cb25 |
Compute (J^T J + mu D^T D) x = -J^T f
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Inputs: f - right hand side vector f
|
|
Packit |
67cb25 |
x - (output) solution vector
|
|
Packit |
67cb25 |
vstate - cholesky workspace
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
cholesky_solve(const gsl_vector * f, gsl_vector *x,
|
|
Packit |
67cb25 |
const void * vtrust_state, void *vstate)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const gsl_multifit_nlinear_trust_state *trust_state =
|
|
Packit |
67cb25 |
(const gsl_multifit_nlinear_trust_state *) vtrust_state;
|
|
Packit |
67cb25 |
cholesky_state_t *state = (cholesky_state_t *) vstate;
|
|
Packit |
67cb25 |
int status;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* compute rhs = -J^T f */
|
|
Packit |
67cb25 |
gsl_blas_dgemv(CblasTrans, -1.0, trust_state->J, f, 0.0, state->rhs);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
status = cholesky_solve_rhs(state->rhs, x, state);
|
|
Packit |
67cb25 |
if (status)
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
cholesky_rcond(double * rcond, void * vstate)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status;
|
|
Packit |
67cb25 |
cholesky_state_t *state = (cholesky_state_t *) vstate;
|
|
Packit |
67cb25 |
double rcond_JTJ;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->mu != 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
* Cholesky decomposition hasn't been computed yet, or was computed
|
|
Packit |
67cb25 |
* with mu > 0 - recompute Cholesky decomposition of J^T J
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* copy lower triangle of JTJ to workspace */
|
|
Packit |
67cb25 |
gsl_matrix_tricpy('L', 1, state->work_JTJ, state->JTJ);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* compute modified Cholesky decomposition */
|
|
Packit |
67cb25 |
status = gsl_linalg_mcholesky_decomp(state->work_JTJ, state->perm, NULL);
|
|
Packit |
67cb25 |
if (status)
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
status = gsl_linalg_mcholesky_rcond(state->work_JTJ, state->perm, &rcond_JTJ, state->work3p);
|
|
Packit |
67cb25 |
if (status == GSL_SUCCESS)
|
|
Packit |
67cb25 |
*rcond = sqrt(rcond_JTJ);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* solve: (J^T J + mu D^T D) x = b */
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
cholesky_solve_rhs(const gsl_vector * b, gsl_vector *x, cholesky_state_t *state)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status;
|
|
Packit |
67cb25 |
gsl_matrix *JTJ = state->work_JTJ;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
status = gsl_linalg_mcholesky_solve(JTJ, state->perm, b, x);
|
|
Packit |
67cb25 |
if (status)
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* A <- A + mu D^T D */
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
cholesky_regularize(const double mu, const gsl_vector * diag, gsl_matrix * A,
|
|
Packit |
67cb25 |
cholesky_state_t * state)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
(void) state;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (mu != 0.0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t i;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 0; i < diag->size; ++i)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
double di = gsl_vector_get(diag, i);
|
|
Packit |
67cb25 |
double *Aii = gsl_matrix_ptr(A, i, i);
|
|
Packit |
67cb25 |
*Aii += mu * di * di;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static const gsl_multifit_nlinear_solver cholesky_type =
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
"cholesky",
|
|
Packit |
67cb25 |
cholesky_alloc,
|
|
Packit |
67cb25 |
cholesky_init,
|
|
Packit |
67cb25 |
cholesky_presolve,
|
|
Packit |
67cb25 |
cholesky_solve,
|
|
Packit |
67cb25 |
cholesky_rcond,
|
|
Packit |
67cb25 |
cholesky_free
|
|
Packit |
67cb25 |
};
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
const gsl_multifit_nlinear_solver *gsl_multifit_nlinear_solver_cholesky = &cholesky_type;
|