table of contents
LMDER_(3) | Library Functions Manual | LMDER_(3) |
NAME¶
lmder_, lmder1_ - minimize the sum of squares of m nonlinear functions, with user supplied JacobianSYNOPSIS¶
include <minpack.h>void
lmder1_ ( void (*fcn) (int *m,
int *n, double *x, double
*fvec, double *fjac,
int *ldfjac, int *iflag),
int *m, int *
n, double *x, double
*fvec, double *fjac, int
*ldfjac,
double *tol, int *info, int
*iwa, double *wa, int
*lwa);
void
lmder_ ( void (*fcn)( int *m,
int *n, double *x, double
*fvec, double *fjac,
int *ldfjac, int *iflag),
int *m, int
*n, double *x, double
*fvec, double *fjac, int
*ldfjac,
double *ftol, double *xtol, double
*gtol, int *maxfev, double
*diag, int *mode,
double *factor, int *nprint, int
*info,
int *nfev, int *njev, int
*ipvt,
double *qtf, double *wa1, double
*wa2, double *wa3, double
*wa4 );
DESCRIPTION¶
The purpose of lmder_ 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 function which calculates the functions and the Jacobian. lmder1_ performs the same function as lmder_ but has a simplified calling sequence. lmstr and lmstr1 also perform the same function but use minimal storage.Language notes¶
These functions are written in FORTRAN. If calling from C, keep these points in mind:- Name mangling.
- With g77 version 2.95 or 3.0, all the function names end in an underscore. This may change with future versions of g77.
- Compile with g77.
- Even if your program is all C code, you should link with g77 so it will pull in the FORTRAN libraries automatically. It's easiest just to use g77 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]
User-supplied Function¶
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,fjac,ldfjac,iflag) integer m,n,iflag double precision x(n),fvec(m),fjac(ldfjac,n) ---------- if iflag = 1 calculate the functions at x and return this vector in fvec. do not alter fjac. if iflag = 2 calculate the jacobian at x and return this matrix in fjac. do not alter fvec. ---------- return end
void fcn(int m, int n, double *x, double *fvec, double *fjac, int *ldfjac, int *iflag) { /* if iflag = 1 calculate the functions at x and return this vector in fvec[0] through fvec[m-1]. do not alter fjac. if iflag = 2 calculate the jacobian at x and return this matrix in fjac. do not alter fvec. */ }
Parameters for both lmder_ and lmder1_¶
m is a positive integer input variable set to the number of functions.t t t
p *(jac *jac)*p = r *r,
Parameters for lmder1_¶
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 = 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.
Parameters for lmder_¶
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.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.
SEE ALSO¶
lmdif(3), lmdif1(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 <jrv@debian.org>, for the Debian GNU/Linux system (but may be used by others).March 8, 2002 | Minpack |