HTML automatically generated with rman
Table of Contents

Name

nllsqfit, nr_nllsqfit - (non)linear least squares fit

Synopsis


int nllsqfit(xdat, xdim, ydat, wdat, ddat, ndat, 
        fpar, epar, mpar, npar, 
        tol, its, lab, f, df)
int nr_nllsqfit(xdat, xdim, ydat, wdat, ddat, ndat, 
        fpar, epar, mpar, npar, 
        tol, its, lab, f, df)
  real *xdat, *ydat, *wdat, *ddat, *fpar, *epar, tol, lab;
  int  xdim, ndat, *mpar, npar, its;
  rproc f;
  iproc df;

Description

nllsqfit is a routine for making a least-squares fit of a function to a set of data points. The method used is described in: Marquardt, J.Soc.Ind.Appl.Math. 11, 431 (1963). (see also Numerical Recipes, Ch. 14) This method is a mixture of the steepest descent method and the Taylor method.

nllsqfit returns number of iterations needed to achieve convergence according to tol. When this number is negative, the fitting was not continued because a fatal error occurred:

    -1     Too many free parameters, maximum "#define MAXPAR 32".
    -2     No free parameters.
    -3     Not enough degrees of freedom.
    -4     Too many iterations to get a solution which satisfies tol.
    -5     Diagonal of matrix contains elements which are zero, or less.
    -6     Determinant of the coefficient matrix is zero.
    -7     Square root of negative number.
A linear fit (lab=0) returns 0.

nr_nllsqfit is a wrapper routine with the same calling sequence, but calls the (NEMO adapted) Numerical Recipes routine mrqmin() and its helper functions.

Parameters

xdat
contains the coordinates of the data points. xdat is two-dimensional: xdat(xdim,ndatf) in FORTRAN notation or xdat[ndat][xdim] in C. sense.
xdim
is the dimension of the fit.
ydat
contains the data points.
wdat
contains the weigths for the data points. Can be a NULL pointer, in which case all weights are equal.
ddat
contains the difference between data and fit. Can be a NULL pointer, in which case no fit differences are returned.
ndat
is the number of data points.
fpar
On input contains initial estimates of the parameters for non-linear fits, on output the fitted parameters.
epar
contains estimates of errors in fitted parameters.
mpar
logical mask telling which parameters are free (mpar[j]=non-zero) and which parameters are fixed (mpar[j]=0).
npar
number of parameters (free+fixed).
tol
relative tolerance. nllsqfit stops when successive iterations fail to produce a decrement in reduced chi-squared less than tol. If tol is less than the minimum tolerance possible, tol will be set to this value. This means that maximum accuracy can be obtained by setting tol=0.0.
its
maximum number of iterations.
lab
mixing parameter, lab determines the initial weight of steepest descent method relative to the Taylor method. lab should be a small value (i.e. 0.01). lab can only be zero when the partial derivatives are independent of the parameters. In fact in this case lab should be exactly equal to zero, in which case the fit can (or is assumed to) be done linear.
f    
external function, must return a real value, see below.
df    
external function, returns the partial deriviates to the fitted parameters, see below

Notes

The following routines have to be defined by the user:
      real func(xdat, fpar, npar)
      func          returns the function value of the function to be fitted.
      real xdat[]   (input) coordinate(s) of data point.
      real fpar[]   (input) parameter list.
      int  npar     (input) number of parameters.

void derv(xdat, fpar, dpar, npar)

real xdat[] (input) coordinate(s) of data point.
real fpar[] (input) parameter list.
real dpar[] (output) partial derivatives to the parameters of the
function to be fitted.
int npar (input) number of parameters.

Example

Fitting a straight line y(x) = a * x + b :


    real func(real *xdat, real *fpar, int npar)
    {
        return fpar[0] * (*xdat) + fpar[1];
    }
              
    void derv(real *xdat, real *fpar, real *dpar, int npar)
    {
        dpar[0] = *xdat;
        dpar[1] = 1.0;
    }

Files

nllsqfit.c

See Also

linreg(3NEMO) , matinv(3NEMO) , FIT.DC2(GIPSY)

Author

K.G. Begeman (originally named FIT.SHL, in Sheltran), P.J. Teuben (C)

Copyright

Copyright (c) Kapteyn Laboratorium Groningen 1990; All Rights Reserved.

History


May  7, 1990    Document created(KGB), document refereed(MXV)
Apr 30, 1991    NEMO version written for rotcur, as old PJT
July 23, 1992   manual page written PJT
Aug 20, 1992    turbocharged getvec() considerably  PJT
July 12, 2002    allow ’wdat’ to be a NULL vector if all weights the same    PJT


Table of Contents