XLPack 7.0
XLPack Numerical Library (Excel VBA) Reference Manual
Loading...
Searching...
No Matches

◆ Lmder1_r()

Sub Lmder1_r ( M As  Long,
N As  Long,
X() As  Double,
Fvec() As  Double,
Fjac() As  Double,
Tol As  Double,
Ipvt() As  Long,
Info As  Long,
XX() As  Double,
YY() As  Double,
IRev As  Long,
Optional Info2 As  Long 
)

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

Purpose
This routine 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.

Lmder1_r is equivalent to using Lmder_r with setting FTol = Tol, XTol = Tol, GTol = 0, Maxfev = 100*(n+1), Mode = 1, Factor = 100 and Nprint = 0.
Parameters
[in]MNumber of functions. (M > 0)
[in]NNumber of variables. (0 < N <= M)
[in,out]X()Array X(LX - 1) (LX >= N)
[in] An initial estimate of the solution vector.
[out] IRev = 0: The obtained solution vector.
  IRev = 2: The abscissa where the Jacobian shoule be evaluated.
[out]Fvec()Array Fvec(LFvec - 1) (LFvec >= M)
IRev = 0: The function values evaluated at the solution vector X().
[in,out]Fjac()Array Fjac(LFjac1 - 1, LFjac2 - 1) (LFjac1 >= M, LFjac2 >= N)
[in] IRev = 2: 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]TolRelative error desired in the sum of squares and the approximate solution. (Tol >= 0)
[out]Ipvt()Array Ipvt(LIpvt - 1) (LIpvt >= N)
The pivot indices that define the permutation matrix P.
[out]Info= 0: Successful exit. (Sub-code is set to Info2)
= -1: The argument M had an illegal value. (M < N)
= -2: The argument N had an illegal value. (N <= 0)
= -3: The argument X() is invalid. (Array X() is not big enough)
= -4: The argument Fvec() is invalid. (Array Fvec() is not big enough)
= -5: The argument Fjac() is invalid. (Array Fjac() is not big enough)
= -6: The argument Tol had an illegal value. (Tol < 0)
= -7: The argument Ipiv() is invalid. (Array Ipiv() is not big enough)
= -9: The argument XX() is invalid. (Array XX() is not big enough)
= -10: The argument YY() is invalid. (Array YY() is not big enough)
= 1: Number of function evaluations with IRev = 1 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 - 1) (LXX >= N)
When returned with IRev = 1, XX() contains the abscissa where the function value should be evaluated and given in the next call.
[in]YY()Array YY(LYY - 1) (LYY >= M)
When returned with IRev = 1, 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, set the required values to the specified variables or display the intermediate result as follows and call the routine again.
= 0: Computation finished. See return code in Info.
= 1: User should set the function values at XX() in YY(). Do not alter any variables other than YY().
= 2: User should set the Jacobian at X() in Fjac(). Do not alter any variables other than Fjac().
[out]Info2(Optional)
Sub-code on return with Info = 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.
Reference
netlib/minpack
Example Program
Approximate the following data with model function f(x) = c1*(1 - exp(-c2*x)). Two parameters c1 and c2 are determined by the nonlinear least squares method.
f(x) x
10.07 77.6
29.61 239.9
50.76 434.8
81.78 760.0
The initial estimates of the solution are c1 = 500 and c2 = 0.0001.
Sub FLmder(M As Long, N As Long, X() As Double, Fvec() As Double, Fjac() As Double, IFlag As Long)
Dim Xdata(3) As Double, Ydata(3) As Double, I As Long
Ydata(0) = 10.07: Xdata(0) = 77.6
Ydata(1) = 29.61: Xdata(1) = 239.9
Ydata(2) = 50.76: Xdata(2) = 434.8
Ydata(3) = 81.78: Xdata(3) = 760
If IFlag = 1 Then
For I = 0 To M - 1
Fvec(I) = Ydata(I) - X(0) * (1 - Exp(-Xdata(I) * X(1)))
Next
ElseIf IFlag = 2 Then
For I = 0 To M - 1
Fjac(I, 0) = Exp(-Xdata(I) * X(1)) - 1
Fjac(I, 1) = -Xdata(I) * X(0) * Exp(-X(1) * Xdata(I))
Next
End If
End Sub
Sub Ex_Lmder1_r()
Const M = 4, N = 2
Dim X(N - 1) As Double, Fvec(M - 1) As Double, Fjac(M - 1, N - 1) As Double
Dim Tol As Double, Ipvt(N - 1) As Long, Info As Long
Dim XX(N - 1) As Double, YY(M - 1) As Double, IRev As Long, IFlag As Long
Tol = 0.00000001 '1.0e-8
X(0) = 500: X(1) = 0.0001
IRev = 0
Do
Call Lmder1_r(M, N, X(), Fvec(), Fjac(), Tol, Ipvt(), Info, XX(), YY(), IRev)
If IRev = 1 Or IRev = 2 Then
IFlag = IRev
Call FLmder(M, N, XX(), YY(), Fjac(), IFlag)
End If
Loop While IRev <> 0
Debug.Print "C1, C2 =", X(0), X(1)
Debug.Print "Info =", Info
End Sub
Sub Lmder1_r(M As Long, N As Long, X() As Double, Fvec() As Double, Fjac() As Double, Tol As Double, Ipvt() As Long, Info As Long, XX() As Double, YY() As Double, IRev As Long, Optional Info2 As Long)
Nonlinear least squares approximation by Levenberg-Marquardt method (simple driver) (reverse communic...
Example Results
C1, C2 = 241.084896089973 5.44942234119956E-04
Info = 0