|
Packit |
67cb25 |
/* multifit/gsl_multifit.h
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* Copyright (C) 2000, 2007, 2010 Brian Gough
|
|
Packit |
67cb25 |
* Copyright (C) 2013, 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 |
#ifndef __GSL_MULTIFIT_H__
|
|
Packit |
67cb25 |
#define __GSL_MULTIFIT_H__
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
#include <stdlib.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_types.h>
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
#undef __BEGIN_DECLS
|
|
Packit |
67cb25 |
#undef __END_DECLS
|
|
Packit |
67cb25 |
#ifdef __cplusplus
|
|
Packit |
67cb25 |
# define __BEGIN_DECLS extern "C" {
|
|
Packit |
67cb25 |
# define __END_DECLS }
|
|
Packit |
67cb25 |
#else
|
|
Packit |
67cb25 |
# define __BEGIN_DECLS /* empty */
|
|
Packit |
67cb25 |
# define __END_DECLS /* empty */
|
|
Packit |
67cb25 |
#endif
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
__BEGIN_DECLS
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
typedef struct
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t nmax; /* maximum number of observations */
|
|
Packit |
67cb25 |
size_t pmax; /* maximum number of parameters */
|
|
Packit |
67cb25 |
size_t n; /* number of observations in current SVD decomposition */
|
|
Packit |
67cb25 |
size_t p; /* number of parameters in current SVD decomposition */
|
|
Packit |
67cb25 |
gsl_matrix * A; /* least squares matrix for SVD, n-by-p */
|
|
Packit |
67cb25 |
gsl_matrix * Q;
|
|
Packit |
67cb25 |
gsl_matrix * QSI;
|
|
Packit |
67cb25 |
gsl_vector * S;
|
|
Packit |
67cb25 |
gsl_vector * t;
|
|
Packit |
67cb25 |
gsl_vector * xt;
|
|
Packit |
67cb25 |
gsl_vector * D;
|
|
Packit |
67cb25 |
double rcond; /* reciprocal condition number */
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace *
|
|
Packit |
67cb25 |
gsl_multifit_linear_alloc (const size_t n, const size_t p);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
void
|
|
Packit |
67cb25 |
gsl_multifit_linear_free (gsl_multifit_linear_workspace * w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear (const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_matrix * cov,
|
|
Packit |
67cb25 |
double * chisq,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_tsvd (const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
const double tol,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_matrix * cov,
|
|
Packit |
67cb25 |
double * chisq,
|
|
Packit |
67cb25 |
size_t * rank,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_svd (const gsl_matrix * X,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_bsvd (const gsl_matrix * X,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
size_t
|
|
Packit |
67cb25 |
gsl_multifit_linear_rank(const double tol, const gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_solve (const double lambda,
|
|
Packit |
67cb25 |
const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
double *rnorm,
|
|
Packit |
67cb25 |
double *snorm,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_applyW(const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * w,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_matrix * WX,
|
|
Packit |
67cb25 |
gsl_vector * Wy);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_stdform1 (const gsl_vector * L,
|
|
Packit |
67cb25 |
const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_matrix * Xs,
|
|
Packit |
67cb25 |
gsl_vector * ys,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_wstdform1 (const gsl_vector * L,
|
|
Packit |
67cb25 |
const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * w,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_matrix * Xs,
|
|
Packit |
67cb25 |
gsl_vector * ys,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_L_decomp (gsl_matrix * L, gsl_vector * tau);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_stdform2 (const gsl_matrix * LQR,
|
|
Packit |
67cb25 |
const gsl_vector * Ltau,
|
|
Packit |
67cb25 |
const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_matrix * Xs,
|
|
Packit |
67cb25 |
gsl_vector * ys,
|
|
Packit |
67cb25 |
gsl_matrix * M,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_wstdform2 (const gsl_matrix * LQR,
|
|
Packit |
67cb25 |
const gsl_vector * Ltau,
|
|
Packit |
67cb25 |
const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * w,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_matrix * Xs,
|
|
Packit |
67cb25 |
gsl_vector * ys,
|
|
Packit |
67cb25 |
gsl_matrix * M,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_genform1 (const gsl_vector * L,
|
|
Packit |
67cb25 |
const gsl_vector * cs,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_genform2 (const gsl_matrix * LQR,
|
|
Packit |
67cb25 |
const gsl_vector * Ltau,
|
|
Packit |
67cb25 |
const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
const gsl_vector * cs,
|
|
Packit |
67cb25 |
const gsl_matrix * M,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_wgenform2 (const gsl_matrix * LQR,
|
|
Packit |
67cb25 |
const gsl_vector * Ltau,
|
|
Packit |
67cb25 |
const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * w,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
const gsl_vector * cs,
|
|
Packit |
67cb25 |
const gsl_matrix * M,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_lreg (const double smin, const double smax,
|
|
Packit |
67cb25 |
gsl_vector * reg_param);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_lcurve (const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_vector * reg_param,
|
|
Packit |
67cb25 |
gsl_vector * rho, gsl_vector * eta,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_lcorner(const gsl_vector *rho,
|
|
Packit |
67cb25 |
const gsl_vector *eta,
|
|
Packit |
67cb25 |
size_t *idx);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_lcorner2(const gsl_vector *reg_param,
|
|
Packit |
67cb25 |
const gsl_vector *eta,
|
|
Packit |
67cb25 |
size_t *idx);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_Lk(const size_t p, const size_t k, gsl_matrix *L);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_Lsobolev(const size_t p, const size_t kmax,
|
|
Packit |
67cb25 |
const gsl_vector *alpha, gsl_matrix *L,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace *work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_wlinear (const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * w,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_matrix * cov,
|
|
Packit |
67cb25 |
double * chisq,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_wlinear_tsvd (const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * w,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
const double tol,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_matrix * cov,
|
|
Packit |
67cb25 |
double * chisq,
|
|
Packit |
67cb25 |
size_t * rank,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_wlinear_svd (const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * w,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
double tol,
|
|
Packit |
67cb25 |
size_t * rank,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_matrix * cov,
|
|
Packit |
67cb25 |
double *chisq,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_wlinear_usvd (const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * w,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
double tol,
|
|
Packit |
67cb25 |
size_t * rank,
|
|
Packit |
67cb25 |
gsl_vector * c,
|
|
Packit |
67cb25 |
gsl_matrix * cov,
|
|
Packit |
67cb25 |
double *chisq,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_est (const gsl_vector * x,
|
|
Packit |
67cb25 |
const gsl_vector * c,
|
|
Packit |
67cb25 |
const gsl_matrix * cov, double *y, double *y_err);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double
|
|
Packit |
67cb25 |
gsl_multifit_linear_rcond (const gsl_multifit_linear_workspace * w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_residuals (const gsl_matrix *X, const gsl_vector *y,
|
|
Packit |
67cb25 |
const gsl_vector *c, gsl_vector *r);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* gcv.c */
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_gcv_init(const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_vector * reg_param,
|
|
Packit |
67cb25 |
gsl_vector * UTy,
|
|
Packit |
67cb25 |
double * delta0,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_gcv_curve(const gsl_vector * reg_param,
|
|
Packit |
67cb25 |
const gsl_vector * UTy,
|
|
Packit |
67cb25 |
const double delta0,
|
|
Packit |
67cb25 |
gsl_vector * G,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_gcv_min(const gsl_vector * reg_param,
|
|
Packit |
67cb25 |
const gsl_vector * UTy,
|
|
Packit |
67cb25 |
const gsl_vector * G,
|
|
Packit |
67cb25 |
const double delta0,
|
|
Packit |
67cb25 |
double * lambda,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double
|
|
Packit |
67cb25 |
gsl_multifit_linear_gcv_calc(const double lambda,
|
|
Packit |
67cb25 |
const gsl_vector * UTy,
|
|
Packit |
67cb25 |
const double delta0,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_multifit_linear_gcv(const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_vector * reg_param,
|
|
Packit |
67cb25 |
gsl_vector * G,
|
|
Packit |
67cb25 |
double * lambda,
|
|
Packit |
67cb25 |
double * G_lambda,
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace * work);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
typedef struct
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const char * name; /* method name */
|
|
Packit |
67cb25 |
int (*wfun)(const gsl_vector *r, gsl_vector *w);
|
|
Packit |
67cb25 |
int (*psi_deriv)(const gsl_vector *r, gsl_vector *dpsi);
|
|
Packit |
67cb25 |
double tuning_default; /* default tuning constant */
|
|
Packit |
67cb25 |
} gsl_multifit_robust_type;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
typedef struct
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
double sigma_ols; /* OLS estimate of sigma */
|
|
Packit |
67cb25 |
double sigma_mad; /* MAD estimate of sigma */
|
|
Packit |
67cb25 |
double sigma_rob; /* robust estimate of sigma */
|
|
Packit |
67cb25 |
double sigma; /* final estimate of sigma */
|
|
Packit |
67cb25 |
double Rsq; /* R^2 coefficient of determination */
|
|
Packit |
67cb25 |
double adj_Rsq; /* degree of freedom adjusted R^2 */
|
|
Packit |
67cb25 |
double rmse; /* root mean squared error */
|
|
Packit |
67cb25 |
double sse; /* residual sum of squares */
|
|
Packit |
67cb25 |
size_t dof; /* degrees of freedom */
|
|
Packit |
67cb25 |
size_t numit; /* number of iterations */
|
|
Packit |
67cb25 |
gsl_vector *weights; /* final weights */
|
|
Packit |
67cb25 |
gsl_vector *r; /* final residuals y - X c */
|
|
Packit |
67cb25 |
} gsl_multifit_robust_stats;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
typedef struct
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t n; /* number of observations */
|
|
Packit |
67cb25 |
size_t p; /* number of parameters */
|
|
Packit |
67cb25 |
size_t numit; /* number of iterations */
|
|
Packit |
67cb25 |
size_t maxiter; /* maximum iterations */
|
|
Packit |
67cb25 |
const gsl_multifit_robust_type *type;
|
|
Packit |
67cb25 |
double tune; /* tuning parameter */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector *r; /* residuals at current iteration */
|
|
Packit |
67cb25 |
gsl_vector *weights; /* weights at current iteration */
|
|
Packit |
67cb25 |
gsl_vector *c_prev; /* coefficients from previous iteration */
|
|
Packit |
67cb25 |
gsl_vector *resfac; /* multiplicative factors for residuals */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector *psi; /* psi(r) */
|
|
Packit |
67cb25 |
gsl_vector *dpsi; /* psi'(r) */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_matrix *QSI; /* Q S^{-1} of original matrix X */
|
|
Packit |
67cb25 |
gsl_vector *D; /* balancing parameters of original matrix X */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector *workn; /* workspace of length n */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_multifit_robust_stats stats; /* various statistics */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_multifit_linear_workspace *multifit_p;
|
|
Packit |
67cb25 |
} gsl_multifit_robust_workspace;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* available types */
|
|
Packit |
67cb25 |
GSL_VAR const gsl_multifit_robust_type * gsl_multifit_robust_default;
|
|
Packit |
67cb25 |
GSL_VAR const gsl_multifit_robust_type * gsl_multifit_robust_bisquare;
|
|
Packit |
67cb25 |
GSL_VAR const gsl_multifit_robust_type * gsl_multifit_robust_cauchy;
|
|
Packit |
67cb25 |
GSL_VAR const gsl_multifit_robust_type * gsl_multifit_robust_fair;
|
|
Packit |
67cb25 |
GSL_VAR const gsl_multifit_robust_type * gsl_multifit_robust_huber;
|
|
Packit |
67cb25 |
GSL_VAR const gsl_multifit_robust_type * gsl_multifit_robust_ols;
|
|
Packit |
67cb25 |
GSL_VAR const gsl_multifit_robust_type * gsl_multifit_robust_welsch;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_multifit_robust_workspace *gsl_multifit_robust_alloc(const gsl_multifit_robust_type *T,
|
|
Packit |
67cb25 |
const size_t n, const size_t p);
|
|
Packit |
67cb25 |
void gsl_multifit_robust_free(gsl_multifit_robust_workspace *w);
|
|
Packit |
67cb25 |
int gsl_multifit_robust_tune(const double tune,
|
|
Packit |
67cb25 |
gsl_multifit_robust_workspace *w);
|
|
Packit |
67cb25 |
int gsl_multifit_robust_maxiter(const size_t maxiter,
|
|
Packit |
67cb25 |
gsl_multifit_robust_workspace *w);
|
|
Packit |
67cb25 |
const char *gsl_multifit_robust_name(const gsl_multifit_robust_workspace *w);
|
|
Packit |
67cb25 |
gsl_multifit_robust_stats gsl_multifit_robust_statistics(const gsl_multifit_robust_workspace *w);
|
|
Packit |
67cb25 |
int gsl_multifit_robust_weights(const gsl_vector *r, gsl_vector *wts,
|
|
Packit |
67cb25 |
gsl_multifit_robust_workspace *w);
|
|
Packit |
67cb25 |
int gsl_multifit_robust(const gsl_matrix * X, const gsl_vector * y,
|
|
Packit |
67cb25 |
gsl_vector * c, gsl_matrix *cov,
|
|
Packit |
67cb25 |
gsl_multifit_robust_workspace *w);
|
|
Packit |
67cb25 |
int gsl_multifit_robust_est(const gsl_vector * x, const gsl_vector * c,
|
|
Packit |
67cb25 |
const gsl_matrix * cov, double *y, double *y_err);
|
|
Packit |
67cb25 |
int gsl_multifit_robust_residuals(const gsl_matrix * X,
|
|
Packit |
67cb25 |
const gsl_vector * y,
|
|
Packit |
67cb25 |
const gsl_vector * c, gsl_vector * r,
|
|
Packit |
67cb25 |
gsl_multifit_robust_workspace * w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
__END_DECLS
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
#endif /* __GSL_MULTIFIT_H__ */
|