Description Usage Arguments Details Value References See Also Examples
Non-Monotone spectral approach for Solving Large-Scale Nonlinear Systems of Equations
1 2 3 |
fn |
a function that takes a real vector as argument and returns a real vector of same length (see details). |
par |
A real vector argument to |
method |
An integer (1, 2, or 3) specifying which Barzilai-Borwein steplength to use. The default is 2. See *Details*. |
control |
A list of control parameters. See *Details*. |
quiet |
A logical variable (TRUE/FALSE). If |
alertConvergence |
A logical variable. With the default |
... |
Additional arguments passed to |
The function sane
implements a non-monotone spectral residual method
for finding a root of nonlinear systems. It stands for "spectral approach
for nonlinear equations".
It differs from the function dfsane
in that it requires an
approximation of a directional derivative at every iteration of the merit
function F(x)^t F(x).
R adaptation, with significant modifications, by Ravi Varadhan, Johns Hopkins University (March 25, 2008), from the original FORTRAN code of La Cruz and Raydan (2003).
A major modification in our R adaptation of the original FORTRAN code is the
availability of 3 different options for Barzilai-Borwein (BB) steplengths:
method = 1
is the BB
steplength used in LaCruz and Raydan (2003); method = 2
is equivalent to
the other steplength proposed in Barzilai and Borwein's (1988) original paper.
Finally, method = 3
, is a new steplength, which is equivalent to that
first proposed in Varadhan and Roland (2008) for accelerating the EM algorithm.
In fact, Varadhan and Roland (2008) considered 3 equivalent steplength schemes
in their EM acceleration work. Here, we have chosen method = 2
as the "default" method, as it generally performed better than the other
schemes in our numerical experiments.
Argument control
is a list specifing any changes to default values of
algorithm control parameters. Note that the names of these must be
specified completely. Partial matching will not work.
Argument control
has the following components:
A positive integer, typically between 5-20, that controls the
monotonicity of the algorithm. M=1
would enforce strict monotonicity
in the reduction of L2-norm of fn
, whereas larger values allow for
more non-monotonicity. Global convergence under non-monotonicity is ensured
by enforcing the Grippo-Lampariello-Lucidi condition (Grippo et al. 1986) in a
non-monotone line-search algorithm. Values of M
between 5 to 20 are
generally good, although some problems may require a much larger M.
The default is M = 10
.
The maximum number of iterations. The default is
maxit = 1500
.
The absolute convergence tolerance on the residual L2-norm
of fn
. Convergence is declared
when sqrt(sum(F(x)^2) / npar) < tol.
Default is tol = 1.e-07
.
A logical variable (TRUE/FALSE). If TRUE
, information on
the progress of solving the system is produced.
Default is trace = !quiet
.
An integer that controls the frequency of tracing
when trace=TRUE
. Default is triter=10
, which means that
the L2-norm of fn
is printed at every 10-th iteration.
An integer. Algorithm is terminated when no progress has been
made in reducing the merit function for noimp
consecutive iterations.
Default is noimp=100
.
A logical variable that dictates whether the Nelder-Mead algorithm
in optim
will be called upon to improve user-specified starting value.
Default is NM=FALSE
.
A logical variable that dictates whether the low-memory L-BFGS-B
algorithm in optim
will be called after certain types of unsuccessful
termination of sane
. Default is BFGS=FALSE
.
A list with the following components:
par |
The best set of parameters that solves the nonlinear system. |
residual |
L2-norm of the function evaluated at |
fn.reduction |
Reduction in the L2-norm of the function from the initial L2-norm. |
feval |
Number of times |
iter |
Number of iterations taken by the algorithm. |
convergence |
An integer code indicating type of convergence.
|
message |
A text message explaining which termination criterion was used. |
J Barzilai, and JM Borwein (1988), Two-point step size gradient methods, IMA J Numerical Analysis, 8, 141-148.
L Grippo, F Lampariello, and S Lucidi (1986), A nonmonotone line search technique for Newton's method, SIAM J on Numerical Analysis, 23, 707-716.
W LaCruz, and M Raydan (2003), Nonmonotone spectral methods for large-scale nonlinear systems, Optimization Methods and Software, 18, 583-599.
R Varadhan and C Roland (2008), Simple and globally-convergent methods for accelerating the convergence of any EM algorithm, Scandinavian J Statistics.
R Varadhan and PD Gilbert (2009), BB: An R Package for Solving a Large System of Nonlinear Equations and for Optimizing a High-Dimensional Nonlinear Objective Function, J. Statistical Software, 32:4, http://www.jstatsoft.org/v32/i04/
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 | trigexp <- function(x) {
# Test function No. 12 in the Appendix of LaCruz and Raydan (2003)
n <- length(x)
F <- rep(NA, n)
F[1] <- 3*x[1]^2 + 2*x[2] - 5 + sin(x[1] - x[2]) * sin(x[1] + x[2])
tn1 <- 2:(n-1)
F[tn1] <- -x[tn1-1] * exp(x[tn1-1] - x[tn1]) + x[tn1] * ( 4 + 3*x[tn1]^2) +
2 * x[tn1 + 1] + sin(x[tn1] - x[tn1 + 1]) * sin(x[tn1] + x[tn1 + 1]) - 8
F[n] <- -x[n-1] * exp(x[n-1] - x[n]) + 4*x[n] - 3
F
}
p0 <- rnorm(50)
sane(par=p0, fn=trigexp)
sane(par=p0, fn=trigexp, method=1)
######################################
brent <- function(x) {
n <- length(x)
tnm1 <- 2:(n-1)
F <- rep(NA, n)
F[1] <- 3 * x[1] * (x[2] - 2*x[1]) + (x[2]^2)/4
F[tnm1] <- 3 * x[tnm1] * (x[tnm1+1] - 2 * x[tnm1] + x[tnm1-1]) +
((x[tnm1+1] - x[tnm1-1])^2) / 4
F[n] <- 3 * x[n] * (20 - 2 * x[n] + x[n-1]) + ((20 - x[n-1])^2) / 4
F
}
p0 <- sort(runif(50, 0, 10))
sane(par=p0, fn=brent, control=list(trace=FALSE))
sane(par=p0, fn=brent, control=list(M=200, trace=FALSE))
|
Iteration: 0 ||F(x0)||: 134.0099
iteration: 10 ||F(xn)|| = 0.5632473
iteration: 20 ||F(xn)|| = 0.0002056251
$par
[1] 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1
[39] 1 1 1 1 1 1 1 1 1 1 1 1
$residual
[1] 2.374329e-08
$fn.reduction
[1] 947.5931
$feval
[1] 59
$iter
[1] 29
$convergence
[1] 0
$message
[1] "Successful convergence"
Iteration: 0 ||F(x0)||: 134.0099
iteration: 10 ||F(xn)|| = 0.4085091
iteration: 20 ||F(xn)|| = 0.0001269341
$par
[1] 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1
[39] 1 1 1 1 1 1 1 1 1 1 1 1
$residual
[1] 3.062008e-08
$fn.reduction
[1] 947.5931
$feval
[1] 59
$iter
[1] 29
$convergence
[1] 0
$message
[1] "Successful convergence"
$par
[1] 0.9697797 1.6932030 2.3260891 2.9062578 3.4501941 3.9669560
[7] 4.4622035 4.9397769 5.4024368 5.8522534 6.2908316 6.7194479
[13] 7.1391395 7.5507642 7.9550420 8.3525849 8.7439189 9.1295003
[19] 9.5097285 9.8849556 10.2554941 10.6216229 10.9835925 11.3416289
[25] 11.6959369 12.0467028 12.3940970 12.7382754 13.0793817 13.4175481
[31] 13.7528971 14.0855424 14.4155896 14.7431369 15.0682763 15.3910938
[37] 15.7116700 16.0300806 16.3463968 16.6606857 16.9730107 17.2834317
[43] 17.5920052 17.8987850 18.2038221 18.5071648 18.8088591 19.1089490
[49] 19.4074762 19.7044805
$residual
[1] 3.701602e-08
$fn.reduction
[1] 310.1166
$feval
[1] 2060
$iter
[1] 858
$convergence
[1] 0
$message
[1] "Successful convergence"
$par
[1] 0.9697795 1.6932027 2.3260887 2.9062574 3.4501936 3.9669555
[7] 4.4622029 4.9397763 5.4024361 5.8522527 6.2908309 6.7194471
[13] 7.1391387 7.5507634 7.9550412 8.3525841 8.7439180 9.1294995
[19] 9.5097277 9.8849548 10.2554933 10.6216221 10.9835917 11.3416281
[25] 11.6959361 12.0467021 12.3940963 12.7382747 13.0793810 13.4175474
[31] 13.7528965 14.0855418 14.4155890 14.7431363 15.0682758 15.3910933
[37] 15.7116696 16.0300802 16.3463964 16.6606854 16.9730104 17.2834314
[43] 17.5920049 17.8987848 18.2038219 18.5071646 18.8088590 19.1089489
[49] 19.4074761 19.7044805
$residual
[1] 7.156854e-08
$fn.reduction
[1] 310.1166
$feval
[1] 1400
$iter
[1] 699
$convergence
[1] 0
$message
[1] "Successful convergence"
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.