XLPack 7.0
XLPack Numerical Library (C API) Reference Manual
Loading...
Searching...
No Matches

◆ lmder_r()

void lmder_r ( int  m,
int  n,
double  x[],
double  fvec[],
int  ldfjac,
double  fjac[],
double  ftol,
double  xtol,
double  gtol,
int  maxfev,
double  diag[],
int  mode,
double  factor,
int  nprint,
int *  nfev,
int *  njev,
int  ipvt[],
double  work[],
int  lwork,
int *  info,
double  xx[],
double  yy[],
int *  irev 
)

Nonlinear least squares approximation by Levenberg-Marquardt method (reverse communication version)

Purpose
lmder_r minimizes the sum of the squares of m nonlinear functions in n variables by a modification of the Levenberg-Marquardt algorithm.
minimize the sum of fi(x1, x2, ..., xn)^2 (sum for i = 1 to m)
The user must provide the calculated function values and the Jacobian according to irev.
Parameters
[in]mNumber of functions. (m > 0)
[in]nNumber of variables. (0 < n <= m)
[in,out]x[]Array x[lx] (lx >= n)
[in] An initial estimate of the solution vector.
[out] irev = 0: The obtained solution vector.
  irev = 30: The abscissa where the Jacobian shoule be evaluated.
  irev = 50, 51: Recent approximation of the solution vector.
[out]fvec[]Array fvec[lfvec] (lfvec >= m)
irev = 0: The function values evaluated at the solution vector x[].
irev = 50, 51: The function values evaluated at the recent approximation of the solution vector.
[in]ldfjacLeading dimension of two dimensional array fjac[]. (ldfjac >= m)
[in,out]fjac[][]Array fjac[lfjac][ldfjac] (lfjac >= n)
[in] irev = 30: Jacobian(dfi/dxj) at x[] should be given in fjac[][] in the next call.
[out] irev = 0: The upper n x n submatrix contains an upper triangular matrix R with diagonal elements of nonincreasing magnitude such that
    P^T * (J^T * J)*P = R^T * R
  where P is a permutation matrix and J is the final calculated Jacobian.
[in]ftolRelative error desired in the sum of squares. Termination occurs when both the actual and predicted relative reductions in the sum of squares are at most ftol. (ftol >= 0)
[in]xtolRelative error desired in the approximate solution. Termination occurs when the relative error between two consecutive iterates is at most xtol. (xtol >= 0)
[in]gtolOrthogonality desired between the function vector and the columns of the Jacobian. Termination occurs when the cosine of the angle between fvec[] and any column of the jacobian is at most gtol in absolute value. (gtol >= 0)
[in]maxfevTermination occurs when the number of function evaluations with irev = 1 to 3 has reached this value (maxfev > 0)
[in,out]diag[]Array diag[ldiag] (ldiag >= n)
[in] If mode = 2, diag must contain positive entries that serve as multiplicative scale factors for the variables. (diag[i] > 0)
[out] If mode = 1, diag[] is set by the subroutine.
[in]mode= 1: The variables will be automatically scaled by the subroutine.
= 2: The scaling is specified by the input diag[].
(For other values, mode = 1 is assumed.)
[in]factorUsed in determining the initial step bound. This bound is set to the product of factor and the Euclidean norm of diag*x if nonzero, or else to factor itself. In most cases factor should lie in the interval (0.1, 100). 100 is a generally recommended value. (factor > 0)
[in]nprint> 0: Returns with irev = 50 or 51 at the beginning of the first iteration and every nprint iterations thereafter and at the end of last iteration for printing intermediate result.
<= 0: Does not return for printing intermediate result.
[out]nfevNumber of function evaluations with irev = 1 or 2.
[out]njevNumber of Jacobian evaluations with irev = 30.
[out]ipvt[]Array ipvt[lipvt] (lipvt >= n)
The pivot indices that define the permutation matrix P.
[out]work[]Array work[lwork]
Work array.
On return with info = 0, sub-code is set to work[0].
= 1: Both actual and predicted relative reductions in the sum of squares are at most ftol.
= 2: Relative error between two consecutive iterates is at most xtol.
= 3: Both of above are satisfied.
= 4: The cosine of the angle between fvec and any column of the Jacobian is at most gtol in absolute value.
[in]lworkThe length of work[]. (lwork >= 4*n + m)
[out]info= 0: Successful exit (sub-code is set to work[0])
= -1: The argument m had an illegal value (m < n)
= -2: The argument n had an illegal value (n <= 0)
= -5: The argument ldfjac had an illegal value (ldfjac < m)
= -7: The argument ftol had an illegal value (ftol < 0)
= -8: The argument xtol had an illegal value (xtol < 0)
= -9: The argument gtol had an illegal value (gtol < 0)
= -10: The argument maxfev had an illegal value (maxfev <= 0)
= -11: The argument diag had an illegal value (diag[i] <= 0 when mode = 2)
= -13: The argument factor had an illegal value (factor <= 0)
= -19: The argument lwork had an illegal value (lwork too small)
= 1: Number of function evaluations with irev=1 or 2 has reached maxfev
= 2: ftol is too small. No further reduction in the sum of squares is possible
= 3: xtol is too small. No further improvement in the approximate solution x is possible
= 4: gtol is too small. fvec is orthogonal to the columns of the Jacobian to machine precision
[out]xx[]Array xx[lxx] (lxx >= n)
When returned with irev = 1 or 2, xx[] contains the abscissa where the function value should be evaluated and given in the next call.
[in]yy[]Array yy[lyy] (lyy >= m)
When returned with irev = 1 or 2, the function value fi(xx[]) (i = 1 to m) should be given in yy[] in the next call.
[in,out]irevControl variable for reverse communication.
[in] Before first call, irev should be initialized to zero. On succeeding calls, irev should not be altered.
[out] If irev is not zero, complete the following process and call this routine again.
= 0: Computation finished. See return code in info.
= 1 or 2: User should set the function values at xx[] in yy[]. Do not alter any variables other than yy[].
= 30: User should set the Jacobian at x[] in fjac[][]. Do not alter any variables other than fjac[][].
= 50 or 51: Display the intermediate result (x[], fvec[], etc.) (in the case of nprint > 0). Do not alter any variables.
Reference
netlib/minpack