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

◆ Lmdif()

Sub Lmdif ( F As  LongPtr,
M As  Long,
N As  Long,
X() As  Double,
Fvec() As  Double,
Fjac() As  Double,
FTol As  Double,
XTol As  Double,
GTol As  Double,
Diag() As  Double,
Mode As  Long,
Ipvt() As  Long,
Info As  Long,
Optional Info2 As  Long,
Optional Maxfev As  Long = 0,
Optional Epsfcn As  Double = 0,
Optional Factor As  Double = 0,
Optional Nprint As  Long = 0,
Optional Nfev As  Long 
)

Nonlinear least squares approximation by Levenberg-Marquardt method (Jacobian not required)

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 a subroutine which calculates the functions. Since the Jacobian is calculated by a forward-difference approximation within the routine, the user is not required to provide the Jacobian.
Parameters
[in]FUser supplied subroutine which calculates the functions fi(x) defined as follows.
Sub F(M As Long, N As Long, X() As Double, Fvec() As Double, IFlag As Long)
If IFlag = 1 or 2:
Calculate the function values fi(x) at X() and return in Fvec(i-1) (i = 1 to M). Other variables should not be changed.
If IFlag = 0:
If Nprint is set to positive value, F is called every Nprint iterations with IFlag = 0 for printing intermediate result X() and Fvec(). Any variable should not be changed.
End Sub
The value of IFlag should not be changed unless the user wants to terminate the execution. In this case, set IFlag to a negative integer.
[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] The obtained solution vector.
[out]Fvec()Array Fvec(LFvec - 1) (LFvec >= M)
The function values evaluated at the solution vector X().
[out]Fjac()Array Fjac(LFjac1 - 1, LFjac2 - 1) (LFjac1 >= M, LFjac2 >= N)
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,out]Diag()Array Diag(LD - 1) (LD >= 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.)
[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)
= -2: The argument M had an illegal value. (M < N)
= -3: The argument N had an illegal value. (N <= 0)
= -4: The argument X() is invalid. (Array X() is not big enough)
= -5: The argument Fvec() is invalid. (Array Fvec() is not big enough)
= -6: The argument Fjac() is invalid. (Array Fjac() is not big enough)
= -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 Diag() is invalid. (Array Diag() is not big enough) or had an illegal value (Diag(i) < 0 when Mode = 2)
= -12: The argument Ipiv() is invalid. (Array Ipiv() is not big enough)
= 1: Number of calls to F with IFlag = 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.
= 5: User imposed termination (returned from F with IFlag < 0).
[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.
= 4: The cosine of the angle between Fvec and any column of the Jacobian is at most GTol in absolute value.
[in]Maxfev(Optional)
Termination occurs when the number of calls to F with IFlag = 1 or 2 has reached this value. (default = 200*(N + 1))
(if Maxfev <= 0, the default value will be used.)
[in]Epsfcn(Optional)
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. (default = 0)
[in]Factor(Optional)
Used 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. (default = 100)
(if Factor <= 0, the default value will be used)
[in]Nprint(Optional)
> 0: F is called with IFlag = 0 at the beginning of the first iteration and every Nprint iterations thereafter and at the end of last iteration for printing intermediate result.
<= 0: No special calls of F for printing intermediate results are made.
(default = 0)
[out]Nfev(Optional)
Number of function evaluations (number of calls to F with IFlag = 1 or 2).
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 FLmdif(M As Long, N As Long, X() As Double, Fvec() 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 Or IFlag = 2 Then
For I = 0 To M - 1
Fvec(I) = Ydata(I) - X(0) * (1 - Exp(-Xdata(I) * X(1)))
Next
End If
End Sub
Sub Ex_Lmdif()
Const M = 4, N = 2
Dim X(N - 1) As Double, Fvec(M - 1) As Double, Fjac(M - 1, N - 1) As Double
Dim FTol As Double, XTol As Double, GTol As Double, Diag(N - 1) As Double
Dim Mode As Long, Ipvt(N - 1) As Long, Info As Long
FTol = 0.00000001: XTol = 0.00000001: GTol = 0.00000001 '1.0e-8
X(0) = 500: X(1) = 0.0001
Call Lmdif(AddressOf FLmdif, M, N, X(), Fvec(), Fjac(), FTol, XTol, GTol, Diag(), Mode, Ipvt(), Info)
Debug.Print "C1, C2 =", X(0), X(1)
Debug.Print "Info =", Info
End Sub
Sub Lmdif(F As LongPtr, M As Long, N As Long, X() As Double, Fvec() As Double, Fjac() As Double, FTol As Double, XTol As Double, GTol As Double, Diag() As Double, Mode As Long, Ipvt() As Long, Info As Long, Optional Info2 As Long, Optional Maxfev As Long=0, Optional Epsfcn As Double=0, Optional Factor As Double=0, Optional Nprint As Long=0, Optional Nfev As Long)
Nonlinear least squares approximation by Levenberg-Marquardt method (Jacobian not required)
Example Results
C1, C2 = 241.084897132837 5.44942231312974E-04
Info = 0