marqLevAlg: A parallelized general-purpose optimization based on...

View source: R/marqLevAlg.R

marqLevAlgR Documentation

A parallelized general-purpose optimization based on Marquardt-Levenberg algorithm

Description

This algorithm provides a numerical solution to the problem of unconstrained local optimization. This is more efficient than the Gauss-Newton-like algorithm when starting from points very far from the final maximum. A new convergence test is implemented (RDM) in addition to the usual stopping criterion : stopping rule is when the gradients are small enough in the parameters metric (GH-1G).

Usage

marqLevAlg(
  b,
  m = FALSE,
  fn,
  gr = NULL,
  hess = NULL,
  maxiter = 500,
  epsa = 1e-04,
  epsb = 1e-04,
  epsd = 1e-04,
  partialH = NULL,
  digits = 8,
  print.info = FALSE,
  blinding = TRUE,
  multipleTry = 25,
  nproc = 1,
  clustertype = NULL,
  file = "",
  .packages = NULL,
  minimize = TRUE,
  ...
)

mla(
  b,
  m = FALSE,
  fn,
  gr = NULL,
  hess = NULL,
  maxiter = 500,
  epsa = 1e-04,
  epsb = 1e-04,
  epsd = 1e-04,
  partialH = NULL,
  digits = 8,
  print.info = FALSE,
  blinding = TRUE,
  multipleTry = 25,
  nproc = 1,
  clustertype = NULL,
  file = "",
  .packages = NULL,
  minimize = TRUE,
  ...
)

Arguments

b

an optional vector containing the initial values for the parameters. Default is 0.1 for every parameter.

m

number of parameters. Optional if b is specified.

fn

the function to be optimized, with first argument the vector of parameters over which optimization is to take place (argument b). It should return a scalar result.

gr

a function to return the gradient value for a specific point. If missing, finite-difference approximation will be used.

hess

a function to return the hessian matrix for a specific point. If missing, finite-difference approximation will be used.

maxiter

optional maximum number of iterations for the marqLevAlg iterative algorithm. Default is 500.

epsa

optional threshold for the convergence criterion based on the parameter stability. Default is 0.0001.

epsb

optional threshold for the convergence criterion based on the objective function stability. Default is 0.0001.

epsd

optional threshold for the relative distance to maximum. This criterion has the nice interpretation of estimating the ratio of the approximation error over the statistical error, thus it can be used for stopping the iterative process whathever the problem. Default is 0.0001.

partialH

optional vector giving the indexes of the parameters to be dropped from the Hessian matrix to define the relative distance to maximum/minimum. If specified, this option will only be considered at iterations where the two first convergence criteria are satisfied (epsa and epsb) and if the total Hessian is not invertible. By default, no partial Hessian is defined.

digits

number of digits to print in outputs. Default value is 8.

print.info

logical indicating if information about computation should be reported at each iteration. Default value is FALSE.

blinding

logical. Equals to TRUE if the algorithm is allowed to go on in case of an infinite or not definite value of function. Default value is FALSE.

multipleTry

integer, different from 1 if the algorithm is allowed to go for the first iteration in case of an infinite or not definite value of gradients or hessian. As many tries as requested in multipleTry will be done by changing the starting point of the algorithm. Default value is 25.

nproc

number of processors for parallel computing

clustertype

one of the supported types from makeCluster

file

optional character giving the name of the file where the outputs of each iteration should be written (if print.info=TRUE).

.packages

for parallel setting only, packages used in the fn function

minimize

logical indicating if the fn function should be minimized or maximized. By default minimize=TRUE, the function is minimized.

...

other arguments of the fn, gr and hess functions

Details

Convergence criteria are very strict as they are based on derivatives of the objective function in addition to the parameter and objective function stability. In some cases, the program may not converge and reach the maximum number of iterations fixed at 500. In this case, the user should check that parameter estimates at the last iteration are not on the boundaries of the parameter space. If the parameters are on the boundaries of the parameter space, the identifiability of the model should be assessed. If not, the program should be run again with other initial values, with a higher maximum number of iterations or less strict convergence tolerances. An alternative is to remove some parameters from the Hessian matrix.

Value

cl

summary of the call to the function marqLevAlg.

ni

number of marqLevAlg iterations before reaching stopping criterion.

istop

status of convergence: =1 if the convergence criteria were satisfied, =2 if the maximum number of iterations was reached, =3 if convergence criteria with partial Hessian matrix were satisfied, =4 if the algorithm encountered a problem in the function computation.

v

if istop=1 or istop=3, vector containing the upper triangle matrix of variance-covariance estimates at the stopping point. Otherwise v contains the second derivatives of the fn function with respect to the parameters.

grad

vector containing the gradient at the stopping point.

fn.value

function evaluation at the stopping point.

b

stopping point value.

ca

convergence criteria for parameters stabilisation.

cb

convergence criteria for function stabilisation.

rdm

convergence criteria on the relative distance to minimum (or maximum).

time

a running time.

Author(s)

Melanie Prague, Viviane Philipps, Cecile Proust-Lima, Boris Hejblum, Daniel Commenges, Amadou Diakite

References

marqLevAlg Algorithm

Donald W. marquardt An algorithm for Least-Squares Estimation of Nonlinear Parameters. Journal of the Society for Industrial and Applied Mathematics, Vol. 11, No. 2. (Jun, 1963), pp. 431-441.

Convergence criteria : Relative distance to Minimim (or Maximum)

Commenges D. Jacqmin-Gadda H. Proust C. Guedj J. A Newton-like algorithm for likelihood maximization the robust-variance scoring algorithm arxiv:math/0610402v2 (2006)

Examples



### example 1
### initial values
b <- c(8,9)
### your function
f1 <- function(b){
	return(-4*(b[1]-5)^2-(b[2]-6)^2)
}
### gradient
g1 <- function(b){
     return(c(-8*(b[1]-5),-2*(b[2]-6)))
}
## Call
test1 <- mla(b=b, fn=f1, minimize=FALSE)

## Not run: 
microbenchmark::microbenchmark(mla(b=b, fn=f1, minimize=FALSE),
                              mla(b=b, fn=f1, minimize=FALSE, nproc=2),
                              mla(b=b, fn=f1, gr=g1, minimize=FALSE),
                              mla(b=b, fn=f1, gr=g1, minimize=FALSE, nproc=2),
                              times=10)
        
## End(Not run)



### example 2
## initial values
b <- c(3,-1,0,1)
## your function
f2 <- function(b){
	return(-((b[1]+10*b[2])^2+5*(b[3]-b[4])^2+(b[2]-2*b[3])^4+10*(b[1]-b[4])^4))
}
## Call
test2 <- mla(b=b, fn=f2, minimize=FALSE)
test2$b

test2_par <- mla(b=b, fn=f2, minimize=FALSE, nproc=2)
test2_par$b

## Not run: 
microbenchmark::microbenchmark(mla(b=b, fn=f2, minimize=FALSE),
                              mla(b=b, fn=f2, minimize=FALSE, nproc=2),
                              times=10)
        
## End(Not run)


## Not run: 
### example 3 : a linear mixed model
## the log-likelihood is implemented in the loglikLMM function
## the gradient is implemented in the gradLMM function

## data
Y <- dataEx$Y
X <- as.matrix(cbind(1,dataEx[,c("t","X1","X3")],dataEx$t*dataEx$X1))
ni <- as.numeric(table(dataEx$i))

## initial values
binit <- c(0,0,0,0,0,1,1)

## estimation in sequential mode, with numeric derivatives
estim <- marqLevAlg(b=binit, fn=loglikLMM, minimize=FALSE, X=X, Y=Y, ni=ni)
## estimation in parallel mode, with numeric derivatives
estim2 <- marqLevAlg(b=binit, fn=loglikLMM, minimize=FALSE, X=X, Y=Y, ni=ni, 
nproc=2, clustertype="FORK")
## estimation in sequential mode, with analytic gradient
estim3 <- marqLevAlg(b=binit, fn=loglikLMM, gr=gradLMM, minimize=FALSE, X=X, Y=Y, ni=ni)
## estimation in parallel mode, with analytic gradient
estim4 <- marqLevAlg(b=binit, fn=loglikLMM, gr=gradLMM, minimize=FALSE, X=X, Y=Y, ni=ni, 
nproc=2, clustertype="FORK")

## End(Not run)

marqLevAlg documentation built on March 31, 2023, 6:33 p.m.