|
|
◆ n2f()
| void n2f |
( |
int |
m, |
|
|
int |
n, |
|
|
double |
x[], |
|
|
void(*)(int, int, double *, int *, double *) |
calcr, |
|
|
void(*)(int, double *, int, int, int, double) |
itsum, |
|
|
double |
v[], |
|
|
int |
lv, |
|
|
int |
iv[], |
|
|
int |
liv, |
|
|
int * |
info |
|
) |
| |
Nonlinear least squares approximation (adaptive algorithm) (Jacobian not required)
- Purpose
- n2f minimizes the sum of the squares of m nonlinear functions in n variables by the adaptive algorithm which combines and augments a Gauss-Newton, Levenberg-Marquardt and other techniques for better convergence.
minimize the sum of fi(x1, x2, ..., xn)^2 (sum for i = 1 to m)
The user must provide a subroutine which calculates the function values. Since the Jacobian is calculated by finite difference approximation within the routine, the user is not required to calculate the Jacobian.
- Parameters
-
| [in] | m | Number of data. (m > 0) |
| [in] | n | Number of parameters. (0 < n <= m) |
| [in,out] | x[] | Array x[lx] (lx >= n)
[in] An initial estimate of the solution vector.
[out] The obtained solution vector. |
| [in] | calcr | The user supplied subroutine which calculates the function values fi(x) defined as follows. void calcr(int m, int n, double x[], int *nf, double r[])
{
Computes residual value fi from given x[] and returns in r[i-1] (i = 1 to m)
}
nf is invocation counter. For same x[] value, same nf will be passed to calcr. If given x[] is out of bounds, nf should be set to 0. The other variables than r[] and nf should not be altered. |
| [in] | itsum | The user supplied subroutine to print the intermediate results defined as follows. If iv[18] = 1, the subroutine is called after every iteration. void itsum(int n, double x[], int niter, int nf, int ng, double fval)
{
Output the following information in desired format.
n: Number of variables
x[]: Current approximation of the solution vector
niter: Iteration counter
nf: Number of function calls of calcr
ng: Number of function calls of calcj
fval: Function value at x[]
}
Argument values should not be altered. |
| [in,out] | v[] | Array v[lv]
The array for the floating point parameters and work space.
[in]
v[25] (tuner1): Parameter to check for false convergence (0 <= tuner1 <= 0.5) (default = 0.1)
v[30] (atol): Absolute function convergence tolerance (0 <= atol) (default = max(1e-20, eps^2), where eps is the machine epsilon)
v[31] (rtol): Relative function convergence tolerance (eps <= rtol <= 0.1) (default = max(1e-10, eps^(2/3))
v[32] (xctol): X-convergence tolerance (0 <= xctol <= 1) (default = eps^(1/2))
v[33] (xftol): False convergence tolerance (0 <= xftol <= 1) (default = 100*eps)
v[34] (lmax0): Maximum 2-norm allowed for scaled very first step (0 < lmax0) (default = 1)
v[35] (lmaxs), v[36] (sctol): Parameters to test for singular convergence (0 < lmaxs, 0 <= sctol <= 1) (default: lmaxs = 1, sctol = max(1e-10, eps^(2/3))).
If the function reduction predicted for a step of length bounded by lmaxs is less than sctol*abs(f), returns with info = 7 (f is the function value at the start of the current iteration).
v[37] (dinit): If nonnegative, all components of the scale vector d[] is initialized to this value (-10 <= dinit) (default = 0)
v[38] (dtol): Tolerance for adaptive scaling (0 =< dtol) (default = 1e-6)
v[39] (d0): Initial value for adaptive scaling (0 <= d0) (default = 1)
v[40] (dfac): Factor for adaptive scaling (0 <= dfac <= 1) (default = 0.6)
A scale factor d[i] is chosen by adaptive scaling so that d[i]*x[i] has about the same magnitude for all i.
Let d1[i] = max(||Ji||, dfac*d[i]) where ||Ji|| is the 2-norm of the i-th column of Jacobian matrix, then d[i] is chosen as follows.
if d1[i] >= dtol then d[i] = d1[i]
if d1[i] < dtol then d[i] = d0
v[41] (dltfdc): Step size used to compute a covariance matrix is chosen as follows (in the case that covreq = -1 or -2)
dltfdc*max(|x[i]|, 1/d[i]) (eps <= dltfdc <= 1) (default = eps^(1/3))
v[42] (dltfdj): Step size used to compute finite difference Jacobian approximation is chosen as follows
dltfdj*max(|x[i]|, 1/d[i]) (eps <= dltfdj <= 1) (default = eps^(1/2))
v[43] (delta0): Step size used to compute a covariance matrix is chosen as follows (in the case that covreq = 1 or 2)
delta0*max(|x[i]|, 1/d[i])*sign(x[i]) (eps <= delta0 <= 1) (default = eps^(1/2))
[out]
v[0] (dgnorm): 2-norm of diag(d)^(-1)*g (g is most recently computed gradient)
v[1] (dstnorm): 2-norm of diag(d)^(-1)*s (s is the current step)
v[9] (f): Current function value (half the sum of squares)
v[12] (f0): Function value at the start of the current iteration |
| [in] | lv | The length of array v[]. (lv >= 105 + n*(m + 2*n + 17) + 2*m) |
| [in,out] | iv[] | Array iv[liv]
The array for the integer parameters and work space.
[in]
iv[0]: If = 0, all parameters in v[] and iv[] will be initialized to the default values before starting calculation.
If = 12, v[] and iv[] are assumed to have already been set by the user and will not be initialized by the routine. Since the subroutine ivset assigns iv[0] = 12, user can first call ivset to set default values to v[] and iv[], and then change some necessary entries to non-default values and start calculation with such non-default parameters.
iv[14] (covreq): Which covariance matrix of the form is to be computed (default = 1)
= 1, -1: sigma * H^(-1) * (J^T * J) * H^(-1)
= 2, -2: sigma * H^(-1)
= 3, -3: sigma * (J^T * J)
where, sigma = rssq / max(1, m-n) (rqqs is residual sum of squares), H is the Hessian matrix and J is the Jacobian matrix.
The lower triangle of the covariance matrix is stored rowwise in v[] starting at v[iv[25]-1].
The sign shows which user's routines are to be used.
< 0: Compute using only calcr
> 0: Compute using both calcr and calcj (if calcj is not available, only calcr is to be used)
iv[15] (dtype): Choice of adaptive scaling (default = 1)
= 0: Disable adaptive scaling (scale factor = 1)
= 1: Enable adaptive scaling during all iterations
= 2: Enable adaptive scaling during the first iteration and scale factor is left unchanged thereafter
iv[16] (mxfcal): Maximum number of function evaluations allowed (default = 200)
iv[17] (mxiter): Maximum number of iterations allowed (default = 150)
iv[18] (outlev): Controls the print of intermediate results (default = 0)
!= 1: itsum is not referenced
= 1: itsum is called after every iteration
iv[24] (inits): The regression routines use a "secant update" to obtain an approximation S to part of Hessian matrix. Usually initial S matrix can be set to all zeros, but occasionally it is useful to initialize S to other values (default = 0)
= 0: S matrix to initialize to all zeros
= 3: To use differences of function values to estimate the initial S
= 4: To use differences of gradients to estimate the initial S
iv[56] (rdreq): Controls computing of the covariance matrix and the regression diagnostics (default = 0)
= 0: Neither to be computed
= 1: Only the covariance matrix to be computed
= 2: Only the regression diagnostics to be computed
= 3: Both to be computed
[out] Output parameters
iv[0]: Return code
= 3: x-convergence
= 4: Relative function convergence
= 5: Both x- and relative function convergence
= 6: Absolute function convergence
= 7: Singular convergence. The Hessian near the current iterate appears to be singular
= 8: False convergence. The iterates appear to be converging to a noncritical point
= 9: Function evaluation limit reached
= 10: Iteration limit reached
= 11: Stopx returned true (external interrupt)
= 12: iv[] and v[] have been allocated and initialized
= 13: f(x) cannot be computed at the initial x[]
= 14: Bad parameters passed to assess (which should not occur)
= 15: liv was too small
= 16: lv was too small
= 17: Restart attempted with m or n changed
= 18: iv[24] is out of range.
= 19...45: v[iv[0]] is out of range
= 50: iv[0] is out of range
= 87...(86+m) = jtol[iv[0]-86] is not positive
iv[5] (nfcall): Number of calls so far made on calcr including those used in computing the covariance matrix but not including those used in computing the Jacobian approximation
iv[25] (covmat): Whether a covariance matrix was computed
= -2, -1: Failed to compute a covariance matrix
= 0: A covariance matrix was not computed
> 0: The lower triangle of the covariance matrix is stored rowwise in v[] starting at v[iv[25]-1])
iv[29] (ngcall): Number of calls so far made on calcr used in computing the Jacobian approximation
iv[30] (niter): Number of iterations performed
iv[43] (lastiv): The least acceptable value of liv
iv[44] (lastv): The least acceptable value of lv
iv[51] (nfcov): Number of calls made on calcr when trying to compute the covariance matrix
iv[52] (ngcov): Number of calls made on calcr used in computing the Jacobian approximation when trying to compute the covariance matrix
iv[60] (r): The residual vector r corresponding to x is stored in v[] starting at v[iv[60]-1]
iv[66] (regd): Whether the regression diagnostics were computed
= -2, -1: Failed to compute the regression diagnostics
= 0: No regression diagnostics computation attempted
> 0: The regression diagnostics were stored in v[] starting at v[iv[66]-1]) |
| [in] | liv | The length of array iv[]. (liv >= 82 + n) |
| [out] | info | = 0: Successful exit (iv[0] = 3 to 6)
= -1: The argument m had an illegal value (m < 1) (iv[0] = 66)
= -2: The argument n had an illegal value (n < 1) (iv[0] = 66)
= -6: The argument v[18], ... or v[49] had an illegal value (out of range) (iv[0] = 19 to 50)
= -7: The argument lv had an illegal value (lv too small) (iv[0] = 16)
= -8: The argument iv[0] had an illegal value (iv[0] out of range) (iv[0] = 80)
= -9: The argument liv had an illegal value (liv too small) (iv[0] = 15)
= 7: Singular convergence (the Hessian near the current iterate appears to be singular)
= 8: False convergence (the iterates appear to be converging to a noncritical point)
= 9: Function evaluation limit reached
= 10: Iteration limit reached
= 14: iv[] and v[] have been allocated (normal exit after a call with iv[0]=13)
= 17: Restart attempted with m or n changed
= 63: f(x) cannot be computed at the initial x
= 65: The gradient could not be computed at x |
- Reference
- netlib/port
|