XLPack 6.0
Excel VBA Numerical Library Reference Manual
Loading...
Searching...
No Matches

◆ Mnh_r()

Sub Mnh_r ( N As  Long,
X() As  Double,
Info As  Long,
YY As  Double,
YYp() As  Double,
YYpd() As  Double,
IRev As  Long,
Optional Iout As  Long = 0,
Optional Info2 As  Long,
Optional NFcall As  Long,
Optional NGcall As  Long,
Optional Niter As  Long,
Optional Fval As  Double,
Optional NFGcal As  Long,
Optional Rtol As  Double = -1,
Optional Atol As  Double = -1,
Optional MaxFcall As  Long = 0,
Optional MaxIter As  Long = 0,
Optional Dtype As  Long = 0,
Optional Dfac As  Double = -1,
Optional Dtol As  Double = 0,
Optional D0 As  Double = 0,
Optional Tuner1 As  Double = -1,
Optional Xctol As  Double = -1,
Optional Xftol As  Double = -1,
Optional Lmax0 As  Double = -1,
Optional Lmaxs As  Double = -1,
Optional Sctol As  Double = -1,
Optional Bias As  Double = -1 
)

Minimum of a multivariable nonlinear function (trust region method) (gradient and Hessian computed analytically) (reverse communication version)

Purpose
This routine finds the local minimum point (xs1, xs2, ..., xsn) of general nonlinear function f(x1, x2, ..., xn) (a twice continuously differentiable real-valued function)

The user-supplied analytic gradient and Hessian values of the objective function are used. Steps are computed by the double dogleg trust region method.

Mnh_r is the reverse communication version of Mnh.
Parameters
[in]NThe number of variables. (N > 0)
[in,out]X()Array X(LX - 1) (LX >= N)
[in] Initial approximation of the solution vector.
[out] IRev = 0: Obtained solution vector.
  IRev = 1, 2: The abscissa where the function value or derivatives and Hessian shoule be evaluated.
  IRev = 3: Recent approximation of the solution vector.
[out]Info= 0: Successful exit. (See sub-code in Info2)
= -1: The argument N had an illegal value. (N < 1)
= -2: The argument X() is invalid.
= -5: The argument YYp() is invalid.
= 7: Singular convergence. (Hessian near the current iterate appears to be singular)
= 8: False convergence. (Iterate appears to be converging to a noncritical point. Tolerances may be too small)
= 9: Function evaluation limit reached.
= 10: Iteration limit reached.
= 63: F(X) cannot be computed at the initial X.
= 65: The gradient could not be computed at X.
[in,out]YY[in] If IRev = 1, the function value at X() should be given in YY in the next call.
[out] If IRev = 3, YY contains the function value at current X() for printing.
[in,out]YYp()Array YYp(LYYp - 1) (LYYp >= N)
[in] If IRev = 30, the derivatives at X() should be given in YYp() in the next call.
[out] If IRev = 50, YYp() contains the derivatives at current X() for printing.
[in,out]YYpd()Array YYpd(LYYpd - 1) (LYYpd >= N(N+1)/2)
[in] If IRev = 2, the lower triangle of the Hessian at X() should be given in YYpd() in packed form in the next call.
[out] If IRev = 3, YYpd() contains the lower triangle of the Hessian at current X() for printing.
[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 value at X() in YY. Do not alter any variables other than YY.
= 2: User should set the derivatives at X() in YYp(), and the lower triangle of the Hessian in YYpd(). Do not alter any variables other than YYp() and YYpd().
= 3: To be returned with IRev = 3 on each iteration if Iout = 1. Output the intermediate result (X(), NFcall, NGcall, Niter, Fval, etc.). Do not alter any variables.
[in]Iout(Optional)
Whether the intermediate result output is necessary. (default = 0)
= 0: Output is not necessary (never return with IRev = 3).
= 1: Return on each iteration with IRev = 3 to output the intermediate results.
(If other value is specified , Iout = 0 will be assumed)
[out]Info2(Optional)
Sub-code for Info = 0.
= 1: X convergence.
= 2: Relative function convergence.
= 3: Both X and relative function convergence.
= 4: Absolute function convergence.
[out]NFcall(Optional)
Number of function evaluations of F.
[out]NGcall(Optional)
Number of function evaluations of GH.
[out]Niter(Optional)
Number of iterations.
[out]Fval(Optional)
Function value at the obtained solution vector X().
[out]NFGcal(Optional)
Invocation counter Nf which is used for subroutines F and G in Mnh.
[in]Rtol(Optional)
Relative function convergence tolerance. (Eps <= Rtol <= 0.1) (default = 1e-10) (Eps: machine epsilon)
(If Rtol < Eps or Rtol > 0.1, the default value will be used)
[in]Atol(Optional)
Absolute function convergence tolerance. (default = 1e-20)
(If Atol < 0, the default value will be used)
[in]MaxFcall(Optional)
Maximum number of function evaluations of F. (default = 200)
(If MaxFcall <= 0, the default value will be used)
[in]MaxIter(Optional)
Maximum number of iterations. (default = 150)
(If MaxIter <= 0, the default value will be used)
[in]Dtype(Optional)
Choice of adaptive scaling. (Dtype = 0, 1 or 2) (default = 0)
= 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.
(For other values, the default value will be used)
[in]Dfac(Optional)
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(Sqr(Abs(H(i,i))), Dfac*D(i)) where H(i,i) is the diagonal element of Hessian matrix, then D(i) is chosen as follows.
  if D1(i) >= Dtol: D(i) = D1(i)
  if D1(i) < Dtol: D(i) = max(D0, Dtol)
(If Dfac < 0 or Dfac > 1, the default value will be used)
[in]Dtol(Optional)
Tolerance for adaptive scaling. (Dtol > 0) (default = 1.0e-6)
(If Dtol <= 0, the default value will be used)
[in]D0(Optional)
Initial value for adaptive scaling. (D0 > 0) (default = 1)
(If D0 <= 0, the default value will be used)
[in]Tuner1(Optional)
Parameter to check for false convergence. (0 <= Tuner1 <= 0.5) (default = 0.1)
(If Tuner1 < 0 or Tuner1 > 0.5, the default value will be used)
[in]Xctol(Optional)
X convergence tolerance. (0 <= Xctol <= 1) (default = Eps^(1/2))
(If Xctol < 0 or Xctol > 1, the default value will be used)
[in]Xftol(Optional)
False convergence tolerance. (0 <= Xftol <= 1) (default = 100*Eps)
(If Xftol < 0 or Xftol > 1, the default value will be used)
[in]Lmax0(Optional)
Maximum 2-norm allowed for scaled very first step. (Lmax0 > 0) (default = 1)
(If Lmax0 <= 0, the default value will be used)
[in]Lmaxs(Optional)
[in]Sctol(Optional)
Lmaxs and Sctol are the singular convergence test parameters. (Lmaxs > 0, 0 <= Sctol <= 1) (default: Lmaxs = 1, Sctol = 1e-10)
To test if the function reduction predicted for a step of length bounded by Lmaxs is at most Sctol*abs(f) (f is the function value at the start of the current iteration).
(If Lmaxs <= 0, the default value will be used)
(If Sctol < 0 or Sctol > 1, the default value will be used)
[in]Bias(Optional)
The bias parameter used in the dogleg trust region method. (0 <= Bias <= 1) (default = 0.8)
(If Bias < 0 or Bias > 1, the default value will be used)
Reference
netlib/port
Example Program
Find the minimum point of the following function (Rosenbrock function).
f(x1, x2) = 100(x2 - x1^2)^2 + (1 - x1)^2
The initial approximation is (x1, x2) = (-1.2, 1).
Sub Ex_Mnh_r()
Const N = 2
Dim X(N - 1) As Double, Info As Long
Dim YY As Double, YYp(N - 1) As Double, YYpd(N * (N + 1) / 2) As Double, IRev As Long
X(0) = -1.2: X(1) = 1
IRev = 0
Do
Call Mnh_r(N, X(), Info, YY, YYp(), YYpd(), IRev)
If IRev = 1 Then
YY = 100 * (X(1) - X(0) ^ 2) ^ 2 + (1 - X(0)) ^ 2
ElseIf IRev = 2 Then
YYp(0) = -400 * X(0) * (X(1) - X(0) ^ 2) + 2 * X(0) - 2
YYp(1) = 200 * X(1) - 200 * X(0) ^ 2
YYpd(0) = 1200 * X(0) ^ 2 - 400 * X(1) + 2
YYpd(1) = -400 * X(0)
YYpd(2) = 200
End If
Loop While IRev <> 0
Debug.Print "X1, X2 =", X(0), X(1)
Debug.Print "Info =", Info
End Sub
Example Results
X1, X2 = 0.999999999999989 0.999999999999978
Info = 0