SYNOPSIS

include <minpack.h> void lmdif1_ ( void (*fcn)(int *m, int *n, double *x, double *fvec, int *iflag),

int *m, int * n, double *x, double *fvec,

double *tol, int *info, int *iwa, double *wa, int *lwa);

void lmdif_ ( void (*fcn)(int *m, int *n, double *x, double *fvec, int *iflag),

int *m, int *n, double *x, double *fvec,

double *ftol, double *xtol, double *gtol, int *maxfev, double *epsfcn, double *diag, int *mode, double *factor, int *nprint, int *info, int *nfev, double *fjac,

int *ldfjac, int *ipvt, double *qtf,

double *wa1, double *wa2, double *wa3, double *wa4 );

DESCRIPTION

The purpose of lmdif_ is to minimize the sum of the squares of m nonlinear functions in n variables by a modification of the Levenberg-Marquardt algorithm. The user must provide a subroutine which calculates the functions. The Jacobian is then calculated by a forward-difference approximation.

lmdif1_ serves the same purpose but has a simplified calling sequence.

Language notes

These functions are written in FORTRAN. If calling from C, keep these points in mind:

Name mangling.

With gfortran, all the function names end in an underscore.

Compile with gfortran.

Even if your program is all C code, you should link with gfortran so it will pull in the FORTRAN libraries automatically. It's easiest just to use gfortran to do all the compiling. (It handles C just fine.)

Call by reference.

All function parameters must be pointers.

Column-major arrays.

Suppose a function returns an array with 5 rows and 3 columns in an array z and in the call you have declared a leading dimension of 7. The FORTRAN and equivalent C references are:

	z(1,1)		z[0]
	z(2,1)		z[1]
	z(5,1)		z[4]
	z(1,2)		z[7]
	z(1,3)		z[14]
	z(i,j)		z[(i-1) + (j-1)*7]

fcn is the name of the user-supplied subroutine which calculates the functions. In FORTRAN, fcn must be declared in an external statement in the user calling program, and should be written as follows:

  subroutine fcn(m,n,x,fvec,iflag)
  integer m,n,iflag
  double precision x(n),fvec(m)
  ----------
  calculate the functions at x and
  return this vector in fvec.
  ----------
  return
  end

In C, fcn should be written as follows:

  void fcn(int *m, int *n, double *x, double *fvec, int *iflag)
  {
  	  /* calculate the functions at x and return
  	     the values in fvec[0] through fvec[m-1] */
  }

The value of iflag should not be changed by fcn unless the user wants to terminate execution of lmdif_ (or lmdif1_). In this case set iflag to a negative integer.

Parameters for both \fBlmdif_\fP and \fBlmdif1_\fP

m is a positive integer input variable set to the number of functions.

n is a positive integer input variable set to the number of variables. n must not exceed m.

x is an array of length n. On input x must contain an initial estimate of the solution vector. On output x contains the final estimate of the solution vector.

fvec is an output array of length m which contains the functions evaluated at the output x.

Parameters for \fBlmdif1_\fP

tol is a nonnegative input variable. Termination occurs when the algorithm estimates either that the relative error in the sum of squares is at most tol or that the relative error between x and the solution is at most tol.

info is an integer output variable. if the user has terminated execution, info is set to the (negative) value of iflag. see description of fcn. otherwise, info is set as follows.

  info = 0  improper input parameters.
  info = 1  algorithm estimates that the relative error

in the sum of squares is at most tol.

  info = 2  algorithm estimates that the relative error

between x and the solution is at most tol.

  info = 3  conditions for info = 1 and info = 2 both hold.
  info = 4  fvec is orthogonal to the columns of the

Jacobian to machine precision.

  info = 5  number of calls to fcn has reached or

exceeded 200*(n+1).

  info = 6  tol is too small. no further reduction in

the sum of squares is possible.

  info = 7  tol is too small. no further improvement in

the approximate solution x is possible.

iwa is an integer work array of length n.

wa is a work array of length lwa.

lwa is an integer input variable not less than m*n + 5*n + m.

Parameters for \fBlmdif_\fP

ftol is a nonnegative input variable. Termination occurs when both the actual and predicted relative reductions in the sum of squares are at most ftol. Therefore, ftol measures the relative error desired in the sum of squares.

xtol is a nonnegative input variable. Termination occurs when the relative error between two consecutive iterates is at most xtol. Therefore, xtol measures the relative error desired in the approximate solution.

gtol is a nonnegative input variable. Termination occurs when the cosine of the angle between fvec and any column of the Jacobian is at most gtol in absolute value. Therefore, gtol measures the orthogonality desired between the function vector and the columns of the Jacobian.

maxfev is a positive integer input variable. Termination occurs when the number of calls to fcn is at least maxfev by the end of an iteration.

epsfcn is an input variable used in determining a suitable step length for the forward-difference approximation. This approximation assumes that the relative errors in the functions are of the order of epsfcn. If epsfcn is less than the machine precision, it is assumed that the relative errors in the functions are of the order of the machine precision.

diag is an array of length n. If mode = 1 (see below), diag is internally set. If mode = 2, diag must contain positive entries that serve as multiplicative scale factors for the variables.

mode is an integer input variable. If mode = 1, the variables will be scaled internally. If mode = 2, the scaling is specified by the input diag. Other values of mode are equivalent to mode = 1.

factor is a positive input variable used in determining the initial step bound. This bound is set to the product of factor and the euclidean norm of diag*x if the latter is nonzero, or else to factor itself. In most cases factor should lie in the interval (.1,100.). 100. is a generally recommended value.

nprint is an integer input variable that enables controlled printing of iterates if it is positive. In this case, fcn is called with iflag = 0 at the beginning of the first iteration and every nprint iterations thereafter and immediately prior to return, with x and fvec available for printing. If nprint is not positive, no special calls of fcn with iflag = 0 are made.

info is an integer output variable. If the user has terminated execution, info is set to the (negative) value of iflag. See description of fcn. Otherwise, info is set as follows.

  info = 0  improper input parameters.
  info = 1  both actual and predicted relative reductions

in the sum of squares are at most ftol.

  info = 2  relative error between two consecutive iterates

is at most xtol.

  info = 3  conditions for info = 1 and info = 2 both hold.
  info = 4  the cosine of the angle between fvec and any

column of the Jacobian is at most gtol in absolute value.

  info = 5  number of calls to fcn has reached or

exceeded maxfev.

  info = 6  ftol is too small. No further reduction in

the sum of squares is possible.

  info = 7  xtol is too small. No further improvement in

the approximate solution x is possible.

  info = 8 gtol is too small. fvec is orthogonal to

the columns of the Jacobian to machine precision.

nfev is an integer output variable set to the number of calls to fcn.

fjac is an output m by n array. The upper n by n submatrix of fjac contains an upper triangular matrix r with diagonal elements of nonincreasing magnitude such that

         t     t           t
        p *(jac *jac)*p = r *r,

where p is a permutation matrix and jac is the final calculated Jacobian. column j of p is column ipvt(j) (see below) of the identity matrix. The lower trapezoidal part of fjac contains information generated during the computation of r.

ldfjac is a positive integer input variable not less than m which specifies the leading dimension of the array fjac.

ipvt is an integer output array of length n. ipvt defines a permutation matrix p such that jac*p = q*r, where jac is the final calculated Jacobian, q is orthogonal (not stored), and r is upper triangular with diagonal elements of nonincreasing magnitude. Column j of p is column ipvt(j) of the identity matrix.

qtf is an output array of length n which contains the first n elements of the vector (q transpose)*fvec.

wa1, wa2, and wa3 are work arrays of length n.

wa4 is a work array of length m.

RELATED TO lmdif1_…

lmder(3), lmder1(3), lmstr(3), lmstr1(3).

AUTHORS

Jorge More', Burt Garbow, and Ken Hillstrom at Argonne National Laboratory. This manual page was written by Jim Van Zandt <[email protected]>, for the Debian GNU/Linux system (but may be used by others).